3 #ifndef DUNE_OVLPISTLSOLVERBACKEND_HH
4 #define DUNE_OVLPISTLSOLVERBACKEND_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>
40 template<
class CC,
class M,
class X,
class Y>
42 :
public Dune::AssembledLinearOperator<M,X,Y>
52 enum {
category=Dune::SolverCategory::overlapping};
59 virtual void apply (
const domain_type& x, range_type& y)
const
67 virtual void applyscaleadd (field_type alpha,
const domain_type& x, range_type& y)
const
87 template<
class GFS,
class X>
89 :
public Dune::ScalarProduct<X>
97 enum {
category=Dune::SolverCategory::overlapping};
102 : gfs(gfs_), helper(helper_)
110 virtual field_type
dot (
const X& x,
const X& y)
113 field_type sum = helper.disjointDot(x,y);
116 return gfs.gridView().comm().sum(sum);
122 virtual double norm (
const X& x)
124 return sqrt(static_cast<double>(this->
dot(x,x)));
133 template<
class CC,
class GFS,
class P>
135 :
public Dune::Preconditioner<Dune::PDELab::Backend::Vector<GFS,typename P::domain_type::field_type>,
136 Dune::PDELab::Backend::Vector<GFS,typename P::range_type::field_type>>
153 : gfs(gfs_), prec(prec_), cc(cc_), helper(helper_)
173 if (gfs.gridView().comm().size()>1)
174 gfs.gridView().communicate(adddh,Dune::All_All_Interface,Dune::ForwardCommunication);
195 template<
class GFS,
class M,
class X,
class Y>
196 class SuperLUSubdomainSolver :
public Dune::Preconditioner<X,Y>
198 typedef Backend::Native<M> ISTLM;
202 typedef X domain_type;
204 typedef Y range_type;
206 typedef typename X::ElementType field_type;
212 category=Dune::SolverCategory::overlapping
221 SuperLUSubdomainSolver (
const GFS& gfs_,
const M& A_)
222 : gfs(gfs_), solver(Backend::
native(A_),false)
228 virtual void pre (X& x, Y& b) {}
233 virtual void apply (X& v,
const Y& d)
235 Dune::InverseOperatorResult stat;
238 if (gfs.gridView().comm().size()>1)
240 AddDataHandle<GFS,X> adddh(gfs,v);
241 gfs.gridView().communicate(adddh,Dune::All_All_Interface,Dune::ForwardCommunication);
248 virtual void post (X& x) {}
252 Dune::SuperLU<ISTLM> solver;
256 template<
class GFS,
class M,
class X,
class Y>
257 class RestrictedSuperLUSubdomainSolver :
public Dune::Preconditioner<X,Y>
259 typedef typename M::BaseT ISTLM;
263 typedef X domain_type;
265 typedef Y range_type;
267 typedef typename X::ElementType field_type;
273 category=Dune::SolverCategory::overlapping
283 RestrictedSuperLUSubdomainSolver (
const GFS& gfs_,
const M& A_,
284 const istl::ParallelHelper<GFS>& helper_)
285 : gfs(gfs_), solver(Backend::
native(A_),false), helper(helper_)
291 virtual void pre (X& x, Y& b) {}
296 virtual void apply (X& v,
const Y& d)
299 Dune::InverseOperatorResult stat;
302 if (gfs.gridView().comm().size()>1)
304 helper.maskForeignDOFs(
native(v));
305 AddDataHandle<GFS,X> adddh(gfs,v);
306 gfs.gridView().communicate(adddh,Dune::InteriorBorder_All_Interface,Dune::ForwardCommunication);
313 virtual void post (X& x) {}
317 Dune::SuperLU<ISTLM> solver;
318 const istl::ParallelHelper<GFS>& helper;
322 template<
typename GFS>
327 : gfs(gfs_), helper(gfs_)
335 typename X::ElementType
dot (
const X& x,
const X& y)
const
338 typename X::ElementType sum = helper.disjointDot(x,y);
341 return gfs.gridView().comm().sum(sum);
348 typename X::ElementType
norm (
const X& x)
const
350 return sqrt(static_cast<double>(this->
dot(x,x)));
370 template<
typename GFS,
typename X>
372 :
public ScalarProduct<X>
377 : implementation(implementation_)
380 virtual typename X::BaseT::field_type
dot(
const X& x,
const X& y)
382 return implementation.dot(x,y);
385 virtual typename X::BaseT::field_type
norm (
const X& x)
387 return sqrt(static_cast<double>(this->
dot(x,x)));
394 template<
class GFS,
class C,
395 template<
class,
class,
class,
int>
class Preconditioner,
396 template<
class>
class Solver>
410 int steps_=5,
int verbose_=1)
421 template<
class M,
class V,
class W>
422 void apply(M& A, V& z, W& r,
typename V::ElementType reduction)
430 typedef Preconditioner<
436 SeqPrec seqprec(
native(A),steps,1.0);
440 if (gfs.gridView().comm().rank()==0) verb=verbose;
441 Solver<V> solver(pop,psp,wprec,reduction,maxiter,verb);
442 Dune::InverseOperatorResult stat;
443 solver.apply(z,r,stat);
459 template<
class GFS,
class C,
460 template<
class>
class Solver>
483 template<
class M,
class V,
class W>
484 void apply(M& A, V& z, W& r,
typename V::ElementType reduction)
498 SeqPrec seqprec(
native(A),1.0);
502 if (gfs.gridView().comm().rank()==0) verb=verbose;
503 Solver<V> solver(pop,psp,wprec,reduction,maxiter,verb);
504 Dune::InverseOperatorResult stat;
505 solver.apply(z,r,stat);
521 template<
class GFS,
class C,
522 template<
class>
class Solver>
546 template<
class M,
class V,
class W>
547 void apply(M& A, V& z, W& r,
typename V::ElementType reduction)
561 SeqPrec seqprec(
native(A),n,1.0);
565 if (gfs.gridView().comm().rank()==0) verb=verbose;
566 Solver<V> solver(pop,psp,wprec,reduction,maxiter,verb);
567 Dune::InverseOperatorResult stat;
568 solver.apply(z,r,stat);
592 template<
class GFS,
class CC>
606 int steps=5,
int verbose=1)
615 template<
class GFS,
class CC>
636 template<
class GFS,
class CC>
658 template<
class GFS,
class CC>
672 int steps=5,
int verbose=1)
682 template<
class GFS,
class CC>
706 template<
class M,
class V,
class W>
707 void apply(M& A, V& z, W& r,
typename V::ElementType reduction)
721 SeqPrec seqprec(
native(A),1.0);
725 if (gfs.gridView().comm().rank()==0) verb=verbose;
726 RestartedGMResSolver<V> solver(pop,psp,wprec,reduction,restart,maxiter,verb);
727 Dune::InverseOperatorResult stat;
728 solver.apply(z,r,stat);
747 template<
class GFS,
class C,
template<
typename>
class Solver>
771 template<
class M,
class V,
class W>
772 void apply(M& A, V& z, W& r,
typename V::ElementType reduction)
779 typedef SuperLUSubdomainSolver<GFS,M,V,W> PREC;
782 if (gfs.gridView().comm().rank()==0) verb=verbose;
783 Solver<V> solver(pop,psp,prec,reduction,maxiter,verb);
784 Dune::InverseOperatorResult stat;
785 solver.apply(z,r,stat);
792 std::cout <<
"No superLU support, please install and configure it." << std::endl;
810 template<
class GFS,
class CC>
834 template<
class GFS,
class CC>
848 unsigned maxiter_=5000,
880 typename V::ElementType
norm(
const V& v)
const
884 "ISTLBackend_OVLP_ExplicitDiagonal::norm() should not be "
885 "neccessary, so we skipped the implementation. If you have a "
886 "scenario where you need it, please implement it or report back to "
897 template<
class M,
class V,
class W>
898 void apply(M& A, V& z, W& r,
typename W::ElementType reduction)
910 if (gfs.gridView().comm().size()>1)
913 gfs.gridView().communicate(copydh,Dune::InteriorBorder_All_Interface,Dune::ForwardCommunication);
927 template<
class GO,
int s,
template<
class,
class,
class,
int>
class Preconditioner,
928 template<
class>
class Solver>
931 typedef typename GO::Traits::TrialGridFunctionSpace GFS;
933 typedef typename GO::Traits::Jacobian M;
935 typedef typename GO::Traits::Domain V;
939 typedef Preconditioner<MatrixType,VectorType,VectorType,1> Smoother;
940 typedef Dune::BlockPreconditioner<VectorType,VectorType,Comm,Smoother> ParSmoother;
941 typedef Dune::OverlappingSchwarzOperator<MatrixType,VectorType,VectorType,Comm> Operator;
943 typedef Preconditioner<MatrixType,VectorType,VectorType,1> ParSmoother;
944 typedef Dune::MatrixAdapter<MatrixType,VectorType,VectorType> Operator;
946 typedef typename Dune::Amg::SmootherTraits<ParSmoother>::Arguments SmootherArgs;
947 typedef Dune::Amg::AMG<Operator,VectorType,ParSmoother,Comm> AMG;
949 typedef typename V::ElementType RF;
960 int verbose_=1,
bool reuse_=
false,
961 bool usesuperlu_=
true)
962 : gfs(gfs_), phelper(gfs,verbose_), maxiter(maxiter_), params(15,2000),
963 verbose(verbose_), reuse(reuse_), firstapply(true),
964 usesuperlu(usesuperlu_)
966 params.setDefaultValuesIsotropic(GFS::Traits::GridViewType::Traits::Grid::dimension);
967 params.setDebugLevel(verbose_);
969 if (gfs.gridView().comm().rank() == 0 && usesuperlu ==
true)
971 std::cout <<
"WARNING: You are using AMG without SuperLU!"
972 <<
" Please consider installing SuperLU,"
973 <<
" or set the usesuperlu flag to false"
974 <<
" to suppress this warning." << std::endl;
988 void setparams(Parameters params_) DUNE_DEPRECATED_MSG(
"setparams() is deprecated, use setParameters() instead")
1009 typename V::ElementType
norm (
const V& v)
const
1012 PSP psp(gfs,phelper);
1023 void apply(M& A, V& z, V& r,
typename V::ElementType reduction)
1026 Comm oocc(gfs.gridView().comm());
1028 typedef Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<MatrixType,
1029 Dune::Amg::FirstDiagonal> > Criterion;
1032 Operator oop(mat, oocc);
1033 Dune::OverlappingSchwarzScalarProduct<VectorType,Comm> sp(oocc);
1036 Dune::SeqScalarProduct<VectorType> sp;
1038 SmootherArgs smootherArgs;
1039 smootherArgs.iterations = 1;
1040 smootherArgs.relaxationFactor = 1;
1041 Criterion criterion(params);
1046 if (gfs.gridView().comm().rank()==0) verb=verbose;
1048 if (reuse==
false || firstapply==
true){
1049 amg.reset(
new AMG(oop, criterion, smootherArgs, oocc));
1051 stats.
tsetup = watch.elapsed();
1052 stats.
levels = amg->maxlevels();
1056 Solver<VectorType> solver(oop,sp,*amg,RF(reduction),maxiter,verb);
1057 Dune::InverseOperatorResult stat;
1060 stats.
tsolve= watch.elapsed();
1086 shared_ptr<AMG> amg;
1099 template<
class GO,
int s=96>
1103 typedef typename GO::Traits::TrialGridFunctionSpace GFS;
1115 int verbose_=1,
bool reuse_=
false,
1116 bool usesuperlu_=
true)
1118 (gfs_, maxiter_, verbose_, reuse_, usesuperlu_)
1128 template<
class GO,
int s=96>
1132 typedef typename GO::Traits::TrialGridFunctionSpace GFS;
1144 int verbose_=1,
bool reuse_=
false,
1145 bool usesuperlu_=
true)
1147 (gfs_, maxiter_, verbose_, reuse_, usesuperlu_)
1157 template<
class GO,
int s=96>
1161 typedef typename GO::Traits::TrialGridFunctionSpace GFS;
1173 int verbose_=1,
bool reuse_=
false,
1174 bool usesuperlu_=
true)
1176 (gfs_, maxiter_, verbose_, reuse_, usesuperlu_)
virtual X::BaseT::field_type dot(const X &x, const X &y)
Definition: istl/ovlpistlsolverbackend.hh:380
M matrix_type
export types
Definition: istl/ovlpistlsolverbackend.hh:46
ISTLBackend_AMG(const GFS &gfs_, unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Definition: istl/ovlpistlsolverbackend.hh:959
void apply(M &A, V &z, W &r, typename V::ElementType reduction)
solve the given linear system
Definition: istl/ovlpistlsolverbackend.hh:772
ISTLBackend_CG_AMG_SSOR(const GFS &gfs_, unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Constructor.
Definition: istl/ovlpistlsolverbackend.hh:1114
X::ElementType norm(const X &x) const
Norm of a right-hand side vector. The vector must be consistent on the interior+border partition...
Definition: istl/ovlpistlsolverbackend.hh:348
ISTLBackend_OVLP_GMRES_ILU0(const GFS &gfs_, const CC &cc_, unsigned maxiter_=5000, int verbose_=1, int restart_=20)
make a linear solver object
Definition: istl/ovlpistlsolverbackend.hh:694
Solver to be used for explicit time-steppers with (block-)diagonal mass matrix.
Definition: istl/ovlpistlsolverbackend.hh:859
Definition: istl/ovlpistlsolverbackend.hh:323
void set_constrained_dofs(const CG &cg, typename XG::ElementType x, XG &xg)
construct constraints from given boundary condition function
Definition: constraints.hh:803
Dune::PDELab::Backend::Vector< GFS, typename P::domain_type::field_type > domain_type
The domain type of the preconditioner.
Definition: istl/ovlpistlsolverbackend.hh:140
void apply(M &A, V &z, W &r, typename W::ElementType reduction)
solve the given linear system
Definition: istl/ovlpistlsolverbackend.hh:898
ISTLBackend_BCGS_AMG_SSOR(const GFS &gfs_, unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Constructor.
Definition: istl/ovlpistlsolverbackend.hh:1143
Definition: istl/ovlpistlsolverbackend.hh:375
The category the preconditioner is part of.
Definition: istl/ovlpistlsolverbackend.hh:147
double tsetup
The time needed for building the AMG hierarchy (coarsening).
Definition: istl/seqistlsolverbackend.hh:553
ISTLBackend_OVLP_SuperLU_Base(const GFS &gfs_, const C &c_, unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/ovlpistlsolverbackend.hh:759
Definition: istl/ovlpistlsolverbackend.hh:523
Dune::Amg::Parameters Parameters
Parameters object to customize matrix hierachy building.
Definition: istl/ovlpistlsolverbackend.hh:956
Definition: istl/ovlpistlsolverbackend.hh:134
void setParameters(const Parameters ¶ms_)
set AMG parameters
Definition: istl/ovlpistlsolverbackend.hh:983
const Parameters & parameters() const
Get the parameters describing the behaviuour of AMG.
Definition: istl/ovlpistlsolverbackend.hh:1000
Overlapping parallel BiCGStab solver preconditioned with AMG smoothed by SSOR.
Definition: istl/ovlpistlsolverbackend.hh:1129
Definition: istl/ovlpistlsolverbackend.hh:52
typename impl::BackendVectorSelector< GridFunctionSpace, FieldType >::Type Vector
alias of the return type of BackendVectorSelector
Definition: backend/interface.hh:112
Overlapping parallel BiCGStab solver with ILU0 preconditioner.
Definition: istl/ovlpistlsolverbackend.hh:616
RFType reduction
Definition: solver.hh:35
Definition: parallelhelper.hh:45
double tprepare
The needed for computing the parallel information and for adapting the linear system.
Definition: istl/seqistlsolverbackend.hh:547
ISTLBackend_OVLP_BCGS_ILUn(const GFS &gfs, const CC &cc, int n=1, unsigned maxiter=5000, int verbose=1)
make a linear solver object
Definition: istl/ovlpistlsolverbackend.hh:649
Overlapping parallel CG solver with SuperLU preconditioner.
Definition: istl/ovlpistlsolverbackend.hh:835
X domain_type
export types
Definition: istl/ovlpistlsolverbackend.hh:93
unsigned int iterations
Definition: solver.hh:33
static const unsigned int value
Definition: gridfunctionspace/tags.hh:175
ISTLBackend_OVLP_ILUn_Base(const GFS &gfs_, const C &c_, int n_=1, unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/ovlpistlsolverbackend.hh:535
virtual void post(domain_type &x)
Clean up.
Definition: istl/ovlpistlsolverbackend.hh:180
virtual double norm(const X &x)
Norm of a right-hand side vector. The vector must be consistent on the interior+border partition...
Definition: istl/ovlpistlsolverbackend.hh:122
OVLPScalarProduct(const OVLPScalarProductImplementation< GFS > &implementation_)
Definition: istl/ovlpistlsolverbackend.hh:376
void apply(M &A, V &z, V &r, typename V::ElementType reduction)
solve the given linear system
Definition: istl/ovlpistlsolverbackend.hh:1023
Definition: istl/ovlpistlsolverbackend.hh:88
Y range_type
Definition: istl/ovlpistlsolverbackend.hh:48
Dune::PDELab::Backend::Vector< GFS, typename P::range_type::field_type > range_type
The range type of the preconditioner.
Definition: istl/ovlpistlsolverbackend.hh:142
X::ElementType dot(const X &x, const X &y) const
Dot product of two vectors. It is assumed that the vectors are consistent on the interior+border part...
Definition: istl/ovlpistlsolverbackend.hh:335
Definition: istl/ovlpistlsolverbackend.hh:371
Definition: istl/ovlpistlsolverbackend.hh:97
virtual void apply(const domain_type &x, range_type &y) const
apply operator to x:
Definition: istl/ovlpistlsolverbackend.hh:59
bool directCoarseLevelSolver
True if a direct solver was used on the coarset level.
Definition: istl/seqistlsolverbackend.hh:557
ISTLBackend_OVLP_BCGS_ILU0(const GFS &gfs, const CC &cc, unsigned maxiter=5000, int verbose=1)
make a linear solver object
Definition: istl/ovlpistlsolverbackend.hh:627
std::enable_if< std::is_base_of< impl::WrapperBase, T >::value, Native< T > & >::type native(T &t)
Definition: backend/interface.hh:198
virtual void applyscaleadd(field_type alpha, const domain_type &x, range_type &y) const
apply operator to x, scale and add:
Definition: istl/ovlpistlsolverbackend.hh:67
bool converged
Definition: solver.hh:32
Definition: istl/ovlpistlsolverbackend.hh:461
virtual void pre(domain_type &x, range_type &b)
Prepare the preconditioner.
Definition: istl/ovlpistlsolverbackend.hh:159
Overlapping parallel BiCGStab solver preconditioned with AMG smoothed by ILU0.
Definition: istl/ovlpistlsolverbackend.hh:1158
RFType conv_rate
Definition: solver.hh:36
double elapsed
Definition: solver.hh:34
V::ElementType norm(const V &v) const
compute global norm of a vector
Definition: istl/ovlpistlsolverbackend.hh:1009
Definition: istl/ovlpistlsolverbackend.hh:929
void apply(M &A, V &z, W &r, typename V::ElementType reduction)
solve the given linear system
Definition: istl/ovlpistlsolverbackend.hh:484
X::ElementType field_type
Definition: istl/ovlpistlsolverbackend.hh:94
virtual const M & getmat() const
get matrix via *
Definition: istl/ovlpistlsolverbackend.hh:75
Definition: genericdatahandle.hh:622
OverlappingScalarProduct(const GFS &gfs_, const istl::ParallelHelper< GFS > &helper_)
Constructor needs to know the grid function space.
Definition: istl/ovlpistlsolverbackend.hh:101
Definition: istl/ovlpistlsolverbackend.hh:397
virtual void apply(domain_type &v, const range_type &d)
Apply the preconditioner.
Definition: istl/ovlpistlsolverbackend.hh:167
Definition: adaptivity.hh:27
V::ElementType norm(const V &v) const
compute global norm of a vector
Definition: istl/ovlpistlsolverbackend.hh:880
void apply(M &A, V &z, W &r, typename V::ElementType reduction)
solve the given linear system
Definition: istl/ovlpistlsolverbackend.hh:422
Dune::Amg::SequentialInformation type
Definition: parallelhelper.hh:421
Overlapping parallel CGS solver with SSOR preconditioner.
Definition: istl/ovlpistlsolverbackend.hh:659
void setparams(Parameters params_)
Definition: istl/ovlpistlsolverbackend.hh:988
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
Overlapping parallel BiCGStab solver with SSOR preconditioner.
Definition: istl/ovlpistlsolverbackend.hh:593
void createIndexSetAndProjectForAMG(MatrixType &m, Comm &c)
Makes the matrix consistent and creates the parallel information for AMG.
const istl::ParallelHelper< GFS > & parallelHelper() const
Definition: istl/ovlpistlsolverbackend.hh:353
void apply(M &A, V &z, W &r, typename V::ElementType reduction)
solve the given linear system
Definition: istl/ovlpistlsolverbackend.hh:707
void apply(M &A, V &z, W &r, typename V::ElementType reduction)
solve the given linear system
Definition: istl/ovlpistlsolverbackend.hh:547
virtual field_type dot(const X &x, const X &y)
Dot product of two vectors. It is assumed that the vectors are consistent on the interior+border part...
Definition: istl/ovlpistlsolverbackend.hh:110
Definition: genericdatahandle.hh:685
ISTLBackend_OVLP_CG_SSORk(const GFS &gfs, const CC &cc, unsigned maxiter=5000, int steps=5, int verbose=1)
make a linear solver object
Definition: istl/ovlpistlsolverbackend.hh:671
ISTLBackend_OVLP_Base(const GFS &gfs_, const C &c_, unsigned maxiter_=5000, int steps_=5, int verbose_=1)
make a linear solver object
Definition: istl/ovlpistlsolverbackend.hh:409
ISTLBackend_OVLP_CG_SuperLU(const GFS &gfs_, const CC &cc_, unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/ovlpistlsolverbackend.hh:847
Overlapping parallel restarted GMRes solver with ILU0 preconditioner.
Definition: istl/ovlpistlsolverbackend.hh:683
OVLPScalarProductImplementation(const GFS &gfs_)
Definition: istl/ovlpistlsolverbackend.hh:326
Definition: istl/ovlpistlsolverbackend.hh:41
OverlappingOperator(const CC &cc_, const M &A)
Definition: istl/ovlpistlsolverbackend.hh:54
int levels
the number of levels in the AMG hierarchy.
Definition: istl/seqistlsolverbackend.hh:549
const std::string s
Definition: function.hh:1102
virtual X::BaseT::field_type norm(const X &x)
Definition: istl/ovlpistlsolverbackend.hh:385
Definition: istl/ovlpistlsolverbackend.hh:748
X::ElementType field_type
Definition: istl/ovlpistlsolverbackend.hh:49
Class providing some statistics of the AMG solver.
Definition: istl/seqistlsolverbackend.hh:541
Dune::PDELab::LinearSolverResult< double > res
Definition: solver.hh:52
const ISTLAMGStatistics & statistics() const
Get statistics of the AMG solver (no of levels, timings).
Definition: istl/ovlpistlsolverbackend.hh:1072
ISTLBackend_OVLP_BCGS_SuperLU(const GFS &gfs_, const CC &cc_, unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/ovlpistlsolverbackend.hh:823
ISTLBackend_OVLP_ExplicitDiagonal(const ISTLBackend_OVLP_ExplicitDiagonal &other_)
Definition: istl/ovlpistlsolverbackend.hh:871
ISTLBackend_OVLP_ILU0_Base(const GFS &gfs_, const C &c_, unsigned maxiter_=5000, int verbose_=1)
make a linear solver object
Definition: istl/ovlpistlsolverbackend.hh:472
Overlapping parallel BiCGStab solver with SuperLU preconditioner.
Definition: istl/ovlpistlsolverbackend.hh:811
Overlapping parallel BiCGStab solver with ILU0 preconditioner.
Definition: istl/ovlpistlsolverbackend.hh:637
istl::ParallelHelper< GFS > & parallelHelper()
Definition: istl/ovlpistlsolverbackend.hh:359
double tsolve
The time spent in solving the system (without building the hierarchy.
Definition: istl/seqistlsolverbackend.hh:551
ISTLBackend_OVLP_BCGS_SSORk(const GFS &gfs, const CC &cc, unsigned maxiter=5000, int steps=5, int verbose=1)
make a linear solver object
Definition: istl/ovlpistlsolverbackend.hh:605
ISTLBackend_BCGS_AMG_ILU0(const GFS &gfs_, unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Constructor.
Definition: istl/ovlpistlsolverbackend.hh:1172
OverlappingWrappedPreconditioner(const GFS &gfs_, P &prec_, const CC &cc_, const istl::ParallelHelper< GFS > &helper_)
Constructor.
Definition: istl/ovlpistlsolverbackend.hh:151
X domain_type
Definition: istl/ovlpistlsolverbackend.hh:47
Overlapping parallel conjugate gradient solver preconditioned with AMG smoothed by SSOR...
Definition: istl/ovlpistlsolverbackend.hh:1100
ISTLBackend_OVLP_ExplicitDiagonal(const GFS &gfs_)
make a linear solver object
Definition: istl/ovlpistlsolverbackend.hh:867