dune-istl  2.3.1
solvers.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 
4 #ifndef DUNE_SOLVERS_HH
5 #define DUNE_SOLVERS_HH
6 
7 #include <cmath>
8 #include <complex>
9 #include <iostream>
10 #include <iomanip>
11 #include <string>
12 #include <vector>
13 
14 #include "istlexception.hh"
15 #include "operators.hh"
16 #include "scalarproducts.hh"
17 #include "solver.hh"
18 #include "preconditioner.hh"
19 #include <dune/common/timer.hh>
20 #include <dune/common/ftraits.hh>
21 #include <dune/common/shared_ptr.hh>
22 #include <dune/common/static_assert.hh>
23 
24 namespace Dune {
42  //=====================================================================
43  // Implementation of this interface
44  //=====================================================================
45 
54  template<class X>
55  class LoopSolver : public InverseOperator<X,X> {
56  public:
58  typedef X domain_type;
60  typedef X range_type;
62  typedef typename X::field_type field_type;
63 
83  template<class L, class P>
84  LoopSolver (L& op, P& prec,
85  double reduction, int maxit, int verbose) :
86  ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
87  {
88  dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
89  "L and P have to have the same category!");
90  dune_static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
91  "L has to be sequential!");
92  }
93 
114  template<class L, class S, class P>
115  LoopSolver (L& op, S& sp, P& prec,
116  double reduction, int maxit, int verbose) :
117  _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
118  {
119  dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
120  "L and P must have the same category!");
121  dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category),
122  "L and S must have the same category!");
123  }
124 
125 
127  virtual void apply (X& x, X& b, InverseOperatorResult& res)
128  {
129  // clear solver statistics
130  res.clear();
131 
132  // start a timer
133  Timer watch;
134 
135  // prepare preconditioner
136  _prec.pre(x,b);
137 
138  // overwrite b with defect
139  _op.applyscaleadd(-1,x,b);
140 
141  // compute norm, \todo parallelization
142  double def0 = _sp.norm(b);
143 
144  // printing
145  if (_verbose>0)
146  {
147  std::cout << "=== LoopSolver" << std::endl;
148  if (_verbose>1)
149  {
150  this->printHeader(std::cout);
151  this->printOutput(std::cout,0,def0);
152  }
153  }
154 
155  // allocate correction vector
156  X v(x);
157 
158  // iteration loop
159  int i=1; double def=def0;
160  for ( ; i<=_maxit; i++ )
161  {
162  v = 0; // clear correction
163  _prec.apply(v,b); // apply preconditioner
164  x += v; // update solution
165  _op.applyscaleadd(-1,v,b); // update defect
166  double defnew=_sp.norm(b); // comp defect norm
167  if (_verbose>1) // print
168  this->printOutput(std::cout,i,defnew,def);
169  //std::cout << i << " " << defnew << " " << defnew/def << std::endl;
170  def = defnew; // update norm
171  if (def<def0*_reduction || def<1E-30) // convergence check
172  {
173  res.converged = true;
174  break;
175  }
176  }
177 
178  //correct i which is wrong if convergence was not achieved.
179  i=std::min(_maxit,i);
180 
181  // print
182  if (_verbose==1)
183  this->printOutput(std::cout,i,def);
184 
185  // postprocess preconditioner
186  _prec.post(x);
187 
188  // fill statistics
189  res.iterations = i;
190  res.reduction = def/def0;
191  res.conv_rate = pow(res.reduction,1.0/i);
192  res.elapsed = watch.elapsed();
193 
194  // final print
195  if (_verbose>0)
196  {
197  std::cout << "=== rate=" << res.conv_rate
198  << ", T=" << res.elapsed
199  << ", TIT=" << res.elapsed/i
200  << ", IT=" << i << std::endl;
201  }
202  }
203 
205  virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
206  {
207  std::swap(_reduction,reduction);
208  (*this).apply(x,b,res);
209  std::swap(_reduction,reduction);
210  }
211 
212  private:
214  LinearOperator<X,X>& _op;
215  Preconditioner<X,X>& _prec;
216  ScalarProduct<X>& _sp;
217  double _reduction;
218  int _maxit;
219  int _verbose;
220  };
221 
222 
223  // all these solvers are taken from the SUMO library
225  template<class X>
226  class GradientSolver : public InverseOperator<X,X> {
227  public:
229  typedef X domain_type;
231  typedef X range_type;
233  typedef typename X::field_type field_type;
234 
240  template<class L, class P>
241  GradientSolver (L& op, P& prec,
242  double reduction, int maxit, int verbose) :
243  ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
244  {
245  dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
246  "L and P have to have the same category!");
247  dune_static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
248  "L has to be sequential!");
249  }
255  template<class L, class S, class P>
256  GradientSolver (L& op, S& sp, P& prec,
257  double reduction, int maxit, int verbose) :
258  _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
259  {
260  dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
261  "L and P have to have the same category!");
262  dune_static_assert(static_cast<int>(L::category) == static_cast<int>(S::category),
263  "L and S have to have the same category!");
264  }
265 
271  virtual void apply (X& x, X& b, InverseOperatorResult& res)
272  {
273  res.clear(); // clear solver statistics
274  Timer watch; // start a timer
275  _prec.pre(x,b); // prepare preconditioner
276  _op.applyscaleadd(-1,x,b); // overwrite b with defect
277 
278  X p(x); // create local vectors
279  X q(b);
280 
281  double def0 = _sp.norm(b); // compute norm
282 
283  if (_verbose>0) // printing
284  {
285  std::cout << "=== GradientSolver" << std::endl;
286  if (_verbose>1)
287  {
288  this->printHeader(std::cout);
289  this->printOutput(std::cout,0,def0);
290  }
291  }
292 
293  int i=1; double def=def0; // loop variables
294  field_type lambda;
295  for ( ; i<=_maxit; i++ )
296  {
297  p = 0; // clear correction
298  _prec.apply(p,b); // apply preconditioner
299  _op.apply(p,q); // q=Ap
300  lambda = _sp.dot(p,b)/_sp.dot(q,p); // minimization
301  x.axpy(lambda,p); // update solution
302  b.axpy(-lambda,q); // update defect
303 
304  double defnew=_sp.norm(b); // comp defect norm
305  if (_verbose>1) // print
306  this->printOutput(std::cout,i,defnew,def);
307 
308  def = defnew; // update norm
309  if (def<def0*_reduction || def<1E-30) // convergence check
310  {
311  res.converged = true;
312  break;
313  }
314  }
315 
316  //correct i which is wrong if convergence was not achieved.
317  i=std::min(_maxit,i);
318 
319  if (_verbose==1) // printing for non verbose
320  this->printOutput(std::cout,i,def);
321 
322  _prec.post(x); // postprocess preconditioner
323  res.iterations = i; // fill statistics
324  res.reduction = def/def0;
325  res.conv_rate = pow(res.reduction,1.0/i);
326  res.elapsed = watch.elapsed();
327  if (_verbose>0) // final print
328  std::cout << "=== rate=" << res.conv_rate
329  << ", T=" << res.elapsed
330  << ", TIT=" << res.elapsed/i
331  << ", IT=" << i << std::endl;
332  }
333 
339  virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
340  {
341  std::swap(_reduction,reduction);
342  (*this).apply(x,b,res);
343  std::swap(_reduction,reduction);
344  }
345 
346  private:
348  LinearOperator<X,X>& _op;
349  Preconditioner<X,X>& _prec;
350  ScalarProduct<X>& _sp;
351  double _reduction;
352  int _maxit;
353  int _verbose;
354  };
355 
356 
357 
359  template<class X>
360  class CGSolver : public InverseOperator<X,X> {
361  public:
363  typedef X domain_type;
365  typedef X range_type;
367  typedef typename X::field_type field_type;
368 
374  template<class L, class P>
375  CGSolver (L& op, P& prec, double reduction, int maxit, int verbose) :
376  ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
377  {
378  dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
379  "L and P must have the same category!");
380  dune_static_assert( static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
381  "L must be sequential!");
382  }
388  template<class L, class S, class P>
389  CGSolver (L& op, S& sp, P& prec, double reduction, int maxit, int verbose) :
390  _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
391  {
392  dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
393  "L and P must have the same category!");
394  dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category),
395  "L and S must have the same category!");
396  }
397 
403  virtual void apply (X& x, X& b, InverseOperatorResult& res)
404  {
405  res.clear(); // clear solver statistics
406  Timer watch; // start a timer
407  _prec.pre(x,b); // prepare preconditioner
408  _op.applyscaleadd(-1,x,b); // overwrite b with defect
409 
410  X p(x); // the search direction
411  X q(x); // a temporary vector
412 
413  double def0 = _sp.norm(b); // compute norm
414  if (def0<1E-30) // convergence check
415  {
416  res.converged = true;
417  res.iterations = 0; // fill statistics
418  res.reduction = 0;
419  res.conv_rate = 0;
420  res.elapsed=0;
421  if (_verbose>0) // final print
422  std::cout << "=== rate=" << res.conv_rate
423  << ", T=" << res.elapsed << ", TIT=" << res.elapsed
424  << ", IT=0" << std::endl;
425  return;
426  }
427 
428  if (_verbose>0) // printing
429  {
430  std::cout << "=== CGSolver" << std::endl;
431  if (_verbose>1) {
432  this->printHeader(std::cout);
433  this->printOutput(std::cout,0,def0);
434  }
435  }
436 
437  // some local variables
438  double def=def0; // loop variables
439  field_type rho,rholast,lambda,alpha,beta;
440 
441  // determine initial search direction
442  p = 0; // clear correction
443  _prec.apply(p,b); // apply preconditioner
444  rholast = _sp.dot(p,b); // orthogonalization
445 
446  // the loop
447  int i=1;
448  for ( ; i<=_maxit; i++ )
449  {
450  // minimize in given search direction p
451  _op.apply(p,q); // q=Ap
452  alpha = _sp.dot(p,q); // scalar product
453  lambda = rholast/alpha; // minimization
454  x.axpy(lambda,p); // update solution
455  b.axpy(-lambda,q); // update defect
456 
457  // convergence test
458  double defnew=_sp.norm(b); // comp defect norm
459 
460  if (_verbose>1) // print
461  this->printOutput(std::cout,i,defnew,def);
462 
463  def = defnew; // update norm
464  if (def<def0*_reduction || def<1E-30) // convergence check
465  {
466  res.converged = true;
467  break;
468  }
469 
470  // determine new search direction
471  q = 0; // clear correction
472  _prec.apply(q,b); // apply preconditioner
473  rho = _sp.dot(q,b); // orthogonalization
474  beta = rho/rholast; // scaling factor
475  p *= beta; // scale old search direction
476  p += q; // orthogonalization with correction
477  rholast = rho; // remember rho for recurrence
478  }
479 
480  //correct i which is wrong if convergence was not achieved.
481  i=std::min(_maxit,i);
482 
483  if (_verbose==1) // printing for non verbose
484  this->printOutput(std::cout,i,def);
485 
486  _prec.post(x); // postprocess preconditioner
487  res.iterations = i; // fill statistics
488  res.reduction = def/def0;
489  res.conv_rate = pow(res.reduction,1.0/i);
490  res.elapsed = watch.elapsed();
491 
492  if (_verbose>0) // final print
493  {
494  std::cout << "=== rate=" << res.conv_rate
495  << ", T=" << res.elapsed
496  << ", TIT=" << res.elapsed/i
497  << ", IT=" << i << std::endl;
498  }
499  }
500 
506  virtual void apply (X& x, X& b, double reduction,
508  {
509  std::swap(_reduction,reduction);
510  (*this).apply(x,b,res);
511  std::swap(_reduction,reduction);
512  }
513 
514  private:
516  LinearOperator<X,X>& _op;
517  Preconditioner<X,X>& _prec;
518  ScalarProduct<X>& _sp;
519  double _reduction;
520  int _maxit;
521  int _verbose;
522  };
523 
524 
525  // Ronald Kriemanns BiCG-STAB implementation from Sumo
527  template<class X>
528  class BiCGSTABSolver : public InverseOperator<X,X> {
529  public:
531  typedef X domain_type;
533  typedef X range_type;
535  typedef typename X::field_type field_type;
537  typedef typename FieldTraits<field_type>::real_type real_type;
538 
544  template<class L, class P>
545  BiCGSTABSolver (L& op, P& prec,
546  double reduction, int maxit, int verbose) :
547  ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
548  {
549  dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category), "L and P must be of the same category!");
550  dune_static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential), "L must be sequential!");
551  }
557  template<class L, class S, class P>
558  BiCGSTABSolver (L& op, S& sp, P& prec,
559  double reduction, int maxit, int verbose) :
560  _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
561  {
562  dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
563  "L and P must have the same category!");
564  dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category),
565  "L and S must have the same category!");
566  }
567 
573  virtual void apply (X& x, X& b, InverseOperatorResult& res)
574  {
575  const double EPSILON=1e-80;
576  double it;
577  field_type rho, rho_new, alpha, beta, h, omega;
578  real_type norm, norm_old, norm_0;
579 
580  //
581  // get vectors and matrix
582  //
583  X& r=b;
584  X p(x);
585  X v(x);
586  X t(x);
587  X y(x);
588  X rt(x);
589 
590  //
591  // begin iteration
592  //
593 
594  // r = r - Ax; rt = r
595  res.clear(); // clear solver statistics
596  Timer watch; // start a timer
597  _prec.pre(x,r); // prepare preconditioner
598  _op.applyscaleadd(-1,x,r); // overwrite b with defect
599 
600  rt=r;
601 
602  norm = norm_old = norm_0 = _sp.norm(r);
603 
604  p=0;
605  v=0;
606 
607  rho = 1;
608  alpha = 1;
609  omega = 1;
610 
611  if (_verbose>0) // printing
612  {
613  std::cout << "=== BiCGSTABSolver" << std::endl;
614  if (_verbose>1)
615  {
616  this->printHeader(std::cout);
617  this->printOutput(std::cout,0,norm_0);
618  //std::cout << " Iter Defect Rate" << std::endl;
619  //std::cout << " 0" << std::setw(14) << norm_0 << std::endl;
620  }
621  }
622 
623  if ( norm < (_reduction * norm_0) || norm<1E-30)
624  {
625  res.converged = 1;
626  _prec.post(x); // postprocess preconditioner
627  res.iterations = 0; // fill statistics
628  res.reduction = 0;
629  res.conv_rate = 0;
630  res.elapsed = watch.elapsed();
631  return;
632  }
633 
634  //
635  // iteration
636  //
637 
638  for (it = 0.5; it < _maxit; it+=.5)
639  {
640  //
641  // preprocess, set vecsizes etc.
642  //
643 
644  // rho_new = < rt , r >
645  rho_new = _sp.dot(rt,r);
646 
647  // look if breakdown occured
648  if (std::abs(rho) <= EPSILON)
649  DUNE_THROW(ISTLError,"breakdown in BiCGSTAB - rho "
650  << rho << " <= EPSILON " << EPSILON
651  << " after " << it << " iterations");
652  if (std::abs(omega) <= EPSILON)
653  DUNE_THROW(ISTLError,"breakdown in BiCGSTAB - omega "
654  << omega << " <= EPSILON " << EPSILON
655  << " after " << it << " iterations");
656 
657 
658  if (it<1)
659  p = r;
660  else
661  {
662  beta = ( rho_new / rho ) * ( alpha / omega );
663  p.axpy(-omega,v); // p = r + beta (p - omega*v)
664  p *= beta;
665  p += r;
666  }
667 
668  // y = W^-1 * p
669  y = 0;
670  _prec.apply(y,p); // apply preconditioner
671 
672  // v = A * y
673  _op.apply(y,v);
674 
675  // alpha = rho_new / < rt, v >
676  h = _sp.dot(rt,v);
677 
678  if ( std::abs(h) < EPSILON )
679  DUNE_THROW(ISTLError,"h=0 in BiCGSTAB");
680 
681  alpha = rho_new / h;
682 
683  // apply first correction to x
684  // x <- x + alpha y
685  x.axpy(alpha,y);
686 
687  // r = r - alpha*v
688  r.axpy(-alpha,v);
689 
690  //
691  // test stop criteria
692  //
693 
694  norm = _sp.norm(r);
695 
696  if (_verbose>1) // print
697  {
698  this->printOutput(std::cout,it,norm,norm_old);
699  }
700 
701  if ( norm < (_reduction * norm_0) )
702  {
703  res.converged = 1;
704  break;
705  }
706  it+=.5;
707 
708  norm_old = norm;
709 
710  // y = W^-1 * r
711  y = 0;
712  _prec.apply(y,r);
713 
714  // t = A * y
715  _op.apply(y,t);
716 
717  // omega = < t, r > / < t, t >
718  omega = _sp.dot(t,r)/_sp.dot(t,t);
719 
720  // apply second correction to x
721  // x <- x + omega y
722  x.axpy(omega,y);
723 
724  // r = s - omega*t (remember : r = s)
725  r.axpy(-omega,t);
726 
727  rho = rho_new;
728 
729  //
730  // test stop criteria
731  //
732 
733  norm = _sp.norm(r);
734 
735  if (_verbose > 1) // print
736  {
737  this->printOutput(std::cout,it,norm,norm_old);
738  }
739 
740  if ( norm < (_reduction * norm_0) || norm<1E-30)
741  {
742  res.converged = 1;
743  break;
744  }
745 
746  norm_old = norm;
747  } // end for
748 
749  //correct i which is wrong if convergence was not achieved.
750  it=std::min((double)_maxit,it);
751 
752  if (_verbose==1) // printing for non verbose
753  this->printOutput(std::cout,it,norm);
754 
755  _prec.post(x); // postprocess preconditioner
756  res.iterations = static_cast<int>(std::ceil(it)); // fill statistics
757  res.reduction = norm/norm_0;
758  res.conv_rate = pow(res.reduction,1.0/it);
759  res.elapsed = watch.elapsed();
760  if (_verbose>0) // final print
761  std::cout << "=== rate=" << res.conv_rate
762  << ", T=" << res.elapsed
763  << ", TIT=" << res.elapsed/it
764  << ", IT=" << it << std::endl;
765  }
766 
772  virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
773  {
774  std::swap(_reduction,reduction);
775  (*this).apply(x,b,res);
776  std::swap(_reduction,reduction);
777  }
778 
779  private:
781  LinearOperator<X,X>& _op;
782  Preconditioner<X,X>& _prec;
783  ScalarProduct<X>& _sp;
784  double _reduction;
785  int _maxit;
786  int _verbose;
787  };
788 
795  template<class X>
796  class MINRESSolver : public InverseOperator<X,X> {
797  public:
799  typedef X domain_type;
801  typedef X range_type;
803  typedef typename X::field_type field_type;
805  typedef typename FieldTraits<field_type>::real_type real_type;
806 
812  template<class L, class P>
813  MINRESSolver (L& op, P& prec, double reduction, int maxit, int verbose) :
814  ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
815  {
816  dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
817  "L and P must have the same category!");
818  dune_static_assert( static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
819  "L must be sequential!");
820  }
826  template<class L, class S, class P>
827  MINRESSolver (L& op, S& sp, P& prec, double reduction, int maxit, int verbose) :
828  _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
829  {
830  dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
831  "L and P must have the same category!");
832  dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category),
833  "L and S must have the same category!");
834  }
835 
841  virtual void apply (X& x, X& b, InverseOperatorResult& res)
842  {
843  res.clear(); // clear solver statistics
844  Timer watch; // start a timer
845  _prec.pre(x,b); // prepare preconditioner
846  _op.applyscaleadd(-1,x,b); // overwrite b with defect/residual
847 
848  real_type def0 = _sp.norm(b); // compute residual norm
849 
850  if (def0<1E-30) // convergence check
851  {
852  res.converged = true;
853  res.iterations = 0; // fill statistics
854  res.reduction = 0;
855  res.conv_rate = 0;
856  res.elapsed=0;
857  if (_verbose>0) // final print
858  std::cout << "=== rate=" << res.conv_rate << ", T=" << res.elapsed << ", TIT=" << res.elapsed << ", IT=0" << std::endl;
859  return;
860  }
861 
862  if (_verbose>0) // printing
863  {
864  std::cout << "=== MINRESSolver" << std::endl;
865  if (_verbose>1) {
866  this->printHeader(std::cout);
867  this->printOutput(std::cout,0,def0);
868  }
869  }
870 
871  // some local variables
872  real_type def=def0; // the defect/residual norm
873  field_type alpha, // recurrence coefficients as computed in the Lanczos alg making up the matrix T
874  c[2]={0.0, 0.0}, // diagonal entry of Givens rotation
875  s[2]={0.0, 0.0}; // off-diagonal entries of Givens rotation
876  real_type beta;
877 
878  field_type T[3]={0.0, 0.0, 0.0}; // recurrence coefficients (column k of Matrix T)
879 
880  X z(b.size()), // some temporary vectors
881  dummy(b.size());
882 
883  field_type xi[2]={1.0, 0.0};
884 
885  // initialize
886  z = 0.0; // clear correction
887 
888  _prec.apply(z,b); // apply preconditioner z=M^-1*b
889 
890  beta = std::sqrt(std::abs(_sp.dot(z,b)));
891  real_type beta0 = beta;
892 
893  X p[3]; // the search directions
894  X q[3]; // Orthonormal basis vectors (in unpreconditioned case)
895 
896  q[0].resize(b.size());
897  q[1].resize(b.size());
898  q[2].resize(b.size());
899  q[0] = 0.0;
900  q[1] = b;
901  q[1] /= beta;
902  q[2] = 0.0;
903 
904  p[0].resize(b.size());
905  p[1].resize(b.size());
906  p[2].resize(b.size());
907  p[0] = 0.0;
908  p[1] = 0.0;
909  p[2] = 0.0;
910 
911 
912  z /= beta; // this is w_current
913 
914  // the loop
915  int i=1;
916  for ( ; i<=_maxit; i++)
917  {
918  dummy = z; // remember z_old for the computation of the search direction p in the next iteration
919 
920  int i1 = i%3,
921  i0 = (i1+2)%3,
922  i2 = (i1+1)%3;
923 
924  // Symmetrically Preconditioned Lanczos (Greenbaum p.121)
925  _op.apply(z,q[i2]); // q[i2] = Az
926  q[i2].axpy(-beta, q[i0]);
927  alpha = _sp.dot(q[i2],z);
928  q[i2].axpy(-alpha, q[i1]);
929 
930  z=0.0;
931  _prec.apply(z,q[i2]);
932 
933  beta = std::sqrt(std::abs(_sp.dot(q[i2],z)));
934 
935  q[i2] /= beta;
936  z /= beta;
937 
938  // QR Factorization of recurrence coefficient matrix
939  // apply previous Givens rotations to last column of T
940  T[1] = T[2];
941  if (i>2)
942  {
943  T[0] = s[i%2]*T[1];
944  T[1] = c[i%2]*T[1];
945  }
946  if (i>1)
947  {
948  T[2] = c[(i+1)%2]*alpha - s[(i+1)%2]*T[1];
949  T[1] = c[(i+1)%2]*T[1] + s[(i+1)%2]*alpha;
950  }
951  else
952  T[2] = alpha;
953 
954  // recompute c, s -> current Givens rotation \TODO use BLAS-routine drotg instead for greater robustness
955  // cblas_drotg (a, b, c, s);
956  c[i%2] = 1.0/std::sqrt(T[2]*T[2] + beta*beta);
957  s[i%2] = beta*c[i%2];
958  c[i%2] *= T[2];
959 
960  // apply current Givens rotation to T eliminating the last entry...
961  T[2] = c[i%2]*T[2] + s[i%2]*beta;
962 
963  // ...and to xi, the right hand side of the least squares problem min_y||beta*xi-T*y||
964  xi[i%2] = -s[i%2]*xi[(i+1)%2];
965  xi[(i+1)%2] *= c[i%2];
966 
967  // compute correction direction
968  p[i2] = dummy;
969  p[i2].axpy(-T[1],p[i1]);
970  p[i2].axpy(-T[0],p[i0]);
971  p[i2] /= T[2];
972 
973  // apply correction/update solution
974  x.axpy(beta0*xi[(i+1)%2], p[i2]);
975 
976  // remember beta_old
977  T[2] = beta;
978 
979  // update residual - not necessary if in the preconditioned case we are content with the residual norm of the
980  // preconditioned system as convergence test
981  // _op.apply(p[i2],dummy);
982  // b.axpy(-beta0*xi[(i+1)%2],dummy);
983 
984  // convergence test
985  real_type defnew = std::abs(beta0*xi[i%2]); // the last entry the QR-transformed least squares RHS is the new residual norm
986 
987  if (_verbose>1) // print
988  this->printOutput(std::cout,i,defnew,def);
989 
990  def = defnew; // update norm
991  if (def<def0*_reduction || def<1E-30 || i==_maxit) // convergence check
992  {
993  res.converged = true;
994  break;
995  }
996  }
997 
998  //correct i which is wrong if convergence was not achieved.
999  i=std::min(_maxit,i);
1000 
1001  if (_verbose==1) // printing for non verbose
1002  this->printOutput(std::cout,i,def);
1003 
1004  _prec.post(x); // postprocess preconditioner
1005  res.iterations = i; // fill statistics
1006  res.reduction = def/def0;
1007  res.conv_rate = pow(res.reduction,1.0/i);
1008  res.elapsed = watch.elapsed();
1009 
1010  if (_verbose>0) // final print
1011  {
1012  std::cout << "=== rate=" << res.conv_rate
1013  << ", T=" << res.elapsed
1014  << ", TIT=" << res.elapsed/i
1015  << ", IT=" << i << std::endl;
1016  }
1017 
1018  }
1019 
1025  virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
1026  {
1027  std::swap(_reduction,reduction);
1028  (*this).apply(x,b,res);
1029  std::swap(_reduction,reduction);
1030  }
1031 
1032  private:
1033  SeqScalarProduct<X> ssp;
1034  LinearOperator<X,X>& _op;
1035  Preconditioner<X,X>& _prec;
1036  ScalarProduct<X>& _sp;
1037  double _reduction;
1038  int _maxit;
1039  int _verbose;
1040  };
1041 
1053  template<class X, class Y=X, class F = Y>
1055  {
1056  public:
1058  typedef X domain_type;
1060  typedef Y range_type;
1062  typedef typename X::field_type field_type;
1064  typedef typename FieldTraits<field_type>::real_type real_type;
1066  typedef F basis_type;
1067 
1075  template<class L, class P>
1076  RestartedGMResSolver (L& op, P& prec, double reduction, int restart, int maxit, int verbose, bool recalc_defect = false) :
1077  _A_(op), _M(prec),
1078  ssp(), _sp(ssp), _restart(restart),
1079  _reduction(reduction), _maxit(maxit), _verbose(verbose),
1080  _recalc_defect(recalc_defect)
1081  {
1082  dune_static_assert(static_cast<int>(P::category) == static_cast<int>(L::category),
1083  "P and L must be the same category!");
1084  dune_static_assert( static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
1085  "L must be sequential!");
1086  }
1087 
1095  template<class L, class S, class P>
1096  RestartedGMResSolver (L& op, S& sp, P& prec, double reduction, int restart, int maxit, int verbose, bool recalc_defect = false) :
1097  _A_(op), _M(prec),
1098  _sp(sp), _restart(restart),
1099  _reduction(reduction), _maxit(maxit), _verbose(verbose),
1100  _recalc_defect(recalc_defect)
1101  {
1102  dune_static_assert(static_cast<int>(P::category) == static_cast<int>(L::category),
1103  "P and L must have the same category!");
1104  dune_static_assert(static_cast<int>(P::category) == static_cast<int>(S::category),
1105  "P and S must have the same category!");
1106  }
1107 
1109  virtual void apply (X& x, X& b, InverseOperatorResult& res)
1110  {
1111  apply(x,b,_reduction,res);
1112  }
1113 
1119  virtual void apply (X& x, Y& b, double reduction, InverseOperatorResult& res)
1120  {
1121  int m =_restart;
1122  real_type norm;
1123  real_type norm_old = 0.0;
1124  real_type norm_0;
1125  real_type beta;
1126  int i, j = 1, k;
1127  std::vector<field_type> s(m+1), cs(m), sn(m);
1128  // helper vector
1129  X w(b);
1130  std::vector< std::vector<field_type> > H(m+1,s);
1131  std::vector<F> v(m+1,b);
1132 
1133  // start timer
1134  Timer watch; // start a timer
1135 
1136  // clear solver statistics
1137  res.clear();
1138  _M.pre(x,b);
1139  if (_recalc_defect)
1140  {
1141  // norm_0 = norm(M^-1 b)
1142  w = 0.0; _M.apply(w,b); // w = M^-1 b
1143  norm_0 = _sp.norm(w); // use defect of preconditioned residual
1144  // r = _M.solve(b - A * x);
1145  w = b;
1146  _A_.applyscaleadd(-1,x, /* => */ w); // w = b - Ax;
1147  v[0] = 0.0; _M.apply(v[0],w); // r = M^-1 w
1148  beta = _sp.norm(v[0]);
1149  }
1150  else
1151  {
1152  // norm_0 = norm(b-Ax)
1153  _A_.applyscaleadd(-1,x, /* => */ b); // b = b - Ax;
1154  v[0] = 0.0; _M.apply(v[0],b); // r = M^-1 b
1155  beta = _sp.norm(v[0]);
1156  norm_0 = beta; // use defect of preconditioned residual
1157  }
1158 
1159  // avoid division by zero
1160  if (norm_0 == 0.0)
1161  norm_0 = 1.0;
1162  norm = norm_old = beta;
1163 
1164  // print header
1165  if (_verbose > 0)
1166  {
1167  std::cout << "=== RestartedGMResSolver" << std::endl;
1168  if (_verbose > 1)
1169  {
1170  this->printHeader(std::cout);
1171  this->printOutput(std::cout,0,norm_0);
1172  }
1173  }
1174 
1175  // check convergence
1176  if (norm <= reduction * norm_0) {
1177  _M.post(x); // postprocess preconditioner
1178  res.converged = true;
1179  if (_verbose > 0) // final print
1180  print_result(res);
1181  return;
1182  }
1183 
1184  while (j <= _maxit && res.converged != true) {
1185  v[0] *= (1.0 / beta);
1186  for (i=1; i<=m; i++) s[i] = 0.0;
1187  s[0] = beta;
1188  int end=std::min(m, _maxit-j+1);
1189  for (i = 0; i < end && res.converged != true; i++, j++) {
1190  w = 0.0;
1191  v[i+1] = 0.0; // use v[i+1] as temporary vector
1192  _A_.apply(v[i], /* => */ v[i+1]);
1193  _M.apply(w, v[i+1]);
1194  for (k = 0; k <= i; k++) {
1195  H[k][i] = _sp.dot(w, v[k]);
1196  // w -= H[k][i] * v[k];
1197  w.axpy(-H[k][i], v[k]);
1198  }
1199  H[i+1][i] = _sp.norm(w);
1200  if (H[i+1][i] == 0.0)
1201  DUNE_THROW(ISTLError,"breakdown in GMRes - |w| "
1202  << " == 0.0 after " << j << " iterations");
1203  // v[i+1] = w * (1.0 / H[i+1][i]);
1204  v[i+1] = w; v[i+1] *= (1.0 / H[i+1][i]);
1205 
1206  for (k = 0; k < i; k++)
1207  applyPlaneRotation(H[k][i], H[k+1][i], cs[k], sn[k]);
1208 
1209  generatePlaneRotation(H[i][i], H[i+1][i], cs[i], sn[i]);
1210  applyPlaneRotation(H[i][i], H[i+1][i], cs[i], sn[i]);
1211  applyPlaneRotation(s[i], s[i+1], cs[i], sn[i]);
1212 
1213  norm = std::abs(s[i+1]);
1214 
1215  if (_verbose > 1) // print
1216  {
1217  this->printOutput(std::cout,j,norm,norm_old);
1218  }
1219 
1220  norm_old = norm;
1221 
1222  if (norm < reduction * norm_0) {
1223  res.converged = true;
1224  }
1225  }
1226 
1227  if (_recalc_defect)
1228  {
1229  // update x
1230  update(x, i - 1, H, s, v);
1231 
1232  // update residuum
1233  // r = M^-1 (b - A * x);
1234  w = b; _A_.applyscaleadd(-1,x, /* => */ w);
1235  _M.apply(v[0], w);
1236  beta = _sp.norm(v[0]);
1237  norm = beta;
1238  }
1239  else
1240  {
1241  // calc update vector
1242  w = 0;
1243  update(w, i - 1, H, s, v);
1244 
1245  // update x
1246  x += w;
1247 
1248  // r = M^-1 (b - A * x);
1249  // update defect
1250  _A_.applyscaleadd(-1,w, /* => */ b);
1251  // r = M^-1 (b - A * x);
1252  v[0] = 0.0; _M.apply(v[0],b); // r = M^-1 b
1253  beta = _sp.norm(v[0]);
1254  norm = beta;
1255 
1256  res.converged = false;
1257  }
1258 
1259  norm_old = norm;
1260 
1261  if (norm < reduction * norm_0) {
1262  // fill statistics
1263  res.converged = true;
1264  }
1265 
1266  if (res.converged != true && _verbose > 0)
1267  std::cout << "=== GMRes::restart\n";
1268  }
1269 
1270  _M.post(x); // postprocess preconditioner
1271 
1272  res.iterations = j;
1273  res.reduction = norm / norm_0;
1274  res.conv_rate = pow(res.reduction,1.0/j);
1275  res.elapsed = watch.elapsed();
1276 
1277  if (_verbose>0)
1278  print_result(res);
1279  }
1280  private:
1281 
1282  void
1283  print_result (const InverseOperatorResult & res) const
1284  {
1285  int j = res.iterations>0 ? res.iterations : 1;
1286  std::cout << "=== rate=" << res.conv_rate
1287  << ", T=" << res.elapsed
1288  << ", TIT=" << res.elapsed/j
1289  << ", IT=" << res.iterations
1290  << std::endl;
1291  }
1292 
1293  static void
1294  update(X &x, int k,
1295  const std::vector< std::vector<field_type> > & h,
1296  const std::vector<field_type> & s, const std::vector<F> &v)
1297  {
1298  std::vector<field_type> y(s);
1299 
1300  // Backsolve:
1301  for (int i = k; i >= 0; i--) {
1302  y[i] /= h[i][i];
1303  for (int j = i - 1; j >= 0; j--)
1304  y[j] -= h[j][i] * y[i];
1305  }
1306 
1307  for (int j = 0; j <= k; j++)
1308  // x += v[j] * y[j];
1309  x.axpy(y[j],v[j]);
1310  }
1311 
1312  void
1313  generatePlaneRotation(field_type &dx, field_type &dy, field_type &cs, field_type &sn)
1314  {
1315  if (dy == 0.0) {
1316  cs = 1.0;
1317  sn = 0.0;
1318  } else if (std::abs(dy) > std::abs(dx)) {
1319  field_type temp = dx / dy;
1320  sn = 1.0 / std::sqrt( 1.0 + temp*temp );
1321  cs = temp * sn;
1322  } else {
1323  field_type temp = dy / dx;
1324  cs = 1.0 / std::sqrt( 1.0 + temp*temp );
1325  sn = temp * cs;
1326  }
1327  }
1328 
1329 
1330  void
1331  applyPlaneRotation(field_type &dx, field_type &dy, field_type &cs, field_type &sn)
1332  {
1333  field_type temp = cs * dx + sn * dy;
1334  dy = -sn * dx + cs * dy;
1335  dx = temp;
1336  }
1337 
1338  LinearOperator<X,X>& _A_;
1339  Preconditioner<X,X>& _M;
1340  SeqScalarProduct<X> ssp;
1341  ScalarProduct<X>& _sp;
1342  int _restart;
1343  double _reduction;
1344  int _maxit;
1345  int _verbose;
1346  bool _recalc_defect;
1347  };
1348 
1349 
1363  template<class X>
1365  {
1366  public:
1368  typedef X domain_type;
1370  typedef X range_type;
1372  typedef typename X::field_type field_type;
1373 
1381  template<class L, class P>
1382  GeneralizedPCGSolver (L& op, P& prec, double reduction, int maxit, int verbose,
1383  int restart=10) :
1384  ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit),
1385  _verbose(verbose), _restart(std::min(maxit,restart))
1386  {
1387  dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
1388  "L and P have to have the same category!");
1389  dune_static_assert(static_cast<int>(L::category) ==
1390  static_cast<int>(SolverCategory::sequential),
1391  "L has to be sequential!");
1392  }
1400  template<class L, class P, class S>
1401  GeneralizedPCGSolver (L& op, S& sp, P& prec,
1402  double reduction, int maxit, int verbose, int restart=10) :
1403  _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose),
1404  _restart(std::min(maxit,restart))
1405  {
1406  dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
1407  "L and P must have the same category!");
1408  dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category),
1409  "L and S must have the same category!");
1410  }
1416  virtual void apply (X& x, X& b, InverseOperatorResult& res)
1417  {
1418  res.clear(); // clear solver statistics
1419  Timer watch; // start a timer
1420  _prec.pre(x,b); // prepare preconditioner
1421  _op.applyscaleadd(-1,x,b); // overwrite b with defect
1422 
1423  std::vector<shared_ptr<X> > p(_restart);
1424  std::vector<typename X::field_type> pp(_restart);
1425  X q(x); // a temporary vector
1426  X prec_res(x); // a temporary vector for preconditioner output
1427 
1428  p[0].reset(new X(x));
1429 
1430  double def0 = _sp.norm(b); // compute norm
1431  if (def0<1E-30) // convergence check
1432  {
1433  res.converged = true;
1434  res.iterations = 0; // fill statistics
1435  res.reduction = 0;
1436  res.conv_rate = 0;
1437  res.elapsed=0;
1438  if (_verbose>0) // final print
1439  std::cout << "=== rate=" << res.conv_rate
1440  << ", T=" << res.elapsed << ", TIT=" << res.elapsed
1441  << ", IT=0" << std::endl;
1442  return;
1443  }
1444 
1445  if (_verbose>0) // printing
1446  {
1447  std::cout << "=== GeneralizedPCGSolver" << std::endl;
1448  if (_verbose>1) {
1449  this->printHeader(std::cout);
1450  this->printOutput(std::cout,0,def0);
1451  }
1452  }
1453  // some local variables
1454  double def=def0; // loop variables
1455  field_type rho, lambda;
1456 
1457  int i=0;
1458  int ii=0;
1459  // determine initial search direction
1460  *(p[0]) = 0; // clear correction
1461  _prec.apply(*(p[0]),b); // apply preconditioner
1462  rho = _sp.dot(*(p[0]),b); // orthogonalization
1463  _op.apply(*(p[0]),q); // q=Ap
1464  pp[0] = _sp.dot(*(p[0]),q); // scalar product
1465  lambda = rho/pp[0]; // minimization
1466  x.axpy(lambda,*(p[0])); // update solution
1467  b.axpy(-lambda,q); // update defect
1468 
1469  // convergence test
1470  double defnew=_sp.norm(b); // comp defect norm
1471  if (_verbose>1) // print
1472  this->printOutput(std::cout,++i,defnew,def);
1473  def = defnew; // update norm
1474  if (def<def0*_reduction || def<1E-30) // convergence check
1475  {
1476  res.converged = true;
1477  if (_verbose>0) // final print
1478  {
1479  std::cout << "=== rate=" << res.conv_rate
1480  << ", T=" << res.elapsed
1481  << ", TIT=" << res.elapsed
1482  << ", IT=" << 1 << std::endl;
1483  }
1484  return;
1485  }
1486 
1487  while(i<_maxit) {
1488  // the loop
1489  int end=std::min(_restart, _maxit-i+1);
1490  for (ii=1; ii<end; ++ii )
1491  {
1492  //std::cout<<" ii="<<ii<<" i="<<i<<std::endl;
1493  // compute next conjugate direction
1494  prec_res = 0; // clear correction
1495  _prec.apply(prec_res,b); // apply preconditioner
1496 
1497  p[ii].reset(new X(prec_res));
1498  _op.apply(prec_res, q);
1499 
1500  for(int j=0; j<ii; ++j) {
1501  rho =_sp.dot(q,*(p[j]))/pp[j];
1502  p[ii]->axpy(-rho, *(p[j]));
1503  }
1504 
1505  // minimize in given search direction
1506  _op.apply(*(p[ii]),q); // q=Ap
1507  pp[ii] = _sp.dot(*(p[ii]),q); // scalar product
1508  rho = _sp.dot(*(p[ii]),b); // orthogonalization
1509  lambda = rho/pp[ii]; // minimization
1510  x.axpy(lambda,*(p[ii])); // update solution
1511  b.axpy(-lambda,q); // update defect
1512 
1513  // convergence test
1514  double defnew=_sp.norm(b); // comp defect norm
1515 
1516  if (_verbose>1) // print
1517  this->printOutput(std::cout,++i,defnew,def);
1518 
1519  def = defnew; // update norm
1520  if (def<def0*_reduction || def<1E-30) // convergence check
1521  {
1522  res.converged = true;
1523  break;
1524  }
1525  }
1526  if(res.converged)
1527  break;
1528  if(end==_restart) {
1529  *(p[0])=*(p[_restart-1]);
1530  pp[0]=pp[_restart-1];
1531  }
1532  }
1533 
1534  // postprocess preconditioner
1535  _prec.post(x);
1536 
1537  // fill statistics
1538  res.iterations = i;
1539  res.reduction = def/def0;
1540  res.conv_rate = pow(res.reduction,1.0/i);
1541  res.elapsed = watch.elapsed();
1542 
1543  if (_verbose>0) // final print
1544  {
1545  std::cout << "=== rate=" << res.conv_rate
1546  << ", T=" << res.elapsed
1547  << ", TIT=" << res.elapsed/i
1548  << ", IT=" << i+1 << std::endl;
1549  }
1550  }
1551 
1557  virtual void apply (X& x, X& b, double reduction,
1558  InverseOperatorResult& res)
1559  {
1560  std::swap(_reduction,reduction);
1561  (*this).apply(x,b,res);
1562  std::swap(_reduction,reduction);
1563  }
1564  private:
1565  SeqScalarProduct<X> ssp;
1566  LinearOperator<X,X>& _op;
1567  Preconditioner<X,X>& _prec;
1568  ScalarProduct<X>& _sp;
1569  double _reduction;
1570  int _maxit;
1571  int _verbose;
1572  int _restart;
1573  };
1574 
1577 } // end namespace
1578 
1579 #endif
FieldTraits< field_type >::real_type real_type
The real type of the field type (is the same of using real numbers, but differs for std::complex) ...
Definition: solvers.hh:537
X domain_type
The domain type of the operator to be inverted.
Definition: solvers.hh:1368
void clear()
Resets all data.
Definition: solver.hh:40
RestartedGMResSolver(L &op, S &sp, P &prec, double reduction, int restart, int maxit, int verbose, bool recalc_defect=false)
Set up solver.
Definition: solvers.hh:1096
Default implementation for the scalar case.
Definition: scalarproducts.hh:94
virtual void apply(X &x, X &b, double reduction, InverseOperatorResult &res)
Apply inverse operator with given reduction factor.
Definition: solvers.hh:1557
Minimal Residual Method (MINRES)
Definition: solvers.hh:796
X range_type
The range type of the operator to be inverted.
Definition: solvers.hh:801
Generalized preconditioned conjugate gradient solver.
Definition: solvers.hh:1364
void printOutput(std::ostream &s, const double iter, const DataType &norm, const DataType &norm_old) const
helper function for printing solver output
Definition: solver.hh:130
CGSolver(L &op, P &prec, double reduction, int maxit, int verbose)
Set up conjugate gradient solver.
Definition: solvers.hh:375
Base class for scalar product and norm computation.
Definition: scalarproducts.hh:43
virtual void apply(X &x, X &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: solvers.hh:205
MINRESSolver(L &op, P &prec, double reduction, int maxit, int verbose)
Set up MINRES solver.
Definition: solvers.hh:813
virtual void apply(const X &x, Y &y) const =0
apply operator to x: The input vector is consistent and the output must also be consistent on the in...
X::field_type field_type
The field type of the operator to be inverted.
Definition: solvers.hh:803
Abstract base class for all solvers.
Definition: solver.hh:79
double reduction
Reduction achieved: .
Definition: solver.hh:53
X range_type
The range type of the operator that we do the inverse for.
Definition: solvers.hh:231
bool converged
True if convergence criterion has been met.
Definition: solver.hh:56
X domain_type
The domain type of the operator to be inverted.
Definition: solvers.hh:799
X::field_type field_type
The field type of the operator to be inverted.
Definition: solvers.hh:1372
X::field_type field_type
The field type of the operator that we do the inverse for.
Definition: solvers.hh:233
virtual void apply(X &x, X &b, double reduction, InverseOperatorResult &res)
Apply inverse operator with given reduction factor.
Definition: solvers.hh:1025
X::field_type field_type
The field type of the operator to be inverted.
Definition: solvers.hh:1062
CGSolver(L &op, S &sp, P &prec, double reduction, int maxit, int verbose)
Set up conjugate gradient solver.
Definition: solvers.hh:389
X domain_type
The domain type of the operator to be inverted.
Definition: solvers.hh:531
virtual void apply(X &x, X &b, InverseOperatorResult &res)
Apply inverse operator.
Definition: solvers.hh:573
virtual void apply(X &x, X &b, double reduction, InverseOperatorResult &res)
Apply inverse operator with given reduction factor.
Definition: solvers.hh:506
X range_type
The range type of the operator to be inverted.
Definition: solvers.hh:365
X::field_type field_type
The field type of the operator that we do the inverse for.
Definition: solvers.hh:62
FieldTraits< field_type >::real_type real_type
The real type of the field type (is the same of using real numbers, but differs for std::complex) ...
Definition: solvers.hh:1064
X domain_type
The domain type of the operator that we do the inverse for.
Definition: solvers.hh:58
virtual void apply(X &x, X &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: solvers.hh:1109
GeneralizedPCGSolver(L &op, S &sp, P &prec, double reduction, int maxit, int verbose, int restart=10)
Set up nonlinear preconditioned conjugate gradient solver.
Definition: solvers.hh:1401
Y range_type
The range type of the operator to be inverted.
Definition: solvers.hh:1060
int iterations
Number of iterations.
Definition: solver.hh:50
X range_type
The range type of the operator that we do the inverse for.
Definition: solvers.hh:60
LoopSolver(L &op, S &sp, P &prec, double reduction, int maxit, int verbose)
Set up loop solver.
Definition: solvers.hh:115
virtual void apply(X &x, Y &b, double reduction, InverseOperatorResult &res)
Apply inverse operator.
Definition: solvers.hh:1119
Bi-conjugate Gradient Stabilized (BiCG-STAB)
Definition: solvers.hh:528
Define general, extensible interface for inverse operators.
Category for sequential solvers.
Definition: solvercategory.hh:22
gradient method
Definition: solvers.hh:226
virtual void apply(X &x, X &b, double reduction, InverseOperatorResult &res)
Apply inverse operator with given reduction factor.
Definition: solvers.hh:339
virtual void pre(X &x, Y &b)=0
Prepare the preconditioner.
FieldTraits< field_type >::real_type real_type
The real type of the field type (is the same of using real numbers, but differs for std::complex) ...
Definition: solvers.hh:805
BiCGSTABSolver(L &op, P &prec, double reduction, int maxit, int verbose)
Set up solver.
Definition: solvers.hh:545
Statistics about the application of an inverse operator.
Definition: solver.hh:31
Define general, extensible interface for operators. The available implementation wraps a matrix...
MINRESSolver(L &op, S &sp, P &prec, double reduction, int maxit, int verbose)
Set up MINRES solver.
Definition: solvers.hh:827
Preconditioned loop solver.
Definition: solvers.hh:55
conjugate gradient method
Definition: solvers.hh:360
double elapsed
Elapsed time in seconds.
Definition: solver.hh:62
X::field_type field_type
The field type of the operator to be inverted.
Definition: solvers.hh:367
X domain_type
The domain type of the operator that we do the inverse for.
Definition: solvers.hh:229
virtual void apply(X &x, X &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: solvers.hh:127
virtual void apply(X &x, X &b, InverseOperatorResult &res)
Apply inverse operator.
Definition: solvers.hh:271
implements the Generalized Minimal Residual (GMRes) method
Definition: solvers.hh:1054
virtual void applyscaleadd(field_type alpha, const X &x, Y &y) const =0
apply operator to x, scale and add:
X domain_type
The domain type of the operator to be inverted.
Definition: solvers.hh:363
derive error class from the base class in common
Definition: istlexception.hh:16
LoopSolver(L &op, P &prec, double reduction, int maxit, int verbose)
Set up Loop solver.
Definition: solvers.hh:84
virtual void apply(X &x, X &b, InverseOperatorResult &res)
Apply inverse operator.
Definition: solvers.hh:841
X range_type
The range type of the operator to be inverted.
Definition: solvers.hh:1370
virtual void apply(X &x, X &b, double reduction, InverseOperatorResult &res)
Apply inverse operator with given reduction factor.
Definition: solvers.hh:772
double conv_rate
Convergence rate (average reduction per step)
Definition: solver.hh:59
BiCGSTABSolver(L &op, S &sp, P &prec, double reduction, int maxit, int verbose)
Set up solver.
Definition: solvers.hh:558
GradientSolver(L &op, S &sp, P &prec, double reduction, int maxit, int verbose)
Set up solver.
Definition: solvers.hh:256
virtual void apply(X &x, X &b, InverseOperatorResult &res)
Apply inverse operator.
Definition: solvers.hh:1416
void printHeader(std::ostream &s) const
helper function for printing header of solver output
Definition: solver.hh:121
X domain_type
The domain type of the operator to be inverted.
Definition: solvers.hh:1058
virtual void post(X &x)=0
Clean up.
Define base class for scalar product and norm.
virtual void apply(X &x, X &b, InverseOperatorResult &res)
Apply inverse operator.
Definition: solvers.hh:403
X range_type
The range type of the operator to be inverted.
Definition: solvers.hh:533
GeneralizedPCGSolver(L &op, P &prec, double reduction, int maxit, int verbose, int restart=10)
Set up nonlinear preconditioned conjugate gradient solver.
Definition: solvers.hh:1382
RestartedGMResSolver(L &op, P &prec, double reduction, int restart, int maxit, int verbose, bool recalc_defect=false)
Set up solver.
Definition: solvers.hh:1076
F basis_type
The field type of the basis vectors.
Definition: solvers.hh:1066
X::field_type field_type
The field type of the operator to be inverted.
Definition: solvers.hh:535
virtual void apply(X &v, const Y &d)=0
Apply one step of the preconditioner to the system A(v)=d.
GradientSolver(L &op, P &prec, double reduction, int maxit, int verbose)
Set up solver.
Definition: solvers.hh:241
Definition: example.cc:34