MFEM
v4.0
Finite element discretization library
|
#include <bilininteg.hpp>
Public Member Functions | |
DiffusionIntegrator () | |
Construct a diffusion integrator with coefficient Q = 1. More... | |
DiffusionIntegrator (Coefficient &q) | |
Construct a diffusion integrator with a scalar coefficient q. More... | |
DiffusionIntegrator (MatrixCoefficient &q) | |
Construct a diffusion integrator with a matrix 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 | AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect) |
Perform the local action of the BilinearFormIntegrator. 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 void | AssemblePA (const FiniteElementSpace &) |
Method defining partial assembly. More... | |
virtual void | AddMultPA (const Vector &, Vector &) const |
Method for partially assembled action. More... | |
Public Member Functions inherited from mfem::BilinearFormIntegrator | |
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 | 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 | ~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 () |
Static Public Member Functions | |
static const IntegrationRule & | GetRule (const FiniteElement &trial_fe, const FiniteElement &test_fe) |
Protected Attributes | |
Coefficient * | Q |
MatrixCoefficient * | MQ |
Protected Attributes inherited from mfem::NonlinearFormIntegrator | |
const IntegrationRule * | IntRule |
Additional Inherited Members | |
Protected Member Functions inherited from mfem::BilinearFormIntegrator | |
BilinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Protected Member Functions inherited from mfem::NonlinearFormIntegrator | |
NonlinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Class for integrating the bilinear form a(u,v) := (Q grad u, grad v) where Q can be a scalar or a matrix coefficient.
Definition at line 1663 of file bilininteg.hpp.
|
inline |
Construct a diffusion integrator with coefficient Q = 1.
Definition at line 1684 of file bilininteg.hpp.
|
inline |
Construct a diffusion integrator with a scalar coefficient q.
Definition at line 1687 of file bilininteg.hpp.
|
inline |
Construct a diffusion integrator with a matrix coefficient q.
Definition at line 1691 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 1122 of file bilininteg_diffusion.cpp.
|
virtual |
Given a particular Finite Element computes the element stiffness matrix elmat.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 384 of file bilininteg.cpp.
|
virtual |
Given a trial and test Finite Element computes the element stiffness matrix elmat.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 434 of file bilininteg.cpp.
|
virtual |
Perform the local action of the BilinearFormIntegrator.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 493 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 192 of file bilininteg_diffusion.cpp.
|
virtual |
Virtual method required for Zienkiewicz-Zhu type error estimators.
The purpose of the method is to compute a local "flux" finite element function given a local finite element solution. The "flux" function has to be computed in terms of its coefficients (represented by the Vector flux) which multiply the basis functions defined by the FiniteElement fluxelem. Typically, the "flux" function will have more than one component and consequently flux should be store the coefficients of all components: first all coefficient for component 0, then all coefficients for component 1, etc. What the "flux" function represents depends on the specific integrator. For example, in the case of DiffusionIntegrator, the flux is the gradient of the solution multiplied by the diffusion coefficient.
[in] | el | FiniteElement of the solution. |
[in] | Trans | The ElementTransformation describing the physical position of the mesh element. |
[in] | u | Solution coefficients representing the expansion of the solution function in the basis of el. |
[in] | fluxelem | FiniteElement of the "flux". |
[out] | flux | "Flux" coefficients representing the expansion of the "flux" function in the basis of fluxelem. The size of flux as a Vector has to be set by this method, e.g. using Vector::SetSize(). |
[in] | with_coef | If zero (the default value is 1) the implementation of the method may choose not to scale the "flux" function by any coefficients describing the integrator. |
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 549 of file bilininteg.cpp.
|
virtual |
Virtual method required for Zienkiewicz-Zhu type error estimators.
The purpose of this method is to compute a local number that measures the energy of a given "flux" function (see ComputeElementFlux() for a description of the "flux" function). Typically, the energy of a "flux" function should be equal to a_local(u,u), if the "flux" is defined from a solution u; here a_local(.,.) denotes the element-local bilinear form represented by the integrator.
[in] | fluxelem | FiniteElement of the "flux". |
[in] | Trans | The ElementTransformation describing the physical position of the mesh element. |
[in] | flux | "Flux" coefficients representing the expansion of the "flux" function in the basis of fluxelem. |
[out] | d_energy | If not NULL, the given Vector should be set to represent directional energy split that can be used for anisotropic error estimation. |
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 607 of file bilininteg.cpp.
|
static |
Definition at line 673 of file bilininteg.cpp.
|
protected |
Definition at line 1667 of file bilininteg.hpp.
|
protected |
Definition at line 1666 of file bilininteg.hpp.