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/static_assert.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
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 double 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 double 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
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 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
00211 "L and P have to have the same category!");
00212 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
00213 "L has to be sequential!");
00214 }
00215
00236 template<class L, class S, class P>
00237 LoopSolver (L& op, S& sp, P& prec,
00238 double reduction, int maxit, int verbose) :
00239 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00240 {
00241 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
00242 "L and P must have the same category!");
00243 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category),
00244 "L and S must have the same category!");
00245 }
00246
00247
00249 virtual void apply (X& x, X& b, InverseOperatorResult& res)
00250 {
00251
00252 res.clear();
00253
00254
00255 Timer watch;
00256
00257
00258 _prec.pre(x,b);
00259
00260
00261 _op.applyscaleadd(-1,x,b);
00262
00263
00264 double def0 = _sp.norm(b);
00265
00266
00267 if (_verbose>0)
00268 {
00269 std::cout << "=== LoopSolver" << std::endl;
00270 if (_verbose>1)
00271 {
00272 this->printHeader(std::cout);
00273 this->printOutput(std::cout,0,def0);
00274 }
00275 }
00276
00277
00278 X v(x);
00279
00280
00281 int i=1; double def=def0;
00282 for ( ; i<=_maxit; i++ )
00283 {
00284 v = 0;
00285 _prec.apply(v,b);
00286 x += v;
00287 _op.applyscaleadd(-1,v,b);
00288 double defnew=_sp.norm(b);
00289 if (_verbose>1)
00290 this->printOutput(std::cout,i,defnew,def);
00291
00292 def = defnew;
00293 if (def<def0*_reduction || def<1E-30)
00294 {
00295 res.converged = true;
00296 break;
00297 }
00298 }
00299
00300
00301 if (_verbose==1)
00302 this->printOutput(std::cout,i,def);
00303
00304
00305 _prec.post(x);
00306
00307
00308 res.iterations = i;
00309 res.reduction = def/def0;
00310 res.conv_rate = pow(res.reduction,1.0/i);
00311 res.elapsed = watch.elapsed();
00312
00313
00314 if (_verbose>0)
00315 {
00316 std::cout << "=== rate=" << res.conv_rate
00317 << ", T=" << res.elapsed
00318 << ", TIT=" << res.elapsed/i
00319 << ", IT=" << i << std::endl;
00320 }
00321 }
00322
00324 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
00325 {
00326 _reduction = reduction;
00327 (*this).apply(x,b,res);
00328 }
00329
00330 private:
00331 SeqScalarProduct<X> ssp;
00332 LinearOperator<X,X>& _op;
00333 Preconditioner<X,X>& _prec;
00334 ScalarProduct<X>& _sp;
00335 double _reduction;
00336 int _maxit;
00337 int _verbose;
00338 };
00339
00340
00341
00343 template<class X>
00344 class GradientSolver : public InverseOperator<X,X> {
00345 public:
00347 typedef X domain_type;
00349 typedef X range_type;
00351 typedef typename X::field_type field_type;
00352
00358 template<class L, class P>
00359 GradientSolver (L& op, P& prec,
00360 double reduction, int maxit, int verbose) :
00361 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00362 {
00363 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
00364 "L and P have to have the same category!");
00365 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
00366 "L has to be sequential!");
00367 }
00373 template<class L, class S, class P>
00374 GradientSolver (L& op, S& sp, P& prec,
00375 double reduction, int maxit, int verbose) :
00376 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00377 {
00378 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
00379 "L and P have to have the same category!");
00380 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(S::category),
00381 "L and S have to have the same category!");
00382 }
00383
00389 virtual void apply (X& x, X& b, InverseOperatorResult& res)
00390 {
00391 res.clear();
00392 Timer watch;
00393 _prec.pre(x,b);
00394 _op.applyscaleadd(-1,x,b);
00395
00396 X p(x);
00397 X q(b);
00398
00399 double def0 = _sp.norm(b);
00400
00401 if (_verbose>0)
00402 {
00403 std::cout << "=== GradientSolver" << std::endl;
00404 if (_verbose>1)
00405 {
00406 this->printHeader(std::cout);
00407 this->printOutput(std::cout,0,def0);
00408 }
00409 }
00410
00411 int i=1; double def=def0;
00412 field_type lambda;
00413 for ( ; i<=_maxit; i++ )
00414 {
00415 p = 0;
00416 _prec.apply(p,b);
00417 _op.apply(p,q);
00418 lambda = _sp.dot(p,b)/_sp.dot(q,p);
00419 x.axpy(lambda,p);
00420 b.axpy(-lambda,q);
00421
00422 double defnew=_sp.norm(b);
00423 if (_verbose>1)
00424 this->printOutput(std::cout,i,defnew,def);
00425
00426 def = defnew;
00427 if (def<def0*_reduction || def<1E-30)
00428 {
00429 res.converged = true;
00430 break;
00431 }
00432 }
00433
00434 if (_verbose==1)
00435 this->printOutput(std::cout,i,def);
00436
00437 _prec.post(x);
00438 res.iterations = i;
00439 res.reduction = def/def0;
00440 res.conv_rate = pow(res.reduction,1.0/i);
00441 res.elapsed = watch.elapsed();
00442 if (_verbose>0)
00443 std::cout << "=== rate=" << res.conv_rate
00444 << ", T=" << res.elapsed
00445 << ", TIT=" << res.elapsed/i
00446 << ", IT=" << i << std::endl;
00447 }
00448
00454 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
00455 {
00456 _reduction = reduction;
00457 (*this).apply(x,b,res);
00458 }
00459
00460 private:
00461 SeqScalarProduct<X> ssp;
00462 LinearOperator<X,X>& _op;
00463 Preconditioner<X,X>& _prec;
00464 ScalarProduct<X>& _sp;
00465 double _reduction;
00466 int _maxit;
00467 int _verbose;
00468 };
00469
00470
00471
00473 template<class X>
00474 class CGSolver : public InverseOperator<X,X> {
00475 public:
00477 typedef X domain_type;
00479 typedef X range_type;
00481 typedef typename X::field_type field_type;
00482
00488 template<class L, class P>
00489 CGSolver (L& op, P& prec, double reduction, int maxit, int verbose) :
00490 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00491 {
00492 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
00493 "L and P must have the same category!");
00494 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
00495 "L must be sequential!");
00496 }
00502 template<class L, class S, class P>
00503 CGSolver (L& op, S& sp, P& prec, double reduction, int maxit, int verbose) :
00504 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00505 {
00506 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
00507 "L and P must have the same category!");
00508 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category),
00509 "L and S must have the same category!");
00510 }
00511
00517 virtual void apply (X& x, X& b, InverseOperatorResult& res)
00518 {
00519 res.clear();
00520 Timer watch;
00521 _prec.pre(x,b);
00522 _op.applyscaleadd(-1,x,b);
00523
00524 X p(x);
00525 X q(x);
00526
00527 double def0 = _sp.norm(b);
00528 if (def0<1E-30)
00529 {
00530 res.converged = true;
00531 res.iterations = 0;
00532 res.reduction = 0;
00533 res.conv_rate = 0;
00534 res.elapsed=0;
00535 if (_verbose>0)
00536 std::cout << "=== rate=" << res.conv_rate
00537 << ", T=" << res.elapsed << ", TIT=" << res.elapsed
00538 << ", IT=0" << std::endl;
00539 return;
00540 }
00541
00542 if (_verbose>0)
00543 {
00544 std::cout << "=== CGSolver" << std::endl;
00545 if (_verbose>1) {
00546 this->printHeader(std::cout);
00547 this->printOutput(std::cout,0,def0);
00548 }
00549 }
00550
00551
00552 double def=def0;
00553 field_type rho,rholast,lambda,alpha,beta;
00554
00555
00556 p = 0;
00557 _prec.apply(p,b);
00558 rholast = _sp.dot(p,b);
00559
00560
00561 int i=1;
00562 for ( ; i<=_maxit; i++ )
00563 {
00564
00565 _op.apply(p,q);
00566 alpha = _sp.dot(p,q);
00567 lambda = rholast/alpha;
00568 x.axpy(lambda,p);
00569 b.axpy(-lambda,q);
00570
00571
00572 double defnew=_sp.norm(b);
00573
00574 if (_verbose>1)
00575 this->printOutput(std::cout,i,defnew,def);
00576
00577 def = defnew;
00578 if (def<def0*_reduction || def<1E-30)
00579 {
00580 res.converged = true;
00581 break;
00582 }
00583
00584
00585 q = 0;
00586 _prec.apply(q,b);
00587 rho = _sp.dot(q,b);
00588 beta = rho/rholast;
00589 p *= beta;
00590 p += q;
00591 rholast = rho;
00592 }
00593
00594 if (_verbose==1)
00595 this->printOutput(std::cout,i,def);
00596
00597 _prec.post(x);
00598 res.iterations = i;
00599 res.reduction = def/def0;
00600 res.conv_rate = pow(res.reduction,1.0/i);
00601 res.elapsed = watch.elapsed();
00602
00603 if (_verbose>0)
00604 {
00605 std::cout << "=== rate=" << res.conv_rate
00606 << ", T=" << res.elapsed
00607 << ", TIT=" << res.elapsed/i
00608 << ", IT=" << i << std::endl;
00609 }
00610 }
00611
00617 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
00618 {
00619 _reduction = reduction;
00620 (*this).apply(x,b,res);
00621 }
00622
00623 private:
00624 SeqScalarProduct<X> ssp;
00625 LinearOperator<X,X>& _op;
00626 Preconditioner<X,X>& _prec;
00627 ScalarProduct<X>& _sp;
00628 double _reduction;
00629 int _maxit;
00630 int _verbose;
00631 };
00632
00633
00634
00636 template<class X>
00637 class BiCGSTABSolver : public InverseOperator<X,X> {
00638 public:
00640 typedef X domain_type;
00642 typedef X range_type;
00644 typedef typename X::field_type field_type;
00645
00651 template<class L, class P>
00652 BiCGSTABSolver (L& op, P& prec,
00653 double reduction, int maxit, int verbose) :
00654 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00655 {
00656 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category), "L and P must be of the same category!");
00657 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential), "L must be sequential!");
00658 }
00664 template<class L, class S, class P>
00665 BiCGSTABSolver (L& op, S& sp, P& prec,
00666 double reduction, int maxit, int verbose) :
00667 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00668 {
00669 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
00670 "L and P must have the same category!");
00671 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category),
00672 "L and S must have the same category!");
00673 }
00674
00680 virtual void apply (X& x, X& b, InverseOperatorResult& res)
00681 {
00682 const double EPSILON=1e-80;
00683
00684 double it;
00685 field_type rho, rho_new, alpha, beta, h, omega;
00686 field_type norm, norm_old, norm_0;
00687
00688
00689
00690
00691 X& r=b;
00692 X p(x);
00693 X v(x);
00694 X t(x);
00695 X y(x);
00696 X rt(x);
00697
00698
00699
00700
00701
00702
00703 res.clear();
00704 Timer watch;
00705 _op.applyscaleadd(-1,x,r);
00706 _prec.pre(x,r);
00707
00708 rt=r;
00709
00710 norm = norm_old = norm_0 = _sp.norm(r);
00711
00712 p=0;
00713 v=0;
00714
00715 rho = 1;
00716 alpha = 1;
00717 omega = 1;
00718
00719 if (_verbose>0)
00720 {
00721 std::cout << "=== BiCGSTABSolver" << std::endl;
00722 if (_verbose>1)
00723 {
00724 this->printHeader(std::cout);
00725 this->printOutput(std::cout,0,norm_0);
00726
00727
00728 }
00729 }
00730
00731 if ( norm < (_reduction * norm_0) || norm<1E-30)
00732 {
00733 res.converged = 1;
00734 _prec.post(x);
00735 res.iterations = 0;
00736 res.reduction = 0;
00737 res.conv_rate = 0;
00738 res.elapsed = watch.elapsed();
00739 return;
00740 }
00741
00742
00743
00744
00745
00746 for (it = 0.5; it < _maxit; it+=.5)
00747 {
00748
00749
00750
00751
00752
00753 rho_new = _sp.dot(rt,r);
00754
00755
00756 if (std::abs(rho) <= EPSILON)
00757 DUNE_THROW(ISTLError,"breakdown in BiCGSTAB - rho "
00758 << rho << " <= EPSILON " << EPSILON
00759 << " after " << it << " iterations");
00760 if (std::abs(omega) <= EPSILON)
00761 DUNE_THROW(ISTLError,"breakdown in BiCGSTAB - omega "
00762 << omega << " <= EPSILON " << EPSILON
00763 << " after " << it << " iterations");
00764
00765
00766 if (it<1)
00767 p = r;
00768 else
00769 {
00770 beta = ( rho_new / rho ) * ( alpha / omega );
00771 p.axpy(-omega,v);
00772 p *= beta;
00773 p += r;
00774 }
00775
00776
00777 y = 0;
00778 _prec.apply(y,p);
00779
00780
00781 _op.apply(y,v);
00782
00783
00784 h = _sp.dot(rt,v);
00785
00786 if ( std::abs(h) < EPSILON )
00787 DUNE_THROW(ISTLError,"h=0 in BiCGSTAB");
00788
00789 alpha = rho_new / h;
00790
00791
00792
00793 x.axpy(alpha,y);
00794
00795
00796 r.axpy(-alpha,v);
00797
00798
00799
00800
00801
00802 norm = _sp.norm(r);
00803
00804 if (_verbose>1)
00805 {
00806 this->printOutput(std::cout,it,norm,norm_old);
00807 }
00808
00809 if ( norm < (_reduction * norm_0) )
00810 {
00811 res.converged = 1;
00812 break;
00813 }
00814 it+=.5;
00815
00816 norm_old = norm;
00817
00818
00819 y = 0;
00820 _prec.apply(y,r);
00821
00822
00823 _op.apply(y,t);
00824
00825
00826 omega = _sp.dot(t,r)/_sp.dot(t,t);
00827
00828
00829
00830 x.axpy(omega,y);
00831
00832
00833 r.axpy(-omega,t);
00834
00835 rho = rho_new;
00836
00837
00838
00839
00840
00841 norm = _sp.norm(r);
00842
00843 if (_verbose > 1)
00844 {
00845 this->printOutput(std::cout,it,norm,norm_old);
00846 }
00847
00848 if ( norm < (_reduction * norm_0) || norm<1E-30)
00849 {
00850 res.converged = 1;
00851 break;
00852 }
00853
00854 norm_old = norm;
00855 }
00856
00857 if (_verbose==1)
00858 this->printOutput(std::cout,it,norm);
00859
00860 _prec.post(x);
00861 res.iterations = static_cast<int>(std::ceil(it));
00862 res.reduction = norm/norm_0;
00863 res.conv_rate = pow(res.reduction,1.0/it);
00864 res.elapsed = watch.elapsed();
00865 if (_verbose>0)
00866 std::cout << "=== rate=" << res.conv_rate
00867 << ", T=" << res.elapsed
00868 << ", TIT=" << res.elapsed/it
00869 << ", IT=" << it << std::endl;
00870 }
00871
00877 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
00878 {
00879 _reduction = reduction;
00880 (*this).apply(x,b,res);
00881 }
00882
00883 private:
00884 SeqScalarProduct<X> ssp;
00885 LinearOperator<X,X>& _op;
00886 Preconditioner<X,X>& _prec;
00887 ScalarProduct<X>& _sp;
00888 double _reduction;
00889 int _maxit;
00890 int _verbose;
00891 };
00892
00899 template<class X>
00900 class MINRESSolver : public InverseOperator<X,X> {
00901 public:
00903 typedef X domain_type;
00905 typedef X range_type;
00907 typedef typename X::field_type field_type;
00908
00914 template<class L, class P>
00915 MINRESSolver (L& op, P& prec, double reduction, int maxit, int verbose) :
00916 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00917 {
00918 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
00919 "L and P must have the same category!");
00920 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
00921 "L must be sequential!");
00922 }
00928 template<class L, class S, class P>
00929 MINRESSolver (L& op, S& sp, P& prec, double reduction, int maxit, int verbose) :
00930 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00931 {
00932 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
00933 "L and P must have the same category!");
00934 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category),
00935 "L and S must have the same category!");
00936 }
00937
00943 virtual void apply (X& x, X& b, InverseOperatorResult& res)
00944 {
00945 res.clear();
00946 Timer watch;
00947 _prec.pre(x,b);
00948 _op.applyscaleadd(-1,x,b);
00949
00950 double def0 = _sp.norm(b);
00951
00952 if (def0<1E-30)
00953 {
00954 res.converged = true;
00955 res.iterations = 0;
00956 res.reduction = 0;
00957 res.conv_rate = 0;
00958 res.elapsed=0;
00959 if (_verbose>0)
00960 std::cout << "=== rate=" << res.conv_rate << ", T=" << res.elapsed << ", TIT=" << res.elapsed << ", IT=0" << std::endl;
00961 return;
00962 }
00963
00964 if (_verbose>0)
00965 {
00966 std::cout << "=== MINRESSolver" << std::endl;
00967 if (_verbose>1) {
00968 this->printHeader(std::cout);
00969 this->printOutput(std::cout,0,def0);
00970 }
00971 }
00972
00973
00974 double def=def0;
00975 field_type alpha,
00976 beta,
00977 c[2]={0.0, 0.0},
00978 s[2]={0.0, 0.0};
00979
00980 field_type T[3]={0.0, 0.0, 0.0};
00981
00982 X z(b.size()),
00983 dummy(b.size());
00984
00985 field_type xi[2]={1.0, 0.0};
00986
00987
00988 z = 0.0;
00989
00990 _prec.apply(z,b);
00991
00992 beta = sqrt(fabs(_sp.dot(z,b)));
00993 double beta0 = beta;
00994
00995 X p[3];
00996 X q[3];
00997
00998 q[0].resize(b.size());
00999 q[1].resize(b.size());
01000 q[2].resize(b.size());
01001 q[0] = 0.0;
01002 q[1] = b;
01003 q[1] /= beta;
01004 q[2] = 0.0;
01005
01006 p[0].resize(b.size());
01007 p[1].resize(b.size());
01008 p[2].resize(b.size());
01009 p[0] = 0.0;
01010 p[1] = 0.0;
01011 p[2] = 0.0;
01012
01013
01014 z /= beta;
01015
01016
01017 int i=1;
01018 for ( ; i<=_maxit; i++)
01019 {
01020 dummy = z;
01021
01022 int i1 = i%3,
01023 i0 = (i1+2)%3,
01024 i2 = (i1+1)%3;
01025
01026
01027 _op.apply(z,q[i2]);
01028 q[i2].axpy(-beta, q[i0]);
01029 alpha = _sp.dot(q[i2],z);
01030 q[i2].axpy(-alpha, q[i1]);
01031
01032 z=0.0;
01033 _prec.apply(z,q[i2]);
01034
01035 beta = sqrt(fabs(_sp.dot(q[i2],z)));
01036
01037 q[i2] /= beta;
01038 z /= beta;
01039
01040
01041
01042 T[1] = T[2];
01043 if (i>2)
01044 {
01045 T[0] = s[i%2]*T[1];
01046 T[1] = c[i%2]*T[1];
01047 }
01048 if (i>1)
01049 {
01050 T[2] = c[(i+1)%2]*alpha - s[(i+1)%2]*T[1];
01051 T[1] = c[(i+1)%2]*T[1] + s[(i+1)%2]*alpha;
01052 }
01053 else
01054 T[2] = alpha;
01055
01056
01057
01058 c[i%2] = 1.0/sqrt(T[2]*T[2] + beta*beta);
01059 s[i%2] = beta*c[i%2];
01060 c[i%2] *= T[2];
01061
01062
01063 T[2] = c[i%2]*T[2] + s[i%2]*beta;
01064
01065
01066 xi[i%2] = -s[i%2]*xi[(i+1)%2];
01067 xi[(i+1)%2] *= c[i%2];
01068
01069
01070 p[i2] = dummy;
01071 p[i2].axpy(-T[1],p[i1]);
01072 p[i2].axpy(-T[0],p[i0]);
01073 p[i2] /= T[2];
01074
01075
01076 x.axpy(beta0*xi[(i+1)%2], p[i2]);
01077
01078
01079 T[2] = beta;
01080
01081
01082
01083
01084
01085
01086
01087
01088 double defnew = fabs(beta0*xi[i%2]);
01089
01090 if (_verbose>1)
01091 this->printOutput(std::cout,i,defnew,def);
01092
01093 def = defnew;
01094 if (def<def0*_reduction || def<1E-30 || i==_maxit)
01095 {
01096 res.converged = true;
01097 break;
01098 }
01099 }
01100
01101 if (_verbose==1)
01102 this->printOutput(std::cout,i,def);
01103
01104 _prec.post(x);
01105 res.iterations = i;
01106 res.reduction = def/def0;
01107 res.conv_rate = pow(res.reduction,1.0/i);
01108 res.elapsed = watch.elapsed();
01109
01110 if (_verbose>0)
01111 {
01112 std::cout << "=== rate=" << res.conv_rate
01113 << ", T=" << res.elapsed
01114 << ", TIT=" << res.elapsed/i
01115 << ", IT=" << i << std::endl;
01116 }
01117
01118 }
01119
01125 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
01126 {
01127 _reduction = reduction;
01128 (*this).apply(x,b,res);
01129 }
01130
01131 private:
01132 SeqScalarProduct<X> ssp;
01133 LinearOperator<X,X>& _op;
01134 Preconditioner<X,X>& _prec;
01135 ScalarProduct<X>& _sp;
01136 double _reduction;
01137 int _maxit;
01138 int _verbose;
01139 };
01140
01152 template<class X, class Y=X, class F = Y>
01153 class RestartedGMResSolver : public InverseOperator<X,Y>
01154 {
01155 public:
01157 typedef X domain_type;
01159 typedef Y range_type;
01161 typedef typename X::field_type field_type;
01163 typedef F basis_type;
01164
01171 template<class L, class P>
01172 RestartedGMResSolver (L& A, P& M, double reduction, int restart, int maxit, int verbose, bool recalc_defect = false) :
01173 _A_(A), _M(M),
01174 ssp(), _sp(ssp), _restart(restart),
01175 _reduction(reduction), _maxit(maxit), _verbose(verbose),
01176 _recalc_defect(recalc_defect)
01177 {
01178 dune_static_assert(static_cast<int>(P::category) == static_cast<int>(SolverCategory::sequential),
01179 "P must be sequential!");
01180 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
01181 "L must be sequential!");
01182 }
01183
01190 template<class L, class S, class P>
01191 RestartedGMResSolver (L& A, S& sp, P& M, double reduction, int restart, int maxit, int verbose, bool recalc_defect = false) :
01192 _A_(A), _M(M),
01193 _sp(sp), _restart(restart),
01194 _reduction(reduction), _maxit(maxit), _verbose(verbose),
01195 _recalc_defect(recalc_defect)
01196 {
01197 dune_static_assert(static_cast<int>(P::category) == static_cast<int>(SolverCategory::sequential),
01198 "P must be sequential!");
01199 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
01200 "L must be sequential!");
01201 }
01202
01204 virtual void apply (X& x, X& b, InverseOperatorResult& res)
01205 {
01206 apply(x,b,_reduction,res);
01207 };
01208
01214 virtual void apply (X& x, Y& b, double reduction, InverseOperatorResult& res)
01215 {
01216 int m = _restart;
01217 field_type norm;
01218 field_type norm_old = 0.0;
01219 field_type norm_0;
01220 field_type beta;
01221 int i, j = 1, k;
01222 std::vector<field_type> s(m+1), cs(m), sn(m);
01223
01224 X w(b.size());
01225 std::vector< std::vector<field_type> > H(m+1,s);
01226 std::vector<F> v(m+1,b);
01227
01228
01229 Timer watch;
01230
01231
01232 res.clear();
01233
01234 _M.pre(x,b);
01235
01236 if (_recalc_defect)
01237 {
01238
01239 w = 0.0; _M.apply(w,b);
01240 norm_0 = _sp.norm(w);
01241
01242 w = b;
01243 _A_.applyscaleadd(-1,x, w);
01244 v[0] = 0.0; _M.apply(v[0],w);
01245 beta = _sp.norm(v[0]);
01246 }
01247 else
01248 {
01249
01250 w = 0.0; _M.apply(w,b);
01251 norm_0 = _sp.norm(w);
01252
01253 _A_.applyscaleadd(-1,x, b);
01254 v[0] = 0.0; _M.apply(v[0],b);
01255 beta = _sp.norm(v[0]);
01256 }
01257
01258
01259 if (norm_0 == 0.0)
01260 norm_0 = 1.0;
01261
01262 norm = norm_old = _sp.norm(v[0]);
01263
01264
01265 if (_verbose > 0)
01266 {
01267 std::cout << "=== RestartedGMResSolver" << std::endl;
01268 if (_verbose > 1)
01269 {
01270 this->printHeader(std::cout);
01271 this->printOutput(std::cout,0,norm);
01272 }
01273 }
01274
01275
01276 if (norm <= _reduction * norm_0) {
01277 _M.post(x);
01278 res.converged = true;
01279 if (_verbose > 0)
01280 print_result(res);
01281 return;
01282 }
01283
01284 while (j <= _maxit && res.converged != true) {
01285 v[0] *= (1.0 / beta);
01286 for (i=1; i<=m; i++) s[i] = 0.0;
01287 s[0] = beta;
01288
01289 for (i = 0; i < m && j <= _maxit && res.converged != true; i++, j++) {
01290 w = 0.0;
01291 v[i+1] = 0.0;
01292 _A_.apply(v[i], v[i+1]);
01293 _M.apply(w, v[i+1]);
01294 for (k = 0; k <= i; k++) {
01295 H[k][i] = _sp.dot(w, v[k]);
01296
01297 w.axpy(-H[k][i], v[k]);
01298 }
01299 H[i+1][i] = _sp.norm(w);
01300 if (H[i+1][i] == 0.0)
01301 DUNE_THROW(ISTLError,"breakdown in GMRes - |w| "
01302 << w << " == 0.0 after " << j << " iterations");
01303
01304 v[i+1] = w; v[i+1] *= (1.0 / H[i+1][i]);
01305
01306 for (k = 0; k < i; k++)
01307 applyPlaneRotation(H[k][i], H[k+1][i], cs[k], sn[k]);
01308
01309 generatePlaneRotation(H[i][i], H[i+1][i], cs[i], sn[i]);
01310 applyPlaneRotation(H[i][i], H[i+1][i], cs[i], sn[i]);
01311 applyPlaneRotation(s[i], s[i+1], cs[i], sn[i]);
01312
01313 norm = std::abs(s[i+1]);
01314
01315 if (_verbose > 1)
01316 {
01317 this->printOutput(std::cout,j,norm,norm_old);
01318 }
01319
01320 norm_old = norm;
01321
01322 if (norm < _reduction * norm_0) {
01323 res.converged = true;
01324 }
01325 }
01326
01327 if (_recalc_defect)
01328 {
01329
01330 update(x, i - 1, H, s, v);
01331
01332
01333
01334 w = b; _A_.applyscaleadd(-1,x, w);
01335 _M.apply(v[0], w);
01336 beta = _sp.norm(v[0]);
01337 norm = beta;
01338 }
01339 else
01340 {
01341
01342 w = 0;
01343 update(w, i - 1, H, s, v);
01344
01345
01346 x += w;
01347
01348
01349
01350 _A_.applyscaleadd(-1,w, b);
01351
01352 v[0] = 0.0; _M.apply(v[0],b);
01353 beta = _sp.norm(v[0]);
01354 norm = beta;
01355
01356 res.converged = false;
01357 }
01358
01359 if (_verbose > 1)
01360 {
01361 this->printOutput(std::cout,j,norm,norm_old);
01362 }
01363
01364 norm_old = norm;
01365
01366 if (norm < _reduction * norm_0) {
01367
01368 res.converged = true;
01369 }
01370
01371 if (res.converged != true && _verbose > 0)
01372 std::cout << "=== GMRes::restart\n";
01373 }
01374
01375 _M.post(x);
01376
01377 res.iterations = j;
01378 res.reduction = norm / norm_0;
01379 res.conv_rate = pow(res.reduction,1.0/j);
01380 res.elapsed = watch.elapsed();
01381
01382 if (_verbose>0)
01383 print_result(res);
01384 }
01385 private:
01386
01387 void
01388 print_result (const InverseOperatorResult & res) const
01389 {
01390 int j = res.iterations>0?res.iterations:1;
01391 std::cout << "=== rate=" << res.conv_rate
01392 << ", T=" << res.elapsed
01393 << ", TIT=" << res.elapsed/j
01394 << ", IT=" << res.iterations
01395 << std::endl;
01396 }
01397
01398 static void
01399 update(X &x, int k,
01400 std::vector< std::vector<field_type> > & h,
01401 std::vector<field_type> & s, std::vector<F> v)
01402 {
01403 std::vector<field_type> y(s);
01404
01405
01406 for (int i = k; i >= 0; i--) {
01407 y[i] /= h[i][i];
01408 for (int j = i - 1; j >= 0; j--)
01409 y[j] -= h[j][i] * y[i];
01410 }
01411
01412 for (int j = 0; j <= k; j++)
01413
01414 x.axpy(y[j],v[j]);
01415 }
01416
01417 void
01418 generatePlaneRotation(field_type &dx, field_type &dy, field_type &cs, field_type &sn)
01419 {
01420 if (dy == 0.0) {
01421 cs = 1.0;
01422 sn = 0.0;
01423 } else if (std::abs(dy) > std::abs(dx)) {
01424 field_type temp = dx / dy;
01425 sn = 1.0 / std::sqrt( 1.0 + temp*temp );
01426 cs = temp * sn;
01427 } else {
01428 field_type temp = dy / dx;
01429 cs = 1.0 / std::sqrt( 1.0 + temp*temp );
01430 sn = temp * cs;
01431 }
01432 }
01433
01434
01435 void
01436 applyPlaneRotation(field_type &dx, field_type &dy, field_type &cs, field_type &sn)
01437 {
01438 field_type temp = cs * dx + sn * dy;
01439 dy = -sn * dx + cs * dy;
01440 dx = temp;
01441 }
01442
01443 LinearOperator<X,X>& _A_;
01444 Preconditioner<X,X>& _M;
01445 SeqScalarProduct<X> ssp;
01446 ScalarProduct<X>& _sp;
01447 int _restart;
01448 double _reduction;
01449 int _maxit;
01450 int _verbose;
01451 bool _recalc_defect;
01452 };
01453
01456 }
01457
01458 #endif