solvers.hh

Go to the documentation of this file.
00001 #ifndef DUNE_SOLVERS_HH
00002 #define DUNE_SOLVERS_HH
00003 
00004 #include<cmath>
00005 #include<complex>
00006 #include<iostream>
00007 #include<iomanip>
00008 #include<string>
00009 
00010 #include "istlexception.hh"
00011 #include "operators.hh"
00012 #include "preconditioners.hh"
00013 #include "scalarproducts.hh"
00014 #include <dune/common/timer.hh>
00015 #include <dune/common/helpertemplates.hh>
00016 
00017 namespace Dune {
00042   struct InverseOperatorResult
00043   {
00045   InverseOperatorResult ()
00046   {
00047     clear();
00048   }
00049 
00051   void clear ()
00052   {
00053     iterations = 0;
00054     reduction = 0;
00055     converged = false;
00056     conv_rate = 1;
00057     elapsed = 0;
00058   }
00059 
00061   int iterations; 
00062 
00064   double reduction;
00065 
00067   bool converged; 
00068 
00070   double conv_rate;
00071 
00073   double elapsed;
00074   };
00075 
00076 
00077   //=====================================================================
00089   template<class X, class Y>
00090   class InverseOperator {
00091   public:
00093     typedef X domain_type;
00094     
00096     typedef Y range_type;
00097     
00099     typedef typename X::field_type field_type;
00100 
00110     virtual void apply (X& x, Y& b, InverseOperatorResult& res) = 0;
00111 
00122     virtual void apply (X& x, Y& b, double reduction, InverseOperatorResult& res) = 0;
00123 
00125     virtual ~InverseOperator () {}
00126 
00127   protected: 
00128     // spacing values 
00129     enum { iterationSpacing = 5 , normSpacing = 16 }; 
00130 
00132     void printHeader(std::ostream& s) const 
00133     {
00134       s << std::setw(iterationSpacing)  << " Iter";
00135       s << std::setw(normSpacing) << "Defect";
00136       s << std::setw(normSpacing) << "Rate" << std::endl;
00137     }
00138 
00140     template <class DataType> 
00141     void printOutput(std::ostream& s, 
00142                      const int iter,
00143                      const DataType& norm,
00144                      const DataType& norm_old) const 
00145     {
00146       const DataType rate = norm/norm_old;
00147       s << std::setw(iterationSpacing)  << iter << " ";
00148       s << std::setw(normSpacing) << norm << " ";
00149       s << std::setw(normSpacing) << rate << std::endl;
00150     }
00151 
00153     template <class DataType> 
00154     void printOutput(std::ostream& s, 
00155                      const int iter,
00156                      const DataType& norm) const 
00157     {
00158       s << std::setw(iterationSpacing)  << iter << " ";
00159       s << std::setw(normSpacing) << norm << std::endl;
00160     }
00161   };
00162 
00163 
00164   //=====================================================================
00165   // Implementation of this interface
00166   //=====================================================================
00167 
00176   template<class X>
00177   class LoopSolver : public InverseOperator<X,X> {
00178   public:
00180     typedef X domain_type;
00182     typedef X range_type;
00184     typedef typename X::field_type field_type;
00185 
00205     template<class L, class P>
00206     LoopSolver (L& op, P& prec,
00207     double reduction, int maxit, int verbose) : 
00208       ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00209     {
00210       IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00211       IsTrue< static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential) >::yes();
00212     }
00213 
00234     template<class L, class S, class P>
00235     LoopSolver (L& op, S& sp, P& prec,
00236         double reduction, int maxit, int verbose) : 
00237       _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00238     {
00239       IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00240       IsTrue< static_cast<int>(L::category) == static_cast<int>(S::category) >::yes();
00241     }
00242 
00243 
00245   virtual void apply (X& x, X& b, InverseOperatorResult& res)
00246   {
00247     // clear solver statistics
00248     res.clear();
00249 
00250     // start a timer
00251     Timer watch;
00252 
00253     // prepare preconditioner
00254     _prec.pre(x,b);
00255 
00256     // overwrite b with defect
00257     _op.applyscaleadd(-1,x,b);
00258 
00259     // compute norm, \todo parallelization
00260     double def0 = _sp.norm(b);
00261 
00262     // printing
00263     if (_verbose>0)
00264     {
00265       std::cout << "=== LoopSolver" << std::endl;
00266       if (_verbose>1) 
00267       {
00268         this->printHeader(std::cout);
00269         this->printOutput(std::cout,0,def0);
00270       }
00271     }
00272 
00273     // allocate correction vector
00274     X v(x);
00275 
00276     // iteration loop
00277     int i=1; double def=def0;
00278     for ( ; i<=_maxit; i++ )  
00279     {     
00280       v = 0;                      // clear correction
00281       _prec.apply(v,b);           // apply preconditioner
00282       x += v;                     // update solution
00283       _op.applyscaleadd(-1,v,b);  // update defect
00284       double defnew=_sp.norm(b);// comp defect norm
00285       if (_verbose>1)             // print
00286         this->printOutput(std::cout,i,defnew,def);
00287                       //std::cout << i << " " << defnew << " " << defnew/def << std::endl;
00288       def = defnew;               // update norm
00289       if (def<def0*_reduction || def<1E-30)    // convergence check 
00290       {
00291         res.converged  = true;
00292         break;
00293       }
00294     }
00295 
00296     // print
00297     if (_verbose==1)
00298        this->printOutput(std::cout,i,def);
00299   
00300     // postprocess preconditioner
00301     _prec.post(x);
00302 
00303     // fill statistics
00304     res.iterations = i;
00305     res.reduction = def/def0;
00306     res.conv_rate  = pow(res.reduction,1.0/i);
00307     res.elapsed = watch.elapsed();
00308 
00309     // final print
00310     if (_verbose>0)
00311     {
00312       std::cout << "=== rate=" << res.conv_rate
00313                 << ", T=" << res.elapsed
00314                 << ", TIT=" << res.elapsed/i
00315                 << ", IT=" << i << std::endl;
00316     }
00317   }
00318 
00320   virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
00321   {
00322     _reduction = reduction;
00323     (*this).apply(x,b,res);
00324   }
00325 
00326   private:
00327   SeqScalarProduct<X> ssp;
00328   LinearOperator<X,X>& _op;
00329   Preconditioner<X,X>& _prec;
00330   ScalarProduct<X>& _sp;
00331   double _reduction;
00332   int _maxit;
00333   int _verbose;   
00334   };
00335 
00336 
00337   // all these solvers are taken from the SUMO library
00339   template<class X>
00340   class GradientSolver : public InverseOperator<X,X> {
00341   public:
00343     typedef X domain_type;
00345     typedef X range_type;
00347     typedef typename X::field_type field_type;
00348 
00354     template<class L, class P>
00355     GradientSolver (L& op, P& prec,
00356         double reduction, int maxit, int verbose) : 
00357       ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00358     {
00359       IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00360       IsTrue< static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential) >::yes();
00361     }
00367     template<class L, class S, class P>
00368     GradientSolver (L& op, S& sp, P& prec,
00369         double reduction, int maxit, int verbose) : 
00370       _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00371     {
00372       IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00373       IsTrue< static_cast<int>(L::category) == static_cast<int>(S::category) >::yes();
00374     }
00375 
00381     virtual void apply (X& x, X& b, InverseOperatorResult& res)
00382     {
00383       res.clear();                  // clear solver statistics
00384       Timer watch;                // start a timer
00385       _prec.pre(x,b);             // prepare preconditioner
00386       _op.applyscaleadd(-1,x,b);  // overwrite b with defect
00387       
00388       X p(x);                     // create local vectors
00389       X q(b);
00390   
00391       double def0 = _sp.norm(b);// compute norm
00392       
00393       if (_verbose>0)             // printing
00394       {
00395         std::cout << "=== GradientSolver" << std::endl;
00396         if (_verbose>1) 
00397         {
00398           this->printHeader(std::cout);
00399           this->printOutput(std::cout,0,def0);
00400         }
00401       }
00402       
00403       int i=1; double def=def0;   // loop variables
00404       field_type lambda;
00405       for ( ; i<=_maxit; i++ )  
00406       {
00407         p = 0;                      // clear correction
00408         _prec.apply(p,b);           // apply preconditioner
00409         _op.apply(p,q);             // q=Ap
00410         lambda = _sp.dot(p,b)/_sp.dot(q,p);// minimization
00411         x.axpy(lambda,p);           // update solution
00412         b.axpy(-lambda,q);          // update defect
00413         
00414         double defnew=_sp.norm(b);// comp defect norm
00415         if (_verbose>1)             // print
00416           this->printOutput(std::cout,i,defnew,def);
00417 
00418         def = defnew;               // update norm
00419         if (def<def0*_reduction || def<1E-30)    // convergence check 
00420         {
00421           res.converged  = true;
00422           break;
00423         }
00424       }
00425           
00426       if (_verbose==1)                // printing for non verbose
00427         this->printOutput(std::cout,i,def);
00428 
00429       _prec.post(x);                  // postprocess preconditioner
00430       res.iterations = i;               // fill statistics
00431       res.reduction = def/def0;
00432       res.conv_rate  = pow(res.reduction,1.0/i);
00433       res.elapsed = watch.elapsed();
00434       if (_verbose>0)                 // final print 
00435           std::cout << "=== rate=" << res.conv_rate
00436                     << ", T=" << res.elapsed
00437                     << ", TIT=" << res.elapsed/i
00438                     << ", IT=" << i << std::endl;
00439     }
00440 
00446     virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
00447     {
00448       _reduction = reduction;
00449       (*this).apply(x,b,res);
00450     }
00451 
00452   private:
00453   SeqScalarProduct<X> ssp;
00454   LinearOperator<X,X>& _op;
00455   Preconditioner<X,X>& _prec;
00456   ScalarProduct<X>& _sp;
00457   double _reduction;
00458   int _maxit;
00459   int _verbose;   
00460   };
00461 
00462 
00463 
00465   template<class X>
00466   class CGSolver : public InverseOperator<X,X> {
00467   public:
00469     typedef X domain_type;
00471     typedef X range_type;
00473     typedef typename X::field_type field_type;
00474 
00480     template<class L, class P>
00481     CGSolver (L& op, P& prec, double reduction, int maxit, int verbose) : 
00482       ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00483     {
00484       IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00485       IsTrue< static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential) >::yes();
00486     }
00492     template<class L, class S, class P>
00493     CGSolver (L& op, S& sp, P& prec, double reduction, int maxit, int verbose) : 
00494       _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00495     {
00496       IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00497       IsTrue< static_cast<int>(L::category) == static_cast<int>(S::category) >::yes();
00498     }
00499 
00505     virtual void apply (X& x, X& b, InverseOperatorResult& res)
00506     {
00507       res.clear();                  // clear solver statistics
00508       Timer watch;                // start a timer
00509       _prec.pre(x,b);             // prepare preconditioner
00510       _op.applyscaleadd(-1,x,b);  // overwrite b with defect
00511     
00512       X p(x);              // the search direction
00513       X q(x);              // a temporary vector
00514     
00515       double def0 = _sp.norm(b);// compute norm
00516       if (def0<1E-30)    // convergence check 
00517       {
00518         res.converged  = true;
00519         res.iterations = 0;               // fill statistics
00520         res.reduction = 0;
00521         res.conv_rate  = 0;
00522         res.elapsed=0;
00523         if (_verbose>0)                 // final print 
00524                         std::cout << "=== rate=" << res.conv_rate
00525                                   << ", T=" << res.elapsed << ", TIT=" << res.elapsed
00526                                   << ", IT=0" << std::endl;
00527         return;
00528       }
00529 
00530       if (_verbose>0)             // printing
00531       {
00532         std::cout << "=== CGSolver" << std::endl;
00533         if (_verbose>1) {
00534           this->printHeader(std::cout);
00535           this->printOutput(std::cout,0,def0);
00536         }
00537       }
00538 
00539       // some local variables
00540       double def=def0;   // loop variables
00541       field_type rho,rholast,lambda,alpha,beta;
00542 
00543       // determine initial search direction
00544       p = 0;                          // clear correction
00545       _prec.apply(p,b);               // apply preconditioner
00546       rholast = _sp.dot(p,b);         // orthogonalization
00547 
00548       // the loop
00549       int i=1; 
00550       for ( ; i<=_maxit; i++ )  
00551       {
00552         // minimize in given search direction p
00553         _op.apply(p,q);             // q=Ap
00554         alpha = _sp.dot(p,q);       // scalar product
00555         lambda = rholast/alpha;     // minimization
00556         x.axpy(lambda,p);           // update solution
00557         b.axpy(-lambda,q);          // update defect
00558 
00559         // convergence test
00560         double defnew=_sp.norm(b);// comp defect norm
00561 
00562         if (_verbose>1)             // print
00563           this->printOutput(std::cout,i,defnew,def);
00564 
00565         def = defnew;               // update norm
00566         if (def<def0*_reduction || def<1E-30 || i==_maxit)    // convergence check  
00567         {
00568           res.converged  = true;
00569           break;
00570         }
00571 
00572         // determine new search direction
00573         q = 0;                      // clear correction
00574         _prec.apply(q,b);           // apply preconditioner
00575         rho = _sp.dot(q,b);         // orthogonalization
00576         beta = rho/rholast;         // scaling factor
00577         p *= beta;                  // scale old search direction
00578         p += q;                     // orthogonalization with correction
00579         rholast = rho;              // remember rho for recurrence
00580       }
00581 
00582       if (_verbose==1)                // printing for non verbose
00583         this->printOutput(std::cout,i,def);
00584 
00585       _prec.post(x);                  // postprocess preconditioner
00586       res.iterations = i;               // fill statistics
00587       res.reduction = def/def0;
00588       res.conv_rate  = pow(res.reduction,1.0/i);
00589       res.elapsed = watch.elapsed();
00590 
00591       if (_verbose>0)                 // final print 
00592       {
00593         std::cout << "=== rate=" << res.conv_rate
00594                   << ", T=" << res.elapsed
00595                   << ", TIT=" << res.elapsed/i
00596                   << ", IT=" << i << std::endl;
00597       }
00598   }
00599 
00605   virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
00606   {
00607     _reduction = reduction;
00608     (*this).apply(x,b,res);
00609   }
00610 
00611   private:
00612   SeqScalarProduct<X> ssp;
00613   LinearOperator<X,X>& _op;
00614   Preconditioner<X,X>& _prec;
00615   ScalarProduct<X>& _sp;
00616   double _reduction;
00617   int _maxit;
00618   int _verbose;   
00619   };
00620 
00621 
00622   // Ronald Kriemanns BiCG-STAB implementation from Sumo
00624   template<class X>
00625   class BiCGSTABSolver : public InverseOperator<X,X> {
00626   public:
00628     typedef X domain_type;
00630     typedef X range_type;
00632     typedef typename X::field_type field_type;
00633 
00639     template<class L, class P>
00640     BiCGSTABSolver (L& op, P& prec,
00641         double reduction, int maxit, int verbose) : 
00642       ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00643     {
00644       IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00645       IsTrue< static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential) >::yes();
00646     }
00652     template<class L, class S, class P>
00653     BiCGSTABSolver (L& op, S& sp, P& prec,
00654         double reduction, int maxit, int verbose) : 
00655       _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00656     {
00657       IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00658       IsTrue< static_cast<int>(L::category) == static_cast<int>(S::category) >::yes();
00659     }
00660 
00666     virtual void apply (X& x, X& b, InverseOperatorResult& res)
00667     {
00668     const double EPSILON=1e-80;
00669     
00670     int                  it;
00671     field_type           rho, rho_new, alpha, beta, h, omega;
00672     field_type           norm, norm_old, norm_0;
00673     
00674     //
00675     // get vectors and matrix
00676     //
00677     X& r=b;
00678     X p(x);
00679     X v(x);
00680     X t(x);
00681     X y(x);
00682     X rt(x);
00683 
00684     //
00685     // begin iteration
00686     //
00687     
00688     // r = r - Ax; rt = r
00689       res.clear();                  // clear solver statistics
00690     Timer watch;                // start a timer
00691     _op.applyscaleadd(-1,x,r);  // overwrite b with defect
00692     _prec.pre(x,r);             // prepare preconditioner
00693 
00694     rt=r;
00695 
00696     norm = norm_old = norm_0 = _sp.norm(r);
00697     
00698     p=0;
00699     v=0;
00700 
00701     rho   = 1;
00702     alpha = 1;
00703     omega = 1;
00704 
00705     if (_verbose>0)             // printing
00706     {
00707       std::cout << "=== BiCGSTABSolver" << std::endl;
00708       if (_verbose>1) 
00709       {
00710         this->printHeader(std::cout);
00711         this->printOutput(std::cout,0,norm_0);
00712         //std::cout << " Iter       Defect         Rate" << std::endl;
00713         //std::cout << "    0" << std::setw(14) << norm_0 << std::endl;
00714       }
00715     }
00716 
00717     if ( norm < (_reduction * norm_0)  || norm<1E-30)
00718     {
00719       res.converged = 1;
00720       _prec.post(x);                  // postprocess preconditioner
00721       res.iterations = 0;             // fill statistics
00722       res.reduction = 0;
00723       res.conv_rate  = 0;
00724       res.elapsed = watch.elapsed();
00725       return;
00726     }
00727 
00728     //
00729     // iteration
00730     //
00731     
00732     it = 0;
00733 
00734     while ( true )
00735     {
00736       //
00737       // preprocess, set vecsizes etc.
00738       //
00739         
00740       // rho_new = < rt , r >
00741       rho_new = _sp.dot(rt,r);
00742 
00743       // look if breakdown occured
00744       if (std::abs(rho) <= EPSILON)
00745             DUNE_THROW(ISTLError,"breakdown in BiCGSTAB - rho "
00746                        << rho << " <= EPSILON " << EPSILON
00747                        << " after " << it << " iterations");
00748           if (std::abs(omega) <= EPSILON)
00749             DUNE_THROW(ISTLError,"breakdown in BiCGSTAB - omega "
00750                        << omega << " <= EPSILON " << EPSILON
00751                        << " after " << it << " iterations");
00752 
00753         
00754       if (it==0)
00755       p = r;
00756       else
00757       {
00758         beta = ( rho_new / rho ) * ( alpha / omega );       
00759         p.axpy(-omega,v); // p = r + beta (p - omega*v)
00760         p *= beta;
00761         p += r;
00762       }
00763 
00764       // y = W^-1 * p
00765       y = 0;
00766       _prec.apply(y,p);           // apply preconditioner
00767 
00768       // v = A * y
00769       _op.apply(y,v);
00770 
00771       // alpha = rho_new / < rt, v >
00772       h = _sp.dot(rt,v);
00773 
00774       if ( std::abs(h) < EPSILON )
00775       DUNE_THROW(ISTLError,"h=0 in BiCGSTAB");
00776       
00777       alpha = rho_new / h;
00778         
00779       // apply first correction to x
00780       // x <- x + alpha y
00781       x.axpy(alpha,y);
00782 
00783       // r = r - alpha*v
00784       r.axpy(-alpha,v);
00785         
00786       //
00787       // test stop criteria
00788       //
00789 
00790       it++;
00791       norm = _sp.norm(r);
00792 
00793       if (_verbose>1) // print
00794       {
00795         this->printOutput(std::cout,it,norm,norm_old);
00796       }
00797         
00798       if ( norm < (_reduction * norm_0) )
00799       {
00800         res.converged = 1;
00801         break;
00802       }
00803 
00804       if (it >= _maxit)
00805             break;
00806         
00807       norm_old = norm;
00808 
00809       // y = W^-1 * r
00810       y = 0;
00811       _prec.apply(y,r);
00812 
00813       // t = A * y
00814       _op.apply(y,t);
00815       
00816       // omega = < t, r > / < t, t >
00817       omega = _sp.dot(t,r)/_sp.dot(t,t);
00818         
00819       // apply second correction to x
00820       // x <- x + omega y
00821       x.axpy(omega,y);
00822         
00823       // r = s - omega*t (remember : r = s)
00824       r.axpy(-omega,t);
00825         
00826       rho = rho_new;
00827 
00828       //
00829       // test stop criteria
00830       //
00831 
00832       it++;
00833         
00834       norm = _sp.norm(r);
00835 
00836       if (_verbose > 1)             // print
00837       {
00838         this->printOutput(std::cout,it,norm,norm_old);
00839       }
00840         
00841       if ( norm < (_reduction * norm_0)  || norm<1E-30)
00842       {
00843         res.converged = 1;
00844         break;
00845       }
00846 
00847       if (it >= _maxit)
00848             break;
00849         
00850       norm_old = norm;
00851     }// while
00852 
00853     if (_verbose==1)                // printing for non verbose
00854       this->printOutput(std::cout,it,norm);
00855 
00856     _prec.post(x);                  // postprocess preconditioner
00857     res.iterations = it;              // fill statistics
00858     res.reduction = norm/norm_0;
00859     res.conv_rate  = pow(res.reduction,1.0/it);
00860     res.elapsed = watch.elapsed();
00861     if (_verbose>0)                 // final print 
00862               std::cout << "=== rate=" << res.conv_rate
00863                         << ", T=" << res.elapsed
00864                         << ", TIT=" << res.elapsed/it
00865                         << ", IT=" << it << std::endl;
00866     }
00867 
00873     virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
00874     {
00875       _reduction = reduction;
00876       (*this).apply(x,b,res);
00877     }
00878 
00879   private:
00880     SeqScalarProduct<X> ssp;
00881     LinearOperator<X,X>& _op;
00882     Preconditioner<X,X>& _prec;
00883     ScalarProduct<X>& _sp;
00884     double _reduction;
00885     int _maxit;
00886     int _verbose;   
00887   };
00888 
00889 
00892 } // end namespace
00893 
00894 #endif

Generated on 6 Nov 2008 with Doxygen (ver 1.5.6) [logfile].