11 #ifndef OOPS_ASSIMILATION_MINRES_H_
12 #define OOPS_ASSIMILATION_MINRES_H_
18 #include "oops/util/dot_product.h"
19 #include "oops/util/formats.h"
20 #include "oops/util/Logger.h"
65 template <
typename VECTOR,
typename AMATRIX,
typename PMATRIX>
66 double MINRES(VECTOR & x,
const VECTOR & b,
67 const AMATRIX & A,
const PMATRIX & precond,
68 const int maxiter,
const double tolerance) {
74 double xnrm2 = dot_product(x, x);
86 precond.multiply(r, y);
88 double normReduction = 1.0;
89 double ynrm2 = sqrt(dot_product(y, y));
90 double beta = sqrt(dot_product(y, r));
100 Log::info() << std::endl;
101 for (
int jiter = 0; jiter < maxiter; ++jiter) {
102 Log::info() <<
" MINRES Starting Iteration " << jiter+1 << std::endl;
109 y.axpy(-beta/oldb, r1);
112 double alpha = dot_product(y, v);
113 y.axpy(-alpha/beta, r2);
117 precond.multiply(r2, y);
120 beta = sqrt(dot_product(y, r2));
123 double oldeps = epsln;
124 double delta = cs*dbar + sn*alpha;
125 double gbar = sn*dbar - cs*alpha;
130 double gamma = sqrt(gbar*gbar + beta*beta);
134 double phi = cs*phibar;
138 double denom = 1/
gamma;
143 work.axpy(-oldeps*denom, work1);
144 work.axpy(-delta*denom, work2);
147 normReduction = phibar/ynrm2;
149 Log::info() <<
"MINRES end of iteration " << jiter+1 <<
". PNorm reduction= "
150 << util::full_precision(normReduction) << std::endl << std::endl;
152 if (normReduction <= tolerance) {
153 Log::info() <<
"MINRES: Achieved required reduction in residual norm." << std::endl;
159 Log::info() <<
"MINRES: end" << std::endl;
161 return normReduction;
166 #endif // OOPS_ASSIMILATION_MINRES_H_