MFEM  v4.6.0
Finite element discretization library
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
mfem::DGMassInverse Class Reference

Solver for the discontinuous Galerkin mass matrix. More...

#include <dgmassinv.hpp>

Inheritance diagram for mfem::DGMassInverse:
[legend]
Collaboration diagram for mfem::DGMassInverse:
[legend]

Public Member Functions

 DGMassInverse (FiniteElementSpace &fes_, int btype=BasisType::GaussLegendre)
 Construct the DG inverse mass operator for fes_. More...
 
 DGMassInverse (FiniteElementSpace &fes_, Coefficient &coeff, int btype=BasisType::GaussLegendre)
 Construct the DG inverse mass operator for fes_ with Coefficient coeff. More...
 
 DGMassInverse (FiniteElementSpace &fes_, Coefficient &coeff, const IntegrationRule &ir, int btype=BasisType::GaussLegendre)
 Construct the DG inverse mass operator for fes_ with Coefficient coeff and IntegrationRule ir. More...
 
 DGMassInverse (FiniteElementSpace &fes_, const IntegrationRule &ir, int btype=BasisType::GaussLegendre)
 Construct the DG inverse mass operator for fes_ with IntegrationRule ir. More...
 
void Mult (const Vector &b, Vector &u) const
 Solve the system M b = u. More...
 
void MultTranspose (const Vector &b, Vector &u) const
 Same as Mult() since the mass matrix is symmetric. More...
 
void SetOperator (const Operator &op)
 Not implemented. Aborts. More...
 
void SetRelTol (const double rel_tol_)
 Set the relative tolerance. More...
 
void SetAbsTol (const double abs_tol_)
 Set the absolute tolerance. More...
 
void SetMaxIter (const double max_iter_)
 Set the maximum number of iterations. More...
 
void Update ()
 Recompute operator and preconditioner (when coefficient or mesh changes). More...
 
 ~DGMassInverse ()
 
template<int DIM, int D1D = 0, int Q1D = 0>
void DGMassCGIteration (const Vector &b_, Vector &u_) const
 Solve the system M b = u. Not part of the public interface. More...
 
- Public Member Functions inherited from mfem::Solver
 Solver (int s=0, bool iter_mode=false)
 Initialize a square Solver with size s. More...
 
 Solver (int h, int w, bool iter_mode=false)
 Initialize a Solver with height h and width w. More...
 
- Public Member Functions inherited from mfem::Operator
void InitTVectors (const Operator *Po, const Operator *Ri, const Operator *Pi, Vector &x, Vector &b, Vector &X, Vector &B) const
 Initializes memory for true vectors of linear system. More...
 
 Operator (int s=0)
 Construct a square Operator with given size s (default 0). More...
 
 Operator (int h, int w)
 Construct an Operator with the given height (output size) and width (input size). More...
 
int Height () const
 Get the height (size of output) of the Operator. Synonym with NumRows(). More...
 
int NumRows () const
 Get the number of rows (size of output) of the Operator. Synonym with Height(). More...
 
int Width () const
 Get the width (size of input) of the Operator. Synonym with NumCols(). More...
 
int NumCols () const
 Get the number of columns (size of input) of the Operator. Synonym with Width(). More...
 
virtual MemoryClass GetMemoryClass () const
 Return the MemoryClass preferred by the Operator. More...
 
virtual void AddMult (const Vector &x, Vector &y, const double a=1.0) const
 Operator application: y+=A(x) (default) or y+=a*A(x). More...
 
virtual void AddMultTranspose (const Vector &x, Vector &y, const double a=1.0) const
 Operator transpose application: y+=A^t(x) (default) or y+=a*A^t(x). More...
 
virtual void ArrayMult (const Array< const Vector *> &X, Array< Vector *> &Y) const
 Operator application on a matrix: Y=A(X). More...
 
virtual void ArrayMultTranspose (const Array< const Vector *> &X, Array< Vector *> &Y) const
 Action of the transpose operator on a matrix: Y=A^t(X). More...
 
virtual void ArrayAddMult (const Array< const Vector *> &X, Array< Vector *> &Y, const double a=1.0) const
 Operator application on a matrix: Y+=A(X) (default) or Y+=a*A(X). More...
 
virtual void ArrayAddMultTranspose (const Array< const Vector *> &X, Array< Vector *> &Y, const double a=1.0) const
 Operator transpose application on a matrix: Y+=A^t(X) (default) or Y+=a*A^t(X). More...
 
virtual OperatorGetGradient (const Vector &x) const
 Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate an error. More...
 
virtual void AssembleDiagonal (Vector &diag) const
 Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operators. In some cases, only an approximation of the diagonal is computed. More...
 
virtual const OperatorGetProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity. More...
 
virtual const OperatorGetRestriction () const
 Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More...
 
virtual const OperatorGetOutputProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity. More...
 
virtual const OperatorGetOutputRestrictionTranspose () const
 Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type operators. More...
 
virtual const OperatorGetOutputRestriction () const
 Restriction operator from output vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More...
 
void FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0)
 Form a constrained linear system using a matrix-free approach. More...
 
void FormRectangularLinearSystem (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B)
 Form a column-constrained linear system using a matrix-free approach. More...
 
virtual void RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x)
 Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear system obtained from Operator::FormLinearSystem() or Operator::FormRectangularLinearSystem(). More...
 
void FormSystemOperator (const Array< int > &ess_tdof_list, Operator *&A)
 Return in A a parallel (on truedofs) version of this square operator. More...
 
void FormRectangularSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Operator *&A)
 Return in A a parallel (on truedofs) version of this rectangular operator (including constraints). More...
 
void FormDiscreteOperator (Operator *&A)
 Return in A a parallel (on truedofs) version of this rectangular operator. More...
 
void PrintMatlab (std::ostream &out, int n, int m=0) const
 Prints operator with input size n and output size m in Matlab format. More...
 
virtual void PrintMatlab (std::ostream &out) const
 Prints operator in Matlab format. More...
 
virtual ~Operator ()
 Virtual destructor. More...
 
Type GetType () const
 Return the type ID of the Operator class. More...
 

Protected Member Functions

 DGMassInverse (FiniteElementSpace &fes_, Coefficient *coeff, const IntegrationRule *ir, int btype)
 Protected constructor, used internally. More...
 
- Protected Member Functions inherited from mfem::Operator
void FormConstrainedSystemOperator (const Array< int > &ess_tdof_list, ConstrainedOperator *&Aout)
 see FormSystemOperator() More...
 
void FormRectangularConstrainedSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, RectangularConstrainedOperator *&Aout)
 see FormRectangularSystemOperator() More...
 
OperatorSetupRAP (const Operator *Pi, const Operator *Po)
 Returns RAP Operator of this, using input/output Prolongation matrices Pi corresponds to "P", Po corresponds to "Rt". More...
 

Protected Attributes

DG_FECollection fec
 FE collection in requested basis. More...
 
FiniteElementSpace fes
 FE space in requested basis. More...
 
const DofToQuadd2q
 Change of basis. Not owned. More...
 
Array< double > B_
 Inverse of change of basis. More...
 
Array< double > Bt_
 Inverse of change of basis, transposed. More...
 
class BilinearFormM
 Mass bilinear form, owned. More...
 
class MassIntegratorm
 Mass integrator, owned by the form M. More...
 
Vector diag_inv
 Jacobi preconditioner. More...
 
double rel_tol = 1e-12
 Relative CG tolerance. More...
 
double abs_tol = 1e-12
 Absolute CG tolerance. More...
 
int max_iter = 100
 Maximum number of CG iterations;. More...
 
Intermediate vectors needed for CG three-term recurrence.
Vector r_
 
Vector d_
 
Vector z_
 
Vector b2_
 
- Protected Attributes inherited from mfem::Operator
int height
 Dimension of the output / number of rows in the matrix. More...
 
int width
 Dimension of the input / number of columns in the matrix. More...
 

Additional Inherited Members

- Public Types inherited from mfem::Operator
enum  DiagonalPolicy { DIAG_ZERO, DIAG_ONE, DIAG_KEEP }
 Defines operator diagonal policy upon elimination of rows and/or columns. More...
 
enum  Type {
  ANY_TYPE, MFEM_SPARSEMAT, Hypre_ParCSR, PETSC_MATAIJ,
  PETSC_MATIS, PETSC_MATSHELL, PETSC_MATNEST, PETSC_MATHYPRE,
  PETSC_MATGENERIC, Complex_Operator, MFEM_ComplexSparseMat, Complex_Hypre_ParCSR,
  Complex_DenseMat, MFEM_Block_Matrix, MFEM_Block_Operator
}
 Enumeration defining IDs for some classes derived from Operator. More...
 
- Public Attributes inherited from mfem::Solver
bool iterative_mode
 If true, use the second argument of Mult() as an initial guess. More...
 

Detailed Description

Solver for the discontinuous Galerkin mass matrix.

This class performs a local (diagonally preconditioned) conjugate gradient iteration for each element. Optionally, a change of basis is performed to iterate on a better-conditioned system. This class fully supports execution on device (GPU).

Definition at line 27 of file dgmassinv.hpp.

Constructor & Destructor Documentation

◆ DGMassInverse() [1/5]

mfem::DGMassInverse::DGMassInverse ( FiniteElementSpace fes_,
Coefficient coeff,
const IntegrationRule ir,
int  btype 
)
protected

Protected constructor, used internally.

Custom coefficient and integration rule are used if coeff and ir are non-NULL.

Definition at line 20 of file dgmassinv.cpp.

◆ DGMassInverse() [2/5]

mfem::DGMassInverse::DGMassInverse ( FiniteElementSpace fes_,
int  btype = BasisType::GaussLegendre 
)

Construct the DG inverse mass operator for fes_.

The basis type btype determines which basis should be used internally in the solver. This does not have to be the same basis as fes_. The best choice is typically BasisType::GaussLegendre because it is well-preconditioned by its diagonal.

The solution and right-hand side used for the solver are not affected by this basis (they correspond to the basis of fes_). btype is only used internally, and only has an effect on the convergence rate.

Definition at line 92 of file dgmassinv.cpp.

◆ DGMassInverse() [3/5]

mfem::DGMassInverse::DGMassInverse ( FiniteElementSpace fes_,
Coefficient coeff,
int  btype = BasisType::GaussLegendre 
)

Construct the DG inverse mass operator for fes_ with Coefficient coeff.

See also
DGMassInverse(FiniteElementSpace&, int) for information about btype.

Definition at line 80 of file dgmassinv.cpp.

◆ DGMassInverse() [4/5]

mfem::DGMassInverse::DGMassInverse ( FiniteElementSpace fes_,
Coefficient coeff,
const IntegrationRule ir,
int  btype = BasisType::GaussLegendre 
)

Construct the DG inverse mass operator for fes_ with Coefficient coeff and IntegrationRule ir.

See also
DGMassInverse(FiniteElementSpace&, int) for information about btype.

Definition at line 84 of file dgmassinv.cpp.

◆ DGMassInverse() [5/5]

mfem::DGMassInverse::DGMassInverse ( FiniteElementSpace fes_,
const IntegrationRule ir,
int  btype = BasisType::GaussLegendre 
)

Construct the DG inverse mass operator for fes_ with IntegrationRule ir.

See also
DGMassInverse(FiniteElementSpace&, int) for information about btype.

Definition at line 88 of file dgmassinv.cpp.

◆ ~DGMassInverse()

mfem::DGMassInverse::~DGMassInverse ( )

Definition at line 113 of file dgmassinv.cpp.

Member Function Documentation

◆ DGMassCGIteration()

template<int DIM, int D1D, int Q1D>
void mfem::DGMassInverse::DGMassCGIteration ( const Vector b_,
Vector u_ 
) const

Solve the system M b = u. Not part of the public interface.

Note
This member function must be public because it defines an extended lambda used in an mfem::forall kernel (nvcc limitation)

Definition at line 119 of file dgmassinv.cpp.

◆ Mult()

void mfem::DGMassInverse::Mult ( const Vector b,
Vector u 
) const
virtual

Solve the system M b = u.

If iterative_mode is true, u is used as an initial guess.

Implements mfem::Operator.

Definition at line 263 of file dgmassinv.cpp.

◆ MultTranspose()

void mfem::DGMassInverse::MultTranspose ( const Vector b,
Vector u 
) const
inlinevirtual

Same as Mult() since the mass matrix is symmetric.

Reimplemented from mfem::Operator.

Definition at line 91 of file dgmassinv.hpp.

◆ SetAbsTol()

void mfem::DGMassInverse::SetAbsTol ( const double  abs_tol_)

Set the absolute tolerance.

Definition at line 102 of file dgmassinv.cpp.

◆ SetMaxIter()

void mfem::DGMassInverse::SetMaxIter ( const double  max_iter_)

Set the maximum number of iterations.

Definition at line 104 of file dgmassinv.cpp.

◆ SetOperator()

void mfem::DGMassInverse::SetOperator ( const Operator op)
virtual

Not implemented. Aborts.

Implements mfem::Solver.

Definition at line 95 of file dgmassinv.cpp.

◆ SetRelTol()

void mfem::DGMassInverse::SetRelTol ( const double  rel_tol_)

Set the relative tolerance.

Definition at line 100 of file dgmassinv.cpp.

◆ Update()

void mfem::DGMassInverse::Update ( )

Recompute operator and preconditioner (when coefficient or mesh changes).

Definition at line 106 of file dgmassinv.cpp.

Member Data Documentation

◆ abs_tol

double mfem::DGMassInverse::abs_tol = 1e-12
protected

Absolute CG tolerance.

Definition at line 39 of file dgmassinv.hpp.

◆ b2_

Vector mfem::DGMassInverse::b2_
mutableprotected

Definition at line 44 of file dgmassinv.hpp.

◆ B_

Array<double> mfem::DGMassInverse::B_
protected

Inverse of change of basis.

Definition at line 33 of file dgmassinv.hpp.

◆ Bt_

Array<double> mfem::DGMassInverse::Bt_
protected

Inverse of change of basis, transposed.

Definition at line 34 of file dgmassinv.hpp.

◆ d2q

const DofToQuad* mfem::DGMassInverse::d2q
protected

Change of basis. Not owned.

Definition at line 32 of file dgmassinv.hpp.

◆ d_

Vector mfem::DGMassInverse::d_
mutableprotected

Definition at line 44 of file dgmassinv.hpp.

◆ diag_inv

Vector mfem::DGMassInverse::diag_inv
protected

Jacobi preconditioner.

Definition at line 37 of file dgmassinv.hpp.

◆ fec

DG_FECollection mfem::DGMassInverse::fec
protected

FE collection in requested basis.

Definition at line 30 of file dgmassinv.hpp.

◆ fes

FiniteElementSpace mfem::DGMassInverse::fes
protected

FE space in requested basis.

Definition at line 31 of file dgmassinv.hpp.

◆ M

class BilinearForm* mfem::DGMassInverse::M
protected

Mass bilinear form, owned.

Definition at line 35 of file dgmassinv.hpp.

◆ m

class MassIntegrator* mfem::DGMassInverse::m
protected

Mass integrator, owned by the form M.

Definition at line 36 of file dgmassinv.hpp.

◆ max_iter

int mfem::DGMassInverse::max_iter = 100
protected

Maximum number of CG iterations;.

Definition at line 40 of file dgmassinv.hpp.

◆ r_

Vector mfem::DGMassInverse::r_
mutableprotected

Definition at line 44 of file dgmassinv.hpp.

◆ rel_tol

double mfem::DGMassInverse::rel_tol = 1e-12
protected

Relative CG tolerance.

Definition at line 38 of file dgmassinv.hpp.

◆ z_

Vector mfem::DGMassInverse::z_
mutableprotected

Definition at line 44 of file dgmassinv.hpp.


The documentation for this class was generated from the following files: