chrono::fea::ChMaterial3DStressParallel Class Reference

Description

Class for 3D stress materials A B that are composed by two generic ChMaterial3DStress in parallel, ex.

an elastic law and a viscous law, that contribute in parallel to the total stress. That is, stress is composed additively:
S_total = S_a + S_b. Btw also tangent modulus is C_total = C_a + C_b, and D_total = D_a + D_b, where the contributions a and b are computed by the two materials in parallel.

#include <ChMaterial3DStressParallel.h>

Inheritance diagram for chrono::fea::ChMaterial3DStressParallel:
Collaboration diagram for chrono::fea::ChMaterial3DStressParallel:

Public Member Functions

void SetMaterialA (std::shared_ptr< ChMaterial3DStress > mA)
 Set the first material (A) in the parallel composition.
 
void SetMaterialB (std::shared_ptr< ChMaterial3DStress > mB)
 Set the second material (B) in the parallel composition.
 
std::shared_ptr< ChMaterial3DStressGetMaterialA () const
 Get the first material (A) in the parallel composition.
 
std::shared_ptr< ChMaterial3DStressGetMaterialB () const
 Get the second material (B) in the parallel composition.
 
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 finite strain, passed as 3x3 deformation gradient tensor F. More...
 
virtual void ComputeTangentModulus (ChMatrixNM< double, 6, 6 > &C, ChMatrixNM< double, 6, 6 > *D, const ChMatrix33d &F, const ChMatrix33d *l, ChFieldData *data_per_point, ChElementData *data_per_element) override
 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 override
 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...
 
virtual std::unique_ptr< ChFieldDataCreateMaterialPointData () const override
 A compound material might generate custom data per material point because materials A or B needs it. 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.
 

Protected Attributes

std::shared_ptr< ChMaterial3DStressmaterial_A
 
std::shared_ptr< ChMaterial3DStressmaterial_B
 
- Protected Attributes inherited from chrono::fea::ChMaterial3DDensity
double m_density
 

Member Function Documentation

◆ ComputeStress()

virtual void chrono::fea::ChMaterial3DStressParallel::ComputeStress ( ChStressTensor<> &  S_stress,
const ChMatrix33d F,
const ChMatrix33d l,
ChFieldData data_per_point,
ChElementData data_per_element 
)
inlineoverridevirtual

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

Assuming stress is 2nd Piola-Kirchhoff tensor "S_stress", in Voigt notation.

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

Implements chrono::fea::ChMaterial3DStress.

◆ ComputeTangentModulus()

virtual void chrono::fea::ChMaterial3DStressParallel::ComputeTangentModulus ( ChMatrixNM< double, 6, 6 > &  C,
ChMatrixNM< double, 6, 6 > *  D,
const ChMatrix33d F,
const ChMatrix33d l,
ChFieldData data_per_point,
ChElementData data_per_element 
)
inlineoverridevirtual

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.

Parameters
Coutput C tangent modulus, as dS=C*dE
Doutput D tangent modulus, as dS=C*d(dE/dt) (might be nullptr if IsSpatialVelocityGradientNeeded() is false)
Fcurrent 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

Implements chrono::fea::ChMaterial3DStress.

◆ ComputeUpdateEndStep()

virtual void chrono::fea::ChMaterial3DStressParallel::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 from chrono::fea::ChMaterial3DStress.

◆ CreateMaterialPointData()

virtual std::unique_ptr<ChFieldData> chrono::fea::ChMaterial3DStressParallel::CreateMaterialPointData ( ) const
inlineoverridevirtual

A compound material might generate custom data per material point because materials A or B needs it.

In this case, it will create a custom data that contains the union of the custom data of A and B, and redirect to them the calls to update this data.

Reimplemented from chrono::fea::ChMaterial.

◆ IsSpatialVelocityGradientNeeded()

virtual bool chrono::fea::ChMaterial3DStressParallel::IsSpatialVelocityGradientNeeded ( ) const
inlineoverridevirtual

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(...)

Implements chrono::fea::ChMaterial3DStress.


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