v0.15.5
Loading...
Searching...
No Matches
SnesCtx.cpp
Go to the documentation of this file.
1/**
2 * @brief Context for PETSc SNES, i.e. nonlinear solver
3 * @author Anonymous contributor under MIT license
4 *
5 */
6
7namespace MoFEM {
8
10
15
16 SnesCtxImpl(Interface &m_field, std::string problem_name);
17
18 virtual ~SnesCtxImpl() = default;
19
21
23
25
27
29
33
35
39
43
44 MoFEMErrorCode copyLoops(const SnesCtxImpl &snes_ctx);
45
47
49
51
53
54 friend PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx);
55
56 friend PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx);
57
59 MatAssemblyType type);
61
63 MatAssemblyType type);
64
66
67 friend MoFEMErrorCode MoFEMSNESMonitorFields(SNES snes, PetscInt its,
68 PetscReal fgnorm,
69 SnesCtx *snes_ctx);
70
71 friend MoFEMErrorCode MoFEMSNESMonitorEnergy(SNES snes, PetscInt its,
72 PetscReal fgnorm,
73 SnesCtx *snes_ctx);
74
75 friend PetscErrorCode SnesLoadTangent(SNES snes, Vec u, Vec F, void *ctx);
76
77protected:
78 PetscLogEvent MOFEM_EVENT_SnesRhs; ///< Log events to assemble residual
79 PetscLogEvent MOFEM_EVENT_SnesMat; ///< Log events to assemble tangent matrix
80 PetscLogEvent MOFEM_EVENT_SnesLoadTangent; ///< Log events to assemble tangent
81 ///< of load vector
82
83 MoFEM::Interface &mField; ///< database Interface
84 moab::Interface &moab; ///< moab Interface
85
86 std::string problemName; ///< problem name
87
88 FEMethodsSequence loopsOperator; ///< Sequence of finite elements instances
89 ///< assembling tangent matrix
90 FEMethodsSequence loopsRhs; ///< Sequence of finite elements instances
91 ///< assembling residual vector
93 loopsLoadTangent; ///< Sequence of finite elements instances assembling
94 ///< tangent of load vector in arc length control
95
96 BasicMethodsSequence preProcessOperator; ///< Sequence of methods run before
97 ///< tangent matrix is assembled
98 BasicMethodsSequence postProcessOperator; ///< Sequence of methods run after
99 ///< tangent matrix is assembled
101 preProcessLoadTangent; ///< Sequence of methods run before
102 ///< tangent of load vector in arc
103 ///< length control is assembled
105 preProcessRhs; ///< Sequence of methods run before residual is assembled
107 postProcessRhs; ///< Sequence of methods run after residual is assembled
109 postProcessLoadTangent; ///< Sequence of methods run after tangent of load
110 ///< vector in arc length control is assembled
111
113 MF_EXIST; ///< If set to MF_EXIST check if element exist, default MF_EXIST
114 bool zeroPreCondMatrixB = true; ///< If true zero matrix, otherwise user need
115 ///< to do it, default true
116 bool vErify = false; ///< If true verify vector
117 MatAssemblyType typeOfAssembly =
118 MAT_FINAL_ASSEMBLY; ///< type of assembly at the end
119
120 HookFunction loadTangentHook; ///< Hook function called at the end of load
121 ///< tangent evaluation
122
123 HookFunction rhsHook; ///< Hook function called at the end of residual
124 ///< evaluation
125};
126
127SnesCtx::SnesCtx(Interface &m_field, std::string problem_name)
128 : snesCtxImpl(boost::movelib::make_unique<SnesCtx::SnesCtxImpl>(
129 m_field, problem_name)) {}
130
131SnesCtx::~SnesCtx() = default;
132
134 return snesCtxImpl->getSetOperators();
135}
136
138 return snesCtxImpl->getComputeRhs();
139}
140
142 return snesCtxImpl->getPreProcComputeRhs();
143}
144
146 return snesCtxImpl->getPostProcComputeRhs();
147}
148
150 return snesCtxImpl->getPreProcSetOperators();
151}
152
154 return snesCtxImpl->getPostProcSetOperators();
155}
156
158 return snesCtxImpl->getLoadTangent();
159}
160
162 return snesCtxImpl->getPreProcLoadTangent();
163}
164
166 return snesCtxImpl->getPostProcLoadTangent();
167}
168
170 return snesCtxImpl->copyLoops(*snes_ctx.snesCtxImpl);
171}
172
174
176
178 return snesCtxImpl->getLoadTangentHook();
179}
180
182 return snesCtxImpl->getRhsHook();
183}
184
185SnesCtx::SnesCtxImpl::SnesCtxImpl(Interface &m_field, std::string problem_name)
186 : mField(m_field), moab(m_field.get_moab()), problemName(problem_name) {
187 PetscLogEventRegister("LoopSNESRhs", 0, &MOFEM_EVENT_SnesRhs);
188 PetscLogEventRegister("LoopSNESMat", 0, &MOFEM_EVENT_SnesMat);
189 PetscLogEventRegister("LoopSNESLoadTangent", 0, &MOFEM_EVENT_SnesLoadTangent);
190 if (!LogManager::checkIfChannelExist("SNES_WORLD")) {
191 auto core_log = logging::core::get();
192 core_log->add_sink(
194 LogManager::setLog("SNES_WORLD");
195 MOFEM_LOG_TAG("SNES_WORLD", "SNES");
196 }
197}
198
201 loopsOperator = snes_ctx.loopsOperator;
202 loopsRhs = snes_ctx.loopsRhs;
203 preProcessOperator = snes_ctx.preProcessOperator;
204 postProcessOperator = snes_ctx.postProcessOperator;
205 preProcessRhs = snes_ctx.preProcessRhs;
206 postProcessRhs = snes_ctx.postProcessRhs;
207 loopsLoadTangent = snes_ctx.loopsLoadTangent;
208 preProcessLoadTangent = snes_ctx.preProcessLoadTangent;
209 postProcessLoadTangent = snes_ctx.postProcessLoadTangent;
211}
212
215 loopsOperator.clear();
216 loopsRhs.clear();
217 preProcessOperator.clear();
218 postProcessOperator.clear();
219 preProcessRhs.clear();
220 postProcessRhs.clear();
221 loopsLoadTangent.clear();
222 preProcessLoadTangent.clear();
223 postProcessLoadTangent.clear();
225}
226
227PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx) {
228 auto *snes_ctx = static_cast<SnesCtx *>(ctx)->snesCtxImpl.get();
230 PetscLogEventBegin(snes_ctx->MOFEM_EVENT_SnesRhs, 0, 0, 0, 0);
231 Vec dx;
232 CHKERR SNESGetSolutionUpdate(snes, &dx);
233 CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
234 CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
235 if (snes_ctx->vErify) {
236 // Verify finite elements, check for not a number
237 CHKERR VecAssemblyBegin(f);
238 CHKERR VecAssemblyEnd(f);
239 MPI_Comm comm = PetscObjectComm((PetscObject)f);
240 PetscSynchronizedPrintf(comm, "SNES Verify x\n");
241 const Problem *prb_ptr;
242 CHKERR snes_ctx->mField.get_problem(snes_ctx->problemName, &prb_ptr);
243 CHKERR snes_ctx->mField.getInterface<Tools>()->checkVectorForNotANumber(
244 prb_ptr, COL, x);
245 }
246 CHKERR snes_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
247 snes_ctx->problemName, COL, x, INSERT_VALUES, SCATTER_REVERSE);
248
249 auto zero_ghost_vec = [](Vec g) {
251 Vec l;
252 CHKERR VecGhostGetLocalForm(g, &l);
253 double *a;
254 CHKERR VecGetArray(l, &a);
255 int s;
256 CHKERR VecGetLocalSize(l, &s);
257 for (int i = 0; i != s; ++i)
258 a[i] = 0;
259 CHKERR VecRestoreArray(l, &a);
260 CHKERR VecGhostRestoreLocalForm(g, &l);
262 };
263 CHKERR zero_ghost_vec(f);
264
265 auto vec_assemble_switch = boost::movelib::make_unique<bool>(true);
266 auto cache_ptr = boost::make_shared<CacheTuple>();
267 CHKERR snes_ctx->mField.cache_problem_entities(snes_ctx->problemName,
268 cache_ptr);
269
270 auto set = [&](auto &fe) {
272 fe.snes = snes;
273 fe.snes_x = x;
274 fe.snes_dx = dx;
275 fe.snes_f = f;
277 fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
279
280 CHKERR SNESGetKSP(snes, &fe.ksp);
281
282 // If SNES is of type SNESNEWTONAL then get load parameter, will act as a
283 // pseudo time
284#if PETSC_VERSION_GE(3, 22, 0)
285 PetscBool same;
286 CHKERR PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONAL, &same);
287 if (same) {
288 CHKERR SNESNewtonALGetLoadParameter(snes, &fe.ts_t);
289 fe.data_ctx |= PetscData::CtxSetTime;
290 }
291#endif
292
293 fe.cacheWeakPtr = cache_ptr;
295 };
296
297 auto unset = [&](auto &fe) {
298 fe.snes_ctx = SnesMethod::CTX_SNESNONE;
299 fe.ksp_ctx = KspMethod::CTX_KSPNONE;
300 fe.data_ctx = PetscData::CtxSetNone;
301 };
302
303 for (auto &bit : snes_ctx->preProcessRhs) {
304 bit->vecAssembleSwitch = boost::move(vec_assemble_switch);
305 CHKERR set(*bit);
306 CHKERR snes_ctx->mField.problem_basic_method_preProcess(
307 snes_ctx->problemName, *bit);
308 unset(*bit);
309 vec_assemble_switch = boost::move(bit->vecAssembleSwitch);
310 }
311
312 for (auto &lit : snes_ctx->loopsRhs) {
313 lit.second->vecAssembleSwitch = boost::move(vec_assemble_switch);
314 CHKERR set(*lit.second);
315 CHKERR snes_ctx->mField.loop_finite_elements(
316 snes_ctx->problemName, lit.first, *lit.second, nullptr, snes_ctx->bH,
317 cache_ptr);
318 unset(*lit.second);
319 if (snes_ctx->vErify) {
320 // Verify finite elements, check for not a number
321 CHKERR VecAssemblyBegin(f);
322 CHKERR VecAssemblyEnd(f);
323 MPI_Comm comm = PetscObjectComm((PetscObject)f);
324 PetscSynchronizedPrintf(comm, "SNES Verify f FE < %s >\n",
325 lit.first.c_str());
326 const Problem *prb_ptr;
327 CHKERR snes_ctx->mField.get_problem(snes_ctx->problemName, &prb_ptr);
328 CHKERR snes_ctx->mField.getInterface<Tools>()->checkVectorForNotANumber(
329 prb_ptr, ROW, f);
330 }
331
332 vec_assemble_switch = boost::move(lit.second->vecAssembleSwitch);
333 }
334
335 for (auto &bit : snes_ctx->postProcessRhs) {
336 bit->vecAssembleSwitch = boost::move(vec_assemble_switch);
337 CHKERR set(*bit);
338 CHKERR snes_ctx->mField.problem_basic_method_postProcess(
339 snes_ctx->problemName, *bit);
340 unset(*bit);
341 vec_assemble_switch = boost::move(bit->vecAssembleSwitch);
342 }
343
344 if (*(vec_assemble_switch)) {
345 CHKERR VecGhostUpdateBegin(f, ADD_VALUES, SCATTER_REVERSE);
346 CHKERR VecGhostUpdateEnd(f, ADD_VALUES, SCATTER_REVERSE);
347 CHKERR VecAssemblyBegin(f);
348 CHKERR VecAssemblyEnd(f);
349 }
350
351 if (snes_ctx->rhsHook) {
352 CHKERR snes_ctx->rhsHook(snes, x, f, cache_ptr, ctx);
353 }
354
355 PetscLogEventEnd(snes_ctx->MOFEM_EVENT_SnesRhs, 0, 0, 0, 0);
357}
358
359PetscErrorCode SnesLoadTangent(SNES snes, Vec x, Vec f, void *ctx) {
360 auto *snes_ctx = static_cast<SnesCtx *>(ctx)->snesCtxImpl.get();
362 PetscLogEventBegin(snes_ctx->MOFEM_EVENT_SnesLoadTangent, 0, 0, 0, 0);
363 Vec dx;
364 CHKERR SNESGetSolutionUpdate(snes, &dx);
365 CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
366 CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
367 if (snes_ctx->vErify) {
368 // Verify finite elements, check for not a number
369 CHKERR VecAssemblyBegin(f);
370 CHKERR VecAssemblyEnd(f);
371 MPI_Comm comm = PetscObjectComm((PetscObject)f);
372 PetscSynchronizedPrintf(comm, "SNES Verify x\n");
373 const Problem *prb_ptr;
374 CHKERR snes_ctx->mField.get_problem(snes_ctx->problemName, &prb_ptr);
375 CHKERR snes_ctx->mField.getInterface<Tools>()->checkVectorForNotANumber(
376 prb_ptr, COL, x);
377 }
378 CHKERR snes_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
379 snes_ctx->problemName, COL, x, INSERT_VALUES, SCATTER_REVERSE);
380
381 auto zero_ghost_vec = [](Vec g) {
383 Vec l;
384 CHKERR VecGhostGetLocalForm(g, &l);
385 double *a;
386 CHKERR VecGetArray(l, &a);
387 int s;
388 CHKERR VecGetLocalSize(l, &s);
389 for (int i = 0; i != s; ++i)
390 a[i] = 0;
391 CHKERR VecRestoreArray(l, &a);
392 CHKERR VecGhostRestoreLocalForm(g, &l);
394 };
395 CHKERR zero_ghost_vec(f);
396
397 auto vec_assemble_switch = boost::movelib::make_unique<bool>(true);
398 auto cache_ptr = boost::make_shared<CacheTuple>();
399 CHKERR snes_ctx->mField.cache_problem_entities(snes_ctx->problemName,
400 cache_ptr);
401
402 auto set = [&](auto &fe) {
404 fe.snes = snes;
405 fe.snes_x = x;
406 fe.snes_dx = dx;
407 fe.snes_f = f;
409 fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
411
412 CHKERR SNESGetKSP(snes, &fe.ksp);
413
414 // If SNES is of type SNESNEWTONAL then get load parameter, will act as a
415 // pseudo time
416#if PETSC_VERSION_GE(3, 22, 0)
417 PetscBool same;
418 CHKERR PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONAL, &same);
419 if (same) {
420 CHKERR SNESNewtonALGetLoadParameter(snes, &fe.ts_t);
421 fe.data_ctx |= PetscData::CtxSetTime;
422 }
423#endif
424
425 fe.cacheWeakPtr = cache_ptr;
427 };
428
429 auto unset = [&](auto &fe) {
430 fe.snes_ctx = SnesMethod::CTX_SNESNONE;
431 fe.ksp_ctx = KspMethod::CTX_KSPNONE;
432 fe.data_ctx = PetscData::CtxSetNone;
433 };
434
435 for (auto &bit : snes_ctx->preProcessLoadTangent) {
436 bit->vecAssembleSwitch = boost::move(vec_assemble_switch);
437 CHKERR set(*bit);
438 CHKERR snes_ctx->mField.problem_basic_method_preProcess(
439 snes_ctx->problemName, *bit);
440 unset(*bit);
441 vec_assemble_switch = boost::move(bit->vecAssembleSwitch);
442 }
443
444 for (auto &lit : snes_ctx->loopsLoadTangent) {
445 lit.second->vecAssembleSwitch = boost::move(vec_assemble_switch);
446 CHKERR set(*lit.second);
447 CHKERR snes_ctx->mField.loop_finite_elements(
448 snes_ctx->problemName, lit.first, *lit.second, nullptr, snes_ctx->bH,
449 cache_ptr);
450 unset(*lit.second);
451 if (snes_ctx->vErify) {
452 // Verify finite elements, check for not a number
453 CHKERR VecAssemblyBegin(f);
454 CHKERR VecAssemblyEnd(f);
455 MPI_Comm comm = PetscObjectComm((PetscObject)f);
456 PetscSynchronizedPrintf(comm, "SNES Verify f FE < %s >\n",
457 lit.first.c_str());
458 const Problem *prb_ptr;
459 CHKERR snes_ctx->mField.get_problem(snes_ctx->problemName, &prb_ptr);
460 CHKERR snes_ctx->mField.getInterface<Tools>()->checkVectorForNotANumber(
461 prb_ptr, ROW, f);
462 }
463
464 vec_assemble_switch = boost::move(lit.second->vecAssembleSwitch);
465 }
466
467 for (auto &bit : snes_ctx->postProcessLoadTangent) {
468 bit->vecAssembleSwitch = boost::move(vec_assemble_switch);
469 CHKERR set(*bit);
470 CHKERR snes_ctx->mField.problem_basic_method_postProcess(
471 snes_ctx->problemName, *bit);
472 unset(*bit);
473 vec_assemble_switch = boost::move(bit->vecAssembleSwitch);
474 }
475
476 if (*(vec_assemble_switch)) {
477 CHKERR VecGhostUpdateBegin(f, ADD_VALUES, SCATTER_REVERSE);
478 CHKERR VecGhostUpdateEnd(f, ADD_VALUES, SCATTER_REVERSE);
479 CHKERR VecAssemblyBegin(f);
480 CHKERR VecAssemblyEnd(f);
481 }
482
483 if (snes_ctx->loadTangentHook) {
484 CHKERR snes_ctx->loadTangentHook(snes, x, f, cache_ptr, ctx);
485 }
486
487 PetscLogEventEnd(snes_ctx->MOFEM_EVENT_SnesLoadTangent, 0, 0, 0, 0);
489}
490
491PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx) {
492 auto *snes_ctx = static_cast<SnesCtx *>(ctx)->snesCtxImpl.get();
494 PetscLogEventBegin(snes_ctx->MOFEM_EVENT_SnesMat, 0, 0, 0, 0);
495
496 auto mat_assemble_switch = boost::movelib::make_unique<bool>(true);
497 auto cache_ptr = boost::make_shared<CacheTuple>();
498 CHKERR snes_ctx->mField.cache_problem_entities(snes_ctx->problemName,
499 cache_ptr);
500 Vec dx;
501 CHKERR SNESGetSolutionUpdate(snes, &dx);
502
503 auto set = [&](auto &fe) {
505 fe.snes = snes;
506 fe.snes_x = x;
507 fe.snes_dx = dx;
508 fe.snes_A = A;
509 fe.snes_B = B;
511 fe.ksp_ctx = KspMethod::CTX_OPERATORS;
514
515 CHKERR SNESGetKSP(snes, &fe.ksp);
516
517 fe.cacheWeakPtr = cache_ptr;
519 };
520
521 auto unset = [&](auto &fe) {
522 fe.snes_ctx = SnesMethod::CTX_SNESNONE;
523 fe.ksp_ctx = KspMethod::CTX_KSPNONE;
524 fe.data_ctx = PetscData::CtxSetNone;
525 };
526
527 CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
528 CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
529 CHKERR snes_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
530 snes_ctx->problemName, COL, x, INSERT_VALUES, SCATTER_REVERSE);
531
532 if (*mat_assemble_switch) {
533 if (snes_ctx->zeroPreCondMatrixB)
534 CHKERR MatZeroEntries(B);
535 }
536 for (auto &bit : snes_ctx->preProcessOperator) {
537 bit->matAssembleSwitch = boost::move(mat_assemble_switch);
538 CHKERR set(*bit);
539 CHKERR snes_ctx->mField.problem_basic_method_preProcess(
540 snes_ctx->problemName, *bit);
541 unset(*bit);
542 mat_assemble_switch = boost::move(bit->matAssembleSwitch);
543 }
544
545 for (auto &lit : snes_ctx->loopsOperator) {
546 lit.second->matAssembleSwitch = boost::move(mat_assemble_switch);
547 CHKERR set(*lit.second);
548 CHKERR snes_ctx->mField.loop_finite_elements(
549 snes_ctx->problemName, lit.first, *(lit.second), nullptr, snes_ctx->bH,
550 cache_ptr);
551 unset(*lit.second);
552 mat_assemble_switch = boost::move(lit.second->matAssembleSwitch);
553 }
554
555 for (auto &bit : snes_ctx->postProcessOperator) {
556 bit->matAssembleSwitch = boost::move(mat_assemble_switch);
557 CHKERR set(*bit);
558 CHKERR snes_ctx->mField.problem_basic_method_postProcess(
559 snes_ctx->problemName, *bit);
560 unset(*bit);
561 mat_assemble_switch = boost::move(bit->matAssembleSwitch);
562 }
563
564 if (*mat_assemble_switch) {
565 CHKERR MatAssemblyBegin(B, snes_ctx->typeOfAssembly);
566 CHKERR MatAssemblyEnd(B, snes_ctx->typeOfAssembly);
567 }
568 PetscLogEventEnd(snes_ctx->MOFEM_EVENT_SnesMat, 0, 0, 0, 0);
570}
571
572MoFEMErrorCode SnesMoFEMSetAssemblyType(SNES snes, MatAssemblyType type) {
573
574 auto get_snes_ctx = [](SNES snes) {
575 SnesCtx *snes_ctx;
576 CHK_THROW_MESSAGE(SNESGetApplicationContext(snes, &snes_ctx),
577 "Cannot get SNES context");
578 return snes_ctx->snesCtxImpl.get();
579 };
580
582 get_snes_ctx(snes)->typeOfAssembly = type;
584}
585
587
588 auto get_snes_ctx = [](SNES snes) {
589 SnesCtx *snes_ctx;
590 CHK_THROW_MESSAGE(SNESGetApplicationContext(snes, &snes_ctx),
591 "Cannot get SNES context");
592 return snes_ctx->snesCtxImpl.get();
593 };
594
596 get_snes_ctx(snes)->bH = bh;
598}
599
600MoFEMErrorCode MoFEMSNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm,
601 SnesCtx *ctx) {
603 auto snes_ctx = ctx->snesCtxImpl.get();
604 auto &m_field = snes_ctx->mField;
605 auto problem_ptr = m_field.get_problem(snes_ctx->problemName);
606 auto fields_ptr = m_field.get_fields();
607 auto dofs = problem_ptr->numeredRowDofsPtr;
608
609 std::vector<double> lnorms(fields_ptr->size(), 0),
610 norms(fields_ptr->size(), 0);
611
612 Vec res;
613 CHKERR SNESGetFunction(snes, &res, NULL, NULL);
614
615 const double *r;
616 CHKERR VecGetArrayRead(res, &r);
617 {
618 int f = 0;
619 for (auto fi : *fields_ptr) {
620 const auto lo_uid = FieldEntity::getLoBitNumberUId(fi->bitNumber);
621 const auto hi_uid = FieldEntity::getHiBitNumberUId(fi->bitNumber);
622 const auto hi = dofs->get<Unique_mi_tag>().upper_bound(hi_uid);
623 for (auto lo = dofs->get<Unique_mi_tag>().lower_bound(lo_uid); lo != hi;
624 ++lo) {
625 const DofIdx loc_idx = (*lo)->getPetscLocalDofIdx();
626 if (loc_idx >= 0 && loc_idx < problem_ptr->nbLocDofsRow) {
627 lnorms[f] += PetscRealPart(PetscSqr(r[loc_idx]));
628 }
629 }
630 ++f;
631 }
632 }
633 CHKERR VecRestoreArrayRead(res, &r);
634
635 MPIU_Allreduce(&*lnorms.begin(), &*norms.begin(), lnorms.size(), MPIU_REAL,
636 MPIU_SUM, PetscObjectComm((PetscObject)snes));
637
638 {
639 MOFEM_LOG("SNES_WORLD", Sev::inform)
640 << its << " Function norm " << boost::format("%14.12e") % (double)fgnorm
641 << " [";
642 int f = 0;
643 for (auto fi : *fields_ptr) {
644 MOFEM_LOG("SNES_WORLD", Sev::inform)
645 << its << "\t Field norm "
646 << boost::format("%14.12e") % (double)PetscSqrtReal(norms[f]) << " "
647 << fi->getName();
648 ++f;
649 }
650 MOFEM_LOG("SNES_WORLD", Sev::inform) << its << "]";
651 }
652
654}
655
656MoFEMErrorCode MoFEMSNESMonitorEnergy(SNES snes, PetscInt its, PetscReal fgnorm,
657 SnesCtx *ctx) {
659 auto snes_ctx = ctx->snesCtxImpl.get();
660 auto &m_field = snes_ctx->mField;
661 auto problem_ptr = m_field.get_problem(snes_ctx->problemName);
662 auto fields_ptr = m_field.get_fields();
663 auto dofs = problem_ptr->numeredRowDofsPtr;
664
665 std::vector<double> lnorms(3 * fields_ptr->size(), 0),
666 norms(3 * fields_ptr->size(), 0);
667
668 Vec x_update, res;
669 if (its == 0)
670 CHKERR SNESGetSolution(snes, &x_update);
671 else
672 CHKERR SNESGetSolutionUpdate(snes, &x_update);
673 CHKERR SNESGetFunction(snes, &res, PETSC_NULLPTR, PETSC_NULLPTR);
674
675 const double *x, *r;
676 ;
677 CHKERR VecGetArrayRead(res, &r);
678 CHKERR VecGetArrayRead(x_update, &x);
679
680 {
681 int f = 0;
682 for (auto fi : *fields_ptr) {
683 const auto lo_uid = FieldEntity::getLoBitNumberUId(fi->bitNumber);
684 const auto hi_uid = FieldEntity::getHiBitNumberUId(fi->bitNumber);
685 const auto hi = dofs->get<Unique_mi_tag>().upper_bound(hi_uid);
686 for (auto lo = dofs->get<Unique_mi_tag>().lower_bound(lo_uid); lo != hi;
687 ++lo) {
688 const DofIdx loc_idx = (*lo)->getPetscLocalDofIdx();
689 if (loc_idx >= 0 && loc_idx < problem_ptr->nbLocDofsRow) {
690 lnorms[3 * f] += std::abs(x[loc_idx] * r[loc_idx]);
691 lnorms[3 * f + 1] += r[loc_idx] * r[loc_idx];
692 lnorms[3 * f + 2] += x[loc_idx] * x[loc_idx];
693 }
694 }
695 ++f;
696 }
697 }
698 CHKERR VecRestoreArrayRead(res, &r);
699 CHKERR VecRestoreArrayRead(x_update, &x);
700
701 MPIU_Allreduce(&*lnorms.begin(), &*norms.begin(), lnorms.size(), MPIU_REAL,
702 MPIU_SUM, PetscObjectComm((PetscObject)snes));
703
704 double energy_norm = 0, residual_norm = 0, increment_norm = 0;
705 for (int f = 0; f != fields_ptr->size(); ++f) {
706 energy_norm += norms[3 * f];
707 residual_norm += norms[3 * f + 1];
708 increment_norm += norms[3 * f + 2];
709 }
710 energy_norm = std::sqrt(energy_norm);
711 residual_norm = std::sqrt(residual_norm);
712 increment_norm = std::sqrt(increment_norm);
713
714 {
715 MOFEM_LOG("SNES_WORLD", Sev::inform)
716 << its << " Function norm Energy/Residual/Increment"
717 << boost::format(" %14.12e ") % energy_norm
718 << boost::format(" %14.12e ") % residual_norm
719 << boost::format(" %14.12e ") % increment_norm << "[";
720 int f = 0;
721 for (auto fi : *fields_ptr) {
722 MOFEM_LOG("SNES_WORLD", Sev::inform)
723 << its << "\t Energy/Residual/Increment norm "
724 << boost::format("%14.12e") % std::sqrt(norms[3 * f]) << "\t"
725 << boost::format("%14.12e") % std::sqrt(norms[3 * f + 1]) << " "
726 << boost::format("%14.12e") % std::sqrt(norms[3 * f + 2]) << "\t "
727 << fi->getName();
728 ++f;
729 }
730 MOFEM_LOG("SNES_WORLD", Sev::inform) << its << " ]";
731 }
732
734}
735
736} // namespace MoFEM
constexpr double a
@ COL
@ ROW
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
@ MF_EXIST
#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()
#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.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ F
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
auto bit
set bit
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'l', 3 > l
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
int DofIdx
Index of DOF.
Definition Types.hpp:18
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
std::deque< BasicMethodPtr > BasicMethodsSequence
Definition AuxPETSc.hpp:58
std::deque< PairNameFEMethodPtr > FEMethodsSequence
Definition AuxPETSc.hpp:57
constexpr AssemblyType A
constexpr double g
Deprecated interface functions.
static UId getHiBitNumberUId(const FieldBitNumber bit_number)
static UId getLoBitNumberUId(const FieldBitNumber bit_number)
@ CTX_KSPNONE
No specific KSP context.
@ CTX_SETFUNCTION
Setting up the linear system function.
@ CTX_OPERATORS
Setting up linear operators.
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 bool checkIfChannelExist(const std::string channel)
Check if channel exist.
static constexpr Switches CtxSetA
Jacobian matrix switch.
static constexpr Switches CtxSetX
Solution vector switch.
static constexpr Switches CtxSetNone
No data switch.
static constexpr Switches CtxSetF
Residual vector switch.
static constexpr Switches CtxSetDX
Solution increment switch.
static constexpr Switches CtxSetB
Preconditioner matrix switch.
static constexpr Switches CtxSetTime
Time value switch.
keeps basic data about problem
BasicMethodsSequence preProcessLoadTangent
Definition SnesCtx.cpp:101
FEMethodsSequence & getComputeRhs()
Definition SnesCtx.cpp:22
MatAssemblyType typeOfAssembly
type of assembly at the end
Definition SnesCtx.cpp:117
friend PetscErrorCode SnesLoadTangent(SNES snes, Vec u, Vec F, void *ctx)
This function calls finite element pipeline to compute tangent of of load vector in arc length contro...
Definition SnesCtx.cpp:359
FEMethodsSequence & getLoadTangent()
Definition SnesCtx.cpp:34
BasicMethodsSequence & getPreProcLoadTangent()
Definition SnesCtx.cpp:36
friend MoFEMErrorCode SnesMoFEMSetBehavior(SNES snes, MoFEMTypes bh)
Set behavior if finite element in sequence does not exist.
Definition SnesCtx.cpp:586
PetscLogEvent MOFEM_EVENT_SnesMat
Log events to assemble tangent matrix.
Definition SnesCtx.cpp:79
BasicMethodsSequence preProcessOperator
Definition SnesCtx.cpp:96
friend MoFEMErrorCode MoFEMSNESMonitorEnergy(SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *snes_ctx)
Sens monitor printing residual field by field.
Definition SnesCtx.cpp:656
FEMethodsSequence loopsRhs
Definition SnesCtx.cpp:90
BasicMethodsSequence & getPostProcComputeRhs()
Definition SnesCtx.cpp:26
moab::Interface & moab
moab Interface
Definition SnesCtx.cpp:84
BasicMethodsSequence postProcessLoadTangent
Definition SnesCtx.cpp:109
MoFEM::BasicMethodsSequence BasicMethodsSequence
Definition SnesCtx.cpp:13
std::string problemName
problem name
Definition SnesCtx.cpp:86
friend MoFEMErrorCode SNESMoFEMSetAssmblyType(SNES snes, MatAssemblyType type)
virtual ~SnesCtxImpl()=default
MoFEMErrorCode copyLoops(const SnesCtxImpl &snes_ctx)
Definition SnesCtx.cpp:199
friend PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx)
This is MoFEM implementation for the right hand side (residual vector) evaluation in SNES solver.
Definition SnesCtx.cpp:227
friend PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx)
This is MoFEM implementation for the left hand side (tangent matrix) evaluation in SNES solver.
Definition SnesCtx.cpp:491
HookFunction & getRhsHook()
Definition SnesCtx.cpp:50
FEMethodsSequence loopsOperator
Definition SnesCtx.cpp:88
MoFEM::Interface & getMField()
Definition SnesCtx.cpp:52
SnesCtx::HookFunction HookFunction
Definition SnesCtx.cpp:14
BasicMethodsSequence & getPreProcComputeRhs()
Definition SnesCtx.cpp:24
FEMethodsSequence loopsLoadTangent
Definition SnesCtx.cpp:93
FEMethodsSequence & getSetOperators()
Definition SnesCtx.cpp:20
MoFEM::Interface & mField
database Interface
Definition SnesCtx.cpp:83
BasicMethodsSequence postProcessRhs
Sequence of methods run after residual is assembled.
Definition SnesCtx.cpp:107
BasicMethodsSequence & getPostProcSetOperators()
Definition SnesCtx.cpp:30
HookFunction & getLoadTangentHook()
Definition SnesCtx.cpp:48
MoFEMTypes bH
If set to MF_EXIST check if element exist, default MF_EXIST.
Definition SnesCtx.cpp:112
PetscLogEvent MOFEM_EVENT_SnesLoadTangent
Definition SnesCtx.cpp:80
SnesCtxImpl(Interface &m_field, std::string problem_name)
Definition SnesCtx.cpp:185
BasicMethodsSequence preProcessRhs
Sequence of methods run before residual is assembled.
Definition SnesCtx.cpp:105
MoFEMErrorCode clearLoops()
Definition SnesCtx.cpp:213
BasicMethodsSequence & getPostProcLoadTangent()
Definition SnesCtx.cpp:40
BasicMethodsSequence postProcessOperator
Definition SnesCtx.cpp:98
MoFEM::FEMethodsSequence FEMethodsSequence
Definition SnesCtx.cpp:12
PetscLogEvent MOFEM_EVENT_SnesRhs
Log events to assemble residual.
Definition SnesCtx.cpp:78
friend MoFEMErrorCode MoFEMSNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *snes_ctx)
Sens monitor printing residual field by field.
Definition SnesCtx.cpp:600
bool vErify
If true verify vector.
Definition SnesCtx.cpp:116
BasicMethodsSequence & getPreProcSetOperators()
Definition SnesCtx.cpp:28
friend MoFEMErrorCode SnesMoFEMSetAssemblyType(SNES snes, MatAssemblyType type)
Set assembly type at the end of SnesMat.
Definition SnesCtx.cpp:572
Interface for nonlinear (SNES) solver.
Definition SnesCtx.hpp:15
MoFEM::BasicMethodsSequence BasicMethodsSequence
Definition SnesCtx.hpp:19
friend PetscErrorCode SnesLoadTangent(SNES snes, Vec u, Vec F, void *ctx)
This function calls finite element pipeline to compute tangent of of load vector in arc length contro...
Definition SnesCtx.cpp:359
FEMethodsSequence & getSetOperators()
Definition SnesCtx.cpp:133
friend MoFEMErrorCode SnesMoFEMSetBehavior(SNES snes, MoFEMTypes bh)
Set behavior if finite element in sequence does not exist.
Definition SnesCtx.cpp:586
friend MoFEMErrorCode MoFEMSNESMonitorEnergy(SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *snes_ctx)
Sens monitor printing residual field by field.
Definition SnesCtx.cpp:656
SnesCtx(Interface &m_field, std::string problem_name)
Definition SnesCtx.cpp:127
BasicMethodsSequence & getPreProcLoadTangent()
Get the BasicMethod sequence for preprocessing of FunctionFn.
Definition SnesCtx.cpp:161
BasicMethodsSequence & getPreProcSetOperators()
Definition SnesCtx.cpp:149
MoFEM::Interface & getMField()
Get mField reference.
Definition SnesCtx.cpp:175
MoFEM::FEMethodsSequence FEMethodsSequence
Definition SnesCtx.hpp:18
BasicMethodsSequence & getPostProcLoadTangent()
Get the BasicMethod sequence for postprocessing of FunctionFn.
Definition SnesCtx.cpp:165
MoFEMErrorCode copyLoops(const SnesCtx &snes_ctx)
Copy sequences from other SNES context.
Definition SnesCtx.cpp:169
boost::function< MoFEMErrorCode(SNES snes, Vec x, Vec F, boost::shared_ptr< CacheTuple > cache_ptr, void *ctx)> HookFunction
Definition SnesCtx.hpp:22
friend PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx)
This is MoFEM implementation for the right hand side (residual vector) evaluation in SNES solver.
Definition SnesCtx.cpp:227
friend PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx)
This is MoFEM implementation for the left hand side (tangent matrix) evaluation in SNES solver.
Definition SnesCtx.cpp:491
BasicMethodsSequence & getPostProcComputeRhs()
Definition SnesCtx.cpp:145
HookFunction & getLoadTangentHook()
Get the Load Tangent Hook function.
Definition SnesCtx.cpp:177
BasicMethodsSequence & getPostProcSetOperators()
Definition SnesCtx.cpp:153
boost::movelib::unique_ptr< SnesCtxImpl > snesCtxImpl
Definition SnesCtx.hpp:155
BasicMethodsSequence & getPreProcComputeRhs()
Definition SnesCtx.cpp:141
FEMethodsSequence & getLoadTangent()
Get the finite element pipeline for FunctionFn, that calculate tangent of function load for arc lengt...
Definition SnesCtx.cpp:157
FEMethodsSequence & getComputeRhs()
Definition SnesCtx.cpp:137
HookFunction & getRhsHook()
Get the Right Hand Side Hook function.
Definition SnesCtx.cpp:181
virtual ~SnesCtx()
friend MoFEMErrorCode MoFEMSNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *snes_ctx)
Sens monitor printing residual field by field.
Definition SnesCtx.cpp:600
MoFEMErrorCode clearLoops()
Clear loops.
Definition SnesCtx.cpp:173
friend MoFEMErrorCode SnesMoFEMSetAssemblyType(SNES snes, MatAssemblyType type)
Set assembly type at the end of SnesMat.
Definition SnesCtx.cpp:572
@ CTX_SNESSETJACOBIAN
Setting up Jacobian matrix computation.
@ CTX_SNESSETFUNCTION
Setting up nonlinear function evaluation.
@ CTX_SNESNONE
No specific SNES context.
Auxiliary tools.
Definition Tools.hpp:19
Vector manager is used to create vectors \mofem_vectors.