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