MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Protected Attributes | List of all members
mfem::ElasticityIntegrator Class Reference

#include <bilininteg.hpp>

Inheritance diagram for mfem::ElasticityIntegrator:
[legend]
Collaboration diagram for mfem::ElasticityIntegrator:
[legend]

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, bool with_coef=true)
 
virtual double ComputeFluxEnergy (const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
 
- Public Member Functions inherited from mfem::BilinearFormIntegrator
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
virtual void AssemblePAInteriorFaces (const FiniteElementSpace &fes)
 
virtual void AssemblePABoundaryFaces (const FiniteElementSpace &fes)
 
virtual void AssembleDiagonalPA (Vector &diag)
 Assemble diagonal and add it to Vector diag. More...
 
virtual void AssembleDiagonalPA_ADAt (const Vector &D, Vector &diag)
 Assemble diagonal of ADA^T (A is this integrator) and add it to diag. 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 AssembleEA (const FiniteElementSpace &fes, Vector &emat, const bool add=true)
 Method defining element assembly. More...
 
virtual void AssembleMF (const FiniteElementSpace &fes)
 Method defining matrix-free assembly. More...
 
virtual void AddMultMF (const Vector &x, Vector &y) const
 
virtual void AddMultTransposeMF (const Vector &x, Vector &y) const
 
virtual void AssembleDiagonalMF (Vector &diag)
 Assemble diagonal and add it to Vector diag. More...
 
virtual void AssembleEAInteriorFaces (const FiniteElementSpace &fes, Vector &ea_data_int, Vector &ea_data_ext, const bool add=true)
 
virtual void AssembleEABoundaryFaces (const FiniteElementSpace &fes, Vector &ea_data_bdr, const bool add=true)
 
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. Note that the default implementation in the base class is general but not efficient. More...
 
virtual void AssembleFaceVector (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect)
 Perform the local action of the BilinearFormIntegrator resulting from a face integral term. Note that the default implementation in the base class is general but not efficient. 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 ()
 
- Public Member Functions inherited from mfem::NonlinearFormIntegrator
virtual 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 &ir)
 Prescribe a fixed IntegrationRule to use. More...
 
void SetPAMemoryType (MemoryType mt)
 
const IntegrationRuleGetIntegrationRule () const
 Get the integration rule of the integrator (possibly NULL). More...
 
virtual double GetElementEnergy (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
 Compute the local energy. More...
 
virtual void AssembleGradPA (const Vector &x, const FiniteElementSpace &fes)
 Prepare the integrator for partial assembly (PA) gradient evaluations on the given FE space fes at the state x. More...
 
virtual double GetLocalStateEnergyPA (const Vector &x) const
 Compute the local (to the MPI rank) energy with partial assembly. More...
 
virtual void AddMultGradPA (const Vector &x, Vector &y) const
 Method for partially assembled gradient action. More...
 
virtual void AssembleGradDiagonalPA (Vector &diag) const
 Method for computing the diagonal of the gradient with partial assembly. More...
 
virtual bool SupportsCeed () const
 Indicates whether this integrator can use a Ceed backend. More...
 
ceed::Operator & GetCeedOp ()
 
virtual ~NonlinearFormIntegrator ()
 

Protected Attributes

double q_lambda
 
double q_mu
 
Coefficientlambda
 
Coefficientmu
 
- Protected Attributes inherited from mfem::NonlinearFormIntegrator
const IntegrationRuleIntRule
 
ceed::Operator * ceedOp
 
MemoryType pa_mt = MemoryType::DEFAULT
 

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)
 

Detailed Description

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 2846 of file bilininteg.hpp.

Constructor & Destructor Documentation

mfem::ElasticityIntegrator::ElasticityIntegrator ( Coefficient l,
Coefficient m 
)
inline

Definition at line 2860 of file bilininteg.hpp.

mfem::ElasticityIntegrator::ElasticityIntegrator ( Coefficient m,
double  q_l,
double  q_m 
)
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 2864 of file bilininteg.hpp.

Member Function Documentation

void mfem::ElasticityIntegrator::AssembleElementMatrix ( const FiniteElement el,
ElementTransformation Trans,
DenseMatrix elmat 
)
virtual

Given a particular Finite Element computes the element matrix elmat.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 2760 of file bilininteg.cpp.

void mfem::ElasticityIntegrator::ComputeElementFlux ( const FiniteElement el,
ElementTransformation Trans,
Vector u,
const FiniteElement fluxelem,
Vector flux,
bool  with_coef = true 
)
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 2842 of file bilininteg.cpp.

double mfem::ElasticityIntegrator::ComputeFluxEnergy ( const FiniteElement fluxelem,
ElementTransformation Trans,
Vector flux,
Vector d_energy = NULL 
)
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 2918 of file bilininteg.cpp.

Member Data Documentation

Coefficient* mfem::ElasticityIntegrator::lambda
protected

Definition at line 2850 of file bilininteg.hpp.

Coefficient * mfem::ElasticityIntegrator::mu
protected

Definition at line 2850 of file bilininteg.hpp.

double mfem::ElasticityIntegrator::q_lambda
protected

Definition at line 2849 of file bilininteg.hpp.

double mfem::ElasticityIntegrator::q_mu
protected

Definition at line 2849 of file bilininteg.hpp.


The documentation for this class was generated from the following files: