Dune Core Modules (unstable)

solvers.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4// vi: set et ts=4 sw=2 sts=2:
5
6#ifndef DUNE_ISTL_SOLVERS_HH
7#define DUNE_ISTL_SOLVERS_HH
8
9#include <array>
10#include <cmath>
11#include <complex>
12#include <iostream>
13#include <memory>
14#include <type_traits>
15#include <vector>
16
18#include <dune/common/math.hh>
21#include <dune/common/std/type_traits.hh>
22#include <dune/common/timer.hh>
23
24#include <dune/istl/allocator.hh>
26#include <dune/istl/eigenvalue/arpackpp.hh>
27#include <dune/istl/istlexception.hh>
29#include <dune/istl/preconditioner.hh>
31#include <dune/istl/solver.hh>
32#include <dune/istl/solverregistry.hh>
33
34namespace Dune {
46 //=====================================================================
47 // Implementation of this interface
48 //=====================================================================
49
58 template<class X>
59 class LoopSolver : public IterativeSolver<X,X> {
60 public:
65
66 // copy base class constructors
68
69 // don't shadow four-argument version of apply defined in the base class
70 using IterativeSolver<X,X>::apply;
71
73 void apply (X& x, X& b, InverseOperatorResult& res) override
74 {
75 Iteration iteration(*this, res);
76 _prec->pre(x,b);
77
78 // overwrite b with defect
79 _op->applyscaleadd(-1,x,b);
80
81 // compute norm, \todo parallelization
82 real_type def = _sp->norm(b);
83 if(iteration.step(0, def)){
84 _prec->post(x);
85 return;
86 }
87 // prepare preconditioner
88
89 // allocate correction vector
90 X v(x);
91
92 // iteration loop
93 int i=1;
94 for ( ; i<=_maxit; i++ )
95 {
96 v = 0; // clear correction
97 _prec->apply(v,b); // apply preconditioner
98 x += v; // update solution
99 _op->applyscaleadd(-1,v,b); // update defect
100 def=_sp->norm(b); // comp defect norm
101 if(iteration.step(i, def))
102 break;
103 }
104
105 // postprocess preconditioner
106 _prec->post(x);
107 }
108
109 protected:
110 using IterativeSolver<X,X>::_op;
111 using IterativeSolver<X,X>::_prec;
112 using IterativeSolver<X,X>::_sp;
113 using IterativeSolver<X,X>::_reduction;
114 using IterativeSolver<X,X>::_maxit;
115 using IterativeSolver<X,X>::_verbose;
116 using Iteration = typename IterativeSolver<X,X>::template Iteration<unsigned int>;
117 };
118 DUNE_REGISTER_SOLVER("loopsolver", defaultIterativeSolverCreator<Dune::LoopSolver>());
119
120
121 // all these solvers are taken from the SUMO library
123 template<class X>
124 class GradientSolver : public IterativeSolver<X,X> {
125 public:
129 using typename IterativeSolver<X,X>::real_type;
130
131 // copy base class constructors
133
134 // don't shadow four-argument version of apply defined in the base class
135 using IterativeSolver<X,X>::apply;
136
142 void apply (X& x, X& b, InverseOperatorResult& res) override
143 {
144 Iteration iteration(*this, res);
145 _prec->pre(x,b); // prepare preconditioner
146
147 _op->applyscaleadd(-1,x,b); // overwrite b with defec
148
149 real_type def = _sp->norm(b); // compute norm
150 if(iteration.step(0, def)){
151 _prec->post(x);
152 return;
153 }
154
155 X p(x); // create local vectors
156 X q(b);
157
158 int i=1; // loop variables
159 field_type lambda;
160 for ( ; i<=_maxit; i++ )
161 {
162 p = 0; // clear correction
163 _prec->apply(p,b); // apply preconditioner
164 _op->apply(p,q); // q=Ap
165 auto alpha = _sp->dot(q,p);
166 lambda = Simd::cond(def==field_type(0.),
167 field_type(0.), // no need for minimization if def is already 0
168 field_type(_sp->dot(p,b)/alpha)); // minimization
169 x.axpy(lambda,p); // update solution
170 b.axpy(-lambda,q); // update defect
171
172 def =_sp->norm(b); // comp defect norm
173 if(iteration.step(i, def))
174 break;
175 }
176 // postprocess preconditioner
177 _prec->post(x);
178 }
179
180 protected:
181 using IterativeSolver<X,X>::_op;
182 using IterativeSolver<X,X>::_prec;
183 using IterativeSolver<X,X>::_sp;
184 using IterativeSolver<X,X>::_reduction;
185 using IterativeSolver<X,X>::_maxit;
186 using IterativeSolver<X,X>::_verbose;
187 using Iteration = typename IterativeSolver<X,X>::template Iteration<unsigned int>;
188 };
189 DUNE_REGISTER_SOLVER("gradientsolver", defaultIterativeSolverCreator<Dune::GradientSolver>());
190
192 template<class X>
193 class CGSolver : public IterativeSolver<X,X> {
194 public:
198 using typename IterativeSolver<X,X>::real_type;
199
200 // copy base class constructors
202
203 private:
205
206 protected:
207
208 static constexpr bool enableConditionEstimate = (std::is_same_v<field_type,float> || std::is_same_v<field_type,double>);
209
210 public:
211
212 // don't shadow four-argument version of apply defined in the base class
213 using IterativeSolver<X,X>::apply;
214
223 scalar_real_type reduction, int maxit, int verbose, bool condition_estimate) : IterativeSolver<X,X>(op, prec, reduction, maxit, verbose),
224 condition_estimate_(condition_estimate)
225 {
226 if (condition_estimate && !enableConditionEstimate) {
227 condition_estimate_ = false;
228 std::cerr << "WARNING: Condition estimate was disabled. It is only available for double and float field types!" << std::endl;
229 }
230 }
231
240 scalar_real_type reduction, int maxit, int verbose, bool condition_estimate) : IterativeSolver<X,X>(op, sp, prec, reduction, maxit, verbose),
241 condition_estimate_(condition_estimate)
242 {
243 if (condition_estimate && !(std::is_same<field_type,float>::value || std::is_same<field_type,double>::value)) {
244 condition_estimate_ = false;
245 std::cerr << "WARNING: Condition estimate was disabled. It is only available for double and float field types!" << std::endl;
246 }
247 }
248
257 CGSolver (std::shared_ptr<const LinearOperator<X,X>> op, std::shared_ptr<ScalarProduct<X>> sp,
258 std::shared_ptr<Preconditioner<X,X>> prec,
259 scalar_real_type reduction, int maxit, int verbose, bool condition_estimate)
260 : IterativeSolver<X,X>(op, sp, prec, reduction, maxit, verbose),
261 condition_estimate_(condition_estimate)
262 {
263 if (condition_estimate && !(std::is_same<field_type,float>::value || std::is_same<field_type,double>::value)) {
264 condition_estimate_ = false;
265 std::cerr << "WARNING: Condition estimate was disabled. It is only available for double and float field types!" << std::endl;
266 }
267 }
268
280 void apply (X& x, X& b, InverseOperatorResult& res) override
281 {
282 Iteration iteration(*this,res);
283 _prec->pre(x,b); // prepare preconditioner
284
285 _op->applyscaleadd(-1,x,b); // overwrite b with defect
286
287 real_type def = _sp->norm(b); // compute norm
288 if(iteration.step(0, def)){
289 _prec->post(x);
290 return;
291 }
292
293 X p(x); // the search direction
294 X q(x); // a temporary vector
295
296 // Remember lambda and beta values for condition estimate
297 std::vector<real_type> lambdas(0);
298 std::vector<real_type> betas(0);
299
300 // some local variables
301 field_type rho,rholast,lambda,alpha,beta;
302
303 // determine initial search direction
304 p = 0; // clear correction
305 _prec->apply(p,b); // apply preconditioner
306 rholast = _sp->dot(p,b); // orthogonalization
307
308 // the loop
309 int i=1;
310 for ( ; i<=_maxit; i++ )
311 {
312 // minimize in given search direction p
313 _op->apply(p,q); // q=Ap
314 alpha = _sp->dot(p,q); // scalar product
315 lambda = Simd::cond(def==field_type(0.), field_type(0.), field_type(rholast/alpha)); // minimization
316 if constexpr (enableConditionEstimate)
317 if (condition_estimate_)
318 lambdas.push_back(std::real(lambda));
319 x.axpy(lambda,p); // update solution
320 b.axpy(-lambda,q); // update defect
321
322 // convergence test
323 def=_sp->norm(b); // comp defect norm
324 if(iteration.step(i, def))
325 break;
326
327 // determine new search direction
328 q = 0; // clear correction
329 _prec->apply(q,b); // apply preconditioner
330 rho = _sp->dot(q,b); // orthogonalization
331 beta = Simd::cond(def==field_type(0.), field_type(0.), field_type(rho/rholast)); // scaling factor
332 if constexpr (enableConditionEstimate)
333 if (condition_estimate_)
334 betas.push_back(std::real(beta));
335 p *= beta; // scale old search direction
336 p += q; // orthogonalization with correction
337 rholast = rho; // remember rho for recurrence
338 }
339
340 _prec->post(x); // postprocess preconditioner
341
342 if (condition_estimate_) {
343#if HAVE_ARPACKPP
344 if constexpr (enableConditionEstimate) {
345 using std::sqrt;
346
347 // Build T matrix which has extreme eigenvalues approximating
348 // those of the original system
349 // (see Y. Saad, Iterative methods for sparse linear systems)
350
352
353 for (auto row = T.createbegin(); row != T.createend(); ++row) {
354 if (row.index() > 0)
355 row.insert(row.index()-1);
356 row.insert(row.index());
357 if (row.index() < T.N() - 1)
358 row.insert(row.index()+1);
359 }
360 for (int row = 0; row < i; ++row) {
361 if (row > 0) {
362 T[row][row-1] = sqrt(betas[row-1]) / lambdas[row-1];
363 }
364
365 T[row][row] = 1.0 / lambdas[row];
366 if (row > 0) {
367 T[row][row] += betas[row-1] / lambdas[row-1];
368 }
369
370 if (row < i - 1) {
371 T[row][row+1] = sqrt(betas[row]) / lambdas[row];
372 }
373 }
374
375 // Compute largest and smallest eigenvalue of T matrix and return as estimate
377
378 real_type eps = 0.0;
379 COND_VEC eigv;
380 real_type min_eigv, max_eigv;
381 arpack.computeSymMinMagnitude (eps, eigv, min_eigv);
382 arpack.computeSymMaxMagnitude (eps, eigv, max_eigv);
383
384 res.condition_estimate = max_eigv / min_eigv;
385
386 if (this->_verbose > 0) {
387 std::cout << "Min eigv estimate: " << Simd::io(min_eigv) << '\n';
388 std::cout << "Max eigv estimate: " << Simd::io(max_eigv) << '\n';
389 std::cout << "Condition estimate: "
390 << Simd::io(max_eigv / min_eigv) << std::endl;
391 }
392 }
393#else
394 std::cerr << "WARNING: Condition estimate was requested. This requires ARPACK, but ARPACK was not found!" << std::endl;
395#endif
396 }
397 }
398
399 private:
400 bool condition_estimate_ = false;
401
402 // Matrix and vector types used for condition estimate
405
406 protected:
407 using IterativeSolver<X,X>::_op;
408 using IterativeSolver<X,X>::_prec;
409 using IterativeSolver<X,X>::_sp;
410 using IterativeSolver<X,X>::_reduction;
411 using IterativeSolver<X,X>::_maxit;
412 using IterativeSolver<X,X>::_verbose;
413 using Iteration = typename IterativeSolver<X,X>::template Iteration<unsigned int>;
414 };
415 DUNE_REGISTER_SOLVER("cgsolver", defaultIterativeSolverCreator<Dune::CGSolver>());
416
417 // Ronald Kriemanns BiCG-STAB implementation from Sumo
419 template<class X>
420 class BiCGSTABSolver : public IterativeSolver<X,X> {
421 public:
425 using typename IterativeSolver<X,X>::real_type;
426
427 // copy base class constructors
429
430 // don't shadow four-argument version of apply defined in the base class
431 using IterativeSolver<X,X>::apply;
432
440 void apply (X& x, X& b, InverseOperatorResult& res) override
441 {
442 using std::abs;
443 const Simd::Scalar<real_type> EPSILON=1e-80;
444 using std::abs;
445 double it;
446 field_type rho, rho_new, alpha, beta, h, omega;
447 real_type norm;
448
449 //
450 // get vectors and matrix
451 //
452 X& r=b;
453 X p(x);
454 X v(x);
455 X t(x);
456 X y(x);
457 X rt(x);
458
459 //
460 // begin iteration
461 //
462
463 // r = r - Ax; rt = r
464 Iteration<double> iteration(*this,res);
465 _prec->pre(x,r); // prepare preconditioner
466
467 _op->applyscaleadd(-1,x,r); // overwrite b with defect
468
469 rt=r;
470
471 norm = _sp->norm(r);
472 if(iteration.step(0, norm)){
473 _prec->post(x);
474 return;
475 }
476 p=0;
477 v=0;
478
479 rho = 1;
480 alpha = 1;
481 omega = 1;
482
483 //
484 // iteration
485 //
486
487 for (it = 0.5; it < _maxit; it+=.5)
488 {
489 //
490 // preprocess, set vecsizes etc.
491 //
492
493 // rho_new = < rt , r >
494 rho_new = _sp->dot(rt,r);
495
496 // look if breakdown occurred
497 if (Simd::allTrue(abs(rho) <= EPSILON))
498 DUNE_THROW(SolverAbort,"breakdown in BiCGSTAB - rho "
499 << Simd::io(rho) << " <= EPSILON " << EPSILON
500 << " after " << it << " iterations");
501 if (Simd::allTrue(abs(omega) <= EPSILON))
502 DUNE_THROW(SolverAbort,"breakdown in BiCGSTAB - omega "
503 << Simd::io(omega) << " <= EPSILON " << EPSILON
504 << " after " << it << " iterations");
505
506
507 if (it<1)
508 p = r;
509 else
510 {
511 beta = Simd::cond(norm==field_type(0.),
512 field_type(0.), // no need for orthogonalization if norm is already 0
513 field_type(( rho_new / rho ) * ( alpha / omega )));
514 p.axpy(-omega,v); // p = r + beta (p - omega*v)
515 p *= beta;
516 p += r;
517 }
518
519 // y = W^-1 * p
520 y = 0;
521 _prec->apply(y,p); // apply preconditioner
522
523 // v = A * y
524 _op->apply(y,v);
525
526 // alpha = rho_new / < rt, v >
527 h = _sp->dot(rt,v);
528
529 if ( Simd::allTrue(abs(h) < EPSILON) )
530 DUNE_THROW(SolverAbort,"abs(h) < EPSILON in BiCGSTAB - abs(h) "
531 << Simd::io(abs(h)) << " < EPSILON " << EPSILON
532 << " after " << it << " iterations");
533
534 alpha = Simd::cond(norm==field_type(0.),
535 field_type(0.),
536 field_type(rho_new / h));
537
538 // apply first correction to x
539 // x <- x + alpha y
540 x.axpy(alpha,y);
541
542 // r = r - alpha*v
543 r.axpy(-alpha,v);
544
545 //
546 // test stop criteria
547 //
548
549 norm = _sp->norm(r);
550 if(iteration.step(it, norm)){
551 break;
552 }
553
554 it+=.5;
555
556 // y = W^-1 * r
557 y = 0;
558 _prec->apply(y,r);
559
560 // t = A * y
561 _op->apply(y,t);
562
563 // omega = < t, r > / < t, t >
564 h = _sp->dot(t,t);
565 omega = Simd::cond(norm==field_type(0.),
566 field_type(0.),
567 field_type(_sp->dot(t,r)/h));
568
569 // apply second correction to x
570 // x <- x + omega y
571 x.axpy(omega,y);
572
573 // r = s - omega*t (remember : r = s)
574 r.axpy(-omega,t);
575
576 rho = rho_new;
577
578 //
579 // test stop criteria
580 //
581
582 norm = _sp->norm(r);
583 if(iteration.step(it, norm)){
584 break;
585 }
586 } // end for
587
588 _prec->post(x); // postprocess preconditioner
589 }
590
591 protected:
592 using IterativeSolver<X,X>::_op;
593 using IterativeSolver<X,X>::_prec;
594 using IterativeSolver<X,X>::_sp;
595 using IterativeSolver<X,X>::_reduction;
596 using IterativeSolver<X,X>::_maxit;
597 using IterativeSolver<X,X>::_verbose;
598 template<class CountType>
599 using Iteration = typename IterativeSolver<X,X>::template Iteration<CountType>;
600 };
601 DUNE_REGISTER_SOLVER("bicgstabsolver", defaultIterativeSolverCreator<Dune::BiCGSTABSolver>());
602
609 template<class X>
610 class MINRESSolver : public IterativeSolver<X,X> {
611 public:
615 using typename IterativeSolver<X,X>::real_type;
616
617 // copy base class constructors
619
620 // don't shadow four-argument version of apply defined in the base class
621 using IterativeSolver<X,X>::apply;
622
628 void apply (X& x, X& b, InverseOperatorResult& res) override
629 {
630 using std::sqrt;
631 using std::abs;
632 Iteration iteration(*this, res);
633 // prepare preconditioner
634 _prec->pre(x,b);
635
636 // overwrite rhs with defect
637 _op->applyscaleadd(-1.0,x,b); // b -= Ax
638
639 // some temporary vectors
640 X z(b), dummy(b);
641 z = 0.0;
642
643 // calculate preconditioned defect
644 _prec->apply(z,b); // r = W^-1 (b - Ax)
645 real_type def = _sp->norm(z);
646 if (iteration.step(0, def)){
647 _prec->post(x);
648 return;
649 }
650
651 // recurrence coefficients as computed in Lanczos algorithm
652 field_type alpha, beta;
653 // diagonal entries of givens rotation
654 std::array<real_type,2> c{{0.0,0.0}};
655 // off-diagonal entries of givens rotation
656 std::array<field_type,2> s{{0.0,0.0}};
657
658 // recurrence coefficients (column k of tridiag matrix T_k)
659 std::array<field_type,3> T{{0.0,0.0,0.0}};
660
661 // the rhs vector of the min problem
662 std::array<field_type,2> xi{{1.0,0.0}};
663
664 // beta is real and positive in exact arithmetic
665 // since it is the norm of the basis vectors (in unpreconditioned case)
666 beta = sqrt(_sp->dot(b,z));
667 field_type beta0 = beta;
668
669 // the search directions
670 std::array<X,3> p{{b,b,b}};
671 p[0] = 0.0;
672 p[1] = 0.0;
673 p[2] = 0.0;
674
675 // orthonormal basis vectors (in unpreconditioned case)
676 std::array<X,3> q{{b,b,b}};
677 q[0] = 0.0;
678 q[1] *= Simd::cond(def==field_type(0.),
679 field_type(0.),
680 field_type(real_type(1.0)/beta));
681 q[2] = 0.0;
682
683 z *= Simd::cond(def==field_type(0.),
684 field_type(0.),
685 field_type(real_type(1.0)/beta));
686
687 // the loop
688 int i = 1;
689 for( ; i<=_maxit; i++) {
690
691 dummy = z;
692 int i1 = i%3,
693 i0 = (i1+2)%3,
694 i2 = (i1+1)%3;
695
696 // symmetrically preconditioned Lanczos algorithm (see Greenbaum p.121)
697 _op->apply(z,q[i2]); // q[i2] = Az
698 q[i2].axpy(-beta,q[i0]);
699 // alpha is real since it is the diagonal entry of the hermitian tridiagonal matrix
700 // from the Lanczos Algorithm
701 // so the order in the scalar product doesn't matter even for the complex case
702 alpha = _sp->dot(z,q[i2]);
703 q[i2].axpy(-alpha,q[i1]);
704
705 z = 0.0;
706 _prec->apply(z,q[i2]);
707
708 // beta is real and positive in exact arithmetic
709 // since it is the norm of the basis vectors (in unpreconditioned case)
710 beta = sqrt(_sp->dot(q[i2],z));
711
712 q[i2] *= Simd::cond(def==field_type(0.),
713 field_type(0.),
714 field_type(real_type(1.0)/beta));
715 z *= Simd::cond(def==field_type(0.),
716 field_type(0.),
717 field_type(real_type(1.0)/beta));
718
719 // QR Factorization of recurrence coefficient matrix
720 // apply previous givens rotations to last column of T
721 T[1] = T[2];
722 if(i>2) {
723 T[0] = s[i%2]*T[1];
724 T[1] = c[i%2]*T[1];
725 }
726 if(i>1) {
727 T[2] = c[(i+1)%2]*alpha - s[(i+1)%2]*T[1];
728 T[1] = c[(i+1)%2]*T[1] + s[(i+1)%2]*alpha;
729 }
730 else
731 T[2] = alpha;
732
733 // update QR factorization
734 generateGivensRotation(T[2],beta,c[i%2],s[i%2]);
735 // to last column of T_k
736 T[2] = c[i%2]*T[2] + s[i%2]*beta;
737 // and to the rhs xi of the min problem
738 xi[i%2] = -s[i%2]*xi[(i+1)%2];
739 xi[(i+1)%2] *= c[i%2];
740
741 // compute correction direction
742 p[i2] = dummy;
743 p[i2].axpy(-T[1],p[i1]);
744 p[i2].axpy(-T[0],p[i0]);
745 p[i2] *= real_type(1.0)/T[2];
746
747 // apply correction/update solution
748 x.axpy(beta0*xi[(i+1)%2],p[i2]);
749
750 // remember beta_old
751 T[2] = beta;
752
753 // check for convergence
754 // the last entry in the rhs of the min-problem is the residual
755 def = abs(beta0*xi[i%2]);
756 if(iteration.step(i, def)){
757 break;
758 }
759 } // end for
760
761 // postprocess preconditioner
762 _prec->post(x);
763 }
764
765 private:
766
767 void generateGivensRotation(field_type &dx, field_type &dy, real_type &cs, field_type &sn)
768 {
769 using std::sqrt;
770 using std::abs;
771 using std::max;
772 using std::min;
773 const real_type eps = 1e-15;
774 real_type norm_dx = abs(dx);
775 real_type norm_dy = abs(dy);
776 real_type norm_max = max(norm_dx, norm_dy);
777 real_type norm_min = min(norm_dx, norm_dy);
778 real_type temp = norm_min/norm_max;
779 // we rewrite the code in a vectorizable fashion
780 cs = Simd::cond(norm_dy < eps,
781 real_type(1.0),
782 real_type(Simd::cond(norm_dx < eps,
783 real_type(0.0),
784 real_type(Simd::cond(norm_dy > norm_dx,
785 real_type(1.0)/sqrt(real_type(1.0) + temp*temp)*temp,
786 real_type(1.0)/sqrt(real_type(1.0) + temp*temp)
787 )))));
788 sn = Simd::cond(norm_dy < eps,
789 field_type(0.0),
790 field_type(Simd::cond(norm_dx < eps,
791 field_type(1.0),
792 field_type(Simd::cond(norm_dy > norm_dx,
793 // dy and dx are real in exact arithmetic
794 // thus dx*dy is real so we can explicitly enforce it
795 field_type(1.0)/sqrt(real_type(1.0) + temp*temp)*dx*dy/norm_dx/norm_dy,
796 // dy and dx is real in exact arithmetic
797 // so we don't have to conjugate both of them
798 field_type(1.0)/sqrt(real_type(1.0) + temp*temp)*dy/dx
799 )))));
800 }
801
802 protected:
803 using IterativeSolver<X,X>::_op;
804 using IterativeSolver<X,X>::_prec;
805 using IterativeSolver<X,X>::_sp;
806 using IterativeSolver<X,X>::_reduction;
807 using IterativeSolver<X,X>::_maxit;
808 using IterativeSolver<X,X>::_verbose;
809 using Iteration = typename IterativeSolver<X,X>::template Iteration<unsigned int>;
810 };
811 DUNE_REGISTER_SOLVER("minressolver", defaultIterativeSolverCreator<Dune::MINRESSolver>());
812
826 template<class X, class Y=X, class F = Y>
828 {
829 public:
833 using typename IterativeSolver<X,Y>::real_type;
834
835 protected:
837
839 using fAlloc = ReboundAllocatorType<X,field_type>;
841 using rAlloc = ReboundAllocatorType<X,real_type>;
842
843 public:
844
851 RestartedGMResSolver (const LinearOperator<X,Y>& op, Preconditioner<X,Y>& prec, scalar_real_type reduction, int restart, int maxit, int verbose) :
852 IterativeSolver<X,Y>::IterativeSolver(op,prec,reduction,maxit,verbose),
853 _restart(restart)
854 {}
855
862 RestartedGMResSolver (const LinearOperator<X,Y>& op, const ScalarProduct<X>& sp, Preconditioner<X,Y>& prec, scalar_real_type reduction, int restart, int maxit, int verbose) :
863 IterativeSolver<X,Y>::IterativeSolver(op,sp,prec,reduction,maxit,verbose),
864 _restart(restart)
865 {}
866
879 RestartedGMResSolver (std::shared_ptr<const LinearOperator<X,Y> > op, std::shared_ptr<Preconditioner<X,X> > prec, const ParameterTree& configuration) :
880 IterativeSolver<X,Y>::IterativeSolver(op,prec,configuration),
881 _restart(configuration.get<int>("restart"))
882 {}
883
884 RestartedGMResSolver (std::shared_ptr<const LinearOperator<X,Y> > op, std::shared_ptr<const ScalarProduct<X> > sp, std::shared_ptr<Preconditioner<X,X> > prec, const ParameterTree& configuration) :
885 IterativeSolver<X,Y>::IterativeSolver(op,sp,prec,configuration),
886 _restart(configuration.get<int>("restart"))
887 {}
888
895 RestartedGMResSolver (std::shared_ptr<const LinearOperator<X,Y>> op,
896 std::shared_ptr<const ScalarProduct<X>> sp,
897 std::shared_ptr<Preconditioner<X,Y>> prec,
898 scalar_real_type reduction, int restart, int maxit, int verbose) :
899 IterativeSolver<X,Y>::IterativeSolver(op,sp,prec,reduction,maxit,verbose),
900 _restart(restart)
901 {}
902
903 // don't shadow four-argument version of apply defined in the base class
904 using IterativeSolver<X,Y>::apply;
905
914 void apply (X& x, Y& b, InverseOperatorResult& res) override
915 {
916 using std::abs;
917 const Simd::Scalar<real_type> EPSILON = 1e-80;
918 const int m = _restart;
919 real_type norm = 0.0;
920 int j = 1;
921 std::vector<field_type,fAlloc> s(m+1), sn(m);
922 std::vector<real_type,rAlloc> cs(m);
923 // need copy of rhs if GMRes has to be restarted
924 Y b2(b);
925 // helper vector
926 Y w(b);
927 std::vector< std::vector<field_type,fAlloc> > H(m+1,s);
928 std::vector<F> v(m+1,b);
929
930 Iteration iteration(*this,res);
931
932 // clear solver statistics and set res.converged to false
933 _prec->pre(x,b);
934
935 // calculate defect and overwrite rhs with it
936 _op->applyscaleadd(-1.0,x,b); // b -= Ax
937 // calculate preconditioned defect
938 v[0] = 0.0; _prec->apply(v[0],b); // r = W^-1 b
939 norm = _sp->norm(v[0]);
940 if(iteration.step(0, norm)){
941 _prec->post(x);
942 return;
943 }
944
945 while(j <= _maxit && res.converged != true) {
946
947 int i = 0;
948 v[0] *= Simd::cond(norm==real_type(0.),
949 real_type(0.),
950 real_type(real_type(1.0)/norm));
951 s[0] = norm;
952 for(i=1; i<m+1; i++)
953 s[i] = 0.0;
954
955 for(i=0; i < m && j <= _maxit && res.converged != true; i++, j++) {
956 w = 0.0;
957 // use v[i+1] as temporary vector
958 v[i+1] = 0.0;
959 // do Arnoldi algorithm
960 _op->apply(v[i],v[i+1]);
961 _prec->apply(w,v[i+1]);
962 for(int k=0; k<i+1; k++) {
963 // notice that _sp->dot(v[k],w) = v[k]\adjoint w
964 // so one has to pay attention to the order
965 // in the scalar product for the complex case
966 // doing the modified Gram-Schmidt algorithm
967 H[k][i] = _sp->dot(v[k],w);
968 // w -= H[k][i] * v[k]
969 w.axpy(-H[k][i],v[k]);
970 }
971 H[i+1][i] = _sp->norm(w);
972 if(Simd::allTrue(abs(H[i+1][i]) < EPSILON))
974 "breakdown in GMRes - |w| == 0.0 after " << j << " iterations");
975
976 // normalize new vector
977 v[i+1] = w;
978 v[i+1] *= Simd::cond(norm==real_type(0.),
979 field_type(0.),
980 field_type(real_type(1.0)/H[i+1][i]));
981
982 // update QR factorization
983 for(int k=0; k<i; k++)
984 applyPlaneRotation(H[k][i],H[k+1][i],cs[k],sn[k]);
985
986 // compute new givens rotation
987 generatePlaneRotation(H[i][i],H[i+1][i],cs[i],sn[i]);
988 // finish updating QR factorization
989 applyPlaneRotation(H[i][i],H[i+1][i],cs[i],sn[i]);
990 applyPlaneRotation(s[i],s[i+1],cs[i],sn[i]);
991
992 // norm of the defect is the last component the vector s
993 norm = abs(s[i+1]);
994
995 iteration.step(j, norm);
996
997 } // end for
998
999 // calculate update vector
1000 w = 0.0;
1001 update(w,i,H,s,v);
1002 // and current iterate
1003 x += w;
1004
1005 // restart GMRes if convergence was not achieved,
1006 // i.e. linear defect has not reached desired reduction
1007 // and if j < _maxit (do not restart on last iteration)
1008 if( res.converged != true && j < _maxit ) {
1009
1010 if(_verbose > 0)
1011 std::cout << "=== GMRes::restart" << std::endl;
1012 // get saved rhs
1013 b = b2;
1014 // calculate new defect
1015 _op->applyscaleadd(-1.0,x,b); // b -= Ax;
1016 // calculate preconditioned defect
1017 v[0] = 0.0;
1018 _prec->apply(v[0],b);
1019 norm = _sp->norm(v[0]);
1020 }
1021
1022 } //end while
1023
1024 // postprocess preconditioner
1025 _prec->post(x);
1026 }
1027
1028 protected :
1029
1030 void update(X& w, int i,
1031 const std::vector<std::vector<field_type,fAlloc> >& H,
1032 const std::vector<field_type,fAlloc>& s,
1033 const std::vector<X>& v) {
1034 // solution vector of the upper triangular system
1035 std::vector<field_type,fAlloc> y(s);
1036
1037 // backsolve
1038 for(int a=i-1; a>=0; a--) {
1039 field_type rhs(s[a]);
1040 for(int b=a+1; b<i; b++)
1041 rhs -= H[a][b]*y[b];
1042 y[a] = Simd::cond(rhs==field_type(0.),
1043 field_type(0.),
1044 field_type(rhs/H[a][a]));
1045
1046 // compute update on the fly
1047 // w += y[a]*v[a]
1048 w.axpy(y[a],v[a]);
1049 }
1050 }
1051
1052 template<typename T>
1053 typename std::enable_if<std::is_same<field_type,real_type>::value,T>::type conjugate(const T& t) {
1054 return t;
1055 }
1056
1057 template<typename T>
1058 typename std::enable_if<!std::is_same<field_type,real_type>::value,T>::type conjugate(const T& t) {
1059 using std::conj;
1060 return conj(t);
1061 }
1062
1063 void
1064 generatePlaneRotation(field_type &dx, field_type &dy, real_type &cs, field_type &sn)
1065 {
1066 using std::sqrt;
1067 using std::abs;
1068 using std::max;
1069 using std::min;
1070 const real_type eps = 1e-15;
1071 real_type norm_dx = abs(dx);
1072 real_type norm_dy = abs(dy);
1073 real_type norm_max = max(norm_dx, norm_dy);
1074 real_type norm_min = min(norm_dx, norm_dy);
1075 real_type temp = norm_min/norm_max;
1076 // we rewrite the code in a vectorizable fashion
1077 cs = Simd::cond(norm_dy < eps,
1078 real_type(1.0),
1079 Simd::cond(norm_dx < eps,
1080 real_type(0.0),
1081 real_type(Simd::cond(norm_dy > norm_dx,
1082 real_type(1.0)/sqrt(real_type(1.0) + temp*temp)*temp,
1083 real_type(1.0)/sqrt(real_type(1.0) + temp*temp)
1084 ))));
1085 sn = Simd::cond(norm_dy < eps,
1086 field_type(0.0),
1087 Simd::cond(norm_dx < eps,
1088 field_type(1.0),
1089 field_type(Simd::cond(norm_dy > norm_dx,
1090 field_type(1.0)/sqrt(real_type(1.0) + temp*temp)*dx*conjugate(dy)/norm_dx/norm_dy,
1091 field_type(1.0)/sqrt(real_type(1.0) + temp*temp)*conjugate(dy/dx)
1092 ))));
1093 }
1094
1095
1096 void
1097 applyPlaneRotation(field_type &dx, field_type &dy, real_type &cs, field_type &sn)
1098 {
1099 field_type temp = cs * dx + sn * dy;
1100 dy = -conjugate(sn) * dx + cs * dy;
1101 dx = temp;
1102 }
1103
1104 using IterativeSolver<X,Y>::_op;
1105 using IterativeSolver<X,Y>::_prec;
1106 using IterativeSolver<X,Y>::_sp;
1107 using IterativeSolver<X,Y>::_reduction;
1108 using IterativeSolver<X,Y>::_maxit;
1109 using IterativeSolver<X,Y>::_verbose;
1110 using Iteration = typename IterativeSolver<X,X>::template Iteration<unsigned int>;
1111 int _restart;
1112 };
1113 DUNE_REGISTER_SOLVER("restartedgmressolver", defaultIterativeSolverCreator<Dune::RestartedGMResSolver>());
1114
1128 template<class X, class Y=X, class F = Y>
1130 {
1131 public:
1136
1137 private:
1139
1141 using fAlloc = typename RestartedGMResSolver<X,Y>::fAlloc;
1143 using rAlloc = typename RestartedGMResSolver<X,Y>::rAlloc;
1144
1145 public:
1146 // copy base class constructors
1148
1149 // don't shadow four-argument version of apply defined in the base class
1150 using RestartedGMResSolver<X,Y>::apply;
1151
1160 void apply (X& x, Y& b, InverseOperatorResult& res) override
1161 {
1162 using std::abs;
1163 const Simd::Scalar<real_type> EPSILON = 1e-80;
1164 const int m = _restart;
1165 real_type norm = 0.0;
1166 int i, j = 1, k;
1167 std::vector<field_type,fAlloc> s(m+1), sn(m);
1168 std::vector<real_type,rAlloc> cs(m);
1169 // helper vector
1170 Y tmp(b);
1171 std::vector< std::vector<field_type,fAlloc> > H(m+1,s);
1172 std::vector<F> v(m+1,b);
1173 std::vector<X> w(m+1,b);
1174
1175 Iteration iteration(*this,res);
1176 // setup preconditioner if it does something in pre
1177
1178 // calculate residual and overwrite a copy of the rhs with it
1179 _prec->pre(x, b);
1180 v[0] = b;
1181 _op->applyscaleadd(-1.0, x, v[0]); // b -= Ax
1182
1183 norm = _sp->norm(v[0]); // the residual norm
1184 if(iteration.step(0, norm)){
1185 _prec->post(x);
1186 return;
1187 }
1188
1189 // start iterations
1190 res.converged = false;;
1191 while(j <= _maxit && res.converged != true)
1192 {
1193 v[0] *= (1.0 / norm);
1194 s[0] = norm;
1195 for(i=1; i<m+1; ++i)
1196 s[i] = 0.0;
1197
1198 // inner loop
1199 for(i=0; i < m && j <= _maxit && res.converged != true; i++, j++)
1200 {
1201 w[i] = 0.0;
1202 // compute wi = M^-1*vi (also called zi)
1203 _prec->apply(w[i], v[i]);
1204 // compute vi = A*wi
1205 // use v[i+1] as temporary vector for w
1206 _op->apply(w[i], v[i+1]);
1207 // do Arnoldi algorithm
1208 for(int kk=0; kk<i+1; kk++)
1209 {
1210 // notice that _sp->dot(v[k],v[i+1]) = v[k]\adjoint v[i+1]
1211 // so one has to pay attention to the order
1212 // in the scalar product for the complex case
1213 // doing the modified Gram-Schmidt algorithm
1214 H[kk][i] = _sp->dot(v[kk],v[i+1]);
1215 // w -= H[k][i] * v[kk]
1216 v[i+1].axpy(-H[kk][i], v[kk]);
1217 }
1218 H[i+1][i] = _sp->norm(v[i+1]);
1219 if(Simd::allTrue(abs(H[i+1][i]) < EPSILON))
1220 DUNE_THROW(SolverAbort, "breakdown in fGMRes - |w| (-> "
1221 << w[i] << ") == 0.0 after "
1222 << j << " iterations");
1223
1224 // v[i+1] = w*1/H[i+1][i]
1225 v[i+1] *= real_type(1.0)/H[i+1][i];
1226
1227 // update QR factorization
1228 for(k=0; k<i; k++)
1229 this->applyPlaneRotation(H[k][i],H[k+1][i],cs[k],sn[k]);
1230
1231 // compute new givens rotation
1232 this->generatePlaneRotation(H[i][i],H[i+1][i],cs[i],sn[i]);
1233
1234 // finish updating QR factorization
1235 this->applyPlaneRotation(H[i][i],H[i+1][i],cs[i],sn[i]);
1236 this->applyPlaneRotation(s[i],s[i+1],cs[i],sn[i]);
1237
1238 // norm of the residual is the last component of vector s
1239 using std::abs;
1240 norm = abs(s[i+1]);
1241 iteration.step(j, norm);
1242 } // end inner for loop
1243
1244 // calculate update vector
1245 tmp = 0.0;
1246 this->update(tmp, i, H, s, w);
1247 // and update current iterate
1248 x += tmp;
1249
1250 // restart fGMRes if convergence was not achieved,
1251 // i.e. linear residual has not reached desired reduction
1252 // and if still j < _maxit (do not restart on last iteration)
1253 if( res.converged != true && j < _maxit)
1254 {
1255 if (_verbose > 0)
1256 std::cout << "=== fGMRes::restart" << std::endl;
1257 // get rhs
1258 v[0] = b;
1259 // calculate new defect
1260 _op->applyscaleadd(-1.0, x,v[0]); // b -= Ax;
1261 // calculate preconditioned defect
1262 norm = _sp->norm(v[0]); // update the residual norm
1263 }
1264
1265 } // end outer while loop
1266
1267 // post-process preconditioner
1268 _prec->post(x);
1269 }
1270
1271private:
1272 using RestartedGMResSolver<X,Y>::_op;
1273 using RestartedGMResSolver<X,Y>::_prec;
1274 using RestartedGMResSolver<X,Y>::_sp;
1275 using RestartedGMResSolver<X,Y>::_reduction;
1276 using RestartedGMResSolver<X,Y>::_maxit;
1277 using RestartedGMResSolver<X,Y>::_verbose;
1278 using RestartedGMResSolver<X,Y>::_restart;
1279 using Iteration = typename IterativeSolver<X,X>::template Iteration<unsigned int>;
1280 };
1281 DUNE_REGISTER_SOLVER("restartedflexiblegmressolver", defaultIterativeSolverCreator<Dune::RestartedFlexibleGMResSolver>());
1282
1296 template<class X>
1298 {
1299 public:
1300 using typename IterativeSolver<X,X>::domain_type;
1301 using typename IterativeSolver<X,X>::range_type;
1302 using typename IterativeSolver<X,X>::field_type;
1303 using typename IterativeSolver<X,X>::real_type;
1304
1305 private:
1307
1309 using fAlloc = ReboundAllocatorType<X,field_type>;
1310
1311 public:
1312
1313 // don't shadow four-argument version of apply defined in the base class
1314 using IterativeSolver<X,X>::apply;
1315
1322 GeneralizedPCGSolver (const LinearOperator<X,X>& op, Preconditioner<X,X>& prec, scalar_real_type reduction, int maxit, int verbose, int restart = 10) :
1323 IterativeSolver<X,X>::IterativeSolver(op,prec,reduction,maxit,verbose),
1324 _restart(restart)
1325 {}
1326
1334 GeneralizedPCGSolver (const LinearOperator<X,X>& op, const ScalarProduct<X>& sp, Preconditioner<X,X>& prec, scalar_real_type reduction, int maxit, int verbose, int restart = 10) :
1335 IterativeSolver<X,X>::IterativeSolver(op,sp,prec,reduction,maxit,verbose),
1336 _restart(restart)
1337 {}
1338
1339
1352 GeneralizedPCGSolver (std::shared_ptr<const LinearOperator<X,X> > op, std::shared_ptr<Preconditioner<X,X> > prec, const ParameterTree& configuration) :
1353 IterativeSolver<X,X>::IterativeSolver(op,prec,configuration),
1354 _restart(configuration.get<int>("restart"))
1355 {}
1356
1357 GeneralizedPCGSolver (std::shared_ptr<const LinearOperator<X,X> > op, std::shared_ptr<const ScalarProduct<X> > sp, std::shared_ptr<Preconditioner<X,X> > prec, const ParameterTree& configuration) :
1358 IterativeSolver<X,X>::IterativeSolver(op,sp,prec,configuration),
1359 _restart(configuration.get<int>("restart"))
1360 {}
1367 GeneralizedPCGSolver (std::shared_ptr<const LinearOperator<X,X>> op,
1368 std::shared_ptr<const ScalarProduct<X>> sp,
1369 std::shared_ptr<Preconditioner<X,X>> prec,
1370 scalar_real_type reduction, int maxit, int verbose,
1371 int restart = 10) :
1372 IterativeSolver<X,X>::IterativeSolver(op,sp,prec,reduction,maxit,verbose),
1373 _restart(restart)
1374 {}
1375
1381 void apply (X& x, X& b, InverseOperatorResult& res) override
1382 {
1383 Iteration iteration(*this, res);
1384 _prec->pre(x,b); // prepare preconditioner
1385 _op->applyscaleadd(-1,x,b); // overwrite b with defect
1386
1387 std::vector<std::shared_ptr<X> > p(_restart);
1388 std::vector<field_type,fAlloc> pp(_restart);
1389 X q(x); // a temporary vector
1390 X prec_res(x); // a temporary vector for preconditioner output
1391
1392 p[0].reset(new X(x));
1393
1394 real_type def = _sp->norm(b); // compute norm
1395 if(iteration.step(0, def)){
1396 _prec->post(x);
1397 return;
1398 }
1399 // some local variables
1400 field_type rho, lambda;
1401
1402 int i=0;
1403 // determine initial search direction
1404 *(p[0]) = 0; // clear correction
1405 _prec->apply(*(p[0]),b); // apply preconditioner
1406 rho = _sp->dot(*(p[0]),b); // orthogonalization
1407 _op->apply(*(p[0]),q); // q=Ap
1408 pp[0] = _sp->dot(*(p[0]),q); // scalar product
1409 lambda = rho/pp[0]; // minimization
1410 x.axpy(lambda,*(p[0])); // update solution
1411 b.axpy(-lambda,q); // update defect
1412
1413 // convergence test
1414 def=_sp->norm(b); // comp defect norm
1415 ++i;
1416 if(iteration.step(i, def)){
1417 _prec->post(x);
1418 return;
1419 }
1420
1421 while(i<_maxit) {
1422 // the loop
1423 int end=std::min(_restart, _maxit-i+1);
1424 for (int ii = 1; ii < end; ++ii)
1425 {
1426 //std::cout<<" ii="<<ii<<" i="<<i<<std::endl;
1427 // compute next conjugate direction
1428 prec_res = 0; // clear correction
1429 _prec->apply(prec_res,b); // apply preconditioner
1430
1431 p[ii].reset(new X(prec_res));
1432 _op->apply(prec_res, q);
1433
1434 for(int j=0; j<ii; ++j) {
1435 rho =_sp->dot(q,*(p[j]))/pp[j];
1436 p[ii]->axpy(-rho, *(p[j]));
1437 }
1438
1439 // minimize in given search direction
1440 _op->apply(*(p[ii]),q); // q=Ap
1441 pp[ii] = _sp->dot(*(p[ii]),q); // scalar product
1442 rho = _sp->dot(*(p[ii]),b); // orthogonalization
1443 lambda = rho/pp[ii]; // minimization
1444 x.axpy(lambda,*(p[ii])); // update solution
1445 b.axpy(-lambda,q); // update defect
1446
1447 // convergence test
1448 def = _sp->norm(b); // comp defect norm
1449
1450 ++i;
1451 iteration.step(i, def);
1452 }
1453 if(res.converged)
1454 break;
1455 if(end==_restart) {
1456 *(p[0])=*(p[_restart-1]);
1457 pp[0]=pp[_restart-1];
1458 }
1459 }
1460
1461 // postprocess preconditioner
1462 _prec->post(x);
1463
1464 }
1465
1466 private:
1467 using IterativeSolver<X,X>::_op;
1468 using IterativeSolver<X,X>::_prec;
1469 using IterativeSolver<X,X>::_sp;
1470 using IterativeSolver<X,X>::_reduction;
1471 using IterativeSolver<X,X>::_maxit;
1472 using IterativeSolver<X,X>::_verbose;
1473 using Iteration = typename IterativeSolver<X,X>::template Iteration<unsigned int>;
1474 int _restart;
1475 };
1476 DUNE_REGISTER_SOLVER("generalizedpcgsolver", defaultIterativeSolverCreator<Dune::GeneralizedPCGSolver>());
1477
1489 template<class X>
1491 public:
1492 using typename IterativeSolver<X,X>::domain_type;
1493 using typename IterativeSolver<X,X>::range_type;
1494 using typename IterativeSolver<X,X>::field_type;
1495 using typename IterativeSolver<X,X>::real_type;
1496
1497 private:
1499
1500 public:
1501 // don't shadow four-argument version of apply defined in the base class
1502 using IterativeSolver<X,X>::apply;
1509 scalar_real_type reduction, int maxit, int verbose, int mmax = 10) : IterativeSolver<X,X>(op, prec, reduction, maxit, verbose), _mmax(mmax)
1510 {
1511 }
1512
1519 scalar_real_type reduction, int maxit, int verbose, int mmax = 10) : IterativeSolver<X,X>(op, sp, prec, reduction, maxit, verbose), _mmax(mmax)
1520 {
1521 }
1522
1528 RestartedFCGSolver (std::shared_ptr<const LinearOperator<X,X>> op,
1529 std::shared_ptr<const ScalarProduct<X>> sp,
1530 std::shared_ptr<Preconditioner<X,X>> prec,
1531 scalar_real_type reduction, int maxit, int verbose,
1532 int mmax = 10)
1533 : IterativeSolver<X,X>(op, sp, prec, reduction, maxit, verbose), _mmax(mmax)
1534 {}
1535
1548 RestartedFCGSolver (std::shared_ptr<const LinearOperator<X,X>> op,
1549 std::shared_ptr<Preconditioner<X,X>> prec,
1550 const ParameterTree& configuration)
1551 : IterativeSolver<X,X>(op, prec, configuration), _mmax(configuration.get("mmax", 10))
1552 {}
1553
1554 RestartedFCGSolver (std::shared_ptr<const LinearOperator<X,X>> op,
1555 std::shared_ptr<const ScalarProduct<X>> sp,
1556 std::shared_ptr<Preconditioner<X,X>> prec,
1557 const ParameterTree& config)
1558 : IterativeSolver<X,X>(op, sp, prec, config), _mmax(config.get("mmax", 10))
1559 {}
1560
1573 void apply (X& x, X& b, InverseOperatorResult& res) override
1574 {
1575 using rAlloc = ReboundAllocatorType<X,field_type>;
1576 res.clear();
1577 Iteration iteration(*this,res);
1578 _prec->pre(x,b); // prepare preconditioner
1579 _op->applyscaleadd(-1,x,b); // overwrite b with defect
1580
1581 //arrays for interim values:
1582 std::vector<X> d(_mmax+1, x); // array for directions
1583 std::vector<X> Ad(_mmax+1, x); // array for Ad[i]
1584 std::vector<field_type,rAlloc> ddotAd(_mmax+1,0); // array for <d[i],Ad[i]>
1585 X w(x);
1586
1587 real_type def = _sp->norm(b); // compute norm
1588 if(iteration.step(0, def)){
1589 _prec->post(x);
1590 return;
1591 }
1592
1593 // some local variables
1594 field_type alpha;
1595
1596 // the loop
1597 int i=1;
1598 int i_bounded=0;
1599 while(i<=_maxit && !res.converged) {
1600 for (; i_bounded <= _mmax && i<= _maxit; i_bounded++) {
1601 d[i_bounded] = 0; // reset search direction
1602 _prec->apply(d[i_bounded], b); // apply preconditioner
1603 w = d[i_bounded]; // copy of current d[i]
1604 // orthogonalization with previous directions
1605 orthogonalizations(i_bounded,Ad,w,ddotAd,d);
1606
1607 //saving interim values for future calculating
1608 _op->apply(d[i_bounded], Ad[i_bounded]); // save Ad[i]
1609 ddotAd[i_bounded]=_sp->dot(d[i_bounded],Ad[i_bounded]); // save <d[i],Ad[i]>
1610 alpha = _sp->dot(d[i_bounded], b)/ddotAd[i_bounded]; // <d[i],b>/<d[i],Ad[i]>
1611
1612 //update solution and defect
1613 x.axpy(alpha, d[i_bounded]);
1614 b.axpy(-alpha, Ad[i_bounded]);
1615
1616 // convergence test
1617 def = _sp->norm(b); // comp defect norm
1618
1619 iteration.step(i, def);
1620 i++;
1621 }
1622 //restart: exchange first and last stored values
1623 cycle(Ad,d,ddotAd,i_bounded);
1624 }
1625
1626 //correct i which is wrong if convergence was not achieved.
1627 i=std::min(_maxit,i);
1628
1629 _prec->post(x); // postprocess preconditioner
1630 }
1631
1632 private:
1633 //This function is called every iteration to orthogonalize against the last search directions
1634 virtual void orthogonalizations(const int& i_bounded,const std::vector<X>& Ad, const X& w, const std::vector<field_type,ReboundAllocatorType<X,field_type>>& ddotAd,std::vector<X>& d) {
1635 // The RestartedFCGSolver uses only values with lower array index;
1636 for (int k = 0; k < i_bounded; k++) {
1637 d[i_bounded].axpy(-_sp->dot(Ad[k], w) / ddotAd[k], d[k]); // d[i] -= <<Ad[k],w>/<d[k],Ad[k]>>d[k]
1638 }
1639 }
1640
1641 // This function is called every mmax iterations to handle limited array sizes.
1642 virtual void cycle(std::vector<X>& Ad,std::vector<X>& d,std::vector<field_type,ReboundAllocatorType<X,field_type> >& ddotAd,int& i_bounded) {
1643 // Reset loop index and exchange the first and last arrays
1644 i_bounded = 1;
1645 std::swap(Ad[0], Ad[_mmax]);
1646 std::swap(d[0], d[_mmax]);
1647 std::swap(ddotAd[0], ddotAd[_mmax]);
1648 }
1649
1650 protected:
1651 int _mmax;
1652 using IterativeSolver<X,X>::_op;
1653 using IterativeSolver<X,X>::_prec;
1654 using IterativeSolver<X,X>::_sp;
1655 using IterativeSolver<X,X>::_reduction;
1656 using IterativeSolver<X,X>::_maxit;
1657 using IterativeSolver<X,X>::_verbose;
1658 using Iteration = typename IterativeSolver<X,X>::template Iteration<unsigned int>;
1659 };
1660 DUNE_REGISTER_SOLVER("restartedfcgsolver", defaultIterativeSolverCreator<Dune::RestartedFCGSolver>());
1661
1668 template<class X>
1670 public:
1672 using typename RestartedFCGSolver<X>::range_type;
1673 using typename RestartedFCGSolver<X>::field_type;
1674 using typename RestartedFCGSolver<X>::real_type;
1675
1676 // copy base class constructors
1678
1679 // don't shadow four-argument version of apply defined in the base class
1681
1682 // just a minor part of the RestartedFCGSolver apply method will be modified
1683 void apply (X& x, X& b, InverseOperatorResult& res) override
1684 {
1685 // reset limiter of orthogonalization loop
1686 _k_limit = 0;
1687 this->RestartedFCGSolver<X>::apply(x,b,res);
1688 };
1689
1690 private:
1691 // This function is called every iteration to orthogonalize against the last search directions.
1692 void orthogonalizations(const int& i_bounded,const std::vector<X>& Ad, const X& w, const std::vector<field_type,ReboundAllocatorType<X,field_type>>& ddotAd,std::vector<X>& d) override {
1693 // This FCGSolver uses values with higher array indexes too, if existent.
1694 for (int k = 0; k < _k_limit; k++) {
1695 if(i_bounded!=k)
1696 d[i_bounded].axpy(-_sp->dot(Ad[k], w) / ddotAd[k], d[k]); // d[i] -= <<Ad[k],w>/<d[k],Ad[k]>>d[k]
1697 }
1698 // The loop limit increase, if array is not completely filled.
1699 if(_k_limit<=i_bounded)
1700 _k_limit++;
1701
1702 };
1703
1704 // This function is called every mmax iterations to handle limited array sizes.
1705 void cycle(std::vector<X>& Ad, [[maybe_unused]] std::vector<X>& d, [[maybe_unused]] std::vector<field_type,ReboundAllocatorType<X,field_type> >& ddotAd,int& i_bounded) override {
1706 // Only the loop index i_bounded return to 0, if it reached mmax.
1707 i_bounded = 0;
1708 // Now all arrays are filled and the loop in void orthogonalizations can use the whole arrays.
1709 _k_limit = Ad.size();
1710 };
1711
1712 int _k_limit = 0;
1713
1714 protected:
1715 using RestartedFCGSolver<X>::_mmax;
1716 using RestartedFCGSolver<X>::_op;
1717 using RestartedFCGSolver<X>::_prec;
1718 using RestartedFCGSolver<X>::_sp;
1719 using RestartedFCGSolver<X>::_reduction;
1720 using RestartedFCGSolver<X>::_maxit;
1721 using RestartedFCGSolver<X>::_verbose;
1722 };
1723 DUNE_REGISTER_SOLVER("completefcgsolver", defaultIterativeSolverCreator<Dune::CompleteFCGSolver>());
1725} // end namespace
1726
1727#endif
Implementation of the BCRSMatrix class.
Wrapper to use a range of ARPACK++ eigenvalue solvers.
Definition: arpackpp.hh:245
void computeSymMaxMagnitude(const Real &epsilon, BlockVector &x, Real &lambda) const
Assume the matrix to be square, symmetric and perform IRLM to compute an approximation lambda of its ...
Definition: arpackpp.hh:289
void computeSymMinMagnitude(const Real &epsilon, BlockVector &x, Real &lambda) const
Assume the matrix to be square, symmetric and perform IRLM to compute an approximation lambda of its ...
Definition: arpackpp.hh:391
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:468
@ row_wise
Build in a row-wise manner.
Definition: bcrsmatrix.hh:519
CreateIterator createend()
get create iterator pointing to one after the last block
Definition: bcrsmatrix.hh:1107
CreateIterator createbegin()
get initial create iterator
Definition: bcrsmatrix.hh:1101
size_type N() const
number of rows (counted in blocks)
Definition: bcrsmatrix.hh:2008
Bi-conjugate Gradient Stabilized (BiCG-STAB)
Definition: solvers.hh:420
void apply(X &x, X &b, InverseOperatorResult &res) override
Apply inverse operator.
Definition: solvers.hh:440
A vector of blocks with memory management.
Definition: bvector.hh:392
conjugate gradient method
Definition: solvers.hh:193
CGSolver(std::shared_ptr< const LinearOperator< X, X > > op, std::shared_ptr< ScalarProduct< X > > sp, std::shared_ptr< Preconditioner< X, X > > prec, scalar_real_type reduction, int maxit, int verbose, bool condition_estimate)
Constructor to initialize a CG solver.
Definition: solvers.hh:257
void apply(X &x, X &b, InverseOperatorResult &res) override
Apply inverse operator.
Definition: solvers.hh:280
CGSolver(const LinearOperator< X, X > &op, const ScalarProduct< X > &sp, Preconditioner< X, X > &prec, scalar_real_type reduction, int maxit, int verbose, bool condition_estimate)
Constructor to initialize a CG solver.
Definition: solvers.hh:239
CGSolver(const LinearOperator< X, X > &op, Preconditioner< X, X > &prec, scalar_real_type reduction, int maxit, int verbose, bool condition_estimate)
Constructor to initialize a CG solver.
Definition: solvers.hh:222
Complete flexible conjugate gradient method.
Definition: solvers.hh:1669
void apply(X &x, X &b, InverseOperatorResult &res) override
Apply inverse operator.
Definition: solvers.hh:1683
Generalized preconditioned conjugate gradient solver.
Definition: solvers.hh:1298
GeneralizedPCGSolver(const LinearOperator< X, X > &op, const ScalarProduct< X > &sp, Preconditioner< X, X > &prec, scalar_real_type reduction, int maxit, int verbose, int restart=10)
Set up nonlinear preconditioned conjugate gradient solver.
Definition: solvers.hh:1334
GeneralizedPCGSolver(const LinearOperator< X, X > &op, Preconditioner< X, X > &prec, scalar_real_type reduction, int maxit, int verbose, int restart=10)
Set up nonlinear preconditioned conjugate gradient solver.
Definition: solvers.hh:1322
GeneralizedPCGSolver(std::shared_ptr< const LinearOperator< X, X > > op, std::shared_ptr< Preconditioner< X, X > > prec, const ParameterTree &configuration)
Set up nonlinear preconditioned conjugate gradient solver.
Definition: solvers.hh:1352
GeneralizedPCGSolver(std::shared_ptr< const LinearOperator< X, X > > op, std::shared_ptr< const ScalarProduct< X > > sp, std::shared_ptr< Preconditioner< X, X > > prec, scalar_real_type reduction, int maxit, int verbose, int restart=10)
Set up nonlinear preconditioned conjugate gradient solver.
Definition: solvers.hh:1367
void apply(X &x, X &b, InverseOperatorResult &res) override
Apply inverse operator.
Definition: solvers.hh:1381
gradient method
Definition: solvers.hh:124
void apply(X &x, X &b, InverseOperatorResult &res) override
Apply inverse operator.
Definition: solvers.hh:142
Simd::Scalar< real_type > scalar_real_type
scalar type underlying the field_type
Definition: solver.hh:116
Y range_type
Type of the range of the operator to be inverted.
Definition: solver.hh:107
X domain_type
Type of the domain of the operator to be inverted.
Definition: solver.hh:104
X::field_type field_type
The field type of the operator.
Definition: solver.hh:110
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: solver.hh:113
Base class for all implementations of iterative solvers.
Definition: solver.hh:205
IterativeSolver(const LinearOperator< X, X > &op, Preconditioner< X, X > &prec, scalar_real_type reduction, int maxit, int verbose)
General constructor to initialize an iterative solver.
Definition: solver.hh:232
Preconditioned loop solver.
Definition: solvers.hh:59
void apply(X &x, X &b, InverseOperatorResult &res) override
Apply inverse operator,.
Definition: solvers.hh:73
Minimal Residual Method (MINRES)
Definition: solvers.hh:610
void apply(X &x, X &b, InverseOperatorResult &res) override
Apply inverse operator.
Definition: solvers.hh:628
Hierarchical structure of string parameters.
Definition: parametertree.hh:37
Accelerated flexible conjugate gradient method.
Definition: solvers.hh:1490
RestartedFCGSolver(const LinearOperator< X, X > &op, Preconditioner< X, X > &prec, scalar_real_type reduction, int maxit, int verbose, int mmax=10)
Constructor to initialize a RestartedFCG solver.
Definition: solvers.hh:1508
RestartedFCGSolver(std::shared_ptr< const LinearOperator< X, X > > op, std::shared_ptr< const ScalarProduct< X > > sp, std::shared_ptr< Preconditioner< X, X > > prec, scalar_real_type reduction, int maxit, int verbose, int mmax=10)
Constructor to initialize a RestartedFCG solver.
Definition: solvers.hh:1528
void apply(X &x, X &b, InverseOperatorResult &res) override
Apply inverse operator.
Definition: solvers.hh:1573
RestartedFCGSolver(std::shared_ptr< const LinearOperator< X, X > > op, std::shared_ptr< Preconditioner< X, X > > prec, const ParameterTree &configuration)
Constructor to initialize a RestartedFCG solver.
Definition: solvers.hh:1548
RestartedFCGSolver(const LinearOperator< X, X > &op, const ScalarProduct< X > &sp, Preconditioner< X, X > &prec, scalar_real_type reduction, int maxit, int verbose, int mmax=10)
Constructor to initialize a RestartedFCG solver.
Definition: solvers.hh:1518
implements the Flexible Generalized Minimal Residual (FGMRes) method (right preconditioned)
Definition: solvers.hh:1130
void apply(X &x, Y &b, InverseOperatorResult &res) override
Apply inverse operator.
Definition: solvers.hh:1160
implements the Generalized Minimal Residual (GMRes) method
Definition: solvers.hh:828
RestartedGMResSolver(std::shared_ptr< const LinearOperator< X, Y > > op, std::shared_ptr< Preconditioner< X, X > > prec, const ParameterTree &configuration)
Constructor.
Definition: solvers.hh:879
void apply(X &x, Y &b, InverseOperatorResult &res) override
Apply inverse operator.
Definition: solvers.hh:914
ReboundAllocatorType< X, field_type > fAlloc
field_type Allocator retrieved from domain type
Definition: solvers.hh:839
ReboundAllocatorType< X, real_type > rAlloc
real_type Allocator retrieved from domain type
Definition: solvers.hh:841
RestartedGMResSolver(const LinearOperator< X, Y > &op, Preconditioner< X, Y > &prec, scalar_real_type reduction, int restart, int maxit, int verbose)
Set up RestartedGMResSolver solver.
Definition: solvers.hh:851
RestartedGMResSolver(const LinearOperator< X, Y > &op, const ScalarProduct< X > &sp, Preconditioner< X, Y > &prec, scalar_real_type reduction, int restart, int maxit, int verbose)
Set up RestartedGMResSolver solver.
Definition: solvers.hh:862
RestartedGMResSolver(std::shared_ptr< const LinearOperator< X, Y > > op, std::shared_ptr< const ScalarProduct< X > > sp, std::shared_ptr< Preconditioner< X, Y > > prec, scalar_real_type reduction, int restart, int maxit, int verbose)
Set up RestartedGMResSolver solver.
Definition: solvers.hh:895
Base class for scalar product and norm computation.
Definition: scalarproducts.hh:52
Thrown when a solver aborts due to some problem.
Definition: istlexception.hh:46
A few common exception classes.
IO interface of the SIMD abstraction.
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:489
constexpr auto min
Function object that returns the smaller of the given values.
Definition: hybridutilities.hh:511
V cond(M &&mask, const V &ifTrue, const V &ifFalse)
Like the ?: operator.
Definition: interface.hh:386
auto io(const V &v)
construct a stream inserter
Definition: io.hh:106
bool allTrue(const Mask &mask)
Whether all entries are true
Definition: interface.hh:439
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition: interface.hh:235
Some useful basic math stuff.
Dune namespace
Definition: alignedallocator.hh:13
constexpr auto get(std::integer_sequence< T, II... >, std::integral_constant< std::size_t, pos >={})
Return the entry at position pos of the given sequence.
Definition: integersequence.hh:22
Define general, extensible interface for operators. The available implementation wraps a matrix.
Define base class for scalar product and norm.
Include file for users of the SIMD abstraction layer.
Define general, extensible interface for inverse operators.
Statistics about the application of an inverse operator.
Definition: solver.hh:50
double condition_estimate
Estimate of condition number.
Definition: solver.hh:81
void clear()
Resets all data.
Definition: solver.hh:58
bool converged
True if convergence criterion has been met.
Definition: solver.hh:75
A simple timing class.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Apr 2, 22:34, 2026)