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

#include <bilinearform.hpp>

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

Public Member Functions

 MixedBilinearForm (FiniteElementSpace *tr_fes, FiniteElementSpace *te_fes)
 Construct a MixedBilinearForm on the given trial, tr_fes, and test, te_fes, FiniteElementSpaces.
 
 MixedBilinearForm (FiniteElementSpace *tr_fes, FiniteElementSpace *te_fes, MixedBilinearForm *mbf)
 Create a MixedBilinearForm on the given trial, tr_fes, and test, te_fes, FiniteElementSpaces, using the same integrators as the MixedBilinearForm mbf.
 
virtual real_tElem (int i, int j)
 Returns a reference to: \( M_{ij} \).
 
virtual const real_tElem (int i, int j) const
 Returns a reference to: \( M_{ij} \).
 
virtual void Mult (const Vector &x, Vector &y) const
 Matrix multiplication: \( y = M 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 \).
 
virtual void MultTranspose (const Vector &x, Vector &y) const
 Matrix transpose vector multiplication: \( y = M^T 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 \).
 
virtual MatrixInverseInverse () const
 Returns a pointer to (approximation) of the matrix inverse: \( M^{-1} \) (currently unimplemented and returns NULL)
 
virtual void Finalize (int skip_zeros=1)
 Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY.
 
void GetBlocks (Array2D< SparseMatrix * > &blocks) const
 Extract the associated matrix as SparseMatrix blocks. The number of block rows and columns is given by the vector dimensions (vdim) of the test and trial spaces, respectively.
 
const SparseMatrixSpMat () const
 Returns a const reference to the sparse matrix: \( M \).
 
SparseMatrixSpMat ()
 Returns a reference to the sparse matrix: \( M \).
 
SparseMatrixLoseMat ()
 Nullifies the internal matrix \( M \) and returns a pointer to it. Used for transferring ownership.
 
void AddDomainIntegrator (BilinearFormIntegrator *bfi)
 Adds a domain integrator. Assumes ownership of bfi.
 
void AddDomainIntegrator (BilinearFormIntegrator *bfi, Array< int > &elem_marker)
 Adds a domain integrator. Assumes ownership of bfi.
 
void AddBoundaryIntegrator (BilinearFormIntegrator *bfi)
 Adds a boundary integrator. Assumes ownership of bfi.
 
void AddBoundaryIntegrator (BilinearFormIntegrator *bfi, Array< int > &bdr_marker)
 Adds a boundary integrator. Assumes ownership of bfi.
 
void AddTraceFaceIntegrator (BilinearFormIntegrator *bfi)
 Add a trace face integrator. Assumes ownership of bfi.
 
void AddBdrTraceFaceIntegrator (BilinearFormIntegrator *bfi)
 Adds a boundary trace face integrator. Assumes ownership of bfi.
 
void AddBdrTraceFaceIntegrator (BilinearFormIntegrator *bfi, Array< int > &bdr_marker)
 Adds a boundary trace face integrator. Assumes ownership of bfi.
 
Array< BilinearFormIntegrator * > * GetDBFI ()
 Access all integrators added with AddDomainIntegrator().
 
Array< Array< int > * > * GetDBFI_Marker ()
 Access all domain 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 integrators added with AddBoundaryIntegrator().
 
Array< Array< int > * > * GetBBFI_Marker ()
 Access all boundary markers added with AddBoundaryIntegrator().
 
Array< BilinearFormIntegrator * > * GetTFBFI ()
 Access all integrators added with AddTraceFaceIntegrator().
 
Array< BilinearFormIntegrator * > * GetBTFBFI ()
 Access all integrators added with AddBdrTraceFaceIntegrator().
 
Array< Array< int > * > * GetBTFBFI_Marker ()
 Access all boundary markers added with AddBdrTraceFaceIntegrator()
 
void operator= (const real_t a)
 Sets all sparse values of \( M \) to a.
 
void SetAssemblyLevel (AssemblyLevel assembly_level)
 Set the desired assembly level. The default is AssemblyLevel::LEGACY.
 
void Assemble (int skip_zeros=1)
 
void AssembleDiagonal_ADAt (const Vector &D, Vector &diag) const
 Assemble the diagonal of ADA^T into diag, where A is this mixed bilinear form and D is a diagonal.
 
virtual const OperatorGetProlongation () const
 Get the input finite element space prolongation matrix.
 
virtual const OperatorGetRestriction () const
 Get the input finite element space restriction matrix.
 
virtual const OperatorGetOutputProlongation () const
 Get the test finite element space prolongation matrix.
 
virtual const OperatorGetOutputRestriction () const
 Get the test finite element space restriction matrix.
 
void ConformingAssemble ()
 For partially conforming trial and/or test FE spaces, complete the assembly process by performing \( P2^t A P1 \) where \( A \) is the internal sparse matrix; \( P1 \) and \( P2 \) are the conforming prolongation matrices of the trial and test FE spaces, respectively. After this call the MixedBilinearForm becomes an operator on the conforming FE spaces.
 
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 > &trial_vdofs, Array< int > &test_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 > &trial_vdofs, Array< int > &test_vdofs, int skip_zeros=1)
 Assemble the given boundary element matrix.
 
void EliminateTrialDofs (const Array< int > &bdr_attr_is_ess, const Vector &sol, Vector &rhs)
 Eliminate essential boundary DOFs from the columns of the system.
 
void EliminateEssentialBCFromTrialDofs (const Array< int > &marked_vdofs, const Vector &sol, Vector &rhs)
 Eliminate the list of DOFs from the columns of the system.
 
virtual void EliminateTestDofs (const Array< int > &bdr_attr_is_ess)
 Eliminate essential boundary DOFs from the rows of the system.
 
virtual void FormRectangularSystemMatrix (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)
 Return in A that is column-constrained.
 
template<typename OpType >
void FormRectangularSystemMatrix (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OpType &A)
 Form the column-constrained linear system matrix A.
 
virtual void FormRectangularLinearSystem (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)
 Form the linear system A X = B, corresponding to this mixed bilinear form and the linear form b(.).
 
template<typename OpType >
void FormRectangularLinearSystem (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OpType &A, Vector &X, Vector &B)
 Form the linear system A X = B, corresponding to this bilinear form and the linear form b(.).
 
void Update ()
 Must be called after making changes to trial_fes or test_fes.
 
FiniteElementSpaceTrialFESpace ()
 Return the trial FE space associated with the BilinearForm.
 
const FiniteElementSpaceTrialFESpace () const
 Read-only access to the associated trial FiniteElementSpace.
 
FiniteElementSpaceTestFESpace ()
 Return the test FE space associated with the BilinearForm.
 
const FiniteElementSpaceTestFESpace () const
 Read-only access to the associated test FiniteElementSpace.
 
virtual ~MixedBilinearForm ()
 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.
 
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.
 
virtual const OperatorGetOutputRestrictionTranspose () const
 Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type operators.
 
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.
 
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().
 
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 Attributes

SparseMatrixmat
 Owned.
 
SparseMatrixmat_e
 Owned.
 
FiniteElementSpacetrial_fes
 Not owned.
 
FiniteElementSpacetest_fes
 Not owned.
 
AssemblyLevel assembly
 The form assembly level (full, partial, etc.)
 
MixedBilinearFormExtensionext
 
int extern_bfs
 Indicates the BilinearFormIntegrators stored in MixedBilinearForm::domain_integs, MixedBilinearForm::boundary_integs, MixedBilinearForm::trace_face_integs and MixedBilinearForm::boundary_trace_face_integs are owned by another MixedBilinearForm.
 
Array< BilinearFormIntegrator * > domain_integs
 Domain integrators.
 
Array< Array< int > * > domain_integs_marker
 Entries are not owned.
 
Array< BilinearFormIntegrator * > boundary_integs
 Boundary integrators.
 
Array< Array< int > * > boundary_integs_marker
 Entries are not owned.
 
Array< BilinearFormIntegrator * > trace_face_integs
 Trace face (skeleton) integrators.
 
Array< BilinearFormIntegrator * > boundary_trace_face_integs
 Boundary trace face (skeleton) integrators.
 
Array< Array< int > * > boundary_trace_face_integs_marker
 Entries are not owned.
 
DenseMatrix elemmat
 
Array< int > trial_vdofs
 
Array< int > test_vdofs
 
- 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...
 
- 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".
 

Detailed Description

Class for assembling of bilinear forms a(u,v) defined on different trial and test spaces. The assembled matrix M is such that

a(u,v) = V^t M U

where U and V are the vectors representing the functions u and v, respectively. The first argument, u, of a(,) is in the trial space and the second argument, v, is in the test space. Thus,

# of rows of M = dimension of the test space and
# of cols of M = dimension of the trial space.

Both trial and test spaces should be defined on the same mesh.

Definition at line 733 of file bilinearform.hpp.

Constructor & Destructor Documentation

◆ MixedBilinearForm() [1/2]

mfem::MixedBilinearForm::MixedBilinearForm ( FiniteElementSpace * tr_fes,
FiniteElementSpace * te_fes )

Construct a MixedBilinearForm on the given trial, tr_fes, and test, te_fes, FiniteElementSpaces.

The pointers tr_fes and te_fes are not owned by the newly constructed object.

Definition at line 1205 of file bilinearform.cpp.

◆ MixedBilinearForm() [2/2]

mfem::MixedBilinearForm::MixedBilinearForm ( FiniteElementSpace * tr_fes,
FiniteElementSpace * te_fes,
MixedBilinearForm * mbf )

Create a MixedBilinearForm on the given trial, tr_fes, and test, te_fes, FiniteElementSpaces, using the same integrators as the MixedBilinearForm mbf.

The pointers tr_fes and te_fes are not owned by the newly constructed object.

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

Definition at line 1218 of file bilinearform.cpp.

◆ ~MixedBilinearForm()

mfem::MixedBilinearForm::~MixedBilinearForm ( )
virtual

Deletes internal matrices, bilinear integrators, and the BilinearFormExtension.

Definition at line 1926 of file bilinearform.cpp.

Member Function Documentation

◆ AddBdrTraceFaceIntegrator() [1/2]

void mfem::MixedBilinearForm::AddBdrTraceFaceIntegrator ( BilinearFormIntegrator * bfi)

Adds a boundary trace face integrator. Assumes ownership of bfi.

Definition at line 1390 of file bilinearform.cpp.

◆ AddBdrTraceFaceIntegrator() [2/2]

void mfem::MixedBilinearForm::AddBdrTraceFaceIntegrator ( BilinearFormIntegrator * bfi,
Array< int > & bdr_marker )

Adds a boundary trace face integrator. Assumes ownership of bfi.

Definition at line 1397 of file bilinearform.cpp.

◆ AddBoundaryIntegrator() [1/2]

void mfem::MixedBilinearForm::AddBoundaryIntegrator ( BilinearFormIntegrator * bfi)

Adds a boundary integrator. Assumes ownership of bfi.

Definition at line 1372 of file bilinearform.cpp.

◆ AddBoundaryIntegrator() [2/2]

void mfem::MixedBilinearForm::AddBoundaryIntegrator ( BilinearFormIntegrator * bfi,
Array< int > & bdr_marker )

Adds a boundary integrator. Assumes ownership of bfi.

Definition at line 1378 of file bilinearform.cpp.

◆ AddDomainIntegrator() [1/2]

void mfem::MixedBilinearForm::AddDomainIntegrator ( BilinearFormIntegrator * bfi)

Adds a domain integrator. Assumes ownership of bfi.

Definition at line 1359 of file bilinearform.cpp.

◆ AddDomainIntegrator() [2/2]

void mfem::MixedBilinearForm::AddDomainIntegrator ( BilinearFormIntegrator * bfi,
Array< int > & elem_marker )

Adds a domain integrator. Assumes ownership of bfi.

Definition at line 1365 of file bilinearform.cpp.

◆ AddMult()

void mfem::MixedBilinearForm::AddMult ( const Vector & x,
Vector & y,
const real_t a = 1.0 ) const
virtual

Add the matrix vector multiple to a vector: \( y += a M x \).

Reimplemented from mfem::Operator.

Definition at line 1293 of file bilinearform.cpp.

◆ AddMultTranspose()

void mfem::MixedBilinearForm::AddMultTranspose ( const Vector & x,
Vector & y,
const real_t a = 1.0 ) const
virtual

Add the matrix transpose vector multiplication: \( y += a M^T x \).

Reimplemented from mfem::Operator.

Definition at line 1312 of file bilinearform.cpp.

◆ AddTraceFaceIntegrator()

void mfem::MixedBilinearForm::AddTraceFaceIntegrator ( BilinearFormIntegrator * bfi)

Add a trace face integrator. Assumes ownership of bfi.

This type of integrator assembles terms over all faces of the mesh using the face FE from the trial space and the two adjacent volume FEs from the test space.

Definition at line 1385 of file bilinearform.cpp.

◆ Assemble()

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

Definition at line 1404 of file bilinearform.cpp.

◆ AssembleBdrElementMatrix() [1/2]

void mfem::MixedBilinearForm::AssembleBdrElementMatrix ( int i,
const DenseMatrix & elmat,
Array< int > & trial_vdofs,
Array< int > & test_vdofs,
int skip_zeros = 1 )

Assemble the given boundary element matrix.

The boundary element matrix elmat is assembled for the boundary element i, i.e. added to the system matrix. The vdofs of the element are returned in trial_vdofs and test_vdofs. The flag skip_zeros skips the zero elements of the matrix, unless they are breaking the symmetry of the system matrix.

Definition at line 1770 of file bilinearform.cpp.

◆ AssembleBdrElementMatrix() [2/2]

void mfem::MixedBilinearForm::AssembleBdrElementMatrix ( int i,
const DenseMatrix & elmat,
int skip_zeros = 1 )

Assemble the given boundary element matrix.

The boundary element matrix elmat is assembled for the boundary element i, i.e. added to the system matrix. The flag skip_zeros skips the zero elements of the matrix, unless they are breaking the symmetry of the system matrix.

Definition at line 1764 of file bilinearform.cpp.

◆ AssembleDiagonal_ADAt()

void mfem::MixedBilinearForm::AssembleDiagonal_ADAt ( const Vector & D,
Vector & diag ) const

Assemble the diagonal of ADA^T into diag, where A is this mixed bilinear form and D is a diagonal.

Definition at line 1614 of file bilinearform.cpp.

◆ AssembleElementMatrix() [1/2]

void mfem::MixedBilinearForm::AssembleElementMatrix ( int i,
const DenseMatrix & elmat,
Array< int > & trial_vdofs,
Array< int > & test_vdofs,
int skip_zeros = 1 )

Assemble the given element matrix.

The element matrix elmat is assembled for the element i, i.e. added to the system matrix. The vdofs of the element are returned in trial_vdofs and test_vdofs. The flag skip_zeros skips the zero elements of the matrix, unless they are breaking the symmetry of the system matrix.

Definition at line 1751 of file bilinearform.cpp.

◆ AssembleElementMatrix() [2/2]

void mfem::MixedBilinearForm::AssembleElementMatrix ( int i,
const DenseMatrix & elmat,
int skip_zeros = 1 )

Assemble the given element matrix.

The element matrix elmat is assembled for the element i, i.e. added to the system matrix. The flag skip_zeros skips the zero elements of the matrix, unless they are breaking the symmetry of the system matrix.

Definition at line 1745 of file bilinearform.cpp.

◆ ComputeBdrElementMatrix()

void mfem::MixedBilinearForm::ComputeBdrElementMatrix ( int i,
DenseMatrix & elmat )

Compute the boundary element matrix of the given boundary element.

Definition at line 1720 of file bilinearform.cpp.

◆ ComputeElementMatrix()

void mfem::MixedBilinearForm::ComputeElementMatrix ( int i,
DenseMatrix & elmat )

Compute the element matrix of the given element.

Definition at line 1695 of file bilinearform.cpp.

◆ ConformingAssemble()

void mfem::MixedBilinearForm::ConformingAssemble ( )

For partially conforming trial and/or test FE spaces, complete the assembly process by performing \( P2^t A P1 \) where \( A \) is the internal sparse matrix; \( P1 \) and \( P2 \) are the conforming prolongation matrices of the trial and test FE spaces, respectively. After this call the MixedBilinearForm becomes an operator on the conforming FE spaces.

Definition at line 1662 of file bilinearform.cpp.

◆ Elem() [1/2]

real_t & mfem::MixedBilinearForm::Elem ( int i,
int j )
virtual

Returns a reference to: \( M_{ij} \).

Implements mfem::Matrix.

Definition at line 1277 of file bilinearform.cpp.

◆ Elem() [2/2]

const real_t & mfem::MixedBilinearForm::Elem ( int i,
int j ) const
virtual

Returns a reference to: \( M_{ij} \).

Implements mfem::Matrix.

Definition at line 1282 of file bilinearform.cpp.

◆ EliminateEssentialBCFromTrialDofs()

void mfem::MixedBilinearForm::EliminateEssentialBCFromTrialDofs ( const Array< int > & marked_vdofs,
const Vector & sol,
Vector & rhs )

Eliminate the list of DOFs from the columns of the system.

marked_vdofs is the of colunm numbers that will be eliminated. All entries in the columns will be set to 0.0 through elimination.

Definition at line 1806 of file bilinearform.cpp.

◆ EliminateTestDofs()

void mfem::MixedBilinearForm::EliminateTestDofs ( const Array< int > & bdr_attr_is_ess)
virtual

Eliminate essential boundary DOFs from the rows of the system.

The array bdr_attr_is_ess marks boundary attributes that constitute the essential part of the boundary. All entries in the rows will be set to 0.0 through elimination.

Definition at line 1812 of file bilinearform.cpp.

◆ EliminateTrialDofs()

void mfem::MixedBilinearForm::EliminateTrialDofs ( const Array< int > & bdr_attr_is_ess,
const Vector & sol,
Vector & rhs )

Eliminate essential boundary DOFs from the columns of the system.

The array bdr_attr_is_ess marks boundary attributes that constitute the essential part of the boundary. All entries in the columns will be set to 0.0 through elimination.

Definition at line 1783 of file bilinearform.cpp.

◆ Finalize()

void mfem::MixedBilinearForm::Finalize ( int skip_zeros = 1)
virtual

Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY.

Reimplemented from mfem::Matrix.

Definition at line 1339 of file bilinearform.cpp.

◆ FormRectangularLinearSystem() [1/2]

void mfem::MixedBilinearForm::FormRectangularLinearSystem ( const Array< int > & trial_tdof_list,
const Array< int > & test_tdof_list,
Vector & x,
Vector & b,
OperatorHandle & A,
Vector & X,
Vector & B )
virtual

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

Return in A a reference to the system matrix that is column-constrained. The reference will be invalidated when SetOperatorType(), Update(), or the destructor is called.

Reimplemented in mfem::ParMixedBilinearForm.

Definition at line 1885 of file bilinearform.cpp.

◆ FormRectangularLinearSystem() [2/2]

template<typename OpType >
void mfem::MixedBilinearForm::FormRectangularLinearSystem ( const Array< int > & trial_tdof_list,
const Array< int > & test_tdof_list,
Vector & x,
Vector & b,
OpType & A,
Vector & X,
Vector & B )
inline

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

Version of the method FormRectangularLinearSystem() 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 1058 of file bilinearform.hpp.

◆ FormRectangularSystemMatrix() [1/2]

void mfem::MixedBilinearForm::FormRectangularSystemMatrix ( const Array< int > & trial_tdof_list,
const Array< int > & test_tdof_list,
OperatorHandle & A )
virtual

Return in A that is column-constrained.

This returns the same operator as FormRectangularLinearSystem(), but does without the transformations of the right-hand side.

Reimplemented in mfem::ParDiscreteLinearOperator, and mfem::ParMixedBilinearForm.

Definition at line 1832 of file bilinearform.cpp.

◆ FormRectangularSystemMatrix() [2/2]

template<typename OpType >
void mfem::MixedBilinearForm::FormRectangularSystemMatrix ( const Array< int > & trial_tdof_list,
const Array< int > & test_tdof_list,
OpType & A )
inline

Form the column-constrained linear system matrix A.

Version of the method FormRectangularSystemMatrix() 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 1027 of file bilinearform.hpp.

◆ GetBBFI()

Array< BilinearFormIntegrator * > * mfem::MixedBilinearForm::GetBBFI ( )
inline

Access all integrators added with AddBoundaryIntegrator().

Definition at line 887 of file bilinearform.hpp.

◆ GetBBFI_Marker()

Array< Array< int > * > * mfem::MixedBilinearForm::GetBBFI_Marker ( )
inline

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.

Definition at line 893 of file bilinearform.hpp.

◆ GetBlocks()

void mfem::MixedBilinearForm::GetBlocks ( Array2D< SparseMatrix * > & blocks) const

Extract the associated matrix as SparseMatrix blocks. The number of block rows and columns is given by the vector dimensions (vdim) of the test and trial spaces, respectively.

Definition at line 1347 of file bilinearform.cpp.

◆ GetBTFBFI()

Array< BilinearFormIntegrator * > * mfem::MixedBilinearForm::GetBTFBFI ( )
inline

Access all integrators added with AddBdrTraceFaceIntegrator().

Definition at line 899 of file bilinearform.hpp.

◆ GetBTFBFI_Marker()

Array< Array< int > * > * mfem::MixedBilinearForm::GetBTFBFI_Marker ( )
inline

Access all boundary markers added with AddBdrTraceFaceIntegrator()

If no marker was specified when the integrator was added, the corresponding pointer (to Array<int>) will be NULL.

Definition at line 906 of file bilinearform.hpp.

◆ GetDBFI()

Array< BilinearFormIntegrator * > * mfem::MixedBilinearForm::GetDBFI ( )
inline

Access all integrators added with AddDomainIntegrator().

Definition at line 880 of file bilinearform.hpp.

◆ GetDBFI_Marker()

Array< Array< int > * > * mfem::MixedBilinearForm::GetDBFI_Marker ( )
inline

Access all domain markers added with AddDomainIntegrator(). If no marker was specified when the integrator was added, the corresponding pointer (to Array<int>) will be NULL.

Definition at line 884 of file bilinearform.hpp.

◆ GetOutputProlongation()

virtual const Operator * mfem::MixedBilinearForm::GetOutputProlongation ( ) const
inlinevirtual

Get the test finite element space prolongation matrix.

Reimplemented from mfem::Operator.

Definition at line 931 of file bilinearform.hpp.

◆ GetOutputRestriction()

virtual const Operator * mfem::MixedBilinearForm::GetOutputRestriction ( ) const
inlinevirtual

Get the test finite element space restriction matrix.

Reimplemented from mfem::Operator.

Definition at line 935 of file bilinearform.hpp.

◆ GetProlongation()

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

Get the input finite element space prolongation matrix.

Reimplemented from mfem::Operator.

Definition at line 923 of file bilinearform.hpp.

◆ GetRestriction()

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

Get the input finite element space restriction matrix.

Reimplemented from mfem::Operator.

Definition at line 927 of file bilinearform.hpp.

◆ GetTFBFI()

Array< BilinearFormIntegrator * > * mfem::MixedBilinearForm::GetTFBFI ( )
inline

Access all integrators added with AddTraceFaceIntegrator().

Definition at line 896 of file bilinearform.hpp.

◆ Inverse()

MatrixInverse * mfem::MixedBilinearForm::Inverse ( ) const
virtual

Returns a pointer to (approximation) of the matrix inverse: \( M^{-1} \) (currently unimplemented and returns NULL)

Implements mfem::Matrix.

Definition at line 1325 of file bilinearform.cpp.

◆ LoseMat()

SparseMatrix * mfem::MixedBilinearForm::LoseMat ( )
inline

Nullifies the internal matrix \( M \) and returns a pointer to it. Used for transferring ownership.

Definition at line 849 of file bilinearform.hpp.

◆ Mult()

void mfem::MixedBilinearForm::Mult ( const Vector & x,
Vector & y ) const
virtual

Matrix multiplication: \( y = M x \).

Implements mfem::Operator.

Definition at line 1287 of file bilinearform.cpp.

◆ MultTranspose()

void mfem::MixedBilinearForm::MultTranspose ( const Vector & x,
Vector & y ) const
virtual

Matrix transpose vector multiplication: \( y = M^T x \).

Reimplemented from mfem::Operator.

Definition at line 1306 of file bilinearform.cpp.

◆ operator=()

void mfem::MixedBilinearForm::operator= ( const real_t a)
inline

Sets all sparse values of \( M \) to a.

Definition at line 910 of file bilinearform.hpp.

◆ SetAssemblyLevel()

void mfem::MixedBilinearForm::SetAssemblyLevel ( AssemblyLevel assembly_level)

Set the desired assembly level. The default is AssemblyLevel::LEGACY.

This method must be called before assembly. See AssemblyLevel

Definition at line 1246 of file bilinearform.cpp.

◆ SpMat() [1/2]

SparseMatrix & mfem::MixedBilinearForm::SpMat ( )
inline

Returns a reference to the sparse matrix: \( M \).

Definition at line 845 of file bilinearform.hpp.

◆ SpMat() [2/2]

const SparseMatrix & mfem::MixedBilinearForm::SpMat ( ) const
inline

Returns a const reference to the sparse matrix: \( M \).

This will segfault if the usual sparse mat is not defined like when static condensation is being used or AllocMat() has not yet been called.

Definition at line 842 of file bilinearform.hpp.

◆ TestFESpace() [1/2]

FiniteElementSpace * mfem::MixedBilinearForm::TestFESpace ( )
inline

Return the test FE space associated with the BilinearForm.

Definition at line 1081 of file bilinearform.hpp.

◆ TestFESpace() [2/2]

const FiniteElementSpace * mfem::MixedBilinearForm::TestFESpace ( ) const
inline

Read-only access to the associated test FiniteElementSpace.

Definition at line 1084 of file bilinearform.hpp.

◆ TrialFESpace() [1/2]

FiniteElementSpace * mfem::MixedBilinearForm::TrialFESpace ( )
inline

Return the trial FE space associated with the BilinearForm.

Definition at line 1075 of file bilinearform.hpp.

◆ TrialFESpace() [2/2]

const FiniteElementSpace * mfem::MixedBilinearForm::TrialFESpace ( ) const
inline

Read-only access to the associated trial FiniteElementSpace.

Definition at line 1078 of file bilinearform.hpp.

◆ Update()

void mfem::MixedBilinearForm::Update ( )

Must be called after making changes to trial_fes or test_fes.

Definition at line 1915 of file bilinearform.cpp.

Member Data Documentation

◆ assembly

AssemblyLevel mfem::MixedBilinearForm::assembly
protected

The form assembly level (full, partial, etc.)

Definition at line 743 of file bilinearform.hpp.

◆ boundary_integs

Array<BilinearFormIntegrator*> mfem::MixedBilinearForm::boundary_integs
protected

Boundary integrators.

Definition at line 762 of file bilinearform.hpp.

◆ boundary_integs_marker

Array<Array<int>*> mfem::MixedBilinearForm::boundary_integs_marker
protected

Entries are not owned.

Definition at line 764 of file bilinearform.hpp.

◆ boundary_trace_face_integs

Array<BilinearFormIntegrator*> mfem::MixedBilinearForm::boundary_trace_face_integs
protected

Boundary trace face (skeleton) integrators.

Definition at line 770 of file bilinearform.hpp.

◆ boundary_trace_face_integs_marker

Array<Array<int>*> mfem::MixedBilinearForm::boundary_trace_face_integs_marker
protected

Entries are not owned.

Definition at line 772 of file bilinearform.hpp.

◆ domain_integs

Array<BilinearFormIntegrator*> mfem::MixedBilinearForm::domain_integs
protected

Domain integrators.

Definition at line 757 of file bilinearform.hpp.

◆ domain_integs_marker

Array<Array<int>*> mfem::MixedBilinearForm::domain_integs_marker
protected

Entries are not owned.

Definition at line 759 of file bilinearform.hpp.

◆ elemmat

DenseMatrix mfem::MixedBilinearForm::elemmat
protected

Definition at line 774 of file bilinearform.hpp.

◆ ext

MixedBilinearFormExtension* mfem::MixedBilinearForm::ext
protected

Extension for supporting Full Assembly (FA), Element Assembly (EA), Partial Assembly (PA), or Matrix Free assembly (MF).

Definition at line 747 of file bilinearform.hpp.

◆ extern_bfs

int mfem::MixedBilinearForm::extern_bfs
protected

◆ mat

SparseMatrix* mfem::MixedBilinearForm::mat
protected

Owned.

Definition at line 736 of file bilinearform.hpp.

◆ mat_e

SparseMatrix* mfem::MixedBilinearForm::mat_e
protected

Owned.

Definition at line 737 of file bilinearform.hpp.

◆ test_fes

FiniteElementSpace * mfem::MixedBilinearForm::test_fes
protected

Not owned.

Definition at line 740 of file bilinearform.hpp.

◆ test_vdofs

Array<int> mfem::MixedBilinearForm::test_vdofs
protected

Definition at line 775 of file bilinearform.hpp.

◆ trace_face_integs

Array<BilinearFormIntegrator*> mfem::MixedBilinearForm::trace_face_integs
protected

Trace face (skeleton) integrators.

Definition at line 767 of file bilinearform.hpp.

◆ trial_fes

FiniteElementSpace* mfem::MixedBilinearForm::trial_fes
protected

Not owned.

Definition at line 739 of file bilinearform.hpp.

◆ trial_vdofs

Array<int> mfem::MixedBilinearForm::trial_vdofs
protected

Definition at line 775 of file bilinearform.hpp.


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