3 #ifndef DUNE_SEQISTLSOLVERBACKEND_HH
4 #define DUNE_SEQISTLSOLVERBACKEND_HH
6 #include <dune/common/deprecated.hh>
7 #include <dune/common/parallel/mpihelper.hh>
9 #include <dune/istl/owneroverlapcopy.hh>
10 #include <dune/istl/solvercategory.hh>
11 #include <dune/istl/operators.hh>
12 #include <dune/istl/solvers.hh>
13 #include <dune/istl/preconditioners.hh>
14 #include <dune/istl/scalarproducts.hh>
15 #include <dune/istl/paamg/amg.hh>
16 #include <dune/istl/paamg/pinfo.hh>
17 #include <dune/istl/io.hh>
18 #include <dune/istl/superlu.hh>
19 #include <dune/istl/umfpack.hh>
34 template<
typename X,
typename Y,
typename GOS>
42 enum {
category=Dune::SolverCategory::sequential};
48 virtual void apply (
const X& x, Y& y)
const
51 gos.jacobian_apply(x,y);
58 gos.jacobian_apply(x,temp);
71 template<
template<
class,
class,
class,
int>
class Preconditioner,
72 template<
class>
class Solver>
83 : maxiter(maxiter_), verbose(verbose_)
95 template<
class M,
class V,
class W>
96 void apply(M& A, V& z, W& r,
typename W::ElementType reduction)
101 Dune::MatrixAdapter<Native<M>,
103 Native<W>> opa(
native(A));
104 Preconditioner<Native<M>,
107 1> prec(
native(A), 3, 1.0);
108 Solver<Native<V>> solver(opa, prec, reduction, maxiter, verbose);
109 Dune::InverseOperatorResult stat;
123 template<
template<
typename>
class Solver>
134 : maxiter(maxiter_), verbose(verbose_)
143 template<
class M,
class V,
class W>
144 void apply(M& A, V& z, W& r,
typename Dune::template FieldTraits<typename W::ElementType >::real_type reduction)
148 Dune::MatrixAdapter<Native<M>,
150 Native<W>> opa(
native(A));
151 Dune::SeqILU0<Native<M>,
155 Solver<Native<V>> solver(opa, ilu0, reduction, maxiter, verbose);
156 Dune::InverseOperatorResult stat;
169 template<
template<
typename>
class Solver>
181 : n_(n), w_(w), maxiter(maxiter_), verbose(verbose_)
190 template<
class M,
class V,
class W>
191 void apply(M& A, V& z, W& r,
typename W::ElementType reduction)
195 Dune::MatrixAdapter<Native<M>,
199 Dune::SeqILUn<Native<M>,
202 > ilun(
native(A), n_, w_);
203 Solver<Native<V>> solver(opa, ilun, reduction, maxiter, verbose);
204 Dune::InverseOperatorResult stat;
392 #if HAVE_SUPERLU || DOXYGEN
425 template<
class M,
class V,
class W>
426 void apply(M& A, V& z, W& r,
typename W::ElementType reduction)
430 using ISTLM = Native<M>;
431 Dune::SuperLU<ISTLM> solver(
native(A), verbose);
432 Dune::InverseOperatorResult stat;
444 #endif // HAVE_SUPERLU || DOXYGEN
446 #if HAVE_UMFPACK || DOXYGEN
479 template<
class M,
class V,
class W>
480 void apply(M& A, V& z, W& r,
typename W::ElementType reduction)
484 Dune::UMFPack<ISTLM> solver(
native(A), verbose);
485 Dune::InverseOperatorResult stat;
497 #endif // HAVE_UMFPACK || DOXYGEN
516 template<
class M,
class V,
class W>
517 void apply(M& A, V& z, W& r,
typename W::ElementType reduction)
520 Dune::SeqJac<Native<M>,
560 template<
class GO,
template<
class,
class,
class,
int>
class Preconditioner,
template<
class>
class Solver,
561 bool skipBlocksizeCheck =
false>
564 typedef typename GO::Traits::TrialGridFunctionSpace GFS;
565 typedef typename GO::Traits::Jacobian M;
567 typedef typename GO::Traits::Domain V;
569 typedef Preconditioner<MatrixType,VectorType,VectorType,1> Smoother;
570 typedef Dune::MatrixAdapter<MatrixType,VectorType,VectorType> Operator;
571 typedef typename Dune::Amg::SmootherTraits<Smoother>::Arguments SmootherArgs;
572 typedef Dune::Amg::AMG<Operator,VectorType,Smoother> AMG;
573 typedef Dune::Amg::Parameters Parameters;
577 bool reuse_=
false,
bool usesuperlu_=
true)
578 : maxiter(maxiter_), params(15,2000), verbose(verbose_),
579 reuse(reuse_), firstapply(true), usesuperlu(usesuperlu_)
581 params.setDefaultValuesIsotropic(GFS::Traits::GridViewType::Traits::Grid::dimension);
582 params.setDebugLevel(verbose_);
584 if (usesuperlu ==
true)
586 std::cout <<
"WARNING: You are using AMG without SuperLU!"
587 <<
" Please consider installing SuperLU,"
588 <<
" or set the usesuperlu flag to false"
589 <<
" to suppress this warning." << std::endl;
607 typename V::ElementType
norm (
const V& v)
const
619 void apply(M& A, V& z, V& r,
typename V::ElementType reduction)
623 typedef Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<MatrixType,
624 Dune::Amg::FirstDiagonal> > Criterion;
625 SmootherArgs smootherArgs;
626 smootherArgs.iterations = 1;
627 smootherArgs.relaxationFactor = 1;
629 Criterion criterion(params);
632 if (reuse==
false || firstapply==
true){
633 amg.reset(
new AMG(oop, criterion, smootherArgs));
635 stats.
tsetup = watch.elapsed();
636 stats.
levels = amg->maxlevels();
640 Dune::InverseOperatorResult stat;
642 Solver<VectorType> solver(oop,*amg,reduction,maxiter,verbose);
644 stats.
tsolve= watch.elapsed();
669 std::shared_ptr<AMG> amg;
696 bool reuse_=
false,
bool usesuperlu_=
true)
698 (maxiter_, verbose_, reuse_, usesuperlu_)
722 bool reuse_=
false,
bool usesuperlu_=
true)
724 (maxiter_, verbose_, reuse_, usesuperlu_)
748 bool reuse_=
false,
bool usesuperlu_=
true)
750 (maxiter_, verbose_, reuse_, usesuperlu_)
774 bool reuse_=
false,
bool usesuperlu_=
true)
776 (maxiter_, verbose_, reuse_, usesuperlu_)
800 bool reuse_=
false,
bool usesuperlu_=
true)
802 (maxiter_, verbose_, reuse_, usesuperlu_)
823 : restart(restart_), maxiter(maxiter_), verbose(verbose_)
833 template<
class M,
class V,
class W>
834 void apply(M& A, V& z, W& r,
typename Dune::template FieldTraits<typename W::ElementType>::real_type reduction)
848 Dune::RestartedGMResSolver<Native<V>> solver(opa,ilu0,reduction,restart,maxiter,verbose);
849 Dune::InverseOperatorResult stat;
859 int restart, maxiter, verbose;
void apply(M &A, V &z, W &r, typename W::ElementType reduction)
solve the given linear system
Definition: istl/seqistlsolverbackend.hh:96
ISTLBackend_SEQ_LOOP_Jac(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:234
ISTLBackend_SEQ_BCGS_ILU0(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:284
void apply(M &A, V &z, W &r, typename W::ElementType reduction)
solve the given linear system
Definition: istl/seqistlsolverbackend.hh:426
const ISTLAMGStatistics & statistics() const
Get statistics of the AMG solver (no of levels, timings).
Definition: istl/seqistlsolverbackend.hh:657
Backend for sequential BiCGSTAB solver with Jacobi preconditioner.
Definition: istl/seqistlsolverbackend.hh:242
void apply(M &A, V &z, W &r, typename Dune::template FieldTraits< typename W::ElementType >::real_type reduction)
solve the given linear system
Definition: istl/seqistlsolverbackend.hh:144
virtual void apply(const X &x, Y &y) const
Definition: istl/seqistlsolverbackend.hh:48
double tsetup
The time needed for building the AMG hierarchy (coarsening).
Definition: istl/seqistlsolverbackend.hh:553
void setparams(Parameters params_)
set AMG parameters
Definition: istl/seqistlsolverbackend.hh:598
Linear solver backend for Restarted GMRes preconditioned with ILU(0)
Definition: istl/seqistlsolverbackend.hh:811
Sequential BiCGSTAB solver preconditioned with AMG smoothed by SOR.
Definition: istl/seqistlsolverbackend.hh:734
RFType reduction
Definition: solver.hh:35
Backend for sequential loop solver with Jacobi preconditioner.
Definition: istl/seqistlsolverbackend.hh:226
Sequential BiCGStab solver with ILU0 preconditioner.
Definition: istl/seqistlsolverbackend.hh:307
double tprepare
The needed for computing the parallel information and for adapting the linear system.
Definition: istl/seqistlsolverbackend.hh:547
Solver backend using UMFPack as a direct solver.
Definition: istl/seqistlsolverbackend.hh:450
Backend using a MINRes solver preconditioned by SSOR.
Definition: istl/seqistlsolverbackend.hh:362
ISTLBackend_SEQ_LS_AMG_SSOR(unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Constructor.
Definition: istl/seqistlsolverbackend.hh:773
void apply(M &A, V &z, W &r, typename W::ElementType reduction)
solve the given linear system
Definition: istl/seqistlsolverbackend.hh:480
ISTLBackend_SEQ_BCGS_SSOR(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:267
Sequential conjugate gradient solver preconditioned with AMG smoothed by SSOR.
Definition: istl/seqistlsolverbackend.hh:682
unsigned int iterations
Definition: solver.hh:33
ISTLBackend_SEQ_AMG(unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Definition: istl/seqistlsolverbackend.hh:576
ISTLBackend_SEQ_BCGS_AMG_SSOR(unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Constructor.
Definition: istl/seqistlsolverbackend.hh:721
ISTLBackend_SEQ_UMFPack(int maxiter, int verbose_)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:468
Solver to be used for explicit time-steppers with (block-)diagonal mass matrix.
Definition: istl/seqistlsolverbackend.hh:500
void apply(M &A, V &z, W &r, typename Dune::template FieldTraits< typename W::ElementType >::real_type reduction)
solve the given linear system
Definition: istl/seqistlsolverbackend.hh:834
ISTLBackend_SEQ_CG_SSOR(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:354
Definition: istl/seqistlsolverbackend.hh:73
ISTLBackend_SEQ_BCGS_ILUn(int n_, double w_=1.0, unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:319
ISTLBackend_SEQ_LS_AMG_SOR(unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Constructor.
Definition: istl/seqistlsolverbackend.hh:799
Definition: istl/seqistlsolverbackend.hh:562
virtual void applyscaleadd(field_type alpha, const X &x, Y &y) const
Definition: istl/seqistlsolverbackend.hh:54
bool directCoarseLevelSolver
True if a direct solver was used on the coarset level.
Definition: istl/seqistlsolverbackend.hh:557
std::enable_if< std::is_base_of< impl::WrapperBase, T >::value, Native< T > & >::type native(T &t)
Definition: backend/interface.hh:198
Definition: istl/seqistlsolverbackend.hh:42
bool converged
Definition: solver.hh:32
Solver backend using SuperLU as a direct solver.
Definition: istl/seqistlsolverbackend.hh:396
Sequential Loop solver preconditioned with AMG smoothed by SSOR.
Definition: istl/seqistlsolverbackend.hh:760
ISTLBackend_SEQ_CG_Jac(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:387
RFType conv_rate
Definition: solver.hh:36
double elapsed
Definition: solver.hh:34
ISTLBackend_SEQ_CG_ILU0(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:301
ISTLBackend_SEQ_CG_AMG_SSOR(unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Constructor.
Definition: istl/seqistlsolverbackend.hh:695
ISTLBackend_SEQ_UMFPack(int verbose_=1)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:458
ISTLBackend_SEQ_ExplicitDiagonal()
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:506
void apply(M &A, V &z, W &r, typename W::ElementType reduction)
solve the given linear system
Definition: istl/seqistlsolverbackend.hh:517
ISTLBackend_SEQ_SuperLU(int maxiter, int verbose_)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:414
Definition: istl/seqistlsolverbackend.hh:170
Backend for sequential BiCGSTAB solver with ILU0 preconditioner.
Definition: istl/seqistlsolverbackend.hh:275
Definition: adaptivity.hh:27
ISTLBackend_SEQ_GMRES_ILU0(int restart_=200, int maxiter_=5000, int verbose_=1)
make linear solver object
Definition: istl/seqistlsolverbackend.hh:822
typename native_type< T >::type Native
Alias of the native container type associated with T or T itself if it is not a backend wrapper...
Definition: backend/interface.hh:182
Backend for sequential BiCGSTAB solver with SSOR preconditioner.
Definition: istl/seqistlsolverbackend.hh:258
X domain_type
Definition: istl/seqistlsolverbackend.hh:38
Definition: istl/seqistlsolverbackend.hh:35
V::ElementType norm(const V &v) const
compute global norm of a vector
Definition: istl/seqistlsolverbackend.hh:607
ISTLBackend_SEQ_BCGS_Jac(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:250
ISTLBackend_SEQ_Base(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:82
int levels
the number of levels in the AMG hierarchy.
Definition: istl/seqistlsolverbackend.hh:549
Backend for sequential conjugate gradient solver with SSOR preconditioner.
Definition: istl/seqistlsolverbackend.hh:345
int iterations
The number of iterations performed until convergence was reached.
Definition: istl/seqistlsolverbackend.hh:555
ISTLBackend_SEQ_MINRES_SSOR(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:371
ISTLBackend_SEQ_CG_ILUn(int n_, double w_=1.0, unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:337
Sequential congute gradient solver with ILU0 preconditioner.
Definition: istl/seqistlsolverbackend.hh:325
ISTLBackend_SEQ_SuperLU(int verbose_=1)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:404
Class providing some statistics of the AMG solver.
Definition: istl/seqistlsolverbackend.hh:541
ISTLBackend_SEQ_ILU0(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:133
Dune::PDELab::LinearSolverResult< double > res
Definition: solver.hh:52
Definition: istl/seqistlsolverbackend.hh:124
X::field_type field_type
Definition: istl/seqistlsolverbackend.hh:40
Sequential BiCGStab solver preconditioned with AMG smoothed by SSOR.
Definition: istl/seqistlsolverbackend.hh:708
Sequential Loop solver preconditioned with AMG smoothed by SOR.
Definition: istl/seqistlsolverbackend.hh:786
ISTLBackend_SEQ_BCGS_AMG_SOR(unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Constructor.
Definition: istl/seqistlsolverbackend.hh:747
void apply(M &A, V &z, W &r, typename W::ElementType reduction)
solve the given linear system
Definition: istl/seqistlsolverbackend.hh:191
Y range_type
Definition: istl/seqistlsolverbackend.hh:39
VTKWriter & w
Definition: function.hh:1101
double tsolve
The time spent in solving the system (without building the hierarchy.
Definition: istl/seqistlsolverbackend.hh:551
void apply(M &A, V &z, V &r, typename V::ElementType reduction)
solve the given linear system
Definition: istl/seqistlsolverbackend.hh:619
ISTLBackend_SEQ_ILUn(int n, double w, unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/seqistlsolverbackend.hh:180
OnTheFlyOperator(GOS &gos_)
Definition: istl/seqistlsolverbackend.hh:44
Backend for conjugate gradient solver with Jacobi preconditioner.
Definition: istl/seqistlsolverbackend.hh:379
Backend for sequential conjugate gradient solver with ILU0 preconditioner.
Definition: istl/seqistlsolverbackend.hh:292