v0.14.0
TsCtx.cpp
Go to the documentation of this file.
1 
2 
3 namespace MoFEM {
4 
5 TsCtx::TsCtx(MoFEM::Interface &m_field, const std::string &problem_name)
6  : mField(m_field), moab(m_field.get_moab()), problemName(problem_name),
7  bH(MF_EXIST), zeroMatrix(true) {
8  PetscLogEventRegister("LoopTsIFunction", 0, &MOFEM_EVENT_TsCtxIFunction);
9  PetscLogEventRegister("LoopTsIJacobian", 0, &MOFEM_EVENT_TsCtxIJacobian);
10  PetscLogEventRegister("LoopTsRHSFunction", 0, &MOFEM_EVENT_TsCtxRHSFunction);
11  PetscLogEventRegister("LoopTsRHSJacobian", 0, &MOFEM_EVENT_TsCtxRHSJacobian);
12  PetscLogEventRegister("LoopTsMonitor", 0, &MOFEM_EVENT_TsCtxMonitor);
13  PetscLogEventRegister("LoopTsI2Function", 0, &MOFEM_EVENT_TsCtxI2Function);
14  PetscLogEventRegister("LoopTsI2Jacobian", 0, &MOFEM_EVENT_TsCtxI2Jacobian);
15 
16  if (!LogManager::checkIfChannelExist("TSWORLD")) {
17  auto core_log = logging::core::get();
18 
19  core_log->add_sink(
21  core_log->add_sink(
23  core_log->add_sink(
25 
26  LogManager::setLog("TSWORLD");
27  LogManager::setLog("TSSYNC");
28  LogManager::setLog("TSSELF");
29 
30  MOFEM_LOG_TAG("TSWORLD", "TS");
31  MOFEM_LOG_TAG("TSSYNC", "TS");
32  MOFEM_LOG_TAG("TSSELF", "TS");
33  }
34 }
35 
38  loopsIJacobian.clear();
39  loopsIFunction.clear();
40  loopsMonitor.clear();
41  loopsRHSJacobian.clear();
42  loopsRHSFunction.clear();
43  preProcessIJacobian.clear();
44  postProcessIJacobian.clear();
45  preProcessIFunction.clear();
46  postProcessIFunction.clear();
47  preProcessMonitor.clear();
48  postProcessMonitor.clear();
49  preProcessRHSJacobian.clear();
50  preProcessRHSFunction.clear();
51  postProcessRHSJacobian.clear();
52  postProcessRHSFunction.clear();
54 }
55 
56 PetscErrorCode TsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t, Vec F,
57  void *ctx) {
59  TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
60  PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
61  CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
62  CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
63  CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
64  CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
65  CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
66  ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
67 
68  auto zero_ghost_vec = [](Vec g) {
70  Vec l;
71  CHKERR VecGhostGetLocalForm(g, &l);
72  double *a;
73  CHKERR VecGetArray(l, &a);
74  int s;
75  CHKERR VecGetLocalSize(l, &s);
76  for (int i = 0; i != s; ++i)
77  a[i] = 0;
78  CHKERR VecRestoreArray(l, &a);
79  CHKERR VecGhostRestoreLocalForm(g, &l);
81  };
82  CHKERR zero_ghost_vec(F);
83 
84  ts_ctx->vecAssembleSwitch = boost::movelib::make_unique<bool>(true);
85 
86  int step;
87 #if PETSC_VERSION_GE(3, 8, 0)
88  CHKERR TSGetStepNumber(ts, &step);
89 #else
90  CHKERR TSGetTimeStepNumber(ts, &step);
91 #endif
92 
93  auto cache_ptr = boost::make_shared<CacheTuple>();
95 
96  auto set = [&](auto &fe) {
97  fe.ts = ts;
98  fe.ts_u = u;
99  fe.ts_u_t = u_t;
100  fe.ts_F = F;
101  fe.ts_t = t;
102  fe.ts_step = step;
103  fe.ts_ctx = TSMethod::CTX_TSSETIFUNCTION;
104  fe.snes_ctx = SnesMethod::CTX_SNESSETFUNCTION;
105  fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
106  fe.data_ctx = PetscData::CtxSetF | PetscData::CtxSetX |
108 
109  CHKERR TSGetSNES(ts, &fe.snes);
110  CHKERR SNESGetKSP(fe.snes, &fe.ksp);
111 
112  fe.cacheWeakPtr = cache_ptr;
113  CHK_THROW_MESSAGE(TSGetTimeStep(ts, &fe.ts_dt), "get time step failed");
114  };
115 
116  auto unset = [&](auto &fe) {
117  fe.ts_ctx = TSMethod::CTX_TSNONE;
118  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
119  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
120  fe.data_ctx = PetscData::CtxSetNone;
121  };
122 
123  // preprocess
124  for (auto &bit : ts_ctx->preProcessIFunction) {
125  bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
126  set(*bit);
128  *bit);
129  unset(*bit);
130  ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
131  }
132 
133  // fe loops
134  for (auto &lit : ts_ctx->loopsIFunction) {
135  lit.second->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
136  set(*lit.second);
138  *(lit.second), nullptr,
139  ts_ctx->bH, cache_ptr);
140  unset(*lit.second);
141  ts_ctx->vecAssembleSwitch = boost::move(lit.second->vecAssembleSwitch);
142  }
143 
144  // post process
145  for (auto &bit : ts_ctx->postProcessIFunction) {
146  bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
147  set(*bit);
149  *bit);
150  unset(*bit);
151  ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
152  }
153 
154  if (*ts_ctx->vecAssembleSwitch) {
155  CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
156  CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
157  CHKERR VecAssemblyBegin(F);
158  CHKERR VecAssemblyEnd(F);
159  }
160 
161  if(ts_ctx->tsDebugHook) {
162  CHKERR ts_ctx->tsDebugHook(ts, t, u, u_t, PETSC_NULL, F, ctx);
163  }
164 
165  PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
167 }
168 
169 PetscErrorCode TsSetIJacobian(TS ts, PetscReal t, Vec u, Vec u_t, PetscReal a,
170  Mat A, Mat B, void *ctx) {
172 
173  TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
174  PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
175  CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
176  CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
177  CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
178  CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
179  CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
180  ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
181  if (ts_ctx->zeroMatrix) {
182  CHKERR MatZeroEntries(B);
183  }
184  int step;
185 #if PETSC_VERSION_GE(3, 8, 0)
186  CHKERR TSGetStepNumber(ts, &step);
187 #else
188  CHKERR TSGetTimeStepNumber(ts, &step);
189 #endif
190 
192  boost::movelib::make_unique<bool>(ts_ctx->zeroMatrix);
193  auto cache_ptr = boost::make_shared<CacheTuple>();
195 
196  auto set = [&](auto &fe) {
197  fe.ts = ts;
198  fe.ts_u = u;
199  fe.ts_u_t = u_t;
200  fe.ts_A = A;
201  fe.ts_B = B;
202  fe.ts_t = t;
203  fe.ts_a = a;
204  fe.ts_step = step;
205  fe.ts_ctx = TSMethod::CTX_TSSETIJACOBIAN;
206  fe.snes_ctx = SnesMethod::CTX_SNESSETJACOBIAN;
207  fe.ksp_ctx = KspMethod::CTX_OPERATORS;
210 
211  CHKERR TSGetSNES(ts, &fe.snes);
212  CHKERR SNESGetKSP(fe.snes, &fe.ksp);
213 
214  fe.cacheWeakPtr = cache_ptr;
215  CHK_THROW_MESSAGE(TSGetTimeStep(ts, &fe.ts_dt), "get time step failed");
216  };
217 
218  auto unset = [&](auto &fe) {
219  fe.ts_ctx = TSMethod::CTX_TSNONE;
220  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
221  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
222  fe.data_ctx = PetscData::CtxSetNone;
223  };
224 
225  // preproces
226  for (auto &bit : ts_ctx->preProcessIJacobian) {
227  bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
228  set(*bit);
230  *bit);
231  unset(*bit);
232  ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
233  }
234 
235  for (auto &lit : ts_ctx->loopsIJacobian) {
236  lit.second->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
237  set(*lit.second);
239  *(lit.second), nullptr,
240  ts_ctx->bH, cache_ptr);
241  unset(*lit.second);
242  ts_ctx->matAssembleSwitch = boost::move(lit.second->matAssembleSwitch);
243  }
244 
245  // post process
246  for (auto &bit : ts_ctx->postProcessIJacobian) {
247  bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
248  set(*bit);
250  *bit);
251  unset(*bit);
252  ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
253  }
254 
255  if (*(ts_ctx->matAssembleSwitch)) {
256  CHKERR MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
257  CHKERR MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
258  }
259  PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
261 }
262 
263 PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u,
264  void *ctx) {
266  TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
267  PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxMonitor, 0, 0, 0, 0);
268  CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
269  CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
270  CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
271  ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
272 
273  auto cache_ptr = boost::make_shared<CacheTuple>();
275 
276  auto set = [&](auto &fe) {
277  fe.ts = ts;
278  fe.ts_u = u;
279  fe.ts_t = t;
280  fe.ts_step = step;
281  fe.ts_F = PETSC_NULL;
282  fe.ts_ctx = TSMethod::CTX_TSTSMONITORSET;
283  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
284  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
286 
287  CHKERR TSGetSNES(ts, &fe.snes);
288  CHKERR SNESGetKSP(fe.snes, &fe.ksp);
289 
290  fe.cacheWeakPtr = cache_ptr;
291  CHK_THROW_MESSAGE(TSGetTimeStep(ts, &fe.ts_dt), "get time step failed");
292  };
293 
294  auto unset = [&](auto &fe) {
295  fe.ts_ctx = TSMethod::CTX_TSNONE;
296  fe.data_ctx = PetscData::CtxSetNone;
297  };
298 
299  // preprocess
300  for (auto &bit : ts_ctx->preProcessMonitor) {
301  set(*bit);
303  *bit);
304  unset(*bit);
305  }
306 
307  for (auto &lit : ts_ctx->loopsMonitor) {
308  set(*lit.second);
310  *(lit.second), nullptr,
311  ts_ctx->bH, cache_ptr);
312  unset(*lit.second);
313  }
314 
315  // post process
316  for (auto &bit : ts_ctx->postProcessMonitor) {
317  set(*bit);
319  *bit);
320  unset(*bit);
321  }
322 
323  PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxMonitor, 0, 0, 0, 0);
325 }
326 
327 PetscErrorCode TsSetRHSFunction(TS ts, PetscReal t, Vec u, Vec F, void *ctx) {
329  TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
330  PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxRHSFunction, 0, 0, 0, 0);
331  CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
332  CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
333  CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
334  ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
335 
336  auto zero_ghost_vec = [](Vec g) {
338  Vec l;
339  CHKERR VecGhostGetLocalForm(g, &l);
340  double *a;
341  CHKERR VecGetArray(l, &a);
342  int s;
343  CHKERR VecGetLocalSize(l, &s);
344  for (int i = 0; i != s; ++i)
345  a[i] = 0;
346  CHKERR VecRestoreArray(l, &a);
347  CHKERR VecGhostRestoreLocalForm(g, &l);
349  };
350  CHKERR zero_ghost_vec(F);
351 
352  ts_ctx->vecAssembleSwitch = boost::movelib::make_unique<bool>(true);
353  auto cache_ptr = boost::make_shared<CacheTuple>();
355 
356  int step;
357 #if PETSC_VERSION_GE(3, 8, 0)
358  CHKERR TSGetStepNumber(ts, &step);
359 #else
360  CHKERR TSGetTimeStepNumber(ts, &step);
361 #endif
362 
363  auto set = [&](auto &fe) {
364  fe.ts_u = u;
365  fe.ts_F = F;
366  fe.ts_t = t;
367  fe.ts = ts;
368  fe.ts_step = step;
369  fe.ts_ctx = TSMethod::CTX_TSSETRHSFUNCTION;
370  fe.snes_ctx = SnesMethod::CTX_SNESSETFUNCTION;
371  fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
372  fe.data_ctx =
374  fe.cacheWeakPtr = cache_ptr;
375  CHK_THROW_MESSAGE(TSGetTimeStep(ts, &fe.ts_dt), "get time step failed");
376  };
377 
378  auto unset = [&](auto &fe) {
379  fe.ts_ctx = TSMethod::CTX_TSNONE;
380  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
381  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
382  fe.data_ctx = PetscData::CtxSetNone;
383  };
384 
385  for (auto &bit : ts_ctx->preProcessRHSFunction) {
386  bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
387  set(*bit);
389  *bit);
390  unset(*bit);
391  ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
392  }
393 
394  // fe loops
395  for (auto &lit : ts_ctx->loopsRHSFunction) {
396  lit.second->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
397  set(*lit.second);
399  *(lit.second), nullptr,
400  ts_ctx->bH, cache_ptr);
401  unset(*lit.second);
402  ts_ctx->vecAssembleSwitch = boost::move(lit.second->vecAssembleSwitch);
403  }
404 
405  // post process
406  for (auto &bit : ts_ctx->postProcessRHSFunction) {
407  bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
408  set(*bit);
410  *bit);
411  unset(*bit);
412  ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
413  }
414 
415  if (*ts_ctx->vecAssembleSwitch) {
416  CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
417  CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
418  CHKERR VecAssemblyBegin(F);
419  CHKERR VecAssemblyEnd(F);
420  }
421 
422  if (ts_ctx->tsDebugHook) {
423  CHKERR ts_ctx->tsDebugHook(ts, t, u, PETSC_NULL, PETSC_NULL, F, ctx);
424  }
425 
426  PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxRHSFunction, 0, 0, 0, 0);
428 }
429 
430 PetscErrorCode TsSetRHSJacobian(TS ts, PetscReal t, Vec u, Mat A, Mat B,
431  void *ctx) {
433  TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
434  PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxRHSJacobian, 0, 0, 0, 0);
435  CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
436  CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
437  CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
438  ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
439 
440  if (ts_ctx->zeroMatrix) {
441  CHKERR MatZeroEntries(B);
442  }
443 
445  boost::movelib::make_unique<bool>(ts_ctx->zeroMatrix);
446  auto cache_ptr = boost::make_shared<CacheTuple>();
448 
449  int step;
450 #if PETSC_VERSION_GE(3, 8, 0)
451  CHKERR TSGetStepNumber(ts, &step);
452 #else
453  CHKERR TSGetTimeStepNumber(ts, &step);
454 #endif
455 
456  auto set = [&](auto &fe) {
457  fe.ts_u = u;
458  fe.ts_A = A;
459  fe.ts_B = B;
460  fe.ts_t = t;
461  fe.ts_step = step;
462  fe.ts_ctx = TSMethod::CTX_TSSETRHSJACOBIAN;
463  fe.snes_ctx = SnesMethod::CTX_SNESSETJACOBIAN;
464  fe.ksp_ctx = KspMethod::CTX_OPERATORS;
467  fe.ts = ts;
468  fe.cacheWeakPtr = cache_ptr;
469  CHK_THROW_MESSAGE(TSGetTimeStep(ts, &fe.ts_dt), "get time step failed");
470  };
471 
472  auto unset = [&](auto &fe) {
473  fe.ts_ctx = TSMethod::CTX_TSNONE;
474  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
475  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
476  fe.data_ctx = PetscData::CtxSetNone;
477  };
478 
479  // preprocess
480  for (auto &bit : ts_ctx->preProcessRHSJacobian) {
481  bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
482  set(*bit);
484  *bit);
485  unset(*bit);
486  ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
487  }
488 
489  // fe loops
490  for (auto &lit : ts_ctx->loopsRHSJacobian) {
491  lit.second->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
492  set(*lit.second);
494  *(lit.second), nullptr,
495  ts_ctx->bH, cache_ptr);
496  unset(*lit.second);
497  ts_ctx->matAssembleSwitch = boost::move(lit.second->matAssembleSwitch);
498  }
499 
500  // post process
501  for (auto &bit : ts_ctx->postProcessRHSJacobian) {
502  bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
503  set(*bit);
505  *bit);
506  unset(*bit);
507  ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
508  }
509 
510  if (*(ts_ctx->matAssembleSwitch)) {
511  CHKERR MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
512  CHKERR MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
513  }
514 
515  PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxRHSJacobian, 0, 0, 0, 0);
517 }
518 
519 PetscErrorCode TsSetI2Jacobian(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt,
520  PetscReal a, PetscReal aa, Mat A, Mat B,
521  void *ctx) {
523 
524  TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
525  PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxI2Function, 0, 0, 0, 0);
526  CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
527  CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
528  CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
529  CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
530  CHKERR VecGhostUpdateBegin(u_tt, INSERT_VALUES, SCATTER_FORWARD);
531  CHKERR VecGhostUpdateEnd(u_tt, INSERT_VALUES, SCATTER_FORWARD);
532  CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
533  ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
534  if (ts_ctx->zeroMatrix) {
535  CHKERR MatZeroEntries(B);
536  }
537  int step;
538 #if PETSC_VERSION_GE(3, 8, 0)
539  CHKERR TSGetStepNumber(ts, &step);
540 #else
541  CHKERR TSGetTimeStepNumber(ts, &step);
542 #endif
543 
545  boost::movelib::make_unique<bool>(ts_ctx->zeroMatrix);
546  auto cache_ptr = boost::make_shared<CacheTuple>();
548 
549  auto set = [&](auto &fe) {
550  fe.ts_u = u;
551  fe.ts_u_t = u_t;
552  fe.ts_u_tt = u_tt;
553  fe.ts_A = A;
554  fe.ts_B = B;
555  fe.ts_t = t;
556  fe.ts_a = a;
557  fe.ts_aa = aa;
558  fe.ts_step = step;
559 
560  fe.ts_ctx = TSMethod::CTX_TSSETIJACOBIAN;
561  fe.snes_ctx = SnesMethod::CTX_SNESSETJACOBIAN;
562  fe.ksp_ctx = KspMethod::CTX_OPERATORS;
566 
567  CHKERR TSGetSNES(ts, &fe.snes);
568  CHKERR SNESGetKSP(fe.snes, &fe.ksp);
569 
570  fe.ts = ts;
571  fe.cacheWeakPtr = cache_ptr;
572  CHK_THROW_MESSAGE(TSGetTimeStep(ts, &fe.ts_dt), "get time step failed");
573  };
574 
575  auto unset = [&](auto &fe) {
576  fe.ts_ctx = TSMethod::CTX_TSNONE;
577  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
578  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
579  fe.data_ctx = PetscData::CtxSetNone;
580  };
581 
582  // preprocess
583  for (auto &bit : ts_ctx->preProcessIJacobian) {
584  bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
585  set(*bit);
587  *bit);
588  unset(*bit);
589  ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
590  }
591 
592  for (auto &lit : ts_ctx->loopsIJacobian) {
593  lit.second->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
594  set(*lit.second);
596  *(lit.second), nullptr,
597  ts_ctx->bH, cache_ptr);
598  unset(*lit.second);
599  ts_ctx->matAssembleSwitch = boost::move(lit.second->matAssembleSwitch);
600  }
601 
602  // post process
603  for (auto &bit : ts_ctx->postProcessIJacobian) {
604  bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
605  set(*bit);
607  *bit);
608  unset(*bit);
609  ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
610  }
611 
612  if (*(ts_ctx->matAssembleSwitch)) {
613  CHKERR MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
614  CHKERR MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
615  }
616  PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxI2Function, 0, 0, 0, 0);
618 }
619 
620 PetscErrorCode TsSetI2Function(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt,
621  Vec F, void *ctx) {
623  TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
624  PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
625  CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
626  CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
627  CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
628  CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
629  CHKERR VecGhostUpdateBegin(u_tt, INSERT_VALUES, SCATTER_FORWARD);
630  CHKERR VecGhostUpdateEnd(u_tt, INSERT_VALUES, SCATTER_FORWARD);
631  CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
632  ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
633 
634  auto zero_ghost_vec = [](Vec g) {
636  Vec l;
637  CHKERR VecGhostGetLocalForm(g, &l);
638  double *a;
639  CHKERR VecGetArray(l, &a);
640  int s;
641  CHKERR VecGetLocalSize(l, &s);
642  for (int i = 0; i != s; ++i)
643  a[i] = 0;
644  CHKERR VecRestoreArray(l, &a);
645  CHKERR VecGhostRestoreLocalForm(g, &l);
647  };
648  CHKERR zero_ghost_vec(F);
649 
650  ts_ctx->vecAssembleSwitch = boost::movelib::make_unique<bool>(true);
651  auto cache_ptr = boost::make_shared<CacheTuple>();
653 
654  int step;
655 #if PETSC_VERSION_GE(3, 8, 0)
656  CHKERR TSGetStepNumber(ts, &step);
657 #else
658  CHKERR TSGetTimeStepNumber(ts, &step);
659 #endif
660 
661  auto set = [&](auto &fe) {
662  fe.ts_u = u;
663  fe.ts_u_t = u_t;
664  fe.ts_u_tt = u_tt;
665  fe.ts_F = F;
666  fe.ts_t = t;
667  fe.ts_step = step;
668  fe.ts_ctx = TSMethod::CTX_TSSETIFUNCTION;
669  fe.snes_ctx = SnesMethod::CTX_SNESSETFUNCTION;
670  fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
671  fe.data_ctx = PetscData::CtxSetF | PetscData::CtxSetX |
674 
675  CHKERR TSGetSNES(ts, &fe.snes);
676  CHKERR SNESGetKSP(fe.snes, &fe.ksp);
677 
678  fe.ts = ts;
679  fe.cacheWeakPtr = cache_ptr;
680  CHK_THROW_MESSAGE(TSGetTimeStep(ts, &fe.ts_dt), "get time step failed");
681  };
682 
683  auto unset = [&](auto &fe) {
684  fe.ts_ctx = TSMethod::CTX_TSNONE;
685  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
686  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
687  fe.data_ctx = PetscData::CtxSetNone;
688  };
689 
690  // preprocess
691  for (auto &bit : ts_ctx->preProcessIFunction) {
692  bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
693  set(*bit);
695  *bit);
696  unset(*bit);
697  ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
698  }
699 
700  // fe loops
701  for (auto &lit : ts_ctx->loopsIFunction) {
702  lit.second->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
703  set(*lit.second);
705  *(lit.second), nullptr,
706  ts_ctx->bH, cache_ptr);
707  unset(*lit.second);
708  ts_ctx->vecAssembleSwitch = boost::move(lit.second->vecAssembleSwitch);
709  }
710 
711  // post process
712  for (auto &bit : ts_ctx->postProcessIFunction) {
713  bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
714  set(*bit);
716  *bit);
717  unset(*bit);
718  ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
719  }
720 
721  if (*ts_ctx->vecAssembleSwitch) {
722  CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
723  CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
724  CHKERR VecAssemblyBegin(F);
725  CHKERR VecAssemblyEnd(F);
726  }
727 
728  if (ts_ctx->tsDebugHook) {
729  CHKERR ts_ctx->tsDebugHook(ts, t, u, u_t, u_tt, F, ctx);
730  }
731 
732  PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
734 }
735 
736 /** \brief Custom TSAdaptivity in MoFEM
737  *
738  * \code
739  * CHKERR TSAdaptRegister(TSADAPTMOFEM, TSAdaptCreateMoFEM);
740  * TSAdapt adapt;
741  * CHKERR TSGetAdapt(solver, &adapt);
742  * CHKERR TSAdaptSetType(adapt, TSADAPTMOFEM);
743  * \endcode
744  *
745  */
746 struct TSAdaptMoFEM {
747 
748  TSAdaptMoFEM();
749 
750  double alpha; //< step reduction if divergence
751  double gamma; //< adaptivity exponent
752  int desiredIt; //< desired number of iterations
753  PetscBool offApat; //< off adaptivity
754 };
755 
757  : alpha(0.75), gamma(0.5), desiredIt(6), offApat(PETSC_FALSE) {
758  CHKERR PetscOptionsGetScalar("", "-ts_mofem_adapt_alpha", &alpha, PETSC_NULL);
759  CHKERR PetscOptionsGetScalar("", "-ts_mofem_adapt_gamma", &gamma, PETSC_NULL);
760  CHKERR PetscOptionsGetInt("", "-ts_mofem_adapt_desired_it", &desiredIt,
761  PETSC_NULL);
762  CHKERR PetscOptionsGetBool("", "-ts_mofem_adapt_off", &offApat, PETSC_NULL);
763 
764  MOFEM_LOG("TSWORLD", Sev::inform)
765  << "TS adaptivity: alpha = " << alpha << ", gamma = " << gamma
766  << ", desiredIt = " << desiredIt << ", offAdapt = " << offApat;
767 }
768 
769 PetscErrorCode TSAdaptChooseMoFEM(TSAdapt adapt, TS ts, PetscReal h,
770  PetscInt *next_sc, PetscReal *next_h,
771  PetscBool *accept, PetscReal *wlte,
772  PetscReal *wltea, PetscReal *wlter) {
774 
775  auto ts_adapt_mofem = boost::make_shared<TSAdaptMoFEM>();
776 
777  *next_sc = 0; /* Reuse the same order scheme */
778  *wlte = -1; /* Weighted local truncation error was not evaluated */
779  *wltea = -1; /* Weighted absolute local truncation error is not used */
780  *wlter = -1; /* Weighted relative local truncation error is not used */
781 
782  *accept = PETSC_TRUE;
783  *next_h = h; /* Reuse the old step */
784 
785  if (!ts_adapt_mofem->offApat) {
786 
787  SNES snes;
788  CHKERR TSGetSNES(ts, &snes);
789 
790  SNESConvergedReason reason;
791  CHKERR SNESGetConvergedReason(snes, &reason);
792 
793  int it;
794  CHKERR SNESGetIterationNumber(snes, &it);
795 
796  if (reason < 0) {
797  h *= ts_adapt_mofem->alpha;
798  *next_h = h;
799  MOFEM_LOG_C(
800  "TSWORLD", Sev::warning,
801  "\tDiverged set step length: it = %d, h = %3.4g set h = %3.4g \n", it,
802  h, *next_h);
803  } else if (reason > 0) {
804  h *= pow((static_cast<double>(ts_adapt_mofem->desiredIt) /
805  static_cast<double>(it + 1)),
806  static_cast<double>(ts_adapt_mofem->gamma));
807  *next_h = PetscClipInterval(h, adapt->dt_min, adapt->dt_max);
808  MOFEM_LOG_C(
809  "TSWORLD", Sev::inform,
810  "\tConverged set step length: it = %d, h = %3.4g set h = %3.4g \n",
811  it, h, *next_h);
812  }
813  }
814 
816 }
817 
818 PetscErrorCode TSAdaptResetMoFEM(TSAdapt adapt) {
819  PetscFunctionBegin;
820  PetscFunctionReturn(0);
821 }
822 
823 PetscErrorCode TSAdaptDestroyMoFEM(TSAdapt adapt) {
825  CHKERR TSAdaptResetMoFEM(adapt);
827 }
828 
829 PetscErrorCode TSAdaptCreateMoFEM(TSAdapt adapt) {
830  PetscFunctionBegin;
831  adapt->ops->choose = TSAdaptChooseMoFEM;
832  adapt->ops->reset = TSAdaptResetMoFEM;
833  adapt->ops->destroy = TSAdaptDestroyMoFEM;
834  PetscFunctionReturn(0);
835 }
836 
837 } // namespace MoFEM
MoFEM::TsCtx::preProcessRHSJacobian
BasicMethodsSequence preProcessRHSJacobian
Definition: TsCtx.hpp:40
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::TsCtx::clearLoops
MoFEMErrorCode clearLoops()
Clear loops.
Definition: TsCtx.cpp:36
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::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::PetscData::CtxSetB
static constexpr Switches CtxSetB
Definition: LoopMethods.hpp:39
MoFEM::PetscData::CtxSetA
static constexpr Switches CtxSetA
Definition: LoopMethods.hpp:38
MoFEM::TsSetIFunction
PetscErrorCode TsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t, Vec F, void *ctx)
Set IFunction for TS solver.
Definition: TsCtx.cpp:56
MoFEM::TsCtx::MOFEM_EVENT_TsCtxIJacobian
PetscLogEvent MOFEM_EVENT_TsCtxIJacobian
Definition: TsCtx.hpp:215
MoFEM::TSAdaptDestroyMoFEM
PetscErrorCode TSAdaptDestroyMoFEM(TSAdapt adapt)
Definition: TsCtx.cpp:823
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:609
MoFEM::TsCtx::MOFEM_EVENT_TsCtxI2Function
PetscLogEvent MOFEM_EVENT_TsCtxI2Function
Definition: TsCtx.hpp:217
MoFEM::TsCtx::MOFEM_EVENT_TsCtxIFunction
PetscLogEvent MOFEM_EVENT_TsCtxIFunction
Definition: TsCtx.hpp:214
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::TSAdaptMoFEM::alpha
double alpha
Definition: TsCtx.cpp:750
MoFEM::TsCtx::MOFEM_EVENT_TsCtxMonitor
PetscLogEvent MOFEM_EVENT_TsCtxMonitor
Definition: TsCtx.hpp:216
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
MoFEM::TSMethod::CTX_TSTSMONITORSET
@ CTX_TSTSMONITORSET
Definition: LoopMethods.hpp:148
MoFEM::TsSetI2Jacobian
PetscErrorCode TsSetI2Jacobian(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt, PetscReal a, PetscReal aa, Mat A, Mat B, void *ctx)
Calculation Jacobian for second order PDE in time.
Definition: TsCtx.cpp:519
MoFEM::TsCtx::vecAssembleSwitch
boost::movelib::unique_ptr< bool > vecAssembleSwitch
Definition: TsCtx.hpp:220
MoFEM::TSMethod::CTX_TSSETRHSJACOBIAN
@ CTX_TSSETRHSJACOBIAN
Definition: LoopMethods.hpp:145
MoFEM::TSAdaptMoFEM
Custom TSAdaptivity in MoFEM.
Definition: TsCtx.cpp:746
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
ts_ctx
MoFEM::TsCtx * ts_ctx
Definition: level_set.cpp:1932
MoFEM::TsCtx::loopsIJacobian
FEMethodsSequence loopsIJacobian
Definition: TsCtx.hpp:29
MoFEM::TSAdaptMoFEM::offApat
PetscBool offApat
Definition: TsCtx.cpp:753
MoFEM::KspMethod::CTX_SETFUNCTION
@ CTX_SETFUNCTION
Definition: LoopMethods.hpp:76
MoFEM::SnesMethod::CTX_SNESSETJACOBIAN
@ CTX_SNESSETJACOBIAN
Definition: LoopMethods.hpp:110
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::TsCtx::bH
MoFEMTypes bH
If set to MF_EXIST check if element exist.
Definition: TsCtx.hpp:23
MoFEM::SnesMethod::CTX_SNESSETFUNCTION
@ CTX_SNESSETFUNCTION
Definition: LoopMethods.hpp:110
MoFEM::TsCtx::problemName
std::string problemName
Definition: TsCtx.hpp:22
MoFEM::TsCtx::loopsRHSJacobian
FEMethodsSequence loopsRHSJacobian
Definition: TsCtx.hpp:32
MoFEM::TsMonitorSet
PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Set monitor for TS solver.
Definition: TsCtx.cpp:263
MoFEM::TSAdaptResetMoFEM
PetscErrorCode TSAdaptResetMoFEM(TSAdapt adapt)
Definition: TsCtx.cpp:818
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.
MoFEM::TsCtx::postProcessIJacobian
BasicMethodsSequence postProcessIJacobian
Definition: TsCtx.hpp:35
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
a
constexpr double a
Definition: approx_sphere.cpp:30
MoFEM::TsCtx::mField
MoFEM::Interface & mField
Definition: TsCtx.hpp:19
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
MoFEM::TSAdaptMoFEM::gamma
double gamma
Definition: TsCtx.cpp:751
MoFEM::SnesMethod::CTX_SNESNONE
@ CTX_SNESNONE
Definition: LoopMethods.hpp:110
MoFEM::TsCtx::postProcessMonitor
BasicMethodsSequence postProcessMonitor
Definition: TsCtx.hpp:39
bit
auto bit
set bit
Definition: hanging_node_approx.cpp:75
MoFEM::TSMethod::CTX_TSSETIJACOBIAN
@ CTX_TSSETIJACOBIAN
Definition: LoopMethods.hpp:147
MoFEM::TsCtx::loopsMonitor
FEMethodsSequence loopsMonitor
Definition: TsCtx.hpp:31
MoFEM::LogManager::getStrmSync
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
Definition: LogManager.cpp:348
h
double h
Definition: photon_diffusion.cpp:60
MoFEM::PetscData::CtxSetX_T
static constexpr Switches CtxSetX_T
Definition: LoopMethods.hpp:42
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::TsCtx::zeroMatrix
bool zeroMatrix
Definition: TsCtx.hpp:49
MoFEM::TsCtx
Interface for Time Stepping (TS) solver.
Definition: TsCtx.hpp:17
MoFEM::TSMethod::CTX_TSNONE
@ CTX_TSNONE
Definition: LoopMethods.hpp:149
MoFEM::CoreInterface::cache_problem_entities
virtual MoFEMErrorCode cache_problem_entities(const std::string prb_name, CacheTupleWeakPtr cache_ptr)=0
Cache variables.
MoFEM::TsCtx::matAssembleSwitch
boost::movelib::unique_ptr< bool > matAssembleSwitch
Definition: TsCtx.hpp:221
MoFEM::TsCtx::preProcessMonitor
BasicMethodsSequence preProcessMonitor
Definition: TsCtx.hpp:38
MoFEM::TsCtx::postProcessIFunction
BasicMethodsSequence postProcessIFunction
Definition: TsCtx.hpp:37
MoFEM::TsCtx::MOFEM_EVENT_TsCtxRHSFunction
PetscLogEvent MOFEM_EVENT_TsCtxRHSFunction
Definition: TsCtx.hpp:212
MoFEM::TsCtx::preProcessIFunction
BasicMethodsSequence preProcessIFunction
Definition: TsCtx.hpp:36
MoFEM::TsCtx::loopsIFunction
FEMethodsSequence loopsIFunction
Definition: TsCtx.hpp:30
MoFEM::TsSetIJacobian
PetscErrorCode TsSetIJacobian(TS ts, PetscReal t, Vec u, Vec u_t, PetscReal a, Mat A, Mat B, void *ctx)
Set function evaluating jacobian in TS solver.
Definition: TsCtx.cpp:169
MoFEM::TsCtx::preProcessIJacobian
BasicMethodsSequence preProcessIJacobian
Definition: TsCtx.hpp:34
MoFEM::LogManager::getStrmSelf
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
Definition: LogManager.cpp:340
MoFEM::PetscData::CtxSetTime
static constexpr Switches CtxSetTime
Definition: LoopMethods.hpp:44
MoFEM::PetscData::CtxSetNone
static constexpr Switches CtxSetNone
Definition: LoopMethods.hpp:36
MoFEM::TSAdaptMoFEM::TSAdaptMoFEM
TSAdaptMoFEM()
Definition: TsCtx.cpp:756
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
t
constexpr double t
plate stiffness
Definition: plate.cpp:58
MoFEM::TsCtx::loopsRHSFunction
FEMethodsSequence loopsRHSFunction
Definition: TsCtx.hpp:33
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::TsCtx::MOFEM_EVENT_TsCtxRHSJacobian
PetscLogEvent MOFEM_EVENT_TsCtxRHSJacobian
Definition: TsCtx.hpp:213
MoFEM::VecManager
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:23
MoFEM::TSAdaptChooseMoFEM
PetscErrorCode TSAdaptChooseMoFEM(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept, PetscReal *wlte, PetscReal *wltea, PetscReal *wlter)
Definition: TsCtx.cpp:769
MoFEM::PetscData::CtxSetX_TT
static constexpr Switches CtxSetX_TT
Definition: LoopMethods.hpp:43
MoFEM::TsCtx::postProcessRHSFunction
BasicMethodsSequence postProcessRHSFunction
Definition: TsCtx.hpp:43
MoFEM::PetscData::CtxSetX
static constexpr Switches CtxSetX
Definition: LoopMethods.hpp:40
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::TsCtx::MOFEM_EVENT_TsCtxI2Jacobian
PetscLogEvent MOFEM_EVENT_TsCtxI2Jacobian
Definition: TsCtx.hpp:218
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:64
MoFEM::TSAdaptMoFEM::desiredIt
int desiredIt
Definition: TsCtx.cpp:752
MoFEM::TsCtx::postProcessRHSJacobian
BasicMethodsSequence postProcessRHSJacobian
Definition: TsCtx.hpp:42
MoFEM::TsSetRHSJacobian
PetscErrorCode TsSetRHSJacobian(TS ts, PetscReal t, Vec u, Mat A, Mat B, void *ctx)
TS solver function.
Definition: TsCtx.cpp:430
MoFEM::TSMethod::CTX_TSSETIFUNCTION
@ CTX_TSSETIFUNCTION
Definition: LoopMethods.hpp:146
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::PetscOptionsGetScalar
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:162
MoFEM::TsSetRHSFunction
PetscErrorCode TsSetRHSFunction(TS ts, PetscReal t, Vec u, Vec F, void *ctx)
TS solver function.
Definition: TsCtx.cpp:327
MoFEM::LogManager::setLog
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389
MoFEM::TsSetI2Function
PetscErrorCode TsSetI2Function(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt, Vec F, void *ctx)
Calculation the right hand side for second order PDE in time.
Definition: TsCtx.cpp:620
MoFEM::TsCtx::TsCtx
TsCtx(MoFEM::Interface &m_field, const std::string &problem_name)
Definition: TsCtx.cpp:5
MF_EXIST
@ MF_EXIST
Definition: definitions.h:113
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::TsCtx::preProcessRHSFunction
BasicMethodsSequence preProcessRHSFunction
Definition: TsCtx.hpp:41
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::TsCtx::tsDebugHook
boost::function< MoFEMErrorCode(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt, Vec F, void *ctx)> tsDebugHook
Definition: TsCtx.hpp:47
MoFEM::TSMethod::CTX_TSSETRHSFUNCTION
@ CTX_TSSETRHSFUNCTION
Definition: LoopMethods.hpp:144
F
@ F
Definition: free_surface.cpp:396
MoFEM::PetscOptionsGetBool
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:182
MoFEM::TSAdaptCreateMoFEM
PetscErrorCode TSAdaptCreateMoFEM(TSAdapt adapt)
Craete MOFEM adapt.
Definition: TsCtx.cpp:829