v0.14.0
Loading...
Searching...
No Matches
thermo_elastic.cpp File Reference
#include <MoFEM.hpp>
#include <ThermalConvection.hpp>
#include <ThermalRadiation.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 DomainEle = PipelineManager::ElementsAndOpsByDim<SPACE_DIM>::DomainEle
 
using BoundaryEle
 
using PostProcEle = PostProcBrokenMeshInMoab<DomainEle>
 
using SkinPostProcEle = PostProcBrokenMeshInMoab<BoundaryEle>
 
using SideEle = PipelineManager::ElementsAndOpsByDim<SPACE_DIM>::FaceSideEle
 
using AssemblyDomainEleOp
 
using OpKCauchy
 [Linear elastic problem]
 
using OpInternalForceCauchy
 
using OpHdivHdiv
 [Linear elastic problem]
 
using OpHdivT
 Integrate Lhs div of base of flux time base of temperature (FLUX x T) and transpose of it, i.e. (T x FLAX)
 
using OpCapacity
 Integrate Lhs base of temperature times (heat capacity) times base of temperature (T x T)
 
using OpHdivFlux
 Integrating Rhs flux base (1/k) flux (FLUX)
 
using OpHDivTemp
 Integrate Rhs div flux base times temperature (T)
 
using OpBaseDotT
 Integrate Rhs base of temperature time heat capacity times heat rate (T)
 
using OpBaseDivFlux = OpBaseDotT
 Integrate Rhs base of temperature times divergence of flux (T)
 
using DomainNaturalBCRhs
 [Thermal problem]
 
using OpBodyForce
 
using OpHeatSource
 
using DomainNaturalBCLhs
 
using BoundaryNaturalBC
 [Body and heat source]
 
using OpForce = BoundaryNaturalBC::OpFlux<NaturalForceMeshsets, 1, SPACE_DIM>
 
using OpEssentialFluxRhs
 [Natural boundary conditions]
 
using OpEssentialFluxLhs
 
using OpSetTemperatureRhs
 
using OpSetTemperatureLhs
 

Functions

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

Variables

constexpr int SPACE_DIM
 
double default_young_modulus = 1
 [Essential boundary conditions (Least square approach)]
 
double default_poisson_ratio = 0.25
 
double ref_temp = 0.0
 
double init_temp = 0.0
 
PetscBool is_plane_strain = PETSC_FALSE
 
double default_coeff_expansion = 1
 
double default_heat_conductivity
 
double default_heat_capacity = 1
 
int order_temp = 2
 
int order_flux = 3
 
int order_disp = 3
 
int atom_test = 0
 
int save_every = 1
 
PetscBool do_output_domain
 
PetscBool do_output_skin
 
auto save_range
 
static char help [] = "...\n\n"
 [Solve]
 

Macro Definition Documentation

◆ EXECUTABLE_DIMENSION

#define EXECUTABLE_DIMENSION   3

Definition at line 10 of file thermo_elastic.cpp.

Typedef Documentation

◆ AssemblyDomainEleOp

◆ BoundaryEle

◆ BoundaryNaturalBC

Initial value:

[Body and heat source]

[Natural boundary conditions]

Definition at line 110 of file thermo_elastic.cpp.

◆ DomainEle

◆ DomainNaturalBCLhs

◆ DomainNaturalBCRhs

Initial value:

[Thermal problem]

[Body and heat source]

Definition at line 99 of file thermo_elastic.cpp.

◆ OpBaseDivFlux

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

Examples
thermo_elastic.cpp.

Definition at line 94 of file thermo_elastic.cpp.

◆ OpBaseDotT

using OpBaseDotT
Initial value:

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

Examples
thermo_elastic.cpp.

Definition at line 87 of file thermo_elastic.cpp.

◆ OpBodyForce

◆ OpCapacity

using OpCapacity
Initial value:
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMass
[Only used with Hooke equation (linear material model)]
Definition seepage.cpp:56

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

Examples
thermo_elastic.cpp.

Definition at line 66 of file thermo_elastic.cpp.

◆ OpEssentialFluxLhs

Initial value:
GAUSS>::OpEssentialLhs<HeatFluxCubitBcData, 3, SPACE_DIM>

Definition at line 120 of file thermo_elastic.cpp.

◆ OpEssentialFluxRhs

Initial value:
GAUSS>::OpEssentialRhs<HeatFluxCubitBcData, 3, SPACE_DIM>

[Natural boundary conditions]

[Essential boundary conditions (Least square approach)]

Definition at line 117 of file thermo_elastic.cpp.

◆ OpForce

◆ OpHdivFlux

using OpHdivFlux
Initial value:
GAUSS>::OpBaseTimesVector<3, SPACE_DIM, 1>

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

Examples
test_broken_space.cpp, and thermo_elastic.cpp.

Definition at line 72 of file thermo_elastic.cpp.

◆ OpHdivHdiv

using OpHdivHdiv
Initial value:

[Linear elastic problem]

[Thermal problem]

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

Examples
phase.cpp, test_broken_space.cpp, and thermo_elastic.cpp.

Definition at line 50 of file thermo_elastic.cpp.

◆ OpHdivT

using OpHdivT
Initial value:

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 58 of file thermo_elastic.cpp.

◆ OpHDivTemp

using OpHDivTemp
Initial value:

Integrate Rhs div flux base times temperature (T)

Examples
thermo_elastic.cpp.

Definition at line 79 of file thermo_elastic.cpp.

◆ OpHeatSource

using OpHeatSource
Initial value:

Definition at line 103 of file thermo_elastic.cpp.

◆ OpInternalForceCauchy

◆ OpKCauchy

using OpKCauchy
Initial value:

[Linear elastic problem]

Examples
PlasticOps.hpp, and thermo_elastic.cpp.

Definition at line 36 of file thermo_elastic.cpp.

◆ OpSetTemperatureLhs

◆ OpSetTemperatureRhs

◆ PostProcEle

Definition at line 28 of file thermo_elastic.cpp.

◆ SideEle

◆ SkinPostProcEle

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 1543 of file thermo_elastic.cpp.

1543 {
1544
1545 // Initialisation of MoFEM/PETSc and MOAB data structures
1546 const char param_file[] = "param_file.petsc";
1547 MoFEM::Core::Initialize(&argc, &argv, param_file, help);
1548
1549 // Add logging channel for example
1550 auto core_log = logging::core::get();
1551 core_log->add_sink(
1553 LogManager::setLog("ThermoElastic");
1554 MOFEM_LOG_TAG("ThermoElastic", "ThermoElastic");
1555
1556 core_log->add_sink(
1557 LogManager::createSink(LogManager::getStrmSync(), "ThermoElasticSync"));
1558 LogManager::setLog("ThermoElasticSync");
1559 MOFEM_LOG_TAG("ThermoElasticSync", "ThermoElasticSync");
1560
1561 try {
1562
1563 //! [Register MoFEM discrete manager in PETSc]
1564 DMType dm_name = "DMMOFEM";
1565 CHKERR DMRegister_MoFEM(dm_name);
1566 //! [Register MoFEM discrete manager in PETSc
1567
1568 //! [Create MoAB]
1569 moab::Core mb_instance; ///< mesh database
1570 moab::Interface &moab = mb_instance; ///< mesh database interface
1571 //! [Create MoAB]
1572
1573 //! [Create MoFEM]
1574 MoFEM::Core core(moab); ///< finite element database
1575 MoFEM::Interface &m_field = core; ///< finite element database interface
1576 //! [Create MoFEM]
1577
1578 //! [Load mesh]
1579 Simple *simple = m_field.getInterface<Simple>();
1580 CHKERR simple->getOptions();
1581 CHKERR simple->loadFile();
1582 //! [Load mesh]
1583
1584 //! [ThermoElasticProblem]
1585 ThermoElasticProblem ex(m_field);
1586 CHKERR ex.runProblem();
1587 //! [ThermoElasticProblem]
1588 }
1590
1592}
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
#define CATCH_ERRORS
Catch errors.
#define CHKERR
Inline error check.
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
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.
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
Simple interface for fast problem set-up.
Definition Simple.hpp:27
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
static char help[]
[Solve]

Variable Documentation

◆ atom_test

int atom_test = 0
Examples
thermo_elastic.cpp.

Definition at line 143 of file thermo_elastic.cpp.

◆ default_coeff_expansion

double default_coeff_expansion = 1
Examples
thermo_elastic.cpp.

Definition at line 132 of file thermo_elastic.cpp.

◆ default_heat_capacity

double default_heat_capacity = 1
Examples
thermo_elastic.cpp.

Definition at line 136 of file thermo_elastic.cpp.

◆ default_heat_conductivity

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

Definition at line 133 of file thermo_elastic.cpp.

◆ default_poisson_ratio

double default_poisson_ratio = 0.25
Examples
seepage.cpp, and thermo_elastic.cpp.

Definition at line 126 of file thermo_elastic.cpp.

◆ default_young_modulus

double default_young_modulus = 1

[Essential boundary conditions (Least square approach)]

Examples
seepage.cpp, and thermo_elastic.cpp.

Definition at line 125 of file thermo_elastic.cpp.

◆ do_output_domain

PetscBool do_output_domain
Examples
thermo_elastic.cpp.

Definition at line 146 of file thermo_elastic.cpp.

◆ do_output_skin

PetscBool do_output_skin
Examples
thermo_elastic.cpp.

Definition at line 147 of file thermo_elastic.cpp.

◆ help

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

[Solve]

Definition at line 1541 of file thermo_elastic.cpp.

◆ init_temp

double init_temp = 0.0

◆ is_plane_strain

PetscBool is_plane_strain = PETSC_FALSE
Examples
thermo_elastic.cpp.

Definition at line 130 of file thermo_elastic.cpp.

◆ order_disp

int order_disp = 3
Examples
thermo_elastic.cpp.

Definition at line 141 of file thermo_elastic.cpp.

◆ order_flux

int order_flux = 3
Examples
thermo_elastic.cpp.

Definition at line 140 of file thermo_elastic.cpp.

◆ order_temp

int order_temp = 2
Examples
thermo_elastic.cpp.

Definition at line 139 of file thermo_elastic.cpp.

◆ ref_temp

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

Definition at line 127 of file thermo_elastic.cpp.

◆ save_every

int save_every = 1
Examples
thermo_elastic.cpp.

Definition at line 145 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 ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
auto get_temp_meshset_ptr(moab::Interface &moab)
Create smart pointer to temporary meshset.
Examples
EshelbianPlasticity.cpp, free_surface.cpp, level_set.cpp, seepage.cpp, and thermo_elastic.cpp.

Definition at line 157 of file thermo_elastic.cpp.

158 {
160 auto out_meshset = get_temp_meshset_ptr(moab);
161 CHKERR moab.add_entities(*out_meshset, r);
162 CHKERR moab.write_file(name.c_str(), "VTK", "", out_meshset->get_ptr(), 1);
164};

◆ SPACE_DIM

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

Definition at line 19 of file thermo_elastic.cpp.