chrono::fea::ChMaterial3DStressOgden Class Reference

Description

Class for Ogden materials.

A subtype of hyperelastic materials. It is mostly used for plastics, rubbers, organic tissue, etc. It supports very large strains. Note, it depends on pairs of {mu_i, alpha_i} parameters, plus one parameter for bulk modulus. For stability, it must be (mu_i * alpha_i)>0, for all pairs.

#include <ChMaterial3DStressOgden.h>

Inheritance diagram for chrono::fea::ChMaterial3DStressOgden:
Collaboration diagram for chrono::fea::ChMaterial3DStressOgden:

Public Member Functions

std::vector< double > & CoefficientsMu ()
 access coefficients "mu" in Ogden material
 
std::vector< double > & CoefficientsAlpha ()
 access coefficients "alpha" in Ogden material
 
void SetKappa (double mk)
 Set bulk modulus kappa.
 
double GetKappa ()
 Get bulk modulus kappa.
 
void SetAsEquivalentNeoHookean (double E, double nu)
 For testing reasons or other reasons, initialize {alpha, mu} parameters of this material such that, for small deformations, it behaves like the linear elastic material with a given Young modulus and Poisson coeff. More...
 
virtual void ComputeElasticStress (ChStressTensor<> &stress, const ChMatrix33d &C_deformation) override
 Compute elastic stress from elastic strain. More...
 
- 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(...)
 
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::ChMaterial3DStressOgden::ComputeElasticStress ( ChStressTensor<> &  stress,
const ChMatrix33d C_deformation 
)
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.

◆ SetAsEquivalentNeoHookean()

void chrono::fea::ChMaterial3DStressOgden::SetAsEquivalentNeoHookean ( double  E,
double  nu 
)
inline

For testing reasons or other reasons, initialize {alpha, mu} parameters of this material such that, for small deformations, it behaves like the linear elastic material with a given Young modulus and Poisson coeff.

It takes just two couples of {alpha, mu}.

Parameters
EYoung modulus
nuPoisson

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