v0.13.0
Public Types | Public Member Functions | Public Attributes | Private 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]

Public Types

typedef MoFEM::PairNameFEMethodPtr PairNameFEMethodPtr
 
typedef MoFEM::FEMethodsSequence FEMethodsSequence
 
typedef MoFEM::BasicMethodsSequence BasicMethodsSequence
 

Public Member Functions

 SnesCtx (Interface &m_field, const std::string &problem_name)
 
virtual ~SnesCtx ()=default
 
FEMethodsSequenceget_loops_to_do_Mat ()
 
FEMethodsSequenceget_loops_to_do_Rhs ()
 
BasicMethodsSequenceget_preProcess_to_do_Rhs ()
 
BasicMethodsSequenceget_postProcess_to_do_Rhs ()
 
BasicMethodsSequenceget_preProcess_to_do_Mat ()
 
BasicMethodsSequenceget_postProcess_to_do_Mat ()
 
MoFEMErrorCode copyLoops (const SnesCtx &snes_ctx)
 Copy sequences from other SNES contex. More...
 
MoFEMErrorCode clearLoops ()
 Clear loops. More...
 

Public Attributes

MoFEM::InterfacemField
 database Interface More...
 
moab::Interface & moab
 moab Interface More...
 
std::string problemName
 problem name More...
 
MoFEMTypes bH
 If set to MF_EXIST check if element exist, default MF_EXIST. More...
 
bool zeroPreCondMatrixB
 
MatAssemblyType typeOfAssembly
 type of assembly at the end More...
 
bool vErify
 If true verify vector. More...
 
FEMethodsSequence loops_to_do_Mat
 
FEMethodsSequence loops_to_do_Rhs
 
BasicMethodsSequence preProcess_Mat
 
BasicMethodsSequence postProcess_Mat
 
BasicMethodsSequence preProcess_Rhs
 Sequence of methods run before residual is assembled. More...
 
BasicMethodsSequence postProcess_Rhs
 Sequence of methods run after residual is assembled. More...
 

Private Attributes

boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
PetscLogEvent MOFEM_EVENT_SnesRhs
 Log events to assemble residual. More...
 
PetscLogEvent MOFEM_EVENT_SnesMat
 Log events to assemble tangent matrix. More...
 

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. More...
 
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. More...
 
MoFEMErrorCode SNESMoFEMSetAssmblyType (SNES snes, MatAssemblyType type)
 
MoFEMErrorCode SnesMoFEMSetBehavior (SNES snes, MoFEMTypes bh)
 Set behavior if finite element in sequence does not exist. More...
 

Detailed Description

Interface for nonlinear (SNES) solver.

Examples
testing_jacobian_of_hook_element.cpp, and testing_jacobian_of_hook_scaled_with_density_element.cpp.

Definition at line 27 of file SnesCtx.hpp.

Member Typedef Documentation

◆ BasicMethodsSequence

Definition at line 42 of file SnesCtx.hpp.

◆ FEMethodsSequence

Examples
nonlinear_dynamics.cpp.

Definition at line 41 of file SnesCtx.hpp.

◆ PairNameFEMethodPtr

Examples
nonlinear_dynamics.cpp.

Definition at line 40 of file SnesCtx.hpp.

Constructor & Destructor Documentation

◆ SnesCtx()

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

Definition at line 58 of file SnesCtx.hpp.

59  : mField(m_field), moab(m_field.get_moab()), problemName(problem_name),
61  typeOfAssembly(MAT_FINAL_ASSEMBLY), vErify(false) {
62  PetscLogEventRegister("LoopSNESRhs", 0, &MOFEM_EVENT_SnesRhs);
63  PetscLogEventRegister("LoopSNESMat", 0, &MOFEM_EVENT_SnesMat);
64  }
@ MF_EXIST
Definition: definitions.h:113
PetscLogEvent MOFEM_EVENT_SnesMat
Log events to assemble tangent matrix.
Definition: SnesCtx.hpp:137
MoFEMTypes bH
If set to MF_EXIST check if element exist, default MF_EXIST.
Definition: SnesCtx.hpp:34
moab::Interface & moab
moab Interface
Definition: SnesCtx.hpp:30
MatAssemblyType typeOfAssembly
type of assembly at the end
Definition: SnesCtx.hpp:37
bool zeroPreCondMatrixB
Definition: SnesCtx.hpp:35
bool vErify
If true verify vector.
Definition: SnesCtx.hpp:38
std::string problemName
problem name
Definition: SnesCtx.hpp:32
MoFEM::Interface & mField
database Interface
Definition: SnesCtx.hpp:29
PetscLogEvent MOFEM_EVENT_SnesRhs
Log events to assemble residual.
Definition: SnesCtx.hpp:136

◆ ~SnesCtx()

virtual MoFEM::SnesCtx::~SnesCtx ( )
virtualdefault

Member Function Documentation

◆ clearLoops()

MoFEMErrorCode SnesCtx::clearLoops ( )

Clear loops.

Returns
MoFEMErrorCode

Definition at line 28 of file SnesCtx.cpp.

28  {
30  loops_to_do_Mat.clear();
31  loops_to_do_Rhs.clear();
32  preProcess_Mat.clear();
33  postProcess_Mat.clear();
34  preProcess_Rhs.clear();
35  postProcess_Rhs.clear();
37 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
BasicMethodsSequence preProcess_Rhs
Sequence of methods run before residual is assembled.
Definition: SnesCtx.hpp:53
BasicMethodsSequence postProcess_Mat
Definition: SnesCtx.hpp:50
BasicMethodsSequence postProcess_Rhs
Sequence of methods run after residual is assembled.
Definition: SnesCtx.hpp:55
BasicMethodsSequence preProcess_Mat
Definition: SnesCtx.hpp:48
FEMethodsSequence loops_to_do_Mat
Definition: SnesCtx.hpp:44
FEMethodsSequence loops_to_do_Rhs
Definition: SnesCtx.hpp:46

◆ copyLoops()

MoFEMErrorCode SnesCtx::copyLoops ( const SnesCtx snes_ctx)

Copy sequences from other SNES contex.

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

Definition at line 17 of file SnesCtx.cpp.

17  {
19  loops_to_do_Mat = snes_ctx.loops_to_do_Mat;
20  loops_to_do_Rhs = snes_ctx.loops_to_do_Rhs;
21  preProcess_Mat = snes_ctx.preProcess_Mat;
22  postProcess_Mat = snes_ctx.postProcess_Mat;
23  preProcess_Rhs = snes_ctx.preProcess_Rhs;
24  postProcess_Rhs = snes_ctx.postProcess_Rhs;
26 }

◆ get_loops_to_do_Mat()

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

Definition at line 72 of file SnesCtx.hpp.

72 { return loops_to_do_Mat; }

◆ get_loops_to_do_Rhs()

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

Definition at line 77 of file SnesCtx.hpp.

77 { return loops_to_do_Rhs; }

◆ get_postProcess_to_do_Mat()

BasicMethodsSequence& MoFEM::SnesCtx::get_postProcess_to_do_Mat ( )

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 109 of file SnesCtx.hpp.

109 { return postProcess_Mat; }

◆ get_postProcess_to_do_Rhs()

BasicMethodsSequence& MoFEM::SnesCtx::get_postProcess_to_do_Rhs ( )

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 95 of file SnesCtx.hpp.

95 { return postProcess_Rhs; }

◆ get_preProcess_to_do_Mat()

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

Definition at line 100 of file SnesCtx.hpp.

100 { return preProcess_Mat; }

◆ get_preProcess_to_do_Rhs()

BasicMethodsSequence& MoFEM::SnesCtx::get_preProcess_to_do_Rhs ( )

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 86 of file SnesCtx.hpp.

86 { return preProcess_Rhs; }

Friends And Related Function Documentation

◆ 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 148 of file SnesCtx.cpp.

148  {
149  SnesCtx *snes_ctx = (SnesCtx *)ctx;
150  // PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
152  PetscLogEventBegin(snes_ctx->MOFEM_EVENT_SnesMat, 0, 0, 0, 0);
153  if (snes_ctx->zeroPreCondMatrixB)
154  CHKERR MatZeroEntries(B);
155 
156  snes_ctx->matAssembleSwitch = boost::movelib::make_unique<bool>(true);
157 
158  auto set = [&](auto &fe) {
159  fe.snes = snes;
160  fe.snes_x = x;
161  fe.snes_A = A;
162  fe.snes_B = B;
163  fe.snes_ctx = SnesMethod::CTX_SNESSETJACOBIAN;
164  fe.ksp_ctx = KspMethod::CTX_OPERATORS;
166  };
167 
168  auto unset = [&](auto &fe) {
169  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
170  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
171  fe.data_ctx = PetscData::CtxSetNone;
172  };
173 
174 
175  CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
176  CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
177  CHKERR snes_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
178  snes_ctx->problemName, COL, x, INSERT_VALUES, SCATTER_REVERSE);
179  for (auto &bit : snes_ctx->preProcess_Mat) {
180  bit->matAssembleSwitch = boost::move(snes_ctx->matAssembleSwitch);
181  set(*bit);
182  CHKERR snes_ctx->mField.problem_basic_method_preProcess(
183  snes_ctx->problemName, *bit);
184  unset(*bit);
185  snes_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
186  }
187 
188  auto cache_ptr = boost::make_shared<CacheTuple>();
189  CHKERR snes_ctx->mField.cache_problem_entities(snes_ctx->problemName,
190  cache_ptr);
191 
192  for (auto &lit : snes_ctx->loops_to_do_Mat) {
193  lit.second->matAssembleSwitch = boost::move(snes_ctx->matAssembleSwitch);
194  set(*lit.second);
195  CHKERR snes_ctx->mField.loop_finite_elements(
196  snes_ctx->problemName, lit.first, *(lit.second), nullptr, snes_ctx->bH,
197  cache_ptr);
198  unset(*lit.second);
199  snes_ctx->matAssembleSwitch = boost::move(lit.second->matAssembleSwitch);
200  }
201 
202  for (auto &bit : snes_ctx->postProcess_Mat) {
203  bit->matAssembleSwitch = boost::move(snes_ctx->matAssembleSwitch);
204  set(*bit);
205  CHKERR snes_ctx->mField.problem_basic_method_postProcess(
206  snes_ctx->problemName, *bit);
207  unset(*bit);
208  snes_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
209  }
210 
211  if (*snes_ctx->matAssembleSwitch) {
212  CHKERR MatAssemblyBegin(B, snes_ctx->typeOfAssembly);
213  CHKERR MatAssemblyEnd(B, snes_ctx->typeOfAssembly);
214  }
215  PetscLogEventEnd(snes_ctx->MOFEM_EVENT_SnesMat, 0, 0, 0, 0);
217 }
@ COL
Definition: definitions.h:136
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
#define CHKERR
Inline error check.
Definition: definitions.h:548
auto bit
set bit
double A
static constexpr Switches CtxSetA
Definition: LoopMethods.hpp:47
static constexpr Switches CtxSetX
Definition: LoopMethods.hpp:49
static constexpr Switches CtxSetNone
Definition: LoopMethods.hpp:45
static constexpr Switches CtxSetB
Definition: LoopMethods.hpp:48

◆ 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 228 of file SnesCtx.cpp.

228  {
229  SnesCtx *snes_ctx;
231  CHKERR SNESGetApplicationContext(snes, &snes_ctx);
232  snes_ctx->bH = bh;
234 }

◆ 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 39 of file SnesCtx.cpp.

39  {
40  SnesCtx *snes_ctx = (SnesCtx *)ctx;
41  // PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
43  PetscLogEventBegin(snes_ctx->MOFEM_EVENT_SnesRhs, 0, 0, 0, 0);
44  CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
45  CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
46  if (snes_ctx->vErify) {
47  // Verify finite elements, check for not a number
48  CHKERR VecAssemblyBegin(f);
49  CHKERR VecAssemblyEnd(f);
50  MPI_Comm comm = PetscObjectComm((PetscObject)f);
51  PetscSynchronizedPrintf(comm, "SNES Verify x\n");
52  const Problem *prb_ptr;
53  CHKERR snes_ctx->mField.get_problem(snes_ctx->problemName, &prb_ptr);
54  CHKERR snes_ctx->mField.getInterface<Tools>()->checkVectorForNotANumber(
55  prb_ptr, COL, x);
56  }
57  CHKERR snes_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
58  snes_ctx->problemName, COL, x, INSERT_VALUES, SCATTER_REVERSE);
59 
60  auto zero_ghost_vec = [](Vec g) {
62  Vec l;
63  CHKERR VecGhostGetLocalForm(g, &l);
64  double *a;
65  CHKERR VecGetArray(l, &a);
66  int s;
67  CHKERR VecGetLocalSize(l, &s);
68  for (int i = 0; i != s; ++i)
69  a[i] = 0;
70  CHKERR VecRestoreArray(l, &a);
71  CHKERR VecGhostRestoreLocalForm(g, &l);
73  };
74  CHKERR zero_ghost_vec(f);
75 
76  snes_ctx->vecAssembleSwitch = boost::movelib::make_unique<bool>(true);
77 
78  auto set = [&](auto &fe) {
79  fe.snes = snes;
80  fe.snes_x = x;
81  fe.snes_f = f;
82  fe.snes_ctx = SnesMethod::CTX_SNESSETFUNCTION;
83  fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
84  fe.data_ctx = PetscData::CtxSetF | PetscData::CtxSetX;
85  };
86 
87  auto unset = [&](auto &fe) {
88  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
89  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
90  fe.data_ctx = PetscData::CtxSetNone;
91  };
92 
93  for (auto &bit : snes_ctx->preProcess_Rhs) {
94  bit->vecAssembleSwitch = boost::move(snes_ctx->vecAssembleSwitch);
95  set(*bit);
96  CHKERR snes_ctx->mField.problem_basic_method_preProcess(
97  snes_ctx->problemName, *bit);
98  unset(*bit);
99  snes_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
100  }
101 
102  auto cache_ptr = boost::make_shared<CacheTuple>();
103  CHKERR snes_ctx->mField.cache_problem_entities(snes_ctx->problemName,
104  cache_ptr);
105 
106  for (auto &lit : snes_ctx->loops_to_do_Rhs) {
107  lit.second->vecAssembleSwitch = boost::move(snes_ctx->vecAssembleSwitch);
108  set(*lit.second);
109  CHKERR snes_ctx->mField.loop_finite_elements(
110  snes_ctx->problemName, lit.first, *lit.second, nullptr, snes_ctx->bH,
111  cache_ptr);
112  unset(*lit.second);
113  if (snes_ctx->vErify) {
114  // Verify finite elements, check for not a number
115  CHKERR VecAssemblyBegin(f);
116  CHKERR VecAssemblyEnd(f);
117  MPI_Comm comm = PetscObjectComm((PetscObject)f);
118  PetscSynchronizedPrintf(comm, "SNES Verify f FE < %s >\n",
119  lit.first.c_str());
120  const Problem *prb_ptr;
121  CHKERR snes_ctx->mField.get_problem(snes_ctx->problemName, &prb_ptr);
122  CHKERR snes_ctx->mField.getInterface<Tools>()->checkVectorForNotANumber(
123  prb_ptr, ROW, f);
124  }
125 
126  snes_ctx->vecAssembleSwitch = boost::move(lit.second->vecAssembleSwitch);
127  }
128 
129  for (auto &bit : snes_ctx->postProcess_Rhs) {
130  bit->vecAssembleSwitch = boost::move(snes_ctx->vecAssembleSwitch);
131  set(*bit);
132  CHKERR snes_ctx->mField.problem_basic_method_postProcess(
133  snes_ctx->problemName, *bit);
134  unset(*bit);
135  snes_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
136  }
137 
138  if (snes_ctx->vecAssembleSwitch) {
139  CHKERR VecGhostUpdateBegin(f, ADD_VALUES, SCATTER_REVERSE);
140  CHKERR VecGhostUpdateEnd(f, ADD_VALUES, SCATTER_REVERSE);
141  CHKERR VecAssemblyBegin(f);
142  CHKERR VecAssemblyEnd(f);
143  }
144  PetscLogEventEnd(snes_ctx->MOFEM_EVENT_SnesRhs, 0, 0, 0, 0);
146 }
constexpr double a
@ ROW
Definition: definitions.h:136
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'l', 3 > l
const FTensor::Tensor2< T, Dim, Dim > Vec
constexpr double g
static constexpr Switches CtxSetF
Definition: LoopMethods.hpp:46

Member Data Documentation

◆ bH

MoFEMTypes MoFEM::SnesCtx::bH

If set to MF_EXIST check if element exist, default MF_EXIST.

Definition at line 34 of file SnesCtx.hpp.

◆ loops_to_do_Mat

FEMethodsSequence MoFEM::SnesCtx::loops_to_do_Mat

Sequence of finite elements instances assembling tangent matrix

Definition at line 44 of file SnesCtx.hpp.

◆ loops_to_do_Rhs

FEMethodsSequence MoFEM::SnesCtx::loops_to_do_Rhs

Sequence of finite elements instances assembling residual vector

Definition at line 46 of file SnesCtx.hpp.

◆ matAssembleSwitch

boost::movelib::unique_ptr<bool> MoFEM::SnesCtx::matAssembleSwitch
private

Definition at line 135 of file SnesCtx.hpp.

◆ mField

MoFEM::Interface& MoFEM::SnesCtx::mField

database Interface

Definition at line 29 of file SnesCtx.hpp.

◆ moab

moab::Interface& MoFEM::SnesCtx::moab

moab Interface

Definition at line 30 of file SnesCtx.hpp.

◆ MOFEM_EVENT_SnesMat

PetscLogEvent MoFEM::SnesCtx::MOFEM_EVENT_SnesMat
private

Log events to assemble tangent matrix.

Definition at line 137 of file SnesCtx.hpp.

◆ MOFEM_EVENT_SnesRhs

PetscLogEvent MoFEM::SnesCtx::MOFEM_EVENT_SnesRhs
private

Log events to assemble residual.

Definition at line 136 of file SnesCtx.hpp.

◆ postProcess_Mat

BasicMethodsSequence MoFEM::SnesCtx::postProcess_Mat

Sequence of methods run after tangent matrix is assembled

Definition at line 50 of file SnesCtx.hpp.

◆ postProcess_Rhs

BasicMethodsSequence MoFEM::SnesCtx::postProcess_Rhs

Sequence of methods run after residual is assembled.

Definition at line 55 of file SnesCtx.hpp.

◆ preProcess_Mat

BasicMethodsSequence MoFEM::SnesCtx::preProcess_Mat

Sequence of methods run before tangent matrix is assembled

Definition at line 48 of file SnesCtx.hpp.

◆ preProcess_Rhs

BasicMethodsSequence MoFEM::SnesCtx::preProcess_Rhs

Sequence of methods run before residual is assembled.

Definition at line 53 of file SnesCtx.hpp.

◆ problemName

std::string MoFEM::SnesCtx::problemName

problem name

Definition at line 32 of file SnesCtx.hpp.

◆ typeOfAssembly

MatAssemblyType MoFEM::SnesCtx::typeOfAssembly

type of assembly at the end

Definition at line 37 of file SnesCtx.hpp.

◆ vecAssembleSwitch

boost::movelib::unique_ptr<bool> MoFEM::SnesCtx::vecAssembleSwitch
private

Definition at line 134 of file SnesCtx.hpp.

◆ vErify

bool MoFEM::SnesCtx::vErify

If true verify vector.

Definition at line 38 of file SnesCtx.hpp.

◆ zeroPreCondMatrixB

bool MoFEM::SnesCtx::zeroPreCondMatrixB

If true zero matrix, otherwise user need to do it, default true

Definition at line 35 of file SnesCtx.hpp.


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