v0.14.0
Classes | Macros | Functions
PCMGSetUpViaApproxOrders.cpp File Reference

implementation of multi-grid solver for p- adaptivity More...

#include <MoFEM.hpp>
#include <UnknownInterface.hpp>
#include <PCMGSetUpViaApproxOrders.hpp>
#include <petsc-private/petscimpl.h>
#include <petsc-private/dmimpl.h>
#include <petsc-private/vecimpl.h>

Go to the source code of this file.

Classes

struct  PCMGSubMatrixCtx_private
 

Macros

#define PETSC_VERSION_RELEASE   1
 
#define GET_DM_FIELD(DM)
 

Functions

template<InsertMode MODE>
MoFEMErrorCode sub_mat_mult_generic (Mat a, Vec x, Vec f)
 
MoFEMErrorCode sub_mat_mult (Mat a, Vec x, Vec f)
 
MoFEMErrorCode sub_mat_mult_add (Mat a, Vec x, Vec f)
 
MoFEMErrorCode sub_mat_sor (Mat mat, Vec b, PetscReal omega, MatSORType flag, PetscReal shift, PetscInt its, PetscInt lits, Vec x)
 
MoFEMErrorCode DMMGViaApproxOrdersGetCtx (DM dm, DMMGViaApproxOrdersCtx **ctx)
 
MoFEMErrorCode DMMGViaApproxOrdersSetAO (DM dm, AO ao)
 Set DM ordering. More...
 
MoFEMErrorCode DMMGViaApproxOrdersGetCoarseningISSize (DM dm, int *size)
 Gets size of coarseningIS in internal data struture DMMGViaApproxOrdersCtx. More...
 
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. More...
 
MoFEMErrorCode DMMGViaApproxOrdersPopBackCoarseningIS (DM dm)
 Pop is form MG via approximation orders. More...
 
MoFEMErrorCode DMMGViaApproxOrdersClearCoarseningIS (DM dm)
 Clear approximation orders. More...
 
MoFEMErrorCode DMMGViaApproxOrdersReplaceCoarseningIS (DM dm, IS *is_vec, int nb_elems, Mat A, int verb)
 Replace coarsening IS in DM via approximation orders. More...
 
MoFEMErrorCode DMMGViaApproxOrdersGetCtx (DM dm, const DMMGViaApproxOrdersCtx **ctx)
 Get context for DM via approximation orders. More...
 
MoFEMErrorCode DMRegister_MGViaApproxOrders (const char sname[])
 Register DM for Multi-Grid via approximation orders. More...
 
static MoFEMErrorCode ksp_set_operators (KSP ksp, Mat A, Mat B, void *ctx)
 
MoFEMErrorCode DMCreate_MGViaApproxOrders (DM dm)
 Create DM data structure for Multi-Grid via approximation orders. More...
 
MoFEMErrorCode DMCreateMatrix_MGViaApproxOrders (DM dm, Mat *M)
 Create matrix for Multi-Grid via approximation orders. More...
 
MoFEMErrorCode DMCoarsen_MGViaApproxOrders (DM dm, MPI_Comm comm, DM *dmc)
 Coarsen DM. More...
 
MoFEMErrorCode DMCreateInterpolation_MGViaApproxOrders (DM dm1, DM dm2, Mat *mat, Vec *vec)
 Create interpolation matrix between data managers dm1 and dm2. More...
 
MoFEMErrorCode DMCreateGlobalVector_MGViaApproxOrders (DM dm, Vec *g)
 Create global vector for DMGViaApproxOrders. More...
 
MoFEMErrorCode PCMGSetUpViaApproxOrders (PC pc, PCMGSetUpViaApproxOrdersCtx *ctx, int verb)
 Function build MG structure. More...
 

Detailed Description

implementation of multi-grid solver for p- adaptivity

Definition in file PCMGSetUpViaApproxOrders.cpp.

Macro Definition Documentation

◆ GET_DM_FIELD

#define GET_DM_FIELD (   DM)
Value:
auto dm_field = \
static_cast<DMCtx *>(DM->data)->getInterface<DMMGViaApproxOrdersCtx>(); \
NOT_USED(dm_field)

Definition at line 169 of file PCMGSetUpViaApproxOrders.cpp.

◆ PETSC_VERSION_RELEASE

#define PETSC_VERSION_RELEASE   1

Definition at line 16 of file PCMGSetUpViaApproxOrders.cpp.

Function Documentation

◆ DMCreate_MGViaApproxOrders()

MoFEMErrorCode DMCreate_MGViaApproxOrders ( DM  dm)

Create DM data structure for Multi-Grid via approximation orders.

It set data structure and operators needed

Parameters
dmDiscrete manager
Returns
Error code

Definition at line 408 of file PCMGSetUpViaApproxOrders.cpp.

408  {
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 }

◆ DMCreateGlobalVector_MGViaApproxOrders()

MoFEMErrorCode DMCreateGlobalVector_MGViaApproxOrders ( DM  dm,
Vec *  g 
)

Create global vector for DMGViaApproxOrders.

Parameters
dmDistributed mesh data structure
greturned pointer to vector
Returns
Error code

Definition at line 575 of file PCMGSetUpViaApproxOrders.cpp.

575  {
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 }

◆ DMCreateInterpolation_MGViaApproxOrders()

MoFEMErrorCode DMCreateInterpolation_MGViaApproxOrders ( DM  dm1,
DM  dm2,
Mat *  mat,
Vec *  vec 
)

Create interpolation matrix between data managers dm1 and dm2.

Parameters
dm1Distributed mesh data structure
dm2Distributed mesh data structure
matPointer to returned interpolation matrix
vecPointer to scaling vector here returned NULL
Returns
Error code

Definition at line 476 of file PCMGSetUpViaApproxOrders.cpp.

477  {
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 }

◆ DMMGViaApproxOrdersGetCoarseningISSize()

MoFEMErrorCode DMMGViaApproxOrdersGetCoarseningISSize ( DM  dm,
int *  size 
)

Gets size of coarseningIS in internal data struture DMMGViaApproxOrdersCtx.

Parameters
dmDM
sizesize of coarseningIS
Returns
Error code

Definition at line 197 of file PCMGSetUpViaApproxOrders.cpp.

197  {
198  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
200  GET_DM_FIELD(dm);
201  *size = dm_field->coarseningIS.size();
203 }

◆ DMMGViaApproxOrdersGetCtx() [1/2]

MoFEMErrorCode DMMGViaApproxOrdersGetCtx ( DM  dm,
const DMMGViaApproxOrdersCtx **  ctx 
)

Get context for DM via approximation orders.

Parameters
dmthe DM object
ctxdata context
Returns
error code

Definition at line 366 of file PCMGSetUpViaApproxOrders.cpp.

367  {
368  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
370  GET_DM_FIELD(dm);
371  *ctx = dm_field;
373 }

◆ DMMGViaApproxOrdersGetCtx() [2/2]

MoFEMErrorCode DMMGViaApproxOrdersGetCtx ( DM  dm,
DMMGViaApproxOrdersCtx **  ctx 
)

Get DM Ctx

Definition at line 174 of file PCMGSetUpViaApproxOrders.cpp.

174  {
175  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
177  GET_DM_FIELD(dm);
178  *ctx = dm_field;
180 }

◆ DMMGViaApproxOrdersReplaceCoarseningIS()

MoFEMErrorCode DMMGViaApproxOrdersReplaceCoarseningIS ( DM  dm,
IS *  is_vec,
int  nb_elems,
Mat  A,
int  verb = 0 
)

Replace coarsening IS in DM via approximation orders.

Parameters
dmdm
is_vecPointer to vector of is
nb_elemsNumber of elements
AFine matrix
Returns
Error code

Definition at line 284 of file PCMGSetUpViaApproxOrders.cpp.

286  {
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 }

◆ DMMGViaApproxOrdersSetAO()

MoFEMErrorCode DMMGViaApproxOrdersSetAO ( DM  dm,
AO  ao 
)

Set DM ordering.

IS can be given is some other ordering, AO will transform indices from coarseningIS ordering to ordering used to construct fine matrix.

Parameters
dm[description]
ao[description]
Returns
[description]

Definition at line 182 of file PCMGSetUpViaApproxOrders.cpp.

182  {
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 }

◆ ksp_set_operators()

static MoFEMErrorCode ksp_set_operators ( KSP  ksp,
Mat  A,
Mat  B,
void *  ctx 
)
static

Definition at line 403 of file PCMGSetUpViaApproxOrders.cpp.

403  {
406 }

◆ PCMGSetUpViaApproxOrders()

MoFEMErrorCode PCMGSetUpViaApproxOrders ( PC  pc,
PCMGSetUpViaApproxOrdersCtx ctx,
int  verb = 0 
)

Function build MG structure.

Parameters
pcMG pre-conditioner http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCMG.html
ctxMoFEM data structure for MG
verbverbosity level
Returns
error code
Examples
elasticity.cpp.

Definition at line 727 of file PCMGSetUpViaApproxOrders.cpp.

728  {
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 }

◆ sub_mat_mult()

MoFEMErrorCode sub_mat_mult ( Mat  a,
Vec  x,
Vec  f 
)

Definition at line 106 of file PCMGSetUpViaApproxOrders.cpp.

106  {
107  return sub_mat_mult_generic<INSERT_VALUES>(a, x, f);
108 }

◆ sub_mat_mult_add()

MoFEMErrorCode sub_mat_mult_add ( Mat  a,
Vec  x,
Vec  f 
)

Definition at line 110 of file PCMGSetUpViaApproxOrders.cpp.

110  {
111  return sub_mat_mult_generic<ADD_VALUES>(a, x, f);
112 }

◆ sub_mat_mult_generic()

template<InsertMode MODE>
MoFEMErrorCode sub_mat_mult_generic ( Mat  a,
Vec  x,
Vec  f 
)

Definition at line 88 of file PCMGSetUpViaApproxOrders.cpp.

88  {
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 }

◆ sub_mat_sor()

MoFEMErrorCode sub_mat_sor ( Mat  mat,
Vec  b,
PetscReal  omega,
MatSORType  flag,
PetscReal  shift,
PetscInt  its,
PetscInt  lits,
Vec  x 
)

Definition at line 114 of file PCMGSetUpViaApproxOrders.cpp.

116  {
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 }
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
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
DMCoarsen_MGViaApproxOrders
MoFEMErrorCode DMCoarsen_MGViaApproxOrders(DM dm, MPI_Comm comm, DM *dmc)
Coarsen DM.
Definition: PCMGSetUpViaApproxOrders.cpp:461
DMMGViaApproxOrdersCtx
Structure for DM for multi-grid via approximation orders.
Definition: PCMGSetUpViaApproxOrders.hpp:22
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
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::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
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
a
constexpr double a
Definition: approx_sphere.cpp:30
PCMGSubMatrixCtx::F
Vec F
Definition: PCMGSetUpViaApproxOrders.hpp:11
PCMGSubMatrixCtx_private
Definition: PCMGSetUpViaApproxOrders.cpp:48
PCMGSubMatrixCtx_private::MOFEM_EVENT_mult
PetscLogEvent MOFEM_EVENT_mult
Definition: PCMGSetUpViaApproxOrders.cpp:82
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
convert.n
n
Definition: convert.py:82
PCMGSubMatrixCtx::A
Mat A
Definition: PCMGSetUpViaApproxOrders.hpp:10
N
const int N
Definition: speed_test.cpp:3
PCMGSubMatrixCtx::sCat
VecScatter sCat
Definition: PCMGSetUpViaApproxOrders.hpp:13
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
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
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
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