MFEM v4.8.0
Finite element discretization library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mfem::NeoHookeanModel Class Reference

#include <nonlininteg.hpp>

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

Public Member Functions

 NeoHookeanModel (real_t mu_, real_t K_, real_t g_=1.0)
 
 NeoHookeanModel (Coefficient &mu_, Coefficient &K_, Coefficient *g_=NULL)
 
real_t EvalW (const DenseMatrix &J) const override
 Evaluate the strain energy density function, W = W(Jpt).
 
void EvalP (const DenseMatrix &J, DenseMatrix &P) const override
 Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
 
void AssembleH (const DenseMatrix &J, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
 Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the local gradient matrix 'A'.
 
- Public Member Functions inherited from mfem::HyperelasticModel
 HyperelasticModel ()
 
virtual ~HyperelasticModel ()
 
void SetTransformation (ElementTransformation &Ttr_)
 

Protected Member Functions

void EvalCoeffs () const
 

Protected Attributes

real_t mu
 
real_t K
 
real_t g
 
Coefficientc_mu
 
Coefficientc_K
 
Coefficientc_g
 
bool have_coeffs
 
DenseMatrix Z
 
DenseMatrix G
 
DenseMatrix C
 
- Protected Attributes inherited from mfem::HyperelasticModel
ElementTransformationTtr
 

Detailed Description

Neo-Hookean hyperelastic model with a strain energy density function given by the formula: (μ/2)(I¯1dim)+(K/2)(det(J)/g1)2 where J is the deformation gradient and

I¯1=(det(J))2/dimTr(JJt)

. The parameters μ and K are the shear and bulk moduli, respectively, and g is a reference volumetric scaling.

Definition at line 269 of file nonlininteg.hpp.

Constructor & Destructor Documentation

◆ NeoHookeanModel() [1/2]

mfem::NeoHookeanModel::NeoHookeanModel ( real_t mu_,
real_t K_,
real_t g_ = 1.0 )
inline

Definition at line 282 of file nonlininteg.hpp.

◆ NeoHookeanModel() [2/2]

mfem::NeoHookeanModel::NeoHookeanModel ( Coefficient & mu_,
Coefficient & K_,
Coefficient * g_ = NULL )
inline

Definition at line 285 of file nonlininteg.hpp.

Member Function Documentation

◆ AssembleH()

void mfem::NeoHookeanModel::AssembleH ( const DenseMatrix & Jpt,
const DenseMatrix & DS,
const real_t weight,
DenseMatrix & A ) const
overridevirtual

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.

Implements mfem::HyperelasticModel.

Definition at line 318 of file nonlininteg.cpp.

◆ EvalCoeffs()

void mfem::NeoHookeanModel::EvalCoeffs ( ) const
inlineprotected

Definition at line 271 of file nonlininteg.cpp.

◆ EvalP()

void mfem::NeoHookeanModel::EvalP ( const DenseMatrix & Jpt,
DenseMatrix & P ) const
overridevirtual

Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).

Parameters
[in]JptRepresents the target->physical transformation Jacobian matrix.
[out]PThe evaluated 1st Piola-Kirchhoff stress tensor.

Implements mfem::HyperelasticModel.

Definition at line 297 of file nonlininteg.cpp.

◆ EvalW()

real_t mfem::NeoHookeanModel::EvalW ( const DenseMatrix & Jpt) const
overridevirtual

Evaluate the strain energy density function, W = W(Jpt).

Parameters
[in]JptRepresents the target->physical transformation Jacobian matrix.

Implements mfem::HyperelasticModel.

Definition at line 281 of file nonlininteg.cpp.

Member Data Documentation

◆ C

DenseMatrix mfem::NeoHookeanModel::C
protected

Definition at line 277 of file nonlininteg.hpp.

◆ c_g

Coefficient * mfem::NeoHookeanModel::c_g
protected

Definition at line 273 of file nonlininteg.hpp.

◆ c_K

Coefficient * mfem::NeoHookeanModel::c_K
protected

Definition at line 273 of file nonlininteg.hpp.

◆ c_mu

Coefficient* mfem::NeoHookeanModel::c_mu
protected

Definition at line 273 of file nonlininteg.hpp.

◆ G

DenseMatrix mfem::NeoHookeanModel::G
mutableprotected

Definition at line 277 of file nonlininteg.hpp.

◆ g

real_t mfem::NeoHookeanModel::g
protected

Definition at line 272 of file nonlininteg.hpp.

◆ have_coeffs

bool mfem::NeoHookeanModel::have_coeffs
protected

Definition at line 274 of file nonlininteg.hpp.

◆ K

real_t mfem::NeoHookeanModel::K
protected

Definition at line 272 of file nonlininteg.hpp.

◆ mu

real_t mfem::NeoHookeanModel::mu
mutableprotected

Definition at line 272 of file nonlininteg.hpp.

◆ Z

DenseMatrix mfem::NeoHookeanModel::Z
mutableprotected

Definition at line 276 of file nonlininteg.hpp.


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