MFEM
v4.0
Finite element discretization library
|
#include <bilininteg.hpp>
Public Member Functions | |
MassIntegrator (const IntegrationRule *ir=NULL) | |
MassIntegrator (Coefficient &q, const IntegrationRule *ir=NULL) | |
Construct a mass integrator with coefficient q. More... | |
virtual void | AssembleElementMatrix (const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat) |
virtual void | AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat) |
virtual void | AssemblePA (const FiniteElementSpace &) |
Method defining partial assembly. More... | |
virtual void | AddMultPA (const Vector &, Vector &) 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 () |
![]() | |
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 () |
Static Public Member Functions | |
static const IntegrationRule & | GetRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans) |
Protected Attributes | |
Vector | shape |
Vector | te_shape |
Coefficient * | Q |
Vector | pa_data |
const DofToQuad * | maps |
Not owned. More... | |
const GeometricFactors * | geom |
Not owned. More... | |
int | dim |
int | ne |
int | nq |
int | dofs1D |
int | quad1D |
![]() | |
const IntegrationRule * | IntRule |
Additional Inherited Members | |
![]() | |
BilinearFormIntegrator (const IntegrationRule *ir=NULL) | |
![]() | |
NonlinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Class for local mass matrix assembling a(u,v) := (Q u, v)
Definition at line 1729 of file bilininteg.hpp.
|
inline |
Definition at line 1743 of file bilininteg.hpp.
|
inline |
Construct a mass integrator with coefficient q.
Definition at line 1747 of file bilininteg.hpp.
Method for partially assembled action.
Perform the action of integrator on the input x and add the result to the output y. Both x and y are E-vectors, i.e. they represent the element-wise discontinuous version of the FE space.
This method can be called only after the method AssemblePA() has been called.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 817 of file bilininteg_mass.cpp.
|
virtual |
Given a particular Finite Element computes the element mass matrix elmat.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 696 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 728 of file bilininteg.cpp.
|
virtual |
Method defining partial assembly.
The result of the partial assembly is stored internally so that it can be used later in the methods AddMultPA() and AddMultTransposePA().
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 24 of file bilininteg_mass.cpp.
|
static |
Definition at line 766 of file bilininteg.cpp.
|
protected |
Definition at line 1740 of file bilininteg.hpp.
|
protected |
Definition at line 1740 of file bilininteg.hpp.
|
protected |
Not owned.
Definition at line 1739 of file bilininteg.hpp.
|
protected |
Not owned.
Definition at line 1738 of file bilininteg.hpp.
|
protected |
Definition at line 1740 of file bilininteg.hpp.
|
protected |
Definition at line 1740 of file bilininteg.hpp.
|
protected |
Definition at line 1737 of file bilininteg.hpp.
|
protected |
Definition at line 1735 of file bilininteg.hpp.
|
protected |
Definition at line 1740 of file bilininteg.hpp.
|
protected |
Definition at line 1733 of file bilininteg.hpp.
|
protected |
Definition at line 1733 of file bilininteg.hpp.