v0.14.0
PCMGSetUpViaApproxOrders.cpp
Go to the documentation of this file.
1 /** \file PCMGSetUpViaApproxOrders.cpp
2  * \brief implementation of multi-grid solver for p- adaptivity
3  */
4 
5 
6 
7 #include <MoFEM.hpp>
8 
9 #include <UnknownInterface.hpp>
10 
11 using namespace MoFEM;
12 
14 
15 #undef PETSC_VERSION_RELEASE
16 #define PETSC_VERSION_RELEASE 1
17 
18 #if PETSC_VERSION_GE(3, 6, 0)
19 #include <petsc/private/petscimpl.h>
20 #else
21 #include <petsc-private/petscimpl.h>
22 #endif
23 
24 #if PETSC_VERSION_GE(3, 6, 0)
25 #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/
26 // #include <petsc/private/vecimpl.h> /*I "petscdm.h" I*/
27 #else
28 #include <petsc-private/dmimpl.h> /*I "petscdm.h" I*/
29 #include <petsc-private/vecimpl.h> /*I "petscdm.h" I*/
30 #endif
31 
32 PCMGSubMatrixCtx::PCMGSubMatrixCtx(Mat a, IS is) : A(a), iS(is) {
33  // Increase reference of petsc object (works like shared_ptr but unique for
34  // PETSc)
35  ierr = PetscObjectReference((PetscObject)A);
36  CHKERRABORT(PETSC_COMM_WORLD, ierr);
37  ierr = PetscObjectReference((PetscObject)iS);
38  CHKERRABORT(PETSC_COMM_WORLD, ierr);
39 }
40 
42  ierr = MatDestroy(&A);
43  CHKERRABORT(PETSC_COMM_WORLD, ierr);
44  ierr = ISDestroy(&iS);
45  CHKERRABORT(PETSC_COMM_WORLD, ierr);
46 }
47 
50  : PCMGSubMatrixCtx(a, is), isInitisalised(false) {
51  PetscLogEventRegister("PCMGSubMatrixCtx_mult", 0, &MOFEM_EVENT_mult);
52  PetscLogEventRegister("PCMGSubMatrixCtx_sor", 0, &MOFEM_EVENT_sor);
53  }
55  if (isInitisalised) {
56  ierr = VecScatterDestroy(&sCat);
57  CHKERRABORT(PETSC_COMM_WORLD, ierr);
58  ierr = VecDestroy(&X);
59  CHKERRABORT(PETSC_COMM_WORLD, ierr);
60  ierr = VecDestroy(&F);
61  CHKERRABORT(PETSC_COMM_WORLD, ierr);
62  }
63  }
64  template <InsertMode MODE>
65  friend MoFEMErrorCode sub_mat_mult_generic(Mat a, Vec x, Vec f);
66  friend MoFEMErrorCode sub_mat_sor(Mat mat, Vec b, PetscReal omega,
67  MatSORType flag, PetscReal shift,
68  PetscInt its, PetscInt lits, Vec x);
69 
70 public:
73  if (!isInitisalised) {
74  CHKERR MatCreateVecs(A, &X, &F);
75  CHKERR VecScatterCreate(X, iS, x, PETSC_NULL, &sCat);
76  CHKERR VecZeroEntries(X);
77  CHKERR VecZeroEntries(F);
78  isInitisalised = true;
79  }
81  }
82  PetscLogEvent MOFEM_EVENT_mult;
83  PetscLogEvent MOFEM_EVENT_sor;
85 };
86 
87 template <InsertMode MODE>
89  void *void_ctx;
91  CHKERR MatShellGetContext(a, &void_ctx);
93  if (!ctx->isInitisalised) {
94  CHKERR ctx->initData(x);
95  }
96  PetscLogEventBegin(ctx->MOFEM_EVENT_mult, 0, 0, 0, 0);
97  CHKERR VecScatterBegin(ctx->sCat, x, ctx->X, INSERT_VALUES, SCATTER_REVERSE);
98  CHKERR VecScatterEnd(ctx->sCat, x, ctx->X, INSERT_VALUES, SCATTER_REVERSE);
99  CHKERR MatMult(ctx->A, ctx->X, ctx->F);
100  CHKERR VecScatterBegin(ctx->sCat, ctx->F, f, MODE, SCATTER_FORWARD);
101  CHKERR VecScatterEnd(ctx->sCat, ctx->F, f, MODE, SCATTER_FORWARD);
102  PetscLogEventEnd(ctx->MOFEM_EVENT_mult, 0, 0, 0, 0);
104 }
105 
107  return sub_mat_mult_generic<INSERT_VALUES>(a, x, f);
108 }
109 
111  return sub_mat_mult_generic<ADD_VALUES>(a, x, f);
112 }
113 
114 MoFEMErrorCode sub_mat_sor(Mat mat, Vec b, PetscReal omega, MatSORType flag,
115  PetscReal shift, PetscInt its, PetscInt lits,
116  Vec x) {
117  void *void_ctx;
119  CHKERR MatShellGetContext(mat, &void_ctx);
121  if (!ctx->isInitisalised) {
122  CHKERR ctx->initData(x);
123  }
124  PetscLogEventBegin(ctx->MOFEM_EVENT_sor, 0, 0, 0, 0);
125  CHKERR VecScatterBegin(ctx->sCat, b, ctx->X, INSERT_VALUES, SCATTER_REVERSE);
126  CHKERR VecScatterEnd(ctx->sCat, b, ctx->X, INSERT_VALUES, SCATTER_REVERSE);
127  CHKERR MatSOR(ctx->A, ctx->X, omega, flag, shift, its, lits, ctx->F);
128  CHKERR VecScatterBegin(ctx->sCat, ctx->F, x, INSERT_VALUES, SCATTER_FORWARD);
129  CHKERR VecScatterEnd(ctx->sCat, ctx->F, x, INSERT_VALUES, SCATTER_FORWARD);
130  PetscLogEventEnd(ctx->MOFEM_EVENT_sor, 0, 0, 0, 0);
132 }
133 
135  : MoFEM::DMCtxImpl(), aO(PETSC_NULL) {}
136 
139  CHKERRABORT(PETSC_COMM_WORLD, ierr);
140 }
141 
144  for (unsigned int ii = 0; ii < coarseningIS.size(); ii++) {
145  if (coarseningIS[ii])
146  CHKERR ISDestroy(&coarseningIS[ii]);
147  }
148  for (unsigned int ii = 0; ii < kspOperators.size(); ii++) {
149  if (kspOperators[ii])
150  CHKERR MatDestroy(&kspOperators[ii]);
151  }
152  if (aO) {
153  CHKERR AODestroy(&aO);
154  }
155  coarseningIS.clear();
156  kspOperators.clear();
157  shellMatrixCtxPtr.clear();
159 }
160 
162 DMMGViaApproxOrdersCtx::query_interface(boost::typeindex::type_index type_index,
163  MoFEM::UnknownInterface **iface) const {
164  *iface = static_cast<DMMGViaApproxOrdersCtx *>(
165  const_cast<DMMGViaApproxOrdersCtx *>(this));
166  return 0;
167 }
168 
169 #define GET_DM_FIELD(DM) \
170  auto dm_field = \
171  static_cast<DMCtx *>(DM->data)->getInterface<DMMGViaApproxOrdersCtx>(); \
172  NOT_USED(dm_field)
173 
175  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
177  GET_DM_FIELD(dm);
178  *ctx = dm_field;
180 }
181 
183  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
185  GET_DM_FIELD(dm);
186  if (dm_field->aO) {
187  // std::cerr << dm_field->aO << std::endl;
188  CHKERR AODestroy(&dm_field->aO);
189  // std::cerr << "destroy ao when adding\n";
190  }
191  dm_field->aO = ao;
192  CHKERR PetscObjectReference((PetscObject)ao);
193  // std::cerr << "add ao\n";
195 }
196 
198  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
200  GET_DM_FIELD(dm);
201  *size = dm_field->coarseningIS.size();
203 }
204 
206  Mat *subA,
207  bool create_sub_matrix,
208  bool shell_sub_a) {
209  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
211  GET_DM_FIELD(dm);
212  dm_field->coarseningIS.push_back(is);
213  dm_field->shellMatrixCtxPtr.push_back(new PCMGSubMatrixCtx_private(A, is));
214  if (is) {
215  CHKERR PetscObjectReference((PetscObject)is);
216  }
217  if (is) {
218  IS is2 = is;
219  if (dm_field->aO) {
220  CHKERR ISDuplicate(is, &is2);
221  CHKERR ISCopy(is, is2);
222  CHKERR AOApplicationToPetscIS(dm_field->aO, is2);
223  }
224  if (create_sub_matrix) {
225  if (shell_sub_a) {
226  int n, N;
227  CHKERR ISGetSize(is, &N);
228  CHKERR ISGetLocalSize(is, &n);
229  MPI_Comm comm;
230  CHKERR PetscObjectGetComm((PetscObject)A, &comm);
231  CHKERR MatCreateShell(comm, n, n, N, N,
232  &(dm_field->shellMatrixCtxPtr.back()), subA);
233  CHKERR MatShellSetOperation(*subA, MATOP_MULT,
234  (void (*)(void))sub_mat_mult);
235  CHKERR MatShellSetOperation(*subA, MATOP_MULT_ADD,
236  (void (*)(void))sub_mat_mult_add);
237  CHKERR MatShellSetOperation(*subA, MATOP_SOR,
238  (void (*)(void))sub_mat_sor);
239  } else {
240 #if PETSC_VERSION_GE(3, 8, 0)
241  CHKERR MatCreateSubMatrix(A, is2, is2, MAT_INITIAL_MATRIX, subA);
242 #else
243  CHKERR MatGetSubMatrix(A, is2, is2, MAT_INITIAL_MATRIX, subA);
244 #endif
245  }
246  }
247  if (dm_field->aO) {
248  CHKERR ISDestroy(&is2);
249  }
250  dm_field->kspOperators.push_back(*subA);
251  CHKERR PetscObjectReference((PetscObject)(*subA));
252  } else {
253  SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
254  }
255  PetscInfo(dm, "Push back IS to DMMGViaApproxOrders\n");
257 }
258 
260  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
262  GET_DM_FIELD(dm);
263  if (dm_field->coarseningIS.back()) {
264  CHKERR ISDestroy(&dm_field->coarseningIS.back());
265  dm_field->coarseningIS.pop_back();
266  }
267  if (dm_field->kspOperators.back()) {
268  CHKERR MatDestroy(&dm_field->kspOperators.back());
269  }
270  dm_field->kspOperators.pop_back();
271  PetscInfo(dm, "Pop back IS to DMMGViaApproxOrders\n");
273 }
274 
276  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
278  GET_DM_FIELD(dm);
279  CHKERR dm_field->destroyCoarseningIS();
280  PetscInfo(dm, "Clear DMs data structures\n");
282 }
283 
285  int nb_elems, Mat A,
286  int verb) {
287  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
289  GET_DM_FIELD(dm);
290  int nb_no_changed = 0;
291  int nb_replaced = 0;
292  int nb_deleted = 0;
293  int nb_added = 0;
294  std::vector<IS>::iterator it;
295  it = dm_field->coarseningIS.begin();
296  int ii = 0;
297  for (; it != dm_field->coarseningIS.end(); it++, ii++) {
298  if (ii < nb_elems) {
299  PetscBool flg;
300  CHKERR ISEqual(*it, is_vec[ii], &flg);
301  if (!flg) {
302  CHKERR ISDestroy(&*it);
303  CHKERR MatDestroy(&dm_field->kspOperators[ii]);
304  *it = is_vec[ii];
305  CHKERR PetscObjectReference((PetscObject)is_vec[ii]);
306  if (ii < nb_elems - 1) {
307  IS is = is_vec[ii];
308  if (dm_field->aO) {
309  CHKERR ISDuplicate(is_vec[ii], &is);
310  CHKERR ISCopy(is_vec[ii], is);
311  CHKERR AOApplicationToPetscIS(dm_field->aO, is);
312  }
313  Mat subA;
314 #if PETSC_VERSION_GE(3, 8, 0)
315  CHKERR MatCreateSubMatrix(A, is, is, MAT_INITIAL_MATRIX, &subA);
316 #else
317  CHKERR MatGetSubMatrix(A, is, is, MAT_INITIAL_MATRIX, &subA);
318 #endif
319  CHKERR PetscObjectReference((PetscObject)subA);
320  dm_field->kspOperators[ii] = subA;
321  CHKERR MatDestroy(&subA);
322  if (dm_field->aO) {
323  CHKERR ISDestroy(&is);
324  }
325  } else {
326  CHKERR PetscObjectReference((PetscObject)A);
327  dm_field->kspOperators[ii] = A;
328  }
329  nb_replaced++;
330  }
331  } else {
332  nb_no_changed++;
333  continue;
334  }
335  }
336  if (static_cast<int>(dm_field->coarseningIS.size()) < nb_elems) {
337  for (; ii < nb_elems - 1; ii++) {
338  Mat subA;
339  CHKERR DMMGViaApproxOrdersPushBackCoarseningIS(dm, is_vec[ii], A, &subA,
340  true, false);
341  CHKERR MatDestroy(&subA);
342  nb_added++;
343  }
344  CHKERR DMMGViaApproxOrdersPushBackCoarseningIS(dm, is_vec[ii], A, &A, false,
345  false);
346  nb_added++;
347  } else {
348  for (; ii < static_cast<int>(dm_field->coarseningIS.size()); ii++) {
350  nb_deleted++;
351  }
352  }
353  MPI_Comm comm;
354  CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
355  if (verb > 0) {
356  PetscPrintf(comm,
357  "DMMGViaApproxOrders nb_no_changed = %d, nb_replaced = %d, "
358  "nb_added = %d, nb_deleted = %d, size = %d\n",
359  nb_no_changed, nb_replaced, nb_added, nb_deleted,
360  dm_field->coarseningIS.size());
361  }
362  PetscInfo(dm, "Replace IS to DMMGViaApproxOrders\n");
364 }
365 
367  const DMMGViaApproxOrdersCtx **ctx) {
368  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
370  GET_DM_FIELD(dm);
371  *ctx = dm_field;
373 }
374 
377  CHKERR DMRegister(sname, DMCreate_MGViaApproxOrders);
379 }
380 
381 // MoFEMErrorCode DMDestroy_MGViaApproxOrders(DM dm) {
382 // PetscValidHeaderSpecific(dm,DM_CLASSID,1);
383 // MoFEMFunctionBeginHot;
384 // if(!((DMMGViaApproxOrdersCtx*)dm->data)->referenceNumber) {
385 // DMMGViaApproxOrdersCtx *dm_field = (DMMGViaApproxOrdersCtx*)dm->data;
386 // if(dm_field->destroyProblem) {
387 // if(dm_field->mField_ptr->check_problem(dm_field->problemName)) {
388 // dm_field->mField_ptr->delete_problem(dm_field->problemName);
389 // } // else problem has to be deleted by the user
390 // }
391 // cerr << "Destroy " << dm_field->problemName << endl;
392 // delete (DMMGViaApproxOrdersCtx*)dm->data;
393 // } else {
394 // DMMGViaApproxOrdersCtx *dm_field = (DMMGViaApproxOrdersCtx*)dm->data;
395 //
396 // cerr << "Derefrence " << dm_field->problemName << " " <<
397 // ((DMCtx*)dm->data)->referenceNumber << endl;
398 // (((DMMGViaApproxOrdersCtx*)dm->data)->referenceNumber)--;
399 // }
400 // MoFEMFunctionReturnHot(0);
401 // }
402 
403 static MoFEMErrorCode ksp_set_operators(KSP ksp, Mat A, Mat B, void *ctx) {
406 }
407 
409  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
411  if (!dm->data) {
412  dm->data = new DMMGViaApproxOrdersCtx();
413  } else {
414  ((DMCtxImpl *)(dm->data))->incrementReference();
415  }
416  // cerr << "Create " << ((DMCtx*)(dm->data))->referenceNumber << endl;
418  dm->ops->creatematrix = DMCreateMatrix_MGViaApproxOrders;
419  dm->ops->createglobalvector = DMCreateGlobalVector_MGViaApproxOrders;
420  dm->ops->coarsen = DMCoarsen_MGViaApproxOrders;
421  // dm->ops->destroy = DMDestroy_MGViaApproxOrders;
422  dm->ops->createinterpolation = DMCreateInterpolation_MGViaApproxOrders;
423  CHKERR DMKSPSetComputeOperators(dm, ksp_set_operators, NULL);
424  PetscInfo1(dm, "Create DMMGViaApproxOrders reference = %d\n",
425  ((DMCtxImpl *)(dm->data))->useCount());
427 }
428 
430 
431  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
433  GET_DM_FIELD(dm);
434 
435  int leveldown = dm->leveldown;
436 
437  if (dm_field->kspOperators.empty()) {
439  } else {
440  MPI_Comm comm;
441  CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
442  if (dm_field->kspOperators.empty()) {
443  SETERRQ(comm, MOFEM_DATA_INCONSISTENCY,
444  "data inconsistency, operator can not be set");
445  }
446  if (static_cast<int>(dm_field->kspOperators.size()) < leveldown) {
447  SETERRQ(comm, MOFEM_DATA_INCONSISTENCY,
448  "data inconsistency, no IS for that level");
449  }
450  *M = dm_field->kspOperators[dm_field->kspOperators.size() - 1 - leveldown];
451  CHKERR PetscObjectReference((PetscObject)*M);
452  CHKERR MatSetDM(*M, dm);
453  }
454 
455  PetscInfo1(dm, "Create Matrix DMMGViaApproxOrders leveldown = %d\n",
456  leveldown);
457 
459 }
460 
461 MoFEMErrorCode DMCoarsen_MGViaApproxOrders(DM dm, MPI_Comm comm, DM *dmc) {
462  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
464  GET_DM_FIELD(dm);
465  CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
466  CHKERR DMCreate(comm, dmc);
467  (*dmc)->data = dm->data;
468  DMType type;
469  CHKERR DMGetType(dm, &type);
470  CHKERR DMSetType(*dmc, type);
471  CHKERR PetscObjectReference((PetscObject)(*dmc));
472  PetscInfo1(dm, "Coarsen DMMGViaApproxOrders leveldown = %d\n", dm->leveldown);
474 }
475 
477  Vec *vec) {
478  PetscValidHeaderSpecific(dm1, DM_CLASSID, 1);
479  PetscValidHeaderSpecific(dm2, DM_CLASSID, 1);
481 
482  MPI_Comm comm;
483  CHKERR PetscObjectGetComm((PetscObject)dm1, &comm);
484 
485  int m, n, M, N;
486 
487  DM dm_down = dm1;
488  DM dm_up = dm2;
489 
490  int dm_down_leveldown = dm_down->leveldown;
491  int dm_up_leveldown = dm_up->leveldown;
492 
493  PetscInfo2(dm1,
494  "Create interpolation DMMGViaApproxOrders dm1_leveldown = %d "
495  "dm2_leveldown = %d\n",
496  dm_down_leveldown, dm_up_leveldown);
497 
498  IS is_down, is_up;
499  {
500  // Coarser mesh
501  GET_DM_FIELD(dm_down);
502  if (static_cast<int>(dm_field->coarseningIS.size()) < dm_down_leveldown) {
503  SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
504  }
505  is_down = dm_field->coarseningIS[dm_field->coarseningIS.size() - 1 -
506  dm_down_leveldown];
507  CHKERR ISGetSize(is_down, &M);
508  CHKERR ISGetLocalSize(is_down, &m);
509  }
510  {
511  // Finer mesh
512  GET_DM_FIELD(dm_up);
513  if (static_cast<int>(dm_field->coarseningIS.size()) < dm_up_leveldown) {
514  SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
515  }
516  is_up =
517  dm_field
518  ->coarseningIS[dm_field->coarseningIS.size() - 1 - dm_up_leveldown];
519  CHKERR ISGetSize(is_up, &N);
520  CHKERR ISGetLocalSize(is_up, &n);
521  }
522 
523  // is_dow rows
524  // is_up columns
525 
526  CHKERR MatCreate(comm, mat);
527  CHKERR MatSetSizes(*mat, m, n, M, N);
528  CHKERR MatSetType(*mat, MATMPIAIJ);
529  CHKERR MatMPIAIJSetPreallocation(*mat, 1, PETSC_NULL, 0, PETSC_NULL);
530 
531  // get matrix layout
532  PetscLayout rmap, cmap;
533  CHKERR MatGetLayouts(*mat, &rmap, &cmap);
534  int rstart, rend, cstart, cend;
535  CHKERR PetscLayoutGetRange(rmap, &rstart, &rend);
536  CHKERR PetscLayoutGetRange(cmap, &cstart, &cend);
537 
538  // if(verb>0) {
539  // PetscSynchronizedPrintf(comm,"level %d row start %d row end
540  // %d\n",kk,rstart,rend); PetscSynchronizedPrintf(comm,"level %d col start
541  // %d col end %d\n",kk,cstart,cend);
542  // }
543 
544  const int *row_indices_ptr, *col_indices_ptr;
545  CHKERR ISGetIndices(is_down, &row_indices_ptr);
546  CHKERR ISGetIndices(is_up, &col_indices_ptr);
547 
548  map<int, int> idx_map;
549  for (int ii = 0; ii < m; ii++) {
550  idx_map[row_indices_ptr[ii]] = rstart + ii;
551  }
552 
553  CHKERR MatZeroEntries(*mat);
554  // FIXME: Use MatCreateMPIAIJWithArrays and set array directly
555  for (int jj = 0; jj < n; jj++) {
556  map<int, int>::iterator mit = idx_map.find(col_indices_ptr[jj]);
557  if (mit != idx_map.end()) {
558  CHKERR MatSetValue(*mat, mit->second, cstart + jj, 1, INSERT_VALUES);
559  }
560  }
561 
562  CHKERR ISRestoreIndices(is_down, &row_indices_ptr);
563  CHKERR ISRestoreIndices(is_up, &col_indices_ptr);
564 
565  CHKERR MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY);
566  CHKERR MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY);
567 
568  if (vec != NULL) {
569  *vec = PETSC_NULL;
570  }
571 
573 }
574 
576 
577  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
579  int leveldown = dm->leveldown;
580  GET_DM_FIELD(dm);
581  if (dm_field->kspOperators.empty()) {
583  } else {
584 #if PETSC_VERSION_GE(3, 5, 3)
585  CHKERR MatCreateVecs(
586  dm_field->kspOperators[dm_field->kspOperators.size() - 1 - leveldown],
587  g, NULL);
588 #else
589  CHKERR MatGetVecs(
590  dm_field->kspOperators[dm_field->kspOperators.size() - 1 - leveldown],
591  g, NULL);
592 #endif
593  CHKERR VecSetDM(*g, dm);
594  }
595  PetscInfo1(dm, "Create global vector DMMGViaApproxOrders leveldown = %d\n",
596  dm->leveldown);
598 }
599 
602  ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "",
603  "MOFEM Multi-Grid (Orders) pre-conditioner", "none");
604  CHKERRG(ierr);
605 
606  CHKERR PetscOptionsInt("-mofem_mg_levels", "nb levels of multi-grid solver",
607  "", 2, &nbLevels, PETSC_NULL);
608  CHKERR PetscOptionsInt("-mofem_mg_coarse_order",
609  "approximation order of coarse level", "", 2,
610  &coarseOrder, PETSC_NULL);
611  CHKERR PetscOptionsInt("-mofem_mg_order_at_last_level", "order at last level",
612  "", 100, &orderAtLastLevel, PETSC_NULL);
613  CHKERR PetscOptionsInt("-mofem_mg_verbose", "nb levels of multi-grid solver",
614  "", 0, &verboseLevel, PETSC_NULL);
615  PetscBool shell_sub_a = shellSubA ? PETSC_TRUE : PETSC_FALSE;
616  CHKERR PetscOptionsBool("-mofem_mg_shell_a", "use shell matrix as sub matrix",
617  "", shell_sub_a, &shell_sub_a, NULL);
618  shellSubA = (shellSubA == PETSC_TRUE);
619 
620  ierr = PetscOptionsEnd();
621  CHKERRG(ierr);
623 }
624 
626  MoFEM::Interface *m_field_ptr;
627  MoFEM::ISManager *is_manager_ptr;
629  // if is last level, take all remaining orders dofs, if any left
630  CHKERR DMoFEMGetInterfacePtr(dM, &m_field_ptr);
631  CHKERR m_field_ptr->getInterface(is_manager_ptr);
632  const Problem *problem_ptr;
633  CHKERR DMMoFEMGetProblemPtr(dM, &problem_ptr);
634  int order_at_next_level = kk + coarseOrder;
635  if (kk == nbLevels - 1) {
636  int first = problem_ptr->getNumeredRowDofsPtr()
637  ->get<PetscLocalIdx_mi_tag>()
638  .find(0)
639  ->get()
640  ->getPetscGlobalDofIdx();
641  CHKERR ISCreateStride(PETSC_COMM_WORLD, problem_ptr->getNbLocalDofsRow(),
642  first, 1, is);
644  // order_at_next_level = orderAtLastLevel;
645  }
646  string problem_name = problem_ptr->getName();
647  CHKERR is_manager_ptr->isCreateProblemOrder(problem_name, ROW, 0,
648  order_at_next_level, is);
650 }
651 
654  CHKERR ISDestroy(is);
656 }
657 
660  int verb) {
662  verb = verb > verboseLevel ? verb : verboseLevel;
663 
664  MPI_Comm comm;
665  CHKERR PetscObjectGetComm((PetscObject)dM, &comm);
666 
667  if (verb > QUIET) {
668  PetscPrintf(comm, "set MG levels %u\n", nbLevels);
669  }
670 
671  std::vector<IS> is_vec(nbLevels + 1);
672  std::vector<int> is_glob_size(nbLevels + 1), is_loc_size(nbLevels + 1);
673 
674  for (int kk = 0; kk < nbLevels; kk++) {
675 
676  // get indices up to up to give approximation order
677  CHKERR createIsAtLevel(kk, &is_vec[kk]);
678  CHKERR ISGetSize(is_vec[kk], &is_glob_size[kk]);
679  CHKERR ISGetLocalSize(is_vec[kk], &is_loc_size[kk]);
680 
681  if (verb > QUIET) {
682  PetscSynchronizedPrintf(comm,
683  "Nb. dofs at level [ %d ] global %u local %d\n",
684  kk, is_glob_size[kk], is_loc_size[kk]);
685  }
686 
687  // if no dofs on level kk finish here
688  if (is_glob_size[kk] == 0) {
689  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "no dofs at level");
690  }
691  }
692 
693  for (int kk = 0; kk != nbLevels; kk++) {
694  Mat subA;
695  if (kk == nbLevels - 1 && use_mat_a) {
696  subA = A;
698  false, false);
699  } else {
700  if (kk > 0) {
701  // Not coarse level
703  true, shellSubA);
704  } else {
705  // Coarse lave is compressed matrix allowing for factorization when
706  // needed
708  true, false);
709  }
710  if (subA) {
711  CHKERR MatDestroy(&subA);
712  }
713  }
714  }
715 
716  for (unsigned int kk = 0; kk < is_vec.size(); kk++) {
717  CHKERR destroyIsAtLevel(kk, &is_vec[kk]);
718  }
719 
720  if (verb > QUIET) {
721  PetscSynchronizedFlush(comm, PETSC_STDOUT);
722  }
723 
725 }
726 
728  int verb) {
729 
730  PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
732 
733  MPI_Comm comm;
734  CHKERR PetscObjectGetComm((PetscObject)pc, &comm);
735  if (verb > 0) {
736  PetscPrintf(comm, "Start PCMGSetUpViaApproxOrders\n");
737  }
738 
739  CHKERR ctx->getOptions();
740  CHKERR ctx->buildProlongationOperator(true, verb);
741 
742 #if PETSC_VERSION_GE(3, 8, 0)
743  CHKERR PCMGSetGalerkin(pc, PC_MG_GALERKIN_NONE);
744 #else
745  CHKERR PCMGSetGalerkin(pc, PETSC_FALSE);
746 #endif
747 
748  CHKERR PCMGSetLevels(pc, ctx->nbLevels, NULL);
749 
750  if (verb > 0) {
751  PetscPrintf(comm, "End PCMGSetUpViaApproxOrders\n");
752  }
753 
755 }
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
PCMGSubMatrixCtx::~PCMGSubMatrixCtx
virtual ~PCMGSubMatrixCtx()
Definition: PCMGSetUpViaApproxOrders.cpp:41
g
constexpr double g
Definition: shallow_wave.cpp:63
omega
constexpr double omega
Save field DOFS on vertices/tags.
Definition: dynamic_first_order_con_law.cpp:93
PCMGSetUpViaApproxOrdersCtx::getOptions
virtual MoFEMErrorCode getOptions()
get options from line command
Definition: PCMGSetUpViaApproxOrders.cpp:600
PCMGSubMatrixCtx::X
Vec X
Definition: PCMGSetUpViaApproxOrders.hpp:11
DMMGViaApproxOrdersCtx::~DMMGViaApproxOrdersCtx
virtual ~DMMGViaApproxOrdersCtx()
Definition: PCMGSetUpViaApproxOrders.cpp:137
DMMGViaApproxOrdersClearCoarseningIS
MoFEMErrorCode DMMGViaApproxOrdersClearCoarseningIS(DM dm)
Clear approximation orders.
Definition: PCMGSetUpViaApproxOrders.cpp:275
DMMGViaApproxOrdersCtx::shellMatrixCtxPtr
boost::ptr_vector< PCMGSubMatrixCtx > shellMatrixCtxPtr
Shell sub-matrix context.
Definition: PCMGSetUpViaApproxOrders.hpp:36
PCMGSetUpViaApproxOrdersCtx::shellSubA
bool shellSubA
Definition: PCMGSetUpViaApproxOrders.hpp:208
PCMGSubMatrixCtx_private::~PCMGSubMatrixCtx_private
~PCMGSubMatrixCtx_private()
Definition: PCMGSetUpViaApproxOrders.cpp:54
DMCoarsen_MGViaApproxOrders
MoFEMErrorCode DMCoarsen_MGViaApproxOrders(DM dm, MPI_Comm comm, DM *dmc)
Coarsen DM.
Definition: PCMGSetUpViaApproxOrders.cpp:461
DMMGViaApproxOrdersCtx::kspOperators
std::vector< Mat > kspOperators
Get KSP operators.
Definition: PCMGSetUpViaApproxOrders.hpp:34
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
PCMGSetUpViaApproxOrdersCtx
Set data structures of MG pre-conditioner via approximation orders.
Definition: PCMGSetUpViaApproxOrders.hpp:190
DMMGViaApproxOrdersCtx
Structure for DM for multi-grid via approximation orders.
Definition: PCMGSetUpViaApproxOrders.hpp:22
MoFEM::PetscLocalIdx_mi_tag
Definition: TagMultiIndices.hpp:45
MoFEM.hpp
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
DMMGViaApproxOrdersPushBackCoarseningIS
MoFEMErrorCode DMMGViaApproxOrdersPushBackCoarseningIS(DM dm, IS is, Mat A, Mat *subA, bool create_sub_matrix, bool shell_sub_a)
Push back coarsening level to MG via approximation orders.
Definition: PCMGSetUpViaApproxOrders.cpp:205
PCMGSetUpViaApproxOrders.hpp
header of multi-grid solver for p- adaptivity
PCMGSetUpViaApproxOrdersCtx::nbLevels
int nbLevels
number of multi-grid levels
Definition: PCMGSetUpViaApproxOrders.hpp:204
GET_DM_FIELD
#define GET_DM_FIELD(DM)
Definition: PCMGSetUpViaApproxOrders.cpp:169
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
ROW
@ ROW
Definition: definitions.h:123
sub_mat_mult_generic
MoFEMErrorCode sub_mat_mult_generic(Mat a, Vec x, Vec f)
Definition: PCMGSetUpViaApproxOrders.cpp:88
DMCreate_MGViaApproxOrders
MoFEMErrorCode DMCreate_MGViaApproxOrders(DM dm)
Create DM data structure for Multi-Grid via approximation orders.
Definition: PCMGSetUpViaApproxOrders.cpp:408
MoFEM::DMCreateGlobalVector_MoFEM
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVector.
Definition: DMMoFEM.cpp:1167
PlasticOps::M
FTensor::Index< 'M', 3 > M
Definition: PlasticOps.hpp:117
MoFEM::DMCtxImpl
Definition: DMCtxImpl.hpp:10
MoFEM::DMCreateMatrix_MoFEM
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
Definition: DMMoFEM.cpp:1197
DMMGViaApproxOrdersCtx::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, MoFEM::UnknownInterface **iface) const
Definition: PCMGSetUpViaApproxOrders.cpp:162
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
PCMGSetUpViaApproxOrders
MoFEMErrorCode PCMGSetUpViaApproxOrders(PC pc, PCMGSetUpViaApproxOrdersCtx *ctx, int verb)
Function build MG structure.
Definition: PCMGSetUpViaApproxOrders.cpp:727
DMMGViaApproxOrdersCtx::DMMGViaApproxOrdersCtx
DMMGViaApproxOrdersCtx()
Definition: PCMGSetUpViaApproxOrders.cpp:134
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::ISManager
Section manager is used to create indexes and sections.
Definition: ISManager.hpp:23
a
constexpr double a
Definition: approx_sphere.cpp:30
sub_mat_mult
MoFEMErrorCode sub_mat_mult(Mat a, Vec x, Vec f)
Definition: PCMGSetUpViaApproxOrders.cpp:106
DMMGViaApproxOrdersSetAO
MoFEMErrorCode DMMGViaApproxOrdersSetAO(DM dm, AO ao)
Set DM ordering.
Definition: PCMGSetUpViaApproxOrders.cpp:182
convert.type
type
Definition: convert.py:64
PCMGSubMatrixCtx::F
Vec F
Definition: PCMGSetUpViaApproxOrders.hpp:11
PCMGSetUpViaApproxOrdersCtx::dM
DM dM
Distributed mesh manager.
Definition: PCMGSetUpViaApproxOrders.hpp:195
PCMGSubMatrixCtx_private
Definition: PCMGSetUpViaApproxOrders.cpp:48
DMMGViaApproxOrdersCtx::destroyCoarseningIS
MoFEMErrorCode destroyCoarseningIS()
Definition: PCMGSetUpViaApproxOrders.cpp:142
MoFEM::Problem::getNumeredRowDofsPtr
auto & getNumeredRowDofsPtr() const
get access to numeredRowDofsPtr storing DOFs on rows
Definition: ProblemsMultiIndices.hpp:82
UnknownInterface.hpp
MoFEM interface.
PCMGSubMatrixCtx_private::MOFEM_EVENT_mult
PetscLogEvent MOFEM_EVENT_mult
Definition: PCMGSetUpViaApproxOrders.cpp:82
PCMGSubMatrixCtx
Definition: PCMGSetUpViaApproxOrders.hpp:9
DMMGViaApproxOrdersGetCoarseningISSize
MoFEMErrorCode DMMGViaApproxOrdersGetCoarseningISSize(DM dm, int *size)
Gets size of coarseningIS in internal data struture DMMGViaApproxOrdersCtx.
Definition: PCMGSetUpViaApproxOrders.cpp:197
PCMGSetUpViaApproxOrdersCtx::createIsAtLevel
virtual MoFEMErrorCode createIsAtLevel(int kk, IS *is)
Set IS for levels.
Definition: PCMGSetUpViaApproxOrders.cpp:625
PCMGSubMatrixCtx_private::PCMGSubMatrixCtx_private
PCMGSubMatrixCtx_private(Mat a, IS is)
Definition: PCMGSetUpViaApproxOrders.cpp:49
PCMGSubMatrixCtx_private::sub_mat_sor
friend MoFEMErrorCode sub_mat_sor(Mat mat, Vec b, PetscReal omega, MatSORType flag, PetscReal shift, PetscInt its, PetscInt lits, Vec x)
Definition: PCMGSetUpViaApproxOrders.cpp:114
PCMGSubMatrixCtx::iS
IS iS
Definition: PCMGSetUpViaApproxOrders.hpp:12
MoFEM::Problem::getNbLocalDofsRow
DofIdx getNbLocalDofsRow() const
Definition: ProblemsMultiIndices.hpp:378
DMCreateMatrix_MGViaApproxOrders
MoFEMErrorCode DMCreateMatrix_MGViaApproxOrders(DM dm, Mat *M)
Create matrix for Multi-Grid via approximation orders.
Definition: PCMGSetUpViaApproxOrders.cpp:429
DMMGViaApproxOrdersPopBackCoarseningIS
MoFEMErrorCode DMMGViaApproxOrdersPopBackCoarseningIS(DM dm)
Pop is form MG via approximation orders.
Definition: PCMGSetUpViaApproxOrders.cpp:259
MoFEM::DMSetOperators_MoFEM
PetscErrorCode DMSetOperators_MoFEM(DM dm)
Set operators for MoFEM dm.
Definition: DMMoFEM.cpp:49
PCMGSubMatrixCtx_private::MOFEM_EVENT_sor
PetscLogEvent MOFEM_EVENT_sor
Definition: PCMGSetUpViaApproxOrders.cpp:83
PCMGSubMatrixCtx_private::initData
MoFEMErrorCode initData(Vec x)
Definition: PCMGSetUpViaApproxOrders.cpp:71
DMCreateGlobalVector_MGViaApproxOrders
MoFEMErrorCode DMCreateGlobalVector_MGViaApproxOrders(DM dm, Vec *g)
Create global vector for DMGViaApproxOrders.
Definition: PCMGSetUpViaApproxOrders.cpp:575
PCMGSetUpViaApproxOrdersCtx::orderAtLastLevel
int orderAtLastLevel
set maximal evaluated order
Definition: PCMGSetUpViaApproxOrders.hpp:206
convert.n
n
Definition: convert.py:82
PCMGSetUpViaApproxOrdersCtx::destroyIsAtLevel
virtual MoFEMErrorCode destroyIsAtLevel(int kk, IS *is)
Destroy IS if internally created.
Definition: PCMGSetUpViaApproxOrders.cpp:652
PCMGSubMatrixCtx::A
Mat A
Definition: PCMGSetUpViaApproxOrders.hpp:10
N
const int N
Definition: speed_test.cpp:3
MoFEM::UnknownInterface
base class for all interface classes
Definition: UnknownInterface.hpp:34
sub_mat_mult_add
MoFEMErrorCode sub_mat_mult_add(Mat a, Vec x, Vec f)
Definition: PCMGSetUpViaApproxOrders.cpp:110
DMRegister_MGViaApproxOrders
MoFEMErrorCode DMRegister_MGViaApproxOrders(const char sname[])
Register DM for Multi-Grid via approximation orders.
Definition: PCMGSetUpViaApproxOrders.cpp:375
PCMGSubMatrixCtx::sCat
VecScatter sCat
Definition: PCMGSetUpViaApproxOrders.hpp:13
DMMGViaApproxOrdersGetCtx
MoFEMErrorCode DMMGViaApproxOrdersGetCtx(DM dm, DMMGViaApproxOrdersCtx **ctx)
Definition: PCMGSetUpViaApproxOrders.cpp:174
DMMGViaApproxOrdersCtx::coarseningIS
std::vector< IS > coarseningIS
Coarsening IS.
Definition: PCMGSetUpViaApproxOrders.hpp:33
PCMGSetUpViaApproxOrdersCtx::verboseLevel
int verboseLevel
Definition: PCMGSetUpViaApproxOrders.hpp:209
MoFEM::DMoFEMGetInterfacePtr
PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr)
Get pointer to MoFEM::Interface.
Definition: DMMoFEM.cpp:414
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
PCMGSetUpViaApproxOrdersCtx::coarseOrder
int coarseOrder
approximation order of coarse level
Definition: PCMGSetUpViaApproxOrders.hpp:205
PCMGSubMatrixCtx_private::sub_mat_mult_generic
friend MoFEMErrorCode sub_mat_mult_generic(Mat a, Vec x, Vec f)
Definition: PCMGSetUpViaApproxOrders.cpp:88
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::DMMoFEMGetProblemPtr
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
Definition: DMMoFEM.cpp:426
PCMGSubMatrixCtx::PCMGSubMatrixCtx
PCMGSubMatrixCtx(Mat a, IS is)
Definition: PCMGSetUpViaApproxOrders.cpp:32
sub_mat_sor
MoFEMErrorCode sub_mat_sor(Mat mat, Vec b, PetscReal omega, MatSORType flag, PetscReal shift, PetscInt its, PetscInt lits, Vec x)
Definition: PCMGSetUpViaApproxOrders.cpp:114
ksp_set_operators
static MoFEMErrorCode ksp_set_operators(KSP ksp, Mat A, Mat B, void *ctx)
Definition: PCMGSetUpViaApproxOrders.cpp:403
PCMGSetUpViaApproxOrdersCtx::buildProlongationOperator
virtual MoFEMErrorCode buildProlongationOperator(bool use_mat_a, int verb=0)
Set up data structures for MG.
Definition: PCMGSetUpViaApproxOrders.cpp:659
DMCreateInterpolation_MGViaApproxOrders
MoFEMErrorCode DMCreateInterpolation_MGViaApproxOrders(DM dm1, DM dm2, Mat *mat, Vec *vec)
Create interpolation matrix between data managers dm1 and dm2.
Definition: PCMGSetUpViaApproxOrders.cpp:476
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
DMMGViaApproxOrdersReplaceCoarseningIS
MoFEMErrorCode DMMGViaApproxOrdersReplaceCoarseningIS(DM dm, IS *is_vec, int nb_elems, Mat A, int verb)
Replace coarsening IS in DM via approximation orders.
Definition: PCMGSetUpViaApproxOrders.cpp:284
DMMGViaApproxOrdersCtx::aO
AO aO
Definition: PCMGSetUpViaApproxOrders.hpp:32
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
QUIET
@ QUIET
Definition: definitions.h:208
MoFEM::Problem
keeps basic data about problem
Definition: ProblemsMultiIndices.hpp:54
MoFEM::ISManager::isCreateProblemOrder
MoFEMErrorCode isCreateProblemOrder(const std::string problem_name, RowColData rc, int min_order, int max_order, IS *is) const
create IS for given order range (collective)
Definition: ISManager.cpp:204
MoFEM::Problem::getName
auto getName() const
Definition: ProblemsMultiIndices.hpp:372
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483
PCMGSetUpViaApproxOrdersCtx::A
Mat A
Matrix at fine level.
Definition: PCMGSetUpViaApproxOrders.hpp:196
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
PCMGSubMatrixCtx_private::isInitisalised
bool isInitisalised
Definition: PCMGSetUpViaApproxOrders.cpp:84