v0.14.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
 
FEMethodsSequencegetSetOperators ()
 
FEMethodsSequencegetComputeRhs ()
 
BasicMethodsSequencegetPreProcComputeRhs ()
 
BasicMethodsSequencegetPostProcComputeRhs ()
 
BasicMethodsSequencegetPreProcSetOperators ()
 
BasicMethodsSequencegetPostProcSetOperators ()
 
MoFEMErrorCode copyLoops (const SnesCtx &snes_ctx)
 Copy sequences from other SNES contex. More...
 
MoFEMErrorCode clearLoops ()
 Clear loops. More...
 
DEPRECATED FEMethodsSequenceget_loops_to_do_Mat ()
 
DEPRECATED FEMethodsSequenceget_loops_to_do_Rhs ()
 
DEPRECATED BasicMethodsSequenceget_preProcess_to_do_Rhs ()
 
DEPRECATED BasicMethodsSequenceget_postProcess_to_do_Rhs ()
 
DEPRECATED BasicMethodsSequenceget_preProcess_to_do_Mat ()
 
DEPRECATED BasicMethodsSequenceget_postProcess_to_do_Mat ()
 

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
minimal_surface_area.cpp, navier_stokes.cpp, nonlinear_dynamics.cpp, simple_contact.cpp, simple_contact_thermal.cpp, test_jacobian_of_simple_contact_element.cpp, testing_jacobian_of_hook_element.cpp, and testing_jacobian_of_hook_scaled_with_density_element.cpp.

Definition at line 13 of file SnesCtx.hpp.

Member Typedef Documentation

◆ BasicMethodsSequence

Definition at line 28 of file SnesCtx.hpp.

◆ FEMethodsSequence

Definition at line 27 of file SnesCtx.hpp.

◆ PairNameFEMethodPtr

Definition at line 26 of file SnesCtx.hpp.

Constructor & Destructor Documentation

◆ SnesCtx()

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

Definition at line 43 of file SnesCtx.hpp.

44  : mField(m_field), moab(m_field.get_moab()), problemName(problem_name),
46  typeOfAssembly(MAT_FINAL_ASSEMBLY), vErify(false) {
47  PetscLogEventRegister("LoopSNESRhs", 0, &MOFEM_EVENT_SnesRhs);
48  PetscLogEventRegister("LoopSNESMat", 0, &MOFEM_EVENT_SnesMat);
49  if (!LogManager::checkIfChannelExist("SNES_WORLD")) {
50  auto core_log = logging::core::get();
51  core_log->add_sink(
53  LogManager::setLog("SNES_WORLD");
54  MOFEM_LOG_TAG("SNES_WORLD", "SNES");
55  }
56  }

◆ ~SnesCtx()

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

Member Function Documentation

◆ clearLoops()

MoFEMErrorCode SnesCtx::clearLoops ( )

Clear loops.

Returns
MoFEMErrorCode

Definition at line 16 of file SnesCtx.cpp.

16  {
18  loops_to_do_Mat.clear();
19  loops_to_do_Rhs.clear();
20  preProcess_Mat.clear();
21  postProcess_Mat.clear();
22  preProcess_Rhs.clear();
23  postProcess_Rhs.clear();
25 }

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

5  {
7  loops_to_do_Mat = snes_ctx.loops_to_do_Mat;
8  loops_to_do_Rhs = snes_ctx.loops_to_do_Rhs;
9  preProcess_Mat = snes_ctx.preProcess_Mat;
10  postProcess_Mat = snes_ctx.postProcess_Mat;
11  preProcess_Rhs = snes_ctx.preProcess_Rhs;
12  postProcess_Rhs = snes_ctx.postProcess_Rhs;
14 }

◆ get_loops_to_do_Mat()

DEPRECATED FEMethodsSequence& MoFEM::SnesCtx::get_loops_to_do_Mat ( )
inline
Deprecated:
use getSetOperators

Definition at line 125 of file SnesCtx.hpp.

125  {
126  return getSetOperators();
127  }

◆ get_loops_to_do_Rhs()

DEPRECATED FEMethodsSequence& MoFEM::SnesCtx::get_loops_to_do_Rhs ( )
inline
Deprecated:
use getComputeRhs

Definition at line 130 of file SnesCtx.hpp.

130  {
131  return getComputeRhs();
132  }

◆ get_postProcess_to_do_Mat()

DEPRECATED BasicMethodsSequence& MoFEM::SnesCtx::get_postProcess_to_do_Mat ( )
inline
Deprecated:
use getPostProcSetOperators

Definition at line 150 of file SnesCtx.hpp.

150  {
151  return getPostProcSetOperators();
152  }

◆ get_postProcess_to_do_Rhs()

DEPRECATED BasicMethodsSequence& MoFEM::SnesCtx::get_postProcess_to_do_Rhs ( )
inline
Deprecated:
use getPostProcComputeRhs

Definition at line 140 of file SnesCtx.hpp.

140  {
141  return getPostProcComputeRhs();
142  }

◆ get_preProcess_to_do_Mat()

DEPRECATED BasicMethodsSequence& MoFEM::SnesCtx::get_preProcess_to_do_Mat ( )
inline
Deprecated:
use getPreProcSetOperators

Definition at line 145 of file SnesCtx.hpp.

145  {
146  return getPreProcSetOperators();
147  }

◆ get_preProcess_to_do_Rhs()

DEPRECATED BasicMethodsSequence& MoFEM::SnesCtx::get_preProcess_to_do_Rhs ( )
inline
Deprecated:
use getPreProcComputeRhs

Definition at line 135 of file SnesCtx.hpp.

135  {
136  return getPreProcComputeRhs();
137  }

◆ getComputeRhs()

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

Definition at line 69 of file SnesCtx.hpp.

69 { return loops_to_do_Rhs; }

◆ getPostProcComputeRhs()

BasicMethodsSequence& MoFEM::SnesCtx::getPostProcComputeRhs ( )
inline

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

87 { return postProcess_Rhs; }

◆ getPostProcSetOperators()

BasicMethodsSequence& MoFEM::SnesCtx::getPostProcSetOperators ( )
inline

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

101 { return postProcess_Mat; }

◆ getPreProcComputeRhs()

BasicMethodsSequence& MoFEM::SnesCtx::getPreProcComputeRhs ( )
inline

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

78 { return preProcess_Rhs; }

◆ getPreProcSetOperators()

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

Definition at line 92 of file SnesCtx.hpp.

92 { return preProcess_Mat; }

◆ getSetOperators()

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

Definition at line 64 of file SnesCtx.hpp.

64 { return loops_to_do_Mat; }

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  snes_ctx->matAssembleSwitch = boost::movelib::make_unique<bool>(true);
151  auto cache_ptr = boost::make_shared<CacheTuple>();
152  CHKERR snes_ctx->mField.cache_problem_entities(snes_ctx->problemName,
153  cache_ptr);
154  Vec dx;
155  CHKERR SNESGetSolutionUpdate(snes, &dx);
156 
157  auto set = [&](auto &fe) {
158  fe.snes = snes;
159  fe.snes_x = x;
160  fe.snes_dx = dx;
161  fe.snes_A = A;
162  fe.snes_B = B;
163  fe.snes_ctx = SnesMethod::CTX_SNESSETJACOBIAN;
164  fe.ksp_ctx = KspMethod::CTX_OPERATORS;
167 
168  CHKERR SNESGetKSP(snes, &fe.ksp);
169 
170  fe.cacheWeakPtr = cache_ptr;
171  };
172 
173  auto unset = [&](auto &fe) {
174  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
175  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
176  fe.data_ctx = PetscData::CtxSetNone;
177  };
178 
179 
180  CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
181  CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
182  CHKERR snes_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
183  snes_ctx->problemName, COL, x, INSERT_VALUES, SCATTER_REVERSE);
184  for (auto &bit : snes_ctx->preProcess_Mat) {
185  bit->matAssembleSwitch = boost::move(snes_ctx->matAssembleSwitch);
186  set(*bit);
187  CHKERR snes_ctx->mField.problem_basic_method_preProcess(
188  snes_ctx->problemName, *bit);
189  unset(*bit);
190  snes_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
191  }
192 
193 
194 
195  for (auto &lit : snes_ctx->loops_to_do_Mat) {
196  lit.second->matAssembleSwitch = boost::move(snes_ctx->matAssembleSwitch);
197  set(*lit.second);
198  CHKERR snes_ctx->mField.loop_finite_elements(
199  snes_ctx->problemName, lit.first, *(lit.second), nullptr, snes_ctx->bH,
200  cache_ptr);
201  unset(*lit.second);
202  snes_ctx->matAssembleSwitch = boost::move(lit.second->matAssembleSwitch);
203  }
204 
205  for (auto &bit : snes_ctx->postProcess_Mat) {
206  bit->matAssembleSwitch = boost::move(snes_ctx->matAssembleSwitch);
207  set(*bit);
208  CHKERR snes_ctx->mField.problem_basic_method_postProcess(
209  snes_ctx->problemName, *bit);
210  unset(*bit);
211  snes_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
212  }
213 
214  if (*snes_ctx->matAssembleSwitch) {
215  CHKERR MatAssemblyBegin(B, snes_ctx->typeOfAssembly);
216  CHKERR MatAssemblyEnd(B, snes_ctx->typeOfAssembly);
217  }
218  PetscLogEventEnd(snes_ctx->MOFEM_EVENT_SnesMat, 0, 0, 0, 0);
220 }

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

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

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

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

Member Data Documentation

◆ bH

MoFEMTypes MoFEM::SnesCtx::bH

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

Definition at line 20 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 30 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 32 of file SnesCtx.hpp.

◆ matAssembleSwitch

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

Definition at line 156 of file SnesCtx.hpp.

◆ mField

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

database Interface

Definition at line 15 of file SnesCtx.hpp.

◆ moab

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

moab Interface

Definition at line 16 of file SnesCtx.hpp.

◆ MOFEM_EVENT_SnesMat

PetscLogEvent MoFEM::SnesCtx::MOFEM_EVENT_SnesMat
private

Log events to assemble tangent matrix.

Definition at line 158 of file SnesCtx.hpp.

◆ MOFEM_EVENT_SnesRhs

PetscLogEvent MoFEM::SnesCtx::MOFEM_EVENT_SnesRhs
private

Log events to assemble residual.

Definition at line 157 of file SnesCtx.hpp.

◆ postProcess_Mat

BasicMethodsSequence MoFEM::SnesCtx::postProcess_Mat

Sequence of methods run after tangent matrix is assembled

Definition at line 36 of file SnesCtx.hpp.

◆ postProcess_Rhs

BasicMethodsSequence MoFEM::SnesCtx::postProcess_Rhs

Sequence of methods run after residual is assembled.

Definition at line 41 of file SnesCtx.hpp.

◆ preProcess_Mat

BasicMethodsSequence MoFEM::SnesCtx::preProcess_Mat

Sequence of methods run before tangent matrix is assembled

Definition at line 34 of file SnesCtx.hpp.

◆ preProcess_Rhs

BasicMethodsSequence MoFEM::SnesCtx::preProcess_Rhs

Sequence of methods run before residual is assembled.

Definition at line 39 of file SnesCtx.hpp.

◆ problemName

std::string MoFEM::SnesCtx::problemName

problem name

Definition at line 18 of file SnesCtx.hpp.

◆ typeOfAssembly

MatAssemblyType MoFEM::SnesCtx::typeOfAssembly

type of assembly at the end

Definition at line 23 of file SnesCtx.hpp.

◆ vecAssembleSwitch

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

Definition at line 155 of file SnesCtx.hpp.

◆ vErify

bool MoFEM::SnesCtx::vErify

If true verify vector.

Definition at line 24 of file SnesCtx.hpp.

◆ zeroPreCondMatrixB

bool MoFEM::SnesCtx::zeroPreCondMatrixB

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

Definition at line 21 of file SnesCtx.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::SnesCtx::getPostProcComputeRhs
BasicMethodsSequence & getPostProcComputeRhs()
Definition: SnesCtx.hpp:87
MoFEM::KspMethod::CTX_KSPNONE
@ CTX_KSPNONE
Definition: LoopMethods.hpp:76
g
constexpr double g
Definition: shallow_wave.cpp:63
MoFEM::LogManager::checkIfChannelExist
static bool checkIfChannelExist(const std::string channel)
Check if channel exist.
Definition: LogManager.cpp:404
MoFEM::SnesCtx::MOFEM_EVENT_SnesMat
PetscLogEvent MOFEM_EVENT_SnesMat
Log events to assemble tangent matrix.
Definition: SnesCtx.hpp:158
MoFEM::PetscData::CtxSetB
static constexpr Switches CtxSetB
Definition: LoopMethods.hpp:39
MoFEM::PetscData::CtxSetA
static constexpr Switches CtxSetA
Definition: LoopMethods.hpp:38
MoFEM::SnesCtx::vErify
bool vErify
If true verify vector.
Definition: SnesCtx.hpp:24
MoFEM::SnesCtx::getPreProcSetOperators
BasicMethodsSequence & getPreProcSetOperators()
Definition: SnesCtx.hpp:92
MoFEM::SnesCtx::moab
moab::Interface & moab
moab Interface
Definition: SnesCtx.hpp:16
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
MoFEM::LogManager::createSink
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:298
MoFEM::SnesCtx::postProcess_Mat
BasicMethodsSequence postProcess_Mat
Definition: SnesCtx.hpp:36
MoFEM::KspMethod::CTX_SETFUNCTION
@ CTX_SETFUNCTION
Definition: LoopMethods.hpp:76
MoFEM::SnesMethod::CTX_SNESSETJACOBIAN
@ CTX_SNESSETJACOBIAN
Definition: LoopMethods.hpp:110
MoFEM::SnesMethod::CTX_SNESSETFUNCTION
@ CTX_SNESSETFUNCTION
Definition: LoopMethods.hpp:110
MoFEM::SnesCtx::getPostProcSetOperators
BasicMethodsSequence & getPostProcSetOperators()
Definition: SnesCtx.hpp:101
ROW
@ ROW
Definition: definitions.h:136
MoFEM::SnesCtx::loops_to_do_Rhs
FEMethodsSequence loops_to_do_Rhs
Definition: SnesCtx.hpp:32
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::SnesCtx::zeroPreCondMatrixB
bool zeroPreCondMatrixB
Definition: SnesCtx.hpp:21
MoFEM::SnesCtx::mField
MoFEM::Interface & mField
database Interface
Definition: SnesCtx.hpp:15
a
constexpr double a
Definition: approx_sphere.cpp:30
MoFEM::SnesCtx::typeOfAssembly
MatAssemblyType typeOfAssembly
type of assembly at the end
Definition: SnesCtx.hpp:23
MoFEM::SnesMethod::CTX_SNESNONE
@ CTX_SNESNONE
Definition: LoopMethods.hpp:110
bit
auto bit
set bit
Definition: hanging_node_approx.cpp:75
MoFEM::SnesCtx::problemName
std::string problemName
problem name
Definition: SnesCtx.hpp:18
MoFEM::SnesCtx::postProcess_Rhs
BasicMethodsSequence postProcess_Rhs
Sequence of methods run after residual is assembled.
Definition: SnesCtx.hpp:41
COL
@ COL
Definition: definitions.h:136
MoFEM::LogManager::getStrmWorld
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Definition: LogManager.cpp:344
MoFEM::SnesCtx::getComputeRhs
FEMethodsSequence & getComputeRhs()
Definition: SnesCtx.hpp:69
MoFEM::SnesCtx::MOFEM_EVENT_SnesRhs
PetscLogEvent MOFEM_EVENT_SnesRhs
Log events to assemble residual.
Definition: SnesCtx.hpp:157
MoFEM::SnesCtx::preProcess_Rhs
BasicMethodsSequence preProcess_Rhs
Sequence of methods run before residual is assembled.
Definition: SnesCtx.hpp:39
MoFEM::PetscData::CtxSetNone
static constexpr Switches CtxSetNone
Definition: LoopMethods.hpp:36
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::SnesCtx::getSetOperators
FEMethodsSequence & getSetOperators()
Definition: SnesCtx.hpp:64
MoFEM::PetscData::CtxSetX
static constexpr Switches CtxSetX
Definition: LoopMethods.hpp:40
MoFEM::SnesCtx::loops_to_do_Mat
FEMethodsSequence loops_to_do_Mat
Definition: SnesCtx.hpp:30
MoFEM::SnesCtx::getPreProcComputeRhs
BasicMethodsSequence & getPreProcComputeRhs()
Definition: SnesCtx.hpp:78
SnesCtx
MoFEM::PetscData::CtxSetDX
static constexpr Switches CtxSetDX
Definition: LoopMethods.hpp:41
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
MoFEM::SnesCtx::preProcess_Mat
BasicMethodsSequence preProcess_Mat
Definition: SnesCtx.hpp:34
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:64
MoFEM::SnesCtx::bH
MoFEMTypes bH
If set to MF_EXIST check if element exist, default MF_EXIST.
Definition: SnesCtx.hpp:20
MoFEM::KspMethod::CTX_OPERATORS
@ CTX_OPERATORS
Definition: LoopMethods.hpp:76
MoFEM::PetscData::CtxSetF
static constexpr Switches CtxSetF
Definition: LoopMethods.hpp:37
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::LogManager::setLog
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389
MF_EXIST
@ MF_EXIST
Definition: definitions.h:113
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21