v0.14.0
SnesCtx.cpp
Go to the documentation of this file.
1 
2 
3 namespace MoFEM {
4 
11  preProcess_Rhs = snes_ctx.preProcess_Rhs;
14 }
15 
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 }
26 
27 PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx) {
28  SnesCtx *snes_ctx = (SnesCtx *)ctx;
29  // PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
31  PetscLogEventBegin(snes_ctx->MOFEM_EVENT_SnesRhs, 0, 0, 0, 0);
32  CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
33  CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
34  if (snes_ctx->vErify) {
35  // Verify finite elements, check for not a number
36  CHKERR VecAssemblyBegin(f);
37  CHKERR VecAssemblyEnd(f);
38  MPI_Comm comm = PetscObjectComm((PetscObject)f);
39  PetscSynchronizedPrintf(comm, "SNES Verify x\n");
40  const Problem *prb_ptr;
41  CHKERR snes_ctx->mField.get_problem(snes_ctx->problemName, &prb_ptr);
42  CHKERR snes_ctx->mField.getInterface<Tools>()->checkVectorForNotANumber(
43  prb_ptr, COL, x);
44  }
45  CHKERR snes_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
46  snes_ctx->problemName, COL, x, INSERT_VALUES, SCATTER_REVERSE);
47 
48  auto zero_ghost_vec = [](Vec g) {
50  Vec l;
51  CHKERR VecGhostGetLocalForm(g, &l);
52  double *a;
53  CHKERR VecGetArray(l, &a);
54  int s;
55  CHKERR VecGetLocalSize(l, &s);
56  for (int i = 0; i != s; ++i)
57  a[i] = 0;
58  CHKERR VecRestoreArray(l, &a);
59  CHKERR VecGhostRestoreLocalForm(g, &l);
61  };
62  CHKERR zero_ghost_vec(f);
63 
64  snes_ctx->vecAssembleSwitch = boost::movelib::make_unique<bool>(true);
65  auto cache_ptr = boost::make_shared<CacheTuple>();
66  CHKERR snes_ctx->mField.cache_problem_entities(snes_ctx->problemName,
67  cache_ptr);
68 
69  auto set = [&](auto &fe) {
70  fe.snes = snes;
71  fe.snes_x = x;
72  fe.snes_f = f;
73  fe.snes_ctx = SnesMethod::CTX_SNESSETFUNCTION;
74  fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
75  fe.data_ctx = PetscData::CtxSetF | PetscData::CtxSetX;
76 
77  CHKERR SNESGetKSP(snes, &fe.ksp);
78 
79  fe.cacheWeakPtr = cache_ptr;
80  };
81 
82  auto unset = [&](auto &fe) {
83  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
84  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
85  fe.data_ctx = PetscData::CtxSetNone;
86  };
87 
88  for (auto &bit : snes_ctx->preProcess_Rhs) {
89  bit->vecAssembleSwitch = boost::move(snes_ctx->vecAssembleSwitch);
90  set(*bit);
92  snes_ctx->problemName, *bit);
93  unset(*bit);
94  snes_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
95  }
96 
97  for (auto &lit : snes_ctx->loops_to_do_Rhs) {
98  lit.second->vecAssembleSwitch = boost::move(snes_ctx->vecAssembleSwitch);
99  set(*lit.second);
101  snes_ctx->problemName, lit.first, *lit.second, nullptr, snes_ctx->bH,
102  cache_ptr);
103  unset(*lit.second);
104  if (snes_ctx->vErify) {
105  // Verify finite elements, check for not a number
106  CHKERR VecAssemblyBegin(f);
107  CHKERR VecAssemblyEnd(f);
108  MPI_Comm comm = PetscObjectComm((PetscObject)f);
109  PetscSynchronizedPrintf(comm, "SNES Verify f FE < %s >\n",
110  lit.first.c_str());
111  const Problem *prb_ptr;
112  CHKERR snes_ctx->mField.get_problem(snes_ctx->problemName, &prb_ptr);
113  CHKERR snes_ctx->mField.getInterface<Tools>()->checkVectorForNotANumber(
114  prb_ptr, ROW, f);
115  }
116 
117  snes_ctx->vecAssembleSwitch = boost::move(lit.second->vecAssembleSwitch);
118  }
119 
120  for (auto &bit : snes_ctx->postProcess_Rhs) {
121  bit->vecAssembleSwitch = boost::move(snes_ctx->vecAssembleSwitch);
122  set(*bit);
124  snes_ctx->problemName, *bit);
125  unset(*bit);
126  snes_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
127  }
128 
129  if (*(snes_ctx->vecAssembleSwitch)) {
130  CHKERR VecGhostUpdateBegin(f, ADD_VALUES, SCATTER_REVERSE);
131  CHKERR VecGhostUpdateEnd(f, ADD_VALUES, SCATTER_REVERSE);
132  CHKERR VecAssemblyBegin(f);
133  CHKERR VecAssemblyEnd(f);
134  }
135  PetscLogEventEnd(snes_ctx->MOFEM_EVENT_SnesRhs, 0, 0, 0, 0);
137 }
138 
139 PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx) {
140  SnesCtx *snes_ctx = (SnesCtx *)ctx;
141  // PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
143  PetscLogEventBegin(snes_ctx->MOFEM_EVENT_SnesMat, 0, 0, 0, 0);
144  if (snes_ctx->zeroPreCondMatrixB)
145  CHKERR MatZeroEntries(B);
146 
147  snes_ctx->matAssembleSwitch = boost::movelib::make_unique<bool>(true);
148  auto cache_ptr = boost::make_shared<CacheTuple>();
149  CHKERR snes_ctx->mField.cache_problem_entities(snes_ctx->problemName,
150  cache_ptr);
151 
152  auto set = [&](auto &fe) {
153  fe.snes = snes;
154  fe.snes_x = x;
155  fe.snes_A = A;
156  fe.snes_B = B;
157  fe.snes_ctx = SnesMethod::CTX_SNESSETJACOBIAN;
158  fe.ksp_ctx = KspMethod::CTX_OPERATORS;
160 
161  CHKERR SNESGetKSP(snes, &fe.ksp);
162 
163  fe.cacheWeakPtr = cache_ptr;
164  };
165 
166  auto unset = [&](auto &fe) {
167  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
168  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
169  fe.data_ctx = PetscData::CtxSetNone;
170  };
171 
172 
173  CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
174  CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
175  CHKERR snes_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
176  snes_ctx->problemName, COL, x, INSERT_VALUES, SCATTER_REVERSE);
177  for (auto &bit : snes_ctx->preProcess_Mat) {
178  bit->matAssembleSwitch = boost::move(snes_ctx->matAssembleSwitch);
179  set(*bit);
181  snes_ctx->problemName, *bit);
182  unset(*bit);
183  snes_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
184  }
185 
186 
187 
188  for (auto &lit : snes_ctx->loops_to_do_Mat) {
189  lit.second->matAssembleSwitch = boost::move(snes_ctx->matAssembleSwitch);
190  set(*lit.second);
192  snes_ctx->problemName, lit.first, *(lit.second), nullptr, snes_ctx->bH,
193  cache_ptr);
194  unset(*lit.second);
195  snes_ctx->matAssembleSwitch = boost::move(lit.second->matAssembleSwitch);
196  }
197 
198  for (auto &bit : snes_ctx->postProcess_Mat) {
199  bit->matAssembleSwitch = boost::move(snes_ctx->matAssembleSwitch);
200  set(*bit);
202  snes_ctx->problemName, *bit);
203  unset(*bit);
204  snes_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
205  }
206 
207  if (*snes_ctx->matAssembleSwitch) {
208  CHKERR MatAssemblyBegin(B, snes_ctx->typeOfAssembly);
209  CHKERR MatAssemblyEnd(B, snes_ctx->typeOfAssembly);
210  }
211  PetscLogEventEnd(snes_ctx->MOFEM_EVENT_SnesMat, 0, 0, 0, 0);
213 }
214 
215 MoFEMErrorCode SnesMoFEMSetAssemblyType(SNES snes, MatAssemblyType type) {
216  SnesCtx *snes_ctx;
217  // PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
219  CHKERR SNESGetApplicationContext(snes, &snes_ctx);
220  snes_ctx->typeOfAssembly = type;
222 }
223 
225  SnesCtx *snes_ctx;
227  CHKERR SNESGetApplicationContext(snes, &snes_ctx);
228  snes_ctx->bH = bh;
230 }
231 
232 MoFEMErrorCode MoFEMSNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm,
233  SnesCtx *snes_ctx) {
235  auto &m_field = snes_ctx->mField;
236  auto problem_ptr = m_field.get_problem(snes_ctx->problemName);
237  auto fields_ptr = m_field.get_fields();
238  auto dofs = problem_ptr->numeredRowDofsPtr;
239 
240  std::vector<double> lnorms(fields_ptr->size(), 0),
241  norms(fields_ptr->size(), 0);
242 
243  Vec res;
244  CHKERR SNESGetFunction(snes, &res, NULL, NULL);
245 
246  const double *r;
247  CHKERR VecGetArrayRead(res, &r);
248  {
249  int f = 0;
250  for (auto fi : *fields_ptr) {
251  const auto lo_uid = FieldEntity::getLoBitNumberUId(fi->bitNumber);
252  const auto hi_uid = FieldEntity::getHiBitNumberUId(fi->bitNumber);
253  const auto hi = dofs->get<Unique_mi_tag>().upper_bound(hi_uid);
254  for (auto lo = dofs->get<Unique_mi_tag>().lower_bound(lo_uid); lo != hi;
255  ++lo) {
256  const DofIdx loc_idx = (*lo)->getPetscLocalDofIdx();
257  if (loc_idx >= 0 && loc_idx < problem_ptr->nbLocDofsRow) {
258  lnorms[f] += PetscRealPart(PetscSqr(r[loc_idx]));
259  }
260  }
261  ++f;
262  }
263  }
264  CHKERR VecRestoreArrayRead(res, &r);
265 
266  MPIU_Allreduce(&*lnorms.begin(), &*norms.begin(), lnorms.size(), MPIU_REAL,
267  MPIU_SUM, PetscObjectComm((PetscObject)snes));
268 
269  std::stringstream s;
270  int tl;
271  CHKERR PetscObjectGetTabLevel((PetscObject)snes, &tl);
272  for (auto t = 0; t != tl; ++t)
273  s << " ";
274  s << its << " Function norm " << boost::format("%14.12e") % (double)fgnorm
275  << " [";
276  {
277  int f = 0;
278  for (auto fi : *fields_ptr) {
279  if (f > 0)
280  s << ", ";
281  s << boost::format("%14.12e") % (double)PetscSqrtReal(norms[f]);
282  ++f;
283  }
284  s << "]";
285  }
286 
287  MOFEM_LOG("SNES_WORLD", Sev::inform) << s.str();
288 
290 }
291 
292 } // namespace MoFEM
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::KspMethod::CTX_KSPNONE
@ CTX_KSPNONE
Definition: LoopMethods.hpp:73
g
constexpr double g
Definition: shallow_wave.cpp:63
MoFEM::CoreInterface::problem_basic_method_postProcess
virtual MoFEMErrorCode problem_basic_method_postProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
MoFEM::CoreInterface::loop_finite_elements
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
MoFEM::SnesCtx::MOFEM_EVENT_SnesMat
PetscLogEvent MOFEM_EVENT_SnesMat
Log events to assemble tangent matrix.
Definition: SnesCtx.hpp:158
MoFEM::SnesCtx::vecAssembleSwitch
boost::movelib::unique_ptr< bool > vecAssembleSwitch
Definition: SnesCtx.hpp:155
MoFEM::PetscData::CtxSetB
static constexpr Switches CtxSetB
Definition: LoopMethods.hpp:38
MoFEM::Tools
Auxiliary tools.
Definition: Tools.hpp:19
MoFEM::PetscData::CtxSetA
static constexpr Switches CtxSetA
Definition: LoopMethods.hpp:37
MoFEM::SnesCtx::vErify
bool vErify
If true verify vector.
Definition: SnesCtx.hpp:24
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::CoreInterface::get_problem
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
MoFEM::SnesCtx::postProcess_Mat
BasicMethodsSequence postProcess_Mat
Definition: SnesCtx.hpp:36
MoFEM::KspMethod::CTX_SETFUNCTION
@ CTX_SETFUNCTION
Definition: LoopMethods.hpp:73
MoFEM::SnesMethod::CTX_SNESSETJACOBIAN
@ CTX_SNESSETJACOBIAN
Definition: LoopMethods.hpp:107
sdf.r
int r
Definition: sdf.py:8
MoFEM::SnesMethod::CTX_SNESSETFUNCTION
@ CTX_SNESSETFUNCTION
Definition: LoopMethods.hpp:107
ROW
@ ROW
Definition: definitions.h:136
MoFEM::SnesCtx::loops_to_do_Rhs
FEMethodsSequence loops_to_do_Rhs
Definition: SnesCtx.hpp:32
MoFEM::CoreInterface::problem_basic_method_preProcess
virtual MoFEMErrorCode problem_basic_method_preProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::SnesCtx::copyLoops
MoFEMErrorCode copyLoops(const SnesCtx &snes_ctx)
Copy sequences from other SNES contex.
Definition: SnesCtx.cpp:5
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::SnesCtx::zeroPreCondMatrixB
bool zeroPreCondMatrixB
Definition: SnesCtx.hpp:21
MoFEM::SnesCtx::mField
MoFEM::Interface & mField
database Interface
Definition: SnesCtx.hpp:15
MoFEM::SnesCtx::clearLoops
MoFEMErrorCode clearLoops()
Clear loops.
Definition: SnesCtx.cpp:16
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:107
double
bit
auto bit
set bit
Definition: hanging_node_approx.cpp:75
convert.type
type
Definition: convert.py:64
MoFEM::SnesCtx::problemName
std::string problemName
problem name
Definition: SnesCtx.hpp:18
MoFEM::FieldEntity::getLoBitNumberUId
static UId getLoBitNumberUId(const FieldBitNumber bit_number)
Definition: FieldEntsMultiIndices.hpp:222
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::SnesCtx::MOFEM_EVENT_SnesRhs
PetscLogEvent MOFEM_EVENT_SnesRhs
Log events to assemble residual.
Definition: SnesCtx.hpp:157
MoFEM::CoreInterface::cache_problem_entities
virtual MoFEMErrorCode cache_problem_entities(const std::string prb_name, CacheTupleWeakPtr cache_ptr)=0
Cache variables.
MoFEM::SnesMoFEMSetBehavior
MoFEMErrorCode SnesMoFEMSetBehavior(SNES snes, MoFEMTypes bh)
Set behavior if finite element in sequence does not exist.
Definition: SnesCtx.cpp:224
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:35
t
constexpr double t
plate stiffness
Definition: plate.cpp:58
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::VecManager
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:23
MoFEM::Problem::numeredRowDofsPtr
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredRowDofsPtr
store DOFs on rows for this problem
Definition: ProblemsMultiIndices.hpp:73
MoFEM::PetscData::CtxSetX
static constexpr Switches CtxSetX
Definition: LoopMethods.hpp:39
MoFEM::SnesCtx::loops_to_do_Mat
FEMethodsSequence loops_to_do_Mat
Definition: SnesCtx.hpp:30
MoFEM::SnesRhs
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:27
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::SnesCtx::matAssembleSwitch
boost::movelib::unique_ptr< bool > matAssembleSwitch
Definition: SnesCtx.hpp:156
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:66
MoFEM::SnesCtx::bH
MoFEMTypes bH
If set to MF_EXIST check if element exist, default MF_EXIST.
Definition: SnesCtx.hpp:20
MoFEM::MoFEMSNESMonitorFields
MoFEMErrorCode MoFEMSNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *snes_ctx)
Sens monitor printing residual field by field.
Definition: SnesCtx.cpp:232
MoFEMTypes
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
Definition: definitions.h:110
MoFEM::FieldEntity::getHiBitNumberUId
static UId getHiBitNumberUId(const FieldBitNumber bit_number)
Definition: FieldEntsMultiIndices.hpp:228
MoFEM::KspMethod::CTX_OPERATORS
@ CTX_OPERATORS
Definition: LoopMethods.hpp:73
MoFEM::PetscData::CtxSetF
static constexpr Switches CtxSetF
Definition: LoopMethods.hpp:36
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::SnesCtx
Interface for nonlinear (SNES) solver.
Definition: SnesCtx.hpp:13
MoFEM::Unique_mi_tag
Definition: TagMultiIndices.hpp:18
MoFEM::Problem
keeps basic data about problem
Definition: ProblemsMultiIndices.hpp:54
MoFEM::SnesMat
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:139
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
MoFEM::SnesMoFEMSetAssemblyType
MoFEMErrorCode SnesMoFEMSetAssemblyType(SNES snes, MatAssemblyType type)
Set assembly type at the end of SnesMat.
Definition: SnesCtx.cpp:215
MoFEM::Types::DofIdx
int DofIdx
Index of DOF.
Definition: Types.hpp:18