v0.8.19
Public Types | Public Member Functions | Public 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

MoFEMErrorCode copyLoops (const SnesCtx &snes_ctx)
 Copy sequences from other SNES contex. More...
 
MoFEMErrorCode clearLoops ()
 
 SnesCtx (Interface &m_field, const std::string &problem_name)
 
virtual ~SnesCtx ()
 
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 ()
 

Public Attributes

DEPRECATED typedef MoFEM::PairNameFEMethodPtr loop_pair_type
 
DEPRECATED typedef MoFEM::FEMethodsSequence loops_to_do_type
 
DEPRECATED typedef MoFEM::BasicMethodsSequence basic_method_to_do
 
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...
 
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:
minimal_surface_area.cpp, and testing_jacobian_of_hook_element.cpp.

Definition at line 27 of file SnesCtx.hpp.

Member Typedef Documentation

◆ BasicMethodsSequence

Definition at line 51 of file SnesCtx.hpp.

◆ FEMethodsSequence

Definition at line 50 of file SnesCtx.hpp.

◆ PairNameFEMethodPtr

Definition at line 49 of file SnesCtx.hpp.

Constructor & Destructor Documentation

◆ SnesCtx()

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

Definition at line 96 of file SnesCtx.hpp.

97  : mField(m_field), moab(m_field.get_moab()), problemName(problem_name),
99  typeOfAssembly(MAT_FINAL_ASSEMBLY), vErify(false) {
100  PetscLogEventRegister("LoopSNESRhs", 0, &MOFEM_EVENT_SnesRhs);
101  PetscLogEventRegister("LoopSNESMat", 0, &MOFEM_EVENT_SnesMat);
102  }
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
PetscLogEvent MOFEM_EVENT_SnesMat
Log events to assemble tangent matrix.
Definition: SnesCtx.hpp:94
PetscLogEvent MOFEM_EVENT_SnesRhs
Log events to assemble residual.
Definition: SnesCtx.hpp:93
MoFEM::Interface & mField
database Interface
Definition: SnesCtx.hpp:29
std::string problemName
problem name
Definition: SnesCtx.hpp:32
bool vErify
If true verify vector.
Definition: SnesCtx.hpp:38
MoFEMTypes bH
If set to MF_EXIST check if element exist, default MF_EXIST.
Definition: SnesCtx.hpp:34

◆ ~SnesCtx()

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

Definition at line 104 of file SnesCtx.hpp.

104 {}

Member Function Documentation

◆ clearLoops()

MoFEMErrorCode MoFEM::SnesCtx::clearLoops ( )

Definition at line 82 of file SnesCtx.hpp.

82  {
84  loops_to_do_Mat.clear();
85  loops_to_do_Rhs.clear();
86  preProcess_Mat.clear();
87  postProcess_Mat.clear();
88  preProcess_Rhs.clear();
89  postProcess_Rhs.clear();
91  }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:499
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:506
BasicMethodsSequence postProcess_Rhs
Sequence of methods run after residual is assembled.
Definition: SnesCtx.hpp:64
FEMethodsSequence loops_to_do_Rhs
Definition: SnesCtx.hpp:55
FEMethodsSequence loops_to_do_Mat
Definition: SnesCtx.hpp:53
BasicMethodsSequence postProcess_Mat
Definition: SnesCtx.hpp:59
BasicMethodsSequence preProcess_Rhs
Sequence of methods run before residual is assembled.
Definition: SnesCtx.hpp:62
BasicMethodsSequence preProcess_Mat
Definition: SnesCtx.hpp:57

◆ copyLoops()

MoFEMErrorCode MoFEM::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 71 of file SnesCtx.hpp.

71  {
73  loops_to_do_Mat = snes_ctx.loops_to_do_Mat;
74  loops_to_do_Rhs = snes_ctx.loops_to_do_Rhs;
75  preProcess_Mat = snes_ctx.preProcess_Mat;
76  postProcess_Mat = snes_ctx.postProcess_Mat;
77  preProcess_Rhs = snes_ctx.preProcess_Rhs;
78  postProcess_Rhs = snes_ctx.postProcess_Rhs;
80  }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:499
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:506
BasicMethodsSequence postProcess_Rhs
Sequence of methods run after residual is assembled.
Definition: SnesCtx.hpp:64
FEMethodsSequence loops_to_do_Rhs
Definition: SnesCtx.hpp:55
FEMethodsSequence loops_to_do_Mat
Definition: SnesCtx.hpp:53
BasicMethodsSequence postProcess_Mat
Definition: SnesCtx.hpp:59
BasicMethodsSequence preProcess_Rhs
Sequence of methods run before residual is assembled.
Definition: SnesCtx.hpp:62
BasicMethodsSequence preProcess_Mat
Definition: SnesCtx.hpp:57

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

110 { return loops_to_do_Mat; }
FEMethodsSequence loops_to_do_Mat
Definition: SnesCtx.hpp:53

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

115 { return loops_to_do_Rhs; }
FEMethodsSequence loops_to_do_Rhs
Definition: SnesCtx.hpp:55

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

147 { return postProcess_Mat; }
BasicMethodsSequence postProcess_Mat
Definition: SnesCtx.hpp:59

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

133 { return postProcess_Rhs; }
BasicMethodsSequence postProcess_Rhs
Sequence of methods run after residual is assembled.
Definition: SnesCtx.hpp:64

◆ get_preProcess_to_do_Mat()

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

Definition at line 138 of file SnesCtx.hpp.

138 { return preProcess_Mat; }
BasicMethodsSequence preProcess_Mat
Definition: SnesCtx.hpp:57

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

124 { return preProcess_Rhs; }
BasicMethodsSequence preProcess_Rhs
Sequence of methods run before residual is assembled.
Definition: SnesCtx.hpp:62

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

142  {
143  SnesCtx *snes_ctx = (SnesCtx *)ctx;
144  // PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
146  PetscLogEventBegin(snes_ctx->MOFEM_EVENT_SnesMat, 0, 0, 0, 0);
147  if (snes_ctx->zeroPreCondMatrixB)
148  CHKERR MatZeroEntries(B);
149 
150  CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
151  CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
152  CHKERR snes_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
153  snes_ctx->problemName, COL, x, INSERT_VALUES, SCATTER_REVERSE);
154  for (auto &bit : snes_ctx->preProcess_Mat) {
155  CHKERR bit->setSnes(snes);
156  bit->snes_x = x;
157  bit->snes_A = A;
158  bit->snes_B = B;
159  CHKERR bit->setSnesCtx(SnesMethod::CTX_SNESSETJACOBIAN);
160  CHKERR snes_ctx->mField.problem_basic_method_preProcess(
161  snes_ctx->problemName, *bit);
162  CHKERR bit->setSnesCtx(SnesMethod::CTX_SNESNONE);
163  }
164  for (auto &lit : snes_ctx->loops_to_do_Mat) {
165  CHKERR lit.second->setSnesCtx(SnesMethod::CTX_SNESSETJACOBIAN);
166  CHKERR lit.second->setSnes(snes);
167  lit.second->snes_x = x;
168  lit.second->snes_A = A;
169  lit.second->snes_B = B;
170  CHKERR snes_ctx->mField.loop_finite_elements(
171  snes_ctx->problemName, lit.first, *(lit.second), snes_ctx->bH);
172  CHKERR lit.second->setSnesCtx(SnesMethod::CTX_SNESNONE);
173  }
174  for (auto &bit : snes_ctx->postProcess_Mat) {
175  CHKERR bit->setSnes(snes);
176  bit->snes_x = x;
177  bit->snes_A = A;
178  bit->snes_B = B;
179  CHKERR bit->setSnesCtx(SnesMethod::CTX_SNESSETJACOBIAN);
180  CHKERR snes_ctx->mField.problem_basic_method_postProcess(
181  snes_ctx->problemName, *bit);
182  CHKERR bit->setSnesCtx(SnesMethod::CTX_SNESNONE);
183  }
184  CHKERR MatAssemblyBegin(B, snes_ctx->typeOfAssembly);
185  CHKERR MatAssemblyEnd(B, snes_ctx->typeOfAssembly);
186  PetscLogEventEnd(snes_ctx->MOFEM_EVENT_SnesMat, 0, 0, 0, 0);
188 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:475
#define CHKERR
Inline error check.
Definition: definitions.h:594
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:405

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

199  {
200  SnesCtx *snes_ctx;
202  CHKERR SNESGetApplicationContext(snes, &snes_ctx);
203  snes_ctx->bH = bh;
205 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:475
#define CHKERR
Inline error check.
Definition: definitions.h:594
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:405

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

55  {
56  SnesCtx *snes_ctx = (SnesCtx *)ctx;
57  // PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
59  PetscLogEventBegin(snes_ctx->MOFEM_EVENT_SnesRhs, 0, 0, 0, 0);
60  CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
61  CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
62  if (snes_ctx->vErify) {
63  // Verify finite elements, check for not a number
64  CHKERR VecAssemblyBegin(f);
65  CHKERR VecAssemblyEnd(f);
66  MPI_Comm comm = PetscObjectComm((PetscObject)f);
67  PetscSynchronizedPrintf(comm, "SNES Verify x\n");
68  const Problem *prb_ptr;
69  CHKERR snes_ctx->mField.get_problem(snes_ctx->problemName, &prb_ptr);
70  CHKERR snes_ctx->mField.getInterface<Tools>()->checkVectorForNotANumber(
71  prb_ptr, COL, x);
72  }
73  CHKERR snes_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
74  snes_ctx->problemName, COL, x, INSERT_VALUES, SCATTER_REVERSE);
75 
76  auto zero_ghost_vec = [](Vec g) {
78  Vec l;
79  CHKERR VecGhostGetLocalForm(g, &l);
80  double *a;
81  CHKERR VecGetArray(l, &a);
82  int s;
83  CHKERR VecGetLocalSize(l, &s);
84  for (int i = 0; i != s; ++i)
85  a[i] = 0;
86  CHKERR VecRestoreArray(l, &a);
87  CHKERR VecGhostRestoreLocalForm(g, &l);
89  };
90  CHKERR zero_ghost_vec(f);
91 
92  for (auto &bit : snes_ctx->preProcess_Rhs) {
93  CHKERR bit->setSnes(snes);
94  bit->snes_x = x;
95  bit->snes_f = f;
96  CHKERR bit->setSnesCtx(SnesMethod::CTX_SNESSETFUNCTION);
97  CHKERR snes_ctx->mField.problem_basic_method_preProcess(
98  snes_ctx->problemName, *bit);
99  CHKERR bit->setSnesCtx(SnesMethod::CTX_SNESNONE);
100  }
101 
102  for (auto &lit : snes_ctx->loops_to_do_Rhs) {
103  CHKERR lit.second->setSnesCtx(SnesMethod::CTX_SNESSETFUNCTION);
104  CHKERR lit.second->setSnes(snes);
105  lit.second->snes_x = x;
106  lit.second->snes_f = f;
107  CHKERR snes_ctx->mField.loop_finite_elements(
108  snes_ctx->problemName, lit.first, *lit.second, snes_ctx->bH);
109  CHKERR lit.second->setSnesCtx(SnesMethod::CTX_SNESNONE);
110  if (snes_ctx->vErify) {
111  // Verify finite elements, check for not a number
112  CHKERR VecAssemblyBegin(f);
113  CHKERR VecAssemblyEnd(f);
114  MPI_Comm comm = PetscObjectComm((PetscObject)f);
115  PetscSynchronizedPrintf(comm, "SNES Verify f FE < %s >\n",
116  lit.first.c_str());
117  const Problem *prb_ptr;
118  CHKERR snes_ctx->mField.get_problem(snes_ctx->problemName, &prb_ptr);
119  CHKERR snes_ctx->mField.getInterface<Tools>()->checkVectorForNotANumber(
120  prb_ptr, ROW, f);
121  }
122  }
123 
124  for (auto &bit : snes_ctx->postProcess_Rhs) {
125  CHKERR bit->setSnes(snes);
126  bit->snes_x = x;
127  bit->snes_f = f;
128  CHKERR bit->setSnesCtx(SnesMethod::CTX_SNESSETFUNCTION);
129  CHKERR snes_ctx->mField.problem_basic_method_postProcess(
130  snes_ctx->problemName, *bit);
131  CHKERR bit->setSnesCtx(SnesMethod::CTX_SNESNONE);
132  }
133 
134  CHKERR VecGhostUpdateBegin(f, ADD_VALUES, SCATTER_REVERSE);
135  CHKERR VecGhostUpdateEnd(f, ADD_VALUES, SCATTER_REVERSE);
136  CHKERR VecAssemblyBegin(f);
137  CHKERR VecAssemblyEnd(f);
138  PetscLogEventEnd(snes_ctx->MOFEM_EVENT_SnesRhs, 0, 0, 0, 0);
140 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:475
#define CHKERR
Inline error check.
Definition: definitions.h:594
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:405

Member Data Documentation

◆ basic_method_to_do

DEPRECATED typedef MoFEM::BasicMethodsSequence MoFEM::SnesCtx::basic_method_to_do
Deprecated:
use BasicMethodsSequence

Definition at line 47 of file SnesCtx.hpp.

◆ 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.

◆ loop_pair_type

DEPRECATED typedef MoFEM::PairNameFEMethodPtr MoFEM::SnesCtx::loop_pair_type
Deprecated:
use PairNameFEMethodPtr

Definition at line 41 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 53 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 55 of file SnesCtx.hpp.

◆ loops_to_do_type

DEPRECATED typedef MoFEM::FEMethodsSequence MoFEM::SnesCtx::loops_to_do_type
Deprecated:
use FEMethodsSequence

Definition at line 44 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

Log events to assemble tangent matrix.

Definition at line 94 of file SnesCtx.hpp.

◆ MOFEM_EVENT_SnesRhs

PetscLogEvent MoFEM::SnesCtx::MOFEM_EVENT_SnesRhs

Log events to assemble residual.

Definition at line 93 of file SnesCtx.hpp.

◆ postProcess_Mat

BasicMethodsSequence MoFEM::SnesCtx::postProcess_Mat

Sequence of methods run after tangent matrix is assembled

Definition at line 59 of file SnesCtx.hpp.

◆ postProcess_Rhs

BasicMethodsSequence MoFEM::SnesCtx::postProcess_Rhs

Sequence of methods run after residual is assembled.

Definition at line 64 of file SnesCtx.hpp.

◆ preProcess_Mat

BasicMethodsSequence MoFEM::SnesCtx::preProcess_Mat

Sequence of methods run before tangent matrix is assembled

Definition at line 57 of file SnesCtx.hpp.

◆ preProcess_Rhs

BasicMethodsSequence MoFEM::SnesCtx::preProcess_Rhs

Sequence of methods run before residual is assembled.

Definition at line 62 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.

◆ 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 file: