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

Class representing the parallel weak formulation. (Convenient for DPG Equations) More...

#include <pweakform.hpp>

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

Public Member Functions

 ParDPGWeakForm ()
 
 ParDPGWeakForm (Array< ParFiniteElementSpace * > &trial_pfes_, Array< FiniteElementCollection * > &fecol_)
 Creates bilinear form associated with FE spaces trial_pfes_. More...
 
void SetParSpaces (Array< ParFiniteElementSpace * > &trial_pfes_, Array< FiniteElementCollection * > &fecol_)
 
void Assemble (int skip_zeros=1)
 Assemble the local matrix. More...
 
void ParallelAssemble (BlockMatrix *mat)
 Returns the matrix assembled on the true dofs, i.e. P^t A P. More...
 
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...
 
void FormSystemMatrix (const Array< int > &ess_tdof_list, OperatorHandle &A)
 Form the linear system matrix A, see FormLinearSystem() for details. More...
 
virtual void RecoverFEMSolution (const Vector &X, Vector &x)
 
virtual void Update ()
 Update the DPGWeakForm after mesh modifications (AMR) More...
 
virtual ~ParDPGWeakForm ()
 Destroys bilinear form. More...
 
- Public Member Functions inherited from mfem::DPGWeakForm
 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...
 
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...
 
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...
 
void SetDiagonalPolicy (Operator::DiagonalPolicy policy)
 Sets diagonal policy used upon construction of the linear system. 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 FillEssTdofLists (const Array< int > &ess_tdof_list)
 
void BuildProlongation ()
 
- Protected Member Functions inherited from mfem::DPGWeakForm
void Init ()
 
void ReleaseInitMemory ()
 
void AllocMat ()
 Allocate appropriate BlockMatrix and assign it to mat. More...
 
void ConformingAssemble ()
 
void ComputeOffsets ()
 

Protected Attributes

Array< ParFiniteElementSpace *> trial_pfes
 
Array< Array< int > * > ess_tdofs
 
BlockOperatorP = nullptr
 
BlockMatrixR = nullptr
 
BlockOperatorp_mat = nullptr
 
BlockOperatorp_mat_e = nullptr
 
- Protected Attributes inherited from mfem::DPGWeakForm
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 parallel weak formulation. (Convenient for DPG Equations)

Definition at line 25 of file pweakform.hpp.

Constructor & Destructor Documentation

◆ ParDPGWeakForm() [1/2]

mfem::ParDPGWeakForm::ParDPGWeakForm ( )
inline

Definition at line 53 of file pweakform.hpp.

◆ ParDPGWeakForm() [2/2]

mfem::ParDPGWeakForm::ParDPGWeakForm ( Array< ParFiniteElementSpace * > &  trial_pfes_,
Array< FiniteElementCollection * > &  fecol_ 
)
inline

Creates bilinear form associated with FE spaces trial_pfes_.

Definition at line 56 of file pweakform.hpp.

◆ ~ParDPGWeakForm()

mfem::ParDPGWeakForm::~ParDPGWeakForm ( )
virtual

Destroys bilinear form.

Definition at line 206 of file pweakform.cpp.

Member Function Documentation

◆ Assemble()

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

Assemble the local matrix.

Definition at line 33 of file pweakform.cpp.

◆ BuildProlongation()

void mfem::ParDPGWeakForm::BuildProlongation ( )
protectedvirtual

Reimplemented from mfem::DPGWeakForm.

Definition at line 84 of file pweakform.cpp.

◆ FillEssTdofLists()

void mfem::ParDPGWeakForm::FillEssTdofLists ( const Array< int > &  ess_tdof_list)
protected

Definition at line 19 of file pweakform.cpp.

◆ FormLinearSystem()

void mfem::ParDPGWeakForm::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 from mfem::DPGWeakForm.

Definition at line 100 of file pweakform.cpp.

◆ FormSystemMatrix()

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

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

Reimplemented from mfem::DPGWeakForm.

Definition at line 140 of file pweakform.cpp.

◆ ParallelAssemble()

void mfem::ParDPGWeakForm::ParallelAssemble ( BlockMatrix mat)

Returns the matrix assembled on the true dofs, i.e. P^t A P.

The returned matrix has to be deleted by the caller.

Definition at line 38 of file pweakform.cpp.

◆ RecoverFEMSolution()

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

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

Reimplemented from mfem::DPGWeakForm.

Definition at line 173 of file pweakform.cpp.

◆ SetParSpaces()

void mfem::ParDPGWeakForm::SetParSpaces ( Array< ParFiniteElementSpace * > &  trial_pfes_,
Array< FiniteElementCollection * > &  fecol_ 
)
inline

Definition at line 63 of file pweakform.hpp.

◆ Update()

void mfem::ParDPGWeakForm::Update ( )
virtual

Update the DPGWeakForm after mesh modifications (AMR)

Reimplemented from mfem::DPGWeakForm.

Definition at line 188 of file pweakform.cpp.

Member Data Documentation

◆ ess_tdofs

Array<Array<int> *> mfem::ParDPGWeakForm::ess_tdofs
protected

Definition at line 33 of file pweakform.hpp.

◆ P

BlockOperator* mfem::ParDPGWeakForm::P = nullptr
protected

Definition at line 40 of file pweakform.hpp.

◆ p_mat

BlockOperator* mfem::ParDPGWeakForm::p_mat = nullptr
protected

Definition at line 44 of file pweakform.hpp.

◆ p_mat_e

BlockOperator* mfem::ParDPGWeakForm::p_mat_e = nullptr
protected

Definition at line 45 of file pweakform.hpp.

◆ R

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

Definition at line 41 of file pweakform.hpp.

◆ trial_pfes

Array<ParFiniteElementSpace * > mfem::ParDPGWeakForm::trial_pfes
protected

Definition at line 30 of file pweakform.hpp.


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