MFEM
v4.0
Finite element discretization library
|
#include <bilininteg.hpp>
Public Member Functions | |
virtual void | AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat) |
virtual void | AssembleElementMatrix (const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &elmat) |
Support for use in BilinearForm. Can be used only when appropriate. More... | |
Public Member Functions inherited from mfem::BilinearFormIntegrator | |
virtual void | AssemblePA (const FiniteElementSpace &fes) |
Method defining partial assembly. More... | |
virtual void | AddMultPA (const Vector &x, Vector &y) const |
Method for partially assembled action. More... | |
virtual void | AddMultTransposePA (const Vector &x, Vector &y) const |
Method for partially assembled transposed action. More... | |
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 | AssembleFaceGrad (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, DenseMatrix &elmat) |
Assemble the local action of the gradient of the NonlinearFormIntegrator resulting from a face integral term. More... | |
virtual void | ComputeElementFlux (const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, int with_coef=1) |
Virtual method required for Zienkiewicz-Zhu type error estimators. More... | |
virtual double | ComputeFluxEnergy (const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL) |
Virtual method required for Zienkiewicz-Zhu type error estimators. More... | |
virtual | ~BilinearFormIntegrator () |
Public Member Functions inherited from mfem::NonlinearFormIntegrator | |
void | SetIntRule (const IntegrationRule *ir) |
Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == NULL). More... | |
void | SetIntegrationRule (const IntegrationRule &irule) |
Prescribe a fixed IntegrationRule to use. More... | |
virtual void | AssembleFaceVector (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect) |
Perform the local action of the NonlinearFormIntegrator resulting from a face integral term. More... | |
virtual double | GetElementEnergy (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun) |
Compute the local energy. More... | |
virtual | ~NonlinearFormIntegrator () |
Protected Member Functions | |
MixedScalarIntegrator () | |
MixedScalarIntegrator (Coefficient &q) | |
virtual bool | VerifyFiniteElementTypes (const FiniteElement &trial_fe, const FiniteElement &test_fe) const |
virtual const char * | FiniteElementTypeFailureMessage () const |
virtual int | GetIntegrationOrder (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans) |
virtual void | CalcTestShape (const FiniteElement &test_fe, ElementTransformation &Trans, Vector &shape) |
virtual void | CalcTrialShape (const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape) |
Protected Member Functions inherited from mfem::BilinearFormIntegrator | |
BilinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Protected Member Functions inherited from mfem::NonlinearFormIntegrator | |
NonlinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Protected Attributes | |
bool | same_calc_shape |
Coefficient * | Q |
Protected Attributes inherited from mfem::NonlinearFormIntegrator | |
const IntegrationRule * | IntRule |
An abstract class for integrating the product of two scalar basis functions with an optional scalar coefficient.
Definition at line 258 of file bilininteg.hpp.
|
inlineprotected |
Definition at line 279 of file bilininteg.hpp.
|
inlineprotected |
Definition at line 280 of file bilininteg.hpp.
|
inlinevirtual |
Support for use in BilinearForm. Can be used only when appropriate.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 268 of file bilininteg.hpp.
|
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 147 of file bilininteg.cpp.
|
inlineprotectedvirtual |
Reimplemented in mfem::MixedScalarWeakCurlIntegrator, mfem::MixedScalarWeakGradientIntegrator, and mfem::MixedScalarWeakDerivativeIntegrator.
Definition at line 302 of file bilininteg.hpp.
|
inlineprotectedvirtual |
Reimplemented in mfem::MixedScalarCurlIntegrator, mfem::MixedScalarDivergenceIntegrator, and mfem::MixedScalarDerivativeIntegrator.
Definition at line 307 of file bilininteg.hpp.
|
inlineprotectedvirtual |
Reimplemented in mfem::MixedScalarWeakCurlIntegrator, mfem::MixedScalarCurlIntegrator, mfem::MixedScalarWeakGradientIntegrator, mfem::MixedScalarDivergenceIntegrator, mfem::MixedScalarWeakDerivativeIntegrator, and mfem::MixedScalarDerivativeIntegrator.
Definition at line 290 of file bilininteg.hpp.
|
inlineprotectedvirtual |
Reimplemented in mfem::MixedScalarCurlIntegrator, and mfem::MixedScalarDivergenceIntegrator.
Definition at line 296 of file bilininteg.hpp.
|
inlineprotectedvirtual |
Reimplemented in mfem::MixedScalarWeakCurlIntegrator, mfem::MixedScalarCurlIntegrator, mfem::MixedScalarWeakGradientIntegrator, mfem::MixedScalarDivergenceIntegrator, mfem::MixedScalarWeakDerivativeIntegrator, and mfem::MixedScalarDerivativeIntegrator.
Definition at line 282 of file bilininteg.hpp.
|
protected |
Definition at line 312 of file bilininteg.hpp.
|
protected |
This parameter can be set by derived methods to enable single shape evaluation in case CalcTestShape() and CalcTrialShape() return the same result if given the same FiniteElement. The default is false.
Definition at line 277 of file bilininteg.hpp.