v0.13.2
Loading...
Searching...
No Matches
Functions | Variables
thermal_steady.cpp File Reference

Example of steady thermal analysis. More...

#include <BasicFiniteElements.hpp>

Go to the source code of this file.

Functions

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

Variables

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

Detailed Description

Example of steady thermal analysis.

TODO:

Todo:
Make it work in distributed meshes with multigird solver. At the moment it is not working efficient as can.

Definition in file thermal_steady.cpp.

Function Documentation

◆ main()

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

Definition at line 21 of file thermal_steady.cpp.

21 {
22
23 const string default_options = "-ksp_type fgmres \n"
24 "-pc_type lu \n"
25 "-pc_factor_mat_solver_type mumps \n"
26 "-mat_mumps_icntl_20 0 \n"
27 "-ksp_monitor \n";
28
29 string param_file = "param_file.petsc";
30 if (!static_cast<bool>(ifstream(param_file))) {
31 std::ofstream file(param_file.c_str(), std::ios::ate);
32 if (file.is_open()) {
33 file << default_options;
34 file.close();
35 }
36 }
37
38 MoFEM::Core::Initialize(&argc, &argv, param_file.c_str(), help);
39
40 try {
41
42 moab::Core mb_instance;
43 moab::Interface &moab = mb_instance;
44 int rank;
45 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
46
47 PetscBool flg = PETSC_TRUE;
48 char mesh_file_name[255];
49 CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
50 mesh_file_name, 255, &flg);
51 if (flg != PETSC_TRUE) {
52 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
53 "*** ERROR -my_file (MESH FILE NEEDED)");
54 }
55
56 const char *option;
57 option = ""; //"PARALLEL=BCAST;";//;DEBUG_IO";
58 CHKERR moab.load_file(mesh_file_name, 0, option);
59
60 // Create MoFEM (Joseph) database
61 MoFEM::Core core(moab);
62 MoFEM::Interface &m_field = core;
63
64 // set entities bit level
65 BitRefLevel bit_level0;
66 bit_level0.set(0);
67 EntityHandle meshset_level0;
68 CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
69 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
70 0, 3, bit_level0);
71
72 // Fields
73 CHKERR m_field.add_field("TEMP", H1, AINSWORTH_LEGENDRE_BASE, 1);
74
75 // Problem
76 CHKERR m_field.add_problem("THERMAL_PROBLEM");
77
78 // set refinement level for problem
79 CHKERR m_field.modify_problem_ref_level_add_bit("THERMAL_PROBLEM",
80 bit_level0);
81
82 // meshset consisting all entities in mesh
83 EntityHandle root_set = moab.get_root_set();
84 // add entities to field
85 CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "TEMP");
86
87 // set app. order
88 // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
89 // (Mark Ainsworth & Joe Coyle)
90 PetscInt order;
91 CHKERR PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-my_order", &order,
92 &flg);
93 if (flg != PETSC_TRUE) {
94 order = 2;
95 }
96
97 CHKERR m_field.set_field_order(root_set, MBTET, "TEMP", order);
98 CHKERR m_field.set_field_order(root_set, MBTRI, "TEMP", order);
99 CHKERR m_field.set_field_order(root_set, MBEDGE, "TEMP", order);
100 CHKERR m_field.set_field_order(root_set, MBVERTEX, "TEMP", 1);
101
102 CHKERR m_field.add_field("MESH_NODE_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE,
103 3);
104 CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET,
105 "MESH_NODE_POSITIONS");
106 CHKERR m_field.set_field_order(0, MBTET, "MESH_NODE_POSITIONS", 2);
107 CHKERR m_field.set_field_order(0, MBTRI, "MESH_NODE_POSITIONS", 2);
108 CHKERR m_field.set_field_order(0, MBEDGE, "MESH_NODE_POSITIONS", 2);
109 CHKERR m_field.set_field_order(0, MBVERTEX, "MESH_NODE_POSITIONS", 1);
110
111 ThermalElement thermal_elements(m_field);
112 CHKERR thermal_elements.addThermalElements("TEMP");
113 CHKERR thermal_elements.addThermalFluxElement("TEMP");
114 CHKERR thermal_elements.addThermalConvectionElement("TEMP");
115
116 CHKERR m_field.modify_problem_add_finite_element("THERMAL_PROBLEM",
117 "THERMAL_FE");
118 CHKERR m_field.modify_problem_add_finite_element("THERMAL_PROBLEM",
119 "THERMAL_FLUX_FE");
120 CHKERR m_field.modify_problem_add_finite_element("THERMAL_PROBLEM",
121 "THERMAL_CONVECTION_FE");
122
123 /****/
124 // build database
125 // build field
126 CHKERR m_field.build_fields();
127 // build finite elemnts
129 // build adjacencies
130 CHKERR m_field.build_adjacencies(bit_level0);
131
132 ProblemsManager *prb_mng_ptr;
133 CHKERR m_field.getInterface(prb_mng_ptr);
134 // build problem
135 CHKERR prb_mng_ptr->buildProblem("THERMAL_PROBLEM", true);
136
137 Projection10NodeCoordsOnField ent_method_material(m_field,
138 "MESH_NODE_POSITIONS");
139 CHKERR m_field.loop_dofs("MESH_NODE_POSITIONS", ent_method_material);
140
141 /****/
142 // mesh partitioning
143 // partition
144 CHKERR prb_mng_ptr->partitionProblem("THERMAL_PROBLEM");
145 CHKERR prb_mng_ptr->partitionFiniteElements("THERMAL_PROBLEM");
146 // what are ghost nodes, see Petsc Manual
147 CHKERR prb_mng_ptr->partitionGhostDofs("THERMAL_PROBLEM");
148
149 Vec F;
150 CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("THERMAL_PROBLEM",
151 ROW, &F);
152 Vec T;
153 CHKERR VecDuplicate(F, &T);
154 Mat A;
156 ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>("THERMAL_PROBLEM", &A);
157
158 DirichletTemperatureBc my_dirichlet_bc(m_field, "TEMP", A, T, F);
159
160 CHKERR addHOOpsVol("MESH_NODE_POSITIONS", thermal_elements.getLoopFeRhs(),
161 true, false, false, false);
162 CHKERR addHOOpsVol("MESH_NODE_POSITIONS", thermal_elements.getLoopFeLhs(),
163 true, false, false, false);
164
165 CHKERR thermal_elements.setThermalFiniteElementRhsOperators("TEMP", F);
166 CHKERR thermal_elements.setThermalFiniteElementLhsOperators("TEMP", A);
167 CHKERR thermal_elements.setThermalFluxFiniteElementRhsOperators("TEMP", F);
168 CHKERR thermal_elements.setThermalConvectionFiniteElementLhsOperators(
169 "TEMP", A);
170 CHKERR thermal_elements.setThermalConvectionFiniteElementRhsOperators(
171 "TEMP", F);
172
173 CHKERR VecZeroEntries(T);
174 CHKERR VecGhostUpdateBegin(T, INSERT_VALUES, SCATTER_FORWARD);
175 CHKERR VecGhostUpdateEnd(T, INSERT_VALUES, SCATTER_FORWARD);
176 CHKERR VecZeroEntries(F);
177 CHKERR VecGhostUpdateBegin(F, INSERT_VALUES, SCATTER_FORWARD);
178 CHKERR VecGhostUpdateEnd(F, INSERT_VALUES, SCATTER_FORWARD);
179 CHKERR MatZeroEntries(A);
180
181 // preproc
182 CHKERR m_field.problem_basic_method_preProcess("THERMAL_PROBLEM",
183 my_dirichlet_bc);
184 CHKERR m_field.getInterface<VecManager>()->setGlobalGhostVector(
185 "THERMAL_PROBLEM", ROW, T, INSERT_VALUES, SCATTER_REVERSE);
186
187 CHKERR m_field.loop_finite_elements("THERMAL_PROBLEM", "THERMAL_FE",
188 thermal_elements.getLoopFeRhs());
189 CHKERR m_field.loop_finite_elements("THERMAL_PROBLEM", "THERMAL_FE",
190 thermal_elements.getLoopFeLhs());
191 CHKERR m_field.loop_finite_elements("THERMAL_PROBLEM", "THERMAL_FLUX_FE",
192 thermal_elements.getLoopFeFlux());
194 "THERMAL_PROBLEM", "THERMAL_CONVECTION_FE",
195 thermal_elements.getLoopFeConvectionRhs());
197 "THERMAL_PROBLEM", "THERMAL_CONVECTION_FE",
198 thermal_elements.getLoopFeConvectionLhs());
199
200 // postproc
201 CHKERR m_field.problem_basic_method_postProcess("THERMAL_PROBLEM",
202 my_dirichlet_bc);
203
204 CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
205 CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
206 CHKERR VecAssemblyBegin(F);
207 CHKERR VecAssemblyEnd(F);
208 CHKERR MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
209 CHKERR MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
210
211 CHKERR VecScale(F, -1);
212
213 // Solver
214 KSP solver;
215 CHKERR KSPCreate(PETSC_COMM_WORLD, &solver);
216 CHKERR KSPSetOperators(solver, A, A);
217 CHKERR KSPSetFromOptions(solver);
218 CHKERR KSPSetUp(solver);
219
220 CHKERR KSPSolve(solver, F, T);
221 CHKERR VecGhostUpdateBegin(T, INSERT_VALUES, SCATTER_FORWARD);
222 CHKERR VecGhostUpdateEnd(T, INSERT_VALUES, SCATTER_FORWARD);
223
224 CHKERR m_field.problem_basic_method_preProcess("THERMAL_PROBLEM",
225 my_dirichlet_bc);
226
227 // Save data on mesh
228 CHKERR m_field.getInterface<VecManager>()->setGlobalGhostVector(
229 "THERMAL_PROBLEM", ROW, T, INSERT_VALUES, SCATTER_REVERSE);
230
231 if (m_field.get_comm_rank() == 0) {
232 CHKERR moab.write_file("solution.h5m");
233 }
234
235 ProjectionFieldOn10NodeTet ent_method_on_10nodeTet(m_field, "TEMP", true,
236 false, "TEMP");
237 CHKERR m_field.loop_dofs("TEMP", ent_method_on_10nodeTet);
238 ent_method_on_10nodeTet.setNodes = false;
239 CHKERR m_field.loop_dofs("TEMP", ent_method_on_10nodeTet);
240
241 if (m_field.get_comm_rank() == 0) {
242 EntityHandle out_meshset;
243 CHKERR moab.create_meshset(MESHSET_SET, out_meshset);
245 "THERMAL_PROBLEM", "THERMAL_FE", out_meshset);
246 CHKERR moab.write_file("out.vtk", "VTK", "", &out_meshset, 1);
247 CHKERR moab.delete_entities(&out_meshset, 1);
248 }
249
250 CHKERR MatDestroy(&A);
251 CHKERR VecDestroy(&F);
252 CHKERR VecDestroy(&T);
253 CHKERR KSPDestroy(&solver);
254 }
256
257 return MoFEM::Core::Finalize();
258}
const std::string default_options
std::string param_file
@ ROW
Definition: definitions.h:123
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
@ H1
continuous field
Definition: definitions.h:85
@ MOFEM_NOT_FOUND
Definition: definitions.h:33
#define CHKERR
Inline error check.
Definition: definitions.h:535
@ F
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
virtual MoFEMErrorCode problem_basic_method_postProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
MoFEMErrorCode partitionProblem(const std::string name, int verb=VERBOSE)
partition problem dofs (collective)
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
virtual MoFEMErrorCode get_problem_finite_elements_entities(const std::string name, const std::string &fe_name, const EntityHandle meshset)=0
add finite elements to the meshset
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
char mesh_file_name[255]
const double T
const FTensor::Tensor2< T, Dim, Dim > Vec
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1, bool hcurl, bool hdiv, bool l2)
constexpr AssemblyType A
Managing BitRefLevels.
virtual MoFEMErrorCode problem_basic_method_preProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
virtual MoFEMErrorCode add_field(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
virtual int get_comm_rank() const =0
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.
Matrix manager is used to build and partition problems.
Problem manager is used to build and partition problems.
Projection of edge entities with one mid-node on hierarchical basis.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:23
structure grouping operators and data used for thermal problems
static char help[]

Variable Documentation

◆ help

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

Definition at line 19 of file thermal_steady.cpp.