MFEM v2.0
Public Member Functions | Protected Attributes
MixedBilinearForm Class Reference

#include <bilinearform.hpp>

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

List of all members.

Public Member Functions

 MixedBilinearForm (FiniteElementSpace *tr_fes, FiniteElementSpace *te_fes)
int Width () const
virtual double & Elem (int i, int j)
 Returns reference to a_{ij}. Index i, j = 0 .. size-1.
virtual const double & Elem (int i, int j) const
 Returns constant reference to a_{ij}. Index i, j = 0 .. size-1.
virtual void Mult (const Vector &x, Vector &y) const
 Operator application.
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.
virtual MatrixInverseInverse () const
 Returns a pointer to (approximation) of the matrix inverse.
virtual void Finalize (int skip_zeros=1)
 Finalizes the matrix initialization.
void GetBlocks (Array2D< SparseMatrix * > &blocks) const
const SparseMatrixSpMat () const
SparseMatrixSpMat ()
void AddDomainIntegrator (BilinearFormIntegrator *bfi)
void AddBoundaryIntegrator (BilinearFormIntegrator *bfi)
void operator= (const double a)
virtual void Assemble (int skip_zeros=1)
void EliminateTrialDofs (Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs)
virtual void EliminateTestDofs (Array< int > &bdr_attr_is_ess)
void Update ()
virtual ~MixedBilinearForm ()

Protected Attributes

SparseMatrixmat
FiniteElementSpacetrial_fes
FiniteElementSpacetest_fes
Array< BilinearFormIntegrator * > dom
Array< BilinearFormIntegrator * > bdr
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 182 of file bilinearform.hpp.


Constructor & Destructor Documentation

MixedBilinearForm::MixedBilinearForm ( FiniteElementSpace tr_fes,
FiniteElementSpace te_fes 
)

Definition at line 407 of file bilinearform.cpp.

References FiniteElementSpace::GetVSize(), mat, test_fes, trial_fes, and width.

MixedBilinearForm::~MixedBilinearForm ( ) [virtual]

Definition at line 564 of file bilinearform.cpp.

References bdr, dom, mat, and Array< T >::Size().


Member Function Documentation

void MixedBilinearForm::AddBoundaryIntegrator ( BilinearFormIntegrator bfi)

Definition at line 471 of file bilinearform.cpp.

References Array< T >::Append(), and bdr.

void MixedBilinearForm::AddDomainIntegrator ( BilinearFormIntegrator bfi)

Definition at line 466 of file bilinearform.cpp.

References Array< T >::Append(), and dom.

Referenced by DiscreteLinearOperator::AddDomainInterpolator().

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

Definition at line 432 of file bilinearform.cpp.

References mat.

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

Definition at line 438 of file bilinearform.cpp.

References mat.

Referenced by MultTranspose().

void MixedBilinearForm::Assemble ( int  skip_zeros = 1) [virtual]

Reimplemented in DiscreteLinearOperator.

Definition at line 476 of file bilinearform.cpp.

References bdr, dom, mat, Array< T >::Size(), Operator::size, test_fes, trial_fes, and width.

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

Returns reference to a_{ij}. Index i, j = 0 .. size-1.

Implements Matrix.

Definition at line 417 of file bilinearform.cpp.

References mat.

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

Returns constant reference to a_{ij}. Index i, j = 0 .. size-1.

Implements Matrix.

Definition at line 422 of file bilinearform.cpp.

References mat.

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

Definition at line 538 of file bilinearform.cpp.

References mat, Array< T >::Size(), and test_fes.

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

Definition at line 517 of file bilinearform.cpp.

References mat, Array< T >::Size(), and trial_fes.

void MixedBilinearForm::Finalize ( int  = 1) [virtual]

Finalizes the matrix initialization.

Reimplemented from Matrix.

Definition at line 449 of file bilinearform.cpp.

References mat.

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

References Ordering::byNODES, SparseMatrix::GetBlocks(), FiniteElementSpace::GetOrdering(), FiniteElementSpace::GetVDim(), mat, mfem_error(), Array2D< T >::SetSize(), test_fes, and trial_fes.

MatrixInverse * MixedBilinearForm::Inverse ( ) const [virtual]

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

Implements Matrix.

Definition at line 444 of file bilinearform.cpp.

References mat.

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

Operator application.

Implements Operator.

Definition at line 427 of file bilinearform.cpp.

References mat.

virtual void MixedBilinearForm::MultTranspose ( const Vector x,
Vector y 
) const [inline, virtual]

Action of the transpose operator.

Reimplemented from Operator.

Definition at line 212 of file bilinearform.hpp.

References AddMultTranspose().

void MixedBilinearForm::operator= ( const double  a) [inline]

Definition at line 231 of file bilinearform.hpp.

References mat.

SparseMatrix& MixedBilinearForm::SpMat ( ) [inline]

Definition at line 225 of file bilinearform.hpp.

References mat.

const SparseMatrix& MixedBilinearForm::SpMat ( ) const [inline]

Definition at line 224 of file bilinearform.hpp.

References mat.

void MixedBilinearForm::Update ( )
int MixedBilinearForm::Width ( ) const [inline]

Definition at line 198 of file bilinearform.hpp.

References width.


Member Data Documentation

Definition at line 190 of file bilinearform.hpp.

Referenced by AddBoundaryIntegrator(), Assemble(), and ~MixedBilinearForm().

int MixedBilinearForm::width [protected]

The documentation for this class was generated from the following files:
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Defines