v0.15.0
Loading...
Searching...
No Matches
thermo_elastic.cpp
Go to the documentation of this file.
1/**
2 * \file thermo_elastic.cpp
3 * \example thermo_elastic.cpp
4 *
5 * Thermo plasticity
6 *
7 */
8
9#ifndef EXECUTABLE_DIMENSION
10#define EXECUTABLE_DIMENSION 3
11#endif
12
13#include <MoFEM.hpp>
15
16using namespace MoFEM;
17
18#include <ThermalConvection.hpp>
19#include <ThermalRadiation.hpp>
20
21constexpr int SPACE_DIM =
22 EXECUTABLE_DIMENSION; //< Space dimension of problem, mesh
23
24constexpr bool IS_LARGE_STRAINS =
25 FINITE_DEFORMATION_FLAG; //< Flag to turn off/on geometric nonlinearities
26
29using DomainEleOp = DomainEle::UserDataOperator;
30using BoundaryEle =
32using BoundaryEleOp = BoundaryEle::UserDataOperator;
40constexpr AssemblyType AT = AssemblyType::PETSC; //< selected assembly type
42 IntegrationType::GAUSS; //< selected integration type
43
44//! [Elastic problem]
45#include <HenckyOps.hpp> // Include Hencky operators
46using namespace HenckyOps;
47//! [Elastic problem]
48
49// Include finite deformation operators
50#include <FiniteThermalOps.hpp> // Include operators for finite strain diffusion problem
51using namespace FiniteThermalOps;
52
53//! [Thermal problem]
54#include <ThermalOps.hpp> // Include operators for thermal problem which are agnostic to small/large strains
55using namespace ThermalOps;
56//! [Thermal problem]
57
58//! [Body and heat source]
61using OpBodyForce =
67//! [Body and heat source]
68
69//! [Natural boundary conditions]
73
74//! [Natural boundary conditions]
75
76//! [Essential boundary conditions (Least square approach)]
79 GAUSS>::OpEssentialRhs<HeatFluxCubitBcData, 3, SPACE_DIM>;
82 GAUSS>::OpEssentialLhs<HeatFluxCubitBcData, 3, SPACE_DIM>;
83//! [Essential boundary conditions (Least square approach)]
85//! [Default input parameters]
86double default_young_modulus = 1;
87double default_poisson_ratio = 0.25;
88double default_ref_temp = 0.0;
89double default_init_temp = 0.0;
93 1; // Force / (time temperature ) or Power /
94 // (length temperature) // Time unit is hour. force unit kN
95double default_heat_capacity = 1; // length^2/(time^2 temperature) // length is
96 // millimeter time is hour
97
98int order_temp = 2; //< default approximation order for the temperature field
99int order_flux = 3; //< default approximation order for heat flux field
100int order_disp = 3; //< default approximation order for the displacement field
101
102int atom_test = 0;
104int save_every = 1; //< Save every n-th step
105PetscBool do_output_domain;
107//! [Default input parameters]
109#include <ThermoElasticOps.hpp> //< additional coupling operators
110using namespace ThermoElasticOps; //< name space of coupling operators
111
117auto save_range = [](moab::Interface &moab, const std::string name,
118 const Range r) {
120 auto out_meshset = get_temp_meshset_ptr(moab);
121 CHKERR moab.add_entities(*out_meshset, r);
122 CHKERR moab.write_file(name.c_str(), "VTK", "", out_meshset->get_ptr(), 1);
124};
125
127
129
131
132private:
134
135 PetscBool doEvalField;
136 std::array<double, SPACE_DIM> fieldEvalCoords;
137 boost::shared_ptr<FieldEvaluatorInterface::SetPtsData> fieldEvalData;
138
139 boost::shared_ptr<VectorDouble> tempFieldPtr;
140 boost::shared_ptr<MatrixDouble> fluxFieldPtr;
141 boost::shared_ptr<MatrixDouble> dispFieldPtr;
142 boost::shared_ptr<MatrixDouble> dispGradPtr;
143 boost::shared_ptr<MatrixDouble> strainFieldPtr;
144 boost::shared_ptr<MatrixDouble> stressFieldPtr;
146 MoFEMErrorCode setupProblem(); ///< add fields
148 getCommandLineParameters(); //< read parameters from command line
149 MoFEMErrorCode bC(); //< read boundary conditions
150 MoFEMErrorCode OPs(); //< add operators to pipeline
151 MoFEMErrorCode tsSolve(); //< time solver
152
154 : public boost::enable_shared_from_this<BlockedThermalParameters> {
155 double heatConductivity;
156 double heatCapacity;
159 return boost::shared_ptr<double>(shared_from_this(), &heatConductivity);
161
162 inline auto getHeatCapacityPtr() {
163 return boost::shared_ptr<double>(shared_from_this(), &heatCapacity);
164 }
165 };
168 : public boost::enable_shared_from_this<BlockedThermoElasticParameters> {
170 double refTemp;
172 inline auto getCoeffExpansionPtr() {
173 return boost::shared_ptr<VectorDouble>(shared_from_this(),
175 }
177 inline auto getRefTempPtr() {
178 return boost::shared_ptr<double>(shared_from_this(), &refTemp);
179 }
180 };
183 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
184 std::string block_name,
185 boost::shared_ptr<BlockedThermalParameters> blockedParamsPtr, Sev sev);
186
188 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
189 std::string block_name,
190 boost::shared_ptr<BlockedThermoElasticParameters> blockedParamsPtr,
191 Sev sev);
192
193 template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
195 MoFEM::Interface &m_field,
196 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
197 std::string field_name,
198 boost::shared_ptr<HenckyOps::CommonData> elastic_common_ptr,
199 boost::shared_ptr<ThermoElasticProblem::BlockedThermalParameters>
200 thermal_common_ptr,
201 boost::shared_ptr<ThermoElasticProblem::BlockedThermoElasticParameters>
202 thermoelastic_common_ptr,
203 Sev sev) {
205
206 using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
207 A>::template LinearForm<I>;
209 auto vec_temp_ptr = boost::make_shared<VectorDouble>();
210 pip.push_back(new OpCalculateScalarFieldValues("T", vec_temp_ptr));
211 auto coeff_expansion_ptr = thermoelastic_common_ptr->getCoeffExpansionPtr();
212 auto ref_temp_ptr = thermoelastic_common_ptr->getRefTempPtr();
213 pip.push_back(
214 new typename H::template OpCalculateHenckyThermalStress<DIM, I, 0>(
215 "U", vec_temp_ptr, elastic_common_ptr, coeff_expansion_ptr,
216 ref_temp_ptr));
218 typename B::template OpGradTimesTensor<1, DIM, DIM>;
219 pip.push_back(new OpInternalForcePiola(
220 "U", elastic_common_ptr->getMatFirstPiolaStress()));
221
223 }
224
225 template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
227 MoFEM::Interface &m_field,
228 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
229 std::string field_name, std::string elastic_block_name,
230 std::string thermal_block_name, std::string thermoelastic_block_name,
231 Sev sev, double scale = 1) {
233
234 auto elastic_common_ptr = commonDataFactory<DIM, I, DomainEleOp>(
235 m_field, pip, field_name, elastic_block_name, sev, scale);
236 auto thermal_common_ptr = boost::make_shared<BlockedThermalParameters>();
237 CHKERR addMatThermalBlockOps(pip, thermal_block_name, thermal_common_ptr,
238 Sev::inform);
239 auto thermoelastic_common_ptr =
240 boost::make_shared<BlockedThermoElasticParameters>();
241 CHKERR addMatThermoElasticBlockOps(pip, thermoelastic_block_name,
242 thermoelastic_common_ptr, Sev::inform);
243 CHKERR opThermoElasticFactoryDomainRhs<DIM, A, I, DomainEleOp>(
244 m_field, pip, field_name, elastic_common_ptr, thermal_common_ptr,
245 thermoelastic_common_ptr, sev);
246
248 }
249
250 template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
252 MoFEM::Interface &m_field,
253 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
254 std::string field_name, std::string coupled_field_name,
255 boost::shared_ptr<HenckyOps::CommonData> elastic_common_ptr,
256 boost::shared_ptr<ThermoElasticProblem::BlockedThermalParameters>
257 thermal_common_ptr,
258 boost::shared_ptr<ThermoElasticProblem::BlockedThermoElasticParameters>
259 thermoelastic_common_ptr,
260 Sev sev) {
262
263 using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
264 A>::template BiLinearForm<I>;
265 using OpKPiola = typename B::template OpGradTensorGrad<1, DIM, DIM, 1>;
266
268 auto vec_temp_ptr = boost::make_shared<VectorDouble>();
269 pip.push_back(new OpCalculateScalarFieldValues("T", vec_temp_ptr));
270 auto coeff_expansion_ptr = thermoelastic_common_ptr->getCoeffExpansionPtr();
271 auto ref_temp_ptr = thermoelastic_common_ptr->getRefTempPtr();
272 pip.push_back(
273 new typename H::template OpCalculateHenckyThermalStress<DIM, I, 0>(
274 "U", vec_temp_ptr, elastic_common_ptr, coeff_expansion_ptr,
275 ref_temp_ptr));
276 pip.push_back(new typename H::template OpHenckyTangent<DIM, I, 0>(
277 field_name, elastic_common_ptr));
278 pip.push_back(new OpKPiola(field_name, field_name,
279 elastic_common_ptr->getMatTangent()));
280 pip.push_back(new typename H::template OpCalculateHenckyThermalStressdT<
281 DIM, I, AssemblyDomainEleOp, 0>(
282 field_name, coupled_field_name, elastic_common_ptr,
283 coeff_expansion_ptr));
284
286 }
287
288 template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
290 MoFEM::Interface &m_field,
291 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
292 std::string field_name, std::string coupled_field_name,
293 std::string elastic_block_name, std::string thermal_block_name,
294 std::string thermoelastic_block_name, Sev sev, double scale = 1) {
296
297 auto elastic_common_ptr = commonDataFactory<DIM, I, DomainEleOp>(
298 m_field, pip, field_name, elastic_block_name, sev, scale);
299 auto thermal_common_ptr = boost::make_shared<BlockedThermalParameters>();
300 CHKERR addMatThermalBlockOps(pip, thermal_block_name, thermal_common_ptr,
301 Sev::inform);
302 auto thermoelastic_common_ptr =
303 boost::make_shared<BlockedThermoElasticParameters>();
304 CHKERR addMatThermoElasticBlockOps(pip, thermoelastic_block_name,
305 thermoelastic_common_ptr, Sev::inform);
306 CHKERR opThermoElasticFactoryDomainLhs<DIM, A, I, DomainEleOp>(
307 m_field, pip, field_name, coupled_field_name, elastic_common_ptr,
308 thermal_common_ptr, thermoelastic_common_ptr, sev);
309
311 }
312};
313
315 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
316 std::string block_name,
317 boost::shared_ptr<BlockedThermalParameters> blockedParamsPtr, Sev sev) {
319
320 struct OpMatThermalBlocks : public DomainEleOp {
321 OpMatThermalBlocks(boost::shared_ptr<double> conductivity_ptr,
322 boost::shared_ptr<double> capacity_ptr,
323 MoFEM::Interface &m_field, Sev sev,
324 std::vector<const CubitMeshSets *> meshset_vec_ptr)
325 : DomainEleOp(NOSPACE, DomainEleOp::OPSPACE),
326 conductivityPtr(conductivity_ptr), capacityPtr(capacity_ptr) {
327 CHK_THROW_MESSAGE(extractThermalBlockData(m_field, meshset_vec_ptr, sev),
328 "Cannot get data from thermal block");
329 }
330
331 MoFEMErrorCode doWork(int side, EntityType type,
334
335 for (auto &b : blockData) {
336
337 if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
338 *conductivityPtr = b.conductivity;
339 *capacityPtr = b.capacity;
341 }
342 }
343
344 *conductivityPtr = default_heat_conductivity;
345 *capacityPtr = default_heat_capacity;
346
348 }
349
350 private:
351 struct BlockData {
352 double conductivity;
353 double capacity;
354 Range blockEnts;
355 };
356
357 std::vector<BlockData> blockData;
358
360 extractThermalBlockData(MoFEM::Interface &m_field,
361 std::vector<const CubitMeshSets *> meshset_vec_ptr,
362 Sev sev) {
364
365 for (auto m : meshset_vec_ptr) {
366 MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Thermal Block") << *m;
367 std::vector<double> block_data;
368 CHKERR m->getAttributes(block_data);
369 if (block_data.size() < 2) {
370 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
371 "Expected that block has at least two attributes");
372 }
373 auto get_block_ents = [&]() {
374 Range ents;
375 CHKERR
376 m_field.get_moab().get_entities_by_handle(m->meshset, ents, true);
377 return ents;
378 };
379
380 MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Thermal Block")
381 << m->getName() << ": conductivity = " << block_data[0]
382 << " capacity = " << block_data[1];
383
384 blockData.push_back({block_data[0], block_data[1], get_block_ents()});
385 }
386 MOFEM_LOG_CHANNEL("WORLD");
388 }
389
390 boost::shared_ptr<double> conductivityPtr;
391 boost::shared_ptr<double> capacityPtr;
392 };
393
394 pipeline.push_back(new OpMatThermalBlocks(
395 blockedParamsPtr->getHeatConductivityPtr(),
396 blockedParamsPtr->getHeatCapacityPtr(), mField, sev,
397
398 // Get blockset using regular expression
400
401 (boost::format("%s(.*)") % block_name).str()
402
403 ))
404
405 ));
406
408}
409
411 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
412 std::string block_name,
413 boost::shared_ptr<BlockedThermoElasticParameters> blockedParamsPtr,
414 Sev sev) {
416
417 struct OpMatThermoElasticBlocks : public DomainEleOp {
418 OpMatThermoElasticBlocks(boost::shared_ptr<VectorDouble> expansion_ptr,
419 boost::shared_ptr<double> ref_temp_ptr,
420 MoFEM::Interface &m_field, Sev sev,
421 std::vector<const CubitMeshSets *> meshset_vec_ptr)
422 : DomainEleOp(NOSPACE, DomainEleOp::OPSPACE),
423 expansionPtr(expansion_ptr), refTempPtr(ref_temp_ptr) {
425 extractThermoElasticBlockData(m_field, meshset_vec_ptr, sev),
426 "Cannot get data from thermoelastic block");
427 }
428
429 MoFEMErrorCode doWork(int side, EntityType type,
432
433 for (auto &b : blockData) {
434
435 if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
436 *expansionPtr = b.expansion;
437 *refTempPtr = b.ref_temp;
439 }
440 }
441
443 *refTempPtr = default_ref_temp;
444
446 }
447
448 private:
449 struct BlockData {
450 double ref_temp;
451 VectorDouble expansion;
452 Range blockEnts;
453 };
454
455 std::vector<BlockData> blockData;
456
457 MoFEMErrorCode extractThermoElasticBlockData(
458 MoFEM::Interface &m_field,
459 std::vector<const CubitMeshSets *> meshset_vec_ptr, Sev sev) {
461
462 for (auto m : meshset_vec_ptr) {
463 MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Thermoelastic Block") << *m;
464 std::vector<double> block_data;
465 CHKERR m->getAttributes(block_data);
466 if (block_data.size() < 2) {
467 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
468 "Expected that block has at least two attributes");
469 }
470 auto get_block_ents = [&]() {
471 Range ents;
472 CHKERR
473 m_field.get_moab().get_entities_by_handle(m->meshset, ents, true);
474 return ents;
475 };
476
477 auto get_expansion = [&]() {
478 VectorDouble expansion(SPACE_DIM, block_data[1]);
479 if (block_data.size() > 2) {
480 expansion[1] = block_data[2];
481 }
482 if (SPACE_DIM == 3 && block_data.size() > 3) {
483 expansion[2] = block_data[3];
484 }
485 return expansion;
486 };
487
488 auto coeff_exp_vec = get_expansion();
489
490 MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Thermoelastic Block")
491 << " ref_temp = " << block_data[0]
492 << " expansion = " << coeff_exp_vec;
493
494 blockData.push_back({block_data[0], coeff_exp_vec, get_block_ents()});
495 }
496 MOFEM_LOG_CHANNEL("WORLD");
498 }
499
500 boost::shared_ptr<VectorDouble> expansionPtr;
501 boost::shared_ptr<double> refTempPtr;
502 };
503
504 pipeline.push_back(new OpMatThermoElasticBlocks(
505 blockedParamsPtr->getCoeffExpansionPtr(),
506 blockedParamsPtr->getRefTempPtr(), mField, sev,
507
508 // Get blockset using regular expression
510
511 (boost::format("%s(.*)") % block_name).str()
512
513 ))
514
515 ));
516
518}
519
520//! [Run problem]
526 CHKERR OPs();
527 CHKERR tsSolve();
529}
530//! [Run problem]
531
532//! [Get command line parameters]
535
536 auto get_command_line_parameters = [&]() {
538
539 CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order_temp,
540 PETSC_NULL);
541 CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order_temp", &order_temp,
542 PETSC_NULL);
544 CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order_flux", &order_flux,
545 PETSC_NULL);
547 CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order_disp", &order_disp,
548 PETSC_NULL);
549 CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-atom_test", &atom_test,
550 PETSC_NULL);
551 CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-save_every", &save_every,
552 PETSC_NULL);
553
554 CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-young_modulus",
555 &default_young_modulus, PETSC_NULL);
556 CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-poisson_ratio",
557 &default_poisson_ratio, PETSC_NULL);
558 CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-coeff_expansion",
559 &default_coeff_expansion, PETSC_NULL);
560 CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-ref_temp", &default_ref_temp,
561 PETSC_NULL);
562 CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-init_temp",
563 &default_init_temp, PETSC_NULL);
564
565 CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-capacity",
566 &default_heat_capacity, PETSC_NULL);
567 CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-conductivity",
568 &default_heat_conductivity, PETSC_NULL);
569
570 if constexpr (SPACE_DIM == 2) {
571 do_output_domain = PETSC_TRUE;
572 do_output_skin = PETSC_FALSE;
573 } else {
574 do_output_domain = PETSC_FALSE;
575 do_output_skin = PETSC_TRUE;
576 }
577
578 CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-output_domain",
579 &do_output_domain, PETSC_NULL);
580 CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-output_skin", &do_output_skin,
581 PETSC_NULL);
582
583 MOFEM_LOG("ThermoElastic", Sev::inform)
584 << "Young's modulus " << default_young_modulus;
585 MOFEM_LOG("ThermoElastic", Sev::inform)
586 << "Poisson's ratio " << default_poisson_ratio;
587 MOFEM_LOG("ThermoElastic", Sev::inform)
588 << "Coeff of expansion " << default_coeff_expansion;
589 MOFEM_LOG("ThermoElastic", Sev::inform)
590 << "Default heat capacity " << default_heat_capacity;
591 MOFEM_LOG("ThermoElastic", Sev::inform)
592 << "Heat conductivity " << default_heat_conductivity;
593
594 MOFEM_LOG("ThermoElastic", Sev::inform)
595 << "Reference temperature " << default_ref_temp;
596 MOFEM_LOG("ThermoElastic", Sev::inform)
597 << "Initial temperature " << default_init_temp;
598
600 };
601
602 CHKERR get_command_line_parameters();
603
605}
606//! [Get command line parameters]
607
608//! [Set up problem]
612 // Add field
614 // Mechanical fields
616 CHKERR simple->addBoundaryField("U", H1, AINSWORTH_LEGENDRE_BASE, SPACE_DIM);
617 // Temperature
618 constexpr auto flux_space = (SPACE_DIM == 2) ? HCURL : HDIV;
619 CHKERR simple->addDomainField("T", L2, base, 1);
620 CHKERR simple->addDomainField("FLUX", flux_space, base, 1);
621 CHKERR simple->addBoundaryField("FLUX", flux_space, base, 1);
622 CHKERR simple->addBoundaryField("TBC", L2, base, 1);
623
624 CHKERR simple->setFieldOrder("U", order_disp);
625 CHKERR simple->setFieldOrder("FLUX", order_flux);
626 CHKERR simple->setFieldOrder("T", order_temp);
627 CHKERR simple->setFieldOrder("TBC", order_temp);
628
629 CHKERR simple->setUp();
630
631 int coords_dim = SPACE_DIM;
632 CHKERR PetscOptionsGetRealArray(NULL, NULL, "-field_eval_coords",
633 fieldEvalCoords.data(), &coords_dim,
634 &doEvalField);
635
636 tempFieldPtr = boost::make_shared<VectorDouble>();
637 fluxFieldPtr = boost::make_shared<MatrixDouble>();
638 dispFieldPtr = boost::make_shared<MatrixDouble>();
639 dispGradPtr = boost::make_shared<MatrixDouble>();
640 strainFieldPtr = boost::make_shared<MatrixDouble>();
641 stressFieldPtr = boost::make_shared<MatrixDouble>();
642
643 if (doEvalField) {
645 mField.getInterface<FieldEvaluatorInterface>()->getData<DomainEle>();
646
647 CHKERR mField.getInterface<FieldEvaluatorInterface>()->buildTree<SPACE_DIM>(
648 fieldEvalData, simple->getDomainFEName());
649
650 fieldEvalData->setEvalPoints(fieldEvalCoords.data(), 1);
651 auto no_rule = [](int, int, int) { return -1; };
652
653 auto field_eval_fe_ptr = fieldEvalData->feMethodPtr.lock();
654 field_eval_fe_ptr->getRuleHook = no_rule;
655
656 auto block_thermal_params = boost::make_shared<BlockedThermalParameters>();
657 auto block_thermoelastic_params =
658 boost::make_shared<BlockedThermoElasticParameters>();
659 auto coeff_expansion_ptr =
660 block_thermoelastic_params->getCoeffExpansionPtr();
661 auto ref_temp_ptr = block_thermoelastic_params->getRefTempPtr();
662
663 CHKERR addMatThermalBlockOps(field_eval_fe_ptr->getOpPtrVector(),
664 "MAT_THERMAL", block_thermal_params,
665 Sev::verbose);
667 field_eval_fe_ptr->getOpPtrVector(), "MAT_THERMOELASTIC",
668 block_thermoelastic_params, Sev::verbose);
669
671 field_eval_fe_ptr->getOpPtrVector(), {H1, HDIV});
672
673 auto hencky_common_data_ptr =
674 HenckyOps::commonDataFactory<SPACE_DIM, IT, DomainEleOp>(
675 mField, field_eval_fe_ptr->getOpPtrVector(), "U", "MAT_ELASTIC",
676 Sev::inform);
677 auto mat_D_ptr = hencky_common_data_ptr->matDPtr;
678 auto dispGradPtr = hencky_common_data_ptr->matGradPtr;
679 auto mat_stress_ptr = boost::make_shared<MatrixDouble>();
680
681 field_eval_fe_ptr->getOpPtrVector().push_back(
683 field_eval_fe_ptr->getOpPtrVector().push_back(
685 field_eval_fe_ptr->getOpPtrVector().push_back(
687 field_eval_fe_ptr->getOpPtrVector().push_back(
689 dispGradPtr));
690
692
693 field_eval_fe_ptr->getOpPtrVector().push_back(
694 new
695 typename H::template OpCalculateHenckyThermalStress<SPACE_DIM, IT, 0>(
696 "U", tempFieldPtr, hencky_common_data_ptr, coeff_expansion_ptr,
697 ref_temp_ptr));
698 if (!IS_LARGE_STRAINS) {
699 field_eval_fe_ptr->getOpPtrVector().push_back(
701 hencky_common_data_ptr->getMatFirstPiolaStress(),
703 field_eval_fe_ptr->getOpPtrVector().push_back(
705 } else {
706 field_eval_fe_ptr->getOpPtrVector().push_back(
707 new typename H::template OpCalculateLogC<SPACE_DIM, IT>(
708 "U", hencky_common_data_ptr));
709 stressFieldPtr = hencky_common_data_ptr->getMatFirstPiolaStress();
710 strainFieldPtr = hencky_common_data_ptr->getMatLogC();
711 };
712 }
713
715}
716//! [Set up problem]
717
718//! [Boundary condition]
721
722 MOFEM_LOG("SYNC", Sev::noisy) << "bC";
724
726 auto bc_mng = mField.getInterface<BcManager>();
727
728 auto get_skin = [&]() {
729 Range body_ents;
730 CHKERR mField.get_moab().get_entities_by_dimension(0, SPACE_DIM, body_ents);
731 Skinner skin(&mField.get_moab());
732 Range skin_ents;
733 CHKERR skin.find_skin(0, body_ents, false, skin_ents);
734 return skin_ents;
735 };
736
737 auto filter_flux_blocks = [&](auto skin, bool temp_bc = false) {
738 auto remove_cubit_blocks = [&](auto c) {
740 for (auto m :
741
742 mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(c)
743
744 ) {
745 Range ents;
746 CHKERR mField.get_moab().get_entities_by_dimension(
747 m->getMeshset(), SPACE_DIM - 1, ents, true);
748 skin = subtract(skin, ents);
749 }
751 };
752
753 auto remove_named_blocks = [&](auto n) {
755 for (auto m : mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(
756 std::regex(
757
758 (boost::format("%s(.*)") % n).str()
759
760 ))
761
762 ) {
763 Range ents;
764 CHKERR mField.get_moab().get_entities_by_dimension(
765 m->getMeshset(), SPACE_DIM - 1, ents, true);
766 skin = subtract(skin, ents);
767 }
769 };
770 if (!temp_bc) {
771 CHK_THROW_MESSAGE(remove_cubit_blocks(NODESET | TEMPERATURESET),
772 "remove_cubit_blocks");
773 CHK_THROW_MESSAGE(remove_named_blocks("TEMPERATURE"),
774 "remove_named_blocks");
775 }
776 CHK_THROW_MESSAGE(remove_cubit_blocks(SIDESET | HEATFLUXSET),
777 "remove_cubit_blocks");
778 CHK_THROW_MESSAGE(remove_named_blocks("HEATFLUX"), "remove_named_blocks");
779 CHK_THROW_MESSAGE(remove_named_blocks("CONVECTION"), "remove_named_blocks");
780 CHK_THROW_MESSAGE(remove_named_blocks("RADIATION"), "remove_named_blocks");
781 return skin;
782 };
783
784 auto filter_true_skin = [&](auto skin) {
785 Range boundary_ents;
786 ParallelComm *pcomm =
787 ParallelComm::get_pcomm(&mField.get_moab(), MYPCOMM_INDEX);
788 CHKERR pcomm->filter_pstatus(skin, PSTATUS_SHARED | PSTATUS_MULTISHARED,
789 PSTATUS_NOT, -1, &boundary_ents);
790 return boundary_ents;
791 };
792
793 auto remove_flux_ents = filter_true_skin(filter_flux_blocks(get_skin()));
794 auto remove_temp_bc_ents =
795 filter_true_skin(filter_flux_blocks(get_skin(), true));
796
797 CHKERR mField.getInterface<CommInterface>()->synchroniseEntities(
798 remove_flux_ents);
799 CHKERR mField.getInterface<CommInterface>()->synchroniseEntities(
800 remove_temp_bc_ents);
801
802 MOFEM_LOG("SYNC", Sev::noisy) << remove_flux_ents << endl;
804
805 MOFEM_LOG("SYNC", Sev::noisy) << remove_temp_bc_ents << endl;
807
808#ifndef NDEBUG
809
812 (boost::format("flux_remove_%d.vtk") % mField.get_comm_rank()).str(),
813 remove_flux_ents);
814
817 (boost::format("temp_bc_remove_%d.vtk") % mField.get_comm_rank()).str(),
818 remove_temp_bc_ents);
819
820#endif
821
822 CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
823 simple->getProblemName(), "FLUX", remove_flux_ents);
824 CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
825 simple->getProblemName(), "TBC", remove_temp_bc_ents);
826
827 auto set_init_temp = [](boost::shared_ptr<FieldEntity> field_entity_ptr) {
828 field_entity_ptr->getEntFieldData()[0] = default_init_temp;
829 return 0;
830 };
831 CHKERR mField.getInterface<FieldBlas>()->fieldLambdaOnEntities(set_init_temp,
832 "T");
833
834 CHKERR bc_mng->removeBlockDOFsOnEntities<DisplacementCubitBcData>(
835 simple->getProblemName(), "U");
836 CHKERR bc_mng->pushMarkDOFsOnEntities<HeatFluxCubitBcData>(
837 simple->getProblemName(), "FLUX", false);
838
840}
841//! [Boundary condition]
842
843//! [Push operators to pipeline]
846
847 MOFEM_LOG("SYNC", Sev::noisy) << "OPs";
849
850 auto pipeline_mng = mField.getInterface<PipelineManager>();
852 auto bc_mng = mField.getInterface<BcManager>();
853
854 auto boundary_marker =
855 bc_mng->getMergedBlocksMarker(vector<string>{"HEATFLUX"});
856
857 auto integration_rule = [](int, int, int approx_order) {
858 return 2 * approx_order;
859 };
860 CHKERR pipeline_mng->setDomainRhsIntegrationRule(integration_rule);
861 CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
862 CHKERR pipeline_mng->setBoundaryRhsIntegrationRule(integration_rule);
863 CHKERR pipeline_mng->setBoundaryLhsIntegrationRule(integration_rule);
864
865 auto block_thermal_params = boost::make_shared<BlockedThermalParameters>();
866 auto heat_conductivity_ptr = block_thermal_params->getHeatConductivityPtr();
867 auto heat_capacity_ptr = block_thermal_params->getHeatCapacityPtr();
868
869 auto block_thermoelastic_params =
870 boost::make_shared<BlockedThermoElasticParameters>();
871 auto coeff_expansion_ptr = block_thermoelastic_params->getCoeffExpansionPtr();
872 auto ref_temp_ptr = block_thermoelastic_params->getRefTempPtr();
873
874 // Default time scaling of BCs and sources, from command line arguments
875 auto time_scale =
876 boost::make_shared<TimeScale>("", false, [](const double) { return 1; });
877 auto def_time_scale = [time_scale](const double time) {
878 return (!time_scale->argFileScale) ? time : 1;
879 };
880 auto def_file_name = [time_scale](const std::string &&name) {
881 return (!time_scale->argFileScale) ? name : "";
882 };
883
884 // Files which define scaling for separate variables, if provided
885 auto time_bodyforce_scaling = boost::make_shared<TimeScale>(
886 def_file_name("bodyforce_scale.txt"), false, def_time_scale);
887 auto time_heatsource_scaling = boost::make_shared<TimeScale>(
888 def_file_name("heatsource_scale.txt"), false, def_time_scale);
889 auto time_temperature_scaling = boost::make_shared<TimeScale>(
890 def_file_name("temperature_bc_scale.txt"), false, def_time_scale);
891 auto time_displacement_scaling = boost::make_shared<TimeScale>(
892 def_file_name("displacement_bc_scale.txt"), false, def_time_scale);
893 auto time_flux_scaling = boost::make_shared<TimeScale>(
894 def_file_name("flux_bc_scale.txt"), false, def_time_scale);
895 auto time_force_scaling = boost::make_shared<TimeScale>(
896 def_file_name("force_bc_scale.txt"), false, def_time_scale);
897
898 auto add_domain_rhs_ops = [&](auto &pipeline) {
900 CHKERR addMatThermalBlockOps(pipeline, "MAT_THERMAL", block_thermal_params,
901 Sev::inform);
902 CHKERR addMatThermoElasticBlockOps(pipeline, "MAT_THERMOELASTIC",
903 block_thermoelastic_params, Sev::inform);
905
906 auto hencky_common_data_ptr =
907 HenckyOps::commonDataFactory<SPACE_DIM, IT, DomainEleOp>(
908 mField, pipeline, "U", "MAT_ELASTIC", Sev::inform);
909 auto mat_D_ptr = hencky_common_data_ptr->matDPtr;
910 auto mat_grad_ptr = hencky_common_data_ptr->matGradPtr;
911 auto mat_strain_ptr = boost::make_shared<MatrixDouble>();
912 auto mat_stress_ptr = boost::make_shared<MatrixDouble>();
913
914 auto vec_temp_ptr = boost::make_shared<VectorDouble>();
915 auto vec_temp_dot_ptr = boost::make_shared<VectorDouble>();
916 auto mat_flux_ptr = boost::make_shared<MatrixDouble>();
917 auto vec_temp_div_ptr = boost::make_shared<VectorDouble>();
918
919 pipeline.push_back(new OpCalculateScalarFieldValues("T", vec_temp_ptr));
920 pipeline.push_back(
921 new OpCalculateScalarFieldValuesDot("T", vec_temp_dot_ptr));
923 "FLUX", vec_temp_div_ptr));
924 pipeline.push_back(
925 new OpCalculateHVecVectorField<3, SPACE_DIM>("FLUX", mat_flux_ptr));
926
927 CHKERR opThermoElasticFactoryDomainRhs<SPACE_DIM, AT, IT, DomainEleOp>(
928 mField, pipeline, "U", "MAT_ELASTIC", "MAT_THERMAL",
929 "MAT_THERMOELASTIC", Sev::inform);
930
931 auto resistance = [heat_conductivity_ptr](const double, const double,
932 const double) {
933 return (1. / (*heat_conductivity_ptr));
934 };
935 // negative value is consequence of symmetric system
936 auto capacity = [heat_capacity_ptr](const double, const double,
937 const double) {
938 return -(*heat_capacity_ptr);
939 };
940 auto unity = [](const double, const double, const double) constexpr {
941 return -1.;
942 };
943 pipeline.push_back(
944 new OpHdivFlux("FLUX", mat_flux_ptr, resistance, mat_grad_ptr));
945 pipeline.push_back(new OpHDivTemp("FLUX", vec_temp_ptr, unity));
946 pipeline.push_back(new OpBaseDivFlux("T", vec_temp_div_ptr, unity));
947 pipeline.push_back(new OpBaseDotT("T", vec_temp_dot_ptr, capacity));
948
949 // Set source terms
951 pipeline, mField, "T", {time_scale, time_heatsource_scaling},
952 "HEAT_SOURCE", Sev::inform);
954 pipeline, mField, "U", {time_scale, time_bodyforce_scaling},
955 "BODY_FORCE", Sev::inform);
957 pipeline, mField, "T", vec_temp_ptr, "SETTEMP", Sev::inform);
958
960 };
961
962 auto add_domain_lhs_ops = [&](auto &pipeline) {
964 CHKERR addMatThermalBlockOps(pipeline, "MAT_THERMAL", block_thermal_params,
965 Sev::verbose);
966 CHKERR addMatThermoElasticBlockOps(pipeline, "MAT_THERMOELASTIC",
967 block_thermoelastic_params,
968 Sev::verbose);
970
971 auto hencky_common_data_ptr =
972 HenckyOps::commonDataFactory<SPACE_DIM, IT, DomainEleOp>(
973 mField, pipeline, "U", "MAT_ELASTIC", Sev::inform, 1);
974 auto mat_D_ptr = hencky_common_data_ptr->matDPtr;
975 auto mat_grad_ptr = hencky_common_data_ptr->matGradPtr;
976
977 auto resistance = [heat_conductivity_ptr](const double, const double,
978 const double) {
979 return (1. / (*heat_conductivity_ptr));
980 };
981 auto capacity = [heat_capacity_ptr](const double, const double,
982 const double) {
983 return -(*heat_capacity_ptr);
984 };
985 pipeline.push_back(
986 new OpHdivHdiv("FLUX", "FLUX", resistance, mat_grad_ptr));
987 pipeline.push_back(
988 new OpHdivT("FLUX", "T", []() constexpr { return -1; }, true));
989
990 auto mat_flux_ptr = boost::make_shared<MatrixDouble>();
991 pipeline.push_back(
992 new OpCalculateHVecVectorField<3, SPACE_DIM>("FLUX", mat_flux_ptr));
993
994 pipeline.push_back(
995 new OpHdivU("FLUX", "U", mat_flux_ptr, resistance, mat_grad_ptr));
996
997 CHKERR opThermoElasticFactoryDomainLhs<SPACE_DIM, AT, IT, DomainEleOp>(
998 mField, pipeline, "U", "T", "MAT_ELASTIC", "MAT_THERMAL",
999 "MAT_THERMOELASTIC", Sev::inform);
1000
1001 auto op_capacity = new OpCapacity("T", "T", capacity);
1002 op_capacity->feScalingFun = [](const FEMethod *fe_ptr) {
1003 return fe_ptr->ts_a;
1004 };
1005 pipeline.push_back(op_capacity);
1006
1007 auto vec_temp_ptr = boost::make_shared<VectorDouble>();
1008 pipeline.push_back(new OpCalculateScalarFieldValues("T", vec_temp_ptr));
1010 pipeline, mField, "T", vec_temp_ptr, "SETTEMP", Sev::verbose);
1011
1013 };
1014
1015 auto add_boundary_rhs_ops = [&](auto &pipeline) {
1017
1019
1021 pipeline, mField, "U", {time_scale, time_force_scaling}, "FORCE",
1022 Sev::inform);
1023
1024 // Temperature BC
1025
1026 using OpTemperatureBC =
1029 pipeline, mField, "FLUX", {time_scale, time_temperature_scaling},
1030 "TEMPERATURE", Sev::inform);
1031
1032 // Note: fluxes for temperature should be aggregated. Here separate is
1033 // NaturalMeshsetType<HEATFLUXSET>, we should add version with BLOCKSET,
1034 // convection, see example tutorials/vec-0/src/NaturalBoundaryBC.hpp.
1035 // and vec-0/elastic.cpp
1036
1037 using OpFluxBC =
1040 pipeline, mField, "TBC", {time_scale, time_flux_scaling}, "FLUX",
1041 Sev::inform);
1042
1044 using OpConvectionRhsBC =
1045 T::OpFlux<ThermoElasticOps::ConvectionBcType<BLOCKSET>, 1, 1>;
1046 using OpRadiationRhsBC =
1047 T::OpFlux<ThermoElasticOps::RadiationBcType<BLOCKSET>, 1, 1>;
1048 auto temp_bc_ptr = boost::make_shared<VectorDouble>();
1049 pipeline.push_back(new OpCalculateScalarFieldValues("TBC", temp_bc_ptr));
1050 T::AddFluxToPipeline<OpConvectionRhsBC>::add(
1051 pipeline, mField, "TBC", temp_bc_ptr, "CONVECTION", Sev::inform);
1052 T::AddFluxToPipeline<OpRadiationRhsBC>::add(
1053 pipeline, mField, "TBC", temp_bc_ptr, "RADIATION", Sev::inform);
1054
1055 auto mat_flux_ptr = boost::make_shared<MatrixDouble>();
1056 pipeline.push_back(
1057 new OpCalculateHVecVectorField<3, SPACE_DIM>("FLUX", mat_flux_ptr));
1058
1059 // This is temporary implementation. It be moved to LinearFormsIntegrators,
1060 // however for hybridised case, what is goal of this changes, such function
1061 // is implemented for fluxes in broken space. Thus ultimately this operator
1062 // would be not needed.
1063
1064 struct OpTBCTimesNormalFlux
1065 : public FormsIntegrators<BoundaryEleOp>::Assembly<PETSC>::OpBase {
1066
1068
1069 OpTBCTimesNormalFlux(const std::string field_name,
1070 boost::shared_ptr<MatrixDouble> vec,
1071 boost::shared_ptr<Range> ents_ptr = nullptr)
1072 : OP(field_name, field_name, OP::OPROW, ents_ptr), sourceVec(vec) {}
1073
1074 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data) {
1077 // get integration weights
1078 auto t_w = OP::getFTensor0IntegrationWeight();
1079 // get base function values on rows
1080 auto t_row_base = row_data.getFTensor0N();
1081 // get normal vector
1082 auto t_normal = OP::getFTensor1NormalsAtGaussPts();
1083 // get vector values
1084 auto t_vec = getFTensor1FromMat<SPACE_DIM, 1>(*sourceVec);
1085 // loop over integration points
1086 for (int gg = 0; gg != OP::nbIntegrationPts; gg++) {
1087 // take into account Jacobian
1088 const double alpha = t_w * (t_vec(i) * t_normal(i));
1089 // loop over rows base functions
1090 int rr = 0;
1091 for (; rr != OP::nbRows; ++rr) {
1092 OP::locF[rr] += alpha * t_row_base;
1093 ++t_row_base;
1094 }
1095 for (; rr < OP::nbRowBaseFunctions; ++rr)
1096 ++t_row_base;
1097 ++t_w; // move to another integration weight
1098 ++t_vec;
1099 ++t_normal;
1100 }
1101 EntityType fe_type = OP::getNumeredEntFiniteElementPtr()->getEntType();
1102 if (fe_type == MBTRI) {
1103 OP::locF /= 2;
1104 }
1106 }
1107
1108 protected:
1109 boost::shared_ptr<MatrixDouble> sourceVec;
1110 };
1111 pipeline.push_back(new OpTBCTimesNormalFlux("TBC", mat_flux_ptr));
1112
1113 struct OpBaseNormalTimesTBC
1114 : public FormsIntegrators<BoundaryEleOp>::Assembly<PETSC>::OpBase {
1115
1117
1118 OpBaseNormalTimesTBC(const std::string field_name,
1119 boost::shared_ptr<VectorDouble> vec,
1120 boost::shared_ptr<Range> ents_ptr = nullptr)
1121 : OP(field_name, field_name, OP::OPROW, ents_ptr), sourceVec(vec) {}
1122
1123 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data) {
1126 // get integration weights
1127 auto t_w = OP::getFTensor0IntegrationWeight();
1128 // get base function values on rows
1129 auto t_row_base = row_data.getFTensor1N<3>();
1130 // get normal vector
1131 auto t_normal = OP::getFTensor1NormalsAtGaussPts();
1132 // get vector values
1133 auto t_vec = getFTensor0FromVec(*sourceVec);
1134 // loop over integration points
1135 for (int gg = 0; gg != OP::nbIntegrationPts; gg++) {
1136 // take into account Jacobian
1137 const double alpha = t_w * t_vec;
1138 // loop over rows base functions
1139 int rr = 0;
1140 for (; rr != OP::nbRows; ++rr) {
1141 OP::locF[rr] += alpha * (t_row_base(i) * t_normal(i));
1142 ++t_row_base;
1143 }
1144 for (; rr < OP::nbRowBaseFunctions; ++rr)
1145 ++t_row_base;
1146 ++t_w; // move to another integration weight
1147 ++t_vec;
1148 ++t_normal;
1149 }
1150 EntityType fe_type = OP::getNumeredEntFiniteElementPtr()->getEntType();
1151 if (fe_type == MBTRI) {
1152 OP::locF /= 2;
1153 }
1155 }
1156
1157 protected:
1158 boost::shared_ptr<VectorDouble> sourceVec;
1159 };
1160
1161 pipeline.push_back(new OpBaseNormalTimesTBC("FLUX", temp_bc_ptr));
1162
1164 };
1165
1166 auto add_boundary_lhs_ops = [&](auto &pipeline) {
1168
1170
1171 using T =
1172 NaturalBC<BoundaryEleOp>::template Assembly<PETSC>::BiLinearForm<GAUSS>;
1173 using OpConvectionLhsBC =
1174 T::OpFlux<ThermoElasticOps::ConvectionBcType<BLOCKSET>, 1, 1>;
1175 using OpRadiationLhsBC =
1176 T::OpFlux<ThermoElasticOps::RadiationBcType<BLOCKSET>, 1, 1>;
1177 auto temp_bc_ptr = boost::make_shared<VectorDouble>();
1178 pipeline.push_back(new OpCalculateScalarFieldValues("TBC", temp_bc_ptr));
1179 T::AddFluxToPipeline<OpConvectionLhsBC>::add(pipeline, mField, "TBC", "TBC",
1180 "CONVECTION", Sev::verbose);
1181 T::AddFluxToPipeline<OpRadiationLhsBC>::add(
1182 pipeline, mField, "TBC", "TBC", temp_bc_ptr, "RADIATION", Sev::verbose);
1183
1184 // This is temporary implementation. It be moved to LinearFormsIntegrators,
1185 // however for hybridised case, what is goal of this changes, such function
1186 // is implemented for fluxes in broken space. Thus ultimately this operator
1187 // would be not needed.
1188
1189 struct OpTBCTimesNormalFlux
1190 : public FormsIntegrators<BoundaryEleOp>::Assembly<PETSC>::OpBase {
1191
1193
1194 OpTBCTimesNormalFlux(const std::string row_field_name,
1195 const std::string col_field_name,
1196 boost::shared_ptr<Range> ents_ptr = nullptr)
1197 : OP(row_field_name, col_field_name, OP::OPROWCOL, ents_ptr) {
1198 this->sYmm = false;
1199 this->assembleTranspose = true;
1200 this->onlyTranspose = false;
1201 }
1202
1203 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data,
1204 EntitiesFieldData::EntData &col_data) {
1206
1208
1209 // get integration weights
1210 auto t_w = OP::getFTensor0IntegrationWeight();
1211 // get base function values on rows
1212 auto t_row_base = row_data.getFTensor0N();
1213 // get normal vector
1214 auto t_normal = OP::getFTensor1NormalsAtGaussPts();
1215 // loop over integration points
1216 for (int gg = 0; gg != OP::nbIntegrationPts; gg++) {
1217 // loop over rows base functions
1218 auto a_mat_ptr = &*OP::locMat.data().begin();
1219 int rr = 0;
1220 for (; rr != OP::nbRows; rr++) {
1221 // take into account Jacobian
1222 const double alpha = t_w * t_row_base;
1223 // get column base functions values at gauss point gg
1224 auto t_col_base = col_data.getFTensor1N<3>(gg, 0);
1225 // loop over columns
1226 for (int cc = 0; cc != OP::nbCols; cc++) {
1227 // calculate element of local matrix
1228 // using L2 norm of normal vector for convenient area calculation
1229 // for quads, tris etc.
1230 *a_mat_ptr += alpha * (t_col_base(i) * t_normal(i));
1231 ++t_col_base;
1232 ++a_mat_ptr;
1233 }
1234 ++t_row_base;
1235 }
1236 for (; rr < OP::nbRowBaseFunctions; ++rr)
1237 ++t_row_base;
1238 ++t_normal;
1239 ++t_w; // move to another integration weight
1240 }
1241 EntityType fe_type = OP::getNumeredEntFiniteElementPtr()->getEntType();
1242 if (fe_type == MBTRI) {
1243 OP::locMat /= 2;
1244 }
1246 }
1247 };
1248
1249 pipeline.push_back(new OpTBCTimesNormalFlux("TBC", "FLUX"));
1250
1252 };
1253
1254 // Set BCs by eliminating degrees of freedom
1255 auto get_bc_hook_rhs = [&]() {
1257 mField, pipeline_mng->getDomainRhsFE(),
1258 {time_scale, time_displacement_scaling});
1259 return hook;
1260 };
1261 auto get_bc_hook_lhs = [&]() {
1263 mField, pipeline_mng->getDomainLhsFE(),
1264 {time_scale, time_displacement_scaling});
1265 return hook;
1266 };
1267
1268 pipeline_mng->getDomainRhsFE()->preProcessHook = get_bc_hook_rhs();
1269 pipeline_mng->getDomainLhsFE()->preProcessHook = get_bc_hook_lhs();
1270
1271 CHKERR add_domain_rhs_ops(pipeline_mng->getOpDomainRhsPipeline());
1272 CHKERR add_domain_lhs_ops(pipeline_mng->getOpDomainLhsPipeline());
1273 CHKERR add_boundary_rhs_ops(pipeline_mng->getOpBoundaryRhsPipeline());
1274 CHKERR add_boundary_lhs_ops(pipeline_mng->getOpBoundaryLhsPipeline());
1275
1277}
1278//! [Push operators to pipeline]
1279
1280//! [Solve]
1283
1286 ISManager *is_manager = mField.getInterface<ISManager>();
1287
1288 auto dm = simple->getDM();
1289 auto solver = pipeline_mng->createTSIM();
1290 auto snes_ctx_ptr = getDMSnesCtx(dm);
1291
1292 auto set_section_monitor = [&](auto solver) {
1294 SNES snes;
1295 CHKERR TSGetSNES(solver, &snes);
1296 CHKERR SNESMonitorSet(snes,
1297 (MoFEMErrorCode(*)(SNES, PetscInt, PetscReal,
1298 void *))MoFEMSNESMonitorFields,
1299 (void *)(snes_ctx_ptr.get()), nullptr);
1301 };
1302
1303 auto create_post_process_elements = [&]() {
1304 auto block_thermal_params = boost::make_shared<BlockedThermalParameters>();
1305
1306 auto block_thermoelastic_params =
1307 boost::make_shared<BlockedThermoElasticParameters>();
1308 auto coeff_expansion_ptr =
1309 block_thermoelastic_params->getCoeffExpansionPtr();
1310 auto ref_temp_ptr = block_thermoelastic_params->getRefTempPtr();
1311
1312 auto u_ptr = boost::make_shared<MatrixDouble>();
1313 auto mat_grad_ptr = boost::make_shared<MatrixDouble>();
1314 auto mat_strain_ptr = boost::make_shared<MatrixDouble>();
1315 auto mat_stress_ptr = boost::make_shared<MatrixDouble>();
1316 auto vec_temp_ptr = boost::make_shared<VectorDouble>();
1317 auto mat_flux_ptr = boost::make_shared<MatrixDouble>();
1318
1319 auto push_domain_ops = [&](auto &pp_fe) {
1321 auto &pip = pp_fe->getOpPtrVector();
1322
1323 CHKERR addMatThermalBlockOps(pip, "MAT_THERMAL", block_thermal_params,
1324 Sev::verbose);
1326 pip, "MAT_THERMOELASTIC", block_thermoelastic_params, Sev::verbose);
1327
1329
1330 pip.push_back(new OpCalculateScalarFieldValues("T", vec_temp_ptr));
1331 pip.push_back(
1332 new OpCalculateHVecVectorField<3, SPACE_DIM>("FLUX", mat_flux_ptr));
1333
1334 pip.push_back(new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
1336 "U", mat_grad_ptr));
1337 auto elastic_common_ptr = commonDataFactory<SPACE_DIM, IT, DomainEleOp>(
1338 mField, pip, "U", "MAT_ELASTIC", Sev::inform);
1340 pip.push_back(
1341 new
1342 typename H::template OpCalculateHenckyThermalStress<SPACE_DIM, IT, 0>(
1343 "U", vec_temp_ptr, elastic_common_ptr, coeff_expansion_ptr,
1344 ref_temp_ptr));
1345 if (!IS_LARGE_STRAINS) {
1346 pip.push_back(new OpSymmetrizeTensor<SPACE_DIM>(
1347 elastic_common_ptr->getMatFirstPiolaStress(), mat_stress_ptr));
1348 pip.push_back(
1349 new OpSymmetrizeTensor<SPACE_DIM>(mat_grad_ptr, mat_strain_ptr));
1350 } else {
1351 mat_stress_ptr = elastic_common_ptr->getMatFirstPiolaStress();
1352 mat_strain_ptr = elastic_common_ptr->getMatLogC();
1353 }
1354
1356 };
1357
1358 auto push_post_proc_ops = [&](auto &pp_fe) {
1360 auto &pip = pp_fe->getOpPtrVector();
1362
1363 if (!IS_LARGE_STRAINS) {
1364 pip.push_back(
1365
1366 new OpPPMap(
1367
1368 pp_fe->getPostProcMesh(), pp_fe->getMapGaussPts(),
1369
1370 {{"T", vec_temp_ptr}},
1371
1372 {{"U", u_ptr}, {"FLUX", mat_flux_ptr}},
1373
1374 {},
1375
1376 {{"CAUCHY", mat_stress_ptr}, {"STRAIN", mat_strain_ptr}}
1377
1378 )
1379
1380 );
1381 } else {
1382 pip.push_back(
1383
1384 new OpPPMap(
1385
1386 pp_fe->getPostProcMesh(), pp_fe->getMapGaussPts(),
1387
1388 {{"T", vec_temp_ptr}},
1389
1390 {{"U", u_ptr}, {"FLUX", mat_flux_ptr}},
1391
1392 {{"PIOLA", mat_stress_ptr}},
1393
1394 {{"HENCKY_STRAIN", mat_strain_ptr}}
1395
1396 )
1397
1398 );
1399 }
1400
1402 };
1403
1404 auto domain_post_proc = [&]() {
1405 if (do_output_domain == PETSC_FALSE)
1406 return boost::shared_ptr<PostProcEle>();
1407 auto pp_fe = boost::make_shared<PostProcEle>(mField);
1408 CHK_MOAB_THROW(push_domain_ops(pp_fe),
1409 "push domain ops to domain element");
1410 CHK_MOAB_THROW(push_post_proc_ops(pp_fe),
1411 "push post proc ops to domain element");
1412 return pp_fe;
1413 };
1414
1415 auto skin_post_proc = [&]() {
1416 if (do_output_skin == PETSC_FALSE)
1417 return boost::shared_ptr<SkinPostProcEle>();
1418 auto pp_fe = boost::make_shared<SkinPostProcEle>(mField);
1419 auto simple = mField.getInterface<Simple>();
1420 auto op_side = new OpLoopSide<SideEle>(mField, simple->getDomainFEName(),
1421 SPACE_DIM, Sev::verbose);
1422 CHK_MOAB_THROW(push_domain_ops(op_side),
1423 "push domain ops to side element");
1424 pp_fe->getOpPtrVector().push_back(op_side);
1425 CHK_MOAB_THROW(push_post_proc_ops(pp_fe),
1426 "push post proc ops to skin element");
1427 return pp_fe;
1428 };
1429
1430 return std::make_pair(domain_post_proc(), skin_post_proc());
1431 };
1432
1433 auto monitor_ptr = boost::make_shared<FEMethod>();
1434
1435 auto set_time_monitor = [&](auto dm, auto solver, auto domain_post_proc_fe,
1436 auto skin_post_proc_fe) {
1438 monitor_ptr->preProcessHook = [&]() {
1440
1441 if (save_every && (monitor_ptr->ts_step % save_every == 0)) {
1442 if (do_output_domain) {
1443 CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(),
1444 domain_post_proc_fe,
1445 monitor_ptr->getCacheWeakPtr());
1446 CHKERR domain_post_proc_fe->writeFile(
1447 "out_" + boost::lexical_cast<std::string>(monitor_ptr->ts_step) +
1448 ".h5m");
1449 }
1450 if (do_output_skin) {
1451 CHKERR DMoFEMLoopFiniteElements(dm, simple->getBoundaryFEName(),
1452 skin_post_proc_fe,
1453 monitor_ptr->getCacheWeakPtr());
1454 CHKERR skin_post_proc_fe->writeFile(
1455 "out_skin_" +
1456 boost::lexical_cast<std::string>(monitor_ptr->ts_step) + ".h5m");
1457 }
1458 }
1459
1460 if (doEvalField) {
1461
1463 ->evalFEAtThePoint<SPACE_DIM>(
1464 fieldEvalCoords.data(), 1e-12, simple->getProblemName(),
1465 simple->getDomainFEName(), fieldEvalData,
1466 mField.get_comm_rank(), mField.get_comm_rank(), nullptr,
1467 MF_EXIST, QUIET);
1468
1469 if (atom_test) {
1470 auto eval_num_vec =
1471 createVectorMPI(mField.get_comm(), PETSC_DECIDE, 1);
1472 CHKERR VecZeroEntries(eval_num_vec);
1473 if (tempFieldPtr->size()) {
1474 CHKERR VecSetValue(eval_num_vec, 0, 1, ADD_VALUES);
1475 }
1476 CHKERR VecAssemblyBegin(eval_num_vec);
1477 CHKERR VecAssemblyEnd(eval_num_vec);
1478
1479 double eval_num;
1480 CHKERR VecSum(eval_num_vec, &eval_num);
1481 if (!(int)eval_num) {
1482 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1483 "atom test %d failed: did not find elements to evaluate "
1484 "the field, check the coordinates",
1485 atom_test);
1486 }
1487 }
1488
1489 if (tempFieldPtr->size()) {
1490 auto t_temp = getFTensor0FromVec(*tempFieldPtr);
1491 MOFEM_LOG("ThermoElasticSync", Sev::inform)
1492 << "Eval point T: " << t_temp;
1493 if (atom_test && fabs(monitor_ptr->ts_t - 10) < 1e-12) {
1494 if (atom_test <= 3 && fabs(t_temp - 554.48) > 1e-2) {
1495 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1496 "atom test %d failed: wrong temperature value",
1497 atom_test);
1498 }
1499 if (atom_test == 4 && fabs(t_temp - 250) > 1e-2) {
1500 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1501 "atom test %d failed: wrong temperature value",
1502 atom_test);
1503 }
1504 if (atom_test == 5 && fabs(t_temp - 1) > 1e-2) {
1505 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1506 "atom test %d failed: wrong temperature value",
1507 atom_test);
1508 }
1509 }
1510 }
1511 if (fluxFieldPtr->size1()) {
1513 auto t_flux = getFTensor1FromMat<SPACE_DIM>(*fluxFieldPtr);
1514 auto flux_mag = sqrt(t_flux(i) * t_flux(i));
1515 MOFEM_LOG("ThermoElasticSync", Sev::inform)
1516 << "Eval point FLUX magnitude: " << flux_mag;
1517 if (atom_test && fabs(monitor_ptr->ts_t - 10) < 1e-12) {
1518 if (atom_test <= 3 && fabs(flux_mag - 27008.0) > 2e1) {
1519 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1520 "atom test %d failed: wrong flux value", atom_test);
1521 }
1522 if (atom_test == 4 && fabs(flux_mag - 150e3) > 1e-1) {
1523 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1524 "atom test %d failed: wrong flux value", atom_test);
1525 }
1526 if (atom_test == 5 && fabs(flux_mag) > 1e-6) {
1527 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1528 "atom test %d failed: wrong flux value", atom_test);
1529 }
1530 }
1531 }
1532 if (dispFieldPtr->size1()) {
1534 auto t_disp = getFTensor1FromMat<SPACE_DIM>(*dispFieldPtr);
1535 auto disp_mag = sqrt(t_disp(i) * t_disp(i));
1536 MOFEM_LOG("ThermoElasticSync", Sev::inform)
1537 << "Eval point U magnitude: " << disp_mag;
1538 if (atom_test && fabs(monitor_ptr->ts_t - 10) < 1e-12) {
1539 if (atom_test == 1 && fabs(disp_mag - 0.00345) > 1e-5) {
1540 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1541 "atom test %d failed: wrong displacement value",
1542 atom_test);
1543 }
1544 if ((atom_test == 2 || atom_test == 3) &&
1545 fabs(disp_mag - 0.00265) > 1e-5) {
1546 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1547 "atom test %d failed: wrong displacement value",
1548 atom_test);
1549 }
1550 if ((atom_test == 5) &&
1551 fabs(t_disp(0) - (std::sqrt(std::exp(0.2)) - 1)) > 1e-5 &&
1552 fabs(t_disp(1) - (std::sqrt(std::exp(0.2)) - 1)) > 1e-5) {
1553 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1554 "atom test %d failed: wrong displacement value",
1555 atom_test);
1556 }
1557 }
1558 }
1559 if (strainFieldPtr->size1()) {
1561 auto t_strain =
1562 getFTensor2SymmetricFromMat<SPACE_DIM>(*strainFieldPtr);
1563 auto t_strain_trace = t_strain(i, i);
1564 if (atom_test && fabs(monitor_ptr->ts_t - 10) < 1e-12) {
1565 if (atom_test == 1 && fabs(t_strain_trace - 0.00679) > 1e-5) {
1566 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1567 "atom test %d failed: wrong strain value", atom_test);
1568 }
1569 if ((atom_test == 2 || atom_test == 3) &&
1570 fabs(t_strain_trace - 0.00522) > 1e-5) {
1571 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1572 "atom test %d failed: wrong strain value", atom_test);
1573 }
1574 if ((atom_test == 5) && fabs(t_strain_trace - 0.2) > 1e-5) {
1575 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1576 "atom test %d failed: wrong strain value", atom_test);
1577 }
1578 }
1579 }
1580 if (stressFieldPtr->size1()) {
1581 double von_mises_stress;
1582 auto getVonMisesStress = [&](auto t_stress) {
1584 von_mises_stress = std::sqrt(
1585 0.5 *
1586 ((t_stress(0, 0) - t_stress(1, 1)) *
1587 (t_stress(0, 0) - t_stress(1, 1)) +
1588 (SPACE_DIM == 3 ? (t_stress(1, 1) - t_stress(2, 2)) *
1589 (t_stress(1, 1) - t_stress(2, 2))
1590 : 0) +
1591 (SPACE_DIM == 3 ? (t_stress(2, 2) - t_stress(0, 0)) *
1592 (t_stress(2, 2) - t_stress(0, 0))
1593 : 0) +
1594 6 * (t_stress(0, 1) * t_stress(0, 1) +
1595 (SPACE_DIM == 3 ? t_stress(1, 2) * t_stress(1, 2) : 0) +
1596 (SPACE_DIM == 3 ? t_stress(2, 0) * t_stress(2, 0) : 0))));
1598 };
1599 if (!IS_LARGE_STRAINS) {
1600 auto t_stress =
1601 getFTensor2SymmetricFromMat<SPACE_DIM>(*stressFieldPtr);
1602 CHKERR getVonMisesStress(t_stress);
1603 } else {
1604 auto t_stress =
1605 getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*stressFieldPtr);
1606 CHKERR getVonMisesStress(t_stress);
1607 }
1608 MOFEM_LOG("ThermoElasticSync", Sev::inform)
1609 << "Eval point von Mises Stress: " << von_mises_stress;
1610 if (atom_test && fabs(monitor_ptr->ts_t - 10) < 1e-12) {
1611 if (atom_test == 1 && fabs(von_mises_stress - 523.0) > 5e-1) {
1612 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1613 "atom test %d failed: wrong von Mises stress value",
1614 atom_test);
1615 }
1616 if (atom_test == 2 && fabs(von_mises_stress - 16.3) > 5e-2) {
1617 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1618 "atom test %d failed: wrong von Mises stress value",
1619 atom_test);
1620 }
1621 if (atom_test == 3 && fabs(von_mises_stress - 14.9) > 5e-2) {
1622 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1623 "atom test %d failed: wrong von Mises stress value",
1624 atom_test);
1625 }
1626 if (atom_test == 5 && fabs(von_mises_stress) > 5e-2) {
1627 SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1628 "atom test %d failed: wrong von Mises stress value",
1629 atom_test);
1630 }
1631 }
1632 }
1633
1635 }
1636
1638 };
1639 auto null = boost::shared_ptr<FEMethod>();
1640 CHKERR DMMoFEMTSSetMonitor(dm, solver, simple->getDomainFEName(), null,
1641 monitor_ptr, null);
1643 };
1644
1645 auto set_fieldsplit_preconditioner = [&](auto solver) {
1647
1648 SNES snes;
1649 CHKERR TSGetSNES(solver, &snes);
1650 KSP ksp;
1651 CHKERR SNESGetKSP(snes, &ksp);
1652 PC pc;
1653 CHKERR KSPGetPC(ksp, &pc);
1654 PetscBool is_pcfs = PETSC_FALSE;
1655 PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &is_pcfs);
1656
1657 // Setup fieldsplit (block) solver - optional: yes/no
1658 if (is_pcfs == PETSC_TRUE) {
1659 auto bc_mng = mField.getInterface<BcManager>();
1660 auto is_mng = mField.getInterface<ISManager>();
1661 auto name_prb = simple->getProblemName();
1662
1663 SmartPetscObj<IS> is_u;
1664 CHKERR is_mng->isCreateProblemFieldAndRank(name_prb, ROW, "U", 0,
1665 SPACE_DIM, is_u);
1666 SmartPetscObj<IS> is_flux;
1667 CHKERR is_mng->isCreateProblemFieldAndRank(name_prb, ROW, "FLUX", 0, 0,
1668 is_flux);
1669 SmartPetscObj<IS> is_T;
1670 CHKERR is_mng->isCreateProblemFieldAndRank(name_prb, ROW, "T", 0, 0,
1671 is_T);
1672 SmartPetscObj<IS> is_TBC;
1673 CHKERR is_mng->isCreateProblemFieldAndRank(name_prb, ROW, "TBC", 0, 0,
1674 is_TBC);
1675 IS is_tmp, is_tmp2;
1676 CHKERR ISExpand(is_T, is_flux, &is_tmp);
1677 CHKERR ISExpand(is_TBC, is_tmp, &is_tmp2);
1678 CHKERR ISDestroy(&is_tmp);
1679 auto is_TFlux = SmartPetscObj<IS>(is_tmp2);
1680
1681 CHKERR ISSort(is_u);
1682 CHKERR ISSort(is_TFlux);
1683 CHKERR PCFieldSplitSetIS(pc, PETSC_NULL, is_TFlux);
1684 CHKERR PCFieldSplitSetIS(pc, PETSC_NULL, is_u);
1685 }
1686
1688 };
1689
1690 auto D = createDMVector(dm);
1691 CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_FORWARD);
1692 CHKERR TSSetSolution(solver, D);
1693 CHKERR TSSetFromOptions(solver);
1694
1695 CHKERR set_section_monitor(solver);
1696 CHKERR set_fieldsplit_preconditioner(solver);
1697
1698 auto [domain_post_proc_fe, skin_post_proc_fe] =
1699 create_post_process_elements();
1700 CHKERR set_time_monitor(dm, solver, domain_post_proc_fe, skin_post_proc_fe);
1701
1702 CHKERR TSSetUp(solver);
1703 CHKERR TSSolve(solver, NULL);
1704
1706}
1707//! [Solve]
1708
1709static char help[] = "...\n\n";
1710
1711int main(int argc, char *argv[]) {
1712
1713 // Initialisation of MoFEM/PETSc and MOAB data structures
1714 const char param_file[] = "param_file.petsc";
1715 MoFEM::Core::Initialize(&argc, &argv, param_file, help);
1716
1717 // Add logging channel for example
1718 auto core_log = logging::core::get();
1719 core_log->add_sink(
1721 LogManager::setLog("ThermoElastic");
1722 MOFEM_LOG_TAG("ThermoElastic", "ThermoElastic");
1723
1724 core_log->add_sink(
1725 LogManager::createSink(LogManager::getStrmSync(), "ThermoElasticSync"));
1726 LogManager::setLog("ThermoElasticSync");
1727 MOFEM_LOG_TAG("ThermoElasticSync", "ThermoElasticSync");
1728
1729 try {
1730
1731 //! [Register MoFEM discrete manager in PETSc]
1732 DMType dm_name = "DMMOFEM";
1733 CHKERR DMRegister_MoFEM(dm_name);
1734 //! [Register MoFEM discrete manager in PETSc
1735
1736 //! [Create MoAB]
1737 moab::Core mb_instance; ///< mesh database
1738 moab::Interface &moab = mb_instance; ///< mesh database interface
1739 //! [Create MoAB]
1740
1741 //! [Create MoFEM]
1742 MoFEM::Core core(moab); ///< finite element database
1743 MoFEM::Interface &m_field = core; ///< finite element database interface
1744 //! [Create MoFEM]
1745
1746 //! [Load mesh]
1747 Simple *simple = m_field.getInterface<Simple>();
1749 CHKERR simple->loadFile();
1750 //! [Load mesh]
1751
1752 //! [ThermoElasticProblem]
1753 ThermoElasticProblem ex(m_field);
1754 CHKERR ex.runProblem();
1755 //! [ThermoElasticProblem]
1756 }
1758
1760}
static auto filter_true_skin(MoFEM::Interface &m_field, Range &&skin)
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
#define MOFEM_LOG_SYNCHRONISE(comm)
Synchronise "SYNC" channel.
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
#define FTENSOR_INDEX(DIM, I)
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
int main()
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
constexpr int SPACE_DIM
[Define dimension]
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
DomainEle::UserDataOperator DomainEleOp
Finire element operator type.
@ QUIET
@ ROW
#define CATCH_ERRORS
Catch errors.
@ MF_EXIST
FieldApproximationBase
approximation base
Definition definitions.h:58
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ L2
field with C-1 continuity
Definition definitions.h:88
@ H1
continuous field
Definition definitions.h:85
@ NOSPACE
Definition definitions.h:83
@ HCURL
field with continuous tangents
Definition definitions.h:86
@ HDIV
field with continuous normal traction
Definition definitions.h:87
#define MYPCOMM_INDEX
default communicator number PCOMM
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
@ TEMPERATURESET
@ HEATFLUXSET
@ NODESET
@ SIDESET
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
auto integration_rule
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:514
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition DMMoFEM.cpp:576
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1234
IntegrationType
Form integrator integration types.
AssemblyType
[Storage and set boundary conditions]
@ GAUSS
Gaussian quadrature integration.
@ PETSC
Standard PETSc assembly.
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#define MOFEM_LOG(channel, severity)
Log.
SeverityLevel
Severity levels.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
MoFEMErrorCode getCubitMeshsetPtr(const int ms_id, const CubitBCType cubit_bc_type, const CubitMeshSets **cubit_meshset_ptr) const
get cubit meshset
FTensor::Index< 'i', SPACE_DIM > i
PetscBool doEvalField
const double c
speed of light (cm/ns)
double D
const double n
refractive index of diffusive medium
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasVector< double > VectorDouble
Definition Types.hpp:68
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, 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 Monitor To TS solver.
Definition DMMoFEM.cpp:1046
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
MoFEMErrorCode MoFEMSNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *ctx)
Sens monitor printing residual field by field.
Definition SnesCtx.cpp:596
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)
PetscErrorCode PetscOptionsGetRealArray(PetscOptions *, const char pre[], const char name[], PetscReal dval[], PetscInt *nmax, PetscBool *set)
auto createVectorMPI(MPI_Comm comm, PetscInt n, PetscInt N)
Create MPI Vector.
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
OpCalculateScalarFieldValuesFromPetscVecImpl< PetscData::CTX_SET_X_T > OpCalculateScalarFieldValuesDot
auto get_temp_meshset_ptr(moab::Interface &moab)
Create smart pointer to temporary meshset.
auto getDMSnesCtx(DM dm)
Get SNES context data structure used by DM.
Definition DMMoFEM.hpp:1262
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMixDivTimesScalar< SPACE_DIM > OpHdivT
Integrate Lhs div of base of flux times base of temperature (FLUX x T) and transpose of it,...
OpHdivHdivImpl< SPACE_DIM, IS_LARGE_STRAINS > OpHdivHdiv
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, 1 > OpCapacity
Integrate Lhs base of temperature times (heat capacity) times base of temperature (T x T)
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesScalar< 1 > OpBaseDotT
Integrate Rhs base of temperature time heat capacity times heat rate (T)
OpBaseDotT OpBaseDivFlux
Integrate Rhs base of temperature times divergence of flux (T)
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpMixDivTimesU< 3, 1, SPACE_DIM > OpHDivTemp
Integrate Rhs div flux base times temperature (T)
OpCalculateQdotQLhs_dU< SPACE_DIM, GAUSS, AssemblyDomainEleOp, IS_LARGE_STRAINS > OpHdivU
Integrate Lhs of flux term coupled to displacement field.
OpHdivFluxImpl< SPACE_DIM, IS_LARGE_STRAINS > OpHdivFlux
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
constexpr IntegrationType I
constexpr AssemblyType A
[Define dimension]
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
constexpr auto field_name
static constexpr int approx_order
OpBaseImpl< PETSC, EdgeEleOp > OpBase
Definition radiation.cpp:29
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, SPACE_DIM, SPACE_DIM > OpInternalForcePiola
Definition seepage.cpp:65
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradTensorGrad< 1, SPACE_DIM, SPACE_DIM, 1 > OpKPiola
[Only used for dynamics]
Definition seepage.cpp:63
FTensor::Index< 'm', 3 > m
Add operators pushing bases from local to physical configuration.
Boundary condition manager for finite element problem setup.
Managing BitRefLevels.
virtual moab::Interface & get_moab()=0
virtual MPI_Comm & get_comm() const =0
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:118
Deprecated interface functions.
Definition of the displacement bc data structure.
Definition BCData.hpp:72
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
Class (Function) to enforce essential constrains.
Definition Essential.hpp:25
Structure for user loop methods on finite elements.
Basic algebra on fields.
Definition FieldBlas.hpp:21
Field evaluator interface.
Definition of the heat flux bc data structure.
Definition BCData.hpp:423
Section manager is used to create indexes and sections.
Definition ISManager.hpp:23
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.
Interface for managing meshsets containing materials and boundary conditions.
Assembly methods.
Definition Natural.hpp:65
Natural boundary conditions.
Definition Natural.hpp:57
Get vector field for H-div approximation.
Calculate divergence of vector field.
Specialization for double precision scalar field values calculation.
Get field gradients at integration pts for scalar field rank 0, i.e. vector field.
Specialization for double precision vector field values calculation.
Enforce essential constrains on lhs.
Definition Essential.hpp:81
Enforce essential constrains on rhs.
Definition Essential.hpp:65
Element used to execute operators on side of the element.
Post post-proc data at points from hash maps.
Operator for symmetrizing tensor fields.
Template struct for dimension-specific finite element types.
PipelineManager interface.
Problem manager is used to build and partition problems.
Simple interface for fast problem set-up.
Definition Simple.hpp:27
MoFEMErrorCode addDomainField(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_ZERO, int verb=-1)
Add field on domain.
Definition Simple.cpp:261
MoFEMErrorCode getOptions()
get options
Definition Simple.cpp:180
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
MoFEMErrorCode getCommandLineParameters()
MoFEMErrorCode runProblem()
std::array< double, SPACE_DIM > fieldEvalCoords
boost::shared_ptr< MatrixDouble > dispFieldPtr
boost::shared_ptr< MatrixDouble > strainFieldPtr
ThermoElasticProblem(MoFEM::Interface &m_field)
std::array< double, 3 > fieldEvalCoords
boost::shared_ptr< VectorDouble > tempFieldPtr
boost::shared_ptr< MatrixDouble > fluxFieldPtr
MoFEMErrorCode opThermoElasticFactoryDomainRhs(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string field_name, boost::shared_ptr< HenckyOps::CommonData > elastic_common_ptr, boost::shared_ptr< ThermoElasticProblem::BlockedThermalParameters > thermal_common_ptr, boost::shared_ptr< ThermoElasticProblem::BlockedThermoElasticParameters > thermoelastic_common_ptr, Sev sev)
MoFEM::Interface & mField
MoFEMErrorCode setupProblem()
add fields
MoFEMErrorCode bC()
boost::shared_ptr< MatrixDouble > stressFieldPtr
MoFEMErrorCode OPs()
MoFEMErrorCode tsSolve()
MoFEMErrorCode opThermoElasticFactoryDomainLhs(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string field_name, std::string coupled_field_name, std::string elastic_block_name, std::string thermal_block_name, std::string thermoelastic_block_name, Sev sev, double scale=1)
MoFEMErrorCode opThermoElasticFactoryDomainLhs(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string field_name, std::string coupled_field_name, boost::shared_ptr< HenckyOps::CommonData > elastic_common_ptr, boost::shared_ptr< ThermoElasticProblem::BlockedThermalParameters > thermal_common_ptr, boost::shared_ptr< ThermoElasticProblem::BlockedThermoElasticParameters > thermoelastic_common_ptr, Sev sev)
MoFEMErrorCode addMatThermalBlockOps(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, std::string block_name, boost::shared_ptr< BlockedThermalParameters > blockedParamsPtr, Sev sev)
MoFEMErrorCode addMatThermoElasticBlockOps(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, std::string block_name, boost::shared_ptr< BlockedThermoElasticParameters > blockedParamsPtr, Sev sev)
boost::shared_ptr< MatrixDouble > dispGradPtr
MoFEMErrorCode opThermoElasticFactoryDomainRhs(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string field_name, std::string elastic_block_name, std::string thermal_block_name, std::string thermoelastic_block_name, Sev sev, double scale=1)
boost::shared_ptr< FieldEvaluatorInterface::SetPtsData > fieldEvalData
ElementsAndOps< SPACE_DIM >::SideEle SideEle
Definition plastic.cpp:61
double scale
Definition plastic.cpp:123
double H
Hardening.
Definition plastic.cpp:128
#define FINITE_DEFORMATION_FLAG
auto save_range
constexpr AssemblyType AT
constexpr IntegrationType IT
static char help[]
[Solve]
int save_every
int order_flux
constexpr bool IS_LARGE_STRAINS
int atom_test
#define EXECUTABLE_DIMENSION
constexpr int SPACE_DIM
double default_ref_temp
double default_heat_capacity
int order_disp
double default_young_modulus
[Essential boundary conditions (Least square approach)]
PetscBool do_output_skin
double default_coeff_expansion
int order_temp
PetscBool do_output_domain
double default_heat_conductivity
double default_init_temp
double default_poisson_ratio