chrono::fea::ChMaterial3DThermalStress Class Reference

Description

Material for coupled thermal and stress problems, ex thermoelasticity.

This is implemented here as a composition of two distinct materials, one for the 3D heat problem, which is independent on stress state and which we recycle from the usual ChMaterial3DThermal, and a one for the 3D stress problem, which could be one of the many ChMaterial3DStress classes, ex. StVenant or Mooney-Rivlin or Ogden etc.

#include <ChMaterial3DThermalStress.h>

Inheritance diagram for chrono::fea::ChMaterial3DThermalStress:
Collaboration diagram for chrono::fea::ChMaterial3DThermalStress:

Public Member Functions

void SetThermalExpansionCoefficient (double mte)
 Set thermal expansion coefficient [1/K]. Ex. aluminium (23*10e-6)/K, steel (11*10e-6)/K.
 
double GetThermalExpansionCoefficient () const
 Get thermal expansion coefficient [1/K].
 
void SetRestTemperature (double mt)
 Set rest temperature [K]. No thermal expanson at point, if point temperature is this.
 
double GetRestTemperature () const
 Get rest temperature [K]. No thermal expanson at point, if point temperature is this.
 
virtual void SetDensity (double md) override
 Set density of material. Both for stress problem and for thermal problem.
 
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
 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
 Implement this because the material might need custom data per material point. More...
 
- Public Member Functions inherited from chrono::fea::ChMaterial3DDensity
 ChMaterial3DDensity (double density=1000)
 
virtual double GetDensity () const
 Get the density of the material, in kg/m^2.
 

Public Attributes

std::shared_ptr< ChMaterial3DStressmaterial_stress
 Material for stress field computation.
 
std::shared_ptr< ChMaterial3DThermalmaterial_thermal
 Material for thermal field computation.
 

Additional Inherited Members

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

Member Function Documentation

◆ ComputeUpdateEndStep()

virtual void chrono::fea::ChMaterial3DThermalStress::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

◆ CreateMaterialPointData()

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

Implement this because the material might need custom data per material point.

Here the user can have plugged a material_stress of variuous type, so fall back to its custom data creation

Reimplemented from chrono::fea::ChMaterial.

◆ IsSpatialVelocityGradientNeeded()

virtual bool chrono::fea::ChMaterial3DThermalStress::IsSpatialVelocityGradientNeeded ( ) const
inlinevirtual

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


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