HalleyChebyAD
 All Classes Files Functions Variables Friends
tminres.hpp
Go to the documentation of this file.
1 // tminres is free software; you can redistribute it and/or modify it under the
2 // terms of the GNU Lesser General Public License (as published by the Free
3 // Software Foundation) version 2.1 dated February 1999.
4 //
5 // Authors:
6 // - Umberto Villa, Emory University - uvilla@emory.edu
7 // - Michael Saunders, Stanford University
8 // - Santiago Akle, Stanford University
9 
95 #include <limits>
96 #include <memory>
97 #include <iostream>
98 #include <iomanip>
99 #include <vector>
100 #include <cmath>
101 #include <algorithm>
102 #include <limits>
103 
104 template< typename Operator, typename Vector, typename Preconditioner>
105 int
106 MINRES(const Operator &A, Vector &x, const Vector &b,
107  const Preconditioner * M, double shift, int &max_iter, double &tol, bool show)
108 {
109 
110  typedef Operator operator_Type;
111  typedef Vector vector_Type;
112 
113  double eps(std::numeric_limits<double>::epsilon());
114 
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 ";
127 
128  if(show)
129  {
130  //std::cout<<std::setfill('-')<<std::setw(80)<<"-"<< "\n";
131  //std::cout<<"| Adapted from minres.m, Stanford University, 06 Jul 2009 |\n";
132  std::cout<<"| Solution of symmetric Ax=b or (A-shift*I)x = b |\n";
133  //std::cout<<std::setfill('-')<<std::setw(80)<<"-"<<"\n";
134  std::cout<<std::setfill(' ');
135  std::cout<<"shift = "<< shift << "; tol = " << tol << "; max iter = " << max_iter<<"\n";
136  }
137 
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);
141  bool done(false);
142 
143  // Step 1
144  /*
145  * Set up y and v for the first Lanczos vector v1.
146  * y = beta1 P' v1, where P = C^(-1).
147  * v is really P'v1
148  */
149 
150  const std::auto_ptr<Vector> r1(x.Clone());
151  const std::auto_ptr<Vector> y(x.Clone());
152 
153  A.Apply(x, *r1);
154  add(*r1, -shift, x, *r1);
155  subtract(b, *r1, *r1);
156 
157  if(M)
158  M->Apply(*r1, *y);
159  else
160  *y = *r1;
161 
162  double beta1(0.0);
163  beta1 = InnerProduct(*r1, *y);
164 
165  // Test for an indefined preconditioner
166  // If b = 0 exactly stop with x = x0.
167 
168  if(beta1 < 0.0)
169  {
170  istop = 9;
171  show = true;
172  done = true;
173  }
174  else
175  {
176  if(beta1 == 0.0)
177  {
178  show = true;
179  done = true;
180  }
181  else
182  beta1 = std::sqrt(beta1); // Normalize y to get v1 later
183  }
184 
185  // TODO: port symmetry checks for A and M
186 
187  // STEP 2
188  /* Initialize other quantities */
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);
196  double z(0.0);
197 
198  const std::auto_ptr<Vector> w(x.Clone()), w1(x.Clone()), w2(x.Clone()), r2(x.Clone());
199  ( *w ) = (* w1 ) = ( *w2 ) = 0.0;
200  *r2 = *r1;
201  const std::auto_ptr<Vector> v(x.Clone());
202 
203  if(show)
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";
210 
211  /* Main Iteration */
212  if(!done)
213  {
214  for(itn = 0; itn < max_iter; ++itn)
215  {
216  // STEP 3
217  /*
218  -----------------------------------------------------------------
219  Obtain quantities for the next Lanczos vector vk+1, k = 1, 2,...
220  The general iteration is similar to the case k = 1 with v0 = 0:
221 
222  p1 = Operator * v1 - beta1 * v0,
223  alpha1 = v1'p1,
224  q2 = p2 - alpha1 * v1,
225  beta2^2 = q2'q2,
226  v2 = (1/beta2) q2.
227 
228  Again, y = betak P vk, where P = C**(-1).
229  .... more description needed.
230  -----------------------------------------------------------------
231  */
232  double s(1./beta); //Normalize previous vector (in y)
233  (*v) = (*y);
234  v->Scale(s); // v = vk if P = I
235 
236  A.Apply(*v,*y);
237  add(*y, -shift, *v, *y);
238  if(itn)
239  add(*y, -beta/oldb, *r1, *y);
240 
241  alpha = InnerProduct(*v, *y); // alphak
242  add(*y, -alpha/beta, *r2, *y);
243  (*r1) = (*r2);
244  (*r2) = (*y);
245  if(M)
246  M->Apply(*r2, *y);
247  else
248  *y = *r2;
249 
250  oldb = beta; //oldb = betak
251  beta = InnerProduct(*r2, *y);
252 
253  if(beta < 0)
254  {
255  istop = 9;
256  break;
257  }
258 
259  beta = sqrt(beta);
260  tnorm2 += alpha*alpha + oldb*oldb + beta*beta;
261 
262  if(itn == 0) //Initialize a few things
263  {
264  if(beta/beta1 < 10.0*eps)
265  istop = 10;
266  }
267 
268  // Apply previous rotation Q_{k-1} to get
269  // [delta_k epsln_{k+1}] = [cs sn] [dbar_k 0]
270  // [gbar_k dbar_{k+1}] [sn -cs] [alpha_k beta_{k+1}].
271  oldeps = epsln;
272  delta = cs*dbar + sn*alpha;
273  gbar = sn*dbar - cs*alpha;
274  epsln = sn*beta;
275  dbar = - cs*beta;
276  double root(sqrt(gbar*gbar + dbar*dbar));
277  Arnorm = phibar * root; // ||Ar_{k-1}||
278 
279  // Compute next plane rotation Q_k
280  gamma = sqrt(gbar*gbar + beta*beta); // gamma_k
281  gamma = std::max(gamma, eps);
282  cs = gbar/gamma; // c_k
283  sn = beta/gamma; // s_k
284  phi = cs*phibar; // phi_k
285  phibar = sn*phibar; // phibar_{k+1}
286 
287 
288  // Update x
289  double denom(1./gamma);
290  (*w1) = (*w2);
291  (*w2) = (*w);
292  add(-oldeps, *w1, -delta, *w2, *w);
293  add(denom, *v, *w, *w);
294  add(x, phi, *w, x);
295 
296  // go round again
297  gmax = std::max(gmax, gamma);
298  gmin = std::min(gmin, gamma);
299  z = rhs1/gamma;
300  rhs1 = rhs2 - delta*z;
301  rhs2 = - epsln*z;
302 
303  // Estimate various norms
304 
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);
311  double diag(gbar);
312  if(0 == diag)
313  diag = epsa;
314 
315  qrnorm = phibar;
316  rnorm = qrnorm;
317  double test1(0.0), test2(0.0);
318  test1 = rnorm/(Anorm*ynorm); // ||r||/(||A|| ||x||)
319  test2 = root/ Anorm; // ||A r_{k-1}|| / (||A|| ||r_{k-1}||)
320 
321  // Estimate cond(A)
322  /*
323  In this version we look at the diagonals of R in the
324  factorization of the lower Hessenberg matrix, Q * H = R,
325  where H is the tridiagonal matrix from Lanczos with one
326  extra row, beta(k+1) e_k^T.
327  */
328  Acond = gmax/gmin;
329 
330  //See if any of the stopping criteria is satisfied
331  if(0 == istop)
332  {
333  double t1(1.0+test1), t2(1.0+test2); //This test work if tol < eps
334  if(t2 <= 1.) istop = 2;
335  if(t1 <= 1.) istop = 1;
336 
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;
342  }
343 
344  if(show)
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;
351 
352  if(0 != istop)
353  break;
354 
355  }
356  }/*End of Main Iteration */
357 
358  // Display final status
359  if(show)
360  {
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(' ');
369  }
370 
371  return istop;
372 }
373 
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