chrono::fea::ChMaterial3DStress Class Referenceabstract

Description

Class for the basic properties of materials in 3D continua where stress can be computed from some deformation measure and some state.


#include <ChMaterial3DStress.h>

Inheritance diagram for chrono::fea::ChMaterial3DStress:
Collaboration diagram for chrono::fea::ChMaterial3DStress:

Public Member Functions

virtual void ComputeStress (ChStressTensor<> &S_stress, const ChMatrix33d &F_def, const ChMatrix33d *l, ChFieldData *data_per_point, ChElementData *data_per_element)=0
 Compute elastic stress from finite strain, passed as 3x3 deformation gradient tensor F_def. More...
 
virtual void ComputeTangentModulus (ChMatrixNM< double, 6, 6 > &C, const ChMatrix33d &F_def, const ChMatrix33d *l, ChFieldData *data_per_point, ChElementData *data_per_element)=0
 Computes the 6x6 tangent modulus for a given strain, assuming the definition dS = C * dE, that maps
variations of S, 2nd Piola-Kirchhoff stress, and of E, Green Lagrange strains, both in Voigt notation. More...
 
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...
 
virtual bool IsSpatialVelocityGradientNeeded () const =0
 Some material need info on the spatial velocity gradient l=\nabla_x v , where the time derivative of the deformation gradient F is dF/dt = l*F. 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

◆ ComputeStress()

virtual void chrono::fea::ChMaterial3DStress::ComputeStress ( ChStressTensor<> &  S_stress,
const ChMatrix33d F_def,
const ChMatrix33d l,
ChFieldData data_per_point,
ChElementData data_per_element 
)
pure virtual

Compute elastic stress from finite strain, passed as 3x3 deformation gradient tensor F_def.

Assuming stress is 2nd Piola-Kirchhoff tensor "S_stress", in Voigt notation. Some materials could also make use of spatial velocity gradient "l", ex. for damping effects, and of the "data_per_point" auxiliary structure, that could contain states like in plasticity.

Parameters
S_stressoutput stress, PK2
F_defcurrent 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

Implemented in chrono::fea::ChMaterial3DStressKelvinVoigt, and chrono::fea::ChMaterial3DHyperelastic.

◆ ComputeTangentModulus()

virtual void chrono::fea::ChMaterial3DStress::ComputeTangentModulus ( ChMatrixNM< double, 6, 6 > &  C,
const ChMatrix33d F_def,
const ChMatrix33d l,
ChFieldData data_per_point,
ChElementData data_per_element 
)
pure virtual

Computes the 6x6 tangent modulus for a given strain, assuming the definition dS = C * dE, that maps
variations of S, 2nd Piola-Kirchhoff stress, and of E, Green Lagrange strains, both in Voigt notation.

Assuming current strain is a 3x3 deformation gradient tensor F_def.

Parameters
Coutput C tangent modulus, as dS=C*dE
F_defcurrent 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

Implemented in chrono::fea::ChMaterial3DStressKelvinVoigt, and chrono::fea::ChMaterial3DHyperelastic.

◆ ComputeUpdateEndStep()

virtual void chrono::fea::ChMaterial3DStress::ComputeUpdateEndStep ( ChFieldData data_per_point,
ChElementData data_per_element,
const double  time 
)
inlinevirtual

Update your own auxiliary data, if any, at the end of time step (ex for plasticity).

This is called at the end of every time step (or nl static step)

Parameters
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

Reimplemented in chrono::fea::ChMaterial3DStressKelvinVoigt.

◆ IsSpatialVelocityGradientNeeded()

virtual bool chrono::fea::ChMaterial3DStress::IsSpatialVelocityGradientNeeded ( ) const
pure virtual

Some material need info on the spatial velocity gradient l=\nabla_x v , where the time derivative of the deformation gradient F is dF/dt = l*F.

Some others, do not need this info. For optimization reason, then, the ChDomainXXYY queries this, and knows if the "l" parameter could be left to null when calling ComputeStress(...)

Implemented in chrono::fea::ChMaterial3DHyperelastic, and chrono::fea::ChMaterial3DStressKelvinVoigt.


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