104 template<
typename Operator,
typename Vector,
typename Preconditioner>
106 MINRES(
const Operator &A, Vector &x,
const Vector &b,
107 const Preconditioner * M,
double shift,
int &max_iter,
double &tol,
bool show)
110 typedef Operator operator_Type;
111 typedef Vector vector_Type;
113 double eps(std::numeric_limits<double>::epsilon());
115 std::vector<std::string> msg(11);
116 msg[0] =
" beta1 = 0. The exact solution is x = 0 ";
117 msg[1] =
" A solution to Ax = b was found, given tol ";
118 msg[2] =
" A least-squares solution was found, given tol ";
119 msg[3] =
" Reasonable accuracy achieved, given eps ";
120 msg[4] =
" x has converged to an eigenvector ";
121 msg[5] =
" acond has exceeded 0.1/eps ";
122 msg[6] =
" The iteration limit was reached ";
123 msg[7] =
" A does not define a symmetric matrix ";
124 msg[8] =
" M does not define a symmetric matrix ";
125 msg[9] =
" M does not define a pos-def preconditioner ";
126 msg[10] =
" beta2 = 0. If M = I, b and x are eigenvectors ";
132 std::cout<<
"| Solution of symmetric Ax=b or (A-shift*I)x = b |\n";
134 std::cout<<std::setfill(
' ');
135 std::cout<<
"shift = "<< shift <<
"; tol = " << tol <<
"; max iter = " << max_iter<<
"\n";
138 int istop(0), itn(0);
139 double Anorm(0.0), Acond(0.0), Arnorm(0.0);
140 double rnorm(0.0), ynorm(0.0);
150 const std::auto_ptr<Vector> r1(x.Clone());
151 const std::auto_ptr<Vector> y(x.Clone());
154 add(*r1, -shift, x, *r1);
155 subtract(b, *r1, *r1);
163 beta1 = InnerProduct(*r1, *y);
182 beta1 = std::sqrt(beta1);
189 double oldb(0.0), beta(beta1), dbar(0.0), epsln(0.0), oldeps(0.0);
190 double qrnorm(beta1), phi(0.0), phibar(beta1), rhs1(beta1);
191 double rhs2(0.0), tnorm2(0.0), ynorm2(0.0);
192 double cs(-1.0), sn(0.0);
193 double gmax(0.0), gmin(std::numeric_limits<double>::max( ));
194 double alpha(0.0), gamma(0.0);
195 double delta(0.0), gbar(0.0);
198 const std::auto_ptr<Vector> w(x.Clone()), w1(x.Clone()), w2(x.Clone()), r2(x.Clone());
199 ( *w ) = (* w1 ) = ( *w2 ) = 0.0;
201 const std::auto_ptr<Vector> v(x.Clone());
204 std::cout<<std::setw(6)<<
"Itn"
205 << std::setw(14) <<
"Compatible"
206 << std::setw(14) <<
"LS"
207 << std::setw(14) <<
"norm(A)"
208 << std::setw(14) <<
"cond(A)"
209 << std::setw(14) <<
"gbar/|A|"<<
"\n";
214 for(itn = 0; itn < max_iter; ++itn)
237 add(*y, -shift, *v, *y);
239 add(*y, -beta/oldb, *r1, *y);
241 alpha = InnerProduct(*v, *y);
242 add(*y, -alpha/beta, *r2, *y);
251 beta = InnerProduct(*r2, *y);
260 tnorm2 += alpha*alpha + oldb*oldb + beta*beta;
264 if(beta/beta1 < 10.0*eps)
272 delta = cs*dbar + sn*alpha;
273 gbar = sn*dbar - cs*alpha;
276 double root(sqrt(gbar*gbar + dbar*dbar));
277 Arnorm = phibar * root;
280 gamma = sqrt(gbar*gbar + beta*beta);
281 gamma = std::max(gamma, eps);
289 double denom(1./gamma);
292 add(-oldeps, *w1, -delta, *w2, *w);
293 add(denom, *v, *w, *w);
297 gmax = std::max(gmax, gamma);
298 gmin = std::min(gmin, gamma);
300 rhs1 = rhs2 - delta*z;
305 Anorm = sqrt(tnorm2);
306 ynorm2 = InnerProduct(x, x);
307 ynorm = sqrt(ynorm2);
308 double epsa(Anorm*eps);
309 double epsx(epsa*ynorm);
310 double epsr(Anorm*ynorm*tol);
317 double test1(0.0), test2(0.0);
318 test1 = rnorm/(Anorm*ynorm);
333 double t1(1.0+test1), t2(1.0+test2);
334 if(t2 <= 1.) istop = 2;
335 if(t1 <= 1.) istop = 1;
337 if(itn >= max_iter-1) istop = 6;
338 if(Acond >= .1/eps) istop = 4;
339 if(epsx >= beta1) istop = 3;
340 if(test2 <= tol ) istop = 2;
341 if(test1 <= tol) istop = 1;
345 std::cout<< std::setw(6) << itn
346 << std::setw(14) << test1
347 << std::setw(14) << test2
348 << std::setw(14) << Anorm
349 << std::setw(14) << Acond
350 << std::setw(14) << gbar/Anorm << std::endl;
361 std::cout << std::setfill(
'-') << std::setw(80) <<
"-" <<
"\n";
362 std::cout << msg[istop] <<
"\n";
363 std::cout <<
" Number of iterations: " << itn <<
"\n";
364 std::cout <<
" Anorm = " << Anorm <<
"\t Acond = " << Acond <<
"\n";
365 std::cout <<
" rnorm = " << rnorm <<
"\t ynorm = " << ynorm <<
"\n";
366 std::cout <<
" Arnorm = " << Arnorm <<
"\n";
367 std::cout << std::setfill(
'-') << std::setw(80) <<
"-" << std::endl;
368 std::cout << std::setfill(
' ');
int MINRES(const Operator &A, Vector &x, const Vector &b, const Preconditioner *M, double shift, int &max_iter, double &tol, bool show)
Preconditioner Minres algorithm.
Definition: tminres.hpp:106