Actual source code: dm.c

  1: #include <petscvec.h>
  2: #include <petsc/private/dmimpl.h>
  3: #include <petsc/private/dmlabelimpl.h>
  4: #include <petsc/private/petscdsimpl.h>
  5: #include <petscdmplex.h>
  6: #include <petscdmceed.h>
  7: #include <petscdmfield.h>
  8: #include <petscsf.h>
  9: #include <petscds.h>

 11: #ifdef PETSC_HAVE_LIBCEED
 12: #include <petscfeceed.h>
 13: #endif

 15: PetscClassId DM_CLASSID;
 16: PetscClassId DMLABEL_CLASSID;
 17: PetscLogEvent DM_Convert, DM_GlobalToLocal, DM_LocalToGlobal, DM_LocalToLocal, DM_LocatePoints, DM_Coarsen, DM_Refine, DM_CreateInterpolation, DM_CreateRestriction, DM_CreateInjection, DM_CreateMatrix, DM_CreateMassMatrix, DM_Load, DM_View, DM_AdaptInterpolator, DM_ProjectFunction;

 19: const char *const DMBoundaryTypes[]          = {"NONE", "GHOSTED", "MIRROR", "PERIODIC", "TWIST", "DMBoundaryType", "DM_BOUNDARY_", NULL};
 20: const char *const DMBoundaryConditionTypes[] = {"INVALID", "ESSENTIAL", "NATURAL", "INVALID", "INVALID", "ESSENTIAL_FIELD", "NATURAL_FIELD", "INVALID", "INVALID", "ESSENTIAL_BD_FIELD", "NATURAL_RIEMANN", "DMBoundaryConditionType", "DM_BC_", NULL};
 21: const char *const DMBlockingTypes[]          = {"TOPOLOGICAL_POINT", "FIELD_NODE", "DMBlockingType", "DM_BLOCKING_", NULL};
 22: const char *const DMPolytopeTypes[] =
 23:   {"vertex",  "segment",      "tensor_segment", "triangle", "quadrilateral",  "tensor_quad",  "tetrahedron", "hexahedron", "triangular_prism", "tensor_triangular_prism", "tensor_quadrilateral_prism", "pyramid", "FV_ghost_cell", "interior_ghost_cell",
 24:    "unknown", "unknown_cell", "unknown_face",   "invalid",  "DMPolytopeType", "DM_POLYTOPE_", NULL};
 25: const char *const DMCopyLabelsModes[] = {"replace", "keep", "fail", "DMCopyLabelsMode", "DM_COPY_LABELS_", NULL};

 27: /*@
 28:   DMCreate - Creates an empty `DM` object. `DM`s are the abstract objects in PETSc that mediate between meshes and discretizations and the
 29:   algebraic solvers, time integrators, and optimization algorithms.

 31:   Collective

 33:   Input Parameter:
 34: . comm - The communicator for the `DM` object

 36:   Output Parameter:
 37: . dm - The `DM` object

 39:   Level: beginner

 41:   Notes:
 42:   See `DMType` for a brief summary of available `DM`.

 44:   The type must then be set with `DMSetType()`. If you never call `DMSetType()` it will generate an
 45:   error when you try to use the dm.

 47: .seealso: [](ch_dmbase), `DM`, `DMSetType()`, `DMType`, `DMDACreate()`, `DMDA`, `DMSLICED`, `DMCOMPOSITE`, `DMPLEX`, `DMMOAB`, `DMNETWORK`
 48: @*/
 49: PetscErrorCode DMCreate(MPI_Comm comm, DM *dm)
 50: {
 51:   DM      v;
 52:   PetscDS ds;

 54:   PetscFunctionBegin;
 55:   PetscAssertPointer(dm, 2);

 57:   PetscCall(DMInitializePackage());
 58:   PetscCall(PetscHeaderCreate(v, DM_CLASSID, "DM", "Distribution Manager", "DM", comm, DMDestroy, DMView));
 59:   ((PetscObject)v)->non_cyclic_references = &DMCountNonCyclicReferences;
 60:   v->setupcalled                          = PETSC_FALSE;
 61:   v->setfromoptionscalled                 = PETSC_FALSE;
 62:   v->ltogmap                              = NULL;
 63:   v->bind_below                           = 0;
 64:   v->bs                                   = 1;
 65:   v->coloringtype                         = IS_COLORING_GLOBAL;
 66:   PetscCall(PetscSFCreate(comm, &v->sf));
 67:   PetscCall(PetscSFCreate(comm, &v->sectionSF));
 68:   v->labels                    = NULL;
 69:   v->adjacency[0]              = PETSC_FALSE;
 70:   v->adjacency[1]              = PETSC_TRUE;
 71:   v->depthLabel                = NULL;
 72:   v->celltypeLabel             = NULL;
 73:   v->localSection              = NULL;
 74:   v->globalSection             = NULL;
 75:   v->defaultConstraint.section = NULL;
 76:   v->defaultConstraint.mat     = NULL;
 77:   v->defaultConstraint.bias    = NULL;
 78:   v->coordinates[0].dim        = PETSC_DEFAULT;
 79:   v->coordinates[1].dim        = PETSC_DEFAULT;
 80:   v->sparseLocalize            = PETSC_TRUE;
 81:   v->dim                       = PETSC_DETERMINE;
 82:   {
 83:     PetscInt i;
 84:     for (i = 0; i < 10; ++i) {
 85:       v->nullspaceConstructors[i]     = NULL;
 86:       v->nearnullspaceConstructors[i] = NULL;
 87:     }
 88:   }
 89:   PetscCall(PetscDSCreate(PETSC_COMM_SELF, &ds));
 90:   PetscCall(DMSetRegionDS(v, NULL, NULL, ds, NULL));
 91:   PetscCall(PetscDSDestroy(&ds));
 92:   PetscCall(PetscHMapAuxCreate(&v->auxData));
 93:   v->dmBC              = NULL;
 94:   v->coarseMesh        = NULL;
 95:   v->outputSequenceNum = -1;
 96:   v->outputSequenceVal = 0.0;
 97:   PetscCall(DMSetVecType(v, VECSTANDARD));
 98:   PetscCall(DMSetMatType(v, MATAIJ));

100:   *dm = v;
101:   PetscFunctionReturn(PETSC_SUCCESS);
102: }

104: /*@
105:   DMClone - Creates a `DM` object with the same topology as the original.

107:   Collective

109:   Input Parameter:
110: . dm - The original `DM` object

112:   Output Parameter:
113: . newdm - The new `DM` object

115:   Level: beginner

117:   Notes:
118:   For some `DM` implementations this is a shallow clone, the result of which may share (reference counted) information with its parent. For example,
119:   `DMClone()` applied to a `DMPLEX` object will result in a new `DMPLEX` that shares the topology with the original `DMPLEX`. It does not
120:   share the `PetscSection` of the original `DM`.

122:   The clone is considered set up if the original has been set up.

124:   Use `DMConvert()` for a general way to create new `DM` from a given `DM`

126: .seealso: [](ch_dmbase), `DM`, `DMDestroy()`, `DMCreate()`, `DMSetType()`, `DMSetLocalSection()`, `DMSetGlobalSection()`, `DMPLEX`, `DMConvert()`
127: @*/
128: PetscErrorCode DMClone(DM dm, DM *newdm)
129: {
130:   PetscSF              sf;
131:   Vec                  coords;
132:   void                *ctx;
133:   MatOrderingType      otype;
134:   DMReorderDefaultFlag flg;
135:   PetscInt             dim, cdim, i;

137:   PetscFunctionBegin;
139:   PetscAssertPointer(newdm, 2);
140:   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), newdm));
141:   PetscCall(DMCopyLabels(dm, *newdm, PETSC_COPY_VALUES, PETSC_TRUE, DM_COPY_LABELS_FAIL));
142:   (*newdm)->leveldown     = dm->leveldown;
143:   (*newdm)->levelup       = dm->levelup;
144:   (*newdm)->prealloc_only = dm->prealloc_only;
145:   (*newdm)->prealloc_skip = dm->prealloc_skip;
146:   PetscCall(PetscFree((*newdm)->vectype));
147:   PetscCall(PetscStrallocpy(dm->vectype, (char **)&(*newdm)->vectype));
148:   PetscCall(PetscFree((*newdm)->mattype));
149:   PetscCall(PetscStrallocpy(dm->mattype, (char **)&(*newdm)->mattype));
150:   PetscCall(DMGetDimension(dm, &dim));
151:   PetscCall(DMSetDimension(*newdm, dim));
152:   PetscTryTypeMethod(dm, clone, newdm);
153:   (*newdm)->setupcalled = dm->setupcalled;
154:   PetscCall(DMGetPointSF(dm, &sf));
155:   PetscCall(DMSetPointSF(*newdm, sf));
156:   PetscCall(DMGetApplicationContext(dm, &ctx));
157:   PetscCall(DMSetApplicationContext(*newdm, ctx));
158:   PetscCall(DMReorderSectionGetDefault(dm, &flg));
159:   PetscCall(DMReorderSectionSetDefault(*newdm, flg));
160:   PetscCall(DMReorderSectionGetType(dm, &otype));
161:   PetscCall(DMReorderSectionSetType(*newdm, otype));
162:   for (i = 0; i < 2; ++i) {
163:     if (dm->coordinates[i].dm) {
164:       DM           ncdm;
165:       PetscSection cs;
166:       PetscInt     pEnd = -1, pEndMax = -1;

168:       PetscCall(DMGetLocalSection(dm->coordinates[i].dm, &cs));
169:       if (cs) PetscCall(PetscSectionGetChart(cs, NULL, &pEnd));
170:       PetscCallMPI(MPIU_Allreduce(&pEnd, &pEndMax, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
171:       if (pEndMax >= 0) {
172:         PetscCall(DMClone(dm->coordinates[i].dm, &ncdm));
173:         PetscCall(DMCopyDisc(dm->coordinates[i].dm, ncdm));
174:         PetscCall(DMSetLocalSection(ncdm, cs));
175:         if (dm->coordinates[i].dm->periodic.setup) {
176:           ncdm->periodic.setup = dm->coordinates[i].dm->periodic.setup;
177:           PetscCall(ncdm->periodic.setup(ncdm));
178:         }
179:         if (i) PetscCall(DMSetCellCoordinateDM(*newdm, ncdm));
180:         else PetscCall(DMSetCoordinateDM(*newdm, ncdm));
181:         PetscCall(DMDestroy(&ncdm));
182:       }
183:     }
184:   }
185:   PetscCall(DMGetCoordinateDim(dm, &cdim));
186:   PetscCall(DMSetCoordinateDim(*newdm, cdim));
187:   PetscCall(DMGetCoordinatesLocal(dm, &coords));
188:   if (coords) {
189:     PetscCall(DMSetCoordinatesLocal(*newdm, coords));
190:   } else {
191:     PetscCall(DMGetCoordinates(dm, &coords));
192:     if (coords) PetscCall(DMSetCoordinates(*newdm, coords));
193:   }
194:   PetscCall(DMGetCellCoordinatesLocal(dm, &coords));
195:   if (coords) {
196:     PetscCall(DMSetCellCoordinatesLocal(*newdm, coords));
197:   } else {
198:     PetscCall(DMGetCellCoordinates(dm, &coords));
199:     if (coords) PetscCall(DMSetCellCoordinates(*newdm, coords));
200:   }
201:   {
202:     const PetscReal *maxCell, *Lstart, *L;

204:     PetscCall(DMGetPeriodicity(dm, &maxCell, &Lstart, &L));
205:     PetscCall(DMSetPeriodicity(*newdm, maxCell, Lstart, L));
206:   }
207:   {
208:     PetscBool useCone, useClosure;

210:     PetscCall(DMGetAdjacency(dm, PETSC_DEFAULT, &useCone, &useClosure));
211:     PetscCall(DMSetAdjacency(*newdm, PETSC_DEFAULT, useCone, useClosure));
212:   }
213:   PetscFunctionReturn(PETSC_SUCCESS);
214: }

216: /*@
217:   DMSetVecType - Sets the type of vector to be created with `DMCreateLocalVector()` and `DMCreateGlobalVector()`

219:   Logically Collective

221:   Input Parameters:
222: + dm    - initial distributed array
223: - ctype - the vector type, for example `VECSTANDARD`, `VECCUDA`, or `VECVIENNACL`

225:   Options Database Key:
226: . -dm_vec_type ctype - the type of vector to create

228:   Level: intermediate

230: .seealso: [](ch_dmbase), `DM`, `DMCreate()`, `DMDestroy()`, `DMDAInterpolationType`, `VecType`, `DMGetVecType()`, `DMSetMatType()`, `DMGetMatType()`,
231:           `VECSTANDARD`, `VECCUDA`, `VECVIENNACL`, `DMCreateLocalVector()`, `DMCreateGlobalVector()`
232: @*/
233: PetscErrorCode DMSetVecType(DM dm, VecType ctype)
234: {
235:   char *tmp;

237:   PetscFunctionBegin;
239:   PetscAssertPointer(ctype, 2);
240:   tmp = (char *)dm->vectype;
241:   PetscCall(PetscStrallocpy(ctype, (char **)&dm->vectype));
242:   PetscCall(PetscFree(tmp));
243:   PetscFunctionReturn(PETSC_SUCCESS);
244: }

246: /*@
247:   DMGetVecType - Gets the type of vector created with `DMCreateLocalVector()` and `DMCreateGlobalVector()`

249:   Logically Collective

251:   Input Parameter:
252: . da - initial distributed array

254:   Output Parameter:
255: . ctype - the vector type

257:   Level: intermediate

259: .seealso: [](ch_dmbase), `DM`, `DMCreate()`, `DMDestroy()`, `DMDAInterpolationType`, `VecType`, `DMSetMatType()`, `DMGetMatType()`, `DMSetVecType()`
260: @*/
261: PetscErrorCode DMGetVecType(DM da, VecType *ctype)
262: {
263:   PetscFunctionBegin;
265:   *ctype = da->vectype;
266:   PetscFunctionReturn(PETSC_SUCCESS);
267: }

269: /*@
270:   VecGetDM - Gets the `DM` defining the data layout of the vector

272:   Not Collective

274:   Input Parameter:
275: . v - The `Vec`

277:   Output Parameter:
278: . dm - The `DM`

280:   Level: intermediate

282:   Note:
283:   A `Vec` may not have a `DM` associated with it.

285: .seealso: [](ch_dmbase), `DM`, `VecSetDM()`, `DMGetLocalVector()`, `DMGetGlobalVector()`, `DMSetVecType()`
286: @*/
287: PetscErrorCode VecGetDM(Vec v, DM *dm)
288: {
289:   PetscFunctionBegin;
291:   PetscAssertPointer(dm, 2);
292:   PetscCall(PetscObjectQuery((PetscObject)v, "__PETSc_dm", (PetscObject *)dm));
293:   PetscFunctionReturn(PETSC_SUCCESS);
294: }

296: /*@
297:   VecSetDM - Sets the `DM` defining the data layout of the vector.

299:   Not Collective

301:   Input Parameters:
302: + v  - The `Vec`
303: - dm - The `DM`

305:   Level: developer

307:   Notes:
308:   This is rarely used, generally one uses `DMGetLocalVector()` or  `DMGetGlobalVector()` to create a vector associated with a given `DM`

310:   This is NOT the same as `DMCreateGlobalVector()` since it does not change the view methods or perform other customization, but merely sets the `DM` member.

312: .seealso: [](ch_dmbase), `DM`, `VecGetDM()`, `DMGetLocalVector()`, `DMGetGlobalVector()`, `DMSetVecType()`
313: @*/
314: PetscErrorCode VecSetDM(Vec v, DM dm)
315: {
316:   PetscFunctionBegin;
319:   PetscCall(PetscObjectCompose((PetscObject)v, "__PETSc_dm", (PetscObject)dm));
320:   PetscFunctionReturn(PETSC_SUCCESS);
321: }

323: /*@
324:   DMSetISColoringType - Sets the type of coloring, `IS_COLORING_GLOBAL` or `IS_COLORING_LOCAL` that is created by the `DM`

326:   Logically Collective

328:   Input Parameters:
329: + dm    - the `DM` context
330: - ctype - the matrix type

332:   Options Database Key:
333: . -dm_is_coloring_type - global or local

335:   Level: intermediate

337: .seealso: [](ch_dmbase), `DM`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMSetMatrixPreallocateOnly()`, `MatType`, `DMGetMatType()`,
338:           `DMGetISColoringType()`, `ISColoringType`, `IS_COLORING_GLOBAL`, `IS_COLORING_LOCAL`
339: @*/
340: PetscErrorCode DMSetISColoringType(DM dm, ISColoringType ctype)
341: {
342:   PetscFunctionBegin;
344:   dm->coloringtype = ctype;
345:   PetscFunctionReturn(PETSC_SUCCESS);
346: }

348: /*@
349:   DMGetISColoringType - Gets the type of coloring, `IS_COLORING_GLOBAL` or `IS_COLORING_LOCAL` that is created by the `DM`

351:   Logically Collective

353:   Input Parameter:
354: . dm - the `DM` context

356:   Output Parameter:
357: . ctype - the matrix type

359:   Options Database Key:
360: . -dm_is_coloring_type - global or local

362:   Level: intermediate

364: .seealso: [](ch_dmbase), `DM`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMSetMatrixPreallocateOnly()`, `MatType`, `DMGetMatType()`,
365:           `ISColoringType`, `IS_COLORING_GLOBAL`, `IS_COLORING_LOCAL`
366: @*/
367: PetscErrorCode DMGetISColoringType(DM dm, ISColoringType *ctype)
368: {
369:   PetscFunctionBegin;
371:   *ctype = dm->coloringtype;
372:   PetscFunctionReturn(PETSC_SUCCESS);
373: }

375: /*@
376:   DMSetMatType - Sets the type of matrix created with `DMCreateMatrix()`

378:   Logically Collective

380:   Input Parameters:
381: + dm    - the `DM` context
382: - ctype - the matrix type, for example `MATMPIAIJ`

384:   Options Database Key:
385: . -dm_mat_type ctype - the type of the matrix to create, for example mpiaij

387:   Level: intermediate

389: .seealso: [](ch_dmbase), `DM`, `MatType`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMSetMatrixPreallocateOnly()`, `DMGetMatType()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
390: @*/
391: PetscErrorCode DMSetMatType(DM dm, MatType ctype)
392: {
393:   char *tmp;

395:   PetscFunctionBegin;
397:   PetscAssertPointer(ctype, 2);
398:   tmp = (char *)dm->mattype;
399:   PetscCall(PetscStrallocpy(ctype, (char **)&dm->mattype));
400:   PetscCall(PetscFree(tmp));
401:   PetscFunctionReturn(PETSC_SUCCESS);
402: }

404: /*@
405:   DMGetMatType - Gets the type of matrix that would be created with `DMCreateMatrix()`

407:   Logically Collective

409:   Input Parameter:
410: . dm - the `DM` context

412:   Output Parameter:
413: . ctype - the matrix type

415:   Level: intermediate

417: .seealso: [](ch_dmbase), `DM`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMSetMatrixPreallocateOnly()`, `MatType`, `DMSetMatType()`
418: @*/
419: PetscErrorCode DMGetMatType(DM dm, MatType *ctype)
420: {
421:   PetscFunctionBegin;
423:   *ctype = dm->mattype;
424:   PetscFunctionReturn(PETSC_SUCCESS);
425: }

427: /*@
428:   MatGetDM - Gets the `DM` defining the data layout of the matrix

430:   Not Collective

432:   Input Parameter:
433: . A - The `Mat`

435:   Output Parameter:
436: . dm - The `DM`

438:   Level: intermediate

440:   Note:
441:   A matrix may not have a `DM` associated with it

443:   Developer Note:
444:   Since the `Mat` class doesn't know about the `DM` class the `DM` object is associated with the `Mat` through a `PetscObjectCompose()` operation

446: .seealso: [](ch_dmbase), `DM`, `MatSetDM()`, `DMCreateMatrix()`, `DMSetMatType()`
447: @*/
448: PetscErrorCode MatGetDM(Mat A, DM *dm)
449: {
450:   PetscFunctionBegin;
452:   PetscAssertPointer(dm, 2);
453:   PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_dm", (PetscObject *)dm));
454:   PetscFunctionReturn(PETSC_SUCCESS);
455: }

457: /*@
458:   MatSetDM - Sets the `DM` defining the data layout of the matrix

460:   Not Collective

462:   Input Parameters:
463: + A  - The `Mat`
464: - dm - The `DM`

466:   Level: developer

468:   Note:
469:   This is rarely used in practice, rather `DMCreateMatrix()` is used to create a matrix associated with a particular `DM`

471:   Developer Note:
472:   Since the `Mat` class doesn't know about the `DM` class the `DM` object is associated with
473:   the `Mat` through a `PetscObjectCompose()` operation

475: .seealso: [](ch_dmbase), `DM`, `MatGetDM()`, `DMCreateMatrix()`, `DMSetMatType()`
476: @*/
477: PetscErrorCode MatSetDM(Mat A, DM dm)
478: {
479:   PetscFunctionBegin;
482:   PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_dm", (PetscObject)dm));
483:   PetscFunctionReturn(PETSC_SUCCESS);
484: }

486: /*@
487:   DMSetOptionsPrefix - Sets the prefix prepended to all option names when searching through the options database

489:   Logically Collective

491:   Input Parameters:
492: + dm     - the `DM` context
493: - prefix - the prefix to prepend

495:   Level: advanced

497:   Note:
498:   A hyphen (-) must NOT be given at the beginning of the prefix name.
499:   The first character of all runtime options is AUTOMATICALLY the hyphen.

501: .seealso: [](ch_dmbase), `DM`, `PetscObjectSetOptionsPrefix()`, `DMSetFromOptions()`
502: @*/
503: PetscErrorCode DMSetOptionsPrefix(DM dm, const char prefix[])
504: {
505:   PetscFunctionBegin;
507:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, prefix));
508:   if (dm->sf) PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm->sf, prefix));
509:   if (dm->sectionSF) PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm->sectionSF, prefix));
510:   PetscFunctionReturn(PETSC_SUCCESS);
511: }

513: /*@
514:   DMAppendOptionsPrefix - Appends an additional string to an already existing prefix used for searching for
515:   `DM` options in the options database.

517:   Logically Collective

519:   Input Parameters:
520: + dm     - the `DM` context
521: - prefix - the string to append to the current prefix

523:   Level: advanced

525:   Note:
526:   If the `DM` does not currently have an options prefix then this value is used alone as the prefix as if `DMSetOptionsPrefix()` had been called.
527:   A hyphen (-) must NOT be given at the beginning of the prefix name.
528:   The first character of all runtime options is AUTOMATICALLY the hyphen.

530: .seealso: [](ch_dmbase), `DM`, `DMSetOptionsPrefix()`, `DMGetOptionsPrefix()`, `PetscObjectAppendOptionsPrefix()`, `DMSetFromOptions()`
531: @*/
532: PetscErrorCode DMAppendOptionsPrefix(DM dm, const char prefix[])
533: {
534:   PetscFunctionBegin;
536:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)dm, prefix));
537:   PetscFunctionReturn(PETSC_SUCCESS);
538: }

540: /*@
541:   DMGetOptionsPrefix - Gets the prefix used for searching for all
542:   DM options in the options database.

544:   Not Collective

546:   Input Parameter:
547: . dm - the `DM` context

549:   Output Parameter:
550: . prefix - pointer to the prefix string used is returned

552:   Level: advanced

554:   Fortran Note:
555:   Pass in a string 'prefix' of
556:   sufficient length to hold the prefix.

558: .seealso: [](ch_dmbase), `DM`, `DMSetOptionsPrefix()`, `DMAppendOptionsPrefix()`, `DMSetFromOptions()`
559: @*/
560: PetscErrorCode DMGetOptionsPrefix(DM dm, const char *prefix[])
561: {
562:   PetscFunctionBegin;
564:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, prefix));
565:   PetscFunctionReturn(PETSC_SUCCESS);
566: }

568: static PetscErrorCode DMCountNonCyclicReferences_Internal(DM dm, PetscBool recurseCoarse, PetscBool recurseFine, PetscInt *ncrefct)
569: {
570:   PetscInt refct = ((PetscObject)dm)->refct;

572:   PetscFunctionBegin;
573:   *ncrefct = 0;
574:   if (dm->coarseMesh && dm->coarseMesh->fineMesh == dm) {
575:     refct--;
576:     if (recurseCoarse) {
577:       PetscInt coarseCount;

579:       PetscCall(DMCountNonCyclicReferences_Internal(dm->coarseMesh, PETSC_TRUE, PETSC_FALSE, &coarseCount));
580:       refct += coarseCount;
581:     }
582:   }
583:   if (dm->fineMesh && dm->fineMesh->coarseMesh == dm) {
584:     refct--;
585:     if (recurseFine) {
586:       PetscInt fineCount;

588:       PetscCall(DMCountNonCyclicReferences_Internal(dm->fineMesh, PETSC_FALSE, PETSC_TRUE, &fineCount));
589:       refct += fineCount;
590:     }
591:   }
592:   *ncrefct = refct;
593:   PetscFunctionReturn(PETSC_SUCCESS);
594: }

596: /* Generic wrapper for DMCountNonCyclicReferences_Internal() */
597: PetscErrorCode DMCountNonCyclicReferences(PetscObject dm, PetscInt *ncrefct)
598: {
599:   PetscFunctionBegin;
600:   PetscCall(DMCountNonCyclicReferences_Internal((DM)dm, PETSC_TRUE, PETSC_TRUE, ncrefct));
601:   PetscFunctionReturn(PETSC_SUCCESS);
602: }

604: PetscErrorCode DMDestroyLabelLinkList_Internal(DM dm)
605: {
606:   DMLabelLink next = dm->labels;

608:   PetscFunctionBegin;
609:   /* destroy the labels */
610:   while (next) {
611:     DMLabelLink tmp = next->next;

613:     if (next->label == dm->depthLabel) dm->depthLabel = NULL;
614:     if (next->label == dm->celltypeLabel) dm->celltypeLabel = NULL;
615:     PetscCall(DMLabelDestroy(&next->label));
616:     PetscCall(PetscFree(next));
617:     next = tmp;
618:   }
619:   dm->labels = NULL;
620:   PetscFunctionReturn(PETSC_SUCCESS);
621: }

623: static PetscErrorCode DMDestroyCoordinates_Private(DMCoordinates *c)
624: {
625:   PetscFunctionBegin;
626:   c->dim = PETSC_DEFAULT;
627:   PetscCall(DMDestroy(&c->dm));
628:   PetscCall(VecDestroy(&c->x));
629:   PetscCall(VecDestroy(&c->xl));
630:   PetscCall(DMFieldDestroy(&c->field));
631:   PetscFunctionReturn(PETSC_SUCCESS);
632: }

634: /*@
635:   DMDestroy - Destroys a `DM`.

637:   Collective

639:   Input Parameter:
640: . dm - the `DM` object to destroy

642:   Level: developer

644: .seealso: [](ch_dmbase), `DM`, `DMCreate()`, `DMType`, `DMSetType()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`
645: @*/
646: PetscErrorCode DMDestroy(DM *dm)
647: {
648:   PetscInt cnt;

650:   PetscFunctionBegin;
651:   if (!*dm) PetscFunctionReturn(PETSC_SUCCESS);

654:   /* count all non-cyclic references in the doubly-linked list of coarse<->fine meshes */
655:   PetscCall(DMCountNonCyclicReferences_Internal(*dm, PETSC_TRUE, PETSC_TRUE, &cnt));
656:   --((PetscObject)*dm)->refct;
657:   if (--cnt > 0) {
658:     *dm = NULL;
659:     PetscFunctionReturn(PETSC_SUCCESS);
660:   }
661:   if (((PetscObject)*dm)->refct < 0) PetscFunctionReturn(PETSC_SUCCESS);
662:   ((PetscObject)*dm)->refct = 0;

664:   PetscCall(DMClearGlobalVectors(*dm));
665:   PetscCall(DMClearLocalVectors(*dm));
666:   PetscCall(DMClearNamedGlobalVectors(*dm));
667:   PetscCall(DMClearNamedLocalVectors(*dm));

669:   /* Destroy the list of hooks */
670:   {
671:     DMCoarsenHookLink link, next;
672:     for (link = (*dm)->coarsenhook; link; link = next) {
673:       next = link->next;
674:       PetscCall(PetscFree(link));
675:     }
676:     (*dm)->coarsenhook = NULL;
677:   }
678:   {
679:     DMRefineHookLink link, next;
680:     for (link = (*dm)->refinehook; link; link = next) {
681:       next = link->next;
682:       PetscCall(PetscFree(link));
683:     }
684:     (*dm)->refinehook = NULL;
685:   }
686:   {
687:     DMSubDomainHookLink link, next;
688:     for (link = (*dm)->subdomainhook; link; link = next) {
689:       next = link->next;
690:       PetscCall(PetscFree(link));
691:     }
692:     (*dm)->subdomainhook = NULL;
693:   }
694:   {
695:     DMGlobalToLocalHookLink link, next;
696:     for (link = (*dm)->gtolhook; link; link = next) {
697:       next = link->next;
698:       PetscCall(PetscFree(link));
699:     }
700:     (*dm)->gtolhook = NULL;
701:   }
702:   {
703:     DMLocalToGlobalHookLink link, next;
704:     for (link = (*dm)->ltoghook; link; link = next) {
705:       next = link->next;
706:       PetscCall(PetscFree(link));
707:     }
708:     (*dm)->ltoghook = NULL;
709:   }
710:   /* Destroy the work arrays */
711:   {
712:     DMWorkLink link, next;
713:     PetscCheck(!(*dm)->workout, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Work array still checked out %p %p", (void *)(*dm)->workout, (*dm)->workout->mem);
714:     for (link = (*dm)->workin; link; link = next) {
715:       next = link->next;
716:       PetscCall(PetscFree(link->mem));
717:       PetscCall(PetscFree(link));
718:     }
719:     (*dm)->workin = NULL;
720:   }
721:   /* destroy the labels */
722:   PetscCall(DMDestroyLabelLinkList_Internal(*dm));
723:   /* destroy the fields */
724:   PetscCall(DMClearFields(*dm));
725:   /* destroy the boundaries */
726:   {
727:     DMBoundary next = (*dm)->boundary;
728:     while (next) {
729:       DMBoundary b = next;

731:       next = b->next;
732:       PetscCall(PetscFree(b));
733:     }
734:   }

736:   PetscCall(PetscObjectDestroy(&(*dm)->dmksp));
737:   PetscCall(PetscObjectDestroy(&(*dm)->dmsnes));
738:   PetscCall(PetscObjectDestroy(&(*dm)->dmts));

740:   if ((*dm)->ctx && (*dm)->ctxdestroy) PetscCall((*(*dm)->ctxdestroy)(&(*dm)->ctx));
741:   PetscCall(MatFDColoringDestroy(&(*dm)->fd));
742:   PetscCall(ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap));
743:   PetscCall(PetscFree((*dm)->vectype));
744:   PetscCall(PetscFree((*dm)->mattype));

746:   PetscCall(PetscSectionDestroy(&(*dm)->localSection));
747:   PetscCall(PetscSectionDestroy(&(*dm)->globalSection));
748:   PetscCall(PetscFree((*dm)->reorderSectionType));
749:   PetscCall(PetscLayoutDestroy(&(*dm)->map));
750:   PetscCall(PetscSectionDestroy(&(*dm)->defaultConstraint.section));
751:   PetscCall(MatDestroy(&(*dm)->defaultConstraint.mat));
752:   PetscCall(PetscSFDestroy(&(*dm)->sf));
753:   PetscCall(PetscSFDestroy(&(*dm)->sectionSF));
754:   if ((*dm)->sfNatural) PetscCall(PetscSFDestroy(&(*dm)->sfNatural));
755:   PetscCall(PetscObjectDereference((PetscObject)(*dm)->sfMigration));
756:   PetscCall(DMClearAuxiliaryVec(*dm));
757:   PetscCall(PetscHMapAuxDestroy(&(*dm)->auxData));
758:   if ((*dm)->coarseMesh && (*dm)->coarseMesh->fineMesh == *dm) PetscCall(DMSetFineDM((*dm)->coarseMesh, NULL));

760:   PetscCall(DMDestroy(&(*dm)->coarseMesh));
761:   if ((*dm)->fineMesh && (*dm)->fineMesh->coarseMesh == *dm) PetscCall(DMSetCoarseDM((*dm)->fineMesh, NULL));
762:   PetscCall(DMDestroy(&(*dm)->fineMesh));
763:   PetscCall(PetscFree((*dm)->Lstart));
764:   PetscCall(PetscFree((*dm)->L));
765:   PetscCall(PetscFree((*dm)->maxCell));
766:   PetscCall(DMDestroyCoordinates_Private(&(*dm)->coordinates[0]));
767:   PetscCall(DMDestroyCoordinates_Private(&(*dm)->coordinates[1]));
768:   if ((*dm)->transformDestroy) PetscCall((*(*dm)->transformDestroy)(*dm, (*dm)->transformCtx));
769:   PetscCall(DMDestroy(&(*dm)->transformDM));
770:   PetscCall(VecDestroy(&(*dm)->transform));
771:   for (PetscInt i = 0; i < (*dm)->periodic.num_affines; i++) {
772:     PetscCall(VecScatterDestroy(&(*dm)->periodic.affine_to_local[i]));
773:     PetscCall(VecDestroy(&(*dm)->periodic.affine[i]));
774:   }
775:   if ((*dm)->periodic.num_affines > 0) PetscCall(PetscFree2((*dm)->periodic.affine_to_local, (*dm)->periodic.affine));

777:   PetscCall(DMClearDS(*dm));
778:   PetscCall(DMDestroy(&(*dm)->dmBC));
779:   /* if memory was published with SAWs then destroy it */
780:   PetscCall(PetscObjectSAWsViewOff((PetscObject)*dm));

782:   PetscTryTypeMethod(*dm, destroy);
783:   PetscCall(DMMonitorCancel(*dm));
784:   PetscCall(DMCeedDestroy(&(*dm)->dmceed));
785: #ifdef PETSC_HAVE_LIBCEED
786:   PetscCallCEED(CeedElemRestrictionDestroy(&(*dm)->ceedERestrict));
787:   PetscCallCEED(CeedDestroy(&(*dm)->ceed));
788: #endif
789:   /* We do not destroy (*dm)->data here so that we can reference count backend objects */
790:   PetscCall(PetscHeaderDestroy(dm));
791:   PetscFunctionReturn(PETSC_SUCCESS);
792: }

794: /*@
795:   DMSetUp - sets up the data structures inside a `DM` object

797:   Collective

799:   Input Parameter:
800: . dm - the `DM` object to setup

802:   Level: intermediate

804:   Note:
805:   This is usually called after various parameter setting operations and `DMSetFromOptions()` are called on the `DM`

807: .seealso: [](ch_dmbase), `DM`, `DMCreate()`, `DMSetType()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`
808: @*/
809: PetscErrorCode DMSetUp(DM dm)
810: {
811:   PetscFunctionBegin;
813:   if (dm->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
814:   PetscTryTypeMethod(dm, setup);
815:   dm->setupcalled = PETSC_TRUE;
816:   PetscFunctionReturn(PETSC_SUCCESS);
817: }

819: /*@
820:   DMSetFromOptions - sets parameters in a `DM` from the options database

822:   Collective

824:   Input Parameter:
825: . dm - the `DM` object to set options for

827:   Options Database Keys:
828: + -dm_preallocate_only                               - Only preallocate the matrix for `DMCreateMatrix()` and `DMCreateMassMatrix()`, but do not fill it with zeros
829: . -dm_vec_type <type>                                - type of vector to create inside `DM`
830: . -dm_mat_type <type>                                - type of matrix to create inside `DM`
831: . -dm_is_coloring_type                               - <global or local>
832: . -dm_bind_below <n>                                 - bind (force execution on CPU) for `Vec` and `Mat` objects with local size (number of vector entries or matrix rows) below n; currently only supported for `DMDA`
833: . -dm_plex_option_phases <ph0_, ph1_, ...>           - List of prefixes for option processing phases
834: . -dm_plex_filename <str>                            - File containing a mesh
835: . -dm_plex_boundary_filename <str>                   - File containing a mesh boundary
836: . -dm_plex_name <str>                                - Name of the mesh in the file
837: . -dm_plex_shape <shape>                             - The domain shape, such as `BOX`, `SPHERE`, etc.
838: . -dm_plex_cell <ct>                                 - Cell shape
839: . -dm_plex_reference_cell_domain <bool>              - Use a reference cell domain
840: . -dm_plex_dim <dim>                                 - Set the topological dimension
841: . -dm_plex_simplex <bool>                            - `PETSC_TRUE` for simplex elements, `PETSC_FALSE` for tensor elements
842: . -dm_plex_interpolate <bool>                        - `PETSC_TRUE` turns on topological interpolation (creating edges and faces)
843: . -dm_plex_orient <bool>                             - `PETSC_TRUE` turns on topological orientation (flipping edges and faces)
844: . -dm_plex_scale <sc>                                - Scale factor for mesh coordinates
845: . -dm_coord_remap <bool>                             - Map coordinates using a function
846: . -dm_coord_map <mapname>                            - Select a builtin coordinate map
847: . -dm_coord_map_params <p0,p1,p2,...>                - Set coordinate mapping parameters
848: . -dm_plex_box_faces <m,n,p>                         - Number of faces along each dimension
849: . -dm_plex_box_lower <x,y,z>                         - Specify lower-left-bottom coordinates for the box
850: . -dm_plex_box_upper <x,y,z>                         - Specify upper-right-top coordinates for the box
851: . -dm_plex_box_bd <bx,by,bz>                         - Specify the `DMBoundaryType` for each direction
852: . -dm_plex_sphere_radius <r>                         - The sphere radius
853: . -dm_plex_ball_radius <r>                           - Radius of the ball
854: . -dm_plex_cylinder_bd <bz>                          - Boundary type in the z direction
855: . -dm_plex_cylinder_num_wedges <n>                   - Number of wedges around the cylinder
856: . -dm_plex_reorder <order>                           - Reorder the mesh using the specified algorithm
857: . -dm_refine_pre <n>                                 - The number of refinements before distribution
858: . -dm_refine_uniform_pre <bool>                      - Flag for uniform refinement before distribution
859: . -dm_refine_volume_limit_pre <v>                    - The maximum cell volume after refinement before distribution
860: . -dm_refine <n>                                     - The number of refinements after distribution
861: . -dm_extrude <l>                                    - Activate extrusion and specify the number of layers to extrude
862: . -dm_plex_transform_extrude_thickness <t>           - The total thickness of extruded layers
863: . -dm_plex_transform_extrude_use_tensor <bool>       - Use tensor cells when extruding
864: . -dm_plex_transform_extrude_symmetric <bool>        - Extrude layers symmetrically about the surface
865: . -dm_plex_transform_extrude_normal <n0,...,nd>      - Specify the extrusion direction
866: . -dm_plex_transform_extrude_thicknesses <t0,...,tl> - Specify thickness of each layer
867: . -dm_plex_create_fv_ghost_cells                     - Flag to create finite volume ghost cells on the boundary
868: . -dm_plex_fv_ghost_cells_label <name>               - Label name for ghost cells boundary
869: . -dm_distribute <bool>                              - Flag to redistribute a mesh among processes
870: . -dm_distribute_overlap <n>                         - The size of the overlap halo
871: . -dm_plex_adj_cone <bool>                           - Set adjacency direction
872: . -dm_plex_adj_closure <bool>                        - Set adjacency size
873: . -dm_plex_use_ceed <bool>                           - Use LibCEED as the FEM backend
874: . -dm_plex_check_symmetry                            - Check that the adjacency information in the mesh is symmetric - `DMPlexCheckSymmetry()`
875: . -dm_plex_check_skeleton                            - Check that each cell has the correct number of vertices (only for homogeneous simplex or tensor meshes) - `DMPlexCheckSkeleton()`
876: . -dm_plex_check_faces                               - Check that the faces of each cell give a vertex order this is consistent with what we expect from the cell type - `DMPlexCheckFaces()`
877: . -dm_plex_check_geometry                            - Check that cells have positive volume - `DMPlexCheckGeometry()`
878: . -dm_plex_check_pointsf                             - Check some necessary conditions for `PointSF` - `DMPlexCheckPointSF()`
879: . -dm_plex_check_interface_cones                     - Check points on inter-partition interfaces have conforming order of cone points - `DMPlexCheckInterfaceCones()`
880: - -dm_plex_check_all                                 - Perform all the checks above

882:   Level: intermediate

884:   Note:
885:   For some `DMType` such as `DMDA` this cannot be called after `DMSetUp()` has been called.

887: .seealso: [](ch_dmbase), `DM`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`,
888:          `DMPlexCheckSymmetry()`, `DMPlexCheckSkeleton()`, `DMPlexCheckFaces()`, `DMPlexCheckGeometry()`, `DMPlexCheckPointSF()`, `DMPlexCheckInterfaceCones()`,
889:          `DMSetOptionsPrefix()`, `DMType`, `DMPLEX`, `DMDA`, `DMSetUp()`
890: @*/
891: PetscErrorCode DMSetFromOptions(DM dm)
892: {
893:   char      typeName[256];
894:   PetscBool flg;

896:   PetscFunctionBegin;
898:   dm->setfromoptionscalled = PETSC_TRUE;
899:   if (dm->sf) PetscCall(PetscSFSetFromOptions(dm->sf));
900:   if (dm->sectionSF) PetscCall(PetscSFSetFromOptions(dm->sectionSF));
901:   if (dm->coordinates[0].dm) PetscCall(DMSetFromOptions(dm->coordinates[0].dm));
902:   PetscObjectOptionsBegin((PetscObject)dm);
903:   PetscCall(PetscOptionsBool("-dm_preallocate_only", "only preallocate matrix, but do not set column indices", "DMSetMatrixPreallocateOnly", dm->prealloc_only, &dm->prealloc_only, NULL));
904:   PetscCall(PetscOptionsFList("-dm_vec_type", "Vector type used for created vectors", "DMSetVecType", VecList, dm->vectype, typeName, 256, &flg));
905:   if (flg) PetscCall(DMSetVecType(dm, typeName));
906:   PetscCall(PetscOptionsFList("-dm_mat_type", "Matrix type used for created matrices", "DMSetMatType", MatList, dm->mattype ? dm->mattype : typeName, typeName, sizeof(typeName), &flg));
907:   if (flg) PetscCall(DMSetMatType(dm, typeName));
908:   PetscCall(PetscOptionsEnum("-dm_blocking_type", "Topological point or field node blocking", "DMSetBlockingType", DMBlockingTypes, (PetscEnum)dm->blocking_type, (PetscEnum *)&dm->blocking_type, NULL));
909:   PetscCall(PetscOptionsEnum("-dm_is_coloring_type", "Global or local coloring of Jacobian", "DMSetISColoringType", ISColoringTypes, (PetscEnum)dm->coloringtype, (PetscEnum *)&dm->coloringtype, NULL));
910:   PetscCall(PetscOptionsInt("-dm_bind_below", "Set the size threshold (in entries) below which the Vec is bound to the CPU", "VecBindToCPU", dm->bind_below, &dm->bind_below, &flg));
911:   PetscCall(PetscOptionsBool("-dm_ignore_perm_output", "Ignore the local section permutation on output", "DMGetOutputDM", dm->ignorePermOutput, &dm->ignorePermOutput, NULL));
912:   PetscTryTypeMethod(dm, setfromoptions, PetscOptionsObject);
913:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
914:   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)dm, PetscOptionsObject));
915:   PetscOptionsEnd();
916:   PetscFunctionReturn(PETSC_SUCCESS);
917: }

919: /*@
920:   DMViewFromOptions - View a `DM` in a particular way based on a request in the options database

922:   Collective

924:   Input Parameters:
925: + dm   - the `DM` object
926: . obj  - optional object that provides the prefix for the options database (if `NULL` then the prefix in obj is used)
927: - name - option string that is used to activate viewing

929:   Level: intermediate

931:   Note:
932:   See `PetscObjectViewFromOptions()` for a list of values that can be provided in the options database to determine how the `DM` is viewed

934: .seealso: [](ch_dmbase), `DM`, `DMView()`, `PetscObjectViewFromOptions()`, `DMCreate()`
935: @*/
936: PetscErrorCode DMViewFromOptions(DM dm, PetscObject obj, const char name[])
937: {
938:   PetscFunctionBegin;
940:   PetscCall(PetscObjectViewFromOptions((PetscObject)dm, obj, name));
941:   PetscFunctionReturn(PETSC_SUCCESS);
942: }

944: /*@
945:   DMView - Views a `DM`. Depending on the `PetscViewer` and its `PetscViewerFormat` it may print some ASCII information about the `DM` to the screen or a file or
946:   save the `DM` in a binary file to be loaded later or create a visualization of the `DM`

948:   Collective

950:   Input Parameters:
951: + dm - the `DM` object to view
952: - v  - the viewer

954:   Level: beginner

956:   Notes:

958:   `PetscViewer` = `PETSCVIEWERHDF5` i.e. HDF5 format can be used with `PETSC_VIEWER_HDF5_PETSC` as the `PetscViewerFormat` to save multiple `DMPLEX`
959:   meshes in a single HDF5 file. This in turn requires one to name the `DMPLEX` object with `PetscObjectSetName()`
960:   before saving it with `DMView()` and before loading it with `DMLoad()` for identification of the mesh object.

962:   `PetscViewer` = `PETSCVIEWEREXODUSII` i.e. ExodusII format assumes that element blocks (mapped to "Cell sets" labels)
963:   consists of sequentially numbered cells.

965:   If `dm` has been distributed, only the part of the `DM` on MPI rank 0 (including "ghost" cells and vertices) will be written.

967:   Only TRI, TET, QUAD, and HEX cells are supported.

969:   `DMPLEX` only represents geometry while most post-processing software expect that a mesh also provides information on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2.
970:   The order of the mesh shall be set using `PetscViewerExodusIISetOrder()`

972:   Variable names can be set and queried using `PetscViewerExodusII[Set/Get][Nodal/Zonal]VariableNames[s]`.

974: .seealso: [](ch_dmbase), `DM`, `PetscViewer`, `PetscViewerFormat`, `PetscViewerSetFormat()`, `DMDestroy()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMLoad()`, `PetscObjectSetName()`
975: @*/
976: PetscErrorCode DMView(DM dm, PetscViewer v)
977: {
978:   PetscBool         isbinary;
979:   PetscMPIInt       size;
980:   PetscViewerFormat format;

982:   PetscFunctionBegin;
984:   if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)dm), &v));
986:   /* Ideally, we would like to have this test on.
987:      However, it currently breaks socket viz via GLVis.
988:      During DMView(parallel_mesh,glvis_viewer), each
989:      process opens a sequential ASCII socket to visualize
990:      the local mesh, and PetscObjectView(dm,local_socket)
991:      is internally called inside VecView_GLVis, incurring
992:      in an error here */
993:   /* PetscCheckSameComm(dm,1,v,2); */
994:   PetscCall(PetscViewerCheckWritable(v));

996:   PetscCall(PetscLogEventBegin(DM_View, v, 0, 0, 0));
997:   PetscCall(PetscViewerGetFormat(v, &format));
998:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
999:   if (size == 1 && format == PETSC_VIEWER_LOAD_BALANCE) PetscFunctionReturn(PETSC_SUCCESS);
1000:   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)dm, v));
1001:   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERBINARY, &isbinary));
1002:   if (isbinary) {
1003:     PetscInt classid = DM_FILE_CLASSID;
1004:     char     type[256];

1006:     PetscCall(PetscViewerBinaryWrite(v, &classid, 1, PETSC_INT));
1007:     PetscCall(PetscStrncpy(type, ((PetscObject)dm)->type_name, sizeof(type)));
1008:     PetscCall(PetscViewerBinaryWrite(v, type, 256, PETSC_CHAR));
1009:   }
1010:   PetscTryTypeMethod(dm, view, v);
1011:   PetscCall(PetscLogEventEnd(DM_View, v, 0, 0, 0));
1012:   PetscFunctionReturn(PETSC_SUCCESS);
1013: }

1015: /*@
1016:   DMCreateGlobalVector - Creates a global vector from a `DM` object. A global vector is a parallel vector that has no duplicate values shared between MPI ranks,
1017:   that is it has no ghost locations.

1019:   Collective

1021:   Input Parameter:
1022: . dm - the `DM` object

1024:   Output Parameter:
1025: . vec - the global vector

1027:   Level: beginner

1029: .seealso: [](ch_dmbase), `DM`, `Vec`, `DMCreateLocalVector()`, `DMGetGlobalVector()`, `DMDestroy()`, `DMView()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`,
1030:          `DMGlobalToLocalBegin()`, `DMGlobalToLocalEnd()`
1031: @*/
1032: PetscErrorCode DMCreateGlobalVector(DM dm, Vec *vec)
1033: {
1034:   PetscFunctionBegin;
1036:   PetscAssertPointer(vec, 2);
1037:   PetscUseTypeMethod(dm, createglobalvector, vec);
1038:   if (PetscDefined(USE_DEBUG)) {
1039:     DM vdm;

1041:     PetscCall(VecGetDM(*vec, &vdm));
1042:     PetscCheck(vdm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DM type '%s' did not attach the DM to the vector", ((PetscObject)dm)->type_name);
1043:   }
1044:   PetscFunctionReturn(PETSC_SUCCESS);
1045: }

1047: /*@
1048:   DMCreateLocalVector - Creates a local vector from a `DM` object.

1050:   Not Collective

1052:   Input Parameter:
1053: . dm - the `DM` object

1055:   Output Parameter:
1056: . vec - the local vector

1058:   Level: beginner

1060:   Note:
1061:   A local vector usually has ghost locations that contain values that are owned by different MPI ranks. A global vector has no ghost locations.

1063: .seealso: [](ch_dmbase), `DM`, `Vec`, `DMCreateGlobalVector()`, `DMGetLocalVector()`, `DMDestroy()`, `DMView()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`
1064:          `DMGlobalToLocalBegin()`, `DMGlobalToLocalEnd()`
1065: @*/
1066: PetscErrorCode DMCreateLocalVector(DM dm, Vec *vec)
1067: {
1068:   PetscFunctionBegin;
1070:   PetscAssertPointer(vec, 2);
1071:   PetscUseTypeMethod(dm, createlocalvector, vec);
1072:   if (PetscDefined(USE_DEBUG)) {
1073:     DM vdm;

1075:     PetscCall(VecGetDM(*vec, &vdm));
1076:     PetscCheck(vdm, PETSC_COMM_SELF, PETSC_ERR_LIB, "DM type '%s' did not attach the DM to the vector", ((PetscObject)dm)->type_name);
1077:   }
1078:   PetscFunctionReturn(PETSC_SUCCESS);
1079: }

1081: /*@
1082:   DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a `DM`.

1084:   Collective

1086:   Input Parameter:
1087: . dm - the `DM` that provides the mapping

1089:   Output Parameter:
1090: . ltog - the mapping

1092:   Level: advanced

1094:   Notes:
1095:   The global to local mapping allows one to set values into the global vector or matrix using `VecSetValuesLocal()` and `MatSetValuesLocal()`

1097:   Vectors obtained with  `DMCreateGlobalVector()` and matrices obtained with `DMCreateMatrix()` already contain the global mapping so you do
1098:   need to use this function with those objects.

1100:   This mapping can then be used by `VecSetLocalToGlobalMapping()` or `MatSetLocalToGlobalMapping()`.

1102: .seealso: [](ch_dmbase), `DM`, `DMCreateLocalVector()`, `DMCreateGlobalVector()`, `VecSetLocalToGlobalMapping()`, `MatSetLocalToGlobalMapping()`,
1103:           `DMCreateMatrix()`
1104: @*/
1105: PetscErrorCode DMGetLocalToGlobalMapping(DM dm, ISLocalToGlobalMapping *ltog)
1106: {
1107:   PetscInt bs = -1, bsLocal[2], bsMinMax[2];

1109:   PetscFunctionBegin;
1111:   PetscAssertPointer(ltog, 2);
1112:   if (!dm->ltogmap) {
1113:     PetscSection section, sectionGlobal;

1115:     PetscCall(DMGetLocalSection(dm, &section));
1116:     if (section) {
1117:       const PetscInt *cdofs;
1118:       PetscInt       *ltog;
1119:       PetscInt        pStart, pEnd, n, p, k, l;

1121:       PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
1122:       PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
1123:       PetscCall(PetscSectionGetStorageSize(section, &n));
1124:       PetscCall(PetscMalloc1(n, &ltog)); /* We want the local+overlap size */
1125:       for (p = pStart, l = 0; p < pEnd; ++p) {
1126:         PetscInt bdof, cdof, dof, off, c, cind;

1128:         /* Should probably use constrained dofs */
1129:         PetscCall(PetscSectionGetDof(section, p, &dof));
1130:         PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
1131:         PetscCall(PetscSectionGetConstraintIndices(section, p, &cdofs));
1132:         PetscCall(PetscSectionGetOffset(sectionGlobal, p, &off));
1133:         /* If you have dofs, and constraints, and they are unequal, we set the blocksize to 1 */
1134:         bdof = cdof && (dof - cdof) ? 1 : dof;
1135:         if (dof) bs = bs < 0 ? bdof : PetscGCD(bs, bdof);

1137:         for (c = 0, cind = 0; c < dof; ++c, ++l) {
1138:           if (cind < cdof && c == cdofs[cind]) {
1139:             ltog[l] = off < 0 ? off - c : -(off + c + 1);
1140:             cind++;
1141:           } else {
1142:             ltog[l] = (off < 0 ? -(off + 1) : off) + c - cind;
1143:           }
1144:         }
1145:       }
1146:       /* Must have same blocksize on all procs (some might have no points) */
1147:       bsLocal[0] = bs < 0 ? PETSC_INT_MAX : bs;
1148:       bsLocal[1] = bs;
1149:       PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm), bsLocal, bsMinMax));
1150:       if (bsMinMax[0] != bsMinMax[1]) {
1151:         bs = 1;
1152:       } else {
1153:         bs = bsMinMax[0];
1154:       }
1155:       bs = bs < 0 ? 1 : bs;
1156:       /* Must reduce indices by blocksize */
1157:       if (bs > 1) {
1158:         for (l = 0, k = 0; l < n; l += bs, ++k) {
1159:           // Integer division of negative values truncates toward zero(!), not toward negative infinity
1160:           ltog[k] = ltog[l] >= 0 ? ltog[l] / bs : -(-(ltog[l] + 1) / bs + 1);
1161:         }
1162:         n /= bs;
1163:       }
1164:       PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), bs, n, ltog, PETSC_OWN_POINTER, &dm->ltogmap));
1165:     } else PetscUseTypeMethod(dm, getlocaltoglobalmapping);
1166:   }
1167:   *ltog = dm->ltogmap;
1168:   PetscFunctionReturn(PETSC_SUCCESS);
1169: }

1171: /*@
1172:   DMGetBlockSize - Gets the inherent block size associated with a `DM`

1174:   Not Collective

1176:   Input Parameter:
1177: . dm - the `DM` with block structure

1179:   Output Parameter:
1180: . bs - the block size, 1 implies no exploitable block structure

1182:   Level: intermediate

1184:   Notes:
1185:   This might be the number of degrees of freedom at each grid point for a structured grid.

1187:   Complex `DM` that represent multiphysics or staggered grids or mixed-methods do not generally have a single inherent block size, but
1188:   rather different locations in the vectors may have a different block size.

1190: .seealso: [](ch_dmbase), `DM`, `ISCreateBlock()`, `VecSetBlockSize()`, `MatSetBlockSize()`, `DMGetLocalToGlobalMapping()`
1191: @*/
1192: PetscErrorCode DMGetBlockSize(DM dm, PetscInt *bs)
1193: {
1194:   PetscFunctionBegin;
1196:   PetscAssertPointer(bs, 2);
1197:   PetscCheck(dm->bs >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DM does not have enough information to provide a block size yet");
1198:   *bs = dm->bs;
1199:   PetscFunctionReturn(PETSC_SUCCESS);
1200: }

1202: /*@
1203:   DMCreateInterpolation - Gets the interpolation matrix between two `DM` objects. The resulting matrix map degrees of freedom in the vector obtained by
1204:   `DMCreateGlobalVector()` on the coarse `DM` to similar vectors on the fine grid `DM`.

1206:   Collective

1208:   Input Parameters:
1209: + dmc - the `DM` object
1210: - dmf - the second, finer `DM` object

1212:   Output Parameters:
1213: + mat - the interpolation
1214: - vec - the scaling (optional, pass `NULL` if not needed), see `DMCreateInterpolationScale()`

1216:   Level: developer

1218:   Notes:
1219:   For `DMDA` objects this only works for "uniform refinement", that is the refined mesh was obtained `DMRefine()` or the coarse mesh was obtained by
1220:   DMCoarsen(). The coordinates set into the `DMDA` are completely ignored in computing the interpolation.

1222:   For `DMDA` objects you can use this interpolation (more precisely the interpolation from the `DMGetCoordinateDM()`) to interpolate the mesh coordinate
1223:   vectors EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic.

1225: .seealso: [](ch_dmbase), `DM`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMRefine()`, `DMCoarsen()`, `DMCreateRestriction()`, `DMCreateInterpolationScale()`
1226: @*/
1227: PetscErrorCode DMCreateInterpolation(DM dmc, DM dmf, Mat *mat, Vec *vec)
1228: {
1229:   PetscFunctionBegin;
1232:   PetscAssertPointer(mat, 3);
1233:   PetscCall(PetscLogEventBegin(DM_CreateInterpolation, dmc, dmf, 0, 0));
1234:   PetscUseTypeMethod(dmc, createinterpolation, dmf, mat, vec);
1235:   PetscCall(PetscLogEventEnd(DM_CreateInterpolation, dmc, dmf, 0, 0));
1236:   PetscFunctionReturn(PETSC_SUCCESS);
1237: }

1239: /*@
1240:   DMCreateInterpolationScale - Forms L = 1/(R*1) where 1 is the vector of all ones, and R is
1241:   the transpose of the interpolation between the `DM`.

1243:   Input Parameters:
1244: + dac - `DM` that defines a coarse mesh
1245: . daf - `DM` that defines a fine mesh
1246: - mat - the restriction (or interpolation operator) from fine to coarse

1248:   Output Parameter:
1249: . scale - the scaled vector

1251:   Level: advanced

1253:   Note:
1254:   xcoarse = diag(L)*R*xfine preserves scale and is thus suitable for state (versus residual)
1255:   restriction. In other words xcoarse is the coarse representation of xfine.

1257:   Developer Note:
1258:   If the fine-scale `DMDA` has the -dm_bind_below option set to true, then `DMCreateInterpolationScale()` calls `MatSetBindingPropagates()`
1259:   on the restriction/interpolation operator to set the bindingpropagates flag to true.

1261: .seealso: [](ch_dmbase), `DM`, `MatRestrict()`, `MatInterpolate()`, `DMCreateInterpolation()`, `DMCreateRestriction()`, `DMCreateGlobalVector()`
1262: @*/
1263: PetscErrorCode DMCreateInterpolationScale(DM dac, DM daf, Mat mat, Vec *scale)
1264: {
1265:   Vec         fine;
1266:   PetscScalar one = 1.0;
1267: #if defined(PETSC_HAVE_CUDA)
1268:   PetscBool bindingpropagates, isbound;
1269: #endif

1271:   PetscFunctionBegin;
1272:   PetscCall(DMCreateGlobalVector(daf, &fine));
1273:   PetscCall(DMCreateGlobalVector(dac, scale));
1274:   PetscCall(VecSet(fine, one));
1275: #if defined(PETSC_HAVE_CUDA)
1276:   /* If the 'fine' Vec is bound to the CPU, it makes sense to bind 'mat' as well.
1277:    * Note that we only do this for the CUDA case, right now, but if we add support for MatMultTranspose() via ViennaCL,
1278:    * we'll need to do it for that case, too.*/
1279:   PetscCall(VecGetBindingPropagates(fine, &bindingpropagates));
1280:   if (bindingpropagates) {
1281:     PetscCall(MatSetBindingPropagates(mat, PETSC_TRUE));
1282:     PetscCall(VecBoundToCPU(fine, &isbound));
1283:     PetscCall(MatBindToCPU(mat, isbound));
1284:   }
1285: #endif
1286:   PetscCall(MatRestrict(mat, fine, *scale));
1287:   PetscCall(VecDestroy(&fine));
1288:   PetscCall(VecReciprocal(*scale));
1289:   PetscFunctionReturn(PETSC_SUCCESS);
1290: }

1292: /*@
1293:   DMCreateRestriction - Gets restriction matrix between two `DM` objects. The resulting matrix map degrees of freedom in the vector obtained by
1294:   `DMCreateGlobalVector()` on the fine `DM` to similar vectors on the coarse grid `DM`.

1296:   Collective

1298:   Input Parameters:
1299: + dmc - the `DM` object
1300: - dmf - the second, finer `DM` object

1302:   Output Parameter:
1303: . mat - the restriction

1305:   Level: developer

1307:   Note:
1308:   This only works for `DMSTAG`. For many situations either the transpose of the operator obtained with `DMCreateInterpolation()` or that
1309:   matrix multiplied by the vector obtained with `DMCreateInterpolationScale()` provides the desired object.

1311: .seealso: [](ch_dmbase), `DM`, `DMRestrict()`, `DMInterpolate()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMRefine()`, `DMCoarsen()`, `DMCreateInterpolation()`
1312: @*/
1313: PetscErrorCode DMCreateRestriction(DM dmc, DM dmf, Mat *mat)
1314: {
1315:   PetscFunctionBegin;
1318:   PetscAssertPointer(mat, 3);
1319:   PetscCall(PetscLogEventBegin(DM_CreateRestriction, dmc, dmf, 0, 0));
1320:   PetscUseTypeMethod(dmc, createrestriction, dmf, mat);
1321:   PetscCall(PetscLogEventEnd(DM_CreateRestriction, dmc, dmf, 0, 0));
1322:   PetscFunctionReturn(PETSC_SUCCESS);
1323: }

1325: /*@
1326:   DMCreateInjection - Gets injection matrix between two `DM` objects.

1328:   Collective

1330:   Input Parameters:
1331: + dac - the `DM` object
1332: - daf - the second, finer `DM` object

1334:   Output Parameter:
1335: . mat - the injection

1337:   Level: developer

1339:   Notes:
1340:   This is an operator that applied to a vector obtained with `DMCreateGlobalVector()` on the
1341:   fine grid maps the values to a vector on the vector on the coarse `DM` by simply selecting
1342:   the values on the coarse grid points. This compares to the operator obtained by
1343:   `DMCreateRestriction()` or the transpose of the operator obtained by
1344:   `DMCreateInterpolation()` that uses a "local weighted average" of the values around the
1345:   coarse grid point as the coarse grid value.

1347:   For `DMDA` objects this only works for "uniform refinement", that is the refined mesh was obtained `DMRefine()` or the coarse mesh was obtained by
1348:   `DMCoarsen()`. The coordinates set into the `DMDA` are completely ignored in computing the injection.

1350: .seealso: [](ch_dmbase), `DM`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMCreateInterpolation()`,
1351:           `DMCreateRestriction()`, `MatRestrict()`, `MatInterpolate()`
1352: @*/
1353: PetscErrorCode DMCreateInjection(DM dac, DM daf, Mat *mat)
1354: {
1355:   PetscFunctionBegin;
1358:   PetscAssertPointer(mat, 3);
1359:   PetscCall(PetscLogEventBegin(DM_CreateInjection, dac, daf, 0, 0));
1360:   PetscUseTypeMethod(dac, createinjection, daf, mat);
1361:   PetscCall(PetscLogEventEnd(DM_CreateInjection, dac, daf, 0, 0));
1362:   PetscFunctionReturn(PETSC_SUCCESS);
1363: }

1365: /*@
1366:   DMCreateMassMatrix - Gets the mass matrix between two `DM` objects, M_ij = \int \phi_i \psi_j where the \phi are Galerkin basis functions for a
1367:   a Galerkin finite element model on the `DM`

1369:   Collective

1371:   Input Parameters:
1372: + dmc - the target `DM` object
1373: - dmf - the source `DM` object, can be `NULL`

1375:   Output Parameter:
1376: . mat - the mass matrix

1378:   Level: developer

1380:   Notes:
1381:   For `DMPLEX` the finite element model for the `DM` must have been already provided.

1383:   if `dmc` is `dmf` or `NULL`, then x^t M x is an approximation to the L2 norm of the vector x which is obtained by `DMCreateGlobalVector()`

1385: .seealso: [](ch_dmbase), `DM`, `DMCreateMassMatrixLumped()`, `DMCreateMatrix()`, `DMRefine()`, `DMCoarsen()`, `DMCreateRestriction()`, `DMCreateInterpolation()`, `DMCreateInjection()`
1386: @*/
1387: PetscErrorCode DMCreateMassMatrix(DM dmc, DM dmf, Mat *mat)
1388: {
1389:   PetscFunctionBegin;
1391:   if (!dmf) dmf = dmc;
1393:   PetscAssertPointer(mat, 3);
1394:   PetscCall(PetscLogEventBegin(DM_CreateMassMatrix, dmc, dmf, 0, 0));
1395:   PetscUseTypeMethod(dmc, createmassmatrix, dmf, mat);
1396:   PetscCall(PetscLogEventEnd(DM_CreateMassMatrix, dmc, dmf, 0, 0));
1397:   PetscFunctionReturn(PETSC_SUCCESS);
1398: }

1400: /*@
1401:   DMCreateMassMatrixLumped - Gets the lumped mass matrix for a given `DM`

1403:   Collective

1405:   Input Parameter:
1406: . dm - the `DM` object

1408:   Output Parameters:
1409: + llm - the local lumped mass matrix, which is a diagonal matrix, represented as a vector
1410: - lm  - the global lumped mass matrix, which is a diagonal matrix, represented as a vector

1412:   Level: developer

1414:   Note:
1415:   See `DMCreateMassMatrix()` for how to create the non-lumped version of the mass matrix.

1417: .seealso: [](ch_dmbase), `DM`, `DMCreateMassMatrix()`, `DMCreateMatrix()`, `DMRefine()`, `DMCoarsen()`, `DMCreateRestriction()`, `DMCreateInterpolation()`, `DMCreateInjection()`
1418: @*/
1419: PetscErrorCode DMCreateMassMatrixLumped(DM dm, Vec *llm, Vec *lm)
1420: {
1421:   PetscFunctionBegin;
1423:   if (llm) PetscAssertPointer(llm, 2);
1424:   if (lm) PetscAssertPointer(lm, 3);
1425:   if (llm || lm) PetscUseTypeMethod(dm, createmassmatrixlumped, llm, lm);
1426:   PetscFunctionReturn(PETSC_SUCCESS);
1427: }

1429: /*@
1430:   DMCreateColoring - Gets coloring of a graph associated with the `DM`. Often the graph represents the operator matrix associated with the discretization
1431:   of a PDE on the `DM`.

1433:   Collective

1435:   Input Parameters:
1436: + dm    - the `DM` object
1437: - ctype - `IS_COLORING_LOCAL` or `IS_COLORING_GLOBAL`

1439:   Output Parameter:
1440: . coloring - the coloring

1442:   Level: developer

1444:   Notes:
1445:   Coloring of matrices can also be computed directly from the sparse matrix nonzero structure via the `MatColoring` object or from the mesh from which the
1446:   matrix comes from (what this function provides). In general using the mesh produces a more optimal coloring (fewer colors).

1448:   This produces a coloring with the distance of 2, see `MatSetColoringDistance()` which can be used for efficiently computing Jacobians with `MatFDColoringCreate()`
1449:   For `DMDA` in three dimensions with periodic boundary conditions the number of grid points in each dimension must be divisible by 2*stencil_width + 1,
1450:   otherwise an error will be generated.

1452: .seealso: [](ch_dmbase), `DM`, `ISColoring`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMSetMatType()`, `MatColoring`, `MatFDColoringCreate()`
1453: @*/
1454: PetscErrorCode DMCreateColoring(DM dm, ISColoringType ctype, ISColoring *coloring)
1455: {
1456:   PetscFunctionBegin;
1458:   PetscAssertPointer(coloring, 3);
1459:   PetscUseTypeMethod(dm, getcoloring, ctype, coloring);
1460:   PetscFunctionReturn(PETSC_SUCCESS);
1461: }

1463: /*@
1464:   DMCreateMatrix - Gets an empty matrix for a `DM` that is most commonly used to store the Jacobian of a discrete PDE operator.

1466:   Collective

1468:   Input Parameter:
1469: . dm - the `DM` object

1471:   Output Parameter:
1472: . mat - the empty Jacobian

1474:   Options Database Key:
1475: . -dm_preallocate_only - Only preallocate the matrix for `DMCreateMatrix()` and `DMCreateMassMatrix()`, but do not fill it with zeros

1477:   Level: beginner

1479:   Notes:
1480:   This properly preallocates the number of nonzeros in the sparse matrix so you
1481:   do not need to do it yourself.

1483:   By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
1484:   the nonzero pattern call `DMSetMatrixPreallocateOnly()`

1486:   For `DMDA`, when you call `MatView()` on this matrix it is displayed using the global natural ordering, NOT in the ordering used
1487:   internally by PETSc.

1489:   For `DMDA`, in general it is easiest to use `MatSetValuesStencil()` or `MatSetValuesLocal()` to put values into the matrix because
1490:   `MatSetValues()` requires the indices for the global numbering for the `DMDA` which is complic`ated to compute

1492: .seealso: [](ch_dmbase), `DM`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMSetMatType()`, `DMCreateMassMatrix()`
1493: @*/
1494: PetscErrorCode DMCreateMatrix(DM dm, Mat *mat)
1495: {
1496:   PetscFunctionBegin;
1498:   PetscAssertPointer(mat, 2);
1499:   PetscCall(MatInitializePackage());
1500:   PetscCall(PetscLogEventBegin(DM_CreateMatrix, 0, 0, 0, 0));
1501:   PetscUseTypeMethod(dm, creatematrix, mat);
1502:   if (PetscDefined(USE_DEBUG)) {
1503:     DM mdm;

1505:     PetscCall(MatGetDM(*mat, &mdm));
1506:     PetscCheck(mdm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DM type '%s' did not attach the DM to the matrix", ((PetscObject)dm)->type_name);
1507:   }
1508:   /* Handle nullspace and near nullspace */
1509:   if (dm->Nf) {
1510:     MatNullSpace nullSpace;
1511:     PetscInt     Nf, f;

1513:     PetscCall(DMGetNumFields(dm, &Nf));
1514:     for (f = 0; f < Nf; ++f) {
1515:       if (dm->nullspaceConstructors[f]) {
1516:         PetscCall((*dm->nullspaceConstructors[f])(dm, f, f, &nullSpace));
1517:         PetscCall(MatSetNullSpace(*mat, nullSpace));
1518:         PetscCall(MatNullSpaceDestroy(&nullSpace));
1519:         break;
1520:       }
1521:     }
1522:     for (f = 0; f < Nf; ++f) {
1523:       if (dm->nearnullspaceConstructors[f]) {
1524:         PetscCall((*dm->nearnullspaceConstructors[f])(dm, f, f, &nullSpace));
1525:         PetscCall(MatSetNearNullSpace(*mat, nullSpace));
1526:         PetscCall(MatNullSpaceDestroy(&nullSpace));
1527:       }
1528:     }
1529:   }
1530:   PetscCall(PetscLogEventEnd(DM_CreateMatrix, 0, 0, 0, 0));
1531:   PetscFunctionReturn(PETSC_SUCCESS);
1532: }

1534: /*@
1535:   DMSetMatrixPreallocateSkip - When `DMCreateMatrix()` is called the matrix sizes and
1536:   `ISLocalToGlobalMapping` will be properly set, but the data structures to store values in the
1537:   matrices will not be preallocated.

1539:   Logically Collective

1541:   Input Parameters:
1542: + dm   - the `DM`
1543: - skip - `PETSC_TRUE` to skip preallocation

1545:   Level: developer

1547:   Note:
1548:   This is most useful to reduce initialization costs when `MatSetPreallocationCOO()` and
1549:   `MatSetValuesCOO()` will be used.

1551: .seealso: [](ch_dmbase), `DM`, `DMCreateMatrix()`, `DMSetMatrixStructureOnly()`, `DMSetMatrixPreallocateOnly()`
1552: @*/
1553: PetscErrorCode DMSetMatrixPreallocateSkip(DM dm, PetscBool skip)
1554: {
1555:   PetscFunctionBegin;
1557:   dm->prealloc_skip = skip;
1558:   PetscFunctionReturn(PETSC_SUCCESS);
1559: }

1561: /*@
1562:   DMSetMatrixPreallocateOnly - When `DMCreateMatrix()` is called the matrix will be properly
1563:   preallocated but the nonzero structure and zero values will not be set.

1565:   Logically Collective

1567:   Input Parameters:
1568: + dm   - the `DM`
1569: - only - `PETSC_TRUE` if only want preallocation

1571:   Options Database Key:
1572: . -dm_preallocate_only - Only preallocate the matrix for `DMCreateMatrix()`, `DMCreateMassMatrix()`, but do not fill it with zeros

1574:   Level: developer

1576: .seealso: [](ch_dmbase), `DM`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMSetMatrixStructureOnly()`, `DMSetMatrixPreallocateSkip()`
1577: @*/
1578: PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only)
1579: {
1580:   PetscFunctionBegin;
1582:   dm->prealloc_only = only;
1583:   PetscFunctionReturn(PETSC_SUCCESS);
1584: }

1586: /*@
1587:   DMSetMatrixStructureOnly - When `DMCreateMatrix()` is called, the matrix structure will be created
1588:   but the array for numerical values will not be allocated.

1590:   Logically Collective

1592:   Input Parameters:
1593: + dm   - the `DM`
1594: - only - `PETSC_TRUE` if you only want matrix structure

1596:   Level: developer

1598: .seealso: [](ch_dmbase), `DM`, `DMCreateMatrix()`, `DMSetMatrixPreallocateOnly()`, `DMSetMatrixPreallocateSkip()`
1599: @*/
1600: PetscErrorCode DMSetMatrixStructureOnly(DM dm, PetscBool only)
1601: {
1602:   PetscFunctionBegin;
1604:   dm->structure_only = only;
1605:   PetscFunctionReturn(PETSC_SUCCESS);
1606: }

1608: /*@
1609:   DMSetBlockingType - set the blocking granularity to be used for variable block size `DMCreateMatrix()` is called

1611:   Logically Collective

1613:   Input Parameters:
1614: + dm    - the `DM`
1615: - btype - block by topological point or field node

1617:   Options Database Key:
1618: . -dm_blocking_type [topological_point, field_node] - use topological point blocking or field node blocking

1620:   Level: advanced

1622: .seealso: [](ch_dmbase), `DM`, `DMCreateMatrix()`, `MatSetVariableBlockSizes()`
1623: @*/
1624: PetscErrorCode DMSetBlockingType(DM dm, DMBlockingType btype)
1625: {
1626:   PetscFunctionBegin;
1628:   dm->blocking_type = btype;
1629:   PetscFunctionReturn(PETSC_SUCCESS);
1630: }

1632: /*@
1633:   DMGetBlockingType - get the blocking granularity to be used for variable block size `DMCreateMatrix()` is called

1635:   Not Collective

1637:   Input Parameter:
1638: . dm - the `DM`

1640:   Output Parameter:
1641: . btype - block by topological point or field node

1643:   Level: advanced

1645: .seealso: [](ch_dmbase), `DM`, `DMCreateMatrix()`, `MatSetVariableBlockSizes()`
1646: @*/
1647: PetscErrorCode DMGetBlockingType(DM dm, DMBlockingType *btype)
1648: {
1649:   PetscFunctionBegin;
1651:   PetscAssertPointer(btype, 2);
1652:   *btype = dm->blocking_type;
1653:   PetscFunctionReturn(PETSC_SUCCESS);
1654: }

1656: /*@C
1657:   DMGetWorkArray - Gets a work array guaranteed to be at least the input size, restore with `DMRestoreWorkArray()`

1659:   Not Collective

1661:   Input Parameters:
1662: + dm    - the `DM` object
1663: . count - The minimum size
1664: - dtype - MPI data type, often `MPIU_REAL`, `MPIU_SCALAR`, or `MPIU_INT`)

1666:   Output Parameter:
1667: . mem - the work array

1669:   Level: developer

1671:   Notes:
1672:   A `DM` may stash the array between instantiations so using this routine may be more efficient than calling `PetscMalloc()`

1674:   The array may contain nonzero values

1676: .seealso: [](ch_dmbase), `DM`, `DMDestroy()`, `DMCreate()`, `DMRestoreWorkArray()`, `PetscMalloc()`
1677: @*/
1678: PetscErrorCode DMGetWorkArray(DM dm, PetscInt count, MPI_Datatype dtype, void *mem)
1679: {
1680:   DMWorkLink  link;
1681:   PetscMPIInt dsize;

1683:   PetscFunctionBegin;
1685:   PetscAssertPointer(mem, 4);
1686:   if (!count) {
1687:     *(void **)mem = NULL;
1688:     PetscFunctionReturn(PETSC_SUCCESS);
1689:   }
1690:   if (dm->workin) {
1691:     link       = dm->workin;
1692:     dm->workin = dm->workin->next;
1693:   } else {
1694:     PetscCall(PetscNew(&link));
1695:   }
1696:   /* Avoid MPI_Type_size for most used datatypes
1697:      Get size directly */
1698:   if (dtype == MPIU_INT) dsize = sizeof(PetscInt);
1699:   else if (dtype == MPIU_REAL) dsize = sizeof(PetscReal);
1700: #if defined(PETSC_USE_64BIT_INDICES)
1701:   else if (dtype == MPI_INT) dsize = sizeof(int);
1702: #endif
1703: #if defined(PETSC_USE_COMPLEX)
1704:   else if (dtype == MPIU_SCALAR) dsize = sizeof(PetscScalar);
1705: #endif
1706:   else PetscCallMPI(MPI_Type_size(dtype, &dsize));

1708:   if (((size_t)dsize * count) > link->bytes) {
1709:     PetscCall(PetscFree(link->mem));
1710:     PetscCall(PetscMalloc(dsize * count, &link->mem));
1711:     link->bytes = dsize * count;
1712:   }
1713:   link->next    = dm->workout;
1714:   dm->workout   = link;
1715:   *(void **)mem = link->mem;
1716:   PetscFunctionReturn(PETSC_SUCCESS);
1717: }

1719: /*@C
1720:   DMRestoreWorkArray - Restores a work array obtained with `DMCreateWorkArray()`

1722:   Not Collective

1724:   Input Parameters:
1725: + dm    - the `DM` object
1726: . count - The minimum size
1727: - dtype - MPI data type, often `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_INT`

1729:   Output Parameter:
1730: . mem - the work array

1732:   Level: developer

1734:   Developer Note:
1735:   count and dtype are ignored, they are only needed for `DMGetWorkArray()`

1737: .seealso: [](ch_dmbase), `DM`, `DMDestroy()`, `DMCreate()`, `DMGetWorkArray()`
1738: @*/
1739: PetscErrorCode DMRestoreWorkArray(DM dm, PetscInt count, MPI_Datatype dtype, void *mem)
1740: {
1741:   DMWorkLink *p, link;

1743:   PetscFunctionBegin;
1745:   PetscAssertPointer(mem, 4);
1746:   (void)count;
1747:   (void)dtype;
1748:   if (!*(void **)mem) PetscFunctionReturn(PETSC_SUCCESS);
1749:   for (p = &dm->workout; (link = *p); p = &link->next) {
1750:     if (link->mem == *(void **)mem) {
1751:       *p            = link->next;
1752:       link->next    = dm->workin;
1753:       dm->workin    = link;
1754:       *(void **)mem = NULL;
1755:       PetscFunctionReturn(PETSC_SUCCESS);
1756:     }
1757:   }
1758:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Array was not checked out");
1759: }

1761: /*@C
1762:   DMSetNullSpaceConstructor - Provide a callback function which constructs the nullspace for a given field, defined with `DMAddField()`, when function spaces
1763:   are joined or split, such as in `DMCreateSubDM()`

1765:   Logically Collective; No Fortran Support

1767:   Input Parameters:
1768: + dm     - The `DM`
1769: . field  - The field number for the nullspace
1770: - nullsp - A callback to create the nullspace

1772:   Calling sequence of `nullsp`:
1773: + dm        - The present `DM`
1774: . origField - The field number given above, in the original `DM`
1775: . field     - The field number in dm
1776: - nullSpace - The nullspace for the given field

1778:   Level: intermediate

1780: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMGetNullSpaceConstructor()`, `DMSetNearNullSpaceConstructor()`, `DMGetNearNullSpaceConstructor()`, `DMCreateSubDM()`, `DMCreateSuperDM()`
1781: @*/
1782: PetscErrorCode DMSetNullSpaceConstructor(DM dm, PetscInt field, PetscErrorCode (*nullsp)(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullSpace))
1783: {
1784:   PetscFunctionBegin;
1786:   PetscCheck(field < 10, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %" PetscInt_FMT " >= 10 fields", field);
1787:   dm->nullspaceConstructors[field] = nullsp;
1788:   PetscFunctionReturn(PETSC_SUCCESS);
1789: }

1791: /*@C
1792:   DMGetNullSpaceConstructor - Return the callback function which constructs the nullspace for a given field, defined with `DMAddField()`

1794:   Not Collective; No Fortran Support

1796:   Input Parameters:
1797: + dm    - The `DM`
1798: - field - The field number for the nullspace

1800:   Output Parameter:
1801: . nullsp - A callback to create the nullspace

1803:   Calling sequence of `nullsp`:
1804: + dm        - The present DM
1805: . origField - The field number given above, in the original DM
1806: . field     - The field number in dm
1807: - nullSpace - The nullspace for the given field

1809:   Level: intermediate

1811: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMGetField()`, `DMSetNullSpaceConstructor()`, `DMSetNearNullSpaceConstructor()`, `DMGetNearNullSpaceConstructor()`, `DMCreateSubDM()`, `DMCreateSuperDM()`
1812: @*/
1813: PetscErrorCode DMGetNullSpaceConstructor(DM dm, PetscInt field, PetscErrorCode (**nullsp)(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullSpace))
1814: {
1815:   PetscFunctionBegin;
1817:   PetscAssertPointer(nullsp, 3);
1818:   PetscCheck(field < 10, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %" PetscInt_FMT " >= 10 fields", field);
1819:   *nullsp = dm->nullspaceConstructors[field];
1820:   PetscFunctionReturn(PETSC_SUCCESS);
1821: }

1823: /*@C
1824:   DMSetNearNullSpaceConstructor - Provide a callback function which constructs the near-nullspace for a given field, defined with `DMAddField()`

1826:   Logically Collective; No Fortran Support

1828:   Input Parameters:
1829: + dm     - The `DM`
1830: . field  - The field number for the nullspace
1831: - nullsp - A callback to create the near-nullspace

1833:   Calling sequence of `nullsp`:
1834: + dm        - The present `DM`
1835: . origField - The field number given above, in the original `DM`
1836: . field     - The field number in dm
1837: - nullSpace - The nullspace for the given field

1839:   Level: intermediate

1841: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMGetNearNullSpaceConstructor()`, `DMSetNullSpaceConstructor()`, `DMGetNullSpaceConstructor()`, `DMCreateSubDM()`, `DMCreateSuperDM()`,
1842:           `MatNullSpace`
1843: @*/
1844: PetscErrorCode DMSetNearNullSpaceConstructor(DM dm, PetscInt field, PetscErrorCode (*nullsp)(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullSpace))
1845: {
1846:   PetscFunctionBegin;
1848:   PetscCheck(field < 10, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %" PetscInt_FMT " >= 10 fields", field);
1849:   dm->nearnullspaceConstructors[field] = nullsp;
1850:   PetscFunctionReturn(PETSC_SUCCESS);
1851: }

1853: /*@C
1854:   DMGetNearNullSpaceConstructor - Return the callback function which constructs the near-nullspace for a given field, defined with `DMAddField()`

1856:   Not Collective; No Fortran Support

1858:   Input Parameters:
1859: + dm    - The `DM`
1860: - field - The field number for the nullspace

1862:   Output Parameter:
1863: . nullsp - A callback to create the near-nullspace

1865:   Calling sequence of `nullsp`:
1866: + dm        - The present `DM`
1867: . origField - The field number given above, in the original `DM`
1868: . field     - The field number in dm
1869: - nullSpace - The nullspace for the given field

1871:   Level: intermediate

1873: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMGetField()`, `DMSetNearNullSpaceConstructor()`, `DMSetNullSpaceConstructor()`, `DMGetNullSpaceConstructor()`, `DMCreateSubDM()`,
1874:           `MatNullSpace`, `DMCreateSuperDM()`
1875: @*/
1876: PetscErrorCode DMGetNearNullSpaceConstructor(DM dm, PetscInt field, PetscErrorCode (**nullsp)(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullSpace))
1877: {
1878:   PetscFunctionBegin;
1880:   PetscAssertPointer(nullsp, 3);
1881:   PetscCheck(field < 10, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %" PetscInt_FMT " >= 10 fields", field);
1882:   *nullsp = dm->nearnullspaceConstructors[field];
1883:   PetscFunctionReturn(PETSC_SUCCESS);
1884: }

1886: /*@C
1887:   DMCreateFieldIS - Creates a set of `IS` objects with the global indices of dofs for each field defined with `DMAddField()`

1889:   Not Collective; No Fortran Support

1891:   Input Parameter:
1892: . dm - the `DM` object

1894:   Output Parameters:
1895: + numFields  - The number of fields (or `NULL` if not requested)
1896: . fieldNames - The number of each field (or `NULL` if not requested)
1897: - fields     - The global indices for each field (or `NULL` if not requested)

1899:   Level: intermediate

1901:   Note:
1902:   The user is responsible for freeing all requested arrays. In particular, every entry of `fieldNames` should be freed with
1903:   `PetscFree()`, every entry of `fields` should be destroyed with `ISDestroy()`, and both arrays should be freed with
1904:   `PetscFree()`.

1906:   Developer Note:
1907:   It is not clear why both this function and `DMCreateFieldDecomposition()` exist. Having two seems redundant and confusing. This function should
1908:   likely be removed.

1910: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMGetField()`, `DMDestroy()`, `DMView()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`,
1911:           `DMCreateFieldDecomposition()`
1912: @*/
1913: PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
1914: {
1915:   PetscSection section, sectionGlobal;

1917:   PetscFunctionBegin;
1919:   if (numFields) {
1920:     PetscAssertPointer(numFields, 2);
1921:     *numFields = 0;
1922:   }
1923:   if (fieldNames) {
1924:     PetscAssertPointer(fieldNames, 3);
1925:     *fieldNames = NULL;
1926:   }
1927:   if (fields) {
1928:     PetscAssertPointer(fields, 4);
1929:     *fields = NULL;
1930:   }
1931:   PetscCall(DMGetLocalSection(dm, &section));
1932:   if (section) {
1933:     PetscInt *fieldSizes, *fieldNc, **fieldIndices;
1934:     PetscInt  nF, f, pStart, pEnd, p;

1936:     PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
1937:     PetscCall(PetscSectionGetNumFields(section, &nF));
1938:     PetscCall(PetscMalloc3(nF, &fieldSizes, nF, &fieldNc, nF, &fieldIndices));
1939:     PetscCall(PetscSectionGetChart(sectionGlobal, &pStart, &pEnd));
1940:     for (f = 0; f < nF; ++f) {
1941:       fieldSizes[f] = 0;
1942:       PetscCall(PetscSectionGetFieldComponents(section, f, &fieldNc[f]));
1943:     }
1944:     for (p = pStart; p < pEnd; ++p) {
1945:       PetscInt gdof;

1947:       PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof));
1948:       if (gdof > 0) {
1949:         for (f = 0; f < nF; ++f) {
1950:           PetscInt fdof, fcdof, fpdof;

1952:           PetscCall(PetscSectionGetFieldDof(section, p, f, &fdof));
1953:           PetscCall(PetscSectionGetFieldConstraintDof(section, p, f, &fcdof));
1954:           fpdof = fdof - fcdof;
1955:           if (fpdof && fpdof != fieldNc[f]) {
1956:             /* Layout does not admit a pointwise block size */
1957:             fieldNc[f] = 1;
1958:           }
1959:           fieldSizes[f] += fpdof;
1960:         }
1961:       }
1962:     }
1963:     for (f = 0; f < nF; ++f) {
1964:       PetscCall(PetscMalloc1(fieldSizes[f], &fieldIndices[f]));
1965:       fieldSizes[f] = 0;
1966:     }
1967:     for (p = pStart; p < pEnd; ++p) {
1968:       PetscInt gdof, goff;

1970:       PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof));
1971:       if (gdof > 0) {
1972:         PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
1973:         for (f = 0; f < nF; ++f) {
1974:           PetscInt fdof, fcdof, fc;

1976:           PetscCall(PetscSectionGetFieldDof(section, p, f, &fdof));
1977:           PetscCall(PetscSectionGetFieldConstraintDof(section, p, f, &fcdof));
1978:           for (fc = 0; fc < fdof - fcdof; ++fc, ++fieldSizes[f]) fieldIndices[f][fieldSizes[f]] = goff++;
1979:         }
1980:       }
1981:     }
1982:     if (numFields) *numFields = nF;
1983:     if (fieldNames) {
1984:       PetscCall(PetscMalloc1(nF, fieldNames));
1985:       for (f = 0; f < nF; ++f) {
1986:         const char *fieldName;

1988:         PetscCall(PetscSectionGetFieldName(section, f, &fieldName));
1989:         PetscCall(PetscStrallocpy(fieldName, &(*fieldNames)[f]));
1990:       }
1991:     }
1992:     if (fields) {
1993:       PetscCall(PetscMalloc1(nF, fields));
1994:       for (f = 0; f < nF; ++f) {
1995:         PetscInt bs, in[2], out[2];

1997:         PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), fieldSizes[f], fieldIndices[f], PETSC_OWN_POINTER, &(*fields)[f]));
1998:         in[0] = -fieldNc[f];
1999:         in[1] = fieldNc[f];
2000:         PetscCallMPI(MPIU_Allreduce(in, out, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
2001:         bs = (-out[0] == out[1]) ? out[1] : 1;
2002:         PetscCall(ISSetBlockSize((*fields)[f], bs));
2003:       }
2004:     }
2005:     PetscCall(PetscFree3(fieldSizes, fieldNc, fieldIndices));
2006:   } else PetscTryTypeMethod(dm, createfieldis, numFields, fieldNames, fields);
2007:   PetscFunctionReturn(PETSC_SUCCESS);
2008: }

2010: /*@C
2011:   DMCreateFieldDecomposition - Returns a list of `IS` objects defining a decomposition of a problem into subproblems
2012:   corresponding to different fields.

2014:   Not Collective; No Fortran Support

2016:   Input Parameter:
2017: . dm - the `DM` object

2019:   Output Parameters:
2020: + len      - The number of fields (or `NULL` if not requested)
2021: . namelist - The name for each field (or `NULL` if not requested)
2022: . islist   - The global indices for each field (or `NULL` if not requested)
2023: - dmlist   - The `DM`s for each field subproblem (or `NULL`, if not requested; if `NULL` is returned, no `DM`s are defined)

2025:   Level: intermediate

2027:   Notes:
2028:   Each `IS` contains the global indices of the dofs of the corresponding field, defined by
2029:   `DMAddField()`. The optional list of `DM`s define the `DM` for each subproblem.

2031:   The same as `DMCreateFieldIS()` but also returns a `DM` for each field.

2033:   The user is responsible for freeing all requested arrays. In particular, every entry of `namelist` should be freed with
2034:   `PetscFree()`, every entry of `islist` should be destroyed with `ISDestroy()`, every entry of `dmlist` should be destroyed with `DMDestroy()`,
2035:   and all of the arrays should be freed with `PetscFree()`.

2037:   Developer Notes:
2038:   It is not clear why this function and `DMCreateFieldIS()` exist. Having two seems redundant and confusing.

2040:   Unlike  `DMRefine()`, `DMCoarsen()`, and `DMCreateDomainDecomposition()` this provides no mechanism to provide hooks that are called after the
2041:   decomposition is computed.

2043: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMCreateFieldIS()`, `DMCreateSubDM()`, `DMCreateDomainDecomposition()`, `DMDestroy()`, `DMView()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMRefine()`, `DMCoarsen()`
2044: @*/
2045: PetscErrorCode DMCreateFieldDecomposition(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
2046: {
2047:   PetscFunctionBegin;
2049:   if (len) {
2050:     PetscAssertPointer(len, 2);
2051:     *len = 0;
2052:   }
2053:   if (namelist) {
2054:     PetscAssertPointer(namelist, 3);
2055:     *namelist = NULL;
2056:   }
2057:   if (islist) {
2058:     PetscAssertPointer(islist, 4);
2059:     *islist = NULL;
2060:   }
2061:   if (dmlist) {
2062:     PetscAssertPointer(dmlist, 5);
2063:     *dmlist = NULL;
2064:   }
2065:   /*
2066:    Is it a good idea to apply the following check across all impls?
2067:    Perhaps some impls can have a well-defined decomposition before DMSetUp?
2068:    This, however, follows the general principle that accessors are not well-behaved until the object is set up.
2069:    */
2070:   PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Decomposition defined only after DMSetUp");
2071:   if (!dm->ops->createfielddecomposition) {
2072:     PetscSection section;
2073:     PetscInt     numFields, f;

2075:     PetscCall(DMGetLocalSection(dm, &section));
2076:     if (section) PetscCall(PetscSectionGetNumFields(section, &numFields));
2077:     if (section && numFields && dm->ops->createsubdm) {
2078:       if (len) *len = numFields;
2079:       if (namelist) PetscCall(PetscMalloc1(numFields, namelist));
2080:       if (islist) PetscCall(PetscMalloc1(numFields, islist));
2081:       if (dmlist) PetscCall(PetscMalloc1(numFields, dmlist));
2082:       for (f = 0; f < numFields; ++f) {
2083:         const char *fieldName;

2085:         PetscCall(DMCreateSubDM(dm, 1, &f, islist ? &(*islist)[f] : NULL, dmlist ? &(*dmlist)[f] : NULL));
2086:         if (namelist) {
2087:           PetscCall(PetscSectionGetFieldName(section, f, &fieldName));
2088:           PetscCall(PetscStrallocpy(fieldName, &(*namelist)[f]));
2089:         }
2090:       }
2091:     } else {
2092:       PetscCall(DMCreateFieldIS(dm, len, namelist, islist));
2093:       /* By default there are no DMs associated with subproblems. */
2094:       if (dmlist) *dmlist = NULL;
2095:     }
2096:   } else PetscUseTypeMethod(dm, createfielddecomposition, len, namelist, islist, dmlist);
2097:   PetscFunctionReturn(PETSC_SUCCESS);
2098: }

2100: /*@
2101:   DMCreateSubDM - Returns an `IS` and `DM` encapsulating a subproblem defined by the fields passed in.
2102:   The fields are defined by `DMCreateFieldIS()`.

2104:   Not collective

2106:   Input Parameters:
2107: + dm        - The `DM` object
2108: . numFields - The number of fields to select
2109: - fields    - The field numbers of the selected fields

2111:   Output Parameters:
2112: + is    - The global indices for all the degrees of freedom in the new sub `DM`, use `NULL` if not needed
2113: - subdm - The `DM` for the subproblem, use `NULL` if not needed

2115:   Level: intermediate

2117:   Note:
2118:   You need to call `DMPlexSetMigrationSF()` on the original `DM` if you want the Global-To-Natural map to be automatically constructed

2120: .seealso: [](ch_dmbase), `DM`, `DMCreateFieldIS()`, `DMCreateFieldDecomposition()`, `DMAddField()`, `DMCreateSuperDM()`, `IS`, `DMPlexSetMigrationSF()`, `DMDestroy()`, `DMView()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`
2121: @*/
2122: PetscErrorCode DMCreateSubDM(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
2123: {
2124:   PetscFunctionBegin;
2126:   PetscAssertPointer(fields, 3);
2127:   if (is) PetscAssertPointer(is, 4);
2128:   if (subdm) PetscAssertPointer(subdm, 5);
2129:   PetscUseTypeMethod(dm, createsubdm, numFields, fields, is, subdm);
2130:   PetscFunctionReturn(PETSC_SUCCESS);
2131: }

2133: /*@C
2134:   DMCreateSuperDM - Returns an arrays of `IS` and `DM` encapsulating a superproblem defined by multiple `DM`s passed in.

2136:   Not collective

2138:   Input Parameters:
2139: + dms - The `DM` objects
2140: - n   - The number of `DM`s

2142:   Output Parameters:
2143: + is      - The global indices for each of subproblem within the super `DM`, or NULL
2144: - superdm - The `DM` for the superproblem

2146:   Level: intermediate

2148:   Note:
2149:   You need to call `DMPlexSetMigrationSF()` on the original `DM` if you want the Global-To-Natural map to be automatically constructed

2151: .seealso: [](ch_dmbase), `DM`, `DMCreateSubDM()`, `DMPlexSetMigrationSF()`, `DMDestroy()`, `DMView()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMCreateFieldIS()`, `DMCreateDomainDecomposition()`
2152: @*/
2153: PetscErrorCode DMCreateSuperDM(DM dms[], PetscInt n, IS *is[], DM *superdm)
2154: {
2155:   PetscInt i;

2157:   PetscFunctionBegin;
2158:   PetscAssertPointer(dms, 1);
2160:   if (is) PetscAssertPointer(is, 3);
2161:   PetscAssertPointer(superdm, 4);
2162:   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of DMs must be nonnegative: %" PetscInt_FMT, n);
2163:   if (n) {
2164:     DM dm = dms[0];
2165:     PetscCheck(dm->ops->createsuperdm, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No method createsuperdm for DM of type %s", ((PetscObject)dm)->type_name);
2166:     PetscCall((*dm->ops->createsuperdm)(dms, n, is, superdm));
2167:   }
2168:   PetscFunctionReturn(PETSC_SUCCESS);
2169: }

2171: /*@C
2172:   DMCreateDomainDecomposition - Returns lists of `IS` objects defining a decomposition of a
2173:   problem into subproblems corresponding to restrictions to pairs of nested subdomains.

2175:   Not Collective

2177:   Input Parameter:
2178: . dm - the `DM` object

2180:   Output Parameters:
2181: + n           - The number of subproblems in the domain decomposition (or `NULL` if not requested)
2182: . namelist    - The name for each subdomain (or `NULL` if not requested)
2183: . innerislist - The global indices for each inner subdomain (or `NULL`, if not requested)
2184: . outerislist - The global indices for each outer subdomain (or `NULL`, if not requested)
2185: - dmlist      - The `DM`s for each subdomain subproblem (or `NULL`, if not requested; if `NULL` is returned, no `DM`s are defined)

2187:   Level: intermediate

2189:   Notes:
2190:   Each `IS` contains the global indices of the dofs of the corresponding subdomains with in the
2191:   dofs of the original `DM`. The inner subdomains conceptually define a nonoverlapping
2192:   covering, while outer subdomains can overlap.

2194:   The optional list of `DM`s define a `DM` for each subproblem.

2196:   The user is responsible for freeing all requested arrays. In particular, every entry of `namelist` should be freed with
2197:   `PetscFree()`, every entry of `innerislist` and `outerislist` should be destroyed with `ISDestroy()`, every entry of `dmlist` should be destroyed with `DMDestroy()`,
2198:   and all of the arrays should be freed with `PetscFree()`.

2200:   Developer Notes:
2201:   The `dmlist` is for the inner subdomains or the outer subdomains or all subdomains?

2203:   The names are inconsistent, the hooks use `DMSubDomainHook` which is nothing like `DMCreateDomainDecomposition()` while `DMRefineHook` is used for `DMRefine()`.

2205: .seealso: [](ch_dmbase), `DM`, `DMCreateFieldDecomposition()`, `DMDestroy()`, `DMCreateDomainDecompositionScatters()`, `DMView()`, `DMCreateInterpolation()`,
2206:           `DMSubDomainHookAdd()`, `DMSubDomainHookRemove()`,`DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMRefine()`, `DMCoarsen()`
2207: @*/
2208: PetscErrorCode DMCreateDomainDecomposition(DM dm, PetscInt *n, char ***namelist, IS **innerislist, IS **outerislist, DM **dmlist)
2209: {
2210:   DMSubDomainHookLink link;
2211:   PetscInt            i, l;

2213:   PetscFunctionBegin;
2215:   if (n) {
2216:     PetscAssertPointer(n, 2);
2217:     *n = 0;
2218:   }
2219:   if (namelist) {
2220:     PetscAssertPointer(namelist, 3);
2221:     *namelist = NULL;
2222:   }
2223:   if (innerislist) {
2224:     PetscAssertPointer(innerislist, 4);
2225:     *innerislist = NULL;
2226:   }
2227:   if (outerislist) {
2228:     PetscAssertPointer(outerislist, 5);
2229:     *outerislist = NULL;
2230:   }
2231:   if (dmlist) {
2232:     PetscAssertPointer(dmlist, 6);
2233:     *dmlist = NULL;
2234:   }
2235:   /*
2236:    Is it a good idea to apply the following check across all impls?
2237:    Perhaps some impls can have a well-defined decomposition before DMSetUp?
2238:    This, however, follows the general principle that accessors are not well-behaved until the object is set up.
2239:    */
2240:   PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Decomposition defined only after DMSetUp");
2241:   if (dm->ops->createdomaindecomposition) {
2242:     PetscUseTypeMethod(dm, createdomaindecomposition, &l, namelist, innerislist, outerislist, dmlist);
2243:     /* copy subdomain hooks and context over to the subdomain DMs */
2244:     if (dmlist && *dmlist) {
2245:       for (i = 0; i < l; i++) {
2246:         for (link = dm->subdomainhook; link; link = link->next) {
2247:           if (link->ddhook) PetscCall((*link->ddhook)(dm, (*dmlist)[i], link->ctx));
2248:         }
2249:         if (dm->ctx) (*dmlist)[i]->ctx = dm->ctx;
2250:       }
2251:     }
2252:     if (n) *n = l;
2253:   }
2254:   PetscFunctionReturn(PETSC_SUCCESS);
2255: }

2257: /*@C
2258:   DMCreateDomainDecompositionScatters - Returns scatters to the subdomain vectors from the global vector for subdomains created with
2259:   `DMCreateDomainDecomposition()`

2261:   Not Collective

2263:   Input Parameters:
2264: + dm     - the `DM` object
2265: . n      - the number of subdomains
2266: - subdms - the local subdomains

2268:   Output Parameters:
2269: + iscat - scatter from global vector to nonoverlapping global vector entries on subdomain
2270: . oscat - scatter from global vector to overlapping global vector entries on subdomain
2271: - gscat - scatter from global vector to local vector on subdomain (fills in ghosts)

2273:   Level: developer

2275:   Note:
2276:   This is an alternative to the iis and ois arguments in `DMCreateDomainDecomposition()` that allow for the solution
2277:   of general nonlinear problems with overlapping subdomain methods.  While merely having index sets that enable subsets
2278:   of the residual equations to be created is fine for linear problems, nonlinear problems require local assembly of
2279:   solution and residual data.

2281:   Developer Note:
2282:   Can the subdms input be anything or are they exactly the `DM` obtained from
2283:   `DMCreateDomainDecomposition()`?

2285: .seealso: [](ch_dmbase), `DM`, `DMCreateDomainDecomposition()`, `DMDestroy()`, `DMView()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMCreateFieldIS()`
2286: @*/
2287: PetscErrorCode DMCreateDomainDecompositionScatters(DM dm, PetscInt n, DM *subdms, VecScatter *iscat[], VecScatter *oscat[], VecScatter *gscat[])
2288: {
2289:   PetscFunctionBegin;
2291:   PetscAssertPointer(subdms, 3);
2292:   PetscUseTypeMethod(dm, createddscatters, n, subdms, iscat, oscat, gscat);
2293:   PetscFunctionReturn(PETSC_SUCCESS);
2294: }

2296: /*@
2297:   DMRefine - Refines a `DM` object using a standard nonadaptive refinement of the underlying mesh

2299:   Collective

2301:   Input Parameters:
2302: + dm   - the `DM` object
2303: - comm - the communicator to contain the new `DM` object (or `MPI_COMM_NULL`)

2305:   Output Parameter:
2306: . dmf - the refined `DM`, or `NULL`

2308:   Options Database Key:
2309: . -dm_plex_cell_refiner <strategy> - chooses the refinement strategy, e.g. regular, tohex

2311:   Level: developer

2313:   Note:
2314:   If no refinement was done, the return value is `NULL`

2316: .seealso: [](ch_dmbase), `DM`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateDomainDecomposition()`,
2317:           `DMRefineHookAdd()`, `DMRefineHookRemove()`
2318: @*/
2319: PetscErrorCode DMRefine(DM dm, MPI_Comm comm, DM *dmf)
2320: {
2321:   DMRefineHookLink link;

2323:   PetscFunctionBegin;
2325:   PetscCall(PetscLogEventBegin(DM_Refine, dm, 0, 0, 0));
2326:   PetscUseTypeMethod(dm, refine, comm, dmf);
2327:   if (*dmf) {
2328:     (*dmf)->ops->creatematrix = dm->ops->creatematrix;

2330:     PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)dm, (PetscObject)*dmf));

2332:     (*dmf)->ctx       = dm->ctx;
2333:     (*dmf)->leveldown = dm->leveldown;
2334:     (*dmf)->levelup   = dm->levelup + 1;

2336:     PetscCall(DMSetMatType(*dmf, dm->mattype));
2337:     for (link = dm->refinehook; link; link = link->next) {
2338:       if (link->refinehook) PetscCall((*link->refinehook)(dm, *dmf, link->ctx));
2339:     }
2340:   }
2341:   PetscCall(PetscLogEventEnd(DM_Refine, dm, 0, 0, 0));
2342:   PetscFunctionReturn(PETSC_SUCCESS);
2343: }

2345: /*@C
2346:   DMRefineHookAdd - adds a callback to be run when interpolating a nonlinear problem to a finer grid

2348:   Logically Collective; No Fortran Support

2350:   Input Parameters:
2351: + coarse     - `DM` on which to run a hook when interpolating to a finer level
2352: . refinehook - function to run when setting up the finer level
2353: . interphook - function to run to update data on finer levels (once per `SNESSolve()`)
2354: - ctx        - [optional] user-defined context for provide data for the hooks (may be `NULL`)

2356:   Calling sequence of `refinehook`:
2357: + coarse - coarse level `DM`
2358: . fine   - fine level `DM` to interpolate problem to
2359: - ctx    - optional user-defined function context

2361:   Calling sequence of `interphook`:
2362: + coarse - coarse level `DM`
2363: . interp - matrix interpolating a coarse-level solution to the finer grid
2364: . fine   - fine level `DM` to update
2365: - ctx    - optional user-defined function context

2367:   Level: advanced

2369:   Notes:
2370:   This function is only needed if auxiliary data that is attached to the `DM`s via, for example, `PetscObjectCompose()`, needs to be
2371:   passed to fine grids while grid sequencing.

2373:   The actual interpolation is done when `DMInterpolate()` is called.

2375:   If this function is called multiple times, the hooks will be run in the order they are added.

2377: .seealso: [](ch_dmbase), `DM`, `DMCoarsenHookAdd()`, `DMInterpolate()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`
2378: @*/
2379: PetscErrorCode DMRefineHookAdd(DM coarse, PetscErrorCode (*refinehook)(DM coarse, DM fine, void *ctx), PetscErrorCode (*interphook)(DM coarse, Mat interp, DM fine, void *ctx), void *ctx)
2380: {
2381:   DMRefineHookLink link, *p;

2383:   PetscFunctionBegin;
2385:   for (p = &coarse->refinehook; *p; p = &(*p)->next) { /* Scan to the end of the current list of hooks */
2386:     if ((*p)->refinehook == refinehook && (*p)->interphook == interphook && (*p)->ctx == ctx) PetscFunctionReturn(PETSC_SUCCESS);
2387:   }
2388:   PetscCall(PetscNew(&link));
2389:   link->refinehook = refinehook;
2390:   link->interphook = interphook;
2391:   link->ctx        = ctx;
2392:   link->next       = NULL;
2393:   *p               = link;
2394:   PetscFunctionReturn(PETSC_SUCCESS);
2395: }

2397: /*@C
2398:   DMRefineHookRemove - remove a callback from the list of hooks, that have been set with `DMRefineHookAdd()`, to be run when interpolating
2399:   a nonlinear problem to a finer grid

2401:   Logically Collective; No Fortran Support

2403:   Input Parameters:
2404: + coarse     - the `DM` on which to run a hook when restricting to a coarser level
2405: . refinehook - function to run when setting up a finer level
2406: . interphook - function to run to update data on finer levels
2407: - ctx        - [optional] user-defined context for provide data for the hooks (may be `NULL`)

2409:   Level: advanced

2411:   Note:
2412:   This function does nothing if the hook is not in the list.

2414: .seealso: [](ch_dmbase), `DM`, `DMRefineHookAdd()`, `DMCoarsenHookRemove()`, `DMInterpolate()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`
2415: @*/
2416: PetscErrorCode DMRefineHookRemove(DM coarse, PetscErrorCode (*refinehook)(DM, DM, void *), PetscErrorCode (*interphook)(DM, Mat, DM, void *), void *ctx)
2417: {
2418:   DMRefineHookLink link, *p;

2420:   PetscFunctionBegin;
2422:   for (p = &coarse->refinehook; *p; p = &(*p)->next) { /* Search the list of current hooks */
2423:     if ((*p)->refinehook == refinehook && (*p)->interphook == interphook && (*p)->ctx == ctx) {
2424:       link = *p;
2425:       *p   = link->next;
2426:       PetscCall(PetscFree(link));
2427:       break;
2428:     }
2429:   }
2430:   PetscFunctionReturn(PETSC_SUCCESS);
2431: }

2433: /*@
2434:   DMInterpolate - interpolates user-defined problem data attached to a `DM` to a finer `DM` by running hooks registered by `DMRefineHookAdd()`

2436:   Collective if any hooks are

2438:   Input Parameters:
2439: + coarse - coarser `DM` to use as a base
2440: . interp - interpolation matrix, apply using `MatInterpolate()`
2441: - fine   - finer `DM` to update

2443:   Level: developer

2445:   Developer Note:
2446:   This routine is called `DMInterpolate()` while the hook is called `DMRefineHookAdd()`. It would be better to have an
2447:   an API with consistent terminology.

2449: .seealso: [](ch_dmbase), `DM`, `DMRefineHookAdd()`, `MatInterpolate()`
2450: @*/
2451: PetscErrorCode DMInterpolate(DM coarse, Mat interp, DM fine)
2452: {
2453:   DMRefineHookLink link;

2455:   PetscFunctionBegin;
2456:   for (link = fine->refinehook; link; link = link->next) {
2457:     if (link->interphook) PetscCall((*link->interphook)(coarse, interp, fine, link->ctx));
2458:   }
2459:   PetscFunctionReturn(PETSC_SUCCESS);
2460: }

2462: /*@
2463:   DMInterpolateSolution - Interpolates a solution from a coarse mesh to a fine mesh.

2465:   Collective

2467:   Input Parameters:
2468: + coarse    - coarse `DM`
2469: . fine      - fine `DM`
2470: . interp    - (optional) the matrix computed by `DMCreateInterpolation()`.  Implementations may not need this, but if it
2471:             is available it can avoid some recomputation.  If it is provided, `MatInterpolate()` will be used if
2472:             the coarse `DM` does not have a specialized implementation.
2473: - coarseSol - solution on the coarse mesh

2475:   Output Parameter:
2476: . fineSol - the interpolation of coarseSol to the fine mesh

2478:   Level: developer

2480:   Note:
2481:   This function exists because the interpolation of a solution vector between meshes is not always a linear
2482:   map.  For example, if a boundary value problem has an inhomogeneous Dirichlet boundary condition that is compressed
2483:   out of the solution vector.  Or if interpolation is inherently a nonlinear operation, such as a method using
2484:   slope-limiting reconstruction.

2486:   Developer Note:
2487:   This doesn't just interpolate "solutions" so its API name is questionable.

2489: .seealso: [](ch_dmbase), `DM`, `DMInterpolate()`, `DMCreateInterpolation()`
2490: @*/
2491: PetscErrorCode DMInterpolateSolution(DM coarse, DM fine, Mat interp, Vec coarseSol, Vec fineSol)
2492: {
2493:   PetscErrorCode (*interpsol)(DM, DM, Mat, Vec, Vec) = NULL;

2495:   PetscFunctionBegin;

2501:   PetscCall(PetscObjectQueryFunction((PetscObject)coarse, "DMInterpolateSolution_C", &interpsol));
2502:   if (interpsol) {
2503:     PetscCall((*interpsol)(coarse, fine, interp, coarseSol, fineSol));
2504:   } else if (interp) {
2505:     PetscCall(MatInterpolate(interp, coarseSol, fineSol));
2506:   } else SETERRQ(PetscObjectComm((PetscObject)coarse), PETSC_ERR_SUP, "DM %s does not implement DMInterpolateSolution()", ((PetscObject)coarse)->type_name);
2507:   PetscFunctionReturn(PETSC_SUCCESS);
2508: }

2510: /*@
2511:   DMGetRefineLevel - Gets the number of refinements that have generated this `DM` from some initial `DM`.

2513:   Not Collective

2515:   Input Parameter:
2516: . dm - the `DM` object

2518:   Output Parameter:
2519: . level - number of refinements

2521:   Level: developer

2523:   Note:
2524:   This can be used, by example, to set the number of coarser levels associated with this `DM` for a multigrid solver.

2526: .seealso: [](ch_dmbase), `DM`, `DMRefine()`, `DMCoarsen()`, `DMGetCoarsenLevel()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`
2527: @*/
2528: PetscErrorCode DMGetRefineLevel(DM dm, PetscInt *level)
2529: {
2530:   PetscFunctionBegin;
2532:   *level = dm->levelup;
2533:   PetscFunctionReturn(PETSC_SUCCESS);
2534: }

2536: /*@
2537:   DMSetRefineLevel - Sets the number of refinements that have generated this `DM`.

2539:   Not Collective

2541:   Input Parameters:
2542: + dm    - the `DM` object
2543: - level - number of refinements

2545:   Level: advanced

2547:   Notes:
2548:   This value is used by `PCMG` to determine how many multigrid levels to use

2550:   The values are usually set automatically by the process that is causing the refinements of an initial `DM` by calling this routine.

2552: .seealso: [](ch_dmbase), `DM`, `DMGetRefineLevel()`, `DMCoarsen()`, `DMGetCoarsenLevel()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`
2553: @*/
2554: PetscErrorCode DMSetRefineLevel(DM dm, PetscInt level)
2555: {
2556:   PetscFunctionBegin;
2558:   dm->levelup = level;
2559:   PetscFunctionReturn(PETSC_SUCCESS);
2560: }

2562: /*@
2563:   DMExtrude - Extrude a `DM` object from a surface

2565:   Collective

2567:   Input Parameters:
2568: + dm     - the `DM` object
2569: - layers - the number of extruded cell layers

2571:   Output Parameter:
2572: . dme - the extruded `DM`, or `NULL`

2574:   Level: developer

2576:   Note:
2577:   If no extrusion was done, the return value is `NULL`

2579: .seealso: [](ch_dmbase), `DM`, `DMRefine()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`
2580: @*/
2581: PetscErrorCode DMExtrude(DM dm, PetscInt layers, DM *dme)
2582: {
2583:   PetscFunctionBegin;
2585:   PetscUseTypeMethod(dm, extrude, layers, dme);
2586:   if (*dme) {
2587:     (*dme)->ops->creatematrix = dm->ops->creatematrix;
2588:     PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)dm, (PetscObject)*dme));
2589:     (*dme)->ctx = dm->ctx;
2590:     PetscCall(DMSetMatType(*dme, dm->mattype));
2591:   }
2592:   PetscFunctionReturn(PETSC_SUCCESS);
2593: }

2595: PetscErrorCode DMGetBasisTransformDM_Internal(DM dm, DM *tdm)
2596: {
2597:   PetscFunctionBegin;
2599:   PetscAssertPointer(tdm, 2);
2600:   *tdm = dm->transformDM;
2601:   PetscFunctionReturn(PETSC_SUCCESS);
2602: }

2604: PetscErrorCode DMGetBasisTransformVec_Internal(DM dm, Vec *tv)
2605: {
2606:   PetscFunctionBegin;
2608:   PetscAssertPointer(tv, 2);
2609:   *tv = dm->transform;
2610:   PetscFunctionReturn(PETSC_SUCCESS);
2611: }

2613: /*@
2614:   DMHasBasisTransform - Whether the `DM` employs a basis transformation from functions in global vectors to functions in local vectors

2616:   Input Parameter:
2617: . dm - The `DM`

2619:   Output Parameter:
2620: . flg - `PETSC_TRUE` if a basis transformation should be done

2622:   Level: developer

2624: .seealso: [](ch_dmbase), `DM`, `DMPlexGlobalToLocalBasis()`, `DMPlexLocalToGlobalBasis()`, `DMPlexCreateBasisRotation()`
2625: @*/
2626: PetscErrorCode DMHasBasisTransform(DM dm, PetscBool *flg)
2627: {
2628:   Vec tv;

2630:   PetscFunctionBegin;
2632:   PetscAssertPointer(flg, 2);
2633:   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
2634:   *flg = tv ? PETSC_TRUE : PETSC_FALSE;
2635:   PetscFunctionReturn(PETSC_SUCCESS);
2636: }

2638: PetscErrorCode DMConstructBasisTransform_Internal(DM dm)
2639: {
2640:   PetscSection s, ts;
2641:   PetscScalar *ta;
2642:   PetscInt     cdim, pStart, pEnd, p, Nf, f, Nc, dof;

2644:   PetscFunctionBegin;
2645:   PetscCall(DMGetCoordinateDim(dm, &cdim));
2646:   PetscCall(DMGetLocalSection(dm, &s));
2647:   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
2648:   PetscCall(PetscSectionGetNumFields(s, &Nf));
2649:   PetscCall(DMClone(dm, &dm->transformDM));
2650:   PetscCall(DMGetLocalSection(dm->transformDM, &ts));
2651:   PetscCall(PetscSectionSetNumFields(ts, Nf));
2652:   PetscCall(PetscSectionSetChart(ts, pStart, pEnd));
2653:   for (f = 0; f < Nf; ++f) {
2654:     PetscCall(PetscSectionGetFieldComponents(s, f, &Nc));
2655:     /* We could start to label fields by their transformation properties */
2656:     if (Nc != cdim) continue;
2657:     for (p = pStart; p < pEnd; ++p) {
2658:       PetscCall(PetscSectionGetFieldDof(s, p, f, &dof));
2659:       if (!dof) continue;
2660:       PetscCall(PetscSectionSetFieldDof(ts, p, f, PetscSqr(cdim)));
2661:       PetscCall(PetscSectionAddDof(ts, p, PetscSqr(cdim)));
2662:     }
2663:   }
2664:   PetscCall(PetscSectionSetUp(ts));
2665:   PetscCall(DMCreateLocalVector(dm->transformDM, &dm->transform));
2666:   PetscCall(VecGetArray(dm->transform, &ta));
2667:   for (p = pStart; p < pEnd; ++p) {
2668:     for (f = 0; f < Nf; ++f) {
2669:       PetscCall(PetscSectionGetFieldDof(ts, p, f, &dof));
2670:       if (dof) {
2671:         PetscReal          x[3] = {0.0, 0.0, 0.0};
2672:         PetscScalar       *tva;
2673:         const PetscScalar *A;

2675:         /* TODO Get quadrature point for this dual basis vector for coordinate */
2676:         PetscCall((*dm->transformGetMatrix)(dm, x, PETSC_TRUE, &A, dm->transformCtx));
2677:         PetscCall(DMPlexPointLocalFieldRef(dm->transformDM, p, f, ta, (void *)&tva));
2678:         PetscCall(PetscArraycpy(tva, A, PetscSqr(cdim)));
2679:       }
2680:     }
2681:   }
2682:   PetscCall(VecRestoreArray(dm->transform, &ta));
2683:   PetscFunctionReturn(PETSC_SUCCESS);
2684: }

2686: PetscErrorCode DMCopyTransform(DM dm, DM newdm)
2687: {
2688:   PetscFunctionBegin;
2691:   newdm->transformCtx       = dm->transformCtx;
2692:   newdm->transformSetUp     = dm->transformSetUp;
2693:   newdm->transformDestroy   = NULL;
2694:   newdm->transformGetMatrix = dm->transformGetMatrix;
2695:   if (newdm->transformSetUp) PetscCall(DMConstructBasisTransform_Internal(newdm));
2696:   PetscFunctionReturn(PETSC_SUCCESS);
2697: }

2699: /*@C
2700:   DMGlobalToLocalHookAdd - adds a callback to be run when `DMGlobalToLocal()` is called

2702:   Logically Collective

2704:   Input Parameters:
2705: + dm        - the `DM`
2706: . beginhook - function to run at the beginning of `DMGlobalToLocalBegin()`
2707: . endhook   - function to run after `DMGlobalToLocalEnd()` has completed
2708: - ctx       - [optional] user-defined context for provide data for the hooks (may be `NULL`)

2710:   Calling sequence of `beginhook`:
2711: + dm   - global `DM`
2712: . g    - global vector
2713: . mode - mode
2714: . l    - local vector
2715: - ctx  - optional user-defined function context

2717:   Calling sequence of `endhook`:
2718: + dm   - global `DM`
2719: . g    - global vector
2720: . mode - mode
2721: . l    - local vector
2722: - ctx  - optional user-defined function context

2724:   Level: advanced

2726:   Note:
2727:   The hook may be used to provide, for example, values that represent boundary conditions in the local vectors that do not exist on the global vector.

2729: .seealso: [](ch_dmbase), `DM`, `DMGlobalToLocal()`, `DMRefineHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`
2730: @*/
2731: PetscErrorCode DMGlobalToLocalHookAdd(DM dm, PetscErrorCode (*beginhook)(DM dm, Vec g, InsertMode mode, Vec l, void *ctx), PetscErrorCode (*endhook)(DM dm, Vec g, InsertMode mode, Vec l, void *ctx), void *ctx)
2732: {
2733:   DMGlobalToLocalHookLink link, *p;

2735:   PetscFunctionBegin;
2737:   for (p = &dm->gtolhook; *p; p = &(*p)->next) { } /* Scan to the end of the current list of hooks */
2738:   PetscCall(PetscNew(&link));
2739:   link->beginhook = beginhook;
2740:   link->endhook   = endhook;
2741:   link->ctx       = ctx;
2742:   link->next      = NULL;
2743:   *p              = link;
2744:   PetscFunctionReturn(PETSC_SUCCESS);
2745: }

2747: static PetscErrorCode DMGlobalToLocalHook_Constraints(DM dm, Vec g, InsertMode mode, Vec l, void *ctx)
2748: {
2749:   Mat          cMat;
2750:   Vec          cVec, cBias;
2751:   PetscSection section, cSec;
2752:   PetscInt     pStart, pEnd, p, dof;

2754:   PetscFunctionBegin;
2755:   (void)g;
2756:   (void)ctx;
2758:   PetscCall(DMGetDefaultConstraints(dm, &cSec, &cMat, &cBias));
2759:   if (cMat && (mode == INSERT_VALUES || mode == INSERT_ALL_VALUES || mode == INSERT_BC_VALUES)) {
2760:     PetscInt nRows;

2762:     PetscCall(MatGetSize(cMat, &nRows, NULL));
2763:     if (nRows <= 0) PetscFunctionReturn(PETSC_SUCCESS);
2764:     PetscCall(DMGetLocalSection(dm, &section));
2765:     PetscCall(MatCreateVecs(cMat, NULL, &cVec));
2766:     PetscCall(MatMult(cMat, l, cVec));
2767:     if (cBias) PetscCall(VecAXPY(cVec, 1., cBias));
2768:     PetscCall(PetscSectionGetChart(cSec, &pStart, &pEnd));
2769:     for (p = pStart; p < pEnd; p++) {
2770:       PetscCall(PetscSectionGetDof(cSec, p, &dof));
2771:       if (dof) {
2772:         PetscScalar *vals;
2773:         PetscCall(VecGetValuesSection(cVec, cSec, p, &vals));
2774:         PetscCall(VecSetValuesSection(l, section, p, vals, INSERT_ALL_VALUES));
2775:       }
2776:     }
2777:     PetscCall(VecDestroy(&cVec));
2778:   }
2779:   PetscFunctionReturn(PETSC_SUCCESS);
2780: }

2782: /*@
2783:   DMGlobalToLocal - update local vectors from global vector

2785:   Neighbor-wise Collective

2787:   Input Parameters:
2788: + dm   - the `DM` object
2789: . g    - the global vector
2790: . mode - `INSERT_VALUES` or `ADD_VALUES`
2791: - l    - the local vector

2793:   Level: beginner

2795:   Notes:
2796:   The communication involved in this update can be overlapped with computation by instead using
2797:   `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()`.

2799:   `DMGlobalToLocalHookAdd()` may be used to provide additional operations that are performed during the update process.

2801: .seealso: [](ch_dmbase), `DM`, `DMGlobalToLocalHookAdd()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`,
2802:           `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToGlobal()`, `DMLocalToGlobalEnd()`,
2803:           `DMGlobalToLocalBegin()` `DMGlobalToLocalEnd()`
2804: @*/
2805: PetscErrorCode DMGlobalToLocal(DM dm, Vec g, InsertMode mode, Vec l)
2806: {
2807:   PetscFunctionBegin;
2808:   PetscCall(DMGlobalToLocalBegin(dm, g, mode, l));
2809:   PetscCall(DMGlobalToLocalEnd(dm, g, mode, l));
2810:   PetscFunctionReturn(PETSC_SUCCESS);
2811: }

2813: /*@
2814:   DMGlobalToLocalBegin - Begins updating local vectors from global vector

2816:   Neighbor-wise Collective

2818:   Input Parameters:
2819: + dm   - the `DM` object
2820: . g    - the global vector
2821: . mode - `INSERT_VALUES` or `ADD_VALUES`
2822: - l    - the local vector

2824:   Level: intermediate

2826:   Notes:
2827:   The operation is completed with `DMGlobalToLocalEnd()`

2829:   One can perform local computations between the `DMGlobalToLocalBegin()` and  `DMGlobalToLocalEnd()` to overlap communication and computation

2831:   `DMGlobalToLocal()` is a short form of  `DMGlobalToLocalBegin()` and  `DMGlobalToLocalEnd()`

2833:   `DMGlobalToLocalHookAdd()` may be used to provide additional operations that are performed during the update process.

2835: .seealso: [](ch_dmbase), `DM`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocal()`, `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToGlobal()`, `DMLocalToGlobalEnd()`
2836: @*/
2837: PetscErrorCode DMGlobalToLocalBegin(DM dm, Vec g, InsertMode mode, Vec l)
2838: {
2839:   PetscSF                 sf;
2840:   DMGlobalToLocalHookLink link;

2842:   PetscFunctionBegin;
2844:   for (link = dm->gtolhook; link; link = link->next) {
2845:     if (link->beginhook) PetscCall((*link->beginhook)(dm, g, mode, l, link->ctx));
2846:   }
2847:   PetscCall(DMGetSectionSF(dm, &sf));
2848:   if (sf) {
2849:     const PetscScalar *gArray;
2850:     PetscScalar       *lArray;
2851:     PetscMemType       lmtype, gmtype;

2853:     PetscCheck(mode != ADD_VALUES, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %d", (int)mode);
2854:     PetscCall(VecGetArrayAndMemType(l, &lArray, &lmtype));
2855:     PetscCall(VecGetArrayReadAndMemType(g, &gArray, &gmtype));
2856:     PetscCall(PetscSFBcastWithMemTypeBegin(sf, MPIU_SCALAR, gmtype, gArray, lmtype, lArray, MPI_REPLACE));
2857:     PetscCall(VecRestoreArrayAndMemType(l, &lArray));
2858:     PetscCall(VecRestoreArrayReadAndMemType(g, &gArray));
2859:   } else {
2860:     PetscUseTypeMethod(dm, globaltolocalbegin, g, mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode), l);
2861:   }
2862:   PetscFunctionReturn(PETSC_SUCCESS);
2863: }

2865: /*@
2866:   DMGlobalToLocalEnd - Ends updating local vectors from global vector

2868:   Neighbor-wise Collective

2870:   Input Parameters:
2871: + dm   - the `DM` object
2872: . g    - the global vector
2873: . mode - `INSERT_VALUES` or `ADD_VALUES`
2874: - l    - the local vector

2876:   Level: intermediate

2878:   Note:
2879:   See `DMGlobalToLocalBegin()` for details.

2881: .seealso: [](ch_dmbase), `DM`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocal()`, `DMLocalToGlobalBegin()`, `DMLocalToGlobal()`, `DMLocalToGlobalEnd()`
2882: @*/
2883: PetscErrorCode DMGlobalToLocalEnd(DM dm, Vec g, InsertMode mode, Vec l)
2884: {
2885:   PetscSF                 sf;
2886:   const PetscScalar      *gArray;
2887:   PetscScalar            *lArray;
2888:   PetscBool               transform;
2889:   DMGlobalToLocalHookLink link;
2890:   PetscMemType            lmtype, gmtype;

2892:   PetscFunctionBegin;
2894:   PetscCall(DMGetSectionSF(dm, &sf));
2895:   PetscCall(DMHasBasisTransform(dm, &transform));
2896:   if (sf) {
2897:     PetscCheck(mode != ADD_VALUES, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %d", (int)mode);

2899:     PetscCall(VecGetArrayAndMemType(l, &lArray, &lmtype));
2900:     PetscCall(VecGetArrayReadAndMemType(g, &gArray, &gmtype));
2901:     PetscCall(PetscSFBcastEnd(sf, MPIU_SCALAR, gArray, lArray, MPI_REPLACE));
2902:     PetscCall(VecRestoreArrayAndMemType(l, &lArray));
2903:     PetscCall(VecRestoreArrayReadAndMemType(g, &gArray));
2904:     if (transform) PetscCall(DMPlexGlobalToLocalBasis(dm, l));
2905:   } else {
2906:     PetscUseTypeMethod(dm, globaltolocalend, g, mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode), l);
2907:   }
2908:   PetscCall(DMGlobalToLocalHook_Constraints(dm, g, mode, l, NULL));
2909:   for (link = dm->gtolhook; link; link = link->next) {
2910:     if (link->endhook) PetscCall((*link->endhook)(dm, g, mode, l, link->ctx));
2911:   }
2912:   PetscFunctionReturn(PETSC_SUCCESS);
2913: }

2915: /*@C
2916:   DMLocalToGlobalHookAdd - adds a callback to be run when a local to global is called

2918:   Logically Collective

2920:   Input Parameters:
2921: + dm        - the `DM`
2922: . beginhook - function to run at the beginning of `DMLocalToGlobalBegin()`
2923: . endhook   - function to run after `DMLocalToGlobalEnd()` has completed
2924: - ctx       - [optional] user-defined context for provide data for the hooks (may be `NULL`)

2926:   Calling sequence of `beginhook`:
2927: + global - global `DM`
2928: . l      - local vector
2929: . mode   - mode
2930: . g      - global vector
2931: - ctx    - optional user-defined function context

2933:   Calling sequence of `endhook`:
2934: + global - global `DM`
2935: . l      - local vector
2936: . mode   - mode
2937: . g      - global vector
2938: - ctx    - optional user-defined function context

2940:   Level: advanced

2942: .seealso: [](ch_dmbase), `DM`, `DMLocalToGlobal()`, `DMRefineHookAdd()`, `DMGlobalToLocalHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`
2943: @*/
2944: PetscErrorCode DMLocalToGlobalHookAdd(DM dm, PetscErrorCode (*beginhook)(DM global, Vec l, InsertMode mode, Vec g, void *ctx), PetscErrorCode (*endhook)(DM global, Vec l, InsertMode mode, Vec g, void *ctx), void *ctx)
2945: {
2946:   DMLocalToGlobalHookLink link, *p;

2948:   PetscFunctionBegin;
2950:   for (p = &dm->ltoghook; *p; p = &(*p)->next) { } /* Scan to the end of the current list of hooks */
2951:   PetscCall(PetscNew(&link));
2952:   link->beginhook = beginhook;
2953:   link->endhook   = endhook;
2954:   link->ctx       = ctx;
2955:   link->next      = NULL;
2956:   *p              = link;
2957:   PetscFunctionReturn(PETSC_SUCCESS);
2958: }

2960: static PetscErrorCode DMLocalToGlobalHook_Constraints(DM dm, Vec l, InsertMode mode, Vec g, void *ctx)
2961: {
2962:   PetscFunctionBegin;
2963:   (void)g;
2964:   (void)ctx;
2966:   if (mode == ADD_VALUES || mode == ADD_ALL_VALUES || mode == ADD_BC_VALUES) {
2967:     Mat          cMat;
2968:     Vec          cVec;
2969:     PetscInt     nRows;
2970:     PetscSection section, cSec;
2971:     PetscInt     pStart, pEnd, p, dof;

2973:     PetscCall(DMGetDefaultConstraints(dm, &cSec, &cMat, NULL));
2974:     if (!cMat) PetscFunctionReturn(PETSC_SUCCESS);

2976:     PetscCall(MatGetSize(cMat, &nRows, NULL));
2977:     if (nRows <= 0) PetscFunctionReturn(PETSC_SUCCESS);
2978:     PetscCall(DMGetLocalSection(dm, &section));
2979:     PetscCall(MatCreateVecs(cMat, NULL, &cVec));
2980:     PetscCall(PetscSectionGetChart(cSec, &pStart, &pEnd));
2981:     for (p = pStart; p < pEnd; p++) {
2982:       PetscCall(PetscSectionGetDof(cSec, p, &dof));
2983:       if (dof) {
2984:         PetscInt     d;
2985:         PetscScalar *vals;
2986:         PetscCall(VecGetValuesSection(l, section, p, &vals));
2987:         PetscCall(VecSetValuesSection(cVec, cSec, p, vals, mode));
2988:         /* for this to be the true transpose, we have to zero the values that
2989:          * we just extracted */
2990:         for (d = 0; d < dof; d++) vals[d] = 0.;
2991:       }
2992:     }
2993:     PetscCall(MatMultTransposeAdd(cMat, cVec, l, l));
2994:     PetscCall(VecDestroy(&cVec));
2995:   }
2996:   PetscFunctionReturn(PETSC_SUCCESS);
2997: }
2998: /*@
2999:   DMLocalToGlobal - updates global vectors from local vectors

3001:   Neighbor-wise Collective

3003:   Input Parameters:
3004: + dm   - the `DM` object
3005: . l    - the local vector
3006: . mode - if `INSERT_VALUES` then no parallel communication is used, if `ADD_VALUES` then all ghost points from the same base point accumulate into that base point.
3007: - g    - the global vector

3009:   Level: beginner

3011:   Notes:
3012:   The communication involved in this update can be overlapped with computation by using
3013:   `DMLocalToGlobalBegin()` and `DMLocalToGlobalEnd()`.

3015:   In the `ADD_VALUES` case you normally would zero the receiving vector before beginning this operation.

3017:   `INSERT_VALUES` is not supported for `DMDA`; in that case simply compute the values directly into a global vector instead of a local one.

3019:   Use `DMLocalToGlobalHookAdd()` to add additional operations that are performed on the data during the update process

3021: .seealso: [](ch_dmbase), `DM`, `DMLocalToGlobalBegin()`, `DMLocalToGlobalEnd()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocal()`, `DMGlobalToLocalEnd()`, `DMGlobalToLocalBegin()`, `DMLocalToGlobalHookAdd()`, `DMGlobaToLocallHookAdd()`
3022: @*/
3023: PetscErrorCode DMLocalToGlobal(DM dm, Vec l, InsertMode mode, Vec g)
3024: {
3025:   PetscFunctionBegin;
3026:   PetscCall(DMLocalToGlobalBegin(dm, l, mode, g));
3027:   PetscCall(DMLocalToGlobalEnd(dm, l, mode, g));
3028:   PetscFunctionReturn(PETSC_SUCCESS);
3029: }

3031: /*@
3032:   DMLocalToGlobalBegin - begins updating global vectors from local vectors

3034:   Neighbor-wise Collective

3036:   Input Parameters:
3037: + dm   - the `DM` object
3038: . l    - the local vector
3039: . mode - if `INSERT_VALUES` then no parallel communication is used, if `ADD_VALUES` then all ghost points from the same base point accumulate into that base point.
3040: - g    - the global vector

3042:   Level: intermediate

3044:   Notes:
3045:   In the `ADD_VALUES` case you normally would zero the receiving vector before beginning this operation.

3047:   `INSERT_VALUES is` not supported for `DMDA`, in that case simply compute the values directly into a global vector instead of a local one.

3049:   Use `DMLocalToGlobalEnd()` to complete the communication process.

3051:   `DMLocalToGlobal()` is a short form of  `DMLocalToGlobalBegin()` and  `DMLocalToGlobalEnd()`

3053:   `DMLocalToGlobalHookAdd()` may be used to provide additional operations that are performed during the update process.

3055: .seealso: [](ch_dmbase), `DM`, `DMLocalToGlobal()`, `DMLocalToGlobalEnd()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocal()`, `DMGlobalToLocalEnd()`, `DMGlobalToLocalBegin()`
3056: @*/
3057: PetscErrorCode DMLocalToGlobalBegin(DM dm, Vec l, InsertMode mode, Vec g)
3058: {
3059:   PetscSF                 sf;
3060:   PetscSection            s, gs;
3061:   DMLocalToGlobalHookLink link;
3062:   Vec                     tmpl;
3063:   const PetscScalar      *lArray;
3064:   PetscScalar            *gArray;
3065:   PetscBool               isInsert, transform, l_inplace = PETSC_FALSE, g_inplace = PETSC_FALSE;
3066:   PetscMemType            lmtype = PETSC_MEMTYPE_HOST, gmtype = PETSC_MEMTYPE_HOST;

3068:   PetscFunctionBegin;
3070:   for (link = dm->ltoghook; link; link = link->next) {
3071:     if (link->beginhook) PetscCall((*link->beginhook)(dm, l, mode, g, link->ctx));
3072:   }
3073:   PetscCall(DMLocalToGlobalHook_Constraints(dm, l, mode, g, NULL));
3074:   PetscCall(DMGetSectionSF(dm, &sf));
3075:   PetscCall(DMGetLocalSection(dm, &s));
3076:   switch (mode) {
3077:   case INSERT_VALUES:
3078:   case INSERT_ALL_VALUES:
3079:   case INSERT_BC_VALUES:
3080:     isInsert = PETSC_TRUE;
3081:     break;
3082:   case ADD_VALUES:
3083:   case ADD_ALL_VALUES:
3084:   case ADD_BC_VALUES:
3085:     isInsert = PETSC_FALSE;
3086:     break;
3087:   default:
3088:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %d", mode);
3089:   }
3090:   if ((sf && !isInsert) || (s && isInsert)) {
3091:     PetscCall(DMHasBasisTransform(dm, &transform));
3092:     if (transform) {
3093:       PetscCall(DMGetNamedLocalVector(dm, "__petsc_dm_transform_local_copy", &tmpl));
3094:       PetscCall(VecCopy(l, tmpl));
3095:       PetscCall(DMPlexLocalToGlobalBasis(dm, tmpl));
3096:       PetscCall(VecGetArrayRead(tmpl, &lArray));
3097:     } else if (isInsert) {
3098:       PetscCall(VecGetArrayRead(l, &lArray));
3099:     } else {
3100:       PetscCall(VecGetArrayReadAndMemType(l, &lArray, &lmtype));
3101:       l_inplace = PETSC_TRUE;
3102:     }
3103:     if (s && isInsert) {
3104:       PetscCall(VecGetArray(g, &gArray));
3105:     } else {
3106:       PetscCall(VecGetArrayAndMemType(g, &gArray, &gmtype));
3107:       g_inplace = PETSC_TRUE;
3108:     }
3109:     if (sf && !isInsert) {
3110:       PetscCall(PetscSFReduceWithMemTypeBegin(sf, MPIU_SCALAR, lmtype, lArray, gmtype, gArray, MPIU_SUM));
3111:     } else if (s && isInsert) {
3112:       PetscInt gStart, pStart, pEnd, p;

3114:       PetscCall(DMGetGlobalSection(dm, &gs));
3115:       PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
3116:       PetscCall(VecGetOwnershipRange(g, &gStart, NULL));
3117:       for (p = pStart; p < pEnd; ++p) {
3118:         PetscInt dof, gdof, cdof, gcdof, off, goff, d, e;

3120:         PetscCall(PetscSectionGetDof(s, p, &dof));
3121:         PetscCall(PetscSectionGetDof(gs, p, &gdof));
3122:         PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
3123:         PetscCall(PetscSectionGetConstraintDof(gs, p, &gcdof));
3124:         PetscCall(PetscSectionGetOffset(s, p, &off));
3125:         PetscCall(PetscSectionGetOffset(gs, p, &goff));
3126:         /* Ignore off-process data and points with no global data */
3127:         if (!gdof || goff < 0) continue;
3128:         PetscCheck(dof == gdof, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Inconsistent sizes, p: %" PetscInt_FMT " dof: %" PetscInt_FMT " gdof: %" PetscInt_FMT " cdof: %" PetscInt_FMT " gcdof: %" PetscInt_FMT, p, dof, gdof, cdof, gcdof);
3129:         /* If no constraints are enforced in the global vector */
3130:         if (!gcdof) {
3131:           for (d = 0; d < dof; ++d) gArray[goff - gStart + d] = lArray[off + d];
3132:           /* If constraints are enforced in the global vector */
3133:         } else if (cdof == gcdof) {
3134:           const PetscInt *cdofs;
3135:           PetscInt        cind = 0;

3137:           PetscCall(PetscSectionGetConstraintIndices(s, p, &cdofs));
3138:           for (d = 0, e = 0; d < dof; ++d) {
3139:             if ((cind < cdof) && (d == cdofs[cind])) {
3140:               ++cind;
3141:               continue;
3142:             }
3143:             gArray[goff - gStart + e++] = lArray[off + d];
3144:           }
3145:         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Inconsistent sizes, p: %" PetscInt_FMT " dof: %" PetscInt_FMT " gdof: %" PetscInt_FMT " cdof: %" PetscInt_FMT " gcdof: %" PetscInt_FMT, p, dof, gdof, cdof, gcdof);
3146:       }
3147:     }
3148:     if (g_inplace) {
3149:       PetscCall(VecRestoreArrayAndMemType(g, &gArray));
3150:     } else {
3151:       PetscCall(VecRestoreArray(g, &gArray));
3152:     }
3153:     if (transform) {
3154:       PetscCall(VecRestoreArrayRead(tmpl, &lArray));
3155:       PetscCall(DMRestoreNamedLocalVector(dm, "__petsc_dm_transform_local_copy", &tmpl));
3156:     } else if (l_inplace) {
3157:       PetscCall(VecRestoreArrayReadAndMemType(l, &lArray));
3158:     } else {
3159:       PetscCall(VecRestoreArrayRead(l, &lArray));
3160:     }
3161:   } else {
3162:     PetscUseTypeMethod(dm, localtoglobalbegin, l, mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode), g);
3163:   }
3164:   PetscFunctionReturn(PETSC_SUCCESS);
3165: }

3167: /*@
3168:   DMLocalToGlobalEnd - updates global vectors from local vectors

3170:   Neighbor-wise Collective

3172:   Input Parameters:
3173: + dm   - the `DM` object
3174: . l    - the local vector
3175: . mode - `INSERT_VALUES` or `ADD_VALUES`
3176: - g    - the global vector

3178:   Level: intermediate

3180:   Note:
3181:   See `DMLocalToGlobalBegin()` for full details

3183: .seealso: [](ch_dmbase), `DM`, `DMLocalToGlobalBegin()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocalEnd()`
3184: @*/
3185: PetscErrorCode DMLocalToGlobalEnd(DM dm, Vec l, InsertMode mode, Vec g)
3186: {
3187:   PetscSF                 sf;
3188:   PetscSection            s;
3189:   DMLocalToGlobalHookLink link;
3190:   PetscBool               isInsert, transform;

3192:   PetscFunctionBegin;
3194:   PetscCall(DMGetSectionSF(dm, &sf));
3195:   PetscCall(DMGetLocalSection(dm, &s));
3196:   switch (mode) {
3197:   case INSERT_VALUES:
3198:   case INSERT_ALL_VALUES:
3199:     isInsert = PETSC_TRUE;
3200:     break;
3201:   case ADD_VALUES:
3202:   case ADD_ALL_VALUES:
3203:     isInsert = PETSC_FALSE;
3204:     break;
3205:   default:
3206:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %d", mode);
3207:   }
3208:   if (sf && !isInsert) {
3209:     const PetscScalar *lArray;
3210:     PetscScalar       *gArray;
3211:     Vec                tmpl;

3213:     PetscCall(DMHasBasisTransform(dm, &transform));
3214:     if (transform) {
3215:       PetscCall(DMGetNamedLocalVector(dm, "__petsc_dm_transform_local_copy", &tmpl));
3216:       PetscCall(VecGetArrayRead(tmpl, &lArray));
3217:     } else {
3218:       PetscCall(VecGetArrayReadAndMemType(l, &lArray, NULL));
3219:     }
3220:     PetscCall(VecGetArrayAndMemType(g, &gArray, NULL));
3221:     PetscCall(PetscSFReduceEnd(sf, MPIU_SCALAR, lArray, gArray, MPIU_SUM));
3222:     if (transform) {
3223:       PetscCall(VecRestoreArrayRead(tmpl, &lArray));
3224:       PetscCall(DMRestoreNamedLocalVector(dm, "__petsc_dm_transform_local_copy", &tmpl));
3225:     } else {
3226:       PetscCall(VecRestoreArrayReadAndMemType(l, &lArray));
3227:     }
3228:     PetscCall(VecRestoreArrayAndMemType(g, &gArray));
3229:   } else if (s && isInsert) {
3230:   } else {
3231:     PetscUseTypeMethod(dm, localtoglobalend, l, mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode), g);
3232:   }
3233:   for (link = dm->ltoghook; link; link = link->next) {
3234:     if (link->endhook) PetscCall((*link->endhook)(dm, g, mode, l, link->ctx));
3235:   }
3236:   PetscFunctionReturn(PETSC_SUCCESS);
3237: }

3239: /*@
3240:   DMLocalToLocalBegin - Begins the process of mapping values from a local vector (that include
3241:   ghost points that contain irrelevant values) to another local vector where the ghost points
3242:   in the second are set correctly from values on other MPI ranks.

3244:   Neighbor-wise Collective

3246:   Input Parameters:
3247: + dm   - the `DM` object
3248: . g    - the original local vector
3249: - mode - one of `INSERT_VALUES` or `ADD_VALUES`

3251:   Output Parameter:
3252: . l - the local vector with correct ghost values

3254:   Level: intermediate

3256:   Note:
3257:   Must be followed by `DMLocalToLocalEnd()`.

3259: .seealso: [](ch_dmbase), `DM`, `DMLocalToLocalEnd()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateLocalVector()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`
3260: @*/
3261: PetscErrorCode DMLocalToLocalBegin(DM dm, Vec g, InsertMode mode, Vec l)
3262: {
3263:   PetscFunctionBegin;
3267:   PetscUseTypeMethod(dm, localtolocalbegin, g, mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode), l);
3268:   PetscFunctionReturn(PETSC_SUCCESS);
3269: }

3271: /*@
3272:   DMLocalToLocalEnd - Maps from a local vector to another local vector where the ghost
3273:   points in the second are set correctly. Must be preceded by `DMLocalToLocalBegin()`.

3275:   Neighbor-wise Collective

3277:   Input Parameters:
3278: + dm   - the `DM` object
3279: . g    - the original local vector
3280: - mode - one of `INSERT_VALUES` or `ADD_VALUES`

3282:   Output Parameter:
3283: . l - the local vector with correct ghost values

3285:   Level: intermediate

3287: .seealso: [](ch_dmbase), `DM`, `DMLocalToLocalBegin()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateLocalVector()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`
3288: @*/
3289: PetscErrorCode DMLocalToLocalEnd(DM dm, Vec g, InsertMode mode, Vec l)
3290: {
3291:   PetscFunctionBegin;
3295:   PetscUseTypeMethod(dm, localtolocalend, g, mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode), l);
3296:   PetscFunctionReturn(PETSC_SUCCESS);
3297: }

3299: /*@
3300:   DMCoarsen - Coarsens a `DM` object using a standard, non-adaptive coarsening of the underlying mesh

3302:   Collective

3304:   Input Parameters:
3305: + dm   - the `DM` object
3306: - comm - the communicator to contain the new `DM` object (or `MPI_COMM_NULL`)

3308:   Output Parameter:
3309: . dmc - the coarsened `DM`

3311:   Level: developer

3313: .seealso: [](ch_dmbase), `DM`, `DMRefine()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateDomainDecomposition()`,
3314:           `DMCoarsenHookAdd()`, `DMCoarsenHookRemove()`
3315: @*/
3316: PetscErrorCode DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
3317: {
3318:   DMCoarsenHookLink link;

3320:   PetscFunctionBegin;
3322:   PetscCall(PetscLogEventBegin(DM_Coarsen, dm, 0, 0, 0));
3323:   PetscUseTypeMethod(dm, coarsen, comm, dmc);
3324:   if (*dmc) {
3325:     (*dmc)->bind_below = dm->bind_below; /* Propagate this from parent DM; otherwise -dm_bind_below will be useless for multigrid cases. */
3326:     PetscCall(DMSetCoarseDM(dm, *dmc));
3327:     (*dmc)->ops->creatematrix = dm->ops->creatematrix;
3328:     PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)dm, (PetscObject)*dmc));
3329:     (*dmc)->ctx       = dm->ctx;
3330:     (*dmc)->levelup   = dm->levelup;
3331:     (*dmc)->leveldown = dm->leveldown + 1;
3332:     PetscCall(DMSetMatType(*dmc, dm->mattype));
3333:     for (link = dm->coarsenhook; link; link = link->next) {
3334:       if (link->coarsenhook) PetscCall((*link->coarsenhook)(dm, *dmc, link->ctx));
3335:     }
3336:   }
3337:   PetscCall(PetscLogEventEnd(DM_Coarsen, dm, 0, 0, 0));
3338:   PetscCheck(*dmc, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "NULL coarse mesh produced");
3339:   PetscFunctionReturn(PETSC_SUCCESS);
3340: }

3342: /*@C
3343:   DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid

3345:   Logically Collective; No Fortran Support

3347:   Input Parameters:
3348: + fine         - `DM` on which to run a hook when restricting to a coarser level
3349: . coarsenhook  - function to run when setting up a coarser level
3350: . restricthook - function to run to update data on coarser levels (called once per `SNESSolve()`)
3351: - ctx          - [optional] user-defined context for provide data for the hooks (may be `NULL`)

3353:   Calling sequence of `coarsenhook`:
3354: + fine   - fine level `DM`
3355: . coarse - coarse level `DM` to restrict problem to
3356: - ctx    - optional user-defined function context

3358:   Calling sequence of `restricthook`:
3359: + fine      - fine level `DM`
3360: . mrestrict - matrix restricting a fine-level solution to the coarse grid, usually the transpose of the interpolation
3361: . rscale    - scaling vector for restriction
3362: . inject    - matrix restricting by injection
3363: . coarse    - coarse level DM to update
3364: - ctx       - optional user-defined function context

3366:   Level: advanced

3368:   Notes:
3369:   This function is only needed if auxiliary data, attached to the `DM` with `PetscObjectCompose()`, needs to be set up or passed from the fine `DM` to the coarse `DM`.

3371:   If this function is called multiple times, the hooks will be run in the order they are added.

3373:   In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
3374:   extract the finest level information from its context (instead of from the `SNES`).

3376:   The hooks are automatically called by `DMRestrict()`

3378: .seealso: [](ch_dmbase), `DM`, `DMCoarsenHookRemove()`, `DMRefineHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`
3379: @*/
3380: PetscErrorCode DMCoarsenHookAdd(DM fine, PetscErrorCode (*coarsenhook)(DM fine, DM coarse, void *ctx), PetscErrorCode (*restricthook)(DM fine, Mat mrestrict, Vec rscale, Mat inject, DM coarse, void *ctx), void *ctx)
3381: {
3382:   DMCoarsenHookLink link, *p;

3384:   PetscFunctionBegin;
3386:   for (p = &fine->coarsenhook; *p; p = &(*p)->next) { /* Scan to the end of the current list of hooks */
3387:     if ((*p)->coarsenhook == coarsenhook && (*p)->restricthook == restricthook && (*p)->ctx == ctx) PetscFunctionReturn(PETSC_SUCCESS);
3388:   }
3389:   PetscCall(PetscNew(&link));
3390:   link->coarsenhook  = coarsenhook;
3391:   link->restricthook = restricthook;
3392:   link->ctx          = ctx;
3393:   link->next         = NULL;
3394:   *p                 = link;
3395:   PetscFunctionReturn(PETSC_SUCCESS);
3396: }

3398: /*@C
3399:   DMCoarsenHookRemove - remove a callback set with `DMCoarsenHookAdd()`

3401:   Logically Collective; No Fortran Support

3403:   Input Parameters:
3404: + fine         - `DM` on which to run a hook when restricting to a coarser level
3405: . coarsenhook  - function to run when setting up a coarser level
3406: . restricthook - function to run to update data on coarser levels
3407: - ctx          - [optional] user-defined context for provide data for the hooks (may be `NULL`)

3409:   Level: advanced

3411:   Notes:
3412:   This function does nothing if the `coarsenhook` is not in the list.

3414:   See `DMCoarsenHookAdd()` for the calling sequence of `coarsenhook` and `restricthook`

3416: .seealso: [](ch_dmbase), `DM`, `DMCoarsenHookAdd()`, `DMRefineHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`
3417: @*/
3418: PetscErrorCode DMCoarsenHookRemove(DM fine, PetscErrorCode (*coarsenhook)(DM, DM, void *), PetscErrorCode (*restricthook)(DM, Mat, Vec, Mat, DM, void *), void *ctx)
3419: {
3420:   DMCoarsenHookLink link, *p;

3422:   PetscFunctionBegin;
3424:   for (p = &fine->coarsenhook; *p; p = &(*p)->next) { /* Search the list of current hooks */
3425:     if ((*p)->coarsenhook == coarsenhook && (*p)->restricthook == restricthook && (*p)->ctx == ctx) {
3426:       link = *p;
3427:       *p   = link->next;
3428:       PetscCall(PetscFree(link));
3429:       break;
3430:     }
3431:   }
3432:   PetscFunctionReturn(PETSC_SUCCESS);
3433: }

3435: /*@
3436:   DMRestrict - restricts user-defined problem data to a coarser `DM` by running hooks registered by `DMCoarsenHookAdd()`

3438:   Collective if any hooks are

3440:   Input Parameters:
3441: + fine    - finer `DM` from which the data is obtained
3442: . restrct - restriction matrix, apply using `MatRestrict()`, usually the transpose of the interpolation
3443: . rscale  - scaling vector for restriction
3444: . inject  - injection matrix, also use `MatRestrict()`
3445: - coarse  - coarser `DM` to update

3447:   Level: developer

3449:   Developer Note:
3450:   Though this routine is called `DMRestrict()` the hooks are added with `DMCoarsenHookAdd()`, a consistent terminology would be better

3452: .seealso: [](ch_dmbase), `DM`, `DMCoarsenHookAdd()`, `MatRestrict()`, `DMInterpolate()`, `DMRefineHookAdd()`
3453: @*/
3454: PetscErrorCode DMRestrict(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse)
3455: {
3456:   DMCoarsenHookLink link;

3458:   PetscFunctionBegin;
3459:   for (link = fine->coarsenhook; link; link = link->next) {
3460:     if (link->restricthook) PetscCall((*link->restricthook)(fine, restrct, rscale, inject, coarse, link->ctx));
3461:   }
3462:   PetscFunctionReturn(PETSC_SUCCESS);
3463: }

3465: /*@C
3466:   DMSubDomainHookAdd - adds a callback to be run when restricting a problem to subdomain `DM`s with `DMCreateDomainDecomposition()`

3468:   Logically Collective; No Fortran Support

3470:   Input Parameters:
3471: + global       - global `DM`
3472: . ddhook       - function to run to pass data to the decomposition `DM` upon its creation
3473: . restricthook - function to run to update data on block solve (at the beginning of the block solve)
3474: - ctx          - [optional] user-defined context for provide data for the hooks (may be `NULL`)

3476:   Calling sequence of `ddhook`:
3477: + global - global `DM`
3478: . block  - subdomain `DM`
3479: - ctx    - optional user-defined function context

3481:   Calling sequence of `restricthook`:
3482: + global - global `DM`
3483: . out    - scatter to the outer (with ghost and overlap points) sub vector
3484: . in     - scatter to sub vector values only owned locally
3485: . block  - subdomain `DM`
3486: - ctx    - optional user-defined function context

3488:   Level: advanced

3490:   Notes:
3491:   This function can be used if auxiliary data needs to be set up on subdomain `DM`s.

3493:   If this function is called multiple times, the hooks will be run in the order they are added.

3495:   In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
3496:   extract the global information from its context (instead of from the `SNES`).

3498:   Developer Note:
3499:   It is unclear what "block solve" means within the definition of `restricthook`

3501: .seealso: [](ch_dmbase), `DM`, `DMSubDomainHookRemove()`, `DMRefineHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`, `DMCreateDomainDecomposition()`
3502: @*/
3503: PetscErrorCode DMSubDomainHookAdd(DM global, PetscErrorCode (*ddhook)(DM global, DM block, void *ctx), PetscErrorCode (*restricthook)(DM global, VecScatter out, VecScatter in, DM block, void *ctx), void *ctx)
3504: {
3505:   DMSubDomainHookLink link, *p;

3507:   PetscFunctionBegin;
3509:   for (p = &global->subdomainhook; *p; p = &(*p)->next) { /* Scan to the end of the current list of hooks */
3510:     if ((*p)->ddhook == ddhook && (*p)->restricthook == restricthook && (*p)->ctx == ctx) PetscFunctionReturn(PETSC_SUCCESS);
3511:   }
3512:   PetscCall(PetscNew(&link));
3513:   link->restricthook = restricthook;
3514:   link->ddhook       = ddhook;
3515:   link->ctx          = ctx;
3516:   link->next         = NULL;
3517:   *p                 = link;
3518:   PetscFunctionReturn(PETSC_SUCCESS);
3519: }

3521: /*@C
3522:   DMSubDomainHookRemove - remove a callback from the list to be run when restricting a problem to subdomain `DM`s with `DMCreateDomainDecomposition()`

3524:   Logically Collective; No Fortran Support

3526:   Input Parameters:
3527: + global       - global `DM`
3528: . ddhook       - function to run to pass data to the decomposition `DM` upon its creation
3529: . restricthook - function to run to update data on block solve (at the beginning of the block solve)
3530: - ctx          - [optional] user-defined context for provide data for the hooks (may be `NULL`)

3532:   Level: advanced

3534:   Note:
3535:   See `DMSubDomainHookAdd()` for the calling sequences of `ddhook` and `restricthook`

3537: .seealso: [](ch_dmbase), `DM`, `DMSubDomainHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`,
3538:           `DMCreateDomainDecomposition()`
3539: @*/
3540: PetscErrorCode DMSubDomainHookRemove(DM global, PetscErrorCode (*ddhook)(DM, DM, void *), PetscErrorCode (*restricthook)(DM, VecScatter, VecScatter, DM, void *), void *ctx)
3541: {
3542:   DMSubDomainHookLink link, *p;

3544:   PetscFunctionBegin;
3546:   for (p = &global->subdomainhook; *p; p = &(*p)->next) { /* Search the list of current hooks */
3547:     if ((*p)->ddhook == ddhook && (*p)->restricthook == restricthook && (*p)->ctx == ctx) {
3548:       link = *p;
3549:       *p   = link->next;
3550:       PetscCall(PetscFree(link));
3551:       break;
3552:     }
3553:   }
3554:   PetscFunctionReturn(PETSC_SUCCESS);
3555: }

3557: /*@
3558:   DMSubDomainRestrict - restricts user-defined problem data to a subdomain `DM` by running hooks registered by `DMSubDomainHookAdd()`

3560:   Collective if any hooks are

3562:   Input Parameters:
3563: + global   - The global `DM` to use as a base
3564: . oscatter - The scatter from domain global vector filling subdomain global vector with overlap
3565: . gscatter - The scatter from domain global vector filling subdomain local vector with ghosts
3566: - subdm    - The subdomain `DM` to update

3568:   Level: developer

3570: .seealso: [](ch_dmbase), `DM`, `DMCoarsenHookAdd()`, `MatRestrict()`, `DMCreateDomainDecomposition()`
3571: @*/
3572: PetscErrorCode DMSubDomainRestrict(DM global, VecScatter oscatter, VecScatter gscatter, DM subdm)
3573: {
3574:   DMSubDomainHookLink link;

3576:   PetscFunctionBegin;
3577:   for (link = global->subdomainhook; link; link = link->next) {
3578:     if (link->restricthook) PetscCall((*link->restricthook)(global, oscatter, gscatter, subdm, link->ctx));
3579:   }
3580:   PetscFunctionReturn(PETSC_SUCCESS);
3581: }

3583: /*@
3584:   DMGetCoarsenLevel - Gets the number of coarsenings that have generated this `DM`.

3586:   Not Collective

3588:   Input Parameter:
3589: . dm - the `DM` object

3591:   Output Parameter:
3592: . level - number of coarsenings

3594:   Level: developer

3596: .seealso: [](ch_dmbase), `DM`, `DMCoarsen()`, `DMSetCoarsenLevel()`, `DMGetRefineLevel()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`
3597: @*/
3598: PetscErrorCode DMGetCoarsenLevel(DM dm, PetscInt *level)
3599: {
3600:   PetscFunctionBegin;
3602:   PetscAssertPointer(level, 2);
3603:   *level = dm->leveldown;
3604:   PetscFunctionReturn(PETSC_SUCCESS);
3605: }

3607: /*@
3608:   DMSetCoarsenLevel - Sets the number of coarsenings that have generated this `DM`.

3610:   Collective

3612:   Input Parameters:
3613: + dm    - the `DM` object
3614: - level - number of coarsenings

3616:   Level: developer

3618:   Note:
3619:   This is rarely used directly, the information is automatically set when a `DM` is created with `DMCoarsen()`

3621: .seealso: [](ch_dmbase), `DM`, `DMCoarsen()`, `DMGetCoarsenLevel()`, `DMGetRefineLevel()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`
3622: @*/
3623: PetscErrorCode DMSetCoarsenLevel(DM dm, PetscInt level)
3624: {
3625:   PetscFunctionBegin;
3627:   dm->leveldown = level;
3628:   PetscFunctionReturn(PETSC_SUCCESS);
3629: }

3631: /*@
3632:   DMRefineHierarchy - Refines a `DM` object, all levels at once

3634:   Collective

3636:   Input Parameters:
3637: + dm      - the `DM` object
3638: - nlevels - the number of levels of refinement

3640:   Output Parameter:
3641: . dmf - the refined `DM` hierarchy

3643:   Level: developer

3645: .seealso: [](ch_dmbase), `DM`, `DMCoarsen()`, `DMCoarsenHierarchy()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`
3646: @*/
3647: PetscErrorCode DMRefineHierarchy(DM dm, PetscInt nlevels, DM dmf[])
3648: {
3649:   PetscFunctionBegin;
3651:   PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative");
3652:   if (nlevels == 0) PetscFunctionReturn(PETSC_SUCCESS);
3653:   PetscAssertPointer(dmf, 3);
3654:   if (dm->ops->refine && !dm->ops->refinehierarchy) {
3655:     PetscInt i;

3657:     PetscCall(DMRefine(dm, PetscObjectComm((PetscObject)dm), &dmf[0]));
3658:     for (i = 1; i < nlevels; i++) PetscCall(DMRefine(dmf[i - 1], PetscObjectComm((PetscObject)dm), &dmf[i]));
3659:   } else PetscUseTypeMethod(dm, refinehierarchy, nlevels, dmf);
3660:   PetscFunctionReturn(PETSC_SUCCESS);
3661: }

3663: /*@
3664:   DMCoarsenHierarchy - Coarsens a `DM` object, all levels at once

3666:   Collective

3668:   Input Parameters:
3669: + dm      - the `DM` object
3670: - nlevels - the number of levels of coarsening

3672:   Output Parameter:
3673: . dmc - the coarsened `DM` hierarchy

3675:   Level: developer

3677: .seealso: [](ch_dmbase), `DM`, `DMCoarsen()`, `DMRefineHierarchy()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`
3678: @*/
3679: PetscErrorCode DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
3680: {
3681:   PetscFunctionBegin;
3683:   PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative");
3684:   if (nlevels == 0) PetscFunctionReturn(PETSC_SUCCESS);
3685:   PetscAssertPointer(dmc, 3);
3686:   if (dm->ops->coarsen && !dm->ops->coarsenhierarchy) {
3687:     PetscInt i;

3689:     PetscCall(DMCoarsen(dm, PetscObjectComm((PetscObject)dm), &dmc[0]));
3690:     for (i = 1; i < nlevels; i++) PetscCall(DMCoarsen(dmc[i - 1], PetscObjectComm((PetscObject)dm), &dmc[i]));
3691:   } else PetscUseTypeMethod(dm, coarsenhierarchy, nlevels, dmc);
3692:   PetscFunctionReturn(PETSC_SUCCESS);
3693: }

3695: /*@C
3696:   DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the `DM` is destroyed

3698:   Logically Collective if the function is collective

3700:   Input Parameters:
3701: + dm      - the `DM` object
3702: - destroy - the destroy function, see `PetscCtxDestroyFn` for the calling sequence

3704:   Level: intermediate

3706: .seealso: [](ch_dmbase), `DM`, `DMSetApplicationContext()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`,
3707:           `DMGetApplicationContext()`, `PetscCtxDestroyFn`
3708: @*/
3709: PetscErrorCode DMSetApplicationContextDestroy(DM dm, PetscCtxDestroyFn *destroy)
3710: {
3711:   PetscFunctionBegin;
3713:   dm->ctxdestroy = destroy;
3714:   PetscFunctionReturn(PETSC_SUCCESS);
3715: }

3717: /*@
3718:   DMSetApplicationContext - Set a user context into a `DM` object

3720:   Not Collective

3722:   Input Parameters:
3723: + dm  - the `DM` object
3724: - ctx - the user context

3726:   Level: intermediate

3728:   Notes:
3729:   A user context is a way to pass problem specific information that is accessible whenever the `DM` is available
3730:   In a multilevel solver, the user context is shared by all the `DM` in the hierarchy; it is thus not advisable
3731:   to store objects that represent discretized quantities inside the context.

3733: .seealso: [](ch_dmbase), `DM`, `DMGetApplicationContext()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`
3734: @*/
3735: PetscErrorCode DMSetApplicationContext(DM dm, void *ctx)
3736: {
3737:   PetscFunctionBegin;
3739:   dm->ctx = ctx;
3740:   PetscFunctionReturn(PETSC_SUCCESS);
3741: }

3743: /*@
3744:   DMGetApplicationContext - Gets a user context from a `DM` object

3746:   Not Collective

3748:   Input Parameter:
3749: . dm - the `DM` object

3751:   Output Parameter:
3752: . ctx - the user context

3754:   Level: intermediate

3756:   Note:
3757:   A user context is a way to pass problem specific information that is accessible whenever the `DM` is available

3759: .seealso: [](ch_dmbase), `DM`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`
3760: @*/
3761: PetscErrorCode DMGetApplicationContext(DM dm, void *ctx)
3762: {
3763:   PetscFunctionBegin;
3765:   *(void **)ctx = dm->ctx;
3766:   PetscFunctionReturn(PETSC_SUCCESS);
3767: }

3769: /*@C
3770:   DMSetVariableBounds - sets a function to compute the lower and upper bound vectors for `SNESVI`.

3772:   Logically Collective

3774:   Input Parameters:
3775: + dm - the DM object
3776: - f  - the function that computes variable bounds used by SNESVI (use `NULL` to cancel a previous function that was set)

3778:   Level: intermediate

3780: .seealso: [](ch_dmbase), `DM`, `DMComputeVariableBounds()`, `DMHasVariableBounds()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMGetApplicationContext()`,
3781:          `DMSetJacobian()`
3782: @*/
3783: PetscErrorCode DMSetVariableBounds(DM dm, PetscErrorCode (*f)(DM, Vec, Vec))
3784: {
3785:   PetscFunctionBegin;
3787:   dm->ops->computevariablebounds = f;
3788:   PetscFunctionReturn(PETSC_SUCCESS);
3789: }

3791: /*@
3792:   DMHasVariableBounds - does the `DM` object have a variable bounds function?

3794:   Not Collective

3796:   Input Parameter:
3797: . dm - the `DM` object to destroy

3799:   Output Parameter:
3800: . flg - `PETSC_TRUE` if the variable bounds function exists

3802:   Level: developer

3804: .seealso: [](ch_dmbase), `DM`, `DMComputeVariableBounds()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMGetApplicationContext()`
3805: @*/
3806: PetscErrorCode DMHasVariableBounds(DM dm, PetscBool *flg)
3807: {
3808:   PetscFunctionBegin;
3810:   PetscAssertPointer(flg, 2);
3811:   *flg = (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE;
3812:   PetscFunctionReturn(PETSC_SUCCESS);
3813: }

3815: /*@
3816:   DMComputeVariableBounds - compute variable bounds used by `SNESVI`.

3818:   Logically Collective

3820:   Input Parameter:
3821: . dm - the `DM` object

3823:   Output Parameters:
3824: + xl - lower bound
3825: - xu - upper bound

3827:   Level: advanced

3829:   Note:
3830:   This is generally not called by users. It calls the function provided by the user with DMSetVariableBounds()

3832: .seealso: [](ch_dmbase), `DM`, `DMHasVariableBounds()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMGetApplicationContext()`
3833: @*/
3834: PetscErrorCode DMComputeVariableBounds(DM dm, Vec xl, Vec xu)
3835: {
3836:   PetscFunctionBegin;
3840:   PetscUseTypeMethod(dm, computevariablebounds, xl, xu);
3841:   PetscFunctionReturn(PETSC_SUCCESS);
3842: }

3844: /*@
3845:   DMHasColoring - does the `DM` object have a method of providing a coloring?

3847:   Not Collective

3849:   Input Parameter:
3850: . dm - the DM object

3852:   Output Parameter:
3853: . flg - `PETSC_TRUE` if the `DM` has facilities for `DMCreateColoring()`.

3855:   Level: developer

3857: .seealso: [](ch_dmbase), `DM`, `DMCreateColoring()`
3858: @*/
3859: PetscErrorCode DMHasColoring(DM dm, PetscBool *flg)
3860: {
3861:   PetscFunctionBegin;
3863:   PetscAssertPointer(flg, 2);
3864:   *flg = (dm->ops->getcoloring) ? PETSC_TRUE : PETSC_FALSE;
3865:   PetscFunctionReturn(PETSC_SUCCESS);
3866: }

3868: /*@
3869:   DMHasCreateRestriction - does the `DM` object have a method of providing a restriction?

3871:   Not Collective

3873:   Input Parameter:
3874: . dm - the `DM` object

3876:   Output Parameter:
3877: . flg - `PETSC_TRUE` if the `DM` has facilities for `DMCreateRestriction()`.

3879:   Level: developer

3881: .seealso: [](ch_dmbase), `DM`, `DMCreateRestriction()`, `DMHasCreateInterpolation()`, `DMHasCreateInjection()`
3882: @*/
3883: PetscErrorCode DMHasCreateRestriction(DM dm, PetscBool *flg)
3884: {
3885:   PetscFunctionBegin;
3887:   PetscAssertPointer(flg, 2);
3888:   *flg = (dm->ops->createrestriction) ? PETSC_TRUE : PETSC_FALSE;
3889:   PetscFunctionReturn(PETSC_SUCCESS);
3890: }

3892: /*@
3893:   DMHasCreateInjection - does the `DM` object have a method of providing an injection?

3895:   Not Collective

3897:   Input Parameter:
3898: . dm - the `DM` object

3900:   Output Parameter:
3901: . flg - `PETSC_TRUE` if the `DM` has facilities for `DMCreateInjection()`.

3903:   Level: developer

3905: .seealso: [](ch_dmbase), `DM`, `DMCreateInjection()`, `DMHasCreateRestriction()`, `DMHasCreateInterpolation()`
3906: @*/
3907: PetscErrorCode DMHasCreateInjection(DM dm, PetscBool *flg)
3908: {
3909:   PetscFunctionBegin;
3911:   PetscAssertPointer(flg, 2);
3912:   if (dm->ops->hascreateinjection) PetscUseTypeMethod(dm, hascreateinjection, flg);
3913:   else *flg = (dm->ops->createinjection) ? PETSC_TRUE : PETSC_FALSE;
3914:   PetscFunctionReturn(PETSC_SUCCESS);
3915: }

3917: PetscFunctionList DMList              = NULL;
3918: PetscBool         DMRegisterAllCalled = PETSC_FALSE;

3920: /*@
3921:   DMSetType - Builds a `DM`, for a particular `DM` implementation.

3923:   Collective

3925:   Input Parameters:
3926: + dm     - The `DM` object
3927: - method - The name of the `DMType`, for example `DMDA`, `DMPLEX`

3929:   Options Database Key:
3930: . -dm_type <type> - Sets the `DM` type; use -help for a list of available types

3932:   Level: intermediate

3934:   Note:
3935:   Of the `DM` is constructed by directly calling a function to construct a particular `DM`, for example, `DMDACreate2d()` or `DMPlexCreateBoxMesh()`

3937: .seealso: [](ch_dmbase), `DM`, `DMType`, `DMDA`, `DMPLEX`, `DMGetType()`, `DMCreate()`, `DMDACreate2d()`
3938: @*/
3939: PetscErrorCode DMSetType(DM dm, DMType method)
3940: {
3941:   PetscErrorCode (*r)(DM);
3942:   PetscBool match;

3944:   PetscFunctionBegin;
3946:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, method, &match));
3947:   if (match) PetscFunctionReturn(PETSC_SUCCESS);

3949:   PetscCall(DMRegisterAll());
3950:   PetscCall(PetscFunctionListFind(DMList, method, &r));
3951:   PetscCheck(r, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);

3953:   PetscTryTypeMethod(dm, destroy);
3954:   PetscCall(PetscMemzero(dm->ops, sizeof(*dm->ops)));
3955:   PetscCall(PetscObjectChangeTypeName((PetscObject)dm, method));
3956:   PetscCall((*r)(dm));
3957:   PetscFunctionReturn(PETSC_SUCCESS);
3958: }

3960: /*@
3961:   DMGetType - Gets the `DM` type name (as a string) from the `DM`.

3963:   Not Collective

3965:   Input Parameter:
3966: . dm - The `DM`

3968:   Output Parameter:
3969: . type - The `DMType` name

3971:   Level: intermediate

3973: .seealso: [](ch_dmbase), `DM`, `DMType`, `DMDA`, `DMPLEX`, `DMSetType()`, `DMCreate()`
3974: @*/
3975: PetscErrorCode DMGetType(DM dm, DMType *type)
3976: {
3977:   PetscFunctionBegin;
3979:   PetscAssertPointer(type, 2);
3980:   PetscCall(DMRegisterAll());
3981:   *type = ((PetscObject)dm)->type_name;
3982:   PetscFunctionReturn(PETSC_SUCCESS);
3983: }

3985: /*@
3986:   DMConvert - Converts a `DM` to another `DM`, either of the same or different type.

3988:   Collective

3990:   Input Parameters:
3991: + dm      - the `DM`
3992: - newtype - new `DM` type (use "same" for the same type)

3994:   Output Parameter:
3995: . M - pointer to new `DM`

3997:   Level: intermediate

3999:   Note:
4000:   Cannot be used to convert a sequential `DM` to a parallel or a parallel to sequential,
4001:   the MPI communicator of the generated `DM` is always the same as the communicator
4002:   of the input `DM`.

4004: .seealso: [](ch_dmbase), `DM`, `DMSetType()`, `DMCreate()`, `DMClone()`
4005: @*/
4006: PetscErrorCode DMConvert(DM dm, DMType newtype, DM *M)
4007: {
4008:   DM        B;
4009:   char      convname[256];
4010:   PetscBool sametype /*, issame */;

4012:   PetscFunctionBegin;
4015:   PetscAssertPointer(M, 3);
4016:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, newtype, &sametype));
4017:   /* PetscCall(PetscStrcmp(newtype, "same", &issame)); */
4018:   if (sametype) {
4019:     *M = dm;
4020:     PetscCall(PetscObjectReference((PetscObject)dm));
4021:     PetscFunctionReturn(PETSC_SUCCESS);
4022:   } else {
4023:     PetscErrorCode (*conv)(DM, DMType, DM *) = NULL;

4025:     /*
4026:        Order of precedence:
4027:        1) See if a specialized converter is known to the current DM.
4028:        2) See if a specialized converter is known to the desired DM class.
4029:        3) See if a good general converter is registered for the desired class
4030:        4) See if a good general converter is known for the current matrix.
4031:        5) Use a really basic converter.
4032:     */

4034:     /* 1) See if a specialized converter is known to the current DM and the desired class */
4035:     PetscCall(PetscStrncpy(convname, "DMConvert_", sizeof(convname)));
4036:     PetscCall(PetscStrlcat(convname, ((PetscObject)dm)->type_name, sizeof(convname)));
4037:     PetscCall(PetscStrlcat(convname, "_", sizeof(convname)));
4038:     PetscCall(PetscStrlcat(convname, newtype, sizeof(convname)));
4039:     PetscCall(PetscStrlcat(convname, "_C", sizeof(convname)));
4040:     PetscCall(PetscObjectQueryFunction((PetscObject)dm, convname, &conv));
4041:     if (conv) goto foundconv;

4043:     /* 2)  See if a specialized converter is known to the desired DM class. */
4044:     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &B));
4045:     PetscCall(DMSetType(B, newtype));
4046:     PetscCall(PetscStrncpy(convname, "DMConvert_", sizeof(convname)));
4047:     PetscCall(PetscStrlcat(convname, ((PetscObject)dm)->type_name, sizeof(convname)));
4048:     PetscCall(PetscStrlcat(convname, "_", sizeof(convname)));
4049:     PetscCall(PetscStrlcat(convname, newtype, sizeof(convname)));
4050:     PetscCall(PetscStrlcat(convname, "_C", sizeof(convname)));
4051:     PetscCall(PetscObjectQueryFunction((PetscObject)B, convname, &conv));
4052:     if (conv) {
4053:       PetscCall(DMDestroy(&B));
4054:       goto foundconv;
4055:     }

4057: #if 0
4058:     /* 3) See if a good general converter is registered for the desired class */
4059:     conv = B->ops->convertfrom;
4060:     PetscCall(DMDestroy(&B));
4061:     if (conv) goto foundconv;

4063:     /* 4) See if a good general converter is known for the current matrix */
4064:     if (dm->ops->convert) {
4065:       conv = dm->ops->convert;
4066:     }
4067:     if (conv) goto foundconv;
4068: #endif

4070:     /* 5) Use a really basic converter. */
4071:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject)dm)->type_name, newtype);

4073:   foundconv:
4074:     PetscCall(PetscLogEventBegin(DM_Convert, dm, 0, 0, 0));
4075:     PetscCall((*conv)(dm, newtype, M));
4076:     /* Things that are independent of DM type: We should consult DMClone() here */
4077:     {
4078:       const PetscReal *maxCell, *Lstart, *L;

4080:       PetscCall(DMGetPeriodicity(dm, &maxCell, &Lstart, &L));
4081:       PetscCall(DMSetPeriodicity(*M, maxCell, Lstart, L));
4082:       (*M)->prealloc_only = dm->prealloc_only;
4083:       PetscCall(PetscFree((*M)->vectype));
4084:       PetscCall(PetscStrallocpy(dm->vectype, (char **)&(*M)->vectype));
4085:       PetscCall(PetscFree((*M)->mattype));
4086:       PetscCall(PetscStrallocpy(dm->mattype, (char **)&(*M)->mattype));
4087:     }
4088:     PetscCall(PetscLogEventEnd(DM_Convert, dm, 0, 0, 0));
4089:   }
4090:   PetscCall(PetscObjectStateIncrease((PetscObject)*M));
4091:   PetscFunctionReturn(PETSC_SUCCESS);
4092: }

4094: /*--------------------------------------------------------------------------------------------------------------------*/

4096: /*@C
4097:   DMRegister -  Adds a new `DM` type implementation

4099:   Not Collective, No Fortran Support

4101:   Input Parameters:
4102: + sname    - The name of a new user-defined creation routine
4103: - function - The creation routine itself

4105:   Level: advanced

4107:   Note:
4108:   `DMRegister()` may be called multiple times to add several user-defined `DM`s

4110:   Example Usage:
4111: .vb
4112:     DMRegister("my_da", MyDMCreate);
4113: .ve

4115:   Then, your `DM` type can be chosen with the procedural interface via
4116: .vb
4117:     DMCreate(MPI_Comm, DM *);
4118:     DMSetType(DM,"my_da");
4119: .ve
4120:   or at runtime via the option
4121: .vb
4122:     -da_type my_da
4123: .ve

4125: .seealso: [](ch_dmbase), `DM`, `DMType`, `DMSetType()`, `DMRegisterAll()`, `DMRegisterDestroy()`
4126: @*/
4127: PetscErrorCode DMRegister(const char sname[], PetscErrorCode (*function)(DM))
4128: {
4129:   PetscFunctionBegin;
4130:   PetscCall(DMInitializePackage());
4131:   PetscCall(PetscFunctionListAdd(&DMList, sname, function));
4132:   PetscFunctionReturn(PETSC_SUCCESS);
4133: }

4135: /*@
4136:   DMLoad - Loads a DM that has been stored in binary  with `DMView()`.

4138:   Collective

4140:   Input Parameters:
4141: + newdm  - the newly loaded `DM`, this needs to have been created with `DMCreate()` or
4142:            some related function before a call to `DMLoad()`.
4143: - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()` or
4144:            `PETSCVIEWERHDF5` file viewer, obtained from `PetscViewerHDF5Open()`

4146:   Level: intermediate

4148:   Notes:
4149:   The type is determined by the data in the file, any type set into the DM before this call is ignored.

4151:   Using `PETSCVIEWERHDF5` type with `PETSC_VIEWER_HDF5_PETSC` format, one can save multiple `DMPLEX`
4152:   meshes in a single HDF5 file. This in turn requires one to name the `DMPLEX` object with `PetscObjectSetName()`
4153:   before saving it with `DMView()` and before loading it with `DMLoad()` for identification of the mesh object.

4155: .seealso: [](ch_dmbase), `DM`, `PetscViewerBinaryOpen()`, `DMView()`, `MatLoad()`, `VecLoad()`
4156: @*/
4157: PetscErrorCode DMLoad(DM newdm, PetscViewer viewer)
4158: {
4159:   PetscBool isbinary, ishdf5;

4161:   PetscFunctionBegin;
4164:   PetscCall(PetscViewerCheckReadable(viewer));
4165:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
4166:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
4167:   PetscCall(PetscLogEventBegin(DM_Load, viewer, 0, 0, 0));
4168:   if (isbinary) {
4169:     PetscInt classid;
4170:     char     type[256];

4172:     PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT));
4173:     PetscCheck(classid == DM_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not DM next in file, classid found %" PetscInt_FMT, classid);
4174:     PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
4175:     PetscCall(DMSetType(newdm, type));
4176:     PetscTryTypeMethod(newdm, load, viewer);
4177:   } else if (ishdf5) {
4178:     PetscTryTypeMethod(newdm, load, viewer);
4179:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen() or PetscViewerHDF5Open()");
4180:   PetscCall(PetscLogEventEnd(DM_Load, viewer, 0, 0, 0));
4181:   PetscFunctionReturn(PETSC_SUCCESS);
4182: }

4184: /* FEM Support */

4186: PetscErrorCode DMPrintCellIndices(PetscInt c, const char name[], PetscInt len, const PetscInt x[])
4187: {
4188:   PetscInt f;

4190:   PetscFunctionBegin;
4191:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " Element %s\n", c, name));
4192:   for (f = 0; f < len; ++f) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  | %" PetscInt_FMT " |\n", x[f]));
4193:   PetscFunctionReturn(PETSC_SUCCESS);
4194: }

4196: PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[])
4197: {
4198:   PetscInt f;

4200:   PetscFunctionBegin;
4201:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " Element %s\n", c, name));
4202:   for (f = 0; f < len; ++f) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  | %g |\n", (double)PetscRealPart(x[f])));
4203:   PetscFunctionReturn(PETSC_SUCCESS);
4204: }

4206: PetscErrorCode DMPrintCellVectorReal(PetscInt c, const char name[], PetscInt len, const PetscReal x[])
4207: {
4208:   PetscInt f;

4210:   PetscFunctionBegin;
4211:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " Element %s\n", c, name));
4212:   for (f = 0; f < len; ++f) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  | %g |\n", (double)x[f]));
4213:   PetscFunctionReturn(PETSC_SUCCESS);
4214: }

4216: PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[])
4217: {
4218:   PetscInt f, g;

4220:   PetscFunctionBegin;
4221:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " Element %s\n", c, name));
4222:   for (f = 0; f < rows; ++f) {
4223:     PetscCall(PetscPrintf(PETSC_COMM_SELF, "  |"));
4224:     for (g = 0; g < cols; ++g) PetscCall(PetscPrintf(PETSC_COMM_SELF, " % 9.5g", (double)PetscRealPart(A[f * cols + g])));
4225:     PetscCall(PetscPrintf(PETSC_COMM_SELF, " |\n"));
4226:   }
4227:   PetscFunctionReturn(PETSC_SUCCESS);
4228: }

4230: PetscErrorCode DMPrintLocalVec(DM dm, const char name[], PetscReal tol, Vec X)
4231: {
4232:   PetscInt           localSize, bs;
4233:   PetscMPIInt        size;
4234:   Vec                x, xglob;
4235:   const PetscScalar *xarray;

4237:   PetscFunctionBegin;
4238:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
4239:   PetscCall(VecDuplicate(X, &x));
4240:   PetscCall(VecCopy(X, x));
4241:   PetscCall(VecFilter(x, tol));
4242:   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "%s:\n", name));
4243:   if (size > 1) {
4244:     PetscCall(VecGetLocalSize(x, &localSize));
4245:     PetscCall(VecGetArrayRead(x, &xarray));
4246:     PetscCall(VecGetBlockSize(x, &bs));
4247:     PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)dm), bs, localSize, PETSC_DETERMINE, xarray, &xglob));
4248:   } else {
4249:     xglob = x;
4250:   }
4251:   PetscCall(VecView(xglob, PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm))));
4252:   if (size > 1) {
4253:     PetscCall(VecDestroy(&xglob));
4254:     PetscCall(VecRestoreArrayRead(x, &xarray));
4255:   }
4256:   PetscCall(VecDestroy(&x));
4257:   PetscFunctionReturn(PETSC_SUCCESS);
4258: }

4260: /*@
4261:   DMGetLocalSection - Get the `PetscSection` encoding the local data layout for the `DM`.

4263:   Input Parameter:
4264: . dm - The `DM`

4266:   Output Parameter:
4267: . section - The `PetscSection`

4269:   Options Database Key:
4270: . -dm_petscsection_view - View the section created by the `DM`

4272:   Level: intermediate

4274:   Note:
4275:   This gets a borrowed reference, so the user should not destroy this `PetscSection`.

4277: .seealso: [](ch_dmbase), `DM`, `DMSetLocalSection()`, `DMGetGlobalSection()`
4278: @*/
4279: PetscErrorCode DMGetLocalSection(DM dm, PetscSection *section)
4280: {
4281:   PetscFunctionBegin;
4283:   PetscAssertPointer(section, 2);
4284:   if (!dm->localSection && dm->ops->createlocalsection) {
4285:     PetscInt d;

4287:     if (dm->setfromoptionscalled) {
4288:       PetscObject       obj = (PetscObject)dm;
4289:       PetscViewer       viewer;
4290:       PetscViewerFormat format;
4291:       PetscBool         flg;

4293:       PetscCall(PetscOptionsCreateViewer(PetscObjectComm(obj), obj->options, obj->prefix, "-dm_petscds_view", &viewer, &format, &flg));
4294:       if (flg) PetscCall(PetscViewerPushFormat(viewer, format));
4295:       for (d = 0; d < dm->Nds; ++d) {
4296:         PetscCall(PetscDSSetFromOptions(dm->probs[d].ds));
4297:         if (flg) PetscCall(PetscDSView(dm->probs[d].ds, viewer));
4298:       }
4299:       if (flg) {
4300:         PetscCall(PetscViewerFlush(viewer));
4301:         PetscCall(PetscViewerPopFormat(viewer));
4302:         PetscCall(PetscViewerDestroy(&viewer));
4303:       }
4304:     }
4305:     PetscUseTypeMethod(dm, createlocalsection);
4306:     if (dm->localSection) PetscCall(PetscObjectViewFromOptions((PetscObject)dm->localSection, NULL, "-dm_petscsection_view"));
4307:   }
4308:   *section = dm->localSection;
4309:   PetscFunctionReturn(PETSC_SUCCESS);
4310: }

4312: /*@
4313:   DMSetLocalSection - Set the `PetscSection` encoding the local data layout for the `DM`.

4315:   Input Parameters:
4316: + dm      - The `DM`
4317: - section - The `PetscSection`

4319:   Level: intermediate

4321:   Note:
4322:   Any existing Section will be destroyed

4324: .seealso: [](ch_dmbase), `DM`, `PetscSection`, `DMGetLocalSection()`, `DMSetGlobalSection()`
4325: @*/
4326: PetscErrorCode DMSetLocalSection(DM dm, PetscSection section)
4327: {
4328:   PetscInt numFields = 0;
4329:   PetscInt f;

4331:   PetscFunctionBegin;
4334:   PetscCall(PetscObjectReference((PetscObject)section));
4335:   PetscCall(PetscSectionDestroy(&dm->localSection));
4336:   dm->localSection = section;
4337:   if (section) PetscCall(PetscSectionGetNumFields(dm->localSection, &numFields));
4338:   if (numFields) {
4339:     PetscCall(DMSetNumFields(dm, numFields));
4340:     for (f = 0; f < numFields; ++f) {
4341:       PetscObject disc;
4342:       const char *name;

4344:       PetscCall(PetscSectionGetFieldName(dm->localSection, f, &name));
4345:       PetscCall(DMGetField(dm, f, NULL, &disc));
4346:       PetscCall(PetscObjectSetName(disc, name));
4347:     }
4348:   }
4349:   /* The global section and the SectionSF will be rebuilt
4350:      in the next call to DMGetGlobalSection() and DMGetSectionSF(). */
4351:   PetscCall(PetscSectionDestroy(&dm->globalSection));
4352:   PetscCall(PetscSFDestroy(&dm->sectionSF));
4353:   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &dm->sectionSF));

4355:   /* Clear scratch vectors */
4356:   PetscCall(DMClearGlobalVectors(dm));
4357:   PetscCall(DMClearLocalVectors(dm));
4358:   PetscCall(DMClearNamedGlobalVectors(dm));
4359:   PetscCall(DMClearNamedLocalVectors(dm));
4360:   PetscFunctionReturn(PETSC_SUCCESS);
4361: }

4363: /*@C
4364:   DMCreateSectionPermutation - Create a permutation of the `PetscSection` chart and optionally a block structure.

4366:   Input Parameter:
4367: . dm - The `DM`

4369:   Output Parameters:
4370: + perm        - A permutation of the mesh points in the chart
4371: - blockStarts - A high bit is set for the point that begins every block, or `NULL` for default blocking

4373:   Level: developer

4375: .seealso: [](ch_dmbase), `DM`, `PetscSection`, `DMGetLocalSection()`, `DMGetGlobalSection()`
4376: @*/
4377: PetscErrorCode DMCreateSectionPermutation(DM dm, IS *perm, PetscBT *blockStarts)
4378: {
4379:   PetscFunctionBegin;
4380:   *perm        = NULL;
4381:   *blockStarts = NULL;
4382:   PetscTryTypeMethod(dm, createsectionpermutation, perm, blockStarts);
4383:   PetscFunctionReturn(PETSC_SUCCESS);
4384: }

4386: /*@
4387:   DMGetDefaultConstraints - Get the `PetscSection` and `Mat` that specify the local constraint interpolation. See `DMSetDefaultConstraints()` for a description of the purpose of constraint interpolation.

4389:   not Collective

4391:   Input Parameter:
4392: . dm - The `DM`

4394:   Output Parameters:
4395: + section - The `PetscSection` describing the range of the constraint matrix: relates rows of the constraint matrix to dofs of the default section.  Returns `NULL` if there are no local constraints.
4396: . mat     - The `Mat` that interpolates local constraints: its width should be the layout size of the default section.  Returns `NULL` if there are no local constraints.
4397: - bias    - Vector containing bias to be added to constrained dofs

4399:   Level: advanced

4401:   Note:
4402:   This gets borrowed references, so the user should not destroy the `PetscSection`, `Mat`, or `Vec`.

4404: .seealso: [](ch_dmbase), `DM`, `DMSetDefaultConstraints()`
4405: @*/
4406: PetscErrorCode DMGetDefaultConstraints(DM dm, PetscSection *section, Mat *mat, Vec *bias)
4407: {
4408:   PetscFunctionBegin;
4410:   if (!dm->defaultConstraint.section && !dm->defaultConstraint.mat && dm->ops->createdefaultconstraints) PetscUseTypeMethod(dm, createdefaultconstraints);
4411:   if (section) *section = dm->defaultConstraint.section;
4412:   if (mat) *mat = dm->defaultConstraint.mat;
4413:   if (bias) *bias = dm->defaultConstraint.bias;
4414:   PetscFunctionReturn(PETSC_SUCCESS);
4415: }

4417: /*@
4418:   DMSetDefaultConstraints - Set the `PetscSection` and `Mat` that specify the local constraint interpolation.

4420:   Collective

4422:   Input Parameters:
4423: + dm      - The `DM`
4424: . section - The `PetscSection` describing the range of the constraint matrix: relates rows of the constraint matrix to dofs of the default section.  Must have a local communicator (`PETSC_COMM_SELF` or derivative).
4425: . mat     - The `Mat` that interpolates local constraints: its width should be the layout size of the default section:  `NULL` indicates no constraints.  Must have a local communicator (`PETSC_COMM_SELF` or derivative).
4426: - bias    - A bias vector to be added to constrained values in the local vector.  `NULL` indicates no bias.  Must have a local communicator (`PETSC_COMM_SELF` or derivative).

4428:   Level: advanced

4430:   Notes:
4431:   If a constraint matrix is specified, then it is applied during `DMGlobalToLocalEnd()` when mode is `INSERT_VALUES`, `INSERT_BC_VALUES`, or `INSERT_ALL_VALUES`.  Without a constraint matrix, the local vector l returned by `DMGlobalToLocalEnd()` contains values that have been scattered from a global vector without modification; with a constraint matrix A, l is modified by computing c = A * l + bias, l[s[i]] = c[i], where the scatter s is defined by the `PetscSection` returned by `DMGetDefaultConstraints()`.

4433:   If a constraint matrix is specified, then its adjoint is applied during `DMLocalToGlobalBegin()` when mode is `ADD_VALUES`, `ADD_BC_VALUES`, or `ADD_ALL_VALUES`.  Without a constraint matrix, the local vector l is accumulated into a global vector without modification; with a constraint matrix A, l is first modified by computing c[i] = l[s[i]], l[s[i]] = 0, l = l + A'*c, which is the adjoint of the operation described above.  Any bias, if specified, is ignored when accumulating.

4435:   This increments the references of the `PetscSection`, `Mat`, and `Vec`, so they user can destroy them.

4437: .seealso: [](ch_dmbase), `DM`, `DMGetDefaultConstraints()`
4438: @*/
4439: PetscErrorCode DMSetDefaultConstraints(DM dm, PetscSection section, Mat mat, Vec bias)
4440: {
4441:   PetscMPIInt result;

4443:   PetscFunctionBegin;
4445:   if (section) {
4447:     PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, PetscObjectComm((PetscObject)section), &result));
4448:     PetscCheck(result == MPI_CONGRUENT || result == MPI_IDENT, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMECOMM, "constraint section must have local communicator");
4449:   }
4450:   if (mat) {
4452:     PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, PetscObjectComm((PetscObject)mat), &result));
4453:     PetscCheck(result == MPI_CONGRUENT || result == MPI_IDENT, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMECOMM, "constraint matrix must have local communicator");
4454:   }
4455:   if (bias) {
4457:     PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, PetscObjectComm((PetscObject)bias), &result));
4458:     PetscCheck(result == MPI_CONGRUENT || result == MPI_IDENT, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMECOMM, "constraint bias must have local communicator");
4459:   }
4460:   PetscCall(PetscObjectReference((PetscObject)section));
4461:   PetscCall(PetscSectionDestroy(&dm->defaultConstraint.section));
4462:   dm->defaultConstraint.section = section;
4463:   PetscCall(PetscObjectReference((PetscObject)mat));
4464:   PetscCall(MatDestroy(&dm->defaultConstraint.mat));
4465:   dm->defaultConstraint.mat = mat;
4466:   PetscCall(PetscObjectReference((PetscObject)bias));
4467:   PetscCall(VecDestroy(&dm->defaultConstraint.bias));
4468:   dm->defaultConstraint.bias = bias;
4469:   PetscFunctionReturn(PETSC_SUCCESS);
4470: }

4472: #if defined(PETSC_USE_DEBUG)
4473: /*
4474:   DMDefaultSectionCheckConsistency - Check the consistentcy of the global and local sections. Generates and error if they are not consistent.

4476:   Input Parameters:
4477: + dm - The `DM`
4478: . localSection - `PetscSection` describing the local data layout
4479: - globalSection - `PetscSection` describing the global data layout

4481:   Level: intermediate

4483: .seealso: [](ch_dmbase), `DM`, `DMGetSectionSF()`, `DMSetSectionSF()`
4484: */
4485: static PetscErrorCode DMDefaultSectionCheckConsistency_Internal(DM dm, PetscSection localSection, PetscSection globalSection)
4486: {
4487:   MPI_Comm        comm;
4488:   PetscLayout     layout;
4489:   const PetscInt *ranges;
4490:   PetscInt        pStart, pEnd, p, nroots;
4491:   PetscMPIInt     size, rank;
4492:   PetscBool       valid = PETSC_TRUE, gvalid;

4494:   PetscFunctionBegin;
4495:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
4497:   PetscCallMPI(MPI_Comm_size(comm, &size));
4498:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
4499:   PetscCall(PetscSectionGetChart(globalSection, &pStart, &pEnd));
4500:   PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &nroots));
4501:   PetscCall(PetscLayoutCreate(comm, &layout));
4502:   PetscCall(PetscLayoutSetBlockSize(layout, 1));
4503:   PetscCall(PetscLayoutSetLocalSize(layout, nroots));
4504:   PetscCall(PetscLayoutSetUp(layout));
4505:   PetscCall(PetscLayoutGetRanges(layout, &ranges));
4506:   for (p = pStart; p < pEnd; ++p) {
4507:     PetscInt dof, cdof, off, gdof, gcdof, goff, gsize, d;

4509:     PetscCall(PetscSectionGetDof(localSection, p, &dof));
4510:     PetscCall(PetscSectionGetOffset(localSection, p, &off));
4511:     PetscCall(PetscSectionGetConstraintDof(localSection, p, &cdof));
4512:     PetscCall(PetscSectionGetDof(globalSection, p, &gdof));
4513:     PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof));
4514:     PetscCall(PetscSectionGetOffset(globalSection, p, &goff));
4515:     if (!gdof) continue; /* Censored point */
4516:     if ((gdof < 0 ? -(gdof + 1) : gdof) != dof) {
4517:       PetscCall(PetscSynchronizedPrintf(comm, "[%d]Global dof %" PetscInt_FMT " for point %" PetscInt_FMT " not equal to local dof %" PetscInt_FMT "\n", rank, gdof, p, dof));
4518:       valid = PETSC_FALSE;
4519:     }
4520:     if (gcdof && (gcdof != cdof)) {
4521:       PetscCall(PetscSynchronizedPrintf(comm, "[%d]Global constraints %" PetscInt_FMT " for point %" PetscInt_FMT " not equal to local constraints %" PetscInt_FMT "\n", rank, gcdof, p, cdof));
4522:       valid = PETSC_FALSE;
4523:     }
4524:     if (gdof < 0) {
4525:       gsize = gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof;
4526:       for (d = 0; d < gsize; ++d) {
4527:         PetscInt offset = -(goff + 1) + d, r;

4529:         PetscCall(PetscFindInt(offset, size + 1, ranges, &r));
4530:         if (r < 0) r = -(r + 2);
4531:         if ((r < 0) || (r >= size)) {
4532:           PetscCall(PetscSynchronizedPrintf(comm, "[%d]Point %" PetscInt_FMT " mapped to invalid process %" PetscInt_FMT " (%" PetscInt_FMT ", %" PetscInt_FMT ")\n", rank, p, r, gdof, goff));
4533:           valid = PETSC_FALSE;
4534:           break;
4535:         }
4536:       }
4537:     }
4538:   }
4539:   PetscCall(PetscLayoutDestroy(&layout));
4540:   PetscCall(PetscSynchronizedFlush(comm, NULL));
4541:   PetscCallMPI(MPIU_Allreduce(&valid, &gvalid, 1, MPIU_BOOL, MPI_LAND, comm));
4542:   if (!gvalid) {
4543:     PetscCall(DMView(dm, NULL));
4544:     SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Inconsistent local and global sections");
4545:   }
4546:   PetscFunctionReturn(PETSC_SUCCESS);
4547: }
4548: #endif

4550: PetscErrorCode DMGetIsoperiodicPointSF_Internal(DM dm, PetscSF *sf)
4551: {
4552:   PetscErrorCode (*f)(DM, PetscSF *);

4554:   PetscFunctionBegin;
4556:   PetscAssertPointer(sf, 2);
4557:   PetscCall(PetscObjectQueryFunction((PetscObject)dm, "DMGetIsoperiodicPointSF_C", &f));
4558:   if (f) PetscCall(f(dm, sf));
4559:   else *sf = dm->sf;
4560:   PetscFunctionReturn(PETSC_SUCCESS);
4561: }

4563: /*@
4564:   DMGetGlobalSection - Get the `PetscSection` encoding the global data layout for the `DM`.

4566:   Collective

4568:   Input Parameter:
4569: . dm - The `DM`

4571:   Output Parameter:
4572: . section - The `PetscSection`

4574:   Level: intermediate

4576:   Note:
4577:   This gets a borrowed reference, so the user should not destroy this `PetscSection`.

4579: .seealso: [](ch_dmbase), `DM`, `DMSetLocalSection()`, `DMGetLocalSection()`
4580: @*/
4581: PetscErrorCode DMGetGlobalSection(DM dm, PetscSection *section)
4582: {
4583:   PetscFunctionBegin;
4585:   PetscAssertPointer(section, 2);
4586:   if (!dm->globalSection) {
4587:     PetscSection s;
4588:     PetscSF      sf;

4590:     PetscCall(DMGetLocalSection(dm, &s));
4591:     PetscCheck(s, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM must have a default PetscSection in order to create a global PetscSection");
4592:     PetscCheck(dm->sf, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM must have a point PetscSF in order to create a global PetscSection");
4593:     PetscCall(DMGetIsoperiodicPointSF_Internal(dm, &sf));
4594:     PetscCall(PetscSectionCreateGlobalSection(s, sf, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &dm->globalSection));
4595:     PetscCall(PetscLayoutDestroy(&dm->map));
4596:     PetscCall(PetscSectionGetValueLayout(PetscObjectComm((PetscObject)dm), dm->globalSection, &dm->map));
4597:     PetscCall(PetscSectionViewFromOptions(dm->globalSection, NULL, "-global_section_view"));
4598:   }
4599:   *section = dm->globalSection;
4600:   PetscFunctionReturn(PETSC_SUCCESS);
4601: }

4603: /*@
4604:   DMSetGlobalSection - Set the `PetscSection` encoding the global data layout for the `DM`.

4606:   Input Parameters:
4607: + dm      - The `DM`
4608: - section - The PetscSection, or `NULL`

4610:   Level: intermediate

4612:   Note:
4613:   Any existing `PetscSection` will be destroyed

4615: .seealso: [](ch_dmbase), `DM`, `DMGetGlobalSection()`, `DMSetLocalSection()`
4616: @*/
4617: PetscErrorCode DMSetGlobalSection(DM dm, PetscSection section)
4618: {
4619:   PetscFunctionBegin;
4622:   PetscCall(PetscObjectReference((PetscObject)section));
4623:   PetscCall(PetscSectionDestroy(&dm->globalSection));
4624:   dm->globalSection = section;
4625: #if defined(PETSC_USE_DEBUG)
4626:   if (section) PetscCall(DMDefaultSectionCheckConsistency_Internal(dm, dm->localSection, section));
4627: #endif
4628:   /* Clear global scratch vectors and sectionSF */
4629:   PetscCall(PetscSFDestroy(&dm->sectionSF));
4630:   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &dm->sectionSF));
4631:   PetscCall(DMClearGlobalVectors(dm));
4632:   PetscCall(DMClearNamedGlobalVectors(dm));
4633:   PetscFunctionReturn(PETSC_SUCCESS);
4634: }

4636: /*@
4637:   DMGetSectionSF - Get the `PetscSF` encoding the parallel dof overlap for the `DM`. If it has not been set,
4638:   it is created from the default `PetscSection` layouts in the `DM`.

4640:   Input Parameter:
4641: . dm - The `DM`

4643:   Output Parameter:
4644: . sf - The `PetscSF`

4646:   Level: intermediate

4648:   Note:
4649:   This gets a borrowed reference, so the user should not destroy this `PetscSF`.

4651: .seealso: [](ch_dmbase), `DM`, `DMSetSectionSF()`, `DMCreateSectionSF()`
4652: @*/
4653: PetscErrorCode DMGetSectionSF(DM dm, PetscSF *sf)
4654: {
4655:   PetscInt nroots;

4657:   PetscFunctionBegin;
4659:   PetscAssertPointer(sf, 2);
4660:   if (!dm->sectionSF) PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &dm->sectionSF));
4661:   PetscCall(PetscSFGetGraph(dm->sectionSF, &nroots, NULL, NULL, NULL));
4662:   if (nroots < 0) {
4663:     PetscSection section, gSection;

4665:     PetscCall(DMGetLocalSection(dm, &section));
4666:     if (section) {
4667:       PetscCall(DMGetGlobalSection(dm, &gSection));
4668:       PetscCall(DMCreateSectionSF(dm, section, gSection));
4669:     } else {
4670:       *sf = NULL;
4671:       PetscFunctionReturn(PETSC_SUCCESS);
4672:     }
4673:   }
4674:   *sf = dm->sectionSF;
4675:   PetscFunctionReturn(PETSC_SUCCESS);
4676: }

4678: /*@
4679:   DMSetSectionSF - Set the `PetscSF` encoding the parallel dof overlap for the `DM`

4681:   Input Parameters:
4682: + dm - The `DM`
4683: - sf - The `PetscSF`

4685:   Level: intermediate

4687:   Note:
4688:   Any previous `PetscSF` is destroyed

4690: .seealso: [](ch_dmbase), `DM`, `DMGetSectionSF()`, `DMCreateSectionSF()`
4691: @*/
4692: PetscErrorCode DMSetSectionSF(DM dm, PetscSF sf)
4693: {
4694:   PetscFunctionBegin;
4697:   PetscCall(PetscObjectReference((PetscObject)sf));
4698:   PetscCall(PetscSFDestroy(&dm->sectionSF));
4699:   dm->sectionSF = sf;
4700:   PetscFunctionReturn(PETSC_SUCCESS);
4701: }

4703: /*@
4704:   DMCreateSectionSF - Create the `PetscSF` encoding the parallel dof overlap for the `DM` based upon the `PetscSection`s
4705:   describing the data layout.

4707:   Input Parameters:
4708: + dm            - The `DM`
4709: . localSection  - `PetscSection` describing the local data layout
4710: - globalSection - `PetscSection` describing the global data layout

4712:   Level: developer

4714:   Note:
4715:   One usually uses `DMGetSectionSF()` to obtain the `PetscSF`

4717:   Developer Note:
4718:   Since this routine has for arguments the two sections from the `DM` and puts the resulting `PetscSF`
4719:   directly into the `DM`, perhaps this function should not take the local and global sections as
4720:   input and should just obtain them from the `DM`? Plus PETSc creation functions return the thing
4721:   they create, this returns nothing

4723: .seealso: [](ch_dmbase), `DM`, `DMGetSectionSF()`, `DMSetSectionSF()`, `DMGetLocalSection()`, `DMGetGlobalSection()`
4724: @*/
4725: PetscErrorCode DMCreateSectionSF(DM dm, PetscSection localSection, PetscSection globalSection)
4726: {
4727:   PetscFunctionBegin;
4729:   PetscCall(PetscSFSetGraphSection(dm->sectionSF, localSection, globalSection));
4730:   PetscFunctionReturn(PETSC_SUCCESS);
4731: }

4733: /*@
4734:   DMGetPointSF - Get the `PetscSF` encoding the parallel section point overlap for the `DM`.

4736:   Not collective but the resulting `PetscSF` is collective

4738:   Input Parameter:
4739: . dm - The `DM`

4741:   Output Parameter:
4742: . sf - The `PetscSF`

4744:   Level: intermediate

4746:   Note:
4747:   This gets a borrowed reference, so the user should not destroy this `PetscSF`.

4749: .seealso: [](ch_dmbase), `DM`, `DMSetPointSF()`, `DMGetSectionSF()`, `DMSetSectionSF()`, `DMCreateSectionSF()`
4750: @*/
4751: PetscErrorCode DMGetPointSF(DM dm, PetscSF *sf)
4752: {
4753:   PetscFunctionBegin;
4755:   PetscAssertPointer(sf, 2);
4756:   *sf = dm->sf;
4757:   PetscFunctionReturn(PETSC_SUCCESS);
4758: }

4760: /*@
4761:   DMSetPointSF - Set the `PetscSF` encoding the parallel section point overlap for the `DM`.

4763:   Collective

4765:   Input Parameters:
4766: + dm - The `DM`
4767: - sf - The `PetscSF`

4769:   Level: intermediate

4771: .seealso: [](ch_dmbase), `DM`, `DMGetPointSF()`, `DMGetSectionSF()`, `DMSetSectionSF()`, `DMCreateSectionSF()`
4772: @*/
4773: PetscErrorCode DMSetPointSF(DM dm, PetscSF sf)
4774: {
4775:   PetscFunctionBegin;
4778:   PetscCall(PetscObjectReference((PetscObject)sf));
4779:   PetscCall(PetscSFDestroy(&dm->sf));
4780:   dm->sf = sf;
4781:   PetscFunctionReturn(PETSC_SUCCESS);
4782: }

4784: /*@
4785:   DMGetNaturalSF - Get the `PetscSF` encoding the map back to the original mesh ordering

4787:   Input Parameter:
4788: . dm - The `DM`

4790:   Output Parameter:
4791: . sf - The `PetscSF`

4793:   Level: intermediate

4795:   Note:
4796:   This gets a borrowed reference, so the user should not destroy this `PetscSF`.

4798: .seealso: [](ch_dmbase), `DM`, `DMSetNaturalSF()`, `DMSetUseNatural()`, `DMGetUseNatural()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexDistribute()`
4799: @*/
4800: PetscErrorCode DMGetNaturalSF(DM dm, PetscSF *sf)
4801: {
4802:   PetscFunctionBegin;
4804:   PetscAssertPointer(sf, 2);
4805:   *sf = dm->sfNatural;
4806:   PetscFunctionReturn(PETSC_SUCCESS);
4807: }

4809: /*@
4810:   DMSetNaturalSF - Set the PetscSF encoding the map back to the original mesh ordering

4812:   Input Parameters:
4813: + dm - The DM
4814: - sf - The PetscSF

4816:   Level: intermediate

4818: .seealso: [](ch_dmbase), `DM`, `DMGetNaturalSF()`, `DMSetUseNatural()`, `DMGetUseNatural()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexDistribute()`
4819: @*/
4820: PetscErrorCode DMSetNaturalSF(DM dm, PetscSF sf)
4821: {
4822:   PetscFunctionBegin;
4825:   PetscCall(PetscObjectReference((PetscObject)sf));
4826:   PetscCall(PetscSFDestroy(&dm->sfNatural));
4827:   dm->sfNatural = sf;
4828:   PetscFunctionReturn(PETSC_SUCCESS);
4829: }

4831: static PetscErrorCode DMSetDefaultAdjacency_Private(DM dm, PetscInt f, PetscObject disc)
4832: {
4833:   PetscClassId id;

4835:   PetscFunctionBegin;
4836:   PetscCall(PetscObjectGetClassId(disc, &id));
4837:   if (id == PETSCFE_CLASSID) {
4838:     PetscCall(DMSetAdjacency(dm, f, PETSC_FALSE, PETSC_TRUE));
4839:   } else if (id == PETSCFV_CLASSID) {
4840:     PetscCall(DMSetAdjacency(dm, f, PETSC_TRUE, PETSC_FALSE));
4841:   } else {
4842:     PetscCall(DMSetAdjacency(dm, f, PETSC_FALSE, PETSC_TRUE));
4843:   }
4844:   PetscFunctionReturn(PETSC_SUCCESS);
4845: }

4847: static PetscErrorCode DMFieldEnlarge_Static(DM dm, PetscInt NfNew)
4848: {
4849:   RegionField *tmpr;
4850:   PetscInt     Nf = dm->Nf, f;

4852:   PetscFunctionBegin;
4853:   if (Nf >= NfNew) PetscFunctionReturn(PETSC_SUCCESS);
4854:   PetscCall(PetscMalloc1(NfNew, &tmpr));
4855:   for (f = 0; f < Nf; ++f) tmpr[f] = dm->fields[f];
4856:   for (f = Nf; f < NfNew; ++f) {
4857:     tmpr[f].disc        = NULL;
4858:     tmpr[f].label       = NULL;
4859:     tmpr[f].avoidTensor = PETSC_FALSE;
4860:   }
4861:   PetscCall(PetscFree(dm->fields));
4862:   dm->Nf     = NfNew;
4863:   dm->fields = tmpr;
4864:   PetscFunctionReturn(PETSC_SUCCESS);
4865: }

4867: /*@
4868:   DMClearFields - Remove all fields from the `DM`

4870:   Logically Collective

4872:   Input Parameter:
4873: . dm - The `DM`

4875:   Level: intermediate

4877: .seealso: [](ch_dmbase), `DM`, `DMGetNumFields()`, `DMSetNumFields()`, `DMSetField()`
4878: @*/
4879: PetscErrorCode DMClearFields(DM dm)
4880: {
4881:   PetscInt f;

4883:   PetscFunctionBegin;
4885:   for (f = 0; f < dm->Nf; ++f) {
4886:     PetscCall(PetscObjectDestroy(&dm->fields[f].disc));
4887:     PetscCall(DMLabelDestroy(&dm->fields[f].label));
4888:   }
4889:   PetscCall(PetscFree(dm->fields));
4890:   dm->fields = NULL;
4891:   dm->Nf     = 0;
4892:   PetscFunctionReturn(PETSC_SUCCESS);
4893: }

4895: /*@
4896:   DMGetNumFields - Get the number of fields in the `DM`

4898:   Not Collective

4900:   Input Parameter:
4901: . dm - The `DM`

4903:   Output Parameter:
4904: . numFields - The number of fields

4906:   Level: intermediate

4908: .seealso: [](ch_dmbase), `DM`, `DMSetNumFields()`, `DMSetField()`
4909: @*/
4910: PetscErrorCode DMGetNumFields(DM dm, PetscInt *numFields)
4911: {
4912:   PetscFunctionBegin;
4914:   PetscAssertPointer(numFields, 2);
4915:   *numFields = dm->Nf;
4916:   PetscFunctionReturn(PETSC_SUCCESS);
4917: }

4919: /*@
4920:   DMSetNumFields - Set the number of fields in the `DM`

4922:   Logically Collective

4924:   Input Parameters:
4925: + dm        - The `DM`
4926: - numFields - The number of fields

4928:   Level: intermediate

4930: .seealso: [](ch_dmbase), `DM`, `DMGetNumFields()`, `DMSetField()`
4931: @*/
4932: PetscErrorCode DMSetNumFields(DM dm, PetscInt numFields)
4933: {
4934:   PetscInt Nf, f;

4936:   PetscFunctionBegin;
4938:   PetscCall(DMGetNumFields(dm, &Nf));
4939:   for (f = Nf; f < numFields; ++f) {
4940:     PetscContainer obj;

4942:     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)dm), &obj));
4943:     PetscCall(DMAddField(dm, NULL, (PetscObject)obj));
4944:     PetscCall(PetscContainerDestroy(&obj));
4945:   }
4946:   PetscFunctionReturn(PETSC_SUCCESS);
4947: }

4949: /*@
4950:   DMGetField - Return the `DMLabel` and discretization object for a given `DM` field

4952:   Not Collective

4954:   Input Parameters:
4955: + dm - The `DM`
4956: - f  - The field number

4958:   Output Parameters:
4959: + label - The label indicating the support of the field, or `NULL` for the entire mesh (pass in `NULL` if not needed)
4960: - disc  - The discretization object (pass in `NULL` if not needed)

4962:   Level: intermediate

4964: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMSetField()`
4965: @*/
4966: PetscErrorCode DMGetField(DM dm, PetscInt f, DMLabel *label, PetscObject *disc)
4967: {
4968:   PetscFunctionBegin;
4970:   PetscAssertPointer(disc, 4);
4971:   PetscCheck((f >= 0) && (f < dm->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, dm->Nf);
4972:   if (label) *label = dm->fields[f].label;
4973:   if (disc) *disc = dm->fields[f].disc;
4974:   PetscFunctionReturn(PETSC_SUCCESS);
4975: }

4977: /* Does not clear the DS */
4978: PetscErrorCode DMSetField_Internal(DM dm, PetscInt f, DMLabel label, PetscObject disc)
4979: {
4980:   PetscFunctionBegin;
4981:   PetscCall(DMFieldEnlarge_Static(dm, f + 1));
4982:   PetscCall(DMLabelDestroy(&dm->fields[f].label));
4983:   PetscCall(PetscObjectDestroy(&dm->fields[f].disc));
4984:   dm->fields[f].label = label;
4985:   dm->fields[f].disc  = disc;
4986:   PetscCall(PetscObjectReference((PetscObject)label));
4987:   PetscCall(PetscObjectReference(disc));
4988:   PetscFunctionReturn(PETSC_SUCCESS);
4989: }

4991: /*@
4992:   DMSetField - Set the discretization object for a given `DM` field. Usually one would call `DMAddField()` which automatically handles
4993:   the field numbering.

4995:   Logically Collective

4997:   Input Parameters:
4998: + dm    - The `DM`
4999: . f     - The field number
5000: . label - The label indicating the support of the field, or `NULL` for the entire mesh
5001: - disc  - The discretization object

5003:   Level: intermediate

5005: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMGetField()`
5006: @*/
5007: PetscErrorCode DMSetField(DM dm, PetscInt f, DMLabel label, PetscObject disc)
5008: {
5009:   PetscFunctionBegin;
5013:   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
5014:   PetscCall(DMSetField_Internal(dm, f, label, disc));
5015:   PetscCall(DMSetDefaultAdjacency_Private(dm, f, disc));
5016:   PetscCall(DMClearDS(dm));
5017:   PetscFunctionReturn(PETSC_SUCCESS);
5018: }

5020: /*@
5021:   DMAddField - Add a field to a `DM` object. A field is a function space defined by of a set of discretization points (geometric entities)
5022:   and a discretization object that defines the function space associated with those points.

5024:   Logically Collective

5026:   Input Parameters:
5027: + dm    - The `DM`
5028: . label - The label indicating the support of the field, or `NULL` for the entire mesh
5029: - disc  - The discretization object

5031:   Level: intermediate

5033:   Notes:
5034:   The label already exists or will be added to the `DM` with `DMSetLabel()`.

5036:   For example, a piecewise continuous pressure field can be defined by coefficients at the cell centers of a mesh and piecewise constant functions
5037:   within each cell. Thus a specific function in the space is defined by the combination of a `Vec` containing the coefficients, a `DM` defining the
5038:   geometry entities, a `DMLabel` indicating a subset of those geometric entities, and a discretization object, such as a `PetscFE`.

5040: .seealso: [](ch_dmbase), `DM`, `DMSetLabel()`, `DMSetField()`, `DMGetField()`, `PetscFE`
5041: @*/
5042: PetscErrorCode DMAddField(DM dm, DMLabel label, PetscObject disc)
5043: {
5044:   PetscInt Nf = dm->Nf;

5046:   PetscFunctionBegin;
5050:   PetscCall(DMFieldEnlarge_Static(dm, Nf + 1));
5051:   dm->fields[Nf].label = label;
5052:   dm->fields[Nf].disc  = disc;
5053:   PetscCall(PetscObjectReference((PetscObject)label));
5054:   PetscCall(PetscObjectReference(disc));
5055:   PetscCall(DMSetDefaultAdjacency_Private(dm, Nf, disc));
5056:   PetscCall(DMClearDS(dm));
5057:   PetscFunctionReturn(PETSC_SUCCESS);
5058: }

5060: /*@
5061:   DMSetFieldAvoidTensor - Set flag to avoid defining the field on tensor cells

5063:   Logically Collective

5065:   Input Parameters:
5066: + dm          - The `DM`
5067: . f           - The field index
5068: - avoidTensor - `PETSC_TRUE` to skip defining the field on tensor cells

5070:   Level: intermediate

5072: .seealso: [](ch_dmbase), `DM`, `DMGetFieldAvoidTensor()`, `DMSetField()`, `DMGetField()`
5073: @*/
5074: PetscErrorCode DMSetFieldAvoidTensor(DM dm, PetscInt f, PetscBool avoidTensor)
5075: {
5076:   PetscFunctionBegin;
5077:   PetscCheck((f >= 0) && (f < dm->Nf), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", f, dm->Nf);
5078:   dm->fields[f].avoidTensor = avoidTensor;
5079:   PetscFunctionReturn(PETSC_SUCCESS);
5080: }

5082: /*@
5083:   DMGetFieldAvoidTensor - Get flag to avoid defining the field on tensor cells

5085:   Not Collective

5087:   Input Parameters:
5088: + dm - The `DM`
5089: - f  - The field index

5091:   Output Parameter:
5092: . avoidTensor - The flag to avoid defining the field on tensor cells

5094:   Level: intermediate

5096: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMSetField()`, `DMGetField()`, `DMSetFieldAvoidTensor()`
5097: @*/
5098: PetscErrorCode DMGetFieldAvoidTensor(DM dm, PetscInt f, PetscBool *avoidTensor)
5099: {
5100:   PetscFunctionBegin;
5101:   PetscCheck((f >= 0) && (f < dm->Nf), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", f, dm->Nf);
5102:   *avoidTensor = dm->fields[f].avoidTensor;
5103:   PetscFunctionReturn(PETSC_SUCCESS);
5104: }

5106: /*@
5107:   DMCopyFields - Copy the discretizations for the `DM` into another `DM`

5109:   Collective

5111:   Input Parameters:
5112: + dm        - The `DM`
5113: . minDegree - Minimum degree for a discretization, or `PETSC_DETERMINE` for no limit
5114: - maxDegree - Maximum degree for a discretization, or `PETSC_DETERMINE` for no limit

5116:   Output Parameter:
5117: . newdm - The `DM`

5119:   Level: advanced

5121: .seealso: [](ch_dmbase), `DM`, `DMGetField()`, `DMSetField()`, `DMAddField()`, `DMCopyDS()`, `DMGetDS()`, `DMGetCellDS()`
5122: @*/
5123: PetscErrorCode DMCopyFields(DM dm, PetscInt minDegree, PetscInt maxDegree, DM newdm)
5124: {
5125:   PetscInt Nf, f;

5127:   PetscFunctionBegin;
5128:   if (dm == newdm) PetscFunctionReturn(PETSC_SUCCESS);
5129:   PetscCall(DMGetNumFields(dm, &Nf));
5130:   PetscCall(DMClearFields(newdm));
5131:   for (f = 0; f < Nf; ++f) {
5132:     DMLabel      label;
5133:     PetscObject  field;
5134:     PetscClassId id;
5135:     PetscBool    useCone, useClosure;

5137:     PetscCall(DMGetField(dm, f, &label, &field));
5138:     PetscCall(PetscObjectGetClassId(field, &id));
5139:     if (id == PETSCFE_CLASSID) {
5140:       PetscFE newfe;

5142:       PetscCall(PetscFELimitDegree((PetscFE)field, minDegree, maxDegree, &newfe));
5143:       PetscCall(DMSetField(newdm, f, label, (PetscObject)newfe));
5144:       PetscCall(PetscFEDestroy(&newfe));
5145:     } else {
5146:       PetscCall(DMSetField(newdm, f, label, field));
5147:     }
5148:     PetscCall(DMGetAdjacency(dm, f, &useCone, &useClosure));
5149:     PetscCall(DMSetAdjacency(newdm, f, useCone, useClosure));
5150:   }
5151:   PetscFunctionReturn(PETSC_SUCCESS);
5152: }

5154: /*@
5155:   DMGetAdjacency - Returns the flags for determining variable influence

5157:   Not Collective

5159:   Input Parameters:
5160: + dm - The `DM` object
5161: - f  - The field number, or `PETSC_DEFAULT` for the default adjacency

5163:   Output Parameters:
5164: + useCone    - Flag for variable influence starting with the cone operation
5165: - useClosure - Flag for variable influence using transitive closure

5167:   Level: developer

5169:   Notes:
5170: .vb
5171:      FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
5172:      FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
5173:      FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
5174: .ve
5175:   Further explanation can be found in the User's Manual Section on the Influence of Variables on One Another.

5177: .seealso: [](ch_dmbase), `DM`, `DMSetAdjacency()`, `DMGetField()`, `DMSetField()`
5178: @*/
5179: PetscErrorCode DMGetAdjacency(DM dm, PetscInt f, PetscBool *useCone, PetscBool *useClosure)
5180: {
5181:   PetscFunctionBegin;
5183:   if (useCone) PetscAssertPointer(useCone, 3);
5184:   if (useClosure) PetscAssertPointer(useClosure, 4);
5185:   if (f < 0) {
5186:     if (useCone) *useCone = dm->adjacency[0];
5187:     if (useClosure) *useClosure = dm->adjacency[1];
5188:   } else {
5189:     PetscInt Nf;

5191:     PetscCall(DMGetNumFields(dm, &Nf));
5192:     PetscCheck(f < Nf, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, Nf);
5193:     if (useCone) *useCone = dm->fields[f].adjacency[0];
5194:     if (useClosure) *useClosure = dm->fields[f].adjacency[1];
5195:   }
5196:   PetscFunctionReturn(PETSC_SUCCESS);
5197: }

5199: /*@
5200:   DMSetAdjacency - Set the flags for determining variable influence

5202:   Not Collective

5204:   Input Parameters:
5205: + dm         - The `DM` object
5206: . f          - The field number
5207: . useCone    - Flag for variable influence starting with the cone operation
5208: - useClosure - Flag for variable influence using transitive closure

5210:   Level: developer

5212:   Notes:
5213: .vb
5214:      FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
5215:      FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
5216:      FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
5217: .ve
5218:   Further explanation can be found in the User's Manual Section on the Influence of Variables on One Another.

5220: .seealso: [](ch_dmbase), `DM`, `DMGetAdjacency()`, `DMGetField()`, `DMSetField()`
5221: @*/
5222: PetscErrorCode DMSetAdjacency(DM dm, PetscInt f, PetscBool useCone, PetscBool useClosure)
5223: {
5224:   PetscFunctionBegin;
5226:   if (f < 0) {
5227:     dm->adjacency[0] = useCone;
5228:     dm->adjacency[1] = useClosure;
5229:   } else {
5230:     PetscInt Nf;

5232:     PetscCall(DMGetNumFields(dm, &Nf));
5233:     PetscCheck(f < Nf, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, Nf);
5234:     dm->fields[f].adjacency[0] = useCone;
5235:     dm->fields[f].adjacency[1] = useClosure;
5236:   }
5237:   PetscFunctionReturn(PETSC_SUCCESS);
5238: }

5240: /*@
5241:   DMGetBasicAdjacency - Returns the flags for determining variable influence, using either the default or field 0 if it is defined

5243:   Not collective

5245:   Input Parameter:
5246: . dm - The `DM` object

5248:   Output Parameters:
5249: + useCone    - Flag for variable influence starting with the cone operation
5250: - useClosure - Flag for variable influence using transitive closure

5252:   Level: developer

5254:   Notes:
5255: .vb
5256:      FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
5257:      FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
5258:      FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
5259: .ve

5261: .seealso: [](ch_dmbase), `DM`, `DMSetBasicAdjacency()`, `DMGetField()`, `DMSetField()`
5262: @*/
5263: PetscErrorCode DMGetBasicAdjacency(DM dm, PetscBool *useCone, PetscBool *useClosure)
5264: {
5265:   PetscInt Nf;

5267:   PetscFunctionBegin;
5269:   if (useCone) PetscAssertPointer(useCone, 2);
5270:   if (useClosure) PetscAssertPointer(useClosure, 3);
5271:   PetscCall(DMGetNumFields(dm, &Nf));
5272:   if (!Nf) {
5273:     PetscCall(DMGetAdjacency(dm, PETSC_DEFAULT, useCone, useClosure));
5274:   } else {
5275:     PetscCall(DMGetAdjacency(dm, 0, useCone, useClosure));
5276:   }
5277:   PetscFunctionReturn(PETSC_SUCCESS);
5278: }

5280: /*@
5281:   DMSetBasicAdjacency - Set the flags for determining variable influence, using either the default or field 0 if it is defined

5283:   Not Collective

5285:   Input Parameters:
5286: + dm         - The `DM` object
5287: . useCone    - Flag for variable influence starting with the cone operation
5288: - useClosure - Flag for variable influence using transitive closure

5290:   Level: developer

5292:   Notes:
5293: .vb
5294:      FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
5295:      FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
5296:      FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
5297: .ve

5299: .seealso: [](ch_dmbase), `DM`, `DMGetBasicAdjacency()`, `DMGetField()`, `DMSetField()`
5300: @*/
5301: PetscErrorCode DMSetBasicAdjacency(DM dm, PetscBool useCone, PetscBool useClosure)
5302: {
5303:   PetscInt Nf;

5305:   PetscFunctionBegin;
5307:   PetscCall(DMGetNumFields(dm, &Nf));
5308:   if (!Nf) {
5309:     PetscCall(DMSetAdjacency(dm, PETSC_DEFAULT, useCone, useClosure));
5310:   } else {
5311:     PetscCall(DMSetAdjacency(dm, 0, useCone, useClosure));
5312:   }
5313:   PetscFunctionReturn(PETSC_SUCCESS);
5314: }

5316: PetscErrorCode DMCompleteBCLabels_Internal(DM dm)
5317: {
5318:   DM           plex;
5319:   DMLabel     *labels, *glabels;
5320:   const char **names;
5321:   char        *sendNames, *recvNames;
5322:   PetscInt     Nds, s, maxLabels = 0, maxLen = 0, gmaxLen, Nl = 0, gNl, l, gl, m;
5323:   size_t       len;
5324:   MPI_Comm     comm;
5325:   PetscMPIInt  rank, size, p, *counts, *displs;

5327:   PetscFunctionBegin;
5328:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
5329:   PetscCallMPI(MPI_Comm_size(comm, &size));
5330:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5331:   PetscCall(DMGetNumDS(dm, &Nds));
5332:   for (s = 0; s < Nds; ++s) {
5333:     PetscDS  dsBC;
5334:     PetscInt numBd;

5336:     PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC, NULL));
5337:     PetscCall(PetscDSGetNumBoundary(dsBC, &numBd));
5338:     maxLabels += numBd;
5339:   }
5340:   PetscCall(PetscCalloc1(maxLabels, &labels));
5341:   /* Get list of labels to be completed */
5342:   for (s = 0; s < Nds; ++s) {
5343:     PetscDS  dsBC;
5344:     PetscInt numBd, bd;

5346:     PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC, NULL));
5347:     PetscCall(PetscDSGetNumBoundary(dsBC, &numBd));
5348:     for (bd = 0; bd < numBd; ++bd) {
5349:       DMLabel      label;
5350:       PetscInt     field;
5351:       PetscObject  obj;
5352:       PetscClassId id;

5354:       PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, NULL, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL));
5355:       PetscCall(DMGetField(dm, field, NULL, &obj));
5356:       PetscCall(PetscObjectGetClassId(obj, &id));
5357:       if (!(id == PETSCFE_CLASSID) || !label) continue;
5358:       for (l = 0; l < Nl; ++l)
5359:         if (labels[l] == label) break;
5360:       if (l == Nl) labels[Nl++] = label;
5361:     }
5362:   }
5363:   /* Get label names */
5364:   PetscCall(PetscMalloc1(Nl, &names));
5365:   for (l = 0; l < Nl; ++l) PetscCall(PetscObjectGetName((PetscObject)labels[l], &names[l]));
5366:   for (l = 0; l < Nl; ++l) {
5367:     PetscCall(PetscStrlen(names[l], &len));
5368:     maxLen = PetscMax(maxLen, (PetscInt)len + 2);
5369:   }
5370:   PetscCall(PetscFree(labels));
5371:   PetscCallMPI(MPIU_Allreduce(&maxLen, &gmaxLen, 1, MPIU_INT, MPI_MAX, comm));
5372:   PetscCall(PetscCalloc1(Nl * gmaxLen, &sendNames));
5373:   for (l = 0; l < Nl; ++l) PetscCall(PetscStrncpy(&sendNames[gmaxLen * l], names[l], gmaxLen));
5374:   PetscCall(PetscFree(names));
5375:   /* Put all names on all processes */
5376:   PetscCall(PetscCalloc2(size, &counts, size + 1, &displs));
5377:   PetscCallMPI(MPI_Allgather(&Nl, 1, MPI_INT, counts, 1, MPI_INT, comm));
5378:   for (p = 0; p < size; ++p) displs[p + 1] = displs[p] + counts[p];
5379:   gNl = displs[size];
5380:   for (p = 0; p < size; ++p) {
5381:     counts[p] *= gmaxLen;
5382:     displs[p] *= gmaxLen;
5383:   }
5384:   PetscCall(PetscCalloc2(gNl * gmaxLen, &recvNames, gNl, &glabels));
5385:   PetscCallMPI(MPI_Allgatherv(sendNames, counts[rank], MPI_CHAR, recvNames, counts, displs, MPI_CHAR, comm));
5386:   PetscCall(PetscFree2(counts, displs));
5387:   PetscCall(PetscFree(sendNames));
5388:   for (l = 0, gl = 0; l < gNl; ++l) {
5389:     PetscCall(DMGetLabel(dm, &recvNames[l * gmaxLen], &glabels[gl]));
5390:     PetscCheck(glabels[gl], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Label %s missing on rank %d", &recvNames[l * gmaxLen], rank);
5391:     for (m = 0; m < gl; ++m)
5392:       if (glabels[m] == glabels[gl]) goto next_label;
5393:     PetscCall(DMConvert(dm, DMPLEX, &plex));
5394:     PetscCall(DMPlexLabelComplete(plex, glabels[gl]));
5395:     PetscCall(DMDestroy(&plex));
5396:     ++gl;
5397:   next_label:
5398:     continue;
5399:   }
5400:   PetscCall(PetscFree2(recvNames, glabels));
5401:   PetscFunctionReturn(PETSC_SUCCESS);
5402: }

5404: static PetscErrorCode DMDSEnlarge_Static(DM dm, PetscInt NdsNew)
5405: {
5406:   DMSpace *tmpd;
5407:   PetscInt Nds = dm->Nds, s;

5409:   PetscFunctionBegin;
5410:   if (Nds >= NdsNew) PetscFunctionReturn(PETSC_SUCCESS);
5411:   PetscCall(PetscMalloc1(NdsNew, &tmpd));
5412:   for (s = 0; s < Nds; ++s) tmpd[s] = dm->probs[s];
5413:   for (s = Nds; s < NdsNew; ++s) {
5414:     tmpd[s].ds     = NULL;
5415:     tmpd[s].label  = NULL;
5416:     tmpd[s].fields = NULL;
5417:   }
5418:   PetscCall(PetscFree(dm->probs));
5419:   dm->Nds   = NdsNew;
5420:   dm->probs = tmpd;
5421:   PetscFunctionReturn(PETSC_SUCCESS);
5422: }

5424: /*@
5425:   DMGetNumDS - Get the number of discrete systems in the `DM`

5427:   Not Collective

5429:   Input Parameter:
5430: . dm - The `DM`

5432:   Output Parameter:
5433: . Nds - The number of `PetscDS` objects

5435:   Level: intermediate

5437: .seealso: [](ch_dmbase), `DM`, `DMGetDS()`, `DMGetCellDS()`
5438: @*/
5439: PetscErrorCode DMGetNumDS(DM dm, PetscInt *Nds)
5440: {
5441:   PetscFunctionBegin;
5443:   PetscAssertPointer(Nds, 2);
5444:   *Nds = dm->Nds;
5445:   PetscFunctionReturn(PETSC_SUCCESS);
5446: }

5448: /*@
5449:   DMClearDS - Remove all discrete systems from the `DM`

5451:   Logically Collective

5453:   Input Parameter:
5454: . dm - The `DM`

5456:   Level: intermediate

5458: .seealso: [](ch_dmbase), `DM`, `DMGetNumDS()`, `DMGetDS()`, `DMSetField()`
5459: @*/
5460: PetscErrorCode DMClearDS(DM dm)
5461: {
5462:   PetscInt s;

5464:   PetscFunctionBegin;
5466:   for (s = 0; s < dm->Nds; ++s) {
5467:     PetscCall(PetscDSDestroy(&dm->probs[s].ds));
5468:     PetscCall(PetscDSDestroy(&dm->probs[s].dsIn));
5469:     PetscCall(DMLabelDestroy(&dm->probs[s].label));
5470:     PetscCall(ISDestroy(&dm->probs[s].fields));
5471:   }
5472:   PetscCall(PetscFree(dm->probs));
5473:   dm->probs = NULL;
5474:   dm->Nds   = 0;
5475:   PetscFunctionReturn(PETSC_SUCCESS);
5476: }

5478: /*@
5479:   DMGetDS - Get the default `PetscDS`

5481:   Not Collective

5483:   Input Parameter:
5484: . dm - The `DM`

5486:   Output Parameter:
5487: . ds - The default `PetscDS`

5489:   Level: intermediate

5491: .seealso: [](ch_dmbase), `DM`, `DMGetCellDS()`, `DMGetRegionDS()`
5492: @*/
5493: PetscErrorCode DMGetDS(DM dm, PetscDS *ds)
5494: {
5495:   PetscFunctionBeginHot;
5497:   PetscAssertPointer(ds, 2);
5498:   PetscCheck(dm->Nds > 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Need to call DMCreateDS() before calling DMGetDS()");
5499:   *ds = dm->probs[0].ds;
5500:   PetscFunctionReturn(PETSC_SUCCESS);
5501: }

5503: /*@
5504:   DMGetCellDS - Get the `PetscDS` defined on a given cell

5506:   Not Collective

5508:   Input Parameters:
5509: + dm    - The `DM`
5510: - point - Cell for the `PetscDS`

5512:   Output Parameters:
5513: + ds   - The `PetscDS` defined on the given cell
5514: - dsIn - The `PetscDS` for input on the given cell, or NULL if the same ds

5516:   Level: developer

5518: .seealso: [](ch_dmbase), `DM`, `DMGetDS()`, `DMSetRegionDS()`
5519: @*/
5520: PetscErrorCode DMGetCellDS(DM dm, PetscInt point, PetscDS *ds, PetscDS *dsIn)
5521: {
5522:   PetscDS  dsDef = NULL;
5523:   PetscInt s;

5525:   PetscFunctionBeginHot;
5527:   if (ds) PetscAssertPointer(ds, 3);
5528:   if (dsIn) PetscAssertPointer(dsIn, 4);
5529:   PetscCheck(point >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Mesh point cannot be negative: %" PetscInt_FMT, point);
5530:   if (ds) *ds = NULL;
5531:   if (dsIn) *dsIn = NULL;
5532:   for (s = 0; s < dm->Nds; ++s) {
5533:     PetscInt val;

5535:     if (!dm->probs[s].label) {
5536:       dsDef = dm->probs[s].ds;
5537:     } else {
5538:       PetscCall(DMLabelGetValue(dm->probs[s].label, point, &val));
5539:       if (val >= 0) {
5540:         if (ds) *ds = dm->probs[s].ds;
5541:         if (dsIn) *dsIn = dm->probs[s].dsIn;
5542:         break;
5543:       }
5544:     }
5545:   }
5546:   if (ds && !*ds) *ds = dsDef;
5547:   PetscFunctionReturn(PETSC_SUCCESS);
5548: }

5550: /*@
5551:   DMGetRegionDS - Get the `PetscDS` for a given mesh region, defined by a `DMLabel`

5553:   Not Collective

5555:   Input Parameters:
5556: + dm    - The `DM`
5557: - label - The `DMLabel` defining the mesh region, or `NULL` for the entire mesh

5559:   Output Parameters:
5560: + fields - The `IS` containing the `DM` field numbers for the fields in this `PetscDS`, or `NULL`
5561: . ds     - The `PetscDS` defined on the given region, or `NULL`
5562: - dsIn   - The `PetscDS` for input in the given region, or `NULL`

5564:   Level: advanced

5566:   Note:
5567:   If a non-`NULL` label is given, but there is no `PetscDS` on that specific label,
5568:   the `PetscDS` for the full domain (if present) is returned. Returns with
5569:   fields = `NULL` and ds = `NULL` if there is no `PetscDS` for the full domain.

5571: .seealso: [](ch_dmbase), `DM`, `DMGetRegionNumDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5572: @*/
5573: PetscErrorCode DMGetRegionDS(DM dm, DMLabel label, IS *fields, PetscDS *ds, PetscDS *dsIn)
5574: {
5575:   PetscInt Nds = dm->Nds, s;

5577:   PetscFunctionBegin;
5580:   if (fields) {
5581:     PetscAssertPointer(fields, 3);
5582:     *fields = NULL;
5583:   }
5584:   if (ds) {
5585:     PetscAssertPointer(ds, 4);
5586:     *ds = NULL;
5587:   }
5588:   if (dsIn) {
5589:     PetscAssertPointer(dsIn, 5);
5590:     *dsIn = NULL;
5591:   }
5592:   for (s = 0; s < Nds; ++s) {
5593:     if (dm->probs[s].label == label || !dm->probs[s].label) {
5594:       if (fields) *fields = dm->probs[s].fields;
5595:       if (ds) *ds = dm->probs[s].ds;
5596:       if (dsIn) *dsIn = dm->probs[s].dsIn;
5597:       if (dm->probs[s].label) PetscFunctionReturn(PETSC_SUCCESS);
5598:     }
5599:   }
5600:   PetscFunctionReturn(PETSC_SUCCESS);
5601: }

5603: /*@
5604:   DMSetRegionDS - Set the `PetscDS` for a given mesh region, defined by a `DMLabel`

5606:   Collective

5608:   Input Parameters:
5609: + dm     - The `DM`
5610: . label  - The `DMLabel` defining the mesh region, or `NULL` for the entire mesh
5611: . fields - The `IS` containing the `DM` field numbers for the fields in this `PetscDS`, or `NULL` for all fields
5612: . ds     - The `PetscDS` defined on the given region
5613: - dsIn   - The `PetscDS` for input on the given cell, or `NULL` if it is the same `PetscDS`

5615:   Level: advanced

5617:   Note:
5618:   If the label has a `PetscDS` defined, it will be replaced. Otherwise, it will be added to the `DM`. If the `PetscDS` is replaced,
5619:   the fields argument is ignored.

5621: .seealso: [](ch_dmbase), `DM`, `DMGetRegionDS()`, `DMSetRegionNumDS()`, `DMGetDS()`, `DMGetCellDS()`
5622: @*/
5623: PetscErrorCode DMSetRegionDS(DM dm, DMLabel label, IS fields, PetscDS ds, PetscDS dsIn)
5624: {
5625:   PetscInt Nds = dm->Nds, s;

5627:   PetscFunctionBegin;
5633:   for (s = 0; s < Nds; ++s) {
5634:     if (dm->probs[s].label == label) {
5635:       PetscCall(PetscDSDestroy(&dm->probs[s].ds));
5636:       PetscCall(PetscDSDestroy(&dm->probs[s].dsIn));
5637:       dm->probs[s].ds   = ds;
5638:       dm->probs[s].dsIn = dsIn;
5639:       PetscFunctionReturn(PETSC_SUCCESS);
5640:     }
5641:   }
5642:   PetscCall(DMDSEnlarge_Static(dm, Nds + 1));
5643:   PetscCall(PetscObjectReference((PetscObject)label));
5644:   PetscCall(PetscObjectReference((PetscObject)fields));
5645:   PetscCall(PetscObjectReference((PetscObject)ds));
5646:   PetscCall(PetscObjectReference((PetscObject)dsIn));
5647:   if (!label) {
5648:     /* Put the NULL label at the front, so it is returned as the default */
5649:     for (s = Nds - 1; s >= 0; --s) dm->probs[s + 1] = dm->probs[s];
5650:     Nds = 0;
5651:   }
5652:   dm->probs[Nds].label  = label;
5653:   dm->probs[Nds].fields = fields;
5654:   dm->probs[Nds].ds     = ds;
5655:   dm->probs[Nds].dsIn   = dsIn;
5656:   PetscFunctionReturn(PETSC_SUCCESS);
5657: }

5659: /*@
5660:   DMGetRegionNumDS - Get the `PetscDS` for a given mesh region, defined by the region number

5662:   Not Collective

5664:   Input Parameters:
5665: + dm  - The `DM`
5666: - num - The region number, in [0, Nds)

5668:   Output Parameters:
5669: + label  - The region label, or `NULL`
5670: . fields - The `IS` containing the `DM` field numbers for the fields in this `PetscDS`, or `NULL`
5671: . ds     - The `PetscDS` defined on the given region, or `NULL`
5672: - dsIn   - The `PetscDS` for input in the given region, or `NULL`

5674:   Level: advanced

5676: .seealso: [](ch_dmbase), `DM`, `DMGetRegionDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5677: @*/
5678: PetscErrorCode DMGetRegionNumDS(DM dm, PetscInt num, DMLabel *label, IS *fields, PetscDS *ds, PetscDS *dsIn)
5679: {
5680:   PetscInt Nds;

5682:   PetscFunctionBegin;
5684:   PetscCall(DMGetNumDS(dm, &Nds));
5685:   PetscCheck((num >= 0) && (num < Nds), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Region number %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", num, Nds);
5686:   if (label) {
5687:     PetscAssertPointer(label, 3);
5688:     *label = dm->probs[num].label;
5689:   }
5690:   if (fields) {
5691:     PetscAssertPointer(fields, 4);
5692:     *fields = dm->probs[num].fields;
5693:   }
5694:   if (ds) {
5695:     PetscAssertPointer(ds, 5);
5696:     *ds = dm->probs[num].ds;
5697:   }
5698:   if (dsIn) {
5699:     PetscAssertPointer(dsIn, 6);
5700:     *dsIn = dm->probs[num].dsIn;
5701:   }
5702:   PetscFunctionReturn(PETSC_SUCCESS);
5703: }

5705: /*@
5706:   DMSetRegionNumDS - Set the `PetscDS` for a given mesh region, defined by the region number

5708:   Not Collective

5710:   Input Parameters:
5711: + dm     - The `DM`
5712: . num    - The region number, in [0, Nds)
5713: . label  - The region label, or `NULL`
5714: . fields - The `IS` containing the `DM` field numbers for the fields in this `PetscDS`, or `NULL` to prevent setting
5715: . ds     - The `PetscDS` defined on the given region, or `NULL` to prevent setting
5716: - dsIn   - The `PetscDS` for input on the given cell, or `NULL` if it is the same `PetscDS`

5718:   Level: advanced

5720: .seealso: [](ch_dmbase), `DM`, `DMGetRegionDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5721: @*/
5722: PetscErrorCode DMSetRegionNumDS(DM dm, PetscInt num, DMLabel label, IS fields, PetscDS ds, PetscDS dsIn)
5723: {
5724:   PetscInt Nds;

5726:   PetscFunctionBegin;
5729:   PetscCall(DMGetNumDS(dm, &Nds));
5730:   PetscCheck((num >= 0) && (num < Nds), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Region number %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", num, Nds);
5731:   PetscCall(PetscObjectReference((PetscObject)label));
5732:   PetscCall(DMLabelDestroy(&dm->probs[num].label));
5733:   dm->probs[num].label = label;
5734:   if (fields) {
5736:     PetscCall(PetscObjectReference((PetscObject)fields));
5737:     PetscCall(ISDestroy(&dm->probs[num].fields));
5738:     dm->probs[num].fields = fields;
5739:   }
5740:   if (ds) {
5742:     PetscCall(PetscObjectReference((PetscObject)ds));
5743:     PetscCall(PetscDSDestroy(&dm->probs[num].ds));
5744:     dm->probs[num].ds = ds;
5745:   }
5746:   if (dsIn) {
5748:     PetscCall(PetscObjectReference((PetscObject)dsIn));
5749:     PetscCall(PetscDSDestroy(&dm->probs[num].dsIn));
5750:     dm->probs[num].dsIn = dsIn;
5751:   }
5752:   PetscFunctionReturn(PETSC_SUCCESS);
5753: }

5755: /*@
5756:   DMFindRegionNum - Find the region number for a given `PetscDS`, or -1 if it is not found.

5758:   Not Collective

5760:   Input Parameters:
5761: + dm - The `DM`
5762: - ds - The `PetscDS` defined on the given region

5764:   Output Parameter:
5765: . num - The region number, in [0, Nds), or -1 if not found

5767:   Level: advanced

5769: .seealso: [](ch_dmbase), `DM`, `DMGetRegionNumDS()`, `DMGetRegionDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5770: @*/
5771: PetscErrorCode DMFindRegionNum(DM dm, PetscDS ds, PetscInt *num)
5772: {
5773:   PetscInt Nds, n;

5775:   PetscFunctionBegin;
5778:   PetscAssertPointer(num, 3);
5779:   PetscCall(DMGetNumDS(dm, &Nds));
5780:   for (n = 0; n < Nds; ++n)
5781:     if (ds == dm->probs[n].ds) break;
5782:   if (n >= Nds) *num = -1;
5783:   else *num = n;
5784:   PetscFunctionReturn(PETSC_SUCCESS);
5785: }

5787: /*@
5788:   DMCreateFEDefault - Create a `PetscFE` based on the celltype for the mesh

5790:   Not Collective

5792:   Input Parameters:
5793: + dm     - The `DM`
5794: . Nc     - The number of components for the field
5795: . prefix - The options prefix for the output `PetscFE`, or `NULL`
5796: - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree

5798:   Output Parameter:
5799: . fem - The `PetscFE`

5801:   Level: intermediate

5803:   Note:
5804:   This is a convenience method that just calls `PetscFECreateByCell()` underneath.

5806: .seealso: [](ch_dmbase), `DM`, `PetscFECreateByCell()`, `DMAddField()`, `DMCreateDS()`, `DMGetCellDS()`, `DMGetRegionDS()`
5807: @*/
5808: PetscErrorCode DMCreateFEDefault(DM dm, PetscInt Nc, const char prefix[], PetscInt qorder, PetscFE *fem)
5809: {
5810:   DMPolytopeType ct;
5811:   PetscInt       dim, cStart;

5813:   PetscFunctionBegin;
5816:   if (prefix) PetscAssertPointer(prefix, 3);
5818:   PetscAssertPointer(fem, 5);
5819:   PetscCall(DMGetDimension(dm, &dim));
5820:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
5821:   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
5822:   PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, Nc, ct, prefix, qorder, fem));
5823:   PetscFunctionReturn(PETSC_SUCCESS);
5824: }

5826: /*@
5827:   DMCreateDS - Create the discrete systems for the `DM` based upon the fields added to the `DM`

5829:   Collective

5831:   Input Parameter:
5832: . dm - The `DM`

5834:   Options Database Key:
5835: . -dm_petscds_view - View all the `PetscDS` objects in this `DM`

5837:   Level: intermediate

5839: .seealso: [](ch_dmbase), `DM`, `DMSetField`, `DMAddField()`, `DMGetDS()`, `DMGetCellDS()`, `DMGetRegionDS()`, `DMSetRegionDS()`
5840: @*/
5841: PetscErrorCode DMCreateDS(DM dm)
5842: {
5843:   MPI_Comm  comm;
5844:   PetscDS   dsDef;
5845:   DMLabel  *labelSet;
5846:   PetscInt  dE, Nf = dm->Nf, f, s, Nl, l, Ndef, k;
5847:   PetscBool doSetup = PETSC_TRUE, flg;

5849:   PetscFunctionBegin;
5851:   if (!dm->fields) PetscFunctionReturn(PETSC_SUCCESS);
5852:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
5853:   PetscCall(DMGetCoordinateDim(dm, &dE));
5854:   /* Determine how many regions we have */
5855:   PetscCall(PetscMalloc1(Nf, &labelSet));
5856:   Nl   = 0;
5857:   Ndef = 0;
5858:   for (f = 0; f < Nf; ++f) {
5859:     DMLabel  label = dm->fields[f].label;
5860:     PetscInt l;

5862: #ifdef PETSC_HAVE_LIBCEED
5863:     /* Move CEED context to discretizations */
5864:     {
5865:       PetscClassId id;

5867:       PetscCall(PetscObjectGetClassId(dm->fields[f].disc, &id));
5868:       if (id == PETSCFE_CLASSID) {
5869:         Ceed ceed;

5871:         PetscCall(DMGetCeed(dm, &ceed));
5872:         PetscCall(PetscFESetCeed((PetscFE)dm->fields[f].disc, ceed));
5873:       }
5874:     }
5875: #endif
5876:     if (!label) {
5877:       ++Ndef;
5878:       continue;
5879:     }
5880:     for (l = 0; l < Nl; ++l)
5881:       if (label == labelSet[l]) break;
5882:     if (l < Nl) continue;
5883:     labelSet[Nl++] = label;
5884:   }
5885:   /* Create default DS if there are no labels to intersect with */
5886:   PetscCall(DMGetRegionDS(dm, NULL, NULL, &dsDef, NULL));
5887:   if (!dsDef && Ndef && !Nl) {
5888:     IS        fields;
5889:     PetscInt *fld, nf;

5891:     for (f = 0, nf = 0; f < Nf; ++f)
5892:       if (!dm->fields[f].label) ++nf;
5893:     PetscCheck(nf, comm, PETSC_ERR_PLIB, "All fields have labels, but we are trying to create a default DS");
5894:     PetscCall(PetscMalloc1(nf, &fld));
5895:     for (f = 0, nf = 0; f < Nf; ++f)
5896:       if (!dm->fields[f].label) fld[nf++] = f;
5897:     PetscCall(ISCreate(PETSC_COMM_SELF, &fields));
5898:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fields, "dm_fields_"));
5899:     PetscCall(ISSetType(fields, ISGENERAL));
5900:     PetscCall(ISGeneralSetIndices(fields, nf, fld, PETSC_OWN_POINTER));

5902:     PetscCall(PetscDSCreate(PETSC_COMM_SELF, &dsDef));
5903:     PetscCall(DMSetRegionDS(dm, NULL, fields, dsDef, NULL));
5904:     PetscCall(PetscDSDestroy(&dsDef));
5905:     PetscCall(ISDestroy(&fields));
5906:   }
5907:   PetscCall(DMGetRegionDS(dm, NULL, NULL, &dsDef, NULL));
5908:   if (dsDef) PetscCall(PetscDSSetCoordinateDimension(dsDef, dE));
5909:   /* Intersect labels with default fields */
5910:   if (Ndef && Nl) {
5911:     DM              plex;
5912:     DMLabel         cellLabel;
5913:     IS              fieldIS, allcellIS, defcellIS = NULL;
5914:     PetscInt       *fields;
5915:     const PetscInt *cells;
5916:     PetscInt        depth, nf = 0, n, c;

5918:     PetscCall(DMConvert(dm, DMPLEX, &plex));
5919:     PetscCall(DMPlexGetDepth(plex, &depth));
5920:     PetscCall(DMGetStratumIS(plex, "dim", depth, &allcellIS));
5921:     if (!allcellIS) PetscCall(DMGetStratumIS(plex, "depth", depth, &allcellIS));
5922:     /* TODO This looks like it only works for one label */
5923:     for (l = 0; l < Nl; ++l) {
5924:       DMLabel label = labelSet[l];
5925:       IS      pointIS;

5927:       PetscCall(ISDestroy(&defcellIS));
5928:       PetscCall(DMLabelGetStratumIS(label, 1, &pointIS));
5929:       PetscCall(ISDifference(allcellIS, pointIS, &defcellIS));
5930:       PetscCall(ISDestroy(&pointIS));
5931:     }
5932:     PetscCall(ISDestroy(&allcellIS));

5934:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "defaultCells", &cellLabel));
5935:     PetscCall(ISGetLocalSize(defcellIS, &n));
5936:     PetscCall(ISGetIndices(defcellIS, &cells));
5937:     for (c = 0; c < n; ++c) PetscCall(DMLabelSetValue(cellLabel, cells[c], 1));
5938:     PetscCall(ISRestoreIndices(defcellIS, &cells));
5939:     PetscCall(ISDestroy(&defcellIS));
5940:     PetscCall(DMPlexLabelComplete(plex, cellLabel));

5942:     PetscCall(PetscMalloc1(Ndef, &fields));
5943:     for (f = 0; f < Nf; ++f)
5944:       if (!dm->fields[f].label) fields[nf++] = f;
5945:     PetscCall(ISCreate(PETSC_COMM_SELF, &fieldIS));
5946:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fieldIS, "dm_fields_"));
5947:     PetscCall(ISSetType(fieldIS, ISGENERAL));
5948:     PetscCall(ISGeneralSetIndices(fieldIS, nf, fields, PETSC_OWN_POINTER));

5950:     PetscCall(PetscDSCreate(PETSC_COMM_SELF, &dsDef));
5951:     PetscCall(DMSetRegionDS(dm, cellLabel, fieldIS, dsDef, NULL));
5952:     PetscCall(PetscDSSetCoordinateDimension(dsDef, dE));
5953:     PetscCall(DMLabelDestroy(&cellLabel));
5954:     PetscCall(PetscDSDestroy(&dsDef));
5955:     PetscCall(ISDestroy(&fieldIS));
5956:     PetscCall(DMDestroy(&plex));
5957:   }
5958:   /* Create label DSes
5959:      - WE ONLY SUPPORT IDENTICAL OR DISJOINT LABELS
5960:   */
5961:   /* TODO Should check that labels are disjoint */
5962:   for (l = 0; l < Nl; ++l) {
5963:     DMLabel   label = labelSet[l];
5964:     PetscDS   ds, dsIn = NULL;
5965:     IS        fields;
5966:     PetscInt *fld, nf;

5968:     PetscCall(PetscDSCreate(PETSC_COMM_SELF, &ds));
5969:     for (f = 0, nf = 0; f < Nf; ++f)
5970:       if (label == dm->fields[f].label || !dm->fields[f].label) ++nf;
5971:     PetscCall(PetscMalloc1(nf, &fld));
5972:     for (f = 0, nf = 0; f < Nf; ++f)
5973:       if (label == dm->fields[f].label || !dm->fields[f].label) fld[nf++] = f;
5974:     PetscCall(ISCreate(PETSC_COMM_SELF, &fields));
5975:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fields, "dm_fields_"));
5976:     PetscCall(ISSetType(fields, ISGENERAL));
5977:     PetscCall(ISGeneralSetIndices(fields, nf, fld, PETSC_OWN_POINTER));
5978:     PetscCall(PetscDSSetCoordinateDimension(ds, dE));
5979:     {
5980:       DMPolytopeType ct;
5981:       PetscInt       lStart, lEnd;
5982:       PetscBool      isCohesiveLocal = PETSC_FALSE, isCohesive;

5984:       PetscCall(DMLabelGetBounds(label, &lStart, &lEnd));
5985:       if (lStart >= 0) {
5986:         PetscCall(DMPlexGetCellType(dm, lStart, &ct));
5987:         switch (ct) {
5988:         case DM_POLYTOPE_POINT_PRISM_TENSOR:
5989:         case DM_POLYTOPE_SEG_PRISM_TENSOR:
5990:         case DM_POLYTOPE_TRI_PRISM_TENSOR:
5991:         case DM_POLYTOPE_QUAD_PRISM_TENSOR:
5992:           isCohesiveLocal = PETSC_TRUE;
5993:           break;
5994:         default:
5995:           break;
5996:         }
5997:       }
5998:       PetscCallMPI(MPIU_Allreduce(&isCohesiveLocal, &isCohesive, 1, MPIU_BOOL, MPI_LOR, comm));
5999:       if (isCohesive) {
6000:         PetscCall(PetscDSCreate(PETSC_COMM_SELF, &dsIn));
6001:         PetscCall(PetscDSSetCoordinateDimension(dsIn, dE));
6002:       }
6003:       for (f = 0, nf = 0; f < Nf; ++f) {
6004:         if (label == dm->fields[f].label || !dm->fields[f].label) {
6005:           if (label == dm->fields[f].label) {
6006:             PetscCall(PetscDSSetDiscretization(ds, nf, NULL));
6007:             PetscCall(PetscDSSetCohesive(ds, nf, isCohesive));
6008:             if (dsIn) {
6009:               PetscCall(PetscDSSetDiscretization(dsIn, nf, NULL));
6010:               PetscCall(PetscDSSetCohesive(dsIn, nf, isCohesive));
6011:             }
6012:           }
6013:           ++nf;
6014:         }
6015:       }
6016:     }
6017:     PetscCall(DMSetRegionDS(dm, label, fields, ds, dsIn));
6018:     PetscCall(ISDestroy(&fields));
6019:     PetscCall(PetscDSDestroy(&ds));
6020:     PetscCall(PetscDSDestroy(&dsIn));
6021:   }
6022:   PetscCall(PetscFree(labelSet));
6023:   /* Set fields in DSes */
6024:   for (s = 0; s < dm->Nds; ++s) {
6025:     PetscDS         ds     = dm->probs[s].ds;
6026:     PetscDS         dsIn   = dm->probs[s].dsIn;
6027:     IS              fields = dm->probs[s].fields;
6028:     const PetscInt *fld;
6029:     PetscInt        nf, dsnf;
6030:     PetscBool       isCohesive;

6032:     PetscCall(PetscDSGetNumFields(ds, &dsnf));
6033:     PetscCall(PetscDSIsCohesive(ds, &isCohesive));
6034:     PetscCall(ISGetLocalSize(fields, &nf));
6035:     PetscCall(ISGetIndices(fields, &fld));
6036:     for (f = 0; f < nf; ++f) {
6037:       PetscObject  disc = dm->fields[fld[f]].disc;
6038:       PetscBool    isCohesiveField;
6039:       PetscClassId id;

6041:       /* Handle DS with no fields */
6042:       if (dsnf) PetscCall(PetscDSGetCohesive(ds, f, &isCohesiveField));
6043:       /* If this is a cohesive cell, then regular fields need the lower dimensional discretization */
6044:       if (isCohesive) {
6045:         if (!isCohesiveField) {
6046:           PetscObject bdDisc;

6048:           PetscCall(PetscFEGetHeightSubspace((PetscFE)disc, 1, (PetscFE *)&bdDisc));
6049:           PetscCall(PetscDSSetDiscretization(ds, f, bdDisc));
6050:           PetscCall(PetscDSSetDiscretization(dsIn, f, disc));
6051:         } else {
6052:           PetscCall(PetscDSSetDiscretization(ds, f, disc));
6053:           PetscCall(PetscDSSetDiscretization(dsIn, f, disc));
6054:         }
6055:       } else {
6056:         PetscCall(PetscDSSetDiscretization(ds, f, disc));
6057:       }
6058:       /* We allow people to have placeholder fields and construct the Section by hand */
6059:       PetscCall(PetscObjectGetClassId(disc, &id));
6060:       if ((id != PETSCFE_CLASSID) && (id != PETSCFV_CLASSID)) doSetup = PETSC_FALSE;
6061:     }
6062:     PetscCall(ISRestoreIndices(fields, &fld));
6063:   }
6064:   /* Allow k-jet tabulation */
6065:   PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)dm)->prefix, "-dm_ds_jet_degree", &k, &flg));
6066:   if (flg) {
6067:     for (s = 0; s < dm->Nds; ++s) {
6068:       PetscDS  ds   = dm->probs[s].ds;
6069:       PetscDS  dsIn = dm->probs[s].dsIn;
6070:       PetscInt Nf, f;

6072:       PetscCall(PetscDSGetNumFields(ds, &Nf));
6073:       for (f = 0; f < Nf; ++f) {
6074:         PetscCall(PetscDSSetJetDegree(ds, f, k));
6075:         if (dsIn) PetscCall(PetscDSSetJetDegree(dsIn, f, k));
6076:       }
6077:     }
6078:   }
6079:   /* Setup DSes */
6080:   if (doSetup) {
6081:     for (s = 0; s < dm->Nds; ++s) {
6082:       if (dm->setfromoptionscalled) {
6083:         PetscCall(PetscDSSetFromOptions(dm->probs[s].ds));
6084:         if (dm->probs[s].dsIn) PetscCall(PetscDSSetFromOptions(dm->probs[s].dsIn));
6085:       }
6086:       PetscCall(PetscDSSetUp(dm->probs[s].ds));
6087:       if (dm->probs[s].dsIn) PetscCall(PetscDSSetUp(dm->probs[s].dsIn));
6088:     }
6089:   }
6090:   PetscFunctionReturn(PETSC_SUCCESS);
6091: }

6093: /*@
6094:   DMUseTensorOrder - Use a tensor product closure ordering for the default section

6096:   Input Parameters:
6097: + dm     - The DM
6098: - tensor - Flag for tensor order

6100:   Level: developer

6102: .seealso: `DMPlexSetClosurePermutationTensor()`, `PetscSectionResetClosurePermutation()`
6103: @*/
6104: PetscErrorCode DMUseTensorOrder(DM dm, PetscBool tensor)
6105: {
6106:   PetscInt  Nf;
6107:   PetscBool reorder = PETSC_TRUE, isPlex;

6109:   PetscFunctionBegin;
6110:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex));
6111:   PetscCall(DMGetNumFields(dm, &Nf));
6112:   for (PetscInt f = 0; f < Nf; ++f) {
6113:     PetscObject  obj;
6114:     PetscClassId id;

6116:     PetscCall(DMGetField(dm, f, NULL, &obj));
6117:     PetscCall(PetscObjectGetClassId(obj, &id));
6118:     if (id == PETSCFE_CLASSID) {
6119:       PetscSpace sp;
6120:       PetscBool  tensor;

6122:       PetscCall(PetscFEGetBasisSpace((PetscFE)obj, &sp));
6123:       PetscCall(PetscSpacePolynomialGetTensor(sp, &tensor));
6124:       reorder = reorder && tensor ? PETSC_TRUE : PETSC_FALSE;
6125:     } else reorder = PETSC_FALSE;
6126:   }
6127:   if (tensor) {
6128:     if (reorder && isPlex) PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
6129:   } else {
6130:     PetscSection s;

6132:     PetscCall(DMGetLocalSection(dm, &s));
6133:     if (s) PetscCall(PetscSectionResetClosurePermutation(s));
6134:   }
6135:   PetscFunctionReturn(PETSC_SUCCESS);
6136: }

6138: /*@
6139:   DMComputeExactSolution - Compute the exact solution for a given `DM`, using the `PetscDS` information.

6141:   Collective

6143:   Input Parameters:
6144: + dm   - The `DM`
6145: - time - The time

6147:   Output Parameters:
6148: + u   - The vector will be filled with exact solution values, or `NULL`
6149: - u_t - The vector will be filled with the time derivative of exact solution values, or `NULL`

6151:   Level: developer

6153:   Note:
6154:   The user must call `PetscDSSetExactSolution()` before using this routine

6156: .seealso: [](ch_dmbase), `DM`, `PetscDSSetExactSolution()`
6157: @*/
6158: PetscErrorCode DMComputeExactSolution(DM dm, PetscReal time, Vec u, Vec u_t)
6159: {
6160:   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
6161:   void   **ectxs;
6162:   Vec      locu, locu_t;
6163:   PetscInt Nf, Nds, s;

6165:   PetscFunctionBegin;
6167:   if (u) {
6169:     PetscCall(DMGetLocalVector(dm, &locu));
6170:     PetscCall(VecSet(locu, 0.));
6171:   }
6172:   if (u_t) {
6174:     PetscCall(DMGetLocalVector(dm, &locu_t));
6175:     PetscCall(VecSet(locu_t, 0.));
6176:   }
6177:   PetscCall(DMGetNumFields(dm, &Nf));
6178:   PetscCall(PetscMalloc2(Nf, &exacts, Nf, &ectxs));
6179:   PetscCall(DMGetNumDS(dm, &Nds));
6180:   for (s = 0; s < Nds; ++s) {
6181:     PetscDS         ds;
6182:     DMLabel         label;
6183:     IS              fieldIS;
6184:     const PetscInt *fields, id = 1;
6185:     PetscInt        dsNf, f;

6187:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL));
6188:     PetscCall(PetscDSGetNumFields(ds, &dsNf));
6189:     PetscCall(ISGetIndices(fieldIS, &fields));
6190:     PetscCall(PetscArrayzero(exacts, Nf));
6191:     PetscCall(PetscArrayzero(ectxs, Nf));
6192:     if (u) {
6193:       for (f = 0; f < dsNf; ++f) PetscCall(PetscDSGetExactSolution(ds, fields[f], &exacts[fields[f]], &ectxs[fields[f]]));
6194:       if (label) PetscCall(DMProjectFunctionLabelLocal(dm, time, label, 1, &id, 0, NULL, exacts, ectxs, INSERT_ALL_VALUES, locu));
6195:       else PetscCall(DMProjectFunctionLocal(dm, time, exacts, ectxs, INSERT_ALL_VALUES, locu));
6196:     }
6197:     if (u_t) {
6198:       PetscCall(PetscArrayzero(exacts, Nf));
6199:       PetscCall(PetscArrayzero(ectxs, Nf));
6200:       for (f = 0; f < dsNf; ++f) PetscCall(PetscDSGetExactSolutionTimeDerivative(ds, fields[f], &exacts[fields[f]], &ectxs[fields[f]]));
6201:       if (label) PetscCall(DMProjectFunctionLabelLocal(dm, time, label, 1, &id, 0, NULL, exacts, ectxs, INSERT_ALL_VALUES, locu_t));
6202:       else PetscCall(DMProjectFunctionLocal(dm, time, exacts, ectxs, INSERT_ALL_VALUES, locu_t));
6203:     }
6204:     PetscCall(ISRestoreIndices(fieldIS, &fields));
6205:   }
6206:   if (u) {
6207:     PetscCall(PetscObjectSetName((PetscObject)u, "Exact Solution"));
6208:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)u, "exact_"));
6209:   }
6210:   if (u_t) {
6211:     PetscCall(PetscObjectSetName((PetscObject)u, "Exact Solution Time Derivative"));
6212:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)u_t, "exact_t_"));
6213:   }
6214:   PetscCall(PetscFree2(exacts, ectxs));
6215:   if (u) {
6216:     PetscCall(DMLocalToGlobalBegin(dm, locu, INSERT_ALL_VALUES, u));
6217:     PetscCall(DMLocalToGlobalEnd(dm, locu, INSERT_ALL_VALUES, u));
6218:     PetscCall(DMRestoreLocalVector(dm, &locu));
6219:   }
6220:   if (u_t) {
6221:     PetscCall(DMLocalToGlobalBegin(dm, locu_t, INSERT_ALL_VALUES, u_t));
6222:     PetscCall(DMLocalToGlobalEnd(dm, locu_t, INSERT_ALL_VALUES, u_t));
6223:     PetscCall(DMRestoreLocalVector(dm, &locu_t));
6224:   }
6225:   PetscFunctionReturn(PETSC_SUCCESS);
6226: }

6228: static PetscErrorCode DMTransferDS_Internal(DM dm, DMLabel label, IS fields, PetscInt minDegree, PetscInt maxDegree, PetscDS ds, PetscDS dsIn)
6229: {
6230:   PetscDS dsNew, dsInNew = NULL;

6232:   PetscFunctionBegin;
6233:   PetscCall(PetscDSCreate(PetscObjectComm((PetscObject)ds), &dsNew));
6234:   PetscCall(PetscDSCopy(ds, minDegree, maxDegree, dm, dsNew));
6235:   if (dsIn) {
6236:     PetscCall(PetscDSCreate(PetscObjectComm((PetscObject)dsIn), &dsInNew));
6237:     PetscCall(PetscDSCopy(dsIn, minDegree, maxDegree, dm, dsInNew));
6238:   }
6239:   PetscCall(DMSetRegionDS(dm, label, fields, dsNew, dsInNew));
6240:   PetscCall(PetscDSDestroy(&dsNew));
6241:   PetscCall(PetscDSDestroy(&dsInNew));
6242:   PetscFunctionReturn(PETSC_SUCCESS);
6243: }

6245: /*@
6246:   DMCopyDS - Copy the discrete systems for the `DM` into another `DM`

6248:   Collective

6250:   Input Parameters:
6251: + dm        - The `DM`
6252: . minDegree - Minimum degree for a discretization, or `PETSC_DETERMINE` for no limit
6253: - maxDegree - Maximum degree for a discretization, or `PETSC_DETERMINE` for no limit

6255:   Output Parameter:
6256: . newdm - The `DM`

6258:   Level: advanced

6260: .seealso: [](ch_dmbase), `DM`, `DMCopyFields()`, `DMAddField()`, `DMGetDS()`, `DMGetCellDS()`, `DMGetRegionDS()`, `DMSetRegionDS()`
6261: @*/
6262: PetscErrorCode DMCopyDS(DM dm, PetscInt minDegree, PetscInt maxDegree, DM newdm)
6263: {
6264:   PetscInt Nds, s;

6266:   PetscFunctionBegin;
6267:   if (dm == newdm) PetscFunctionReturn(PETSC_SUCCESS);
6268:   PetscCall(DMGetNumDS(dm, &Nds));
6269:   PetscCall(DMClearDS(newdm));
6270:   for (s = 0; s < Nds; ++s) {
6271:     DMLabel  label;
6272:     IS       fields;
6273:     PetscDS  ds, dsIn, newds;
6274:     PetscInt Nbd, bd;

6276:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fields, &ds, &dsIn));
6277:     /* TODO: We need to change all keys from labels in the old DM to labels in the new DM */
6278:     PetscCall(DMTransferDS_Internal(newdm, label, fields, minDegree, maxDegree, ds, dsIn));
6279:     /* Complete new labels in the new DS */
6280:     PetscCall(DMGetRegionDS(newdm, label, NULL, &newds, NULL));
6281:     PetscCall(PetscDSGetNumBoundary(newds, &Nbd));
6282:     for (bd = 0; bd < Nbd; ++bd) {
6283:       PetscWeakForm wf;
6284:       DMLabel       label;
6285:       PetscInt      field;

6287:       PetscCall(PetscDSGetBoundary(newds, bd, &wf, NULL, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL));
6288:       PetscCall(PetscWeakFormReplaceLabel(wf, label));
6289:     }
6290:   }
6291:   PetscCall(DMCompleteBCLabels_Internal(newdm));
6292:   PetscFunctionReturn(PETSC_SUCCESS);
6293: }

6295: /*@
6296:   DMCopyDisc - Copy the fields and discrete systems for the `DM` into another `DM`

6298:   Collective

6300:   Input Parameter:
6301: . dm - The `DM`

6303:   Output Parameter:
6304: . newdm - The `DM`

6306:   Level: advanced

6308:   Developer Note:
6309:   Really ugly name, nothing in PETSc is called a `Disc` plus it is an ugly abbreviation

6311: .seealso: [](ch_dmbase), `DM`, `DMCopyFields()`, `DMCopyDS()`
6312: @*/
6313: PetscErrorCode DMCopyDisc(DM dm, DM newdm)
6314: {
6315:   PetscFunctionBegin;
6316:   PetscCall(DMCopyFields(dm, PETSC_DETERMINE, PETSC_DETERMINE, newdm));
6317:   PetscCall(DMCopyDS(dm, PETSC_DETERMINE, PETSC_DETERMINE, newdm));
6318:   PetscFunctionReturn(PETSC_SUCCESS);
6319: }

6321: /*@
6322:   DMGetDimension - Return the topological dimension of the `DM`

6324:   Not Collective

6326:   Input Parameter:
6327: . dm - The `DM`

6329:   Output Parameter:
6330: . dim - The topological dimension

6332:   Level: beginner

6334: .seealso: [](ch_dmbase), `DM`, `DMSetDimension()`, `DMCreate()`
6335: @*/
6336: PetscErrorCode DMGetDimension(DM dm, PetscInt *dim)
6337: {
6338:   PetscFunctionBegin;
6340:   PetscAssertPointer(dim, 2);
6341:   *dim = dm->dim;
6342:   PetscFunctionReturn(PETSC_SUCCESS);
6343: }

6345: /*@
6346:   DMSetDimension - Set the topological dimension of the `DM`

6348:   Collective

6350:   Input Parameters:
6351: + dm  - The `DM`
6352: - dim - The topological dimension

6354:   Level: beginner

6356: .seealso: [](ch_dmbase), `DM`, `DMGetDimension()`, `DMCreate()`
6357: @*/
6358: PetscErrorCode DMSetDimension(DM dm, PetscInt dim)
6359: {
6360:   PetscDS  ds;
6361:   PetscInt Nds, n;

6363:   PetscFunctionBegin;
6366:   dm->dim = dim;
6367:   if (dm->dim >= 0) {
6368:     PetscCall(DMGetNumDS(dm, &Nds));
6369:     for (n = 0; n < Nds; ++n) {
6370:       PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds, NULL));
6371:       if (ds->dimEmbed < 0) PetscCall(PetscDSSetCoordinateDimension(ds, dim));
6372:     }
6373:   }
6374:   PetscFunctionReturn(PETSC_SUCCESS);
6375: }

6377: /*@
6378:   DMGetDimPoints - Get the half-open interval for all points of a given dimension

6380:   Collective

6382:   Input Parameters:
6383: + dm  - the `DM`
6384: - dim - the dimension

6386:   Output Parameters:
6387: + pStart - The first point of the given dimension
6388: - pEnd   - The first point following points of the given dimension

6390:   Level: intermediate

6392:   Note:
6393:   The points are vertices in the Hasse diagram encoding the topology. This is explained in
6394:   https://arxiv.org/abs/0908.4427. If no points exist of this dimension in the storage scheme,
6395:   then the interval is empty.

6397: .seealso: [](ch_dmbase), `DM`, `DMPLEX`, `DMPlexGetDepthStratum()`, `DMPlexGetHeightStratum()`
6398: @*/
6399: PetscErrorCode DMGetDimPoints(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
6400: {
6401:   PetscInt d;

6403:   PetscFunctionBegin;
6405:   PetscCall(DMGetDimension(dm, &d));
6406:   PetscCheck((dim >= 0) && (dim <= d), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim);
6407:   PetscUseTypeMethod(dm, getdimpoints, dim, pStart, pEnd);
6408:   PetscFunctionReturn(PETSC_SUCCESS);
6409: }

6411: /*@
6412:   DMGetOutputDM - Retrieve the `DM` associated with the layout for output

6414:   Collective

6416:   Input Parameter:
6417: . dm - The original `DM`

6419:   Output Parameter:
6420: . odm - The `DM` which provides the layout for output

6422:   Level: intermediate

6424:   Note:
6425:   In some situations the vector obtained with `DMCreateGlobalVector()` excludes points for degrees of freedom that are associated with fixed (Dirichelet) boundary
6426:   conditions since the algebraic solver does not solve for those variables. The output `DM` includes these excluded points and its global vector contains the
6427:   locations for those dof so that they can be output to a file or other viewer along with the unconstrained dof.

6429: .seealso: [](ch_dmbase), `DM`, `VecView()`, `DMGetGlobalSection()`, `DMCreateGlobalVector()`, `PetscSectionHasConstraints()`, `DMSetGlobalSection()`
6430: @*/
6431: PetscErrorCode DMGetOutputDM(DM dm, DM *odm)
6432: {
6433:   PetscSection section;
6434:   IS           perm;
6435:   PetscBool    hasConstraints, newDM, gnewDM;

6437:   PetscFunctionBegin;
6439:   PetscAssertPointer(odm, 2);
6440:   PetscCall(DMGetLocalSection(dm, &section));
6441:   PetscCall(PetscSectionHasConstraints(section, &hasConstraints));
6442:   PetscCall(PetscSectionGetPermutation(section, &perm));
6443:   newDM = hasConstraints || perm ? PETSC_TRUE : PETSC_FALSE;
6444:   PetscCallMPI(MPIU_Allreduce(&newDM, &gnewDM, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)dm)));
6445:   if (!gnewDM) {
6446:     *odm = dm;
6447:     PetscFunctionReturn(PETSC_SUCCESS);
6448:   }
6449:   if (!dm->dmBC) {
6450:     PetscSection newSection, gsection;
6451:     PetscSF      sf;
6452:     PetscBool    usePerm = dm->ignorePermOutput ? PETSC_FALSE : PETSC_TRUE;

6454:     PetscCall(DMClone(dm, &dm->dmBC));
6455:     PetscCall(DMCopyDisc(dm, dm->dmBC));
6456:     PetscCall(PetscSectionClone(section, &newSection));
6457:     PetscCall(DMSetLocalSection(dm->dmBC, newSection));
6458:     PetscCall(PetscSectionDestroy(&newSection));
6459:     PetscCall(DMGetPointSF(dm->dmBC, &sf));
6460:     PetscCall(PetscSectionCreateGlobalSection(section, sf, usePerm, PETSC_TRUE, PETSC_FALSE, &gsection));
6461:     PetscCall(DMSetGlobalSection(dm->dmBC, gsection));
6462:     PetscCall(PetscSectionDestroy(&gsection));
6463:   }
6464:   *odm = dm->dmBC;
6465:   PetscFunctionReturn(PETSC_SUCCESS);
6466: }

6468: /*@
6469:   DMGetOutputSequenceNumber - Retrieve the sequence number/value for output

6471:   Input Parameter:
6472: . dm - The original `DM`

6474:   Output Parameters:
6475: + num - The output sequence number
6476: - val - The output sequence value

6478:   Level: intermediate

6480:   Note:
6481:   This is intended for output that should appear in sequence, for instance
6482:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6484:   Developer Note:
6485:   The `DM` serves as a convenient place to store the current iteration value. The iteration is not
6486:   not directly related to the `DM`.

6488: .seealso: [](ch_dmbase), `DM`, `VecView()`
6489: @*/
6490: PetscErrorCode DMGetOutputSequenceNumber(DM dm, PetscInt *num, PetscReal *val)
6491: {
6492:   PetscFunctionBegin;
6494:   if (num) {
6495:     PetscAssertPointer(num, 2);
6496:     *num = dm->outputSequenceNum;
6497:   }
6498:   if (val) {
6499:     PetscAssertPointer(val, 3);
6500:     *val = dm->outputSequenceVal;
6501:   }
6502:   PetscFunctionReturn(PETSC_SUCCESS);
6503: }

6505: /*@
6506:   DMSetOutputSequenceNumber - Set the sequence number/value for output

6508:   Input Parameters:
6509: + dm  - The original `DM`
6510: . num - The output sequence number
6511: - val - The output sequence value

6513:   Level: intermediate

6515:   Note:
6516:   This is intended for output that should appear in sequence, for instance
6517:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6519: .seealso: [](ch_dmbase), `DM`, `VecView()`
6520: @*/
6521: PetscErrorCode DMSetOutputSequenceNumber(DM dm, PetscInt num, PetscReal val)
6522: {
6523:   PetscFunctionBegin;
6525:   dm->outputSequenceNum = num;
6526:   dm->outputSequenceVal = val;
6527:   PetscFunctionReturn(PETSC_SUCCESS);
6528: }

6530: /*@
6531:   DMOutputSequenceLoad - Retrieve the sequence value from a `PetscViewer`

6533:   Input Parameters:
6534: + dm     - The original `DM`
6535: . viewer - The `PetscViewer` to get it from
6536: . name   - The sequence name
6537: - num    - The output sequence number

6539:   Output Parameter:
6540: . val - The output sequence value

6542:   Level: intermediate

6544:   Note:
6545:   This is intended for output that should appear in sequence, for instance
6546:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6548:   Developer Note:
6549:   It is unclear at the user API level why a `DM` is needed as input

6551: .seealso: [](ch_dmbase), `DM`, `DMGetOutputSequenceNumber()`, `DMSetOutputSequenceNumber()`, `VecView()`
6552: @*/
6553: PetscErrorCode DMOutputSequenceLoad(DM dm, PetscViewer viewer, const char name[], PetscInt num, PetscReal *val)
6554: {
6555:   PetscBool ishdf5;

6557:   PetscFunctionBegin;
6560:   PetscAssertPointer(name, 3);
6561:   PetscAssertPointer(val, 5);
6562:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
6563:   if (ishdf5) {
6564: #if defined(PETSC_HAVE_HDF5)
6565:     PetscScalar value;

6567:     PetscCall(DMSequenceLoad_HDF5_Internal(dm, name, num, &value, viewer));
6568:     *val = PetscRealPart(value);
6569: #endif
6570:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerHDF5Open()");
6571:   PetscFunctionReturn(PETSC_SUCCESS);
6572: }

6574: /*@
6575:   DMGetOutputSequenceLength - Retrieve the number of sequence values from a `PetscViewer`

6577:   Input Parameters:
6578: + dm     - The original `DM`
6579: . viewer - The `PetscViewer` to get it from
6580: - name   - The sequence name

6582:   Output Parameter:
6583: . len - The length of the output sequence

6585:   Level: intermediate

6587:   Note:
6588:   This is intended for output that should appear in sequence, for instance
6589:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6591:   Developer Note:
6592:   It is unclear at the user API level why a `DM` is needed as input

6594: .seealso: [](ch_dmbase), `DM`, `DMGetOutputSequenceNumber()`, `DMSetOutputSequenceNumber()`, `VecView()`
6595: @*/
6596: PetscErrorCode DMGetOutputSequenceLength(DM dm, PetscViewer viewer, const char name[], PetscInt *len)
6597: {
6598:   PetscBool ishdf5;

6600:   PetscFunctionBegin;
6603:   PetscAssertPointer(name, 3);
6604:   PetscAssertPointer(len, 4);
6605:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
6606:   if (ishdf5) {
6607: #if defined(PETSC_HAVE_HDF5)
6608:     PetscCall(DMSequenceGetLength_HDF5_Internal(dm, name, len, viewer));
6609: #endif
6610:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerHDF5Open()");
6611:   PetscFunctionReturn(PETSC_SUCCESS);
6612: }

6614: /*@
6615:   DMGetUseNatural - Get the flag for creating a mapping to the natural order when a `DM` is (re)distributed in parallel

6617:   Not Collective

6619:   Input Parameter:
6620: . dm - The `DM`

6622:   Output Parameter:
6623: . useNatural - `PETSC_TRUE` to build the mapping to a natural order during distribution

6625:   Level: beginner

6627: .seealso: [](ch_dmbase), `DM`, `DMSetUseNatural()`, `DMCreate()`
6628: @*/
6629: PetscErrorCode DMGetUseNatural(DM dm, PetscBool *useNatural)
6630: {
6631:   PetscFunctionBegin;
6633:   PetscAssertPointer(useNatural, 2);
6634:   *useNatural = dm->useNatural;
6635:   PetscFunctionReturn(PETSC_SUCCESS);
6636: }

6638: /*@
6639:   DMSetUseNatural - Set the flag for creating a mapping to the natural order when a `DM` is (re)distributed in parallel

6641:   Collective

6643:   Input Parameters:
6644: + dm         - The `DM`
6645: - useNatural - `PETSC_TRUE` to build the mapping to a natural order during distribution

6647:   Level: beginner

6649:   Note:
6650:   This also causes the map to be build after `DMCreateSubDM()` and `DMCreateSuperDM()`

6652: .seealso: [](ch_dmbase), `DM`, `DMGetUseNatural()`, `DMCreate()`, `DMPlexDistribute()`, `DMCreateSubDM()`, `DMCreateSuperDM()`
6653: @*/
6654: PetscErrorCode DMSetUseNatural(DM dm, PetscBool useNatural)
6655: {
6656:   PetscFunctionBegin;
6659:   dm->useNatural = useNatural;
6660:   PetscFunctionReturn(PETSC_SUCCESS);
6661: }

6663: /*@
6664:   DMCreateLabel - Create a label of the given name if it does not already exist in the `DM`

6666:   Not Collective

6668:   Input Parameters:
6669: + dm   - The `DM` object
6670: - name - The label name

6672:   Level: intermediate

6674: .seealso: [](ch_dmbase), `DM`, `DMLabelCreate()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6675: @*/
6676: PetscErrorCode DMCreateLabel(DM dm, const char name[])
6677: {
6678:   PetscBool flg;
6679:   DMLabel   label;

6681:   PetscFunctionBegin;
6683:   PetscAssertPointer(name, 2);
6684:   PetscCall(DMHasLabel(dm, name, &flg));
6685:   if (!flg) {
6686:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, &label));
6687:     PetscCall(DMAddLabel(dm, label));
6688:     PetscCall(DMLabelDestroy(&label));
6689:   }
6690:   PetscFunctionReturn(PETSC_SUCCESS);
6691: }

6693: /*@
6694:   DMCreateLabelAtIndex - Create a label of the given name at the given index. If it already exists in the `DM`, move it to this index.

6696:   Not Collective

6698:   Input Parameters:
6699: + dm   - The `DM` object
6700: . l    - The index for the label
6701: - name - The label name

6703:   Level: intermediate

6705: .seealso: [](ch_dmbase), `DM`, `DMCreateLabel()`, `DMLabelCreate()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6706: @*/
6707: PetscErrorCode DMCreateLabelAtIndex(DM dm, PetscInt l, const char name[])
6708: {
6709:   DMLabelLink orig, prev = NULL;
6710:   DMLabel     label;
6711:   PetscInt    Nl, m;
6712:   PetscBool   flg, match;
6713:   const char *lname;

6715:   PetscFunctionBegin;
6717:   PetscAssertPointer(name, 3);
6718:   PetscCall(DMHasLabel(dm, name, &flg));
6719:   if (!flg) {
6720:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, &label));
6721:     PetscCall(DMAddLabel(dm, label));
6722:     PetscCall(DMLabelDestroy(&label));
6723:   }
6724:   PetscCall(DMGetNumLabels(dm, &Nl));
6725:   PetscCheck(l < Nl, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label index %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", l, Nl);
6726:   for (m = 0, orig = dm->labels; m < Nl; ++m, prev = orig, orig = orig->next) {
6727:     PetscCall(PetscObjectGetName((PetscObject)orig->label, &lname));
6728:     PetscCall(PetscStrcmp(name, lname, &match));
6729:     if (match) break;
6730:   }
6731:   if (m == l) PetscFunctionReturn(PETSC_SUCCESS);
6732:   if (!m) dm->labels = orig->next;
6733:   else prev->next = orig->next;
6734:   if (!l) {
6735:     orig->next = dm->labels;
6736:     dm->labels = orig;
6737:   } else {
6738:     for (m = 0, prev = dm->labels; m < l - 1; ++m, prev = prev->next);
6739:     orig->next = prev->next;
6740:     prev->next = orig;
6741:   }
6742:   PetscFunctionReturn(PETSC_SUCCESS);
6743: }

6745: /*@
6746:   DMGetLabelValue - Get the value in a `DMLabel` for the given point, with -1 as the default

6748:   Not Collective

6750:   Input Parameters:
6751: + dm    - The `DM` object
6752: . name  - The label name
6753: - point - The mesh point

6755:   Output Parameter:
6756: . value - The label value for this point, or -1 if the point is not in the label

6758:   Level: beginner

6760: .seealso: [](ch_dmbase), `DM`, `DMLabelGetValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6761: @*/
6762: PetscErrorCode DMGetLabelValue(DM dm, const char name[], PetscInt point, PetscInt *value)
6763: {
6764:   DMLabel label;

6766:   PetscFunctionBegin;
6768:   PetscAssertPointer(name, 2);
6769:   PetscCall(DMGetLabel(dm, name, &label));
6770:   PetscCheck(label, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No label named %s was found", name);
6771:   PetscCall(DMLabelGetValue(label, point, value));
6772:   PetscFunctionReturn(PETSC_SUCCESS);
6773: }

6775: /*@
6776:   DMSetLabelValue - Add a point to a `DMLabel` with given value

6778:   Not Collective

6780:   Input Parameters:
6781: + dm    - The `DM` object
6782: . name  - The label name
6783: . point - The mesh point
6784: - value - The label value for this point

6786:   Output Parameter:

6788:   Level: beginner

6790: .seealso: [](ch_dmbase), `DM`, `DMLabelSetValue()`, `DMGetStratumIS()`, `DMClearLabelValue()`
6791: @*/
6792: PetscErrorCode DMSetLabelValue(DM dm, const char name[], PetscInt point, PetscInt value)
6793: {
6794:   DMLabel label;

6796:   PetscFunctionBegin;
6798:   PetscAssertPointer(name, 2);
6799:   PetscCall(DMGetLabel(dm, name, &label));
6800:   if (!label) {
6801:     PetscCall(DMCreateLabel(dm, name));
6802:     PetscCall(DMGetLabel(dm, name, &label));
6803:   }
6804:   PetscCall(DMLabelSetValue(label, point, value));
6805:   PetscFunctionReturn(PETSC_SUCCESS);
6806: }

6808: /*@
6809:   DMClearLabelValue - Remove a point from a `DMLabel` with given value

6811:   Not Collective

6813:   Input Parameters:
6814: + dm    - The `DM` object
6815: . name  - The label name
6816: . point - The mesh point
6817: - value - The label value for this point

6819:   Level: beginner

6821: .seealso: [](ch_dmbase), `DM`, `DMLabelClearValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6822: @*/
6823: PetscErrorCode DMClearLabelValue(DM dm, const char name[], PetscInt point, PetscInt value)
6824: {
6825:   DMLabel label;

6827:   PetscFunctionBegin;
6829:   PetscAssertPointer(name, 2);
6830:   PetscCall(DMGetLabel(dm, name, &label));
6831:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6832:   PetscCall(DMLabelClearValue(label, point, value));
6833:   PetscFunctionReturn(PETSC_SUCCESS);
6834: }

6836: /*@
6837:   DMGetLabelSize - Get the value of `DMLabelGetNumValues()` of a `DMLabel` in the `DM`

6839:   Not Collective

6841:   Input Parameters:
6842: + dm   - The `DM` object
6843: - name - The label name

6845:   Output Parameter:
6846: . size - The number of different integer ids, or 0 if the label does not exist

6848:   Level: beginner

6850:   Developer Note:
6851:   This should be renamed to something like `DMGetLabelNumValues()` or removed.

6853: .seealso: [](ch_dmbase), `DM`, `DMLabelGetNumValues()`, `DMSetLabelValue()`, `DMGetLabel()`
6854: @*/
6855: PetscErrorCode DMGetLabelSize(DM dm, const char name[], PetscInt *size)
6856: {
6857:   DMLabel label;

6859:   PetscFunctionBegin;
6861:   PetscAssertPointer(name, 2);
6862:   PetscAssertPointer(size, 3);
6863:   PetscCall(DMGetLabel(dm, name, &label));
6864:   *size = 0;
6865:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6866:   PetscCall(DMLabelGetNumValues(label, size));
6867:   PetscFunctionReturn(PETSC_SUCCESS);
6868: }

6870: /*@
6871:   DMGetLabelIdIS - Get the `DMLabelGetValueIS()` from a `DMLabel` in the `DM`

6873:   Not Collective

6875:   Input Parameters:
6876: + dm   - The `DM` object
6877: - name - The label name

6879:   Output Parameter:
6880: . ids - The integer ids, or `NULL` if the label does not exist

6882:   Level: beginner

6884: .seealso: [](ch_dmbase), `DM`, `DMLabelGetValueIS()`, `DMGetLabelSize()`
6885: @*/
6886: PetscErrorCode DMGetLabelIdIS(DM dm, const char name[], IS *ids)
6887: {
6888:   DMLabel label;

6890:   PetscFunctionBegin;
6892:   PetscAssertPointer(name, 2);
6893:   PetscAssertPointer(ids, 3);
6894:   PetscCall(DMGetLabel(dm, name, &label));
6895:   *ids = NULL;
6896:   if (label) {
6897:     PetscCall(DMLabelGetValueIS(label, ids));
6898:   } else {
6899:     /* returning an empty IS */
6900:     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_USE_POINTER, ids));
6901:   }
6902:   PetscFunctionReturn(PETSC_SUCCESS);
6903: }

6905: /*@
6906:   DMGetStratumSize - Get the number of points in a label stratum

6908:   Not Collective

6910:   Input Parameters:
6911: + dm    - The `DM` object
6912: . name  - The label name of the stratum
6913: - value - The stratum value

6915:   Output Parameter:
6916: . size - The number of points, also called the stratum size

6918:   Level: beginner

6920: .seealso: [](ch_dmbase), `DM`, `DMLabelGetStratumSize()`, `DMGetLabelSize()`, `DMGetLabelIds()`
6921: @*/
6922: PetscErrorCode DMGetStratumSize(DM dm, const char name[], PetscInt value, PetscInt *size)
6923: {
6924:   DMLabel label;

6926:   PetscFunctionBegin;
6928:   PetscAssertPointer(name, 2);
6929:   PetscAssertPointer(size, 4);
6930:   PetscCall(DMGetLabel(dm, name, &label));
6931:   *size = 0;
6932:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6933:   PetscCall(DMLabelGetStratumSize(label, value, size));
6934:   PetscFunctionReturn(PETSC_SUCCESS);
6935: }

6937: /*@
6938:   DMGetStratumIS - Get the points in a label stratum

6940:   Not Collective

6942:   Input Parameters:
6943: + dm    - The `DM` object
6944: . name  - The label name
6945: - value - The stratum value

6947:   Output Parameter:
6948: . points - The stratum points, or `NULL` if the label does not exist or does not have that value

6950:   Level: beginner

6952: .seealso: [](ch_dmbase), `DM`, `DMLabelGetStratumIS()`, `DMGetStratumSize()`
6953: @*/
6954: PetscErrorCode DMGetStratumIS(DM dm, const char name[], PetscInt value, IS *points)
6955: {
6956:   DMLabel label;

6958:   PetscFunctionBegin;
6960:   PetscAssertPointer(name, 2);
6961:   PetscAssertPointer(points, 4);
6962:   PetscCall(DMGetLabel(dm, name, &label));
6963:   *points = NULL;
6964:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6965:   PetscCall(DMLabelGetStratumIS(label, value, points));
6966:   PetscFunctionReturn(PETSC_SUCCESS);
6967: }

6969: /*@
6970:   DMSetStratumIS - Set the points in a label stratum

6972:   Not Collective

6974:   Input Parameters:
6975: + dm     - The `DM` object
6976: . name   - The label name
6977: . value  - The stratum value
6978: - points - The stratum points

6980:   Level: beginner

6982: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMClearLabelStratum()`, `DMLabelClearStratum()`, `DMLabelSetStratumIS()`, `DMGetStratumSize()`
6983: @*/
6984: PetscErrorCode DMSetStratumIS(DM dm, const char name[], PetscInt value, IS points)
6985: {
6986:   DMLabel label;

6988:   PetscFunctionBegin;
6990:   PetscAssertPointer(name, 2);
6992:   PetscCall(DMGetLabel(dm, name, &label));
6993:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6994:   PetscCall(DMLabelSetStratumIS(label, value, points));
6995:   PetscFunctionReturn(PETSC_SUCCESS);
6996: }

6998: /*@
6999:   DMClearLabelStratum - Remove all points from a stratum from a `DMLabel`

7001:   Not Collective

7003:   Input Parameters:
7004: + dm    - The `DM` object
7005: . name  - The label name
7006: - value - The label value for this point

7008:   Output Parameter:

7010:   Level: beginner

7012: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMLabelClearStratum()`, `DMSetLabelValue()`, `DMGetStratumIS()`, `DMClearLabelValue()`
7013: @*/
7014: PetscErrorCode DMClearLabelStratum(DM dm, const char name[], PetscInt value)
7015: {
7016:   DMLabel label;

7018:   PetscFunctionBegin;
7020:   PetscAssertPointer(name, 2);
7021:   PetscCall(DMGetLabel(dm, name, &label));
7022:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
7023:   PetscCall(DMLabelClearStratum(label, value));
7024:   PetscFunctionReturn(PETSC_SUCCESS);
7025: }

7027: /*@
7028:   DMGetNumLabels - Return the number of labels defined by on the `DM`

7030:   Not Collective

7032:   Input Parameter:
7033: . dm - The `DM` object

7035:   Output Parameter:
7036: . numLabels - the number of Labels

7038:   Level: intermediate

7040: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetLabelByNum()`, `DMGetLabelName()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7041: @*/
7042: PetscErrorCode DMGetNumLabels(DM dm, PetscInt *numLabels)
7043: {
7044:   DMLabelLink next = dm->labels;
7045:   PetscInt    n    = 0;

7047:   PetscFunctionBegin;
7049:   PetscAssertPointer(numLabels, 2);
7050:   while (next) {
7051:     ++n;
7052:     next = next->next;
7053:   }
7054:   *numLabels = n;
7055:   PetscFunctionReturn(PETSC_SUCCESS);
7056: }

7058: /*@
7059:   DMGetLabelName - Return the name of nth label

7061:   Not Collective

7063:   Input Parameters:
7064: + dm - The `DM` object
7065: - n  - the label number

7067:   Output Parameter:
7068: . name - the label name

7070:   Level: intermediate

7072:   Developer Note:
7073:   Some of the functions that appropriate on labels using their number have the suffix ByNum, others do not.

7075: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetLabelByNum()`, `DMGetLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7076: @*/
7077: PetscErrorCode DMGetLabelName(DM dm, PetscInt n, const char *name[])
7078: {
7079:   DMLabelLink next = dm->labels;
7080:   PetscInt    l    = 0;

7082:   PetscFunctionBegin;
7084:   PetscAssertPointer(name, 3);
7085:   while (next) {
7086:     if (l == n) {
7087:       PetscCall(PetscObjectGetName((PetscObject)next->label, name));
7088:       PetscFunctionReturn(PETSC_SUCCESS);
7089:     }
7090:     ++l;
7091:     next = next->next;
7092:   }
7093:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %" PetscInt_FMT " does not exist in this DM", n);
7094: }

7096: /*@
7097:   DMHasLabel - Determine whether the `DM` has a label of a given name

7099:   Not Collective

7101:   Input Parameters:
7102: + dm   - The `DM` object
7103: - name - The label name

7105:   Output Parameter:
7106: . hasLabel - `PETSC_TRUE` if the label is present

7108:   Level: intermediate

7110: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetLabel()`, `DMGetLabelByNum()`, `DMCreateLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7111: @*/
7112: PetscErrorCode DMHasLabel(DM dm, const char name[], PetscBool *hasLabel)
7113: {
7114:   DMLabelLink next = dm->labels;
7115:   const char *lname;

7117:   PetscFunctionBegin;
7119:   PetscAssertPointer(name, 2);
7120:   PetscAssertPointer(hasLabel, 3);
7121:   *hasLabel = PETSC_FALSE;
7122:   while (next) {
7123:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7124:     PetscCall(PetscStrcmp(name, lname, hasLabel));
7125:     if (*hasLabel) break;
7126:     next = next->next;
7127:   }
7128:   PetscFunctionReturn(PETSC_SUCCESS);
7129: }

7131: // PetscClangLinter pragma ignore: -fdoc-section-header-unknown
7132: /*@
7133:   DMGetLabel - Return the label of a given name, or `NULL`, from a `DM`

7135:   Not Collective

7137:   Input Parameters:
7138: + dm   - The `DM` object
7139: - name - The label name

7141:   Output Parameter:
7142: . label - The `DMLabel`, or `NULL` if the label is absent

7144:   Default labels in a `DMPLEX`:
7145: + "depth"       - Holds the depth (co-dimension) of each mesh point
7146: . "celltype"    - Holds the topological type of each cell
7147: . "ghost"       - If the DM is distributed with overlap, this marks the cells and faces in the overlap
7148: . "Cell Sets"   - Mirrors the cell sets defined by GMsh and ExodusII
7149: . "Face Sets"   - Mirrors the face sets defined by GMsh and ExodusII
7150: - "Vertex Sets" - Mirrors the vertex sets defined by GMsh

7152:   Level: intermediate

7154: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMHasLabel()`, `DMGetLabelByNum()`, `DMAddLabel()`, `DMCreateLabel()`, `DMPlexGetDepthLabel()`, `DMPlexGetCellType()`
7155: @*/
7156: PetscErrorCode DMGetLabel(DM dm, const char name[], DMLabel *label)
7157: {
7158:   DMLabelLink next = dm->labels;
7159:   PetscBool   hasLabel;
7160:   const char *lname;

7162:   PetscFunctionBegin;
7164:   PetscAssertPointer(name, 2);
7165:   PetscAssertPointer(label, 3);
7166:   *label = NULL;
7167:   while (next) {
7168:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7169:     PetscCall(PetscStrcmp(name, lname, &hasLabel));
7170:     if (hasLabel) {
7171:       *label = next->label;
7172:       break;
7173:     }
7174:     next = next->next;
7175:   }
7176:   PetscFunctionReturn(PETSC_SUCCESS);
7177: }

7179: /*@
7180:   DMGetLabelByNum - Return the nth label on a `DM`

7182:   Not Collective

7184:   Input Parameters:
7185: + dm - The `DM` object
7186: - n  - the label number

7188:   Output Parameter:
7189: . label - the label

7191:   Level: intermediate

7193: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMAddLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7194: @*/
7195: PetscErrorCode DMGetLabelByNum(DM dm, PetscInt n, DMLabel *label)
7196: {
7197:   DMLabelLink next = dm->labels;
7198:   PetscInt    l    = 0;

7200:   PetscFunctionBegin;
7202:   PetscAssertPointer(label, 3);
7203:   while (next) {
7204:     if (l == n) {
7205:       *label = next->label;
7206:       PetscFunctionReturn(PETSC_SUCCESS);
7207:     }
7208:     ++l;
7209:     next = next->next;
7210:   }
7211:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %" PetscInt_FMT " does not exist in this DM", n);
7212: }

7214: /*@
7215:   DMAddLabel - Add the label to this `DM`

7217:   Not Collective

7219:   Input Parameters:
7220: + dm    - The `DM` object
7221: - label - The `DMLabel`

7223:   Level: developer

7225: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7226: @*/
7227: PetscErrorCode DMAddLabel(DM dm, DMLabel label)
7228: {
7229:   DMLabelLink l, *p, tmpLabel;
7230:   PetscBool   hasLabel;
7231:   const char *lname;
7232:   PetscBool   flg;

7234:   PetscFunctionBegin;
7236:   PetscCall(PetscObjectGetName((PetscObject)label, &lname));
7237:   PetscCall(DMHasLabel(dm, lname, &hasLabel));
7238:   PetscCheck(!hasLabel, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %s already exists in this DM", lname);
7239:   PetscCall(PetscCalloc1(1, &tmpLabel));
7240:   tmpLabel->label  = label;
7241:   tmpLabel->output = PETSC_TRUE;
7242:   for (p = &dm->labels; (l = *p); p = &l->next) { }
7243:   *p = tmpLabel;
7244:   PetscCall(PetscObjectReference((PetscObject)label));
7245:   PetscCall(PetscStrcmp(lname, "depth", &flg));
7246:   if (flg) dm->depthLabel = label;
7247:   PetscCall(PetscStrcmp(lname, "celltype", &flg));
7248:   if (flg) dm->celltypeLabel = label;
7249:   PetscFunctionReturn(PETSC_SUCCESS);
7250: }

7252: // PetscClangLinter pragma ignore: -fdoc-section-header-unknown
7253: /*@
7254:   DMSetLabel - Replaces the label of a given name, or ignores it if the name is not present

7256:   Not Collective

7258:   Input Parameters:
7259: + dm    - The `DM` object
7260: - label - The `DMLabel`, having the same name, to substitute

7262:   Default labels in a `DMPLEX`:
7263: + "depth"       - Holds the depth (co-dimension) of each mesh point
7264: . "celltype"    - Holds the topological type of each cell
7265: . "ghost"       - If the DM is distributed with overlap, this marks the cells and faces in the overlap
7266: . "Cell Sets"   - Mirrors the cell sets defined by GMsh and ExodusII
7267: . "Face Sets"   - Mirrors the face sets defined by GMsh and ExodusII
7268: - "Vertex Sets" - Mirrors the vertex sets defined by GMsh

7270:   Level: intermediate

7272: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMPlexGetDepthLabel()`, `DMPlexGetCellType()`
7273: @*/
7274: PetscErrorCode DMSetLabel(DM dm, DMLabel label)
7275: {
7276:   DMLabelLink next = dm->labels;
7277:   PetscBool   hasLabel, flg;
7278:   const char *name, *lname;

7280:   PetscFunctionBegin;
7283:   PetscCall(PetscObjectGetName((PetscObject)label, &name));
7284:   while (next) {
7285:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7286:     PetscCall(PetscStrcmp(name, lname, &hasLabel));
7287:     if (hasLabel) {
7288:       PetscCall(PetscObjectReference((PetscObject)label));
7289:       PetscCall(PetscStrcmp(lname, "depth", &flg));
7290:       if (flg) dm->depthLabel = label;
7291:       PetscCall(PetscStrcmp(lname, "celltype", &flg));
7292:       if (flg) dm->celltypeLabel = label;
7293:       PetscCall(DMLabelDestroy(&next->label));
7294:       next->label = label;
7295:       break;
7296:     }
7297:     next = next->next;
7298:   }
7299:   PetscFunctionReturn(PETSC_SUCCESS);
7300: }

7302: /*@
7303:   DMRemoveLabel - Remove the label given by name from this `DM`

7305:   Not Collective

7307:   Input Parameters:
7308: + dm   - The `DM` object
7309: - name - The label name

7311:   Output Parameter:
7312: . label - The `DMLabel`, or `NULL` if the label is absent. Pass in `NULL` to call `DMLabelDestroy()` on the label, otherwise the
7313:           caller is responsible for calling `DMLabelDestroy()`.

7315:   Level: developer

7317: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMLabelDestroy()`, `DMRemoveLabelBySelf()`
7318: @*/
7319: PetscErrorCode DMRemoveLabel(DM dm, const char name[], DMLabel *label)
7320: {
7321:   DMLabelLink link, *pnext;
7322:   PetscBool   hasLabel;
7323:   const char *lname;

7325:   PetscFunctionBegin;
7327:   PetscAssertPointer(name, 2);
7328:   if (label) {
7329:     PetscAssertPointer(label, 3);
7330:     *label = NULL;
7331:   }
7332:   for (pnext = &dm->labels; (link = *pnext); pnext = &link->next) {
7333:     PetscCall(PetscObjectGetName((PetscObject)link->label, &lname));
7334:     PetscCall(PetscStrcmp(name, lname, &hasLabel));
7335:     if (hasLabel) {
7336:       *pnext = link->next; /* Remove from list */
7337:       PetscCall(PetscStrcmp(name, "depth", &hasLabel));
7338:       if (hasLabel) dm->depthLabel = NULL;
7339:       PetscCall(PetscStrcmp(name, "celltype", &hasLabel));
7340:       if (hasLabel) dm->celltypeLabel = NULL;
7341:       if (label) *label = link->label;
7342:       else PetscCall(DMLabelDestroy(&link->label));
7343:       PetscCall(PetscFree(link));
7344:       break;
7345:     }
7346:   }
7347:   PetscFunctionReturn(PETSC_SUCCESS);
7348: }

7350: /*@
7351:   DMRemoveLabelBySelf - Remove the label from this `DM`

7353:   Not Collective

7355:   Input Parameters:
7356: + dm           - The `DM` object
7357: . label        - The `DMLabel` to be removed from the `DM`
7358: - failNotFound - Should it fail if the label is not found in the `DM`?

7360:   Level: developer

7362:   Note:
7363:   Only exactly the same instance is removed if found, name match is ignored.
7364:   If the `DM` has an exclusive reference to the label, the label gets destroyed and
7365:   *label nullified.

7367: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabel()` `DMGetLabelValue()`, `DMSetLabelValue()`, `DMLabelDestroy()`, `DMRemoveLabel()`
7368: @*/
7369: PetscErrorCode DMRemoveLabelBySelf(DM dm, DMLabel *label, PetscBool failNotFound)
7370: {
7371:   DMLabelLink link, *pnext;
7372:   PetscBool   hasLabel = PETSC_FALSE;

7374:   PetscFunctionBegin;
7376:   PetscAssertPointer(label, 2);
7377:   if (!*label && !failNotFound) PetscFunctionReturn(PETSC_SUCCESS);
7380:   for (pnext = &dm->labels; (link = *pnext); pnext = &link->next) {
7381:     if (*label == link->label) {
7382:       hasLabel = PETSC_TRUE;
7383:       *pnext   = link->next; /* Remove from list */
7384:       if (*label == dm->depthLabel) dm->depthLabel = NULL;
7385:       if (*label == dm->celltypeLabel) dm->celltypeLabel = NULL;
7386:       if (((PetscObject)link->label)->refct < 2) *label = NULL; /* nullify if exclusive reference */
7387:       PetscCall(DMLabelDestroy(&link->label));
7388:       PetscCall(PetscFree(link));
7389:       break;
7390:     }
7391:   }
7392:   PetscCheck(hasLabel || !failNotFound, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Given label not found in DM");
7393:   PetscFunctionReturn(PETSC_SUCCESS);
7394: }

7396: /*@
7397:   DMGetLabelOutput - Get the output flag for a given label

7399:   Not Collective

7401:   Input Parameters:
7402: + dm   - The `DM` object
7403: - name - The label name

7405:   Output Parameter:
7406: . output - The flag for output

7408:   Level: developer

7410: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMSetLabelOutput()`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7411: @*/
7412: PetscErrorCode DMGetLabelOutput(DM dm, const char name[], PetscBool *output)
7413: {
7414:   DMLabelLink next = dm->labels;
7415:   const char *lname;

7417:   PetscFunctionBegin;
7419:   PetscAssertPointer(name, 2);
7420:   PetscAssertPointer(output, 3);
7421:   while (next) {
7422:     PetscBool flg;

7424:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7425:     PetscCall(PetscStrcmp(name, lname, &flg));
7426:     if (flg) {
7427:       *output = next->output;
7428:       PetscFunctionReturn(PETSC_SUCCESS);
7429:     }
7430:     next = next->next;
7431:   }
7432:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No label named %s was present in this dm", name);
7433: }

7435: /*@
7436:   DMSetLabelOutput - Set if a given label should be saved to a `PetscViewer` in calls to `DMView()`

7438:   Not Collective

7440:   Input Parameters:
7441: + dm     - The `DM` object
7442: . name   - The label name
7443: - output - `PETSC_TRUE` to save the label to the viewer

7445:   Level: developer

7447: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetOutputFlag()`, `DMGetLabelOutput()`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7448: @*/
7449: PetscErrorCode DMSetLabelOutput(DM dm, const char name[], PetscBool output)
7450: {
7451:   DMLabelLink next = dm->labels;
7452:   const char *lname;

7454:   PetscFunctionBegin;
7456:   PetscAssertPointer(name, 2);
7457:   while (next) {
7458:     PetscBool flg;

7460:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7461:     PetscCall(PetscStrcmp(name, lname, &flg));
7462:     if (flg) {
7463:       next->output = output;
7464:       PetscFunctionReturn(PETSC_SUCCESS);
7465:     }
7466:     next = next->next;
7467:   }
7468:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No label named %s was present in this dm", name);
7469: }

7471: /*@
7472:   DMCopyLabels - Copy labels from one `DM` mesh to another `DM` with a superset of the points

7474:   Collective

7476:   Input Parameters:
7477: + dmA   - The `DM` object with initial labels
7478: . dmB   - The `DM` object to which labels are copied
7479: . mode  - Copy labels by pointers (`PETSC_OWN_POINTER`) or duplicate them (`PETSC_COPY_VALUES`)
7480: . all   - Copy all labels including "depth", "dim", and "celltype" (`PETSC_TRUE`) which are otherwise ignored (`PETSC_FALSE`)
7481: - emode - How to behave when a `DMLabel` in the source and destination `DM`s with the same name is encountered (see `DMCopyLabelsMode`)

7483:   Level: intermediate

7485:   Note:
7486:   This is typically used when interpolating or otherwise adding to a mesh, or testing.

7488: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMAddLabel()`, `DMCopyLabelsMode`
7489: @*/
7490: PetscErrorCode DMCopyLabels(DM dmA, DM dmB, PetscCopyMode mode, PetscBool all, DMCopyLabelsMode emode)
7491: {
7492:   DMLabel     label, labelNew, labelOld;
7493:   const char *name;
7494:   PetscBool   flg;
7495:   DMLabelLink link;

7497:   PetscFunctionBegin;
7502:   PetscCheck(mode != PETSC_USE_POINTER, PetscObjectComm((PetscObject)dmA), PETSC_ERR_SUP, "PETSC_USE_POINTER not supported for objects");
7503:   if (dmA == dmB) PetscFunctionReturn(PETSC_SUCCESS);
7504:   for (link = dmA->labels; link; link = link->next) {
7505:     label = link->label;
7506:     PetscCall(PetscObjectGetName((PetscObject)label, &name));
7507:     if (!all) {
7508:       PetscCall(PetscStrcmp(name, "depth", &flg));
7509:       if (flg) continue;
7510:       PetscCall(PetscStrcmp(name, "dim", &flg));
7511:       if (flg) continue;
7512:       PetscCall(PetscStrcmp(name, "celltype", &flg));
7513:       if (flg) continue;
7514:     }
7515:     PetscCall(DMGetLabel(dmB, name, &labelOld));
7516:     if (labelOld) {
7517:       switch (emode) {
7518:       case DM_COPY_LABELS_KEEP:
7519:         continue;
7520:       case DM_COPY_LABELS_REPLACE:
7521:         PetscCall(DMRemoveLabelBySelf(dmB, &labelOld, PETSC_TRUE));
7522:         break;
7523:       case DM_COPY_LABELS_FAIL:
7524:         SETERRQ(PetscObjectComm((PetscObject)dmA), PETSC_ERR_ARG_OUTOFRANGE, "Label %s already exists in destination DM", name);
7525:       default:
7526:         SETERRQ(PetscObjectComm((PetscObject)dmA), PETSC_ERR_ARG_OUTOFRANGE, "Unhandled DMCopyLabelsMode %d", (int)emode);
7527:       }
7528:     }
7529:     if (mode == PETSC_COPY_VALUES) {
7530:       PetscCall(DMLabelDuplicate(label, &labelNew));
7531:     } else {
7532:       labelNew = label;
7533:     }
7534:     PetscCall(DMAddLabel(dmB, labelNew));
7535:     if (mode == PETSC_COPY_VALUES) PetscCall(DMLabelDestroy(&labelNew));
7536:   }
7537:   PetscFunctionReturn(PETSC_SUCCESS);
7538: }

7540: /*@C
7541:   DMCompareLabels - Compare labels between two `DM` objects

7543:   Collective; No Fortran Support

7545:   Input Parameters:
7546: + dm0 - First `DM` object
7547: - dm1 - Second `DM` object

7549:   Output Parameters:
7550: + equal   - (Optional) Flag whether labels of dm0 and dm1 are the same
7551: - message - (Optional) Message describing the difference, or `NULL` if there is no difference

7553:   Level: intermediate

7555:   Notes:
7556:   The output flag equal will be the same on all processes.

7558:   If equal is passed as `NULL` and difference is found, an error is thrown on all processes.

7560:   Make sure to pass equal is `NULL` on all processes or none of them.

7562:   The output message is set independently on each rank.

7564:   message must be freed with `PetscFree()`

7566:   If message is passed as `NULL` and a difference is found, the difference description is printed to stderr in synchronized manner.

7568:   Make sure to pass message as `NULL` on all processes or no processes.

7570:   Labels are matched by name. If the number of labels and their names are equal,
7571:   `DMLabelCompare()` is used to compare each pair of labels with the same name.

7573:   Developer Note:
7574:   Can automatically generate the Fortran stub because `message` must be freed with `PetscFree()`

7576: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMAddLabel()`, `DMCopyLabelsMode`, `DMLabelCompare()`
7577: @*/
7578: PetscErrorCode DMCompareLabels(DM dm0, DM dm1, PetscBool *equal, char **message)
7579: {
7580:   PetscInt    n, i;
7581:   char        msg[PETSC_MAX_PATH_LEN] = "";
7582:   PetscBool   eq;
7583:   MPI_Comm    comm;
7584:   PetscMPIInt rank;

7586:   PetscFunctionBegin;
7589:   PetscCheckSameComm(dm0, 1, dm1, 2);
7590:   if (equal) PetscAssertPointer(equal, 3);
7591:   if (message) PetscAssertPointer(message, 4);
7592:   PetscCall(PetscObjectGetComm((PetscObject)dm0, &comm));
7593:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7594:   {
7595:     PetscInt n1;

7597:     PetscCall(DMGetNumLabels(dm0, &n));
7598:     PetscCall(DMGetNumLabels(dm1, &n1));
7599:     eq = (PetscBool)(n == n1);
7600:     if (!eq) PetscCall(PetscSNPrintf(msg, sizeof(msg), "Number of labels in dm0 = %" PetscInt_FMT " != %" PetscInt_FMT " = Number of labels in dm1", n, n1));
7601:     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPIU_BOOL, MPI_LAND, comm));
7602:     if (!eq) goto finish;
7603:   }
7604:   for (i = 0; i < n; i++) {
7605:     DMLabel     l0, l1;
7606:     const char *name;
7607:     char       *msgInner;

7609:     /* Ignore label order */
7610:     PetscCall(DMGetLabelByNum(dm0, i, &l0));
7611:     PetscCall(PetscObjectGetName((PetscObject)l0, &name));
7612:     PetscCall(DMGetLabel(dm1, name, &l1));
7613:     if (!l1) {
7614:       PetscCall(PetscSNPrintf(msg, sizeof(msg), "Label \"%s\" (#%" PetscInt_FMT " in dm0) not found in dm1", name, i));
7615:       eq = PETSC_FALSE;
7616:       break;
7617:     }
7618:     PetscCall(DMLabelCompare(comm, l0, l1, &eq, &msgInner));
7619:     PetscCall(PetscStrncpy(msg, msgInner, sizeof(msg)));
7620:     PetscCall(PetscFree(msgInner));
7621:     if (!eq) break;
7622:   }
7623:   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPIU_BOOL, MPI_LAND, comm));
7624: finish:
7625:   /* If message output arg not set, print to stderr */
7626:   if (message) {
7627:     *message = NULL;
7628:     if (msg[0]) PetscCall(PetscStrallocpy(msg, message));
7629:   } else {
7630:     if (msg[0]) PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] %s\n", rank, msg));
7631:     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR));
7632:   }
7633:   /* If same output arg not ser and labels are not equal, throw error */
7634:   if (equal) *equal = eq;
7635:   else PetscCheck(eq, comm, PETSC_ERR_ARG_INCOMP, "DMLabels are not the same in dm0 and dm1");
7636:   PetscFunctionReturn(PETSC_SUCCESS);
7637: }

7639: PetscErrorCode DMSetLabelValue_Fast(DM dm, DMLabel *label, const char name[], PetscInt point, PetscInt value)
7640: {
7641:   PetscFunctionBegin;
7642:   PetscAssertPointer(label, 2);
7643:   if (!*label) {
7644:     PetscCall(DMCreateLabel(dm, name));
7645:     PetscCall(DMGetLabel(dm, name, label));
7646:   }
7647:   PetscCall(DMLabelSetValue(*label, point, value));
7648:   PetscFunctionReturn(PETSC_SUCCESS);
7649: }

7651: /*
7652:   Many mesh programs, such as Triangle and TetGen, allow only a single label for each mesh point. Therefore, we would
7653:   like to encode all label IDs using a single, universal label. We can do this by assigning an integer to every
7654:   (label, id) pair in the DM.

7656:   However, a mesh point can have multiple labels, so we must separate all these values. We will assign a bit range to
7657:   each label.
7658: */
7659: PetscErrorCode DMUniversalLabelCreate(DM dm, DMUniversalLabel *universal)
7660: {
7661:   DMUniversalLabel ul;
7662:   PetscBool       *active;
7663:   PetscInt         pStart, pEnd, p, Nl, l, m;

7665:   PetscFunctionBegin;
7666:   PetscCall(PetscMalloc1(1, &ul));
7667:   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "universal", &ul->label));
7668:   PetscCall(DMGetNumLabels(dm, &Nl));
7669:   PetscCall(PetscCalloc1(Nl, &active));
7670:   ul->Nl = 0;
7671:   for (l = 0; l < Nl; ++l) {
7672:     PetscBool   isdepth, iscelltype;
7673:     const char *name;

7675:     PetscCall(DMGetLabelName(dm, l, &name));
7676:     PetscCall(PetscStrncmp(name, "depth", 6, &isdepth));
7677:     PetscCall(PetscStrncmp(name, "celltype", 9, &iscelltype));
7678:     active[l] = !(isdepth || iscelltype) ? PETSC_TRUE : PETSC_FALSE;
7679:     if (active[l]) ++ul->Nl;
7680:   }
7681:   PetscCall(PetscCalloc5(ul->Nl, &ul->names, ul->Nl, &ul->indices, ul->Nl + 1, &ul->offsets, ul->Nl + 1, &ul->bits, ul->Nl, &ul->masks));
7682:   ul->Nv = 0;
7683:   for (l = 0, m = 0; l < Nl; ++l) {
7684:     DMLabel     label;
7685:     PetscInt    nv;
7686:     const char *name;

7688:     if (!active[l]) continue;
7689:     PetscCall(DMGetLabelName(dm, l, &name));
7690:     PetscCall(DMGetLabelByNum(dm, l, &label));
7691:     PetscCall(DMLabelGetNumValues(label, &nv));
7692:     PetscCall(PetscStrallocpy(name, &ul->names[m]));
7693:     ul->indices[m] = l;
7694:     ul->Nv += nv;
7695:     ul->offsets[m + 1] = nv;
7696:     ul->bits[m + 1]    = PetscCeilReal(PetscLog2Real(nv + 1));
7697:     ++m;
7698:   }
7699:   for (l = 1; l <= ul->Nl; ++l) {
7700:     ul->offsets[l] = ul->offsets[l - 1] + ul->offsets[l];
7701:     ul->bits[l]    = ul->bits[l - 1] + ul->bits[l];
7702:   }
7703:   for (l = 0; l < ul->Nl; ++l) {
7704:     PetscInt b;

7706:     ul->masks[l] = 0;
7707:     for (b = ul->bits[l]; b < ul->bits[l + 1]; ++b) ul->masks[l] |= 1 << b;
7708:   }
7709:   PetscCall(PetscMalloc1(ul->Nv, &ul->values));
7710:   for (l = 0, m = 0; l < Nl; ++l) {
7711:     DMLabel         label;
7712:     IS              valueIS;
7713:     const PetscInt *varr;
7714:     PetscInt        nv, v;

7716:     if (!active[l]) continue;
7717:     PetscCall(DMGetLabelByNum(dm, l, &label));
7718:     PetscCall(DMLabelGetNumValues(label, &nv));
7719:     PetscCall(DMLabelGetValueIS(label, &valueIS));
7720:     PetscCall(ISGetIndices(valueIS, &varr));
7721:     for (v = 0; v < nv; ++v) ul->values[ul->offsets[m] + v] = varr[v];
7722:     PetscCall(ISRestoreIndices(valueIS, &varr));
7723:     PetscCall(ISDestroy(&valueIS));
7724:     PetscCall(PetscSortInt(nv, &ul->values[ul->offsets[m]]));
7725:     ++m;
7726:   }
7727:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
7728:   for (p = pStart; p < pEnd; ++p) {
7729:     PetscInt  uval   = 0;
7730:     PetscBool marked = PETSC_FALSE;

7732:     for (l = 0, m = 0; l < Nl; ++l) {
7733:       DMLabel  label;
7734:       PetscInt val, defval, loc, nv;

7736:       if (!active[l]) continue;
7737:       PetscCall(DMGetLabelByNum(dm, l, &label));
7738:       PetscCall(DMLabelGetValue(label, p, &val));
7739:       PetscCall(DMLabelGetDefaultValue(label, &defval));
7740:       if (val == defval) {
7741:         ++m;
7742:         continue;
7743:       }
7744:       nv     = ul->offsets[m + 1] - ul->offsets[m];
7745:       marked = PETSC_TRUE;
7746:       PetscCall(PetscFindInt(val, nv, &ul->values[ul->offsets[m]], &loc));
7747:       PetscCheck(loc >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Label value %" PetscInt_FMT " not found in compression array", val);
7748:       uval += (loc + 1) << ul->bits[m];
7749:       ++m;
7750:     }
7751:     if (marked) PetscCall(DMLabelSetValue(ul->label, p, uval));
7752:   }
7753:   PetscCall(PetscFree(active));
7754:   *universal = ul;
7755:   PetscFunctionReturn(PETSC_SUCCESS);
7756: }

7758: PetscErrorCode DMUniversalLabelDestroy(DMUniversalLabel *universal)
7759: {
7760:   PetscInt l;

7762:   PetscFunctionBegin;
7763:   for (l = 0; l < (*universal)->Nl; ++l) PetscCall(PetscFree((*universal)->names[l]));
7764:   PetscCall(DMLabelDestroy(&(*universal)->label));
7765:   PetscCall(PetscFree5((*universal)->names, (*universal)->indices, (*universal)->offsets, (*universal)->bits, (*universal)->masks));
7766:   PetscCall(PetscFree((*universal)->values));
7767:   PetscCall(PetscFree(*universal));
7768:   *universal = NULL;
7769:   PetscFunctionReturn(PETSC_SUCCESS);
7770: }

7772: PetscErrorCode DMUniversalLabelGetLabel(DMUniversalLabel ul, DMLabel *ulabel)
7773: {
7774:   PetscFunctionBegin;
7775:   PetscAssertPointer(ulabel, 2);
7776:   *ulabel = ul->label;
7777:   PetscFunctionReturn(PETSC_SUCCESS);
7778: }

7780: PetscErrorCode DMUniversalLabelCreateLabels(DMUniversalLabel ul, PetscBool preserveOrder, DM dm)
7781: {
7782:   PetscInt Nl = ul->Nl, l;

7784:   PetscFunctionBegin;
7786:   for (l = 0; l < Nl; ++l) {
7787:     if (preserveOrder) PetscCall(DMCreateLabelAtIndex(dm, ul->indices[l], ul->names[l]));
7788:     else PetscCall(DMCreateLabel(dm, ul->names[l]));
7789:   }
7790:   if (preserveOrder) {
7791:     for (l = 0; l < ul->Nl; ++l) {
7792:       const char *name;
7793:       PetscBool   match;

7795:       PetscCall(DMGetLabelName(dm, ul->indices[l], &name));
7796:       PetscCall(PetscStrcmp(name, ul->names[l], &match));
7797:       PetscCheck(match, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Label %" PetscInt_FMT " name %s does not match new name %s", l, name, ul->names[l]);
7798:     }
7799:   }
7800:   PetscFunctionReturn(PETSC_SUCCESS);
7801: }

7803: PetscErrorCode DMUniversalLabelSetLabelValue(DMUniversalLabel ul, DM dm, PetscBool useIndex, PetscInt p, PetscInt value)
7804: {
7805:   PetscInt l;

7807:   PetscFunctionBegin;
7808:   for (l = 0; l < ul->Nl; ++l) {
7809:     DMLabel  label;
7810:     PetscInt lval = (value & ul->masks[l]) >> ul->bits[l];

7812:     if (lval) {
7813:       if (useIndex) PetscCall(DMGetLabelByNum(dm, ul->indices[l], &label));
7814:       else PetscCall(DMGetLabel(dm, ul->names[l], &label));
7815:       PetscCall(DMLabelSetValue(label, p, ul->values[ul->offsets[l] + lval - 1]));
7816:     }
7817:   }
7818:   PetscFunctionReturn(PETSC_SUCCESS);
7819: }

7821: /*@
7822:   DMGetCoarseDM - Get the coarse `DM`from which this `DM` was obtained by refinement

7824:   Not Collective

7826:   Input Parameter:
7827: . dm - The `DM` object

7829:   Output Parameter:
7830: . cdm - The coarse `DM`

7832:   Level: intermediate

7834: .seealso: [](ch_dmbase), `DM`, `DMSetCoarseDM()`, `DMCoarsen()`
7835: @*/
7836: PetscErrorCode DMGetCoarseDM(DM dm, DM *cdm)
7837: {
7838:   PetscFunctionBegin;
7840:   PetscAssertPointer(cdm, 2);
7841:   *cdm = dm->coarseMesh;
7842:   PetscFunctionReturn(PETSC_SUCCESS);
7843: }

7845: /*@
7846:   DMSetCoarseDM - Set the coarse `DM` from which this `DM` was obtained by refinement

7848:   Input Parameters:
7849: + dm  - The `DM` object
7850: - cdm - The coarse `DM`

7852:   Level: intermediate

7854:   Note:
7855:   Normally this is set automatically by `DMRefine()`

7857: .seealso: [](ch_dmbase), `DM`, `DMGetCoarseDM()`, `DMCoarsen()`, `DMSetRefine()`, `DMSetFineDM()`
7858: @*/
7859: PetscErrorCode DMSetCoarseDM(DM dm, DM cdm)
7860: {
7861:   PetscFunctionBegin;
7864:   if (dm == cdm) cdm = NULL;
7865:   PetscCall(PetscObjectReference((PetscObject)cdm));
7866:   PetscCall(DMDestroy(&dm->coarseMesh));
7867:   dm->coarseMesh = cdm;
7868:   PetscFunctionReturn(PETSC_SUCCESS);
7869: }

7871: /*@
7872:   DMGetFineDM - Get the fine mesh from which this `DM` was obtained by coarsening

7874:   Input Parameter:
7875: . dm - The `DM` object

7877:   Output Parameter:
7878: . fdm - The fine `DM`

7880:   Level: intermediate

7882: .seealso: [](ch_dmbase), `DM`, `DMSetFineDM()`, `DMCoarsen()`, `DMRefine()`
7883: @*/
7884: PetscErrorCode DMGetFineDM(DM dm, DM *fdm)
7885: {
7886:   PetscFunctionBegin;
7888:   PetscAssertPointer(fdm, 2);
7889:   *fdm = dm->fineMesh;
7890:   PetscFunctionReturn(PETSC_SUCCESS);
7891: }

7893: /*@
7894:   DMSetFineDM - Set the fine mesh from which this was obtained by coarsening

7896:   Input Parameters:
7897: + dm  - The `DM` object
7898: - fdm - The fine `DM`

7900:   Level: developer

7902:   Note:
7903:   Normally this is set automatically by `DMCoarsen()`

7905: .seealso: [](ch_dmbase), `DM`, `DMGetFineDM()`, `DMCoarsen()`, `DMRefine()`
7906: @*/
7907: PetscErrorCode DMSetFineDM(DM dm, DM fdm)
7908: {
7909:   PetscFunctionBegin;
7912:   if (dm == fdm) fdm = NULL;
7913:   PetscCall(PetscObjectReference((PetscObject)fdm));
7914:   PetscCall(DMDestroy(&dm->fineMesh));
7915:   dm->fineMesh = fdm;
7916:   PetscFunctionReturn(PETSC_SUCCESS);
7917: }

7919: /*@C
7920:   DMAddBoundary - Add a boundary condition to a model represented by a `DM`

7922:   Collective

7924:   Input Parameters:
7925: + dm       - The `DM`, with a `PetscDS` that matches the problem being constrained
7926: . type     - The type of condition, e.g. `DM_BC_ESSENTIAL_ANALYTIC`, `DM_BC_ESSENTIAL_FIELD` (Dirichlet), or `DM_BC_NATURAL` (Neumann)
7927: . name     - The BC name
7928: . label    - The label defining constrained points
7929: . Nv       - The number of `DMLabel` values for constrained points
7930: . values   - An array of values for constrained points
7931: . field    - The field to constrain
7932: . Nc       - The number of constrained field components (0 will constrain all fields)
7933: . comps    - An array of constrained component numbers
7934: . bcFunc   - A pointwise function giving boundary values
7935: . bcFunc_t - A pointwise function giving the time deriative of the boundary values, or NULL
7936: - ctx      - An optional user context for bcFunc

7938:   Output Parameter:
7939: . bd - (Optional) Boundary number

7941:   Options Database Keys:
7942: + -bc_<boundary name> <num>      - Overrides the boundary ids
7943: - -bc_<boundary name>_comp <num> - Overrides the boundary components

7945:   Level: intermediate

7947:   Notes:
7948:   Both bcFunc and bcFunc_t will depend on the boundary condition type. If the type if `DM_BC_ESSENTIAL`, then the calling sequence is\:
7949: .vb
7950:  void bcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[])
7951: .ve

7953:   If the type is `DM_BC_ESSENTIAL_FIELD` or other _FIELD value, then the calling sequence is\:

7955: .vb
7956:   void bcFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
7957:               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
7958:               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
7959:               PetscReal time, const PetscReal x[], PetscScalar bcval[])
7960: .ve
7961: + dim - the spatial dimension
7962: . Nf - the number of fields
7963: . uOff - the offset into u[] and u_t[] for each field
7964: . uOff_x - the offset into u_x[] for each field
7965: . u - each field evaluated at the current point
7966: . u_t - the time derivative of each field evaluated at the current point
7967: . u_x - the gradient of each field evaluated at the current point
7968: . aOff - the offset into a[] and a_t[] for each auxiliary field
7969: . aOff_x - the offset into a_x[] for each auxiliary field
7970: . a - each auxiliary field evaluated at the current point
7971: . a_t - the time derivative of each auxiliary field evaluated at the current point
7972: . a_x - the gradient of auxiliary each field evaluated at the current point
7973: . t - current time
7974: . x - coordinates of the current point
7975: . numConstants - number of constant parameters
7976: . constants - constant parameters
7977: - bcval - output values at the current point

7979: .seealso: [](ch_dmbase), `DM`, `DSGetBoundary()`, `PetscDSAddBoundary()`
7980: @*/
7981: PetscErrorCode DMAddBoundary(DM dm, DMBoundaryConditionType type, const char name[], DMLabel label, PetscInt Nv, const PetscInt values[], PetscInt field, PetscInt Nc, const PetscInt comps[], void (*bcFunc)(void), void (*bcFunc_t)(void), void *ctx, PetscInt *bd)
7982: {
7983:   PetscDS ds;

7985:   PetscFunctionBegin;
7992:   PetscCheck(!dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot add boundary to DM after creating local section");
7993:   PetscCall(DMGetDS(dm, &ds));
7994:   /* Complete label */
7995:   if (label) {
7996:     PetscObject  obj;
7997:     PetscClassId id;

7999:     PetscCall(DMGetField(dm, field, NULL, &obj));
8000:     PetscCall(PetscObjectGetClassId(obj, &id));
8001:     if (id == PETSCFE_CLASSID) {
8002:       DM plex;

8004:       PetscCall(DMConvert(dm, DMPLEX, &plex));
8005:       if (plex) PetscCall(DMPlexLabelComplete(plex, label));
8006:       PetscCall(DMDestroy(&plex));
8007:     }
8008:   }
8009:   PetscCall(PetscDSAddBoundary(ds, type, name, label, Nv, values, field, Nc, comps, bcFunc, bcFunc_t, ctx, bd));
8010:   PetscFunctionReturn(PETSC_SUCCESS);
8011: }

8013: /* TODO Remove this since now the structures are the same */
8014: static PetscErrorCode DMPopulateBoundary(DM dm)
8015: {
8016:   PetscDS     ds;
8017:   DMBoundary *lastnext;
8018:   DSBoundary  dsbound;

8020:   PetscFunctionBegin;
8021:   PetscCall(DMGetDS(dm, &ds));
8022:   dsbound = ds->boundary;
8023:   if (dm->boundary) {
8024:     DMBoundary next = dm->boundary;

8026:     /* quick check to see if the PetscDS has changed */
8027:     if (next->dsboundary == dsbound) PetscFunctionReturn(PETSC_SUCCESS);
8028:     /* the PetscDS has changed: tear down and rebuild */
8029:     while (next) {
8030:       DMBoundary b = next;

8032:       next = b->next;
8033:       PetscCall(PetscFree(b));
8034:     }
8035:     dm->boundary = NULL;
8036:   }

8038:   lastnext = &dm->boundary;
8039:   while (dsbound) {
8040:     DMBoundary dmbound;

8042:     PetscCall(PetscNew(&dmbound));
8043:     dmbound->dsboundary = dsbound;
8044:     dmbound->label      = dsbound->label;
8045:     /* push on the back instead of the front so that it is in the same order as in the PetscDS */
8046:     *lastnext = dmbound;
8047:     lastnext  = &dmbound->next;
8048:     dsbound   = dsbound->next;
8049:   }
8050:   PetscFunctionReturn(PETSC_SUCCESS);
8051: }

8053: /* TODO: missing manual page */
8054: PetscErrorCode DMIsBoundaryPoint(DM dm, PetscInt point, PetscBool *isBd)
8055: {
8056:   DMBoundary b;

8058:   PetscFunctionBegin;
8060:   PetscAssertPointer(isBd, 3);
8061:   *isBd = PETSC_FALSE;
8062:   PetscCall(DMPopulateBoundary(dm));
8063:   b = dm->boundary;
8064:   while (b && !*isBd) {
8065:     DMLabel    label = b->label;
8066:     DSBoundary dsb   = b->dsboundary;
8067:     PetscInt   i;

8069:     if (label) {
8070:       for (i = 0; i < dsb->Nv && !*isBd; ++i) PetscCall(DMLabelStratumHasPoint(label, dsb->values[i], point, isBd));
8071:     }
8072:     b = b->next;
8073:   }
8074:   PetscFunctionReturn(PETSC_SUCCESS);
8075: }

8077: /*@C
8078:   DMProjectFunction - This projects the given function into the function space provided by a `DM`, putting the coefficients in a global vector.

8080:   Collective

8082:   Input Parameters:
8083: + dm    - The `DM`
8084: . time  - The time
8085: . funcs - The coordinate functions to evaluate, one per field
8086: . ctxs  - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
8087: - mode  - The insertion mode for values

8089:   Output Parameter:
8090: . X - vector

8092:   Calling sequence of `funcs`:
8093: + dim  - The spatial dimension
8094: . time - The time at which to sample
8095: . x    - The coordinates
8096: . Nc   - The number of components
8097: . u    - The output field values
8098: - ctx  - optional user-defined function context

8100:   Level: developer

8102:   Developer Notes:
8103:   This API is specific to only particular usage of `DM`

8105:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8107: .seealso: [](ch_dmbase), `DM`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
8108: @*/
8109: PetscErrorCode DMProjectFunction(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx), void **ctxs, InsertMode mode, Vec X)
8110: {
8111:   Vec localX;

8113:   PetscFunctionBegin;
8115:   PetscCall(PetscLogEventBegin(DM_ProjectFunction, dm, X, 0, 0));
8116:   PetscCall(DMGetLocalVector(dm, &localX));
8117:   PetscCall(VecSet(localX, 0.));
8118:   PetscCall(DMProjectFunctionLocal(dm, time, funcs, ctxs, mode, localX));
8119:   PetscCall(DMLocalToGlobalBegin(dm, localX, mode, X));
8120:   PetscCall(DMLocalToGlobalEnd(dm, localX, mode, X));
8121:   PetscCall(DMRestoreLocalVector(dm, &localX));
8122:   PetscCall(PetscLogEventEnd(DM_ProjectFunction, dm, X, 0, 0));
8123:   PetscFunctionReturn(PETSC_SUCCESS);
8124: }

8126: /*@C
8127:   DMProjectFunctionLocal - This projects the given function into the function space provided by a `DM`, putting the coefficients in a local vector.

8129:   Not Collective

8131:   Input Parameters:
8132: + dm    - The `DM`
8133: . time  - The time
8134: . funcs - The coordinate functions to evaluate, one per field
8135: . ctxs  - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
8136: - mode  - The insertion mode for values

8138:   Output Parameter:
8139: . localX - vector

8141:   Calling sequence of `funcs`:
8142: + dim  - The spatial dimension
8143: . time - The current timestep
8144: . x    - The coordinates
8145: . Nc   - The number of components
8146: . u    - The output field values
8147: - ctx  - optional user-defined function context

8149:   Level: developer

8151:   Developer Notes:
8152:   This API is specific to only particular usage of `DM`

8154:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8156: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
8157: @*/
8158: PetscErrorCode DMProjectFunctionLocal(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx), void **ctxs, InsertMode mode, Vec localX)
8159: {
8160:   PetscFunctionBegin;
8163:   PetscUseTypeMethod(dm, projectfunctionlocal, time, funcs, ctxs, mode, localX);
8164:   PetscFunctionReturn(PETSC_SUCCESS);
8165: }

8167: /*@C
8168:   DMProjectFunctionLabel - This projects the given function into the function space provided by the `DM`, putting the coefficients in a global vector, setting values only for points in the given label.

8170:   Collective

8172:   Input Parameters:
8173: + dm     - The `DM`
8174: . time   - The time
8175: . numIds - The number of ids
8176: . ids    - The ids
8177: . Nc     - The number of components
8178: . comps  - The components
8179: . label  - The `DMLabel` selecting the portion of the mesh for projection
8180: . funcs  - The coordinate functions to evaluate, one per field
8181: . ctxs   - Optional array of contexts to pass to each coordinate function.  ctxs may be null.
8182: - mode   - The insertion mode for values

8184:   Output Parameter:
8185: . X - vector

8187:   Calling sequence of `funcs`:
8188: + dim  - The spatial dimension
8189: . time - The current timestep
8190: . x    - The coordinates
8191: . Nc   - The number of components
8192: . u    - The output field values
8193: - ctx  - optional user-defined function context

8195:   Level: developer

8197:   Developer Notes:
8198:   This API is specific to only particular usage of `DM`

8200:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8202: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabelLocal()`, `DMComputeL2Diff()`
8203: @*/
8204: PetscErrorCode DMProjectFunctionLabel(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Nc, const PetscInt comps[], PetscErrorCode (**funcs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx), void **ctxs, InsertMode mode, Vec X)
8205: {
8206:   Vec localX;

8208:   PetscFunctionBegin;
8210:   PetscCall(DMGetLocalVector(dm, &localX));
8211:   PetscCall(VecSet(localX, 0.));
8212:   PetscCall(DMProjectFunctionLabelLocal(dm, time, label, numIds, ids, Nc, comps, funcs, ctxs, mode, localX));
8213:   PetscCall(DMLocalToGlobalBegin(dm, localX, mode, X));
8214:   PetscCall(DMLocalToGlobalEnd(dm, localX, mode, X));
8215:   PetscCall(DMRestoreLocalVector(dm, &localX));
8216:   PetscFunctionReturn(PETSC_SUCCESS);
8217: }

8219: /*@C
8220:   DMProjectFunctionLabelLocal - This projects the given function into the function space provided by the `DM`, putting the coefficients in a local vector, setting values only for points in the given label.

8222:   Not Collective

8224:   Input Parameters:
8225: + dm     - The `DM`
8226: . time   - The time
8227: . label  - The `DMLabel` selecting the portion of the mesh for projection
8228: . numIds - The number of ids
8229: . ids    - The ids
8230: . Nc     - The number of components
8231: . comps  - The components
8232: . funcs  - The coordinate functions to evaluate, one per field
8233: . ctxs   - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
8234: - mode   - The insertion mode for values

8236:   Output Parameter:
8237: . localX - vector

8239:   Calling sequence of `funcs`:
8240: + dim  - The spatial dimension
8241: . time - The current time
8242: . x    - The coordinates
8243: . Nc   - The number of components
8244: . u    - The output field values
8245: - ctx  - optional user-defined function context

8247:   Level: developer

8249:   Developer Notes:
8250:   This API is specific to only particular usage of `DM`

8252:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8254: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
8255: @*/
8256: PetscErrorCode DMProjectFunctionLabelLocal(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Nc, const PetscInt comps[], PetscErrorCode (**funcs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx), void **ctxs, InsertMode mode, Vec localX)
8257: {
8258:   PetscFunctionBegin;
8261:   PetscUseTypeMethod(dm, projectfunctionlabellocal, time, label, numIds, ids, Nc, comps, funcs, ctxs, mode, localX);
8262:   PetscFunctionReturn(PETSC_SUCCESS);
8263: }

8265: /*@C
8266:   DMProjectFieldLocal - This projects the given function of the input fields into the function space provided by the `DM`, putting the coefficients in a local vector.

8268:   Not Collective

8270:   Input Parameters:
8271: + dm     - The `DM`
8272: . time   - The time
8273: . localU - The input field vector; may be `NULL` if projection is defined purely by coordinates
8274: . funcs  - The functions to evaluate, one per field
8275: - mode   - The insertion mode for values

8277:   Output Parameter:
8278: . localX - The output vector

8280:   Calling sequence of `funcs`:
8281: + dim          - The spatial dimension
8282: . Nf           - The number of input fields
8283: . NfAux        - The number of input auxiliary fields
8284: . uOff         - The offset of each field in u[]
8285: . uOff_x       - The offset of each field in u_x[]
8286: . u            - The field values at this point in space
8287: . u_t          - The field time derivative at this point in space (or NULL)
8288: . u_x          - The field derivatives at this point in space
8289: . aOff         - The offset of each auxiliary field in u[]
8290: . aOff_x       - The offset of each auxiliary field in u_x[]
8291: . a            - The auxiliary field values at this point in space
8292: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8293: . a_x          - The auxiliary field derivatives at this point in space
8294: . t            - The current time
8295: . x            - The coordinates of this point
8296: . numConstants - The number of constants
8297: . constants    - The value of each constant
8298: - f            - The value of the function at this point in space

8300:   Level: intermediate

8302:   Note:
8303:   There are three different `DM`s that potentially interact in this function. The output `DM`, dm, specifies the layout of the values calculates by funcs.
8304:   The input `DM`, attached to U, may be different. For example, you can input the solution over the full domain, but output over a piece of the boundary, or
8305:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8306:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8308:   Developer Notes:
8309:   This API is specific to only particular usage of `DM`

8311:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8313: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabelLocal()`,
8314: `DMProjectFunction()`, `DMComputeL2Diff()`
8315: @*/
8316: PetscErrorCode DMProjectFieldLocal(DM dm, PetscReal time, Vec localU, void (**funcs)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]), InsertMode mode, Vec localX)
8317: {
8318:   PetscFunctionBegin;
8322:   PetscUseTypeMethod(dm, projectfieldlocal, time, localU, funcs, mode, localX);
8323:   PetscFunctionReturn(PETSC_SUCCESS);
8324: }

8326: /*@C
8327:   DMProjectFieldLabelLocal - This projects the given function of the input fields into the function space provided, putting the coefficients in a local vector, calculating only over the portion of the domain specified by the label.

8329:   Not Collective

8331:   Input Parameters:
8332: + dm     - The `DM`
8333: . time   - The time
8334: . label  - The `DMLabel` marking the portion of the domain to output
8335: . numIds - The number of label ids to use
8336: . ids    - The label ids to use for marking
8337: . Nc     - The number of components to set in the output, or `PETSC_DETERMINE` for all components
8338: . comps  - The components to set in the output, or `NULL` for all components
8339: . localU - The input field vector
8340: . funcs  - The functions to evaluate, one per field
8341: - mode   - The insertion mode for values

8343:   Output Parameter:
8344: . localX - The output vector

8346:   Calling sequence of `funcs`:
8347: + dim          - The spatial dimension
8348: . Nf           - The number of input fields
8349: . NfAux        - The number of input auxiliary fields
8350: . uOff         - The offset of each field in u[]
8351: . uOff_x       - The offset of each field in u_x[]
8352: . u            - The field values at this point in space
8353: . u_t          - The field time derivative at this point in space (or NULL)
8354: . u_x          - The field derivatives at this point in space
8355: . aOff         - The offset of each auxiliary field in u[]
8356: . aOff_x       - The offset of each auxiliary field in u_x[]
8357: . a            - The auxiliary field values at this point in space
8358: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8359: . a_x          - The auxiliary field derivatives at this point in space
8360: . t            - The current time
8361: . x            - The coordinates of this point
8362: . numConstants - The number of constants
8363: . constants    - The value of each constant
8364: - f            - The value of the function at this point in space

8366:   Level: intermediate

8368:   Note:
8369:   There are three different `DM`s that potentially interact in this function. The output `DM`, dm, specifies the layout of the values calculates by funcs.
8370:   The input `DM`, attached to localU, may be different. For example, you can input the solution over the full domain, but output over a piece of the boundary, or
8371:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8372:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8374:   Developer Notes:
8375:   This API is specific to only particular usage of `DM`

8377:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8379: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabel()`, `DMProjectFunction()`, `DMComputeL2Diff()`
8380: @*/
8381: PetscErrorCode DMProjectFieldLabelLocal(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Nc, const PetscInt comps[], Vec localU, void (**funcs)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]), InsertMode mode, Vec localX)
8382: {
8383:   PetscFunctionBegin;
8387:   PetscUseTypeMethod(dm, projectfieldlabellocal, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX);
8388:   PetscFunctionReturn(PETSC_SUCCESS);
8389: }

8391: /*@C
8392:   DMProjectFieldLabel - This projects the given function of the input fields into the function space provided, putting the coefficients in a global vector, calculating only over the portion of the domain specified by the label.

8394:   Not Collective

8396:   Input Parameters:
8397: + dm     - The `DM`
8398: . time   - The time
8399: . label  - The `DMLabel` marking the portion of the domain to output
8400: . numIds - The number of label ids to use
8401: . ids    - The label ids to use for marking
8402: . Nc     - The number of components to set in the output, or `PETSC_DETERMINE` for all components
8403: . comps  - The components to set in the output, or `NULL` for all components
8404: . U      - The input field vector
8405: . funcs  - The functions to evaluate, one per field
8406: - mode   - The insertion mode for values

8408:   Output Parameter:
8409: . X - The output vector

8411:   Calling sequence of `funcs`:
8412: + dim          - The spatial dimension
8413: . Nf           - The number of input fields
8414: . NfAux        - The number of input auxiliary fields
8415: . uOff         - The offset of each field in u[]
8416: . uOff_x       - The offset of each field in u_x[]
8417: . u            - The field values at this point in space
8418: . u_t          - The field time derivative at this point in space (or NULL)
8419: . u_x          - The field derivatives at this point in space
8420: . aOff         - The offset of each auxiliary field in u[]
8421: . aOff_x       - The offset of each auxiliary field in u_x[]
8422: . a            - The auxiliary field values at this point in space
8423: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8424: . a_x          - The auxiliary field derivatives at this point in space
8425: . t            - The current time
8426: . x            - The coordinates of this point
8427: . numConstants - The number of constants
8428: . constants    - The value of each constant
8429: - f            - The value of the function at this point in space

8431:   Level: intermediate

8433:   Note:
8434:   There are three different `DM`s that potentially interact in this function. The output `DM`, dm, specifies the layout of the values calculates by funcs.
8435:   The input `DM`, attached to U, may be different. For example, you can input the solution over the full domain, but output over a piece of the boundary, or
8436:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8437:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8439:   Developer Notes:
8440:   This API is specific to only particular usage of `DM`

8442:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8444: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabelLocal()`, `DMProjectFunction()`, `DMComputeL2Diff()`
8445: @*/
8446: PetscErrorCode DMProjectFieldLabel(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Nc, const PetscInt comps[], Vec U, void (**funcs)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]), InsertMode mode, Vec X)
8447: {
8448:   DM  dmIn;
8449:   Vec localU, localX;

8451:   PetscFunctionBegin;
8453:   PetscCall(VecGetDM(U, &dmIn));
8454:   PetscCall(DMGetLocalVector(dmIn, &localU));
8455:   PetscCall(DMGetLocalVector(dm, &localX));
8456:   PetscCall(VecSet(localX, 0.));
8457:   PetscCall(DMGlobalToLocalBegin(dmIn, U, mode, localU));
8458:   PetscCall(DMGlobalToLocalEnd(dmIn, U, mode, localU));
8459:   PetscCall(DMProjectFieldLabelLocal(dm, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX));
8460:   PetscCall(DMLocalToGlobalBegin(dm, localX, mode, X));
8461:   PetscCall(DMLocalToGlobalEnd(dm, localX, mode, X));
8462:   PetscCall(DMRestoreLocalVector(dm, &localX));
8463:   PetscCall(DMRestoreLocalVector(dmIn, &localU));
8464:   PetscFunctionReturn(PETSC_SUCCESS);
8465: }

8467: /*@C
8468:   DMProjectBdFieldLabelLocal - This projects the given function of the input fields into the function space provided, putting the coefficients in a local vector, calculating only over the portion of the domain boundary specified by the label.

8470:   Not Collective

8472:   Input Parameters:
8473: + dm     - The `DM`
8474: . time   - The time
8475: . label  - The `DMLabel` marking the portion of the domain boundary to output
8476: . numIds - The number of label ids to use
8477: . ids    - The label ids to use for marking
8478: . Nc     - The number of components to set in the output, or `PETSC_DETERMINE` for all components
8479: . comps  - The components to set in the output, or `NULL` for all components
8480: . localU - The input field vector
8481: . funcs  - The functions to evaluate, one per field
8482: - mode   - The insertion mode for values

8484:   Output Parameter:
8485: . localX - The output vector

8487:   Calling sequence of `funcs`:
8488: + dim          - The spatial dimension
8489: . Nf           - The number of input fields
8490: . NfAux        - The number of input auxiliary fields
8491: . uOff         - The offset of each field in u[]
8492: . uOff_x       - The offset of each field in u_x[]
8493: . u            - The field values at this point in space
8494: . u_t          - The field time derivative at this point in space (or NULL)
8495: . u_x          - The field derivatives at this point in space
8496: . aOff         - The offset of each auxiliary field in u[]
8497: . aOff_x       - The offset of each auxiliary field in u_x[]
8498: . a            - The auxiliary field values at this point in space
8499: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8500: . a_x          - The auxiliary field derivatives at this point in space
8501: . t            - The current time
8502: . x            - The coordinates of this point
8503: . n            - The face normal
8504: . numConstants - The number of constants
8505: . constants    - The value of each constant
8506: - f            - The value of the function at this point in space

8508:   Level: intermediate

8510:   Note:
8511:   There are three different `DM`s that potentially interact in this function. The output `DM`, dm, specifies the layout of the values calculates by funcs.
8512:   The input `DM`, attached to U, may be different. For example, you can input the solution over the full domain, but output over a piece of the boundary, or
8513:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8514:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8516:   Developer Notes:
8517:   This API is specific to only particular usage of `DM`

8519:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8521: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabelLocal()`, `DMProjectFunction()`, `DMComputeL2Diff()`
8522: @*/
8523: PetscErrorCode DMProjectBdFieldLabelLocal(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Nc, const PetscInt comps[], Vec localU, void (**funcs)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]), InsertMode mode, Vec localX)
8524: {
8525:   PetscFunctionBegin;
8529:   PetscUseTypeMethod(dm, projectbdfieldlabellocal, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX);
8530:   PetscFunctionReturn(PETSC_SUCCESS);
8531: }

8533: /*@C
8534:   DMComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.

8536:   Collective

8538:   Input Parameters:
8539: + dm    - The `DM`
8540: . time  - The time
8541: . funcs - The functions to evaluate for each field component
8542: . ctxs  - Optional array of contexts to pass to each function, or NULL.
8543: - X     - The coefficient vector u_h, a global vector

8545:   Output Parameter:
8546: . diff - The diff ||u - u_h||_2

8548:   Level: developer

8550:   Developer Notes:
8551:   This API is specific to only particular usage of `DM`

8553:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8555: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMComputeL2FieldDiff()`, `DMComputeL2GradientDiff()`
8556: @*/
8557: PetscErrorCode DMComputeL2Diff(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
8558: {
8559:   PetscFunctionBegin;
8562:   PetscUseTypeMethod(dm, computel2diff, time, funcs, ctxs, X, diff);
8563:   PetscFunctionReturn(PETSC_SUCCESS);
8564: }

8566: /*@C
8567:   DMComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h.

8569:   Collective

8571:   Input Parameters:
8572: + dm    - The `DM`
8573: . time  - The time
8574: . funcs - The gradient functions to evaluate for each field component
8575: . ctxs  - Optional array of contexts to pass to each function, or NULL.
8576: . X     - The coefficient vector u_h, a global vector
8577: - n     - The vector to project along

8579:   Output Parameter:
8580: . diff - The diff ||(grad u - grad u_h) . n||_2

8582:   Level: developer

8584:   Developer Notes:
8585:   This API is specific to only particular usage of `DM`

8587:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8589: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMComputeL2Diff()`, `DMComputeL2FieldDiff()`
8590: @*/
8591: PetscErrorCode DMComputeL2GradientDiff(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
8592: {
8593:   PetscFunctionBegin;
8596:   PetscUseTypeMethod(dm, computel2gradientdiff, time, funcs, ctxs, X, n, diff);
8597:   PetscFunctionReturn(PETSC_SUCCESS);
8598: }

8600: /*@C
8601:   DMComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components.

8603:   Collective

8605:   Input Parameters:
8606: + dm    - The `DM`
8607: . time  - The time
8608: . funcs - The functions to evaluate for each field component
8609: . ctxs  - Optional array of contexts to pass to each function, or NULL.
8610: - X     - The coefficient vector u_h, a global vector

8612:   Output Parameter:
8613: . diff - The array of differences, ||u^f - u^f_h||_2

8615:   Level: developer

8617:   Developer Notes:
8618:   This API is specific to only particular usage of `DM`

8620:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8622: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMComputeL2GradientDiff()`
8623: @*/
8624: PetscErrorCode DMComputeL2FieldDiff(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
8625: {
8626:   PetscFunctionBegin;
8629:   PetscUseTypeMethod(dm, computel2fielddiff, time, funcs, ctxs, X, diff);
8630:   PetscFunctionReturn(PETSC_SUCCESS);
8631: }

8633: /*@C
8634:   DMGetNeighbors - Gets an array containing the MPI ranks of all the processes neighbors

8636:   Not Collective

8638:   Input Parameter:
8639: . dm - The `DM`

8641:   Output Parameters:
8642: + nranks - the number of neighbours
8643: - ranks  - the neighbors ranks

8645:   Level: beginner

8647:   Note:
8648:   Do not free the array, it is freed when the `DM` is destroyed.

8650: .seealso: [](ch_dmbase), `DM`, `DMDAGetNeighbors()`, `PetscSFGetRootRanks()`
8651: @*/
8652: PetscErrorCode DMGetNeighbors(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
8653: {
8654:   PetscFunctionBegin;
8656:   PetscUseTypeMethod(dm, getneighbors, nranks, ranks);
8657:   PetscFunctionReturn(PETSC_SUCCESS);
8658: }

8660: #include <petsc/private/matimpl.h>

8662: /*
8663:     Converts the input vector to a ghosted vector and then calls the standard coloring code.
8664:     This must be a different function because it requires DM which is not defined in the Mat library
8665: */
8666: static PetscErrorCode MatFDColoringApply_AIJDM(Mat J, MatFDColoring coloring, Vec x1, void *sctx)
8667: {
8668:   PetscFunctionBegin;
8669:   if (coloring->ctype == IS_COLORING_LOCAL) {
8670:     Vec x1local;
8671:     DM  dm;
8672:     PetscCall(MatGetDM(J, &dm));
8673:     PetscCheck(dm, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_INCOMP, "IS_COLORING_LOCAL requires a DM");
8674:     PetscCall(DMGetLocalVector(dm, &x1local));
8675:     PetscCall(DMGlobalToLocalBegin(dm, x1, INSERT_VALUES, x1local));
8676:     PetscCall(DMGlobalToLocalEnd(dm, x1, INSERT_VALUES, x1local));
8677:     x1 = x1local;
8678:   }
8679:   PetscCall(MatFDColoringApply_AIJ(J, coloring, x1, sctx));
8680:   if (coloring->ctype == IS_COLORING_LOCAL) {
8681:     DM dm;
8682:     PetscCall(MatGetDM(J, &dm));
8683:     PetscCall(DMRestoreLocalVector(dm, &x1));
8684:   }
8685:   PetscFunctionReturn(PETSC_SUCCESS);
8686: }

8688: /*@
8689:   MatFDColoringUseDM - allows a `MatFDColoring` object to use the `DM` associated with the matrix to compute a `IS_COLORING_LOCAL` coloring

8691:   Input Parameters:
8692: + coloring   - The matrix to get the `DM` from
8693: - fdcoloring - the `MatFDColoring` object

8695:   Level: advanced

8697:   Developer Note:
8698:   This routine exists because the PETSc `Mat` library does not know about the `DM` objects

8700: .seealso: [](ch_dmbase), `DM`, `MatFDColoring`, `MatFDColoringCreate()`, `ISColoringType`
8701: @*/
8702: PetscErrorCode MatFDColoringUseDM(Mat coloring, MatFDColoring fdcoloring)
8703: {
8704:   PetscFunctionBegin;
8705:   coloring->ops->fdcoloringapply = MatFDColoringApply_AIJDM;
8706:   PetscFunctionReturn(PETSC_SUCCESS);
8707: }

8709: /*@
8710:   DMGetCompatibility - determine if two `DM`s are compatible

8712:   Collective

8714:   Input Parameters:
8715: + dm1 - the first `DM`
8716: - dm2 - the second `DM`

8718:   Output Parameters:
8719: + compatible - whether or not the two `DM`s are compatible
8720: - set        - whether or not the compatible value was actually determined and set

8722:   Level: advanced

8724:   Notes:
8725:   Two `DM`s are deemed compatible if they represent the same parallel decomposition
8726:   of the same topology. This implies that the section (field data) on one
8727:   "makes sense" with respect to the topology and parallel decomposition of the other.
8728:   Loosely speaking, compatible `DM`s represent the same domain and parallel
8729:   decomposition, but hold different data.

8731:   Typically, one would confirm compatibility if intending to simultaneously iterate
8732:   over a pair of vectors obtained from different `DM`s.

8734:   For example, two `DMDA` objects are compatible if they have the same local
8735:   and global sizes and the same stencil width. They can have different numbers
8736:   of degrees of freedom per node. Thus, one could use the node numbering from
8737:   either `DM` in bounds for a loop over vectors derived from either `DM`.

8739:   Consider the operation of summing data living on a 2-dof `DMDA` to data living
8740:   on a 1-dof `DMDA`, which should be compatible, as in the following snippet.
8741: .vb
8742:   ...
8743:   PetscCall(DMGetCompatibility(da1,da2,&compatible,&set));
8744:   if (set && compatible)  {
8745:     PetscCall(DMDAVecGetArrayDOF(da1,vec1,&arr1));
8746:     PetscCall(DMDAVecGetArrayDOF(da2,vec2,&arr2));
8747:     PetscCall(DMDAGetCorners(da1,&x,&y,NULL,&m,&n,NULL));
8748:     for (j=y; j<y+n; ++j) {
8749:       for (i=x; i<x+m, ++i) {
8750:         arr1[j][i][0] = arr2[j][i][0] + arr2[j][i][1];
8751:       }
8752:     }
8753:     PetscCall(DMDAVecRestoreArrayDOF(da1,vec1,&arr1));
8754:     PetscCall(DMDAVecRestoreArrayDOF(da2,vec2,&arr2));
8755:   } else {
8756:     SETERRQ(PetscObjectComm((PetscObject)da1,PETSC_ERR_ARG_INCOMP,"DMDA objects incompatible");
8757:   }
8758:   ...
8759: .ve

8761:   Checking compatibility might be expensive for a given implementation of `DM`,
8762:   or might be impossible to unambiguously confirm or deny. For this reason,
8763:   this function may decline to determine compatibility, and hence users should
8764:   always check the "set" output parameter.

8766:   A `DM` is always compatible with itself.

8768:   In the current implementation, `DM`s which live on "unequal" communicators
8769:   (MPI_UNEQUAL in the terminology of MPI_Comm_compare()) are always deemed
8770:   incompatible.

8772:   This function is labeled "Collective," as information about all subdomains
8773:   is required on each rank. However, in `DM` implementations which store all this
8774:   information locally, this function may be merely "Logically Collective".

8776:   Developer Note:
8777:   Compatibility is assumed to be a symmetric concept; `DM` A is compatible with `DM` B
8778:   iff B is compatible with A. Thus, this function checks the implementations
8779:   of both dm and dmc (if they are of different types), attempting to determine
8780:   compatibility. It is left to `DM` implementers to ensure that symmetry is
8781:   preserved. The simplest way to do this is, when implementing type-specific
8782:   logic for this function, is to check for existing logic in the implementation
8783:   of other `DM` types and let *set = PETSC_FALSE if found.

8785: .seealso: [](ch_dmbase), `DM`, `DMDACreateCompatibleDMDA()`, `DMStagCreateCompatibleDMStag()`
8786: @*/
8787: PetscErrorCode DMGetCompatibility(DM dm1, DM dm2, PetscBool *compatible, PetscBool *set)
8788: {
8789:   PetscMPIInt compareResult;
8790:   DMType      type, type2;
8791:   PetscBool   sameType;

8793:   PetscFunctionBegin;

8797:   /* Declare a DM compatible with itself */
8798:   if (dm1 == dm2) {
8799:     *set        = PETSC_TRUE;
8800:     *compatible = PETSC_TRUE;
8801:     PetscFunctionReturn(PETSC_SUCCESS);
8802:   }

8804:   /* Declare a DM incompatible with a DM that lives on an "unequal"
8805:      communicator. Note that this does not preclude compatibility with
8806:      DMs living on "congruent" or "similar" communicators, but this must be
8807:      determined by the implementation-specific logic */
8808:   PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dm1), PetscObjectComm((PetscObject)dm2), &compareResult));
8809:   if (compareResult == MPI_UNEQUAL) {
8810:     *set        = PETSC_TRUE;
8811:     *compatible = PETSC_FALSE;
8812:     PetscFunctionReturn(PETSC_SUCCESS);
8813:   }

8815:   /* Pass to the implementation-specific routine, if one exists. */
8816:   if (dm1->ops->getcompatibility) {
8817:     PetscUseTypeMethod(dm1, getcompatibility, dm2, compatible, set);
8818:     if (*set) PetscFunctionReturn(PETSC_SUCCESS);
8819:   }

8821:   /* If dm1 and dm2 are of different types, then attempt to check compatibility
8822:      with an implementation of this function from dm2 */
8823:   PetscCall(DMGetType(dm1, &type));
8824:   PetscCall(DMGetType(dm2, &type2));
8825:   PetscCall(PetscStrcmp(type, type2, &sameType));
8826:   if (!sameType && dm2->ops->getcompatibility) {
8827:     PetscUseTypeMethod(dm2, getcompatibility, dm1, compatible, set); /* Note argument order */
8828:   } else {
8829:     *set = PETSC_FALSE;
8830:   }
8831:   PetscFunctionReturn(PETSC_SUCCESS);
8832: }

8834: /*@C
8835:   DMMonitorSet - Sets an additional monitor function that is to be used after a solve to monitor discretization performance.

8837:   Logically Collective

8839:   Input Parameters:
8840: + dm             - the `DM`
8841: . f              - the monitor function
8842: . mctx           - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired)
8843: - monitordestroy - [optional] routine that frees monitor context (may be `NULL`), see `PetscCtxDestroyFn` for the calling sequence

8845:   Options Database Key:
8846: . -dm_monitor_cancel - cancels all monitors that have been hardwired into a code by calls to `DMMonitorSet()`, but
8847:                        does not cancel those set via the options database.

8849:   Level: intermediate

8851:   Note:
8852:   Several different monitoring routines may be set by calling
8853:   `DMMonitorSet()` multiple times or with `DMMonitorSetFromOptions()`; all will be called in the
8854:   order in which they were set.

8856:   Fortran Note:
8857:   Only a single monitor function can be set for each `DM` object

8859:   Developer Note:
8860:   This API has a generic name but seems specific to a very particular aspect of the use of `DM`

8862: .seealso: [](ch_dmbase), `DM`, `DMMonitorCancel()`, `DMMonitorSetFromOptions()`, `DMMonitor()`, `PetscCtxDestroyFn`
8863: @*/
8864: PetscErrorCode DMMonitorSet(DM dm, PetscErrorCode (*f)(DM, void *), void *mctx, PetscCtxDestroyFn *monitordestroy)
8865: {
8866:   PetscInt m;

8868:   PetscFunctionBegin;
8870:   for (m = 0; m < dm->numbermonitors; ++m) {
8871:     PetscBool identical;

8873:     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))f, mctx, monitordestroy, (PetscErrorCode (*)(void))dm->monitor[m], dm->monitorcontext[m], dm->monitordestroy[m], &identical));
8874:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
8875:   }
8876:   PetscCheck(dm->numbermonitors < MAXDMMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
8877:   dm->monitor[dm->numbermonitors]          = f;
8878:   dm->monitordestroy[dm->numbermonitors]   = monitordestroy;
8879:   dm->monitorcontext[dm->numbermonitors++] = mctx;
8880:   PetscFunctionReturn(PETSC_SUCCESS);
8881: }

8883: /*@
8884:   DMMonitorCancel - Clears all the monitor functions for a `DM` object.

8886:   Logically Collective

8888:   Input Parameter:
8889: . dm - the DM

8891:   Options Database Key:
8892: . -dm_monitor_cancel - cancels all monitors that have been hardwired
8893:   into a code by calls to `DMonitorSet()`, but does not cancel those
8894:   set via the options database

8896:   Level: intermediate

8898:   Note:
8899:   There is no way to clear one specific monitor from a `DM` object.

8901: .seealso: [](ch_dmbase), `DM`, `DMMonitorSet()`, `DMMonitorSetFromOptions()`, `DMMonitor()`
8902: @*/
8903: PetscErrorCode DMMonitorCancel(DM dm)
8904: {
8905:   PetscInt m;

8907:   PetscFunctionBegin;
8909:   for (m = 0; m < dm->numbermonitors; ++m) {
8910:     if (dm->monitordestroy[m]) PetscCall((*dm->monitordestroy[m])(&dm->monitorcontext[m]));
8911:   }
8912:   dm->numbermonitors = 0;
8913:   PetscFunctionReturn(PETSC_SUCCESS);
8914: }

8916: /*@C
8917:   DMMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user

8919:   Collective

8921:   Input Parameters:
8922: + dm           - `DM` object you wish to monitor
8923: . name         - the monitor type one is seeking
8924: . help         - message indicating what monitoring is done
8925: . manual       - manual page for the monitor
8926: . monitor      - the monitor function, this must use a `PetscViewerFormat` as its context
8927: - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the `DM` or `PetscViewer` objects

8929:   Output Parameter:
8930: . flg - Flag set if the monitor was created

8932:   Level: developer

8934: .seealso: [](ch_dmbase), `DM`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
8935:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
8936:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
8937:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
8938:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
8939:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
8940:           `PetscOptionsFList()`, `PetscOptionsEList()`, `DMMonitor()`, `DMMonitorSet()`
8941: @*/
8942: PetscErrorCode DMMonitorSetFromOptions(DM dm, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(DM, void *), PetscErrorCode (*monitorsetup)(DM, PetscViewerAndFormat *), PetscBool *flg)
8943: {
8944:   PetscViewer       viewer;
8945:   PetscViewerFormat format;

8947:   PetscFunctionBegin;
8949:   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)dm), ((PetscObject)dm)->options, ((PetscObject)dm)->prefix, name, &viewer, &format, flg));
8950:   if (*flg) {
8951:     PetscViewerAndFormat *vf;

8953:     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
8954:     PetscCall(PetscViewerDestroy(&viewer));
8955:     if (monitorsetup) PetscCall((*monitorsetup)(dm, vf));
8956:     PetscCall(DMMonitorSet(dm, monitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy));
8957:   }
8958:   PetscFunctionReturn(PETSC_SUCCESS);
8959: }

8961: /*@
8962:   DMMonitor - runs the user provided monitor routines, if they exist

8964:   Collective

8966:   Input Parameter:
8967: . dm - The `DM`

8969:   Level: developer

8971:   Developer Note:
8972:   Note should indicate when during the life of the `DM` the monitor is run. It appears to be
8973:   related to the discretization process seems rather specialized since some `DM` have no
8974:   concept of discretization.

8976: .seealso: [](ch_dmbase), `DM`, `DMMonitorSet()`, `DMMonitorSetFromOptions()`
8977: @*/
8978: PetscErrorCode DMMonitor(DM dm)
8979: {
8980:   PetscInt m;

8982:   PetscFunctionBegin;
8983:   if (!dm) PetscFunctionReturn(PETSC_SUCCESS);
8985:   for (m = 0; m < dm->numbermonitors; ++m) PetscCall((*dm->monitor[m])(dm, dm->monitorcontext[m]));
8986:   PetscFunctionReturn(PETSC_SUCCESS);
8987: }

8989: /*@
8990:   DMComputeError - Computes the error assuming the user has provided the exact solution functions

8992:   Collective

8994:   Input Parameters:
8995: + dm  - The `DM`
8996: - sol - The solution vector

8998:   Input/Output Parameter:
8999: . errors - An array of length Nf, the number of fields, or `NULL` for no output; on output
9000:            contains the error in each field

9002:   Output Parameter:
9003: . errorVec - A vector to hold the cellwise error (may be `NULL`)

9005:   Level: developer

9007:   Note:
9008:   The exact solutions come from the `PetscDS` object, and the time comes from `DMGetOutputSequenceNumber()`.

9010: .seealso: [](ch_dmbase), `DM`, `DMMonitorSet()`, `DMGetRegionNumDS()`, `PetscDSGetExactSolution()`, `DMGetOutputSequenceNumber()`
9011: @*/
9012: PetscErrorCode DMComputeError(DM dm, Vec sol, PetscReal errors[], Vec *errorVec)
9013: {
9014:   PetscErrorCode (**exactSol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
9015:   void    **ctxs;
9016:   PetscReal time;
9017:   PetscInt  Nf, f, Nds, s;

9019:   PetscFunctionBegin;
9020:   PetscCall(DMGetNumFields(dm, &Nf));
9021:   PetscCall(PetscCalloc2(Nf, &exactSol, Nf, &ctxs));
9022:   PetscCall(DMGetNumDS(dm, &Nds));
9023:   for (s = 0; s < Nds; ++s) {
9024:     PetscDS         ds;
9025:     DMLabel         label;
9026:     IS              fieldIS;
9027:     const PetscInt *fields;
9028:     PetscInt        dsNf;

9030:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL));
9031:     PetscCall(PetscDSGetNumFields(ds, &dsNf));
9032:     if (fieldIS) PetscCall(ISGetIndices(fieldIS, &fields));
9033:     for (f = 0; f < dsNf; ++f) {
9034:       const PetscInt field = fields[f];
9035:       PetscCall(PetscDSGetExactSolution(ds, field, &exactSol[field], &ctxs[field]));
9036:     }
9037:     if (fieldIS) PetscCall(ISRestoreIndices(fieldIS, &fields));
9038:   }
9039:   for (f = 0; f < Nf; ++f) PetscCheck(exactSol[f], PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DS must contain exact solution functions in order to calculate error, missing for field %" PetscInt_FMT, f);
9040:   PetscCall(DMGetOutputSequenceNumber(dm, NULL, &time));
9041:   if (errors) PetscCall(DMComputeL2FieldDiff(dm, time, exactSol, ctxs, sol, errors));
9042:   if (errorVec) {
9043:     DM             edm;
9044:     DMPolytopeType ct;
9045:     PetscBool      simplex;
9046:     PetscInt       dim, cStart, Nf;

9048:     PetscCall(DMClone(dm, &edm));
9049:     PetscCall(DMGetDimension(edm, &dim));
9050:     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
9051:     PetscCall(DMPlexGetCellType(dm, cStart, &ct));
9052:     simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
9053:     PetscCall(DMGetNumFields(dm, &Nf));
9054:     for (f = 0; f < Nf; ++f) {
9055:       PetscFE         fe, efe;
9056:       PetscQuadrature q;
9057:       const char     *name;

9059:       PetscCall(DMGetField(dm, f, NULL, (PetscObject *)&fe));
9060:       PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, Nf, simplex, 0, PETSC_DETERMINE, &efe));
9061:       PetscCall(PetscObjectGetName((PetscObject)fe, &name));
9062:       PetscCall(PetscObjectSetName((PetscObject)efe, name));
9063:       PetscCall(PetscFEGetQuadrature(fe, &q));
9064:       PetscCall(PetscFESetQuadrature(efe, q));
9065:       PetscCall(DMSetField(edm, f, NULL, (PetscObject)efe));
9066:       PetscCall(PetscFEDestroy(&efe));
9067:     }
9068:     PetscCall(DMCreateDS(edm));

9070:     PetscCall(DMCreateGlobalVector(edm, errorVec));
9071:     PetscCall(PetscObjectSetName((PetscObject)*errorVec, "Error"));
9072:     PetscCall(DMPlexComputeL2DiffVec(dm, time, exactSol, ctxs, sol, *errorVec));
9073:     PetscCall(DMDestroy(&edm));
9074:   }
9075:   PetscCall(PetscFree2(exactSol, ctxs));
9076:   PetscFunctionReturn(PETSC_SUCCESS);
9077: }

9079: /*@
9080:   DMGetNumAuxiliaryVec - Get the number of auxiliary vectors associated with this `DM`

9082:   Not Collective

9084:   Input Parameter:
9085: . dm - The `DM`

9087:   Output Parameter:
9088: . numAux - The number of auxiliary data vectors

9090:   Level: advanced

9092: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMSetAuxiliaryVec()`, `DMGetAuxiliaryLabels()`, `DMGetAuxiliaryVec()`
9093: @*/
9094: PetscErrorCode DMGetNumAuxiliaryVec(DM dm, PetscInt *numAux)
9095: {
9096:   PetscFunctionBegin;
9098:   PetscCall(PetscHMapAuxGetSize(dm->auxData, numAux));
9099:   PetscFunctionReturn(PETSC_SUCCESS);
9100: }

9102: /*@
9103:   DMGetAuxiliaryVec - Get the auxiliary vector for region specified by the given label and value, and equation part

9105:   Not Collective

9107:   Input Parameters:
9108: + dm    - The `DM`
9109: . label - The `DMLabel`
9110: . value - The label value indicating the region
9111: - part  - The equation part, or 0 if unused

9113:   Output Parameter:
9114: . aux - The `Vec` holding auxiliary field data

9116:   Level: advanced

9118:   Note:
9119:   If no auxiliary vector is found for this (label, value), (NULL, 0, 0) is checked as well.

9121: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMSetAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryLabels()`
9122: @*/
9123: PetscErrorCode DMGetAuxiliaryVec(DM dm, DMLabel label, PetscInt value, PetscInt part, Vec *aux)
9124: {
9125:   PetscHashAuxKey key, wild = {NULL, 0, 0};
9126:   PetscBool       has;

9128:   PetscFunctionBegin;
9131:   key.label = label;
9132:   key.value = value;
9133:   key.part  = part;
9134:   PetscCall(PetscHMapAuxHas(dm->auxData, key, &has));
9135:   if (has) PetscCall(PetscHMapAuxGet(dm->auxData, key, aux));
9136:   else PetscCall(PetscHMapAuxGet(dm->auxData, wild, aux));
9137:   PetscFunctionReturn(PETSC_SUCCESS);
9138: }

9140: /*@
9141:   DMSetAuxiliaryVec - Set an auxiliary vector for region specified by the given label and value, and equation part

9143:   Not Collective because auxiliary vectors are not parallel

9145:   Input Parameters:
9146: + dm    - The `DM`
9147: . label - The `DMLabel`
9148: . value - The label value indicating the region
9149: . part  - The equation part, or 0 if unused
9150: - aux   - The `Vec` holding auxiliary field data

9152:   Level: advanced

9154: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMGetAuxiliaryLabels()`, `DMCopyAuxiliaryVec()`
9155: @*/
9156: PetscErrorCode DMSetAuxiliaryVec(DM dm, DMLabel label, PetscInt value, PetscInt part, Vec aux)
9157: {
9158:   Vec             old;
9159:   PetscHashAuxKey key;

9161:   PetscFunctionBegin;
9164:   key.label = label;
9165:   key.value = value;
9166:   key.part  = part;
9167:   PetscCall(PetscHMapAuxGet(dm->auxData, key, &old));
9168:   PetscCall(PetscObjectReference((PetscObject)aux));
9169:   if (!aux) PetscCall(PetscHMapAuxDel(dm->auxData, key));
9170:   else PetscCall(PetscHMapAuxSet(dm->auxData, key, aux));
9171:   PetscCall(VecDestroy(&old));
9172:   PetscFunctionReturn(PETSC_SUCCESS);
9173: }

9175: /*@
9176:   DMGetAuxiliaryLabels - Get the labels, values, and parts for all auxiliary vectors in this `DM`

9178:   Not Collective

9180:   Input Parameter:
9181: . dm - The `DM`

9183:   Output Parameters:
9184: + labels - The `DMLabel`s for each `Vec`
9185: . values - The label values for each `Vec`
9186: - parts  - The equation parts for each `Vec`

9188:   Level: advanced

9190:   Note:
9191:   The arrays passed in must be at least as large as `DMGetNumAuxiliaryVec()`.

9193: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMSetAuxiliaryVec()`, `DMCopyAuxiliaryVec()`
9194: @*/
9195: PetscErrorCode DMGetAuxiliaryLabels(DM dm, DMLabel labels[], PetscInt values[], PetscInt parts[])
9196: {
9197:   PetscHashAuxKey *keys;
9198:   PetscInt         n, i, off = 0;

9200:   PetscFunctionBegin;
9202:   PetscAssertPointer(labels, 2);
9203:   PetscAssertPointer(values, 3);
9204:   PetscAssertPointer(parts, 4);
9205:   PetscCall(DMGetNumAuxiliaryVec(dm, &n));
9206:   PetscCall(PetscMalloc1(n, &keys));
9207:   PetscCall(PetscHMapAuxGetKeys(dm->auxData, &off, keys));
9208:   for (i = 0; i < n; ++i) {
9209:     labels[i] = keys[i].label;
9210:     values[i] = keys[i].value;
9211:     parts[i]  = keys[i].part;
9212:   }
9213:   PetscCall(PetscFree(keys));
9214:   PetscFunctionReturn(PETSC_SUCCESS);
9215: }

9217: /*@
9218:   DMCopyAuxiliaryVec - Copy the auxiliary vector data on a `DM` to a new `DM`

9220:   Not Collective

9222:   Input Parameter:
9223: . dm - The `DM`

9225:   Output Parameter:
9226: . dmNew - The new `DM`, now with the same auxiliary data

9228:   Level: advanced

9230:   Note:
9231:   This is a shallow copy of the auxiliary vectors

9233: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMSetAuxiliaryVec()`
9234: @*/
9235: PetscErrorCode DMCopyAuxiliaryVec(DM dm, DM dmNew)
9236: {
9237:   PetscFunctionBegin;
9240:   if (dm == dmNew) PetscFunctionReturn(PETSC_SUCCESS);
9241:   PetscCall(DMClearAuxiliaryVec(dmNew));

9243:   PetscCall(PetscHMapAuxDestroy(&dmNew->auxData));
9244:   PetscCall(PetscHMapAuxDuplicate(dm->auxData, &dmNew->auxData));
9245:   {
9246:     Vec     *auxData;
9247:     PetscInt n, i, off = 0;

9249:     PetscCall(PetscHMapAuxGetSize(dmNew->auxData, &n));
9250:     PetscCall(PetscMalloc1(n, &auxData));
9251:     PetscCall(PetscHMapAuxGetVals(dmNew->auxData, &off, auxData));
9252:     for (i = 0; i < n; ++i) PetscCall(PetscObjectReference((PetscObject)auxData[i]));
9253:     PetscCall(PetscFree(auxData));
9254:   }
9255:   PetscFunctionReturn(PETSC_SUCCESS);
9256: }

9258: /*@
9259:   DMClearAuxiliaryVec - Destroys the auxiliary vector information and creates a new empty one

9261:   Not Collective

9263:   Input Parameter:
9264: . dm - The `DM`

9266:   Level: advanced

9268: .seealso: [](ch_dmbase), `DM`, `DMCopyAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMSetAuxiliaryVec()`
9269: @*/
9270: PetscErrorCode DMClearAuxiliaryVec(DM dm)
9271: {
9272:   Vec     *auxData;
9273:   PetscInt n, i, off = 0;

9275:   PetscFunctionBegin;
9276:   PetscCall(PetscHMapAuxGetSize(dm->auxData, &n));
9277:   PetscCall(PetscMalloc1(n, &auxData));
9278:   PetscCall(PetscHMapAuxGetVals(dm->auxData, &off, auxData));
9279:   for (i = 0; i < n; ++i) PetscCall(VecDestroy(&auxData[i]));
9280:   PetscCall(PetscFree(auxData));
9281:   PetscCall(PetscHMapAuxDestroy(&dm->auxData));
9282:   PetscCall(PetscHMapAuxCreate(&dm->auxData));
9283:   PetscFunctionReturn(PETSC_SUCCESS);
9284: }

9286: /*@
9287:   DMPolytopeMatchOrientation - Determine an orientation (transformation) that takes the source face arrangement to the target face arrangement

9289:   Not Collective

9291:   Input Parameters:
9292: + ct         - The `DMPolytopeType`
9293: . sourceCone - The source arrangement of faces
9294: - targetCone - The target arrangement of faces

9296:   Output Parameters:
9297: + ornt  - The orientation (transformation) which will take the source arrangement to the target arrangement
9298: - found - Flag indicating that a suitable orientation was found

9300:   Level: advanced

9302:   Note:
9303:   An arrangement is a face order combined with an orientation for each face

9305:   Each orientation (transformation) is labeled with an integer from negative `DMPolytopeTypeGetNumArrangements(ct)`/2 to `DMPolytopeTypeGetNumArrangements(ct)`/2
9306:   that labels each arrangement (face ordering plus orientation for each face).

9308:   See `DMPolytopeMatchVertexOrientation()` to find a new vertex orientation that takes the source vertex arrangement to the target vertex arrangement

9310: .seealso: [](ch_dmbase), `DM`, `DMPolytopeGetOrientation()`, `DMPolytopeMatchVertexOrientation()`, `DMPolytopeGetVertexOrientation()`
9311: @*/
9312: PetscErrorCode DMPolytopeMatchOrientation(DMPolytopeType ct, const PetscInt sourceCone[], const PetscInt targetCone[], PetscInt *ornt, PetscBool *found)
9313: {
9314:   const PetscInt cS = DMPolytopeTypeGetConeSize(ct);
9315:   const PetscInt nO = DMPolytopeTypeGetNumArrangements(ct) / 2;
9316:   PetscInt       o, c;

9318:   PetscFunctionBegin;
9319:   if (!nO) {
9320:     *ornt  = 0;
9321:     *found = PETSC_TRUE;
9322:     PetscFunctionReturn(PETSC_SUCCESS);
9323:   }
9324:   for (o = -nO; o < nO; ++o) {
9325:     const PetscInt *arr = DMPolytopeTypeGetArrangement(ct, o);

9327:     for (c = 0; c < cS; ++c)
9328:       if (sourceCone[arr[c * 2]] != targetCone[c]) break;
9329:     if (c == cS) {
9330:       *ornt = o;
9331:       break;
9332:     }
9333:   }
9334:   *found = o == nO ? PETSC_FALSE : PETSC_TRUE;
9335:   PetscFunctionReturn(PETSC_SUCCESS);
9336: }

9338: /*@
9339:   DMPolytopeGetOrientation - Determine an orientation (transformation) that takes the source face arrangement to the target face arrangement

9341:   Not Collective

9343:   Input Parameters:
9344: + ct         - The `DMPolytopeType`
9345: . sourceCone - The source arrangement of faces
9346: - targetCone - The target arrangement of faces

9348:   Output Parameter:
9349: . ornt - The orientation (transformation) which will take the source arrangement to the target arrangement

9351:   Level: advanced

9353:   Note:
9354:   This function is the same as `DMPolytopeMatchOrientation()` except it will generate an error if no suitable orientation can be found.

9356:   Developer Note:
9357:   It is unclear why this function needs to exist since one can simply call `DMPolytopeMatchOrientation()` and error if none is found

9359: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMPolytopeMatchOrientation()`, `DMPolytopeGetVertexOrientation()`, `DMPolytopeMatchVertexOrientation()`
9360: @*/
9361: PetscErrorCode DMPolytopeGetOrientation(DMPolytopeType ct, const PetscInt sourceCone[], const PetscInt targetCone[], PetscInt *ornt)
9362: {
9363:   PetscBool found;

9365:   PetscFunctionBegin;
9366:   PetscCall(DMPolytopeMatchOrientation(ct, sourceCone, targetCone, ornt, &found));
9367:   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find orientation for %s", DMPolytopeTypes[ct]);
9368:   PetscFunctionReturn(PETSC_SUCCESS);
9369: }

9371: /*@
9372:   DMPolytopeMatchVertexOrientation - Determine an orientation (transformation) that takes the source vertex arrangement to the target vertex arrangement

9374:   Not Collective

9376:   Input Parameters:
9377: + ct         - The `DMPolytopeType`
9378: . sourceVert - The source arrangement of vertices
9379: - targetVert - The target arrangement of vertices

9381:   Output Parameters:
9382: + ornt  - The orientation (transformation) which will take the source arrangement to the target arrangement
9383: - found - Flag indicating that a suitable orientation was found

9385:   Level: advanced

9387:   Notes:
9388:   An arrangement is a vertex order

9390:   Each orientation (transformation) is labeled with an integer from negative `DMPolytopeTypeGetNumArrangements(ct)`/2 to `DMPolytopeTypeGetNumArrangements(ct)`/2
9391:   that labels each arrangement (vertex ordering).

9393:   See `DMPolytopeMatchOrientation()` to find a new face orientation that takes the source face arrangement to the target face arrangement

9395: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMPolytopeGetOrientation()`, `DMPolytopeMatchOrientation()`, `DMPolytopeTypeGetNumVertices()`, `DMPolytopeTypeGetVertexArrangement()`
9396: @*/
9397: PetscErrorCode DMPolytopeMatchVertexOrientation(DMPolytopeType ct, const PetscInt sourceVert[], const PetscInt targetVert[], PetscInt *ornt, PetscBool *found)
9398: {
9399:   const PetscInt cS = DMPolytopeTypeGetNumVertices(ct);
9400:   const PetscInt nO = DMPolytopeTypeGetNumArrangements(ct) / 2;
9401:   PetscInt       o, c;

9403:   PetscFunctionBegin;
9404:   if (!nO) {
9405:     *ornt  = 0;
9406:     *found = PETSC_TRUE;
9407:     PetscFunctionReturn(PETSC_SUCCESS);
9408:   }
9409:   for (o = -nO; o < nO; ++o) {
9410:     const PetscInt *arr = DMPolytopeTypeGetVertexArrangement(ct, o);

9412:     for (c = 0; c < cS; ++c)
9413:       if (sourceVert[arr[c]] != targetVert[c]) break;
9414:     if (c == cS) {
9415:       *ornt = o;
9416:       break;
9417:     }
9418:   }
9419:   *found = o == nO ? PETSC_FALSE : PETSC_TRUE;
9420:   PetscFunctionReturn(PETSC_SUCCESS);
9421: }

9423: /*@
9424:   DMPolytopeGetVertexOrientation - Determine an orientation (transformation) that takes the source vertex arrangement to the target vertex arrangement

9426:   Not Collective

9428:   Input Parameters:
9429: + ct         - The `DMPolytopeType`
9430: . sourceCone - The source arrangement of vertices
9431: - targetCone - The target arrangement of vertices

9433:   Output Parameter:
9434: . ornt - The orientation (transformation) which will take the source arrangement to the target arrangement

9436:   Level: advanced

9438:   Note:
9439:   This function is the same as `DMPolytopeMatchVertexOrientation()` except it errors if not orientation is possible.

9441:   Developer Note:
9442:   It is unclear why this function needs to exist since one can simply call `DMPolytopeMatchVertexOrientation()` and error if none is found

9444: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMPolytopeMatchVertexOrientation()`, `DMPolytopeGetOrientation()`
9445: @*/
9446: PetscErrorCode DMPolytopeGetVertexOrientation(DMPolytopeType ct, const PetscInt sourceCone[], const PetscInt targetCone[], PetscInt *ornt)
9447: {
9448:   PetscBool found;

9450:   PetscFunctionBegin;
9451:   PetscCall(DMPolytopeMatchVertexOrientation(ct, sourceCone, targetCone, ornt, &found));
9452:   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find orientation for %s", DMPolytopeTypes[ct]);
9453:   PetscFunctionReturn(PETSC_SUCCESS);
9454: }

9456: /*@
9457:   DMPolytopeInCellTest - Check whether a point lies inside the reference cell of given type

9459:   Not Collective

9461:   Input Parameters:
9462: + ct    - The `DMPolytopeType`
9463: - point - Coordinates of the point

9465:   Output Parameter:
9466: . inside - Flag indicating whether the point is inside the reference cell of given type

9468:   Level: advanced

9470: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMLocatePoints()`
9471: @*/
9472: PetscErrorCode DMPolytopeInCellTest(DMPolytopeType ct, const PetscReal point[], PetscBool *inside)
9473: {
9474:   PetscReal sum = 0.0;
9475:   PetscInt  d;

9477:   PetscFunctionBegin;
9478:   *inside = PETSC_TRUE;
9479:   switch (ct) {
9480:   case DM_POLYTOPE_TRIANGLE:
9481:   case DM_POLYTOPE_TETRAHEDRON:
9482:     for (d = 0; d < DMPolytopeTypeGetDim(ct); ++d) {
9483:       if (point[d] < -1.0) {
9484:         *inside = PETSC_FALSE;
9485:         break;
9486:       }
9487:       sum += point[d];
9488:     }
9489:     if (sum > PETSC_SMALL) {
9490:       *inside = PETSC_FALSE;
9491:       break;
9492:     }
9493:     break;
9494:   case DM_POLYTOPE_QUADRILATERAL:
9495:   case DM_POLYTOPE_HEXAHEDRON:
9496:     for (d = 0; d < DMPolytopeTypeGetDim(ct); ++d)
9497:       if (PetscAbsReal(point[d]) > 1. + PETSC_SMALL) {
9498:         *inside = PETSC_FALSE;
9499:         break;
9500:       }
9501:     break;
9502:   default:
9503:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported polytope type %s", DMPolytopeTypes[ct]);
9504:   }
9505:   PetscFunctionReturn(PETSC_SUCCESS);
9506: }

9508: /*@
9509:   DMReorderSectionSetDefault - Set flag indicating whether the local section should be reordered by default

9511:   Logically collective

9513:   Input Parameters:
9514: + dm      - The DM
9515: - reorder - Flag for reordering

9517:   Level: intermediate

9519: .seealso: `DMReorderSectionGetDefault()`
9520: @*/
9521: PetscErrorCode DMReorderSectionSetDefault(DM dm, DMReorderDefaultFlag reorder)
9522: {
9523:   PetscFunctionBegin;
9525:   PetscTryMethod(dm, "DMReorderSectionSetDefault_C", (DM, DMReorderDefaultFlag), (dm, reorder));
9526:   PetscFunctionReturn(PETSC_SUCCESS);
9527: }

9529: /*@
9530:   DMReorderSectionGetDefault - Get flag indicating whether the local section should be reordered by default

9532:   Not collective

9534:   Input Parameter:
9535: . dm - The DM

9537:   Output Parameter:
9538: . reorder - Flag for reordering

9540:   Level: intermediate

9542: .seealso: `DMReorderSetDefault()`
9543: @*/
9544: PetscErrorCode DMReorderSectionGetDefault(DM dm, DMReorderDefaultFlag *reorder)
9545: {
9546:   PetscFunctionBegin;
9548:   PetscAssertPointer(reorder, 2);
9549:   *reorder = DM_REORDER_DEFAULT_NOTSET;
9550:   PetscTryMethod(dm, "DMReorderSectionGetDefault_C", (DM, DMReorderDefaultFlag *), (dm, reorder));
9551:   PetscFunctionReturn(PETSC_SUCCESS);
9552: }

9554: /*@
9555:   DMReorderSectionSetType - Set the type of local section reordering

9557:   Logically collective

9559:   Input Parameters:
9560: + dm      - The DM
9561: - reorder - The reordering method

9563:   Level: intermediate

9565: .seealso: `DMReorderSectionGetType()`, `DMReorderSectionSetDefault()`
9566: @*/
9567: PetscErrorCode DMReorderSectionSetType(DM dm, MatOrderingType reorder)
9568: {
9569:   PetscFunctionBegin;
9571:   PetscTryMethod(dm, "DMReorderSectionSetType_C", (DM, MatOrderingType), (dm, reorder));
9572:   PetscFunctionReturn(PETSC_SUCCESS);
9573: }

9575: /*@
9576:   DMReorderSectionGetType - Get the reordering type for the local section

9578:   Not collective

9580:   Input Parameter:
9581: . dm - The DM

9583:   Output Parameter:
9584: . reorder - The reordering method

9586:   Level: intermediate

9588: .seealso: `DMReorderSetDefault()`, `DMReorderSectionGetDefault()`
9589: @*/
9590: PetscErrorCode DMReorderSectionGetType(DM dm, MatOrderingType *reorder)
9591: {
9592:   PetscFunctionBegin;
9594:   PetscAssertPointer(reorder, 2);
9595:   *reorder = NULL;
9596:   PetscTryMethod(dm, "DMReorderSectionGetType_C", (DM, MatOrderingType *), (dm, reorder));
9597:   PetscFunctionReturn(PETSC_SUCCESS);
9598: }