MFEM  v3.0
 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 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 martix. More...
 
SparseMatrixSpMat ()
 
SparseMatrixLoseMat ()
 
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 ConformingAssemble ()
 
void ConformingAssemble (GridFunction &sol, LinearForm &rhs)
 
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 EliminateEssentialBC (Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs, int d=0)
 
void EliminateVDofs (Array< int > &vdofs, Vector &sol, Vector &rhs, int d=0)
 Here, vdofs is a list of DOFs. More...
 
void EliminateVDofs (Array< int > &vdofs, int d=0)
 
void EliminateVDofsInRHS (Array< int > &vdofs, const Vector &x, Vector &b)
 
double FullInnerProduct (const Vector &x, const Vector &y) const
 
void EliminateEssentialBC (Array< int > &bdr_attr_is_ess, 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 Update (FiniteElementSpace *nfes=NULL)
 
FiniteElementSpaceGetFES ()
 
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)
 Prints operator with input size n and output size m in matlab format. More...
 
virtual ~Operator ()
 

Protected Member Functions

void AllocMat ()
 
 BilinearForm ()
 

Protected Attributes

SparseMatrixmat
 Sparse matrix to be associated with the form. More...
 
SparseMatrixmat_e
 
FiniteElementSpacefes
 FE space on which the form lives. More...
 
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
 
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 27 of file bilinearform.hpp.

Constructor & Destructor Documentation

mfem::BilinearForm::BilinearForm ( )
inlineprotected

Definition at line 63 of file bilinearform.hpp.

mfem::BilinearForm::BilinearForm ( FiniteElementSpace f)

Creates bilinear form associated with FE space *f.

Definition at line 62 of file bilinearform.cpp.

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

Definition at line 72 of file bilinearform.cpp.

mfem::BilinearForm::~BilinearForm ( )
virtual

Destroys bilinear form.

Definition at line 473 of file bilinearform.cpp.

Member Function Documentation

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

Adds new boundary Face Integrator.

Definition at line 149 of file bilinearform.cpp.

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

Adds new Boundary Integrator.

Definition at line 139 of file bilinearform.cpp.

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

Adds new Domain Integrator.

Definition at line 134 of file bilinearform.cpp.

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

Adds new interior Face Integrator.

Definition at line 144 of file bilinearform.cpp.

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

Definition at line 108 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 84 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 191 of file bilinearform.cpp.

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

Definition at line 182 of file bilinearform.cpp.

void mfem::BilinearForm::ComputeElementMatrices ( )

Compute and store internally all element matrices.

Definition at line 333 of file bilinearform.cpp.

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

Definition at line 154 of file bilinearform.cpp.

void mfem::BilinearForm::ConformingAssemble ( )

For partially conforming FE spaces, complete the assembly process by performing A := P^t A P where A is the internal sparse matrix and P is the conforming prolongation of the FE space. After this call the BilinearForm becomes an operator on the conforming FE space.

Definition at line 299 of file bilinearform.cpp.

void mfem::BilinearForm::ConformingAssemble ( GridFunction sol,
LinearForm rhs 
)
inline

A shortcut for converting the whole linear system to conforming DOFs.

Definition at line 153 of file bilinearform.hpp.

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

Returns reference to a_{ij}.

Implements mfem::Matrix.

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

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

If d == 0 the diagonal at the ess. b.c. is set to 1.0, otherwise leave it the same.

Definition at line 373 of file bilinearform.cpp.

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

Definition at line 424 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 439 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 451 of file bilinearform.cpp.

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

Here, vdofs is a list of DOFs.

Definition at line 389 of file bilinearform.cpp.

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

Eliminate the given vdofs storing the eliminated part internally; vdofs is a list of DOFs.

Definition at line 402 of file bilinearform.cpp.

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

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

Definition at line 417 of file bilinearform.cpp.

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

Finalizes the matrix initialization.

Reimplemented from mfem::Matrix.

Definition at line 127 of file bilinearform.cpp.

void mfem::BilinearForm::FreeElementMatrices ( )
inline

Free the memory used by the element matrices.

Definition at line 164 of file bilinearform.hpp.

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

Definition at line 111 of file bilinearform.hpp.

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

Definition at line 187 of file bilinearform.hpp.

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

Definition at line 105 of file bilinearform.hpp.

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

Definition at line 88 of file bilinearform.hpp.

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

Definition at line 92 of file bilinearform.hpp.

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

Definition at line 86 of file bilinearform.hpp.

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

Definition at line 90 of file bilinearform.hpp.

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

Definition at line 203 of file bilinearform.hpp.

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

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

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

Definition at line 126 of file bilinearform.hpp.

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

Matrix vector multiplication.

Implements mfem::Operator.

Definition at line 117 of file bilinearform.cpp.

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

Definition at line 94 of file bilinearform.hpp.

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

Definition at line 140 of file bilinearform.hpp.

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

Get the size of the BilinearForm as a square matrix.

Definition at line 74 of file bilinearform.hpp.

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

Returns a reference to the sparse martix.

Definition at line 124 of file bilinearform.hpp.

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

Definition at line 125 of file bilinearform.hpp.

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

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

Member Data Documentation

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

Set of Boundary Integrators to be applied.

Definition at line 45 of file bilinearform.hpp.

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

Set of boundary face Integrators to be applied.

Definition at line 51 of file bilinearform.hpp.

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

Set of Domain Integrators to be applied.

Definition at line 42 of file bilinearform.hpp.

DenseTensor* mfem::BilinearForm::element_matrices
protected

Definition at line 56 of file bilinearform.hpp.

DenseMatrix mfem::BilinearForm::elemmat
protected

Definition at line 53 of file bilinearform.hpp.

int mfem::BilinearForm::extern_bfs
protected

Definition at line 39 of file bilinearform.hpp.

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

Set of interior face Integrators to be applied.

Definition at line 48 of file bilinearform.hpp.

FiniteElementSpace* mfem::BilinearForm::fes
protected

FE space on which the form lives.

Definition at line 37 of file bilinearform.hpp.

SparseMatrix* mfem::BilinearForm::mat
protected

Sparse matrix to be associated with the form.

Definition at line 31 of file bilinearform.hpp.

SparseMatrix* mfem::BilinearForm::mat_e
protected

Definition at line 34 of file bilinearform.hpp.

int mfem::BilinearForm::precompute_sparsity
protected

Definition at line 58 of file bilinearform.hpp.

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

Definition at line 54 of file bilinearform.hpp.


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