11 #ifndef OOPS_ASSIMILATION_SPECTRALLMP_H_
12 #define OOPS_ASSIMILATION_SPECTRALLMP_H_
21 #include "eckit/config/LocalConfiguration.h"
23 #include "oops/util/dot_product.h"
24 #include "oops/util/Logger.h"
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> &);
67 std::vector<std::unique_ptr<VECTOR>>
X_;
68 std::vector<std::unique_ptr<VECTOR>>
U_;
73 std::vector<std::unique_ptr<VECTOR>>
Y_;
74 std::vector<std::unique_ptr<VECTOR>>
S_;
83 template<
typename VECTOR>
85 : maxpairs_(0), useoldpairs_(false), RitzPrecond_(false), maxouter_(0), update_(0)
88 Log::info() <<
"SpectralLMP: maxouter = " <<
maxouter_ << std::endl;
90 if (conf.has(
"preconditioner")) {
91 const eckit::LocalConfiguration precond(conf,
"preconditioner");
92 if (precond.has(
"maxpairs"))
maxpairs_ = precond.getInt(
"maxpairs");
95 Log::info() <<
"SpectralLMP: maxpairs_ = " <<
maxpairs_ << std::endl;
98 if (precond.has(
"ritz")) ritz = precond.getString(
"ritz");
100 Log::info() <<
"SpectralLMP: RitzPrecond_ = " <<
RitzPrecond_ << std::endl;
103 if (precond.has(
"useoldpairs")) old = precond.getString(
"useoldpairs");
105 Log::info() <<
"SpectralLMP: useoldpairs_ = " <<
useoldpairs_ << std::endl;
112 template<
typename VECTOR>
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) {
124 usedpairIndx_.clear();
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);
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;
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];
155 if (std::abs(erritz) <= 0.001*lambda) {
156 convIndx.push_back(jiter);
157 eigvals_.push_back(lambda);
158 omega_.push_back(erritz/lambda);
161 if (convIndx.size() > maxpairs_) {
162 convIndx.erase(convIndx.begin());
164 if (eigvals_.size() > maxpairs_) {
165 eigvals_.erase(eigvals_.begin());
166 omega_.erase(omega_.begin());
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;
179 reduc -= usedpairIndx_[kiter];
180 usedpairIndx_[kiter] = 0;
183 unsigned counter = 0;
184 for (
unsigned kiter = 0; kiter < usedpairIndx_.size(); ++kiter) {
185 counter += usedpairIndx_[kiter];
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;
200 if (X_.size() != 0) {
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));
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]);
213 X_.emplace_back(std::move(ww));
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]);
221 U_.emplace_back(std::move(ww));
225 Zlast_.push_back(*Zl[nvec-1]);
226 Zhlast_.push_back(*Zhl[nvec-1]);
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]);
239 kk += usedpairIndx_[kiter];
240 Y_.emplace_back(std::move(ww));
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]);
253 kk += usedpairIndx_[kiter];
254 S_.emplace_back(std::move(ww));
270 template<
typename VECTOR>
274 for (
unsigned iiter = 0; iiter < eigvals_.size(); ++iiter) {
275 double zeval = std::min(10.0, eigvals_[iiter]);
277 double zz = 1.0/zeval - 1.0;
278 b.axpy(zz*dot_product(
a, *X_[iiter]), *U_[iiter]);
281 if (RitzPrecond_ && eigvals_.size() != 0) {
283 std::vector<double> sta;
285 for (
unsigned iiter = 0; iiter < S_.size(); ++iiter) {
286 sta.push_back(dot_product(*S_[iiter],
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]);
298 for (
unsigned kiter = 0; kiter < S_.size(); ++kiter) {
299 b.axpy(-sta[kiter], Zhlast_[zcount[kiter]]);
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)