DUNE PDELab (2.8)
seq_amg_dg_backend.hh
29 class SeqDGAMGPrec : public Dune::Preconditioner<typename DGPrec::domain_type,typename DGPrec::range_type>
145 template<class DGGO, class CGGFS, class TransferLOP, template<class,class,class,int> class DGPrec, template<class> class Solver>
166 typedef Dune::PDELab::GridOperator<CGGFS,GFS,CGTODGLOP,MBE,field_type,field_type,field_type,CC,CC> PGO;
204 ISTLBackend_SEQ_AMG_4_DG(DGGO& dggo_, CGGFS& cggfs_, unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
233 if (verbose>0) std::cout << "allocated prolongation matrix of size " << pmatrix.N() << " x " << pmatrix.M() << std::endl;
238 ISTLBackend_SEQ_AMG_4_DG(DGGO& dggo_, CGGFS& cggfs_, const ParameterTree& params)//unsigned maxiter_=5000, int verbose_=1, bool usesuperlu_=true)
247 , low_order_space_entries_per_row(params.get<std::size_t>("low_order_space.entries_per_row",StaticPower<3,GFS::Traits::GridView::dimension>::power))
267 if (verbose>0) std::cout << "allocated prolongation matrix of size " << pmatrix.N() << " x " << pmatrix.M() << std::endl;
339 void apply (M& A, V& z, V& r, typename Dune::template FieldTraits<typename V::ElementType >::real_type reduction)
365 typedef Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<CGMatrix,Dune::Amg::FirstDiagonal> > Criterion;
396 if (verbose>0) std::cout << "=== Hybrid total solve time " << amg_solve_time+amg_setup_time+triple_product_time << " s" << std::endl;
The criterion describing the stop criteria for the coarsening process.
Definition: matrixhierarchy.hh:281
void jacobian(const Domain &x, Jacobian &a) const
Assembler jacobian.
Definition: gridoperator.hh:180
Dune::PDELab::Backend::Matrix< MBE, Domain, Range, field_type > Jacobian
The type of the jacobian.
Definition: gridoperator.hh:47
Definition: seq_amg_dg_backend.hh:147
void setNoDGPreSmoothSteps(int n1_)
set number of presmoothing steps on the DG level
Definition: seq_amg_dg_backend.hh:321
void setNoDGPostSmoothSteps(int n2_)
set number of postsmoothing steps on the DG level
Definition: seq_amg_dg_backend.hh:327
void setParameters(const Parameters &amg_parameters_)
set AMG parameters
Definition: seq_amg_dg_backend.hh:285
void setDGSmootherRelaxation(double relaxation_)
set number of presmoothing steps on the DG level
Definition: seq_amg_dg_backend.hh:315
void setReuse(bool reuse_)
Set whether the AMG should be reused again during call to apply().
Definition: seq_amg_dg_backend.hh:303
const Parameters & parameters() const
Get the parameters describing the behaviuour of AMG.
Definition: seq_amg_dg_backend.hh:297
bool getReuse() const
Return whether the AMG is reused during call to apply()
Definition: seq_amg_dg_backend.hh:309
V::ElementType norm(const V &v) const
compute global norm of a vector
Definition: seq_amg_dg_backend.hh:276
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:339
Definition: seq_amg_dg_backend.hh:30
SeqDGAMGPrec(DGMatrix &dgmatrix_, DGPrec &dgprec_, CGPrec &cgprec_, P &p_, int n1_, int n2_)
Constructor.
Definition: seq_amg_dg_backend.hh:56
virtual void pre(X &x, Y &b) override
Prepare the preconditioner.
Definition: seq_amg_dg_backend.hh:67
virtual void apply(X &x, const Y &b) override
Apply the precondioner.
Definition: seq_amg_dg_backend.hh:82
virtual void post(X &x) override
Clean up.
Definition: seq_amg_dg_backend.hh:126
SolverCategory::Category category() const override
Category of the preconditioner (see SolverCategory::Category)
Definition: seq_amg_dg_backend.hh:44
std::string get(const std::string &key, const std::string &defaultValue) const
get value as string
Definition: parametertree.cc:183
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:30
void reset() noexcept
Reset timer while keeping the running/stopped state.
Definition: timer.hh:55
double elapsed() const noexcept
Get elapsed user-time from last reset until now/last stop in seconds.
Definition: timer.hh:75
Describes the parallel communication interface class for MessageBuffers and DataHandles.
Various implementations of the power function for run-time and static arguments.
void setDebugLevel(int level)
Set the debugging level.
Definition: parameters.hh:399
void setDefaultValuesIsotropic(std::size_t dim, std::size_t diameter=2)
Sets reasonable default values for an isotropic problem.
Definition: parameters.hh:107
int iterations
The numbe of iterations to perform.
Definition: smoother.hh:45
void transposeMatMultMat(BCRSMatrix< FieldMatrix< T, n, m >, A > &res, const BCRSMatrix< FieldMatrix< T, k, n >, A1 > &mat, const BCRSMatrix< FieldMatrix< T, k, m >, A2 > &matt, bool tryHard=false)
Calculate product of a transposed sparse matrix with another sparse matrices ( ).
Definition: matrixmatrix.hh:588
void matMultMat(BCRSMatrix< FieldMatrix< T, n, m >, A > &res, const BCRSMatrix< FieldMatrix< T, n, k >, A1 > &mat, const BCRSMatrix< FieldMatrix< T, k, m >, A2 > &matt, bool tryHard=false)
Calculate product of two sparse matrices ( ).
Definition: matrixmatrix.hh:573
provides functions for sparse matrix matrix multiplication.
A hierarchical structure of string parameters.
The default class for the smoother arguments.
Definition: smoother.hh:36
Statistics about the application of an inverse operator.
Definition: solver.hh:46
double conv_rate
Convergence rate (average reduction per step)
Definition: solver.hh:74
bool converged
True if convergence criterion has been met.
Definition: solver.hh:71
Helper TMP to get the result type of a sparse matrix matrix multiplication ( )
Definition: matrixmatrix.hh:508
Backend using (possibly nested) ISTL BCRSMatrices.
Definition: bcrsmatrixbackend.hh:188
@ sequential
Category for sequential solvers.
Definition: solvercategory.hh:23
Helper TMP to get the result type of a sparse matrix matrix multiplication ( )
Definition: matrixmatrix.hh:533
|
Legal Statements / Impressum |
Hosted by TU Dresden |
generated with Hugo v0.111.3
(Dec 21, 23:30, 2024)