DUNE PDELab (git)

Dune::PDELab::ISTLBackend_OVLP_AMG_4_DG< DGGO, DGCC, CGGFS, CGCC, TransferLOP, DGPrec, Solver, s > Class Template Reference

#include <dune/pdelab/backend/istl/ovlp_amg_dg_backend.hh>

Public Member Functions

void setParameters (const Parameters &amg_parameters_)
 set AMG parameters More...
 
const Parametersparameters () const
 Get the parameters describing the behaviuour of AMG. More...
 
void setReuse (bool reuse_)
 Set whether the AMG should be reused again during call to apply().
 
bool getReuse () const
 Return whether the AMG is reused during call to apply()
 
 ISTLBackend_OVLP_AMG_4_DG (DGGO &dggo_, const DGCC &dgcc_, CGGFS &cggfs_, const CGCC &cgcc_, unsigned maxiter_=5000, int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
 
 ISTLBackend_OVLP_AMG_4_DG (DGGO &dggo_, const DGCC &dgcc_, CGGFS &cggfs_, const CGCC &cgcc_, const ParameterTree &params)
 
void setDGSmootherRelaxation (double relaxation_)
 set number of presmoothing steps on the DG level
 
void setNoDGPreSmoothSteps (int n1_)
 set number of presmoothing steps on the DG level
 
void setNoDGPostSmoothSteps (int n2_)
 set number of postsmoothing steps on the DG level
 
void apply (M &A, V &z, V &r, typename Dune::template FieldTraits< typename V::ElementType >::real_type reduction)
 solve the given linear system More...
 
template<typename X >
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 partition.
 
template<typename X >
Dune::template FieldTraits< typenameX::ElementType >::real_type norm (const X &x) const
 Norm of a right-hand side vector. The vector must be consistent on the interior+border partition.
 
const Dune::PDELab::LinearSolverResult< double > & result () const
 Return access to result data.
 

Detailed Description

template<class DGGO, class DGCC, class CGGFS, class CGCC, class TransferLOP, template< class, class, class, int > class DGPrec, template< class > class Solver, int s = 96>
class Dune::PDELab::ISTLBackend_OVLP_AMG_4_DG< DGGO, DGCC, CGGFS, CGCC, TransferLOP, DGPrec, Solver, s >

Overlapping solver backend for using AMG for DG in PDELab

The template parameters are: DGGO GridOperator for DG discretization, allows access to matrix, vector and grid function space DGCC constraints container for DG problem CGGFS grid function space for CG subspace CGCC constraints container for CG problem TransferLOP local operator to assemble prolongation from CGGFS to DGGFS DGPrec preconditioner for DG problem Solver solver to be used on the complete problem int s size of global index to be used in AMG

Note: The subspace matrix is calculated by a triple matrix product with the fine space DG matrix passed into the apply method. This only works if this DG matrix does not contain P0ParallelConstraints. In order to achieve this behaviour you should use this solver backend with a StationaryLinearProblem or Newton solver that was created using a grid operator with empty constranints. See the overlapping test case in

dune-pdelab/dune/pdelab/test-dg-amg.cc

for an example.

Constructor & Destructor Documentation

◆ ISTLBackend_OVLP_AMG_4_DG() [1/2]

template<class DGGO , class DGCC , class CGGFS , class CGCC , class TransferLOP , template< class, class, class, int > class DGPrec, template< class > class Solver, int s = 96>
Dune::PDELab::ISTLBackend_OVLP_AMG_4_DG< DGGO, DGCC, CGGFS, CGCC, TransferLOP, DGPrec, Solver, s >::ISTLBackend_OVLP_AMG_4_DG ( DGGO &  dggo_,
const DGCC &  dgcc_,
CGGFS &  cggfs_,
const CGCC &  cgcc_,
unsigned  maxiter_ = 5000,
int  verbose_ = 1,
bool  reuse_ = false,
bool  usesuperlu_ = true 
)
inline

◆ ISTLBackend_OVLP_AMG_4_DG() [2/2]

template<class DGGO , class DGCC , class CGGFS , class CGCC , class TransferLOP , template< class, class, class, int > class DGPrec, template< class > class Solver, int s = 96>
Dune::PDELab::ISTLBackend_OVLP_AMG_4_DG< DGGO, DGCC, CGGFS, CGCC, TransferLOP, DGPrec, Solver, s >::ISTLBackend_OVLP_AMG_4_DG ( DGGO &  dggo_,
const DGCC &  dgcc_,
CGGFS &  cggfs_,
const CGCC &  cgcc_,
const ParameterTree params 
)
inline

Member Function Documentation

◆ apply()

template<class DGGO , class DGCC , class CGGFS , class CGCC , class TransferLOP , template< class, class, class, int > class DGPrec, template< class > class Solver, int s = 96>
void Dune::PDELab::ISTLBackend_OVLP_AMG_4_DG< DGGO, DGCC, CGGFS, CGCC, TransferLOP, DGPrec, Solver, s >::apply ( M &  A,
V &  z,
V &  r,
typename Dune::template FieldTraits< typename V::ElementType >::real_type  reduction 
)
inline

solve the given linear system

Parameters
[in]Athe given matrix
[out]zthe solution vector to be computed
[in]rright hand side
[in]reductionto be achieved

◆ parameters()

template<class DGGO , class DGCC , class CGGFS , class CGCC , class TransferLOP , template< class, class, class, int > class DGPrec, template< class > class Solver, int s = 96>
const Parameters & Dune::PDELab::ISTLBackend_OVLP_AMG_4_DG< DGGO, DGCC, CGGFS, CGCC, TransferLOP, DGPrec, Solver, s >::parameters ( ) const
inline

Get the parameters describing the behaviuour of AMG.

The returned object can be adjusted to ones needs and then can be reset using setParameters.

Returns
The object holding the parameters of AMG.

◆ setParameters()

template<class DGGO , class DGCC , class CGGFS , class CGCC , class TransferLOP , template< class, class, class, int > class DGPrec, template< class > class Solver, int s = 96>
void Dune::PDELab::ISTLBackend_OVLP_AMG_4_DG< DGGO, DGCC, CGGFS, CGCC, TransferLOP, DGPrec, Solver, s >::setParameters ( const Parameters amg_parameters_)
inline

set AMG parameters

Parameters
[in]amg_parameters_a parameter object of Type Dune::Amg::Parameters

The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 23, 23:29, 2024)