MFEM  v3.2
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 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
 
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 AddInteriorFaceIntegrator (BilinearFormIntegrator *bfi)
 Adds new interior Face Integrator. More...
 
void AddBdrFaceIntegrator (BilinearFormIntegrator *bfi)
 Adds new boundary Face Integrator. More...
 
void operator= (const double a)
 
void Assemble (int skip_zeros=1)
 Assembles the form i.e. sums over all domain/bdr integrators. More...
 
void FormLinearSystem (Array< int > &ess_tdof_list, Vector &x, Vector &b, SparseMatrix &A, Vector &X, Vector &B, int copy_interior=0)
 
void RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x)
 
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 (Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs, int d=0)
 
void EliminateEssentialBC (Array< int > &bdr_attr_is_ess, int d=0)
 
void EliminateEssentialBCDiag (Array< int > &bdr_attr_is_ess, double value)
 Perform elimination and set the diagonal entry to the given value. More...
 
void EliminateVDofs (Array< int > &vdofs, Vector &sol, Vector &rhs, int d=0)
 Eliminate the given vdofs. NOTE: here, vdofs is a list of DOFs. More...
 
void EliminateVDofs (Array< int > &vdofs, int d=0)
 
void EliminateEssentialBCFromDofs (Array< int > &ess_dofs, Vector &sol, Vector &rhs, int d=0)
 
void EliminateEssentialBCFromDofs (Array< int > &ess_dofs, int d=0)
 
void EliminateEssentialBCFromDofsDiag (Array< int > &ess_dofs, double value)
 Perform elimination and set the diagonal entry to the given value. More...
 
void EliminateVDofsInRHS (Array< int > &vdofs, const Vector &x, Vector &b)
 
double FullInnerProduct (const Vector &x, const Vector &y) const
 
virtual void Update (FiniteElementSpace *nfes=NULL)
 
FiniteElementSpaceGetFES ()
 Return the FE space associated with the BilinearForm. 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=std::cout, 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)
 
int Height () const
 Get the height (size of output) of the Operator. Synonym with NumRows. More...
 
int NumRows () const
 
int Width () const
 Get the width (size of input) of the Operator. Synonym with NumCols. More...
 
int NumCols () const
 
virtual void MultTranspose (const Vector &x, Vector &y) const
 Action of the transpose operator. More...
 
virtual OperatorGetGradient (const Vector &x) const
 Evaluate the gradient operator at the point x. 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 ()
 

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< BilinearFormIntegrator * > fbfi
 Set of interior face Integrators to be applied. More...
 
Array< BilinearFormIntegrator * > bfbfi
 Set of boundary face Integrators to be applied. More...
 
DenseMatrix elemmat
 
Array< int > vdofs
 
DenseTensorelement_matrices
 
StaticCondensationstatic_cond
 
Hybridizationhybridization
 
int precompute_sparsity
 
- Protected Attributes inherited from mfem::Operator
int height
 
int width
 

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

mfem::BilinearForm::~BilinearForm ( )
virtual

Destroys bilinear form.

Definition at line 807 of file bilinearform.cpp.

Member Function Documentation

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

Adds new boundary Face Integrator.

Definition at line 194 of file bilinearform.cpp.

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

Adds new Boundary Integrator.

Definition at line 184 of file bilinearform.cpp.

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

Adds new Domain Integrator.

Definition at line 179 of file bilinearform.cpp.

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

Adds new interior Face Integrator.

Definition at line 189 of file bilinearform.cpp.

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

Definition at line 144 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 120 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 271 of file bilinearform.cpp.

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

Definition at line 249 of file bilinearform.cpp.

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

Definition at line 227 of file bilinearform.cpp.

void mfem::BilinearForm::ComputeElementMatrices ( )

Compute and store internally all element matrices.

Definition at line 593 of file bilinearform.cpp.

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

Definition at line 199 of file bilinearform.cpp.

void mfem::BilinearForm::ConformingAssemble ( )
protected

Definition at line 415 of file bilinearform.cpp.

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

Returns reference to a_{ij}.

Implements mfem::Matrix.

Definition at line 151 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 156 of file bilinearform.cpp.

void mfem::BilinearForm::EliminateEssentialBC ( Array< int > &  bdr_attr_is_ess,
Vector sol,
Vector rhs,
int  d = 0 
)

Eliminate essential boundary DOFs from the system. The array 'bdr_attr_is_ess' marks boundary attributes that constitute the essential part of the boundary. If d == 0, the diagonal at the essential DOFs is set to 1.0, otherwise it is left the same.

Definition at line 635 of file bilinearform.cpp.

void mfem::BilinearForm::EliminateEssentialBC ( Array< int > &  bdr_attr_is_ess,
int  d = 0 
)

Definition at line 652 of file bilinearform.cpp.

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

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

Definition at line 668 of file bilinearform.cpp.

void mfem::BilinearForm::EliminateEssentialBCFromDofs ( Array< int > &  ess_dofs,
Vector sol,
Vector rhs,
int  d = 0 
)

Similar to EliminateVDofs but here ess_dofs is a marker (boolean) array on all vdofs (ess_dofs[i] < 0 is true).

Definition at line 723 of file bilinearform.cpp.

void mfem::BilinearForm::EliminateEssentialBCFromDofs ( Array< int > &  ess_dofs,
int  d = 0 
)

Similar to EliminateVDofs but here ess_dofs is a marker (boolean) array on all vdofs (ess_dofs[i] < 0 is true).

Definition at line 737 of file bilinearform.cpp.

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

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

Definition at line 748 of file bilinearform.cpp.

void mfem::BilinearForm::EliminateVDofs ( Array< int > &  vdofs,
Vector sol,
Vector rhs,
int  d = 0 
)

Eliminate the given vdofs. NOTE: here, vdofs is a list of DOFs.

Definition at line 685 of file bilinearform.cpp.

void mfem::BilinearForm::EliminateVDofs ( Array< int > &  vdofs,
int  d = 0 
)

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

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

Use the stored eliminated part of the matrix (see EliminateVDofs) to modify r.h.s.; vdofs is a list of DOFs (non-directional, i.e. >= 0).

Definition at line 760 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 141 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 124 of file bilinearform.cpp.

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

Finalizes the matrix initialization.

Reimplemented from mfem::Matrix.

Definition at line 171 of file bilinearform.cpp.

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

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

void mfem::BilinearForm::FreeElementMatrices ( )
inline

Free the memory used by the element matrices.

Definition at line 240 of file bilinearform.hpp.

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

Definition at line 147 of file bilinearform.hpp.

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

Definition at line 284 of file bilinearform.hpp.

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

Definition at line 141 of file bilinearform.hpp.

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

Definition at line 124 of file bilinearform.hpp.

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

Definition at line 128 of file bilinearform.hpp.

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

Definition at line 122 of file bilinearform.hpp.

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

Definition at line 126 of file bilinearform.hpp.

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

Return the FE space associated with the BilinearForm.

Definition at line 290 of file bilinearform.hpp.

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

Definition at line 150 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 166 of file bilinearform.cpp.

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

Definition at line 170 of file bilinearform.hpp.

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

Matrix vector multiplication.

Implements mfem::Operator.

Definition at line 161 of file bilinearform.cpp.

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

Definition at line 130 of file bilinearform.hpp.

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

Definition at line 196 of file bilinearform.hpp.

void mfem::BilinearForm::RecoverFEMSolution ( const Vector X,
const Vector b,
Vector x 
)

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.

Definition at line 545 of file bilinearform.cpp.

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

Return the trace FE space associated with static condensation.

Definition at line 102 of file bilinearform.hpp.

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

Get the size of the BilinearForm as a square matrix.

Definition at line 89 of file bilinearform.hpp.

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

Returns a reference to the sparse matrix.

Definition at line 160 of file bilinearform.hpp.

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

Definition at line 165 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 173 of file bilinearform.hpp.

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

Definition at line 178 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 99 of file bilinearform.hpp.

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

Reimplemented in mfem::ParBilinearForm.

Definition at line 767 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 115 of file bilinearform.hpp.

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<BilinearFormIntegrator*> mfem::BilinearForm::bfbfi
protected

Set of boundary face Integrators to be applied.

Definition at line 57 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.

DenseTensor* mfem::BilinearForm::element_matrices
protected

Definition at line 62 of file bilinearform.hpp.

DenseMatrix mfem::BilinearForm::elemmat
protected

Definition at line 59 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 54 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 65 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 67 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 64 of file bilinearform.hpp.

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

Definition at line 60 of file bilinearform.hpp.


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