MFEM
v4.0
Finite element discretization library
|
#include <bilininteg.hpp>
Public Member Functions | |
ElasticityIntegrator (Coefficient &l, Coefficient &m) | |
ElasticityIntegrator (Coefficient &m, double q_l, double q_m) | |
virtual void | AssembleElementMatrix (const FiniteElement &, ElementTransformation &, DenseMatrix &) |
Given a particular Finite Element computes the element matrix elmat. 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 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 | AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat) |
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 | ~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 () |
Protected Attributes | |
double | q_lambda |
double | q_mu |
Coefficient * | lambda |
Coefficient * | mu |
![]() | |
const IntegrationRule * | IntRule |
Additional Inherited Members | |
![]() | |
BilinearFormIntegrator (const IntegrationRule *ir=NULL) | |
![]() | |
NonlinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Integrator for the linear elasticity form: a(u,v) = (lambda div(u), div(v)) + (2 mu e(u), e(v)), where e(v) = (1/2) (grad(v) + grad(v)^T). This is a 'Vector' integrator, i.e. defined for FE spaces using multiple copies of a scalar FE space.
Definition at line 2158 of file bilininteg.hpp.
|
inline |
Definition at line 2172 of file bilininteg.hpp.
|
inline |
With this constructor lambda = q_l * m and mu = q_m * m; if dim * q_l + 2 * q_m = 0 then trace(sigma) = 0.
Definition at line 2176 of file bilininteg.hpp.
|
virtual |
Given a particular Finite Element computes the element matrix elmat.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 2144 of file bilininteg.cpp.
|
virtual |
Compute the stress corresponding to the local displacement u and interpolate it at the nodes of the given fluxelem. Only the symmetric part of the stress is stored, so that the size of flux is equal to the number of DOFs in fluxelem times dim*(dim+1)/2. In 2D, the order of the stress components is: s_xx, s_yy, s_xy. In 3D, it is: s_xx, s_yy, s_zz, s_xy, s_xz, s_yz. In other words, flux is the local vector for a FE space with dim*(dim+1)/2 vector components, based on the finite element fluxelem.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 2226 of file bilininteg.cpp.
|
virtual |
Compute the element energy (integral of the strain energy density) corresponding to the stress represented by flux which is a vector of coefficients multiplying the basis functions defined by fluxelem. In other words, flux is the local vector for a FE space with dim*(dim+1)/2 vector components, based on the finite element fluxelem. The number of components, dim*(dim+1)/2 is such that it represents the symmetric part of the (symmetric) stress tensor. The order of the components is: s_xx, s_yy, s_xy in 2D, and s_xx, s_yy, s_zz, s_xy, s_xz, s_yz in 3D.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 2302 of file bilininteg.cpp.
|
protected |
Definition at line 2162 of file bilininteg.hpp.
|
protected |
Definition at line 2162 of file bilininteg.hpp.
|
protected |
Definition at line 2161 of file bilininteg.hpp.
|
protected |
Definition at line 2161 of file bilininteg.hpp.