MFEM  v3.4
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
mfem::BilinearForm Class Reference

#include <bilinearform.hpp>

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

Public Member Functions

 BilinearForm (FiniteElementSpace *f)
 Creates bilinear form associated with FE space *f. More...
 
 BilinearForm (FiniteElementSpace *f, BilinearForm *bf, int ps=0)
 
int Size () const
 Get the size of the BilinearForm as a square matrix. More...
 
void EnableStaticCondensation ()
 
bool StaticCondensationIsEnabled () const
 
FiniteElementSpaceSCFESpace () const
 Return the trace FE space associated with static condensation. More...
 
void EnableHybridization (FiniteElementSpace *constr_space, BilinearFormIntegrator *constr_integ, const Array< int > &ess_tdof_list)
 
void UsePrecomputedSparsity (int ps=1)
 
void UseSparsity (int *I, int *J, bool isSorted)
 Use the given CSR sparsity pattern to allocate the internal SparseMatrix. More...
 
void UseSparsity (SparseMatrix &A)
 Use the sparsity of A to allocate the internal SparseMatrix. More...
 
void AllocateMatrix ()
 
Array< BilinearFormIntegrator * > * GetDBFI ()
 
Array< BilinearFormIntegrator * > * GetBBFI ()
 
Array< BilinearFormIntegrator * > * GetFBFI ()
 
Array< BilinearFormIntegrator * > * GetBFBFI ()
 
const double & operator() (int i, int j)
 
virtual double & Elem (int i, int j)
 Returns reference to a_{ij}. More...
 
virtual const double & Elem (int i, int j) const
 Returns constant reference to a_{ij}. More...
 
virtual void Mult (const Vector &x, Vector &y) const
 Matrix vector multiplication. More...
 
void FullMult (const Vector &x, Vector &y) const
 
virtual void AddMult (const Vector &x, Vector &y, const double a=1.0) const
 
void FullAddMult (const Vector &x, Vector &y) const
 
virtual void AddMultTranspose (const Vector &x, Vector &y, const double a=1.0) const
 
void FullAddMultTranspose (const Vector &x, Vector &y) const
 
virtual void MultTranspose (const Vector &x, Vector &y) const
 Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an error. More...
 
double InnerProduct (const Vector &x, const Vector &y) const
 
virtual MatrixInverseInverse () const
 Returns a pointer to (approximation) of the matrix inverse. More...
 
virtual void Finalize (int skip_zeros=1)
 Finalizes the matrix initialization. More...
 
const SparseMatrixSpMat () const
 Returns a reference to the sparse matrix. More...
 
SparseMatrixSpMat ()
 
SparseMatrixLoseMat ()
 
const SparseMatrixSpMatElim () const
 Returns a reference to the sparse matrix of eliminated b.c. More...
 
SparseMatrixSpMatElim ()
 
void AddDomainIntegrator (BilinearFormIntegrator *bfi)
 Adds new Domain Integrator. More...
 
void AddBoundaryIntegrator (BilinearFormIntegrator *bfi)
 Adds new Boundary Integrator. More...
 
void AddBoundaryIntegrator (BilinearFormIntegrator *bfi, Array< int > &bdr_marker)
 Adds new Boundary Integrator, restricted to specific boundary attributes. More...
 
void AddInteriorFaceIntegrator (BilinearFormIntegrator *bfi)
 Adds new interior Face Integrator. More...
 
void AddBdrFaceIntegrator (BilinearFormIntegrator *bfi)
 Adds new boundary Face Integrator. More...
 
void AddBdrFaceIntegrator (BilinearFormIntegrator *bfi, Array< int > &bdr_marker)
 Adds new boundary Face Integrator, restricted to specific boundary attributes. More...
 
void operator= (const double a)
 
void Assemble (int skip_zeros=1)
 Assembles the form i.e. sums over all domain/bdr integrators. More...
 
virtual const OperatorGetProlongation () const
 Get the finite element space prolongation matrix. More...
 
virtual const OperatorGetRestriction () const
 Get the finite element space restriction matrix. More...
 
void FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, SparseMatrix &A, Vector &X, Vector &B, int copy_interior=0)
 Form a linear system, A X = B. More...
 
void FormSystemMatrix (const Array< int > &ess_tdof_list, SparseMatrix &A)
 Form the linear system matrix A, see FormLinearSystem() for details. More...
 
virtual void RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x)
 Recover the solution of a linear system formed with FormLinearSystem(). More...
 
void ComputeElementMatrices ()
 Compute and store internally all element matrices. More...
 
void FreeElementMatrices ()
 Free the memory used by the element matrices. More...
 
void ComputeElementMatrix (int i, DenseMatrix &elmat)
 
void AssembleElementMatrix (int i, const DenseMatrix &elmat, Array< int > &vdofs, int skip_zeros=1)
 
void AssembleBdrElementMatrix (int i, const DenseMatrix &elmat, Array< int > &vdofs, int skip_zeros=1)
 
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. More...
 
void EliminateEssentialBC (const Array< int > &bdr_attr_is_ess, DiagonalPolicy dpolicy=DIAG_ONE)
 Eliminate essential boundary DOFs from the system matrix. More...
 
void EliminateEssentialBCDiag (const Array< int > &bdr_attr_is_ess, double value)
 Perform elimination and set the diagonal entry to the given value. More...
 
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. More...
 
void EliminateVDofs (const Array< int > &vdofs, DiagonalPolicy dpolicy=DIAG_ONE)
 Eliminate the given vdofs, storing the eliminated part internally. More...
 
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). More...
 
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). More...
 
void EliminateEssentialBCFromDofsDiag (const Array< int > &ess_dofs, double value)
 Perform elimination and set the diagonal entry to the given value. 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 (non-directional, i.e. >= 0). More...
 
double FullInnerProduct (const Vector &x, const Vector &y) const
 
virtual void Update (FiniteElementSpace *nfes=NULL)
 
FiniteElementSpaceGetFES ()
 (DEPRECATED) Return the FE space associated with the BilinearForm. More...
 
FiniteElementSpaceFESpace ()
 Return the FE space associated with the BilinearForm. More...
 
const FiniteElementSpaceFESpace () const
 Read-only access to the associated FiniteElementSpace. More...
 
void SetDiagonalPolicy (DiagonalPolicy policy)
 Sets diagonal policy used upon construction of the linear system. More...
 
virtual ~BilinearForm ()
 Destroys bilinear form. More...
 
- Public Member Functions inherited from mfem::Matrix
 Matrix (int s)
 Creates a square matrix of size s. More...
 
 Matrix (int h, int w)
 Creates a matrix of the given height and width. More...
 
virtual void Print (std::ostream &out=mfem::out, int width_=4) const
 Prints matrix to stream out. More...
 
virtual ~Matrix ()
 Destroys matrix. More...
 
- Public Member Functions inherited from mfem::Operator
 Operator (int s=0)
 Construct a square Operator with given size s (default 0). More...
 
 Operator (int h, int w)
 Construct an Operator with the given height (output size) and width (input size). More...
 
int Height () const
 Get the height (size of output) of the Operator. Synonym with NumRows(). More...
 
int NumRows () const
 Get the number of rows (size of output) of the Operator. Synonym with Height(). More...
 
int Width () const
 Get the width (size of input) of the Operator. Synonym with NumCols(). More...
 
int NumCols () const
 Get the number of columns (size of input) of the Operator. Synonym with Width(). More...
 
virtual OperatorGetGradient (const Vector &x) const
 Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate an error. More...
 
void FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0)
 Form a constrained linear system using a matrix-free approach. More...
 
void PrintMatlab (std::ostream &out, int n=0, int m=0) const
 Prints operator with input size n and output size m in Matlab format. More...
 
virtual ~Operator ()
 Virtual destructor. More...
 
Type GetType () const
 Return the type ID of the Operator class. More...
 

Protected Member Functions

void AllocMat ()
 
void ConformingAssemble ()
 
 BilinearForm ()
 

Protected Attributes

SparseMatrixmat
 Sparse matrix to be associated with the form. More...
 
SparseMatrixmat_e
 Matrix used to eliminate b.c. More...
 
FiniteElementSpacefes
 FE space on which the form lives. More...
 
long sequence
 
int extern_bfs
 
Array< BilinearFormIntegrator * > dbfi
 Set of Domain Integrators to be applied. More...
 
Array< BilinearFormIntegrator * > bbfi
 Set of Boundary Integrators to be applied. More...
 
Array< Array< int > * > bbfi_marker
 
Array< BilinearFormIntegrator * > fbfi
 Set of interior face Integrators to be applied. More...
 
Array< BilinearFormIntegrator * > bfbfi
 Set of boundary face Integrators to be applied. More...
 
Array< Array< int > * > bfbfi_marker
 
DenseMatrix elemmat
 
Array< int > vdofs
 
DenseTensorelement_matrices
 
StaticCondensationstatic_cond
 
Hybridizationhybridization
 
DiagonalPolicy diag_policy
 
int precompute_sparsity
 
- Protected Attributes inherited from mfem::Operator
int height
 Dimension of the output / number of rows in the matrix. More...
 
int width
 Dimension of the input / number of columns in the matrix. More...
 

Additional Inherited Members

- Public Types inherited from mfem::Matrix
enum  DiagonalPolicy { DIAG_ZERO, DIAG_ONE, DIAG_KEEP }
 
- Public Types inherited from mfem::Operator
enum  Type {
  ANY_TYPE, MFEM_SPARSEMAT, Hypre_ParCSR, PETSC_MATAIJ,
  PETSC_MATIS, PETSC_MATSHELL, PETSC_MATNEST, PETSC_MATHYPRE,
  PETSC_MATGENERIC
}
 Enumeration defining IDs for some classes derived from Operator. More...
 

Detailed Description

Class for bilinear form - "Matrix" with associated FE space and BLFIntegrators.

Definition at line 29 of file bilinearform.hpp.

Constructor & Destructor Documentation

mfem::BilinearForm::BilinearForm ( )
inlineprotected

Definition at line 83 of file bilinearform.hpp.

mfem::BilinearForm::BilinearForm ( FiniteElementSpace f)

Creates bilinear form associated with FE space *f.

Definition at line 65 of file bilinearform.cpp.

mfem::BilinearForm::BilinearForm ( FiniteElementSpace f,
BilinearForm bf,
int  ps = 0 
)

Definition at line 79 of file bilinearform.cpp.

mfem::BilinearForm::~BilinearForm ( )
virtual

Destroys bilinear form.

Definition at line 918 of file bilinearform.cpp.

Member Function Documentation

void mfem::BilinearForm::AddBdrFaceIntegrator ( BilinearFormIntegrator bfi)

Adds new boundary Face Integrator.

Definition at line 225 of file bilinearform.cpp.

void mfem::BilinearForm::AddBdrFaceIntegrator ( BilinearFormIntegrator bfi,
Array< int > &  bdr_marker 
)

Adds new boundary Face Integrator, restricted to specific boundary attributes.

Definition at line 231 of file bilinearform.cpp.

void mfem::BilinearForm::AddBoundaryIntegrator ( BilinearFormIntegrator bfi)

Adds new Boundary Integrator.

Definition at line 207 of file bilinearform.cpp.

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

Adds new Boundary Integrator, restricted to specific boundary attributes.

Definition at line 213 of file bilinearform.cpp.

void mfem::BilinearForm::AddDomainIntegrator ( BilinearFormIntegrator bfi)

Adds new Domain Integrator.

Definition at line 202 of file bilinearform.cpp.

void mfem::BilinearForm::AddInteriorFaceIntegrator ( BilinearFormIntegrator bfi)

Adds new interior Face Integrator.

Definition at line 220 of file bilinearform.cpp.

virtual void mfem::BilinearForm::AddMult ( const Vector x,
Vector y,
const double  a = 1.0 
) const
inlinevirtual

Definition at line 168 of file bilinearform.hpp.

virtual void mfem::BilinearForm::AddMultTranspose ( const Vector x,
Vector y,
const double  a = 1.0 
) const
inlinevirtual

Definition at line 174 of file bilinearform.hpp.

void mfem::BilinearForm::AllocateMatrix ( )
inline

Pre-allocate the internal SparseMatrix before assembly. If the flag 'precompute sparsity' is set, the matrix is allocated in CSR format (i.e. finalized) and the entries are initialized with zeros.

Definition at line 144 of file bilinearform.hpp.

void mfem::BilinearForm::AllocMat ( )
protected

Definition at line 20 of file bilinearform.cpp.

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

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

Definition at line 310 of file bilinearform.cpp.

void mfem::BilinearForm::AssembleBdrElementMatrix ( int  i,
const DenseMatrix elmat,
Array< int > &  vdofs,
int  skip_zeros = 1 
)

Definition at line 288 of file bilinearform.cpp.

void mfem::BilinearForm::AssembleElementMatrix ( int  i,
const DenseMatrix elmat,
Array< int > &  vdofs,
int  skip_zeros = 1 
)

Definition at line 266 of file bilinearform.cpp.

void mfem::BilinearForm::ComputeElementMatrices ( )

Compute and store internally all element matrices.

Definition at line 694 of file bilinearform.cpp.

void mfem::BilinearForm::ComputeElementMatrix ( int  i,
DenseMatrix elmat 
)

Definition at line 238 of file bilinearform.cpp.

void mfem::BilinearForm::ConformingAssemble ( )
protected

Definition at line 508 of file bilinearform.cpp.

double & mfem::BilinearForm::Elem ( int  i,
int  j 
)
virtual

Returns reference to a_{ij}.

Implements mfem::Matrix.

Definition at line 179 of file bilinearform.cpp.

const double & mfem::BilinearForm::Elem ( int  i,
int  j 
) const
virtual

Returns constant reference to a_{ij}.

Implements mfem::Matrix.

Definition at line 184 of file bilinearform.cpp.

void mfem::BilinearForm::EliminateEssentialBC ( const Array< int > &  bdr_attr_is_ess,
const Vector sol,
Vector rhs,
DiagonalPolicy  dpolicy = DIAG_ONE 
)

Eliminate essential boundary DOFs from the system.

The array bdr_attr_is_ess marks boundary attributes that constitute the essential part of the boundary. By default, the diagonal at the essential DOFs is set to 1.0. This behavior is controlled by the argument dpolicy.

Definition at line 736 of file bilinearform.cpp.

void mfem::BilinearForm::EliminateEssentialBC ( const Array< int > &  bdr_attr_is_ess,
DiagonalPolicy  dpolicy = DIAG_ONE 
)

Eliminate essential boundary DOFs from the system matrix.

Definition at line 753 of file bilinearform.cpp.

void mfem::BilinearForm::EliminateEssentialBCDiag ( const Array< int > &  bdr_attr_is_ess,
double  value 
)

Perform elimination and set the diagonal entry to the given value.

Definition at line 770 of file bilinearform.cpp.

void mfem::BilinearForm::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).

Definition at line 827 of file bilinearform.cpp.

void mfem::BilinearForm::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).

Definition at line 842 of file bilinearform.cpp.

void mfem::BilinearForm::EliminateEssentialBCFromDofsDiag ( const Array< int > &  ess_dofs,
double  value 
)

Perform elimination and set the diagonal entry to the given value.

Definition at line 854 of file bilinearform.cpp.

void mfem::BilinearForm::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.

Definition at line 787 of file bilinearform.cpp.

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

Eliminate the given vdofs, storing the eliminated part internally.

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 805 of file bilinearform.cpp.

void mfem::BilinearForm::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).

Definition at line 866 of file bilinearform.cpp.

void mfem::BilinearForm::EnableHybridization ( FiniteElementSpace constr_space,
BilinearFormIntegrator constr_integ,
const Array< int > &  ess_tdof_list 
)

Enable hybridization; for details see the description for class Hybridization in fem/hybridization.hpp. This method should be called before assembly.

Definition at line 143 of file bilinearform.cpp.

void mfem::BilinearForm::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.

Definition at line 126 of file bilinearform.cpp.

FiniteElementSpace* mfem::BilinearForm::FESpace ( )
inline

Return the FE space associated with the BilinearForm.

Definition at line 367 of file bilinearform.hpp.

const FiniteElementSpace* mfem::BilinearForm::FESpace ( ) const
inline

Read-only access to the associated FiniteElementSpace.

Definition at line 369 of file bilinearform.hpp.

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

Finalizes the matrix initialization.

Reimplemented from mfem::Matrix.

Definition at line 194 of file bilinearform.cpp.

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

Form a linear system, A X = B.

Form the linear system A X = B, corresponding to the current bilinear form and b(.), by applying any necessary transformations such as: eliminating boundary conditions; applying conforming constraints for non-conforming AMR; 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.

Definition at line 543 of file bilinearform.cpp.

void mfem::BilinearForm::FormSystemMatrix ( const Array< int > &  ess_tdof_list,
SparseMatrix A 
)

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

Definition at line 609 of file bilinearform.cpp.

void mfem::BilinearForm::FreeElementMatrices ( )
inline

Free the memory used by the element matrices.

Definition at line 299 of file bilinearform.hpp.

void mfem::BilinearForm::FullAddMult ( const Vector x,
Vector y 
) const
inline

Definition at line 171 of file bilinearform.hpp.

void mfem::BilinearForm::FullAddMultTranspose ( const Vector x,
Vector y 
) const
inline

Definition at line 178 of file bilinearform.hpp.

double mfem::BilinearForm::FullInnerProduct ( const Vector x,
const Vector y 
) const
inline

Definition at line 357 of file bilinearform.hpp.

void mfem::BilinearForm::FullMult ( const Vector x,
Vector y 
) const
inline

Definition at line 165 of file bilinearform.hpp.

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

Definition at line 148 of file bilinearform.hpp.

Array<BilinearFormIntegrator*>* mfem::BilinearForm::GetBFBFI ( )
inline

Definition at line 152 of file bilinearform.hpp.

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

Definition at line 146 of file bilinearform.hpp.

Array<BilinearFormIntegrator*>* mfem::BilinearForm::GetFBFI ( )
inline

Definition at line 150 of file bilinearform.hpp.

FiniteElementSpace* mfem::BilinearForm::GetFES ( )
inline

(DEPRECATED) Return the FE space associated with the BilinearForm.

Deprecated:
Use FESpace() instead.

Definition at line 364 of file bilinearform.hpp.

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

Get the finite element space prolongation matrix.

Reimplemented from mfem::Operator.

Reimplemented in mfem::ParBilinearForm.

Definition at line 250 of file bilinearform.hpp.

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

Get the finite element space restriction matrix.

Reimplemented from mfem::Operator.

Reimplemented in mfem::ParBilinearForm.

Definition at line 253 of file bilinearform.hpp.

double mfem::BilinearForm::InnerProduct ( const Vector x,
const Vector y 
) const
inline

Definition at line 184 of file bilinearform.hpp.

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

Returns a pointer to (approximation) of the matrix inverse.

Implements mfem::Matrix.

Definition at line 189 of file bilinearform.cpp.

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

Definition at line 204 of file bilinearform.hpp.

virtual void mfem::BilinearForm::Mult ( const Vector x,
Vector y 
) const
inlinevirtual

Matrix vector multiplication.

Implements mfem::Operator.

Definition at line 163 of file bilinearform.hpp.

virtual void mfem::BilinearForm::MultTranspose ( const Vector x,
Vector y 
) const
inlinevirtual

Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an error.

Reimplemented from mfem::Operator.

Definition at line 181 of file bilinearform.hpp.

const double& mfem::BilinearForm::operator() ( int  i,
int  j 
)
inline

Definition at line 154 of file bilinearform.hpp.

void mfem::BilinearForm::operator= ( const double  a)
inline

Definition at line 240 of file bilinearform.hpp.

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

Reimplemented in mfem::ParBilinearForm.

Definition at line 646 of file bilinearform.cpp.

FiniteElementSpace* mfem::BilinearForm::SCFESpace ( ) const
inline

Return the trace FE space associated with static condensation.

Definition at line 112 of file bilinearform.hpp.

void mfem::BilinearForm::SetDiagonalPolicy ( DiagonalPolicy  policy)

Sets diagonal policy used upon construction of the linear system.

Definition at line 913 of file bilinearform.cpp.

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

Get the size of the BilinearForm as a square matrix.

Definition at line 99 of file bilinearform.hpp.

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

Returns a reference to the sparse matrix.

Definition at line 194 of file bilinearform.hpp.

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

Definition at line 199 of file bilinearform.hpp.

const SparseMatrix& mfem::BilinearForm::SpMatElim ( ) const
inline

Returns a reference to the sparse matrix of eliminated b.c.

Definition at line 207 of file bilinearform.hpp.

SparseMatrix& mfem::BilinearForm::SpMatElim ( )
inline

Definition at line 212 of file bilinearform.hpp.

bool mfem::BilinearForm::StaticCondensationIsEnabled ( ) const
inline

Check if static condensation was actually enabled by a previous call to EnableStaticCondensation().

Definition at line 109 of file bilinearform.hpp.

void mfem::BilinearForm::Update ( FiniteElementSpace nfes = NULL)
virtual

Reimplemented in mfem::ParBilinearForm.

Definition at line 873 of file bilinearform.cpp.

void mfem::BilinearForm::UsePrecomputedSparsity ( int  ps = 1)
inline

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.

Definition at line 125 of file bilinearform.hpp.

void mfem::BilinearForm::UseSparsity ( int *  I,
int *  J,
bool  isSorted 
)

Use the given CSR sparsity pattern to allocate the internal SparseMatrix.

  • The I and J arrays must define a square graph with size equal to GetVSize() of the associated FiniteElementSpace.
  • This method should be called after enabling static condensation or hybridization, if used.
  • In the case of static condensation, I and J are not used.
  • The ownership of the arrays I and J remains with the caller.

Definition at line 153 of file bilinearform.cpp.

void mfem::BilinearForm::UseSparsity ( SparseMatrix A)

Use the sparsity of A to allocate the internal SparseMatrix.

Definition at line 169 of file bilinearform.cpp.

Member Data Documentation

Array<BilinearFormIntegrator*> mfem::BilinearForm::bbfi
protected

Set of Boundary Integrators to be applied.

Definition at line 51 of file bilinearform.hpp.

Array<Array<int>*> mfem::BilinearForm::bbfi_marker
protected

Definition at line 52 of file bilinearform.hpp.

Array<BilinearFormIntegrator*> mfem::BilinearForm::bfbfi
protected

Set of boundary face Integrators to be applied.

Definition at line 58 of file bilinearform.hpp.

Array<Array<int>*> mfem::BilinearForm::bfbfi_marker
protected

Definition at line 59 of file bilinearform.hpp.

Array<BilinearFormIntegrator*> mfem::BilinearForm::dbfi
protected

Set of Domain Integrators to be applied.

Definition at line 48 of file bilinearform.hpp.

DiagonalPolicy mfem::BilinearForm::diag_policy
protected

This member allows one to specify what should be done to the diagonal matrix entries and corresponding RHS values upon elimination of the constrained DoFs.

Definition at line 74 of file bilinearform.hpp.

DenseTensor* mfem::BilinearForm::element_matrices
protected

Definition at line 64 of file bilinearform.hpp.

DenseMatrix mfem::BilinearForm::elemmat
protected

Definition at line 61 of file bilinearform.hpp.

int mfem::BilinearForm::extern_bfs
protected

Definition at line 45 of file bilinearform.hpp.

Array<BilinearFormIntegrator*> mfem::BilinearForm::fbfi
protected

Set of interior face Integrators to be applied.

Definition at line 55 of file bilinearform.hpp.

FiniteElementSpace* mfem::BilinearForm::fes
protected

FE space on which the form lives.

Definition at line 39 of file bilinearform.hpp.

Hybridization* mfem::BilinearForm::hybridization
protected

Definition at line 67 of file bilinearform.hpp.

SparseMatrix* mfem::BilinearForm::mat
protected

Sparse matrix to be associated with the form.

Definition at line 33 of file bilinearform.hpp.

SparseMatrix* mfem::BilinearForm::mat_e
protected

Matrix used to eliminate b.c.

Definition at line 36 of file bilinearform.hpp.

int mfem::BilinearForm::precompute_sparsity
protected

Definition at line 76 of file bilinearform.hpp.

long mfem::BilinearForm::sequence
protected

Indicates the Mesh::sequence corresponding to the current state of the BilinearForm.

Definition at line 43 of file bilinearform.hpp.

StaticCondensation* mfem::BilinearForm::static_cond
protected

Definition at line 66 of file bilinearform.hpp.

Array<int> mfem::BilinearForm::vdofs
protected

Definition at line 62 of file bilinearform.hpp.


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