v0.9.1
TsCtx.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 TsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t, Vec F,
18  void *ctx) {
20  TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
21  PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
22  CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
23  CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
24  CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
25  CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
26  CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
27  ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
28 
29  auto zero_ghost_vec = [](Vec g) {
31  Vec l;
32  CHKERR VecGhostGetLocalForm(g, &l);
33  double *a;
34  CHKERR VecGetArray(l, &a);
35  int s;
36  CHKERR VecGetLocalSize(l, &s);
37  for (int i = 0; i != s; ++i)
38  a[i] = 0;
39  CHKERR VecRestoreArray(l, &a);
40  CHKERR VecGhostRestoreLocalForm(g, &l);
42  };
43  CHKERR zero_ghost_vec(F);
44 
45  ts_ctx->vecAssembleSwitch = boost::movelib::make_unique<bool>(true);
46 
47  int step;
48 #if PETSC_VERSION_GE(3, 8, 0)
49  CHKERR TSGetStepNumber(ts, &step);
50 #else
51  CHKERR TSGetTimeStepNumber(ts, &step);
52 #endif
53 
54  auto set = [&](auto &fe) {
55  fe.ts = ts;
56  fe.ts_u = u;
57  fe.ts_u_t = u_t;
58  fe.ts_F = F;
59  fe.ts_t = t;
60  fe.ts_step = step;
61  fe.ts_ctx = TSMethod::CTX_TSSETIFUNCTION;
62  fe.snes_ctx = SnesMethod::CTX_SNESSETFUNCTION;
63  fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
64  fe.data_ctx = PetscData::CtxSetF | PetscData::CtxSetX |
66  };
67 
68  auto unset = [&](auto &fe) {
69  fe.ts_ctx = TSMethod::CTX_TSNONE;
70  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
71  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
72  fe.data_ctx = PetscData::CtxSetNone;
73  };
74 
75  // preprocess
76  for (auto &bit : ts_ctx->preProcess_IFunction) {
77  bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
78  set(*bit);
80  *bit);
81  unset(*bit);
82  ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
83  }
84 
85  // fe loops
86  for (auto &lit : ts_ctx->loops_to_do_IFunction) {
87  lit.second->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
88  set(*lit.second);
90  ts_ctx->problemName, lit.first, *(lit.second), nullptr, ts_ctx->bH);
91  unset(*lit.second);
92  ts_ctx->vecAssembleSwitch = boost::move(lit.second->vecAssembleSwitch);
93  }
94 
95  // post process
96  for (auto &bit : ts_ctx->postProcess_IFunction) {
97  bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
98  set(*bit);
100  *bit);
101  unset(*bit);
102  ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
103  }
104 
105  if (*ts_ctx->vecAssembleSwitch) {
106  CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
107  CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
108  CHKERR VecAssemblyBegin(F);
109  CHKERR VecAssemblyEnd(F);
110  }
111 
112  PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
114 }
115 
116 PetscErrorCode TsSetIJacobian(TS ts, PetscReal t, Vec u, Vec u_t, PetscReal a,
117  Mat A, Mat B, void *ctx) {
119 
120  TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
121  PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
122  CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
123  CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
124  CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
125  CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
126  CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
127  ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
128  if (ts_ctx->zeroMatrix) {
129  CHKERR MatZeroEntries(B);
130  }
131  int step;
132 #if PETSC_VERSION_GE(3, 8, 0)
133  CHKERR TSGetStepNumber(ts, &step);
134 #else
135  CHKERR TSGetTimeStepNumber(ts, &step);
136 #endif
137 
138  ts_ctx->matAssembleSwitch =
139  boost::movelib::make_unique<bool>(ts_ctx->zeroMatrix);
140 
141  auto set = [&](auto &fe) {
142  fe.ts = ts;
143  fe.ts_u = u;
144  fe.ts_u_t = u_t;
145  fe.ts_A = A;
146  fe.ts_B = B;
147  fe.ts_t = t;
148  fe.ts_a = a;
149  fe.ts_step = step;
150  fe.ts_ctx = TSMethod::CTX_TSSETIJACOBIAN;
151  fe.snes_ctx = SnesMethod::CTX_SNESSETJACOBIAN;
152  fe.ksp_ctx = KspMethod::CTX_OPERATORS;
153  fe.data_ctx = PetscData::CtxSetA | PetscData::CtxSetB |
156  };
157 
158  auto unset = [&](auto &fe) {
159  fe.ts_ctx = TSMethod::CTX_TSNONE;
160  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
161  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
162  fe.data_ctx = PetscData::CtxSetNone;
163  };
164 
165  // preproces
166  for (auto &bit : ts_ctx->preProcess_IJacobian) {
167  bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
168  set(*bit);
170  *bit);
171  unset(*bit);
172  ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
173  }
174 
175  for (auto &lit : ts_ctx->loops_to_do_IJacobian) {
176  lit.second->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
177  set(*lit.second);
179  ts_ctx->problemName, lit.first, *(lit.second), nullptr, ts_ctx->bH);
180  unset(*lit.second);
181  ts_ctx->matAssembleSwitch = boost::move(lit.second->matAssembleSwitch);
182  }
183 
184  // post process
185  for (auto &bit : ts_ctx->postProcess_IJacobian) {
186  bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
187  set(*bit);
189  *bit);
190  unset(*bit);
191  ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
192  }
193 
194  if (ts_ctx->matAssembleSwitch) {
195  CHKERR MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
196  CHKERR MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
197  }
198  PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
200 }
201 
202 PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u,
203  void *ctx) {
205  TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
206  PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxMonitor, 0, 0, 0, 0);
207  CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
208  CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
209  CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
210  ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
211 
212  auto set = [&](auto &fe) {
213  fe.ts = ts;
214  fe.ts_u = u;
215  fe.ts_t = t;
216  fe.ts_step = step;
217  fe.ts_F = PETSC_NULL;
218  fe.ts_ctx = TSMethod::CTX_TSTSMONITORSET;
219  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
220  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
222  };
223 
224  auto unset = [&](auto &fe) {
225  fe.ts_ctx = TSMethod::CTX_TSNONE;
226  fe.data_ctx = PetscData::CtxSetNone;
227  };
228 
229  // preproces
230  for (auto &bit : ts_ctx->preProcess_Monitor) {
231  set(*bit);
233  *bit);
234  unset(*bit);
235  }
236 
237  for (auto &lit : ts_ctx->loops_to_do_Monitor) {
238  set(*lit.second);
240  ts_ctx->problemName, lit.first, *(lit.second), nullptr, ts_ctx->bH);
241  unset(*lit.second);
242  }
243 
244  // post process
245  for (auto &bit : ts_ctx->postProcess_Monitor) {
246  set(*bit);
248  *bit);
249  unset(*bit);
250  }
251 
252  PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxMonitor, 0, 0, 0, 0);
254 }
255 
256 PetscErrorCode TsSetRHSFunction(TS ts, PetscReal t, Vec u, Vec F, void *ctx) {
258  TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
259  PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxRHSFunction, 0, 0, 0, 0);
260  CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
261  CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
262  CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
263  ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
264 
265  auto zero_ghost_vec = [](Vec g) {
267  Vec l;
268  CHKERR VecGhostGetLocalForm(g, &l);
269  double *a;
270  CHKERR VecGetArray(l, &a);
271  int s;
272  CHKERR VecGetLocalSize(l, &s);
273  for (int i = 0; i != s; ++i)
274  a[i] = 0;
275  CHKERR VecRestoreArray(l, &a);
276  CHKERR VecGhostRestoreLocalForm(g, &l);
278  };
279  CHKERR zero_ghost_vec(F);
280 
281  ts_ctx->vecAssembleSwitch = boost::movelib::make_unique<bool>(true);
282 
283  int step;
284 #if PETSC_VERSION_GE(3, 8, 0)
285  CHKERR TSGetStepNumber(ts, &step);
286 #else
287  CHKERR TSGetTimeStepNumber(ts, &step);
288 #endif
289 
290  auto set = [&](auto &fe) {
291  fe.ts_u = u;
292  fe.ts_F = F;
293  fe.ts_t = t;
294  fe.ts = ts;
295  fe.ts_step = step;
296  fe.ts_ctx = TSMethod::CTX_TSSETRHSFUNCTION;
297  fe.snes_ctx = SnesMethod::CTX_SNESSETFUNCTION;
298  fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
299  fe.data_ctx =
301  };
302 
303  auto unset = [&](auto &fe) {
304  fe.ts_ctx = TSMethod::CTX_TSNONE;
305  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
306  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
307  fe.data_ctx = PetscData::CtxSetNone;
308  };
309 
310  for (auto &bit : ts_ctx->preProcess_RHSJacobian) {
311  bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
312  set(*bit);
314  *bit);
315  unset(*bit);
316  ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
317  }
318 
319  // fe loops
320  for (auto &lit : ts_ctx->loops_to_do_RHSFunction) {
321  lit.second->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
322  set(*lit.second);
324  ts_ctx->problemName, lit.first, *(lit.second), nullptr, ts_ctx->bH);
325  unset(*lit.second);
326  ts_ctx->vecAssembleSwitch = boost::move(lit.second->vecAssembleSwitch);
327  }
328 
329  // post process
330  for (auto &bit : ts_ctx->postProcess_RHSJacobian) {
331  bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
332  set(*bit);
334  *bit);
335  unset(*bit);
336  ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
337  }
338 
339  if (*ts_ctx->vecAssembleSwitch) {
340  CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
341  CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
342  CHKERR VecAssemblyBegin(F);
343  CHKERR VecAssemblyEnd(F);
344  }
345 
346  PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxRHSFunction, 0, 0, 0, 0);
348 }
349 
350 PetscErrorCode TsSetRHSJacobian(TS ts, PetscReal t, Vec u, Mat A, Mat B,
351  void *ctx) {
353  TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
354  PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxRHSJacobian, 0, 0, 0, 0);
355  CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
356  CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
357  CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
358  ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
359 
360  if (ts_ctx->zeroMatrix) {
361  CHKERR MatZeroEntries(B);
362  }
363 
364  ts_ctx->matAssembleSwitch =
365  boost::movelib::make_unique<bool>(ts_ctx->zeroMatrix);
366 
367  int step;
368 #if PETSC_VERSION_GE(3, 8, 0)
369  CHKERR TSGetStepNumber(ts, &step);
370 #else
371  CHKERR TSGetTimeStepNumber(ts, &step);
372 #endif
373 
374  auto set = [&](auto &fe) {
375  fe.ts_u = u;
376  fe.ts_A = A;
377  fe.ts_B = B;
378  fe.ts_t = t;
379  fe.ts_step = step;
380  fe.ts_ctx = TSMethod::CTX_TSSETRHSJACOBIAN;
381  fe.snes_ctx = SnesMethod::CTX_SNESSETJACOBIAN;
382  fe.ksp_ctx = KspMethod::CTX_OPERATORS;
383  fe.data_ctx = PetscData::CtxSetA | PetscData::CtxSetB |
385  fe.ts = ts;
386  };
387 
388  auto unset = [&](auto &fe) {
389  fe.ts_ctx = TSMethod::CTX_TSNONE;
390  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
391  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
392  fe.data_ctx = PetscData::CtxSetNone;
393  };
394 
395  // preprocess
396  for (auto &bit : ts_ctx->preProcess_RHSJacobian) {
397  bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
398  set(*bit);
400  *bit);
401  unset(*bit);
402  ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
403  }
404 
405  // fe loops
406  for (auto &lit : ts_ctx->loops_to_do_RHSJacobian) {
407  lit.second->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
408  set(*lit.second);
410  ts_ctx->problemName, lit.first, *(lit.second), nullptr, ts_ctx->bH);
411  unset(*lit.second);
412  ts_ctx->matAssembleSwitch = boost::move(lit.second->matAssembleSwitch);
413  }
414 
415  // post process
416  for (auto &bit : ts_ctx->postProcess_RHSJacobian) {
417  bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
418  set(*bit);
420  *bit);
421  unset(*bit);
422  ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
423  }
424 
425  if (ts_ctx->matAssembleSwitch) {
426  CHKERR MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
427  CHKERR MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
428  }
429 
430  PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxRHSJacobian, 0, 0, 0, 0);
432 }
433 
434 PetscErrorCode TsSetI2Jacobian(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt,
435  PetscReal v, PetscReal a, Mat A, Mat B,
436  void *ctx) {
438 
439  TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
440  PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxI2Function, 0, 0, 0, 0);
441  CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
442  CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
443  CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
444  CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
445  CHKERR VecGhostUpdateBegin(u_tt, INSERT_VALUES, SCATTER_FORWARD);
446  CHKERR VecGhostUpdateEnd(u_tt, INSERT_VALUES, SCATTER_FORWARD);
447  CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
448  ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
449  if (ts_ctx->zeroMatrix) {
450  CHKERR MatZeroEntries(B);
451  }
452  int step;
453 #if PETSC_VERSION_GE(3, 8, 0)
454  CHKERR TSGetStepNumber(ts, &step);
455 #else
456  CHKERR TSGetTimeStepNumber(ts, &step);
457 #endif
458 
459  ts_ctx->matAssembleSwitch =
460  boost::movelib::make_unique<bool>(ts_ctx->zeroMatrix);
461 
462  auto set = [&](auto &fe) {
463  fe.ts_u = u;
464  fe.ts_u_t = u_t;
465  fe.ts_u_tt = u_tt;
466  fe.ts_A = A;
467  fe.ts_B = B;
468  fe.ts_t = t;
469  fe.ts_a = a;
470  fe.ts_step = step;
471 
472  fe.ts_ctx = TSMethod::CTX_TSSETIJACOBIAN;
473  fe.snes_ctx = SnesMethod::CTX_SNESSETJACOBIAN;
474  fe.ksp_ctx = KspMethod::CTX_OPERATORS;
475  fe.data_ctx = PetscData::CtxSetA | PetscData::CtxSetB |
478  fe.ts = ts;
479  };
480 
481  auto unset = [&](auto &fe) {
482  fe.ts_ctx = TSMethod::CTX_TSNONE;
483  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
484  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
485  fe.data_ctx = PetscData::CtxSetNone;
486  };
487 
488  // preproces
489  for (auto &bit : ts_ctx->preProcess_IJacobian) {
490  bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
491  set(*bit);
493  *bit);
494  unset(*bit);
495  ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
496  }
497 
498  for (auto &lit : ts_ctx->loops_to_do_IJacobian) {
499  lit.second->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
500  set(*lit.second);
502  ts_ctx->problemName, lit.first, *(lit.second), nullptr, ts_ctx->bH);
503  unset(*lit.second);
504  ts_ctx->matAssembleSwitch = boost::move(lit.second->matAssembleSwitch);
505  }
506 
507  // post process
508  for (auto &bit : ts_ctx->postProcess_IJacobian) {
509  bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
510  set(*bit);
512  *bit);
513  unset(*bit);
514  ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
515  }
516 
517  if (ts_ctx->matAssembleSwitch) {
518  CHKERR MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
519  CHKERR MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
520  }
521  PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxI2Function, 0, 0, 0, 0);
523 }
524 
525 PetscErrorCode TsSetI2Function(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt,
526  Vec F, void *ctx) {
528  TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
529  PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
530  CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
531  CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
532  CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
533  CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
534  CHKERR VecGhostUpdateBegin(u_tt, INSERT_VALUES, SCATTER_FORWARD);
535  CHKERR VecGhostUpdateEnd(u_tt, INSERT_VALUES, SCATTER_FORWARD);
536  CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
537  ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
538 
539  auto zero_ghost_vec = [](Vec g) {
541  Vec l;
542  CHKERR VecGhostGetLocalForm(g, &l);
543  double *a;
544  CHKERR VecGetArray(l, &a);
545  int s;
546  CHKERR VecGetLocalSize(l, &s);
547  for (int i = 0; i != s; ++i)
548  a[i] = 0;
549  CHKERR VecRestoreArray(l, &a);
550  CHKERR VecGhostRestoreLocalForm(g, &l);
552  };
553  CHKERR zero_ghost_vec(F);
554 
555  ts_ctx->vecAssembleSwitch = boost::movelib::make_unique<bool>(true);
556 
557  int step;
558 #if PETSC_VERSION_GE(3, 8, 0)
559  CHKERR TSGetStepNumber(ts, &step);
560 #else
561  CHKERR TSGetTimeStepNumber(ts, &step);
562 #endif
563 
564  auto set = [&](auto &fe) {
565  fe.ts_u = u;
566  fe.ts_u_t = u_t;
567  fe.ts_u_tt = u_tt;
568  fe.ts_F = F;
569  fe.ts_t = t;
570  fe.ts_step = step;
571  fe.ts_ctx = TSMethod::CTX_TSSETIFUNCTION;
572  fe.snes_ctx = SnesMethod::CTX_SNESSETFUNCTION;
573  fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
574  fe.data_ctx = PetscData::CtxSetF | PetscData::CtxSetX |
577  fe.ts = ts;
578  };
579 
580  auto unset = [&](auto &fe) {
581  fe.ts_ctx = TSMethod::CTX_TSNONE;
582  fe.snes_ctx = SnesMethod::CTX_SNESNONE;
583  fe.ksp_ctx = KspMethod::CTX_KSPNONE;
584  fe.data_ctx = PetscData::CtxSetNone;
585  };
586 
587  // preprocess
588  for (auto &bit : ts_ctx->preProcess_IFunction) {
589  bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
590  set(*bit);
592  *bit);
593  unset(*bit);
594  ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
595  }
596 
597  // fe loops
598  for (auto &lit : ts_ctx->loops_to_do_IFunction) {
599  lit.second->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
600  set(*lit.second);
602  ts_ctx->problemName, lit.first, *(lit.second), nullptr, ts_ctx->bH);
603  unset(*lit.second);
604  ts_ctx->vecAssembleSwitch = boost::move(lit.second->vecAssembleSwitch);
605  }
606 
607  // post process
608  for (auto &bit : ts_ctx->postProcess_IFunction) {
609  bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
610  set(*bit);
612  *bit);
613  unset(*bit);
614  ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
615  }
616 
617  if (*ts_ctx->vecAssembleSwitch) {
618  CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
619  CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
620  CHKERR VecAssemblyBegin(F);
621  CHKERR VecAssemblyEnd(F);
622  }
623 
624  PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
626 }
627 
628 } // namespace MoFEM
PetscErrorCode TsSetI2Jacobian(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt, PetscReal v, PetscReal a, Mat A, Mat B, void *ctx)
Calculation Jaconian for second order PDE in time.
Definition: TsCtx.cpp:434
FEMethodsSequence loops_to_do_Monitor
Definition: TsCtx.hpp:50
static constexpr Switches CtxSetX_T
Definition: LoopMethods.hpp:66
PetscLogEvent MOFEM_EVENT_TsCtxIFunction
Definition: TsCtx.hpp:248
bool zeroMatrix
Definition: TsCtx.hpp:64
BasicMethodsSequence preProcess_RHSJacobian
Definition: TsCtx.hpp:59
static constexpr Switches CtxSetB
Definition: LoopMethods.hpp:64
Interface for Time Stepping (TS) solver.
Definition: TsCtx.hpp:27
MoFEMTypes bH
If set to MF_EXIST check if element exist.
Definition: TsCtx.hpp:33
BasicMethodsSequence preProcess_IJacobian
Definition: TsCtx.hpp:53
PetscErrorCode TsSetRHSFunction(TS ts, PetscReal t, Vec u, Vec F, void *ctx)
TS solver function.
Definition: TsCtx.cpp:256
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:482
PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Set monitor for TS solver.
Definition: TsCtx.cpp:202
std::string problemName
Definition: TsCtx.hpp:32
boost::movelib::unique_ptr< bool > vecAssembleSwitch
Definition: TsCtx.hpp:254
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:21
PetscLogEvent MOFEM_EVENT_TsCtxRHSFunction
Definition: TsCtx.hpp:246
MoFEMErrorCode getInterface(const MOFEMuuid &uuid, IFACE *&iface) const
Get interface by uuid and return reference to pointer of interface.
static constexpr Switches CtxSetTime
Definition: LoopMethods.hpp:68
PetscErrorCode TsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t, Vec F, void *ctx)
Set IFunction for TS solver.
Definition: TsCtx.cpp:17
FEMethodsSequence loops_to_do_RHSJacobian
Definition: TsCtx.hpp:51
BasicMethodsSequence postProcess_IFunction
Definition: TsCtx.hpp:56
BasicMethodsSequence preProcess_IFunction
Definition: TsCtx.hpp:55
BasicMethodsSequence postProcess_IJacobian
Definition: TsCtx.hpp:54
static constexpr Switches CtxSetF
Definition: LoopMethods.hpp:62
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:36
FEMethodsSequence loops_to_do_IJacobian
Definition: TsCtx.hpp:48
static constexpr Switches CtxSetX_TT
Definition: LoopMethods.hpp:67
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:601
PetscLogEvent MOFEM_EVENT_TsCtxI2Function
Definition: TsCtx.hpp:251
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:525
PetscLogEvent MOFEM_EVENT_TsCtxMonitor
Definition: TsCtx.hpp:250
boost::movelib::unique_ptr< bool > matAssembleSwitch
Definition: TsCtx.hpp:255
static constexpr Switches CtxSetA
Definition: LoopMethods.hpp:63
PetscErrorCode TsSetRHSJacobian(TS ts, PetscReal t, Vec u, Mat A, Mat B, void *ctx)
TS solver function.
Definition: TsCtx.cpp:350
FEMethodsSequence loops_to_do_IFunction
Definition: TsCtx.hpp:49
MoFEM::Interface & mField
Definition: TsCtx.hpp:29
BasicMethodsSequence preProcess_Monitor
Definition: TsCtx.hpp:57
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:412
FEMethodsSequence loops_to_do_RHSFunction
Definition: TsCtx.hpp:52
PetscLogEvent MOFEM_EVENT_TsCtxRHSJacobian
Definition: TsCtx.hpp:247
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
FTensor::Index< 'l', 2 > l
Definition: ContactOps.hpp:29
BasicMethodsSequence postProcess_Monitor
Definition: TsCtx.hpp:58
BasicMethodsSequence postProcess_RHSJacobian
Definition: TsCtx.hpp:61
PetscErrorCode TsSetIJacobian(TS ts, PetscReal t, Vec u, Vec u_t, PetscReal a, Mat A, Mat B, void *ctx)
Set function evaluating jacobina in TS solver.
Definition: TsCtx.cpp:116
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