Loading [MathJax]/extensions/tex2jax.js

Description

Class for all elastic materials that can undergo plastic flow.

Defines simply some interface functions.

#include <ChContinuumMaterial.h>

Inheritance diagram for chrono::fea::ChContinuumElastoplastic:
Collaboration diagram for chrono::fea::ChContinuumElastoplastic:

Public Member Functions

 ChContinuumElastoplastic (double young=10000000, double poisson=0.4, double density=1000)
 
virtual double ComputeYieldFunction (const ChStressTensor<> &stress) const =0
 Return a scalar value that is 0 on the yield surface, <0 inside (elastic), >0 outside (incompatible->plastic flow)
 
virtual void ComputePlasticStrainFlow (ChStrainTensor<> &plasticstrainflow, const ChStrainTensor<> &totstrain) const =0
 Compute plastic strain flow (flow derivative dE_plast/dt) from strain, according to VonMises strain yield theory.
 
virtual void ComputeReturnMapping (ChStrainTensor<> &plasticstrainflow, const ChStrainTensor<> &incrementstrain, const ChStrainTensor<> &lastelasticstrain, const ChStrainTensor<> &lastplasticstrain) const =0
 Correct the strain-stress by enforcing that elastic stress must remain on the yield surface, computing a plastic flow to be added to plastic strain while integrating.
 
virtual void SetPlasticFlowRate (double flow_rate)=0
 Set the plastic flow rate, i.e. More...
 
virtual double GetPlasticFlowRate () const =0
 Set the plastic flow rate.
 
virtual void ArchiveOut (ChArchiveOut &archive_out) override
 
virtual void ArchiveIn (ChArchiveIn &archive_in) override
 
- Public Member Functions inherited from chrono::fea::ChContinuumElastic
 ChContinuumElastic (double young=10000000, double poisson=0.4, double density=1000)
 
 ChContinuumElastic (const ChContinuumElastic &other)
 
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 )
 
void ComputeStressStrainMatrix ()
 Computes Elasticity matrix and stores the value in this->StressStrainMatrix Note: is performed every time you change a material parameter.
 
ChMatrixDynamicGetStressStrainMatrix ()
 Get the Elasticity matrix.
 
void ComputeElasticStress (ChStressTensor<> &stress, const ChStrainTensor<> &strain) const
 Compute elastic stress from elastic strain (using column tensors, in Voight notation)
 
void ComputeElasticStrain (ChStrainTensor<> &strain, const ChStressTensor<> &stress) const
 Compute elastic strain from elastic stress (using column tensors, in Voight notation)
 
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::ChContinuumMaterial
 ChContinuumMaterial (double density=1000)
 
 ChContinuumMaterial (const ChContinuumMaterial &other)
 
void SetDensity (double density)
 Set the density of the material, in kg/m^2.
 
double GetDensity () const
 Get the density of the material, in kg/m^2.
 

Additional Inherited Members

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

Member Function Documentation

◆ SetPlasticFlowRate()

virtual void chrono::fea::ChContinuumElastoplastic::SetPlasticFlowRate ( double  flow_rate)
pure virtual

Set the plastic flow rate, i.e.

the 'creep' speed. The lower the value, the slower the plastic flow during dynamic simulations, with delayed plasticity.

Implemented in chrono::fea::ChContinuumDruckerPrager, and chrono::fea::ChContinuumPlasticVonMises.


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