![]() |
v0.14.0 |
Files | |
file | ThermalElement.cpp |
file | ThermalElement.hpp |
Operators and data structures for thermal analysis. | |
file | ThermalStressElement.hpp |
Implementation of thermal stresses element. | |
file | thermal_steady.cpp |
Example of steady thermal analysis. | |
file | thermal_unsteady.cpp |
Example of thermal unsteady analyze. | |
Classes | |
struct | ThermalElement |
structure grouping operators and data used for thermal problems More... | |
struct | ThermalElement::MyVolumeFE |
definition of volume element More... | |
struct | ThermalElement::MyTriFE |
define surface element More... | |
struct | ThermalElement::BlockData |
data for calculation heat conductivity and heat capacity elements More... | |
struct | ThermalElement::FluxData |
data for calculation heat flux More... | |
struct | ThermalElement::ConvectionData |
data for convection More... | |
struct | ThermalElement::RadiationData |
data for radiation More... | |
struct | ThermalElement::CommonData |
common data used by volume elements More... | |
struct | ThermalElement::OpGetGradAtGaussPts |
operator to calculate temperature gradient at Gauss points More... | |
struct | ThermalElement::OpGetFieldAtGaussPts< OP > |
operator to calculate temperature and rate of temperature at Gauss points More... | |
struct | ThermalElement::OpGetTetTemperatureAtGaussPts |
operator to calculate temperature at Gauss pts More... | |
struct | ThermalElement::OpGetTriTemperatureAtGaussPts |
operator to calculate temperature at Gauss pts More... | |
struct | ThermalElement::OpGetTetRateAtGaussPts |
operator to calculate temperature rate at Gauss pts More... | |
struct | ThermalElement::OpThermalRhs |
struct | ThermalElement::OpThermalLhs |
struct | ThermalElement::OpHeatCapacityRhs |
operator to calculate right hand side of heat capacity terms More... | |
struct | ThermalElement::OpHeatCapacityLhs |
operator to calculate left hand side of heat capacity terms More... | |
struct | ThermalElement::OpHeatFlux |
operator for calculate heat flux and assemble to right hand side More... | |
struct | ThermalElement::OpRadiationLhs |
struct | ThermalElement::OpRadiationRhs |
operator to calculate radiation therms on body surface and assemble to rhs of transient equations(Residual Vector) More... | |
struct | ThermalElement::OpConvectionRhs |
operator to calculate convection therms on body surface and assemble to rhs of equations More... | |
struct | ThermalElement::OpConvectionLhs |
struct | ThermalElement::UpdateAndControl |
this calass is to control time stepping More... | |
struct | ThermalElement::TimeSeriesMonitor |
TS monitore it records temperature at time steps. More... | |
struct | ThermalStressElement |
Implentation of thermal stress element. More... | |
struct | ThermalStressElement::MyVolumeFE |
struct | ThermalStressElement::BlockData |
struct | ThermalStressElement::CommonData |
struct | ThermalStressElement::OpGetTemperatureAtGaussPts |
struct | ThermalStressElement::OpThermalStressRhs |
Functions | |
MoFEMErrorCode | ThermalElement::addThermalElements (const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS") |
add thermal element on tets More... | |
MoFEMErrorCode | ThermalElement::addThermalFluxElement (const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS") |
add heat flux element More... | |
MoFEMErrorCode | ThermalElement::addThermalConvectionElement (const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS") |
add convection element More... | |
MoFEMErrorCode | ThermalElement::addThermalRadiationElement (const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS") |
add Non-linear Radiation element More... | |
MoFEMErrorCode | ThermalElement::setThermalFiniteElementRhsOperators (string field_name, Vec &F) |
this function is used in case of stationary problem to set elements for rhs More... | |
MoFEMErrorCode | ThermalElement::setThermalFiniteElementLhsOperators (string field_name, Mat A) |
this function is used in case of stationary heat conductivity problem for lhs More... | |
MoFEMErrorCode | ThermalElement::setThermalFluxFiniteElementRhsOperators (string field_name, Vec &F, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS") |
this function is used in case of stationary problem for heat flux terms More... | |
MoFEMErrorCode | ThermalElement::setTimeSteppingProblem (string field_name, string rate_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS") |
set up operators for unsteady heat flux; convection; radiation problem More... | |
MoFEMErrorCode | ThermalElement::setTimeSteppingProblem (TsCtx &ts_ctx, string field_name, string rate_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS") |
set up operators for unsteady heat flux; convection; radiation problem More... | |
MoFEMErrorCode ThermalElement::addThermalConvectionElement | ( | const std::string | field_name, |
const std::string | mesh_nodals_positions = "MESH_NODE_POSITIONS" |
||
) |
add convection element
It get data from convection set and define element in moab. Alternatively uses block set with name including substring CONVECTION.
field | name |
name | of mesh nodal positions (if not defined nodal coordinates are used) |
Definition at line 621 of file ThermalElement.cpp.
MoFEMErrorCode ThermalElement::addThermalElements | ( | const std::string | field_name, |
const std::string | mesh_nodals_positions = "MESH_NODE_POSITIONS" |
||
) |
add thermal element on tets
It get data from block set and define element in moab w
field | name |
name | of mesh nodal positions (if not defined nodal coordinates are used) |
Definition at line 527 of file ThermalElement.cpp.
MoFEMErrorCode ThermalElement::addThermalFluxElement | ( | const std::string | field_name, |
const std::string | mesh_nodals_positions = "MESH_NODE_POSITIONS" |
||
) |
add heat flux element
It get data from heat flux set and define element in moab. Alternatively uses block set with name including substring HEAT_FLUX.
field | name |
name | of mesh nodal positions (if not defined nodal coordinates are used) |
Definition at line 574 of file ThermalElement.cpp.
MoFEMErrorCode ThermalElement::addThermalRadiationElement | ( | const std::string | field_name, |
const std::string | mesh_nodals_positions = "MESH_NODE_POSITIONS" |
||
) |
add Non-linear Radiation element
It get data from Radiation set and define element in moab. Alternatively uses block set with name including substring RADIATION.
field | name |
name | of mesh nodal positions (if not defined nodal coordinates are used) |
Definition at line 659 of file ThermalElement.cpp.
MoFEMErrorCode ThermalElement::setThermalFiniteElementLhsOperators | ( | string | field_name, |
Mat | A | ||
) |
this function is used in case of stationary heat conductivity problem for lhs
Definition at line 713 of file ThermalElement.cpp.
MoFEMErrorCode ThermalElement::setThermalFiniteElementRhsOperators | ( | string | field_name, |
Vec & | F | ||
) |
this function is used in case of stationary problem to set elements for rhs
Definition at line 699 of file ThermalElement.cpp.
MoFEMErrorCode ThermalElement::setThermalFluxFiniteElementRhsOperators | ( | string | field_name, |
Vec & | F, | ||
const std::string | mesh_nodals_positions = "MESH_NODE_POSITIONS" |
||
) |
this function is used in case of stationary problem for heat flux terms
Definition at line 724 of file ThermalElement.cpp.
MoFEMErrorCode ThermalElement::setTimeSteppingProblem | ( | string | field_name, |
string | rate_name, | ||
const std::string | mesh_nodals_positions = "MESH_NODE_POSITIONS" |
||
) |
set up operators for unsteady heat flux; convection; radiation problem
Definition at line 774 of file ThermalElement.cpp.
MoFEMErrorCode ThermalElement::setTimeSteppingProblem | ( | TsCtx & | ts_ctx, |
string | field_name, | ||
string | rate_name, | ||
const std::string | mesh_nodals_positions = "MESH_NODE_POSITIONS" |
||
) |
set up operators for unsteady heat flux; convection; radiation problem
Definition at line 859 of file ThermalElement.cpp.