11 #ifndef OOPS_ASSIMILATION_MINRES_H_
12 #define OOPS_ASSIMILATION_MINRES_H_
19 #include "oops/util/dot_product.h"
20 #include "oops/util/formats.h"
21 #include "oops/util/Logger.h"
66 template <
typename VECTOR,
typename AMATRIX,
typename PMATRIX>
67 double MINRES(VECTOR & x,
const VECTOR & b,
68 const AMATRIX & A,
const PMATRIX & precond,
69 const int maxiter,
const double tolerance) {
75 double xnrm2 = dot_product(x, x);
87 precond.multiply(r, y);
89 double normReduction = 1.0;
90 double ynrm2 = sqrt(dot_product(y, y));
91 double beta = sqrt(dot_product(y, r));
101 Log::info() << std::endl;
102 for (
int jiter = 0; jiter < maxiter; ++jiter) {
103 Log::info() <<
" MINRES Starting Iteration " << jiter+1 << std::endl;
110 y.axpy(-beta/oldb, r1);
113 double alpha = dot_product(y, v);
114 y.axpy(-alpha/beta, r2);
118 precond.multiply(r2, y);
121 beta = sqrt(dot_product(y, r2));
124 double oldeps = epsln;
125 double delta = cs*dbar + sn*alpha;
126 double gbar = sn*dbar - cs*alpha;
131 double gamma = sqrt(gbar*gbar + beta*beta);
135 double phi = cs*phibar;
139 double denom = 1/
gamma;
144 work.axpy(-oldeps*denom, work1);
145 work.axpy(-delta*denom, work2);
148 normReduction = phibar/ynrm2;
150 Log::info() <<
"MINRES end of iteration " << jiter+1 << std::endl;
153 if (normReduction <= tolerance) {
154 Log::info() <<
"MINRES: Achieved required reduction in residual norm." << std::endl;
160 Log::info() <<
"MINRES: end" << std::endl;
162 return normReduction;
The namespace for the main oops code.
double MINRES(VECTOR &x, const VECTOR &b, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance)
void printNormReduction(int iteration, const double &grad, const double &norm)