11 #ifndef OOPS_ASSIMILATION_FULLGMRES_H_
12 #define OOPS_ASSIMILATION_FULLGMRES_H_
21 #include "oops/util/dot_product.h"
22 #include "oops/util/formats.h"
23 #include "oops/util/Logger.h"
67 template <
typename VECTOR,
typename AMATRIX,
typename PMATRIX>
68 double FullGMRES(VECTOR & xx,
const VECTOR & bb,
const AMATRIX & A,
69 const PMATRIX & precond,
const int maxiter,
70 const double tolerance, std::vector<VECTOR> & pqVEC,
71 std::vector<VECTOR> & xyVEC) {
72 std::vector<VECTOR> VV;
74 VV.reserve(maxiter+1);
76 std::vector< std::vector<double> > H;
77 std::vector<double> cs(maxiter+1, 0);
78 std::vector<double> sn(maxiter+1, 0);
79 std::vector<double> ss;
80 std::vector<double> yy(maxiter+1, 0);
86 double xnrm2 = dot_product(xx, xx);
92 precond.multiply(rr, zz);
94 double znrm2 = sqrt(dot_product(zz, zz));
95 double normReduction = 1.0;
103 for (
int ii = 0; ii <= maxiter-1; ii++) {
104 H[ii].resize(maxiter + 1);
105 for (
int jj = 0; jj <= maxiter; jj++) {
115 Log::info() << std::endl;
116 for (jiter = 0; jiter < maxiter; ++jiter) {
117 Log::info() <<
" FullGMRES Starting Iteration " << jiter+1 << std::endl;
119 A.multiply(VV[jiter], zz);
120 precond.multiply(zz, ww);
122 xyVEC.push_back(VV[jiter]);
126 double avnrm2 = sqrt(dot_product(ww, ww));
129 for (
int jj = 0; jj <= jiter; ++jj) {
130 H[jiter][jj] = dot_product(VV[jj], ww);
131 ww.axpy(-H[jiter][jj], VV[jj]);
134 H[jiter][jiter+1] = sqrt(dot_product(ww, ww));
136 double av2nrm2 = H[jiter][jiter+1];
139 if (avnrm2 + 0.001*av2nrm2 == avnrm2) {
140 for (
int jj = 0; jj <= jiter; ++jj) {
141 double hr = dot_product(VV[jj], ww);
143 ww.axpy(-hr, VV[jj]);
145 H[jiter][jiter+1] = sqrt(dot_product(ww, ww));
149 if (H[jiter][jiter+1] != 0.0) {
150 ww *= 1/H[jiter][jiter+1];
157 for (
int jj = 0; jj < jiter; ++jj) {
158 double temp = cs[jj]*H[jiter][jj] + sn[jj]*H[jiter][jj+1];
159 H[jiter][jj+1] = -sn[jj]*H[jiter][jj] + cs[jj]*H[jiter][jj+1];
165 rotmat(H[jiter][jiter], H[jiter][jiter+1], cs[jiter], sn[jiter]);
167 H[jiter][jiter] = cs[jiter]*H[jiter][jiter] + sn[jiter]*H[jiter][jiter+1];
168 H[jiter][jiter+1] = 0.0;
170 double temp = cs[jiter]*ss[jiter];
171 ss.push_back(-sn[jiter]*ss[jiter]);
174 normReduction = std::abs(ss[jiter+1])/znrm2;
175 Log::info() <<
"FullGMRES end of iteration " << jiter+1 << std::endl;
178 if (normReduction <= tolerance) {
179 Log::info() <<
"FullGMRES: Achieved required reduction in presidual norm." << std::endl;
187 for (
int jj = 0; jj < jiter; ++jj) {
188 xx.axpy(yy[jj], VV[jj]);
191 return normReduction;
Solves the upper triangular system sol = H \ rhs.
The namespace for the main oops code.
double FullGMRES(VECTOR &xx, const VECTOR &bb, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance, std::vector< VECTOR > &pqVEC, std::vector< VECTOR > &xyVEC)
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)
Compute the Givens rotation matrix parameters.