v0.15.0
Loading...
Searching...
No Matches
ExampleArcLength.hpp
Go to the documentation of this file.
1/**
2 * @file ExampleArcLength.hpp
3 * @example ExampleArcLength.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
12
14
17
18 MoFEMErrorCode runProblem();
19
20protected:
21 MoFEMErrorCode ArcLengthSolve();
22 MoFEMErrorCode checkResults();
23};
24
25//! [Run problem]
28
37
38 ExampleTimeScale::snesPtr.reset();
39
41}
42//! [Run problem]
43
44//! [Arc Length Solve]
47 auto *simple = mField.getInterface<Simple>();
48 auto *pipeline_mng = mField.getInterface<PipelineManager>();
49
50 auto dm = simple->getDM();
51 auto snes = pipeline_mng->createSNES();
52 ExampleTimeScale::snesPtr = snes;
53
54 PetscBool post_proc_vol;
55 PetscBool post_proc_skin;
56
57 if constexpr (SPACE_DIM == 2) {
58 post_proc_vol = PETSC_TRUE;
59 post_proc_skin = PETSC_FALSE;
60 } else {
61 post_proc_vol = PETSC_FALSE;
62 post_proc_skin = PETSC_TRUE;
63 }
64 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-post_proc_vol",
65 &post_proc_vol, PETSC_NULLPTR);
66 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-post_proc_skin",
67 &post_proc_skin, PETSC_NULLPTR);
68
69 // Setup postprocessing
70 auto create_post_proc_fe = [dm, this, simple, post_proc_vol,
71 post_proc_skin]() {
72 auto post_proc_ele_domain = [dm, this](auto &pip_domain) {
74 "GEOMETRY");
75 auto common_ptr = commonDataFactory<SPACE_DIM, GAUSS, DomainEleOp>(
76 mField, pip_domain, "U", "MAT_ELASTIC", Sev::inform);
77 return common_ptr;
78 };
79
80 auto post_proc_map = [this](auto &pip, auto u_ptr, auto common_ptr) {
82
83 using OpPPMap = OpPostProcMapInMoab<SPACE_DIM, SPACE_DIM>;
84
85 auto X_ptr = boost::make_shared<MatrixDouble>();
86
87 pip->getOpPtrVector().push_back(
88 new OpCalculateVectorFieldValues<SPACE_DIM>("GEOMETRY", X_ptr));
89
90 pip->getOpPtrVector().push_back(
91
92 new OpPPMap(pip->getPostProcMesh(), pip->getMapGaussPts(), {},
93 {{"U", u_ptr}, {"GEOMETRY", X_ptr}},
94 {{"GRAD", common_ptr->matGradPtr},
95 {"FIRST_PIOLA", common_ptr->getMatFirstPiolaStress()}},
96 {}));
98 };
99
100 auto push_post_proc_bdy = [dm, this, simple, post_proc_skin,
101 &post_proc_ele_domain,
102 &post_proc_map](auto &pip_bdy) {
103 if (post_proc_skin == PETSC_FALSE)
104 return boost::shared_ptr<PostProcEleBdy>();
105
106 auto u_ptr = boost::make_shared<MatrixDouble>();
107 pip_bdy->getOpPtrVector().push_back(
108 new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
109 auto op_loop_side =
110 new OpLoopSide<SideEle>(mField, simple->getDomainFEName(), SPACE_DIM);
111 auto common_ptr = post_proc_ele_domain(op_loop_side->getOpPtrVector());
112 pip_bdy->getOpPtrVector().push_back(op_loop_side);
113
114 CHKERR post_proc_map(pip_bdy, u_ptr, common_ptr);
115
116 return pip_bdy;
117 };
118
119 auto push_post_proc_domain = [dm, this, simple, post_proc_vol,
120 &post_proc_ele_domain,
121 &post_proc_map](auto &pip_domain) {
122 if (post_proc_vol == PETSC_FALSE)
123 return boost::shared_ptr<PostProcEleDomain>();
124
125 auto u_ptr = boost::make_shared<MatrixDouble>();
126 pip_domain->getOpPtrVector().push_back(
127 new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
128 auto common_ptr = post_proc_ele_domain(pip_domain->getOpPtrVector());
129
130 CHKERR post_proc_map(pip_domain, u_ptr, common_ptr);
131
132 return pip_domain;
133 };
134
135 auto post_proc_fe_domain = boost::make_shared<PostProcEleDomain>(mField);
136 auto post_proc_fe_bdy = boost::make_shared<PostProcEleBdy>(mField);
137
138 return std::make_pair(push_post_proc_domain(post_proc_fe_domain),
139 push_post_proc_bdy(post_proc_fe_bdy));
140 };
141
142 auto add_arc_length_load_tangent_operator = [&](auto snes) {
144 auto bdy_pipeline = boost::make_shared<BoundaryEle>(mField);
145 auto domain_pipeline = boost::make_shared<DomainEle>(mField);
146
147 auto integration_rule = [](int, int, int approx_order) {
148 return 2 * (approx_order - 1);
149 };
150 bdy_pipeline->getRuleHook = integration_rule;
151 domain_pipeline->getRuleHook = integration_rule;
152
153 struct TimeOne : TimeScale {
154 double getScale(const double) override { return -1.; }
155 };
156 auto time_one = boost::make_shared<TimeOne>();
157
158 // The trick is here, that is assumed that load is applied linearly
159 // proportional to the load factor (time in this case). So derivative of the
160 // load with respect to the load factor is the same as the load at time
161 // = 1.0
162
163 CHKERR BoundaryNaturalBC::AddFluxToPipeline<OpForce>::add(
164 bdy_pipeline->getOpPtrVector(), mField, "U", {time_one}, "FORCE",
165 Sev::inform);
166 //! [Define gravity vector]
167 CHKERR DomainNaturalBC::AddFluxToPipeline<OpBodyForce>::add(
168 domain_pipeline->getOpPtrVector(), mField, "U", {time_one},
169 "BODY_FORCE", Sev::inform);
170
171 auto snes_ctx_ptr = getDMSnesCtx(simple->getDM());
172 snes_ctx_ptr->getLoadTangent().push_back(
173 {simple->getBoundaryFEName(), bdy_pipeline});
174 snes_ctx_ptr->getLoadTangent().push_back(
175 {simple->getDomainFEName(), domain_pipeline});
176
177 auto time_scale = boost::make_shared<ExampleTimeScale>();
178
179 auto pre_proc_ptr = boost::make_shared<FEMethod>();
180 auto get_bc_hook_rhs = [this, pre_proc_ptr, time_scale]() {
182 CHKERR EssentialPreProc<DisplacementCubitBcData>(mField, pre_proc_ptr,
183 {time_scale}, false)();
185 };
186 pre_proc_ptr->preProcessHook = get_bc_hook_rhs;
187 snes_ctx_ptr->getPreProcLoadTangent().push_back(pre_proc_ptr);
188 auto post_proc_rhs_ptr = boost::make_shared<FEMethod>();
189 auto get_post_proc_hook_rhs = [this, post_proc_rhs_ptr]() {
191 CHKERR EssentialPostProcRhs<DisplacementCubitBcData>(
192 mField, post_proc_rhs_ptr, 0.)();
194 };
195 post_proc_rhs_ptr->postProcessHook = get_post_proc_hook_rhs;
196 snes_ctx_ptr->getPostProcLoadTangent().push_back(post_proc_rhs_ptr);
197
198 CHKERR SNESNewtonALSetFunction(snes, SnesLoadTangent, snes_ctx_ptr.get());
199
201 };
202
203 auto add_extra_finite_elements_to_solver_pipelines = [&]() {
205
206 auto pre_proc_ptr = boost::make_shared<FEMethod>();
207 auto post_proc_rhs_ptr = boost::make_shared<FEMethod>();
208 auto post_proc_lhs_ptr = boost::make_shared<FEMethod>();
209
210 auto time_scale = boost::make_shared<ExampleTimeScale>();
211
212 auto get_bc_hook_rhs = [this, pre_proc_ptr, time_scale]() {
214 CHKERR EssentialPreProc<DisplacementCubitBcData>(mField, pre_proc_ptr,
215 {time_scale}, false)();
217 };
218
219 pre_proc_ptr->preProcessHook = get_bc_hook_rhs;
220
221 auto get_post_proc_hook_rhs = [this, post_proc_rhs_ptr]() {
223 CHKERR EssentialPreProcReaction<DisplacementCubitBcData>(
224 mField, post_proc_rhs_ptr, nullptr, Sev::verbose)();
225 CHKERR EssentialPostProcRhs<DisplacementCubitBcData>(
226 mField, post_proc_rhs_ptr, 1.)();
228 };
229 auto get_post_proc_hook_lhs = [this, post_proc_lhs_ptr]() {
231 CHKERR EssentialPostProcLhs<DisplacementCubitBcData>(
232 mField, post_proc_lhs_ptr, 1.)();
234 };
235 post_proc_rhs_ptr->postProcessHook = get_post_proc_hook_rhs;
236 post_proc_lhs_ptr->postProcessHook = get_post_proc_hook_lhs;
237
238 // This is low level pushing finite elements (pipelines) to solver
239 auto snes_ctx_ptr = getDMSnesCtx(simple->getDM());
240 snes_ctx_ptr->getPreProcComputeRhs().push_front(pre_proc_ptr);
241 snes_ctx_ptr->getPreProcSetOperators().push_front(pre_proc_ptr);
242 snes_ctx_ptr->getPostProcComputeRhs().push_back(post_proc_rhs_ptr);
243 snes_ctx_ptr->getPostProcSetOperators().push_back(post_proc_lhs_ptr);
245 };
246
247 auto set_log_and_post_hook = [this](auto post_proc_fe, auto post_proc_bdy) {
249
250 auto make_vtk = [this, post_proc_fe, post_proc_bdy](auto iter,
251 auto cache_ptr) {
253 auto simple = mField.getInterface<Simple>();
254 auto dm = simple->getDM();
255 if (post_proc_fe) {
256 CHKERR DMoFEMLoopFiniteElements(dm, "dFE", post_proc_fe, cache_ptr);
257 CHKERR post_proc_fe->writeFile(
258 "out_step_" + boost::lexical_cast<std::string>(iter) + ".h5m");
259 }
260 if (post_proc_bdy) {
261 CHKERR DMoFEMLoopFiniteElements(dm, "bFE", post_proc_bdy, cache_ptr);
262 CHKERR post_proc_bdy->writeFile(
263 "out_skin_" + boost::lexical_cast<std::string>(iter) + ".h5m");
264 }
266 };
267
268 auto hook = [this, make_vtk](SNES snes, Vec x, Vec F,
269 boost::shared_ptr<CacheTuple> cache_ptr,
270 void *ctx) -> MoFEMErrorCode {
272
273 double atol, rtol, stol;
274 int maxit, maxf;
275 CHKERR SNESGetTolerances(snes, &atol, &rtol, &stol, &maxit, &maxf);
276
277 double lambda;
278 CHKERR SNESNewtonALGetLoadParameter(snes, &lambda);
279
280 double nrm;
281 CHKERR VecNorm(F, NORM_2, &nrm);
282
283 int iter;
284 CHKERR SNESGetIterationNumber(snes, &iter);
285 SNESConvergedReason reason;
286 CHKERR SNESGetConvergedReason(snes, &reason);
287 MOFEM_LOG_C("EXAMPLE", Sev::inform,
288 "SNES iteration %d, reason %d nrm2 %g lambda %g", iter,
289 reason, nrm, lambda);
290
291 CHKERR make_vtk(iter, cache_ptr);
292
294 };
295
296 auto simple = mField.getInterface<Simple>();
297 auto snes_ctx_ptr = getDMSnesCtx(simple->getDM());
298 snes_ctx_ptr->getRhsHook() = hook;
299
301 };
302
303 // Set arc length solver
304 auto D = createDMVector(simple->getDM());
305 CHKERR SNESSetType(snes, SNESNEWTONLS);
306 CHKERR SNESSetFromOptions(snes);
307
308 auto B = createDMMatrix(dm);
309 CHKERR SNESSetJacobian(snes, B, B, PETSC_NULLPTR, PETSC_NULLPTR);
310 // Add extra finite elements to SNES solver pipelines to resolve essential
311 // boundary conditions
312 CHKERR add_extra_finite_elements_to_solver_pipelines();
313 CHKERR add_arc_length_load_tangent_operator(snes);
314
315 auto [post_proc_fe, post_proc_bdy] = create_post_proc_fe();
316 CHKERR set_log_and_post_hook(post_proc_fe, post_proc_bdy);
317
318 CHKERR SNESSolve(snes, PETSC_NULLPTR, D);
319
321}
322//! [Arc Length Solve]
323
324//! [Check]
327 PetscInt test_nb = 0;
328 PetscBool test_flg = PETSC_FALSE;
329 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-test", &test_nb, &test_flg);
330
331 if (test_flg) {
332 auto simple = mField.getInterface<Simple>();
333 auto T = createDMVector(simple->getDM());
334 CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
335 SCATTER_FORWARD);
336 double nrm2;
337 CHKERR VecNorm(T, NORM_2, &nrm2);
338 MOFEM_LOG("EXAMPLE", Sev::verbose) << "Regression norm " << nrm2;
339 double regression_value = 0;
340 switch (test_nb) {
341 case 1:
342 regression_value = 1.215134e-04;
343 break;
344 default:
345 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Wrong test number.");
346 break;
347 }
348 if (fabs(nrm2 - regression_value) > 1e-2)
349 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
350 "Regression test field; wrong norm value. %6.4e != %6.4e", nrm2,
351 regression_value);
352 }
354}
355//! [Check]
#define MOFEM_LOG_C(channel, severity, format,...)
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
constexpr int SPACE_DIM
#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
#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 ...
@ F
auto integration_rule
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:1102
auto createDMMatrix(DM dm)
Get smart matrix from DM.
Definition DMMoFEM.hpp:1059
#define MOFEM_LOG(channel, severity)
Log.
static double lambda
double D
auto getDMSnesCtx(DM dm)
Get SNES context data structure used by DM.
Definition DMMoFEM.hpp:1130
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
static constexpr int approx_order
MoFEMErrorCode checkResults()
[Arc Length Solve]
MoFEMErrorCode ArcLengthSolve()
[Run problem]
ExampleArcLength(MoFEM::Interface &m_field)
MoFEMErrorCode runProblem()
[Run problem]
MoFEMErrorCode readMesh()
[Run problem]
MoFEMErrorCode gettingNorms()
[TS Solve]
MoFEMErrorCode outputResults()
[Getting norms]
MoFEMErrorCode assembleSystem()
[Boundary condition]
MoFEMErrorCode setupProblem()
[Read mesh]
MoFEMErrorCode boundaryCondition()
[Set up problem]
Deprecated interface functions.
Post post-proc data at points from hash maps.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.