OOPS
oops/assimilation/SpectralLMP.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_SPECTRALLMP_H_
12 #define OOPS_ASSIMILATION_SPECTRALLMP_H_
13 
14 #include <algorithm>
15 #include <cmath>
16 #include <memory>
17 #include <string>
18 #include <utility>
19 #include <vector>
20 
21 #include "eckit/config/LocalConfiguration.h"
23 #include "oops/util/dot_product.h"
24 #include "oops/util/Logger.h"
25 
26 namespace oops {
27 
28  /* The SPECTRAL preconditioner matrix \f$ P \f$ for the system matrix AB = (I + Ht Rinv H B)
29  * It is defined as \f$ P_(k+1) = P_k + U (D^{-1} - I_l) X' \f$
30  * where
31  * \f$ D \f$ is a l-by-l diagonal matrix whose diagonal entries are \f$ lambda_i \f$,
32  * \f$ X = [x_1, x_2, ..., x_l] \f$ and \f$ U = [u_1, u_2, ..., u_l] \f$
33  * are column matrices.
34  *
35  * The pairs \f$ (lambda_i,x_i) \f$ define the Ritz pairs of the matrix
36  * \f$ PA \f$ with respect to the subspace \f$ K(PA,Pr0) \f$ and \f$ U = Binv X \f$.
37  * Note that \f$ r0 \f$ is the initial residual.
38  *
39  * Note that if RitzPrecond_ is active, it applies the RITZ LMP. For details,
40  * we refer to S. Gratton, A. Sartenaer and J. Tshimanga, SIAM J. Optim.,21(3),912–935,2011
41  * and S. Gurol, PhD Manuscript, 2013.
42  */
43  /*!
44  * The solvers represent matrices as objects that implement a "multiply"
45  * method. This class defines objects that apply the preconditioner matrix \f$ P \f$.
46  */
47 
48 // -----------------------------------------------------------------------------
49 
50 template<typename VECTOR> class SpectralLMP {
51  public:
52  explicit SpectralLMP(const eckit::Configuration &);
54 
55  void update(std::vector<std::unique_ptr<VECTOR>> &, std::vector<std::unique_ptr<VECTOR>> &,
56  std::vector<std::unique_ptr<VECTOR>> &, std::vector<double> &, std::vector<double> &);
57 
58  void multiply(const VECTOR &, VECTOR &) const;
59 
60  private:
61  unsigned maxpairs_;
64  int maxouter_;
65  int update_;
66 
67  std::vector<std::unique_ptr<VECTOR>> X_;
68  std::vector<std::unique_ptr<VECTOR>> U_;
69  std::vector<double> eigvals_;
70  std::vector<double> omega_;
71 
72  // For RitzPrecond
73  std::vector<std::unique_ptr<VECTOR>> Y_;
74  std::vector<std::unique_ptr<VECTOR>> S_;
75  std::vector<VECTOR> Zlast_;
76  std::vector<VECTOR> Zhlast_;
77  std::vector<unsigned> usedpairIndx_;
78  std::vector<unsigned> zcount;
79 };
80 
81 // =============================================================================
82 
83 template<typename VECTOR>
84 SpectralLMP<VECTOR>::SpectralLMP(const eckit::Configuration & conf)
85  : maxpairs_(0), useoldpairs_(false), RitzPrecond_(false), maxouter_(0), update_(0)
86 {
87  maxouter_ = conf.getInt("nouter");
88  Log::info() << "SpectralLMP: maxouter = " << maxouter_ << std::endl;
89 
90  if (conf.has("preconditioner")) {
91  const eckit::LocalConfiguration precond(conf, "preconditioner");
92  if (precond.has("maxpairs")) maxpairs_ = precond.getInt("maxpairs");
93 
94  if (maxpairs_ > 0) {
95  Log::info() << "SpectralLMP: maxpairs_ = " << maxpairs_ << std::endl;
96 
97  std::string ritz;
98  if (precond.has("ritz")) ritz = precond.getString("ritz");
99  RitzPrecond_ = (ritz == "on" || ritz == "true");
100  Log::info() << "SpectralLMP: RitzPrecond_ = " << RitzPrecond_ << std::endl;
101 
102  std::string old;
103  if (precond.has("useoldpairs")) old = precond.getString("useoldpairs");
104  useoldpairs_ = (old == "on" || old == "true");
105  Log::info() << "SpectralLMP: useoldpairs_ = " << useoldpairs_ << std::endl;
106  }
107  }
108 }
109 
110 // -----------------------------------------------------------------------------
111 
112 template<typename VECTOR>
113 void SpectralLMP<VECTOR>::update(std::vector<std::unique_ptr<VECTOR>> & Zv,
114  std::vector<std::unique_ptr<VECTOR>> & Zhl,
115  std::vector<std::unique_ptr<VECTOR>> & Zl,
116  std::vector<double> & alphas,
117  std::vector<double> & betas) {
118 // If useoldpairs = false, use only current information
119  if (!useoldpairs_) {
120  eigvals_.clear();
121  omega_.clear();
122  X_.clear();
123  U_.clear();
124  usedpairIndx_.clear();
125  Zlast_.clear();
126  Zhlast_.clear();
127  }
128 
129  ++update_;
130  Log::debug() << "SpectralLMP size Zl = " << Zl.size() << std::endl;
131  Log::debug() << "SpectralLMP size Zhl = " << Zhl.size() << std::endl;
132  Log::debug() << "SpectralLMP size alpha = " << alphas.size() << std::endl;
133  Log::debug() << "SpectralLMP size beta = " << betas.size() << std::endl;
134  const unsigned nvec = Zl.size();
135  ASSERT(Zhl.size() == nvec);
136 
137  if (nvec > 0) {
138  Log::info() << std::endl;
139  Log::info() << "SpectralLMP: update " << update_ << ", max = " << maxouter_-1 << std::endl;
140  ASSERT(alphas.size() == nvec-1);
141  ASSERT(betas.size() == nvec-1);
142  std::vector<double> evals;
143  std::vector< std::vector<double> > evecs;
144 
145 // Compute spectrum of tri-diagonal matrix
146  TriDiagSpectrum(alphas, betas, evals, evecs);
147 
148 // Determine the converged eigenvalues
149  std::vector<unsigned> convIndx;
150  unsigned oldpairIndx = 0;
151  for (unsigned jiter = 0; jiter < nvec-1; ++jiter) {
152  double erritz = evecs[jiter][nvec-2]*betas[nvec-2];
153  double lambda = evals[jiter];
154 // Add new converged eigenvalues
155  if (std::abs(erritz) <= 0.001*lambda) {
156  convIndx.push_back(jiter);
157  eigvals_.push_back(lambda);
158  omega_.push_back(erritz/lambda);
159 // If maximum number of pairs are exceeded, remove old information from
160 // the beginning
161  if (convIndx.size() > maxpairs_) {
162  convIndx.erase(convIndx.begin());
163  }
164  if (eigvals_.size() > maxpairs_) {
165  eigvals_.erase(eigvals_.begin());
166  omega_.erase(omega_.begin());
167  oldpairIndx += 1;
168  }
169  }
170  }
171 // Keep information on how many pairs are used at each minimization loop
172  usedpairIndx_.push_back(convIndx.size());
173  unsigned reduc = oldpairIndx;
174  for (unsigned kiter = 0; kiter < usedpairIndx_.size()-1; ++kiter) {
175  if (reduc <= usedpairIndx_[kiter]) {
176  usedpairIndx_[kiter] -= reduc;
177  break;
178  } else {
179  reduc -= usedpairIndx_[kiter];
180  usedpairIndx_[kiter] = 0;
181  }
182  }
183  unsigned counter = 0;
184  for (unsigned kiter = 0; kiter < usedpairIndx_.size(); ++kiter) {
185  counter += usedpairIndx_[kiter];
186  }
187  ASSERT(counter <= maxpairs_);
188  Log::info() << "SpectralLMP: Number of inner iterations : " << Zl.size()-1 << std::endl;
189  Log::info() << "SpectralLMP: Number of maximum pairs : " << maxpairs_ << std::endl;
190  Log::info() << "SpectralLMP: Number of used total pairs : " << counter << std::endl;
191  Log::info() << "SpectralLMP: Number of new converged pairs : " << convIndx.size() << std::endl;
192  Log::info() << "SpectralLMP: Number of removed old pairs : " << oldpairIndx << std::endl;
193  for (unsigned kiter = 0; kiter < usedpairIndx_.size(); ++kiter) {
194  Log::info() << "SpectralLMP: Number of used pairs from outer loop " << kiter + 1
195  << " : " << usedpairIndx_[kiter] << std::endl;
196  }
197 
198 // Calculate/Update the matrix X and U
199 
200  if (X_.size() != 0) {
201 // Remove the first oldpairIndx elements
202  unsigned minsize = std::min(oldpairIndx, maxpairs_);
203  unsigned xsize = X_.size();
204  X_.erase(X_.begin(), X_.begin() + std::min(xsize, minsize));
205  U_.erase(U_.begin(), U_.begin() + std::min(xsize, minsize));
206  }
207  for (unsigned jiter = 0; jiter < convIndx.size(); ++jiter) {
208  std::unique_ptr<VECTOR> ww(new VECTOR(*Zl[0], false));
209  for (unsigned iiter = 0; iiter < nvec-1; ++iiter) {
210  ww->axpy(evecs[convIndx[jiter]][iiter], *Zl[iiter]);
211  }
212 // Add new information
213  X_.emplace_back(std::move(ww));
214  }
215  for (unsigned jiter = 0; jiter < convIndx.size(); ++jiter) {
216  std::unique_ptr<VECTOR> ww(new VECTOR(*Zl[0], false));
217  for (unsigned iiter = 0; iiter < nvec-1; ++iiter) {
218  ww->axpy(evecs[convIndx[jiter]][iiter], *Zhl[iiter]);
219  }
220 // Add new information
221  U_.emplace_back(std::move(ww));
222  }
223 
224  if (RitzPrecond_) {
225  Zlast_.push_back(*Zl[nvec-1]);
226  Zhlast_.push_back(*Zhl[nvec-1]);
227 
228 // Calculate the matrix Y = [U1*omega1, ..., Uk*omegak]
229  Y_.clear();
230  zcount.clear();
231  unsigned kk = 0;
232  for (unsigned kiter = 0; kiter < usedpairIndx_.size(); ++kiter) {
233  if (usedpairIndx_[kiter] != 0) {
234  zcount.push_back(kiter);
235  std::unique_ptr<VECTOR> ww(new VECTOR(*Zl[0], false));
236  for (unsigned jiter = 0; jiter < usedpairIndx_[kiter]; ++jiter) {
237  ww->axpy(omega_[jiter + kk], *U_[jiter + kk]);
238  }
239  kk += usedpairIndx_[kiter];
240  Y_.emplace_back(std::move(ww));
241  }
242  }
243 
244 // Calculate the matrix S = [X1*omega1, ..., Xk*omegak]
245  S_.clear();
246  kk = 0;
247  for (unsigned kiter = 0; kiter < usedpairIndx_.size(); ++kiter) {
248  if (usedpairIndx_[kiter] != 0) {
249  std::unique_ptr<VECTOR> ww(new VECTOR(*Zl[0], false));
250  for (unsigned jiter = 0; jiter < usedpairIndx_[kiter]; ++jiter) {
251  ww->axpy(omega_[jiter + kk], *X_[jiter + kk]);
252  }
253  kk += usedpairIndx_[kiter];
254  S_.emplace_back(std::move(ww));
255  }
256  }
257  }
258 
259 // Release remaining memory
260  Zv.clear();
261  Zl.clear();
262  Zhl.clear();
263  alphas.clear();
264  betas.clear();
265  }
266 }
267 
268 // -----------------------------------------------------------------------------
269 
270 template<typename VECTOR>
271 void SpectralLMP<VECTOR>::multiply(const VECTOR & a, VECTOR & b) const {
272  b = a; // P_0 = I
273 
274  for (unsigned iiter = 0; iiter < eigvals_.size(); ++iiter) {
275  double zeval = std::min(10.0, eigvals_[iiter]);
276 // double zeval = eigvals_[iiter];
277  double zz = 1.0/zeval - 1.0;
278  b.axpy(zz*dot_product(a, *X_[iiter]), *U_[iiter]);
279  }
280 
281  if (RitzPrecond_ && eigvals_.size() != 0) {
282 // Sk'*a
283  std::vector<double> sta;
284  sta.clear();
285  for (unsigned iiter = 0; iiter < S_.size(); ++iiter) {
286  sta.push_back(dot_product(*S_[iiter], a));
287  }
288 
289 // Yk (sta(k) - Zlast' a)
290  for (unsigned kiter = 0; kiter < S_.size(); ++kiter) {
291  for (unsigned iiter = 0; iiter < eigvals_.size(); ++iiter) {
292  double wxap = sta[kiter] - dot_product(a, Zlast_[zcount[kiter]]);
293  b.axpy(wxap, *Y_[kiter]);
294  }
295  }
296 
297 // -Zhlast sta
298  for (unsigned kiter = 0; kiter < S_.size(); ++kiter) {
299  b.axpy(-sta[kiter], Zhlast_[zcount[kiter]]);
300  }
301  }
302 }
303 
304 // -----------------------------------------------------------------------------
305 
306 } // namespace oops
307 
308 #endif // OOPS_ASSIMILATION_SPECTRALLMP_H_
The solvers represent matrices as objects that implement a "multiply" method.
std::vector< std::unique_ptr< VECTOR > > U_
void multiply(const VECTOR &, VECTOR &) const
SpectralLMP(const eckit::Configuration &)
std::vector< std::unique_ptr< VECTOR > > Y_
std::vector< std::unique_ptr< VECTOR > > X_
std::vector< double > eigvals_
std::vector< std::unique_ptr< VECTOR > > S_
std::vector< double > omega_
std::vector< unsigned > usedpairIndx_
std::vector< unsigned > zcount
void update(std::vector< std::unique_ptr< VECTOR >> &, std::vector< std::unique_ptr< VECTOR >> &, std::vector< std::unique_ptr< VECTOR >> &, std::vector< double > &, std::vector< double > &)
std::vector< VECTOR > Zhlast_
std::vector< VECTOR > Zlast_
The namespace for the main oops code.
void TriDiagSpectrum(const std::vector< double > &diag, const std::vector< double > &sub, std::vector< double > &evals, std::vector< std::vector< double > > &evecs)