MFEM  v3.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Protected Attributes | List of all members
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)
 
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
 Operator application. More...
 
virtual void AddMult (const Vector &x, Vector &y, const double a=1.0) const
 
virtual void AddMultTranspose (const Vector &x, Vector &y, const double a=1.0) const
 
virtual void MultTranspose (const Vector &x, Vector &y) const
 Action of the transpose operator. More...
 
virtual MatrixInverseInverse () const
 Returns a pointer to (an approximation) of the matrix inverse. More...
 
virtual void Finalize (int skip_zeros=1)
 Finalizes the matrix initialization. More...
 
void GetBlocks (Array2D< SparseMatrix * > &blocks) const
 
const SparseMatrixSpMat () const
 
SparseMatrixSpMat ()
 
SparseMatrixLoseMat ()
 
void AddDomainIntegrator (BilinearFormIntegrator *bfi)
 
void AddBoundaryIntegrator (BilinearFormIntegrator *bfi)
 
void AddTraceFaceIntegrator (BilinearFormIntegrator *bfi)
 
Array< BilinearFormIntegrator * > * GetDBFI ()
 
Array< BilinearFormIntegrator * > * GetBBFI ()
 
Array< BilinearFormIntegrator * > * GetTFBFI ()
 
void operator= (const double a)
 
void Assemble (int skip_zeros=1)
 
void ConformingAssemble ()
 
void EliminateTrialDofs (Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs)
 
void EliminateEssentialBCFromTrialDofs (Array< int > &marked_vdofs, Vector &sol, Vector &rhs)
 
virtual void EliminateTestDofs (Array< int > &bdr_attr_is_ess)
 
void Update ()
 
virtual ~MixedBilinearForm ()
 
- 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 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 Attributes

SparseMatrixmat
 
FiniteElementSpacetrial_fes
 
FiniteElementSpacetest_fes
 
Array< BilinearFormIntegrator * > dom
 
Array< BilinearFormIntegrator * > bdr
 
Array< BilinearFormIntegrator * > skt
 
- Protected Attributes inherited from mfem::Operator
int height
 
int width
 

Detailed Description

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

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

where U and V are the vectors representing the function 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 A = dimension of the test space and

of cols of A = dimension of the trial space.

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

Definition at line 294 of file bilinearform.hpp.

Constructor & Destructor Documentation

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

Definition at line 773 of file bilinearform.cpp.

mfem::MixedBilinearForm::~MixedBilinearForm ( )
virtual

Definition at line 1014 of file bilinearform.cpp.

Member Function Documentation

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

Definition at line 836 of file bilinearform.cpp.

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

Definition at line 831 of file bilinearform.cpp.

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

Definition at line 797 of file bilinearform.cpp.

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

Definition at line 803 of file bilinearform.cpp.

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

Add a trace face integrator. 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 841 of file bilinearform.cpp.

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

Definition at line 846 of file bilinearform.cpp.

void mfem::MixedBilinearForm::ConformingAssemble ( )

For partially conforming trial and/or test FE spaces, complete the assembly process by performing A := 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 931 of file bilinearform.cpp.

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

Returns reference to a_{ij}.

Implements mfem::Matrix.

Definition at line 782 of file bilinearform.cpp.

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

Returns constant reference to a_{ij}.

Implements mfem::Matrix.

Definition at line 787 of file bilinearform.cpp.

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

Definition at line 980 of file bilinearform.cpp.

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

Definition at line 986 of file bilinearform.cpp.

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

Definition at line 957 of file bilinearform.cpp.

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

Finalizes the matrix initialization.

Reimplemented from mfem::Matrix.

Definition at line 814 of file bilinearform.cpp.

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

Definition at line 348 of file bilinearform.hpp.

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

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

Definition at line 346 of file bilinearform.hpp.

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

Definition at line 350 of file bilinearform.hpp.

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

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

Implements mfem::Matrix.

Definition at line 809 of file bilinearform.cpp.

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

Definition at line 335 of file bilinearform.hpp.

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

Operator application.

Implements mfem::Operator.

Definition at line 792 of file bilinearform.cpp.

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

Action of the transpose operator.

Reimplemented from mfem::Operator.

Definition at line 321 of file bilinearform.hpp.

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

Definition at line 352 of file bilinearform.hpp.

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

Definition at line 333 of file bilinearform.hpp.

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

Definition at line 334 of file bilinearform.hpp.

void mfem::MixedBilinearForm::Update ( )

Definition at line 1006 of file bilinearform.cpp.

Member Data Documentation

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

Definition at line 302 of file bilinearform.hpp.

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

Definition at line 301 of file bilinearform.hpp.

SparseMatrix* mfem::MixedBilinearForm::mat
protected

Definition at line 297 of file bilinearform.hpp.

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

Definition at line 303 of file bilinearform.hpp.

FiniteElementSpace * mfem::MixedBilinearForm::test_fes
protected

Definition at line 299 of file bilinearform.hpp.

FiniteElementSpace* mfem::MixedBilinearForm::trial_fes
protected

Definition at line 299 of file bilinearform.hpp.


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