00001 #ifndef DUNE_SOLVERS_HH
00002 #define DUNE_SOLVERS_HH
00003
00004 #include<math.h>
00005 #include<complex>
00006 #include<iostream>
00007 #include<iomanip>
00008 #include<string>
00009 #include<stdio.h>
00010
00011 #include "istlexception.hh"
00012 #include "operators.hh"
00013 #include "preconditioners.hh"
00014 #include "scalarproducts.hh"
00015 #include <dune/common/timer.hh>
00016 #include <dune/common/helpertemplates.hh>
00017
00018 namespace Dune {
00043 struct InverseOperatorResult
00044 {
00046 InverseOperatorResult ()
00047 {
00048 clear();
00049 }
00050
00052 void clear ()
00053 {
00054 iterations = 0;
00055 reduction = 0;
00056 converged = false;
00057 conv_rate = 1;
00058 elapsed = 0;
00059 }
00060
00062 int iterations;
00063
00065 double reduction;
00066
00068 bool converged;
00069
00071 double conv_rate;
00072
00074 double elapsed;
00075 };
00076
00077
00078
00090 template<class X, class Y>
00091 class InverseOperator {
00092 public:
00094 typedef X domain_type;
00095
00097 typedef Y range_type;
00098
00100 typedef typename X::field_type field_type;
00101
00111 virtual void apply (X& x, Y& b, InverseOperatorResult& res) = 0;
00112
00123 virtual void apply (X& x, Y& b, double reduction, InverseOperatorResult& res) = 0;
00124
00126 virtual ~InverseOperator () {}
00127 };
00128
00129
00130
00131
00132
00133
00142 template<class X>
00143 class LoopSolver : public InverseOperator<X,X> {
00144 public:
00146 typedef X domain_type;
00148 typedef X range_type;
00150 typedef typename X::field_type field_type;
00151
00171 template<class L, class P>
00172 LoopSolver (L& op, P& prec,
00173 double reduction, int maxit, int verbose) :
00174 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00175 {
00176 IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00177 IsTrue< static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential) >::yes();
00178 }
00179
00200 template<class L, class S, class P>
00201 LoopSolver (L& op, S& sp, P& prec,
00202 double reduction, int maxit, int verbose) :
00203 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00204 {
00205 IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00206 IsTrue< static_cast<int>(L::category) == static_cast<int>(S::category) >::yes();
00207 }
00208
00209
00211 virtual void apply (X& x, X& b, InverseOperatorResult& res)
00212 {
00213
00214 res.clear();
00215
00216
00217 Timer watch;
00218
00219
00220 _prec.pre(x,b);
00221
00222
00223 _op.applyscaleadd(-1,x,b);
00224
00225
00226 double def0 = _sp.norm(b);
00227
00228
00229 if (_verbose>0)
00230 {
00231 printf("=== LoopSolver\n");
00232 if (_verbose>1) printf(" Iter Defect Rate\n");
00233 if (_verbose>1) printf("%5d %12.4E\n",0,def0);
00234 }
00235
00236
00237 X v(x);
00238
00239
00240 int i=1; double def=def0;
00241 for ( ; i<=_maxit; i++ )
00242 {
00243 v = 0;
00244 _prec.apply(v,b);
00245 x += v;
00246 _op.applyscaleadd(-1,v,b);
00247 double defnew=_sp.norm(b);
00248 if (_verbose>1)
00249 printf("%5d %12.4E %12.4g\n",i,defnew,defnew/def);
00250 def = defnew;
00251 if (def<def0*_reduction || def<1E-30)
00252 {
00253 res.converged = true;
00254 break;
00255 }
00256 }
00257
00258
00259 if (_verbose==1)
00260 printf("%5d %12.4E\n",i,def);
00261
00262
00263 _prec.post(x);
00264
00265
00266 res.iterations = i;
00267 res.reduction = def/def0;
00268 res.conv_rate = pow(res.reduction,1.0/i);
00269 res.elapsed = watch.elapsed();
00270
00271
00272 if (_verbose>0)
00273 printf("=== rate=%g, T=%g, TIT=%g, IT=%d\n",res.conv_rate,res.elapsed,res.elapsed/i,i);
00274 }
00275
00277 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
00278 {
00279 _reduction = reduction;
00280 (*this).apply(x,b,res);
00281 }
00282
00283 private:
00284 SeqScalarProduct<X> ssp;
00285 LinearOperator<X,X>& _op;
00286 Preconditioner<X,X>& _prec;
00287 ScalarProduct<X>& _sp;
00288 double _reduction;
00289 int _maxit;
00290 int _verbose;
00291 };
00292
00293
00294
00296 template<class X>
00297 class GradientSolver : public InverseOperator<X,X> {
00298 public:
00300 typedef X domain_type;
00302 typedef X range_type;
00304 typedef typename X::field_type field_type;
00305
00311 template<class L, class P>
00312 GradientSolver (L& op, P& prec,
00313 double reduction, int maxit, int verbose) :
00314 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00315 {
00316 IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00317 IsTrue< static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential) >::yes();
00318 }
00324 template<class L, class S, class P>
00325 GradientSolver (L& op, S& sp, P& prec,
00326 double reduction, int maxit, int verbose) :
00327 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00328 {
00329 IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00330 IsTrue< static_cast<int>(L::category) == static_cast<int>(S::category) >::yes();
00331 }
00332
00338 virtual void apply (X& x, X& b, InverseOperatorResult& res)
00339 {
00340 res.clear();
00341 Timer watch;
00342 _prec.pre(x,b);
00343 _op.applyscaleadd(-1,x,b);
00344
00345 X p(x);
00346 X q(b);
00347
00348 double def0 = _sp.norm(b);
00349
00350 if (_verbose>0)
00351 {
00352 printf("=== GradientSolver\n");
00353 if (_verbose>1) printf(" Iter Defect Rate\n");
00354 if (_verbose>1) printf("%5d %12.4E\n",0,def0);
00355 }
00356
00357 int i=1; double def=def0;
00358 field_type lambda;
00359 for ( ; i<=_maxit; i++ )
00360 {
00361 p = 0;
00362 _prec.apply(p,b);
00363 _op.apply(p,q);
00364 lambda = _sp.dot(p,b)/_sp.dot(q,p);
00365 x.axpy(lambda,p);
00366 b.axpy(-lambda,q);
00367
00368 double defnew=_sp.norm(b);
00369 if (_verbose>1)
00370 printf("%5d %12.4E %12.4g\n",i,defnew,defnew/def);
00371 def = defnew;
00372 if (def<def0*_reduction || def<1E-30)
00373 {
00374 res.converged = true;
00375 break;
00376 }
00377 }
00378
00379 if (_verbose==1)
00380 printf("%5d %12.4E\n",i,def);
00381 _prec.post(x);
00382 res.iterations = i;
00383 res.reduction = def/def0;
00384 res.conv_rate = pow(res.reduction,1.0/i);
00385 res.elapsed = watch.elapsed();
00386 if (_verbose>0)
00387 printf("=== rate=%g, T=%g, TIT=%g, IT=%d\n",res.conv_rate,res.elapsed,res.elapsed/i,i);
00388 }
00389
00395 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
00396 {
00397 _reduction = reduction;
00398 (*this).apply(x,b,res);
00399 }
00400
00401 private:
00402 SeqScalarProduct<X> ssp;
00403 LinearOperator<X,X>& _op;
00404 Preconditioner<X,X>& _prec;
00405 ScalarProduct<X>& _sp;
00406 double _reduction;
00407 int _maxit;
00408 int _verbose;
00409 };
00410
00411
00412
00414 template<class X>
00415 class CGSolver : public InverseOperator<X,X> {
00416 public:
00418 typedef X domain_type;
00420 typedef X range_type;
00422 typedef typename X::field_type field_type;
00423
00429 template<class L, class P>
00430 CGSolver (L& op, P& prec, double reduction, int maxit, int verbose) :
00431 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00432 {
00433 IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00434 IsTrue< static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential) >::yes();
00435 }
00441 template<class L, class S, class P>
00442 CGSolver (L& op, S& sp, P& prec, double reduction, int maxit, int verbose) :
00443 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00444 {
00445 IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00446 IsTrue< static_cast<int>(L::category) == static_cast<int>(S::category) >::yes();
00447 }
00448
00454 virtual void apply (X& x, X& b, InverseOperatorResult& res)
00455 {
00456 res.clear();
00457 Timer watch;
00458 _prec.pre(x,b);
00459 _op.applyscaleadd(-1,x,b);
00460
00461 X p(x);
00462 X q(x);
00463
00464 double def0 = _sp.norm(b);
00465 if (def0<1E-30)
00466 {
00467 res.converged = true;
00468 res.iterations = 0;
00469 res.reduction = 0;
00470 res.conv_rate = 0;
00471 res.elapsed=0;
00472 if (_verbose>0)
00473 printf("=== rate=%g, T=%g, TIT=%g, IT=%d\n",res.conv_rate,res.elapsed,res.elapsed,0);
00474 return;
00475 }
00476
00477 if (_verbose>0)
00478 {
00479 printf("=== CGSolver\n");
00480 if (_verbose>1) printf(" Iter Defect Rate\n");
00481 if (_verbose>1) printf("%5d %12.4E\n",0,def0);
00482 }
00483
00484
00485 double def=def0;
00486 field_type rho,rholast,lambda,alpha,beta;
00487
00488
00489 p = 0;
00490 _prec.apply(p,b);
00491 rholast = _sp.dot(p,b);
00492
00493
00494 int i=1;
00495 for ( ; i<=_maxit; i++ )
00496 {
00497
00498 _op.apply(p,q);
00499 alpha = _sp.dot(p,q);
00500 lambda = rholast/alpha;
00501 x.axpy(lambda,p);
00502 b.axpy(-lambda,q);
00503
00504
00505 double defnew=_sp.norm(b);
00506 if (_verbose>1)
00507 printf("%5d %12.4E %12.4E\n",i,defnew,defnew/def);
00508 def = defnew;
00509 if (def<def0*_reduction || def<1E-30 || i==_maxit)
00510 {
00511 res.converged = true;
00512 break;
00513 }
00514
00515
00516 q = 0;
00517 _prec.apply(q,b);
00518 rho = _sp.dot(q,b);
00519 beta = rho/rholast;
00520 p *= beta;
00521 p += q;
00522 rholast = rho;
00523 }
00524
00525 if (_verbose==1)
00526 printf("%5d %12.4E\n",i,def);
00527 _prec.post(x);
00528 res.iterations = i;
00529 res.reduction = def/def0;
00530 res.conv_rate = pow(res.reduction,1.0/i);
00531 res.elapsed = watch.elapsed();
00532 if (_verbose>0)
00533 printf("=== rate=%g, T=%g, TIT=%g, IT=%d\n",res.conv_rate,res.elapsed,res.elapsed/i,i);
00534 }
00535
00541 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
00542 {
00543 _reduction = reduction;
00544 (*this).apply(x,b,res);
00545 }
00546
00547 private:
00548 SeqScalarProduct<X> ssp;
00549 LinearOperator<X,X>& _op;
00550 Preconditioner<X,X>& _prec;
00551 ScalarProduct<X>& _sp;
00552 double _reduction;
00553 int _maxit;
00554 int _verbose;
00555 };
00556
00557
00558
00560 template<class X>
00561 class BiCGSTABSolver : public InverseOperator<X,X> {
00562 public:
00564 typedef X domain_type;
00566 typedef X range_type;
00568 typedef typename X::field_type field_type;
00569
00575 template<class L, class P>
00576 BiCGSTABSolver (L& op, P& prec,
00577 double reduction, int maxit, int verbose) :
00578 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00579 {
00580 IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00581 IsTrue< static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential) >::yes();
00582 }
00588 template<class L, class S, class P>
00589 BiCGSTABSolver (L& op, S& sp, P& prec,
00590 double reduction, int maxit, int verbose) :
00591 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
00592 {
00593 IsTrue< static_cast<int>(L::category) == static_cast<int>(P::category) >::yes();
00594 IsTrue< static_cast<int>(L::category) == static_cast<int>(S::category) >::yes();
00595 }
00596
00602 virtual void apply (X& x, X& b, InverseOperatorResult& res)
00603 {
00604 const double EPSILON=1e-80;
00605
00606 int it;
00607 field_type rho, rho_new, alpha, beta, h, omega;
00608 field_type norm, norm_old, norm_0;
00609
00610
00611
00612
00613 X& r=b;
00614 X p(x);
00615 X v(x);
00616 X t(x);
00617 X y(x);
00618 X rt(x);
00619
00620
00621
00622
00623
00624
00625 res.clear();
00626 Timer watch;
00627 _op.applyscaleadd(-1,x,r);
00628 _prec.pre(x,r);
00629
00630 rt=r;
00631
00632 norm = norm_old = norm_0 = _sp.norm(r);
00633
00634 p=0;
00635 v=0;
00636
00637 rho = 1;
00638 alpha = 1;
00639 omega = 1;
00640
00641 if (_verbose>0)
00642 {
00643 printf("=== BiCGSTABSolver\n");
00644 if (_verbose>1) printf(" Iter Defect Rate\n");
00645 if (_verbose>1) printf("%5d %12.4E\n",0,norm_0);
00646 }
00647
00648 if ( norm < (_reduction * norm_0) || norm<1E-30)
00649 {
00650 res.converged = 1;
00651 _prec.post(x);
00652 res.iterations = 0;
00653 res.reduction = 0;
00654 res.conv_rate = 0;
00655 res.elapsed = watch.elapsed();
00656 return;
00657 }
00658
00659
00660
00661
00662
00663 it = 0;
00664
00665 while ( true )
00666 {
00667
00668
00669
00670
00671
00672 rho_new = _sp.dot(rt,r);
00673
00674
00675 if (std::abs(rho) <= EPSILON)
00676 DUNE_THROW(ISTLError,"breakdown in BiCGSTAB - rho "
00677 << rho << " <= EPSILON " << EPSILON
00678 << " after " << it << " iterations");
00679 if (std::abs(omega) <= EPSILON)
00680 DUNE_THROW(ISTLError,"breakdown in BiCGSTAB - omega "
00681 << omega << " <= EPSILON " << EPSILON
00682 << " after " << it << " iterations");
00683
00684
00685 if (it==0)
00686 p = r;
00687 else
00688 {
00689 beta = ( rho_new / rho ) * ( alpha / omega );
00690 p.axpy(-omega,v);
00691 p *= beta;
00692 p += r;
00693 }
00694
00695
00696 y = 0;
00697 _prec.apply(y,p);
00698
00699
00700 _op.apply(y,v);
00701
00702
00703 h = _sp.dot(rt,v);
00704
00705 if ( std::abs(h) < EPSILON )
00706 DUNE_THROW(ISTLError,"h=0 in BiCGSTAB");
00707
00708 alpha = rho_new / h;
00709
00710
00711
00712 x.axpy(alpha,y);
00713
00714
00715 r.axpy(-alpha,v);
00716
00717
00718
00719
00720
00721 it++;
00722 norm = _sp.norm(r);
00723
00724 if (_verbose>1)
00725 printf("%5d %12.4E %12.4g\n",it,norm,norm/norm_old);
00726
00727 if ( norm < (_reduction * norm_0) )
00728 {
00729 res.converged = 1;
00730 break;
00731 }
00732
00733 if (it >= _maxit)
00734 break;
00735
00736 norm_old = norm;
00737
00738
00739 y = 0;
00740 _prec.apply(y,r);
00741
00742
00743 _op.apply(y,t);
00744
00745
00746 omega = _sp.dot(t,r)/_sp.dot(t,t);
00747
00748
00749
00750 x.axpy(omega,y);
00751
00752
00753 r.axpy(-omega,t);
00754
00755 rho = rho_new;
00756
00757
00758
00759
00760
00761 it++;
00762
00763 norm = _sp.norm(r);
00764
00765 if (_verbose>1)
00766 printf("%5d %12.4E %12.4g\n",it,norm,norm/norm_old);
00767
00768 if ( norm < (_reduction * norm_0) || norm<1E-30)
00769 {
00770 res.converged = 1;
00771 break;
00772 }
00773
00774 if (it >= _maxit)
00775 break;
00776
00777 norm_old = norm;
00778 }
00779
00780 if (_verbose==1)
00781 printf("%5d %12.4E\n",it,norm);
00782 _prec.post(x);
00783 res.iterations = it;
00784 res.reduction = norm/norm_0;
00785 res.conv_rate = pow(res.reduction,1.0/it);
00786 res.elapsed = watch.elapsed();
00787 if (_verbose>0)
00788 printf("=== rate=%g, T=%g, TIT=%g, IT=%d\n",res.conv_rate,res.elapsed,res.elapsed/it,it);
00789 }
00790
00796 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
00797 {
00798 _reduction = reduction;
00799 (*this).apply(x,b,res);
00800 }
00801
00802 private:
00803 SeqScalarProduct<X> ssp;
00804 LinearOperator<X,X>& _op;
00805 Preconditioner<X,X>& _prec;
00806 ScalarProduct<X>& _sp;
00807 double _reduction;
00808 int _maxit;
00809 int _verbose;
00810 };
00811
00812
00815 }
00816
00817 #endif