MFEM
v3.3.2
Finite element discretization library
|
#include <bilininteg.hpp>
Inherits mfem::BilinearFormIntegrator.
Inherited by mfem::MixedDirectionalDerivativeIntegrator, mfem::MixedDivGradIntegrator, mfem::MixedDotProductIntegrator, mfem::MixedGradDivIntegrator, mfem::MixedScalarCrossCurlIntegrator, mfem::MixedScalarCrossGradIntegrator, mfem::MixedScalarCrossProductIntegrator, mfem::MixedScalarWeakCrossProductIntegrator, mfem::MixedScalarWeakCurlCrossIntegrator, mfem::MixedScalarWeakDivergenceIntegrator, mfem::MixedVectorDivergenceIntegrator, mfem::MixedVectorProductIntegrator, and mfem::MixedWeakGradDotIntegrator.
Public Member Functions | |
virtual void | AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat) |
Public Member Functions inherited from mfem::BilinearFormIntegrator | |
virtual void | AssembleElementMatrix (const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat) |
Given a particular Finite Element computes the element matrix elmat. 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 double | ComputeFluxEnergy (const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL) |
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 | |
MixedScalarVectorIntegrator (VectorCoefficient &vq, bool _transpose=false, bool _cross_2d=false) | |
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 | CalcVShape (const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape) |
virtual void | CalcShape (const FiniteElement &scalar_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) | |
Additional Inherited Members | |
Protected Attributes inherited from mfem::NonlinearFormIntegrator | |
const IntegrationRule * | IntRule |
An abstract class for integrating the product of a scalar basis function and the inner product of a vector basis function with a vector coefficient. In 2D the inner product can be replaced with a cross product.
Definition at line 323 of file bilininteg.hpp.
|
inlineprotected |
Definition at line 334 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 297 of file bilininteg.cpp.
|
inlineprotectedvirtual |
Reimplemented in mfem::MixedDivGradIntegrator, mfem::MixedGradDivIntegrator, mfem::MixedScalarWeakCrossProductIntegrator, mfem::MixedScalarCrossCurlIntegrator, mfem::MixedScalarWeakCurlCrossIntegrator, mfem::MixedWeakGradDotIntegrator, and mfem::MixedVectorDivergenceIntegrator.
Definition at line 378 of file bilininteg.hpp.
|
inlineprotectedvirtual |
Reimplemented in mfem::MixedScalarWeakDivergenceIntegrator, mfem::MixedDivGradIntegrator, mfem::MixedGradDivIntegrator, mfem::MixedDirectionalDerivativeIntegrator, and mfem::MixedScalarCrossGradIntegrator.
Definition at line 373 of file bilininteg.hpp.
|
inlineprotectedvirtual |
Reimplemented in mfem::MixedScalarWeakDivergenceIntegrator, mfem::MixedDivGradIntegrator, mfem::MixedGradDivIntegrator, mfem::MixedDirectionalDerivativeIntegrator, mfem::MixedScalarWeakCrossProductIntegrator, mfem::MixedScalarCrossProductIntegrator, mfem::MixedScalarCrossGradIntegrator, mfem::MixedScalarCrossCurlIntegrator, mfem::MixedScalarWeakCurlCrossIntegrator, mfem::MixedWeakGradDotIntegrator, mfem::MixedDotProductIntegrator, and mfem::MixedVectorDivergenceIntegrator.
Definition at line 351 of file bilininteg.hpp.
|
inlineprotectedvirtual |
Reimplemented in mfem::MixedVectorDivergenceIntegrator.
Definition at line 367 of file bilininteg.hpp.
|
inlineprotectedvirtual |
Reimplemented in mfem::MixedScalarWeakDivergenceIntegrator, mfem::MixedDivGradIntegrator, mfem::MixedGradDivIntegrator, mfem::MixedDirectionalDerivativeIntegrator, mfem::MixedScalarWeakCrossProductIntegrator, mfem::MixedScalarCrossProductIntegrator, mfem::MixedScalarCrossGradIntegrator, mfem::MixedScalarCrossCurlIntegrator, mfem::MixedScalarWeakCurlCrossIntegrator, mfem::MixedWeakGradDotIntegrator, mfem::MixedDotProductIntegrator, and mfem::MixedVectorDivergenceIntegrator.
Definition at line 338 of file bilininteg.hpp.