v0.15.5
Loading...
Searching...
No Matches
NonlinearElasticExample.hpp
Go to the documentation of this file.
1/**
2 * @file NonlinearElasticExample.hpp
3 * @example mofem/tutorials/vec-2_nonlinear_elasticity/src/NonlinearElasticExample.hpp
4 * @date 2025-10-5
5 *
6 * @copyright Anonymous authors (c) 2025 under the MIT license (see LICENSE for
7 * details)
8 */
9
10#ifdef WITH_ADOL_C
11 #include <AdolCOps.hpp>
12 #include <AdolCElastic.hpp>
13 #include <AdolCHuHu.hpp>
14#endif
15
17
19
20 MoFEMErrorCode runProblem();
21
22protected:
24 boost::shared_ptr<DomainEle> reactionFe;
25 boost::shared_ptr<PostProcEleDomain> postProcDomainFe;
26 boost::shared_ptr<PostProcEleBdy> postProcBdyFe;
27 PetscBool useAdolcMaterial = PETSC_FALSE;
28
29 MoFEMErrorCode readMesh();
30 MoFEMErrorCode setupProblem();
31 MoFEMErrorCode boundaryCondition();
32 MoFEMErrorCode assembleSystemHencky();
33 MoFEMErrorCode assembleSystemAdolc();
34 MoFEMErrorCode opTest();
35 MoFEMErrorCode TsSolve();
36 MoFEMErrorCode gettingNorms();
37 MoFEMErrorCode outputResults();
38 MoFEMErrorCode checkResults();
39
40 std::vector<Tag>
41 listTagsToTransfer; ///< list of tags to transfer to postprocessor
42
43#ifdef WITH_ADOL_C
44 static inline constexpr int modelType =
45 (SPACE_DIM == 3) ? AdolCOps::MODEL_3D : AdolCOps::MODEL_2D_PLANE_STRAIN;
46 boost::shared_ptr<AdolCOps::PhysicalEquations> physicalEquationsPtr;
47 boost::shared_ptr<AdolCOps::PhysicalEquations> physicalHuHuPtr;
48#endif
49
50 static auto get_skin(MoFEM::Interface &m_field, Range body_ents) {
51 Skinner skin(&m_field.get_moab());
52 Range skin_ents;
53 CHK_MOAB_THROW(skin.find_skin(0, body_ents, false, skin_ents), "find_skin");
54 return skin_ents;
55 };
56
57 static auto filter_true_skin(MoFEM::Interface &m_field, Range &&skin) {
58 Range boundary_ents;
59 ParallelComm *pcomm =
60 ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
61 CHK_MOAB_THROW(pcomm->filter_pstatus(skin,
62 PSTATUS_SHARED | PSTATUS_MULTISHARED,
63 PSTATUS_NOT, -1, &boundary_ents),
64 "filter_pstatus");
65 return boundary_ents;
66 };
67};
68
69//! [Run problem]
75
76 PetscBool test_op = PETSC_FALSE;
77 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-test_op", &test_op,
78 PETSC_NULLPTR);
79
80 if (useAdolcMaterial == PETSC_TRUE) {
82 if (test_op == PETSC_TRUE) {
83 CHKERR opTest();
84 }
89 } else {
91 if (test_op == PETSC_TRUE) {
92 CHKERR opTest();
93 }
98 }
99
100 CHKERR TsSolve();
101
102 if (useAdolcMaterial == PETSC_TRUE) {
103 } else {
105 }
109}
110//! [Run problem]
111
112//! [Read mesh]
115 auto simple = mField.getInterface<Simple>();
116 CHKERR simple->getOptions();
117 CHKERR simple->loadFile();
119}
120//! [Read mesh]
121
122//! [Set up problem]
125 Simple *simple = mField.getInterface<Simple>();
126
127 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOPT };
128 const char *list_bases[LASBASETOPT] = {"ainsworth", "demkowicz"};
129 PetscInt choice_base_value = AINSWORTH;
130 CHKERR PetscOptionsGetEList(PETSC_NULLPTR, NULL, "-base", list_bases,
131 LASBASETOPT, &choice_base_value, PETSC_NULLPTR);
132
134 switch (choice_base_value) {
135 case AINSWORTH:
137 MOFEM_LOG("WORLD", Sev::inform)
138 << "Set AINSWORTH_LEGENDRE_BASE for displacements";
139 break;
140 case DEMKOWICZ:
142 MOFEM_LOG("WORLD", Sev::inform)
143 << "Set DEMKOWICZ_JACOBI_BASE for displacements";
144 break;
145 default:
146 base = LASTBASE;
147 break;
148 }
149
150 // Add field
151 CHKERR simple->addDomainField("U", H1, base, SPACE_DIM);
152 CHKERR simple->addBoundaryField("U", H1, base, SPACE_DIM);
153 CHKERR simple->addDataField("GEOMETRY", H1, base, SPACE_DIM);
154 int order = 2;
155 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-order", &order, PETSC_NULLPTR);
156 CHKERR simple->setFieldOrder("U", order);
157 CHKERR simple->setFieldOrder("GEOMETRY", 2);
158 CHKERR simple->setUp();
159
160 auto project_ho_geometry = [&]() {
161 Projection10NodeCoordsOnField ent_method(mField, "GEOMETRY");
162 return mField.loop_dofs("GEOMETRY", ent_method);
163 };
164 CHKERR project_ho_geometry();
165
166 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-use_adolc_material",
167 &useAdolcMaterial, PETSC_NULLPTR);
168#ifndef WITH_ADOL_C
169 if (useAdolcMaterial == PETSC_TRUE) {
170 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
171 "ADOL-C support is not enabled. Please reconfigure with "
172 "-DWITH_ADOL_C=ON and recompile to use AdolC material model.");
173 }
174#endif
175 if (useAdolcMaterial == PETSC_TRUE) {
176 MOFEM_LOG("WORLD", Sev::inform) << "Use ADOL-C material model";
177 } else {
178 MOFEM_LOG("WORLD", Sev::inform) << "Use Hencky material model";
179 }
180
181 auto tag_meshest = [&]() {
183
184 auto set_block = [&](auto name, int dim) {
185 std::map<int, Range> map;
186 auto set_tag_impl = [&](auto name) {
188 auto mesh_mng = mField.getInterface<MeshsetsManager>();
189 auto bcs = mesh_mng->getCubitMeshsetPtr(
190
191 std::regex((boost::format("%s(.*)") % name).str())
192
193 );
194 std::map<int, int> ids_map;
195 int bit_id = 1;
196 for (auto bc : bcs) {
197 int id = bc->getMeshsetId();
198 ids_map[id] = bit_id;
199 bit_id <<= 1;
200 }
201 for (auto bc : bcs) {
202 Range r;
203 CHKERR bc->getMeshsetIdEntitiesByDimension(mField.get_moab(), dim, r,
204 true);
205 map[ids_map[bc->getMeshsetId()]] = r;
206 MOFEM_LOG("EXAMPLE", Sev::inform)
207 << "Block " << name << " id " << bc->getMeshsetId() << " : "
208 << ids_map[bc->getMeshsetId()] << " has " << r.size()
209 << " entities";
210 }
212 };
213
214 CHKERR set_tag_impl(name);
215
216 return std::make_pair(name, map);
217 };
218
219 auto set_skin = [&](auto &&map) {
220 for (auto &m : map.second) {
221 auto s = filter_true_skin(mField, get_skin(mField, m.second));
222 m.second.swap(s);
223 MOFEM_LOG("EXAMPLE", Sev::inform)
224 << "Skin for block " << map.first << " id " << m.first << " has "
225 << m.second.size() << " entities";
226 }
227 return map;
228 };
229
230 auto set_tag = [&](auto &&map) {
231 Tag th;
232 auto name = map.first;
233 int def_val[] = {-1};
234 CHK_MOAB_THROW(mField.get_moab().tag_get_handle(
235 name, 1, MB_TYPE_INTEGER, th,
236 MB_TAG_SPARSE | MB_TAG_CREAT, def_val),
237 "create tag");
238 for (auto &m : map.second) {
239 int id = m.first;
240 for (auto ent : m.second) {
241 int current_id;
243 mField.get_moab().tag_get_data(th, &ent, 1, &current_id),
244 "get tag data");
245 if (current_id != -1) {
246 id |= current_id;
247 }
248 CHK_MOAB_THROW(mField.get_moab().tag_set_data(th, &ent, 1, &id),
249 "set tag data");
250 }
251 }
252 return th;
253 };
254
255 if (SPACE_DIM == 3) {
256 listTagsToTransfer.push_back(set_tag(set_block("MAT_", 3)));
257 listTagsToTransfer.push_back(set_tag(set_skin(set_block("MAT_", 3))));
258 } else {
259 listTagsToTransfer.push_back(set_tag(set_block("MAT_", 2)));
260 }
261
263 };
264
265 CHKERR tag_meshest();
266
268}
269//! [Set up problem]
270
271//! [Boundary condition]
272
275 auto *pipeline_mng = mField.getInterface<PipelineManager>();
276 auto simple = mField.getInterface<Simple>();
277 auto bc_mng = mField.getInterface<BcManager>();
278 auto time_scale = boost::make_shared<ExampleTimeScale>();
279
280 auto integration_rule = [](int, int, int approx_order) {
281 return 2 * (approx_order - 1) + 1;
282 };
283
284 CHKERR pipeline_mng->setDomainRhsIntegrationRule(integration_rule);
285 CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
286 CHKERR pipeline_mng->setBoundaryRhsIntegrationRule(integration_rule);
287
288 reactionFe = boost::make_shared<DomainEle>(mField);
289 reactionFe->getRuleHook = integration_rule;
290
292 pipeline_mng->getOpBoundaryRhsPipeline(), mField, "U", {time_scale},
293 "FORCE", "PRESSURE", Sev::inform);
294
295 //! [Define gravity vector]
297 pipeline_mng->getOpDomainRhsPipeline(), mField, "U", {time_scale},
298 "BODY_FORCE", Sev::inform);
299
300 // Essential BC
301 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_X",
302 "U", 0, 0);
303 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_Y",
304 "U", 1, 1);
305 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_Z",
306 "U", 2, 2);
307 CHKERR bc_mng->pushMarkDOFsOnEntities<DisplacementCubitBcData>(
308 simple->getProblemName(), "U");
309
311}
312//! [Boundary condition]
313
314//! [Push operators to pipeline]
317 auto *pipeline_mng = mField.getInterface<PipelineManager>();
318
319 auto add_domain_ops_lhs = [&](auto &pip) {
322 "GEOMETRY");
323 CHKERR
324 HenckyOps::opFactoryDomainLhs<SPACE_DIM, PETSC, GAUSS, DomainEleOp>(
325 mField, pip, "U", "MAT_ELASTIC", Sev::inform);
327 };
328
329 auto add_domain_ops_rhs = [&](auto &pip) {
332 "GEOMETRY");
333 CHKERR
334 HenckyOps::opFactoryDomainRhs<SPACE_DIM, PETSC, GAUSS, DomainEleOp>(
335 mField, pip, "U", "MAT_ELASTIC", Sev::inform);
337 };
338
339 CHKERR add_domain_ops_lhs(pipeline_mng->getOpDomainLhsPipeline());
340 CHKERR add_domain_ops_rhs(pipeline_mng->getOpDomainRhsPipeline());
341
342 // push operators to reaction pipeline
343 CHKERR add_domain_ops_rhs(reactionFe->getOpPtrVector());
344 reactionFe->postProcessHook =
345 EssentialPreProcReaction<DisplacementCubitBcData>(mField, reactionFe);
346
347 PetscBool post_proc_vol;
348 PetscBool post_proc_skin;
349
350 if constexpr (SPACE_DIM == 2) {
351 post_proc_vol = PETSC_TRUE;
352 post_proc_skin = PETSC_FALSE;
353 } else {
354 post_proc_vol = PETSC_FALSE;
355 post_proc_skin = PETSC_TRUE;
356 }
357
358 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-post_proc_vol",
359 &post_proc_vol, PETSC_NULLPTR);
360 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-post_proc_skin",
361 &post_proc_skin, PETSC_NULLPTR);
362
363 // Setup postprocessing
364 auto create_post_proc_fe = [&]() {
365 auto post_proc_ele_domain = [this](auto &pip_domain) {
367 "GEOMETRY");
368 auto common_ptr = commonDataFactory<SPACE_DIM, GAUSS, DomainEleOp>(
369 mField, pip_domain, "U", "MAT_ELASTIC", Sev::inform);
370 return common_ptr;
371 };
372
373 auto post_proc_map = [&](auto &pip, auto u_ptr, auto common_ptr) {
375
376 using OpPPMap = OpPostProcMapInMoab<SPACE_DIM, SPACE_DIM>;
377
378 pip->getOpPtrVector().push_back(
379
380 new OpPPMap(pip->getPostProcMesh(), pip->getMapGaussPts(), {},
381 {{"U", u_ptr}},
382 {{"GRAD", common_ptr->matGradPtr},
383 {"FIRST_PIOLA", common_ptr->getMatFirstPiolaStress()}},
384 {}));
386 };
387
388 auto push_post_proc_bdy = [&](auto &pip_bdy) {
389 if (post_proc_skin == PETSC_FALSE)
390 return boost::shared_ptr<PostProcEleBdy>();
391 auto simple = mField.getInterface<Simple>();
392 auto domain_fe_name = simple->getDomainFEName();
393 auto u_ptr = boost::make_shared<MatrixDouble>();
394 pip_bdy->getOpPtrVector().push_back(
395 new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
396 auto op_loop_side =
397 new OpLoopSide<SideEle>(mField, domain_fe_name, SPACE_DIM);
398 auto common_ptr = post_proc_ele_domain(op_loop_side->getOpPtrVector());
399 pip_bdy->getOpPtrVector().push_back(op_loop_side);
400 CHKERR post_proc_map(pip_bdy, u_ptr, common_ptr);
401 return pip_bdy;
402 };
403
404 auto push_post_proc_domain = [&](auto &pip_domain) {
405 if (post_proc_vol == PETSC_FALSE)
406 return boost::shared_ptr<PostProcEleDomain>();
407
408 auto u_ptr = boost::make_shared<MatrixDouble>();
409 pip_domain->getOpPtrVector().push_back(
410 new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
411 auto common_ptr = post_proc_ele_domain(pip_domain->getOpPtrVector());
412
413 CHKERR post_proc_map(pip_domain, u_ptr, common_ptr);
414
415 return pip_domain;
416 };
417
418 auto post_proc_fe_domain = boost::make_shared<PostProcEleDomain>(mField);
419 auto post_proc_fe_bdy = boost::make_shared<PostProcEleBdy>(mField);
420
421 return std::make_pair(push_post_proc_domain(post_proc_fe_domain),
422 push_post_proc_bdy(post_proc_fe_bdy));
423 };
424
425 auto post_proc_pair = create_post_proc_fe();
426 postProcDomainFe = post_proc_pair.first;
427 postProcBdyFe = post_proc_pair.second;
428
430}
431//! [Push operators to pipeline]
432
435
436 // get operators tester
437 auto simple = mField.getInterface<Simple>();
438 auto opt = mField.getInterface<OperatorsTester>(); // get interface to
439 // OperatorsTester
440 auto pip = mField.getInterface<PipelineManager>(); // get interface to
441 // pipeline manager
442
443 constexpr double eps = 1e-9;
444
445 auto x = opt->setRandomFields(simple->getDM(), {
446
447 {"U", {-1e-6, 1e-6}}
448
449 });
450
451 auto dot_x = opt->setRandomFields(simple->getDM(), {
452
453 {"U", {-1, 1}}
454
455 });
456
457 auto diff_x = opt->setRandomFields(simple->getDM(), {
458
459 {"U", {-1, 1}}
460
461 });
462
463 auto test_domain_ops = [&](auto fe_name, auto lhs_pipeline,
464 auto rhs_pipeline) {
466
467 auto diff_res = opt->checkCentralFiniteDifference(
468 simple->getDM(), fe_name, rhs_pipeline, lhs_pipeline, x, dot_x,
469 SmartPetscObj<Vec>(), diff_x, 0, 0.5, eps);
470
471 // Calculate norm of difference between directional derivative calculated
472 // from finite difference, and tangent matrix.
473 double fnorm;
474 CHKERR VecNorm(diff_res, NORM_2, &fnorm);
475 MOFEM_LOG_C("EXAMPLE", Sev::inform,
476 "Test consistency of tangent matrix %3.4e", fnorm);
477
478 constexpr double err = 1e-5;
479 if (fnorm > err)
480 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
481 "Norm of directional derivative too large err = %3.4e", fnorm);
482
484 };
485
486 MOFEM_LOG("EXAMPLE", Sev::inform)
487 << "Test operators with finite difference of directional "
488 "derivative";
489 CHKERR test_domain_ops(simple->getDomainFEName(), pip->getDomainLhsFE(),
490 pip->getDomainRhsFE());
491
493}
494
497#ifdef WITH_ADOL_C
501 modelType>(
504 }
505 if (!physicalHuHuPtr) {
506 const int tag_hu_hu = AdolCOps::AdolCTagsRegistry::setTagName("HuHU");
508 AdolCOps::createAdolCPhysicalEquationsPtr<AdolCOps::HUHU, modelType>(
509 AdolCOps::createADolCDataPtr(), tag_hu_hu);
510 }
511
512 auto *simple = mField.getInterface<Simple>();
513 auto *pipeline_mng = mField.getInterface<PipelineManager>();
514
515 CHKERR physicalEquationsPtr->getOptions(&mField);
516 CHKERR physicalEquationsPtr->recordTape();
517 CHKERR physicalHuHuPtr->getOptions(&mField);
518 CHKERR physicalHuHuPtr->recordTape();
519
520 auto add_domain_ops_lhs = [&](auto &pip) {
523 "GEOMETRY");
525 template opLhsFactory<PETSC, GAUSS, DomainEleOp>(mField, pip, "U",
528 template opLhsFactory<PETSC, GAUSS, DomainEleOp>(
529 mField, pip, simple->getDomainFEName(), "U", physicalHuHuPtr);
531 };
532
533 auto add_domain_ops_rhs = [&](auto &pip) {
536 "GEOMETRY");
538 template opRhsFactory<PETSC, GAUSS, DomainEleOp>(mField, pip, "U",
541 template opRhsFactory<PETSC, GAUSS, DomainEleOp>(
542 mField, pip, simple->getDomainFEName(), "U", physicalHuHuPtr);
544 };
545
546 CHKERR add_domain_ops_lhs(pipeline_mng->getOpDomainLhsPipeline());
547 CHKERR add_domain_ops_rhs(pipeline_mng->getOpDomainRhsPipeline());
548
549 CHKERR add_domain_ops_rhs(reactionFe->getOpPtrVector());
550 reactionFe->postProcessHook =
551 EssentialPreProcReaction<DisplacementCubitBcData>(mField, reactionFe);
552
553 PetscBool post_proc_vol;
554 PetscBool post_proc_skin;
555
556 if constexpr (SPACE_DIM == 2) {
557 post_proc_vol = PETSC_TRUE;
558 post_proc_skin = PETSC_FALSE;
559 } else {
560 post_proc_vol = PETSC_FALSE;
561 post_proc_skin = PETSC_TRUE;
562 }
563
564 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-post_proc_vol",
565 &post_proc_vol, PETSC_NULLPTR);
566 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-post_proc_skin",
567 &post_proc_skin, PETSC_NULLPTR);
568
569 // Setup postprocessing
570 auto create_post_proc_fe = [&]() {
571 auto post_proc_ele_domain = [&](auto &pip_domain) {
574 "GEOMETRY");
578 };
579
580 auto post_proc_map = [&](auto &pip, auto u_ptr) {
582 using OpPPMap = OpPostProcMapInMoab<SPACE_DIM, SPACE_DIM>;
583 auto grad_ptr =
584 physicalEquationsPtr->adolcDataPtr->getCommonDataPtr("grad");
585 auto first_piola_ptr =
586 physicalEquationsPtr->adolcDataPtr->getCommonDataPtr("P");
587 pip->getOpPtrVector().push_back(new OpPPMap(
588 pip->getPostProcMesh(), pip->getMapGaussPts(), {}, {{"U", u_ptr}},
589 {{"GRAD", grad_ptr}, {"FIRST_PIOLA", first_piola_ptr}}, {}));
591 };
592
593 auto push_post_proc_bdy = [&](auto &pip_bdy) {
594 if (post_proc_skin == PETSC_FALSE)
595 return boost::shared_ptr<PostProcEleBdy>();
596 auto simple = mField.getInterface<Simple>();
597 auto domain_fe_name = simple->getDomainFEName();
598 auto u_ptr = boost::make_shared<MatrixDouble>();
599 pip_bdy->getOpPtrVector().push_back(
600 new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
601 auto op_loop_side =
602 new OpLoopSide<SideEle>(mField, domain_fe_name, SPACE_DIM);
603 CHKERR post_proc_ele_domain(op_loop_side->getOpPtrVector());
604 pip_bdy->getOpPtrVector().push_back(op_loop_side);
605 CHKERR post_proc_map(pip_bdy, u_ptr);
606 return pip_bdy;
607 };
608
609 auto push_post_proc_domain = [&](auto &pip_domain) {
610 if (post_proc_vol == PETSC_FALSE)
611 return boost::shared_ptr<PostProcEleDomain>();
612 auto u_ptr = boost::make_shared<MatrixDouble>();
613 pip_domain->getOpPtrVector().push_back(
614 new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
615 CHKERR post_proc_ele_domain(pip_domain->getOpPtrVector());
616 CHKERR post_proc_map(pip_domain, u_ptr);
617
618 return pip_domain;
619 };
620
621 auto post_proc_fe_domain = boost::make_shared<PostProcEleDomain>(mField);
622 auto post_proc_fe_bdy = boost::make_shared<PostProcEleBdy>(mField);
623
624 return std::make_pair(push_post_proc_domain(post_proc_fe_domain),
625 push_post_proc_bdy(post_proc_fe_bdy));
626 };
627
628 auto post_proc_pair = create_post_proc_fe();
629 postProcDomainFe = post_proc_pair.first;
630 postProcBdyFe = post_proc_pair.second;
631
632#else
633 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
634 "ADOL-C support is not enabled. Please reconfigure with "
635 "-DWITH_ADOL_C=ON and recompile to use AdolC material model.");
636
637#endif // WITH_ADOL
639}
640
641//! [TS Solve]
644 auto *simple = mField.getInterface<Simple>();
645 auto *pipeline_mng = mField.getInterface<PipelineManager>();
646
647 auto dm = simple->getDM();
648 auto ts = pipeline_mng->createTSIM();
649
650 auto add_extra_finite_elements_to_solver_pipelines = [&]() {
652
653 auto pre_proc_ptr = boost::make_shared<FEMethod>();
654 auto post_proc_rhs_ptr = boost::make_shared<FEMethod>();
655 auto post_proc_lhs_ptr = boost::make_shared<FEMethod>();
656
657 auto time_scale = boost::make_shared<ExampleTimeScale>();
658
659 auto get_bc_hook_rhs = [this, pre_proc_ptr, time_scale]() {
661 CHKERR EssentialPreProc<DisplacementCubitBcData>(mField, pre_proc_ptr,
662 {time_scale}, false)();
664 };
665
666 pre_proc_ptr->preProcessHook = get_bc_hook_rhs;
667
668 auto get_post_proc_hook_rhs = [this, post_proc_rhs_ptr]() {
670 CHKERR EssentialPreProcReaction<DisplacementCubitBcData>(
671 mField, post_proc_rhs_ptr, nullptr, Sev::verbose)();
672 CHKERR EssentialPostProcRhs<DisplacementCubitBcData>(
673 mField, post_proc_rhs_ptr, 1.)();
675 };
676 auto get_post_proc_hook_lhs = [this, post_proc_lhs_ptr]() {
678 CHKERR EssentialPostProcLhs<DisplacementCubitBcData>(
679 mField, post_proc_lhs_ptr, 1.)();
681 };
682 post_proc_rhs_ptr->postProcessHook = get_post_proc_hook_rhs;
683 post_proc_lhs_ptr->postProcessHook = get_post_proc_hook_lhs;
684
685 // This is low level pushing finite elements (pipelines) to solver
686 auto ts_ctx_ptr = getDMTsCtx(simple->getDM());
687 ts_ctx_ptr->getPreProcessIFunction().push_front(pre_proc_ptr);
688 ts_ctx_ptr->getPreProcessIJacobian().push_front(pre_proc_ptr);
689 ts_ctx_ptr->getPostProcessIFunction().push_back(post_proc_rhs_ptr);
690 ts_ctx_ptr->getPostProcessIJacobian().push_back(post_proc_lhs_ptr);
692 };
693
694 // Add extra finite elements to SNES solver pipelines to resolve essential
695 // boundary conditions
696 CHKERR add_extra_finite_elements_to_solver_pipelines();
697
698 auto create_monitor_fe = [dm](auto &&post_proc_fe, auto &&reactionFe) {
699 return boost::make_shared<Monitor>(dm.get(), post_proc_fe, reactionFe);
700 };
701
702 // Set monitor which postprocessing results and saves them to the hard drive
703 boost::shared_ptr<FEMethod> null_fe;
705 CHKERR postProcDomainFe->setTagsToTransfer(
706 std::vector<Tag>(listTagsToTransfer));
707 if (postProcBdyFe)
708 CHKERR postProcBdyFe->setTagsToTransfer(
709 std::vector<Tag>(listTagsToTransfer));
710 auto monitor_ptr = create_monitor_fe(
711 std::make_pair(postProcDomainFe, postProcBdyFe), reactionFe);
712 CHKERR DMMoFEMTSSetMonitor(dm, ts, simple->getDomainFEName(), null_fe,
713 null_fe, monitor_ptr);
714
715 // Set time solver
716 double ftime = 1;
717 CHKERR TSSetMaxTime(ts, ftime);
718 CHKERR TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);
719
720 auto B = createDMMatrix(dm);
721 CHKERR TSSetI2Jacobian(ts, B, B, PETSC_NULLPTR, PETSC_NULLPTR);
722 auto D = createDMVector(simple->getDM());
723 CHKERR TSSetSolution(ts, D);
724 CHKERR TSSetFromOptions(ts);
725
726 CHKERR TSSolve(ts, NULL);
727 CHKERR TSGetTime(ts, &ftime);
728
729 PetscInt steps, snesfails, rejects, nonlinits, linits;
730 CHKERR TSGetStepNumber(ts, &steps);
731 CHKERR TSGetSNESFailures(ts, &snesfails);
732 CHKERR TSGetStepRejections(ts, &rejects);
733 CHKERR TSGetSNESIterations(ts, &nonlinits);
734 CHKERR TSGetKSPIterations(ts, &linits);
735 MOFEM_LOG_C("EXAMPLE", Sev::inform,
736 "steps %d (%d rejected, %d SNES fails), ftime %g, nonlinits "
737 "%d, linits %d",
738 steps, rejects, snesfails, ftime, nonlinits, linits);
739
741}
742//! [TS Solve]
743
744//! [Getting norms]
747
748 auto simple = mField.getInterface<Simple>();
749 auto dm = simple->getDM();
750
751 auto T = createDMVector(simple->getDM());
752 CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
753 SCATTER_FORWARD);
754 double nrm2;
755 CHKERR VecNorm(T, NORM_2, &nrm2);
756 MOFEM_LOG("EXAMPLE", Sev::inform) << "Solution norm " << nrm2;
757
758 auto post_proc_norm_fe = boost::make_shared<DomainEle>(mField);
759
760 auto post_proc_norm_rule_hook = [](int, int, int p) -> int { return 2 * p; };
761 post_proc_norm_fe->getRuleHook = post_proc_norm_rule_hook;
762
764 post_proc_norm_fe->getOpPtrVector(), {H1}, "GEOMETRY");
765
766 enum NORMS { U_NORM_L2 = 0, PIOLA_NORM, LAST_NORM };
767 auto norms_vec =
768 createVectorMPI(mField.get_comm(),
769 (mField.get_comm_rank() == 0) ? LAST_NORM : 0, LAST_NORM);
770 CHKERR VecZeroEntries(norms_vec);
771
772 auto u_ptr = boost::make_shared<MatrixDouble>();
773 post_proc_norm_fe->getOpPtrVector().push_back(
774 new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
775
776 post_proc_norm_fe->getOpPtrVector().push_back(
777 new OpCalcNormL2Tensor1<SPACE_DIM>(u_ptr, norms_vec, U_NORM_L2));
778
779 if (useAdolcMaterial == PETSC_TRUE) {
780#ifdef WITH_ADOL_C
782 opPostProcFactory(mField, post_proc_norm_fe->getOpPtrVector(), "U",
784 auto m_P = physicalEquationsPtr->adolcDataPtr->getCommonDataPtr("P");
785 post_proc_norm_fe->getOpPtrVector().push_back(
786 new OpCalcNormL2Tensor2<SPACE_DIM, SPACE_DIM>(m_P, norms_vec,
787 PIOLA_NORM));
788#else
789 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
790 "ADOL-C support is not enabled. Please reconfigure with "
791 "-DWITH_ADOL_C=ON and recompile to use AdolC material model.");
792#endif
793 } else {
794 auto common_ptr = commonDataFactory<SPACE_DIM, GAUSS, DomainEleOp>(
795 mField, post_proc_norm_fe->getOpPtrVector(), "U", "MAT_ELASTIC",
796 Sev::inform);
797 post_proc_norm_fe->getOpPtrVector().push_back(
798 new OpCalcNormL2Tensor2<SPACE_DIM, SPACE_DIM>(
799 common_ptr->getMatFirstPiolaStress(), norms_vec, PIOLA_NORM));
800 }
801
802 CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(),
803 post_proc_norm_fe);
804
805 CHKERR VecAssemblyBegin(norms_vec);
806 CHKERR VecAssemblyEnd(norms_vec);
807
808 MOFEM_LOG_CHANNEL("SELF"); // Clear channel from old tags
809 if (mField.get_comm_rank() == 0) {
810 const double *norms;
811 CHKERR VecGetArrayRead(norms_vec, &norms);
812 MOFEM_TAG_AND_LOG("SELF", Sev::inform, "example")
813 << "norm_u: " << std::scientific << std::sqrt(norms[U_NORM_L2]);
814 MOFEM_TAG_AND_LOG("SELF", Sev::inform, "example")
815 << "norm_piola: " << std::scientific << std::sqrt(norms[PIOLA_NORM]);
816 CHKERR VecRestoreArrayRead(norms_vec, &norms);
817 }
818
820}
821//! [Getting norms]
822
823//! [Postprocessing results]
828//! [Postprocessing results]
829
830//! [Check]
833 PetscInt test_nb = 0;
834 PetscBool test_flg = PETSC_FALSE;
835 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-test", &test_nb, &test_flg);
836
837 if (test_flg) {
838 auto simple = mField.getInterface<Simple>();
839 auto T = createDMVector(simple->getDM());
840 CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
841 SCATTER_FORWARD);
842 double nrm2;
843 CHKERR VecNorm(T, NORM_2, &nrm2);
844 MOFEM_LOG("EXAMPLE", Sev::verbose) << "Regression norm " << nrm2;
845 double regression_value = 0;
846 switch (test_nb) {
847 case 1:
848 regression_value = 1.02789;
849 break;
850 case 2:
851 regression_value = 1.8841e+00;
852 break;
853 case 3:
854 regression_value = 1.8841e+00;
855 break;
856 default:
857 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Wrong test number.");
858 break;
859 }
860 if (fabs(nrm2 - regression_value) > 1e-2)
861 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
862 "Regression test field; wrong norm value. %6.4e != %6.4e", nrm2,
863 regression_value);
864 }
866}
867//! [Check]
#define MOFEM_LOG_C(channel, severity, format,...)
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
static const double eps
constexpr int SPACE_DIM
FieldApproximationBase
approximation base
Definition definitions.h:58
@ LASTBASE
Definition definitions.h:69
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ H1
continuous field
Definition definitions.h:85
#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.
@ 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 ...
constexpr int order
auto integration_rule
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
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.
double D
boost::shared_ptr< ADolCData > createADolCDataPtr()
Definition AdolOps.cpp:245
boost::shared_ptr< PhysicalEquations > createAdolCPhysicalEquationsPtr(boost::shared_ptr< ADolCData > adolc_data_ptr, int tag)
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
static constexpr int approx_order
FTensor::Index< 'm', 3 > m
static auto setTagName(std::string name, int tag=-1)
Definition AdolCOps.hpp:16
virtual moab::Interface & get_moab()=0
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
Deprecated interface functions.
Post post-proc data at points from hash maps.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
NonlinearElasticExample(MoFEM::Interface &m_field)
MoFEMErrorCode runProblem()
[Run problem]
boost::shared_ptr< PostProcEleBdy > postProcBdyFe
boost::shared_ptr< AdolCOps::PhysicalEquations > physicalEquationsPtr
boost::shared_ptr< DomainEle > reactionFe
std::vector< Tag > listTagsToTransfer
list of tags to transfer to postprocessor
MoFEMErrorCode setupProblem()
[Read mesh]
MoFEMErrorCode readMesh()
[Run problem]
MoFEMErrorCode assembleSystemHencky()
[Boundary condition]
MoFEMErrorCode boundaryCondition()
[Set up problem]
MoFEMErrorCode opTest()
[Push operators to pipeline]
boost::shared_ptr< PostProcEleDomain > postProcDomainFe
boost::shared_ptr< AdolCOps::PhysicalEquations > physicalHuHuPtr
MoFEMErrorCode TsSolve()
[TS Solve]
MoFEMErrorCode outputResults()
[Getting norms]
static auto get_skin(MoFEM::Interface &m_field, Range body_ents)
static auto filter_true_skin(MoFEM::Interface &m_field, Range &&skin)
MoFEMErrorCode checkResults()
[Postprocessing results]
MoFEMErrorCode gettingNorms()
[TS Solve]