chrono::fea::ChMaterial3DHyperelastic Class Referenceabstract

Description

Class for hyperelastic materials, that is materials with a constitutive equation of the type S = f(C) , where S is the 2nd Piola-Kirchhoff stress and C is the right Cauchy-Green deformation tensor C=F^T*F (in material space) as in E = 1/2(C -I) where E the Green Lagrange strain tensor.

Or in general, materials for which the stress–strain relationship derives from a strain energy density function.

#include <ChMaterial3DHyperelastic.h>

Inheritance diagram for chrono::fea::ChMaterial3DHyperelastic:
Collaboration diagram for chrono::fea::ChMaterial3DHyperelastic:

Public Member Functions

virtual void ComputeStress (ChStressTensor<> &S_stress, const ChMatrix33d &F, const ChMatrix33d *l, ChFieldData *data_per_point, ChElementData *data_per_element) override
 Implement interface to lower level stress material. More...
 
virtual void ComputeTangentModulus (ChMatrixNM< double, 6, 6 > &C, const ChMatrix33d &F, const ChMatrix33d *l, ChFieldData *data_per_point, ChElementData *data_per_element) override
 Implement interface to lower level stress material. More...
 
virtual bool IsSpatialVelocityGradientNeeded () const override
 Hyperelastic materials do not need info on the spatial velocity gradient l=\nabla_x v , Returning false from this means that the ChDomainXXYY can know that the computation of the "l" parameter could be skipped left to null when calling ComputeStress(...)
 
virtual void ComputeElasticStress (ChStressTensor<> &S_stress, const ChMatrix33d &C_def)=0
 Compute elastic stress from the right Cauchy-Green deformation tensor C_def = F^T*F. More...
 
virtual void ComputeElasticTangentModulus (ChMatrixNM< double, 6, 6 > &C_tangent, const ChMatrix33d &C_def)
 Computes the tangent modulus for a given deformation. More...
 
- Public Member Functions inherited from chrono::fea::ChMaterial3DStress
virtual void ComputeUpdateEndStep (ChFieldData *data_per_point, ChElementData *data_per_element, const double time)
 Update your own auxiliary data, if any, at the end of time step (ex for plasticity). More...
 
- Public Member Functions inherited from chrono::fea::ChMaterial3DDensity
 ChMaterial3DDensity (double density=1000)
 
virtual void SetDensity (double density)
 Set the density of the material, in kg/m^2.
 
virtual double GetDensity () const
 Get the density of the material, in kg/m^2.
 
- Public Member Functions inherited from chrono::fea::ChMaterial
virtual std::unique_ptr< ChFieldDataCreateMaterialPointData () const
 Implement this if the material needs custom data per material point, returning a std::make_unique<ChFieldDataCustom>() where ChFieldDataCustom is your custom class with aditional states/properties per-point.
 

Additional Inherited Members

- Protected Attributes inherited from chrono::fea::ChMaterial3DDensity
double m_density
 

Member Function Documentation

◆ ComputeElasticStress()

virtual void chrono::fea::ChMaterial3DHyperelastic::ComputeElasticStress ( ChStressTensor<> &  S_stress,
const ChMatrix33d C_def 
)
pure virtual

Compute elastic stress from the right Cauchy-Green deformation tensor C_def = F^T*F.

Return stress as 2nd Piola-Kirchhoff S tensor, in Voigt notation. Children materials MUST implement this.

Implemented in chrono::fea::ChMaterial3DStressStVenant, chrono::fea::ChMaterial3DStressNeoHookean, and chrono::fea::ChMaterial3DStressOgden.

◆ ComputeElasticTangentModulus()

virtual void chrono::fea::ChMaterial3DHyperelastic::ComputeElasticTangentModulus ( ChMatrixNM< double, 6, 6 > &  C_tangent,
const ChMatrix33d C_def 
)
inlinevirtual

Computes the tangent modulus for a given deformation.

It takes the right Cauchy-Green deformation tensor C_def, and returns the tangent modulus in 6x6 matrix C_tangent, for dS = C_tangent * dE where S is 2nd Piola-Kirchhoff stress, and E is Green Lagrange strain. A fallback implementation is provided here, that computes C via numerical differentiation by calling 6 times the ComputeElasticStress() function. However, if you are able to compute a custom analytical expression of C_tangent, and the analytical expression is fast to compute, then you SHOULD override this.

Reimplemented in chrono::fea::ChMaterial3DStressNeoHookean, and chrono::fea::ChMaterial3DStressStVenant.

◆ ComputeStress()

virtual void chrono::fea::ChMaterial3DHyperelastic::ComputeStress ( ChStressTensor<> &  S_stress,
const ChMatrix33d F,
const ChMatrix33d l,
ChFieldData data_per_point,
ChElementData data_per_element 
)
inlineoverridevirtual

Implement interface to lower level stress material.

Parameters
S_stressoutput stress, PK2
Fcurrent deformation gradient tensor
lcurrent spatial velocity gradient (might be nullptr if IsSpatialVelocityGradientNeeded() is false)
data_per_pointpointer to auxiliary data (ex states), if any, per quadrature point
data_per_elementpointer to auxiliary data (ex states), if any, per element

Implements chrono::fea::ChMaterial3DStress.

◆ ComputeTangentModulus()

virtual void chrono::fea::ChMaterial3DHyperelastic::ComputeTangentModulus ( ChMatrixNM< double, 6, 6 > &  C,
const ChMatrix33d F,
const ChMatrix33d l,
ChFieldData data_per_point,
ChElementData data_per_element 
)
inlineoverridevirtual

Implement interface to lower level stress material.

Parameters
Coutput C tangent modulus, dP=C*dE
Fcurrent deformation gradient tensor
lcurrent spatial velocity gradient (might be nullptr if IsSpatialVelocityGradientNeeded() is false)
data_per_pointpointer to auxiliary data (ex states), if any, per quadrature point
data_per_elementpointer to auxiliary data (ex states), if any, per element

Implements chrono::fea::ChMaterial3DStress.


The documentation for this class was generated from the following file:
  • /builds/uwsbel/chrono/src/chrono/fea/ChMaterial3DHyperelastic.h