FV3 Bundle
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 <string>
17 #include <vector>
18 
19 #include <boost/ptr_container/ptr_vector.hpp>
20 #include <boost/scoped_ptr.hpp>
21 
22 #include "eckit/config/LocalConfiguration.h"
24 #include "oops/util/dot_product.h"
25 #include "oops/util/Logger.h"
26 
27 namespace oops {
28 
29  /* The SPECTRAL preconditioner matrix \f$ P \f$ for the system matrix AB = (I + Ht Rinv H B)
30  * It is defined as \f$ P_(k+1) = P_k + U (D^{-1} - I_l) X' \f$
31  * where
32  * \f$ D \f$ is a l-by-l diagonal matrix whose diagonal entries are \f$ lambda_i \f$,
33  * \f$ X = [x_1, x_2, ..., x_l] \f$ and \f$ U = [u_1, u_2, ..., u_l] \f$
34  * are column matrices.
35  *
36  * The pairs \f$ (lambda_i,x_i) \f$ define the Ritz pairs of the matrix
37  * \f$ PA \f$ with respect to the subspace \f$ K(PA,Pr0) \f$ and \f$ U = Binv X \f$.
38  * Note that \f$ r0 \f$ is the initial residual.
39  *
40  * Note that if RitzPrecond_ is active, it applies the RITZ LMP. For details,
41  * we refer to S. Gratton, A. Sartenaer and J. Tshimanga, SIAM J. Optim.,21(3),912–935,2011
42  * and S. Gurol, PhD Manuscript, 2013.
43  */
44  /*!
45  * The solvers represent matrices as objects that implement a "multiply"
46  * method. This class defines objects that apply the preconditioner matrix \f$ P \f$.
47  */
48 
49 // -----------------------------------------------------------------------------
50 
51 template<typename VECTOR> class SpectralLMP {
52  public:
53  explicit SpectralLMP(const eckit::Configuration &);
55 
56  void update(boost::ptr_vector<VECTOR> &, boost::ptr_vector<VECTOR> &,
57  boost::ptr_vector<VECTOR> &, std::vector<double> &, std::vector<double> &);
58 
59  void multiply(const VECTOR &, VECTOR &) const;
60 
61  private:
62  unsigned maxpairs_;
65  int maxouter_;
66  int update_;
67 
68  boost::ptr_vector<VECTOR> X_;
69  boost::ptr_vector<VECTOR> U_;
70  std::vector<double> eigvals_;
71  std::vector<double> omega_;
72 
73  // For RitzPrecond
74  boost::ptr_vector<VECTOR> Y_;
75  boost::ptr_vector<VECTOR> S_;
76  std::vector<VECTOR> Zlast_;
77  std::vector<VECTOR> Zhlast_;
78  std::vector<unsigned> usedpairIndx_;
79  std::vector<unsigned> zcount;
80 };
81 
82 // =============================================================================
83 
84 template<typename VECTOR>
85 SpectralLMP<VECTOR>::SpectralLMP(const eckit::Configuration & conf)
86  : maxpairs_(0), useoldpairs_(false), RitzPrecond_(false), maxouter_(0), update_(0)
87 {
88  maxouter_ = conf.getInt("nouter");
89  Log::info() << "SpectralLMP: maxouter = " << maxouter_ << std::endl;
90 
91  if (conf.has("preconditioner")) {
92  const eckit::LocalConfiguration precond(conf, "preconditioner");
93  if (precond.has("maxpairs")) maxpairs_ = precond.getInt("maxpairs");
94 
95  if (maxpairs_ > 0) {
96  Log::info() << "SpectralLMP: maxpairs_ = " << maxpairs_ << std::endl;
97 
98  std::string ritz;
99  if (precond.has("ritz")) ritz = precond.getString("ritz");
100  RitzPrecond_ = (ritz == "on" || ritz == "true");
101  Log::info() << "SpectralLMP: RitzPrecond_ = " << RitzPrecond_ << std::endl;
102 
103  std::string old;
104  if (precond.has("useoldpairs")) old = precond.getString("useoldpairs");
105  useoldpairs_ = (old == "on" || old == "true");
106  Log::info() << "SpectralLMP: useoldpairs_ = " << useoldpairs_ << std::endl;
107  }
108  }
109 }
110 
111 // -----------------------------------------------------------------------------
112 
113 template<typename VECTOR>
114 void SpectralLMP<VECTOR>::update(boost::ptr_vector<VECTOR> & Zv,
115  boost::ptr_vector<VECTOR> & Zhl,
116  boost::ptr_vector<VECTOR> & Zl,
117  std::vector<double> & alphas,
118  std::vector<double> & betas) {
119 // If useoldpairs = false, use only current information
120  if (!useoldpairs_) {
121  eigvals_.clear();
122  omega_.clear();
123  X_.clear();
124  U_.clear();
125  usedpairIndx_.clear();
126  Zlast_.clear();
127  Zhlast_.clear();
128  }
129 
130  ++update_;
131  Log::debug() << "SpectralLMP size Zl = " << Zl.size() << std::endl;
132  Log::debug() << "SpectralLMP size Zhl = " << Zhl.size() << std::endl;
133  Log::debug() << "SpectralLMP size alpha = " << alphas.size() << std::endl;
134  Log::debug() << "SpectralLMP size beta = " << betas.size() << std::endl;
135  const unsigned nvec = Zl.size();
136  ASSERT(Zhl.size() == nvec);
137 
138  if (nvec > 0) {
139  Log::info() << std::endl;
140  Log::info() << "SpectralLMP: update " << update_ << ", max = " << maxouter_-1 << std::endl;
141  ASSERT(alphas.size() == nvec-1);
142  ASSERT(betas.size() == nvec-1);
143  std::vector<double> evals;
144  std::vector< std::vector<double> > evecs;
145 
146 // Compute spectrum of tri-diagonal matrix
147  TriDiagSpectrum(alphas, betas, evals, evecs);
148 
149 // Determine the converged eigenvalues
150  std::vector<unsigned> convIndx;
151  unsigned oldpairIndx = 0;
152  for (unsigned jiter = 0; jiter < nvec-1; ++jiter) {
153  double erritz = evecs[jiter][nvec-2]*betas[nvec-2];
154  double lambda = evals[jiter];
155 // Add new converged eigenvalues
156  if (std::abs(erritz) <= 0.001*lambda) {
157  convIndx.push_back(jiter);
158  eigvals_.push_back(lambda);
159  omega_.push_back(erritz/lambda);
160 // If maximum number of pairs are exceeded, remove old information from
161 // the beginning
162  if (convIndx.size() > maxpairs_) {
163  convIndx.erase(convIndx.begin());
164  }
165  if (eigvals_.size() > maxpairs_) {
166  eigvals_.erase(eigvals_.begin());
167  omega_.erase(omega_.begin());
168  oldpairIndx += 1;
169  }
170  }
171  }
172 // Keep information on how many pairs are used at each minimization loop
173  usedpairIndx_.push_back(convIndx.size());
174  unsigned reduc = oldpairIndx;
175  for (unsigned kiter = 0; kiter < usedpairIndx_.size()-1; ++kiter) {
176  if (reduc <= usedpairIndx_[kiter]) {
177  usedpairIndx_[kiter] -= reduc;
178  break;
179  } else {
180  reduc -= usedpairIndx_[kiter];
181  usedpairIndx_[kiter] = 0;
182  }
183  }
184  unsigned counter = 0;
185  for (unsigned kiter = 0; kiter < usedpairIndx_.size(); ++kiter) {
186  counter += usedpairIndx_[kiter];
187  }
188  ASSERT(counter <= maxpairs_);
189  Log::info() << "SpectralLMP: Number of inner iterations : " << Zl.size()-1 << std::endl;
190  Log::info() << "SpectralLMP: Number of maximum pairs : " << maxpairs_ << std::endl;
191  Log::info() << "SpectralLMP: Number of used total pairs : " << counter << std::endl;
192  Log::info() << "SpectralLMP: Number of new converged pairs : " << convIndx.size() << std::endl;
193  Log::info() << "SpectralLMP: Number of removed old pairs : " << oldpairIndx << std::endl;
194  for (unsigned kiter = 0; kiter < usedpairIndx_.size(); ++kiter) {
195  Log::info() << "SpectralLMP: Number of used pairs from outer loop " << kiter + 1
196  << " : " << usedpairIndx_[kiter] << std::endl;
197  }
198 
199 // Calculate/Update the matrix X and U
200 
201  if (X_.size() != 0) {
202 // Remove the first oldpairIndx elements
203  unsigned minsize = std::min(oldpairIndx, maxpairs_);
204  unsigned xsize = X_.size();
205  X_.erase(X_.begin(), X_.begin() + std::min(xsize, minsize));
206  U_.erase(U_.begin(), U_.begin() + std::min(xsize, minsize));
207  }
208  for (unsigned jiter = 0; jiter < convIndx.size(); ++jiter) {
209  VECTOR * ww = new VECTOR(Zl[0], false);
210  for (unsigned iiter = 0; iiter < nvec-1; ++iiter) {
211  ww->axpy(evecs[convIndx[jiter]][iiter], Zl[iiter]);
212  }
213 // Add new information
214  X_.push_back(ww);
215  }
216  for (unsigned jiter = 0; jiter < convIndx.size(); ++jiter) {
217  VECTOR * ww = new VECTOR(Zl[0], false);
218  for (unsigned iiter = 0; iiter < nvec-1; ++iiter) {
219  ww->axpy(evecs[convIndx[jiter]][iiter], Zhl[iiter]);
220  }
221 // Add new information
222  U_.push_back(ww);
223  }
224 
225  if (RitzPrecond_) {
226  Zlast_.push_back(Zl[nvec-1]);
227  Zhlast_.push_back(Zhl[nvec-1]);
228 
229 // Calculate the matrix Y = [U1*omega1, ..., Uk*omegak]
230  Y_.clear();
231  zcount.clear();
232  unsigned kk = 0;
233  for (unsigned kiter = 0; kiter < usedpairIndx_.size(); ++kiter) {
234  if (usedpairIndx_[kiter] != 0) {
235  zcount.push_back(kiter);
236  VECTOR * ww = new VECTOR(Zl[0], false);
237  for (unsigned jiter = 0; jiter < usedpairIndx_[kiter]; ++jiter) {
238  ww->axpy(omega_[jiter + kk], U_[jiter + kk]);
239  }
240  kk += usedpairIndx_[kiter];
241  Y_.push_back(ww);
242  }
243  }
244 
245 // Calculate the matrix S = [X1*omega1, ..., Xk*omegak]
246  S_.clear();
247  kk = 0;
248  for (unsigned kiter = 0; kiter < usedpairIndx_.size(); ++kiter) {
249  if (usedpairIndx_[kiter] != 0) {
250  VECTOR * ww = new VECTOR(Zl[0], false);
251  for (unsigned jiter = 0; jiter < usedpairIndx_[kiter]; ++jiter) {
252  ww->axpy(omega_[jiter + kk], X_[jiter + kk]);
253  }
254  kk += usedpairIndx_[kiter];
255  S_.push_back(ww);
256  }
257  }
258  }
259 
260 // Release remaining memory
261  Zv.clear();
262  Zl.clear();
263  Zhl.clear();
264  alphas.clear();
265  betas.clear();
266  }
267 }
268 
269 // -----------------------------------------------------------------------------
270 
271 template<typename VECTOR>
272 void SpectralLMP<VECTOR>::multiply(const VECTOR & a, VECTOR & b) const {
273  b = a; // P_0 = I
274 
275  for (unsigned iiter = 0; iiter < eigvals_.size(); ++iiter) {
276  double zeval = std::min(10.0, eigvals_[iiter]);
277 // double zeval = eigvals_[iiter];
278  double zz = 1.0/zeval - 1.0;
279  b.axpy(zz*dot_product(a, X_[iiter]), U_[iiter]);
280  }
281 
282  if (RitzPrecond_ && eigvals_.size() != 0) {
283 // Sk'*a
284  std::vector<double> sta;
285  sta.clear();
286  for (unsigned iiter = 0; iiter < S_.size(); ++iiter) {
287  sta.push_back(dot_product(S_[iiter], a));
288  }
289 
290 // Yk (sta(k) - Zlast' a)
291  for (unsigned kiter = 0; kiter < S_.size(); ++kiter) {
292  for (unsigned iiter = 0; iiter < eigvals_.size(); ++iiter) {
293  double wxap = sta[kiter] - dot_product(a, Zlast_[zcount[kiter]]);
294  b.axpy(wxap, Y_[kiter]);
295  }
296  }
297 
298 // -Zhlast sta
299  for (unsigned kiter = 0; kiter < S_.size(); ++kiter) {
300  b.axpy(-sta[kiter], Zhlast_[zcount[kiter]]);
301  }
302  }
303 }
304 
305 // -----------------------------------------------------------------------------
306 
307 } // namespace oops
308 
309 #endif // OOPS_ASSIMILATION_SPECTRALLMP_H_
boost::ptr_vector< VECTOR > U_
Definition: SpectralLMP.h:69
void multiply(const VECTOR &, VECTOR &) const
Definition: SpectralLMP.h:272
std::vector< double > omega_
Definition: SpectralLMP.h:71
std::vector< VECTOR > Zlast_
Definition: SpectralLMP.h:76
Definition: conf.py:1
std::vector< double > eigvals_
Definition: SpectralLMP.h:70
std::vector< unsigned > zcount
Definition: SpectralLMP.h:79
boost::ptr_vector< VECTOR > X_
Definition: SpectralLMP.h:68
boost::ptr_vector< VECTOR > Y_
Definition: SpectralLMP.h:74
integer(long), parameter false
std::vector< unsigned > usedpairIndx_
Definition: SpectralLMP.h:78
boost::ptr_vector< VECTOR > S_
Definition: SpectralLMP.h:75
The solvers represent matrices as objects that implement a "multiply" method.
Definition: SpectralLMP.h:51
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)
logical debug
Definition: mpp.F90:1297
subroutine, public multiply(self, geom, xctl, xmod)
SpectralLMP(const eckit::Configuration &)
Definition: SpectralLMP.h:85
void update(boost::ptr_vector< VECTOR > &, boost::ptr_vector< VECTOR > &, boost::ptr_vector< VECTOR > &, std::vector< double > &, std::vector< double > &)
Definition: SpectralLMP.h:114
unsigned maxpairs_
Definition: SpectralLMP.h:62
subroutine, public info(self)
std::vector< VECTOR > Zhlast_
Definition: SpectralLMP.h:77
************************************************************************GNU Lesser General Public License **This file is part of the GFDL Flexible Modeling System(FMS). ! *! *FMS is free software without even the implied warranty of MERCHANTABILITY or *FITNESS FOR A PARTICULAR PURPOSE See the GNU General Public License *for more details **You should have received a copy of the GNU Lesser General Public *License along with FMS If see< http:! ***********************************************************************!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! !! MPP_TRANSMIT !! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! subroutine MPP_TRANSMIT_(put_data, put_len, to_pe, get_data, get_len, from_pe, block, tag, recv_request, send_request)!a message-passing routine intended to be reminiscent equally of both MPI and SHMEM!put_data and get_data are contiguous MPP_TYPE_ arrays!at each call, your put_data array is put to to_pe 's get_data! your get_data array is got from from_pe 's put_data!i.e we assume that typically(e.g updating halo regions) each PE performs a put _and_ a get!special PE designations:! NULL_PE:to disable a put or a get(e.g at boundaries)! ANY_PE:if remote PE for the put or get is to be unspecific! ALL_PES:broadcast and collect operations(collect not yet implemented)!ideally we would not pass length, but this f77-style call performs better(arrays passed by address, not descriptor)!further, this permits< length > contiguous words from an array of any rank to be passed(avoiding f90 rank conformance check)!caller is responsible for completion checks(mpp_sync_self) before and after integer, intent(in) ::put_len, to_pe, get_len, from_pe MPP_TYPE_, intent(in) ::put_data(*) MPP_TYPE_, intent(out) ::get_data(*) logical, intent(in), optional ::block integer, intent(in), optional ::tag integer, intent(out), optional ::recv_request, send_request logical ::block_comm integer ::i MPP_TYPE_, allocatable, save ::local_data(:) !local copy used by non-parallel code(no SHMEM or MPI) integer ::comm_tag integer ::rsize if(.NOT.module_is_initialized) call mpp_error(FATAL, 'MPP_TRANSMIT:You must first call mpp_init.') if(to_pe.EQ.NULL_PE .AND. from_pe.EQ.NULL_PE) return block_comm=.true. if(PRESENT(block)) block_comm=block if(debug) then call SYSTEM_CLOCK(tick) write(stdout_unit,'(a, i18, a, i6, a, 2i6, 2i8)')&'T=', tick, ' PE=', pe, ' MPP_TRANSMIT begin:to_pe, from_pe, put_len, get_len=', to_pe, from_pe, put_len, get_len end if comm_tag=DEFAULT_TAG if(present(tag)) comm_tag=tag!do put first and then get if(to_pe.GE.0 .AND. to_pe.LT.npes) then!use non-blocking sends if(debug .and.(current_clock.NE.0)) call SYSTEM_CLOCK(start_tick)!z1l:truly non-blocking send.! if(request(to_pe).NE.MPI_REQUEST_NULL) then !only one message from pe-> to_pe in queue *PE waiting for to_pe ! call error else get_len so only do gets but you cannot have a pure get with MPI call a get means do a wait to ensure put on remote PE is complete error call increase mpp_nml request_multiply call MPP_TRANSMIT get_len end if return end subroutine MPP_TRANSMIT_ ! MPP_BROADCAST ! subroutine but that doesn t allow !broadcast to a subset of PEs This version and mpp_transmit will remain !backward compatible intent(inout) a
#define min(a, b)
Definition: mosaic_util.h:32