MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::HyperelasticModel Class Referenceabstract

Abstract class for hyperelastic models. More...

#include <nonlininteg.hpp>

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

Public Member Functions

 HyperelasticModel ()
 
virtual ~HyperelasticModel ()
 
void SetTransformation (ElementTransformation &Ttr_)
 
virtual real_t EvalW (const DenseMatrix &Jpt) const =0
 Evaluate the strain energy density function, W = W(Jpt).
 
virtual void EvalP (const DenseMatrix &Jpt, DenseMatrix &P) const =0
 Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
 
virtual void AssembleH (const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const =0
 Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the local gradient matrix 'A'.
 

Protected Attributes

ElementTransformationTtr
 

Detailed Description

Abstract class for hyperelastic models.

Definition at line 215 of file nonlininteg.hpp.

Constructor & Destructor Documentation

◆ HyperelasticModel()

mfem::HyperelasticModel::HyperelasticModel ( )
inline

Definition at line 222 of file nonlininteg.hpp.

◆ ~HyperelasticModel()

virtual mfem::HyperelasticModel::~HyperelasticModel ( )
inlinevirtual

Definition at line 223 of file nonlininteg.hpp.

Member Function Documentation

◆ AssembleH()

virtual void mfem::HyperelasticModel::AssembleH ( const DenseMatrix & Jpt,
const DenseMatrix & DS,
const real_t weight,
DenseMatrix & A ) const
pure virtual

Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the local gradient matrix 'A'.

Parameters
[in]JptRepresents the target->physical transformation Jacobian matrix.
[in]DSGradient of the basis matrix (dof x dim).
[in]weightQuadrature weight coefficient for the point.
[in,out]ALocal gradient matrix where the contribution from this point will be added.

Computes weight * d(dW_dxi)_d(xj) at the current point, for all i and j, where x1 ... xn are the FE dofs. This function is usually defined using the matrix invariants and their derivatives.

Implemented in mfem::InverseHarmonicModel, mfem::NeoHookeanModel, mfem::TMOP_AMetric_011, mfem::TMOP_AMetric_014a, mfem::TMOP_AMetric_036, mfem::TMOP_AMetric_107a, mfem::TMOP_Combo_QualityMetric, mfem::TMOP_Metric_001, mfem::TMOP_Metric_002, mfem::TMOP_Metric_004, mfem::TMOP_Metric_007, mfem::TMOP_Metric_009, mfem::TMOP_Metric_014, mfem::TMOP_Metric_022, mfem::TMOP_Metric_050, mfem::TMOP_Metric_055, mfem::TMOP_Metric_056, mfem::TMOP_Metric_058, mfem::TMOP_Metric_077, mfem::TMOP_Metric_085, mfem::TMOP_Metric_098, mfem::TMOP_Metric_211, mfem::TMOP_Metric_252, mfem::TMOP_Metric_301, mfem::TMOP_Metric_302, mfem::TMOP_Metric_303, mfem::TMOP_Metric_304, mfem::TMOP_Metric_311, mfem::TMOP_Metric_313, mfem::TMOP_Metric_315, mfem::TMOP_Metric_316, mfem::TMOP_Metric_318, mfem::TMOP_Metric_321, mfem::TMOP_Metric_322, mfem::TMOP_Metric_323, mfem::TMOP_Metric_352, mfem::TMOP_Metric_360, mfem::TMOP_Metric_aspratio2D, mfem::TMOP_Metric_aspratio3D, mfem::TMOP_Metric_skew2D, mfem::TMOP_Metric_skew3D, mfem::TMOP_QualityMetric, and mfem::TMOP_WorstCaseUntangleOptimizer_Metric.

◆ EvalP()

◆ EvalW()

◆ SetTransformation()

void mfem::HyperelasticModel::SetTransformation ( ElementTransformation & Ttr_)
inline

A reference-element to target-element transformation that can be used to evaluate Coefficients.

Note
It is assumed that Ttr_.SetIntPoint() is already called for the point of interest.

Definition at line 229 of file nonlininteg.hpp.

Member Data Documentation

◆ Ttr

ElementTransformation* mfem::HyperelasticModel::Ttr
protected

Reference-element to target-element transformation.

Definition at line 218 of file nonlininteg.hpp.


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