Description
Class for properties of a linear thermal material in finite element problems.
This is the simplest case: it contains constant properties (c specific heat, k thermal conductivity, rho density) for the linear thermal problems PDEs of the type (c*rho) dT/dt + div (-[k]*grad T) = q_source
#include <ChMaterial3DThermalLinear.h>


Public Member Functions | |
| void | SetThermalConductivity (double mk) |
| Sets the k conductivity constant of the material, expressed in watts per meter kelvin [ W/(m K) ]. More... | |
| double | GetThermalConductivity () const |
| Gets the k conductivity constant of the material, expressed in watts per meter kelvin (W/(m K)). | |
| void | SetSpecificHeatCapacity (double mc) |
| Sets the c mass-specific heat capacity of the material, expressed as Joule per kg Kelvin [ J / (kg K) ]. | |
| double | GetSpecificHeatCapacity () const |
| Sets the c mass-specific heat capacity of the material, expressed as Joule per kg Kelvin [ J / (kg K) ]. | |
| ChMatrixDynamic & | GetConductivityMatrix () |
| Get the k conductivity matrix. More... | |
| virtual double | Get_DtMultiplier () override |
| override base: (the dT/dt term has multiplier rho*c with rho=density, c=heat capacity) | |
| virtual void | ComputeHeatFlux (ChVector3d &q_flux, const ChVector3d grad_T, const double T, ChFieldData *data_per_point, ChElementData *data_per_element) override |
| Compute heat flux vector q_flux from a temperature gradient grad_T q_flux = f(grad_T) Because of linear model, this is q_flux = -[k]*grad_T with [k] conductivity matrix. More... | |
| virtual void | ComputeTangentModulus (ChMatrix33d &tangent_matrix, const ChVector3d grad_T, const double T, ChFieldData *data_per_point, ChElementData *data_per_element) override |
| Computes the tangent conductivity [C] in the linearization of q_flux = f(grad_T) about temperature: δq_flux = [C] δgrad_T Because of linear model, this is [C] = [k] conductivity matrix. More... | |
| virtual void | ComputeDtMultiplier (double &Dt_multiplier, const double T, ChFieldData *data_per_point, ChElementData *data_per_element) override |
| Computes the Z multiplier in Z*dT/dt + div q_flux = q_source: Because of linear model, this is Z=(m_density * c_mass_specific_heat_capacity) More... | |
Public Member Functions inherited from chrono::fea::ChMaterialPoisson | |
| ChMatrixDynamic & | ConstitutiveMatrix () |
| The constitutive matrix [C] to compute the bilinear form in the weak formulation. 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::ChMaterialPoisson | |
| ChMatrixDynamic | constitutiveMatrix |
Protected Attributes inherited from chrono::fea::ChMaterial3DDensity | |
| double | m_density |
Member Function Documentation
◆ ComputeDtMultiplier()
|
inlineoverridevirtual |
Computes the Z multiplier in Z*dT/dt + div q_flux = q_source: Because of linear model, this is Z=(m_density * c_mass_specific_heat_capacity)
- Parameters
-
Dt_multiplier the Z term in Z*dT/dt + div(q_flux) = q_source T current temperature (not used here) data_per_point pointer to auxiliary data (ex states), if any, per quadrature point (not used here) data_per_element pointer to auxiliary data (ex states), if any, per element (not used here)
Implements chrono::fea::ChMaterial3DThermal.
◆ ComputeHeatFlux()
|
inlineoverridevirtual |
Compute heat flux vector q_flux from a temperature gradient grad_T q_flux = f(grad_T) Because of linear model, this is q_flux = -[k]*grad_T with [k] conductivity matrix.
- Parameters
-
q_flux output stress, PK2 grad_T current temperature gradient T current temperature (not used here) data_per_point pointer to auxiliary data (ex states), if any, per quadrature point (not used here) data_per_element pointer to auxiliary data (ex states), if any, per element (not used here)
Implements chrono::fea::ChMaterial3DThermal.
◆ ComputeTangentModulus()
|
inlineoverridevirtual |
Computes the tangent conductivity [C] in the linearization of q_flux = f(grad_T) about temperature: δq_flux = [C] δgrad_T Because of linear model, this is [C] = [k] conductivity matrix.
- Parameters
-
tangent_matrix output [C] tangent as δq_flux = [C] δgrad_T grad_T current temperature gradient (not used here) T current temperature (not used here) data_per_point pointer to auxiliary data (ex states), if any, per quadrature point (not used here) data_per_element pointer to auxiliary data (ex states), if any, per element (not used here)
Implements chrono::fea::ChMaterial3DThermal.
◆ GetConductivityMatrix()
|
inline |
Get the k conductivity matrix.
It is also the constitutive matrix of the Poisson problem. You can modify its values, for setting material property.
◆ SetThermalConductivity()
|
inline |
Sets the k conductivity constant of the material, expressed in watts per meter kelvin [ W/(m K) ].
Ex. (approx.): water = 0.6, aluminium = 200, steel = 50, plastics=0.9-0.2 Sets the conductivity matrix as isotropic (diagonal k)
The documentation for this class was generated from the following file:
- /builds/uwsbel/chrono/src/chrono/fea/ChMaterial3DThermalLinear.h
Public Member Functions inherited from