11 #ifndef OOPS_ASSIMILATION_PLANCZOS_H_
12 #define OOPS_ASSIMILATION_PLANCZOS_H_
19 #include "oops/util/dot_product.h"
20 #include "oops/util/formats.h"
21 #include "oops/util/Logger.h"
62 template <
typename VECTOR,
typename AMATRIX,
typename PMATRIX>
63 double PLanczos(VECTOR & xx,
const VECTOR & bb,
64 const AMATRIX & A,
const PMATRIX & precond,
65 const int maxiter,
const double tolerance) {
69 std::vector<VECTOR> vVEC;
70 std::vector<VECTOR> zVEC;
72 vVEC.reserve(maxiter+1);
73 zVEC.reserve(maxiter+1);
75 Log::debug() <<
"PLanczos xx = " << xx << std::endl;
76 Log::debug() <<
"PLanczos bb = " << bb << std::endl;
78 std::vector<double> alphas;
79 std::vector<double> betas;
80 std::vector<double> dd;
81 std::vector<double> yy;
85 double xnrm2 = dot_product(xx, xx);
90 Log::debug() <<
"PLanczos rr = " << rr << std::endl;
93 precond.multiply(rr, zz);
94 Log::debug() <<
"PLanczos zz = " << zz << std::endl;
96 double normReduction = 1.0;
97 double beta0 = sqrt(dot_product(rr, zz));
99 Log::debug() <<
"PLanczos beta0 = " << beta0 << std::endl;
108 Log::debug() <<
"PLanczos vv = " << vv << std::endl;
109 Log::debug() <<
"PLanczos zz = " << zz << std::endl;
112 Log::info() << std::endl;
113 while (jiter < maxiter) {
114 Log::info() <<
" PLanczos Starting Iteration " << jiter+1 << std::endl;
118 Log::debug() <<
"PLanczos ww = " << ww << std::endl;
119 if (jiter > 0) ww.axpy(-beta, vVEC[jiter-1]);
121 double alpha = dot_product(zz, ww);
122 Log::debug() <<
"PLanczos alpha = " << alpha << std::endl;
125 Log::debug() <<
"PLanczos ww = " << ww << std::endl;
128 for (
int iiter = 0; iiter < jiter; ++iiter) {
129 double proj = dot_product(ww, zVEC[iiter]);
130 ww.axpy(-proj, vVEC[iiter]);
132 Log::debug() <<
"PLanczos ww = " << ww << std::endl;
134 precond.multiply(ww, zz);
135 Log::debug() <<
"PLanczos zz = " << zz << std::endl;
137 beta = sqrt(dot_product(zz, ww));
138 Log::debug() <<
"PLanczos beta = " << beta << std::endl;
143 Log::debug() <<
"PLanczos vv = " << vv << std::endl;
144 Log::debug() <<
"PLanczos zz = " << zz << std::endl;
149 alphas.push_back(alpha);
152 yy.push_back(beta0/alpha);
156 dd.push_back(beta0*dot_product(zVEC[0], vv));
161 double rznorm = beta*std::abs(yy[jiter]);
162 Log::debug() <<
"PLanczos rznorm = " << rznorm << std::endl;
164 normReduction = rznorm/beta0;
166 betas.push_back(beta);
168 Log::info() <<
"PLanczos end of iteration " << jiter+1 << std::endl;
173 if (normReduction < tolerance) {
174 Log::info() <<
"PLanczos: Achieved required reduction in residual norm." << std::endl;
178 Log::debug() <<
"PLanczos jiter = " << jiter << std::endl;
181 for (
int iiter = 0; iiter < jiter; ++iiter) {
182 xx.axpy(yy[iiter], zVEC[iiter]);
185 Log::info() <<
"PLanczos: end" << std::endl;
187 return normReduction;
The namespace for the main oops code.
void TriDiagSolve(const std::vector< double > &diag, const std::vector< double > &sub, const std::vector< double > &rhs, std::vector< double > &sol)
double PLanczos(VECTOR &xx, const VECTOR &bb, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance)
void printNormReduction(int iteration, const double &grad, const double &norm)