Description
Class for linear damping (Newton damping in lagrangian solid, as S=D:E_dot where E_dot is the time derivative of Green Lagrange strain and D is a constant tensor build from viscosity coefficients, or more in detail: S = lambda tr(E_dot) I + 2μ E_dot ).
This is a bit different from the ChMaterial3DStressViscoNewton, which provides damping forces proportional to the spatial velocity, as in CFD fluids. But this is faster to compute. It is assumed to be used in parallel with some elastic stress model like ChMaterial3DStressStVenant. If used alone, it provides a fluid-like viscous damping.
#include <ChMaterial3DStressViscoLinear.h>


Public Member Functions | |
| virtual void | SetDeviatoricDamping (double mdamp) |
| Set the damping for the deviatoric effect. More... | |
| virtual double | GetDeviatoricDamping () const |
| virtual void | SetVolumetricDamping (double mdamp) |
| Set the damping for the volumetric effect. More... | |
| virtual double | GetVolumetricDamping () const |
| virtual void | ComputeStress (ChStressTensor<> &S_stress, const ChMatrix33d &F, const ChMatrix33d *l, ChFieldData *data_per_point, ChElementData *data_per_element) override |
| Compute elastic stress from spatial velocity gradient "l". More... | |
| virtual void | ComputeTangentModulus (ChMatrixNM< double, 6, 6 > &C, ChMatrixNM< double, 6, 6 > *D, const ChMatrix33d &F_def, const ChMatrix33d *l, ChFieldData *data_per_point, ChElementData *data_per_element) override |
| Computes the tangent modulus. More... | |
| virtual bool | IsSpatialVelocityGradientNeeded () const |
| This material need info on the spatial velocity gradient l=\nabla_x v. | |
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
◆ ComputeStress()
|
inlineoverridevirtual |
Compute elastic stress from spatial velocity gradient "l".
- Parameters
-
S_stress output stress, PK2 F current deformation gradient tensor l current spatial velocity gradient (might be nullptr if IsSpatialVelocityGradientNeeded() is false) data_per_point pointer to auxiliary data (ex states), if any, per quadrature point data_per_element pointer to auxiliary data (ex states), if any, per element
Implements chrono::fea::ChMaterial3DStress.
◆ ComputeTangentModulus()
|
inlineoverridevirtual |
Computes the tangent modulus.
- Parameters
-
C output C tangent modulus, as dS=C*dE D output D tangent modulus, as dS=C*d(E_dot) (maybe nullptr if IsSpatialVelocityGradientNeeded() is false) F_def current deformation gradient tensor l current spatial velocity gradient (might be nullptr if IsSpatialVelocityGradientNeeded() is false) data_per_point pointer to auxiliary data (ex states), if any, per quadrature point data_per_element pointer to auxiliary data (ex states), if any, per element
Implements chrono::fea::ChMaterial3DStress.
◆ SetDeviatoricDamping()
|
inlinevirtual |
Set the damping for the deviatoric effect.
This corresponds to the viscous damping μ in Newtonian fluids.
◆ SetVolumetricDamping()
|
inlinevirtual |
Set the damping for the volumetric effect.
If the material is incompressible or nearly incompressible, this term has no effect.
The documentation for this class was generated from the following file:
- /builds/uwsbel/chrono/src/chrono/fea/multiphysics/ChMaterial3DStressViscoLinear.h
Public Member Functions inherited from