MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::HyperelasticNLFIntegrator Class Reference

#include <nonlininteg.hpp>

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

Public Member Functions

 HyperelasticNLFIntegrator (HyperelasticModel *m)
 
virtual real_t GetElementEnergy (const FiniteElement &el, ElementTransformation &Ttr, const Vector &elfun)
 Computes the integral of W(Jacobian(Trt)) over a target zone.
 
virtual void AssembleElementVector (const FiniteElement &el, ElementTransformation &Ttr, const Vector &elfun, Vector &elvect)
 Perform the local action of the NonlinearFormIntegrator.
 
virtual void AssembleElementGrad (const FiniteElement &el, ElementTransformation &Ttr, const Vector &elfun, DenseMatrix &elmat)
 Assemble the local gradient matrix.
 
- 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).
 
void SetIntegrationMode (Mode m)
 
void SetNURBSPatchIntRule (NURBSMeshRules *pr)
 For patchwise integration, SetNURBSPatchIntRule must be called.
 
bool HasNURBSPatchIntRule () const
 
bool Patchwise () const
 
void SetIntegrationRule (const IntegrationRule &ir)
 Prescribe a fixed IntegrationRule to use.
 
void SetPAMemoryType (MemoryType mt)
 
const IntegrationRuleGetIntegrationRule () const
 Get the integration rule of the integrator (possibly NULL).
 
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.
 
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.
 
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly.
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
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.
 
virtual real_t GetLocalStateEnergyPA (const Vector &x) const
 Compute the local (to the MPI rank) energy with partial assembly.
 
virtual void AddMultPA (const Vector &x, Vector &y) const
 Method for partially assembled action.
 
virtual void AddMultGradPA (const Vector &x, Vector &y) const
 Method for partially assembled gradient action.
 
virtual void AssembleGradDiagonalPA (Vector &diag) const
 Method for computing the diagonal of the gradient with partial assembly.
 
virtual bool SupportsCeed () const
 Indicates whether this integrator can use a Ceed backend.
 
virtual void AssembleMF (const FiniteElementSpace &fes)
 Method defining fully unassembled operator.
 
virtual void AddMultMF (const Vector &x, Vector &y) const
 
ceed::OperatorGetCeedOp ()
 
virtual ~NonlinearFormIntegrator ()
 

Additional Inherited Members

- Public Types inherited from mfem::NonlinearFormIntegrator
enum  Mode { ELEMENTWISE = 0 , PATCHWISE = 1 , PATCHWISE_REDUCED = 2 }
 
- Protected Member Functions inherited from mfem::NonlinearFormIntegrator
 NonlinearFormIntegrator (const IntegrationRule *ir=NULL)
 
- Protected Attributes inherited from mfem::NonlinearFormIntegrator
const IntegrationRuleIntRule
 
Mode integrationMode = Mode::ELEMENTWISE
 
NURBSMeshRulespatchRules = nullptr
 
ceed::OperatorceedOp
 
MemoryType pa_mt = MemoryType::DEFAULT
 

Detailed Description

Hyperelastic integrator for any given HyperelasticModel.

Represents \( \int W(Jpt) dx \) over a target zone, where W is the model's strain energy density function, and Jpt is the Jacobian of the target->physical coordinates transformation. The target configuration is given by the current mesh at the time of the evaluation of the integrator.

Definition at line 320 of file nonlininteg.hpp.

Constructor & Destructor Documentation

◆ HyperelasticNLFIntegrator()

mfem::HyperelasticNLFIntegrator::HyperelasticNLFIntegrator ( HyperelasticModel * m)
inline
Parameters
[in]mHyperelasticModel that will be integrated.

Definition at line 339 of file nonlininteg.hpp.

Member Function Documentation

◆ AssembleElementGrad()

void mfem::HyperelasticNLFIntegrator::AssembleElementGrad ( const FiniteElement & el,
ElementTransformation & Tr,
const Vector & elfun,
DenseMatrix & elmat )
virtual

Assemble the local gradient matrix.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 464 of file nonlininteg.cpp.

◆ AssembleElementVector()

void mfem::HyperelasticNLFIntegrator::AssembleElementVector ( const FiniteElement & el,
ElementTransformation & Tr,
const Vector & elfun,
Vector & elvect )
virtual

Perform the local action of the NonlinearFormIntegrator.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 424 of file nonlininteg.cpp.

◆ GetElementEnergy()

real_t mfem::HyperelasticNLFIntegrator::GetElementEnergy ( const FiniteElement & el,
ElementTransformation & Ttr,
const Vector & elfun )
virtual

Computes the integral of W(Jacobian(Trt)) over a target zone.

Parameters
[in]elType of FiniteElement.
[in]TtrRepresents ref->target coordinates transformation.
[in]elfunPhysical coordinates of the zone.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 387 of file nonlininteg.cpp.


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