v0.13.2
Loading...
Searching...
No Matches
thermal_unsteady.cpp
Go to the documentation of this file.
1/** \file thermal_unsteady.cpp
2 \ingroup mofem_thermal_elem
3 \brief Example of thermal unsteady analyze.
4
5 TODO:
6 \todo Make it work in distributed meshes with multigird solver. At the moment
7 it is not working efficient as can.
8*/
9
10
11
13using namespace MoFEM;
14
15#ifdef __GROUND_SURFACE_TEMPERATURE_HPP
16
18 #include <GroundSurfaceTemperature.hpp>
19
20 #include <time.h>
21 extern "C" {
22 #include <spa.h>
23 }
24 #include <CrudeClimateModel.hpp>
25
26#endif // __GROUND_SURFACE_TEMPERATURE_HPP
27
30
31static char help[] =
32 "-my_file mesh file\n"
33 "-order set approx. order to all blocks\n"
34 "-my_block_config set block data\n"
35 "-my_ground_analysis_data data for crude climate model\n"
36 "\n";
37
38struct BlockOptionData {
39 int oRder;
41 double cApacity;
42 double initTemp;
44 oRder(-1),
45 cOnductivity(-1),
46 cApacity(-1),
47 initTemp(0) {}
48};
49
50struct MonitorPostProc : public FEMethod {
51
55
56 bool iNit;
57 int pRT;
58
60 : FEMethod(), mField(m_field), postProc(m_field), skinPostProc(m_field),
61 iNit(false) {
62
63 PetscBool flg = PETSC_TRUE;
64 CHKERR PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-my_output_prt", &pRT,
65 &flg);
66 CHKERRABORT(PETSC_COMM_WORLD, ierr);
67 if (flg != PETSC_TRUE) {
68 pRT = 1;
69 }
70 }
71
75 }
76
80 }
81
84
85 if (!iNit) {
86 CHKERR addHOOpsVol("MESH_NODE_POSITIONS", postProc, true, false, false,
87 false);
92 CHKERR postProc.addFieldValuesPostProc("MESH_NODE_POSITIONS");
93
94 CHKERR addHOOpsFace3D("MESH_NODE_POSITIONS", skinPostProc, false, false);
97
98 iNit = true;
99 }
100 int step;
101 CHKERR TSGetTimeStepNumber(ts, &step);
102
103 if (pRT && (step) % pRT == 0) {
104 // CHKERR mField.loop_finite_elements("DMTHERMAL","THERMAL_FE",postProc);
105 // std::ostringstream sss;
106 // sss << "out_thermal_" << step << ".h5m";
107 // CHKERR postProc.writeFile(sss.str().c_str());
108 CHKERR mField.loop_finite_elements("DMTHERMAL", "POST_PROC_SKIN",
110 std::ostringstream sss;
111 sss << "out_skin_" << step << ".h5m";
112 CHKERR skinPostProc.writeFile(sss.str().c_str());
113 }
115 }
116};
117
118int main(int argc, char *argv[]) {
119
120 const string default_options = "-ksp_type fgmres \n"
121 "-pc_type lu \n"
122 "-pc_factor_mat_solver_type mumps \n"
123 "-mat_mumps_icntl_20 0 \n"
124 "-ksp_monitor \n"
125 "-snes_type newtonls \n"
126 "-snes_linesearch_type basic \n"
127 "-snes_max_it 100 \n"
128 "-snes_atol 1e-8 \n"
129 "-snes_rtol 1e-8 \n"
130 "-snes_monitor \n"
131 "-ts_monitor \n"
132 "-ts_type beuler \n"
133 "-ts_exact_final_time stepover \n";
134
135 string param_file = "param_file.petsc";
136 if (!static_cast<bool>(ifstream(param_file))) {
137 std::ofstream file(param_file.c_str(), std::ios::ate);
138 if (file.is_open()) {
139 file << default_options;
140 file.close();
141 }
142 }
143
144 MoFEM::Core::Initialize(&argc, &argv, param_file.c_str(), help);
145
146 auto core_log = logging::core::get();
147 core_log->add_sink(
149 LogManager::setLog("THERMALSYNC");
150 MOFEM_LOG_TAG("THERMALSYNC", "thermal");
151
152 try {
153
154 PetscBool flg = PETSC_TRUE;
155 char mesh_file_name[255];
156 CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
157 mesh_file_name, 255, &flg);
158 if(flg != PETSC_TRUE) {
159 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
160 "*** ERROR -my_file (MESH FILE NEEDED)");
161 }
162
163 char time_data_file_for_ground_surface[255];
164 PetscBool ground_temperature_analysis = PETSC_FALSE;
165 CHKERR PetscOptionsGetString(PETSC_NULL,PETSC_NULL,"-my_ground_analysis_data",
166 time_data_file_for_ground_surface,255,&ground_temperature_analysis);
167 if(ground_temperature_analysis) {
168#ifndef WITH_ADOL_C
169 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_INSTALLED,
170 "*** ERROR to do ground thermal analysis MoFEM need to be compiled "
171 "with ADOL-C");
172#endif // WITH_ADOL_C
173 }
174
175 //create MoAB database
176 moab::Core mb_instance;
177 moab::Interface& moab = mb_instance;
178 const char *option;
179 option = "";
180 CHKERR moab.load_file(mesh_file_name, 0, option);
181 //create MoFEM database
182 MoFEM::Core core(moab);
183 MoFEM::Interface& m_field = core;
184
185 DMType dm_name = "DMTHERMAL";
186 CHKERR DMRegister_MoFEM(dm_name);
187 // create dm instance
188 DM dm;
189 CHKERR DMCreate(PETSC_COMM_WORLD, &dm);
190 CHKERR DMSetType(dm, dm_name);
191
192 //set entities bit level (this allow to set refinement levels for h-adaptivity)
193 //only one level is used in this example
194 BitRefLevel bit_level0;
195 bit_level0.set(0);
196 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(0, 3,
197 bit_level0);
198
199 //Fields H1 space rank 1
200 CHKERR m_field.add_field("TEMP", H1, AINSWORTH_LEGENDRE_BASE, 1,
201 MB_TAG_SPARSE, MF_ZERO);
202 CHKERR m_field.add_field("TEMP_RATE", H1, AINSWORTH_LEGENDRE_BASE, 1,
203 MB_TAG_SPARSE, MF_ZERO);
204
205 //Add field H1 space rank 3 to approximate geometry using hierarchical basis
206 //For 10 node tets, before use, geometry is projected on that field (see below)
207 CHKERR m_field.add_field(
208 "MESH_NODE_POSITIONS",H1,AINSWORTH_LEGENDRE_BASE,3,MB_TAG_SPARSE,MF_ZERO
209 );
210
211 //meshset consisting all entities in mesh
212 EntityHandle root_set = moab.get_root_set();
213 //add entities to field (root_mesh, i.e. on all mesh etities fields are approx.)
214 CHKERR m_field.add_ents_to_field_by_type(root_set,MBTET,"TEMP");
215 CHKERR m_field.add_ents_to_field_by_type(root_set,MBTET,"TEMP_RATE");
216
217 int order;
218 CHKERR PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-my_order", &order, &flg);
219
220 if (flg != PETSC_TRUE) {
221 order = 1;
222 }
223 // set app. order
224 // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
225 // (Mark Ainsworth & Joe Coyle) for simplicity of example to all entities is
226 // applied the same order
227 CHKERR m_field.set_field_order(root_set,MBTET,"TEMP",order);
228 CHKERR m_field.set_field_order(root_set,MBTRI,"TEMP",order);
229 CHKERR m_field.set_field_order(root_set,MBEDGE,"TEMP",order);
230 CHKERR m_field.set_field_order(root_set,MBVERTEX,"TEMP",1);
231
232 CHKERR m_field.set_field_order(root_set,MBTET,"TEMP_RATE",order);
233 CHKERR m_field.set_field_order(root_set,MBTRI,"TEMP_RATE",order);
234 CHKERR m_field.set_field_order(root_set,MBEDGE,"TEMP_RATE",order);
235 CHKERR m_field.set_field_order(root_set,MBVERTEX,"TEMP_RATE",1);
236
237 //geometry approximation is set to 2nd oreder
238 CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET,
239 "MESH_NODE_POSITIONS");
240 CHKERR m_field.set_field_order(0,MBTET,"MESH_NODE_POSITIONS",2);
241 CHKERR m_field.set_field_order(0,MBTRI,"MESH_NODE_POSITIONS",2);
242 CHKERR m_field.set_field_order(0,MBEDGE,"MESH_NODE_POSITIONS",2);
243 CHKERR m_field.set_field_order(0,MBVERTEX,"MESH_NODE_POSITIONS",1);
244
245 // configure blocks by parsing config file
246 // it allow to set approximation order for each block independently
247 PetscBool block_config;
248 char block_config_file[255];
249 CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_block_config",
250 block_config_file, 255, &block_config);
251 std::map<int,BlockOptionData> block_data;
252 bool solar_radiation = false;
253 if (block_config) {
254 try {
255 ifstream ini_file(block_config_file);
256 // std::cerr << block_config_file << std::endl;
257 po::variables_map vm;
258 po::options_description config_file_options;
260
261 std::ostringstream str_order;
262 str_order << "block_" << it->getMeshsetId() << ".temperature_order";
263 config_file_options.add_options()(
264 str_order.str().c_str(),
265 po::value<int>(&block_data[it->getMeshsetId()].oRder)
266 ->default_value(order));
267
268 std::ostringstream str_cond;
269 str_cond << "block_" << it->getMeshsetId() << ".heat_conductivity";
270 config_file_options.add_options()(
271 str_cond.str().c_str(),
272 po::value<double>(&block_data[it->getMeshsetId()].cOnductivity)
273 ->default_value(-1));
274
275 std::ostringstream str_capa;
276 str_capa << "block_" << it->getMeshsetId() << ".heat_capacity";
277 config_file_options.add_options()(
278 str_capa.str().c_str(),
279 po::value<double>(&block_data[it->getMeshsetId()].cApacity)
280 ->default_value(-1));
281
282 std::ostringstream str_init_temp;
283 str_init_temp << "block_" << it->getMeshsetId()
284 << ".initial_temperature";
285 config_file_options.add_options()(
286 str_init_temp.str().c_str(),
287 po::value<double>(&block_data[it->getMeshsetId()].initTemp)
288 ->default_value(0));
289 }
290 config_file_options.add_options()(
291 "climate_model.solar_radiation",
292 po::value<bool>(&solar_radiation)->default_value(false));
293
294 po::parsed_options parsed =
295 parse_config_file(ini_file, config_file_options, true);
296 store(parsed,vm);
297 po::notify(vm);
298
300 if (block_data[it->getMeshsetId()].oRder == -1)
301 continue;
302 if (block_data[it->getMeshsetId()].oRder == order)
303 continue;
304 PetscPrintf(PETSC_COMM_WORLD, "Set block %d oRder to %d\n",
305 it->getMeshsetId(), block_data[it->getMeshsetId()].oRder);
306 Range block_ents;
307 CHKERR moab.get_entities_by_handle(it->meshset, block_ents, true);
308 Range ents_to_set_order;
309 CHKERR moab.get_adjacencies(block_ents, 3, false, ents_to_set_order,
310 moab::Interface::UNION);
311 ents_to_set_order = ents_to_set_order.subset_by_type(MBTET);
312 CHKERR moab.get_adjacencies(block_ents, 2, false, ents_to_set_order,
313 moab::Interface::UNION);
314 CHKERR moab.get_adjacencies(block_ents, 1, false, ents_to_set_order,
315 moab::Interface::UNION);
316 CHKERR m_field.set_field_order(ents_to_set_order, "TEMP",
317 block_data[it->getMeshsetId()].oRder);
318 CHKERR m_field.set_field_order(ents_to_set_order, "TEMP_RATE",
319 block_data[it->getMeshsetId()].oRder);
320 }
321 std::vector<std::string> additional_parameters;
322 additional_parameters =
323 collect_unrecognized(parsed.options, po::include_positional);
324 for (std::vector<std::string>::iterator vit =
325 additional_parameters.begin();
326 vit != additional_parameters.end(); vit++) {
327 CHKERR PetscPrintf(PETSC_COMM_WORLD,
328 "** WARRING Unrecognised option %s\n", vit->c_str());
329 }
330
331 } catch (const std::exception& ex) {
332 std::ostringstream ss;
333 ss << ex.what() << std::endl;
334 SETERRQ(PETSC_COMM_SELF,MOFEM_STD_EXCEPTION_THROW,ss.str().c_str());
335 }
336 }
337
338 // this default class to calculate thermal elements
339 ThermalElement thermal_elements(m_field);
340 CHKERR thermal_elements.addThermalElements("TEMP");
341 CHKERR thermal_elements.addThermalFluxElement("TEMP");
342 CHKERR thermal_elements.addThermalConvectionElement("TEMP");
343 CHKERR thermal_elements.addThermalRadiationElement("TEMP");
344 // add rate of temperature to data field of finite element
345 CHKERR m_field.modify_finite_element_add_field_data("THERMAL_FE",
346 "TEMP_RATE");
347 // and temperature element default element operators at integration (gauss)
348 // points
349 CHKERR addHOOpsVol("MESH_NODE_POSITIONS", thermal_elements.getLoopFeRhs(), true,
350 false, false, false);
351 CHKERR addHOOpsVol("MESH_NODE_POSITIONS", thermal_elements.getLoopFeLhs(), true,
352 false, false, false);
353 CHKERR thermal_elements.setTimeSteppingProblem("TEMP", "TEMP_RATE");
354
355 // set block material data from option file
356 std::map<int, ThermalElement::BlockData>::iterator mit;
357 mit = thermal_elements.setOfBlocks.begin();
358 for (; mit != thermal_elements.setOfBlocks.end(); mit++) {
359 // std::cerr << mit->first << std::endl;
360 // std::cerr << block_data[mit->first].cOnductivity << " " <<
361 // block_data[mit->first].cApacity << std::endl;
362 if (block_data[mit->first].cOnductivity != -1) {
363 PetscPrintf(PETSC_COMM_WORLD, "Set block %d heat conductivity to %3.2e\n",
364 mit->first, block_data[mit->first].cOnductivity);
365 for (int dd = 0; dd < 3; dd++) {
366 mit->second.cOnductivity_mat(dd, dd) =
367 block_data[mit->first].cOnductivity;
368 }
369 }
370 if (block_data[mit->first].cApacity != -1) {
371 PetscPrintf(PETSC_COMM_WORLD, "Set block %d heat capacity to %3.2e\n",
372 mit->first, block_data[mit->first].cApacity);
373 mit->second.cApacity = block_data[mit->first].cApacity;
374 }
375 }
376
377#ifdef __GROUND_SURFACE_TEMPERATURE_HPP
378 GroundSurfaceTemperature ground_surface(m_field);
379 CrudeClimateModel time_data(time_data_file_for_ground_surface);
380 GroundSurfaceTemperature::PreProcess exectuteGenericClimateModel(&time_data);
381 if (ground_temperature_analysis) {
382 CHKERR ground_surface.addSurfaces("TEMP");
383 CHKERR ground_surface.setOperators(&time_data, "TEMP");
384 }
385#endif //__GROUND_SURFACE_TEMPERATURE_HPP
386
387 //build database, i.e. declare dofs, elements and adjacencies
388
389 // build field
390 CHKERR m_field.build_fields();
391 // project 10 node tet approximation of geometry on hierarchical basis
392 Projection10NodeCoordsOnField ent_method_material(m_field,
393 "MESH_NODE_POSITIONS");
394 CHKERR m_field.loop_dofs("MESH_NODE_POSITIONS", ent_method_material);
395
396 // set initial temperature from Cubit blocksets
397 mit = thermal_elements.setOfBlocks.begin();
398 for (; mit != thermal_elements.setOfBlocks.end(); mit++) {
399 if (mit->second.initTemp != 0) {
400 Range vertices;
401 CHKERR moab.get_connectivity(mit->second.tEts, vertices, true);
402 CHKERR m_field.getInterface<FieldBlas>()->setField(
403 mit->second.initTemp, MBVERTEX, vertices, "TEMP");
404 }
405 }
406
408 if (std::regex_match(it->getName(), std::regex("INT_THERMAL(.*)"))) {
409 std::vector<double> data;
410 CHKERR it->getAttributes(data);
411 if (data.size() != 1)
412 SETERRQ(PETSC_COMM_SELF, 1, "Data inconsistency");
413 Range block_ents, block_verts;
414 CHKERR moab.get_entities_by_handle(it->getMeshset(), block_ents, true);
415 CHKERR moab.get_connectivity(block_ents, block_verts, true);
416 CHKERR m_field.getInterface<FieldBlas>()->setField(data[0], MBVERTEX,
417 block_verts, "TEMP");
418 }
419 }
420
422 if (block_data[it->getMeshsetId()].initTemp != 0) {
423 Range block_ents;
424 CHKERR moab.get_entities_by_handle(it->meshset, block_ents, true);
425 Range vertices;
426 CHKERR moab.get_connectivity(block_ents, vertices, true);
427 CHKERR m_field.getInterface<FieldBlas>()->setField(
428 block_data[it->getMeshsetId()].initTemp, MBVERTEX, vertices, "TEMP");
429 }
430 }
431
432 MPI_Comm moab_comm_world;
433 MPI_Comm_dup(PETSC_COMM_WORLD, &moab_comm_world);
434 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
435 if (pcomm == NULL)
436 pcomm = new ParallelComm(&moab, moab_comm_world);
437
438 PetscBool is_partitioned = PETSC_FALSE;
439 CHKERR PetscOptionsGetBool(PETSC_NULL, PETSC_NULL, "-dm_is_partitioned",
440 &is_partitioned, PETSC_NULL);
441
442 Range thermal_element_ents;
443 CHKERR m_field.get_finite_element_entities_by_dimension("THERMAL_FE", 3,
444 thermal_element_ents);
445 Skinner skin(&m_field.get_moab());
446 Range skin_faces; // skin faces from 3d ents
447 CHKERR skin.find_skin(0, thermal_element_ents, false, skin_faces);
448 Range proc_skin;
449 if (is_partitioned) {
450 CHKERR pcomm->filter_pstatus(skin_faces,
451 PSTATUS_SHARED | PSTATUS_MULTISHARED,
452 PSTATUS_NOT, -1, &proc_skin);
453 } else {
454 proc_skin = skin_faces;
455 }
456 CHKERR m_field.add_finite_element("POST_PROC_SKIN");
457 CHKERR m_field.modify_finite_element_add_field_row("POST_PROC_SKIN", "TEMP");
458 CHKERR m_field.modify_finite_element_add_field_col("POST_PROC_SKIN", "TEMP");
459 CHKERR m_field.modify_finite_element_add_field_data("POST_PROC_SKIN", "TEMP");
460 CHKERR m_field.modify_finite_element_add_field_data("POST_PROC_SKIN",
461 "MESH_NODE_POSITIONS");
462 CHKERR m_field.add_ents_to_finite_element_by_dim(proc_skin, 2,
463 "POST_PROC_SKIN");
464
465 // build finite elemnts
467 // build adjacencies
468 CHKERR m_field.build_adjacencies(bit_level0);
469
470 // delete old temperature recorded series
471 SeriesRecorder *recorder_ptr;
472 CHKERR m_field.getInterface(recorder_ptr);
473 if (recorder_ptr->check_series("THEMP_SERIES")) {
474 /*for(_IT_SERIES_STEPS_BY_NAME_FOR_LOOP_(recorder_ptr,"THEMP_SERIES",sit)) {
475 CHKERR
476 recorder_ptr->load_series_data("THEMP_SERIES",sit->get_step_number());
477 }*/
478 CHKERR recorder_ptr->delete_recorder_series("THEMP_SERIES");
479 }
480
481 std::vector<std::array<double, 3>> eval_points;
482 eval_points.resize(0);
483 PetscBool eval_points_flg = PETSC_FALSE;
484 char eval_points_file[255];
485 CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_eval_points_file",
486 eval_points_file, 255, &eval_points_flg);
487 if (eval_points_flg) {
488 std::ifstream in_file(eval_points_file, std::ios::in);
489 if (!in_file.is_open()) {
490 SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Cannot open file %s",
491 eval_points_file);
492 }
493 double x, y, z;
494 while (in_file >> x >> y >> z) {
495 eval_points.push_back({x, y, z});
496 }
497 }
498
499 // set dm data structure which created mofem data structures
500 CHKERR DMMoFEMCreateMoFEM(dm, &m_field, dm_name, bit_level0);
501 CHKERR DMSetFromOptions(dm);
502 // add elements to dm
503 CHKERR DMMoFEMAddElement(dm, "THERMAL_FE");
504 CHKERR DMMoFEMAddElement(dm, "THERMAL_FLUX_FE");
505 CHKERR DMMoFEMAddElement(dm, "THERMAL_CONVECTION_FE");
506 CHKERR DMMoFEMAddElement(dm, "THERMAL_RADIATION_FE");
507 CHKERR DMMoFEMAddElement(dm, "POST_PROC_SKIN");
508#ifdef __GROUND_SURFACE_TEMPERATURE_HPP
509 if (ground_temperature_analysis) {
510 CHKERR DMMoFEMAddElement(dm, "GROUND_SURFACE_FE");
511 }
512#endif //__GROUND_SURFACE_TEMPERATURE_HPP
513
514 CHKERR DMSetUp(dm);
515
516 // create matrices
517 Vec T, F;
519 CHKERR VecDuplicate(T, &F);
520 Mat A;
522
523 DirichletTemperatureBc dirichlet_bc(m_field, "TEMP", A, T, F);
524 ThermalElement::UpdateAndControl update_velocities(m_field, "TEMP",
525 "TEMP_RATE");
526 ThermalElement::TimeSeriesMonitor monitor(m_field, "THEMP_SERIES", "TEMP",
527 eval_points);
528 MonitorPostProc post_proc(m_field);
529
530 // Initialize data with values save of on the field
531 CHKERR VecZeroEntries(T);
532 CHKERR DMoFEMMeshToLocalVector(dm, T, INSERT_VALUES, SCATTER_FORWARD);
533 CHKERR DMoFEMPreProcessFiniteElements(dm, &dirichlet_bc);
534 CHKERR DMoFEMMeshToGlobalVector(dm, T, INSERT_VALUES, SCATTER_REVERSE);
535
536 // preprocess
537 CHKERR DMMoFEMTSSetIFunction(dm, DM_NO_ELEMENT, NULL, &update_velocities,
538 NULL);
539 CHKERR DMMoFEMTSSetIFunction(dm, DM_NO_ELEMENT, NULL, &dirichlet_bc, NULL);
540 CHKERR DMMoFEMTSSetIJacobian(dm, DM_NO_ELEMENT, NULL, &dirichlet_bc, NULL);
541#ifdef __GROUND_SURFACE_TEMPERATURE_HPP
543 &exectuteGenericClimateModel, NULL);
544 { // add preprocessor, calculating angle on which sun ray on the surface
545 if (solar_radiation) {
546 boost::ptr_vector<
548 hi_it;
549 it = ground_surface.preProcessShade.begin();
550 hi_it = ground_surface.preProcessShade.end();
551 for (; it != hi_it; it++) {
552 CHKERR DMMoFEMTSSetIFunction(dm, DM_NO_ELEMENT, NULL, &*it, NULL);
553 }
554 }
555 }
556#endif //__GROUND_SURFACE_TEMPERATURE_HPP
557
558 // loops rhs
559 CHKERR DMMoFEMTSSetIFunction(dm, "THERMAL_FE", &thermal_elements.feRhs, NULL,
560 NULL);
561 CHKERR DMMoFEMTSSetIFunction(dm, "THERMAL_FLUX_FE", &thermal_elements.feFlux,
562 NULL, NULL);
563 CHKERR DMMoFEMTSSetIFunction(dm, "THERMAL_CONVECTION_FE",
564 &thermal_elements.feConvectionRhs, NULL, NULL);
565 CHKERR DMMoFEMTSSetIFunction(dm, "THERMAL_RADIATION_FE",
566 &thermal_elements.feRadiationRhs, NULL, NULL);
567#ifdef __GROUND_SURFACE_TEMPERATURE_HPP
568 if (ground_temperature_analysis) {
569 CHKERR DMMoFEMTSSetIFunction(dm, "GROUND_SURFACE_FE",
570 &ground_surface.getFeGroundSurfaceRhs(), NULL,
571 NULL);
572 }
573#endif //__GROUND_SURFACE_TEMPERATURE_HPP
574
575 // loops lhs
576 CHKERR DMMoFEMTSSetIJacobian(dm, "THERMAL_FE", &thermal_elements.feLhs, NULL,
577 NULL);
578 CHKERR DMMoFEMTSSetIJacobian(dm, "THERMAL_CONVECTION_FE",
579 &thermal_elements.feConvectionLhs, NULL, NULL);
580 CHKERR DMMoFEMTSSetIJacobian(dm, "THERMAL_RADIATION_FE",
581 &thermal_elements.feRadiationLhs, NULL, NULL);
582#ifdef __GROUND_SURFACE_TEMPERATURE_HPP
583 if (ground_temperature_analysis) {
584 CHKERR DMMoFEMTSSetIJacobian(dm, "GROUND_SURFACE_FE",
585 &ground_surface.getFeGroundSurfaceLhs(), NULL,
586 NULL);
587 }
588#endif //__GROUND_SURFACE_TEMPERATURE_HPP
589
590 //postprocess
591 CHKERR DMMoFEMTSSetIFunction(dm,DM_NO_ELEMENT,NULL,NULL,&dirichlet_bc);
592 CHKERR DMMoFEMTSSetIJacobian(dm,DM_NO_ELEMENT,NULL,NULL,&dirichlet_bc);
593
594 TsCtx *ts_ctx;
596 //add monitor operator
597 ts_ctx->getPostProcessMonitor().push_back(&monitor);
598 ts_ctx->getPostProcessMonitor().push_back(&post_proc);
599
600 //create time solver
601 TS ts;
602 CHKERR TSCreate(PETSC_COMM_WORLD,&ts);
603 CHKERR TSSetType(ts,TSBEULER);
604
605 CHKERR TSSetIFunction(ts,F,PETSC_NULL,PETSC_NULL);
606 CHKERR TSSetIJacobian(ts,A,A,PETSC_NULL,PETSC_NULL);
607 //add monitor to TS solver
608 CHKERR TSMonitorSet(ts,TsMonitorSet,ts_ctx,PETSC_NULL); // !!!
609
610 CHKERR recorder_ptr->add_series_recorder("THEMP_SERIES");
611 //start to record
612 CHKERR recorder_ptr->initialize_series_recorder("THEMP_SERIES");
613
614 double ftime = 1;
615 CHKERR TSSetDuration(ts,PETSC_DEFAULT,ftime);
616 CHKERR TSSetFromOptions(ts);
617 CHKERR TSSetDM(ts,dm);
618
619 CHKERR TSSolve(ts,T);
620 CHKERR TSGetTime(ts,&ftime);
621
622 //end recoder
623 CHKERR recorder_ptr->finalize_series_recorder("THEMP_SERIES");
624
625 PetscInt steps,snesfails,rejects,nonlinits,linits;
626 CHKERR TSGetTimeStepNumber(ts,&steps);
627 CHKERR TSGetSNESFailures(ts,&snesfails);
628 CHKERR TSGetStepRejections(ts,&rejects);
629 CHKERR TSGetSNESIterations(ts,&nonlinits);
630 CHKERR TSGetKSPIterations(ts,&linits);
631
632 PetscPrintf(PETSC_COMM_WORLD,
633 "steps %D (%D rejected, %D SNES fails), ftime %g, nonlinits %D, "
634 "linits %D\n",
635 steps, rejects, snesfails, ftime, nonlinits, linits);
636
637 // save solution, if boundary conditions are defined you can use that file in
638 // mechanical problem to calculate thermal stresses
639 PetscBool save_solution = PETSC_TRUE;
640 CHKERR PetscOptionsGetBool(PETSC_NULL, PETSC_NULL, "-my_save_solution",
641 &save_solution, PETSC_NULL);
642
643 if (save_solution) {
644 if (is_partitioned) {
645 CHKERR moab.write_file("solution.h5m");
646 } else {
647 if (m_field.get_comm_rank() == 0) {
648 CHKERR moab.write_file("solution.h5m");
649 }
650 }
651 }
652
653 CHKERR MatDestroy(&A);
654 CHKERR VecDestroy(&T);
655 CHKERR VecDestroy(&F);
656
657 CHKERR TSDestroy(&ts);
658
659 }
661
662 return MoFEM::Core::Finalize();
663}
#define DM_NO_ELEMENT
Definition: DMMoFEM.hpp:10
const std::string default_options
std::string param_file
int main()
Definition: adol-c_atom.cpp:46
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
@ MF_ZERO
Definition: definitions.h:98
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
@ H1
continuous field
Definition: definitions.h:85
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ BLOCKSET
Definition: definitions.h:148
@ MOFEM_NOT_FOUND
Definition: definitions.h:33
@ MOFEM_NOT_INSTALLED
Definition: definitions.h:37
@ MOFEM_STD_EXCEPTION_THROW
Definition: definitions.h:39
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
@ F
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition: DMMoFEM.cpp:483
PetscErrorCode DMMoFEMTSSetIFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set TS implicit function evaluation function
Definition: DMMoFEM.cpp:786
PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[], const MoFEM::BitRefLevel bit_level, const MoFEM::BitRefLevel bit_mask=MoFEM::BitRefLevel().set())
Must be called by user to set MoFEM data structures.
Definition: DMMoFEM.cpp:118
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMoFEM.cpp:509
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
Definition: DMMoFEM.cpp:1183
PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx)
get MoFEM::TsCtx data structure
Definition: DMMoFEM.cpp:1128
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS Jacobian evaluation function
Definition: DMMoFEM.cpp:839
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVector.
Definition: DMMoFEM.cpp:1153
PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode)
set ghosted vector values on all existing mesh entities
Definition: DMMoFEM.cpp:521
PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition: DMMoFEM.cpp:532
virtual MoFEMErrorCode add_ents_to_finite_element_by_dim(const EntityHandle entities, const int dim, const std::string &name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
virtual MoFEMErrorCode get_finite_element_entities_by_dimension(const std::string name, int dim, Range &ents) const =0
get entities in the finite element by dimension
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.
MoFEMErrorCode addFieldValuesGradientPostProc(const std::string field_name, Vec v=PETSC_NULL)
Add operator to post-process L2 or H1 field gradient.
MoFEMErrorCode addFieldValuesPostProc(const std::string field_name, Vec v=PETSC_NULL)
Add operator to post-process L2, H1, Hdiv, Hcurl field value.
MoFEMErrorCode writeFile(const std::string file_name, const char *file_type="MOAB", const char *file_options="PARALLEL=WRITE_PART")
wrote results in (MOAB) format, use "file_name.h5m"
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
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 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.
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
virtual MoFEMErrorCode finalize_series_recorder(const std::string &serie_name)
virtual bool check_series(const std::string &name) const
check if series is in database
virtual MoFEMErrorCode initialize_series_recorder(const std::string &serie_name)
virtual MoFEMErrorCode add_series_recorder(const std::string &series_name)
virtual MoFEMErrorCode delete_recorder_series(const std::string &series_name)
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 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
MoFEM::TsCtx * ts_ctx
Definition: level_set.cpp:1884
char mesh_file_name[255]
const double T
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Set monitor for TS solver.
Definition: TsCtx.cpp:221
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
MoFEMErrorCode addHOOpsFace3D(const std::string field, E &e, bool hcurl, bool hdiv)
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
PetscErrorCode monitor(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx)
Implementation of ground surface temperature.
MoFEMErrorCode setOperators(GenericClimateModel *time_data_ptr, string field_name, const string mesh_nodals_positions="MESH_NODE_POSITIONS")
boost::ptr_vector< SolarRadiationPreProcessor > preProcessShade
MoFEMErrorCode addSurfaces(const string field_name, const string mesh_nodals_positions="MESH_NODE_POSITIONS")
Managing BitRefLevels.
virtual moab::Interface & get_moab()=0
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.
structure for User Loop Methods on finite elements
Basic algebra on fields.
Definition: FieldBlas.hpp:21
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 > getStrmSync()
Get the strm sync object.
Definition: LogManager.cpp:348
Projection of edge entities with one mid-node on hierarchical basis.
TS ts
time solver
Interface for Time Stepping (TS) solver.
Definition: TsCtx.hpp:15
BasicMethodsSequence & getPostProcessMonitor()
Get the postProcess to do Monitor object.
Definition: TsCtx.hpp:154
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
MoFEMErrorCode postProcess()
function is run at the end of loop
MoFEMErrorCode operator()()
function is run for every finite element
MoFEMErrorCode preProcess()
function is run at the beginning of loop
MoFEM::Interface & mField
bool iNit
int pRT
PostProcFaceOnRefinedMesh skinPostProc
MonitorPostProc(MoFEM::Interface &m_field)
PostProcVolumeOnRefinedMesh postProc
int * step
Postprocess on face.
MoFEMErrorCode generateReferenceElementMesh()
MoFEMErrorCode generateReferenceElementMesh()
Generate reference mesh on single element.
Post processing.
TS monitore it records temperature at time steps.
this calass is to control time stepping
structure grouping operators and data used for thermal problems
MyTriFE feRadiationLhs
MyTriFE feConvectionLhs
MyVolumeFE feRhs
cauclate right hand side for tetrahedral elements
MyVolumeFE & getLoopFeLhs()
get lhs volume element
MyTriFE feConvectionRhs
MyVolumeFE & getLoopFeRhs()
get rhs volume element
MyVolumeFE feLhs
MyTriFE feRadiationRhs
std::map< int, BlockData > setOfBlocks
maps block set id with appropriate BlockData
static char help[]