- Home
- About DUNE
- Download
- Documentation
- Community
- Development
00001 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*- 00002 // vi: set ts=4 sw=2 et sts=2: 00003 00004 #ifndef DUNE_SOLVERS_HH 00005 #define DUNE_SOLVERS_HH 00006 00007 #include<cmath> 00008 #include<complex> 00009 #include<iostream> 00010 #include<iomanip> 00011 #include<string> 00012 00013 #include "istlexception.hh" 00014 #include "operators.hh" 00015 #include "preconditioners.hh" 00016 #include "scalarproducts.hh" 00017 #include <dune/common/timer.hh> 00018 #include <dune/common/static_assert.hh> 00019 00020 namespace Dune { 00045 struct InverseOperatorResult 00046 { 00048 InverseOperatorResult () 00049 { 00050 clear(); 00051 } 00052 00054 void clear () 00055 { 00056 iterations = 0; 00057 reduction = 0; 00058 converged = false; 00059 conv_rate = 1; 00060 elapsed = 0; 00061 } 00062 00064 int iterations; 00065 00067 double reduction; 00068 00070 bool converged; 00071 00073 double conv_rate; 00074 00076 double elapsed; 00077 }; 00078 00079 00080 //===================================================================== 00092 template<class X, class Y> 00093 class InverseOperator { 00094 public: 00096 typedef X domain_type; 00097 00099 typedef Y range_type; 00100 00102 typedef typename X::field_type field_type; 00103 00113 virtual void apply (X& x, Y& b, InverseOperatorResult& res) = 0; 00114 00125 virtual void apply (X& x, Y& b, double reduction, InverseOperatorResult& res) = 0; 00126 00128 virtual ~InverseOperator () {} 00129 00130 protected: 00131 // spacing values 00132 enum { iterationSpacing = 5 , normSpacing = 16 }; 00133 00135 void printHeader(std::ostream& s) const 00136 { 00137 s << std::setw(iterationSpacing) << " Iter"; 00138 s << std::setw(normSpacing) << "Defect"; 00139 s << std::setw(normSpacing) << "Rate" << std::endl; 00140 } 00141 00143 template <class DataType> 00144 void printOutput(std::ostream& s, 00145 const double iter, 00146 const DataType& norm, 00147 const DataType& norm_old) const 00148 { 00149 const DataType rate = norm/norm_old; 00150 s << std::setw(iterationSpacing) << iter << " "; 00151 s << std::setw(normSpacing) << norm << " "; 00152 s << std::setw(normSpacing) << rate << std::endl; 00153 } 00154 00156 template <class DataType> 00157 void printOutput(std::ostream& s, 00158 const double iter, 00159 const DataType& norm) const 00160 { 00161 s << std::setw(iterationSpacing) << iter << " "; 00162 s << std::setw(normSpacing) << norm << std::endl; 00163 } 00164 }; 00165 00166 00167 //===================================================================== 00168 // Implementation of this interface 00169 //===================================================================== 00170 00179 template<class X> 00180 class LoopSolver : public InverseOperator<X,X> { 00181 public: 00183 typedef X domain_type; 00185 typedef X range_type; 00187 typedef typename X::field_type field_type; 00188 00208 template<class L, class P> 00209 LoopSolver (L& op, P& prec, 00210 double reduction, int maxit, int verbose) : 00211 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose) 00212 { 00213 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category), 00214 "L and P have to have the same category!"); 00215 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential), 00216 "L has to be sequential!"); 00217 } 00218 00239 template<class L, class S, class P> 00240 LoopSolver (L& op, S& sp, P& prec, 00241 double reduction, int maxit, int verbose) : 00242 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose) 00243 { 00244 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category), 00245 "L and P must have the same category!"); 00246 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category), 00247 "L and S must have the same category!"); 00248 } 00249 00250 00252 virtual void apply (X& x, X& b, InverseOperatorResult& res) 00253 { 00254 // clear solver statistics 00255 res.clear(); 00256 00257 // start a timer 00258 Timer watch; 00259 00260 // prepare preconditioner 00261 _prec.pre(x,b); 00262 00263 // overwrite b with defect 00264 _op.applyscaleadd(-1,x,b); 00265 00266 // compute norm, \todo parallelization 00267 double def0 = _sp.norm(b); 00268 00269 // printing 00270 if (_verbose>0) 00271 { 00272 std::cout << "=== LoopSolver" << std::endl; 00273 if (_verbose>1) 00274 { 00275 this->printHeader(std::cout); 00276 this->printOutput(std::cout,0,def0); 00277 } 00278 } 00279 00280 // allocate correction vector 00281 X v(x); 00282 00283 // iteration loop 00284 int i=1; double def=def0; 00285 for ( ; i<=_maxit; i++ ) 00286 { 00287 v = 0; // clear correction 00288 _prec.apply(v,b); // apply preconditioner 00289 x += v; // update solution 00290 _op.applyscaleadd(-1,v,b); // update defect 00291 double defnew=_sp.norm(b); // comp defect norm 00292 if (_verbose>1) // print 00293 this->printOutput(std::cout,i,defnew,def); 00294 //std::cout << i << " " << defnew << " " << defnew/def << std::endl; 00295 def = defnew; // update norm 00296 if (def<def0*_reduction || def<1E-30) // convergence check 00297 { 00298 res.converged = true; 00299 break; 00300 } 00301 } 00302 00303 // print 00304 if (_verbose==1) 00305 this->printOutput(std::cout,i,def); 00306 00307 // postprocess preconditioner 00308 _prec.post(x); 00309 00310 // fill statistics 00311 res.iterations = i; 00312 res.reduction = def/def0; 00313 res.conv_rate = pow(res.reduction,1.0/i); 00314 res.elapsed = watch.elapsed(); 00315 00316 // final print 00317 if (_verbose>0) 00318 { 00319 std::cout << "=== rate=" << res.conv_rate 00320 << ", T=" << res.elapsed 00321 << ", TIT=" << res.elapsed/i 00322 << ", IT=" << i << std::endl; 00323 } 00324 } 00325 00327 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res) 00328 { 00329 _reduction = reduction; 00330 (*this).apply(x,b,res); 00331 } 00332 00333 private: 00334 SeqScalarProduct<X> ssp; 00335 LinearOperator<X,X>& _op; 00336 Preconditioner<X,X>& _prec; 00337 ScalarProduct<X>& _sp; 00338 double _reduction; 00339 int _maxit; 00340 int _verbose; 00341 }; 00342 00343 00344 // all these solvers are taken from the SUMO library 00346 template<class X> 00347 class GradientSolver : public InverseOperator<X,X> { 00348 public: 00350 typedef X domain_type; 00352 typedef X range_type; 00354 typedef typename X::field_type field_type; 00355 00361 template<class L, class P> 00362 GradientSolver (L& op, P& prec, 00363 double reduction, int maxit, int verbose) : 00364 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose) 00365 { 00366 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category), 00367 "L and P have to have the same category!"); 00368 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential), 00369 "L has to be sequential!"); 00370 } 00376 template<class L, class S, class P> 00377 GradientSolver (L& op, S& sp, P& prec, 00378 double reduction, int maxit, int verbose) : 00379 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose) 00380 { 00381 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category), 00382 "L and P have to have the same category!"); 00383 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(S::category), 00384 "L and S have to have the same category!"); 00385 } 00386 00392 virtual void apply (X& x, X& b, InverseOperatorResult& res) 00393 { 00394 res.clear(); // clear solver statistics 00395 Timer watch; // start a timer 00396 _prec.pre(x,b); // prepare preconditioner 00397 _op.applyscaleadd(-1,x,b); // overwrite b with defect 00398 00399 X p(x); // create local vectors 00400 X q(b); 00401 00402 double def0 = _sp.norm(b);// compute norm 00403 00404 if (_verbose>0) // printing 00405 { 00406 std::cout << "=== GradientSolver" << std::endl; 00407 if (_verbose>1) 00408 { 00409 this->printHeader(std::cout); 00410 this->printOutput(std::cout,0,def0); 00411 } 00412 } 00413 00414 int i=1; double def=def0; // loop variables 00415 field_type lambda; 00416 for ( ; i<=_maxit; i++ ) 00417 { 00418 p = 0; // clear correction 00419 _prec.apply(p,b); // apply preconditioner 00420 _op.apply(p,q); // q=Ap 00421 lambda = _sp.dot(p,b)/_sp.dot(q,p);// minimization 00422 x.axpy(lambda,p); // update solution 00423 b.axpy(-lambda,q); // update defect 00424 00425 double defnew=_sp.norm(b);// comp defect norm 00426 if (_verbose>1) // print 00427 this->printOutput(std::cout,i,defnew,def); 00428 00429 def = defnew; // update norm 00430 if (def<def0*_reduction || def<1E-30) // convergence check 00431 { 00432 res.converged = true; 00433 break; 00434 } 00435 } 00436 00437 if (_verbose==1) // printing for non verbose 00438 this->printOutput(std::cout,i,def); 00439 00440 _prec.post(x); // postprocess preconditioner 00441 res.iterations = i; // fill statistics 00442 res.reduction = def/def0; 00443 res.conv_rate = pow(res.reduction,1.0/i); 00444 res.elapsed = watch.elapsed(); 00445 if (_verbose>0) // final print 00446 std::cout << "=== rate=" << res.conv_rate 00447 << ", T=" << res.elapsed 00448 << ", TIT=" << res.elapsed/i 00449 << ", IT=" << i << std::endl; 00450 } 00451 00457 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res) 00458 { 00459 _reduction = reduction; 00460 (*this).apply(x,b,res); 00461 } 00462 00463 private: 00464 SeqScalarProduct<X> ssp; 00465 LinearOperator<X,X>& _op; 00466 Preconditioner<X,X>& _prec; 00467 ScalarProduct<X>& _sp; 00468 double _reduction; 00469 int _maxit; 00470 int _verbose; 00471 }; 00472 00473 00474 00476 template<class X> 00477 class CGSolver : public InverseOperator<X,X> { 00478 public: 00480 typedef X domain_type; 00482 typedef X range_type; 00484 typedef typename X::field_type field_type; 00485 00491 template<class L, class P> 00492 CGSolver (L& op, P& prec, double reduction, int maxit, int verbose) : 00493 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose) 00494 { 00495 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category), 00496 "L and P must have the same category!"); 00497 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential), 00498 "L must be sequential!"); 00499 } 00505 template<class L, class S, class P> 00506 CGSolver (L& op, S& sp, P& prec, double reduction, int maxit, int verbose) : 00507 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose) 00508 { 00509 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category), 00510 "L and P must have the same category!"); 00511 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category), 00512 "L and S must have the same category!"); 00513 } 00514 00520 virtual void apply (X& x, X& b, InverseOperatorResult& res) 00521 { 00522 res.clear(); // clear solver statistics 00523 Timer watch; // start a timer 00524 _prec.pre(x,b); // prepare preconditioner 00525 _op.applyscaleadd(-1,x,b); // overwrite b with defect 00526 00527 X p(x); // the search direction 00528 X q(x); // a temporary vector 00529 00530 double def0 = _sp.norm(b);// compute norm 00531 if (def0<1E-30) // convergence check 00532 { 00533 res.converged = true; 00534 res.iterations = 0; // fill statistics 00535 res.reduction = 0; 00536 res.conv_rate = 0; 00537 res.elapsed=0; 00538 if (_verbose>0) // final print 00539 std::cout << "=== rate=" << res.conv_rate 00540 << ", T=" << res.elapsed << ", TIT=" << res.elapsed 00541 << ", IT=0" << std::endl; 00542 return; 00543 } 00544 00545 if (_verbose>0) // printing 00546 { 00547 std::cout << "=== CGSolver" << std::endl; 00548 if (_verbose>1) { 00549 this->printHeader(std::cout); 00550 this->printOutput(std::cout,0,def0); 00551 } 00552 } 00553 00554 // some local variables 00555 double def=def0; // loop variables 00556 field_type rho,rholast,lambda,alpha,beta; 00557 00558 // determine initial search direction 00559 p = 0; // clear correction 00560 _prec.apply(p,b); // apply preconditioner 00561 rholast = _sp.dot(p,b); // orthogonalization 00562 00563 // the loop 00564 int i=1; 00565 for ( ; i<=_maxit; i++ ) 00566 { 00567 // minimize in given search direction p 00568 _op.apply(p,q); // q=Ap 00569 alpha = _sp.dot(p,q); // scalar product 00570 lambda = rholast/alpha; // minimization 00571 x.axpy(lambda,p); // update solution 00572 b.axpy(-lambda,q); // update defect 00573 00574 // convergence test 00575 double defnew=_sp.norm(b);// comp defect norm 00576 00577 if (_verbose>1) // print 00578 this->printOutput(std::cout,i,defnew,def); 00579 00580 def = defnew; // update norm 00581 if (def<def0*_reduction || def<1E-30) // convergence check 00582 { 00583 res.converged = true; 00584 break; 00585 } 00586 00587 // determine new search direction 00588 q = 0; // clear correction 00589 _prec.apply(q,b); // apply preconditioner 00590 rho = _sp.dot(q,b); // orthogonalization 00591 beta = rho/rholast; // scaling factor 00592 p *= beta; // scale old search direction 00593 p += q; // orthogonalization with correction 00594 rholast = rho; // remember rho for recurrence 00595 } 00596 00597 if (_verbose==1) // printing for non verbose 00598 this->printOutput(std::cout,i,def); 00599 00600 _prec.post(x); // postprocess preconditioner 00601 res.iterations = i; // fill statistics 00602 res.reduction = def/def0; 00603 res.conv_rate = pow(res.reduction,1.0/i); 00604 res.elapsed = watch.elapsed(); 00605 00606 if (_verbose>0) // final print 00607 { 00608 std::cout << "=== rate=" << res.conv_rate 00609 << ", T=" << res.elapsed 00610 << ", TIT=" << res.elapsed/i 00611 << ", IT=" << i << std::endl; 00612 } 00613 } 00614 00620 virtual void apply (X& x, X& b, double reduction, 00621 InverseOperatorResult& res) 00622 { 00623 _reduction = reduction; 00624 (*this).apply(x,b,res); 00625 } 00626 00627 private: 00628 SeqScalarProduct<X> ssp; 00629 LinearOperator<X,X>& _op; 00630 Preconditioner<X,X>& _prec; 00631 ScalarProduct<X>& _sp; 00632 double _reduction; 00633 int _maxit; 00634 int _verbose; 00635 }; 00636 00637 00638 // Ronald Kriemanns BiCG-STAB implementation from Sumo 00640 template<class X> 00641 class BiCGSTABSolver : public InverseOperator<X,X> { 00642 public: 00644 typedef X domain_type; 00646 typedef X range_type; 00648 typedef typename X::field_type field_type; 00649 00655 template<class L, class P> 00656 BiCGSTABSolver (L& op, P& prec, 00657 double reduction, int maxit, int verbose) : 00658 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose) 00659 { 00660 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category), "L and P must be of the same category!"); 00661 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential), "L must be sequential!"); 00662 } 00668 template<class L, class S, class P> 00669 BiCGSTABSolver (L& op, S& sp, P& prec, 00670 double reduction, int maxit, int verbose) : 00671 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose) 00672 { 00673 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category), 00674 "L and P must have the same category!"); 00675 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category), 00676 "L and S must have the same category!"); 00677 } 00678 00684 virtual void apply (X& x, X& b, InverseOperatorResult& res) 00685 { 00686 const double EPSILON=1e-80; 00687 00688 double it; 00689 field_type rho, rho_new, alpha, beta, h, omega; 00690 field_type norm, norm_old, norm_0; 00691 00692 // 00693 // get vectors and matrix 00694 // 00695 X& r=b; 00696 X p(x); 00697 X v(x); 00698 X t(x); 00699 X y(x); 00700 X rt(x); 00701 00702 // 00703 // begin iteration 00704 // 00705 00706 // r = r - Ax; rt = r 00707 res.clear(); // clear solver statistics 00708 Timer watch; // start a timer 00709 _prec.pre(x,r); // prepare preconditioner 00710 _op.applyscaleadd(-1,x,r); // overwrite b with defect 00711 00712 rt=r; 00713 00714 norm = norm_old = norm_0 = _sp.norm(r); 00715 00716 p=0; 00717 v=0; 00718 00719 rho = 1; 00720 alpha = 1; 00721 omega = 1; 00722 00723 if (_verbose>0) // printing 00724 { 00725 std::cout << "=== BiCGSTABSolver" << std::endl; 00726 if (_verbose>1) 00727 { 00728 this->printHeader(std::cout); 00729 this->printOutput(std::cout,0,norm_0); 00730 //std::cout << " Iter Defect Rate" << std::endl; 00731 //std::cout << " 0" << std::setw(14) << norm_0 << std::endl; 00732 } 00733 } 00734 00735 if ( norm < (_reduction * norm_0) || norm<1E-30) 00736 { 00737 res.converged = 1; 00738 _prec.post(x); // postprocess preconditioner 00739 res.iterations = 0; // fill statistics 00740 res.reduction = 0; 00741 res.conv_rate = 0; 00742 res.elapsed = watch.elapsed(); 00743 return; 00744 } 00745 00746 // 00747 // iteration 00748 // 00749 00750 for (it = 0.5; it < _maxit; it+=.5) 00751 { 00752 // 00753 // preprocess, set vecsizes etc. 00754 // 00755 00756 // rho_new = < rt , r > 00757 rho_new = _sp.dot(rt,r); 00758 00759 // look if breakdown occured 00760 if (std::abs(rho) <= EPSILON) 00761 DUNE_THROW(ISTLError,"breakdown in BiCGSTAB - rho " 00762 << rho << " <= EPSILON " << EPSILON 00763 << " after " << it << " iterations"); 00764 if (std::abs(omega) <= EPSILON) 00765 DUNE_THROW(ISTLError,"breakdown in BiCGSTAB - omega " 00766 << omega << " <= EPSILON " << EPSILON 00767 << " after " << it << " iterations"); 00768 00769 00770 if (it<1) 00771 p = r; 00772 else 00773 { 00774 beta = ( rho_new / rho ) * ( alpha / omega ); 00775 p.axpy(-omega,v); // p = r + beta (p - omega*v) 00776 p *= beta; 00777 p += r; 00778 } 00779 00780 // y = W^-1 * p 00781 y = 0; 00782 _prec.apply(y,p); // apply preconditioner 00783 00784 // v = A * y 00785 _op.apply(y,v); 00786 00787 // alpha = rho_new / < rt, v > 00788 h = _sp.dot(rt,v); 00789 00790 if ( std::abs(h) < EPSILON ) 00791 DUNE_THROW(ISTLError,"h=0 in BiCGSTAB"); 00792 00793 alpha = rho_new / h; 00794 00795 // apply first correction to x 00796 // x <- x + alpha y 00797 x.axpy(alpha,y); 00798 00799 // r = r - alpha*v 00800 r.axpy(-alpha,v); 00801 00802 // 00803 // test stop criteria 00804 // 00805 00806 norm = _sp.norm(r); 00807 00808 if (_verbose>1) // print 00809 { 00810 this->printOutput(std::cout,it,norm,norm_old); 00811 } 00812 00813 if ( norm < (_reduction * norm_0) ) 00814 { 00815 res.converged = 1; 00816 break; 00817 } 00818 it+=.5; 00819 00820 norm_old = norm; 00821 00822 // y = W^-1 * r 00823 y = 0; 00824 _prec.apply(y,r); 00825 00826 // t = A * y 00827 _op.apply(y,t); 00828 00829 // omega = < t, r > / < t, t > 00830 omega = _sp.dot(t,r)/_sp.dot(t,t); 00831 00832 // apply second correction to x 00833 // x <- x + omega y 00834 x.axpy(omega,y); 00835 00836 // r = s - omega*t (remember : r = s) 00837 r.axpy(-omega,t); 00838 00839 rho = rho_new; 00840 00841 // 00842 // test stop criteria 00843 // 00844 00845 norm = _sp.norm(r); 00846 00847 if (_verbose > 1) // print 00848 { 00849 this->printOutput(std::cout,it,norm,norm_old); 00850 } 00851 00852 if ( norm < (_reduction * norm_0) || norm<1E-30) 00853 { 00854 res.converged = 1; 00855 break; 00856 } 00857 00858 norm_old = norm; 00859 } // end for 00860 00861 if (_verbose==1) // printing for non verbose 00862 this->printOutput(std::cout,it,norm); 00863 00864 _prec.post(x); // postprocess preconditioner 00865 res.iterations = static_cast<int>(std::ceil(it)); // fill statistics 00866 res.reduction = norm/norm_0; 00867 res.conv_rate = pow(res.reduction,1.0/it); 00868 res.elapsed = watch.elapsed(); 00869 if (_verbose>0) // final print 00870 std::cout << "=== rate=" << res.conv_rate 00871 << ", T=" << res.elapsed 00872 << ", TIT=" << res.elapsed/it 00873 << ", IT=" << it << std::endl; 00874 } 00875 00881 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res) 00882 { 00883 _reduction = reduction; 00884 (*this).apply(x,b,res); 00885 } 00886 00887 private: 00888 SeqScalarProduct<X> ssp; 00889 LinearOperator<X,X>& _op; 00890 Preconditioner<X,X>& _prec; 00891 ScalarProduct<X>& _sp; 00892 double _reduction; 00893 int _maxit; 00894 int _verbose; 00895 }; 00896 00903 template<class X> 00904 class MINRESSolver : public InverseOperator<X,X> { 00905 public: 00907 typedef X domain_type; 00909 typedef X range_type; 00911 typedef typename X::field_type field_type; 00912 00918 template<class L, class P> 00919 MINRESSolver (L& op, P& prec, double reduction, int maxit, int verbose) : 00920 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose) 00921 { 00922 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category), 00923 "L and P must have the same category!"); 00924 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential), 00925 "L must be sequential!"); 00926 } 00932 template<class L, class S, class P> 00933 MINRESSolver (L& op, S& sp, P& prec, double reduction, int maxit, int verbose) : 00934 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose) 00935 { 00936 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category), 00937 "L and P must have the same category!"); 00938 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category), 00939 "L and S must have the same category!"); 00940 } 00941 00947 virtual void apply (X& x, X& b, InverseOperatorResult& res) 00948 { 00949 res.clear(); // clear solver statistics 00950 Timer watch; // start a timer 00951 _prec.pre(x,b); // prepare preconditioner 00952 _op.applyscaleadd(-1,x,b); // overwrite b with defect/residual 00953 00954 double def0 = _sp.norm(b); // compute residual norm 00955 00956 if (def0<1E-30) // convergence check 00957 { 00958 res.converged = true; 00959 res.iterations = 0; // fill statistics 00960 res.reduction = 0; 00961 res.conv_rate = 0; 00962 res.elapsed=0; 00963 if (_verbose>0) // final print 00964 std::cout << "=== rate=" << res.conv_rate << ", T=" << res.elapsed << ", TIT=" << res.elapsed << ", IT=0" << std::endl; 00965 return; 00966 } 00967 00968 if (_verbose>0) // printing 00969 { 00970 std::cout << "=== MINRESSolver" << std::endl; 00971 if (_verbose>1) { 00972 this->printHeader(std::cout); 00973 this->printOutput(std::cout,0,def0); 00974 } 00975 } 00976 00977 // some local variables 00978 double def=def0; // the defect/residual norm 00979 field_type alpha, // recurrence coefficients as computed in the Lanczos alg making up the matrix T 00980 beta, // 00981 c[2]={0.0, 0.0}, // diagonal entry of Givens rotation 00982 s[2]={0.0, 0.0}; // off-diagonal entries of Givens rotation 00983 00984 field_type T[3]={0.0, 0.0, 0.0}; // recurrence coefficients (column k of Matrix T) 00985 00986 X z(b.size()), // some temporary vectors 00987 dummy(b.size()); 00988 00989 field_type xi[2]={1.0, 0.0}; 00990 00991 // initialize 00992 z = 0.0; // clear correction 00993 00994 _prec.apply(z,b); // apply preconditioner z=M^-1*b 00995 00996 beta = sqrt(fabs(_sp.dot(z,b))); 00997 double beta0 = beta; 00998 00999 X p[3]; // the search directions 01000 X q[3]; // Orthonormal basis vectors (in unpreconditioned case) 01001 01002 q[0].resize(b.size()); 01003 q[1].resize(b.size()); 01004 q[2].resize(b.size()); 01005 q[0] = 0.0; 01006 q[1] = b; 01007 q[1] /= beta; 01008 q[2] = 0.0; 01009 01010 p[0].resize(b.size()); 01011 p[1].resize(b.size()); 01012 p[2].resize(b.size()); 01013 p[0] = 0.0; 01014 p[1] = 0.0; 01015 p[2] = 0.0; 01016 01017 01018 z /= beta; // this is w_current 01019 01020 // the loop 01021 int i=1; 01022 for ( ; i<=_maxit; i++) 01023 { 01024 dummy = z; // remember z_old for the computation of the search direction p in the next iteration 01025 01026 int i1 = i%3, 01027 i0 = (i1+2)%3, 01028 i2 = (i1+1)%3; 01029 01030 // Symmetrically Preconditioned Lanczos (Greenbaum p.121) 01031 _op.apply(z,q[i2]); // q[i2] = Az 01032 q[i2].axpy(-beta, q[i0]); 01033 alpha = _sp.dot(q[i2],z); 01034 q[i2].axpy(-alpha, q[i1]); 01035 01036 z=0.0; 01037 _prec.apply(z,q[i2]); 01038 01039 beta = sqrt(fabs(_sp.dot(q[i2],z))); 01040 01041 q[i2] /= beta; 01042 z /= beta; 01043 01044 // QR Factorization of recurrence coefficient matrix 01045 // apply previous Givens rotations to last column of T 01046 T[1] = T[2]; 01047 if (i>2) 01048 { 01049 T[0] = s[i%2]*T[1]; 01050 T[1] = c[i%2]*T[1]; 01051 } 01052 if (i>1) 01053 { 01054 T[2] = c[(i+1)%2]*alpha - s[(i+1)%2]*T[1]; 01055 T[1] = c[(i+1)%2]*T[1] + s[(i+1)%2]*alpha; 01056 } 01057 else 01058 T[2] = alpha; 01059 01060 // recompute c, s -> current Givens rotation \TODO use BLAS-routine drotg instead for greater robustness 01061 // cblas_drotg (a, b, c, s); 01062 c[i%2] = 1.0/sqrt(T[2]*T[2] + beta*beta); 01063 s[i%2] = beta*c[i%2]; 01064 c[i%2] *= T[2]; 01065 01066 // apply current Givens rotation to T eliminating the last entry... 01067 T[2] = c[i%2]*T[2] + s[i%2]*beta; 01068 01069 // ...and to xi, the right hand side of the least squares problem min_y||beta*xi-T*y|| 01070 xi[i%2] = -s[i%2]*xi[(i+1)%2]; 01071 xi[(i+1)%2] *= c[i%2]; 01072 01073 // compute correction direction 01074 p[i2] = dummy; 01075 p[i2].axpy(-T[1],p[i1]); 01076 p[i2].axpy(-T[0],p[i0]); 01077 p[i2] /= T[2]; 01078 01079 // apply correction/update solution 01080 x.axpy(beta0*xi[(i+1)%2], p[i2]); 01081 01082 // remember beta_old 01083 T[2] = beta; 01084 01085 // update residual - not necessary if in the preconditioned case we are content with the residual norm of the 01086 // preconditioned system as convergence test 01087 // _op.apply(p[i2],dummy); 01088 // b.axpy(-beta0*xi[(i+1)%2],dummy); 01089 01090 // convergence test 01091 // double defnew=_sp.norm(b); // residual norm of original system 01092 double defnew = fabs(beta0*xi[i%2]); // the last entry the QR-transformed least squares RHS is the new residual norm 01093 01094 if (_verbose>1) // print 01095 this->printOutput(std::cout,i,defnew,def); 01096 01097 def = defnew; // update norm 01098 if (def<def0*_reduction || def<1E-30 || i==_maxit) // convergence check 01099 { 01100 res.converged = true; 01101 break; 01102 } 01103 } 01104 01105 if (_verbose==1) // printing for non verbose 01106 this->printOutput(std::cout,i,def); 01107 01108 _prec.post(x); // postprocess preconditioner 01109 res.iterations = i; // fill statistics 01110 res.reduction = def/def0; 01111 res.conv_rate = pow(res.reduction,1.0/i); 01112 res.elapsed = watch.elapsed(); 01113 01114 if (_verbose>0) // final print 01115 { 01116 std::cout << "=== rate=" << res.conv_rate 01117 << ", T=" << res.elapsed 01118 << ", TIT=" << res.elapsed/i 01119 << ", IT=" << i << std::endl; 01120 } 01121 01122 } 01123 01129 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res) 01130 { 01131 _reduction = reduction; 01132 (*this).apply(x,b,res); 01133 } 01134 01135 private: 01136 SeqScalarProduct<X> ssp; 01137 LinearOperator<X,X>& _op; 01138 Preconditioner<X,X>& _prec; 01139 ScalarProduct<X>& _sp; 01140 double _reduction; 01141 int _maxit; 01142 int _verbose; 01143 }; 01144 01156 template<class X, class Y=X, class F = Y> 01157 class RestartedGMResSolver : public InverseOperator<X,Y> 01158 { 01159 public: 01161 typedef X domain_type; 01163 typedef Y range_type; 01165 typedef typename X::field_type field_type; 01167 typedef F basis_type; 01168 01176 template<class L, class P> 01177 RestartedGMResSolver (L& op, P& prec, double reduction, int restart, int maxit, int verbose, bool recalc_defect = false) : 01178 _A_(op), _M(prec), 01179 ssp(), _sp(ssp), _restart(restart), 01180 _reduction(reduction), _maxit(maxit), _verbose(verbose), 01181 _recalc_defect(recalc_defect) 01182 { 01183 dune_static_assert(static_cast<int>(P::category) == static_cast<int>(L::category), 01184 "P and L must be the same category!"); 01185 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential), 01186 "L must be sequential!"); 01187 } 01188 01196 template<class L, class S, class P> 01197 RestartedGMResSolver (L& op, S& sp, P& prec, double reduction, int restart, int maxit, int verbose, bool recalc_defect = false) : 01198 _A_(op), _M(prec), 01199 _sp(sp), _restart(restart), 01200 _reduction(reduction), _maxit(maxit), _verbose(verbose), 01201 _recalc_defect(recalc_defect) 01202 { 01203 dune_static_assert(static_cast<int>(P::category) == static_cast<int>(L::category), 01204 "P and L must have the same category!"); 01205 dune_static_assert(static_cast<int>(P::category) == static_cast<int>(S::category), 01206 "P and S must have the same category!"); 01207 } 01208 01210 virtual void apply (X& x, X& b, InverseOperatorResult& res) 01211 { 01212 apply(x,b,_reduction,res); 01213 }; 01214 01220 virtual void apply (X& x, Y& b, double reduction, InverseOperatorResult& res) 01221 { 01222 int m = _restart; 01223 field_type norm; 01224 field_type norm_old = 0.0; 01225 field_type norm_0; 01226 field_type beta; 01227 int i, j = 1, k; 01228 std::vector<field_type> s(m+1), cs(m), sn(m); 01229 // helper vector 01230 X w(b); 01231 std::vector< std::vector<field_type> > H(m+1,s); 01232 std::vector<F> v(m+1,b); 01233 01234 // start timer 01235 Timer watch; // start a timer 01236 01237 // clear solver statistics 01238 res.clear(); 01239 01240 _M.pre(x,b); 01241 01242 if (_recalc_defect) 01243 { 01244 // norm_0 = norm(M^-1 b) 01245 w = 0.0; _M.apply(w,b); // w = M^-1 b 01246 norm_0 = _sp.norm(w); 01247 // r = _M.solve(b - A * x); 01248 w = b; 01249 _A_.applyscaleadd(-1,x, /* => */ w); // w = b - Ax; 01250 v[0] = 0.0; _M.apply(v[0],w); // r = M^-1 w 01251 beta = _sp.norm(v[0]); 01252 } 01253 else 01254 { 01255 // norm_0 = norm(M^-1 b) 01256 w = 0.0; _M.apply(w,b); // w = M^-1 b 01257 norm_0 = _sp.norm(w); 01258 // r = _M.solve(b - A * x); 01259 _A_.applyscaleadd(-1,x, /* => */ b); // b = b - Ax; 01260 v[0] = 0.0; _M.apply(v[0],b); // r = M^-1 b 01261 beta = _sp.norm(v[0]); 01262 } 01263 01264 // avoid division by zero 01265 if (norm_0 == 0.0) 01266 norm_0 = 1.0; 01267 01268 norm = norm_old = _sp.norm(v[0]); 01269 01270 // print header 01271 if (_verbose > 0) 01272 { 01273 std::cout << "=== RestartedGMResSolver" << std::endl; 01274 if (_verbose > 1) 01275 { 01276 this->printHeader(std::cout); 01277 this->printOutput(std::cout,0,norm_0); 01278 this->printOutput(std::cout,0,norm); 01279 } 01280 } 01281 01282 // check convergence 01283 if (norm <= reduction * norm_0) { 01284 _M.post(x); // postprocess preconditioner 01285 res.converged = true; 01286 if (_verbose > 0) // final print 01287 print_result(res); 01288 return; 01289 } 01290 01291 while (j <= _maxit && res.converged != true) { 01292 v[0] *= (1.0 / beta); 01293 for (i=1; i<=m; i++) s[i] = 0.0; 01294 s[0] = beta; 01295 01296 for (i = 0; i < m && j <= _maxit && res.converged != true; i++, j++) { 01297 w = 0.0; 01298 v[i+1] = 0.0; // use v[i+1] as temporary vector 01299 _A_.apply(v[i], /* => */ v[i+1]); 01300 _M.apply(w, v[i+1]); 01301 for (k = 0; k <= i; k++) { 01302 H[k][i] = _sp.dot(w, v[k]); 01303 // w -= H[k][i] * v[k]; 01304 w.axpy(-H[k][i], v[k]); 01305 } 01306 H[i+1][i] = _sp.norm(w); 01307 if (H[i+1][i] == 0.0) 01308 DUNE_THROW(ISTLError,"breakdown in GMRes - |w| " 01309 << w << " == 0.0 after " << j << " iterations"); 01310 // v[i+1] = w * (1.0 / H[i+1][i]); 01311 v[i+1] = w; v[i+1] *= (1.0 / H[i+1][i]); 01312 01313 for (k = 0; k < i; k++) 01314 applyPlaneRotation(H[k][i], H[k+1][i], cs[k], sn[k]); 01315 01316 generatePlaneRotation(H[i][i], H[i+1][i], cs[i], sn[i]); 01317 applyPlaneRotation(H[i][i], H[i+1][i], cs[i], sn[i]); 01318 applyPlaneRotation(s[i], s[i+1], cs[i], sn[i]); 01319 01320 norm = std::abs(s[i+1]); 01321 01322 if (_verbose > 1) // print 01323 { 01324 this->printOutput(std::cout,j,norm,norm_old); 01325 } 01326 01327 norm_old = norm; 01328 01329 if (norm < reduction * norm_0) { 01330 res.converged = true; 01331 } 01332 } 01333 01334 if (_recalc_defect) 01335 { 01336 // update x 01337 update(x, i - 1, H, s, v); 01338 01339 // update residuum 01340 // r = M^-1 (b - A * x); 01341 w = b; _A_.applyscaleadd(-1,x, /* => */ w); 01342 _M.apply(v[0], w); 01343 beta = _sp.norm(v[0]); 01344 norm = beta; 01345 } 01346 else 01347 { 01348 // calc update vector 01349 w = 0; 01350 update(w, i - 1, H, s, v); 01351 01352 // update x 01353 x += w; 01354 01355 // r = M^-1 (b - A * x); 01356 // update defect 01357 _A_.applyscaleadd(-1,w, /* => */ b); 01358 // r = M^-1 (b - A * x); 01359 v[0] = 0.0; _M.apply(v[0],b); // r = M^-1 b 01360 beta = _sp.norm(v[0]); 01361 norm = beta; 01362 01363 res.converged = false; 01364 } 01365 01366 if (_verbose > 1) // print 01367 { 01368 this->printOutput(std::cout,j,norm,norm_old); 01369 } 01370 01371 norm_old = norm; 01372 01373 if (norm < reduction * norm_0) { 01374 // fill statistics 01375 res.converged = true; 01376 } 01377 01378 if (res.converged != true && _verbose > 0) 01379 std::cout << "=== GMRes::restart\n"; 01380 } 01381 01382 _M.post(x); // postprocess preconditioner 01383 01384 res.iterations = j; 01385 res.reduction = norm / norm_0; 01386 res.conv_rate = pow(res.reduction,1.0/j); 01387 res.elapsed = watch.elapsed(); 01388 01389 if (_verbose>0) 01390 print_result(res); 01391 } 01392 private: 01393 01394 void 01395 print_result (const InverseOperatorResult & res) const 01396 { 01397 int j = res.iterations>0?res.iterations:1; 01398 std::cout << "=== rate=" << res.conv_rate 01399 << ", T=" << res.elapsed 01400 << ", TIT=" << res.elapsed/j 01401 << ", IT=" << res.iterations 01402 << std::endl; 01403 } 01404 01405 static void 01406 update(X &x, int k, 01407 std::vector< std::vector<field_type> > & h, 01408 std::vector<field_type> & s, std::vector<F> v) 01409 { 01410 std::vector<field_type> y(s); 01411 01412 // Backsolve: 01413 for (int i = k; i >= 0; i--) { 01414 y[i] /= h[i][i]; 01415 for (int j = i - 1; j >= 0; j--) 01416 y[j] -= h[j][i] * y[i]; 01417 } 01418 01419 for (int j = 0; j <= k; j++) 01420 // x += v[j] * y[j]; 01421 x.axpy(y[j],v[j]); 01422 } 01423 01424 void 01425 generatePlaneRotation(field_type &dx, field_type &dy, field_type &cs, field_type &sn) 01426 { 01427 if (dy == 0.0) { 01428 cs = 1.0; 01429 sn = 0.0; 01430 } else if (std::abs(dy) > std::abs(dx)) { 01431 field_type temp = dx / dy; 01432 sn = 1.0 / std::sqrt( 1.0 + temp*temp ); 01433 cs = temp * sn; 01434 } else { 01435 field_type temp = dy / dx; 01436 cs = 1.0 / std::sqrt( 1.0 + temp*temp ); 01437 sn = temp * cs; 01438 } 01439 } 01440 01441 01442 void 01443 applyPlaneRotation(field_type &dx, field_type &dy, field_type &cs, field_type &sn) 01444 { 01445 field_type temp = cs * dx + sn * dy; 01446 dy = -sn * dx + cs * dy; 01447 dx = temp; 01448 } 01449 01450 LinearOperator<X,X>& _A_; 01451 Preconditioner<X,X>& _M; 01452 SeqScalarProduct<X> ssp; 01453 ScalarProduct<X>& _sp; 01454 int _restart; 01455 double _reduction; 01456 int _maxit; 01457 int _verbose; 01458 bool _recalc_defect; 01459 }; 01460 01463 } // end namespace 01464 01465 #endif
Generated on Fri Apr 29 2011 with Doxygen (ver 1.7.1) [doxygen-log,error-log].