MFEM
v3.1
Finite element discretization library
|
#include <bilininteg.hpp>
Public Member Functions | |
VectorMassIntegrator () | |
Construct an integrator with coefficient 1.0. More... | |
VectorMassIntegrator (Coefficient &q, int qo=0) | |
VectorMassIntegrator (Coefficient &q, const IntegrationRule *ir) | |
VectorMassIntegrator (VectorCoefficient &q, int qo=0) | |
Construct an integrator with diagonal coefficient q. More... | |
VectorMassIntegrator (MatrixCoefficient &q, int qo=0) | |
Construct an integrator with matrix coefficient q. More... | |
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) |
Public Member Functions inherited from mfem::BilinearFormIntegrator | |
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 () |
Additional Inherited Members | |
Protected Member Functions inherited from mfem::BilinearFormIntegrator | |
BilinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Protected Attributes inherited from mfem::BilinearFormIntegrator | |
const IntegrationRule * | IntRule |
Class for integrating the bilinear form a(u,v) := (Q u, v), where u=(u1,...,un) and v=(v1,...,vn); ui and vi are defined by scalar FE through standard transformation.
Definition at line 288 of file bilininteg.hpp.
|
inline |
Construct an integrator with coefficient 1.0.
Definition at line 302 of file bilininteg.hpp.
|
inline |
Construct an integrator with scalar coefficient q. If possible, save memory by using a scalar integrator since the resulting matrix is block diagonal with the same diagonal block repeated.
Definition at line 308 of file bilininteg.hpp.
|
inline |
Definition at line 310 of file bilininteg.hpp.
|
inline |
Construct an integrator with diagonal coefficient q.
Definition at line 314 of file bilininteg.hpp.
|
inline |
Construct an integrator with matrix coefficient q.
Definition at line 317 of file bilininteg.hpp.
|
virtual |
Given a particular Finite Element computes the element matrix elmat.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 670 of file bilininteg.cpp.
|
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 from mfem::BilinearFormIntegrator.
Definition at line 751 of file bilininteg.cpp.