v0.15.0
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ElasticExample Struct Reference

#include "tutorials/vec-0/src/ElasticExample.hpp"

Inheritance diagram for ElasticExample:
[legend]
Collaboration diagram for ElasticExample:
[legend]

Public Member Functions

 ElasticExample (MoFEM::Interface &m_field)
 
MoFEMErrorCode runProblem ()
 [Run problem]
 

Protected Member Functions

MoFEMErrorCode readMesh ()
 [Run problem]
 
MoFEMErrorCode setupProblem ()
 [Read mesh]
 
MoFEMErrorCode boundaryCondition ()
 [Set up problem]
 
MoFEMErrorCode assembleSystem ()
 [Boundary condition]
 
MoFEMErrorCode solveSystem ()
 [Solve]
 
MoFEMErrorCode outputResults ()
 [Solve]
 
MoFEMErrorCode checkResults ()
 [Postprocess results]
 
virtual MoFEMErrorCode kspSetUpAndSolve (SmartPetscObj< KSP > solver)
 [Push operators to pipeline]
 

Protected Attributes

MoFEM::InterfacemField
 
boost::shared_ptr< MatrixDouble > vectorFieldPtr = nullptr
 

Detailed Description

Examples
mofem/tutorials/vec-0/elastic.cpp, and mofem/tutorials/vec-10/schur_elastic.cpp.

Definition at line 17 of file ElasticExample.hpp.

Constructor & Destructor Documentation

◆ ElasticExample()

ElasticExample::ElasticExample ( MoFEM::Interface m_field)
inline

Definition at line 19 of file ElasticExample.hpp.

19: mField(m_field) {}
MoFEM::Interface & mField

Member Function Documentation

◆ assembleSystem()

MoFEMErrorCode ElasticExample::assembleSystem ( )
protected

[Boundary condition]

[Push operators to pipeline]

[Integration rule]

[Integration rule]

[Push domain stiffness matrix to pipeline]

[Push domain stiffness matrix to pipeline]

[Push Body forces]

[Push Body forces]

[Push natural boundary conditions]

[Push natural boundary conditions]

Examples
mofem/tutorials/vec-10/schur_elastic.cpp.

Definition at line 137 of file ElasticExample.hpp.

137 {
139 auto pip = mField.getInterface<PipelineManager>();
140
141 //! [Integration rule]
142 auto integration_rule = [](int, int, int approx_order) {
143 return 2 * approx_order + 1;
144 };
145
146 auto integration_rule_bc = [](int, int, int approx_order) {
147 return 2 * approx_order + 1;
148 };
149
150 CHKERR pip->setDomainRhsIntegrationRule(integration_rule);
151 CHKERR pip->setDomainLhsIntegrationRule(integration_rule);
152 CHKERR pip->setBoundaryRhsIntegrationRule(integration_rule_bc);
153 CHKERR pip->setBoundaryLhsIntegrationRule(integration_rule_bc);
154 //! [Integration rule]
155
157 pip->getOpDomainLhsPipeline(), {H1}, "GEOMETRY");
159 pip->getOpDomainRhsPipeline(), {H1}, "GEOMETRY");
161 pip->getOpBoundaryRhsPipeline(), {NOSPACE}, "GEOMETRY");
163 pip->getOpBoundaryLhsPipeline(), {NOSPACE}, "GEOMETRY");
164
165 //! [Push domain stiffness matrix to pipeline]
166 // Add LHS operator for elasticity (stiffness matrix)
167 CHKERR HookeOps::opFactoryDomainLhs<SPACE_DIM, A, I, DomainEleOp>(
168 mField, pip->getOpDomainLhsPipeline(), "U", "MAT_ELASTIC", Sev::verbose);
169 //! [Push domain stiffness matrix to pipeline]
170
171 // Add RHS operator for internal forces
172 CHKERR HookeOps::opFactoryDomainRhs<SPACE_DIM, A, I, DomainEleOp>(
173 mField, pip->getOpDomainRhsPipeline(), "U", "MAT_ELASTIC", Sev::verbose);
174
175 //! [Push Body forces]
176 CHKERR DomainRhsBCs::AddFluxToPipeline<OpDomainRhsBCs>::add(
177 pip->getOpDomainRhsPipeline(), mField, "U", Sev::inform);
178 //! [Push Body forces]
179
180 //! [Push natural boundary conditions]
181 // Add force boundary condition
182 CHKERR BoundaryRhsBCs::AddFluxToPipeline<OpBoundaryRhsBCs>::add(
183 pip->getOpBoundaryRhsPipeline(), mField, "U", -1, Sev::inform);
184 // Add case for mix type of BCs
185 CHKERR BoundaryLhsBCs::AddFluxToPipeline<OpBoundaryLhsBCs>::add(
186 pip->getOpBoundaryLhsPipeline(), mField, "U", Sev::verbose);
187 //! [Push natural boundary conditions]
189}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
auto integration_rule
static constexpr int approx_order
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

◆ boundaryCondition()

MoFEMErrorCode ElasticExample::boundaryCondition ( )
protected

[Set up problem]

[Boundary condition]

Examples
mofem/tutorials/vec-10/schur_elastic.cpp.

Definition at line 116 of file ElasticExample.hpp.

116 {
118 auto simple = mField.getInterface<Simple>();
119 auto bc_mng = mField.getInterface<BcManager>();
120
121 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_X",
122 "U", 0, 0);
123 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_Y",
124 "U", 1, 1);
125 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_Z",
126 "U", 2, 2);
127 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(),
128 "REMOVE_ALL", "U", 0, 3);
129 CHKERR bc_mng->pushMarkDOFsOnEntities<DisplacementCubitBcData>(
130 simple->getProblemName(), "U");
131
133}
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69

◆ checkResults()

MoFEMErrorCode ElasticExample::checkResults ( )
protected

[Postprocess results]

[Check]

Examples
mofem/tutorials/vec-10/schur_elastic.cpp.

Definition at line 503 of file ElasticExample.hpp.

503 {
504 MOFEM_LOG_CHANNEL("WORLD");
505 auto simple = mField.getInterface<Simple>();
506 auto pip = mField.getInterface<PipelineManager>();
508 pip->getDomainRhsFE().reset();
509 pip->getDomainLhsFE().reset();
510 pip->getBoundaryRhsFE().reset();
511 pip->getBoundaryLhsFE().reset();
512
513 auto integration_rule = [](int, int, int p_data) { return 2 * p_data + 1; };
514 CHKERR pip->setDomainRhsIntegrationRule(integration_rule);
515 CHKERR pip->setBoundaryRhsIntegrationRule(integration_rule);
516
518 pip->getOpDomainRhsPipeline(), {H1}, "GEOMETRY");
520 pip->getOpBoundaryRhsPipeline(), {}, "GEOMETRY");
521
522 // Add RHS operators for internal forces
523 CHKERR HookeOps::opFactoryDomainRhs<SPACE_DIM, A, I, DomainEleOp>(
524 mField, pip->getOpDomainRhsPipeline(), "U", "MAT_ELASTIC", Sev::verbose);
525
526 CHKERR DomainRhsBCs::AddFluxToPipeline<OpDomainRhsBCs>::add(
527 pip->getOpDomainRhsPipeline(), mField, "U", Sev::verbose);
528
529 CHKERR BoundaryRhsBCs::AddFluxToPipeline<OpBoundaryRhsBCs>::add(
530 pip->getOpBoundaryRhsPipeline(), mField, "U", -1, Sev::verbose);
531
532 auto dm = simple->getDM();
533 auto res = createDMVector(dm);
534 CHKERR VecSetDM(res, PETSC_NULLPTR);
535
536 pip->getDomainRhsFE()->f = res;
537 pip->getBoundaryRhsFE()->f = res;
538
539 CHKERR VecZeroEntries(res);
540
541 CHKERR pip->loopFiniteElements();
542 // CHKERR mField.getInterface<FieldBlas>()->fieldScale(-1, "U");
543
544 CHKERR VecGhostUpdateBegin(res, ADD_VALUES, SCATTER_REVERSE);
545 CHKERR VecGhostUpdateEnd(res, ADD_VALUES, SCATTER_REVERSE);
546 CHKERR VecAssemblyBegin(res);
547 CHKERR VecAssemblyEnd(res);
548
549 auto zero_residual_at_constrains = [&]() {
551 auto fe_post_proc_ptr = boost::make_shared<FEMethod>();
552 auto get_post_proc_hook_rhs = [&]() {
554 CHKERR EssentialPreProcReaction<DisplacementCubitBcData>(
555 mField, fe_post_proc_ptr, res)();
556 CHKERR EssentialPostProcRhs<DisplacementCubitBcData>(
557 mField, fe_post_proc_ptr, 0, res)();
559 };
560 fe_post_proc_ptr->postProcessHook = get_post_proc_hook_rhs;
561 CHKERR DMoFEMPostProcessFiniteElements(dm, fe_post_proc_ptr.get());
563 };
564
565 CHKERR zero_residual_at_constrains();
566
567 double nrm2;
568 CHKERR VecNorm(res, NORM_2, &nrm2);
569 MOFEM_LOG_CHANNEL("WORLD");
570 MOFEM_LOG_C("WORLD", Sev::inform, "residual = %3.4e\n", nrm2);
571
572 int test = 0;
573 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-test", &test, PETSC_NULLPTR);
574 if (test > 0) {
575 auto post_proc_residual = [&](auto dm, auto f_res, auto out_name) {
577 auto post_proc_fe =
578 boost::make_shared<PostProcBrokenMeshInMoab<DomainEle>>(mField);
579 using OpPPMap = OpPostProcMapInMoab<SPACE_DIM, SPACE_DIM>;
580 auto u_vec = boost::make_shared<MatrixDouble>();
581 post_proc_fe->getOpPtrVector().push_back(
582 new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_vec, f_res));
583 post_proc_fe->getOpPtrVector().push_back(
584
585 new OpPPMap(
586
587 post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
588
589 {},
590
591 {{"RES", u_vec}},
592
593 {}, {})
594
595 );
596
597 CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(),
598 post_proc_fe);
599 post_proc_fe->writeFile(out_name);
601 };
602
603 CHKERR post_proc_residual(simple->getDM(), res, "res.h5m");
604
605 constexpr double eps = 1e-8;
606 if (nrm2 > eps)
607 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
608 "Residual is not zero");
609 }
610 if (test == 2) {
611 if (!vectorFieldPtr || vectorFieldPtr->size1() == 0) {
612 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
613 "atom test %d failed: Field Evaluator did not provide result",
614 test);
615 }
616 auto t_disp = getFTensor1FromMat<SPACE_DIM>(*vectorFieldPtr);
617 double Ux_ref = 0.46;
618 double Uy_ref = -0.03;
619 constexpr double eps = 1e-8;
620 if (fabs(t_disp(0) - Ux_ref) > eps || fabs(t_disp(1) - Uy_ref) > eps) {
621 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
622 "atom test %d failed: Ux_ref = %3.6e, computed = %3.6e, Uy_ref "
623 "= %3.6e, computed = %3.6e",
624 test, Ux_ref, t_disp(0), Uy_ref, t_disp(1));
625 }
626 }
628}
#define MOFEM_LOG_C(channel, severity, format,...)
static const double eps
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition DMMoFEM.cpp:546
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
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
boost::shared_ptr< MatrixDouble > vectorFieldPtr

◆ kspSetUpAndSolve()

MoFEMErrorCode ElasticExample::kspSetUpAndSolve ( SmartPetscObj< KSP >  solver)
protectedvirtual

[Push operators to pipeline]

Reimplemented in ElasticSchurExample.

Definition at line 193 of file ElasticExample.hpp.

193 {
195
196 MOFEM_LOG_CHANNEL("TIMER");
197 MOFEM_LOG_TAG("TIMER", "timer");
198
199 BOOST_LOG_SCOPED_THREAD_ATTR("Timeline", attrs::timer());
200 MOFEM_LOG("TIMER", Sev::inform) << "KSPSetUp";
201 CHKERR KSPSetUp(solver);
202 MOFEM_LOG("TIMER", Sev::inform) << "KSPSetUp <= Done";
203
204 DM dm;
205 CHKERR KSPGetDM(solver, &dm);
206 auto D = createDMVector(dm);
207 auto F = vectorDuplicate(D);
208
209 MOFEM_LOG("TIMER", Sev::inform) << "KSPSolve";
210 CHKERR KSPSolve(solver, F, D);
211 MOFEM_LOG("TIMER", Sev::inform) << "KSPSolve <= Done";
212
213 CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
214 CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
215 CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
216
218};
@ F
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
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
double D
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.

◆ outputResults()

MoFEMErrorCode ElasticExample::outputResults ( )
protected

[Solve]

[Postprocess results]

[Postprocess clean]

[Postprocess clean]

[Postprocess initialise]

[Postprocess initialise]

Examples
mofem/tutorials/vec-10/schur_elastic.cpp.

Definition at line 331 of file ElasticExample.hpp.

331 {
333 auto simple = mField.getInterface<Simple>();
334 auto pip = mField.getInterface<PipelineManager>();
335 auto det_ptr = boost::make_shared<VectorDouble>();
336 auto jac_ptr = boost::make_shared<MatrixDouble>();
337 auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
338 //! [Postprocess clean]
339 pip->getDomainRhsFE().reset();
340 pip->getDomainLhsFE().reset();
341 pip->getBoundaryRhsFE().reset();
342 pip->getBoundaryLhsFE().reset();
343 //! [Postprocess clean]
344
345 //! [Postprocess initialise]
346 auto post_proc_mesh = boost::make_shared<moab::Core>();
347 auto post_proc_begin = boost::make_shared<PostProcBrokenMeshInMoabBaseBegin>(
348 mField, post_proc_mesh);
349 auto post_proc_end = boost::make_shared<PostProcBrokenMeshInMoabBaseEnd>(
350 mField, post_proc_mesh);
351 //! [Postprocess initialise]
352
353 // lamdaa function to calculate dispalcement, strain and stress
354 auto calculate_stress_ops = [&](auto &pip) {
355 // Add H1 geometry ops
357
358 // Use HookeOps commonDataFactory to get matStrainPtr and matCauchyStressPtr
359 auto common_ptr = HookeOps::commonDataFactory<SPACE_DIM, GAUSS, DomainEle>(
360 mField, pip, "U", "MAT_ELASTIC", Sev::verbose);
361
362 // Store U and GEOMETRY values if needed
363 auto u_ptr = boost::make_shared<MatrixDouble>();
364 pip.push_back(new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
365 auto x_ptr = boost::make_shared<MatrixDouble>();
366 pip.push_back(
367 new OpCalculateVectorFieldValues<SPACE_DIM>("GEOMETRY", x_ptr));
368
369 // Return what you need: displacements, coordinates, strain, stress
370 return boost::make_tuple(u_ptr, x_ptr, common_ptr->getMatStrain(),
371 common_ptr->getMatCauchyStress());
372 };
373
374 auto get_tag_id_on_pmesh = [&](bool post_proc_skin) {
375 int def_val_int = 0;
376 Tag tag_mat;
377 CHKERR mField.get_moab().tag_get_handle(
378 "MAT_ELASTIC", 1, MB_TYPE_INTEGER, tag_mat,
379 MB_TAG_CREAT | MB_TAG_SPARSE, &def_val_int);
380 auto meshset_vec_ptr =
381 mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(
382 std::regex((boost::format("%s(.*)") % "MAT_ELASTIC").str()));
383
384 Range skin_ents;
385 std::unique_ptr<Skinner> skin_ptr;
386 if (post_proc_skin) {
387 skin_ptr = std::make_unique<Skinner>(&mField.get_moab());
388 auto boundary_meshset = simple->getBoundaryMeshSet();
389 CHKERR mField.get_moab().get_entities_by_handle(boundary_meshset,
390 skin_ents, true);
391 }
392
393 for (auto m : meshset_vec_ptr) {
394 Range ents_3d;
395 CHKERR mField.get_moab().get_entities_by_handle(m->getMeshset(), ents_3d,
396 true);
397 int const id = m->getMeshsetId();
398 ents_3d = ents_3d.subset_by_dimension(SPACE_DIM);
399 if (post_proc_skin) {
400 Range skin_faces;
401 CHKERR skin_ptr->find_skin(0, ents_3d, false, skin_faces);
402 ents_3d = intersect(skin_ents, skin_faces);
403 }
404 CHKERR mField.get_moab().tag_clear_data(tag_mat, ents_3d, &id);
405 }
406
407 return tag_mat;
408 };
409
410 auto post_proc_domain = [&](auto post_proc_mesh) {
411 auto post_proc_fe =
412 boost::make_shared<PostProcEleDomain>(mField, post_proc_mesh);
413 using OpPPMap = OpPostProcMapInMoab<SPACE_DIM, SPACE_DIM>;
414
415 auto [u_ptr, x_ptr, mat_strain_ptr, mat_stress_ptr] =
416 calculate_stress_ops(post_proc_fe->getOpPtrVector());
417
418 post_proc_fe->getOpPtrVector().push_back(
419
420 new OpPPMap(
421
422 post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
423
424 {},
425
426 {{"U", u_ptr}, {"GEOMETRY", x_ptr}},
427
428 {},
429
430 {{"STRAIN", mat_strain_ptr}, {"STRESS", mat_stress_ptr}}
431
432 )
433
434 );
435
436 post_proc_fe->setTagsToTransfer({get_tag_id_on_pmesh(false)});
437 return post_proc_fe;
438 };
439
440 auto post_proc_boundary = [&](auto post_proc_mesh) {
441 auto post_proc_fe =
442 boost::make_shared<PostProcEleBdy>(mField, post_proc_mesh);
444 post_proc_fe->getOpPtrVector(), {}, "GEOMETRY");
445 auto op_loop_side =
446 new OpLoopSide<SideEle>(mField, simple->getDomainFEName(), SPACE_DIM);
447 // push ops to side element, through op_loop_side operator
448 auto [u_ptr, x_ptr, mat_strain_ptr, mat_stress_ptr] =
449 calculate_stress_ops(op_loop_side->getOpPtrVector());
450 post_proc_fe->getOpPtrVector().push_back(op_loop_side);
451 auto mat_traction_ptr = boost::make_shared<MatrixDouble>();
452 post_proc_fe->getOpPtrVector().push_back(
453 new ElasticOps::OpCalculateTraction(mat_stress_ptr, mat_traction_ptr));
454
455 using OpPPMap = OpPostProcMapInMoab<SPACE_DIM, SPACE_DIM>;
456
457 post_proc_fe->getOpPtrVector().push_back(
458
459 new OpPPMap(
460
461 post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
462
463 {},
464
465 {{"U", u_ptr}, {"GEOMETRY", x_ptr}, {"T", mat_traction_ptr}},
466
467 {},
468
469 {{"STRAIN", mat_strain_ptr}, {"STRESS", mat_stress_ptr}}
470
471 )
472
473 );
474
475 post_proc_fe->setTagsToTransfer({get_tag_id_on_pmesh(true)});
476 return post_proc_fe;
477 };
478
479 PetscBool post_proc_skin_only = PETSC_FALSE;
480 if (SPACE_DIM == 3) {
481 post_proc_skin_only = PETSC_TRUE;
482 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-post_proc_skin_only",
483 &post_proc_skin_only, PETSC_NULLPTR);
484 }
485 if (post_proc_skin_only == PETSC_FALSE) {
486 pip->getDomainRhsFE() = post_proc_domain(post_proc_mesh);
487 } else {
488 pip->getBoundaryRhsFE() = post_proc_boundary(post_proc_mesh);
489 }
490
492 post_proc_begin->getFEMethod());
493 CHKERR pip->loopFiniteElements();
495 post_proc_end->getFEMethod());
496
497 CHKERR post_proc_end->writeFile("out_elastic.h5m");
499}
constexpr int SPACE_DIM
@ H1
continuous field
Definition definitions.h:85
PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition DMMoFEM.cpp:536
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
FTensor::Index< 'm', 3 > m
virtual moab::Interface & get_moab()=0

◆ readMesh()

MoFEMErrorCode ElasticExample::readMesh ( )
protected

[Run problem]

[Read mesh]

Examples
mofem/tutorials/vec-10/schur_elastic.cpp.

Definition at line 54 of file ElasticExample.hpp.

54 {
56 auto simple = mField.getInterface<Simple>();
57 CHKERR simple->getOptions();
58 CHKERR simple->loadFile();
59 // Add meshsets if config file provided
60 CHKERR mField.getInterface<MeshsetsManager>()->setMeshsetFromFile();
62}

◆ runProblem()

MoFEMErrorCode ElasticExample::runProblem ( )

[Run problem]

Definition at line 40 of file ElasticExample.hpp.

40 {
50}
MoFEMErrorCode solveSystem()
[Solve]
MoFEMErrorCode outputResults()
[Solve]
MoFEMErrorCode assembleSystem()
[Boundary condition]
MoFEMErrorCode setupProblem()
[Read mesh]
MoFEMErrorCode boundaryCondition()
[Set up problem]
MoFEMErrorCode checkResults()
[Postprocess results]
MoFEMErrorCode readMesh()
[Run problem]

◆ setupProblem()

MoFEMErrorCode ElasticExample::setupProblem ( )
protected

[Read mesh]

[Set up problem]

Examples
mofem/tutorials/vec-10/schur_elastic.cpp.

Definition at line 66 of file ElasticExample.hpp.

66 {
68 Simple *simple = mField.getInterface<Simple>();
69
70 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOPT };
71 const char *list_bases[LASBASETOPT] = {"ainsworth", "demkowicz"};
72 PetscInt choice_base_value = AINSWORTH;
73 CHKERR PetscOptionsGetEList(PETSC_NULLPTR, NULL, "-base", list_bases,
74 LASBASETOPT, &choice_base_value, PETSC_NULLPTR);
75
77 switch (choice_base_value) {
78 case AINSWORTH:
80 MOFEM_LOG("WORLD", Sev::inform)
81 << "Set AINSWORTH_LEGENDRE_BASE for displacements";
82 break;
83 case DEMKOWICZ:
85 MOFEM_LOG("WORLD", Sev::inform)
86 << "Set DEMKOWICZ_JACOBI_BASE for displacements";
87 break;
88 default:
89 base = LASTBASE;
90 break;
91 }
92
93 // Add field
94 CHKERR simple->addDomainField("U", H1, base, SPACE_DIM);
95 CHKERR simple->addBoundaryField("U", H1, base, SPACE_DIM);
96 int order = 3;
97 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-order", &order, PETSC_NULLPTR);
98
99 CHKERR simple->addDataField("GEOMETRY", H1, base, SPACE_DIM);
100
101 CHKERR simple->setFieldOrder("U", order);
102 CHKERR simple->setFieldOrder("GEOMETRY", 2);
103 CHKERR simple->setUp();
104
105 auto project_ho_geometry = [&]() {
106 Projection10NodeCoordsOnField ent_method(mField, "GEOMETRY");
107 return mField.loop_dofs("GEOMETRY", ent_method);
108 };
109 CHKERR project_ho_geometry();
110
112}
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
constexpr int order
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.
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)

◆ solveSystem()

MoFEMErrorCode ElasticExample::solveSystem ( )
protected

[Solve]

Examples
mofem/tutorials/vec-10/schur_elastic.cpp.

Definition at line 222 of file ElasticExample.hpp.

222 {
224
225 auto simple = mField.getInterface<Simple>();
226 auto pip = mField.getInterface<PipelineManager>();
227 auto solver = pip->createKSP();
228 CHKERR KSPSetFromOptions(solver);
229
230 auto set_essential_bc = [this]() {
232 // This is low level pushing finite elements (pipelines) to solver
233 auto simple = mField.getInterface<Simple>();
234 auto dm = simple->getDM();
235 auto ksp_ctx_ptr = getDMKspCtx(dm);
236
237 auto pre_proc_rhs = boost::make_shared<FEMethod>();
238 auto post_proc_rhs = boost::make_shared<FEMethod>();
239 auto post_proc_lhs = boost::make_shared<FEMethod>();
240
241 auto get_pre_proc_hook = [this, pre_proc_rhs]() {
242 return EssentialPreProc<DisplacementCubitBcData>(mField, pre_proc_rhs,
243 {});
244 };
245 pre_proc_rhs->preProcessHook = get_pre_proc_hook();
246
247 auto get_post_proc_hook_rhs = [this, post_proc_rhs]() {
249
250 CHKERR EssentialPostProcRhs<DisplacementCubitBcData>(mField,
251 post_proc_rhs, 1.)();
253 };
254
255 auto get_post_proc_hook_lhs = [this, post_proc_lhs]() {
257
258 CHKERR EssentialPostProcLhs<DisplacementCubitBcData>(mField,
259 post_proc_lhs, 1.)();
261 };
262
263 post_proc_rhs->postProcessHook = get_post_proc_hook_rhs;
264 post_proc_lhs->postProcessHook = get_post_proc_hook_lhs;
265
266 ksp_ctx_ptr->getPreProcComputeRhs().push_front(pre_proc_rhs);
267 ksp_ctx_ptr->getPostProcComputeRhs().push_back(post_proc_rhs);
268 ksp_ctx_ptr->getPostProcSetOperators().push_back(post_proc_lhs);
270 };
271
272 auto evaluate_field_at_the_point = [&]() {
274
275 int coords_dim = 3;
276 std::array<double, 3> field_eval_coords{0.0, 0.0, 0.0};
277 PetscBool do_eval_field = PETSC_FALSE;
278 CHKERR PetscOptionsGetRealArray(NULL, NULL, "-field_eval_coords",
279 field_eval_coords.data(), &coords_dim,
281
282 if (do_eval_field) {
283
284 vectorFieldPtr = boost::make_shared<MatrixDouble>();
285 auto field_eval_data =
286 mField.getInterface<FieldEvaluatorInterface>()->getData<DomainEle>();
287
288 CHKERR mField.getInterface<FieldEvaluatorInterface>()
289 ->buildTree<SPACE_DIM>(field_eval_data, simple->getDomainFEName());
290
291 field_eval_data->setEvalPoints(field_eval_coords.data(), 1);
292 auto no_rule = [](int, int, int) { return -1; };
293 auto field_eval_fe_ptr = field_eval_data->feMethodPtr.lock();
294 field_eval_fe_ptr->getRuleHook = no_rule;
295
296 field_eval_fe_ptr->getOpPtrVector().push_back(
297 new OpCalculateVectorFieldValues<SPACE_DIM>("U", vectorFieldPtr));
298
299 CHKERR mField.getInterface<FieldEvaluatorInterface>()
300 ->evalFEAtThePoint<SPACE_DIM>(
301 field_eval_coords.data(), 1e-12, simple->getProblemName(),
302 simple->getDomainFEName(), field_eval_data,
304 QUIET);
305
306 if (vectorFieldPtr->size1()) {
307 auto t_disp = getFTensor1FromMat<SPACE_DIM>(*vectorFieldPtr);
308 if constexpr (SPACE_DIM == 2)
309 MOFEM_LOG("FieldEvaluator", Sev::inform)
310 << "U_X: " << t_disp(0) << " U_Y: " << t_disp(1);
311 else
312 MOFEM_LOG("FieldEvaluator", Sev::inform)
313 << "U_X: " << t_disp(0) << " U_Y: " << t_disp(1)
314 << " U_Z: " << t_disp(2);
315 }
316
318 }
320 };
321
322 CHKERR set_essential_bc();
323 CHKERR kspSetUpAndSolve(solver);
324 CHKERR evaluate_field_at_the_point();
325
327}
#define MOFEM_LOG_SYNCHRONISE(comm)
Synchronise "SYNC" channel.
@ QUIET
@ MF_EXIST
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
PetscErrorCode PetscOptionsGetRealArray(PetscOptions *, const char pre[], const char name[], PetscReal dval[], PetscInt *nmax, PetscBool *set)
auto getDMKspCtx(DM dm)
Get KSP context data structure used by DM.
Definition DMMoFEM.hpp:1248
PetscBool do_eval_field
Evaluate field.
Definition plastic.cpp:119
virtual MoFEMErrorCode kspSetUpAndSolve(SmartPetscObj< KSP > solver)
[Push operators to pipeline]
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0

Member Data Documentation

◆ mField

MoFEM::Interface& ElasticExample::mField
protected
Examples
mofem/tutorials/vec-10/schur_elastic.cpp.

Definition at line 24 of file ElasticExample.hpp.

◆ vectorFieldPtr

boost::shared_ptr<MatrixDouble> ElasticExample::vectorFieldPtr = nullptr
protected

Definition at line 25 of file ElasticExample.hpp.


The documentation for this struct was generated from the following file: