3 #ifndef DUNE_PDELAB_BACKEND_ISTL_SEQISTLSOLVERBACKEND_HH 4 #define DUNE_PDELAB_BACKEND_ISTL_SEQISTLSOLVERBACKEND_HH 9 #include <dune/common/deprecated.hh> 10 #include <dune/common/parallel/mpihelper.hh> 12 #include <dune/istl/owneroverlapcopy.hh> 13 #include <dune/istl/solvercategory.hh> 14 #include <dune/istl/operators.hh> 15 #include <dune/istl/solvers.hh> 16 #include <dune/istl/preconditioners.hh> 17 #include <dune/istl/scalarproducts.hh> 18 #include <dune/istl/paamg/amg.hh> 19 #include <dune/istl/paamg/pinfo.hh> 20 #include <dune/istl/io.hh> 21 #include <dune/istl/superlu.hh> 22 #include <dune/istl/umfpack.hh> 37 template<
typename X,
typename Y,
typename GOS>
45 enum {
category=Dune::SolverCategory::sequential};
51 virtual void apply (
const X& x, Y& y)
const 54 gos.jacobian_apply(x,y);
61 gos.jacobian_apply(x,temp);
74 template<
template<
class,
class,
class,
int>
class Preconditioner,
75 template<
class>
class Solver>
86 : maxiter(maxiter_), verbose(verbose_)
98 template<
class M,
class V,
class W>
99 void apply(M& A, V& z, W& r,
typename Dune::template FieldTraits<typename W::ElementType >::real_type reduction)
104 Dune::MatrixAdapter<Native<M>,
106 Native<W>> opa(
native(A));
107 Preconditioner<Native<M>,
110 1> prec(
native(A), 3, 1.0);
111 Solver<Native<V>> solver(opa, prec, reduction, maxiter, verbose);
112 Dune::InverseOperatorResult stat;
114 res.converged = stat.converged;
115 res.iterations = stat.iterations;
116 res.elapsed = stat.elapsed;
117 res.reduction = stat.reduction;
118 res.conv_rate = stat.conv_rate;
126 template<
template<
typename>
class Solver>
137 : maxiter(maxiter_), verbose(verbose_)
146 template<
class M,
class V,
class W>
147 void apply(M& A, V& z, W& r,
typename Dune::template FieldTraits<typename W::ElementType >::real_type reduction)
151 Dune::MatrixAdapter<Native<M>,
153 Native<W>> opa(
native(A));
154 Dune::SeqILU0<Native<M>,
158 Solver<Native<V>> solver(opa, ilu0, reduction, maxiter, verbose);
159 Dune::InverseOperatorResult stat;
161 res.converged = stat.converged;
162 res.iterations = stat.iterations;
163 res.elapsed = stat.elapsed;
164 res.reduction = stat.reduction;
165 res.conv_rate = stat.conv_rate;
172 template<
template<
typename>
class Solver>
184 : n_(n), w_(w), maxiter(maxiter_), verbose(verbose_)
193 template<
class M,
class V,
class W>
194 void apply(M& A, V& z, W& r,
typename Dune::template FieldTraits<typename W::ElementType >::real_type reduction)
198 Dune::MatrixAdapter<Native<M>,
202 Dune::SeqILUn<Native<M>,
205 > ilun(
native(A), n_, w_);
206 Solver<Native<V>> solver(opa, ilun, reduction, maxiter, verbose);
207 Dune::InverseOperatorResult stat;
209 res.converged = stat.converged;
210 res.iterations = stat.iterations;
211 res.elapsed = stat.elapsed;
212 res.reduction = stat.reduction;
213 res.conv_rate = stat.conv_rate;
395 #if HAVE_SUPERLU || DOXYGEN 428 template<
class M,
class V,
class W>
429 void apply(M& A, V& z, W& r,
typename Dune::template FieldTraits<typename W::ElementType >::real_type reduction)
433 using ISTLM = Native<M>;
434 Dune::SuperLU<ISTLM> solver(
native(A), verbose);
435 Dune::InverseOperatorResult stat;
437 res.converged = stat.converged;
438 res.iterations = stat.iterations;
439 res.elapsed = stat.elapsed;
440 res.reduction = stat.reduction;
441 res.conv_rate = stat.conv_rate;
447 #endif // HAVE_SUPERLU || DOXYGEN 449 #if HAVE_SUITESPARSE_UMFPACK || DOXYGEN 482 template<
class M,
class V,
class W>
483 void apply(M& A, V& z, W& r,
typename Dune::template FieldTraits<typename W::ElementType >::real_type reduction)
487 Dune::UMFPack<ISTLM> solver(
native(A), verbose);
488 Dune::InverseOperatorResult stat;
490 res.converged = stat.converged;
491 res.iterations = stat.iterations;
492 res.elapsed = stat.elapsed;
493 res.reduction = stat.reduction;
494 res.conv_rate = stat.conv_rate;
500 #endif // HAVE_SUITESPARSE_UMFPACK || DOXYGEN 519 template<
class M,
class V,
class W>
520 void apply(M& A, V& z, W& r,
typename Dune::template FieldTraits<typename W::ElementType >::real_type reduction)
523 Dune::SeqJac<Native<M>,
530 res.converged =
true;
533 res.reduction = reduction;
534 res.conv_rate = reduction;
563 template<
class GO,
template<
class,
class,
class,
int>
class Preconditioner,
template<
class>
class Solver,
564 bool skipBlocksizeCheck =
false>
567 typedef typename GO::Traits::TrialGridFunctionSpace GFS;
568 typedef typename GO::Traits::Jacobian M;
570 typedef typename GO::Traits::Domain V;
572 typedef Preconditioner<MatrixType,VectorType,VectorType,1> Smoother;
573 typedef Dune::MatrixAdapter<MatrixType,VectorType,VectorType> Operator;
574 typedef typename Dune::Amg::SmootherTraits<Smoother>::Arguments SmootherArgs;
575 typedef Dune::Amg::AMG<Operator,VectorType,Smoother> AMG;
576 typedef Dune::Amg::Parameters Parameters;
580 bool reuse_=
false,
bool usesuperlu_=
true)
581 : maxiter(maxiter_), params(15,2000), verbose(verbose_),
582 reuse(reuse_), firstapply(true), usesuperlu(usesuperlu_)
584 params.setDefaultValuesIsotropic(GFS::Traits::GridViewType::Traits::Grid::dimension);
585 params.setDebugLevel(verbose_);
587 if (usesuperlu ==
true)
589 std::cout <<
"WARNING: You are using AMG without SuperLU!" 590 <<
" Please consider installing SuperLU," 591 <<
" or set the usesuperlu flag to false" 592 <<
" to suppress this warning." << std::endl;
610 typename V::ElementType
norm (
const V& v)
const 622 void apply(M& A, V& z, V& r,
typename Dune::template FieldTraits<typename V::ElementType >::real_type reduction)
626 typedef Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<MatrixType,
627 Dune::Amg::FirstDiagonal> > Criterion;
628 SmootherArgs smootherArgs;
629 smootherArgs.iterations = 1;
630 smootherArgs.relaxationFactor = 1;
632 Criterion criterion(params);
634 if (reuse==
false || firstapply==
true){
635 oop.reset(
new Operator(mat));
636 amg.reset(
new AMG(*oop, criterion, smootherArgs));
638 stats.tsetup = watch.elapsed();
639 stats.levels = amg->maxlevels();
640 stats.directCoarseLevelSolver=amg->usesDirectCoarseLevelSolver();
643 Dune::InverseOperatorResult stat;
645 Solver<VectorType> solver(*oop,*amg,reduction,maxiter,verbose);
647 stats.tsolve= watch.elapsed();
648 res.converged = stat.converged;
649 res.iterations = stat.iterations;
650 res.elapsed = stat.elapsed;
651 res.reduction = stat.reduction;
652 res.conv_rate = stat.conv_rate;
672 std::shared_ptr<Operator> oop;
673 std::shared_ptr<AMG> amg;
700 bool reuse_=
false,
bool usesuperlu_=
true)
702 (maxiter_, verbose_, reuse_, usesuperlu_)
726 bool reuse_=
false,
bool usesuperlu_=
true)
728 (maxiter_, verbose_, reuse_, usesuperlu_)
752 bool reuse_=
false,
bool usesuperlu_=
true)
754 (maxiter_, verbose_, reuse_, usesuperlu_)
778 bool reuse_=
false,
bool usesuperlu_=
true)
780 (maxiter_, verbose_, reuse_, usesuperlu_)
804 bool reuse_=
false,
bool usesuperlu_=
true)
806 (maxiter_, verbose_, reuse_, usesuperlu_)
827 : restart(restart_), maxiter(maxiter_), verbose(verbose_)
837 template<
class M,
class V,
class W>
838 void apply(M& A, V& z, W& r,
typename Dune::template FieldTraits<typename W::ElementType>::real_type reduction)
852 Dune::RestartedGMResSolver<Native<V>> solver(opa,ilu0,reduction,restart,maxiter,verbose);
853 Dune::InverseOperatorResult stat;
855 res.converged = stat.converged;
856 res.iterations = stat.iterations;
857 res.elapsed = stat.elapsed;
858 res.reduction = stat.reduction;
859 res.conv_rate = stat.conv_rate;
863 int restart, maxiter, verbose;
872 #endif // DUNE_PDELAB_BACKEND_ISTL_SEQISTLSOLVERBACKEND_HH ISTLBackend_SEQ_SuperLU(int maxiter, int verbose_)
make a linear solver object
Definition: seqistlsolverbackend.hh:417
ISTLBackend_SEQ_AMG(unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Definition: seqistlsolverbackend.hh:579
void apply(M &A, V &z, W &r, typename Dune::template FieldTraits< typename W::ElementType >::real_type reduction)
solve the given linear system
Definition: seqistlsolverbackend.hh:429
ISTLBackend_SEQ_UMFPack(int verbose_=1)
make a linear solver object
Definition: seqistlsolverbackend.hh:461
Backend using a MINRes solver preconditioned by SSOR.
Definition: seqistlsolverbackend.hh:365
void apply(M &A, V &z, V &r, typename Dune::template FieldTraits< typename V::ElementType >::real_type reduction)
solve the given linear system
Definition: seqistlsolverbackend.hh:622
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:176
void apply(M &A, V &z, W &r, typename Dune::template FieldTraits< typename W::ElementType >::real_type reduction)
solve the given linear system
Definition: seqistlsolverbackend.hh:838
Class providing some statistics of the AMG solver.
Definition: seqistlsolverbackend.hh:544
int levels
the number of levels in the AMG hierarchy.
Definition: seqistlsolverbackend.hh:552
Solver backend using SuperLU as a direct solver.
Definition: seqistlsolverbackend.hh:399
ISTLBackend_SEQ_BCGS_AMG_SSOR(unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Constructor.
Definition: seqistlsolverbackend.hh:725
ISTLBackend_SEQ_UMFPack(int maxiter, int verbose_)
make a linear solver object
Definition: seqistlsolverbackend.hh:471
ISTLBackend_SEQ_BCGS_AMG_SOR(unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Constructor.
Definition: seqistlsolverbackend.hh:751
V::ElementType norm(const V &v) const
compute global norm of a vector
Definition: seqistlsolverbackend.hh:610
Sequential Loop solver preconditioned with AMG smoothed by SSOR.
Definition: seqistlsolverbackend.hh:764
void apply(M &A, V &z, W &r, typename Dune::template FieldTraits< typename W::ElementType >::real_type reduction)
solve the given linear system
Definition: seqistlsolverbackend.hh:483
ISTLBackend_SEQ_BCGS_Jac(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: seqistlsolverbackend.hh:253
Y range_type
Definition: seqistlsolverbackend.hh:42
Sequential conjugate gradient solver preconditioned with AMG smoothed by SSOR.
Definition: seqistlsolverbackend.hh:686
ISTLBackend_SEQ_CG_SSOR(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: seqistlsolverbackend.hh:357
ISTLBackend_SEQ_Base(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: seqistlsolverbackend.hh:85
const ISTLAMGStatistics & statistics() const
Get statistics of the AMG solver (no of levels, timings).
Definition: seqistlsolverbackend.hh:660
Definition: seqistlsolverbackend.hh:127
ISTLBackend_SEQ_CG_ILUn(int n_, double w_=1.0, unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: seqistlsolverbackend.hh:340
For backward compatibility – Do not use this!
Definition: adaptivity.hh:27
Linear solver backend for Restarted GMRes preconditioned with ILU(0)
Definition: seqistlsolverbackend.hh:815
ISTLBackend_SEQ_GMRES_ILU0(int restart_=200, int maxiter_=5000, int verbose_=1)
make linear solver object
Definition: seqistlsolverbackend.hh:826
Backend for sequential BiCGSTAB solver with ILU0 preconditioner.
Definition: seqistlsolverbackend.hh:278
Solver to be used for explicit time-steppers with (block-)diagonal mass matrix.
Definition: seqistlsolverbackend.hh:503
double tsetup
The time needed for building the AMG hierarchy (coarsening).
Definition: seqistlsolverbackend.hh:556
Backend for sequential loop solver with Jacobi preconditioner.
Definition: seqistlsolverbackend.hh:229
Backend for conjugate gradient solver with Jacobi preconditioner.
Definition: seqistlsolverbackend.hh:382
void setparams(Parameters params_)
set AMG parameters
Definition: seqistlsolverbackend.hh:601
OnTheFlyOperator(GOS &gos_)
Definition: seqistlsolverbackend.hh:47
ISTLBackend_SEQ_CG_AMG_SSOR(unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Constructor.
Definition: seqistlsolverbackend.hh:699
Sequential BiCGSTAB solver preconditioned with AMG smoothed by SOR.
Definition: seqistlsolverbackend.hh:738
ISTLBackend_SEQ_LS_AMG_SSOR(unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Constructor.
Definition: seqistlsolverbackend.hh:777
Sequential BiCGStab solver with ILU0 preconditioner.
Definition: seqistlsolverbackend.hh:310
Definition: seqistlsolverbackend.hh:565
virtual void applyscaleadd(field_type alpha, const X &x, Y &y) const
Definition: seqistlsolverbackend.hh:57
Solver backend using UMFPack as a direct solver.
Definition: seqistlsolverbackend.hh:453
ISTLBackend_SEQ_LS_AMG_SOR(unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Constructor.
Definition: seqistlsolverbackend.hh:803
ISTLBackend_SEQ_LOOP_Jac(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: seqistlsolverbackend.hh:237
ISTLBackend_SEQ_MINRES_SSOR(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: seqistlsolverbackend.hh:374
VTKWriter & w
Definition: function.hh:1100
ISTLBackend_SEQ_BCGS_SSOR(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: seqistlsolverbackend.hh:270
int iterations
The number of iterations performed until convergence was reached.
Definition: seqistlsolverbackend.hh:558
Definition: seqistlsolverbackend.hh:76
Sequential Loop solver preconditioned with AMG smoothed by SOR.
Definition: seqistlsolverbackend.hh:790
Sequential BiCGStab solver preconditioned with AMG smoothed by SSOR.
Definition: seqistlsolverbackend.hh:712
ISTLBackend_SEQ_CG_Jac(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: seqistlsolverbackend.hh:390
ISTLBackend_SEQ_BCGS_ILUn(int n_, double w_=1.0, unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: seqistlsolverbackend.hh:322
Backend for sequential BiCGSTAB solver with SSOR preconditioner.
Definition: seqistlsolverbackend.hh:261
Definition: seqistlsolverbackend.hh:45
void apply(M &A, V &z, W &r, typename Dune::template FieldTraits< typename W::ElementType >::real_type reduction)
solve the given linear system
Definition: seqistlsolverbackend.hh:520
X domain_type
Definition: seqistlsolverbackend.hh:41
Definition: seqistlsolverbackend.hh:173
Definition: seqistlsolverbackend.hh:38
double tsolve
The time spent in solving the system (without building the hierarchy.
Definition: seqistlsolverbackend.hh:554
void apply(M &A, V &z, W &r, typename Dune::template FieldTraits< typename W::ElementType >::real_type reduction)
solve the given linear system
Definition: seqistlsolverbackend.hh:99
virtual void apply(const X &x, Y &y) const
Definition: seqistlsolverbackend.hh:51
X::field_type field_type
Definition: seqistlsolverbackend.hh:43
ISTLBackend_SEQ_ExplicitDiagonal()
make a linear solver object
Definition: seqistlsolverbackend.hh:509
void apply(M &A, V &z, W &r, typename Dune::template FieldTraits< typename W::ElementType >::real_type reduction)
solve the given linear system
Definition: seqistlsolverbackend.hh:147
ISTLBackend_SEQ_SuperLU(int verbose_=1)
make a linear solver object
Definition: seqistlsolverbackend.hh:407
Backend for sequential conjugate gradient solver with SSOR preconditioner.
Definition: seqistlsolverbackend.hh:348
Backend for sequential conjugate gradient solver with ILU0 preconditioner.
Definition: seqistlsolverbackend.hh:295
ISTLBackend_SEQ_BCGS_ILU0(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: seqistlsolverbackend.hh:287
double tprepare
The needed for computing the parallel information and for adapting the linear system.
Definition: seqistlsolverbackend.hh:550
std::enable_if< std::is_base_of< impl::WrapperBase, T >::value, Native< T > &>::type native(T &t)
Definition: backend/interface.hh:192
Backend for sequential BiCGSTAB solver with Jacobi preconditioner.
Definition: seqistlsolverbackend.hh:245
ISTLBackend_SEQ_CG_ILU0(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: seqistlsolverbackend.hh:304
ISTLBackend_SEQ_ILUn(int n, double w, unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: seqistlsolverbackend.hh:183
ISTLBackend_SEQ_ILU0(unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: seqistlsolverbackend.hh:136
Sequential congute gradient solver with ILU0 preconditioner.
Definition: seqistlsolverbackend.hh:328
void apply(M &A, V &z, W &r, typename Dune::template FieldTraits< typename W::ElementType >::real_type reduction)
solve the given linear system
Definition: seqistlsolverbackend.hh:194
bool directCoarseLevelSolver
True if a direct solver was used on the coarset level.
Definition: seqistlsolverbackend.hh:560