Description
Class for properties of the 3D elasticity from StVenant-Kirchhoff model.
It is the simplest type of hyperelastic material. For S PiolaKirchhoff stress and E Green Lagrange strain tensors, it is S=C:E with C constant 4th order elastic tensor. For small deformations it corresponds to the classical σ=C:ε. In Voigt notation, S=[C]*E, σ=[C]*ε
#include <ChMaterial3DStressStVenant.h>


Public Member Functions | |
| ChMaterial3DStressStVenant (double young=10000000, double poisson=0.4, double density=1000) | |
| void | SetYoungModulus (double E) |
| Set the Young elastic modulus, in Pa (N/m^2), as the ratio of the uniaxial stress over the uniaxial strain, for hookean materials. | |
| double | GetYoungModulus () const |
| Get the Young elastic modulus, in Pa (N/m^2). | |
| void | SetPoissonRatio (double v) |
| Set the Poisson ratio, as v=-transverse_strain/axial_strain, so takes into account the 'squeezing' effect of materials that are pulled. More... | |
| double | GetPoissonRatio () const |
| Get the Poisson ratio, as v=-transverse_strain/axial_strain. | |
| void | SetShearModulus (double G) |
| Set the shear modulus G, in Pa (N/m^2). More... | |
| double | GetShearModulus () const |
| Get the shear modulus G, in Pa (N/m^2) | |
| double | GetLameFirstParam () const |
| Get Lamé first parameter (the second is shear modulus, so GetShearModulus() ) | |
| double | GetBulkModulus () const |
| Get bulk modulus (increase of pressure for decrease of volume), in Pa. | |
| double | GetPWaveModulus () const |
| Get P-wave modulus (if V=speed of propagation of a P-wave, then (M/density)=V^2 ) | |
| virtual void | ComputeElasticStress (ChStressTensor<> &stress, const ChMatrix33d &C_deformation) override |
| Compute elastic stress from elastic strain. More... | |
| virtual void | ComputeElasticTangentModulus (ChMatrixNM< double, 6, 6 > &C, const ChMatrix33d &C_def) override |
| Computes the tangent modulus C. More... | |
| void | SetRayleighDampingAlpha (double alpha) |
| Set the Rayleigh mass-proportional damping factor alpha, to build damping R as R=alpha*M + beta*K. | |
| double | GetRayleighDampingAlpha () const |
| Set the Rayleigh mass-proportional damping factor alpha, in R=alpha*M + beta*K. | |
| void | SetRayleighDampingBeta (double beta) |
| Set the Rayleigh stiffness-proportional damping factor beta, to build damping R as R=alpha*M + beta*K. | |
| double | GetRayleighDampingBeta () const |
| Set the Rayleigh stiffness-proportional damping factor beta, in R=alpha*M + beta*K. | |
Public Member Functions inherited from chrono::fea::ChMaterial3DHyperelastic | |
| 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(...) | |
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< ChFieldData > | CreateMaterialPointData () 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()
|
inlineoverridevirtual |
Compute elastic stress from elastic strain.
Starts computing Green-Lagrange strain E from C_deformation, the right Cauchy-Green deformation. For small strains the Green Lagrange strain in Voigt notation coincides with espilon tensor. Return stress as Piola-Kirchhoff S tensor, in Voigt notation. This is a very simple material, ie. a linear funciton S=C:E with C as 4th order constant tensor, also S=[C]*E with 6x6 C in Voigt notation.
Implements chrono::fea::ChMaterial3DHyperelastic.
◆ ComputeElasticTangentModulus()
|
inlineoverridevirtual |
Computes the tangent modulus C.
It takes the right Cauchy-Green deformation tensor C_def, and returns the tangent modulus in 6x6 matrix C, for dS = C * dE where S is 2nd PK stress, and E is Green Lagrange strain. (The Cauchy-Green deformation "C_def" is not used here, as C in S=C:E is independent on C_def)
Reimplemented from chrono::fea::ChMaterial3DHyperelastic.
◆ SetPoissonRatio()
|
inline |
Set the Poisson ratio, as v=-transverse_strain/axial_strain, so takes into account the 'squeezing' effect of materials that are pulled.
Note: v=0.5 means perfectly incompressible material, that could give problems with some type of solvers. Setting v also changes G.
◆ SetShearModulus()
|
inline |
Set the shear modulus G, in Pa (N/m^2).
Setting G also changes Poisson ratio v.
The documentation for this class was generated from the following file:
- /builds/uwsbel/chrono/src/chrono/fea/ChMaterial3DStressStVenant.h
Public Member Functions inherited from