v0.15.5
Loading...
Searching...
No Matches
Classes | Public Types | Public Member Functions | Protected Attributes | Friends | List of all members
MoFEM::SnesCtx Struct Reference

Interface for nonlinear (SNES) solver. More...

#include "src/petsc/SnesCtx.hpp"

Collaboration diagram for MoFEM::SnesCtx:
[legend]

Classes

struct  SnesCtxImpl
 

Public Types

using PairNameFEMethodPtr = MoFEM::PairNameFEMethodPtr
 
using FEMethodsSequence = MoFEM::FEMethodsSequence
 
using BasicMethodsSequence = MoFEM::BasicMethodsSequence
 
using HookFunction = boost::function< MoFEMErrorCode(SNES snes, Vec x, Vec F, boost::shared_ptr< CacheTuple > cache_ptr, void *ctx)>
 

Public Member Functions

 SnesCtx (Interface &m_field, std::string problem_name)
 
virtual ~SnesCtx ()
 
FEMethodsSequencegetSetOperators ()
 
FEMethodsSequencegetComputeRhs ()
 
BasicMethodsSequencegetPreProcComputeRhs ()
 
BasicMethodsSequencegetPostProcComputeRhs ()
 
BasicMethodsSequencegetPreProcSetOperators ()
 
BasicMethodsSequencegetPostProcSetOperators ()
 
FEMethodsSequencegetLoadTangent ()
 Get the finite element pipeline for FunctionFn, that calculate tangent of function load for arc length control.
 
BasicMethodsSequencegetPreProcLoadTangent ()
 Get the BasicMethod sequence for preprocessing of FunctionFn.
 
BasicMethodsSequencegetPostProcLoadTangent ()
 Get the BasicMethod sequence for postprocessing of FunctionFn.
 
HookFunctiongetLoadTangentHook ()
 Get the Load Tangent Hook function.
 
HookFunctiongetRhsHook ()
 Get the Right Hand Side Hook function.
 
MoFEMErrorCode copyLoops (const SnesCtx &snes_ctx)
 Copy sequences from other SNES context.
 
MoFEMErrorCode clearLoops ()
 Clear loops.
 
MoFEM::InterfacegetMField ()
 Get mField reference.
 

Protected Attributes

boost::movelib::unique_ptr< SnesCtxImplsnesCtxImpl
 

Friends

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.
 
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.
 
MoFEMErrorCode SNESMoFEMSetAssmblyType (SNES snes, MatAssemblyType type)
 
MoFEMErrorCode SnesMoFEMSetBehavior (SNES snes, MoFEMTypes bh)
 Set behavior if finite element in sequence does not exist.
 
MoFEMErrorCode SnesMoFEMSetAssemblyType (SNES snes, MatAssemblyType type)
 Set assembly type at the end of SnesMat.
 
MoFEMErrorCode MoFEMSNESMonitorFields (SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *snes_ctx)
 Sens monitor printing residual field by field.
 
MoFEMErrorCode MoFEMSNESMonitorEnergy (SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *snes_ctx)
 Sens monitor printing residual field by field.
 
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 control.
 

Detailed Description

Interface for nonlinear (SNES) solver.

Examples
mofem/tutorials/cor-10_navier_stokes/navier_stokes.cpp, mofem/users_modules/basic_finite_elements/atom_tests/testing_jacobian_of_hook_element.cpp, mofem/users_modules/basic_finite_elements/atom_tests/testing_jacobian_of_hook_scaled_with_density_element.cpp, mofem/users_modules/basic_finite_elements/nonlinear_elasticity/nonlinear_dynamics.cpp, and mofem/users_modules/minimal_surface_equation/minimal_surface_area.cpp.

Definition at line 15 of file SnesCtx.hpp.

Member Typedef Documentation

◆ BasicMethodsSequence

Definition at line 19 of file SnesCtx.hpp.

◆ FEMethodsSequence

Definition at line 18 of file SnesCtx.hpp.

◆ HookFunction

using MoFEM::SnesCtx::HookFunction = boost::function<MoFEMErrorCode( SNES snes, Vec x, Vec F, boost::shared_ptr<CacheTuple> cache_ptr, void *ctx)>

Definition at line 20 of file SnesCtx.hpp.

◆ PairNameFEMethodPtr

Definition at line 17 of file SnesCtx.hpp.

Constructor & Destructor Documentation

◆ SnesCtx()

MoFEM::SnesCtx::SnesCtx ( Interface m_field,
std::string  problem_name 
)

Definition at line 127 of file SnesCtx.cpp.

128 : snesCtxImpl(boost::movelib::make_unique<SnesCtx::SnesCtxImpl>(
129 m_field, problem_name)) {}
boost::movelib::unique_ptr< SnesCtxImpl > snesCtxImpl
Definition SnesCtx.hpp:155

◆ ~SnesCtx()

MoFEM::SnesCtx::~SnesCtx ( )
virtualdefault

Member Function Documentation

◆ clearLoops()

MoFEMErrorCode MoFEM::SnesCtx::clearLoops ( )

Clear loops.

Returns
MoFEMErrorCode

Definition at line 173 of file SnesCtx.cpp.

173{ return snesCtxImpl->clearLoops(); }

◆ copyLoops()

MoFEMErrorCode MoFEM::SnesCtx::copyLoops ( const SnesCtx snes_ctx)

Copy sequences from other SNES context.

Parameters
snes_ctxSNES context from which Sequence is copied from
Returns
error code

Definition at line 169 of file SnesCtx.cpp.

169 {
170 return snesCtxImpl->copyLoops(*snes_ctx.snesCtxImpl);
171}

◆ getComputeRhs()

FEMethodsSequence & MoFEM::SnesCtx::getComputeRhs ( )
Returns
return vector to vector with FEMethod to calculate residual

Definition at line 137 of file SnesCtx.cpp.

137 {
138 return snesCtxImpl->getComputeRhs();
139}

◆ getLoadTangent()

FEMethodsSequence & MoFEM::SnesCtx::getLoadTangent ( )

Get the finite element pipeline for FunctionFn, that calculate tangent of function load for arc length control.

Returns
FEMethodsSequence&

Definition at line 157 of file SnesCtx.cpp.

157 {
158 return snesCtxImpl->getLoadTangent();
159}

◆ getLoadTangentHook()

SnesCtx::HookFunction & MoFEM::SnesCtx::getLoadTangentHook ( )

Get the Load Tangent Hook function.

Using hook you can debug or post-process while load vector tangent is evaluated

Returns
HookFunction&

Definition at line 177 of file SnesCtx.cpp.

177 {
178 return snesCtxImpl->getLoadTangentHook();
179}

◆ getMField()

MoFEM::Interface & MoFEM::SnesCtx::getMField ( )

Get mField reference.

Definition at line 175 of file SnesCtx.cpp.

175{ return snesCtxImpl->getMField(); }

◆ getPostProcComputeRhs()

BasicMethodsSequence & MoFEM::SnesCtx::getPostProcComputeRhs ( )

The sequence of BasicMethod is executed after residual is calculated. It can be used to setup data structures, e.g. aggregate data from processors or to apply essential boundary conditions.

Returns
reference to BasicMethod for postprocessing

Definition at line 145 of file SnesCtx.cpp.

145 {
146 return snesCtxImpl->getPostProcComputeRhs();
147}

◆ getPostProcLoadTangent()

BasicMethodsSequence & MoFEM::SnesCtx::getPostProcLoadTangent ( )

Get the BasicMethod sequence for postprocessing of FunctionFn.

Returns
BasicMethodsSequence&

Definition at line 165 of file SnesCtx.cpp.

165 {
166 return snesCtxImpl->getPostProcLoadTangent();
167}

◆ getPostProcSetOperators()

BasicMethodsSequence & MoFEM::SnesCtx::getPostProcSetOperators ( )

The sequence of BasicMethod is executed after tangent matrix is calculated. It can be used to setup data structures, e.g. aggregate data from processors or to apply essential boundary conditions.

Returns
reference to BasicMethod for postprocessing

Definition at line 153 of file SnesCtx.cpp.

153 {
154 return snesCtxImpl->getPostProcSetOperators();
155}

◆ getPreProcComputeRhs()

BasicMethodsSequence & MoFEM::SnesCtx::getPreProcComputeRhs ( )

The sequence of BasicMethod is executed before residual is calculated. It can be used to setup data structures, e.g. zero global variable which is integrated in domain, e.g. for calculation of strain energy.

Returns
reference to BasicMethod for preprocessing

Definition at line 141 of file SnesCtx.cpp.

141 {
142 return snesCtxImpl->getPreProcComputeRhs();
143}

◆ getPreProcLoadTangent()

BasicMethodsSequence & MoFEM::SnesCtx::getPreProcLoadTangent ( )

Get the BasicMethod sequence for preprocessing of FunctionFn.

Returns
BasicMethodsSequence&

Definition at line 161 of file SnesCtx.cpp.

161 {
162 return snesCtxImpl->getPreProcLoadTangent();
163}

◆ getPreProcSetOperators()

BasicMethodsSequence & MoFEM::SnesCtx::getPreProcSetOperators ( )
Returns
reference to BasicMethod for preprocessing

Definition at line 149 of file SnesCtx.cpp.

149 {
150 return snesCtxImpl->getPreProcSetOperators();
151}

◆ getRhsHook()

SnesCtx::HookFunction & MoFEM::SnesCtx::getRhsHook ( )

Get the Right Hand Side Hook function.

Using hook you can debug or post-process while right hand side vector is evaluated

Returns
HookFunction&

Definition at line 181 of file SnesCtx.cpp.

181 {
182 return snesCtxImpl->getRhsHook();
183}

◆ getSetOperators()

FEMethodsSequence & MoFEM::SnesCtx::getSetOperators ( )
Returns
return reference to vector with FEMethod to calculate tangent matrix

Definition at line 133 of file SnesCtx.cpp.

133 {
134 return snesCtxImpl->getSetOperators();
135}

Friends And Related Symbol Documentation

◆ MoFEMSNESMonitorEnergy

MoFEMErrorCode MoFEMSNESMonitorEnergy ( SNES  snes,
PetscInt  its,
PetscReal  fgnorm,
SnesCtx snes_ctx 
)
friend

Sens monitor printing residual field by field.

Definition at line 656 of file SnesCtx.cpp.

657 {
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}
#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 MOFEM_LOG(channel, severity)
Log.
const FTensor::Tensor2< T, Dim, Dim > Vec
int DofIdx
Index of DOF.
Definition Types.hpp:18
int r
Definition sdf.py:205
static UId getHiBitNumberUId(const FieldBitNumber bit_number)
static UId getLoBitNumberUId(const FieldBitNumber bit_number)

◆ MoFEMSNESMonitorFields

MoFEMErrorCode MoFEMSNESMonitorFields ( SNES  snes,
PetscInt  its,
PetscReal  fgnorm,
SnesCtx snes_ctx 
)
friend

Sens monitor printing residual field by field.

Definition at line 600 of file SnesCtx.cpp.

601 {
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}

◆ SnesLoadTangent

PetscErrorCode SnesLoadTangent ( SNES  snes,
Vec  u,
Vec  F,
void *  ctx 
)
friend

This function calls finite element pipeline to compute tangent of of load vector in arc length control.

For more information pleas look to PETSc manual, i.e. SNESSetFunction https://petsc.org/main/manualpages/SNES/SNESSetFunction/

Parameters
snesSNES solver
uSolution vector at current iteration
FThe right hand side vector
ctxPointer to context i.e. SnesCtx
Returns
Error code

Definition at line 359 of file SnesCtx.cpp.

359 {
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}
constexpr double a
@ COL
@ ROW
#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 ...
auto bit
set bit
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'l', 3 > l
constexpr double g
@ CTX_KSPNONE
No specific KSP context.
@ CTX_SETFUNCTION
Setting up the linear system function.
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 CtxSetTime
Time value switch.
SnesCtx(Interface &m_field, std::string problem_name)
Definition SnesCtx.cpp:127
@ CTX_SNESSETFUNCTION
Setting up nonlinear function evaluation.
@ CTX_SNESNONE
No specific SNES context.

◆ SnesMat

PetscErrorCode SnesMat ( SNES  snes,
Vec  x,
Mat  A,
Mat  B,
void *  ctx 
)
friend

This is MoFEM implementation for the left hand side (tangent matrix) evaluation in SNES solver.

For more information pleas look to PETSc manual, i.e. SNESSetJacobian http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetJacobian.html#SNESSetJacobian

Parameters
snesSNES solver
xSolution vector at current iteration
ATangent matrix
BPreconditioner tangent matrix
ctxPointer to context i.e. SnesCtx
Returns
Error code

Definition at line 491 of file SnesCtx.cpp.

491 {
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}
constexpr AssemblyType A
@ CTX_OPERATORS
Setting up linear operators.
static constexpr Switches CtxSetA
Jacobian matrix switch.
static constexpr Switches CtxSetB
Preconditioner matrix switch.
@ CTX_SNESSETJACOBIAN
Setting up Jacobian matrix computation.

◆ SnesMoFEMSetAssemblyType

MoFEMErrorCode SnesMoFEMSetAssemblyType ( SNES  snes,
MatAssemblyType  type 
)
friend

Set assembly type at the end of SnesMat.

Note
Note that tangent matrix need have to have final assembly, you would use flush assembly in special case that you call SnesMat form other function set to SNESSetJacobian
Parameters
snes
typetype of assembly, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY
Returns
error code

Definition at line 572 of file SnesCtx.cpp.

572 {
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}
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.

◆ SNESMoFEMSetAssmblyType

MoFEMErrorCode SNESMoFEMSetAssmblyType ( SNES  snes,
MatAssemblyType  type 
)
friend

◆ SnesMoFEMSetBehavior

MoFEMErrorCode SnesMoFEMSetBehavior ( SNES  snes,
MoFEMTypes  bh 
)
friend

Set behavior if finite element in sequence does not exist.

Parameters
snes
bhIf set to MF_EXIST check if element exist, default MF_EXIST. Otherwise set MF_ZERO
Returns
error code

Definition at line 586 of file SnesCtx.cpp.

586 {
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}

◆ SnesRhs

PetscErrorCode SnesRhs ( SNES  snes,
Vec  x,
Vec  f,
void *  ctx 
)
friend

This is MoFEM implementation for the right hand side (residual vector) evaluation in SNES solver.

For more information pleas look to PETSc manual, i.e. SNESSetFunction http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFunction.html

Parameters
snesSNES solver
xSolution vector at current iteration
fThe right hand side vector
ctxPointer to context i.e. SnesCtx
Returns
Error code

Definition at line 227 of file SnesCtx.cpp.

227 {
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}

Member Data Documentation

◆ snesCtxImpl

boost::movelib::unique_ptr<SnesCtxImpl> MoFEM::SnesCtx::snesCtxImpl
protected

Definition at line 155 of file SnesCtx.hpp.


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