1 #ifndef DUNE_PDELAB_BACKEND_ISTL_SEQ_AMG_DG_BACKEND_HH 2 #define DUNE_PDELAB_BACKEND_ISTL_SEQ_AMG_DG_BACKEND_HH 7 #include <dune/common/power.hh> 8 #include <dune/common/parametertree.hh> 10 #include <dune/istl/matrixmatrix.hh> 12 #include <dune/grid/common/datahandleif.hh> 31 template<
class DGMatrix,
class DGPrec,
class CGPrec,
class P>
32 class SeqDGAMGPrec :
public Dune::Preconditioner<typename DGPrec::domain_type,typename DGPrec::range_type>
42 typedef typename DGPrec::domain_type
X;
43 typedef typename DGPrec::range_type
Y;
44 typedef typename CGPrec::domain_type
CGX;
45 typedef typename CGPrec::range_type
CGY;
60 SeqDGAMGPrec (DGMatrix& dgmatrix_, DGPrec& dgprec_, CGPrec& cgprec_, P& p_,
int n1_,
int n2_)
61 : dgmatrix(dgmatrix_), dgprec(dgprec_), cgprec(cgprec_), p(p_), n1(n1_), n2(n2_),
71 virtual void pre (X& x, Y& b)
86 virtual void apply (X& x,
const Y& b)
93 for (
int i=0; i<n1; i++)
108 cgprec.apply(cgv,cgd);
116 for (
int i=0; i<n2; i++)
149 template<
class DGGO,
class CGGFS,
class TransferLOP,
template<
class,
class,
class,
int>
class DGPrec,
template<
class>
class Solver>
153 typedef typename DGGO::Traits::TrialGridFunctionSpace GFS;
156 typedef typename DGGO::Traits::Jacobian M;
157 typedef typename DGGO::Traits::Domain V;
160 typedef typename Vector::field_type field_type;
169 typedef TransferLOP CGTODGLOP;
171 typedef typename PGO::Jacobian PMatrix;
175 typedef typename Dune::TransposedMatMultMatResult<P,Matrix>::type PTADG;
176 typedef typename Dune::MatMultMatResult<PTADG,P>::type ACG;
177 typedef ACG CGMatrix;
180 typedef Dune::MatrixAdapter<CGMatrix,CGVector,CGVector> CGOperator;
181 typedef Dune::SeqSSOR<CGMatrix,CGVector,CGVector,1> Smoother;
182 typedef Dune::Amg::AMG<CGOperator,CGVector,Smoother> AMG;
183 typedef Dune::Amg::Parameters Parameters;
187 std::shared_ptr<CGOperator> cgop;
188 std::shared_ptr<AMG> amg;
189 Parameters amg_parameters;
195 std::size_t low_order_space_entries_per_row;
203 ISTLBackend_SEQ_AMG_4_DG(DGGO& dggo_, CGGFS& cggfs_,
unsigned maxiter_=5000,
int verbose_=1,
bool reuse_=
false,
bool usesuperlu_=
true)
206 , amg_parameters(15,2000)
211 , usesuperlu(usesuperlu_)
212 , low_order_space_entries_per_row(StaticPower<3,GFS::Traits::GridView::dimension>::power)
214 , pgo(cggfs,dggo.trialGridFunctionSpace(),cgtodglop,MBE(low_order_space_entries_per_row))
218 amg_parameters.setDefaultValuesIsotropic(GFS::Traits::GridViewType::Traits::Grid::dimension);
219 amg_parameters.setDebugLevel(verbose_);
221 if (usesuperlu ==
true)
223 std::cout <<
"WARNING: You are using AMG without SuperLU!" 224 <<
" Please consider installing SuperLU," 225 <<
" or set the usesuperlu flag to false" 226 <<
" to suppress this warning." << std::endl;
232 if (verbose>0) std::cout <<
"allocated prolongation matrix of size " << pmatrix.N() <<
" x " << pmatrix.M() << std::endl;
240 , amg_parameters(15,2000)
241 , maxiter(params.get<int>(
"max_iterations",5000))
242 , verbose(params.get<int>(
"verbose",1))
243 , usesuperlu(params.get<bool>(
"use_superlu",true))
244 , low_order_space_entries_per_row(params.get<
std::size_t>(
"low_order_space.entries_per_row",StaticPower<3,GFS::Traits::GridView::dimension>::power))
246 , pgo(cggfs,dggo.trialGridFunctionSpace(),cgtodglop,MBE(low_order_space_entries_per_row))
249 amg_parameters.setDefaultValuesIsotropic(GFS::Traits::GridViewType::Traits::Grid::dimension);
250 amg_parameters.setDebugLevel(params.get<
int>(
"verbose",1));
252 if (usesuperlu ==
true)
254 std::cout <<
"WARNING: You are using AMG without SuperLU!" 255 <<
" Please consider installing SuperLU," 256 <<
" or set the usesuperlu flag to false" 257 <<
" to suppress this warning." << std::endl;
264 if (verbose>0) std::cout <<
"allocated prolongation matrix of size " << pmatrix.N() <<
" x " << pmatrix.M() << std::endl;
273 typename V::ElementType
norm (
const V& v)
const 284 amg_parameters = amg_parameters_;
296 return amg_parameters;
318 void apply (M& A, V& z, V& r,
typename Dune::template FieldTraits<typename V::ElementType >::real_type reduction)
325 double triple_product_time = 0.0;
327 if(reuse ==
false || firstapply ==
true) {
329 Dune::transposeMatMultMat(ptadg,
native(pmatrix),
native(A));
330 Dune::matMultMat(acg,ptadg,
native(pmatrix));
331 triple_product_time = watch.elapsed();
333 std::cout <<
"=== triple matrix product " << triple_product_time <<
" s" << std::endl;
337 std::cout <<
"=== reuse CG matrix, SKIPPING triple matrix product " << std::endl;
340 typedef typename Dune::Amg::SmootherTraits<Smoother>::Arguments SmootherArgs;
341 SmootherArgs smootherArgs;
342 smootherArgs.iterations = 1;
343 smootherArgs.relaxationFactor = 1.0;
344 typedef Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<CGMatrix,Dune::Amg::FirstDiagonal> > Criterion;
345 Criterion criterion(amg_parameters);
349 double amg_setup_time = 0.0;
350 if(reuse ==
false || firstapply ==
true) {
351 cgop.reset(
new CGOperator(acg));
352 amg.reset(
new AMG(*cgop,criterion,smootherArgs));
354 amg_setup_time = watch.elapsed();
356 std::cout <<
"=== AMG setup " <<amg_setup_time <<
" s" << std::endl;
359 std::cout <<
"=== reuse CG matrix, SKIPPING AMG setup " << std::endl;
362 Dune::MatrixAdapter<Matrix,Vector,Vector> op(
native(A));
363 DGPrec<Matrix,Vector,Vector,1> dgprec(
native(A),1,1);
365 HybridPrec hybridprec(
native(A),dgprec,*amg,
native(pmatrix),2,2);
368 Solver<Vector> solver(op,hybridprec,reduction,maxiter,verbose);
371 Dune::InverseOperatorResult stat;
374 double amg_solve_time = watch.elapsed();
375 if (verbose>0) std::cout <<
"=== Hybrid total solve time " << amg_solve_time+amg_setup_time+triple_product_time <<
" s" << std::endl;
376 res.converged = stat.converged;
377 res.iterations = stat.iterations;
378 res.elapsed = amg_solve_time+amg_setup_time+triple_product_time;
379 res.reduction = stat.reduction;
380 res.conv_rate = stat.conv_rate;
386 #endif // DUNE_PDELAB_BACKEND_ISTL_SEQ_AMG_DG_BACKEND_HH ISTLBackend_SEQ_AMG_4_DG(DGGO &dggo_, CGGFS &cggfs_, unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
Definition: seq_amg_dg_backend.hh:203
typename impl::BackendMatrixSelector< Backend, VU, VV, E >::Type Matrix
alias of the return type of BackendMatrixSelector
Definition: backend/interface.hh:127
std::enable_if< std::is_base_of< impl::WrapperBase, T >::value, Native< T > &>::type native(T &t)
Definition: backend/interface.hh:192
Definition: seq_amg_dg_backend.hh:150
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: seq_amg_dg_backend.hh:71
SeqDGAMGPrec(DGMatrix &dgmatrix_, DGPrec &dgprec_, CGPrec &cgprec_, P &p_, int n1_, int n2_)
Constructor.
Definition: seq_amg_dg_backend.hh:60
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
bool getReuse() const
Return whether the AMG is reused during call to apply()
Definition: seq_amg_dg_backend.hh:306
typename impl::BackendVectorSelector< GridFunctionSpace, FieldType >::Type Vector
alias of the return type of BackendVectorSelector
Definition: backend/interface.hh:106
DGPrec::domain_type X
Definition: seq_amg_dg_backend.hh:42
Definition: constraintstransformation.hh:111
For backward compatibility – Do not use this!
Definition: adaptivity.hh:27
void apply(M &A, V &z, V &r, typename Dune::template FieldTraits< typename V::ElementType >::real_type reduction)
solve the given linear system
Definition: seq_amg_dg_backend.hh:318
void setParameters(const Parameters &amg_parameters_)
set AMG parameters
Definition: seq_amg_dg_backend.hh:282
virtual void apply(X &x, const Y &b)
Apply the precondioner.
Definition: seq_amg_dg_backend.hh:86
virtual void post(X &x)
Clean up.
Definition: seq_amg_dg_backend.hh:130
CGPrec::range_type CGY
Definition: seq_amg_dg_backend.hh:45
V::ElementType norm(const V &v) const
compute global norm of a vector
Definition: seq_amg_dg_backend.hh:273
const Parameters & parameters() const
Get the parameters describing the behaviuour of AMG.
Definition: seq_amg_dg_backend.hh:294
Backend using (possibly nested) ISTL BCRSMatrices.
Definition: bcrsmatrixbackend.hh:190
Definition: seq_amg_dg_backend.hh:32
void setReuse(bool reuse_)
Set whether the AMG should be reused again during call to apply().
Definition: seq_amg_dg_backend.hh:300
CGPrec::domain_type CGX
Definition: seq_amg_dg_backend.hh:44
DGPrec::range_type Y
Definition: seq_amg_dg_backend.hh:43
The category the preconditioner is part of.
Definition: seq_amg_dg_backend.hh:50
void jacobian(const Domain &x, Jacobian &a) const
Assembler jacobian.
Definition: gridoperator.hh:183
ISTLBackend_SEQ_AMG_4_DG(DGGO &dggo_, CGGFS &cggfs_, const ParameterTree ¶ms)
Definition: seq_amg_dg_backend.hh:237