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

#include <nonlininteg.hpp>

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

Public Member Functions

virtual real_t EvalW (const DenseMatrix &J) const
 Evaluate the strain energy density function, W = W(Jpt).
 
virtual void EvalP (const DenseMatrix &J, DenseMatrix &P) const
 Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
 
virtual void AssembleH (const DenseMatrix &J, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const
 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 Attributes

DenseMatrix Z
 
DenseMatrix S
 
DenseMatrix G
 
DenseMatrix C
 
- Protected Attributes inherited from mfem::HyperelasticModel
ElementTransformationTtr
 

Detailed Description

Inverse-harmonic hyperelastic model with a strain energy density function given by the formula: W(J) = (1/2) det(J) Tr((J J^t)^{-1}) where J is the deformation gradient.

Definition at line 263 of file nonlininteg.hpp.

Member Function Documentation

◆ AssembleH()

void mfem::InverseHarmonicModel::AssembleH ( const DenseMatrix & Jpt,
const DenseMatrix & DS,
const real_t 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'.

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 196 of file nonlininteg.cpp.

◆ EvalP()

void mfem::InverseHarmonicModel::EvalP ( const DenseMatrix & Jpt,
DenseMatrix & P ) const
virtual

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 177 of file nonlininteg.cpp.

◆ EvalW()

real_t mfem::InverseHarmonicModel::EvalW ( const DenseMatrix & Jpt) const
virtual

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

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

Implements mfem::HyperelasticModel.

Definition at line 170 of file nonlininteg.cpp.

Member Data Documentation

◆ C

DenseMatrix mfem::InverseHarmonicModel::C
protected

Definition at line 267 of file nonlininteg.hpp.

◆ G

DenseMatrix mfem::InverseHarmonicModel::G
mutableprotected

Definition at line 267 of file nonlininteg.hpp.

◆ S

DenseMatrix mfem::InverseHarmonicModel::S
protected

Definition at line 266 of file nonlininteg.hpp.

◆ Z

DenseMatrix mfem::InverseHarmonicModel::Z
mutableprotected

Definition at line 266 of file nonlininteg.hpp.


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