v0.14.0
Loading...
Searching...
No Matches
Classes | Macros | Typedefs | Functions | Variables
thermo_elastic.cpp File Reference
#include <MoFEM.hpp>
#include <ThermoElasticOps.hpp>

Go to the source code of this file.

Classes

struct  ThermoElasticProblem
 
struct  ThermoElasticProblem::BlockedParameters
 

Macros

#define EXECUTABLE_DIMENSION   3
 

Typedefs

using EntData = EntitiesFieldData::EntData
 
using DomainEle = PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::DomainEle
 
using DomainEleOp = DomainEle::UserDataOperator
 
using BoundaryEle = PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::BoundaryEle
 
using BoundaryEleOp = BoundaryEle::UserDataOperator
 
using PostProcEle = PostProcBrokenMeshInMoab< DomainEle >
 
using AssemblyDomainEleOp = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::OpBase
 
using OpKCauchy = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradSymTensorGrad< 1, SPACE_DIM, SPACE_DIM, 0 >
 [Linear elastic problem] More...
 
using OpInternalForceCauchy = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesSymTensor< 1, SPACE_DIM, SPACE_DIM >
 
using OpHdivHdiv = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 3, SPACE_DIM >
 [Linear elastic problem] More...
 
using OpHdivT = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMixDivTimesScalar< SPACE_DIM >
 Integrate Lhs div of base of flux time base of temperature (FLUX x T) and transpose of it, i.e. (T x FLAX) More...
 
using OpCapacity = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, 1 >
 Integrate Lhs base of temperature times (heat capacity) times base of temperature (T x T) More...
 
using OpHdivFlux = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesVector< 3, SPACE_DIM, 1 >
 Integrating Rhs flux base (1/k) flux (FLUX) More...
 
using OpHDivTemp = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpMixDivTimesU< 3, 1, SPACE_DIM >
 Integrate Rhs div flux base times temperature (T) More...
 
using OpBaseDotT = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesScalar< 1 >
 Integrate Rhs base of temperature time heat capacity times heat rate (T) More...
 
using OpBaseDivFlux = OpBaseDotT
 Integrate Rhs base of temperature times divergence of flux (T) More...
 
using DomainNaturalBCRhs = NaturalBC< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >
 [Thermal problem] More...
 
using OpBodyForce = DomainNaturalBCRhs::OpFlux< NaturalMeshsetType< BLOCKSET >, 1, SPACE_DIM >
 
using OpHeatSource = DomainNaturalBCRhs::OpFlux< NaturalMeshsetType< BLOCKSET >, 1, 1 >
 
using DomainNaturalBCLhs = NaturalBC< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >
 
using BoundaryNaturalBC = NaturalBC< BoundaryEleOp >::Assembly< PETSC >::LinearForm< GAUSS >
 [Body and heat source] More...
 
using OpForce = BoundaryNaturalBC::OpFlux< NaturalForceMeshsets, 1, SPACE_DIM >
 
using OpTemperatureBC = BoundaryNaturalBC::OpFlux< NaturalTemperatureMeshsets, 3, SPACE_DIM >
 
using OpEssentialFluxRhs = EssentialBC< BoundaryEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpEssentialRhs< HeatFluxCubitBcData, 3, SPACE_DIM >
 [Natural boundary conditions] More...
 
using OpEssentialFluxLhs = EssentialBC< BoundaryEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpEssentialLhs< HeatFluxCubitBcData, 3, SPACE_DIM >
 
using OpSetTemperatureRhs = DomainNaturalBCRhs::OpFlux< SetTargetTemperature, 1, 1 >
 
using OpSetTemperatureLhs = DomainNaturalBCLhs::OpFlux< SetTargetTemperature, 1, 1 >
 

Functions

int main (int argc, char *argv[])
 

Variables

constexpr int SPACE_DIM
 
double default_young_modulus = 1
 [Essential boundary conditions (Least square approach)] More...
 
double default_poisson_ratio = 0.25
 
double ref_temp = 0.0
 
double default_coeff_expansion = 1
 
double default_heat_conductivity
 
double default_heat_capacity = 1
 
int order = 2
 
auto save_range
 
static char help [] = "...\n\n"
 [Solve] More...
 

Macro Definition Documentation

◆ EXECUTABLE_DIMENSION

#define EXECUTABLE_DIMENSION   3

Definition at line 10 of file thermo_elastic.cpp.

Typedef Documentation

◆ AssemblyDomainEleOp

Definition at line 28 of file thermo_elastic.cpp.

◆ BoundaryEle

Definition at line 23 of file thermo_elastic.cpp.

◆ BoundaryEleOp

Definition at line 25 of file thermo_elastic.cpp.

◆ BoundaryNaturalBC

using BoundaryNaturalBC = NaturalBC<BoundaryEleOp>::Assembly<PETSC>::LinearForm<GAUSS>

[Body and heat source]

[Natural boundary conditions]

Definition at line 106 of file thermo_elastic.cpp.

◆ DomainEle

Definition at line 21 of file thermo_elastic.cpp.

◆ DomainEleOp

Definition at line 22 of file thermo_elastic.cpp.

◆ DomainNaturalBCLhs

using DomainNaturalBCLhs = NaturalBC<DomainEleOp>::Assembly<PETSC>::BiLinearForm<GAUSS>

Definition at line 101 of file thermo_elastic.cpp.

◆ DomainNaturalBCRhs

using DomainNaturalBCRhs = NaturalBC<DomainEleOp>::Assembly<PETSC>::LinearForm<GAUSS>

[Thermal problem]

[Body and heat source]

Definition at line 95 of file thermo_elastic.cpp.

◆ EntData

Definition at line 20 of file thermo_elastic.cpp.

◆ OpBaseDivFlux

Integrate Rhs base of temperature times divergence of flux (T)

Examples
thermo_elastic.cpp.

Definition at line 90 of file thermo_elastic.cpp.

◆ OpBaseDotT

using OpBaseDotT = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::LinearForm< GAUSS>::OpBaseTimesScalar<1>

Integrate Rhs base of temperature time heat capacity times heat rate (T)

Examples
thermo_elastic.cpp.

Definition at line 83 of file thermo_elastic.cpp.

◆ OpBodyForce

Definition at line 97 of file thermo_elastic.cpp.

◆ OpCapacity

using OpCapacity = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::BiLinearForm< GAUSS>::OpMass<1, 1>

Integrate Lhs base of temperature times (heat capacity) times base of temperature (T x T)

Examples
thermo_elastic.cpp.

Definition at line 62 of file thermo_elastic.cpp.

◆ OpEssentialFluxLhs

using OpEssentialFluxLhs = EssentialBC<BoundaryEleOp>::Assembly<PETSC>::BiLinearForm< GAUSS>::OpEssentialLhs<HeatFluxCubitBcData, 3, SPACE_DIM>

Definition at line 117 of file thermo_elastic.cpp.

◆ OpEssentialFluxRhs

using OpEssentialFluxRhs = EssentialBC<BoundaryEleOp>::Assembly<PETSC>::LinearForm< GAUSS>::OpEssentialRhs<HeatFluxCubitBcData, 3, SPACE_DIM>

[Natural boundary conditions]

[Essential boundary conditions (Least square approach)]

Definition at line 114 of file thermo_elastic.cpp.

◆ OpForce

Definition at line 108 of file thermo_elastic.cpp.

◆ OpHdivFlux

using OpHdivFlux = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::LinearForm< GAUSS>::OpBaseTimesVector<3, SPACE_DIM, 1>

Integrating Rhs flux base (1/k) flux (FLUX)

Examples
thermo_elastic.cpp.

Definition at line 68 of file thermo_elastic.cpp.

◆ OpHdivHdiv

[Linear elastic problem]

[Thermal problem]

Integrate Lhs base of flux (1/k) base of flux (FLUX x FLUX)

Examples
phase.cpp, and thermo_elastic.cpp.

Definition at line 46 of file thermo_elastic.cpp.

◆ OpHdivT

using OpHdivT = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::BiLinearForm< GAUSS>::OpMixDivTimesScalar<SPACE_DIM>

Integrate Lhs div of base of flux time base of temperature (FLUX x T) and transpose of it, i.e. (T x FLAX)

Examples
thermo_elastic.cpp.

Definition at line 54 of file thermo_elastic.cpp.

◆ OpHDivTemp

using OpHDivTemp = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::LinearForm< GAUSS>::OpMixDivTimesU<3, 1, SPACE_DIM>

Integrate Rhs div flux base times temperature (T)

Examples
thermo_elastic.cpp.

Definition at line 75 of file thermo_elastic.cpp.

◆ OpHeatSource

Definition at line 99 of file thermo_elastic.cpp.

◆ OpInternalForceCauchy

using OpInternalForceCauchy = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::LinearForm< GAUSS>::OpGradTimesSymTensor<1, SPACE_DIM, SPACE_DIM>
Examples
PlasticOps.hpp, and thermo_elastic.cpp.

Definition at line 35 of file thermo_elastic.cpp.

◆ OpKCauchy

using OpKCauchy = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::BiLinearForm< GAUSS>::OpGradSymTensorGrad<1, SPACE_DIM, SPACE_DIM, 0>

[Linear elastic problem]

Examples
PlasticOps.hpp, and thermo_elastic.cpp.

Definition at line 32 of file thermo_elastic.cpp.

◆ OpSetTemperatureLhs

using OpSetTemperatureLhs = DomainNaturalBCLhs::OpFlux<SetTargetTemperature, 1, 1>

Definition at line 140 of file thermo_elastic.cpp.

◆ OpSetTemperatureRhs

using OpSetTemperatureRhs = DomainNaturalBCRhs::OpFlux<SetTargetTemperature, 1, 1>

Definition at line 138 of file thermo_elastic.cpp.

◆ OpTemperatureBC

Definition at line 109 of file thermo_elastic.cpp.

◆ PostProcEle

Definition at line 26 of file thermo_elastic.cpp.

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

[Register MoFEM discrete manager in PETSc]

[Register MoFEM discrete manager in PETSc

[Create MoAB]

< mesh database

< mesh database interface

[Create MoAB]

[Create MoFEM]

< finite element database

< finite element database interface

[Create MoFEM]

[Load mesh]

[Load mesh]

[ThermoElasticProblem]

[ThermoElasticProblem]

Definition at line 1066 of file thermo_elastic.cpp.

1066 {
1067
1068 // Initialisation of MoFEM/PETSc and MOAB data structures
1069 const char param_file[] = "param_file.petsc";
1070 MoFEM::Core::Initialize(&argc, &argv, param_file, help);
1071
1072 // Add logging channel for example
1073 auto core_log = logging::core::get();
1074 core_log->add_sink(
1076 LogManager::setLog("ThermoElastic");
1077 MOFEM_LOG_TAG("ThermoElastic", "ThermoElastic");
1078
1079 core_log->add_sink(
1080 LogManager::createSink(LogManager::getStrmSync(), "ThermoElasticSync"));
1081 LogManager::setLog("ThermoElasticSync");
1082 MOFEM_LOG_TAG("ThermoElasticSync", "ThermoElasticSync");
1083
1084 try {
1085
1086 //! [Register MoFEM discrete manager in PETSc]
1087 DMType dm_name = "DMMOFEM";
1088 CHKERR DMRegister_MoFEM(dm_name);
1089 //! [Register MoFEM discrete manager in PETSc
1090
1091 //! [Create MoAB]
1092 moab::Core mb_instance; ///< mesh database
1093 moab::Interface &moab = mb_instance; ///< mesh database interface
1094 //! [Create MoAB]
1095
1096 //! [Create MoFEM]
1097 MoFEM::Core core(moab); ///< finite element database
1098 MoFEM::Interface &m_field = core; ///< finite element database interface
1099 //! [Create MoFEM]
1100
1101 //! [Load mesh]
1102 Simple *simple = m_field.getInterface<Simple>();
1103 CHKERR simple->getOptions();
1104 CHKERR simple->loadFile();
1105 //! [Load mesh]
1106
1107 //! [ThermoElasticProblem]
1108 ThermoElasticProblem ex(m_field);
1109 CHKERR ex.runProblem();
1110 //! [ThermoElasticProblem]
1111 }
1113
1115}
std::string param_file
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
#define CHKERR
Inline error check.
Definition: definitions.h:535
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
Core (interface) class.
Definition: Core.hpp:82
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
Deprecated interface functions.
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:298
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Definition: LogManager.cpp:344
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
Definition: LogManager.cpp:348
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
static char help[]
[Solve]

Variable Documentation

◆ default_coeff_expansion

double default_coeff_expansion = 1
Examples
thermo_elastic.cpp.

Definition at line 126 of file thermo_elastic.cpp.

◆ default_heat_capacity

double default_heat_capacity = 1
Examples
thermo_elastic.cpp.

Definition at line 130 of file thermo_elastic.cpp.

◆ default_heat_conductivity

double default_heat_conductivity
Initial value:
=
1
Examples
thermo_elastic.cpp.

Definition at line 127 of file thermo_elastic.cpp.

◆ default_poisson_ratio

double default_poisson_ratio = 0.25
Examples
thermo_elastic.cpp.

Definition at line 123 of file thermo_elastic.cpp.

◆ default_young_modulus

double default_young_modulus = 1

[Essential boundary conditions (Least square approach)]

Examples
thermo_elastic.cpp.

Definition at line 122 of file thermo_elastic.cpp.

◆ help

char help[] = "...\n\n"
static

[Solve]

Definition at line 1064 of file thermo_elastic.cpp.

◆ order

int order = 2
Examples
thermo_elastic.cpp.

Definition at line 133 of file thermo_elastic.cpp.

◆ ref_temp

double ref_temp = 0.0
Examples
ThermoElasticOps.hpp, and thermo_elastic.cpp.

Definition at line 124 of file thermo_elastic.cpp.

◆ save_range

auto save_range
Initial value:
= [](moab::Interface &moab, const std::string name,
const Range r) {
auto out_meshset = get_temp_meshset_ptr(moab);
CHKERR moab.add_entities(*out_meshset, r);
CHKERR moab.write_file(name.c_str(), "VTK", "", out_meshset->get_ptr(), 1);
}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
auto get_temp_meshset_ptr(moab::Interface &moab)
Create smart pointer to temporary meshset.
Definition: Templates.hpp:1906
Examples
free_surface.cpp, level_set.cpp, and thermo_elastic.cpp.

Definition at line 143 of file thermo_elastic.cpp.

◆ SPACE_DIM

constexpr int SPACE_DIM
constexpr
Initial value:
=
#define EXECUTABLE_DIMENSION
Examples
thermo_elastic.cpp.

Definition at line 17 of file thermo_elastic.cpp.