11 #ifndef OOPS_ASSIMILATION_FGMRES_H_
12 #define OOPS_ASSIMILATION_FGMRES_H_
21 #include "oops/util/dot_product.h"
22 #include "oops/util/formats.h"
23 #include "oops/util/Logger.h"
68 template <
typename VECTOR,
typename AMATRIX,
typename PMATRIX>
69 double FGMRES(VECTOR & x,
const VECTOR & b,
70 const AMATRIX & A,
const PMATRIX & precond,
71 const int maxiter,
const double tolerance) {
72 std::vector<VECTOR> V;
73 std::vector<VECTOR> Z;
74 std::vector< std::vector<double> > H;
75 std::vector<double> cs(maxiter+1, 0);
76 std::vector<double> sn(maxiter+1, 0);
77 std::vector<double> s;
78 std::vector<double> y(maxiter+1, 0);
88 double xnrm2 = dot_product(x, x);
94 double bnrm2 = sqrt(dot_product(b, b));
95 double rnrm2 = sqrt(dot_product(r, r));
96 double normReduction = rnrm2 / bnrm2;
100 for (
int ii = 0; ii <= maxiter-1; ii++) {
101 H[ii].resize(maxiter + 1);
102 for (
int jj = 0; jj <= maxiter; jj++) {
113 Log::info() << std::endl;
114 for (jiter = 0; jiter < maxiter; ++jiter) {
115 Log::info() <<
" FGMRES Starting Iteration " << jiter+1 << std::endl;
117 precond.multiply(V[jiter], z);
121 double avnrm2 = sqrt(dot_product(w, w));
124 for (
int jj = 0; jj <= jiter; ++jj) {
125 H[jiter][jj] = dot_product(V[jj], w);
126 w.axpy(-H[jiter][jj], V[jj]);
128 H[jiter][jiter+1] = sqrt(dot_product(w, w));
129 double av2nrm2 = H[jiter][jiter+1];
132 if (avnrm2 + 0.001*av2nrm2 == avnrm2) {
133 for (
int jj = 0; jj <= jiter; ++jj) {
134 double hr = dot_product(V[jj], w);
138 H[jiter][jiter+1] = sqrt(dot_product(w, w));
142 if (H[jiter][jiter+1] != 0) {
143 w *= 1/H[jiter][jiter+1];
149 for (
int jj = 0; jj < jiter; ++jj) {
150 double temp = cs[jj]*H[jiter][jj] + sn[jj]*H[jiter][jj+1];
151 H[jiter][jj+1] = -sn[jj]*H[jiter][jj] + cs[jj]*H[jiter][jj+1];
157 rotmat(H[jiter][jiter], H[jiter][jiter+1], cs[jiter], sn[jiter]);
159 H[jiter][jiter] = cs[jiter]*H[jiter][jiter] + sn[jiter]*H[jiter][jiter+1];
160 H[jiter][jiter+1] = 0.0;
162 double temp = cs[jiter]*s[jiter];
163 s.push_back(-sn[jiter]*s[jiter]);
169 for (
int jj = 0; jj < jiter+1; ++jj) {
170 x.axpy(y[jj], Z[jj]);
173 normReduction = std::abs(s[jiter+1])/bnrm2;
174 Log::info() <<
"FGMRES end of iteration " << jiter+1 << std::endl;
177 if (normReduction <= tolerance) {
178 Log::info() <<
"FGMRES: Achieved required reduction in residual norm." << std::endl;
187 for (
int jj = 0; jj < jiter; ++jj) {
188 x.axpy(y[jj], Z[jj]);
191 Log::info() <<
"FGMRES: end" << std::endl;
193 return normReduction;
Solves the upper triangular system sol = H \ rhs.
The namespace for the main oops code.
void UpTriSolve(const std::vector< std::vector< double > > &H, const std::vector< double > &rhs, std::vector< double > &sol, const int &dim)
void rotmat(const double &a, const double &b, double &c, double &s)
void printNormReduction(int iteration, const double &grad, const double &norm)
double FGMRES(VECTOR &x, const VECTOR &b, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance)
Compute the Givens rotation matrix parameters.