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::BilinearFormIntegrator Class Reference

Abstract base class BilinearFormIntegrator. More...

#include <bilininteg.hpp>

Inherits mfem::NonlinearFormIntegrator.

Inherited by mfem::ConvectionIntegrator, mfem::CurlCurlIntegrator, mfem::DerivativeIntegrator, mfem::DGDiffusionIntegrator, mfem::DGTraceIntegrator, mfem::DiffusionIntegrator, mfem::DiscreteInterpolator, mfem::DivDivIntegrator, mfem::ElasticityIntegrator, mfem::GroupConvectionIntegrator, mfem::InverseIntegrator, mfem::LumpedIntegrator, mfem::MassIntegrator, mfem::NormalTraceJumpIntegrator, mfem::SumIntegrator, mfem::TraceJumpIntegrator, mfem::TransposeIntegrator, mfem::VectorCurlCurlIntegrator, mfem::VectorDiffusionIntegrator, mfem::VectorDivergenceIntegrator, mfem::VectorFECurlIntegrator, mfem::VectorFEDivergenceIntegrator, mfem::VectorFEMassIntegrator, and mfem::VectorMassIntegrator.

Collaboration diagram for mfem::BilinearFormIntegrator:
[legend]

Public Member Functions

virtual void AssembleElementMatrix (const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
 Given a particular Finite Element computes the element matrix elmat. More...
 
virtual void AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
 
virtual void AssembleFaceMatrix (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
virtual void AssembleFaceMatrix (const FiniteElement &trial_face_fe, const FiniteElement &test_fe1, const FiniteElement &test_fe2, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
virtual void AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
 Perform the local action of the BilinearFormIntegrator. More...
 
virtual void AssembleElementGrad (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat)
 Assemble the local gradient matrix. More...
 
virtual void ComputeElementFlux (const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, int with_coef=1)
 
virtual double ComputeFluxEnergy (const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
 
void SetIntRule (const IntegrationRule *ir)
 
virtual ~BilinearFormIntegrator ()
 
- Public Member Functions inherited from mfem::NonlinearFormIntegrator
virtual double GetElementEnergy (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
 Compute the local energy. More...
 
virtual ~NonlinearFormIntegrator ()
 

Protected Member Functions

 BilinearFormIntegrator (const IntegrationRule *ir=NULL)
 

Protected Attributes

const IntegrationRuleIntRule
 

Detailed Description

Abstract base class BilinearFormIntegrator.

Definition at line 22 of file bilininteg.hpp.

Constructor & Destructor Documentation

mfem::BilinearFormIntegrator::BilinearFormIntegrator ( const IntegrationRule ir = NULL)
inlineprotected

Definition at line 27 of file bilininteg.hpp.

virtual mfem::BilinearFormIntegrator::~BilinearFormIntegrator ( )
inlinevirtual

Definition at line 81 of file bilininteg.hpp.

Member Function Documentation

virtual void mfem::BilinearFormIntegrator::AssembleElementGrad ( const FiniteElement el,
ElementTransformation Tr,
const Vector elfun,
DenseMatrix elmat 
)
inlinevirtual

Assemble the local gradient matrix.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 63 of file bilininteg.hpp.

void mfem::BilinearFormIntegrator::AssembleElementMatrix ( const FiniteElement el,
ElementTransformation Trans,
DenseMatrix elmat 
)
virtual
void mfem::BilinearFormIntegrator::AssembleElementMatrix2 ( const FiniteElement trial_fe,
const FiniteElement test_fe,
ElementTransformation Trans,
DenseMatrix elmat 
)
virtual

Compute the local matrix representation of a bilinear form a(u,v) defined on different trial (given by u) and test (given by v) spaces. The rows in the local matrix correspond to the test dofs and the columns – to the trial dofs.

Reimplemented in mfem::NormalInterpolator, mfem::DivergenceInterpolator, mfem::CurlInterpolator, mfem::IdentityInterpolator, mfem::GradientInterpolator, mfem::VectorDivergenceIntegrator, mfem::VectorFEMassIntegrator, mfem::DerivativeIntegrator, mfem::VectorFECurlIntegrator, mfem::VectorFEDivergenceIntegrator, mfem::VectorMassIntegrator, mfem::MassIntegrator, mfem::DiffusionIntegrator, and mfem::TransposeIntegrator.

Definition at line 31 of file bilininteg.cpp.

void mfem::BilinearFormIntegrator::AssembleElementVector ( const FiniteElement el,
ElementTransformation Tr,
const Vector elfun,
Vector elvect 
)
virtual

Perform the local action of the BilinearFormIntegrator.

Implements mfem::NonlinearFormIntegrator.

Reimplemented in mfem::DiffusionIntegrator.

Definition at line 56 of file bilininteg.cpp.

void mfem::BilinearFormIntegrator::AssembleFaceMatrix ( const FiniteElement el1,
const FiniteElement el2,
FaceElementTransformations Trans,
DenseMatrix elmat 
)
virtual
void mfem::BilinearFormIntegrator::AssembleFaceMatrix ( const FiniteElement trial_face_fe,
const FiniteElement test_fe1,
const FiniteElement test_fe2,
FaceElementTransformations Trans,
DenseMatrix elmat 
)
virtual

Abstract method used for assembling TraceFaceIntegrators in a MixedBilinearForm.

Reimplemented in mfem::NormalTraceJumpIntegrator, and mfem::TraceJumpIntegrator.

Definition at line 47 of file bilininteg.cpp.

virtual void mfem::BilinearFormIntegrator::ComputeElementFlux ( const FiniteElement el,
ElementTransformation Trans,
Vector u,
const FiniteElement fluxelem,
Vector flux,
int  with_coef = 1 
)
inlinevirtual

Reimplemented in mfem::CurlCurlIntegrator, and mfem::DiffusionIntegrator.

Definition at line 68 of file bilininteg.hpp.

virtual double mfem::BilinearFormIntegrator::ComputeFluxEnergy ( const FiniteElement fluxelem,
ElementTransformation Trans,
Vector flux,
Vector d_energy = NULL 
)
inlinevirtual

Reimplemented in mfem::CurlCurlIntegrator, and mfem::DiffusionIntegrator.

Definition at line 74 of file bilininteg.hpp.

void mfem::BilinearFormIntegrator::SetIntRule ( const IntegrationRule ir)
inline

Definition at line 79 of file bilininteg.hpp.

Member Data Documentation

const IntegrationRule* mfem::BilinearFormIntegrator::IntRule
protected

Definition at line 25 of file bilininteg.hpp.


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