Dune Core Modules (2.5.0)

solvers.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4#ifndef DUNE_ISTL_SOLVERS_HH
5#define DUNE_ISTL_SOLVERS_HH
6
7#include <cmath>
8#include <complex>
9#include <iostream>
10#include <iomanip>
11#include <memory>
12#include <string>
13#include <vector>
14#include <array>
15
16#include "istlexception.hh"
17#include "operators.hh"
18#include "scalarproducts.hh"
19#include "solver.hh"
20#include "preconditioner.hh"
23#include <dune/common/timer.hh>
26
27namespace Dune {
45 //=====================================================================
46 // Implementation of this interface
47 //=====================================================================
48
57 template<class X>
58 class LoopSolver : public InverseOperator<X,X> {
59 public:
61 typedef X domain_type;
63 typedef X range_type;
65 typedef typename X::field_type field_type;
67 typedef typename FieldTraits<field_type>::real_type real_type;
68
88 template<class L, class P>
89 LoopSolver (L& op, P& prec,
90 real_type reduction, int maxit, int verbose) :
91 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
92 {
93 static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
94 "L and P have to have the same category!");
95 static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
96 "L has to be sequential!");
97 }
98
119 template<class L, class S, class P>
120 LoopSolver (L& op, S& sp, P& prec,
121 real_type reduction, int maxit, int verbose) :
122 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
123 {
124 static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
125 "L and P must have the same category!");
126 static_assert(static_cast<int>(L::category) == static_cast<int>(S::category),
127 "L and S must have the same category!");
128 }
129
130
132 virtual void apply (X& x, X& b, InverseOperatorResult& res)
133 {
134 // clear solver statistics
135 res.clear();
136
137 // start a timer
138 Timer watch;
139
140 // prepare preconditioner
141 _prec.pre(x,b);
142
143 // overwrite b with defect
144 _op.applyscaleadd(-1,x,b);
145
146 // compute norm, \todo parallelization
147 real_type def0 = _sp.norm(b);
148
149 // printing
150 if (_verbose>0)
151 {
152 std::cout << "=== LoopSolver" << std::endl;
153 if (_verbose>1)
154 {
155 this->printHeader(std::cout);
156 this->printOutput(std::cout,0,def0);
157 }
158 }
159
160 // allocate correction vector
161 X v(x);
162
163 // iteration loop
164 int i=1; real_type def=def0;
165 for ( ; i<=_maxit; i++ )
166 {
167 v = 0; // clear correction
168 _prec.apply(v,b); // apply preconditioner
169 x += v; // update solution
170 _op.applyscaleadd(-1,v,b); // update defect
171 real_type defnew=_sp.norm(b); // comp defect norm
172 if (_verbose>1) // print
173 this->printOutput(std::cout,i,defnew,def);
174 //std::cout << i << " " << defnew << " " << defnew/def << std::endl;
175 def = defnew; // update norm
176 if (def<def0*_reduction || def<1E-30) // convergence check
177 {
178 res.converged = true;
179 break;
180 }
181 }
182
183 //correct i which is wrong if convergence was not achieved.
184 i=std::min(_maxit,i);
185
186 // print
187 if (_verbose==1)
188 this->printOutput(std::cout,i,def);
189
190 // postprocess preconditioner
191 _prec.post(x);
192
193 // fill statistics
194 res.iterations = i;
195 res.reduction = def/def0;
196 res.conv_rate = pow(res.reduction,1.0/i);
197 res.elapsed = watch.elapsed();
198
199 // final print
200 if (_verbose>0)
201 {
202 std::cout << "=== rate=" << res.conv_rate
203 << ", T=" << res.elapsed
204 << ", TIT=" << res.elapsed/i
205 << ", IT=" << i << std::endl;
206 }
207 }
208
210 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
211 {
212 real_type saved_reduction = _reduction;
213 _reduction = reduction;
214 (*this).apply(x,b,res);
215 _reduction = saved_reduction;
216 }
217
218 private:
221 Preconditioner<X,X>& _prec;
222 ScalarProduct<X>& _sp;
223 real_type _reduction;
224 int _maxit;
225 int _verbose;
226 };
227
228
229 // all these solvers are taken from the SUMO library
231 template<class X>
232 class GradientSolver : public InverseOperator<X,X> {
233 public:
235 typedef X domain_type;
237 typedef X range_type;
239 typedef typename X::field_type field_type;
241 typedef typename FieldTraits<field_type>::real_type real_type;
242
243
249 template<class L, class P>
250 GradientSolver (L& op, P& prec,
251 real_type reduction, int maxit, int verbose) :
252 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
253 {
254 static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
255 "L and P have to have the same category!");
256 static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
257 "L has to be sequential!");
258 }
264 template<class L, class S, class P>
265 GradientSolver (L& op, S& sp, P& prec,
266 real_type reduction, int maxit, int verbose) :
267 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
268 {
269 static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
270 "L and P have to have the same category!");
271 static_assert(static_cast<int>(L::category) == static_cast<int>(S::category),
272 "L and S have to have the same category!");
273 }
274
280 virtual void apply (X& x, X& b, InverseOperatorResult& res)
281 {
282 res.clear(); // clear solver statistics
283 Timer watch; // start a timer
284 _prec.pre(x,b); // prepare preconditioner
285 _op.applyscaleadd(-1,x,b); // overwrite b with defect
286
287 X p(x); // create local vectors
288 X q(b);
289
290 real_type def0 = _sp.norm(b); // compute norm
291
292 if (_verbose>0) // printing
293 {
294 std::cout << "=== GradientSolver" << std::endl;
295 if (_verbose>1)
296 {
297 this->printHeader(std::cout);
298 this->printOutput(std::cout,0,def0);
299 }
300 }
301
302 int i=1; real_type def=def0; // loop variables
303 field_type lambda;
304 for ( ; i<=_maxit; i++ )
305 {
306 p = 0; // clear correction
307 _prec.apply(p,b); // apply preconditioner
308 _op.apply(p,q); // q=Ap
309 lambda = _sp.dot(p,b)/_sp.dot(q,p); // minimization
310 x.axpy(lambda,p); // update solution
311 b.axpy(-lambda,q); // update defect
312
313 real_type defnew=_sp.norm(b); // comp defect norm
314 if (_verbose>1) // print
315 this->printOutput(std::cout,i,defnew,def);
316
317 def = defnew; // update norm
318 if (def<def0*_reduction || def<1E-30) // convergence check
319 {
320 res.converged = true;
321 break;
322 }
323 }
324
325 //correct i which is wrong if convergence was not achieved.
326 i=std::min(_maxit,i);
327
328 if (_verbose==1) // printing for non verbose
329 this->printOutput(std::cout,i,def);
330
331 _prec.post(x); // postprocess preconditioner
332 res.iterations = i; // fill statistics
333 res.reduction = static_cast<double>(def/def0);
334 res.conv_rate = static_cast<double>(pow(res.reduction,1.0/i));
335 res.elapsed = watch.elapsed();
336 if (_verbose>0) // final print
337 std::cout << "=== rate=" << res.conv_rate
338 << ", T=" << res.elapsed
339 << ", TIT=" << res.elapsed/i
340 << ", IT=" << i << std::endl;
341 }
342
348 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
349 {
350 real_type saved_reduction = _reduction;
351 _reduction = reduction;
352 (*this).apply(x,b,res);
353 _reduction = saved_reduction;
354 }
355
356 private:
359 Preconditioner<X,X>& _prec;
360 ScalarProduct<X>& _sp;
361 real_type _reduction;
362 int _maxit;
363 int _verbose;
364 };
365
366
367
369 template<class X>
370 class CGSolver : public InverseOperator<X,X> {
371 public:
373 typedef X domain_type;
375 typedef X range_type;
377 typedef typename X::field_type field_type;
379 typedef typename FieldTraits<field_type>::real_type real_type;
380
386 template<class L, class P>
387 CGSolver (L& op, P& prec, real_type reduction, int maxit, int verbose) :
388 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
389 {
390 static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
391 "L and P must have the same category!");
392 static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
393 "L must be sequential!");
394 }
400 template<class L, class S, class P>
401 CGSolver (L& op, S& sp, P& prec, real_type reduction, int maxit, int verbose) :
402 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
403 {
404 static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
405 "L and P must have the same category!");
406 static_assert(static_cast<int>(L::category) == static_cast<int>(S::category),
407 "L and S must have the same category!");
408 }
409
421 virtual void apply (X& x, X& b, InverseOperatorResult& res)
422 {
423 using std::isfinite;
424
425 res.clear(); // clear solver statistics
426 Timer watch; // start a timer
427 _prec.pre(x,b); // prepare preconditioner
428 _op.applyscaleadd(-1,x,b); // overwrite b with defect
429
430 X p(x); // the search direction
431 X q(x); // a temporary vector
432
433 real_type def0 = _sp.norm(b); // compute norm
434
435 if (!isfinite(def0)) // check for inf or NaN
436 {
437 if (_verbose>0)
438 std::cout << "=== CGSolver: abort due to infinite or NaN initial defect"
439 << std::endl;
440 DUNE_THROW(SolverAbort, "CGSolver: initial defect=" << def0
441 << " is infinite or NaN");
442 }
443
444 if (def0<1E-30) // convergence check
445 {
446 res.converged = true;
447 res.iterations = 0; // fill statistics
448 res.reduction = 0;
449 res.conv_rate = 0;
450 res.elapsed=0;
451 if (_verbose>0) // final print
452 std::cout << "=== rate=" << res.conv_rate
453 << ", T=" << res.elapsed << ", TIT=" << res.elapsed
454 << ", IT=0" << std::endl;
455 return;
456 }
457
458 if (_verbose>0) // printing
459 {
460 std::cout << "=== CGSolver" << std::endl;
461 if (_verbose>1) {
462 this->printHeader(std::cout);
463 this->printOutput(std::cout,real_type(0),def0);
464 }
465 }
466
467 // some local variables
468 real_type def=def0; // loop variables
469 field_type rho,rholast,lambda,alpha,beta;
470
471 // determine initial search direction
472 p = 0; // clear correction
473 _prec.apply(p,b); // apply preconditioner
474 rholast = _sp.dot(p,b); // orthogonalization
475
476 // the loop
477 int i=1;
478 for ( ; i<=_maxit; i++ )
479 {
480 // minimize in given search direction p
481 _op.apply(p,q); // q=Ap
482 alpha = _sp.dot(p,q); // scalar product
483 lambda = rholast/alpha; // minimization
484 x.axpy(lambda,p); // update solution
485 b.axpy(-lambda,q); // update defect
486
487 // convergence test
488 real_type defnew=_sp.norm(b); // comp defect norm
489
490 if (_verbose>1) // print
491 this->printOutput(std::cout,real_type(i),defnew,def);
492
493 def = defnew; // update norm
494 if (!isfinite(def)) // check for inf or NaN
495 {
496 if (_verbose>0)
497 std::cout << "=== CGSolver: abort due to infinite or NaN defect"
498 << std::endl;
500 "CGSolver: defect=" << def << " is infinite or NaN");
501 }
502
503 if (def<def0*_reduction || def<1E-30) // convergence check
504 {
505 res.converged = true;
506 break;
507 }
508
509 // determine new search direction
510 q = 0; // clear correction
511 _prec.apply(q,b); // apply preconditioner
512 rho = _sp.dot(q,b); // orthogonalization
513 beta = rho/rholast; // scaling factor
514 p *= beta; // scale old search direction
515 p += q; // orthogonalization with correction
516 rholast = rho; // remember rho for recurrence
517 }
518
519 //correct i which is wrong if convergence was not achieved.
520 i=std::min(_maxit,i);
521
522 if (_verbose==1) // printing for non verbose
523 this->printOutput(std::cout,real_type(i),def);
524
525 _prec.post(x); // postprocess preconditioner
526 res.iterations = i; // fill statistics
527 res.reduction = static_cast<double>(def/def0);
528 res.conv_rate = static_cast<double>(pow(res.reduction,1.0/i));
529 res.elapsed = watch.elapsed();
530
531 if (_verbose>0) // final print
532 {
533 std::cout << "=== rate=" << res.conv_rate
534 << ", T=" << res.elapsed
535 << ", TIT=" << res.elapsed/i
536 << ", IT=" << i << std::endl;
537 }
538 }
539
551 virtual void apply (X& x, X& b, double reduction,
553 {
554 real_type saved_reduction = _reduction;
555 _reduction = reduction;
556 (*this).apply(x,b,res);
557 _reduction = saved_reduction;
558 }
559
560 private:
563 Preconditioner<X,X>& _prec;
564 ScalarProduct<X>& _sp;
565 real_type _reduction;
566 int _maxit;
567 int _verbose;
568 };
569
570
571 // Ronald Kriemanns BiCG-STAB implementation from Sumo
573 template<class X>
574 class BiCGSTABSolver : public InverseOperator<X,X> {
575 public:
577 typedef X domain_type;
579 typedef X range_type;
581 typedef typename X::field_type field_type;
583 typedef typename FieldTraits<field_type>::real_type real_type;
584
590 template<class L, class P>
591 BiCGSTABSolver (L& op, P& prec,
592 real_type reduction, int maxit, int verbose) :
593 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
594 {
595 static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
596 "L and P must be of the same category!");
597 static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
598 "L must be sequential!");
599 }
605 template<class L, class S, class P>
606 BiCGSTABSolver (L& op, S& sp, P& prec,
607 real_type reduction, int maxit, int verbose) :
608 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
609 {
610 static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
611 "L and P must have the same category!");
612 static_assert(static_cast<int>(L::category) == static_cast<int>(S::category),
613 "L and S must have the same category!");
614 }
615
623 virtual void apply (X& x, X& b, InverseOperatorResult& res)
624 {
625 using std::abs;
626 const real_type EPSILON=1e-80;
627 double it;
628 field_type rho, rho_new, alpha, beta, h, omega;
629 real_type norm, norm_old, norm_0;
630
631 //
632 // get vectors and matrix
633 //
634 X& r=b;
635 X p(x);
636 X v(x);
637 X t(x);
638 X y(x);
639 X rt(x);
640
641 //
642 // begin iteration
643 //
644
645 // r = r - Ax; rt = r
646 res.clear(); // clear solver statistics
647 Timer watch; // start a timer
648 _prec.pre(x,r); // prepare preconditioner
649 _op.applyscaleadd(-1,x,r); // overwrite b with defect
650
651 rt=r;
652
653 norm = norm_old = norm_0 = _sp.norm(r);
654
655 p=0;
656 v=0;
657
658 rho = 1;
659 alpha = 1;
660 omega = 1;
661
662 if (_verbose>0) // printing
663 {
664 std::cout << "=== BiCGSTABSolver" << std::endl;
665 if (_verbose>1)
666 {
667 this->printHeader(std::cout);
668 this->printOutput(std::cout,0,norm_0);
669 //std::cout << " Iter Defect Rate" << std::endl;
670 //std::cout << " 0" << std::setw(14) << norm_0 << std::endl;
671 }
672 }
673
674 if ( norm < (_reduction * norm_0) || norm<1E-30)
675 {
676 res.converged = 1;
677 _prec.post(x); // postprocess preconditioner
678 res.iterations = 0; // fill statistics
679 res.reduction = 0;
680 res.conv_rate = 0;
681 res.elapsed = watch.elapsed();
682 return;
683 }
684
685 //
686 // iteration
687 //
688
689 for (it = 0.5; it < _maxit; it+=.5)
690 {
691 //
692 // preprocess, set vecsizes etc.
693 //
694
695 // rho_new = < rt , r >
696 rho_new = _sp.dot(rt,r);
697
698 // look if breakdown occurred
699 if (abs(rho) <= EPSILON)
700 DUNE_THROW(SolverAbort,"breakdown in BiCGSTAB - rho "
701 << rho << " <= EPSILON " << EPSILON
702 << " after " << it << " iterations");
703 if (abs(omega) <= EPSILON)
704 DUNE_THROW(SolverAbort,"breakdown in BiCGSTAB - omega "
705 << omega << " <= EPSILON " << EPSILON
706 << " after " << it << " iterations");
707
708
709 if (it<1)
710 p = r;
711 else
712 {
713 beta = ( rho_new / rho ) * ( alpha / omega );
714 p.axpy(-omega,v); // p = r + beta (p - omega*v)
715 p *= beta;
716 p += r;
717 }
718
719 // y = W^-1 * p
720 y = 0;
721 _prec.apply(y,p); // apply preconditioner
722
723 // v = A * y
724 _op.apply(y,v);
725
726 // alpha = rho_new / < rt, v >
727 h = _sp.dot(rt,v);
728
729 if (abs(h) < EPSILON)
730 DUNE_THROW(SolverAbort,"abs(h) < EPSILON in BiCGSTAB - abs(h) "
731 << abs(h) << " < EPSILON " << EPSILON
732 << " after " << it << " iterations");
733
734 alpha = rho_new / h;
735
736 // apply first correction to x
737 // x <- x + alpha y
738 x.axpy(alpha,y);
739
740 // r = r - alpha*v
741 r.axpy(-alpha,v);
742
743 //
744 // test stop criteria
745 //
746
747 norm = _sp.norm(r);
748
749 if (_verbose>1) // print
750 {
751 this->printOutput(std::cout,it,norm,norm_old);
752 }
753
754 if ( norm < (_reduction * norm_0) )
755 {
756 res.converged = 1;
757 break;
758 }
759 it+=.5;
760
761 norm_old = norm;
762
763 // y = W^-1 * r
764 y = 0;
765 _prec.apply(y,r);
766
767 // t = A * y
768 _op.apply(y,t);
769
770 // omega = < t, r > / < t, t >
771 omega = _sp.dot(t,r)/_sp.dot(t,t);
772
773 // apply second correction to x
774 // x <- x + omega y
775 x.axpy(omega,y);
776
777 // r = s - omega*t (remember : r = s)
778 r.axpy(-omega,t);
779
780 rho = rho_new;
781
782 //
783 // test stop criteria
784 //
785
786 norm = _sp.norm(r);
787
788 if (_verbose > 1) // print
789 {
790 this->printOutput(std::cout,it,norm,norm_old);
791 }
792
793 if ( norm < (_reduction * norm_0) || norm<1E-30)
794 {
795 res.converged = 1;
796 break;
797 }
798
799 norm_old = norm;
800 } // end for
801
802 //correct i which is wrong if convergence was not achieved.
803 it=std::min(static_cast<double>(_maxit),it);
804
805 if (_verbose==1) // printing for non verbose
806 this->printOutput(std::cout,it,norm);
807
808 _prec.post(x); // postprocess preconditioner
809 res.iterations = static_cast<int>(std::ceil(it)); // fill statistics
810 res.reduction = static_cast<double>(norm/norm_0);
811 res.conv_rate = static_cast<double>(pow(res.reduction,1.0/it));
812 res.elapsed = watch.elapsed();
813 if (_verbose>0) // final print
814 std::cout << "=== rate=" << res.conv_rate
815 << ", T=" << res.elapsed
816 << ", TIT=" << res.elapsed/it
817 << ", IT=" << it << std::endl;
818 }
819
827 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
828 {
829 real_type saved_reduction = _reduction;
830 _reduction = reduction;
831 (*this).apply(x,b,res);
832 _reduction = saved_reduction;
833 }
834
835 private:
838 Preconditioner<X,X>& _prec;
839 ScalarProduct<X>& _sp;
840 real_type _reduction;
841 int _maxit;
842 int _verbose;
843 };
844
851 template<class X>
852 class MINRESSolver : public InverseOperator<X,X> {
853 public:
855 typedef X domain_type;
857 typedef X range_type;
859 typedef typename X::field_type field_type;
861 typedef typename FieldTraits<field_type>::real_type real_type;
862
868 template<class L, class P>
869 MINRESSolver (L& op, P& prec, real_type reduction, int maxit, int verbose) :
870 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
871 {
872 static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
873 "L and P must have the same category!");
874 static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
875 "L must be sequential!");
876 }
882 template<class L, class S, class P>
883 MINRESSolver (L& op, S& sp, P& prec, real_type reduction, int maxit, int verbose) :
884 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
885 {
886 static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
887 "L and P must have the same category!");
888 static_assert(static_cast<int>(L::category) == static_cast<int>(S::category),
889 "L and S must have the same category!");
890 }
891
897 virtual void apply (X& x, X& b, InverseOperatorResult& res)
898 {
899 using std::sqrt;
900 using std::abs;
901 // clear solver statistics
902 res.clear();
903 // start a timer
904 Dune::Timer watch;
905 watch.reset();
906 // prepare preconditioner
907 _prec.pre(x,b);
908 // overwrite rhs with defect
909 _op.applyscaleadd(-1,x,b);
910
911 // compute residual norm
912 real_type def0 = _sp.norm(b);
913
914 // printing
915 if(_verbose > 0) {
916 std::cout << "=== MINRESSolver" << std::endl;
917 if(_verbose > 1) {
918 this->printHeader(std::cout);
919 this->printOutput(std::cout,0,def0);
920 }
921 }
922
923 // check for convergence
924 if(def0 < 1e-30 ) {
925 res.converged = true;
926 res.iterations = 0;
927 res.reduction = 0;
928 res.conv_rate = 0;
929 res.elapsed = 0.0;
930 // final print
931 if(_verbose > 0)
932 std::cout << "=== rate=" << res.conv_rate
933 << ", T=" << res.elapsed
934 << ", TIT=" << res.elapsed
935 << ", IT=" << res.iterations
936 << std::endl;
937 return;
938 }
939
940 // the defect norm
941 real_type def = def0;
942 // recurrence coefficients as computed in Lanczos algorithm
943 field_type alpha, beta;
944 // diagonal entries of givens rotation
945 std::array<real_type,2> c{{0.0,0.0}};
946 // off-diagonal entries of givens rotation
947 std::array<field_type,2> s{{0.0,0.0}};
948
949 // recurrence coefficients (column k of tridiag matrix T_k)
950 std::array<field_type,3> T{{0.0,0.0,0.0}};
951
952 // the rhs vector of the min problem
953 std::array<field_type,2> xi{{1.0,0.0}};
954
955 // some temporary vectors
956 X z(b), dummy(b);
957
958 // initialize and clear correction
959 z = 0.0;
960 _prec.apply(z,b);
961
962 // beta is real and positive in exact arithmetic
963 // since it is the norm of the basis vectors (in unpreconditioned case)
964 beta = sqrt(_sp.dot(b,z));
965 field_type beta0 = beta;
966
967 // the search directions
968 std::array<X,3> p{{b,b,b}};
969 p[0] = 0.0;
970 p[1] = 0.0;
971 p[2] = 0.0;
972
973 // orthonormal basis vectors (in unpreconditioned case)
974 std::array<X,3> q{{b,b,b}};
975 q[0] = 0.0;
976 q[1] *= 1.0/beta;
977 q[2] = 0.0;
978
979 z *= 1.0/beta;
980
981 // the loop
982 int i = 1;
983 for( ; i<=_maxit; i++) {
984
985 dummy = z;
986 int i1 = i%3,
987 i0 = (i1+2)%3,
988 i2 = (i1+1)%3;
989
990 // symmetrically preconditioned Lanczos algorithm (see Greenbaum p.121)
991 _op.apply(z,q[i2]); // q[i2] = Az
992 q[i2].axpy(-beta,q[i0]);
993 // alpha is real since it is the diagonal entry of the hermitian tridiagonal matrix
994 // from the Lanczos Algorithm
995 // so the order in the scalar product doesn't matter even for the complex case
996 alpha = _sp.dot(z,q[i2]);
997 q[i2].axpy(-alpha,q[i1]);
998
999 z = 0.0;
1000 _prec.apply(z,q[i2]);
1001
1002 // beta is real and positive in exact arithmetic
1003 // since it is the norm of the basis vectors (in unpreconditioned case)
1004 beta = sqrt(_sp.dot(q[i2],z));
1005
1006 q[i2] *= 1.0/beta;
1007 z *= 1.0/beta;
1008
1009 // QR Factorization of recurrence coefficient matrix
1010 // apply previous givens rotations to last column of T
1011 T[1] = T[2];
1012 if(i>2) {
1013 T[0] = s[i%2]*T[1];
1014 T[1] = c[i%2]*T[1];
1015 }
1016 if(i>1) {
1017 T[2] = c[(i+1)%2]*alpha - s[(i+1)%2]*T[1];
1018 T[1] = c[(i+1)%2]*T[1] + s[(i+1)%2]*alpha;
1019 }
1020 else
1021 T[2] = alpha;
1022
1023 // update QR factorization
1024 generateGivensRotation(T[2],beta,c[i%2],s[i%2]);
1025 // to last column of T_k
1026 T[2] = c[i%2]*T[2] + s[i%2]*beta;
1027 // and to the rhs xi of the min problem
1028 xi[i%2] = -s[i%2]*xi[(i+1)%2];
1029 xi[(i+1)%2] *= c[i%2];
1030
1031 // compute correction direction
1032 p[i2] = dummy;
1033 p[i2].axpy(-T[1],p[i1]);
1034 p[i2].axpy(-T[0],p[i0]);
1035 p[i2] *= 1.0/T[2];
1036
1037 // apply correction/update solution
1038 x.axpy(beta0*xi[(i+1)%2],p[i2]);
1039
1040 // remember beta_old
1041 T[2] = beta;
1042
1043 // check for convergence
1044 // the last entry in the rhs of the min-problem is the residual
1045 real_type defnew = abs(beta0*xi[i%2]);
1046
1047 if(_verbose > 1)
1048 this->printOutput(std::cout,i,defnew,def);
1049
1050 def = defnew;
1051 if(def < def0*_reduction || def < 1e-30 || i == _maxit ) {
1052 res.converged = true;
1053 break;
1054 }
1055 } // end for
1056
1057 if(_verbose == 1)
1058 this->printOutput(std::cout,i,def);
1059
1060 // postprocess preconditioner
1061 _prec.post(x);
1062 // fill statistics
1063 res.iterations = i;
1064 res.reduction = static_cast<double>(def/def0);
1065 res.conv_rate = static_cast<double>(pow(res.reduction,1.0/i));
1066 res.elapsed = watch.elapsed();
1067
1068 // final print
1069 if(_verbose > 0) {
1070 std::cout << "=== rate=" << res.conv_rate
1071 << ", T=" << res.elapsed
1072 << ", TIT=" << res.elapsed/i
1073 << ", IT=" << i << std::endl;
1074 }
1075 }
1076
1082 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res)
1083 {
1084 real_type saved_reduction = _reduction;
1085 _reduction = reduction;
1086 (*this).apply(x,b,res);
1087 _reduction = saved_reduction;
1088 }
1089
1090 private:
1091
1092 void generateGivensRotation(field_type &dx, field_type &dy, real_type &cs, field_type &sn)
1093 {
1094 using std::sqrt;
1095 using std::abs;
1096 real_type norm_dx = abs(dx);
1097 real_type norm_dy = abs(dy);
1098 if(norm_dy < 1e-15) {
1099 cs = 1.0;
1100 sn = 0.0;
1101 } else if(norm_dx < 1e-15) {
1102 cs = 0.0;
1103 sn = 1.0;
1104 } else if(norm_dy > norm_dx) {
1105 real_type temp = norm_dx/norm_dy;
1106 cs = 1.0/sqrt(1.0 + temp*temp);
1107 sn = cs;
1108 cs *= temp;
1109 sn *= dx/norm_dx;
1110 // dy is real in exact arithmetic
1111 // so we don't need to conjugate here
1112 sn *= dy/norm_dy;
1113 } else {
1114 real_type temp = norm_dy/norm_dx;
1115 cs = 1.0/sqrt(1.0 + temp*temp);
1116 sn = cs;
1117 sn *= dy/dx;
1118 // dy and dx is real in exact arithmetic
1119 // so we don't have to conjugate both of them
1120 }
1121 }
1122
1123 SeqScalarProduct<X> ssp;
1124 LinearOperator<X,X>& _op;
1125 Preconditioner<X,X>& _prec;
1126 ScalarProduct<X>& _sp;
1127 real_type _reduction;
1128 int _maxit;
1129 int _verbose;
1130 };
1131
1145 template<class X, class Y=X, class F = Y>
1147 {
1148 public:
1150 typedef X domain_type;
1152 typedef Y range_type;
1154 typedef typename X::field_type field_type;
1156 typedef typename FieldTraits<field_type>::real_type real_type;
1158 typedef F basis_type;
1159
1160 template<class L, class P>
1161 DUNE_DEPRECATED_MSG("recalc_defect is a unused parameter! Use RestartedGMResSolver(L& op, P& prec, real_type reduction, int restart, int maxit, int verbose) instead")
1162 RestartedGMResSolver (L& op, P& prec, real_type reduction, int restart, int maxit, int verbose, bool recalc_defect)
1163 : _A(op)
1164 , _W(prec)
1165 , ssp()
1166 , _sp(ssp)
1167 , _restart(restart)
1168 , _reduction(reduction)
1169 , _maxit(maxit)
1170 , _verbose(verbose)
1171 {
1172 static_assert(static_cast<int>(P::category) == static_cast<int>(L::category),
1173 "P and L must be the same category!");
1174 static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
1175 "L must be sequential!");
1176 }
1177
1178
1185 template<class L, class P>
1186 RestartedGMResSolver (L& op, P& prec, real_type reduction, int restart, int maxit, int verbose) :
1187 _A(op), _W(prec),
1188 ssp(), _sp(ssp), _restart(restart),
1189 _reduction(reduction), _maxit(maxit), _verbose(verbose)
1190 {
1191 static_assert(static_cast<int>(P::category) == static_cast<int>(L::category),
1192 "P and L must be the same category!");
1193 static_assert(static_cast<int>(L::category) == static_cast<int>(SolverCategory::sequential),
1194 "L must be sequential!");
1195 }
1196
1197 template<class L, class S, class P>
1198 DUNE_DEPRECATED_MSG("recalc_defect is a unused parameter! Use RestartedGMResSolver(L& op, S& sp, P& prec, real_type reduction, int restart, int maxit, int verbose) instead")
1199 RestartedGMResSolver(L& op, S& sp, P& prec, real_type reduction, int restart, int maxit, int verbose, bool recalc_defect)
1200 : _A(op)
1201 , _W(prec)
1202 , _sp(sp)
1203 , _restart(restart)
1204 , _reduction(reduction)
1205 , _maxit(maxit)
1206 , _verbose(verbose)
1207 {
1208 static_assert(static_cast<int>(P::category) == static_cast<int>(L::category),
1209 " P and L must have the same category!");
1210 static_assert(static_cast<int>(P::category) == static_cast<int>(S::category),
1211 "P and S must have the same category!");
1212 }
1213
1220 template<class L, class S, class P>
1221 RestartedGMResSolver (L& op, S& sp, P& prec, real_type reduction, int restart, int maxit, int verbose) :
1222 _A(op), _W(prec),
1223 _sp(sp), _restart(restart),
1224 _reduction(reduction), _maxit(maxit), _verbose(verbose)
1225 {
1226 static_assert(static_cast<int>(P::category) == static_cast<int>(L::category),
1227 "P and L must have the same category!");
1228 static_assert(static_cast<int>(P::category) == static_cast<int>(S::category),
1229 "P and S must have the same category!");
1230 }
1231
1240 virtual void apply (X& x, Y& b, InverseOperatorResult& res)
1241 {
1242 apply(x,b,_reduction,res);
1243 }
1244
1253 virtual void apply (X& x, Y& b, double reduction, InverseOperatorResult& res)
1254 {
1255 using std::abs;
1256 const real_type EPSILON = 1e-80;
1257 const int m = _restart;
1258 real_type norm, norm_old = 0.0, norm_0;
1259 int j = 1;
1260 std::vector<field_type> s(m+1), sn(m);
1261 std::vector<real_type> cs(m);
1262 // need copy of rhs if GMRes has to be restarted
1263 Y b2(b);
1264 // helper vector
1265 Y w(b);
1266 std::vector< std::vector<field_type> > H(m+1,s);
1267 std::vector<F> v(m+1,b);
1268
1269 // start timer
1270 Dune::Timer watch;
1271 watch.reset();
1272
1273 // clear solver statistics and set res.converged to false
1274 res.clear();
1275 _W.pre(x,b);
1276
1277 // calculate defect and overwrite rhs with it
1278 _A.applyscaleadd(-1.0,x,b); // b -= Ax
1279 // calculate preconditioned defect
1280 v[0] = 0.0; _W.apply(v[0],b); // r = W^-1 b
1281 norm_0 = _sp.norm(v[0]);
1282 norm = norm_0;
1283 norm_old = norm;
1284
1285 // print header
1286 if(_verbose > 0)
1287 {
1288 std::cout << "=== RestartedGMResSolver" << std::endl;
1289 if(_verbose > 1) {
1290 this->printHeader(std::cout);
1291 this->printOutput(std::cout,0,norm_0);
1292 }
1293 }
1294
1295 if(norm_0 < EPSILON) {
1296 _W.post(x);
1297 res.converged = true;
1298 if(_verbose > 0) // final print
1299 print_result(res);
1300 }
1301
1302 while(j <= _maxit && res.converged != true) {
1303
1304 int i = 0;
1305 v[0] *= 1.0/norm;
1306 s[0] = norm;
1307 for(i=1; i<m+1; i++)
1308 s[i] = 0.0;
1309
1310 for(i=0; i < m && j <= _maxit && res.converged != true; i++, j++) {
1311 w = 0.0;
1312 // use v[i+1] as temporary vector
1313 v[i+1] = 0.0;
1314 // do Arnoldi algorithm
1315 _A.apply(v[i],v[i+1]);
1316 _W.apply(w,v[i+1]);
1317 for(int k=0; k<i+1; k++) {
1318 // notice that _sp.dot(v[k],w) = v[k]\adjoint w
1319 // so one has to pay attention to the order
1320 // in the scalar product for the complex case
1321 // doing the modified Gram-Schmidt algorithm
1322 H[k][i] = _sp.dot(v[k],w);
1323 // w -= H[k][i] * v[k]
1324 w.axpy(-H[k][i],v[k]);
1325 }
1326 H[i+1][i] = _sp.norm(w);
1327 if(abs(H[i+1][i]) < EPSILON)
1329 "breakdown in GMRes - |w| == 0.0 after " << j << " iterations");
1330
1331 // normalize new vector
1332 v[i+1] = w; v[i+1] *= 1.0/H[i+1][i];
1333
1334 // update QR factorization
1335 for(int k=0; k<i; k++)
1336 applyPlaneRotation(H[k][i],H[k+1][i],cs[k],sn[k]);
1337
1338 // compute new givens rotation
1339 generatePlaneRotation(H[i][i],H[i+1][i],cs[i],sn[i]);
1340 // finish updating QR factorization
1341 applyPlaneRotation(H[i][i],H[i+1][i],cs[i],sn[i]);
1342 applyPlaneRotation(s[i],s[i+1],cs[i],sn[i]);
1343
1344 // norm of the defect is the last component the vector s
1345 norm = abs(s[i+1]);
1346
1347 // print current iteration statistics
1348 if(_verbose > 1) {
1349 this->printOutput(std::cout,j,norm,norm_old);
1350 }
1351
1352 norm_old = norm;
1353
1354 // check convergence
1355 if(norm < reduction * norm_0)
1356 res.converged = true;
1357
1358 } // end for
1359
1360 // calculate update vector
1361 w = 0.0;
1362 update(w,i,H,s,v);
1363 // and current iterate
1364 x += w;
1365
1366 // restart GMRes if convergence was not achieved,
1367 // i.e. linear defect has not reached desired reduction
1368 // and if j < _maxit
1369 if( res.converged != true && j <= _maxit ) {
1370
1371 if(_verbose > 0)
1372 std::cout << "=== GMRes::restart" << std::endl;
1373 // get saved rhs
1374 b = b2;
1375 // calculate new defect
1376 _A.applyscaleadd(-1.0,x,b); // b -= Ax;
1377 // calculate preconditioned defect
1378 v[0] = 0.0;
1379 _W.apply(v[0],b);
1380 norm = _sp.norm(v[0]);
1381 norm_old = norm;
1382 }
1383
1384 } //end while
1385
1386 // postprocess preconditioner
1387 _W.post(x);
1388
1389 // save solver statistics
1390 res.iterations = j-1; // it has to be j-1!!!
1391 res.reduction = static_cast<double>(norm/norm_0);
1392 res.conv_rate = static_cast<double>(pow(res.reduction,1.0/(j-1)));
1393 res.elapsed = watch.elapsed();
1394
1395 if(_verbose>0)
1396 print_result(res);
1397
1398 }
1399
1400 private :
1401
1402 void print_result(const InverseOperatorResult& res) const {
1403 int k = res.iterations>0 ? res.iterations : 1;
1404 std::cout << "=== rate=" << res.conv_rate
1405 << ", T=" << res.elapsed
1406 << ", TIT=" << res.elapsed/k
1407 << ", IT=" << res.iterations
1408 << std::endl;
1409 }
1410
1411 void update(X& w, int i,
1412 const std::vector<std::vector<field_type> >& H,
1413 const std::vector<field_type>& s,
1414 const std::vector<X>& v) {
1415 // solution vector of the upper triangular system
1416 std::vector<field_type> y(s);
1417
1418 // backsolve
1419 for(int a=i-1; a>=0; a--) {
1420 field_type rhs(s[a]);
1421 for(int b=a+1; b<i; b++)
1422 rhs -= H[a][b]*y[b];
1423 y[a] = rhs/H[a][a];
1424
1425 // compute update on the fly
1426 // w += y[a]*v[a]
1427 w.axpy(y[a],v[a]);
1428 }
1429 }
1430
1431 template<typename T>
1432 typename std::enable_if<std::is_same<field_type,real_type>::value,T>::type conjugate(const T& t) {
1433 return t;
1434 }
1435
1436 template<typename T>
1437 typename std::enable_if<!std::is_same<field_type,real_type>::value,T>::type conjugate(const T& t) {
1438 return conj(t);
1439 }
1440
1441 void
1442 generatePlaneRotation(field_type &dx, field_type &dy, real_type &cs, field_type &sn)
1443 {
1444 using std::sqrt;
1445 using std::abs;
1446 real_type norm_dx = abs(dx);
1447 real_type norm_dy = abs(dy);
1448 if(norm_dy < 1e-15) {
1449 cs = 1.0;
1450 sn = 0.0;
1451 } else if(norm_dx < 1e-15) {
1452 cs = 0.0;
1453 sn = 1.0;
1454 } else if(norm_dy > norm_dx) {
1455 real_type temp = norm_dx/norm_dy;
1456 cs = 1.0/sqrt(1.0 + temp*temp);
1457 sn = cs;
1458 cs *= temp;
1459 sn *= dx/norm_dx;
1460 sn *= conjugate(dy)/norm_dy;
1461 } else {
1462 real_type temp = norm_dy/norm_dx;
1463 cs = 1.0/sqrt(1.0 + temp*temp);
1464 sn = cs;
1465 sn *= conjugate(dy/dx);
1466 }
1467 }
1468
1469
1470 void
1471 applyPlaneRotation(field_type &dx, field_type &dy, real_type &cs, field_type &sn)
1472 {
1473 field_type temp = cs * dx + sn * dy;
1474 dy = -conjugate(sn) * dx + cs * dy;
1475 dx = temp;
1476 }
1477
1478 LinearOperator<X,Y>& _A;
1479 Preconditioner<X,Y>& _W;
1480 SeqScalarProduct<X> ssp;
1481 ScalarProduct<X>& _sp;
1482 int _restart;
1483 real_type _reduction;
1484 int _maxit;
1485 int _verbose;
1486 };
1487
1488
1502 template<class X>
1504 {
1505 public:
1507 typedef X domain_type;
1509 typedef X range_type;
1511 typedef typename X::field_type field_type;
1513 typedef typename FieldTraits<field_type>::real_type real_type;
1514
1522 template<class L, class P>
1523 GeneralizedPCGSolver (L& op, P& prec, real_type reduction, int maxit, int verbose,
1524 int restart=10) :
1525 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit),
1526 _verbose(verbose), _restart(std::min(maxit,restart))
1527 {
1528 static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
1529 "L and P have to have the same category!");
1530 static_assert(static_cast<int>(L::category) ==
1531 static_cast<int>(SolverCategory::sequential),
1532 "L has to be sequential!");
1533 }
1541 template<class L, class P, class S>
1542 GeneralizedPCGSolver (L& op, S& sp, P& prec,
1543 real_type reduction, int maxit, int verbose, int restart=10) :
1544 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose),
1545 _restart(std::min(maxit,restart))
1546 {
1547 static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
1548 "L and P must have the same category!");
1549 static_assert(static_cast<int>(L::category) == static_cast<int>(S::category),
1550 "L and S must have the same category!");
1551 }
1557 virtual void apply (X& x, X& b, InverseOperatorResult& res)
1558 {
1559 res.clear(); // clear solver statistics
1560 Timer watch; // start a timer
1561 _prec.pre(x,b); // prepare preconditioner
1562 _op.applyscaleadd(-1,x,b); // overwrite b with defect
1563
1564 std::vector<std::shared_ptr<X> > p(_restart);
1565 std::vector<typename X::field_type> pp(_restart);
1566 X q(x); // a temporary vector
1567 X prec_res(x); // a temporary vector for preconditioner output
1568
1569 p[0].reset(new X(x));
1570
1571 real_type def0 = _sp.norm(b); // compute norm
1572 if (def0<1E-30) // convergence check
1573 {
1574 res.converged = true;
1575 res.iterations = 0; // fill statistics
1576 res.reduction = 0;
1577 res.conv_rate = 0;
1578 res.elapsed=0;
1579 if (_verbose>0) // final print
1580 std::cout << "=== rate=" << res.conv_rate
1581 << ", T=" << res.elapsed << ", TIT=" << res.elapsed
1582 << ", IT=0" << std::endl;
1583 return;
1584 }
1585
1586 if (_verbose>0) // printing
1587 {
1588 std::cout << "=== GeneralizedPCGSolver" << std::endl;
1589 if (_verbose>1) {
1590 this->printHeader(std::cout);
1591 this->printOutput(std::cout,0,def0);
1592 }
1593 }
1594 // some local variables
1595 real_type def=def0; // loop variables
1596 field_type rho, lambda;
1597
1598 int i=0;
1599 int ii=0;
1600 // determine initial search direction
1601 *(p[0]) = 0; // clear correction
1602 _prec.apply(*(p[0]),b); // apply preconditioner
1603 rho = _sp.dot(*(p[0]),b); // orthogonalization
1604 _op.apply(*(p[0]),q); // q=Ap
1605 pp[0] = _sp.dot(*(p[0]),q); // scalar product
1606 lambda = rho/pp[0]; // minimization
1607 x.axpy(lambda,*(p[0])); // update solution
1608 b.axpy(-lambda,q); // update defect
1609
1610 // convergence test
1611 real_type defnew=_sp.norm(b); // comp defect norm
1612 if (_verbose>1) // print
1613 this->printOutput(std::cout,++i,defnew,def);
1614 def = defnew; // update norm
1615 if (def<def0*_reduction || def<1E-30) // convergence check
1616 {
1617 res.converged = true;
1618 if (_verbose>0) // final print
1619 {
1620 std::cout << "=== rate=" << res.conv_rate
1621 << ", T=" << res.elapsed
1622 << ", TIT=" << res.elapsed
1623 << ", IT=" << 1 << std::endl;
1624 }
1625 return;
1626 }
1627
1628 while(i<_maxit) {
1629 // the loop
1630 int end=std::min(_restart, _maxit-i+1);
1631 for (ii=1; ii<end; ++ii )
1632 {
1633 //std::cout<<" ii="<<ii<<" i="<<i<<std::endl;
1634 // compute next conjugate direction
1635 prec_res = 0; // clear correction
1636 _prec.apply(prec_res,b); // apply preconditioner
1637
1638 p[ii].reset(new X(prec_res));
1639 _op.apply(prec_res, q);
1640
1641 for(int j=0; j<ii; ++j) {
1642 rho =_sp.dot(q,*(p[j]))/pp[j];
1643 p[ii]->axpy(-rho, *(p[j]));
1644 }
1645
1646 // minimize in given search direction
1647 _op.apply(*(p[ii]),q); // q=Ap
1648 pp[ii] = _sp.dot(*(p[ii]),q); // scalar product
1649 rho = _sp.dot(*(p[ii]),b); // orthogonalization
1650 lambda = rho/pp[ii]; // minimization
1651 x.axpy(lambda,*(p[ii])); // update solution
1652 b.axpy(-lambda,q); // update defect
1653
1654 // convergence test
1655 real_type defNew=_sp.norm(b); // comp defect norm
1656
1657 if (_verbose>1) // print
1658 this->printOutput(std::cout,++i,defNew,def);
1659
1660 def = defNew; // update norm
1661 if (def<def0*_reduction || def<1E-30) // convergence check
1662 {
1663 res.converged = true;
1664 break;
1665 }
1666 }
1667 if(res.converged)
1668 break;
1669 if(end==_restart) {
1670 *(p[0])=*(p[_restart-1]);
1671 pp[0]=pp[_restart-1];
1672 }
1673 }
1674
1675 // postprocess preconditioner
1676 _prec.post(x);
1677
1678 // fill statistics
1679 res.iterations = i;
1680 res.reduction = def/def0;
1681 res.conv_rate = pow(res.reduction,1.0/i);
1682 res.elapsed = watch.elapsed();
1683
1684 if (_verbose>0) // final print
1685 {
1686 std::cout << "=== rate=" << res.conv_rate
1687 << ", T=" << res.elapsed
1688 << ", TIT=" << res.elapsed/i
1689 << ", IT=" << i+1 << std::endl;
1690 }
1691 }
1692
1698 virtual void apply (X& x, X& b, double reduction,
1700 {
1701 real_type saved_reduction = _reduction;
1702 _reduction = reduction;
1703 (*this).apply(x,b,res);
1704 _reduction = saved_reduction;
1705 }
1706 private:
1709 Preconditioner<X,X>& _prec;
1710 ScalarProduct<X>& _sp;
1711 real_type _reduction;
1712 int _maxit;
1713 int _verbose;
1714 int _restart;
1715 };
1716
1719} // end namespace
1720
1721#endif
Bi-conjugate Gradient Stabilized (BiCG-STAB)
Definition: solvers.hh:574
BiCGSTABSolver(L &op, S &sp, P &prec, real_type reduction, int maxit, int verbose)
Set up solver.
Definition: solvers.hh:606
X domain_type
The domain type of the operator to be inverted.
Definition: solvers.hh:577
FieldTraits< field_type >::real_type real_type
The real type of the field type (is the same if using real numbers, but differs for std::complex)
Definition: solvers.hh:583
X range_type
The range type of the operator to be inverted.
Definition: solvers.hh:579
virtual void apply(X &x, X &b, InverseOperatorResult &res)
Apply inverse operator.
Definition: solvers.hh:623
X::field_type field_type
The field type of the operator to be inverted.
Definition: solvers.hh:581
virtual void apply(X &x, X &b, double reduction, InverseOperatorResult &res)
Apply inverse operator with given reduction factor.
Definition: solvers.hh:827
BiCGSTABSolver(L &op, P &prec, real_type reduction, int maxit, int verbose)
Set up solver.
Definition: solvers.hh:591
conjugate gradient method
Definition: solvers.hh:370
X range_type
The range type of the operator to be inverted.
Definition: solvers.hh:375
X domain_type
The domain type of the operator to be inverted.
Definition: solvers.hh:373
CGSolver(L &op, S &sp, P &prec, real_type reduction, int maxit, int verbose)
Set up conjugate gradient solver.
Definition: solvers.hh:401
virtual void apply(X &x, X &b, InverseOperatorResult &res)
Apply inverse operator.
Definition: solvers.hh:421
virtual void apply(X &x, X &b, double reduction, InverseOperatorResult &res)
Apply inverse operator with given reduction factor.
Definition: solvers.hh:551
CGSolver(L &op, P &prec, real_type reduction, int maxit, int verbose)
Set up conjugate gradient solver.
Definition: solvers.hh:387
FieldTraits< field_type >::real_type real_type
The real type of the field type (is the same if using real numbers, but differs for std::complex)
Definition: solvers.hh:379
X::field_type field_type
The field type of the operator to be inverted.
Definition: solvers.hh:377
Generalized preconditioned conjugate gradient solver.
Definition: solvers.hh:1504
FieldTraits< field_type >::real_type real_type
The real type of the field type (is the same if using real numbers, but differs for std::complex)
Definition: solvers.hh:1513
GeneralizedPCGSolver(L &op, P &prec, real_type reduction, int maxit, int verbose, int restart=10)
Set up nonlinear preconditioned conjugate gradient solver.
Definition: solvers.hh:1523
virtual void apply(X &x, X &b, double reduction, InverseOperatorResult &res)
Apply inverse operator with given reduction factor.
Definition: solvers.hh:1698
virtual void apply(X &x, X &b, InverseOperatorResult &res)
Apply inverse operator.
Definition: solvers.hh:1557
GeneralizedPCGSolver(L &op, S &sp, P &prec, real_type reduction, int maxit, int verbose, int restart=10)
Set up nonlinear preconditioned conjugate gradient solver.
Definition: solvers.hh:1542
X::field_type field_type
The field type of the operator to be inverted.
Definition: solvers.hh:1511
X range_type
The range type of the operator to be inverted.
Definition: solvers.hh:1509
X domain_type
The domain type of the operator to be inverted.
Definition: solvers.hh:1507
gradient method
Definition: solvers.hh:232
X domain_type
The domain type of the operator that we do the inverse for.
Definition: solvers.hh:235
virtual void apply(X &x, X &b, InverseOperatorResult &res)
Apply inverse operator.
Definition: solvers.hh:280
GradientSolver(L &op, P &prec, real_type reduction, int maxit, int verbose)
Set up solver.
Definition: solvers.hh:250
GradientSolver(L &op, S &sp, P &prec, real_type reduction, int maxit, int verbose)
Set up solver.
Definition: solvers.hh:265
X range_type
The range type of the operator that we do the inverse for.
Definition: solvers.hh:237
X::field_type field_type
The field type of the operator that we do the inverse for.
Definition: solvers.hh:239
FieldTraits< field_type >::real_type real_type
The real type of the field type (is the same if using real numbers, but differs for std::complex)
Definition: solvers.hh:241
virtual void apply(X &x, X &b, double reduction, InverseOperatorResult &res)
Apply inverse operator with given reduction factor.
Definition: solvers.hh:348
Abstract base class for all solvers.
Definition: solver.hh:79
void printHeader(std::ostream &s) const
helper function for printing header of solver output
Definition: solver.hh:127
void printOutput(std::ostream &s, const CountType &iter, const DataType &norm, const DataType &norm_old) const
helper function for printing solver output
Definition: solver.hh:136
virtual void applyscaleadd(field_type alpha, const X &x, Y &y) const =0
apply operator to x, scale and add:
virtual void apply(const X &x, Y &y) const =0
apply operator to x: The input vector is consistent and the output must also be consistent on the in...
Preconditioned loop solver.
Definition: solvers.hh:58
X range_type
The range type of the operator that we do the inverse for.
Definition: solvers.hh:63
X domain_type
The domain type of the operator that we do the inverse for.
Definition: solvers.hh:61
FieldTraits< field_type >::real_type real_type
The real type of the field type (is the same if using real numbers, but differs for std::complex)
Definition: solvers.hh:67
LoopSolver(L &op, P &prec, real_type reduction, int maxit, int verbose)
Set up Loop solver.
Definition: solvers.hh:89
virtual void apply(X &x, X &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: solvers.hh:132
virtual void apply(X &x, X &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: solvers.hh:210
LoopSolver(L &op, S &sp, P &prec, real_type reduction, int maxit, int verbose)
Set up loop solver.
Definition: solvers.hh:120
X::field_type field_type
The field type of the operator that we do the inverse for.
Definition: solvers.hh:65
Minimal Residual Method (MINRES)
Definition: solvers.hh:852
virtual void apply(X &x, X &b, InverseOperatorResult &res)
Apply inverse operator.
Definition: solvers.hh:897
X::field_type field_type
The field type of the operator to be inverted.
Definition: solvers.hh:859
MINRESSolver(L &op, S &sp, P &prec, real_type reduction, int maxit, int verbose)
Set up MINRES solver.
Definition: solvers.hh:883
MINRESSolver(L &op, P &prec, real_type reduction, int maxit, int verbose)
Set up MINRES solver.
Definition: solvers.hh:869
FieldTraits< field_type >::real_type real_type
The real type of the field type (is the same if using real numbers, but differs for std::complex)
Definition: solvers.hh:861
X range_type
The range type of the operator to be inverted.
Definition: solvers.hh:857
virtual void apply(X &x, X &b, double reduction, InverseOperatorResult &res)
Apply inverse operator with given reduction factor.
Definition: solvers.hh:1082
X domain_type
The domain type of the operator to be inverted.
Definition: solvers.hh:855
virtual void post(X &x)=0
Clean up.
virtual void apply(X &v, const Y &d)=0
Apply one step of the preconditioner to the system A(v)=d.
virtual void pre(X &x, Y &b)=0
Prepare the preconditioner.
implements the Generalized Minimal Residual (GMRes) method
Definition: solvers.hh:1147
F basis_type
The field type of the basis vectors.
Definition: solvers.hh:1158
RestartedGMResSolver(L &op, P &prec, real_type reduction, int restart, int maxit, int verbose)
Set up solver.
Definition: solvers.hh:1186
Y range_type
The range type of the operator to be inverted.
Definition: solvers.hh:1152
virtual void apply(X &x, Y &b, double reduction, InverseOperatorResult &res)
Apply inverse operator.
Definition: solvers.hh:1253
X::field_type field_type
The field type of the operator to be inverted.
Definition: solvers.hh:1154
FieldTraits< field_type >::real_type real_type
The real type of the field type (is the same if using real numbers, but differs for std::complex)
Definition: solvers.hh:1156
X domain_type
The domain type of the operator to be inverted.
Definition: solvers.hh:1150
virtual void apply(X &x, Y &b, InverseOperatorResult &res)
Apply inverse operator.
Definition: solvers.hh:1240
RestartedGMResSolver(L &op, S &sp, P &prec, real_type reduction, int restart, int maxit, int verbose)
Set up solver.
Definition: solvers.hh:1221
Base class for scalar product and norm computation.
Definition: scalarproducts.hh:44
Default implementation for the scalar case.
Definition: scalarproducts.hh:96
Thrown when a solver aborts due to some problem.
Definition: istlexception.hh:43
A simple stop watch.
Definition: timer.hh:52
void reset()
Reset timer while keeping the running/stopped state.
Definition: timer.hh:66
double elapsed() const
Get elapsed user-time from last reset until now/last stop in seconds.
Definition: timer.hh:86
Definition of the DUNE_DEPRECATED macro for the case that config.h is not available.
A few common exception classes.
Type traits to determine the type of reals (when working with complex numbers)
struct DUNE_DEPRECATED_MSG("Use <type_traits> instead!") ConstantVolatileTraits
Determines whether a type is const or volatile and provides the unqualified types.
Definition: typetraits.hh:57
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
Dune namespace.
Definition: alignment.hh:11
STL namespace.
Define general, extensible interface for operators. The available implementation wraps a matrix.
Define base class for scalar product and norm.
Define general, extensible interface for inverse operators.
Statistics about the application of an inverse operator.
Definition: solver.hh:32
double elapsed
Elapsed time in seconds.
Definition: solver.hh:62
int iterations
Number of iterations.
Definition: solver.hh:50
double reduction
Reduction achieved: .
Definition: solver.hh:53
void clear()
Resets all data.
Definition: solver.hh:40
double conv_rate
Convergence rate (average reduction per step)
Definition: solver.hh:59
bool converged
True if convergence criterion has been met.
Definition: solver.hh:56
@ sequential
Category for sequential solvers.
Definition: solvercategory.hh:21
A simple timing class.
Traits for type conversions and type information.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)