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

#include <nonlininteg.hpp>

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

Public Member Functions

 NeoHookeanModel (double _mu, double _K, double _g=1.0)
 
 NeoHookeanModel (Coefficient &_mu, Coefficient &_K, Coefficient *_g=NULL)
 
virtual double EvalW (const DenseMatrix &J) const
 Evaluate the strain energy density function, W=W(J). More...
 
virtual void EvalP (const DenseMatrix &J, DenseMatrix &P) const
 Evaluate the 1st Piola-Kirchhoff stress tensor, P=P(J). More...
 
virtual void AssembleH (const DenseMatrix &J, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
 
- Public Member Functions inherited from mfem::HyperelasticModel
 HyperelasticModel ()
 
void SetTransformation (ElementTransformation &_T)
 An element transformation that can be used to evaluate coefficients. More...
 
virtual ~HyperelasticModel ()
 

Protected Member Functions

void EvalCoeffs () const
 

Protected Attributes

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

Detailed Description

Neo-Hookean hyperelastic model with a strain energy density function given by the formula: $(\mu/2)(\bar{I}_1 - dim) + (K/2)(det(J)/g - 1)^2$ where J is the deformation gradient and $\bar{I}_1 = (det(J))^{-2/dim} Tr(J J^t)$. The parameters $\mu$ and K are the shear and bulk moduli, respectively, and g is a reference volumetric scaling.

Definition at line 101 of file nonlininteg.hpp.

Constructor & Destructor Documentation

mfem::NeoHookeanModel::NeoHookeanModel ( double  _mu,
double  _K,
double  _g = 1.0 
)
inline

Definition at line 114 of file nonlininteg.hpp.

mfem::NeoHookeanModel::NeoHookeanModel ( Coefficient _mu,
Coefficient _K,
Coefficient _g = NULL 
)
inline

Definition at line 117 of file nonlininteg.hpp.

Member Function Documentation

void mfem::NeoHookeanModel::AssembleH ( const DenseMatrix J,
const DenseMatrix DS,
const double  weight,
DenseMatrix A 
) const
virtual

Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the local gradient matrix 'A'. 'DS' is the gradient of the basis matrix (dof x dim), and 'weight' is the quadrature weight.

Implements mfem::HyperelasticModel.

Definition at line 182 of file nonlininteg.cpp.

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

Definition at line 135 of file nonlininteg.cpp.

void mfem::NeoHookeanModel::EvalP ( const DenseMatrix J,
DenseMatrix P 
) const
virtual

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

Implements mfem::HyperelasticModel.

Definition at line 161 of file nonlininteg.cpp.

double mfem::NeoHookeanModel::EvalW ( const DenseMatrix J) const
virtual

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

Implements mfem::HyperelasticModel.

Definition at line 145 of file nonlininteg.cpp.

Member Data Documentation

DenseMatrix mfem::NeoHookeanModel::C
mutableprotected

Definition at line 109 of file nonlininteg.hpp.

Coefficient * mfem::NeoHookeanModel::c_g
protected

Definition at line 105 of file nonlininteg.hpp.

Coefficient * mfem::NeoHookeanModel::c_K
protected

Definition at line 105 of file nonlininteg.hpp.

Coefficient* mfem::NeoHookeanModel::c_mu
protected

Definition at line 105 of file nonlininteg.hpp.

double mfem::NeoHookeanModel::g
mutableprotected

Definition at line 104 of file nonlininteg.hpp.

DenseMatrix mfem::NeoHookeanModel::G
mutableprotected

Definition at line 109 of file nonlininteg.hpp.

bool mfem::NeoHookeanModel::have_coeffs
protected

Definition at line 106 of file nonlininteg.hpp.

double mfem::NeoHookeanModel::K
mutableprotected

Definition at line 104 of file nonlininteg.hpp.

double mfem::NeoHookeanModel::mu
mutableprotected

Definition at line 104 of file nonlininteg.hpp.

DenseMatrix mfem::NeoHookeanModel::Z
mutableprotected

Definition at line 108 of file nonlininteg.hpp.


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