v0.9.2
SnesCtx.cpp
Go to the documentation of this file.
1 /* This file is part of MoFEM.
2  * MoFEM is free software: you can redistribute it and/or modify it under
3  * the terms of the GNU Lesser General Public License as published by the
4  * Free Software Foundation, either version 3 of the License, or (at your
5  * option) any later version.
6  *
7  * MoFEM is distributed in the hope that it will be useful, but WITHOUT
8  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
9  * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
10  * License for more details.
11  *
12  * You should have received a copy of the GNU Lesser General Public
13  * License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
14 
15 namespace MoFEM {
16 
17 PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx) {
18  SnesCtx *snes_ctx = (SnesCtx *)ctx;
19  // PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
21  PetscLogEventBegin(snes_ctx->MOFEM_EVENT_SnesRhs, 0, 0, 0, 0);
22  CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
23  CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
24  if (snes_ctx->vErify) {
25  // Verify finite elements, check for not a number
26  CHKERR VecAssemblyBegin(f);
27  CHKERR VecAssemblyEnd(f);
28  MPI_Comm comm = PetscObjectComm((PetscObject)f);
29  PetscSynchronizedPrintf(comm, "SNES Verify x\n");
30  const Problem *prb_ptr;
31  CHKERR snes_ctx->mField.get_problem(snes_ctx->problemName, &prb_ptr);
32  CHKERR snes_ctx->mField.getInterface<Tools>()->checkVectorForNotANumber(
33  prb_ptr, COL, x);
34  }
35  CHKERR snes_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
36  snes_ctx->problemName, COL, x, INSERT_VALUES, SCATTER_REVERSE);
37 
38  auto zero_ghost_vec = [](Vec g) {
40  Vec l;
41  CHKERR VecGhostGetLocalForm(g, &l);
42  double *a;
43  CHKERR VecGetArray(l, &a);
44  int s;
45  CHKERR VecGetLocalSize(l, &s);
46  for (int i = 0; i != s; ++i)
47  a[i] = 0;
48  CHKERR VecRestoreArray(l, &a);
49  CHKERR VecGhostRestoreLocalForm(g, &l);
51  };
52  CHKERR zero_ghost_vec(f);
53 
54  snes_ctx->vecAssembleSwitch = boost::movelib::make_unique<bool>(true);
55 
56  auto set = [&](auto &fe) {
57  fe.snes = snes;
58  fe.snes_x = x;
59  fe.snes_f = f;
60  fe.snes_ctx = SnesMethod::CTX_SNESSETFUNCTION;
61  fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
62  fe.data_ctx = PetscData::CtxSetF | PetscData::CtxSetX;
63  };
64 
65  auto unset = [&](auto &fe) {
66  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
67  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
68  fe.data_ctx = PetscData::CtxSetNone;
69  };
70 
71  for (auto &bit : snes_ctx->preProcess_Rhs) {
72  bit->vecAssembleSwitch = boost::move(snes_ctx->vecAssembleSwitch);
73  set(*bit);
75  snes_ctx->problemName, *bit);
76  unset(*bit);
77  snes_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
78  }
79 
80  for (auto &lit : snes_ctx->loops_to_do_Rhs) {
81  lit.second->vecAssembleSwitch = boost::move(snes_ctx->vecAssembleSwitch);
82  set(*lit.second);
84  snes_ctx->problemName, lit.first, *lit.second, nullptr, snes_ctx->bH);
85  unset(*lit.second);
86  if (snes_ctx->vErify) {
87  // Verify finite elements, check for not a number
88  CHKERR VecAssemblyBegin(f);
89  CHKERR VecAssemblyEnd(f);
90  MPI_Comm comm = PetscObjectComm((PetscObject)f);
91  PetscSynchronizedPrintf(comm, "SNES Verify f FE < %s >\n",
92  lit.first.c_str());
93  const Problem *prb_ptr;
94  CHKERR snes_ctx->mField.get_problem(snes_ctx->problemName, &prb_ptr);
95  CHKERR snes_ctx->mField.getInterface<Tools>()->checkVectorForNotANumber(
96  prb_ptr, ROW, f);
97  }
98 
99  snes_ctx->vecAssembleSwitch = boost::move(lit.second->vecAssembleSwitch);
100  }
101 
102  for (auto &bit : snes_ctx->postProcess_Rhs) {
103  bit->vecAssembleSwitch = boost::move(snes_ctx->vecAssembleSwitch);
104  set(*bit);
106  snes_ctx->problemName, *bit);
107  unset(*bit);
108  snes_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
109  }
110 
111  if (snes_ctx->vecAssembleSwitch) {
112  CHKERR VecGhostUpdateBegin(f, ADD_VALUES, SCATTER_REVERSE);
113  CHKERR VecGhostUpdateEnd(f, ADD_VALUES, SCATTER_REVERSE);
114  CHKERR VecAssemblyBegin(f);
115  CHKERR VecAssemblyEnd(f);
116  }
117  PetscLogEventEnd(snes_ctx->MOFEM_EVENT_SnesRhs, 0, 0, 0, 0);
119 }
120 
121 PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx) {
122  SnesCtx *snes_ctx = (SnesCtx *)ctx;
123  // PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
125  PetscLogEventBegin(snes_ctx->MOFEM_EVENT_SnesMat, 0, 0, 0, 0);
126  if (snes_ctx->zeroPreCondMatrixB)
127  CHKERR MatZeroEntries(B);
128 
129  snes_ctx->matAssembleSwitch = boost::movelib::make_unique<bool>(true);
130 
131  auto set = [&](auto &fe) {
132  fe.snes = snes;
133  fe.snes_x = x;
134  fe.snes_A = A;
135  fe.snes_B = B;
136  fe.snes_ctx = SnesMethod::CTX_SNESSETJACOBIAN;
137  fe.ksp_ctx = KspMethod::CTX_OPERATORS;
138  fe.data_ctx =
140  };
141 
142  auto unset = [&](auto &fe) {
143  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
144  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
145  fe.data_ctx = PetscData::CtxSetNone;
146  };
147 
148  CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
149  CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
150  CHKERR snes_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
151  snes_ctx->problemName, COL, x, INSERT_VALUES, SCATTER_REVERSE);
152  for (auto &bit : snes_ctx->preProcess_Mat) {
153  bit->matAssembleSwitch = boost::move(snes_ctx->matAssembleSwitch);
154  set(*bit);
156  snes_ctx->problemName, *bit);
157  unset(*bit);
158  snes_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
159  }
160 
161  for (auto &lit : snes_ctx->loops_to_do_Mat) {
162  lit.second->matAssembleSwitch = boost::move(snes_ctx->matAssembleSwitch);
163  set(*lit.second);
165  snes_ctx->problemName, lit.first, *(lit.second), nullptr, snes_ctx->bH);
166  unset(*lit.second);
167  snes_ctx->matAssembleSwitch = boost::move(lit.second->matAssembleSwitch);
168  }
169 
170  for (auto &bit : snes_ctx->postProcess_Mat) {
171  bit->matAssembleSwitch = boost::move(snes_ctx->matAssembleSwitch);
172  set(*bit);
174  snes_ctx->problemName, *bit);
175  unset(*bit);
176  snes_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
177  }
178 
179  if (*snes_ctx->matAssembleSwitch) {
180  CHKERR MatAssemblyBegin(B, snes_ctx->typeOfAssembly);
181  CHKERR MatAssemblyEnd(B, snes_ctx->typeOfAssembly);
182  }
183  PetscLogEventEnd(snes_ctx->MOFEM_EVENT_SnesMat, 0, 0, 0, 0);
185 }
186 
187 MoFEMErrorCode SnesMoFEMSetAssemblyType(SNES snes, MatAssemblyType type) {
188  SnesCtx *snes_ctx;
189  // PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
191  CHKERR SNESGetApplicationContext(snes, &snes_ctx);
192  snes_ctx->typeOfAssembly = type;
194 }
195 
197  SnesCtx *snes_ctx;
199  CHKERR SNESGetApplicationContext(snes, &snes_ctx);
200  snes_ctx->bH = bh;
202 }
203 
204 } // namespace MoFEM
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:158
PetscLogEvent MOFEM_EVENT_SnesRhs
Log events to assemble residual.
Definition: SnesCtx.hpp:157
static constexpr Switches CtxSetB
Definition: LoopMethods.hpp:64
boost::movelib::unique_ptr< bool > matAssembleSwitch
Definition: SnesCtx.hpp:156
MoFEM::Interface & mField
database Interface
Definition: SnesCtx.hpp:29
std::string problemName
problem name
Definition: SnesCtx.hpp:32
PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx)
This is MoFEM implementation for the left hand side (tangent matrix) evaluation in SNES solver.
Definition: SnesCtx.cpp:121
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
Auxiliary tools.
Definition: Tools.hpp:30
keeps basic data about problemThis is low level structure with information about problem,...
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:21
MoFEMErrorCode SnesMoFEMSetAssemblyType(SNES snes, MatAssemblyType type)
Set assembly type at the end of SnesMat.
Definition: SnesCtx.cpp:187
BasicMethodsSequence postProcess_Rhs
Sequence of methods run after residual is assembled.
Definition: SnesCtx.hpp:64
MoFEMErrorCode getInterface(const MOFEMuuid &uuid, IFACE *&iface) const
Get interface by uuid and return reference to pointer of interface.
bool vErify
If true verify vector.
Definition: SnesCtx.hpp:38
FEMethodsSequence loops_to_do_Rhs
Definition: SnesCtx.hpp:55
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:67
FEMethodsSequence loops_to_do_Mat
Definition: SnesCtx.hpp:53
Interface for nonlinear (SNES) solver.
Definition: SnesCtx.hpp:27
static constexpr Switches CtxSetF
Definition: LoopMethods.hpp:62
BasicMethodsSequence postProcess_Mat
Definition: SnesCtx.hpp:59
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:36
boost::movelib::unique_ptr< bool > vecAssembleSwitch
Definition: SnesCtx.hpp:155
virtual MoFEMErrorCode problem_basic_method_preProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
#define CHKERR
Inline error check.
Definition: definitions.h:602
BasicMethodsSequence preProcess_Rhs
Sequence of methods run before residual is assembled.
Definition: SnesCtx.hpp:62
static constexpr Switches CtxSetA
Definition: LoopMethods.hpp:63
virtual const Problem * get_problem(const std::string &problem_name) const =0
Get the problem object.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413
MoFEMTypes bH
If set to MF_EXIST check if element exist, default MF_EXIST.
Definition: SnesCtx.hpp:34
PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx)
This is MoFEM implementation for the right hand side (residual vector) evaluation in SNES solver.
Definition: SnesCtx.cpp:17
virtual MoFEMErrorCode problem_basic_method_postProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethodThis function set data about problem, adjacencies and other multi-indices in ...
static constexpr Switches CtxSetX
Definition: LoopMethods.hpp:65
MoFEMErrorCode SnesMoFEMSetBehavior(SNES snes, MoFEMTypes bh)
Set behavior if finite element in sequence does not exist.
Definition: SnesCtx.cpp:196
FTensor::Index< 'l', 2 > l
Definition: ContactOps.hpp:29
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
Definition: definitions.h:189
BasicMethodsSequence preProcess_Mat
Definition: SnesCtx.hpp:57
DEPRECATED MoFEMErrorCode loop_finite_elements(const Problem *problem_ptr, const std::string &fe_name, FEMethod &method, int lower_rank, int upper_rank, MoFEMTypes bh, int verb=DEFAULT_VERBOSITY)
static constexpr Switches CtxSetNone
Definition: LoopMethods.hpp:61
FTensor::Index< 'i', 2 > i
[Common data]
Definition: ContactOps.hpp:26