11 #ifndef OOPS_ASSIMILATION_GMRESR_H_
12 #define OOPS_ASSIMILATION_GMRESR_H_
18 #include "oops/util/dot_product.h"
19 #include "oops/util/formats.h"
20 #include "oops/util/Logger.h"
63 template <
typename VECTOR,
typename AMATRIX,
typename PMATRIX>
64 double GMRESR(VECTOR & xx,
const VECTOR & bb,
65 const AMATRIX & A,
const PMATRIX & precond,
66 const int maxiter,
const double tolerance) {
67 const double stagthresh = 1.0e-3;
68 const double smallres = 1.0e-6;
70 std::vector<VECTOR> cvecs;
71 std::vector<VECTOR> zvecs;
79 double normReduction = 1.0;
80 double rrnorm = sqrt(dot_product(rr, rr));
81 double cdotr = rrnorm;
82 const double rrnorm0 = rrnorm;
84 if (rrnorm > smallres) {
85 Log::info() << std::endl;
86 for (
int jiter = 0; jiter < maxiter; ++jiter) {
87 Log::info() <<
" GMRESR Starting Iteration " << jiter+1 << std::endl;
90 if (std::abs(cdotr) < stagthresh*rrnorm) {
91 Log::info() <<
"GMRESR stagnated. Doing an LSQR step." << std::endl;
94 precond.multiply(rr, zz);
99 for (
int jj = 0; jj < jiter; ++jj) {
100 double alpha = -dot_product(cvecs[jj], cc);
101 cc.axpy(alpha, cvecs[jj]);
102 zz.axpy(alpha, zvecs[jj]);
105 double ccnorm = sqrt(dot_product(cc, cc));
107 cvecs[jiter] *= 1.0/ccnorm;
109 zvecs[jiter] *= 1.0/ccnorm;
111 cdotr = dot_product(cvecs[jiter], rr);
112 rr.axpy(-cdotr, cvecs[jiter]);
113 xx.axpy(cdotr, zvecs[jiter]);
115 rrnorm = sqrt(dot_product(rr, rr));
116 normReduction = rrnorm/rrnorm0;
117 Log::info() <<
"GMRESR end of iteration " << jiter+1 << std::endl;
120 if (normReduction < tolerance) {
121 Log::info() <<
"GMRESR: Achieved required reduction in residual norm." << std::endl;
127 return normReduction;
The namespace for the main oops code.
double GMRESR(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)