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

Class representing the DPG weak formulation. Given the variational formulation a(u,v) = b(v), (or A u = b, where <Au,v> = a(u,v)) this class forms the DPG linear system A^T G^-1 A u = A^T G^-1 b This system results from the minimum residual formulation u = argmin_w ||G^-1(b - Aw)||. Here G is a symmetric positive definite matrix resulting from the discretization of the Riesz operator on the test space. Since the test space is broken (discontinuous), G is defined and inverted element-wise and the assembly of the global system is performed in the same manner as the standard FEM method. Note that DPGWeakForm can handle multiple Finite Element spaces. More...

#include <weakform.hpp>

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

Public Member Functions

 DPGWeakForm ()
 Default constructor. User must call SetSpaces to setup the FE spaces. More...
 
 DPGWeakForm (Array< FiniteElementSpace * > &fes_, Array< FiniteElementCollection *> &fecol_)
 Creates bilinear form associated with FE spaces fes_. More...
 
void SetTestFECollVdim (int test_fec, int vdim)
 
void SetSpaces (Array< FiniteElementSpace * > &fes_, Array< FiniteElementCollection *> &fecol_)
 
int Size () const
 Get the size of the bilinear form of the DPGWeakForm. More...
 
void AllocateMatrix ()
 Pre-allocate the internal BlockMatrix before assembly. More...
 
void Finalize (int skip_zeros=1)
 Finalizes the matrix initialization. More...
 
BlockMatrixBlockMat ()
 Returns a reference to the BlockMatrix: \( M \). More...
 
BlockMatrixBlockMatElim ()
 Returns a reference to the sparse matrix of eliminated b.c.: \( M_e \). More...
 
void AddTrialIntegrator (BilinearFormIntegrator *bfi, int n, int m)
 Adds new Trial Integrator. Assumes ownership of bfi. More...
 
void AddTestIntegrator (BilinearFormIntegrator *bfi, int n, int m)
 Adds new Test Integrator. Assumes ownership of bfi. More...
 
void AddDomainLFIntegrator (LinearFormIntegrator *lfi, int n)
 Adds new Domain LF Integrator. Assumes ownership of bfi. More...
 
void Assemble (int skip_zeros=1)
 Assembles the form i.e. sums over all integrators. More...
 
virtual void FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
 Form the linear system A X = B, corresponding to this DPG weak form. More...
 
template<typename OpType >
void FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, OpType &A, Vector &X, Vector &B, int copy_interior=0)
 Form the linear system A X = B, corresponding to this DPG weak form Version of the method FormLinearSystem() where the system matrix is returned in the variable A, of type OpType, holding a reference to the system matrix (created with the method OpType::MakeRef()). More...
 
virtual void FormSystemMatrix (const Array< int > &ess_tdof_list, OperatorHandle &A)
 Form the linear system matrix A, see FormLinearSystem() for details. More...
 
template<typename OpType >
void FormSystemMatrix (const Array< int > &ess_tdof_list, OpType &A)
 Form the linear system matrix A, see FormLinearSystem() for details. More...
 
void EliminateVDofs (const Array< int > &vdofs, Operator::DiagonalPolicy dpolicy=Operator::DIAG_ONE)
 Eliminate the given vdofs, storing the eliminated part internally in \( M_e \). More...
 
void EliminateVDofsInRHS (const Array< int > &vdofs, const Vector &x, Vector &b)
 Use the stored eliminated part of the matrix (see EliminateVDofs(const Array<int> &, DiagonalPolicy)) to modify the r.h.s. b; vdofs is a list of DOFs. More...
 
virtual void RecoverFEMSolution (const Vector &X, Vector &x)
 Recover the solution of a linear system formed with FormLinearSystem(). More...
 
void SetDiagonalPolicy (Operator::DiagonalPolicy policy)
 Sets diagonal policy used upon construction of the linear system. More...
 
virtual void Update ()
 Update the DPGWeakForm after mesh modifications (AMR) More...
 
void StoreMatrices (bool store_matrices_=true)
 Store internal element matrices used for computation of residual after solve. More...
 
void EnableStaticCondensation ()
 
VectorComputeResidual (const BlockVector &x)
 Compute DPG residual based error estimator. More...
 
virtual ~DPGWeakForm ()
 

Protected Member Functions

void Init ()
 
void ReleaseInitMemory ()
 
void AllocMat ()
 Allocate appropriate BlockMatrix and assign it to mat. More...
 
void ConformingAssemble ()
 
void ComputeOffsets ()
 
virtual void BuildProlongation ()
 

Protected Attributes

BlockStaticCondensationstatic_cond
 Owned. More...
 
bool initialized = false
 
Meshmesh = nullptr
 
int height
 
int width
 
int nblocks
 
Array< int > dof_offsets
 
Array< int > tdof_offsets
 
BlockMatrixmat = nullptr
 Block matrix \( M \) to be associated with the Block bilinear form. Owned. More...
 
BlockVectory = nullptr
 Block vector \( y \) to be associated with the Block linear form. More...
 
BlockMatrixmat_e = nullptr
 Block Matrix \( M_e \) used to store the eliminations from the b.c. Owned. \( M + M_e = M_{original} \). More...
 
Array< FiniteElementSpace *> trial_fes
 Trial FE spaces. More...
 
Array< int > IsTraceFes
 Flags to determine if a FiniteElementSpace is Trace. More...
 
Array< FiniteElementCollection * > test_fecols
 Test FE Collections (Broken) More...
 
Array< int > test_fecols_vdims
 
Array2D< Array< BilinearFormIntegrator *> *> trial_integs
 Set of Trial Integrators to be applied for matrix B. More...
 
Array2D< Array< BilinearFormIntegrator *> *> test_integs
 Set of Test Space (broken) Integrators to be applied for matrix G. More...
 
Array< Array< LinearFormIntegrator *> *> lfis
 Set of Linear Form Integrators to be applied. More...
 
BlockMatrixP = nullptr
 Block Prolongation. More...
 
BlockMatrixR = nullptr
 Block Restriction. More...
 
mfem::Operator::DiagonalPolicy diag_policy
 
bool store_matrices = false
 
Array< DenseMatrix *> Bmat
 
Array< Vector *> fvec
 
Vector residuals
 

Detailed Description

Class representing the DPG weak formulation. Given the variational formulation a(u,v) = b(v), (or A u = b, where <Au,v> = a(u,v)) this class forms the DPG linear system A^T G^-1 A u = A^T G^-1 b This system results from the minimum residual formulation u = argmin_w ||G^-1(b - Aw)||. Here G is a symmetric positive definite matrix resulting from the discretization of the Riesz operator on the test space. Since the test space is broken (discontinuous), G is defined and inverted element-wise and the assembly of the global system is performed in the same manner as the standard FEM method. Note that DPGWeakForm can handle multiple Finite Element spaces.

Definition at line 33 of file weakform.hpp.

Constructor & Destructor Documentation

◆ DPGWeakForm() [1/2]

mfem::DPGWeakForm::DPGWeakForm ( )
inline

Default constructor. User must call SetSpaces to setup the FE spaces.

Definition at line 107 of file weakform.hpp.

◆ DPGWeakForm() [2/2]

mfem::DPGWeakForm::DPGWeakForm ( Array< FiniteElementSpace * > &  fes_,
Array< FiniteElementCollection *> &  fecol_ 
)
inline

Creates bilinear form associated with FE spaces fes_.

Definition at line 114 of file weakform.hpp.

◆ ~DPGWeakForm()

mfem::DPGWeakForm::~DPGWeakForm ( )
virtual

Definition at line 811 of file weakform.cpp.

Member Function Documentation

◆ AddDomainLFIntegrator()

void mfem::DPGWeakForm::AddDomainLFIntegrator ( LinearFormIntegrator lfi,
int  n 
)

Adds new Domain LF Integrator. Assumes ownership of bfi.

Definition at line 125 of file weakform.cpp.

◆ AddTestIntegrator()

void mfem::DPGWeakForm::AddTestIntegrator ( BilinearFormIntegrator bfi,
int  n,
int  m 
)

Adds new Test Integrator. Assumes ownership of bfi.

Adds new Domain BF Integrator. Assumes ownership of bfi.

Definition at line 117 of file weakform.cpp.

◆ AddTrialIntegrator()

void mfem::DPGWeakForm::AddTrialIntegrator ( BilinearFormIntegrator bfi,
int  n,
int  m 
)

Adds new Trial Integrator. Assumes ownership of bfi.

Adds new Domain BF Integrator. Assumes ownership of bfi.

Definition at line 105 of file weakform.cpp.

◆ AllocateMatrix()

void mfem::DPGWeakForm::AllocateMatrix ( )
inline

Pre-allocate the internal BlockMatrix before assembly.

Definition at line 150 of file weakform.hpp.

◆ AllocMat()

void mfem::DPGWeakForm::AllocMat ( )
protected

Allocate appropriate BlockMatrix and assign it to mat.

Definition at line 77 of file weakform.cpp.

◆ Assemble()

void mfem::DPGWeakForm::Assemble ( int  skip_zeros = 1)

Assembles the form i.e. sums over all integrators.

Assembles the form i.e. sums over all domain integrators.

Definition at line 247 of file weakform.cpp.

◆ BlockMat()

BlockMatrix& mfem::DPGWeakForm::BlockMat ( )
inline

Returns a reference to the BlockMatrix: \( M \).

Definition at line 156 of file weakform.hpp.

◆ BlockMatElim()

BlockMatrix& mfem::DPGWeakForm::BlockMatElim ( )
inline

Returns a reference to the sparse matrix of eliminated b.c.: \( M_e \).

Definition at line 163 of file weakform.hpp.

◆ BuildProlongation()

void mfem::DPGWeakForm::BuildProlongation ( )
protectedvirtual

Reimplemented in mfem::ParDPGWeakForm.

Definition at line 133 of file weakform.cpp.

◆ ComputeOffsets()

void mfem::DPGWeakForm::ComputeOffsets ( )
protected

Definition at line 61 of file weakform.cpp.

◆ ComputeResidual()

Vector & mfem::DPGWeakForm::ComputeResidual ( const BlockVector x)

Compute DPG residual based error estimator.

Definition at line 723 of file weakform.cpp.

◆ ConformingAssemble()

void mfem::DPGWeakForm::ConformingAssemble ( )
protected

Definition at line 151 of file weakform.cpp.

◆ EliminateVDofs()

void mfem::DPGWeakForm::EliminateVDofs ( const Array< int > &  vdofs,
Operator::DiagonalPolicy  dpolicy = Operator::DIAG_ONE 
)

Eliminate the given vdofs, storing the eliminated part internally in \( M_e \).

This method works in conjunction with EliminateVDofsInRHS() and allows elimination of boundary conditions in multiple right-hand sides. In this method, vdofs is a list of DOFs.

Definition at line 574 of file weakform.cpp.

◆ EliminateVDofsInRHS()

void mfem::DPGWeakForm::EliminateVDofsInRHS ( const Array< int > &  vdofs,
const Vector x,
Vector b 
)

Use the stored eliminated part of the matrix (see EliminateVDofs(const Array<int> &, DiagonalPolicy)) to modify the r.h.s. b; vdofs is a list of DOFs.

Definition at line 567 of file weakform.cpp.

◆ EnableStaticCondensation()

void mfem::DPGWeakForm::EnableStaticCondensation ( )

Definition at line 717 of file weakform.cpp.

◆ Finalize()

void mfem::DPGWeakForm::Finalize ( int  skip_zeros = 1)

Finalizes the matrix initialization.

Definition at line 97 of file weakform.cpp.

◆ FormLinearSystem() [1/2]

void mfem::DPGWeakForm::FormLinearSystem ( const Array< int > &  ess_tdof_list,
Vector x,
OperatorHandle A,
Vector X,
Vector B,
int  copy_interior = 0 
)
virtual

Form the linear system A X = B, corresponding to this DPG weak form.

This method applies any necessary transformations to the linear system such as: eliminating boundary conditions; applying conforming constraints for non-conforming AMR; static condensation;

The GridFunction-size vector x must contain the essential b.c. The DPGWeakForm must be assembled.

The vector X is initialized with a suitable initial guess: the essential entries of X are set to the corresponding b.c. and all other entries are set to zero (copy_interior == 0) or copied from x (copy_interior != 0).

After solving the linear system, the finite element solution x can be recovered by calling RecoverFEMSolution() (with the same vectors X, and x).

NOTE: If there are no transformations, X simply reuses the data of x.

Reimplemented in mfem::ParDPGWeakForm.

Definition at line 476 of file weakform.cpp.

◆ FormLinearSystem() [2/2]

template<typename OpType >
void mfem::DPGWeakForm::FormLinearSystem ( const Array< int > &  ess_tdof_list,
Vector x,
OpType &  A,
Vector X,
Vector B,
int  copy_interior = 0 
)
inline

Form the linear system A X = B, corresponding to this DPG weak form Version of the method FormLinearSystem() where the system matrix is returned in the variable A, of type OpType, holding a reference to the system matrix (created with the method OpType::MakeRef()).

Definition at line 210 of file weakform.hpp.

◆ FormSystemMatrix() [1/2]

void mfem::DPGWeakForm::FormSystemMatrix ( const Array< int > &  ess_tdof_list,
OperatorHandle A 
)
virtual

Form the linear system matrix A, see FormLinearSystem() for details.

Reimplemented in mfem::ParDPGWeakForm.

Definition at line 531 of file weakform.cpp.

◆ FormSystemMatrix() [2/2]

template<typename OpType >
void mfem::DPGWeakForm::FormSystemMatrix ( const Array< int > &  ess_tdof_list,
OpType &  A 
)
inline

Form the linear system matrix A, see FormLinearSystem() for details.

Version of the method FormSystemMatrix() where the system matrix is returned in the variable A, of type OpType, holding a reference to the system matrix (created with the method OpType::MakeRef()).

Definition at line 230 of file weakform.hpp.

◆ Init()

void mfem::DPGWeakForm::Init ( )
protected

Definition at line 17 of file weakform.cpp.

◆ RecoverFEMSolution()

void mfem::DPGWeakForm::RecoverFEMSolution ( const Vector X,
Vector x 
)
virtual

Recover the solution of a linear system formed with FormLinearSystem().

Call this method after solving a linear system constructed using the FormLinearSystem() method to recover the solution as a GridFunction-size vector in x. Use the same arguments as in the FormLinearSystem() call.

Reimplemented in mfem::ParDPGWeakForm.

Definition at line 598 of file weakform.cpp.

◆ ReleaseInitMemory()

void mfem::DPGWeakForm::ReleaseInitMemory ( )
protected

Definition at line 629 of file weakform.cpp.

◆ SetDiagonalPolicy()

void mfem::DPGWeakForm::SetDiagonalPolicy ( Operator::DiagonalPolicy  policy)
inline

Sets diagonal policy used upon construction of the linear system.

Policies include:

  • DIAG_ZERO (Set the diagonal values to zero)
  • DIAG_ONE (Set the diagonal values to one)
  • DIAG_KEEP (Keep the diagonal values)

Definition at line 264 of file weakform.hpp.

◆ SetSpaces()

void mfem::DPGWeakForm::SetSpaces ( Array< FiniteElementSpace * > &  fes_,
Array< FiniteElementCollection *> &  fecol_ 
)
inline

Definition at line 125 of file weakform.hpp.

◆ SetTestFECollVdim()

void mfem::DPGWeakForm::SetTestFECollVdim ( int  test_fec,
int  vdim 
)
inline

Definition at line 120 of file weakform.hpp.

◆ Size()

int mfem::DPGWeakForm::Size ( ) const
inline

Get the size of the bilinear form of the DPGWeakForm.

Definition at line 147 of file weakform.hpp.

◆ StoreMatrices()

void mfem::DPGWeakForm::StoreMatrices ( bool  store_matrices_ = true)
inline

Store internal element matrices used for computation of residual after solve.

Definition at line 273 of file weakform.hpp.

◆ Update()

void mfem::DPGWeakForm::Update ( )
virtual

Update the DPGWeakForm after mesh modifications (AMR)

Reimplemented in mfem::ParDPGWeakForm.

Definition at line 671 of file weakform.cpp.

Member Data Documentation

◆ Bmat

Array<DenseMatrix * > mfem::DPGWeakForm::Bmat
protected

Store the matrix L^-1 B and Vector L^-1 l where G = L L^t

Definition at line 101 of file weakform.hpp.

◆ diag_policy

mfem::Operator::DiagonalPolicy mfem::DPGWeakForm::diag_policy
protected

Definition at line 83 of file weakform.hpp.

◆ dof_offsets

Array<int> mfem::DPGWeakForm::dof_offsets
protected

Definition at line 45 of file weakform.hpp.

◆ fvec

Array<Vector * > mfem::DPGWeakForm::fvec
protected

Definition at line 102 of file weakform.hpp.

◆ height

int mfem::DPGWeakForm::height
protected

Definition at line 43 of file weakform.hpp.

◆ initialized

bool mfem::DPGWeakForm::initialized = false
protected

Definition at line 40 of file weakform.hpp.

◆ IsTraceFes

Array<int> mfem::DPGWeakForm::IsTraceFes
protected

Flags to determine if a FiniteElementSpace is Trace.

Definition at line 63 of file weakform.hpp.

◆ lfis

Array<Array<LinearFormIntegrator * > * > mfem::DPGWeakForm::lfis
protected

Set of Linear Form Integrators to be applied.

Definition at line 76 of file weakform.hpp.

◆ mat

BlockMatrix* mfem::DPGWeakForm::mat = nullptr
protected

Block matrix \( M \) to be associated with the Block bilinear form. Owned.

Definition at line 49 of file weakform.hpp.

◆ mat_e

BlockMatrix* mfem::DPGWeakForm::mat_e = nullptr
protected

Block Matrix \( M_e \) used to store the eliminations from the b.c. Owned. \( M + M_e = M_{original} \).

Definition at line 57 of file weakform.hpp.

◆ mesh

Mesh* mfem::DPGWeakForm::mesh = nullptr
protected

Definition at line 42 of file weakform.hpp.

◆ nblocks

int mfem::DPGWeakForm::nblocks
protected

Definition at line 44 of file weakform.hpp.

◆ P

BlockMatrix* mfem::DPGWeakForm::P = nullptr
protected

Block Prolongation.

Definition at line 79 of file weakform.hpp.

◆ R

BlockMatrix* mfem::DPGWeakForm::R = nullptr
protected

Block Restriction.

Definition at line 81 of file weakform.hpp.

◆ residuals

Vector mfem::DPGWeakForm::residuals
protected

Definition at line 103 of file weakform.hpp.

◆ static_cond

BlockStaticCondensation* mfem::DPGWeakForm::static_cond
protected

Owned.

Definition at line 38 of file weakform.hpp.

◆ store_matrices

bool mfem::DPGWeakForm::store_matrices = false
protected

Definition at line 97 of file weakform.hpp.

◆ tdof_offsets

Array<int> mfem::DPGWeakForm::tdof_offsets
protected

Definition at line 46 of file weakform.hpp.

◆ test_fecols

Array<FiniteElementCollection *> mfem::DPGWeakForm::test_fecols
protected

Test FE Collections (Broken)

Definition at line 66 of file weakform.hpp.

◆ test_fecols_vdims

Array<int> mfem::DPGWeakForm::test_fecols_vdims
protected

Definition at line 67 of file weakform.hpp.

◆ test_integs

Array2D<Array<BilinearFormIntegrator * > * > mfem::DPGWeakForm::test_integs
protected

Set of Test Space (broken) Integrators to be applied for matrix G.

Definition at line 73 of file weakform.hpp.

◆ trial_fes

Array<FiniteElementSpace * > mfem::DPGWeakForm::trial_fes
protected

Trial FE spaces.

Definition at line 60 of file weakform.hpp.

◆ trial_integs

Array2D<Array<BilinearFormIntegrator * > * > mfem::DPGWeakForm::trial_integs
protected

Set of Trial Integrators to be applied for matrix B.

Definition at line 70 of file weakform.hpp.

◆ width

int mfem::DPGWeakForm::width
protected

Definition at line 43 of file weakform.hpp.

◆ y

BlockVector* mfem::DPGWeakForm::y = nullptr
protected

Block vector \( y \) to be associated with the Block linear form.

Definition at line 52 of file weakform.hpp.


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