v0.15.0
Loading...
Searching...
No Matches
ThermalElement Struct Reference

structure grouping operators and data used for thermal problems More...

#include "users_modules/basic_finite_elements/src/ThermalElement.hpp"

Collaboration diagram for ThermalElement:
[legend]

Classes

struct  BlockData
 data for calculation heat conductivity and heat capacity elements More...
 
struct  CommonData
 common data used by volume elements More...
 
struct  ConvectionData
 data for convection More...
 
struct  FluxData
 data for calculation heat flux More...
 
struct  MyTriFE
 define surface element More...
 
struct  MyVolumeFE
 definition of volume element More...
 
struct  OpConvectionLhs
 
struct  OpConvectionRhs
 operator to calculate convection therms on body surface and assemble to rhs of equations More...
 
struct  OpGetFieldAtGaussPts
 operator to calculate temperature and rate of temperature at Gauss points More...
 
struct  OpGetGradAtGaussPts
 operator to calculate temperature gradient at Gauss points More...
 
struct  OpGetTetRateAtGaussPts
 operator to calculate temperature rate at Gauss pts More...
 
struct  OpGetTetTemperatureAtGaussPts
 operator to calculate temperature at Gauss pts More...
 
struct  OpGetTriTemperatureAtGaussPts
 operator to calculate temperature at Gauss pts More...
 
struct  OpHeatCapacityLhs
 operator to calculate left hand side of heat capacity terms More...
 
struct  OpHeatCapacityRhs
 operator to calculate right hand side of heat capacity terms More...
 
struct  OpHeatFlux
 operator for calculate heat flux and assemble to right hand side More...
 
struct  OpRadiationLhs
 
struct  OpRadiationRhs
 operator to calculate radiation therms on body surface and assemble to rhs of transient equations(Residual Vector) More...
 
struct  OpThermalLhs
 
struct  OpThermalRhs
 
struct  RadiationData
 data for radiation More...
 
struct  TimeSeriesMonitor
 TS monitore it records temperature at time steps. More...
 
struct  UpdateAndControl
 this calass is to control time stepping More...
 

Public Member Functions

MyVolumeFEgetLoopFeRhs ()
 get rhs volume element
 
MyVolumeFEgetLoopFeLhs ()
 get lhs volume element
 
MyTriFEgetLoopFeFlux ()
 
MyTriFEgetLoopFeConvectionRhs ()
 
MyTriFEgetLoopFeConvectionLhs ()
 
MyTriFEgetLoopFeRadiationRhs ()
 
MyTriFEgetLoopFeRadiationLhs ()
 
 ThermalElement (MoFEM::Interface &m_field)
 
MoFEMErrorCode addThermalElements (const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
 add thermal element on tets
 
MoFEMErrorCode addThermalFluxElement (const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
 add heat flux element
 
MoFEMErrorCode addThermalConvectionElement (const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
 add convection element
 
MoFEMErrorCode addThermalRadiationElement (const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
 add Non-linear Radiation element
 
MoFEMErrorCode setThermalFiniteElementRhsOperators (string field_name, Vec &F)
 this function is used in case of stationary problem to set elements for rhs
 
MoFEMErrorCode setThermalFiniteElementLhsOperators (string field_name, Mat A)
 this function is used in case of stationary heat conductivity problem for lhs
 
MoFEMErrorCode 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
 
MoFEMErrorCode setThermalConvectionFiniteElementRhsOperators (string field_name, Vec &F, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
 
MoFEMErrorCode setThermalConvectionFiniteElementLhsOperators (string field_name, Mat A, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
 
MoFEMErrorCode 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
 
MoFEMErrorCode 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
 

Public Attributes

MyVolumeFE feRhs
 cauclate right hand side for tetrahedral elements
 
MyVolumeFE feLhs
 
MyTriFE feFlux
 
MyTriFE feConvectionRhs
 
MyTriFE feConvectionLhs
 
MyTriFE feRadiationRhs
 
MyTriFE feRadiationLhs
 
MoFEM::InterfacemField
 
std::map< int, BlockDatasetOfBlocks
 maps block set id with appropriate BlockData
 
std::map< int, FluxDatasetOfFluxes
 maps side set id with appropriate FluxData
 
std::map< int, ConvectionDatasetOfConvection
 
std::map< int, RadiationDatasetOfRadiation
 
CommonData commonData
 

Detailed Description

structure grouping operators and data used for thermal problems

In order to assemble matrices and right hand vectors, the loops over elements, entities within the element and finally loop over integration points are executed.

Following implementation separate those three types of loops and to each loop attach operator.

Definition at line 27 of file ThermalElement.hpp.

Constructor & Destructor Documentation

◆ ThermalElement()

ThermalElement::ThermalElement ( MoFEM::Interface & m_field)
inline

Definition at line 84 of file ThermalElement.hpp.

85 : feRhs(m_field), feLhs(m_field), feFlux(m_field),
86 feConvectionRhs(m_field), feConvectionLhs(m_field),
87 feRadiationRhs(m_field), feRadiationLhs(m_field), mField(m_field) {}
MyVolumeFE feRhs
cauclate right hand side for tetrahedral elements
MoFEM::Interface & mField

Member Function Documentation

◆ getLoopFeConvectionLhs()

MyTriFE & ThermalElement::getLoopFeConvectionLhs ( )
inline

Definition at line 74 of file ThermalElement.hpp.

74{ return feConvectionLhs; }

◆ getLoopFeConvectionRhs()

MyTriFE & ThermalElement::getLoopFeConvectionRhs ( )
inline

Definition at line 71 of file ThermalElement.hpp.

71 {
72 return feConvectionRhs;
73 } //< get convection element

◆ getLoopFeFlux()

MyTriFE & ThermalElement::getLoopFeFlux ( )
inline

Definition at line 67 of file ThermalElement.hpp.

67{ return feFlux; } //< get heat flux element

◆ getLoopFeLhs()

MyVolumeFE & ThermalElement::getLoopFeLhs ( )
inline

get lhs volume element

Definition at line 54 of file ThermalElement.hpp.

◆ getLoopFeRadiationLhs()

MyTriFE & ThermalElement::getLoopFeRadiationLhs ( )
inline

Definition at line 81 of file ThermalElement.hpp.

81{ return feRadiationLhs; }

◆ getLoopFeRadiationRhs()

MyTriFE & ThermalElement::getLoopFeRadiationRhs ( )
inline

Definition at line 78 of file ThermalElement.hpp.

78 {
79 return feRadiationRhs;
80 } //< get radiation element

◆ getLoopFeRhs()

MyVolumeFE & ThermalElement::getLoopFeRhs ( )
inline

get rhs volume element

Definition at line 52 of file ThermalElement.hpp.

◆ setThermalConvectionFiniteElementLhsOperators()

MoFEMErrorCode ThermalElement::setThermalConvectionFiniteElementLhsOperators ( string field_name,
Mat A,
const std::string mesh_nodals_positions = "MESH_NODE_POSITIONS" )

Definition at line 759 of file ThermalElement.cpp.

760 {
762 bool hoGeometry = false;
763 if (mField.check_field(mesh_nodals_positions)) {
764 hoGeometry = true;
765 }
766 std::map<int, ConvectionData>::iterator sit = setOfConvection.begin();
767 for (; sit != setOfConvection.end(); sit++) {
768 // add finite element
770 new OpConvectionLhs(field_name, A, sit->second, hoGeometry));
771 }
773}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
virtual bool check_field(const std::string &name) const =0
check if field is in database
constexpr auto field_name
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
std::map< int, ConvectionData > setOfConvection

◆ setThermalConvectionFiniteElementRhsOperators()

MoFEMErrorCode ThermalElement::setThermalConvectionFiniteElementRhsOperators ( string field_name,
Vec & F,
const std::string mesh_nodals_positions = "MESH_NODE_POSITIONS" )

Definition at line 741 of file ThermalElement.cpp.

742 {
744 bool hoGeometry = false;
745 if (mField.check_field(mesh_nodals_positions)) {
746 hoGeometry = true;
747 }
748 std::map<int, ConvectionData>::iterator sit = setOfConvection.begin();
749 for (; sit != setOfConvection.end(); sit++) {
750 // add finite element
752 new OpGetTriTemperatureAtGaussPts(field_name, commonData));
753 feConvectionRhs.getOpPtrVector().push_back(new OpConvectionRhs(
754 field_name, F, sit->second, commonData, hoGeometry));
755 }
757}
@ F
CommonData commonData

Member Data Documentation

◆ commonData

CommonData ThermalElement::commonData

Definition at line 152 of file ThermalElement.hpp.

◆ feConvectionLhs

MyTriFE ThermalElement::feConvectionLhs

Definition at line 70 of file ThermalElement.hpp.

◆ feConvectionRhs

MyTriFE ThermalElement::feConvectionRhs

Definition at line 69 of file ThermalElement.hpp.

◆ feFlux

MyTriFE ThermalElement::feFlux

Definition at line 66 of file ThermalElement.hpp.

◆ feLhs

MyVolumeFE ThermalElement::feLhs

Definition at line 53 of file ThermalElement.hpp.

◆ feRadiationLhs

MyTriFE ThermalElement::feRadiationLhs

Definition at line 77 of file ThermalElement.hpp.

◆ feRadiationRhs

MyTriFE ThermalElement::feRadiationRhs

Definition at line 76 of file ThermalElement.hpp.

◆ feRhs

MyVolumeFE ThermalElement::feRhs

cauclate right hand side for tetrahedral elements

Definition at line 51 of file ThermalElement.hpp.

◆ mField

MoFEM::Interface& ThermalElement::mField

Definition at line 83 of file ThermalElement.hpp.

◆ setOfBlocks

std::map<int, BlockData> ThermalElement::setOfBlocks

maps block set id with appropriate BlockData

Definition at line 100 of file ThermalElement.hpp.

◆ setOfConvection

std::map<int, ConvectionData> ThermalElement::setOfConvection

Definition at line 124 of file ThermalElement.hpp.

◆ setOfFluxes

std::map<int, FluxData> ThermalElement::setOfFluxes

maps side set id with appropriate FluxData

Definition at line 111 of file ThermalElement.hpp.

◆ setOfRadiation

std::map<int, RadiationData> ThermalElement::setOfRadiation

Definition at line 139 of file ThermalElement.hpp.


The documentation for this struct was generated from the following files: