OOPS
QNewtonLMP.h
Go to the documentation of this file.
1 /*
2  * (C) Copyright 2009-2016 ECMWF.
3  *
4  * This software is licensed under the terms of the Apache Licence Version 2.0
5  * which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
6  * In applying this licence, ECMWF does not waive the privileges and immunities
7  * granted to it by virtue of its status as an intergovernmental organisation nor
8  * does it submit to any jurisdiction.
9  */
10 
11 #ifndef OOPS_ASSIMILATION_QNEWTONLMP_H_
12 #define OOPS_ASSIMILATION_QNEWTONLMP_H_
13 
14 #include <algorithm>
15 #include <string>
16 #include <vector>
17 
18 #include "eckit/config/LocalConfiguration.h"
19 #include "oops/util/dot_product.h"
20 #include "oops/util/Logger.h"
21 
22 
23 namespace oops {
24 
25  /* The QUASI-NEWTON Limited Memory preconditioner matrix \f$ C \f$
26  * for the system matrix \f$ AB = (I + Ht Rinv H B) \f$
27  *
28  * It is defined as
29  * \f$ C_(k+1) = (I - rho_k ph_k q_k^T) C_k (I - rho_k q_k p_k^T) + rho_k ph_k p_k^T\f$
30  *
31  * For details, we refer to S. Gratton, A. Sartenaer and J. Tshimanga,
32  * SIAM J. Optim.,21(3),912–935,2011 and S. Gurol, PhD Manuscript, 2013.
33  *
34  * The solvers represent matrices as objects that implement a "multiply"
35  * method. This class defines objects that apply the preconditioner matrix \f$ P \f$.
36  */
37 
38 // -----------------------------------------------------------------------------
39 
40 template<typename VECTOR, typename BMATRIX> class QNewtonLMP {
41  public:
42  explicit QNewtonLMP(const eckit::Configuration &);
44 
45  void push(const VECTOR &, const VECTOR &, const VECTOR &, const double &);
46  void update(const BMATRIX & B);
47 
48  void multiply(const VECTOR &, VECTOR &) const;
49 
50  // Matrix-vector products with the transpose of the preconditioner
51  void tmultiply(const VECTOR &, VECTOR &) const;
52 
53  private:
54  unsigned maxpairs_;
55  unsigned maxnewpairs_;
57  int maxouter_;
58  int update_;
59 
60  std::vector<VECTOR> P_;
61  std::vector<VECTOR> Ph_;
62  std::vector<VECTOR> AP_;
63  std::vector<VECTOR> BAP_;
64  std::vector<double> rhos_;
65  std::vector<unsigned> usedpairIndx_;
66 
67  std::vector<VECTOR> savedP_;
68  std::vector<VECTOR> savedPh_;
69  std::vector<VECTOR> savedAP_;
70  std::vector<double> savedrhos_;
71 };
72 
73 // =============================================================================
74 
75 template<typename VECTOR, typename BMATRIX>
76 QNewtonLMP<VECTOR, BMATRIX>::QNewtonLMP(const eckit::Configuration & conf)
77  : maxpairs_(0), maxnewpairs_(0), useoldpairs_(false), maxouter_(0), update_(1)
78 {
79  maxouter_ = conf.getInt("nouter");
80  Log::info() << "QNewtonLMP: maxouter : " << maxouter_ << std::endl;
81 
82  if (conf.has("preconditioner")) {
83  const eckit::LocalConfiguration precond(conf, "preconditioner");
84  if (precond.has("maxpairs")) maxpairs_ = precond.getInt("maxpairs");
85  if (maxpairs_ > 0) {
86  std::string old;
87  if (precond.has("useoldpairs")) old = precond.getString("useoldpairs");
88  useoldpairs_ = (old == "on" || old == "true");
89  if (useoldpairs_ && precond.has("maxnewpairs")) {
90  maxnewpairs_ = precond.getInt("maxnewpairs");
91  maxnewpairs_ = std::min(maxpairs_, maxnewpairs_);
92  } else {
94  }
95  }
96  }
97 }
98 
99 // -----------------------------------------------------------------------------
100 
101 template<typename VECTOR, typename BMATRIX>
102 void QNewtonLMP<VECTOR, BMATRIX>::push(const VECTOR & p, const VECTOR & ph,
103  const VECTOR & ap, const double & rho) {
104  ASSERT(savedP_.size() <= maxnewpairs_);
105  if (maxnewpairs_ > 0 && update_ < maxouter_) {
106  if (savedP_.size() == maxnewpairs_) {
107  savedP_.erase(savedP_.begin());
108  savedPh_.erase(savedPh_.begin());
109  savedAP_.erase(savedAP_.begin());
110  savedrhos_.erase(savedrhos_.begin());
111  }
112  savedP_.push_back(p);
113  savedPh_.push_back(ph);
114  savedAP_.push_back(ap);
115  savedrhos_.push_back(1.0/rho);
116  }
117 }
118 
119 // -----------------------------------------------------------------------------
120 
121 template<typename VECTOR, typename BMATRIX>
122 void QNewtonLMP<VECTOR, BMATRIX>::update(const BMATRIX & Bmat) {
123  Log::debug() << "QNewtonLMP size saved Ph = " << savedPh_.size() << std::endl;
124  Log::debug() << "QNewtonLMP size saved P = " << savedP_.size() << std::endl;
125  Log::debug() << "QNewtonLMP size saved AP = " << savedAP_.size() << std::endl;
126  Log::debug() << "QNewtonLMP size saved rhos = " << savedrhos_.size() << std::endl;
127 
128  const unsigned nvec = savedPh_.size();
129  ASSERT(savedP_.size() == nvec);
130  ASSERT(savedAP_.size() == nvec);
131  ASSERT(savedrhos_.size() == nvec);
132 
133  if (!useoldpairs_ || nvec >= maxpairs_) {
134  Ph_.clear();
135  P_.clear();
136  AP_.clear();
137  BAP_.clear();
138  rhos_.clear();
139  for (unsigned kiter = 0; kiter < usedpairIndx_.size(); ++kiter) {
140  usedpairIndx_[kiter] = 0;
141  }
142  }
143 
144  if (nvec > 0 && update_ < maxouter_) {
145  Log::info() << "QNewtonLMP: update " << update_ << ", max = " << maxouter_-1 << std::endl;
146  unsigned newpairs = std::min(nvec, maxpairs_);
147  unsigned oldpairs = P_.size();
148  unsigned rmpairs = 0;
149 // First remove pairs we no longer need.
150  if (oldpairs + newpairs > maxpairs_) {
151  rmpairs = oldpairs + newpairs - maxpairs_;
152  for (unsigned jv = 0; jv < rmpairs; ++jv) {
153  Ph_.erase(Ph_.begin());
154  P_.erase(P_.begin());
155  AP_.erase(AP_.begin());
156  BAP_.erase(BAP_.begin());
157  rhos_.erase(rhos_.begin());
158  }
159  oldpairs -= rmpairs;
160 // Keep information on how many pairs are used at each minimization loop
161  unsigned removed = rmpairs;
162  for (unsigned kiter = 0; kiter < usedpairIndx_.size(); ++kiter) {
163  unsigned rmiter = std::min(usedpairIndx_[kiter], removed);
164  usedpairIndx_[kiter] -= rmiter;
165  removed -= rmiter;
166  }
167  ASSERT(removed == 0);
168  }
169  ASSERT(P_.size() == oldpairs);
170  ASSERT(oldpairs + newpairs <= maxpairs_);
171 
172 // Add the new pairs.
173  for (unsigned jv = nvec - newpairs; jv < nvec; ++jv) {
174  P_.push_back(savedP_[jv]);
175  Ph_.push_back(savedPh_[jv]);
176  AP_.push_back(savedAP_[jv]);
177  rhos_.push_back(savedrhos_[jv]);
178  // Save B*ap
179  VECTOR ww(savedAP_[jv]);
180  Bmat.multiply(savedAP_[jv], ww);
181  BAP_.push_back(ww);
182  }
183  ASSERT(P_.size() == oldpairs + newpairs);
184 
185  Log::info() << "Number of inner iterations : " << nvec << std::endl;
186  Log::info() << "Number of maximum pairs : " << maxpairs_ << std::endl;
187  Log::info() << "Number of used total pairs : " << oldpairs + newpairs << std::endl;
188  Log::info() << "Number of new pairs : " << newpairs << std::endl;
189  Log::info() << "Number of removed old pairs : " << rmpairs << std::endl;
190  for (unsigned kiter = 0; kiter < usedpairIndx_.size(); ++kiter) {
191  Log::info() << "Number of used pairs from outer loop " << kiter + 1
192  << " : " << usedpairIndx_[kiter];
193  }
194  }
195 
196  ++update_;
197  savedP_.clear();
198  savedPh_.clear();
199  savedAP_.clear();
200  savedrhos_.clear();
201 }
202 
203 // -----------------------------------------------------------------------------
204 
205 template<typename VECTOR, typename BMATRIX>
206 void QNewtonLMP<VECTOR, BMATRIX>::multiply(const VECTOR & a, VECTOR & b) const {
207  b = a;
208  const unsigned nvec = P_.size();
209  if (nvec != 0) {
210  std::vector<double> etas;
211  etas.clear();
212  for (int iiter = nvec-1; iiter >= 0; iiter--) {
213  double eta = dot_product(b, P_[iiter]);
214  eta *= rhos_[iiter];
215  etas.push_back(eta);
216  b.axpy(-eta, AP_[iiter]);
217  }
218  b *= dot_product(AP_[nvec-1], AP_[nvec-1])/dot_product(AP_[nvec-1], Ph_[nvec-1]);
219  for (unsigned iiter = 0; iiter < nvec; ++iiter) {
220  double sigma = dot_product(b, BAP_[iiter]);
221  sigma *= rhos_[iiter];
222  sigma -= etas[iiter];
223  b.axpy(-sigma, Ph_[iiter]);
224  }
225  }
226 }
227 
228 // -----------------------------------------------------------------------------
229 
230 template<typename VECTOR, typename BMATRIX>
231 void QNewtonLMP<VECTOR, BMATRIX>::tmultiply(const VECTOR & a, VECTOR & b) const {
232  b = a;
233  const unsigned nvec = P_.size();
234  if (nvec != 0) {
235  std::vector<double> etas;
236  etas.clear();
237  for (int iiter = nvec-1; iiter >= 0; iiter--) {
238  double eta = dot_product(b, Ph_[iiter]);
239  eta *= rhos_[iiter];
240  etas.push_back(eta);
241  b.axpy(-eta, BAP_[iiter]);
242  }
243  b *= dot_product(AP_[nvec-1], AP_[nvec-1])/dot_product(AP_[nvec-1], Ph_[nvec-1]);
244  for (unsigned iiter = 0; iiter < nvec; ++iiter) {
245  double sigma = dot_product(b, AP_[iiter]);
246  sigma *= rhos_[iiter];
247  sigma -= etas[iiter];
248  b.axpy(-sigma, P_[iiter]);
249  }
250  }
251 }
252 
253 // -----------------------------------------------------------------------------
254 
255 } // namespace oops
256 
257 #endif // OOPS_ASSIMILATION_QNEWTONLMP_H_
unsigned maxpairs_
Definition: QNewtonLMP.h:54
std::vector< unsigned > usedpairIndx_
Definition: QNewtonLMP.h:65
QNewtonLMP(const eckit::Configuration &)
Definition: QNewtonLMP.h:76
std::vector< VECTOR > Ph_
Definition: QNewtonLMP.h:61
void push(const VECTOR &, const VECTOR &, const VECTOR &, const double &)
Definition: QNewtonLMP.h:102
void multiply(const VECTOR &, VECTOR &) const
Definition: QNewtonLMP.h:206
std::vector< double > rhos_
Definition: QNewtonLMP.h:64
std::vector< VECTOR > savedAP_
Definition: QNewtonLMP.h:69
std::vector< VECTOR > BAP_
Definition: QNewtonLMP.h:63
unsigned maxnewpairs_
Definition: QNewtonLMP.h:55
void update(const BMATRIX &B)
Definition: QNewtonLMP.h:122
std::vector< VECTOR > P_
Definition: QNewtonLMP.h:60
std::vector< VECTOR > savedPh_
Definition: QNewtonLMP.h:68
std::vector< VECTOR > AP_
Definition: QNewtonLMP.h:62
void tmultiply(const VECTOR &, VECTOR &) const
Definition: QNewtonLMP.h:231
std::vector< double > savedrhos_
Definition: QNewtonLMP.h:70
std::vector< VECTOR > savedP_
Definition: QNewtonLMP.h:67
The namespace for the main oops code.