Description
Class for properties of a non-linear thermal material in finite element problems.
This a quite pupular type of non-linearity, where the thermal properties (c specific heat, k thermal conductivity) are scalar functions of temperature T, i.e. c=c(T), k=k(T),
assuming isotropic thermal conductivity matrix [k](T)=[I]*k(T), as in: (c(T)*rho) dT/dt - div ([k](T) * grad T) = q_source
#include <ChMaterial3DThermalNonlinear.h>


Public Member Functions | |
| void | SetThermalConductivity (std::shared_ptr< ChFunction > mk) |
| Sets the k conductivity constant of the material as a function of temperature T, k=k(T) where k has the units of watts per meter kelvin [ W/(m K) ]. More... | |
| void | SetThermalConductivity (double mk) |
| Shortcut: sets the k conductivity constant of the material as a constant ChFunction. More... | |
| std::shared_ptr< ChFunction > | GetThermalConductivity () const |
| Gets the k conductivity of the material as a function of temperature T, k=k(T) expressed in watts per meter kelvin (W/(m K)). | |
| void | SetSpecificHeatCapacity (std::shared_ptr< ChFunction > mc) |
| Sets the c mass-specific heat capacity of the material, as a function of temperature T, c=c(T) where c has the units of Joule per kg Kelvin [ J / (kg K) ]. | |
| void | SetSpecificHeatCapacity (double mc) |
| Shortcut: sets the c mass-specific heat capacity of the material, as a function of temperature T, c=c(T) where c has the units of Joule per kg Kelvin [ J / (kg K) ]. | |
| std::shared_ptr< ChFunction > | GetSpecificHeatCapacity () const |
| Gets the c mass-specific heat capacity of the material, as a function of temperature T, c=c(T) where c has the units of Joule per kg Kelvin [ J / (kg K) ]. | |
| ChMatrixDynamic & | GetConductivityMatrix (double T) |
| Get the k conductivity matrix for a given temperature T. It is also the constitutive matrix of the Poisson problem. | |
| 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 nonlinear model with temperature-dependant conductivity k=k(T) this is q_flux = -[k](T) * grad_T with conductivity matrix [k](T) = [I] * k(T) . 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 nonlinear model, this is [C] = [k](T) = [I] * k(T) 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 nonlinear model with temperature-dependant conductivity k=k(T) this is q_flux = -[k](T) * grad_T with conductivity matrix [k](T) = [I] * k(T) .
- 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 nonlinear model, this is [C] = [k](T) = [I] * k(T)
- 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.
◆ SetThermalConductivity() [1/2]
|
inline |
Shortcut: sets the k conductivity constant of the material as a constant ChFunction.
where k has the units of watts per meter kelvin [ W/(m K) ].
◆ SetThermalConductivity() [2/2]
|
inline |
Sets the k conductivity constant of the material as a function of temperature T, k=k(T) where k has the units of watts per meter kelvin [ W/(m K) ].
Ex. (approx.): water = 0.6, aluminium = 200, steel = 50, plastics=0.9-0.2 Note: conductivity matrix [k] is assumed isotropic (diagonal matrix ie. [k]=[I]*k(T) )
The documentation for this class was generated from the following file:
- /builds/uwsbel/chrono/src/chrono/fea/ChMaterial3DThermalNonlinear.h
Public Member Functions inherited from