v0.15.5
Loading...
Searching...
No Matches
ElasticExample.hpp
Go to the documentation of this file.
1/**
2 * @file ElasticExample.hpp
3 *
4 * @brief Implementation of elastic example class
5 *
6 * @copyright Copyright (c) 2025
7 *
8 */
9
10#include <ElasticSpring.hpp>
11#include <FluidLevel.hpp>
12#include <CalculateTraction.hpp>
13#include <NaturalDomainBC.hpp>
14#include <NaturalBoundaryBC.hpp>
15#include <HookeOps.hpp>
16#include <ElasticPostProc.hpp>
17
19
20 ElasticExample(MoFEM::Interface &m_field) : mField(m_field) {}
21
22 MoFEMErrorCode runProblem();
23
24protected:
26 boost::shared_ptr<MatrixDouble> vectorFieldPtr = nullptr;
27
28 MoFEMErrorCode readMesh();
29 MoFEMErrorCode setupProblem();
30 MoFEMErrorCode boundaryCondition();
31 MoFEMErrorCode assembleSystem();
32 MoFEMErrorCode solveSystem();
33 MoFEMErrorCode outputResults();
34 MoFEMErrorCode checkResults();
35
36 // Auxiliary functions
37 virtual MoFEMErrorCode kspSetUpAndSolve(SmartPetscObj<KSP> solver);
38};
39
40//! [Run problem]
52//! [Run problem]
53
54//! [Read mesh]
55MoFEMErrorCode ElasticExample::readMesh() {
57 auto simple = mField.getInterface<Simple>();
58 CHKERR simple->getOptions();
59 CHKERR simple->loadFile();
60 // Add meshsets if config file provided
61 CHKERR mField.getInterface<MeshsetsManager>()->setMeshsetFromFile();
63}
64//! [Read mesh]
65
66//! [Set up problem]
69 Simple *simple = mField.getInterface<Simple>();
70
71 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOPT };
72 const char *list_bases[LASBASETOPT] = {"ainsworth", "demkowicz"};
73 PetscInt choice_base_value = AINSWORTH;
74 CHKERR PetscOptionsGetEList(PETSC_NULLPTR, NULL, "-base", list_bases,
75 LASBASETOPT, &choice_base_value, PETSC_NULLPTR);
76
78 switch (choice_base_value) {
79 case AINSWORTH:
81 MOFEM_LOG("WORLD", Sev::inform)
82 << "Set AINSWORTH_LEGENDRE_BASE for displacements";
83 break;
84 case DEMKOWICZ:
86 MOFEM_LOG("WORLD", Sev::inform)
87 << "Set DEMKOWICZ_JACOBI_BASE for displacements";
88 break;
89 default:
90 base = LASTBASE;
91 break;
92 }
93
94 // Add field
95 CHKERR simple->addDomainField("U", H1, base, SPACE_DIM);
96 CHKERR simple->addBoundaryField("U", H1, base, SPACE_DIM);
97 int order = 3;
98 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-order", &order, PETSC_NULLPTR);
99
100 CHKERR simple->addDataField("GEOMETRY", H1, base, SPACE_DIM);
101
102 CHKERR simple->setFieldOrder("U", order);
103 CHKERR simple->setFieldOrder("GEOMETRY", 2);
104 CHKERR simple->setUp();
105
106 auto project_ho_geometry = [&]() {
107 Projection10NodeCoordsOnField ent_method(mField, "GEOMETRY");
108 return mField.loop_dofs("GEOMETRY", ent_method);
109 };
110 CHKERR project_ho_geometry();
111
113}
114//! [Set up problem]
115
116//! [Boundary condition]
119 auto simple = mField.getInterface<Simple>();
120 auto bc_mng = mField.getInterface<BcManager>();
121
122 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_X",
123 "U", 0, 0);
124 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_Y",
125 "U", 1, 1);
126 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_Z",
127 "U", 2, 2);
128 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(),
129 "REMOVE_ALL", "U", 0, 3);
130 CHKERR bc_mng->pushMarkDOFsOnEntities<DisplacementCubitBcData>(
131 simple->getProblemName(), "U");
132
134}
135//! [Boundary condition]
136
137//! [Push operators to pipeline]
140 auto pip = mField.getInterface<PipelineManager>();
141
142 //! [Integration rule]
143 auto integration_rule = [](int, int, int approx_order) {
144 return 2 * approx_order + 1;
145 };
146
147 auto integration_rule_bc = [](int, int, int approx_order) {
148 return 2 * approx_order + 1;
149 };
150
151 CHKERR pip->setDomainRhsIntegrationRule(integration_rule);
152 CHKERR pip->setDomainLhsIntegrationRule(integration_rule);
153 CHKERR pip->setBoundaryRhsIntegrationRule(integration_rule_bc);
154 CHKERR pip->setBoundaryLhsIntegrationRule(integration_rule_bc);
155 //! [Integration rule]
156
158 pip->getOpDomainLhsPipeline(), {H1}, "GEOMETRY");
160 pip->getOpDomainRhsPipeline(), {H1}, "GEOMETRY");
162 pip->getOpBoundaryRhsPipeline(), {NOSPACE}, "GEOMETRY");
164 pip->getOpBoundaryLhsPipeline(), {NOSPACE}, "GEOMETRY");
165
166 //! [Push domain stiffness matrix to pipeline]
167 // Add LHS operator for elasticity (stiffness matrix)
168 CHKERR HookeOps::opFactoryDomainLhs<SPACE_DIM, A, I, DomainEleOp>(
169 mField, pip->getOpDomainLhsPipeline(), "U", "MAT_ELASTIC", Sev::verbose);
170 //! [Push domain stiffness matrix to pipeline]
171
172 // Add RHS operator for internal forces
173 CHKERR HookeOps::opFactoryDomainRhs<SPACE_DIM, A, I, DomainEleOp>(
174 mField, pip->getOpDomainRhsPipeline(), "U", "MAT_ELASTIC", Sev::verbose);
175
176 //! [Push Body forces]
178 pip->getOpDomainRhsPipeline(), mField, "U", Sev::inform);
179 //! [Push Body forces]
180
181 //! [Push natural boundary conditions]
182 // Add force boundary condition
184 pip->getOpBoundaryRhsPipeline(), mField, "U", -1, Sev::inform);
185 // Add case for mix type of BCs
187 pip->getOpBoundaryLhsPipeline(), mField, "U", Sev::verbose);
188 //! [Push natural boundary conditions]
190}
191//! [Push operators to pipeline]
192
193// ! [KSP set up]
194MoFEMErrorCode ElasticExample::kspSetUpAndSolve(SmartPetscObj<KSP> solver) {
196
197 MOFEM_LOG_CHANNEL("TIMER");
198 MOFEM_LOG_TAG("TIMER", "timer");
199
200 BOOST_LOG_SCOPED_THREAD_ATTR("Timeline", attrs::timer());
201 MOFEM_LOG("TIMER", Sev::inform) << "KSPSetUp";
202 CHKERR KSPSetUp(solver);
203 MOFEM_LOG("TIMER", Sev::inform) << "KSPSetUp <= Done";
204
205 DM dm;
206 CHKERR KSPGetDM(solver, &dm);
207 auto D = createDMVector(dm);
208 auto F = vectorDuplicate(D);
209
210 MOFEM_LOG("TIMER", Sev::inform) << "KSPSolve";
211 CHKERR KSPSolve(solver, F, D);
212 MOFEM_LOG("TIMER", Sev::inform) << "KSPSolve <= Done";
213
214 CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
215 CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
216 CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
217
219};
220// ! [KSP set up]
221
222//! [Solve]
225
226 auto simple = mField.getInterface<Simple>();
227 auto pip = mField.getInterface<PipelineManager>();
228 auto solver = pip->createKSP();
229 CHKERR KSPSetFromOptions(solver);
230
231 auto set_essential_bc = [this]() {
233 // This is low level pushing finite elements (pipelines) to solver
234 auto simple = mField.getInterface<Simple>();
235 auto dm = simple->getDM();
236 auto ksp_ctx_ptr = getDMKspCtx(dm);
237
238 auto pre_proc_rhs = boost::make_shared<FEMethod>();
239 auto post_proc_rhs = boost::make_shared<FEMethod>();
240 auto post_proc_lhs = boost::make_shared<FEMethod>();
241
242 auto get_pre_proc_hook = [this, pre_proc_rhs]() {
243 return EssentialPreProc<DisplacementCubitBcData>(mField, pre_proc_rhs,
244 {});
245 };
246 pre_proc_rhs->preProcessHook = get_pre_proc_hook();
247
248 auto get_post_proc_hook_rhs = [this, post_proc_rhs]() {
250
251 CHKERR EssentialPostProcRhs<DisplacementCubitBcData>(mField,
252 post_proc_rhs, 1.)();
254 };
255
256 auto get_post_proc_hook_lhs = [this, post_proc_lhs]() {
258
259 CHKERR EssentialPostProcLhs<DisplacementCubitBcData>(mField,
260 post_proc_lhs, 1.)();
262 };
263
264 post_proc_rhs->postProcessHook = get_post_proc_hook_rhs;
265 post_proc_lhs->postProcessHook = get_post_proc_hook_lhs;
266
267 ksp_ctx_ptr->getPreProcComputeRhs().push_front(pre_proc_rhs);
268 ksp_ctx_ptr->getPostProcComputeRhs().push_back(post_proc_rhs);
269 ksp_ctx_ptr->getPostProcSetOperators().push_back(post_proc_lhs);
271 };
272
273 auto evaluate_field_at_the_point = [&]() {
275
276 int coords_dim = 3;
277 std::array<double, 3> field_eval_coords{0.0, 0.0, 0.0};
278 PetscBool do_eval_field = PETSC_FALSE;
279 CHKERR PetscOptionsGetRealArray(NULL, NULL, "-field_eval_coords",
280 field_eval_coords.data(), &coords_dim,
282
283 if (do_eval_field) {
284
285 vectorFieldPtr = boost::make_shared<MatrixDouble>();
286 auto field_eval_data =
287 mField.getInterface<FieldEvaluatorInterface>()->getData<DomainEle>();
288
289 CHKERR mField.getInterface<FieldEvaluatorInterface>()
290 ->buildTree<SPACE_DIM>(field_eval_data, simple->getDomainFEName());
291
292 field_eval_data->setEvalPoints(field_eval_coords.data(), 1);
293 auto no_rule = [](int, int, int) { return -1; };
294 auto field_eval_fe_ptr = field_eval_data->feMethodPtr.lock();
295 field_eval_fe_ptr->getRuleHook = no_rule;
296
297 field_eval_fe_ptr->getOpPtrVector().push_back(
298 new OpCalculateVectorFieldValues<SPACE_DIM>("U", vectorFieldPtr));
299
300 CHKERR mField.getInterface<FieldEvaluatorInterface>()
301 ->evalFEAtThePoint<SPACE_DIM>(
302 field_eval_coords.data(), 1e-12, simple->getProblemName(),
303 simple->getDomainFEName(), field_eval_data,
305 QUIET);
306
307 if (vectorFieldPtr->size1()) {
308 auto t_disp = getFTensor1FromMat<SPACE_DIM>(*vectorFieldPtr);
309 if constexpr (SPACE_DIM == 2)
310 MOFEM_LOG("FieldEvaluator", Sev::inform)
311 << "U_X: " << t_disp(0) << " U_Y: " << t_disp(1);
312 else
313 MOFEM_LOG("FieldEvaluator", Sev::inform)
314 << "U_X: " << t_disp(0) << " U_Y: " << t_disp(1)
315 << " U_Z: " << t_disp(2);
316 }
317
319 }
321 };
322
323 CHKERR set_essential_bc();
324 CHKERR kspSetUpAndSolve(solver);
325 CHKERR evaluate_field_at_the_point();
326
328}
329//! [Solve]
330
331//! [Postprocess results]
334 auto simple = mField.getInterface<Simple>();
337 mField, simple->getDM(), simple->getDomainFEName(), "out_elastic.h5m");
339}
340//! [Postprocess results]
341
342//! [Check]
344 MOFEM_LOG_CHANNEL("WORLD");
345 auto simple = mField.getInterface<Simple>();
346 auto pip = mField.getInterface<PipelineManager>();
348 pip->getDomainRhsFE().reset();
349 pip->getDomainLhsFE().reset();
350 pip->getBoundaryRhsFE().reset();
351 pip->getBoundaryLhsFE().reset();
352
353 auto integration_rule = [](int, int, int p_data) { return 2 * p_data + 1; };
354 CHKERR pip->setDomainRhsIntegrationRule(integration_rule);
355 CHKERR pip->setBoundaryRhsIntegrationRule(integration_rule);
356
358 pip->getOpDomainRhsPipeline(), {H1}, "GEOMETRY");
360 pip->getOpBoundaryRhsPipeline(), {}, "GEOMETRY");
361
362 // Add RHS operators for internal forces
363 CHKERR HookeOps::opFactoryDomainRhs<SPACE_DIM, A, I, DomainEleOp>(
364 mField, pip->getOpDomainRhsPipeline(), "U", "MAT_ELASTIC", Sev::verbose);
365
367 pip->getOpDomainRhsPipeline(), mField, "U", Sev::verbose);
368
370 pip->getOpBoundaryRhsPipeline(), mField, "U", -1, Sev::verbose);
371
372 auto dm = simple->getDM();
373 auto res = createDMVector(dm);
374 CHKERR VecSetDM(res, PETSC_NULLPTR);
375
376 pip->getDomainRhsFE()->f = res;
377 pip->getBoundaryRhsFE()->f = res;
378
379 CHKERR VecZeroEntries(res);
380
381 CHKERR pip->loopFiniteElements();
382 // CHKERR mField.getInterface<FieldBlas>()->fieldScale(-1, "U");
383
384 CHKERR VecGhostUpdateBegin(res, ADD_VALUES, SCATTER_REVERSE);
385 CHKERR VecGhostUpdateEnd(res, ADD_VALUES, SCATTER_REVERSE);
386 CHKERR VecAssemblyBegin(res);
387 CHKERR VecAssemblyEnd(res);
388
389 auto zero_residual_at_constrains = [&]() {
391 auto fe_post_proc_ptr = boost::make_shared<FEMethod>();
392 auto get_post_proc_hook_rhs = [&]() {
394 CHKERR EssentialPreProcReaction<DisplacementCubitBcData>(
395 mField, fe_post_proc_ptr, res)();
396 CHKERR EssentialPostProcRhs<DisplacementCubitBcData>(
397 mField, fe_post_proc_ptr, 0, res)();
399 };
400 fe_post_proc_ptr->postProcessHook = get_post_proc_hook_rhs;
401 CHKERR DMoFEMPostProcessFiniteElements(dm, fe_post_proc_ptr.get());
403 };
404
405 CHKERR zero_residual_at_constrains();
406
407 double nrm2;
408 CHKERR VecNorm(res, NORM_2, &nrm2);
409 MOFEM_LOG_CHANNEL("WORLD");
410 MOFEM_LOG_C("WORLD", Sev::inform, "residual = %3.4e\n", nrm2);
411
412 int test = 0;
413 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-test", &test, PETSC_NULLPTR);
414 if (test > 0) {
415 auto post_proc_residual = [&](auto dm, auto f_res, auto out_name) {
417 auto post_proc_fe =
418 boost::make_shared<PostProcBrokenMeshInMoab<DomainEle>>(mField);
419 using OpPPMap = OpPostProcMapInMoab<SPACE_DIM, SPACE_DIM>;
420 auto u_vec = boost::make_shared<MatrixDouble>();
421 post_proc_fe->getOpPtrVector().push_back(
422 new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_vec, f_res));
423 post_proc_fe->getOpPtrVector().push_back(
424
425 new OpPPMap(
426
427 post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
428
429 {},
430
431 {{"RES", u_vec}},
432
433 {}, {})
434
435 );
436
437 CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(),
438 post_proc_fe);
439 post_proc_fe->writeFile(out_name);
441 };
442
443 CHKERR post_proc_residual(simple->getDM(), res, "res.h5m");
444
445 constexpr double eps = 1e-8;
446 if (nrm2 > eps)
447 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
448 "Residual is not zero");
449 }
450 if (test == 2) {
451 if (!vectorFieldPtr || vectorFieldPtr->size1() == 0) {
452 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
453 "atom test %d failed: Field Evaluator did not provide result",
454 test);
455 }
456 auto t_disp = getFTensor1FromMat<SPACE_DIM>(*vectorFieldPtr);
457 double Ux_ref = 0.46;
458 double Uy_ref = -0.03;
459 constexpr double eps = 1e-8;
460 if (fabs(t_disp(0) - Ux_ref) > eps || fabs(t_disp(1) - Uy_ref) > eps) {
461 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
462 "atom test %d failed: Ux_ref = %3.6e, computed = %3.6e, Uy_ref "
463 "= %3.6e, computed = %3.6e",
464 test, Ux_ref, t_disp(0), Uy_ref, t_disp(1));
465 }
466 }
468}
469//! [Check]
#define MOFEM_LOG_SYNCHRONISE(comm)
Synchronise "SYNC" channel.
#define MOFEM_LOG_C(channel, severity, format,...)
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
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
@ QUIET
@ MF_EXIST
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 MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ 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
@ F
auto integration_rule
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
#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
MoFEMErrorCode postProcessElasticResults(MoFEM::Interface &mField, SmartPetscObj< DM > dm, const std::string &domain_fe_name, const std::string &out_file_name, SmartPetscObj< Vec > extra_vector=nullptr, const std::string &extra_vector_name="", const Sev hooke_ops_sev=Sev::verbose)
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
static constexpr int approx_order
MoFEMErrorCode solveSystem()
[Solve]
MoFEMErrorCode outputResults()
[Solve]
MoFEMErrorCode assembleSystem()
[Boundary condition]
MoFEMErrorCode setupProblem()
[Read mesh]
virtual MoFEMErrorCode kspSetUpAndSolve(SmartPetscObj< KSP > solver)
[Push operators to pipeline]
MoFEMErrorCode runProblem()
[Run problem]
MoFEM::Interface & mField
MoFEMErrorCode boundaryCondition()
[Set up problem]
MoFEMErrorCode checkResults()
[Postprocess results]
MoFEMErrorCode readMesh()
[Run problem]
boost::shared_ptr< MatrixDouble > vectorFieldPtr
ElasticExample(MoFEM::Interface &m_field)
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.
PetscBool do_eval_field
Evaluate field.
Definition plastic.cpp:119
ElementsAndOps< SPACE_DIM >::SideEle SideEle
Definition plastic.cpp:61
Calculate traction for linear problem.
Implementation of elastic spring bc.
Natural boundary condition applying pressure from fluid.
Implementation of Hookes operator Hookes for linear elastic problems in MoFEM.
Implementation of natural boundary conditions.
Boundary conditions in domain, i.e. body forces.