MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::ParBilinearForm Class Reference

Class for parallel bilinear form. More...

#include <pbilinearform.hpp>

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

Public Member Functions

 ParBilinearForm (ParFiniteElementSpace *pf)
 Creates parallel bilinear form associated with the FE space *pf.
 
 ParBilinearForm (ParFiniteElementSpace *pf, ParBilinearForm *bf)
 Create a ParBilinearForm on the ParFiniteElementSpace *pf, using the same integrators as the ParBilinearForm *bf.
 
void KeepNbrBlock (bool knb=true)
 
void SetOperatorType (Operator::Type tid)
 Set the operator type id for the parallel matrix/operator when using AssemblyLevel::LEGACY.
 
void Assemble (int skip_zeros=1)
 Assemble the local matrix.
 
virtual void AssembleDiagonal (Vector &diag) const
 Assemble the diagonal of the bilinear form into diag. Note that diag is a true-dof Vector.
 
HypreParMatrixParallelAssemble ()
 Returns the matrix assembled on the true dofs, i.e. P^t A P.
 
HypreParMatrixParallelAssembleElim ()
 Returns the eliminated matrix assembled on the true dofs, i.e. P^t A_e P.
 
HypreParMatrixParallelAssemble (SparseMatrix *m)
 Return the matrix m assembled on the true dofs, i.e. P^t A P.
 
void ParallelRAP (SparseMatrix &loc_A, OperatorHandle &A, bool steal_loc_A=false)
 Compute parallel RAP operator and store it in A as a HypreParMatrix.
 
void ParallelAssemble (OperatorHandle &A)
 Returns the matrix assembled on the true dofs, i.e. A = P^t A_local P, in the format (type id) specified by A.
 
void ParallelAssembleElim (OperatorHandle &A_elim)
 
void ParallelAssemble (OperatorHandle &A, SparseMatrix *A_local)
 
void ParallelEliminateEssentialBC (const Array< int > &bdr_attr_is_ess, HypreParMatrix &A, const HypreParVector &X, HypreParVector &B) const
 Eliminate essential boundary DOFs from a parallel assembled system.
 
HypreParMatrixParallelEliminateEssentialBC (const Array< int > &bdr_attr_is_ess, HypreParMatrix &A) const
 Eliminate essential boundary DOFs from a parallel assembled matrix A.
 
HypreParMatrixParallelEliminateTDofs (const Array< int > &tdofs_list, HypreParMatrix &A) const
 Eliminate essential true DOFs from a parallel assembled matrix A.
 
void TrueAddMult (const Vector &x, Vector &y, const real_t a=1.0) const
 Compute y += a (P^t A P) x, where x and y are vectors on the true dofs.
 
real_t ParInnerProduct (const ParGridFunction &x, const ParGridFunction &y) const
 Compute \( y^T M x \).
 
real_t TrueInnerProduct (const ParGridFunction &x, const ParGridFunction &y) const
 Compute \( y^T M x \) on true dofs (grid function version)
 
real_t TrueInnerProduct (HypreParVector &x, HypreParVector &y) const
 Compute \( y^T M x \) on true dofs (Hypre vector version)
 
real_t TrueInnerProduct (const Vector &x, const Vector &y) const
 Compute \( y^T M x \) on true dofs (true-vector version)
 
ParFiniteElementSpaceParFESpace () const
 Return the parallel FE space associated with the ParBilinearForm.
 
ParFiniteElementSpaceSCParFESpace () const
 Return the parallel trace FE space associated with static condensation.
 
virtual const OperatorGetProlongation () const
 Get the parallel finite element space prolongation matrix.
 
virtual const OperatorGetRestrictionTranspose () const
 Get the transpose of GetRestriction, useful for matrix-free RAP.
 
virtual const OperatorGetRestriction () const
 Get the parallel finite element space restriction matrix.
 
virtual void FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
 Form the linear system A X = B, corresponding to this bilinear form and the linear form b(.).
 
virtual void FormSystemMatrix (const Array< int > &ess_tdof_list, OperatorHandle &A)
 Form the linear system matrix A, see FormLinearSystem() for details.
 
virtual void RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x)
 
virtual void Update (FiniteElementSpace *nfes=NULL)
 Update the FiniteElementSpace and delete all data associated with the old one.
 
void EliminateVDofsInRHS (const Array< int > &vdofs, const Vector &x, Vector &b)
 
virtual ~ParBilinearForm ()
 
template<typename OpType >
void FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, OpType &A, Vector &X, Vector &B, int copy_interior=0)
 Form the linear system A X = B, corresponding to this bilinear form and the linear form b(.).
 
template<typename OpType >
void FormSystemMatrix (const Array< int > &ess_tdof_list, OpType &A)
 Form the linear system matrix A, see FormLinearSystem() for details.
 
- Public Member Functions inherited from mfem::BilinearForm
 BilinearForm (FiniteElementSpace *f)
 Creates bilinear form associated with FE space *f.
 
 BilinearForm (FiniteElementSpace *f, BilinearForm *bf, int ps=0)
 Create a BilinearForm on the FiniteElementSpace f, using the same integrators as the BilinearForm bf.
 
int Size () const
 Get the size of the BilinearForm as a square matrix.
 
void SetAssemblyLevel (AssemblyLevel assembly_level)
 Set the desired assembly level.
 
void EnableSparseMatrixSorting (bool enable_it)
 Force the sparse matrix column indices to be sorted when using AssemblyLevel::FULL.
 
AssemblyLevel GetAssemblyLevel () const
 Returns the assembly level.
 
HybridizationGetHybridization () const
 
void EnableStaticCondensation ()
 Enable the use of static condensation. For details see the description for class StaticCondensation in fem/staticcond.hpp This method should be called before assembly. If the number of unknowns after static condensation is not reduced, it is not enabled.
 
bool StaticCondensationIsEnabled () const
 Check if static condensation was actually enabled by a previous call to EnableStaticCondensation().
 
FiniteElementSpaceSCFESpace () const
 Return the trace FE space associated with static condensation.
 
void EnableHybridization (FiniteElementSpace *constr_space, BilinearFormIntegrator *constr_integ, const Array< int > &ess_tdof_list)
 Enable hybridization.
 
void UsePrecomputedSparsity (int ps=1)
 For scalar FE spaces, precompute the sparsity pattern of the matrix (assuming dense element matrices) based on the types of integrators present in the bilinear form.
 
void UseSparsity (int *I, int *J, bool isSorted)
 Use the given CSR sparsity pattern to allocate the internal SparseMatrix.
 
void UseSparsity (SparseMatrix &A)
 Use the sparsity of A to allocate the internal SparseMatrix.
 
void AllocateMatrix ()
 Pre-allocate the internal SparseMatrix before assembly. If the internal flag precompute_sparsity is set, the matrix is allocated in CSR format (i.e. finalized) and the entries are initialized with zeros.
 
Array< BilinearFormIntegrator * > * GetDBFI ()
 Access all the integrators added with AddDomainIntegrator().
 
Array< Array< int > * > * GetDBFI_Marker ()
 Access all boundary markers added with AddDomainIntegrator(). If no marker was specified when the integrator was added, the corresponding pointer (to Array<int>) will be NULL.
 
Array< BilinearFormIntegrator * > * GetBBFI ()
 Access all the integrators added with AddBoundaryIntegrator().
 
Array< Array< int > * > * GetBBFI_Marker ()
 Access all boundary markers added with AddBoundaryIntegrator(). If no marker was specified when the integrator was added, the corresponding pointer (to Array<int>) will be NULL.
 
Array< BilinearFormIntegrator * > * GetFBFI ()
 Access all integrators added with AddInteriorFaceIntegrator().
 
Array< BilinearFormIntegrator * > * GetBFBFI ()
 Access all integrators added with AddBdrFaceIntegrator().
 
Array< Array< int > * > * GetBFBFI_Marker ()
 Access all boundary markers added with AddBdrFaceIntegrator(). If no marker was specified when the integrator was added, the corresponding pointer (to Array<int>) will be NULL.
 
const real_toperator() (int i, int j)
 Returns a reference to: \( M_{ij} \).
 
virtual real_tElem (int i, int j)
 Returns a reference to: \( M_{ij} \).
 
virtual const real_tElem (int i, int j) const
 Returns constant reference to: \( M_{ij} \).
 
virtual void Mult (const Vector &x, Vector &y) const
 Matrix vector multiplication: \( y = M x \).
 
void FullMult (const Vector &x, Vector &y) const
 Matrix vector multiplication with the original uneliminated matrix. The original matrix is \( M + M_e \) so we have: \( y = M x + M_e x \).
 
virtual void AddMult (const Vector &x, Vector &y, const real_t a=1.0) const
 Add the matrix vector multiple to a vector: \( y += a M x \).
 
void FullAddMult (const Vector &x, Vector &y) const
 Add the original uneliminated matrix vector multiple to a vector. The original matrix is \( M + Me \) so we have: \( y += M x + M_e x \).
 
virtual void AddMultTranspose (const Vector &x, Vector &y, const real_t a=1.0) const
 Add the matrix transpose vector multiplication: \( y += a M^T x \).
 
void FullAddMultTranspose (const Vector &x, Vector &y) const
 Add the original uneliminated matrix transpose vector multiple to a vector. The original matrix is \( M + M_e \) so we have: \( y += M^T x + {M_e}^T x \).
 
virtual void MultTranspose (const Vector &x, Vector &y) const
 Matrix transpose vector multiplication: \( y = M^T x \).
 
real_t InnerProduct (const Vector &x, const Vector &y) const
 Compute \( y^T M x \).
 
virtual MatrixInverseInverse () const
 Returns a pointer to (approximation) of the matrix inverse: \( M^{-1} \) (currently returns NULL)
 
virtual void Finalize (int skip_zeros=1)
 Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY. THe matrix that gets finalized is different if you are using static condensation or hybridization.
 
const SparseMatrixSpMat () const
 Returns a const reference to the sparse matrix: \( M \).
 
SparseMatrixSpMat ()
 Returns a reference to the sparse matrix: \( M \).
 
bool HasSpMat ()
 Returns true if the sparse matrix is not null, false otherwise.
 
SparseMatrixLoseMat ()
 Nullifies the internal matrix \( M \) and returns a pointer to it. Used for transferring ownership.
 
const SparseMatrixSpMatElim () const
 Returns a const reference to the sparse matrix of eliminated b.c.: \( M_e \).
 
SparseMatrixSpMatElim ()
 Returns a reference to the sparse matrix of eliminated b.c.: \( M_e \).
 
bool HasSpMatElim ()
 Returns true if the sparse matrix of eliminated b.c.s is not null, false otherwise.
 
void AddDomainIntegrator (BilinearFormIntegrator *bfi)
 Adds new Domain Integrator. Assumes ownership of bfi.
 
void AddDomainIntegrator (BilinearFormIntegrator *bfi, Array< int > &elem_marker)
 
void AddBoundaryIntegrator (BilinearFormIntegrator *bfi)
 Adds new Boundary Integrator. Assumes ownership of bfi.
 
void AddBoundaryIntegrator (BilinearFormIntegrator *bfi, Array< int > &bdr_marker)
 Adds new Boundary Integrator, restricted to specific boundary attributes.
 
void AddInteriorFaceIntegrator (BilinearFormIntegrator *bfi)
 Adds new interior Face Integrator. Assumes ownership of bfi.
 
void AddBdrFaceIntegrator (BilinearFormIntegrator *bfi)
 Adds new boundary Face Integrator. Assumes ownership of bfi.
 
void AddBdrFaceIntegrator (BilinearFormIntegrator *bfi, Array< int > &bdr_marker)
 Adds new boundary Face Integrator, restricted to specific boundary attributes.
 
void operator= (const real_t a)
 Sets all sparse values of \( M \) and \( M_e \) to 'a'.
 
void Assemble (int skip_zeros=1)
 Assembles the form i.e. sums over all domain/bdr integrators.
 
virtual const OperatorGetOutputProlongation () const
 Get the output finite element space prolongation matrix.
 
virtual const OperatorGetOutputRestrictionTranspose () const
 Returns the output fe space restriction matrix, transposed.
 
virtual const OperatorGetOutputRestriction () const
 Get the output finite element space restriction matrix.
 
void SerialRAP (OperatorHandle &A)
 Compute serial RAP operator and store it in A as a SparseMatrix.
 
template<typename OpType >
void FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, OpType &A, Vector &X, Vector &B, int copy_interior=0)
 Form the linear system A X = B, corresponding to this bilinear form and the linear form b(.).
 
template<typename OpType >
void FormSystemMatrix (const Array< int > &ess_tdof_list, OpType &A)
 Form the linear system matrix A, see FormLinearSystem() for details.
 
void ComputeElementMatrices ()
 Compute and store internally all element matrices.
 
void FreeElementMatrices ()
 Free the memory used by the element matrices.
 
void ComputeElementMatrix (int i, DenseMatrix &elmat)
 Compute the element matrix of the given element.
 
void ComputeBdrElementMatrix (int i, DenseMatrix &elmat)
 Compute the boundary element matrix of the given boundary element.
 
void AssembleElementMatrix (int i, const DenseMatrix &elmat, int skip_zeros=1)
 Assemble the given element matrix.
 
void AssembleElementMatrix (int i, const DenseMatrix &elmat, Array< int > &vdofs, int skip_zeros=1)
 Assemble the given element matrix.
 
void AssembleBdrElementMatrix (int i, const DenseMatrix &elmat, int skip_zeros=1)
 Assemble the given boundary element matrix.
 
void AssembleBdrElementMatrix (int i, const DenseMatrix &elmat, Array< int > &vdofs, int skip_zeros=1)
 Assemble the given boundary element matrix.
 
void EliminateEssentialBC (const Array< int > &bdr_attr_is_ess, const Vector &sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
 Eliminate essential boundary DOFs from the system.
 
void EliminateEssentialBC (const Array< int > &bdr_attr_is_ess, DiagonalPolicy dpolicy=DIAG_ONE)
 Eliminate essential boundary DOFs from the system matrix.
 
void EliminateEssentialBCDiag (const Array< int > &bdr_attr_is_ess, real_t value)
 Perform elimination and set the diagonal entry to the given value.
 
void EliminateVDofs (const Array< int > &vdofs, const Vector &sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
 Eliminate the given vdofs. NOTE: here, vdofs is a list of DOFs.
 
void EliminateVDofs (const Array< int > &vdofs, DiagonalPolicy dpolicy=DIAG_ONE)
 Eliminate the given vdofs, storing the eliminated part internally in \( M_e \).
 
void EliminateEssentialBCFromDofs (const Array< int > &ess_dofs, const Vector &sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
 Similar to EliminateVDofs(const Array<int> &, const Vector &, Vector &, DiagonalPolicy) but here ess_dofs is a marker (boolean) array on all vector-dofs (ess_dofs[i] < 0 is true).
 
void EliminateEssentialBCFromDofs (const Array< int > &ess_dofs, DiagonalPolicy dpolicy=DIAG_ONE)
 Similar to EliminateVDofs(const Array<int> &, DiagonalPolicy) but here ess_dofs is a marker (boolean) array on all vector-dofs (ess_dofs[i] < 0 is true).
 
void EliminateEssentialBCFromDofsDiag (const Array< int > &ess_dofs, real_t value)
 Perform elimination and set the diagonal entry to the given value.
 
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 (non-directional, i.e. >= 0).
 
real_t FullInnerProduct (const Vector &x, const Vector &y) const
 Compute inner product for full uneliminated matrix: \( y^T M x + y^T M_e x \).
 
MFEM_DEPRECATED FiniteElementSpaceGetFES ()
 (DEPRECATED) Return the FE space associated with the BilinearForm.
 
FiniteElementSpaceFESpace ()
 Return the FE space associated with the BilinearForm.
 
const FiniteElementSpaceFESpace () const
 Read-only access to the associated FiniteElementSpace.
 
void SetDiagonalPolicy (DiagonalPolicy policy)
 Sets Operator::DiagonalPolicy used upon construction of the linear system. Policies include:
 
void UseExternalIntegrators ()
 Indicate that integrators are not owned by the BilinearForm.
 
virtual ~BilinearForm ()
 Deletes internal matrices, bilinear integrators, and the BilinearFormExtension.
 
- Public Member Functions inherited from mfem::Matrix
 Matrix (int s)
 Creates a square matrix of size s.
 
 Matrix (int h, int w)
 Creates a matrix of the given height and width.
 
bool IsSquare () const
 Returns whether the matrix is a square matrix.
 
virtual void Print (std::ostream &out=mfem::out, int width_=4) const
 Prints matrix to stream out.
 
virtual ~Matrix ()
 Destroys matrix.
 
- 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.
 
 Operator (int s=0)
 Construct a square Operator with given size s (default 0).
 
 Operator (int h, int w)
 Construct an Operator with the given height (output size) and width (input size).
 
int Height () const
 Get the height (size of output) of the Operator. Synonym with NumRows().
 
int NumRows () const
 Get the number of rows (size of output) of the Operator. Synonym with Height().
 
int Width () const
 Get the width (size of input) of the Operator. Synonym with NumCols().
 
int NumCols () const
 Get the number of columns (size of input) of the Operator. Synonym with Width().
 
virtual MemoryClass GetMemoryClass () const
 Return the MemoryClass preferred by the Operator.
 
virtual void ArrayMult (const Array< const Vector * > &X, Array< Vector * > &Y) const
 Operator application on a matrix: Y=A(X).
 
virtual void ArrayMultTranspose (const Array< const Vector * > &X, Array< Vector * > &Y) const
 Action of the transpose operator on a matrix: Y=A^t(X).
 
virtual void ArrayAddMult (const Array< const Vector * > &X, Array< Vector * > &Y, const real_t a=1.0) const
 Operator application on a matrix: Y+=A(X) (default) or Y+=a*A(X).
 
virtual void ArrayAddMultTranspose (const Array< const Vector * > &X, Array< Vector * > &Y, const real_t a=1.0) const
 Operator transpose application on a matrix: Y+=A^t(X) (default) or Y+=a*A^t(X).
 
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.
 
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.
 
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.
 
void FormSystemOperator (const Array< int > &ess_tdof_list, Operator *&A)
 Return in A a parallel (on truedofs) version of this square operator.
 
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).
 
void FormDiscreteOperator (Operator *&A)
 Return in A a parallel (on truedofs) version of this rectangular operator.
 
void PrintMatlab (std::ostream &out, int n, int m=0) const
 Prints operator with input size n and output size m in Matlab format.
 
virtual void PrintMatlab (std::ostream &out) const
 Prints operator in Matlab format.
 
virtual ~Operator ()
 Virtual destructor.
 
Type GetType () const
 Return the type ID of the Operator class.
 

Protected Member Functions

void pAllocMat ()
 
void AssembleSharedFaces (int skip_zeros=1)
 
- Protected Member Functions inherited from mfem::BilinearForm
void AllocMat ()
 Allocate appropriate SparseMatrix and assign it to mat.
 
void ConformingAssemble ()
 For partially conforming trial and/or test FE spaces, complete the assembly process by performing \( P^t A P \) where \( A \) is the internal sparse matrix and \( P \) is the conforming prolongation matrix of the trial/test FE space. After this call the BilinearForm becomes an operator on the conforming FE space.
 
 BilinearForm ()
 may be used in the construction of derived classes
 
- Protected Member Functions inherited from mfem::Operator
void FormConstrainedSystemOperator (const Array< int > &ess_tdof_list, ConstrainedOperator *&Aout)
 see FormSystemOperator()
 
void FormRectangularConstrainedSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, RectangularConstrainedOperator *&Aout)
 see FormRectangularSystemOperator()
 
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".
 

Protected Attributes

ParFiniteElementSpacepfes
 Points to the same object as fes.
 
Vector Xaux
 Auxiliary vectors used in TrueAddMult(): L-, L-, and T-vector, resp.
 
Vector Yaux
 
Vector Ytmp
 
OperatorHandle p_mat
 
OperatorHandle p_mat_e
 
bool keep_nbr_block
 
- Protected Attributes inherited from mfem::BilinearForm
SparseMatrixmat
 Sparse matrix \( M \) to be associated with the form. Owned.
 
SparseMatrixmat_e
 Sparse Matrix \( M_e \) used to store the eliminations from the b.c. Owned. \( M + M_e = M_{original} \).
 
FiniteElementSpacefes
 FE space on which the form lives. Not owned.
 
AssemblyLevel assembly
 The AssemblyLevel of the form (AssemblyLevel::LEGACY, AssemblyLevel::FULL, AssemblyLevel::ELEMENT, AssemblyLevel::PARTIAL)
 
int batch
 Element batch size used in the form action (1, 8, num_elems, etc.)
 
BilinearFormExtensionext
 Extension for supporting Full Assembly (FA), Element Assembly (EA),Partial Assembly (PA), or Matrix Free assembly (MF).
 
bool sort_sparse_matrix = false
 
long sequence
 Indicates the Mesh::sequence corresponding to the current state of the BilinearForm.
 
int extern_bfs
 Indicates the BilinearFormIntegrators stored in domain_integs, boundary_integs, interior_face_integs, and boundary_face_integs are owned by another BilinearForm.
 
Array< BilinearFormIntegrator * > domain_integs
 Set of Domain Integrators to be applied.
 
Array< Array< int > * > domain_integs_marker
 Entries are not owned.
 
Array< BilinearFormIntegrator * > boundary_integs
 Set of Boundary Integrators to be applied.
 
Array< Array< int > * > boundary_integs_marker
 Entries are not owned.
 
Array< BilinearFormIntegrator * > interior_face_integs
 Set of interior face Integrators to be applied.
 
Array< BilinearFormIntegrator * > boundary_face_integs
 Set of boundary face Integrators to be applied.
 
Array< Array< int > * > boundary_face_integs_marker
 Entries are not owned.
 
DenseMatrix elemmat
 
Array< int > vdofs
 
DenseTensorelement_matrices
 Owned.
 
StaticCondensationstatic_cond
 Owned.
 
Hybridizationhybridization
 Owned.
 
DiagonalPolicy diag_policy
 This data member allows one to specify what should be done to the diagonal matrix entries and corresponding RHS values upon elimination of the constrained DoFs.
 
int precompute_sparsity
 
- Protected Attributes inherited from mfem::Operator
int height
 Dimension of the output / number of rows in the matrix.
 
int width
 Dimension of the input / number of columns in the matrix.
 

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...
 

Detailed Description

Class for parallel bilinear form.

Definition at line 28 of file pbilinearform.hpp.

Constructor & Destructor Documentation

◆ ParBilinearForm() [1/2]

mfem::ParBilinearForm::ParBilinearForm ( ParFiniteElementSpace * pf)
inline

Creates parallel bilinear form associated with the FE space *pf.

The pointer pf is not owned by the newly constructed object.

Definition at line 56 of file pbilinearform.hpp.

◆ ParBilinearForm() [2/2]

mfem::ParBilinearForm::ParBilinearForm ( ParFiniteElementSpace * pf,
ParBilinearForm * bf )
inline

Create a ParBilinearForm on the ParFiniteElementSpace *pf, using the same integrators as the ParBilinearForm *bf.

The pointer pf is not owned by the newly constructed object.

The integrators in bf are copied as pointers and they are not owned by the newly constructed ParBilinearForm.

Definition at line 68 of file pbilinearform.hpp.

◆ ~ParBilinearForm()

virtual mfem::ParBilinearForm::~ParBilinearForm ( )
inlinevirtual

Definition at line 243 of file pbilinearform.hpp.

Member Function Documentation

◆ Assemble()

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

Assemble the local matrix.

Definition at line 265 of file pbilinearform.cpp.

◆ AssembleDiagonal()

void mfem::ParBilinearForm::AssembleDiagonal ( Vector & diag) const
virtual

Assemble the diagonal of the bilinear form into diag. Note that diag is a true-dof Vector.

When the AssemblyLevel is not LEGACY, and the mesh is nonconforming, this method returns |P^T| d_l, where d_l is the local diagonal of the form before applying parallel/conforming assembly, P^T is the transpose of the parallel/conforming prolongation, and |.| denotes the entry-wise absolute value. In general, this is just an approximation of the exact diagonal for this case.

Reimplemented from mfem::BilinearForm.

Definition at line 284 of file pbilinearform.cpp.

◆ AssembleSharedFaces()

void mfem::ParBilinearForm::AssembleSharedFaces ( int skip_zeros = 1)
protected

Definition at line 220 of file pbilinearform.cpp.

◆ EliminateVDofsInRHS()

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

Definition at line 487 of file pbilinearform.cpp.

◆ FormLinearSystem() [1/2]

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

Form the linear system A X = B, corresponding to this bilinear form and the linear form b(.).

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

The GridFunction-size vector x must contain the essential b.c. The BilinearForm and the LinearForm-size vector b must be assembled.

The vector X is initialized with a suitable initial guess: when using hybridization, the vector X is set to zero; otherwise, 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).

This method can be called multiple times (with the same ess_tdof_list array) to initialize different right-hand sides and boundary condition values.

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

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

Reimplemented from mfem::BilinearForm.

Definition at line 437 of file pbilinearform.cpp.

◆ FormLinearSystem() [2/2]

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

Form the linear system A X = B, corresponding to this bilinear form and the linear form b(.).

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()). The reference will be invalidated when SetOperatorType(), Update(), or the destructor is called.

Definition at line 533 of file bilinearform.hpp.

◆ FormSystemMatrix() [1/2]

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

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

Reimplemented from mfem::BilinearForm.

Definition at line 493 of file pbilinearform.cpp.

◆ FormSystemMatrix() [2/2]

template<typename OpType >
void mfem::BilinearForm::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()). The reference will be invalidated when SetOperatorType(), Update(), or the destructor is called.

Definition at line 555 of file bilinearform.hpp.

◆ GetProlongation()

virtual const Operator * mfem::ParBilinearForm::GetProlongation ( ) const
inlinevirtual

Get the parallel finite element space prolongation matrix.

Reimplemented from mfem::BilinearForm.

Definition at line 215 of file pbilinearform.hpp.

◆ GetRestriction()

virtual const Operator * mfem::ParBilinearForm::GetRestriction ( ) const
inlinevirtual

Get the parallel finite element space restriction matrix.

Reimplemented from mfem::BilinearForm.

Definition at line 221 of file pbilinearform.hpp.

◆ GetRestrictionTranspose()

virtual const Operator * mfem::ParBilinearForm::GetRestrictionTranspose ( ) const
inlinevirtual

Get the transpose of GetRestriction, useful for matrix-free RAP.

Definition at line 218 of file pbilinearform.hpp.

◆ KeepNbrBlock()

void mfem::ParBilinearForm::KeepNbrBlock ( bool knb = true)
inline

When set to true and the ParBilinearForm has interior face integrators, the local SparseMatrix will include the rows (in addition to the columns) corresponding to face-neighbor dofs. The default behavior is to disregard those rows. Must be called before the first Assemble call.

Definition at line 77 of file pbilinearform.hpp.

◆ pAllocMat()

void mfem::ParBilinearForm::pAllocMat ( )
protected

Definition at line 22 of file pbilinearform.cpp.

◆ ParallelAssemble() [1/4]

HypreParMatrix * mfem::ParBilinearForm::ParallelAssemble ( )
inline

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 106 of file pbilinearform.hpp.

◆ ParallelAssemble() [2/4]

void mfem::ParBilinearForm::ParallelAssemble ( OperatorHandle & A)
inline

Returns the matrix assembled on the true dofs, i.e. A = P^t A_local P, in the format (type id) specified by A.

Definition at line 129 of file pbilinearform.hpp.

◆ ParallelAssemble() [3/4]

void mfem::ParBilinearForm::ParallelAssemble ( OperatorHandle & A,
SparseMatrix * A_local )

Returns the matrix A_local assembled on the true dofs, i.e. A = P^t A_local P in the format (type id) specified by A.

Definition at line 154 of file pbilinearform.cpp.

◆ ParallelAssemble() [4/4]

HypreParMatrix * mfem::ParBilinearForm::ParallelAssemble ( SparseMatrix * m)

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

The returned matrix has to be deleted by the caller.

Definition at line 212 of file pbilinearform.cpp.

◆ ParallelAssembleElim() [1/2]

HypreParMatrix * mfem::ParBilinearForm::ParallelAssembleElim ( )
inline

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

The returned matrix has to be deleted by the caller.

Definition at line 110 of file pbilinearform.hpp.

◆ ParallelAssembleElim() [2/2]

void mfem::ParBilinearForm::ParallelAssembleElim ( OperatorHandle & A_elim)
inline

Returns the eliminated matrix assembled on the true dofs, i.e. A_elim = P^t A_elim_local P in the format (type id) specified by A.

Definition at line 134 of file pbilinearform.hpp.

◆ ParallelEliminateEssentialBC() [1/2]

HypreParMatrix * mfem::ParBilinearForm::ParallelEliminateEssentialBC ( const Array< int > & bdr_attr_is_ess,
HypreParMatrix & A ) const

Eliminate essential boundary DOFs from a parallel assembled matrix A.

The array bdr_attr_is_ess marks boundary attributes that constitute the essential part of the boundary. The eliminated part is stored in a matrix A_elim such that A_original = A_new + A_elim. Returns a pointer to the newly allocated matrix A_elim which should be deleted by the caller. The matrices A and A_elim can be used to eliminate boundary conditions in multiple right-hand sides, by calling the function EliminateBC() (from hypre.hpp).

Definition at line 336 of file pbilinearform.cpp.

◆ ParallelEliminateEssentialBC() [2/2]

void mfem::ParBilinearForm::ParallelEliminateEssentialBC ( const Array< int > & bdr_attr_is_ess,
HypreParMatrix & A,
const HypreParVector & X,
HypreParVector & B ) const

Eliminate essential boundary DOFs from a parallel assembled system.

The array bdr_attr_is_ess marks boundary attributes that constitute the essential part of the boundary.

Definition at line 323 of file pbilinearform.cpp.

◆ ParallelEliminateTDofs()

HypreParMatrix * mfem::ParBilinearForm::ParallelEliminateTDofs ( const Array< int > & tdofs_list,
HypreParMatrix & A ) const
inline

Eliminate essential true DOFs from a parallel assembled matrix A.

Given a list of essential true dofs and the parallel assembled matrix A, eliminate the true dofs from the matrix, storing the eliminated part in a matrix A_elim such that A_original = A_new + A_elim. Returns a pointer to the newly allocated matrix A_elim which should be deleted by the caller. The matrices A and A_elim can be used to eliminate boundary conditions in multiple right-hand sides, by calling the function EliminateBC() (from hypre.hpp).

Definition at line 168 of file pbilinearform.hpp.

◆ ParallelRAP()

void mfem::ParBilinearForm::ParallelRAP ( SparseMatrix & loc_A,
OperatorHandle & A,
bool steal_loc_A = false )

Compute parallel RAP operator and store it in A as a HypreParMatrix.

Parameters
[in]loc_AThe rank-local SparseMatrix.
[out]AThe OperatorHandle containing the global HypreParMatrix.
[in]steal_loc_AHave the HypreParMatrix in A take ownership of the memory objects in loc_A.

Definition at line 124 of file pbilinearform.cpp.

◆ ParFESpace()

ParFiniteElementSpace * mfem::ParBilinearForm::ParFESpace ( ) const
inline

Return the parallel FE space associated with the ParBilinearForm.

Definition at line 208 of file pbilinearform.hpp.

◆ ParInnerProduct()

real_t mfem::ParBilinearForm::ParInnerProduct ( const ParGridFunction & x,
const ParGridFunction & y ) const

Compute \( y^T M x \).

Warning
The calculation is performed on local dofs, assuming that the local vectors are consistent with the prolongations of the true vectors (see ParGridFunction::Distribute()). If this is not the case, use TrueInnerProduct(const ParGridFunction &, const ParGridFunction &) instead.
Note
It is assumed that the local matrix is assembled and it has not been replaced by the parallel matrix through FormSystemMatrix().
See also
TrueInnerProduct(const ParGridFunction&, const ParGridFunction&)

Definition at line 371 of file pbilinearform.cpp.

◆ RecoverFEMSolution()

void mfem::ParBilinearForm::RecoverFEMSolution ( const Vector & X,
const Vector & b,
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::BilinearForm.

Definition at line 541 of file pbilinearform.cpp.

◆ SCParFESpace()

ParFiniteElementSpace * mfem::ParBilinearForm::SCParFESpace ( ) const
inline

Return the parallel trace FE space associated with static condensation.

Definition at line 211 of file pbilinearform.hpp.

◆ SetOperatorType()

void mfem::ParBilinearForm::SetOperatorType ( Operator::Type tid)
inline

Set the operator type id for the parallel matrix/operator when using AssemblyLevel::LEGACY.

If using static condensation or hybridization, call this method after enabling it.

Definition at line 83 of file pbilinearform.hpp.

◆ TrueAddMult()

void mfem::ParBilinearForm::TrueAddMult ( const Vector & x,
Vector & y,
const real_t a = 1.0 ) const

Compute y += a (P^t A P) x, where x and y are vectors on the true dofs.

Definition at line 347 of file pbilinearform.cpp.

◆ TrueInnerProduct() [1/3]

real_t mfem::ParBilinearForm::TrueInnerProduct ( const ParGridFunction & x,
const ParGridFunction & y ) const

Compute \( y^T M x \) on true dofs (grid function version)

Note
The ParGridFunctions are restricted to the true-vectors for for calculation.
It is assumed that the parallel system matrix is assembled, see FormSystemMatrix().
See also
ParInnerProduct(const ParGridFunction&, const ParGridFunction&)

Definition at line 385 of file pbilinearform.cpp.

◆ TrueInnerProduct() [2/3]

real_t mfem::ParBilinearForm::TrueInnerProduct ( const Vector & x,
const Vector & y ) const

Compute \( y^T M x \) on true dofs (true-vector version)

Note
It is assumed that the parallel system matrix is assembled, see FormSystemMatrix().

Definition at line 424 of file pbilinearform.cpp.

◆ TrueInnerProduct() [3/3]

real_t mfem::ParBilinearForm::TrueInnerProduct ( HypreParVector & x,
HypreParVector & y ) const

Compute \( y^T M x \) on true dofs (Hypre vector version)

Note
It is assumed that the parallel system matrix is assembled, see FormSystemMatrix().

Definition at line 402 of file pbilinearform.cpp.

◆ Update()

void mfem::ParBilinearForm::Update ( FiniteElementSpace * nfes = NULL)
virtual

Update the FiniteElementSpace and delete all data associated with the old one.

Reimplemented from mfem::BilinearForm.

Definition at line 576 of file pbilinearform.cpp.

Member Data Documentation

◆ keep_nbr_block

bool mfem::ParBilinearForm::keep_nbr_block
protected

Definition at line 39 of file pbilinearform.hpp.

◆ p_mat

OperatorHandle mfem::ParBilinearForm::p_mat
protected

Definition at line 37 of file pbilinearform.hpp.

◆ p_mat_e

OperatorHandle mfem::ParBilinearForm::p_mat_e
protected

Definition at line 37 of file pbilinearform.hpp.

◆ pfes

ParFiniteElementSpace* mfem::ParBilinearForm::pfes
protected

Points to the same object as fes.

Definition at line 32 of file pbilinearform.hpp.

◆ Xaux

Vector mfem::ParBilinearForm::Xaux
mutableprotected

Auxiliary vectors used in TrueAddMult(): L-, L-, and T-vector, resp.

Definition at line 35 of file pbilinearform.hpp.

◆ Yaux

Vector mfem::ParBilinearForm::Yaux
protected

Definition at line 35 of file pbilinearform.hpp.

◆ Ytmp

Vector mfem::ParBilinearForm::Ytmp
protected

Definition at line 35 of file pbilinearform.hpp.


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