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: #if !defined(PETSC_HAVE_WINDOWS_COMPILERS)
 16: #include <petsc/private/valgrind/memcheck.h>
 17: #endif

 19: PetscClassId DM_CLASSID;
 20: PetscClassId DMLABEL_CLASSID;
 21: 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_AdaptInterpolator, DM_ProjectFunction;

 23: const char *const DMBoundaryTypes[]          = {"NONE", "GHOSTED", "MIRROR", "PERIODIC", "TWIST", "DMBoundaryType", "DM_BOUNDARY_", NULL};
 24: const char *const DMBoundaryConditionTypes[] = {"INVALID", "ESSENTIAL", "NATURAL", "INVALID", "INVALID", "ESSENTIAL_FIELD", "NATURAL_FIELD", "INVALID", "INVALID", "ESSENTIAL_BD_FIELD", "NATURAL_RIEMANN", "DMBoundaryConditionType", "DM_BC_", NULL};
 25: const char *const DMBlockingTypes[]          = {"TOPOLOGICAL_POINT", "FIELD_NODE", "DMBlockingType", "DM_BLOCKING_", NULL};
 26: const char *const DMPolytopeTypes[] =
 27:   {"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",
 28:    "unknown", "unknown_cell", "unknown_face",   "invalid",  "DMPolytopeType", "DM_POLYTOPE_", NULL};
 29: const char *const DMCopyLabelsModes[] = {"replace", "keep", "fail", "DMCopyLabelsMode", "DM_COPY_LABELS_", NULL};

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

 35:   Collective

 37:   Input Parameter:
 38: . comm - The communicator for the `DM` object

 40:   Output Parameter:
 41: . dm - The `DM` object

 43:   Level: beginner

 45:   Notes:
 46:   See `DMType` for a brief summary of available `DM`.

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

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

 58:   PetscFunctionBegin;
 59:   PetscAssertPointer(dm, 2);

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

104:   *dm = v;
105:   PetscFunctionReturn(PETSC_SUCCESS);
106: }

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

111:   Collective

113:   Input Parameter:
114: . dm - The original `DM` object

116:   Output Parameter:
117: . newdm - The new `DM` object

119:   Level: beginner

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

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

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

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

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

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

208:     PetscCall(DMGetPeriodicity(dm, &maxCell, &Lstart, &L));
209:     PetscCall(DMSetPeriodicity(*newdm, maxCell, Lstart, L));
210:   }
211:   {
212:     PetscBool useCone, useClosure;

214:     PetscCall(DMGetAdjacency(dm, PETSC_DEFAULT, &useCone, &useClosure));
215:     PetscCall(DMSetAdjacency(*newdm, PETSC_DEFAULT, useCone, useClosure));
216:   }
217:   PetscFunctionReturn(PETSC_SUCCESS);
218: }

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

223:   Logically Collective

225:   Input Parameters:
226: + dm    - initial distributed array
227: - ctype - the vector type, for example `VECSTANDARD`, `VECCUDA`, or `VECVIENNACL`

229:   Options Database Key:
230: . -dm_vec_type ctype - the type of vector to create

232:   Level: intermediate

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

241:   PetscFunctionBegin;
243:   PetscAssertPointer(ctype, 2);
244:   tmp = (char *)dm->vectype;
245:   PetscCall(PetscStrallocpy(ctype, (char **)&dm->vectype));
246:   PetscCall(PetscFree(tmp));
247:   PetscFunctionReturn(PETSC_SUCCESS);
248: }

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

253:   Logically Collective

255:   Input Parameter:
256: . da - initial distributed array

258:   Output Parameter:
259: . ctype - the vector type

261:   Level: intermediate

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

273: /*@
274:   VecGetDM - Gets the `DM` defining the data layout of the vector

276:   Not Collective

278:   Input Parameter:
279: . v - The `Vec`

281:   Output Parameter:
282: . dm - The `DM`

284:   Level: intermediate

286:   Note:
287:   A `Vec` may not have a `DM` associated with it.

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

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

303:   Not Collective

305:   Input Parameters:
306: + v  - The `Vec`
307: - dm - The `DM`

309:   Level: developer

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

314:   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.

316: .seealso: [](ch_dmbase), `DM`, `VecGetDM()`, `DMGetLocalVector()`, `DMGetGlobalVector()`, `DMSetVecType()`
317: @*/
318: PetscErrorCode VecSetDM(Vec v, DM dm)
319: {
320:   PetscFunctionBegin;
323:   PetscCall(PetscObjectCompose((PetscObject)v, "__PETSc_dm", (PetscObject)dm));
324:   PetscFunctionReturn(PETSC_SUCCESS);
325: }

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

330:   Logically Collective

332:   Input Parameters:
333: + dm    - the `DM` context
334: - ctype - the matrix type

336:   Options Database Key:
337: . -dm_is_coloring_type - global or local

339:   Level: intermediate

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

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

355:   Logically Collective

357:   Input Parameter:
358: . dm - the `DM` context

360:   Output Parameter:
361: . ctype - the matrix type

363:   Options Database Key:
364: . -dm_is_coloring_type - global or local

366:   Level: intermediate

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

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

382:   Logically Collective

384:   Input Parameters:
385: + dm    - the `DM` context
386: - ctype - the matrix type, for example `MATMPIAIJ`

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

391:   Level: intermediate

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

399:   PetscFunctionBegin;
401:   PetscAssertPointer(ctype, 2);
402:   tmp = (char *)dm->mattype;
403:   PetscCall(PetscStrallocpy(ctype, (char **)&dm->mattype));
404:   PetscCall(PetscFree(tmp));
405:   PetscFunctionReturn(PETSC_SUCCESS);
406: }

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

411:   Logically Collective

413:   Input Parameter:
414: . dm - the `DM` context

416:   Output Parameter:
417: . ctype - the matrix type

419:   Level: intermediate

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

431: /*@
432:   MatGetDM - Gets the `DM` defining the data layout of the matrix

434:   Not Collective

436:   Input Parameter:
437: . A - The `Mat`

439:   Output Parameter:
440: . dm - The `DM`

442:   Level: intermediate

444:   Note:
445:   A matrix may not have a `DM` associated with it

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

450: .seealso: [](ch_dmbase), `DM`, `MatSetDM()`, `DMCreateMatrix()`, `DMSetMatType()`
451: @*/
452: PetscErrorCode MatGetDM(Mat A, DM *dm)
453: {
454:   PetscFunctionBegin;
456:   PetscAssertPointer(dm, 2);
457:   PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_dm", (PetscObject *)dm));
458:   PetscFunctionReturn(PETSC_SUCCESS);
459: }

461: /*@
462:   MatSetDM - Sets the `DM` defining the data layout of the matrix

464:   Not Collective

466:   Input Parameters:
467: + A  - The `Mat`
468: - dm - The `DM`

470:   Level: developer

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

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

479: .seealso: [](ch_dmbase), `DM`, `MatGetDM()`, `DMCreateMatrix()`, `DMSetMatType()`
480: @*/
481: PetscErrorCode MatSetDM(Mat A, DM dm)
482: {
483:   PetscFunctionBegin;
486:   PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_dm", (PetscObject)dm));
487:   PetscFunctionReturn(PETSC_SUCCESS);
488: }

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

493:   Logically Collective

495:   Input Parameters:
496: + dm     - the `DM` context
497: - prefix - the prefix to prepend

499:   Level: advanced

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

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

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

521:   Logically Collective

523:   Input Parameters:
524: + dm     - the `DM` context
525: - prefix - the string to append to the current prefix

527:   Level: advanced

529:   Note:
530:   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.
531:   A hyphen (-) must NOT be given at the beginning of the prefix name.
532:   The first character of all runtime options is AUTOMATICALLY the hyphen.

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

544: /*@
545:   DMGetOptionsPrefix - Gets the prefix used for searching for all
546:   DM options in the options database.

548:   Not Collective

550:   Input Parameter:
551: . dm - the `DM` context

553:   Output Parameter:
554: . prefix - pointer to the prefix string used is returned

556:   Level: advanced

558:   Fortran Note:
559:   Pass in a string 'prefix' of
560:   sufficient length to hold the prefix.

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

572: static PetscErrorCode DMCountNonCyclicReferences_Internal(DM dm, PetscBool recurseCoarse, PetscBool recurseFine, PetscInt *ncrefct)
573: {
574:   PetscInt refct = ((PetscObject)dm)->refct;

576:   PetscFunctionBegin;
577:   *ncrefct = 0;
578:   if (dm->coarseMesh && dm->coarseMesh->fineMesh == dm) {
579:     refct--;
580:     if (recurseCoarse) {
581:       PetscInt coarseCount;

583:       PetscCall(DMCountNonCyclicReferences_Internal(dm->coarseMesh, PETSC_TRUE, PETSC_FALSE, &coarseCount));
584:       refct += coarseCount;
585:     }
586:   }
587:   if (dm->fineMesh && dm->fineMesh->coarseMesh == dm) {
588:     refct--;
589:     if (recurseFine) {
590:       PetscInt fineCount;

592:       PetscCall(DMCountNonCyclicReferences_Internal(dm->fineMesh, PETSC_FALSE, PETSC_TRUE, &fineCount));
593:       refct += fineCount;
594:     }
595:   }
596:   *ncrefct = refct;
597:   PetscFunctionReturn(PETSC_SUCCESS);
598: }

600: /* Generic wrapper for DMCountNonCyclicReferences_Internal() */
601: PetscErrorCode DMCountNonCyclicReferences(PetscObject dm, PetscInt *ncrefct)
602: {
603:   PetscFunctionBegin;
604:   PetscCall(DMCountNonCyclicReferences_Internal((DM)dm, PETSC_TRUE, PETSC_TRUE, ncrefct));
605:   PetscFunctionReturn(PETSC_SUCCESS);
606: }

608: PetscErrorCode DMDestroyLabelLinkList_Internal(DM dm)
609: {
610:   DMLabelLink next = dm->labels;

612:   PetscFunctionBegin;
613:   /* destroy the labels */
614:   while (next) {
615:     DMLabelLink tmp = next->next;

617:     if (next->label == dm->depthLabel) dm->depthLabel = NULL;
618:     if (next->label == dm->celltypeLabel) dm->celltypeLabel = NULL;
619:     PetscCall(DMLabelDestroy(&next->label));
620:     PetscCall(PetscFree(next));
621:     next = tmp;
622:   }
623:   dm->labels = NULL;
624:   PetscFunctionReturn(PETSC_SUCCESS);
625: }

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

638: /*@
639:   DMDestroy - Destroys a `DM`.

641:   Collective

643:   Input Parameter:
644: . dm - the `DM` object to destroy

646:   Level: developer

648: .seealso: [](ch_dmbase), `DM`, `DMCreate()`, `DMType`, `DMSetType()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`
649: @*/
650: PetscErrorCode DMDestroy(DM *dm)
651: {
652:   PetscInt cnt;

654:   PetscFunctionBegin;
655:   if (!*dm) PetscFunctionReturn(PETSC_SUCCESS);

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

668:   PetscCall(DMClearGlobalVectors(*dm));
669:   PetscCall(DMClearLocalVectors(*dm));
670:   PetscCall(DMClearNamedGlobalVectors(*dm));
671:   PetscCall(DMClearNamedLocalVectors(*dm));

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

735:       next = b->next;
736:       PetscCall(PetscFree(b));
737:     }
738:   }

740:   PetscCall(PetscObjectDestroy(&(*dm)->dmksp));
741:   PetscCall(PetscObjectDestroy(&(*dm)->dmsnes));
742:   PetscCall(PetscObjectDestroy(&(*dm)->dmts));

744:   if ((*dm)->ctx && (*dm)->ctxdestroy) PetscCall((*(*dm)->ctxdestroy)(&(*dm)->ctx));
745:   PetscCall(MatFDColoringDestroy(&(*dm)->fd));
746:   PetscCall(ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap));
747:   PetscCall(PetscFree((*dm)->vectype));
748:   PetscCall(PetscFree((*dm)->mattype));

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

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

781:   PetscCall(DMClearDS(*dm));
782:   PetscCall(DMDestroy(&(*dm)->dmBC));
783:   /* if memory was published with SAWs then destroy it */
784:   PetscCall(PetscObjectSAWsViewOff((PetscObject)*dm));

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

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

801:   Collective

803:   Input Parameter:
804: . dm - the `DM` object to setup

806:   Level: intermediate

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

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

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

826:   Collective

828:   Input Parameter:
829: . dm - the `DM` object to set options for

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

884:   Level: intermediate

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

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

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

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

924:   Collective

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

931:   Level: intermediate

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

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

946: /*@
947:   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
948:   save the `DM` in a binary file to be loaded later or create a visualization of the `DM`

950:   Collective

952:   Input Parameters:
953: + dm - the `DM` object to view
954: - v  - the viewer

956:   Level: beginner

958:   Note:
959:   Using `PETSCVIEWERHDF5` type with `PETSC_VIEWER_HDF5_PETSC` as the `PetscViewerFormat` one can save multiple `DMPLEX`
960:   meshes in a single HDF5 file. This in turn requires one to name the `DMPLEX` object with `PetscObjectSetName()`
961:   before saving it with `DMView()` and before loading it with `DMLoad()` for identification of the mesh object.

963: .seealso: [](ch_dmbase), `DM`, `PetscViewer`, `PetscViewerFormat`, `PetscViewerSetFormat()`, `DMDestroy()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMLoad()`, `PetscObjectSetName()`
964: @*/
965: PetscErrorCode DMView(DM dm, PetscViewer v)
966: {
967:   PetscBool         isbinary;
968:   PetscMPIInt       size;
969:   PetscViewerFormat format;

971:   PetscFunctionBegin;
973:   if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)dm), &v));
975:   /* Ideally, we would like to have this test on.
976:      However, it currently breaks socket viz via GLVis.
977:      During DMView(parallel_mesh,glvis_viewer), each
978:      process opens a sequential ASCII socket to visualize
979:      the local mesh, and PetscObjectView(dm,local_socket)
980:      is internally called inside VecView_GLVis, incurring
981:      in an error here */
982:   /* PetscCheckSameComm(dm,1,v,2); */
983:   PetscCall(PetscViewerCheckWritable(v));

985:   PetscCall(PetscViewerGetFormat(v, &format));
986:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
987:   if (size == 1 && format == PETSC_VIEWER_LOAD_BALANCE) PetscFunctionReturn(PETSC_SUCCESS);
988:   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)dm, v));
989:   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERBINARY, &isbinary));
990:   if (isbinary) {
991:     PetscInt classid = DM_FILE_CLASSID;
992:     char     type[256];

994:     PetscCall(PetscViewerBinaryWrite(v, &classid, 1, PETSC_INT));
995:     PetscCall(PetscStrncpy(type, ((PetscObject)dm)->type_name, sizeof(type)));
996:     PetscCall(PetscViewerBinaryWrite(v, type, 256, PETSC_CHAR));
997:   }
998:   PetscTryTypeMethod(dm, view, v);
999:   PetscFunctionReturn(PETSC_SUCCESS);
1000: }

1002: /*@
1003:   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,
1004:   that is it has no ghost locations.

1006:   Collective

1008:   Input Parameter:
1009: . dm - the `DM` object

1011:   Output Parameter:
1012: . vec - the global vector

1014:   Level: beginner

1016: .seealso: [](ch_dmbase), `DM`, `Vec`, `DMCreateLocalVector()`, `DMGetGlobalVector()`, `DMDestroy()`, `DMView()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`,
1017:          `DMGlobalToLocalBegin()`, `DMGlobalToLocalEnd()`
1018: @*/
1019: PetscErrorCode DMCreateGlobalVector(DM dm, Vec *vec)
1020: {
1021:   PetscFunctionBegin;
1023:   PetscAssertPointer(vec, 2);
1024:   PetscUseTypeMethod(dm, createglobalvector, vec);
1025:   if (PetscDefined(USE_DEBUG)) {
1026:     DM vdm;

1028:     PetscCall(VecGetDM(*vec, &vdm));
1029:     PetscCheck(vdm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DM type '%s' did not attach the DM to the vector", ((PetscObject)dm)->type_name);
1030:   }
1031:   PetscFunctionReturn(PETSC_SUCCESS);
1032: }

1034: /*@
1035:   DMCreateLocalVector - Creates a local vector from a `DM` object.

1037:   Not Collective

1039:   Input Parameter:
1040: . dm - the `DM` object

1042:   Output Parameter:
1043: . vec - the local vector

1045:   Level: beginner

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

1050: .seealso: [](ch_dmbase), `DM`, `Vec`, `DMCreateGlobalVector()`, `DMGetLocalVector()`, `DMDestroy()`, `DMView()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`
1051:          `DMGlobalToLocalBegin()`, `DMGlobalToLocalEnd()`
1052: @*/
1053: PetscErrorCode DMCreateLocalVector(DM dm, Vec *vec)
1054: {
1055:   PetscFunctionBegin;
1057:   PetscAssertPointer(vec, 2);
1058:   PetscUseTypeMethod(dm, createlocalvector, vec);
1059:   if (PetscDefined(USE_DEBUG)) {
1060:     DM vdm;

1062:     PetscCall(VecGetDM(*vec, &vdm));
1063:     PetscCheck(vdm, PETSC_COMM_SELF, PETSC_ERR_LIB, "DM type '%s' did not attach the DM to the vector", ((PetscObject)dm)->type_name);
1064:   }
1065:   PetscFunctionReturn(PETSC_SUCCESS);
1066: }

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

1071:   Collective

1073:   Input Parameter:
1074: . dm - the `DM` that provides the mapping

1076:   Output Parameter:
1077: . ltog - the mapping

1079:   Level: advanced

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

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

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

1089: .seealso: [](ch_dmbase), `DM`, `DMCreateLocalVector()`, `DMCreateGlobalVector()`, `VecSetLocalToGlobalMapping()`, `MatSetLocalToGlobalMapping()`,
1090:           `DMCreateMatrix()`
1091: @*/
1092: PetscErrorCode DMGetLocalToGlobalMapping(DM dm, ISLocalToGlobalMapping *ltog)
1093: {
1094:   PetscInt bs = -1, bsLocal[2], bsMinMax[2];

1096:   PetscFunctionBegin;
1098:   PetscAssertPointer(ltog, 2);
1099:   if (!dm->ltogmap) {
1100:     PetscSection section, sectionGlobal;

1102:     PetscCall(DMGetLocalSection(dm, &section));
1103:     if (section) {
1104:       const PetscInt *cdofs;
1105:       PetscInt       *ltog;
1106:       PetscInt        pStart, pEnd, n, p, k, l;

1108:       PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
1109:       PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
1110:       PetscCall(PetscSectionGetStorageSize(section, &n));
1111:       PetscCall(PetscMalloc1(n, &ltog)); /* We want the local+overlap size */
1112:       for (p = pStart, l = 0; p < pEnd; ++p) {
1113:         PetscInt bdof, cdof, dof, off, c, cind;

1115:         /* Should probably use constrained dofs */
1116:         PetscCall(PetscSectionGetDof(section, p, &dof));
1117:         PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
1118:         PetscCall(PetscSectionGetConstraintIndices(section, p, &cdofs));
1119:         PetscCall(PetscSectionGetOffset(sectionGlobal, p, &off));
1120:         /* If you have dofs, and constraints, and they are unequal, we set the blocksize to 1 */
1121:         bdof = cdof && (dof - cdof) ? 1 : dof;
1122:         if (dof) bs = bs < 0 ? bdof : PetscGCD(bs, bdof);

1124:         for (c = 0, cind = 0; c < dof; ++c, ++l) {
1125:           if (cind < cdof && c == cdofs[cind]) {
1126:             ltog[l] = off < 0 ? off - c : -(off + c + 1);
1127:             cind++;
1128:           } else {
1129:             ltog[l] = (off < 0 ? -(off + 1) : off) + c - cind;
1130:           }
1131:         }
1132:       }
1133:       /* Must have same blocksize on all procs (some might have no points) */
1134:       bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs;
1135:       bsLocal[1] = bs;
1136:       PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm), bsLocal, bsMinMax));
1137:       if (bsMinMax[0] != bsMinMax[1]) {
1138:         bs = 1;
1139:       } else {
1140:         bs = bsMinMax[0];
1141:       }
1142:       bs = bs < 0 ? 1 : bs;
1143:       /* Must reduce indices by blocksize */
1144:       if (bs > 1) {
1145:         for (l = 0, k = 0; l < n; l += bs, ++k) {
1146:           // Integer division of negative values truncates toward zero(!), not toward negative infinity
1147:           ltog[k] = ltog[l] >= 0 ? ltog[l] / bs : -(-(ltog[l] + 1) / bs + 1);
1148:         }
1149:         n /= bs;
1150:       }
1151:       PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), bs, n, ltog, PETSC_OWN_POINTER, &dm->ltogmap));
1152:     } else PetscUseTypeMethod(dm, getlocaltoglobalmapping);
1153:   }
1154:   *ltog = dm->ltogmap;
1155:   PetscFunctionReturn(PETSC_SUCCESS);
1156: }

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

1161:   Not Collective

1163:   Input Parameter:
1164: . dm - the `DM` with block structure

1166:   Output Parameter:
1167: . bs - the block size, 1 implies no exploitable block structure

1169:   Level: intermediate

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

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

1177: .seealso: [](ch_dmbase), `DM`, `ISCreateBlock()`, `VecSetBlockSize()`, `MatSetBlockSize()`, `DMGetLocalToGlobalMapping()`
1178: @*/
1179: PetscErrorCode DMGetBlockSize(DM dm, PetscInt *bs)
1180: {
1181:   PetscFunctionBegin;
1183:   PetscAssertPointer(bs, 2);
1184:   PetscCheck(dm->bs >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DM does not have enough information to provide a block size yet");
1185:   *bs = dm->bs;
1186:   PetscFunctionReturn(PETSC_SUCCESS);
1187: }

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

1193:   Collective

1195:   Input Parameters:
1196: + dmc - the `DM` object
1197: - dmf - the second, finer `DM` object

1199:   Output Parameters:
1200: + mat - the interpolation
1201: - vec - the scaling (optional, pass `NULL` if not needed), see `DMCreateInterpolationScale()`

1203:   Level: developer

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

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

1212: .seealso: [](ch_dmbase), `DM`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMRefine()`, `DMCoarsen()`, `DMCreateRestriction()`, `DMCreateInterpolationScale()`
1213: @*/
1214: PetscErrorCode DMCreateInterpolation(DM dmc, DM dmf, Mat *mat, Vec *vec)
1215: {
1216:   PetscFunctionBegin;
1219:   PetscAssertPointer(mat, 3);
1220:   PetscCall(PetscLogEventBegin(DM_CreateInterpolation, dmc, dmf, 0, 0));
1221:   PetscUseTypeMethod(dmc, createinterpolation, dmf, mat, vec);
1222:   PetscCall(PetscLogEventEnd(DM_CreateInterpolation, dmc, dmf, 0, 0));
1223:   PetscFunctionReturn(PETSC_SUCCESS);
1224: }

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

1230:   Input Parameters:
1231: + dac - `DM` that defines a coarse mesh
1232: . daf - `DM` that defines a fine mesh
1233: - mat - the restriction (or interpolation operator) from fine to coarse

1235:   Output Parameter:
1236: . scale - the scaled vector

1238:   Level: advanced

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

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

1248: .seealso: [](ch_dmbase), `DM`, `MatRestrict()`, `MatInterpolate()`, `DMCreateInterpolation()`, `DMCreateRestriction()`, `DMCreateGlobalVector()`
1249: @*/
1250: PetscErrorCode DMCreateInterpolationScale(DM dac, DM daf, Mat mat, Vec *scale)
1251: {
1252:   Vec         fine;
1253:   PetscScalar one = 1.0;
1254: #if defined(PETSC_HAVE_CUDA)
1255:   PetscBool bindingpropagates, isbound;
1256: #endif

1258:   PetscFunctionBegin;
1259:   PetscCall(DMCreateGlobalVector(daf, &fine));
1260:   PetscCall(DMCreateGlobalVector(dac, scale));
1261:   PetscCall(VecSet(fine, one));
1262: #if defined(PETSC_HAVE_CUDA)
1263:   /* If the 'fine' Vec is bound to the CPU, it makes sense to bind 'mat' as well.
1264:    * Note that we only do this for the CUDA case, right now, but if we add support for MatMultTranspose() via ViennaCL,
1265:    * we'll need to do it for that case, too.*/
1266:   PetscCall(VecGetBindingPropagates(fine, &bindingpropagates));
1267:   if (bindingpropagates) {
1268:     PetscCall(MatSetBindingPropagates(mat, PETSC_TRUE));
1269:     PetscCall(VecBoundToCPU(fine, &isbound));
1270:     PetscCall(MatBindToCPU(mat, isbound));
1271:   }
1272: #endif
1273:   PetscCall(MatRestrict(mat, fine, *scale));
1274:   PetscCall(VecDestroy(&fine));
1275:   PetscCall(VecReciprocal(*scale));
1276:   PetscFunctionReturn(PETSC_SUCCESS);
1277: }

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

1283:   Collective

1285:   Input Parameters:
1286: + dmc - the `DM` object
1287: - dmf - the second, finer `DM` object

1289:   Output Parameter:
1290: . mat - the restriction

1292:   Level: developer

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

1298: .seealso: [](ch_dmbase), `DM`, `DMRestrict()`, `DMInterpolate()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMRefine()`, `DMCoarsen()`, `DMCreateInterpolation()`
1299: @*/
1300: PetscErrorCode DMCreateRestriction(DM dmc, DM dmf, Mat *mat)
1301: {
1302:   PetscFunctionBegin;
1305:   PetscAssertPointer(mat, 3);
1306:   PetscCall(PetscLogEventBegin(DM_CreateRestriction, dmc, dmf, 0, 0));
1307:   PetscUseTypeMethod(dmc, createrestriction, dmf, mat);
1308:   PetscCall(PetscLogEventEnd(DM_CreateRestriction, dmc, dmf, 0, 0));
1309:   PetscFunctionReturn(PETSC_SUCCESS);
1310: }

1312: /*@
1313:   DMCreateInjection - Gets injection matrix between two `DM` objects.

1315:   Collective

1317:   Input Parameters:
1318: + dac - the `DM` object
1319: - daf - the second, finer `DM` object

1321:   Output Parameter:
1322: . mat - the injection

1324:   Level: developer

1326:   Notes:
1327:   This is an operator that applied to a vector obtained with `DMCreateGlobalVector()` on the
1328:   fine grid maps the values to a vector on the vector on the coarse `DM` by simply selecting
1329:   the values on the coarse grid points. This compares to the operator obtained by
1330:   `DMCreateRestriction()` or the transpose of the operator obtained by
1331:   `DMCreateInterpolation()` that uses a "local weighted average" of the values around the
1332:   coarse grid point as the coarse grid value.

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

1337: .seealso: [](ch_dmbase), `DM`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMCreateInterpolation()`,
1338:           `DMCreateRestriction()`, `MatRestrict()`, `MatInterpolate()`
1339: @*/
1340: PetscErrorCode DMCreateInjection(DM dac, DM daf, Mat *mat)
1341: {
1342:   PetscFunctionBegin;
1345:   PetscAssertPointer(mat, 3);
1346:   PetscCall(PetscLogEventBegin(DM_CreateInjection, dac, daf, 0, 0));
1347:   PetscUseTypeMethod(dac, createinjection, daf, mat);
1348:   PetscCall(PetscLogEventEnd(DM_CreateInjection, dac, daf, 0, 0));
1349:   PetscFunctionReturn(PETSC_SUCCESS);
1350: }

1352: /*@
1353:   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
1354:   a Galerkin finite element model on the `DM`

1356:   Collective

1358:   Input Parameters:
1359: + dmc - the target `DM` object
1360: - dmf - the source `DM` object

1362:   Output Parameter:
1363: . mat - the mass matrix

1365:   Level: developer

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

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

1372: .seealso: [](ch_dmbase), `DM`, `DMCreateMassMatrixLumped()`, `DMCreateMatrix()`, `DMRefine()`, `DMCoarsen()`, `DMCreateRestriction()`, `DMCreateInterpolation()`, `DMCreateInjection()`
1373: @*/
1374: PetscErrorCode DMCreateMassMatrix(DM dmc, DM dmf, Mat *mat)
1375: {
1376:   PetscFunctionBegin;
1379:   PetscAssertPointer(mat, 3);
1380:   PetscCall(PetscLogEventBegin(DM_CreateMassMatrix, 0, 0, 0, 0));
1381:   PetscUseTypeMethod(dmc, createmassmatrix, dmf, mat);
1382:   PetscCall(PetscLogEventEnd(DM_CreateMassMatrix, 0, 0, 0, 0));
1383:   PetscFunctionReturn(PETSC_SUCCESS);
1384: }

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

1389:   Collective

1391:   Input Parameter:
1392: . dm - the `DM` object

1394:   Output Parameter:
1395: . lm - the lumped mass matrix, which is a diagonal matrix, represented as a vector

1397:   Level: developer

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

1402: .seealso: [](ch_dmbase), `DM`, `DMCreateMassMatrix()`, `DMCreateMatrix()`, `DMRefine()`, `DMCoarsen()`, `DMCreateRestriction()`, `DMCreateInterpolation()`, `DMCreateInjection()`
1403: @*/
1404: PetscErrorCode DMCreateMassMatrixLumped(DM dm, Vec *lm)
1405: {
1406:   PetscFunctionBegin;
1408:   PetscAssertPointer(lm, 2);
1409:   PetscUseTypeMethod(dm, createmassmatrixlumped, lm);
1410:   PetscFunctionReturn(PETSC_SUCCESS);
1411: }

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

1417:   Collective

1419:   Input Parameters:
1420: + dm    - the `DM` object
1421: - ctype - `IS_COLORING_LOCAL` or `IS_COLORING_GLOBAL`

1423:   Output Parameter:
1424: . coloring - the coloring

1426:   Level: developer

1428:   Notes:
1429:   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
1430:   matrix comes from (what this function provides). In general using the mesh produces a more optimal coloring (fewer colors).

1432:   This produces a coloring with the distance of 2, see `MatSetColoringDistance()` which can be used for efficiently computing Jacobians with `MatFDColoringCreate()`
1433:   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,
1434:   otherwise an error will be generated.

1436: .seealso: [](ch_dmbase), `DM`, `ISColoring`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMSetMatType()`, `MatColoring`, `MatFDColoringCreate()`
1437: @*/
1438: PetscErrorCode DMCreateColoring(DM dm, ISColoringType ctype, ISColoring *coloring)
1439: {
1440:   PetscFunctionBegin;
1442:   PetscAssertPointer(coloring, 3);
1443:   PetscUseTypeMethod(dm, getcoloring, ctype, coloring);
1444:   PetscFunctionReturn(PETSC_SUCCESS);
1445: }

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

1450:   Collective

1452:   Input Parameter:
1453: . dm - the `DM` object

1455:   Output Parameter:
1456: . mat - the empty Jacobian

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

1461:   Level: beginner

1463:   Notes:
1464:   This properly preallocates the number of nonzeros in the sparse matrix so you
1465:   do not need to do it yourself.

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

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

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

1476: .seealso: [](ch_dmbase), `DM`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMSetMatType()`, `DMCreateMassMatrix()`
1477: @*/
1478: PetscErrorCode DMCreateMatrix(DM dm, Mat *mat)
1479: {
1480:   PetscFunctionBegin;
1482:   PetscAssertPointer(mat, 2);
1483:   PetscCall(MatInitializePackage());
1484:   PetscCall(PetscLogEventBegin(DM_CreateMatrix, 0, 0, 0, 0));
1485:   PetscUseTypeMethod(dm, creatematrix, mat);
1486:   if (PetscDefined(USE_DEBUG)) {
1487:     DM mdm;

1489:     PetscCall(MatGetDM(*mat, &mdm));
1490:     PetscCheck(mdm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DM type '%s' did not attach the DM to the matrix", ((PetscObject)dm)->type_name);
1491:   }
1492:   /* Handle nullspace and near nullspace */
1493:   if (dm->Nf) {
1494:     MatNullSpace nullSpace;
1495:     PetscInt     Nf, f;

1497:     PetscCall(DMGetNumFields(dm, &Nf));
1498:     for (f = 0; f < Nf; ++f) {
1499:       if (dm->nullspaceConstructors[f]) {
1500:         PetscCall((*dm->nullspaceConstructors[f])(dm, f, f, &nullSpace));
1501:         PetscCall(MatSetNullSpace(*mat, nullSpace));
1502:         PetscCall(MatNullSpaceDestroy(&nullSpace));
1503:         break;
1504:       }
1505:     }
1506:     for (f = 0; f < Nf; ++f) {
1507:       if (dm->nearnullspaceConstructors[f]) {
1508:         PetscCall((*dm->nearnullspaceConstructors[f])(dm, f, f, &nullSpace));
1509:         PetscCall(MatSetNearNullSpace(*mat, nullSpace));
1510:         PetscCall(MatNullSpaceDestroy(&nullSpace));
1511:       }
1512:     }
1513:   }
1514:   PetscCall(PetscLogEventEnd(DM_CreateMatrix, 0, 0, 0, 0));
1515:   PetscFunctionReturn(PETSC_SUCCESS);
1516: }

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

1523:   Logically Collective

1525:   Input Parameters:
1526: + dm   - the `DM`
1527: - skip - `PETSC_TRUE` to skip preallocation

1529:   Level: developer

1531:   Note:
1532:   This is most useful to reduce initialization costs when `MatSetPreallocationCOO()` and
1533:   `MatSetValuesCOO()` will be used.

1535: .seealso: [](ch_dmbase), `DM`, `DMCreateMatrix()`, `DMSetMatrixStructureOnly()`, `DMSetMatrixPreallocateOnly()`
1536: @*/
1537: PetscErrorCode DMSetMatrixPreallocateSkip(DM dm, PetscBool skip)
1538: {
1539:   PetscFunctionBegin;
1541:   dm->prealloc_skip = skip;
1542:   PetscFunctionReturn(PETSC_SUCCESS);
1543: }

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

1549:   Logically Collective

1551:   Input Parameters:
1552: + dm   - the `DM`
1553: - only - `PETSC_TRUE` if only want preallocation

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

1558:   Level: developer

1560: .seealso: [](ch_dmbase), `DM`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMSetMatrixStructureOnly()`, `DMSetMatrixPreallocateSkip()`
1561: @*/
1562: PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only)
1563: {
1564:   PetscFunctionBegin;
1566:   dm->prealloc_only = only;
1567:   PetscFunctionReturn(PETSC_SUCCESS);
1568: }

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

1574:   Logically Collective

1576:   Input Parameters:
1577: + dm   - the `DM`
1578: - only - `PETSC_TRUE` if you only want matrix structure

1580:   Level: developer

1582: .seealso: [](ch_dmbase), `DM`, `DMCreateMatrix()`, `DMSetMatrixPreallocateOnly()`, `DMSetMatrixPreallocateSkip()`
1583: @*/
1584: PetscErrorCode DMSetMatrixStructureOnly(DM dm, PetscBool only)
1585: {
1586:   PetscFunctionBegin;
1588:   dm->structure_only = only;
1589:   PetscFunctionReturn(PETSC_SUCCESS);
1590: }

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

1595:   Logically Collective

1597:   Input Parameters:
1598: + dm    - the `DM`
1599: - btype - block by topological point or field node

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

1604:   Level: advanced

1606: .seealso: [](ch_dmbase), `DM`, `DMCreateMatrix()`, `MatSetVariableBlockSizes()`
1607: @*/
1608: PetscErrorCode DMSetBlockingType(DM dm, DMBlockingType btype)
1609: {
1610:   PetscFunctionBegin;
1612:   dm->blocking_type = btype;
1613:   PetscFunctionReturn(PETSC_SUCCESS);
1614: }

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

1619:   Not Collective

1621:   Input Parameter:
1622: . dm - the `DM`

1624:   Output Parameter:
1625: . btype - block by topological point or field node

1627:   Level: advanced

1629: .seealso: [](ch_dmbase), `DM`, `DMCreateMatrix()`, `MatSetVariableBlockSizes()`
1630: @*/
1631: PetscErrorCode DMGetBlockingType(DM dm, DMBlockingType *btype)
1632: {
1633:   PetscFunctionBegin;
1635:   PetscAssertPointer(btype, 2);
1636:   *btype = dm->blocking_type;
1637:   PetscFunctionReturn(PETSC_SUCCESS);
1638: }

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

1643:   Not Collective

1645:   Input Parameters:
1646: + dm    - the `DM` object
1647: . count - The minimum size
1648: - dtype - MPI data type, often `MPIU_REAL`, `MPIU_SCALAR`, or `MPIU_INT`)

1650:   Output Parameter:
1651: . mem - the work array

1653:   Level: developer

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

1658:   The array may contain nonzero values

1660: .seealso: [](ch_dmbase), `DM`, `DMDestroy()`, `DMCreate()`, `DMRestoreWorkArray()`, `PetscMalloc()`
1661: @*/
1662: PetscErrorCode DMGetWorkArray(DM dm, PetscInt count, MPI_Datatype dtype, void *mem)
1663: {
1664:   DMWorkLink  link;
1665:   PetscMPIInt dsize;

1667:   PetscFunctionBegin;
1669:   PetscAssertPointer(mem, 4);
1670:   if (!count) {
1671:     *(void **)mem = NULL;
1672:     PetscFunctionReturn(PETSC_SUCCESS);
1673:   }
1674:   if (dm->workin) {
1675:     link       = dm->workin;
1676:     dm->workin = dm->workin->next;
1677:   } else {
1678:     PetscCall(PetscNew(&link));
1679:   }
1680:   /* Avoid MPI_Type_size for most used datatypes
1681:      Get size directly */
1682:   if (dtype == MPIU_INT) dsize = sizeof(PetscInt);
1683:   else if (dtype == MPIU_REAL) dsize = sizeof(PetscReal);
1684: #if defined(PETSC_USE_64BIT_INDICES)
1685:   else if (dtype == MPI_INT) dsize = sizeof(int);
1686: #endif
1687: #if defined(PETSC_USE_COMPLEX)
1688:   else if (dtype == MPIU_SCALAR) dsize = sizeof(PetscScalar);
1689: #endif
1690:   else PetscCallMPI(MPI_Type_size(dtype, &dsize));

1692:   if (((size_t)dsize * count) > link->bytes) {
1693:     PetscCall(PetscFree(link->mem));
1694:     PetscCall(PetscMalloc(dsize * count, &link->mem));
1695:     link->bytes = dsize * count;
1696:   }
1697:   link->next  = dm->workout;
1698:   dm->workout = link;
1699: #if defined(__MEMCHECK_H) && (defined(PLAT_amd64_linux) || defined(PLAT_x86_linux) || defined(PLAT_amd64_darwin))
1700:   VALGRIND_MAKE_MEM_NOACCESS((char *)link->mem + (size_t)dsize * count, link->bytes - (size_t)dsize * count);
1701:   VALGRIND_MAKE_MEM_UNDEFINED(link->mem, (size_t)dsize * count);
1702: #endif
1703:   *(void **)mem = link->mem;
1704:   PetscFunctionReturn(PETSC_SUCCESS);
1705: }

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

1710:   Not Collective

1712:   Input Parameters:
1713: + dm    - the `DM` object
1714: . count - The minimum size
1715: - dtype - MPI data type, often `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_INT`

1717:   Output Parameter:
1718: . mem - the work array

1720:   Level: developer

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

1725: .seealso: [](ch_dmbase), `DM`, `DMDestroy()`, `DMCreate()`, `DMGetWorkArray()`
1726: @*/
1727: PetscErrorCode DMRestoreWorkArray(DM dm, PetscInt count, MPI_Datatype dtype, void *mem)
1728: {
1729:   DMWorkLink *p, link;

1731:   PetscFunctionBegin;
1733:   PetscAssertPointer(mem, 4);
1734:   (void)count;
1735:   (void)dtype;
1736:   if (!*(void **)mem) PetscFunctionReturn(PETSC_SUCCESS);
1737:   for (p = &dm->workout; (link = *p); p = &link->next) {
1738:     if (link->mem == *(void **)mem) {
1739:       *p            = link->next;
1740:       link->next    = dm->workin;
1741:       dm->workin    = link;
1742:       *(void **)mem = NULL;
1743:       PetscFunctionReturn(PETSC_SUCCESS);
1744:     }
1745:   }
1746:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Array was not checked out");
1747: }

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

1753:   Logically Collective; No Fortran Support

1755:   Input Parameters:
1756: + dm     - The `DM`
1757: . field  - The field number for the nullspace
1758: - nullsp - A callback to create the nullspace

1760:   Calling sequence of `nullsp`:
1761: + dm        - The present `DM`
1762: . origField - The field number given above, in the original `DM`
1763: . field     - The field number in dm
1764: - nullSpace - The nullspace for the given field

1766:   Level: intermediate

1768: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMGetNullSpaceConstructor()`, `DMSetNearNullSpaceConstructor()`, `DMGetNearNullSpaceConstructor()`, `DMCreateSubDM()`, `DMCreateSuperDM()`
1769: @*/
1770: PetscErrorCode DMSetNullSpaceConstructor(DM dm, PetscInt field, PetscErrorCode (*nullsp)(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullSpace))
1771: {
1772:   PetscFunctionBegin;
1774:   PetscCheck(field < 10, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %" PetscInt_FMT " >= 10 fields", field);
1775:   dm->nullspaceConstructors[field] = nullsp;
1776:   PetscFunctionReturn(PETSC_SUCCESS);
1777: }

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

1782:   Not Collective; No Fortran Support

1784:   Input Parameters:
1785: + dm    - The `DM`
1786: - field - The field number for the nullspace

1788:   Output Parameter:
1789: . nullsp - A callback to create the nullspace

1791:   Calling sequence of `nullsp`:
1792: + dm        - The present DM
1793: . origField - The field number given above, in the original DM
1794: . field     - The field number in dm
1795: - nullSpace - The nullspace for the given field

1797:   Level: intermediate

1799: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMGetField()`, `DMSetNullSpaceConstructor()`, `DMSetNearNullSpaceConstructor()`, `DMGetNearNullSpaceConstructor()`, `DMCreateSubDM()`, `DMCreateSuperDM()`
1800: @*/
1801: PetscErrorCode DMGetNullSpaceConstructor(DM dm, PetscInt field, PetscErrorCode (**nullsp)(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullSpace))
1802: {
1803:   PetscFunctionBegin;
1805:   PetscAssertPointer(nullsp, 3);
1806:   PetscCheck(field < 10, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %" PetscInt_FMT " >= 10 fields", field);
1807:   *nullsp = dm->nullspaceConstructors[field];
1808:   PetscFunctionReturn(PETSC_SUCCESS);
1809: }

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

1814:   Logically Collective; No Fortran Support

1816:   Input Parameters:
1817: + dm     - The `DM`
1818: . field  - The field number for the nullspace
1819: - nullsp - A callback to create the near-nullspace

1821:   Calling sequence of `nullsp`:
1822: + dm        - The present `DM`
1823: . origField - The field number given above, in the original `DM`
1824: . field     - The field number in dm
1825: - nullSpace - The nullspace for the given field

1827:   Level: intermediate

1829: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMGetNearNullSpaceConstructor()`, `DMSetNullSpaceConstructor()`, `DMGetNullSpaceConstructor()`, `DMCreateSubDM()`, `DMCreateSuperDM()`,
1830:           `MatNullSpace`
1831: @*/
1832: PetscErrorCode DMSetNearNullSpaceConstructor(DM dm, PetscInt field, PetscErrorCode (*nullsp)(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullSpace))
1833: {
1834:   PetscFunctionBegin;
1836:   PetscCheck(field < 10, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %" PetscInt_FMT " >= 10 fields", field);
1837:   dm->nearnullspaceConstructors[field] = nullsp;
1838:   PetscFunctionReturn(PETSC_SUCCESS);
1839: }

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

1844:   Not Collective; No Fortran Support

1846:   Input Parameters:
1847: + dm    - The `DM`
1848: - field - The field number for the nullspace

1850:   Output Parameter:
1851: . nullsp - A callback to create the near-nullspace

1853:   Calling sequence of `nullsp`:
1854: + dm        - The present `DM`
1855: . origField - The field number given above, in the original `DM`
1856: . field     - The field number in dm
1857: - nullSpace - The nullspace for the given field

1859:   Level: intermediate

1861: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMGetField()`, `DMSetNearNullSpaceConstructor()`, `DMSetNullSpaceConstructor()`, `DMGetNullSpaceConstructor()`, `DMCreateSubDM()`,
1862:           `MatNullSpace`, `DMCreateSuperDM()`
1863: @*/
1864: PetscErrorCode DMGetNearNullSpaceConstructor(DM dm, PetscInt field, PetscErrorCode (**nullsp)(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullSpace))
1865: {
1866:   PetscFunctionBegin;
1868:   PetscAssertPointer(nullsp, 3);
1869:   PetscCheck(field < 10, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %" PetscInt_FMT " >= 10 fields", field);
1870:   *nullsp = dm->nearnullspaceConstructors[field];
1871:   PetscFunctionReturn(PETSC_SUCCESS);
1872: }

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

1877:   Not Collective; No Fortran Support

1879:   Input Parameter:
1880: . dm - the `DM` object

1882:   Output Parameters:
1883: + numFields  - The number of fields (or `NULL` if not requested)
1884: . fieldNames - The number of each field (or `NULL` if not requested)
1885: - fields     - The global indices for each field (or `NULL` if not requested)

1887:   Level: intermediate

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

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

1898: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMGetField()`, `DMDestroy()`, `DMView()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`,
1899:           `DMCreateFieldDecomposition()`
1900: @*/
1901: PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
1902: {
1903:   PetscSection section, sectionGlobal;

1905:   PetscFunctionBegin;
1907:   if (numFields) {
1908:     PetscAssertPointer(numFields, 2);
1909:     *numFields = 0;
1910:   }
1911:   if (fieldNames) {
1912:     PetscAssertPointer(fieldNames, 3);
1913:     *fieldNames = NULL;
1914:   }
1915:   if (fields) {
1916:     PetscAssertPointer(fields, 4);
1917:     *fields = NULL;
1918:   }
1919:   PetscCall(DMGetLocalSection(dm, &section));
1920:   if (section) {
1921:     PetscInt *fieldSizes, *fieldNc, **fieldIndices;
1922:     PetscInt  nF, f, pStart, pEnd, p;

1924:     PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
1925:     PetscCall(PetscSectionGetNumFields(section, &nF));
1926:     PetscCall(PetscMalloc3(nF, &fieldSizes, nF, &fieldNc, nF, &fieldIndices));
1927:     PetscCall(PetscSectionGetChart(sectionGlobal, &pStart, &pEnd));
1928:     for (f = 0; f < nF; ++f) {
1929:       fieldSizes[f] = 0;
1930:       PetscCall(PetscSectionGetFieldComponents(section, f, &fieldNc[f]));
1931:     }
1932:     for (p = pStart; p < pEnd; ++p) {
1933:       PetscInt gdof;

1935:       PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof));
1936:       if (gdof > 0) {
1937:         for (f = 0; f < nF; ++f) {
1938:           PetscInt fdof, fcdof, fpdof;

1940:           PetscCall(PetscSectionGetFieldDof(section, p, f, &fdof));
1941:           PetscCall(PetscSectionGetFieldConstraintDof(section, p, f, &fcdof));
1942:           fpdof = fdof - fcdof;
1943:           if (fpdof && fpdof != fieldNc[f]) {
1944:             /* Layout does not admit a pointwise block size */
1945:             fieldNc[f] = 1;
1946:           }
1947:           fieldSizes[f] += fpdof;
1948:         }
1949:       }
1950:     }
1951:     for (f = 0; f < nF; ++f) {
1952:       PetscCall(PetscMalloc1(fieldSizes[f], &fieldIndices[f]));
1953:       fieldSizes[f] = 0;
1954:     }
1955:     for (p = pStart; p < pEnd; ++p) {
1956:       PetscInt gdof, goff;

1958:       PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof));
1959:       if (gdof > 0) {
1960:         PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
1961:         for (f = 0; f < nF; ++f) {
1962:           PetscInt fdof, fcdof, fc;

1964:           PetscCall(PetscSectionGetFieldDof(section, p, f, &fdof));
1965:           PetscCall(PetscSectionGetFieldConstraintDof(section, p, f, &fcdof));
1966:           for (fc = 0; fc < fdof - fcdof; ++fc, ++fieldSizes[f]) fieldIndices[f][fieldSizes[f]] = goff++;
1967:         }
1968:       }
1969:     }
1970:     if (numFields) *numFields = nF;
1971:     if (fieldNames) {
1972:       PetscCall(PetscMalloc1(nF, fieldNames));
1973:       for (f = 0; f < nF; ++f) {
1974:         const char *fieldName;

1976:         PetscCall(PetscSectionGetFieldName(section, f, &fieldName));
1977:         PetscCall(PetscStrallocpy(fieldName, (char **)&(*fieldNames)[f]));
1978:       }
1979:     }
1980:     if (fields) {
1981:       PetscCall(PetscMalloc1(nF, fields));
1982:       for (f = 0; f < nF; ++f) {
1983:         PetscInt bs, in[2], out[2];

1985:         PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), fieldSizes[f], fieldIndices[f], PETSC_OWN_POINTER, &(*fields)[f]));
1986:         in[0] = -fieldNc[f];
1987:         in[1] = fieldNc[f];
1988:         PetscCall(MPIU_Allreduce(in, out, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
1989:         bs = (-out[0] == out[1]) ? out[1] : 1;
1990:         PetscCall(ISSetBlockSize((*fields)[f], bs));
1991:       }
1992:     }
1993:     PetscCall(PetscFree3(fieldSizes, fieldNc, fieldIndices));
1994:   } else PetscTryTypeMethod(dm, createfieldis, numFields, fieldNames, fields);
1995:   PetscFunctionReturn(PETSC_SUCCESS);
1996: }

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

2002:   Not Collective; No Fortran Support

2004:   Input Parameter:
2005: . dm - the `DM` object

2007:   Output Parameters:
2008: + len      - The number of fields (or `NULL` if not requested)
2009: . namelist - The name for each field (or `NULL` if not requested)
2010: . islist   - The global indices for each field (or `NULL` if not requested)
2011: - dmlist   - The `DM`s for each field subproblem (or `NULL`, if not requested; if `NULL` is returned, no `DM`s are defined)

2013:   Level: intermediate

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

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

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

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

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

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

2063:     PetscCall(DMGetLocalSection(dm, &section));
2064:     if (section) PetscCall(PetscSectionGetNumFields(section, &numFields));
2065:     if (section && numFields && dm->ops->createsubdm) {
2066:       if (len) *len = numFields;
2067:       if (namelist) PetscCall(PetscMalloc1(numFields, namelist));
2068:       if (islist) PetscCall(PetscMalloc1(numFields, islist));
2069:       if (dmlist) PetscCall(PetscMalloc1(numFields, dmlist));
2070:       for (f = 0; f < numFields; ++f) {
2071:         const char *fieldName;

2073:         PetscCall(DMCreateSubDM(dm, 1, &f, islist ? &(*islist)[f] : NULL, dmlist ? &(*dmlist)[f] : NULL));
2074:         if (namelist) {
2075:           PetscCall(PetscSectionGetFieldName(section, f, &fieldName));
2076:           PetscCall(PetscStrallocpy(fieldName, (char **)&(*namelist)[f]));
2077:         }
2078:       }
2079:     } else {
2080:       PetscCall(DMCreateFieldIS(dm, len, namelist, islist));
2081:       /* By default there are no DMs associated with subproblems. */
2082:       if (dmlist) *dmlist = NULL;
2083:     }
2084:   } else PetscUseTypeMethod(dm, createfielddecomposition, len, namelist, islist, dmlist);
2085:   PetscFunctionReturn(PETSC_SUCCESS);
2086: }

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

2092:   Not collective

2094:   Input Parameters:
2095: + dm        - The `DM` object
2096: . numFields - The number of fields to select
2097: - fields    - The field numbers of the selected fields

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

2103:   Level: intermediate

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

2108: .seealso: [](ch_dmbase), `DM`, `DMCreateFieldIS()`, `DMCreateFieldDecomposition()`, `DMAddField()`, `DMCreateSuperDM()`, `IS`, `DMPlexSetMigrationSF()`, `DMDestroy()`, `DMView()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`
2109: @*/
2110: PetscErrorCode DMCreateSubDM(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
2111: {
2112:   PetscFunctionBegin;
2114:   PetscAssertPointer(fields, 3);
2115:   if (is) PetscAssertPointer(is, 4);
2116:   if (subdm) PetscAssertPointer(subdm, 5);
2117:   PetscUseTypeMethod(dm, createsubdm, numFields, fields, is, subdm);
2118:   PetscFunctionReturn(PETSC_SUCCESS);
2119: }

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

2124:   Not collective

2126:   Input Parameters:
2127: + dms - The `DM` objects
2128: - n   - The number of `DM`s

2130:   Output Parameters:
2131: + is      - The global indices for each of subproblem within the super `DM`, or NULL
2132: - superdm - The `DM` for the superproblem

2134:   Level: intermediate

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

2139: .seealso: [](ch_dmbase), `DM`, `DMCreateSubDM()`, `DMPlexSetMigrationSF()`, `DMDestroy()`, `DMView()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMCreateFieldIS()`, `DMCreateDomainDecomposition()`
2140: @*/
2141: PetscErrorCode DMCreateSuperDM(DM dms[], PetscInt n, IS *is[], DM *superdm)
2142: {
2143:   PetscInt i;

2145:   PetscFunctionBegin;
2146:   PetscAssertPointer(dms, 1);
2148:   if (is) PetscAssertPointer(is, 3);
2149:   PetscAssertPointer(superdm, 4);
2150:   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of DMs must be nonnegative: %" PetscInt_FMT, n);
2151:   if (n) {
2152:     DM dm = dms[0];
2153:     PetscCheck(dm->ops->createsuperdm, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No method createsuperdm for DM of type %s", ((PetscObject)dm)->type_name);
2154:     PetscCall((*dm->ops->createsuperdm)(dms, n, is, superdm));
2155:   }
2156:   PetscFunctionReturn(PETSC_SUCCESS);
2157: }

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

2163:   Not Collective

2165:   Input Parameter:
2166: . dm - the `DM` object

2168:   Output Parameters:
2169: + n           - The number of subproblems in the domain decomposition (or `NULL` if not requested)
2170: . namelist    - The name for each subdomain (or `NULL` if not requested)
2171: . innerislist - The global indices for each inner subdomain (or `NULL`, if not requested)
2172: . outerislist - The global indices for each outer subdomain (or `NULL`, if not requested)
2173: - dmlist      - The `DM`s for each subdomain subproblem (or `NULL`, if not requested; if `NULL` is returned, no `DM`s are defined)

2175:   Level: intermediate

2177:   Notes:
2178:   Each `IS` contains the global indices of the dofs of the corresponding subdomains with in the
2179:   dofs of the original `DM`. The inner subdomains conceptually define a nonoverlapping
2180:   covering, while outer subdomains can overlap.

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

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

2188:   Developer Notes:
2189:   The `dmlist` is for the inner subdomains or the outer subdomains or all subdomains?

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

2193: .seealso: [](ch_dmbase), `DM`, `DMCreateFieldDecomposition()`, `DMDestroy()`, `DMCreateDomainDecompositionScatters()`, `DMView()`, `DMCreateInterpolation()`,
2194:           `DMSubDomainHookAdd()`, `DMSubDomainHookRemove()`,`DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMRefine()`, `DMCoarsen()`
2195: @*/
2196: PetscErrorCode DMCreateDomainDecomposition(DM dm, PetscInt *n, char ***namelist, IS **innerislist, IS **outerislist, DM **dmlist)
2197: {
2198:   DMSubDomainHookLink link;
2199:   PetscInt            i, l;

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

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

2249:   Not Collective

2251:   Input Parameters:
2252: + dm     - the `DM` object
2253: . n      - the number of subdomains
2254: - subdms - the local subdomains

2256:   Output Parameters:
2257: + iscat - scatter from global vector to nonoverlapping global vector entries on subdomain
2258: . oscat - scatter from global vector to overlapping global vector entries on subdomain
2259: - gscat - scatter from global vector to local vector on subdomain (fills in ghosts)

2261:   Level: developer

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

2269:   Developer Note:
2270:   Can the subdms input be anything or are they exactly the `DM` obtained from
2271:   `DMCreateDomainDecomposition()`?

2273: .seealso: [](ch_dmbase), `DM`, `DMCreateDomainDecomposition()`, `DMDestroy()`, `DMView()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMCreateFieldIS()`
2274: @*/
2275: PetscErrorCode DMCreateDomainDecompositionScatters(DM dm, PetscInt n, DM *subdms, VecScatter *iscat[], VecScatter *oscat[], VecScatter *gscat[])
2276: {
2277:   PetscFunctionBegin;
2279:   PetscAssertPointer(subdms, 3);
2280:   PetscUseTypeMethod(dm, createddscatters, n, subdms, iscat, oscat, gscat);
2281:   PetscFunctionReturn(PETSC_SUCCESS);
2282: }

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

2287:   Collective

2289:   Input Parameters:
2290: + dm   - the `DM` object
2291: - comm - the communicator to contain the new `DM` object (or `MPI_COMM_NULL`)

2293:   Output Parameter:
2294: . dmf - the refined `DM`, or `NULL`

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

2299:   Level: developer

2301:   Note:
2302:   If no refinement was done, the return value is `NULL`

2304: .seealso: [](ch_dmbase), `DM`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateDomainDecomposition()`,
2305:           `DMRefineHookAdd()`, `DMRefineHookRemove()`
2306: @*/
2307: PetscErrorCode DMRefine(DM dm, MPI_Comm comm, DM *dmf)
2308: {
2309:   DMRefineHookLink link;

2311:   PetscFunctionBegin;
2313:   PetscCall(PetscLogEventBegin(DM_Refine, dm, 0, 0, 0));
2314:   PetscUseTypeMethod(dm, refine, comm, dmf);
2315:   if (*dmf) {
2316:     (*dmf)->ops->creatematrix = dm->ops->creatematrix;

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

2320:     (*dmf)->ctx       = dm->ctx;
2321:     (*dmf)->leveldown = dm->leveldown;
2322:     (*dmf)->levelup   = dm->levelup + 1;

2324:     PetscCall(DMSetMatType(*dmf, dm->mattype));
2325:     for (link = dm->refinehook; link; link = link->next) {
2326:       if (link->refinehook) PetscCall((*link->refinehook)(dm, *dmf, link->ctx));
2327:     }
2328:   }
2329:   PetscCall(PetscLogEventEnd(DM_Refine, dm, 0, 0, 0));
2330:   PetscFunctionReturn(PETSC_SUCCESS);
2331: }

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

2336:   Logically Collective; No Fortran Support

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

2344:   Calling sequence of `refinehook`:
2345: + coarse - coarse level `DM`
2346: . fine   - fine level `DM` to interpolate problem to
2347: - ctx    - optional user-defined function context

2349:   Calling sequence of `interphook`:
2350: + coarse - coarse level `DM`
2351: . interp - matrix interpolating a coarse-level solution to the finer grid
2352: . fine   - fine level `DM` to update
2353: - ctx    - optional user-defined function context

2355:   Level: advanced

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

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

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

2365: .seealso: [](ch_dmbase), `DM`, `DMCoarsenHookAdd()`, `DMInterpolate()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`
2366: @*/
2367: PetscErrorCode DMRefineHookAdd(DM coarse, PetscErrorCode (*refinehook)(DM coarse, DM fine, void *ctx), PetscErrorCode (*interphook)(DM coarse, Mat interp, DM fine, void *ctx), void *ctx)
2368: {
2369:   DMRefineHookLink link, *p;

2371:   PetscFunctionBegin;
2373:   for (p = &coarse->refinehook; *p; p = &(*p)->next) { /* Scan to the end of the current list of hooks */
2374:     if ((*p)->refinehook == refinehook && (*p)->interphook == interphook && (*p)->ctx == ctx) PetscFunctionReturn(PETSC_SUCCESS);
2375:   }
2376:   PetscCall(PetscNew(&link));
2377:   link->refinehook = refinehook;
2378:   link->interphook = interphook;
2379:   link->ctx        = ctx;
2380:   link->next       = NULL;
2381:   *p               = link;
2382:   PetscFunctionReturn(PETSC_SUCCESS);
2383: }

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

2389:   Logically Collective; No Fortran Support

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

2397:   Level: advanced

2399:   Note:
2400:   This function does nothing if the hook is not in the list.

2402: .seealso: [](ch_dmbase), `DM`, `DMRefineHookAdd()`, `DMCoarsenHookRemove()`, `DMInterpolate()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`
2403: @*/
2404: PetscErrorCode DMRefineHookRemove(DM coarse, PetscErrorCode (*refinehook)(DM, DM, void *), PetscErrorCode (*interphook)(DM, Mat, DM, void *), void *ctx)
2405: {
2406:   DMRefineHookLink link, *p;

2408:   PetscFunctionBegin;
2410:   for (p = &coarse->refinehook; *p; p = &(*p)->next) { /* Search the list of current hooks */
2411:     if ((*p)->refinehook == refinehook && (*p)->interphook == interphook && (*p)->ctx == ctx) {
2412:       link = *p;
2413:       *p   = link->next;
2414:       PetscCall(PetscFree(link));
2415:       break;
2416:     }
2417:   }
2418:   PetscFunctionReturn(PETSC_SUCCESS);
2419: }

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

2424:   Collective if any hooks are

2426:   Input Parameters:
2427: + coarse - coarser `DM` to use as a base
2428: . interp - interpolation matrix, apply using `MatInterpolate()`
2429: - fine   - finer `DM` to update

2431:   Level: developer

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

2437: .seealso: [](ch_dmbase), `DM`, `DMRefineHookAdd()`, `MatInterpolate()`
2438: @*/
2439: PetscErrorCode DMInterpolate(DM coarse, Mat interp, DM fine)
2440: {
2441:   DMRefineHookLink link;

2443:   PetscFunctionBegin;
2444:   for (link = fine->refinehook; link; link = link->next) {
2445:     if (link->interphook) PetscCall((*link->interphook)(coarse, interp, fine, link->ctx));
2446:   }
2447:   PetscFunctionReturn(PETSC_SUCCESS);
2448: }

2450: /*@
2451:   DMInterpolateSolution - Interpolates a solution from a coarse mesh to a fine mesh.

2453:   Collective

2455:   Input Parameters:
2456: + coarse    - coarse `DM`
2457: . fine      - fine `DM`
2458: . interp    - (optional) the matrix computed by `DMCreateInterpolation()`.  Implementations may not need this, but if it
2459:             is available it can avoid some recomputation.  If it is provided, `MatInterpolate()` will be used if
2460:             the coarse `DM` does not have a specialized implementation.
2461: - coarseSol - solution on the coarse mesh

2463:   Output Parameter:
2464: . fineSol - the interpolation of coarseSol to the fine mesh

2466:   Level: developer

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

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

2477: .seealso: [](ch_dmbase), `DM`, `DMInterpolate()`, `DMCreateInterpolation()`
2478: @*/
2479: PetscErrorCode DMInterpolateSolution(DM coarse, DM fine, Mat interp, Vec coarseSol, Vec fineSol)
2480: {
2481:   PetscErrorCode (*interpsol)(DM, DM, Mat, Vec, Vec) = NULL;

2483:   PetscFunctionBegin;

2489:   PetscCall(PetscObjectQueryFunction((PetscObject)coarse, "DMInterpolateSolution_C", &interpsol));
2490:   if (interpsol) {
2491:     PetscCall((*interpsol)(coarse, fine, interp, coarseSol, fineSol));
2492:   } else if (interp) {
2493:     PetscCall(MatInterpolate(interp, coarseSol, fineSol));
2494:   } else SETERRQ(PetscObjectComm((PetscObject)coarse), PETSC_ERR_SUP, "DM %s does not implement DMInterpolateSolution()", ((PetscObject)coarse)->type_name);
2495:   PetscFunctionReturn(PETSC_SUCCESS);
2496: }

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

2501:   Not Collective

2503:   Input Parameter:
2504: . dm - the `DM` object

2506:   Output Parameter:
2507: . level - number of refinements

2509:   Level: developer

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

2514: .seealso: [](ch_dmbase), `DM`, `DMRefine()`, `DMCoarsen()`, `DMGetCoarsenLevel()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`
2515: @*/
2516: PetscErrorCode DMGetRefineLevel(DM dm, PetscInt *level)
2517: {
2518:   PetscFunctionBegin;
2520:   *level = dm->levelup;
2521:   PetscFunctionReturn(PETSC_SUCCESS);
2522: }

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

2527:   Not Collective

2529:   Input Parameters:
2530: + dm    - the `DM` object
2531: - level - number of refinements

2533:   Level: advanced

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

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

2540: .seealso: [](ch_dmbase), `DM`, `DMGetRefineLevel()`, `DMCoarsen()`, `DMGetCoarsenLevel()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`
2541: @*/
2542: PetscErrorCode DMSetRefineLevel(DM dm, PetscInt level)
2543: {
2544:   PetscFunctionBegin;
2546:   dm->levelup = level;
2547:   PetscFunctionReturn(PETSC_SUCCESS);
2548: }

2550: /*@
2551:   DMExtrude - Extrude a `DM` object from a surface

2553:   Collective

2555:   Input Parameters:
2556: + dm     - the `DM` object
2557: - layers - the number of extruded cell layers

2559:   Output Parameter:
2560: . dme - the extruded `DM`, or `NULL`

2562:   Level: developer

2564:   Note:
2565:   If no extrusion was done, the return value is `NULL`

2567: .seealso: [](ch_dmbase), `DM`, `DMRefine()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`
2568: @*/
2569: PetscErrorCode DMExtrude(DM dm, PetscInt layers, DM *dme)
2570: {
2571:   PetscFunctionBegin;
2573:   PetscUseTypeMethod(dm, extrude, layers, dme);
2574:   if (*dme) {
2575:     (*dme)->ops->creatematrix = dm->ops->creatematrix;
2576:     PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)dm, (PetscObject)*dme));
2577:     (*dme)->ctx = dm->ctx;
2578:     PetscCall(DMSetMatType(*dme, dm->mattype));
2579:   }
2580:   PetscFunctionReturn(PETSC_SUCCESS);
2581: }

2583: PetscErrorCode DMGetBasisTransformDM_Internal(DM dm, DM *tdm)
2584: {
2585:   PetscFunctionBegin;
2587:   PetscAssertPointer(tdm, 2);
2588:   *tdm = dm->transformDM;
2589:   PetscFunctionReturn(PETSC_SUCCESS);
2590: }

2592: PetscErrorCode DMGetBasisTransformVec_Internal(DM dm, Vec *tv)
2593: {
2594:   PetscFunctionBegin;
2596:   PetscAssertPointer(tv, 2);
2597:   *tv = dm->transform;
2598:   PetscFunctionReturn(PETSC_SUCCESS);
2599: }

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

2604:   Input Parameter:
2605: . dm - The `DM`

2607:   Output Parameter:
2608: . flg - `PETSC_TRUE` if a basis transformation should be done

2610:   Level: developer

2612: .seealso: [](ch_dmbase), `DM`, `DMPlexGlobalToLocalBasis()`, `DMPlexLocalToGlobalBasis()`, `DMPlexCreateBasisRotation()`
2613: @*/
2614: PetscErrorCode DMHasBasisTransform(DM dm, PetscBool *flg)
2615: {
2616:   Vec tv;

2618:   PetscFunctionBegin;
2620:   PetscAssertPointer(flg, 2);
2621:   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
2622:   *flg = tv ? PETSC_TRUE : PETSC_FALSE;
2623:   PetscFunctionReturn(PETSC_SUCCESS);
2624: }

2626: PetscErrorCode DMConstructBasisTransform_Internal(DM dm)
2627: {
2628:   PetscSection s, ts;
2629:   PetscScalar *ta;
2630:   PetscInt     cdim, pStart, pEnd, p, Nf, f, Nc, dof;

2632:   PetscFunctionBegin;
2633:   PetscCall(DMGetCoordinateDim(dm, &cdim));
2634:   PetscCall(DMGetLocalSection(dm, &s));
2635:   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
2636:   PetscCall(PetscSectionGetNumFields(s, &Nf));
2637:   PetscCall(DMClone(dm, &dm->transformDM));
2638:   PetscCall(DMGetLocalSection(dm->transformDM, &ts));
2639:   PetscCall(PetscSectionSetNumFields(ts, Nf));
2640:   PetscCall(PetscSectionSetChart(ts, pStart, pEnd));
2641:   for (f = 0; f < Nf; ++f) {
2642:     PetscCall(PetscSectionGetFieldComponents(s, f, &Nc));
2643:     /* We could start to label fields by their transformation properties */
2644:     if (Nc != cdim) continue;
2645:     for (p = pStart; p < pEnd; ++p) {
2646:       PetscCall(PetscSectionGetFieldDof(s, p, f, &dof));
2647:       if (!dof) continue;
2648:       PetscCall(PetscSectionSetFieldDof(ts, p, f, PetscSqr(cdim)));
2649:       PetscCall(PetscSectionAddDof(ts, p, PetscSqr(cdim)));
2650:     }
2651:   }
2652:   PetscCall(PetscSectionSetUp(ts));
2653:   PetscCall(DMCreateLocalVector(dm->transformDM, &dm->transform));
2654:   PetscCall(VecGetArray(dm->transform, &ta));
2655:   for (p = pStart; p < pEnd; ++p) {
2656:     for (f = 0; f < Nf; ++f) {
2657:       PetscCall(PetscSectionGetFieldDof(ts, p, f, &dof));
2658:       if (dof) {
2659:         PetscReal          x[3] = {0.0, 0.0, 0.0};
2660:         PetscScalar       *tva;
2661:         const PetscScalar *A;

2663:         /* TODO Get quadrature point for this dual basis vector for coordinate */
2664:         PetscCall((*dm->transformGetMatrix)(dm, x, PETSC_TRUE, &A, dm->transformCtx));
2665:         PetscCall(DMPlexPointLocalFieldRef(dm->transformDM, p, f, ta, (void *)&tva));
2666:         PetscCall(PetscArraycpy(tva, A, PetscSqr(cdim)));
2667:       }
2668:     }
2669:   }
2670:   PetscCall(VecRestoreArray(dm->transform, &ta));
2671:   PetscFunctionReturn(PETSC_SUCCESS);
2672: }

2674: PetscErrorCode DMCopyTransform(DM dm, DM newdm)
2675: {
2676:   PetscFunctionBegin;
2679:   newdm->transformCtx       = dm->transformCtx;
2680:   newdm->transformSetUp     = dm->transformSetUp;
2681:   newdm->transformDestroy   = NULL;
2682:   newdm->transformGetMatrix = dm->transformGetMatrix;
2683:   if (newdm->transformSetUp) PetscCall(DMConstructBasisTransform_Internal(newdm));
2684:   PetscFunctionReturn(PETSC_SUCCESS);
2685: }

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

2690:   Logically Collective

2692:   Input Parameters:
2693: + dm        - the `DM`
2694: . beginhook - function to run at the beginning of `DMGlobalToLocalBegin()`
2695: . endhook   - function to run after `DMGlobalToLocalEnd()` has completed
2696: - ctx       - [optional] user-defined context for provide data for the hooks (may be `NULL`)

2698:   Calling sequence of `beginhook`:
2699: + dm   - global `DM`
2700: . g    - global vector
2701: . mode - mode
2702: . l    - local vector
2703: - ctx  - optional user-defined function context

2705:   Calling sequence of `endhook`:
2706: + dm   - global `DM`
2707: . g    - global vector
2708: . mode - mode
2709: . l    - local vector
2710: - ctx  - optional user-defined function context

2712:   Level: advanced

2714:   Note:
2715:   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.

2717: .seealso: [](ch_dmbase), `DM`, `DMGlobalToLocal()`, `DMRefineHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`
2718: @*/
2719: 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)
2720: {
2721:   DMGlobalToLocalHookLink link, *p;

2723:   PetscFunctionBegin;
2725:   for (p = &dm->gtolhook; *p; p = &(*p)->next) { } /* Scan to the end of the current list of hooks */
2726:   PetscCall(PetscNew(&link));
2727:   link->beginhook = beginhook;
2728:   link->endhook   = endhook;
2729:   link->ctx       = ctx;
2730:   link->next      = NULL;
2731:   *p              = link;
2732:   PetscFunctionReturn(PETSC_SUCCESS);
2733: }

2735: static PetscErrorCode DMGlobalToLocalHook_Constraints(DM dm, Vec g, InsertMode mode, Vec l, void *ctx)
2736: {
2737:   Mat          cMat;
2738:   Vec          cVec, cBias;
2739:   PetscSection section, cSec;
2740:   PetscInt     pStart, pEnd, p, dof;

2742:   PetscFunctionBegin;
2743:   (void)g;
2744:   (void)ctx;
2746:   PetscCall(DMGetDefaultConstraints(dm, &cSec, &cMat, &cBias));
2747:   if (cMat && (mode == INSERT_VALUES || mode == INSERT_ALL_VALUES || mode == INSERT_BC_VALUES)) {
2748:     PetscInt nRows;

2750:     PetscCall(MatGetSize(cMat, &nRows, NULL));
2751:     if (nRows <= 0) PetscFunctionReturn(PETSC_SUCCESS);
2752:     PetscCall(DMGetLocalSection(dm, &section));
2753:     PetscCall(MatCreateVecs(cMat, NULL, &cVec));
2754:     PetscCall(MatMult(cMat, l, cVec));
2755:     if (cBias) PetscCall(VecAXPY(cVec, 1., cBias));
2756:     PetscCall(PetscSectionGetChart(cSec, &pStart, &pEnd));
2757:     for (p = pStart; p < pEnd; p++) {
2758:       PetscCall(PetscSectionGetDof(cSec, p, &dof));
2759:       if (dof) {
2760:         PetscScalar *vals;
2761:         PetscCall(VecGetValuesSection(cVec, cSec, p, &vals));
2762:         PetscCall(VecSetValuesSection(l, section, p, vals, INSERT_ALL_VALUES));
2763:       }
2764:     }
2765:     PetscCall(VecDestroy(&cVec));
2766:   }
2767:   PetscFunctionReturn(PETSC_SUCCESS);
2768: }

2770: /*@
2771:   DMGlobalToLocal - update local vectors from global vector

2773:   Neighbor-wise Collective

2775:   Input Parameters:
2776: + dm   - the `DM` object
2777: . g    - the global vector
2778: . mode - `INSERT_VALUES` or `ADD_VALUES`
2779: - l    - the local vector

2781:   Level: beginner

2783:   Notes:
2784:   The communication involved in this update can be overlapped with computation by instead using
2785:   `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()`.

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

2789: .seealso: [](ch_dmbase), `DM`, `DMGlobalToLocalHookAdd()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`,
2790:           `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToGlobal()`, `DMLocalToGlobalEnd()`,
2791:           `DMGlobalToLocalBegin()` `DMGlobalToLocalEnd()`
2792: @*/
2793: PetscErrorCode DMGlobalToLocal(DM dm, Vec g, InsertMode mode, Vec l)
2794: {
2795:   PetscFunctionBegin;
2796:   PetscCall(DMGlobalToLocalBegin(dm, g, mode, l));
2797:   PetscCall(DMGlobalToLocalEnd(dm, g, mode, l));
2798:   PetscFunctionReturn(PETSC_SUCCESS);
2799: }

2801: /*@
2802:   DMGlobalToLocalBegin - Begins updating local vectors from global vector

2804:   Neighbor-wise Collective

2806:   Input Parameters:
2807: + dm   - the `DM` object
2808: . g    - the global vector
2809: . mode - `INSERT_VALUES` or `ADD_VALUES`
2810: - l    - the local vector

2812:   Level: intermediate

2814:   Notes:
2815:   The operation is completed with `DMGlobalToLocalEnd()`

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

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

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

2823: .seealso: [](ch_dmbase), `DM`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocal()`, `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToGlobal()`, `DMLocalToGlobalEnd()`
2824: @*/
2825: PetscErrorCode DMGlobalToLocalBegin(DM dm, Vec g, InsertMode mode, Vec l)
2826: {
2827:   PetscSF                 sf;
2828:   DMGlobalToLocalHookLink link;

2830:   PetscFunctionBegin;
2832:   for (link = dm->gtolhook; link; link = link->next) {
2833:     if (link->beginhook) PetscCall((*link->beginhook)(dm, g, mode, l, link->ctx));
2834:   }
2835:   PetscCall(DMGetSectionSF(dm, &sf));
2836:   if (sf) {
2837:     const PetscScalar *gArray;
2838:     PetscScalar       *lArray;
2839:     PetscMemType       lmtype, gmtype;

2841:     PetscCheck(mode != ADD_VALUES, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %d", (int)mode);
2842:     PetscCall(VecGetArrayAndMemType(l, &lArray, &lmtype));
2843:     PetscCall(VecGetArrayReadAndMemType(g, &gArray, &gmtype));
2844:     PetscCall(PetscSFBcastWithMemTypeBegin(sf, MPIU_SCALAR, gmtype, gArray, lmtype, lArray, MPI_REPLACE));
2845:     PetscCall(VecRestoreArrayAndMemType(l, &lArray));
2846:     PetscCall(VecRestoreArrayReadAndMemType(g, &gArray));
2847:   } else {
2848:     PetscUseTypeMethod(dm, globaltolocalbegin, g, mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode), l);
2849:   }
2850:   PetscFunctionReturn(PETSC_SUCCESS);
2851: }

2853: /*@
2854:   DMGlobalToLocalEnd - Ends updating local vectors from global vector

2856:   Neighbor-wise Collective

2858:   Input Parameters:
2859: + dm   - the `DM` object
2860: . g    - the global vector
2861: . mode - `INSERT_VALUES` or `ADD_VALUES`
2862: - l    - the local vector

2864:   Level: intermediate

2866:   Note:
2867:   See `DMGlobalToLocalBegin()` for details.

2869: .seealso: [](ch_dmbase), `DM`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocal()`, `DMLocalToGlobalBegin()`, `DMLocalToGlobal()`, `DMLocalToGlobalEnd()`
2870: @*/
2871: PetscErrorCode DMGlobalToLocalEnd(DM dm, Vec g, InsertMode mode, Vec l)
2872: {
2873:   PetscSF                 sf;
2874:   const PetscScalar      *gArray;
2875:   PetscScalar            *lArray;
2876:   PetscBool               transform;
2877:   DMGlobalToLocalHookLink link;
2878:   PetscMemType            lmtype, gmtype;

2880:   PetscFunctionBegin;
2882:   PetscCall(DMGetSectionSF(dm, &sf));
2883:   PetscCall(DMHasBasisTransform(dm, &transform));
2884:   if (sf) {
2885:     PetscCheck(mode != ADD_VALUES, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %d", (int)mode);

2887:     PetscCall(VecGetArrayAndMemType(l, &lArray, &lmtype));
2888:     PetscCall(VecGetArrayReadAndMemType(g, &gArray, &gmtype));
2889:     PetscCall(PetscSFBcastEnd(sf, MPIU_SCALAR, gArray, lArray, MPI_REPLACE));
2890:     PetscCall(VecRestoreArrayAndMemType(l, &lArray));
2891:     PetscCall(VecRestoreArrayReadAndMemType(g, &gArray));
2892:     if (transform) PetscCall(DMPlexGlobalToLocalBasis(dm, l));
2893:   } else {
2894:     PetscUseTypeMethod(dm, globaltolocalend, g, mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode), l);
2895:   }
2896:   PetscCall(DMGlobalToLocalHook_Constraints(dm, g, mode, l, NULL));
2897:   for (link = dm->gtolhook; link; link = link->next) {
2898:     if (link->endhook) PetscCall((*link->endhook)(dm, g, mode, l, link->ctx));
2899:   }
2900:   PetscFunctionReturn(PETSC_SUCCESS);
2901: }

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

2906:   Logically Collective

2908:   Input Parameters:
2909: + dm        - the `DM`
2910: . beginhook - function to run at the beginning of `DMLocalToGlobalBegin()`
2911: . endhook   - function to run after `DMLocalToGlobalEnd()` has completed
2912: - ctx       - [optional] user-defined context for provide data for the hooks (may be `NULL`)

2914:   Calling sequence of `beginhook`:
2915: + global - global `DM`
2916: . l      - local vector
2917: . mode   - mode
2918: . g      - global vector
2919: - ctx    - optional user-defined function context

2921:   Calling sequence of `endhook`:
2922: + global - global `DM`
2923: . l      - local vector
2924: . mode   - mode
2925: . g      - global vector
2926: - ctx    - optional user-defined function context

2928:   Level: advanced

2930: .seealso: [](ch_dmbase), `DM`, `DMLocalToGlobal()`, `DMRefineHookAdd()`, `DMGlobalToLocalHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`
2931: @*/
2932: 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)
2933: {
2934:   DMLocalToGlobalHookLink link, *p;

2936:   PetscFunctionBegin;
2938:   for (p = &dm->ltoghook; *p; p = &(*p)->next) { } /* Scan to the end of the current list of hooks */
2939:   PetscCall(PetscNew(&link));
2940:   link->beginhook = beginhook;
2941:   link->endhook   = endhook;
2942:   link->ctx       = ctx;
2943:   link->next      = NULL;
2944:   *p              = link;
2945:   PetscFunctionReturn(PETSC_SUCCESS);
2946: }

2948: static PetscErrorCode DMLocalToGlobalHook_Constraints(DM dm, Vec l, InsertMode mode, Vec g, void *ctx)
2949: {
2950:   PetscFunctionBegin;
2951:   (void)g;
2952:   (void)ctx;
2954:   if (mode == ADD_VALUES || mode == ADD_ALL_VALUES || mode == ADD_BC_VALUES) {
2955:     Mat          cMat;
2956:     Vec          cVec;
2957:     PetscInt     nRows;
2958:     PetscSection section, cSec;
2959:     PetscInt     pStart, pEnd, p, dof;

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

2964:     PetscCall(MatGetSize(cMat, &nRows, NULL));
2965:     if (nRows <= 0) PetscFunctionReturn(PETSC_SUCCESS);
2966:     PetscCall(DMGetLocalSection(dm, &section));
2967:     PetscCall(MatCreateVecs(cMat, NULL, &cVec));
2968:     PetscCall(PetscSectionGetChart(cSec, &pStart, &pEnd));
2969:     for (p = pStart; p < pEnd; p++) {
2970:       PetscCall(PetscSectionGetDof(cSec, p, &dof));
2971:       if (dof) {
2972:         PetscInt     d;
2973:         PetscScalar *vals;
2974:         PetscCall(VecGetValuesSection(l, section, p, &vals));
2975:         PetscCall(VecSetValuesSection(cVec, cSec, p, vals, mode));
2976:         /* for this to be the true transpose, we have to zero the values that
2977:          * we just extracted */
2978:         for (d = 0; d < dof; d++) vals[d] = 0.;
2979:       }
2980:     }
2981:     PetscCall(MatMultTransposeAdd(cMat, cVec, l, l));
2982:     PetscCall(VecDestroy(&cVec));
2983:   }
2984:   PetscFunctionReturn(PETSC_SUCCESS);
2985: }
2986: /*@
2987:   DMLocalToGlobal - updates global vectors from local vectors

2989:   Neighbor-wise Collective

2991:   Input Parameters:
2992: + dm   - the `DM` object
2993: . l    - the local vector
2994: . 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.
2995: - g    - the global vector

2997:   Level: beginner

2999:   Notes:
3000:   The communication involved in this update can be overlapped with computation by using
3001:   `DMLocalToGlobalBegin()` and `DMLocalToGlobalEnd()`.

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

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

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

3009: .seealso: [](ch_dmbase), `DM`, `DMLocalToGlobalBegin()`, `DMLocalToGlobalEnd()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocal()`, `DMGlobalToLocalEnd()`, `DMGlobalToLocalBegin()`, `DMLocalToGlobalHookAdd()`, `DMGlobaToLocallHookAdd()`
3010: @*/
3011: PetscErrorCode DMLocalToGlobal(DM dm, Vec l, InsertMode mode, Vec g)
3012: {
3013:   PetscFunctionBegin;
3014:   PetscCall(DMLocalToGlobalBegin(dm, l, mode, g));
3015:   PetscCall(DMLocalToGlobalEnd(dm, l, mode, g));
3016:   PetscFunctionReturn(PETSC_SUCCESS);
3017: }

3019: /*@
3020:   DMLocalToGlobalBegin - begins updating global vectors from local vectors

3022:   Neighbor-wise Collective

3024:   Input Parameters:
3025: + dm   - the `DM` object
3026: . l    - the local vector
3027: . 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.
3028: - g    - the global vector

3030:   Level: intermediate

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

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

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

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

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

3043: .seealso: [](ch_dmbase), `DM`, `DMLocalToGlobal()`, `DMLocalToGlobalEnd()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocal()`, `DMGlobalToLocalEnd()`, `DMGlobalToLocalBegin()`
3044: @*/
3045: PetscErrorCode DMLocalToGlobalBegin(DM dm, Vec l, InsertMode mode, Vec g)
3046: {
3047:   PetscSF                 sf;
3048:   PetscSection            s, gs;
3049:   DMLocalToGlobalHookLink link;
3050:   Vec                     tmpl;
3051:   const PetscScalar      *lArray;
3052:   PetscScalar            *gArray;
3053:   PetscBool               isInsert, transform, l_inplace = PETSC_FALSE, g_inplace = PETSC_FALSE;
3054:   PetscMemType            lmtype = PETSC_MEMTYPE_HOST, gmtype = PETSC_MEMTYPE_HOST;

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

3102:       PetscCall(DMGetGlobalSection(dm, &gs));
3103:       PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
3104:       PetscCall(VecGetOwnershipRange(g, &gStart, NULL));
3105:       for (p = pStart; p < pEnd; ++p) {
3106:         PetscInt dof, gdof, cdof, gcdof, off, goff, d, e;

3108:         PetscCall(PetscSectionGetDof(s, p, &dof));
3109:         PetscCall(PetscSectionGetDof(gs, p, &gdof));
3110:         PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
3111:         PetscCall(PetscSectionGetConstraintDof(gs, p, &gcdof));
3112:         PetscCall(PetscSectionGetOffset(s, p, &off));
3113:         PetscCall(PetscSectionGetOffset(gs, p, &goff));
3114:         /* Ignore off-process data and points with no global data */
3115:         if (!gdof || goff < 0) continue;
3116:         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);
3117:         /* If no constraints are enforced in the global vector */
3118:         if (!gcdof) {
3119:           for (d = 0; d < dof; ++d) gArray[goff - gStart + d] = lArray[off + d];
3120:           /* If constraints are enforced in the global vector */
3121:         } else if (cdof == gcdof) {
3122:           const PetscInt *cdofs;
3123:           PetscInt        cind = 0;

3125:           PetscCall(PetscSectionGetConstraintIndices(s, p, &cdofs));
3126:           for (d = 0, e = 0; d < dof; ++d) {
3127:             if ((cind < cdof) && (d == cdofs[cind])) {
3128:               ++cind;
3129:               continue;
3130:             }
3131:             gArray[goff - gStart + e++] = lArray[off + d];
3132:           }
3133:         } 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);
3134:       }
3135:     }
3136:     if (g_inplace) {
3137:       PetscCall(VecRestoreArrayAndMemType(g, &gArray));
3138:     } else {
3139:       PetscCall(VecRestoreArray(g, &gArray));
3140:     }
3141:     if (transform) {
3142:       PetscCall(VecRestoreArrayRead(tmpl, &lArray));
3143:       PetscCall(DMRestoreNamedLocalVector(dm, "__petsc_dm_transform_local_copy", &tmpl));
3144:     } else if (l_inplace) {
3145:       PetscCall(VecRestoreArrayReadAndMemType(l, &lArray));
3146:     } else {
3147:       PetscCall(VecRestoreArrayRead(l, &lArray));
3148:     }
3149:   } else {
3150:     PetscUseTypeMethod(dm, localtoglobalbegin, l, mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode), g);
3151:   }
3152:   PetscFunctionReturn(PETSC_SUCCESS);
3153: }

3155: /*@
3156:   DMLocalToGlobalEnd - updates global vectors from local vectors

3158:   Neighbor-wise Collective

3160:   Input Parameters:
3161: + dm   - the `DM` object
3162: . l    - the local vector
3163: . mode - `INSERT_VALUES` or `ADD_VALUES`
3164: - g    - the global vector

3166:   Level: intermediate

3168:   Note:
3169:   See `DMLocalToGlobalBegin()` for full details

3171: .seealso: [](ch_dmbase), `DM`, `DMLocalToGlobalBegin()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocalEnd()`
3172: @*/
3173: PetscErrorCode DMLocalToGlobalEnd(DM dm, Vec l, InsertMode mode, Vec g)
3174: {
3175:   PetscSF                 sf;
3176:   PetscSection            s;
3177:   DMLocalToGlobalHookLink link;
3178:   PetscBool               isInsert, transform;

3180:   PetscFunctionBegin;
3182:   PetscCall(DMGetSectionSF(dm, &sf));
3183:   PetscCall(DMGetLocalSection(dm, &s));
3184:   switch (mode) {
3185:   case INSERT_VALUES:
3186:   case INSERT_ALL_VALUES:
3187:     isInsert = PETSC_TRUE;
3188:     break;
3189:   case ADD_VALUES:
3190:   case ADD_ALL_VALUES:
3191:     isInsert = PETSC_FALSE;
3192:     break;
3193:   default:
3194:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %d", mode);
3195:   }
3196:   if (sf && !isInsert) {
3197:     const PetscScalar *lArray;
3198:     PetscScalar       *gArray;
3199:     Vec                tmpl;

3201:     PetscCall(DMHasBasisTransform(dm, &transform));
3202:     if (transform) {
3203:       PetscCall(DMGetNamedLocalVector(dm, "__petsc_dm_transform_local_copy", &tmpl));
3204:       PetscCall(VecGetArrayRead(tmpl, &lArray));
3205:     } else {
3206:       PetscCall(VecGetArrayReadAndMemType(l, &lArray, NULL));
3207:     }
3208:     PetscCall(VecGetArrayAndMemType(g, &gArray, NULL));
3209:     PetscCall(PetscSFReduceEnd(sf, MPIU_SCALAR, lArray, gArray, MPIU_SUM));
3210:     if (transform) {
3211:       PetscCall(VecRestoreArrayRead(tmpl, &lArray));
3212:       PetscCall(DMRestoreNamedLocalVector(dm, "__petsc_dm_transform_local_copy", &tmpl));
3213:     } else {
3214:       PetscCall(VecRestoreArrayReadAndMemType(l, &lArray));
3215:     }
3216:     PetscCall(VecRestoreArrayAndMemType(g, &gArray));
3217:   } else if (s && isInsert) {
3218:   } else {
3219:     PetscUseTypeMethod(dm, localtoglobalend, l, mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode), g);
3220:   }
3221:   for (link = dm->ltoghook; link; link = link->next) {
3222:     if (link->endhook) PetscCall((*link->endhook)(dm, g, mode, l, link->ctx));
3223:   }
3224:   PetscFunctionReturn(PETSC_SUCCESS);
3225: }

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

3232:   Neighbor-wise Collective

3234:   Input Parameters:
3235: + dm   - the `DM` object
3236: . g    - the original local vector
3237: - mode - one of `INSERT_VALUES` or `ADD_VALUES`

3239:   Output Parameter:
3240: . l - the local vector with correct ghost values

3242:   Level: intermediate

3244:   Note:
3245:   Must be followed by `DMLocalToLocalEnd()`.

3247: .seealso: [](ch_dmbase), `DM`, `DMLocalToLocalEnd()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateLocalVector()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`
3248: @*/
3249: PetscErrorCode DMLocalToLocalBegin(DM dm, Vec g, InsertMode mode, Vec l)
3250: {
3251:   PetscFunctionBegin;
3255:   PetscUseTypeMethod(dm, localtolocalbegin, g, mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode), l);
3256:   PetscFunctionReturn(PETSC_SUCCESS);
3257: }

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

3263:   Neighbor-wise Collective

3265:   Input Parameters:
3266: + dm   - the `DM` object
3267: . g    - the original local vector
3268: - mode - one of `INSERT_VALUES` or `ADD_VALUES`

3270:   Output Parameter:
3271: . l - the local vector with correct ghost values

3273:   Level: intermediate

3275: .seealso: [](ch_dmbase), `DM`, `DMLocalToLocalBegin()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateLocalVector()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`
3276: @*/
3277: PetscErrorCode DMLocalToLocalEnd(DM dm, Vec g, InsertMode mode, Vec l)
3278: {
3279:   PetscFunctionBegin;
3283:   PetscUseTypeMethod(dm, localtolocalend, g, mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode), l);
3284:   PetscFunctionReturn(PETSC_SUCCESS);
3285: }

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

3290:   Collective

3292:   Input Parameters:
3293: + dm   - the `DM` object
3294: - comm - the communicator to contain the new `DM` object (or `MPI_COMM_NULL`)

3296:   Output Parameter:
3297: . dmc - the coarsened `DM`

3299:   Level: developer

3301: .seealso: [](ch_dmbase), `DM`, `DMRefine()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateDomainDecomposition()`,
3302:           `DMCoarsenHookAdd()`, `DMCoarsenHookRemove()`
3303: @*/
3304: PetscErrorCode DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
3305: {
3306:   DMCoarsenHookLink link;

3308:   PetscFunctionBegin;
3310:   PetscCall(PetscLogEventBegin(DM_Coarsen, dm, 0, 0, 0));
3311:   PetscUseTypeMethod(dm, coarsen, comm, dmc);
3312:   if (*dmc) {
3313:     (*dmc)->bind_below = dm->bind_below; /* Propagate this from parent DM; otherwise -dm_bind_below will be useless for multigrid cases. */
3314:     PetscCall(DMSetCoarseDM(dm, *dmc));
3315:     (*dmc)->ops->creatematrix = dm->ops->creatematrix;
3316:     PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)dm, (PetscObject)*dmc));
3317:     (*dmc)->ctx       = dm->ctx;
3318:     (*dmc)->levelup   = dm->levelup;
3319:     (*dmc)->leveldown = dm->leveldown + 1;
3320:     PetscCall(DMSetMatType(*dmc, dm->mattype));
3321:     for (link = dm->coarsenhook; link; link = link->next) {
3322:       if (link->coarsenhook) PetscCall((*link->coarsenhook)(dm, *dmc, link->ctx));
3323:     }
3324:   }
3325:   PetscCall(PetscLogEventEnd(DM_Coarsen, dm, 0, 0, 0));
3326:   PetscCheck(*dmc, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "NULL coarse mesh produced");
3327:   PetscFunctionReturn(PETSC_SUCCESS);
3328: }

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

3333:   Logically Collective; No Fortran Support

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

3341:   Calling sequence of `coarsenhook`:
3342: + fine   - fine level `DM`
3343: . coarse - coarse level `DM` to restrict problem to
3344: - ctx    - optional user-defined function context

3346:   Calling sequence of `restricthook`:
3347: + fine      - fine level `DM`
3348: . mrestrict - matrix restricting a fine-level solution to the coarse grid, usually the transpose of the interpolation
3349: . rscale    - scaling vector for restriction
3350: . inject    - matrix restricting by injection
3351: . coarse    - coarse level DM to update
3352: - ctx       - optional user-defined function context

3354:   Level: advanced

3356:   Notes:
3357:   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`.

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

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

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

3366: .seealso: [](ch_dmbase), `DM`, `DMCoarsenHookRemove()`, `DMRefineHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`
3367: @*/
3368: 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)
3369: {
3370:   DMCoarsenHookLink link, *p;

3372:   PetscFunctionBegin;
3374:   for (p = &fine->coarsenhook; *p; p = &(*p)->next) { /* Scan to the end of the current list of hooks */
3375:     if ((*p)->coarsenhook == coarsenhook && (*p)->restricthook == restricthook && (*p)->ctx == ctx) PetscFunctionReturn(PETSC_SUCCESS);
3376:   }
3377:   PetscCall(PetscNew(&link));
3378:   link->coarsenhook  = coarsenhook;
3379:   link->restricthook = restricthook;
3380:   link->ctx          = ctx;
3381:   link->next         = NULL;
3382:   *p                 = link;
3383:   PetscFunctionReturn(PETSC_SUCCESS);
3384: }

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

3389:   Logically Collective; No Fortran Support

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

3397:   Level: advanced

3399:   Notes:
3400:   This function does nothing if the `coarsenhook` is not in the list.

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

3404: .seealso: [](ch_dmbase), `DM`, `DMCoarsenHookAdd()`, `DMRefineHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`
3405: @*/
3406: PetscErrorCode DMCoarsenHookRemove(DM fine, PetscErrorCode (*coarsenhook)(DM, DM, void *), PetscErrorCode (*restricthook)(DM, Mat, Vec, Mat, DM, void *), void *ctx)
3407: {
3408:   DMCoarsenHookLink link, *p;

3410:   PetscFunctionBegin;
3412:   for (p = &fine->coarsenhook; *p; p = &(*p)->next) { /* Search the list of current hooks */
3413:     if ((*p)->coarsenhook == coarsenhook && (*p)->restricthook == restricthook && (*p)->ctx == ctx) {
3414:       link = *p;
3415:       *p   = link->next;
3416:       PetscCall(PetscFree(link));
3417:       break;
3418:     }
3419:   }
3420:   PetscFunctionReturn(PETSC_SUCCESS);
3421: }

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

3426:   Collective if any hooks are

3428:   Input Parameters:
3429: + fine    - finer `DM` from which the data is obtained
3430: . restrct - restriction matrix, apply using `MatRestrict()`, usually the transpose of the interpolation
3431: . rscale  - scaling vector for restriction
3432: . inject  - injection matrix, also use `MatRestrict()`
3433: - coarse  - coarser `DM` to update

3435:   Level: developer

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

3440: .seealso: [](ch_dmbase), `DM`, `DMCoarsenHookAdd()`, `MatRestrict()`, `DMInterpolate()`, `DMRefineHookAdd()`
3441: @*/
3442: PetscErrorCode DMRestrict(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse)
3443: {
3444:   DMCoarsenHookLink link;

3446:   PetscFunctionBegin;
3447:   for (link = fine->coarsenhook; link; link = link->next) {
3448:     if (link->restricthook) PetscCall((*link->restricthook)(fine, restrct, rscale, inject, coarse, link->ctx));
3449:   }
3450:   PetscFunctionReturn(PETSC_SUCCESS);
3451: }

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

3456:   Logically Collective; No Fortran Support

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

3464:   Calling sequence of `ddhook`:
3465: + global - global `DM`
3466: . block  - subdomain `DM`
3467: - ctx    - optional user-defined function context

3469:   Calling sequence of `restricthook`:
3470: + global - global `DM`
3471: . out    - scatter to the outer (with ghost and overlap points) sub vector
3472: . in     - scatter to sub vector values only owned locally
3473: . block  - subdomain `DM`
3474: - ctx    - optional user-defined function context

3476:   Level: advanced

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

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

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

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

3489: .seealso: [](ch_dmbase), `DM`, `DMSubDomainHookRemove()`, `DMRefineHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`, `DMCreateDomainDecomposition()`
3490: @*/
3491: 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)
3492: {
3493:   DMSubDomainHookLink link, *p;

3495:   PetscFunctionBegin;
3497:   for (p = &global->subdomainhook; *p; p = &(*p)->next) { /* Scan to the end of the current list of hooks */
3498:     if ((*p)->ddhook == ddhook && (*p)->restricthook == restricthook && (*p)->ctx == ctx) PetscFunctionReturn(PETSC_SUCCESS);
3499:   }
3500:   PetscCall(PetscNew(&link));
3501:   link->restricthook = restricthook;
3502:   link->ddhook       = ddhook;
3503:   link->ctx          = ctx;
3504:   link->next         = NULL;
3505:   *p                 = link;
3506:   PetscFunctionReturn(PETSC_SUCCESS);
3507: }

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

3512:   Logically Collective; No Fortran Support

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

3520:   Level: advanced

3522:   Note:
3523:   See `DMSubDomainHookAdd()` for the calling sequences of `ddhook` and `restricthook`

3525: .seealso: [](ch_dmbase), `DM`, `DMSubDomainHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`,
3526:           `DMCreateDomainDecomposition()`
3527: @*/
3528: PetscErrorCode DMSubDomainHookRemove(DM global, PetscErrorCode (*ddhook)(DM, DM, void *), PetscErrorCode (*restricthook)(DM, VecScatter, VecScatter, DM, void *), void *ctx)
3529: {
3530:   DMSubDomainHookLink link, *p;

3532:   PetscFunctionBegin;
3534:   for (p = &global->subdomainhook; *p; p = &(*p)->next) { /* Search the list of current hooks */
3535:     if ((*p)->ddhook == ddhook && (*p)->restricthook == restricthook && (*p)->ctx == ctx) {
3536:       link = *p;
3537:       *p   = link->next;
3538:       PetscCall(PetscFree(link));
3539:       break;
3540:     }
3541:   }
3542:   PetscFunctionReturn(PETSC_SUCCESS);
3543: }

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

3548:   Collective if any hooks are

3550:   Input Parameters:
3551: + global   - The global `DM` to use as a base
3552: . oscatter - The scatter from domain global vector filling subdomain global vector with overlap
3553: . gscatter - The scatter from domain global vector filling subdomain local vector with ghosts
3554: - subdm    - The subdomain `DM` to update

3556:   Level: developer

3558: .seealso: [](ch_dmbase), `DM`, `DMCoarsenHookAdd()`, `MatRestrict()`, `DMCreateDomainDecomposition()`
3559: @*/
3560: PetscErrorCode DMSubDomainRestrict(DM global, VecScatter oscatter, VecScatter gscatter, DM subdm)
3561: {
3562:   DMSubDomainHookLink link;

3564:   PetscFunctionBegin;
3565:   for (link = global->subdomainhook; link; link = link->next) {
3566:     if (link->restricthook) PetscCall((*link->restricthook)(global, oscatter, gscatter, subdm, link->ctx));
3567:   }
3568:   PetscFunctionReturn(PETSC_SUCCESS);
3569: }

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

3574:   Not Collective

3576:   Input Parameter:
3577: . dm - the `DM` object

3579:   Output Parameter:
3580: . level - number of coarsenings

3582:   Level: developer

3584: .seealso: [](ch_dmbase), `DM`, `DMCoarsen()`, `DMSetCoarsenLevel()`, `DMGetRefineLevel()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`
3585: @*/
3586: PetscErrorCode DMGetCoarsenLevel(DM dm, PetscInt *level)
3587: {
3588:   PetscFunctionBegin;
3590:   PetscAssertPointer(level, 2);
3591:   *level = dm->leveldown;
3592:   PetscFunctionReturn(PETSC_SUCCESS);
3593: }

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

3598:   Collective

3600:   Input Parameters:
3601: + dm    - the `DM` object
3602: - level - number of coarsenings

3604:   Level: developer

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

3609: .seealso: [](ch_dmbase), `DM`, `DMCoarsen()`, `DMGetCoarsenLevel()`, `DMGetRefineLevel()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`
3610: @*/
3611: PetscErrorCode DMSetCoarsenLevel(DM dm, PetscInt level)
3612: {
3613:   PetscFunctionBegin;
3615:   dm->leveldown = level;
3616:   PetscFunctionReturn(PETSC_SUCCESS);
3617: }

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

3622:   Collective

3624:   Input Parameters:
3625: + dm      - the `DM` object
3626: - nlevels - the number of levels of refinement

3628:   Output Parameter:
3629: . dmf - the refined `DM` hierarchy

3631:   Level: developer

3633: .seealso: [](ch_dmbase), `DM`, `DMCoarsen()`, `DMCoarsenHierarchy()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`
3634: @*/
3635: PetscErrorCode DMRefineHierarchy(DM dm, PetscInt nlevels, DM dmf[])
3636: {
3637:   PetscFunctionBegin;
3639:   PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative");
3640:   if (nlevels == 0) PetscFunctionReturn(PETSC_SUCCESS);
3641:   PetscAssertPointer(dmf, 3);
3642:   if (dm->ops->refine && !dm->ops->refinehierarchy) {
3643:     PetscInt i;

3645:     PetscCall(DMRefine(dm, PetscObjectComm((PetscObject)dm), &dmf[0]));
3646:     for (i = 1; i < nlevels; i++) PetscCall(DMRefine(dmf[i - 1], PetscObjectComm((PetscObject)dm), &dmf[i]));
3647:   } else PetscUseTypeMethod(dm, refinehierarchy, nlevels, dmf);
3648:   PetscFunctionReturn(PETSC_SUCCESS);
3649: }

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

3654:   Collective

3656:   Input Parameters:
3657: + dm      - the `DM` object
3658: - nlevels - the number of levels of coarsening

3660:   Output Parameter:
3661: . dmc - the coarsened `DM` hierarchy

3663:   Level: developer

3665: .seealso: [](ch_dmbase), `DM`, `DMCoarsen()`, `DMRefineHierarchy()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`
3666: @*/
3667: PetscErrorCode DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
3668: {
3669:   PetscFunctionBegin;
3671:   PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative");
3672:   if (nlevels == 0) PetscFunctionReturn(PETSC_SUCCESS);
3673:   PetscAssertPointer(dmc, 3);
3674:   if (dm->ops->coarsen && !dm->ops->coarsenhierarchy) {
3675:     PetscInt i;

3677:     PetscCall(DMCoarsen(dm, PetscObjectComm((PetscObject)dm), &dmc[0]));
3678:     for (i = 1; i < nlevels; i++) PetscCall(DMCoarsen(dmc[i - 1], PetscObjectComm((PetscObject)dm), &dmc[i]));
3679:   } else PetscUseTypeMethod(dm, coarsenhierarchy, nlevels, dmc);
3680:   PetscFunctionReturn(PETSC_SUCCESS);
3681: }

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

3686:   Logically Collective if the function is collective

3688:   Input Parameters:
3689: + dm      - the `DM` object
3690: - destroy - the destroy function

3692:   Level: intermediate

3694: .seealso: [](ch_dmbase), `DM`, `DMSetApplicationContext()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMGetApplicationContext()`
3695: @*/
3696: PetscErrorCode DMSetApplicationContextDestroy(DM dm, PetscErrorCode (*destroy)(void **))
3697: {
3698:   PetscFunctionBegin;
3700:   dm->ctxdestroy = destroy;
3701:   PetscFunctionReturn(PETSC_SUCCESS);
3702: }

3704: /*@
3705:   DMSetApplicationContext - Set a user context into a `DM` object

3707:   Not Collective

3709:   Input Parameters:
3710: + dm  - the `DM` object
3711: - ctx - the user context

3713:   Level: intermediate

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

3720: .seealso: [](ch_dmbase), `DM`, `DMGetApplicationContext()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`
3721: @*/
3722: PetscErrorCode DMSetApplicationContext(DM dm, void *ctx)
3723: {
3724:   PetscFunctionBegin;
3726:   dm->ctx = ctx;
3727:   PetscFunctionReturn(PETSC_SUCCESS);
3728: }

3730: /*@
3731:   DMGetApplicationContext - Gets a user context from a `DM` object

3733:   Not Collective

3735:   Input Parameter:
3736: . dm - the `DM` object

3738:   Output Parameter:
3739: . ctx - the user context

3741:   Level: intermediate

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

3746: .seealso: [](ch_dmbase), `DM`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`
3747: @*/
3748: PetscErrorCode DMGetApplicationContext(DM dm, void *ctx)
3749: {
3750:   PetscFunctionBegin;
3752:   *(void **)ctx = dm->ctx;
3753:   PetscFunctionReturn(PETSC_SUCCESS);
3754: }

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

3759:   Logically Collective

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

3765:   Level: intermediate

3767: .seealso: [](ch_dmbase), `DM`, `DMComputeVariableBounds()`, `DMHasVariableBounds()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMGetApplicationContext()`,
3768:          `DMSetJacobian()`
3769: @*/
3770: PetscErrorCode DMSetVariableBounds(DM dm, PetscErrorCode (*f)(DM, Vec, Vec))
3771: {
3772:   PetscFunctionBegin;
3774:   dm->ops->computevariablebounds = f;
3775:   PetscFunctionReturn(PETSC_SUCCESS);
3776: }

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

3781:   Not Collective

3783:   Input Parameter:
3784: . dm - the `DM` object to destroy

3786:   Output Parameter:
3787: . flg - `PETSC_TRUE` if the variable bounds function exists

3789:   Level: developer

3791: .seealso: [](ch_dmbase), `DM`, `DMComputeVariableBounds()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMGetApplicationContext()`
3792: @*/
3793: PetscErrorCode DMHasVariableBounds(DM dm, PetscBool *flg)
3794: {
3795:   PetscFunctionBegin;
3797:   PetscAssertPointer(flg, 2);
3798:   *flg = (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE;
3799:   PetscFunctionReturn(PETSC_SUCCESS);
3800: }

3802: /*@
3803:   DMComputeVariableBounds - compute variable bounds used by `SNESVI`.

3805:   Logically Collective

3807:   Input Parameter:
3808: . dm - the `DM` object

3810:   Output Parameters:
3811: + xl - lower bound
3812: - xu - upper bound

3814:   Level: advanced

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

3819: .seealso: [](ch_dmbase), `DM`, `DMHasVariableBounds()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMGetApplicationContext()`
3820: @*/
3821: PetscErrorCode DMComputeVariableBounds(DM dm, Vec xl, Vec xu)
3822: {
3823:   PetscFunctionBegin;
3827:   PetscUseTypeMethod(dm, computevariablebounds, xl, xu);
3828:   PetscFunctionReturn(PETSC_SUCCESS);
3829: }

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

3834:   Not Collective

3836:   Input Parameter:
3837: . dm - the DM object

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

3842:   Level: developer

3844: .seealso: [](ch_dmbase), `DM`, `DMCreateColoring()`
3845: @*/
3846: PetscErrorCode DMHasColoring(DM dm, PetscBool *flg)
3847: {
3848:   PetscFunctionBegin;
3850:   PetscAssertPointer(flg, 2);
3851:   *flg = (dm->ops->getcoloring) ? PETSC_TRUE : PETSC_FALSE;
3852:   PetscFunctionReturn(PETSC_SUCCESS);
3853: }

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

3858:   Not Collective

3860:   Input Parameter:
3861: . dm - the `DM` object

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

3866:   Level: developer

3868: .seealso: [](ch_dmbase), `DM`, `DMCreateRestriction()`, `DMHasCreateInterpolation()`, `DMHasCreateInjection()`
3869: @*/
3870: PetscErrorCode DMHasCreateRestriction(DM dm, PetscBool *flg)
3871: {
3872:   PetscFunctionBegin;
3874:   PetscAssertPointer(flg, 2);
3875:   *flg = (dm->ops->createrestriction) ? PETSC_TRUE : PETSC_FALSE;
3876:   PetscFunctionReturn(PETSC_SUCCESS);
3877: }

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

3882:   Not Collective

3884:   Input Parameter:
3885: . dm - the `DM` object

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

3890:   Level: developer

3892: .seealso: [](ch_dmbase), `DM`, `DMCreateInjection()`, `DMHasCreateRestriction()`, `DMHasCreateInterpolation()`
3893: @*/
3894: PetscErrorCode DMHasCreateInjection(DM dm, PetscBool *flg)
3895: {
3896:   PetscFunctionBegin;
3898:   PetscAssertPointer(flg, 2);
3899:   if (dm->ops->hascreateinjection) PetscUseTypeMethod(dm, hascreateinjection, flg);
3900:   else *flg = (dm->ops->createinjection) ? PETSC_TRUE : PETSC_FALSE;
3901:   PetscFunctionReturn(PETSC_SUCCESS);
3902: }

3904: PetscFunctionList DMList              = NULL;
3905: PetscBool         DMRegisterAllCalled = PETSC_FALSE;

3907: /*@
3908:   DMSetType - Builds a `DM`, for a particular `DM` implementation.

3910:   Collective

3912:   Input Parameters:
3913: + dm     - The `DM` object
3914: - method - The name of the `DMType`, for example `DMDA`, `DMPLEX`

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

3919:   Level: intermediate

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

3924: .seealso: [](ch_dmbase), `DM`, `DMType`, `DMDA`, `DMPLEX`, `DMGetType()`, `DMCreate()`, `DMDACreate2d()`
3925: @*/
3926: PetscErrorCode DMSetType(DM dm, DMType method)
3927: {
3928:   PetscErrorCode (*r)(DM);
3929:   PetscBool match;

3931:   PetscFunctionBegin;
3933:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, method, &match));
3934:   if (match) PetscFunctionReturn(PETSC_SUCCESS);

3936:   PetscCall(DMRegisterAll());
3937:   PetscCall(PetscFunctionListFind(DMList, method, &r));
3938:   PetscCheck(r, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);

3940:   PetscTryTypeMethod(dm, destroy);
3941:   PetscCall(PetscMemzero(dm->ops, sizeof(*dm->ops)));
3942:   PetscCall(PetscObjectChangeTypeName((PetscObject)dm, method));
3943:   PetscCall((*r)(dm));
3944:   PetscFunctionReturn(PETSC_SUCCESS);
3945: }

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

3950:   Not Collective

3952:   Input Parameter:
3953: . dm - The `DM`

3955:   Output Parameter:
3956: . type - The `DMType` name

3958:   Level: intermediate

3960: .seealso: [](ch_dmbase), `DM`, `DMType`, `DMDA`, `DMPLEX`, `DMSetType()`, `DMCreate()`
3961: @*/
3962: PetscErrorCode DMGetType(DM dm, DMType *type)
3963: {
3964:   PetscFunctionBegin;
3966:   PetscAssertPointer(type, 2);
3967:   PetscCall(DMRegisterAll());
3968:   *type = ((PetscObject)dm)->type_name;
3969:   PetscFunctionReturn(PETSC_SUCCESS);
3970: }

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

3975:   Collective

3977:   Input Parameters:
3978: + dm      - the `DM`
3979: - newtype - new `DM` type (use "same" for the same type)

3981:   Output Parameter:
3982: . M - pointer to new `DM`

3984:   Level: intermediate

3986:   Note:
3987:   Cannot be used to convert a sequential `DM` to a parallel or a parallel to sequential,
3988:   the MPI communicator of the generated `DM` is always the same as the communicator
3989:   of the input `DM`.

3991: .seealso: [](ch_dmbase), `DM`, `DMSetType()`, `DMCreate()`, `DMClone()`
3992: @*/
3993: PetscErrorCode DMConvert(DM dm, DMType newtype, DM *M)
3994: {
3995:   DM        B;
3996:   char      convname[256];
3997:   PetscBool sametype /*, issame */;

3999:   PetscFunctionBegin;
4002:   PetscAssertPointer(M, 3);
4003:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, newtype, &sametype));
4004:   /* PetscCall(PetscStrcmp(newtype, "same", &issame)); */
4005:   if (sametype) {
4006:     *M = dm;
4007:     PetscCall(PetscObjectReference((PetscObject)dm));
4008:     PetscFunctionReturn(PETSC_SUCCESS);
4009:   } else {
4010:     PetscErrorCode (*conv)(DM, DMType, DM *) = NULL;

4012:     /*
4013:        Order of precedence:
4014:        1) See if a specialized converter is known to the current DM.
4015:        2) See if a specialized converter is known to the desired DM class.
4016:        3) See if a good general converter is registered for the desired class
4017:        4) See if a good general converter is known for the current matrix.
4018:        5) Use a really basic converter.
4019:     */

4021:     /* 1) See if a specialized converter is known to the current DM and the desired class */
4022:     PetscCall(PetscStrncpy(convname, "DMConvert_", sizeof(convname)));
4023:     PetscCall(PetscStrlcat(convname, ((PetscObject)dm)->type_name, sizeof(convname)));
4024:     PetscCall(PetscStrlcat(convname, "_", sizeof(convname)));
4025:     PetscCall(PetscStrlcat(convname, newtype, sizeof(convname)));
4026:     PetscCall(PetscStrlcat(convname, "_C", sizeof(convname)));
4027:     PetscCall(PetscObjectQueryFunction((PetscObject)dm, convname, &conv));
4028:     if (conv) goto foundconv;

4030:     /* 2)  See if a specialized converter is known to the desired DM class. */
4031:     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &B));
4032:     PetscCall(DMSetType(B, newtype));
4033:     PetscCall(PetscStrncpy(convname, "DMConvert_", sizeof(convname)));
4034:     PetscCall(PetscStrlcat(convname, ((PetscObject)dm)->type_name, sizeof(convname)));
4035:     PetscCall(PetscStrlcat(convname, "_", sizeof(convname)));
4036:     PetscCall(PetscStrlcat(convname, newtype, sizeof(convname)));
4037:     PetscCall(PetscStrlcat(convname, "_C", sizeof(convname)));
4038:     PetscCall(PetscObjectQueryFunction((PetscObject)B, convname, &conv));
4039:     if (conv) {
4040:       PetscCall(DMDestroy(&B));
4041:       goto foundconv;
4042:     }

4044: #if 0
4045:     /* 3) See if a good general converter is registered for the desired class */
4046:     conv = B->ops->convertfrom;
4047:     PetscCall(DMDestroy(&B));
4048:     if (conv) goto foundconv;

4050:     /* 4) See if a good general converter is known for the current matrix */
4051:     if (dm->ops->convert) {
4052:       conv = dm->ops->convert;
4053:     }
4054:     if (conv) goto foundconv;
4055: #endif

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

4060:   foundconv:
4061:     PetscCall(PetscLogEventBegin(DM_Convert, dm, 0, 0, 0));
4062:     PetscCall((*conv)(dm, newtype, M));
4063:     /* Things that are independent of DM type: We should consult DMClone() here */
4064:     {
4065:       const PetscReal *maxCell, *Lstart, *L;

4067:       PetscCall(DMGetPeriodicity(dm, &maxCell, &Lstart, &L));
4068:       PetscCall(DMSetPeriodicity(*M, maxCell, Lstart, L));
4069:       (*M)->prealloc_only = dm->prealloc_only;
4070:       PetscCall(PetscFree((*M)->vectype));
4071:       PetscCall(PetscStrallocpy(dm->vectype, (char **)&(*M)->vectype));
4072:       PetscCall(PetscFree((*M)->mattype));
4073:       PetscCall(PetscStrallocpy(dm->mattype, (char **)&(*M)->mattype));
4074:     }
4075:     PetscCall(PetscLogEventEnd(DM_Convert, dm, 0, 0, 0));
4076:   }
4077:   PetscCall(PetscObjectStateIncrease((PetscObject)*M));
4078:   PetscFunctionReturn(PETSC_SUCCESS);
4079: }

4081: /*--------------------------------------------------------------------------------------------------------------------*/

4083: /*@C
4084:   DMRegister -  Adds a new `DM` type implementation

4086:   Not Collective, No Fortran Support

4088:   Input Parameters:
4089: + sname    - The name of a new user-defined creation routine
4090: - function - The creation routine itself

4092:   Level: advanced

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

4097:   Example Usage:
4098: .vb
4099:     DMRegister("my_da", MyDMCreate);
4100: .ve

4102:   Then, your `DM` type can be chosen with the procedural interface via
4103: .vb
4104:     DMCreate(MPI_Comm, DM *);
4105:     DMSetType(DM,"my_da");
4106: .ve
4107:   or at runtime via the option
4108: .vb
4109:     -da_type my_da
4110: .ve

4112: .seealso: [](ch_dmbase), `DM`, `DMType`, `DMSetType()`, `DMRegisterAll()`, `DMRegisterDestroy()`
4113: @*/
4114: PetscErrorCode DMRegister(const char sname[], PetscErrorCode (*function)(DM))
4115: {
4116:   PetscFunctionBegin;
4117:   PetscCall(DMInitializePackage());
4118:   PetscCall(PetscFunctionListAdd(&DMList, sname, function));
4119:   PetscFunctionReturn(PETSC_SUCCESS);
4120: }

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

4125:   Collective

4127:   Input Parameters:
4128: + newdm  - the newly loaded `DM`, this needs to have been created with `DMCreate()` or
4129:            some related function before a call to `DMLoad()`.
4130: - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()` or
4131:            `PETSCVIEWERHDF5` file viewer, obtained from `PetscViewerHDF5Open()`

4133:   Level: intermediate

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

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

4142: .seealso: [](ch_dmbase), `DM`, `PetscViewerBinaryOpen()`, `DMView()`, `MatLoad()`, `VecLoad()`
4143: @*/
4144: PetscErrorCode DMLoad(DM newdm, PetscViewer viewer)
4145: {
4146:   PetscBool isbinary, ishdf5;

4148:   PetscFunctionBegin;
4151:   PetscCall(PetscViewerCheckReadable(viewer));
4152:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
4153:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
4154:   PetscCall(PetscLogEventBegin(DM_Load, viewer, 0, 0, 0));
4155:   if (isbinary) {
4156:     PetscInt classid;
4157:     char     type[256];

4159:     PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT));
4160:     PetscCheck(classid == DM_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not DM next in file, classid found %d", (int)classid);
4161:     PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
4162:     PetscCall(DMSetType(newdm, type));
4163:     PetscTryTypeMethod(newdm, load, viewer);
4164:   } else if (ishdf5) {
4165:     PetscTryTypeMethod(newdm, load, viewer);
4166:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen() or PetscViewerHDF5Open()");
4167:   PetscCall(PetscLogEventEnd(DM_Load, viewer, 0, 0, 0));
4168:   PetscFunctionReturn(PETSC_SUCCESS);
4169: }

4171: /******************************** FEM Support **********************************/

4173: PetscErrorCode DMPrintCellIndices(PetscInt c, const char name[], PetscInt len, const PetscInt x[])
4174: {
4175:   PetscInt f;

4177:   PetscFunctionBegin;
4178:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " Element %s\n", c, name));
4179:   for (f = 0; f < len; ++f) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  | %" PetscInt_FMT " |\n", x[f]));
4180:   PetscFunctionReturn(PETSC_SUCCESS);
4181: }

4183: PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[])
4184: {
4185:   PetscInt f;

4187:   PetscFunctionBegin;
4188:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " Element %s\n", c, name));
4189:   for (f = 0; f < len; ++f) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  | %g |\n", (double)PetscRealPart(x[f])));
4190:   PetscFunctionReturn(PETSC_SUCCESS);
4191: }

4193: PetscErrorCode DMPrintCellVectorReal(PetscInt c, const char name[], PetscInt len, const PetscReal x[])
4194: {
4195:   PetscInt f;

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

4203: PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[])
4204: {
4205:   PetscInt f, g;

4207:   PetscFunctionBegin;
4208:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " Element %s\n", c, name));
4209:   for (f = 0; f < rows; ++f) {
4210:     PetscCall(PetscPrintf(PETSC_COMM_SELF, "  |"));
4211:     for (g = 0; g < cols; ++g) PetscCall(PetscPrintf(PETSC_COMM_SELF, " % 9.5g", (double)PetscRealPart(A[f * cols + g])));
4212:     PetscCall(PetscPrintf(PETSC_COMM_SELF, " |\n"));
4213:   }
4214:   PetscFunctionReturn(PETSC_SUCCESS);
4215: }

4217: PetscErrorCode DMPrintLocalVec(DM dm, const char name[], PetscReal tol, Vec X)
4218: {
4219:   PetscInt           localSize, bs;
4220:   PetscMPIInt        size;
4221:   Vec                x, xglob;
4222:   const PetscScalar *xarray;

4224:   PetscFunctionBegin;
4225:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
4226:   PetscCall(VecDuplicate(X, &x));
4227:   PetscCall(VecCopy(X, x));
4228:   PetscCall(VecFilter(x, tol));
4229:   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "%s:\n", name));
4230:   if (size > 1) {
4231:     PetscCall(VecGetLocalSize(x, &localSize));
4232:     PetscCall(VecGetArrayRead(x, &xarray));
4233:     PetscCall(VecGetBlockSize(x, &bs));
4234:     PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)dm), bs, localSize, PETSC_DETERMINE, xarray, &xglob));
4235:   } else {
4236:     xglob = x;
4237:   }
4238:   PetscCall(VecView(xglob, PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm))));
4239:   if (size > 1) {
4240:     PetscCall(VecDestroy(&xglob));
4241:     PetscCall(VecRestoreArrayRead(x, &xarray));
4242:   }
4243:   PetscCall(VecDestroy(&x));
4244:   PetscFunctionReturn(PETSC_SUCCESS);
4245: }

4247: /*@
4248:   DMGetSection - Get the `PetscSection` encoding the local data layout for the `DM`.   This is equivalent to `DMGetLocalSection()`. Deprecated in v3.12

4250:   Input Parameter:
4251: . dm - The `DM`

4253:   Output Parameter:
4254: . section - The `PetscSection`

4256:   Options Database Key:
4257: . -dm_petscsection_view - View the `PetscSection` created by the `DM`

4259:   Level: advanced

4261:   Notes:
4262:   Use `DMGetLocalSection()` in new code.

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

4266: .seealso: [](ch_dmbase), `DM`, `DMGetLocalSection()`, `DMSetLocalSection()`, `DMGetGlobalSection()`
4267: @*/
4268: PetscErrorCode DMGetSection(DM dm, PetscSection *section)
4269: {
4270:   PetscFunctionBegin;
4271:   PetscCall(DMGetLocalSection(dm, section));
4272:   PetscFunctionReturn(PETSC_SUCCESS);
4273: }

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

4278:   Input Parameter:
4279: . dm - The `DM`

4281:   Output Parameter:
4282: . section - The `PetscSection`

4284:   Options Database Key:
4285: . -dm_petscsection_view - View the section created by the `DM`

4287:   Level: intermediate

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

4292: .seealso: [](ch_dmbase), `DM`, `DMSetLocalSection()`, `DMGetGlobalSection()`
4293: @*/
4294: PetscErrorCode DMGetLocalSection(DM dm, PetscSection *section)
4295: {
4296:   PetscFunctionBegin;
4298:   PetscAssertPointer(section, 2);
4299:   if (!dm->localSection && dm->ops->createlocalsection) {
4300:     PetscInt d;

4302:     if (dm->setfromoptionscalled) {
4303:       PetscObject       obj = (PetscObject)dm;
4304:       PetscViewer       viewer;
4305:       PetscViewerFormat format;
4306:       PetscBool         flg;

4308:       PetscCall(PetscOptionsCreateViewer(PetscObjectComm(obj), obj->options, obj->prefix, "-dm_petscds_view", &viewer, &format, &flg));
4309:       if (flg) PetscCall(PetscViewerPushFormat(viewer, format));
4310:       for (d = 0; d < dm->Nds; ++d) {
4311:         PetscCall(PetscDSSetFromOptions(dm->probs[d].ds));
4312:         if (flg) PetscCall(PetscDSView(dm->probs[d].ds, viewer));
4313:       }
4314:       if (flg) {
4315:         PetscCall(PetscViewerFlush(viewer));
4316:         PetscCall(PetscViewerPopFormat(viewer));
4317:         PetscCall(PetscViewerDestroy(&viewer));
4318:       }
4319:     }
4320:     PetscUseTypeMethod(dm, createlocalsection);
4321:     if (dm->localSection) PetscCall(PetscObjectViewFromOptions((PetscObject)dm->localSection, NULL, "-dm_petscsection_view"));
4322:   }
4323:   *section = dm->localSection;
4324:   PetscFunctionReturn(PETSC_SUCCESS);
4325: }

4327: /*@
4328:   DMSetSection - Set the `PetscSection` encoding the local data layout for the `DM`.  This is equivalent to `DMSetLocalSection()`. Deprecated in v3.12

4330:   Input Parameters:
4331: + dm      - The `DM`
4332: - section - The `PetscSection`

4334:   Level: advanced

4336:   Notes:
4337:   Use `DMSetLocalSection()` in new code.

4339:   Any existing `PetscSection` will be destroyed

4341: .seealso: [](ch_dmbase), `DM`, `DMSetLocalSection()`, `DMGetLocalSection()`, `DMSetGlobalSection()`
4342: @*/
4343: PetscErrorCode DMSetSection(DM dm, PetscSection section)
4344: {
4345:   PetscFunctionBegin;
4346:   PetscCall(DMSetLocalSection(dm, section));
4347:   PetscFunctionReturn(PETSC_SUCCESS);
4348: }

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

4353:   Input Parameters:
4354: + dm      - The `DM`
4355: - section - The `PetscSection`

4357:   Level: intermediate

4359:   Note:
4360:   Any existing Section will be destroyed

4362: .seealso: [](ch_dmbase), `DM`, `PetscSection`, `DMGetLocalSection()`, `DMSetGlobalSection()`
4363: @*/
4364: PetscErrorCode DMSetLocalSection(DM dm, PetscSection section)
4365: {
4366:   PetscInt numFields = 0;
4367:   PetscInt f;

4369:   PetscFunctionBegin;
4372:   PetscCall(PetscObjectReference((PetscObject)section));
4373:   PetscCall(PetscSectionDestroy(&dm->localSection));
4374:   dm->localSection = section;
4375:   if (section) PetscCall(PetscSectionGetNumFields(dm->localSection, &numFields));
4376:   if (numFields) {
4377:     PetscCall(DMSetNumFields(dm, numFields));
4378:     for (f = 0; f < numFields; ++f) {
4379:       PetscObject disc;
4380:       const char *name;

4382:       PetscCall(PetscSectionGetFieldName(dm->localSection, f, &name));
4383:       PetscCall(DMGetField(dm, f, NULL, &disc));
4384:       PetscCall(PetscObjectSetName(disc, name));
4385:     }
4386:   }
4387:   /* The global section and the SectionSF will be rebuilt
4388:      in the next call to DMGetGlobalSection() and DMGetSectionSF(). */
4389:   PetscCall(PetscSectionDestroy(&dm->globalSection));
4390:   PetscCall(PetscSFDestroy(&dm->sectionSF));
4391:   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &dm->sectionSF));

4393:   /* Clear scratch vectors */
4394:   PetscCall(DMClearGlobalVectors(dm));
4395:   PetscCall(DMClearLocalVectors(dm));
4396:   PetscCall(DMClearNamedGlobalVectors(dm));
4397:   PetscCall(DMClearNamedLocalVectors(dm));
4398:   PetscFunctionReturn(PETSC_SUCCESS);
4399: }

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

4404:   Input Parameter:
4405: . dm - The `DM`

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

4411:   Level: developer

4413: .seealso: [](ch_dmbase), `DM`, `PetscSection`, `DMGetLocalSection()`, `DMGetGlobalSection()`
4414: @*/
4415: PetscErrorCode DMCreateSectionPermutation(DM dm, IS *perm, PetscBT *blockStarts)
4416: {
4417:   PetscFunctionBegin;
4418:   *perm        = NULL;
4419:   *blockStarts = NULL;
4420:   PetscTryTypeMethod(dm, createsectionpermutation, perm, blockStarts);
4421:   PetscFunctionReturn(PETSC_SUCCESS);
4422: }

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

4427:   not Collective

4429:   Input Parameter:
4430: . dm - The `DM`

4432:   Output Parameters:
4433: + 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.
4434: . 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.
4435: - bias    - Vector containing bias to be added to constrained dofs

4437:   Level: advanced

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

4442: .seealso: [](ch_dmbase), `DM`, `DMSetDefaultConstraints()`
4443: @*/
4444: PetscErrorCode DMGetDefaultConstraints(DM dm, PetscSection *section, Mat *mat, Vec *bias)
4445: {
4446:   PetscFunctionBegin;
4448:   if (!dm->defaultConstraint.section && !dm->defaultConstraint.mat && dm->ops->createdefaultconstraints) PetscUseTypeMethod(dm, createdefaultconstraints);
4449:   if (section) *section = dm->defaultConstraint.section;
4450:   if (mat) *mat = dm->defaultConstraint.mat;
4451:   if (bias) *bias = dm->defaultConstraint.bias;
4452:   PetscFunctionReturn(PETSC_SUCCESS);
4453: }

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

4458:   Collective

4460:   Input Parameters:
4461: + dm      - The `DM`
4462: . 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).
4463: . 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).
4464: - 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).

4466:   Level: advanced

4468:   Notes:
4469:   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()`.

4471:   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.

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

4475: .seealso: [](ch_dmbase), `DM`, `DMGetDefaultConstraints()`
4476: @*/
4477: PetscErrorCode DMSetDefaultConstraints(DM dm, PetscSection section, Mat mat, Vec bias)
4478: {
4479:   PetscMPIInt result;

4481:   PetscFunctionBegin;
4483:   if (section) {
4485:     PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, PetscObjectComm((PetscObject)section), &result));
4486:     PetscCheck(result == MPI_CONGRUENT || result == MPI_IDENT, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMECOMM, "constraint section must have local communicator");
4487:   }
4488:   if (mat) {
4490:     PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, PetscObjectComm((PetscObject)mat), &result));
4491:     PetscCheck(result == MPI_CONGRUENT || result == MPI_IDENT, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMECOMM, "constraint matrix must have local communicator");
4492:   }
4493:   if (bias) {
4495:     PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, PetscObjectComm((PetscObject)bias), &result));
4496:     PetscCheck(result == MPI_CONGRUENT || result == MPI_IDENT, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMECOMM, "constraint bias must have local communicator");
4497:   }
4498:   PetscCall(PetscObjectReference((PetscObject)section));
4499:   PetscCall(PetscSectionDestroy(&dm->defaultConstraint.section));
4500:   dm->defaultConstraint.section = section;
4501:   PetscCall(PetscObjectReference((PetscObject)mat));
4502:   PetscCall(MatDestroy(&dm->defaultConstraint.mat));
4503:   dm->defaultConstraint.mat = mat;
4504:   PetscCall(PetscObjectReference((PetscObject)bias));
4505:   PetscCall(VecDestroy(&dm->defaultConstraint.bias));
4506:   dm->defaultConstraint.bias = bias;
4507:   PetscFunctionReturn(PETSC_SUCCESS);
4508: }

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

4514:   Input Parameters:
4515: + dm - The `DM`
4516: . localSection - `PetscSection` describing the local data layout
4517: - globalSection - `PetscSection` describing the global data layout

4519:   Level: intermediate

4521: .seealso: [](ch_dmbase), `DM`, `DMGetSectionSF()`, `DMSetSectionSF()`
4522: */
4523: static PetscErrorCode DMDefaultSectionCheckConsistency_Internal(DM dm, PetscSection localSection, PetscSection globalSection)
4524: {
4525:   MPI_Comm        comm;
4526:   PetscLayout     layout;
4527:   const PetscInt *ranges;
4528:   PetscInt        pStart, pEnd, p, nroots;
4529:   PetscMPIInt     size, rank;
4530:   PetscBool       valid = PETSC_TRUE, gvalid;

4532:   PetscFunctionBegin;
4533:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
4535:   PetscCallMPI(MPI_Comm_size(comm, &size));
4536:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
4537:   PetscCall(PetscSectionGetChart(globalSection, &pStart, &pEnd));
4538:   PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &nroots));
4539:   PetscCall(PetscLayoutCreate(comm, &layout));
4540:   PetscCall(PetscLayoutSetBlockSize(layout, 1));
4541:   PetscCall(PetscLayoutSetLocalSize(layout, nroots));
4542:   PetscCall(PetscLayoutSetUp(layout));
4543:   PetscCall(PetscLayoutGetRanges(layout, &ranges));
4544:   for (p = pStart; p < pEnd; ++p) {
4545:     PetscInt dof, cdof, off, gdof, gcdof, goff, gsize, d;

4547:     PetscCall(PetscSectionGetDof(localSection, p, &dof));
4548:     PetscCall(PetscSectionGetOffset(localSection, p, &off));
4549:     PetscCall(PetscSectionGetConstraintDof(localSection, p, &cdof));
4550:     PetscCall(PetscSectionGetDof(globalSection, p, &gdof));
4551:     PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof));
4552:     PetscCall(PetscSectionGetOffset(globalSection, p, &goff));
4553:     if (!gdof) continue; /* Censored point */
4554:     if ((gdof < 0 ? -(gdof + 1) : gdof) != dof) {
4555:       PetscCall(PetscSynchronizedPrintf(comm, "[%d]Global dof %" PetscInt_FMT " for point %" PetscInt_FMT " not equal to local dof %" PetscInt_FMT "\n", rank, gdof, p, dof));
4556:       valid = PETSC_FALSE;
4557:     }
4558:     if (gcdof && (gcdof != cdof)) {
4559:       PetscCall(PetscSynchronizedPrintf(comm, "[%d]Global constraints %" PetscInt_FMT " for point %" PetscInt_FMT " not equal to local constraints %" PetscInt_FMT "\n", rank, gcdof, p, cdof));
4560:       valid = PETSC_FALSE;
4561:     }
4562:     if (gdof < 0) {
4563:       gsize = gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof;
4564:       for (d = 0; d < gsize; ++d) {
4565:         PetscInt offset = -(goff + 1) + d, r;

4567:         PetscCall(PetscFindInt(offset, size + 1, ranges, &r));
4568:         if (r < 0) r = -(r + 2);
4569:         if ((r < 0) || (r >= size)) {
4570:           PetscCall(PetscSynchronizedPrintf(comm, "[%d]Point %" PetscInt_FMT " mapped to invalid process %" PetscInt_FMT " (%" PetscInt_FMT ", %" PetscInt_FMT ")\n", rank, p, r, gdof, goff));
4571:           valid = PETSC_FALSE;
4572:           break;
4573:         }
4574:       }
4575:     }
4576:   }
4577:   PetscCall(PetscLayoutDestroy(&layout));
4578:   PetscCall(PetscSynchronizedFlush(comm, NULL));
4579:   PetscCall(MPIU_Allreduce(&valid, &gvalid, 1, MPIU_BOOL, MPI_LAND, comm));
4580:   if (!gvalid) {
4581:     PetscCall(DMView(dm, NULL));
4582:     SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Inconsistent local and global sections");
4583:   }
4584:   PetscFunctionReturn(PETSC_SUCCESS);
4585: }
4586: #endif

4588: static PetscErrorCode DMGetIsoperiodicPointSF_Internal(DM dm, PetscSF *sf)
4589: {
4590:   PetscErrorCode (*f)(DM, PetscSF *);

4592:   PetscFunctionBegin;
4594:   PetscAssertPointer(sf, 2);
4595:   PetscCall(PetscObjectQueryFunction((PetscObject)dm, "DMGetIsoperiodicPointSF_C", &f));
4596:   if (f) PetscCall(f(dm, sf));
4597:   else *sf = dm->sf;
4598:   PetscFunctionReturn(PETSC_SUCCESS);
4599: }

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

4604:   Collective

4606:   Input Parameter:
4607: . dm - The `DM`

4609:   Output Parameter:
4610: . section - The `PetscSection`

4612:   Level: intermediate

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

4617: .seealso: [](ch_dmbase), `DM`, `DMSetLocalSection()`, `DMGetLocalSection()`
4618: @*/
4619: PetscErrorCode DMGetGlobalSection(DM dm, PetscSection *section)
4620: {
4621:   PetscFunctionBegin;
4623:   PetscAssertPointer(section, 2);
4624:   if (!dm->globalSection) {
4625:     PetscSection s;
4626:     PetscSF      sf;

4628:     PetscCall(DMGetLocalSection(dm, &s));
4629:     PetscCheck(s, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM must have a default PetscSection in order to create a global PetscSection");
4630:     PetscCheck(dm->sf, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM must have a point PetscSF in order to create a global PetscSection");
4631:     PetscCall(DMGetIsoperiodicPointSF_Internal(dm, &sf));
4632:     PetscCall(PetscSectionCreateGlobalSection(s, sf, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &dm->globalSection));
4633:     PetscCall(PetscLayoutDestroy(&dm->map));
4634:     PetscCall(PetscSectionGetValueLayout(PetscObjectComm((PetscObject)dm), dm->globalSection, &dm->map));
4635:     PetscCall(PetscSectionViewFromOptions(dm->globalSection, NULL, "-global_section_view"));
4636:   }
4637:   *section = dm->globalSection;
4638:   PetscFunctionReturn(PETSC_SUCCESS);
4639: }

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

4644:   Input Parameters:
4645: + dm      - The `DM`
4646: - section - The PetscSection, or `NULL`

4648:   Level: intermediate

4650:   Note:
4651:   Any existing `PetscSection` will be destroyed

4653: .seealso: [](ch_dmbase), `DM`, `DMGetGlobalSection()`, `DMSetLocalSection()`
4654: @*/
4655: PetscErrorCode DMSetGlobalSection(DM dm, PetscSection section)
4656: {
4657:   PetscFunctionBegin;
4660:   PetscCall(PetscObjectReference((PetscObject)section));
4661:   PetscCall(PetscSectionDestroy(&dm->globalSection));
4662:   dm->globalSection = section;
4663: #if defined(PETSC_USE_DEBUG)
4664:   if (section) PetscCall(DMDefaultSectionCheckConsistency_Internal(dm, dm->localSection, section));
4665: #endif
4666:   /* Clear global scratch vectors and sectionSF */
4667:   PetscCall(PetscSFDestroy(&dm->sectionSF));
4668:   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &dm->sectionSF));
4669:   PetscCall(DMClearGlobalVectors(dm));
4670:   PetscCall(DMClearNamedGlobalVectors(dm));
4671:   PetscFunctionReturn(PETSC_SUCCESS);
4672: }

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

4678:   Input Parameter:
4679: . dm - The `DM`

4681:   Output Parameter:
4682: . sf - The `PetscSF`

4684:   Level: intermediate

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

4689: .seealso: [](ch_dmbase), `DM`, `DMSetSectionSF()`, `DMCreateSectionSF()`
4690: @*/
4691: PetscErrorCode DMGetSectionSF(DM dm, PetscSF *sf)
4692: {
4693:   PetscInt nroots;

4695:   PetscFunctionBegin;
4697:   PetscAssertPointer(sf, 2);
4698:   if (!dm->sectionSF) PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &dm->sectionSF));
4699:   PetscCall(PetscSFGetGraph(dm->sectionSF, &nroots, NULL, NULL, NULL));
4700:   if (nroots < 0) {
4701:     PetscSection section, gSection;

4703:     PetscCall(DMGetLocalSection(dm, &section));
4704:     if (section) {
4705:       PetscCall(DMGetGlobalSection(dm, &gSection));
4706:       PetscCall(DMCreateSectionSF(dm, section, gSection));
4707:     } else {
4708:       *sf = NULL;
4709:       PetscFunctionReturn(PETSC_SUCCESS);
4710:     }
4711:   }
4712:   *sf = dm->sectionSF;
4713:   PetscFunctionReturn(PETSC_SUCCESS);
4714: }

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

4719:   Input Parameters:
4720: + dm - The `DM`
4721: - sf - The `PetscSF`

4723:   Level: intermediate

4725:   Note:
4726:   Any previous `PetscSF` is destroyed

4728: .seealso: [](ch_dmbase), `DM`, `DMGetSectionSF()`, `DMCreateSectionSF()`
4729: @*/
4730: PetscErrorCode DMSetSectionSF(DM dm, PetscSF sf)
4731: {
4732:   PetscFunctionBegin;
4735:   PetscCall(PetscObjectReference((PetscObject)sf));
4736:   PetscCall(PetscSFDestroy(&dm->sectionSF));
4737:   dm->sectionSF = sf;
4738:   PetscFunctionReturn(PETSC_SUCCESS);
4739: }

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

4745:   Input Parameters:
4746: + dm            - The `DM`
4747: . localSection  - `PetscSection` describing the local data layout
4748: - globalSection - `PetscSection` describing the global data layout

4750:   Level: developer

4752:   Note:
4753:   One usually uses `DMGetSectionSF()` to obtain the `PetscSF`

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

4761: .seealso: [](ch_dmbase), `DM`, `DMGetSectionSF()`, `DMSetSectionSF()`, `DMGetLocalSection()`, `DMGetGlobalSection()`
4762: @*/
4763: PetscErrorCode DMCreateSectionSF(DM dm, PetscSection localSection, PetscSection globalSection)
4764: {
4765:   PetscFunctionBegin;
4767:   PetscCall(PetscSFSetGraphSection(dm->sectionSF, localSection, globalSection));
4768:   PetscFunctionReturn(PETSC_SUCCESS);
4769: }

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

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

4776:   Input Parameter:
4777: . dm - The `DM`

4779:   Output Parameter:
4780: . sf - The `PetscSF`

4782:   Level: intermediate

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

4787: .seealso: [](ch_dmbase), `DM`, `DMSetPointSF()`, `DMGetSectionSF()`, `DMSetSectionSF()`, `DMCreateSectionSF()`
4788: @*/
4789: PetscErrorCode DMGetPointSF(DM dm, PetscSF *sf)
4790: {
4791:   PetscFunctionBegin;
4793:   PetscAssertPointer(sf, 2);
4794:   *sf = dm->sf;
4795:   PetscFunctionReturn(PETSC_SUCCESS);
4796: }

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

4801:   Collective

4803:   Input Parameters:
4804: + dm - The `DM`
4805: - sf - The `PetscSF`

4807:   Level: intermediate

4809: .seealso: [](ch_dmbase), `DM`, `DMGetPointSF()`, `DMGetSectionSF()`, `DMSetSectionSF()`, `DMCreateSectionSF()`
4810: @*/
4811: PetscErrorCode DMSetPointSF(DM dm, PetscSF sf)
4812: {
4813:   PetscFunctionBegin;
4816:   PetscCall(PetscObjectReference((PetscObject)sf));
4817:   PetscCall(PetscSFDestroy(&dm->sf));
4818:   dm->sf = sf;
4819:   PetscFunctionReturn(PETSC_SUCCESS);
4820: }

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

4825:   Input Parameter:
4826: . dm - The `DM`

4828:   Output Parameter:
4829: . sf - The `PetscSF`

4831:   Level: intermediate

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

4836: .seealso: [](ch_dmbase), `DM`, `DMSetNaturalSF()`, `DMSetUseNatural()`, `DMGetUseNatural()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexDistribute()`
4837: @*/
4838: PetscErrorCode DMGetNaturalSF(DM dm, PetscSF *sf)
4839: {
4840:   PetscFunctionBegin;
4842:   PetscAssertPointer(sf, 2);
4843:   *sf = dm->sfNatural;
4844:   PetscFunctionReturn(PETSC_SUCCESS);
4845: }

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

4850:   Input Parameters:
4851: + dm - The DM
4852: - sf - The PetscSF

4854:   Level: intermediate

4856: .seealso: [](ch_dmbase), `DM`, `DMGetNaturalSF()`, `DMSetUseNatural()`, `DMGetUseNatural()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexDistribute()`
4857: @*/
4858: PetscErrorCode DMSetNaturalSF(DM dm, PetscSF sf)
4859: {
4860:   PetscFunctionBegin;
4863:   PetscCall(PetscObjectReference((PetscObject)sf));
4864:   PetscCall(PetscSFDestroy(&dm->sfNatural));
4865:   dm->sfNatural = sf;
4866:   PetscFunctionReturn(PETSC_SUCCESS);
4867: }

4869: static PetscErrorCode DMSetDefaultAdjacency_Private(DM dm, PetscInt f, PetscObject disc)
4870: {
4871:   PetscClassId id;

4873:   PetscFunctionBegin;
4874:   PetscCall(PetscObjectGetClassId(disc, &id));
4875:   if (id == PETSCFE_CLASSID) {
4876:     PetscCall(DMSetAdjacency(dm, f, PETSC_FALSE, PETSC_TRUE));
4877:   } else if (id == PETSCFV_CLASSID) {
4878:     PetscCall(DMSetAdjacency(dm, f, PETSC_TRUE, PETSC_FALSE));
4879:   } else {
4880:     PetscCall(DMSetAdjacency(dm, f, PETSC_FALSE, PETSC_TRUE));
4881:   }
4882:   PetscFunctionReturn(PETSC_SUCCESS);
4883: }

4885: static PetscErrorCode DMFieldEnlarge_Static(DM dm, PetscInt NfNew)
4886: {
4887:   RegionField *tmpr;
4888:   PetscInt     Nf = dm->Nf, f;

4890:   PetscFunctionBegin;
4891:   if (Nf >= NfNew) PetscFunctionReturn(PETSC_SUCCESS);
4892:   PetscCall(PetscMalloc1(NfNew, &tmpr));
4893:   for (f = 0; f < Nf; ++f) tmpr[f] = dm->fields[f];
4894:   for (f = Nf; f < NfNew; ++f) {
4895:     tmpr[f].disc        = NULL;
4896:     tmpr[f].label       = NULL;
4897:     tmpr[f].avoidTensor = PETSC_FALSE;
4898:   }
4899:   PetscCall(PetscFree(dm->fields));
4900:   dm->Nf     = NfNew;
4901:   dm->fields = tmpr;
4902:   PetscFunctionReturn(PETSC_SUCCESS);
4903: }

4905: /*@
4906:   DMClearFields - Remove all fields from the `DM`

4908:   Logically Collective

4910:   Input Parameter:
4911: . dm - The `DM`

4913:   Level: intermediate

4915: .seealso: [](ch_dmbase), `DM`, `DMGetNumFields()`, `DMSetNumFields()`, `DMSetField()`
4916: @*/
4917: PetscErrorCode DMClearFields(DM dm)
4918: {
4919:   PetscInt f;

4921:   PetscFunctionBegin;
4923:   for (f = 0; f < dm->Nf; ++f) {
4924:     PetscCall(PetscObjectDestroy(&dm->fields[f].disc));
4925:     PetscCall(DMLabelDestroy(&dm->fields[f].label));
4926:   }
4927:   PetscCall(PetscFree(dm->fields));
4928:   dm->fields = NULL;
4929:   dm->Nf     = 0;
4930:   PetscFunctionReturn(PETSC_SUCCESS);
4931: }

4933: /*@
4934:   DMGetNumFields - Get the number of fields in the `DM`

4936:   Not Collective

4938:   Input Parameter:
4939: . dm - The `DM`

4941:   Output Parameter:
4942: . numFields - The number of fields

4944:   Level: intermediate

4946: .seealso: [](ch_dmbase), `DM`, `DMSetNumFields()`, `DMSetField()`
4947: @*/
4948: PetscErrorCode DMGetNumFields(DM dm, PetscInt *numFields)
4949: {
4950:   PetscFunctionBegin;
4952:   PetscAssertPointer(numFields, 2);
4953:   *numFields = dm->Nf;
4954:   PetscFunctionReturn(PETSC_SUCCESS);
4955: }

4957: /*@
4958:   DMSetNumFields - Set the number of fields in the `DM`

4960:   Logically Collective

4962:   Input Parameters:
4963: + dm        - The `DM`
4964: - numFields - The number of fields

4966:   Level: intermediate

4968: .seealso: [](ch_dmbase), `DM`, `DMGetNumFields()`, `DMSetField()`
4969: @*/
4970: PetscErrorCode DMSetNumFields(DM dm, PetscInt numFields)
4971: {
4972:   PetscInt Nf, f;

4974:   PetscFunctionBegin;
4976:   PetscCall(DMGetNumFields(dm, &Nf));
4977:   for (f = Nf; f < numFields; ++f) {
4978:     PetscContainer obj;

4980:     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)dm), &obj));
4981:     PetscCall(DMAddField(dm, NULL, (PetscObject)obj));
4982:     PetscCall(PetscContainerDestroy(&obj));
4983:   }
4984:   PetscFunctionReturn(PETSC_SUCCESS);
4985: }

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

4990:   Not Collective

4992:   Input Parameters:
4993: + dm - The `DM`
4994: - f  - The field number

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

5000:   Level: intermediate

5002: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMSetField()`
5003: @*/
5004: PetscErrorCode DMGetField(DM dm, PetscInt f, DMLabel *label, PetscObject *disc)
5005: {
5006:   PetscFunctionBegin;
5008:   PetscAssertPointer(disc, 4);
5009:   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);
5010:   if (label) *label = dm->fields[f].label;
5011:   if (disc) *disc = dm->fields[f].disc;
5012:   PetscFunctionReturn(PETSC_SUCCESS);
5013: }

5015: /* Does not clear the DS */
5016: PetscErrorCode DMSetField_Internal(DM dm, PetscInt f, DMLabel label, PetscObject disc)
5017: {
5018:   PetscFunctionBegin;
5019:   PetscCall(DMFieldEnlarge_Static(dm, f + 1));
5020:   PetscCall(DMLabelDestroy(&dm->fields[f].label));
5021:   PetscCall(PetscObjectDestroy(&dm->fields[f].disc));
5022:   dm->fields[f].label = label;
5023:   dm->fields[f].disc  = disc;
5024:   PetscCall(PetscObjectReference((PetscObject)label));
5025:   PetscCall(PetscObjectReference((PetscObject)disc));
5026:   PetscFunctionReturn(PETSC_SUCCESS);
5027: }

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

5033:   Logically Collective

5035:   Input Parameters:
5036: + dm    - The `DM`
5037: . f     - The field number
5038: . label - The label indicating the support of the field, or `NULL` for the entire mesh
5039: - disc  - The discretization object

5041:   Level: intermediate

5043: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMGetField()`
5044: @*/
5045: PetscErrorCode DMSetField(DM dm, PetscInt f, DMLabel label, PetscObject disc)
5046: {
5047:   PetscFunctionBegin;
5051:   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
5052:   PetscCall(DMSetField_Internal(dm, f, label, disc));
5053:   PetscCall(DMSetDefaultAdjacency_Private(dm, f, disc));
5054:   PetscCall(DMClearDS(dm));
5055:   PetscFunctionReturn(PETSC_SUCCESS);
5056: }

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

5062:   Logically Collective

5064:   Input Parameters:
5065: + dm    - The `DM`
5066: . label - The label indicating the support of the field, or `NULL` for the entire mesh
5067: - disc  - The discretization object

5069:   Level: intermediate

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

5074:   For example, a piecewise continuous pressure field can be defined by coefficients at the cell centers of a mesh and piecewise constant functions
5075:   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
5076:   geometry entities, a `DMLabel` indicating a subset of those geometric entities, and a discretization object, such as a `PetscFE`.

5078: .seealso: [](ch_dmbase), `DM`, `DMSetLabel()`, `DMSetField()`, `DMGetField()`, `PetscFE`
5079: @*/
5080: PetscErrorCode DMAddField(DM dm, DMLabel label, PetscObject disc)
5081: {
5082:   PetscInt Nf = dm->Nf;

5084:   PetscFunctionBegin;
5088:   PetscCall(DMFieldEnlarge_Static(dm, Nf + 1));
5089:   dm->fields[Nf].label = label;
5090:   dm->fields[Nf].disc  = disc;
5091:   PetscCall(PetscObjectReference((PetscObject)label));
5092:   PetscCall(PetscObjectReference((PetscObject)disc));
5093:   PetscCall(DMSetDefaultAdjacency_Private(dm, Nf, disc));
5094:   PetscCall(DMClearDS(dm));
5095:   PetscFunctionReturn(PETSC_SUCCESS);
5096: }

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

5101:   Logically Collective

5103:   Input Parameters:
5104: + dm          - The `DM`
5105: . f           - The field index
5106: - avoidTensor - `PETSC_TRUE` to skip defining the field on tensor cells

5108:   Level: intermediate

5110: .seealso: [](ch_dmbase), `DM`, `DMGetFieldAvoidTensor()`, `DMSetField()`, `DMGetField()`
5111: @*/
5112: PetscErrorCode DMSetFieldAvoidTensor(DM dm, PetscInt f, PetscBool avoidTensor)
5113: {
5114:   PetscFunctionBegin;
5115:   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);
5116:   dm->fields[f].avoidTensor = avoidTensor;
5117:   PetscFunctionReturn(PETSC_SUCCESS);
5118: }

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

5123:   Not Collective

5125:   Input Parameters:
5126: + dm - The `DM`
5127: - f  - The field index

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

5132:   Level: intermediate

5134: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMSetField()`, `DMGetField()`, `DMSetFieldAvoidTensor()`
5135: @*/
5136: PetscErrorCode DMGetFieldAvoidTensor(DM dm, PetscInt f, PetscBool *avoidTensor)
5137: {
5138:   PetscFunctionBegin;
5139:   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);
5140:   *avoidTensor = dm->fields[f].avoidTensor;
5141:   PetscFunctionReturn(PETSC_SUCCESS);
5142: }

5144: /*@
5145:   DMCopyFields - Copy the discretizations for the `DM` into another `DM`

5147:   Collective

5149:   Input Parameter:
5150: . dm - The `DM`

5152:   Output Parameter:
5153: . newdm - The `DM`

5155:   Level: advanced

5157: .seealso: [](ch_dmbase), `DM`, `DMGetField()`, `DMSetField()`, `DMAddField()`, `DMCopyDS()`, `DMGetDS()`, `DMGetCellDS()`
5158: @*/
5159: PetscErrorCode DMCopyFields(DM dm, DM newdm)
5160: {
5161:   PetscInt Nf, f;

5163:   PetscFunctionBegin;
5164:   if (dm == newdm) PetscFunctionReturn(PETSC_SUCCESS);
5165:   PetscCall(DMGetNumFields(dm, &Nf));
5166:   PetscCall(DMClearFields(newdm));
5167:   for (f = 0; f < Nf; ++f) {
5168:     DMLabel     label;
5169:     PetscObject field;
5170:     PetscBool   useCone, useClosure;

5172:     PetscCall(DMGetField(dm, f, &label, &field));
5173:     PetscCall(DMSetField(newdm, f, label, field));
5174:     PetscCall(DMGetAdjacency(dm, f, &useCone, &useClosure));
5175:     PetscCall(DMSetAdjacency(newdm, f, useCone, useClosure));
5176:   }
5177:   PetscFunctionReturn(PETSC_SUCCESS);
5178: }

5180: /*@
5181:   DMGetAdjacency - Returns the flags for determining variable influence

5183:   Not Collective

5185:   Input Parameters:
5186: + dm - The `DM` object
5187: - f  - The field number, or `PETSC_DEFAULT` for the default adjacency

5189:   Output Parameters:
5190: + useCone    - Flag for variable influence starting with the cone operation
5191: - useClosure - Flag for variable influence using transitive closure

5193:   Level: developer

5195:   Notes:
5196: .vb
5197:      FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
5198:      FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
5199:      FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
5200: .ve
5201:   Further explanation can be found in the User's Manual Section on the Influence of Variables on One Another.

5203: .seealso: [](ch_dmbase), `DM`, `DMSetAdjacency()`, `DMGetField()`, `DMSetField()`
5204: @*/
5205: PetscErrorCode DMGetAdjacency(DM dm, PetscInt f, PetscBool *useCone, PetscBool *useClosure)
5206: {
5207:   PetscFunctionBegin;
5209:   if (useCone) PetscAssertPointer(useCone, 3);
5210:   if (useClosure) PetscAssertPointer(useClosure, 4);
5211:   if (f < 0) {
5212:     if (useCone) *useCone = dm->adjacency[0];
5213:     if (useClosure) *useClosure = dm->adjacency[1];
5214:   } else {
5215:     PetscInt Nf;

5217:     PetscCall(DMGetNumFields(dm, &Nf));
5218:     PetscCheck(f < Nf, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, Nf);
5219:     if (useCone) *useCone = dm->fields[f].adjacency[0];
5220:     if (useClosure) *useClosure = dm->fields[f].adjacency[1];
5221:   }
5222:   PetscFunctionReturn(PETSC_SUCCESS);
5223: }

5225: /*@
5226:   DMSetAdjacency - Set the flags for determining variable influence

5228:   Not Collective

5230:   Input Parameters:
5231: + dm         - The `DM` object
5232: . f          - The field number
5233: . useCone    - Flag for variable influence starting with the cone operation
5234: - useClosure - Flag for variable influence using transitive closure

5236:   Level: developer

5238:   Notes:
5239: .vb
5240:      FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
5241:      FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
5242:      FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
5243: .ve
5244:   Further explanation can be found in the User's Manual Section on the Influence of Variables on One Another.

5246: .seealso: [](ch_dmbase), `DM`, `DMGetAdjacency()`, `DMGetField()`, `DMSetField()`
5247: @*/
5248: PetscErrorCode DMSetAdjacency(DM dm, PetscInt f, PetscBool useCone, PetscBool useClosure)
5249: {
5250:   PetscFunctionBegin;
5252:   if (f < 0) {
5253:     dm->adjacency[0] = useCone;
5254:     dm->adjacency[1] = useClosure;
5255:   } else {
5256:     PetscInt Nf;

5258:     PetscCall(DMGetNumFields(dm, &Nf));
5259:     PetscCheck(f < Nf, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, Nf);
5260:     dm->fields[f].adjacency[0] = useCone;
5261:     dm->fields[f].adjacency[1] = useClosure;
5262:   }
5263:   PetscFunctionReturn(PETSC_SUCCESS);
5264: }

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

5269:   Not collective

5271:   Input Parameter:
5272: . dm - The `DM` object

5274:   Output Parameters:
5275: + useCone    - Flag for variable influence starting with the cone operation
5276: - useClosure - Flag for variable influence using transitive closure

5278:   Level: developer

5280:   Notes:
5281: .vb
5282:      FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
5283:      FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
5284:      FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
5285: .ve

5287: .seealso: [](ch_dmbase), `DM`, `DMSetBasicAdjacency()`, `DMGetField()`, `DMSetField()`
5288: @*/
5289: PetscErrorCode DMGetBasicAdjacency(DM dm, PetscBool *useCone, PetscBool *useClosure)
5290: {
5291:   PetscInt Nf;

5293:   PetscFunctionBegin;
5295:   if (useCone) PetscAssertPointer(useCone, 2);
5296:   if (useClosure) PetscAssertPointer(useClosure, 3);
5297:   PetscCall(DMGetNumFields(dm, &Nf));
5298:   if (!Nf) {
5299:     PetscCall(DMGetAdjacency(dm, PETSC_DEFAULT, useCone, useClosure));
5300:   } else {
5301:     PetscCall(DMGetAdjacency(dm, 0, useCone, useClosure));
5302:   }
5303:   PetscFunctionReturn(PETSC_SUCCESS);
5304: }

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

5309:   Not Collective

5311:   Input Parameters:
5312: + dm         - The `DM` object
5313: . useCone    - Flag for variable influence starting with the cone operation
5314: - useClosure - Flag for variable influence using transitive closure

5316:   Level: developer

5318:   Notes:
5319: .vb
5320:      FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
5321:      FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
5322:      FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
5323: .ve

5325: .seealso: [](ch_dmbase), `DM`, `DMGetBasicAdjacency()`, `DMGetField()`, `DMSetField()`
5326: @*/
5327: PetscErrorCode DMSetBasicAdjacency(DM dm, PetscBool useCone, PetscBool useClosure)
5328: {
5329:   PetscInt Nf;

5331:   PetscFunctionBegin;
5333:   PetscCall(DMGetNumFields(dm, &Nf));
5334:   if (!Nf) {
5335:     PetscCall(DMSetAdjacency(dm, PETSC_DEFAULT, useCone, useClosure));
5336:   } else {
5337:     PetscCall(DMSetAdjacency(dm, 0, useCone, useClosure));
5338:   }
5339:   PetscFunctionReturn(PETSC_SUCCESS);
5340: }

5342: PetscErrorCode DMCompleteBCLabels_Internal(DM dm)
5343: {
5344:   DM           plex;
5345:   DMLabel     *labels, *glabels;
5346:   const char **names;
5347:   char        *sendNames, *recvNames;
5348:   PetscInt     Nds, s, maxLabels = 0, maxLen = 0, gmaxLen, Nl = 0, gNl, l, gl, m;
5349:   size_t       len;
5350:   MPI_Comm     comm;
5351:   PetscMPIInt  rank, size, p, *counts, *displs;

5353:   PetscFunctionBegin;
5354:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
5355:   PetscCallMPI(MPI_Comm_size(comm, &size));
5356:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5357:   PetscCall(DMGetNumDS(dm, &Nds));
5358:   for (s = 0; s < Nds; ++s) {
5359:     PetscDS  dsBC;
5360:     PetscInt numBd;

5362:     PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC, NULL));
5363:     PetscCall(PetscDSGetNumBoundary(dsBC, &numBd));
5364:     maxLabels += numBd;
5365:   }
5366:   PetscCall(PetscCalloc1(maxLabels, &labels));
5367:   /* Get list of labels to be completed */
5368:   for (s = 0; s < Nds; ++s) {
5369:     PetscDS  dsBC;
5370:     PetscInt numBd, bd;

5372:     PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC, NULL));
5373:     PetscCall(PetscDSGetNumBoundary(dsBC, &numBd));
5374:     for (bd = 0; bd < numBd; ++bd) {
5375:       DMLabel      label;
5376:       PetscInt     field;
5377:       PetscObject  obj;
5378:       PetscClassId id;

5380:       PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, NULL, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL));
5381:       PetscCall(DMGetField(dm, field, NULL, &obj));
5382:       PetscCall(PetscObjectGetClassId(obj, &id));
5383:       if (!(id == PETSCFE_CLASSID) || !label) continue;
5384:       for (l = 0; l < Nl; ++l)
5385:         if (labels[l] == label) break;
5386:       if (l == Nl) labels[Nl++] = label;
5387:     }
5388:   }
5389:   /* Get label names */
5390:   PetscCall(PetscMalloc1(Nl, &names));
5391:   for (l = 0; l < Nl; ++l) PetscCall(PetscObjectGetName((PetscObject)labels[l], &names[l]));
5392:   for (l = 0; l < Nl; ++l) {
5393:     PetscCall(PetscStrlen(names[l], &len));
5394:     maxLen = PetscMax(maxLen, (PetscInt)len + 2);
5395:   }
5396:   PetscCall(PetscFree(labels));
5397:   PetscCall(MPIU_Allreduce(&maxLen, &gmaxLen, 1, MPIU_INT, MPI_MAX, comm));
5398:   PetscCall(PetscCalloc1(Nl * gmaxLen, &sendNames));
5399:   for (l = 0; l < Nl; ++l) PetscCall(PetscStrncpy(&sendNames[gmaxLen * l], names[l], gmaxLen));
5400:   PetscCall(PetscFree(names));
5401:   /* Put all names on all processes */
5402:   PetscCall(PetscCalloc2(size, &counts, size + 1, &displs));
5403:   PetscCallMPI(MPI_Allgather(&Nl, 1, MPI_INT, counts, 1, MPI_INT, comm));
5404:   for (p = 0; p < size; ++p) displs[p + 1] = displs[p] + counts[p];
5405:   gNl = displs[size];
5406:   for (p = 0; p < size; ++p) {
5407:     counts[p] *= gmaxLen;
5408:     displs[p] *= gmaxLen;
5409:   }
5410:   PetscCall(PetscCalloc2(gNl * gmaxLen, &recvNames, gNl, &glabels));
5411:   PetscCallMPI(MPI_Allgatherv(sendNames, counts[rank], MPI_CHAR, recvNames, counts, displs, MPI_CHAR, comm));
5412:   PetscCall(PetscFree2(counts, displs));
5413:   PetscCall(PetscFree(sendNames));
5414:   for (l = 0, gl = 0; l < gNl; ++l) {
5415:     PetscCall(DMGetLabel(dm, &recvNames[l * gmaxLen], &glabels[gl]));
5416:     PetscCheck(glabels[gl], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Label %s missing on rank %d", &recvNames[l * gmaxLen], rank);
5417:     for (m = 0; m < gl; ++m)
5418:       if (glabels[m] == glabels[gl]) continue;
5419:     PetscCall(DMConvert(dm, DMPLEX, &plex));
5420:     PetscCall(DMPlexLabelComplete(plex, glabels[gl]));
5421:     PetscCall(DMDestroy(&plex));
5422:     ++gl;
5423:   }
5424:   PetscCall(PetscFree2(recvNames, glabels));
5425:   PetscFunctionReturn(PETSC_SUCCESS);
5426: }

5428: static PetscErrorCode DMDSEnlarge_Static(DM dm, PetscInt NdsNew)
5429: {
5430:   DMSpace *tmpd;
5431:   PetscInt Nds = dm->Nds, s;

5433:   PetscFunctionBegin;
5434:   if (Nds >= NdsNew) PetscFunctionReturn(PETSC_SUCCESS);
5435:   PetscCall(PetscMalloc1(NdsNew, &tmpd));
5436:   for (s = 0; s < Nds; ++s) tmpd[s] = dm->probs[s];
5437:   for (s = Nds; s < NdsNew; ++s) {
5438:     tmpd[s].ds     = NULL;
5439:     tmpd[s].label  = NULL;
5440:     tmpd[s].fields = NULL;
5441:   }
5442:   PetscCall(PetscFree(dm->probs));
5443:   dm->Nds   = NdsNew;
5444:   dm->probs = tmpd;
5445:   PetscFunctionReturn(PETSC_SUCCESS);
5446: }

5448: /*@
5449:   DMGetNumDS - Get the number of discrete systems in the `DM`

5451:   Not Collective

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

5456:   Output Parameter:
5457: . Nds - The number of `PetscDS` objects

5459:   Level: intermediate

5461: .seealso: [](ch_dmbase), `DM`, `DMGetDS()`, `DMGetCellDS()`
5462: @*/
5463: PetscErrorCode DMGetNumDS(DM dm, PetscInt *Nds)
5464: {
5465:   PetscFunctionBegin;
5467:   PetscAssertPointer(Nds, 2);
5468:   *Nds = dm->Nds;
5469:   PetscFunctionReturn(PETSC_SUCCESS);
5470: }

5472: /*@
5473:   DMClearDS - Remove all discrete systems from the `DM`

5475:   Logically Collective

5477:   Input Parameter:
5478: . dm - The `DM`

5480:   Level: intermediate

5482: .seealso: [](ch_dmbase), `DM`, `DMGetNumDS()`, `DMGetDS()`, `DMSetField()`
5483: @*/
5484: PetscErrorCode DMClearDS(DM dm)
5485: {
5486:   PetscInt s;

5488:   PetscFunctionBegin;
5490:   for (s = 0; s < dm->Nds; ++s) {
5491:     PetscCall(PetscDSDestroy(&dm->probs[s].ds));
5492:     PetscCall(PetscDSDestroy(&dm->probs[s].dsIn));
5493:     PetscCall(DMLabelDestroy(&dm->probs[s].label));
5494:     PetscCall(ISDestroy(&dm->probs[s].fields));
5495:   }
5496:   PetscCall(PetscFree(dm->probs));
5497:   dm->probs = NULL;
5498:   dm->Nds   = 0;
5499:   PetscFunctionReturn(PETSC_SUCCESS);
5500: }

5502: /*@
5503:   DMGetDS - Get the default `PetscDS`

5505:   Not Collective

5507:   Input Parameter:
5508: . dm - The `DM`

5510:   Output Parameter:
5511: . ds - The default `PetscDS`

5513:   Level: intermediate

5515: .seealso: [](ch_dmbase), `DM`, `DMGetCellDS()`, `DMGetRegionDS()`
5516: @*/
5517: PetscErrorCode DMGetDS(DM dm, PetscDS *ds)
5518: {
5519:   PetscFunctionBeginHot;
5521:   PetscAssertPointer(ds, 2);
5522:   PetscCheck(dm->Nds > 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Need to call DMCreateDS() before calling DMGetDS()");
5523:   *ds = dm->probs[0].ds;
5524:   PetscFunctionReturn(PETSC_SUCCESS);
5525: }

5527: /*@
5528:   DMGetCellDS - Get the `PetscDS` defined on a given cell

5530:   Not Collective

5532:   Input Parameters:
5533: + dm    - The `DM`
5534: - point - Cell for the `PetscDS`

5536:   Output Parameters:
5537: + ds   - The `PetscDS` defined on the given cell
5538: - dsIn - The `PetscDS` for input on the given cell, or NULL if the same ds

5540:   Level: developer

5542: .seealso: [](ch_dmbase), `DM`, `DMGetDS()`, `DMSetRegionDS()`
5543: @*/
5544: PetscErrorCode DMGetCellDS(DM dm, PetscInt point, PetscDS *ds, PetscDS *dsIn)
5545: {
5546:   PetscDS  dsDef = NULL;
5547:   PetscInt s;

5549:   PetscFunctionBeginHot;
5551:   if (ds) PetscAssertPointer(ds, 3);
5552:   if (dsIn) PetscAssertPointer(dsIn, 4);
5553:   PetscCheck(point >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Mesh point cannot be negative: %" PetscInt_FMT, point);
5554:   if (ds) *ds = NULL;
5555:   if (dsIn) *dsIn = NULL;
5556:   for (s = 0; s < dm->Nds; ++s) {
5557:     PetscInt val;

5559:     if (!dm->probs[s].label) {
5560:       dsDef = dm->probs[s].ds;
5561:     } else {
5562:       PetscCall(DMLabelGetValue(dm->probs[s].label, point, &val));
5563:       if (val >= 0) {
5564:         if (ds) *ds = dm->probs[s].ds;
5565:         if (dsIn) *dsIn = dm->probs[s].dsIn;
5566:         break;
5567:       }
5568:     }
5569:   }
5570:   if (ds && !*ds) *ds = dsDef;
5571:   PetscFunctionReturn(PETSC_SUCCESS);
5572: }

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

5577:   Not Collective

5579:   Input Parameters:
5580: + dm    - The `DM`
5581: - label - The `DMLabel` defining the mesh region, or `NULL` for the entire mesh

5583:   Output Parameters:
5584: + fields - The `IS` containing the `DM` field numbers for the fields in this `PetscDS`, or `NULL`
5585: . ds     - The `PetscDS` defined on the given region, or `NULL`
5586: - dsIn   - The `PetscDS` for input in the given region, or `NULL`

5588:   Level: advanced

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

5595: .seealso: [](ch_dmbase), `DM`, `DMGetRegionNumDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5596: @*/
5597: PetscErrorCode DMGetRegionDS(DM dm, DMLabel label, IS *fields, PetscDS *ds, PetscDS *dsIn)
5598: {
5599:   PetscInt Nds = dm->Nds, s;

5601:   PetscFunctionBegin;
5604:   if (fields) {
5605:     PetscAssertPointer(fields, 3);
5606:     *fields = NULL;
5607:   }
5608:   if (ds) {
5609:     PetscAssertPointer(ds, 4);
5610:     *ds = NULL;
5611:   }
5612:   if (dsIn) {
5613:     PetscAssertPointer(dsIn, 5);
5614:     *dsIn = NULL;
5615:   }
5616:   for (s = 0; s < Nds; ++s) {
5617:     if (dm->probs[s].label == label || !dm->probs[s].label) {
5618:       if (fields) *fields = dm->probs[s].fields;
5619:       if (ds) *ds = dm->probs[s].ds;
5620:       if (dsIn) *dsIn = dm->probs[s].dsIn;
5621:       if (dm->probs[s].label) PetscFunctionReturn(PETSC_SUCCESS);
5622:     }
5623:   }
5624:   PetscFunctionReturn(PETSC_SUCCESS);
5625: }

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

5630:   Collective

5632:   Input Parameters:
5633: + dm     - The `DM`
5634: . label  - The `DMLabel` defining the mesh region, or `NULL` for the entire mesh
5635: . fields - The `IS` containing the `DM` field numbers for the fields in this `PetscDS`, or `NULL` for all fields
5636: . ds     - The `PetscDS` defined on the given region
5637: - dsIn   - The `PetscDS` for input on the given cell, or `NULL` if it is the same `PetscDS`

5639:   Level: advanced

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

5645: .seealso: [](ch_dmbase), `DM`, `DMGetRegionDS()`, `DMSetRegionNumDS()`, `DMGetDS()`, `DMGetCellDS()`
5646: @*/
5647: PetscErrorCode DMSetRegionDS(DM dm, DMLabel label, IS fields, PetscDS ds, PetscDS dsIn)
5648: {
5649:   PetscInt Nds = dm->Nds, s;

5651:   PetscFunctionBegin;
5657:   for (s = 0; s < Nds; ++s) {
5658:     if (dm->probs[s].label == label) {
5659:       PetscCall(PetscDSDestroy(&dm->probs[s].ds));
5660:       PetscCall(PetscDSDestroy(&dm->probs[s].dsIn));
5661:       dm->probs[s].ds   = ds;
5662:       dm->probs[s].dsIn = dsIn;
5663:       PetscFunctionReturn(PETSC_SUCCESS);
5664:     }
5665:   }
5666:   PetscCall(DMDSEnlarge_Static(dm, Nds + 1));
5667:   PetscCall(PetscObjectReference((PetscObject)label));
5668:   PetscCall(PetscObjectReference((PetscObject)fields));
5669:   PetscCall(PetscObjectReference((PetscObject)ds));
5670:   PetscCall(PetscObjectReference((PetscObject)dsIn));
5671:   if (!label) {
5672:     /* Put the NULL label at the front, so it is returned as the default */
5673:     for (s = Nds - 1; s >= 0; --s) dm->probs[s + 1] = dm->probs[s];
5674:     Nds = 0;
5675:   }
5676:   dm->probs[Nds].label  = label;
5677:   dm->probs[Nds].fields = fields;
5678:   dm->probs[Nds].ds     = ds;
5679:   dm->probs[Nds].dsIn   = dsIn;
5680:   PetscFunctionReturn(PETSC_SUCCESS);
5681: }

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

5686:   Not Collective

5688:   Input Parameters:
5689: + dm  - The `DM`
5690: - num - The region number, in [0, Nds)

5692:   Output Parameters:
5693: + label  - The region label, or `NULL`
5694: . fields - The `IS` containing the `DM` field numbers for the fields in this `PetscDS`, or `NULL`
5695: . ds     - The `PetscDS` defined on the given region, or `NULL`
5696: - dsIn   - The `PetscDS` for input in the given region, or `NULL`

5698:   Level: advanced

5700: .seealso: [](ch_dmbase), `DM`, `DMGetRegionDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5701: @*/
5702: PetscErrorCode DMGetRegionNumDS(DM dm, PetscInt num, DMLabel *label, IS *fields, PetscDS *ds, PetscDS *dsIn)
5703: {
5704:   PetscInt Nds;

5706:   PetscFunctionBegin;
5708:   PetscCall(DMGetNumDS(dm, &Nds));
5709:   PetscCheck((num >= 0) && (num < Nds), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Region number %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", num, Nds);
5710:   if (label) {
5711:     PetscAssertPointer(label, 3);
5712:     *label = dm->probs[num].label;
5713:   }
5714:   if (fields) {
5715:     PetscAssertPointer(fields, 4);
5716:     *fields = dm->probs[num].fields;
5717:   }
5718:   if (ds) {
5719:     PetscAssertPointer(ds, 5);
5720:     *ds = dm->probs[num].ds;
5721:   }
5722:   if (dsIn) {
5723:     PetscAssertPointer(dsIn, 6);
5724:     *dsIn = dm->probs[num].dsIn;
5725:   }
5726:   PetscFunctionReturn(PETSC_SUCCESS);
5727: }

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

5732:   Not Collective

5734:   Input Parameters:
5735: + dm     - The `DM`
5736: . num    - The region number, in [0, Nds)
5737: . label  - The region label, or `NULL`
5738: . fields - The `IS` containing the `DM` field numbers for the fields in this `PetscDS`, or `NULL` to prevent setting
5739: . ds     - The `PetscDS` defined on the given region, or `NULL` to prevent setting
5740: - dsIn   - The `PetscDS` for input on the given cell, or `NULL` if it is the same `PetscDS`

5742:   Level: advanced

5744: .seealso: [](ch_dmbase), `DM`, `DMGetRegionDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5745: @*/
5746: PetscErrorCode DMSetRegionNumDS(DM dm, PetscInt num, DMLabel label, IS fields, PetscDS ds, PetscDS dsIn)
5747: {
5748:   PetscInt Nds;

5750:   PetscFunctionBegin;
5753:   PetscCall(DMGetNumDS(dm, &Nds));
5754:   PetscCheck((num >= 0) && (num < Nds), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Region number %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", num, Nds);
5755:   PetscCall(PetscObjectReference((PetscObject)label));
5756:   PetscCall(DMLabelDestroy(&dm->probs[num].label));
5757:   dm->probs[num].label = label;
5758:   if (fields) {
5760:     PetscCall(PetscObjectReference((PetscObject)fields));
5761:     PetscCall(ISDestroy(&dm->probs[num].fields));
5762:     dm->probs[num].fields = fields;
5763:   }
5764:   if (ds) {
5766:     PetscCall(PetscObjectReference((PetscObject)ds));
5767:     PetscCall(PetscDSDestroy(&dm->probs[num].ds));
5768:     dm->probs[num].ds = ds;
5769:   }
5770:   if (dsIn) {
5772:     PetscCall(PetscObjectReference((PetscObject)dsIn));
5773:     PetscCall(PetscDSDestroy(&dm->probs[num].dsIn));
5774:     dm->probs[num].dsIn = dsIn;
5775:   }
5776:   PetscFunctionReturn(PETSC_SUCCESS);
5777: }

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

5782:   Not Collective

5784:   Input Parameters:
5785: + dm - The `DM`
5786: - ds - The `PetscDS` defined on the given region

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

5791:   Level: advanced

5793: .seealso: [](ch_dmbase), `DM`, `DMGetRegionNumDS()`, `DMGetRegionDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5794: @*/
5795: PetscErrorCode DMFindRegionNum(DM dm, PetscDS ds, PetscInt *num)
5796: {
5797:   PetscInt Nds, n;

5799:   PetscFunctionBegin;
5802:   PetscAssertPointer(num, 3);
5803:   PetscCall(DMGetNumDS(dm, &Nds));
5804:   for (n = 0; n < Nds; ++n)
5805:     if (ds == dm->probs[n].ds) break;
5806:   if (n >= Nds) *num = -1;
5807:   else *num = n;
5808:   PetscFunctionReturn(PETSC_SUCCESS);
5809: }

5811: /*@
5812:   DMCreateFEDefault - Create a `PetscFE` based on the celltype for the mesh

5814:   Not Collective

5816:   Input Parameters:
5817: + dm     - The `DM`
5818: . Nc     - The number of components for the field
5819: . prefix - The options prefix for the output `PetscFE`, or `NULL`
5820: - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree

5822:   Output Parameter:
5823: . fem - The `PetscFE`

5825:   Level: intermediate

5827:   Note:
5828:   This is a convenience method that just calls `PetscFECreateByCell()` underneath.

5830: .seealso: [](ch_dmbase), `DM`, `PetscFECreateByCell()`, `DMAddField()`, `DMCreateDS()`, `DMGetCellDS()`, `DMGetRegionDS()`
5831: @*/
5832: PetscErrorCode DMCreateFEDefault(DM dm, PetscInt Nc, const char prefix[], PetscInt qorder, PetscFE *fem)
5833: {
5834:   DMPolytopeType ct;
5835:   PetscInt       dim, cStart;

5837:   PetscFunctionBegin;
5840:   if (prefix) PetscAssertPointer(prefix, 3);
5842:   PetscAssertPointer(fem, 5);
5843:   PetscCall(DMGetDimension(dm, &dim));
5844:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
5845:   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
5846:   PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, Nc, ct, prefix, qorder, fem));
5847:   PetscFunctionReturn(PETSC_SUCCESS);
5848: }

5850: /*@
5851:   DMCreateDS - Create the discrete systems for the `DM` based upon the fields added to the `DM`

5853:   Collective

5855:   Input Parameter:
5856: . dm - The `DM`

5858:   Options Database Key:
5859: . -dm_petscds_view - View all the `PetscDS` objects in this `DM`

5861:   Level: intermediate

5863: .seealso: [](ch_dmbase), `DM`, `DMSetField`, `DMAddField()`, `DMGetDS()`, `DMGetCellDS()`, `DMGetRegionDS()`, `DMSetRegionDS()`
5864: @*/
5865: PetscErrorCode DMCreateDS(DM dm)
5866: {
5867:   MPI_Comm  comm;
5868:   PetscDS   dsDef;
5869:   DMLabel  *labelSet;
5870:   PetscInt  dE, Nf = dm->Nf, f, s, Nl, l, Ndef, k;
5871:   PetscBool doSetup = PETSC_TRUE, flg;

5873:   PetscFunctionBegin;
5875:   if (!dm->fields) PetscFunctionReturn(PETSC_SUCCESS);
5876:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
5877:   PetscCall(DMGetCoordinateDim(dm, &dE));
5878:   /* Determine how many regions we have */
5879:   PetscCall(PetscMalloc1(Nf, &labelSet));
5880:   Nl   = 0;
5881:   Ndef = 0;
5882:   for (f = 0; f < Nf; ++f) {
5883:     DMLabel  label = dm->fields[f].label;
5884:     PetscInt l;

5886: #ifdef PETSC_HAVE_LIBCEED
5887:     /* Move CEED context to discretizations */
5888:     {
5889:       PetscClassId id;

5891:       PetscCall(PetscObjectGetClassId(dm->fields[f].disc, &id));
5892:       if (id == PETSCFE_CLASSID) {
5893:         Ceed ceed;

5895:         PetscCall(DMGetCeed(dm, &ceed));
5896:         PetscCall(PetscFESetCeed((PetscFE)dm->fields[f].disc, ceed));
5897:       }
5898:     }
5899: #endif
5900:     if (!label) {
5901:       ++Ndef;
5902:       continue;
5903:     }
5904:     for (l = 0; l < Nl; ++l)
5905:       if (label == labelSet[l]) break;
5906:     if (l < Nl) continue;
5907:     labelSet[Nl++] = label;
5908:   }
5909:   /* Create default DS if there are no labels to intersect with */
5910:   PetscCall(DMGetRegionDS(dm, NULL, NULL, &dsDef, NULL));
5911:   if (!dsDef && Ndef && !Nl) {
5912:     IS        fields;
5913:     PetscInt *fld, nf;

5915:     for (f = 0, nf = 0; f < Nf; ++f)
5916:       if (!dm->fields[f].label) ++nf;
5917:     PetscCheck(nf, comm, PETSC_ERR_PLIB, "All fields have labels, but we are trying to create a default DS");
5918:     PetscCall(PetscMalloc1(nf, &fld));
5919:     for (f = 0, nf = 0; f < Nf; ++f)
5920:       if (!dm->fields[f].label) fld[nf++] = f;
5921:     PetscCall(ISCreate(PETSC_COMM_SELF, &fields));
5922:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fields, "dm_fields_"));
5923:     PetscCall(ISSetType(fields, ISGENERAL));
5924:     PetscCall(ISGeneralSetIndices(fields, nf, fld, PETSC_OWN_POINTER));

5926:     PetscCall(PetscDSCreate(PETSC_COMM_SELF, &dsDef));
5927:     PetscCall(DMSetRegionDS(dm, NULL, fields, dsDef, NULL));
5928:     PetscCall(PetscDSDestroy(&dsDef));
5929:     PetscCall(ISDestroy(&fields));
5930:   }
5931:   PetscCall(DMGetRegionDS(dm, NULL, NULL, &dsDef, NULL));
5932:   if (dsDef) PetscCall(PetscDSSetCoordinateDimension(dsDef, dE));
5933:   /* Intersect labels with default fields */
5934:   if (Ndef && Nl) {
5935:     DM              plex;
5936:     DMLabel         cellLabel;
5937:     IS              fieldIS, allcellIS, defcellIS = NULL;
5938:     PetscInt       *fields;
5939:     const PetscInt *cells;
5940:     PetscInt        depth, nf = 0, n, c;

5942:     PetscCall(DMConvert(dm, DMPLEX, &plex));
5943:     PetscCall(DMPlexGetDepth(plex, &depth));
5944:     PetscCall(DMGetStratumIS(plex, "dim", depth, &allcellIS));
5945:     if (!allcellIS) PetscCall(DMGetStratumIS(plex, "depth", depth, &allcellIS));
5946:     /* TODO This looks like it only works for one label */
5947:     for (l = 0; l < Nl; ++l) {
5948:       DMLabel label = labelSet[l];
5949:       IS      pointIS;

5951:       PetscCall(ISDestroy(&defcellIS));
5952:       PetscCall(DMLabelGetStratumIS(label, 1, &pointIS));
5953:       PetscCall(ISDifference(allcellIS, pointIS, &defcellIS));
5954:       PetscCall(ISDestroy(&pointIS));
5955:     }
5956:     PetscCall(ISDestroy(&allcellIS));

5958:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "defaultCells", &cellLabel));
5959:     PetscCall(ISGetLocalSize(defcellIS, &n));
5960:     PetscCall(ISGetIndices(defcellIS, &cells));
5961:     for (c = 0; c < n; ++c) PetscCall(DMLabelSetValue(cellLabel, cells[c], 1));
5962:     PetscCall(ISRestoreIndices(defcellIS, &cells));
5963:     PetscCall(ISDestroy(&defcellIS));
5964:     PetscCall(DMPlexLabelComplete(plex, cellLabel));

5966:     PetscCall(PetscMalloc1(Ndef, &fields));
5967:     for (f = 0; f < Nf; ++f)
5968:       if (!dm->fields[f].label) fields[nf++] = f;
5969:     PetscCall(ISCreate(PETSC_COMM_SELF, &fieldIS));
5970:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fieldIS, "dm_fields_"));
5971:     PetscCall(ISSetType(fieldIS, ISGENERAL));
5972:     PetscCall(ISGeneralSetIndices(fieldIS, nf, fields, PETSC_OWN_POINTER));

5974:     PetscCall(PetscDSCreate(PETSC_COMM_SELF, &dsDef));
5975:     PetscCall(DMSetRegionDS(dm, cellLabel, fieldIS, dsDef, NULL));
5976:     PetscCall(PetscDSSetCoordinateDimension(dsDef, dE));
5977:     PetscCall(DMLabelDestroy(&cellLabel));
5978:     PetscCall(PetscDSDestroy(&dsDef));
5979:     PetscCall(ISDestroy(&fieldIS));
5980:     PetscCall(DMDestroy(&plex));
5981:   }
5982:   /* Create label DSes
5983:      - WE ONLY SUPPORT IDENTICAL OR DISJOINT LABELS
5984:   */
5985:   /* TODO Should check that labels are disjoint */
5986:   for (l = 0; l < Nl; ++l) {
5987:     DMLabel   label = labelSet[l];
5988:     PetscDS   ds, dsIn = NULL;
5989:     IS        fields;
5990:     PetscInt *fld, nf;

5992:     PetscCall(PetscDSCreate(PETSC_COMM_SELF, &ds));
5993:     for (f = 0, nf = 0; f < Nf; ++f)
5994:       if (label == dm->fields[f].label || !dm->fields[f].label) ++nf;
5995:     PetscCall(PetscMalloc1(nf, &fld));
5996:     for (f = 0, nf = 0; f < Nf; ++f)
5997:       if (label == dm->fields[f].label || !dm->fields[f].label) fld[nf++] = f;
5998:     PetscCall(ISCreate(PETSC_COMM_SELF, &fields));
5999:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fields, "dm_fields_"));
6000:     PetscCall(ISSetType(fields, ISGENERAL));
6001:     PetscCall(ISGeneralSetIndices(fields, nf, fld, PETSC_OWN_POINTER));
6002:     PetscCall(PetscDSSetCoordinateDimension(ds, dE));
6003:     {
6004:       DMPolytopeType ct;
6005:       PetscInt       lStart, lEnd;
6006:       PetscBool      isCohesiveLocal = PETSC_FALSE, isCohesive;

6008:       PetscCall(DMLabelGetBounds(label, &lStart, &lEnd));
6009:       if (lStart >= 0) {
6010:         PetscCall(DMPlexGetCellType(dm, lStart, &ct));
6011:         switch (ct) {
6012:         case DM_POLYTOPE_POINT_PRISM_TENSOR:
6013:         case DM_POLYTOPE_SEG_PRISM_TENSOR:
6014:         case DM_POLYTOPE_TRI_PRISM_TENSOR:
6015:         case DM_POLYTOPE_QUAD_PRISM_TENSOR:
6016:           isCohesiveLocal = PETSC_TRUE;
6017:           break;
6018:         default:
6019:           break;
6020:         }
6021:       }
6022:       PetscCall(MPIU_Allreduce(&isCohesiveLocal, &isCohesive, 1, MPIU_BOOL, MPI_LOR, comm));
6023:       if (isCohesive) {
6024:         PetscCall(PetscDSCreate(PETSC_COMM_SELF, &dsIn));
6025:         PetscCall(PetscDSSetCoordinateDimension(dsIn, dE));
6026:       }
6027:       for (f = 0, nf = 0; f < Nf; ++f) {
6028:         if (label == dm->fields[f].label || !dm->fields[f].label) {
6029:           if (label == dm->fields[f].label) {
6030:             PetscCall(PetscDSSetDiscretization(ds, nf, NULL));
6031:             PetscCall(PetscDSSetCohesive(ds, nf, isCohesive));
6032:             if (dsIn) {
6033:               PetscCall(PetscDSSetDiscretization(dsIn, nf, NULL));
6034:               PetscCall(PetscDSSetCohesive(dsIn, nf, isCohesive));
6035:             }
6036:           }
6037:           ++nf;
6038:         }
6039:       }
6040:     }
6041:     PetscCall(DMSetRegionDS(dm, label, fields, ds, dsIn));
6042:     PetscCall(ISDestroy(&fields));
6043:     PetscCall(PetscDSDestroy(&ds));
6044:     PetscCall(PetscDSDestroy(&dsIn));
6045:   }
6046:   PetscCall(PetscFree(labelSet));
6047:   /* Set fields in DSes */
6048:   for (s = 0; s < dm->Nds; ++s) {
6049:     PetscDS         ds     = dm->probs[s].ds;
6050:     PetscDS         dsIn   = dm->probs[s].dsIn;
6051:     IS              fields = dm->probs[s].fields;
6052:     const PetscInt *fld;
6053:     PetscInt        nf, dsnf;
6054:     PetscBool       isCohesive;

6056:     PetscCall(PetscDSGetNumFields(ds, &dsnf));
6057:     PetscCall(PetscDSIsCohesive(ds, &isCohesive));
6058:     PetscCall(ISGetLocalSize(fields, &nf));
6059:     PetscCall(ISGetIndices(fields, &fld));
6060:     for (f = 0; f < nf; ++f) {
6061:       PetscObject  disc = dm->fields[fld[f]].disc;
6062:       PetscBool    isCohesiveField;
6063:       PetscClassId id;

6065:       /* Handle DS with no fields */
6066:       if (dsnf) PetscCall(PetscDSGetCohesive(ds, f, &isCohesiveField));
6067:       /* If this is a cohesive cell, then regular fields need the lower dimensional discretization */
6068:       if (isCohesive) {
6069:         if (!isCohesiveField) {
6070:           PetscObject bdDisc;

6072:           PetscCall(PetscFEGetHeightSubspace((PetscFE)disc, 1, (PetscFE *)&bdDisc));
6073:           PetscCall(PetscDSSetDiscretization(ds, f, bdDisc));
6074:           PetscCall(PetscDSSetDiscretization(dsIn, f, disc));
6075:         } else {
6076:           PetscCall(PetscDSSetDiscretization(ds, f, disc));
6077:           PetscCall(PetscDSSetDiscretization(dsIn, f, disc));
6078:         }
6079:       } else {
6080:         PetscCall(PetscDSSetDiscretization(ds, f, disc));
6081:       }
6082:       /* We allow people to have placeholder fields and construct the Section by hand */
6083:       PetscCall(PetscObjectGetClassId(disc, &id));
6084:       if ((id != PETSCFE_CLASSID) && (id != PETSCFV_CLASSID)) doSetup = PETSC_FALSE;
6085:     }
6086:     PetscCall(ISRestoreIndices(fields, &fld));
6087:   }
6088:   /* Allow k-jet tabulation */
6089:   PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)dm)->prefix, "-dm_ds_jet_degree", &k, &flg));
6090:   if (flg) {
6091:     for (s = 0; s < dm->Nds; ++s) {
6092:       PetscDS  ds   = dm->probs[s].ds;
6093:       PetscDS  dsIn = dm->probs[s].dsIn;
6094:       PetscInt Nf, f;

6096:       PetscCall(PetscDSGetNumFields(ds, &Nf));
6097:       for (f = 0; f < Nf; ++f) {
6098:         PetscCall(PetscDSSetJetDegree(ds, f, k));
6099:         if (dsIn) PetscCall(PetscDSSetJetDegree(dsIn, f, k));
6100:       }
6101:     }
6102:   }
6103:   /* Setup DSes */
6104:   if (doSetup) {
6105:     for (s = 0; s < dm->Nds; ++s) {
6106:       if (dm->setfromoptionscalled) {
6107:         PetscCall(PetscDSSetFromOptions(dm->probs[s].ds));
6108:         if (dm->probs[s].dsIn) PetscCall(PetscDSSetFromOptions(dm->probs[s].dsIn));
6109:       }
6110:       PetscCall(PetscDSSetUp(dm->probs[s].ds));
6111:       if (dm->probs[s].dsIn) PetscCall(PetscDSSetUp(dm->probs[s].dsIn));
6112:     }
6113:   }
6114:   PetscFunctionReturn(PETSC_SUCCESS);
6115: }

6117: /*@
6118:   DMUseTensorOrder - Use a tensor product closure ordering for the default section

6120:   Input Parameters:
6121: + dm     - The DM
6122: - tensor - Flag for tensor order

6124:   Level: developer

6126: .seealso: `DMPlexSetClosurePermutationTensor()`, `PetscSectionResetClosurePermutation()`
6127: @*/
6128: PetscErrorCode DMUseTensorOrder(DM dm, PetscBool tensor)
6129: {
6130:   PetscInt  Nf;
6131:   PetscBool reorder = PETSC_TRUE, isPlex;

6133:   PetscFunctionBegin;
6134:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex));
6135:   PetscCall(DMGetNumFields(dm, &Nf));
6136:   for (PetscInt f = 0; f < Nf; ++f) {
6137:     PetscObject  obj;
6138:     PetscClassId id;

6140:     PetscCall(DMGetField(dm, f, NULL, &obj));
6141:     PetscCall(PetscObjectGetClassId(obj, &id));
6142:     if (id == PETSCFE_CLASSID) {
6143:       PetscSpace sp;
6144:       PetscBool  tensor;

6146:       PetscCall(PetscFEGetBasisSpace((PetscFE)obj, &sp));
6147:       PetscCall(PetscSpacePolynomialGetTensor(sp, &tensor));
6148:       reorder = reorder && tensor ? PETSC_TRUE : PETSC_FALSE;
6149:     } else reorder = PETSC_FALSE;
6150:   }
6151:   if (tensor) {
6152:     if (reorder && isPlex) PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
6153:   } else {
6154:     PetscSection s;

6156:     PetscCall(DMGetLocalSection(dm, &s));
6157:     if (s) PetscCall(PetscSectionResetClosurePermutation(s));
6158:   }
6159:   PetscFunctionReturn(PETSC_SUCCESS);
6160: }

6162: /*@
6163:   DMComputeExactSolution - Compute the exact solution for a given `DM`, using the `PetscDS` information.

6165:   Collective

6167:   Input Parameters:
6168: + dm   - The `DM`
6169: - time - The time

6171:   Output Parameters:
6172: + u   - The vector will be filled with exact solution values, or `NULL`
6173: - u_t - The vector will be filled with the time derivative of exact solution values, or `NULL`

6175:   Level: developer

6177:   Note:
6178:   The user must call `PetscDSSetExactSolution()` before using this routine

6180: .seealso: [](ch_dmbase), `DM`, `PetscDSSetExactSolution()`
6181: @*/
6182: PetscErrorCode DMComputeExactSolution(DM dm, PetscReal time, Vec u, Vec u_t)
6183: {
6184:   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
6185:   void   **ectxs;
6186:   Vec      locu, locu_t;
6187:   PetscInt Nf, Nds, s;

6189:   PetscFunctionBegin;
6191:   if (u) {
6193:     PetscCall(DMGetLocalVector(dm, &locu));
6194:     PetscCall(VecSet(locu, 0.));
6195:   }
6196:   if (u_t) {
6198:     PetscCall(DMGetLocalVector(dm, &locu_t));
6199:     PetscCall(VecSet(locu_t, 0.));
6200:   }
6201:   PetscCall(DMGetNumFields(dm, &Nf));
6202:   PetscCall(PetscMalloc2(Nf, &exacts, Nf, &ectxs));
6203:   PetscCall(DMGetNumDS(dm, &Nds));
6204:   for (s = 0; s < Nds; ++s) {
6205:     PetscDS         ds;
6206:     DMLabel         label;
6207:     IS              fieldIS;
6208:     const PetscInt *fields, id = 1;
6209:     PetscInt        dsNf, f;

6211:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL));
6212:     PetscCall(PetscDSGetNumFields(ds, &dsNf));
6213:     PetscCall(ISGetIndices(fieldIS, &fields));
6214:     PetscCall(PetscArrayzero(exacts, Nf));
6215:     PetscCall(PetscArrayzero(ectxs, Nf));
6216:     if (u) {
6217:       for (f = 0; f < dsNf; ++f) PetscCall(PetscDSGetExactSolution(ds, fields[f], &exacts[fields[f]], &ectxs[fields[f]]));
6218:       if (label) PetscCall(DMProjectFunctionLabelLocal(dm, time, label, 1, &id, 0, NULL, exacts, ectxs, INSERT_ALL_VALUES, locu));
6219:       else PetscCall(DMProjectFunctionLocal(dm, time, exacts, ectxs, INSERT_ALL_VALUES, locu));
6220:     }
6221:     if (u_t) {
6222:       PetscCall(PetscArrayzero(exacts, Nf));
6223:       PetscCall(PetscArrayzero(ectxs, Nf));
6224:       for (f = 0; f < dsNf; ++f) PetscCall(PetscDSGetExactSolutionTimeDerivative(ds, fields[f], &exacts[fields[f]], &ectxs[fields[f]]));
6225:       if (label) PetscCall(DMProjectFunctionLabelLocal(dm, time, label, 1, &id, 0, NULL, exacts, ectxs, INSERT_ALL_VALUES, locu_t));
6226:       else PetscCall(DMProjectFunctionLocal(dm, time, exacts, ectxs, INSERT_ALL_VALUES, locu_t));
6227:     }
6228:     PetscCall(ISRestoreIndices(fieldIS, &fields));
6229:   }
6230:   if (u) {
6231:     PetscCall(PetscObjectSetName((PetscObject)u, "Exact Solution"));
6232:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)u, "exact_"));
6233:   }
6234:   if (u_t) {
6235:     PetscCall(PetscObjectSetName((PetscObject)u, "Exact Solution Time Derivative"));
6236:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)u_t, "exact_t_"));
6237:   }
6238:   PetscCall(PetscFree2(exacts, ectxs));
6239:   if (u) {
6240:     PetscCall(DMLocalToGlobalBegin(dm, locu, INSERT_ALL_VALUES, u));
6241:     PetscCall(DMLocalToGlobalEnd(dm, locu, INSERT_ALL_VALUES, u));
6242:     PetscCall(DMRestoreLocalVector(dm, &locu));
6243:   }
6244:   if (u_t) {
6245:     PetscCall(DMLocalToGlobalBegin(dm, locu_t, INSERT_ALL_VALUES, u_t));
6246:     PetscCall(DMLocalToGlobalEnd(dm, locu_t, INSERT_ALL_VALUES, u_t));
6247:     PetscCall(DMRestoreLocalVector(dm, &locu_t));
6248:   }
6249:   PetscFunctionReturn(PETSC_SUCCESS);
6250: }

6252: static PetscErrorCode DMTransferDS_Internal(DM dm, DMLabel label, IS fields, PetscDS ds, PetscDS dsIn)
6253: {
6254:   PetscDS dsNew, dsInNew = NULL;

6256:   PetscFunctionBegin;
6257:   PetscCall(PetscDSCreate(PetscObjectComm((PetscObject)ds), &dsNew));
6258:   PetscCall(PetscDSCopy(ds, dm, dsNew));
6259:   if (dsIn) {
6260:     PetscCall(PetscDSCreate(PetscObjectComm((PetscObject)dsIn), &dsInNew));
6261:     PetscCall(PetscDSCopy(dsIn, dm, dsInNew));
6262:   }
6263:   PetscCall(DMSetRegionDS(dm, label, fields, dsNew, dsInNew));
6264:   PetscCall(PetscDSDestroy(&dsNew));
6265:   PetscCall(PetscDSDestroy(&dsInNew));
6266:   PetscFunctionReturn(PETSC_SUCCESS);
6267: }

6269: /*@
6270:   DMCopyDS - Copy the discrete systems for the `DM` into another `DM`

6272:   Collective

6274:   Input Parameter:
6275: . dm - The `DM`

6277:   Output Parameter:
6278: . newdm - The `DM`

6280:   Level: advanced

6282: .seealso: [](ch_dmbase), `DM`, `DMCopyFields()`, `DMAddField()`, `DMGetDS()`, `DMGetCellDS()`, `DMGetRegionDS()`, `DMSetRegionDS()`
6283: @*/
6284: PetscErrorCode DMCopyDS(DM dm, DM newdm)
6285: {
6286:   PetscInt Nds, s;

6288:   PetscFunctionBegin;
6289:   if (dm == newdm) PetscFunctionReturn(PETSC_SUCCESS);
6290:   PetscCall(DMGetNumDS(dm, &Nds));
6291:   PetscCall(DMClearDS(newdm));
6292:   for (s = 0; s < Nds; ++s) {
6293:     DMLabel  label;
6294:     IS       fields;
6295:     PetscDS  ds, dsIn, newds;
6296:     PetscInt Nbd, bd;

6298:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fields, &ds, &dsIn));
6299:     /* TODO: We need to change all keys from labels in the old DM to labels in the new DM */
6300:     PetscCall(DMTransferDS_Internal(newdm, label, fields, ds, dsIn));
6301:     /* Complete new labels in the new DS */
6302:     PetscCall(DMGetRegionDS(newdm, label, NULL, &newds, NULL));
6303:     PetscCall(PetscDSGetNumBoundary(newds, &Nbd));
6304:     for (bd = 0; bd < Nbd; ++bd) {
6305:       PetscWeakForm wf;
6306:       DMLabel       label;
6307:       PetscInt      field;

6309:       PetscCall(PetscDSGetBoundary(newds, bd, &wf, NULL, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL));
6310:       PetscCall(PetscWeakFormReplaceLabel(wf, label));
6311:     }
6312:   }
6313:   PetscCall(DMCompleteBCLabels_Internal(newdm));
6314:   PetscFunctionReturn(PETSC_SUCCESS);
6315: }

6317: /*@
6318:   DMCopyDisc - Copy the fields and discrete systems for the `DM` into another `DM`

6320:   Collective

6322:   Input Parameter:
6323: . dm - The `DM`

6325:   Output Parameter:
6326: . newdm - The `DM`

6328:   Level: advanced

6330:   Developer Note:
6331:   Really ugly name, nothing in PETSc is called a `Disc` plus it is an ugly abbreviation

6333: .seealso: [](ch_dmbase), `DM`, `DMCopyFields()`, `DMCopyDS()`
6334: @*/
6335: PetscErrorCode DMCopyDisc(DM dm, DM newdm)
6336: {
6337:   PetscFunctionBegin;
6338:   PetscCall(DMCopyFields(dm, newdm));
6339:   PetscCall(DMCopyDS(dm, newdm));
6340:   PetscFunctionReturn(PETSC_SUCCESS);
6341: }

6343: /*@
6344:   DMGetDimension - Return the topological dimension of the `DM`

6346:   Not Collective

6348:   Input Parameter:
6349: . dm - The `DM`

6351:   Output Parameter:
6352: . dim - The topological dimension

6354:   Level: beginner

6356: .seealso: [](ch_dmbase), `DM`, `DMSetDimension()`, `DMCreate()`
6357: @*/
6358: PetscErrorCode DMGetDimension(DM dm, PetscInt *dim)
6359: {
6360:   PetscFunctionBegin;
6362:   PetscAssertPointer(dim, 2);
6363:   *dim = dm->dim;
6364:   PetscFunctionReturn(PETSC_SUCCESS);
6365: }

6367: /*@
6368:   DMSetDimension - Set the topological dimension of the `DM`

6370:   Collective

6372:   Input Parameters:
6373: + dm  - The `DM`
6374: - dim - The topological dimension

6376:   Level: beginner

6378: .seealso: [](ch_dmbase), `DM`, `DMGetDimension()`, `DMCreate()`
6379: @*/
6380: PetscErrorCode DMSetDimension(DM dm, PetscInt dim)
6381: {
6382:   PetscDS  ds;
6383:   PetscInt Nds, n;

6385:   PetscFunctionBegin;
6388:   dm->dim = dim;
6389:   if (dm->dim >= 0) {
6390:     PetscCall(DMGetNumDS(dm, &Nds));
6391:     for (n = 0; n < Nds; ++n) {
6392:       PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds, NULL));
6393:       if (ds->dimEmbed < 0) PetscCall(PetscDSSetCoordinateDimension(ds, dim));
6394:     }
6395:   }
6396:   PetscFunctionReturn(PETSC_SUCCESS);
6397: }

6399: /*@
6400:   DMGetDimPoints - Get the half-open interval for all points of a given dimension

6402:   Collective

6404:   Input Parameters:
6405: + dm  - the `DM`
6406: - dim - the dimension

6408:   Output Parameters:
6409: + pStart - The first point of the given dimension
6410: - pEnd   - The first point following points of the given dimension

6412:   Level: intermediate

6414:   Note:
6415:   The points are vertices in the Hasse diagram encoding the topology. This is explained in
6416:   https://arxiv.org/abs/0908.4427. If no points exist of this dimension in the storage scheme,
6417:   then the interval is empty.

6419: .seealso: [](ch_dmbase), `DM`, `DMPLEX`, `DMPlexGetDepthStratum()`, `DMPlexGetHeightStratum()`
6420: @*/
6421: PetscErrorCode DMGetDimPoints(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
6422: {
6423:   PetscInt d;

6425:   PetscFunctionBegin;
6427:   PetscCall(DMGetDimension(dm, &d));
6428:   PetscCheck((dim >= 0) && (dim <= d), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim);
6429:   PetscUseTypeMethod(dm, getdimpoints, dim, pStart, pEnd);
6430:   PetscFunctionReturn(PETSC_SUCCESS);
6431: }

6433: /*@
6434:   DMGetOutputDM - Retrieve the `DM` associated with the layout for output

6436:   Collective

6438:   Input Parameter:
6439: . dm - The original `DM`

6441:   Output Parameter:
6442: . odm - The `DM` which provides the layout for output

6444:   Level: intermediate

6446:   Note:
6447:   In some situations the vector obtained with `DMCreateGlobalVector()` excludes points for degrees of freedom that are associated with fixed (Dirichelet) boundary
6448:   conditions since the algebraic solver does not solve for those variables. The output `DM` includes these excluded points and its global vector contains the
6449:   locations for those dof so that they can be output to a file or other viewer along with the unconstrained dof.

6451: .seealso: [](ch_dmbase), `DM`, `VecView()`, `DMGetGlobalSection()`, `DMCreateGlobalVector()`, `PetscSectionHasConstraints()`, `DMSetGlobalSection()`
6452: @*/
6453: PetscErrorCode DMGetOutputDM(DM dm, DM *odm)
6454: {
6455:   PetscSection section;
6456:   IS           perm;
6457:   PetscBool    hasConstraints, newDM, gnewDM;

6459:   PetscFunctionBegin;
6461:   PetscAssertPointer(odm, 2);
6462:   PetscCall(DMGetLocalSection(dm, &section));
6463:   PetscCall(PetscSectionHasConstraints(section, &hasConstraints));
6464:   PetscCall(PetscSectionGetPermutation(section, &perm));
6465:   newDM = hasConstraints || perm ? PETSC_TRUE : PETSC_FALSE;
6466:   PetscCall(MPIU_Allreduce(&newDM, &gnewDM, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)dm)));
6467:   if (!gnewDM) {
6468:     *odm = dm;
6469:     PetscFunctionReturn(PETSC_SUCCESS);
6470:   }
6471:   if (!dm->dmBC) {
6472:     PetscSection newSection, gsection;
6473:     PetscSF      sf;
6474:     PetscBool    usePerm = dm->ignorePermOutput ? PETSC_FALSE : PETSC_TRUE;

6476:     PetscCall(DMClone(dm, &dm->dmBC));
6477:     PetscCall(DMCopyDisc(dm, dm->dmBC));
6478:     PetscCall(PetscSectionClone(section, &newSection));
6479:     PetscCall(DMSetLocalSection(dm->dmBC, newSection));
6480:     PetscCall(PetscSectionDestroy(&newSection));
6481:     PetscCall(DMGetPointSF(dm->dmBC, &sf));
6482:     PetscCall(PetscSectionCreateGlobalSection(section, sf, usePerm, PETSC_TRUE, PETSC_FALSE, &gsection));
6483:     PetscCall(DMSetGlobalSection(dm->dmBC, gsection));
6484:     PetscCall(PetscSectionDestroy(&gsection));
6485:   }
6486:   *odm = dm->dmBC;
6487:   PetscFunctionReturn(PETSC_SUCCESS);
6488: }

6490: /*@
6491:   DMGetOutputSequenceNumber - Retrieve the sequence number/value for output

6493:   Input Parameter:
6494: . dm - The original `DM`

6496:   Output Parameters:
6497: + num - The output sequence number
6498: - val - The output sequence value

6500:   Level: intermediate

6502:   Note:
6503:   This is intended for output that should appear in sequence, for instance
6504:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6506:   Developer Note:
6507:   The `DM` serves as a convenient place to store the current iteration value. The iteration is not
6508:   not directly related to the `DM`.

6510: .seealso: [](ch_dmbase), `DM`, `VecView()`
6511: @*/
6512: PetscErrorCode DMGetOutputSequenceNumber(DM dm, PetscInt *num, PetscReal *val)
6513: {
6514:   PetscFunctionBegin;
6516:   if (num) {
6517:     PetscAssertPointer(num, 2);
6518:     *num = dm->outputSequenceNum;
6519:   }
6520:   if (val) {
6521:     PetscAssertPointer(val, 3);
6522:     *val = dm->outputSequenceVal;
6523:   }
6524:   PetscFunctionReturn(PETSC_SUCCESS);
6525: }

6527: /*@
6528:   DMSetOutputSequenceNumber - Set the sequence number/value for output

6530:   Input Parameters:
6531: + dm  - The original `DM`
6532: . num - The output sequence number
6533: - val - The output sequence value

6535:   Level: intermediate

6537:   Note:
6538:   This is intended for output that should appear in sequence, for instance
6539:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6541: .seealso: [](ch_dmbase), `DM`, `VecView()`
6542: @*/
6543: PetscErrorCode DMSetOutputSequenceNumber(DM dm, PetscInt num, PetscReal val)
6544: {
6545:   PetscFunctionBegin;
6547:   dm->outputSequenceNum = num;
6548:   dm->outputSequenceVal = val;
6549:   PetscFunctionReturn(PETSC_SUCCESS);
6550: }

6552: /*@
6553:   DMOutputSequenceLoad - Retrieve the sequence value from a `PetscViewer`

6555:   Input Parameters:
6556: + dm     - The original `DM`
6557: . viewer - The `PetscViewer` to get it from
6558: . name   - The sequence name
6559: - num    - The output sequence number

6561:   Output Parameter:
6562: . val - The output sequence value

6564:   Level: intermediate

6566:   Note:
6567:   This is intended for output that should appear in sequence, for instance
6568:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6570:   Developer Note:
6571:   It is unclear at the user API level why a `DM` is needed as input

6573: .seealso: [](ch_dmbase), `DM`, `DMGetOutputSequenceNumber()`, `DMSetOutputSequenceNumber()`, `VecView()`
6574: @*/
6575: PetscErrorCode DMOutputSequenceLoad(DM dm, PetscViewer viewer, const char name[], PetscInt num, PetscReal *val)
6576: {
6577:   PetscBool ishdf5;

6579:   PetscFunctionBegin;
6582:   PetscAssertPointer(name, 3);
6583:   PetscAssertPointer(val, 5);
6584:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
6585:   if (ishdf5) {
6586: #if defined(PETSC_HAVE_HDF5)
6587:     PetscScalar value;

6589:     PetscCall(DMSequenceLoad_HDF5_Internal(dm, name, num, &value, viewer));
6590:     *val = PetscRealPart(value);
6591: #endif
6592:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerHDF5Open()");
6593:   PetscFunctionReturn(PETSC_SUCCESS);
6594: }

6596: /*@
6597:   DMGetOutputSequenceLength - Retrieve the number of sequence values from a `PetscViewer`

6599:   Input Parameters:
6600: + dm     - The original `DM`
6601: . viewer - The `PetscViewer` to get it from
6602: - name   - The sequence name

6604:   Output Parameter:
6605: . len - The length of the output sequence

6607:   Level: intermediate

6609:   Note:
6610:   This is intended for output that should appear in sequence, for instance
6611:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6613:   Developer Note:
6614:   It is unclear at the user API level why a `DM` is needed as input

6616: .seealso: [](ch_dmbase), `DM`, `DMGetOutputSequenceNumber()`, `DMSetOutputSequenceNumber()`, `VecView()`
6617: @*/
6618: PetscErrorCode DMGetOutputSequenceLength(DM dm, PetscViewer viewer, const char name[], PetscInt *len)
6619: {
6620:   PetscBool ishdf5;

6622:   PetscFunctionBegin;
6625:   PetscAssertPointer(name, 3);
6626:   PetscAssertPointer(len, 4);
6627:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
6628:   if (ishdf5) {
6629: #if defined(PETSC_HAVE_HDF5)
6630:     PetscCall(DMSequenceGetLength_HDF5_Internal(dm, name, len, viewer));
6631: #endif
6632:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerHDF5Open()");
6633:   PetscFunctionReturn(PETSC_SUCCESS);
6634: }

6636: /*@
6637:   DMGetUseNatural - Get the flag for creating a mapping to the natural order when a `DM` is (re)distributed in parallel

6639:   Not Collective

6641:   Input Parameter:
6642: . dm - The `DM`

6644:   Output Parameter:
6645: . useNatural - `PETSC_TRUE` to build the mapping to a natural order during distribution

6647:   Level: beginner

6649: .seealso: [](ch_dmbase), `DM`, `DMSetUseNatural()`, `DMCreate()`
6650: @*/
6651: PetscErrorCode DMGetUseNatural(DM dm, PetscBool *useNatural)
6652: {
6653:   PetscFunctionBegin;
6655:   PetscAssertPointer(useNatural, 2);
6656:   *useNatural = dm->useNatural;
6657:   PetscFunctionReturn(PETSC_SUCCESS);
6658: }

6660: /*@
6661:   DMSetUseNatural - Set the flag for creating a mapping to the natural order when a `DM` is (re)distributed in parallel

6663:   Collective

6665:   Input Parameters:
6666: + dm         - The `DM`
6667: - useNatural - `PETSC_TRUE` to build the mapping to a natural order during distribution

6669:   Level: beginner

6671:   Note:
6672:   This also causes the map to be build after `DMCreateSubDM()` and `DMCreateSuperDM()`

6674: .seealso: [](ch_dmbase), `DM`, `DMGetUseNatural()`, `DMCreate()`, `DMPlexDistribute()`, `DMCreateSubDM()`, `DMCreateSuperDM()`
6675: @*/
6676: PetscErrorCode DMSetUseNatural(DM dm, PetscBool useNatural)
6677: {
6678:   PetscFunctionBegin;
6681:   dm->useNatural = useNatural;
6682:   PetscFunctionReturn(PETSC_SUCCESS);
6683: }

6685: /*@
6686:   DMCreateLabel - Create a label of the given name if it does not already exist in the `DM`

6688:   Not Collective

6690:   Input Parameters:
6691: + dm   - The `DM` object
6692: - name - The label name

6694:   Level: intermediate

6696: .seealso: [](ch_dmbase), `DM`, `DMLabelCreate()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6697: @*/
6698: PetscErrorCode DMCreateLabel(DM dm, const char name[])
6699: {
6700:   PetscBool flg;
6701:   DMLabel   label;

6703:   PetscFunctionBegin;
6705:   PetscAssertPointer(name, 2);
6706:   PetscCall(DMHasLabel(dm, name, &flg));
6707:   if (!flg) {
6708:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, &label));
6709:     PetscCall(DMAddLabel(dm, label));
6710:     PetscCall(DMLabelDestroy(&label));
6711:   }
6712:   PetscFunctionReturn(PETSC_SUCCESS);
6713: }

6715: /*@
6716:   DMCreateLabelAtIndex - Create a label of the given name at the given index. If it already exists in the `DM`, move it to this index.

6718:   Not Collective

6720:   Input Parameters:
6721: + dm   - The `DM` object
6722: . l    - The index for the label
6723: - name - The label name

6725:   Level: intermediate

6727: .seealso: [](ch_dmbase), `DM`, `DMCreateLabel()`, `DMLabelCreate()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6728: @*/
6729: PetscErrorCode DMCreateLabelAtIndex(DM dm, PetscInt l, const char name[])
6730: {
6731:   DMLabelLink orig, prev = NULL;
6732:   DMLabel     label;
6733:   PetscInt    Nl, m;
6734:   PetscBool   flg, match;
6735:   const char *lname;

6737:   PetscFunctionBegin;
6739:   PetscAssertPointer(name, 3);
6740:   PetscCall(DMHasLabel(dm, name, &flg));
6741:   if (!flg) {
6742:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, &label));
6743:     PetscCall(DMAddLabel(dm, label));
6744:     PetscCall(DMLabelDestroy(&label));
6745:   }
6746:   PetscCall(DMGetNumLabels(dm, &Nl));
6747:   PetscCheck(l < Nl, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label index %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", l, Nl);
6748:   for (m = 0, orig = dm->labels; m < Nl; ++m, prev = orig, orig = orig->next) {
6749:     PetscCall(PetscObjectGetName((PetscObject)orig->label, &lname));
6750:     PetscCall(PetscStrcmp(name, lname, &match));
6751:     if (match) break;
6752:   }
6753:   if (m == l) PetscFunctionReturn(PETSC_SUCCESS);
6754:   if (!m) dm->labels = orig->next;
6755:   else prev->next = orig->next;
6756:   if (!l) {
6757:     orig->next = dm->labels;
6758:     dm->labels = orig;
6759:   } else {
6760:     for (m = 0, prev = dm->labels; m < l - 1; ++m, prev = prev->next);
6761:     orig->next = prev->next;
6762:     prev->next = orig;
6763:   }
6764:   PetscFunctionReturn(PETSC_SUCCESS);
6765: }

6767: /*@
6768:   DMGetLabelValue - Get the value in a `DMLabel` for the given point, with -1 as the default

6770:   Not Collective

6772:   Input Parameters:
6773: + dm    - The `DM` object
6774: . name  - The label name
6775: - point - The mesh point

6777:   Output Parameter:
6778: . value - The label value for this point, or -1 if the point is not in the label

6780:   Level: beginner

6782: .seealso: [](ch_dmbase), `DM`, `DMLabelGetValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6783: @*/
6784: PetscErrorCode DMGetLabelValue(DM dm, const char name[], PetscInt point, PetscInt *value)
6785: {
6786:   DMLabel label;

6788:   PetscFunctionBegin;
6790:   PetscAssertPointer(name, 2);
6791:   PetscCall(DMGetLabel(dm, name, &label));
6792:   PetscCheck(label, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No label named %s was found", name);
6793:   PetscCall(DMLabelGetValue(label, point, value));
6794:   PetscFunctionReturn(PETSC_SUCCESS);
6795: }

6797: /*@
6798:   DMSetLabelValue - Add a point to a `DMLabel` with given value

6800:   Not Collective

6802:   Input Parameters:
6803: + dm    - The `DM` object
6804: . name  - The label name
6805: . point - The mesh point
6806: - value - The label value for this point

6808:   Output Parameter:

6810:   Level: beginner

6812: .seealso: [](ch_dmbase), `DM`, `DMLabelSetValue()`, `DMGetStratumIS()`, `DMClearLabelValue()`
6813: @*/
6814: PetscErrorCode DMSetLabelValue(DM dm, const char name[], PetscInt point, PetscInt value)
6815: {
6816:   DMLabel label;

6818:   PetscFunctionBegin;
6820:   PetscAssertPointer(name, 2);
6821:   PetscCall(DMGetLabel(dm, name, &label));
6822:   if (!label) {
6823:     PetscCall(DMCreateLabel(dm, name));
6824:     PetscCall(DMGetLabel(dm, name, &label));
6825:   }
6826:   PetscCall(DMLabelSetValue(label, point, value));
6827:   PetscFunctionReturn(PETSC_SUCCESS);
6828: }

6830: /*@
6831:   DMClearLabelValue - Remove a point from a `DMLabel` with given value

6833:   Not Collective

6835:   Input Parameters:
6836: + dm    - The `DM` object
6837: . name  - The label name
6838: . point - The mesh point
6839: - value - The label value for this point

6841:   Level: beginner

6843: .seealso: [](ch_dmbase), `DM`, `DMLabelClearValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6844: @*/
6845: PetscErrorCode DMClearLabelValue(DM dm, const char name[], PetscInt point, PetscInt value)
6846: {
6847:   DMLabel label;

6849:   PetscFunctionBegin;
6851:   PetscAssertPointer(name, 2);
6852:   PetscCall(DMGetLabel(dm, name, &label));
6853:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6854:   PetscCall(DMLabelClearValue(label, point, value));
6855:   PetscFunctionReturn(PETSC_SUCCESS);
6856: }

6858: /*@
6859:   DMGetLabelSize - Get the value of `DMLabelGetNumValues()` of a `DMLabel` in the `DM`

6861:   Not Collective

6863:   Input Parameters:
6864: + dm   - The `DM` object
6865: - name - The label name

6867:   Output Parameter:
6868: . size - The number of different integer ids, or 0 if the label does not exist

6870:   Level: beginner

6872:   Developer Note:
6873:   This should be renamed to something like `DMGetLabelNumValues()` or removed.

6875: .seealso: [](ch_dmbase), `DM`, `DMLabelGetNumValues()`, `DMSetLabelValue()`, `DMGetLabel()`
6876: @*/
6877: PetscErrorCode DMGetLabelSize(DM dm, const char name[], PetscInt *size)
6878: {
6879:   DMLabel label;

6881:   PetscFunctionBegin;
6883:   PetscAssertPointer(name, 2);
6884:   PetscAssertPointer(size, 3);
6885:   PetscCall(DMGetLabel(dm, name, &label));
6886:   *size = 0;
6887:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6888:   PetscCall(DMLabelGetNumValues(label, size));
6889:   PetscFunctionReturn(PETSC_SUCCESS);
6890: }

6892: /*@
6893:   DMGetLabelIdIS - Get the `DMLabelGetValueIS()` from a `DMLabel` in the `DM`

6895:   Not Collective

6897:   Input Parameters:
6898: + dm   - The `DM` object
6899: - name - The label name

6901:   Output Parameter:
6902: . ids - The integer ids, or `NULL` if the label does not exist

6904:   Level: beginner

6906: .seealso: [](ch_dmbase), `DM`, `DMLabelGetValueIS()`, `DMGetLabelSize()`
6907: @*/
6908: PetscErrorCode DMGetLabelIdIS(DM dm, const char name[], IS *ids)
6909: {
6910:   DMLabel label;

6912:   PetscFunctionBegin;
6914:   PetscAssertPointer(name, 2);
6915:   PetscAssertPointer(ids, 3);
6916:   PetscCall(DMGetLabel(dm, name, &label));
6917:   *ids = NULL;
6918:   if (label) {
6919:     PetscCall(DMLabelGetValueIS(label, ids));
6920:   } else {
6921:     /* returning an empty IS */
6922:     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_USE_POINTER, ids));
6923:   }
6924:   PetscFunctionReturn(PETSC_SUCCESS);
6925: }

6927: /*@
6928:   DMGetStratumSize - Get the number of points in a label stratum

6930:   Not Collective

6932:   Input Parameters:
6933: + dm    - The `DM` object
6934: . name  - The label name of the stratum
6935: - value - The stratum value

6937:   Output Parameter:
6938: . size - The number of points, also called the stratum size

6940:   Level: beginner

6942: .seealso: [](ch_dmbase), `DM`, `DMLabelGetStratumSize()`, `DMGetLabelSize()`, `DMGetLabelIds()`
6943: @*/
6944: PetscErrorCode DMGetStratumSize(DM dm, const char name[], PetscInt value, PetscInt *size)
6945: {
6946:   DMLabel label;

6948:   PetscFunctionBegin;
6950:   PetscAssertPointer(name, 2);
6951:   PetscAssertPointer(size, 4);
6952:   PetscCall(DMGetLabel(dm, name, &label));
6953:   *size = 0;
6954:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6955:   PetscCall(DMLabelGetStratumSize(label, value, size));
6956:   PetscFunctionReturn(PETSC_SUCCESS);
6957: }

6959: /*@
6960:   DMGetStratumIS - Get the points in a label stratum

6962:   Not Collective

6964:   Input Parameters:
6965: + dm    - The `DM` object
6966: . name  - The label name
6967: - value - The stratum value

6969:   Output Parameter:
6970: . points - The stratum points, or `NULL` if the label does not exist or does not have that value

6972:   Level: beginner

6974: .seealso: [](ch_dmbase), `DM`, `DMLabelGetStratumIS()`, `DMGetStratumSize()`
6975: @*/
6976: PetscErrorCode DMGetStratumIS(DM dm, const char name[], PetscInt value, IS *points)
6977: {
6978:   DMLabel label;

6980:   PetscFunctionBegin;
6982:   PetscAssertPointer(name, 2);
6983:   PetscAssertPointer(points, 4);
6984:   PetscCall(DMGetLabel(dm, name, &label));
6985:   *points = NULL;
6986:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6987:   PetscCall(DMLabelGetStratumIS(label, value, points));
6988:   PetscFunctionReturn(PETSC_SUCCESS);
6989: }

6991: /*@
6992:   DMSetStratumIS - Set the points in a label stratum

6994:   Not Collective

6996:   Input Parameters:
6997: + dm     - The `DM` object
6998: . name   - The label name
6999: . value  - The stratum value
7000: - points - The stratum points

7002:   Level: beginner

7004: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMClearLabelStratum()`, `DMLabelClearStratum()`, `DMLabelSetStratumIS()`, `DMGetStratumSize()`
7005: @*/
7006: PetscErrorCode DMSetStratumIS(DM dm, const char name[], PetscInt value, IS points)
7007: {
7008:   DMLabel label;

7010:   PetscFunctionBegin;
7012:   PetscAssertPointer(name, 2);
7014:   PetscCall(DMGetLabel(dm, name, &label));
7015:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
7016:   PetscCall(DMLabelSetStratumIS(label, value, points));
7017:   PetscFunctionReturn(PETSC_SUCCESS);
7018: }

7020: /*@
7021:   DMClearLabelStratum - Remove all points from a stratum from a `DMLabel`

7023:   Not Collective

7025:   Input Parameters:
7026: + dm    - The `DM` object
7027: . name  - The label name
7028: - value - The label value for this point

7030:   Output Parameter:

7032:   Level: beginner

7034: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMLabelClearStratum()`, `DMSetLabelValue()`, `DMGetStratumIS()`, `DMClearLabelValue()`
7035: @*/
7036: PetscErrorCode DMClearLabelStratum(DM dm, const char name[], PetscInt value)
7037: {
7038:   DMLabel label;

7040:   PetscFunctionBegin;
7042:   PetscAssertPointer(name, 2);
7043:   PetscCall(DMGetLabel(dm, name, &label));
7044:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
7045:   PetscCall(DMLabelClearStratum(label, value));
7046:   PetscFunctionReturn(PETSC_SUCCESS);
7047: }

7049: /*@
7050:   DMGetNumLabels - Return the number of labels defined by on the `DM`

7052:   Not Collective

7054:   Input Parameter:
7055: . dm - The `DM` object

7057:   Output Parameter:
7058: . numLabels - the number of Labels

7060:   Level: intermediate

7062: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetLabelByNum()`, `DMGetLabelName()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7063: @*/
7064: PetscErrorCode DMGetNumLabels(DM dm, PetscInt *numLabels)
7065: {
7066:   DMLabelLink next = dm->labels;
7067:   PetscInt    n    = 0;

7069:   PetscFunctionBegin;
7071:   PetscAssertPointer(numLabels, 2);
7072:   while (next) {
7073:     ++n;
7074:     next = next->next;
7075:   }
7076:   *numLabels = n;
7077:   PetscFunctionReturn(PETSC_SUCCESS);
7078: }

7080: /*@
7081:   DMGetLabelName - Return the name of nth label

7083:   Not Collective

7085:   Input Parameters:
7086: + dm - The `DM` object
7087: - n  - the label number

7089:   Output Parameter:
7090: . name - the label name

7092:   Level: intermediate

7094:   Developer Note:
7095:   Some of the functions that appropriate on labels using their number have the suffix ByNum, others do not.

7097: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetLabelByNum()`, `DMGetLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7098: @*/
7099: PetscErrorCode DMGetLabelName(DM dm, PetscInt n, const char *name[])
7100: {
7101:   DMLabelLink next = dm->labels;
7102:   PetscInt    l    = 0;

7104:   PetscFunctionBegin;
7106:   PetscAssertPointer(name, 3);
7107:   while (next) {
7108:     if (l == n) {
7109:       PetscCall(PetscObjectGetName((PetscObject)next->label, name));
7110:       PetscFunctionReturn(PETSC_SUCCESS);
7111:     }
7112:     ++l;
7113:     next = next->next;
7114:   }
7115:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %" PetscInt_FMT " does not exist in this DM", n);
7116: }

7118: /*@
7119:   DMHasLabel - Determine whether the `DM` has a label of a given name

7121:   Not Collective

7123:   Input Parameters:
7124: + dm   - The `DM` object
7125: - name - The label name

7127:   Output Parameter:
7128: . hasLabel - `PETSC_TRUE` if the label is present

7130:   Level: intermediate

7132: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetLabel()`, `DMGetLabelByNum()`, `DMCreateLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7133: @*/
7134: PetscErrorCode DMHasLabel(DM dm, const char name[], PetscBool *hasLabel)
7135: {
7136:   DMLabelLink next = dm->labels;
7137:   const char *lname;

7139:   PetscFunctionBegin;
7141:   PetscAssertPointer(name, 2);
7142:   PetscAssertPointer(hasLabel, 3);
7143:   *hasLabel = PETSC_FALSE;
7144:   while (next) {
7145:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7146:     PetscCall(PetscStrcmp(name, lname, hasLabel));
7147:     if (*hasLabel) break;
7148:     next = next->next;
7149:   }
7150:   PetscFunctionReturn(PETSC_SUCCESS);
7151: }

7153: // PetscClangLinter pragma ignore: -fdoc-section-header-unknown
7154: /*@
7155:   DMGetLabel - Return the label of a given name, or `NULL`, from a `DM`

7157:   Not Collective

7159:   Input Parameters:
7160: + dm   - The `DM` object
7161: - name - The label name

7163:   Output Parameter:
7164: . label - The `DMLabel`, or `NULL` if the label is absent

7166:   Default labels in a `DMPLEX`:
7167: + "depth"       - Holds the depth (co-dimension) of each mesh point
7168: . "celltype"    - Holds the topological type of each cell
7169: . "ghost"       - If the DM is distributed with overlap, this marks the cells and faces in the overlap
7170: . "Cell Sets"   - Mirrors the cell sets defined by GMsh and ExodusII
7171: . "Face Sets"   - Mirrors the face sets defined by GMsh and ExodusII
7172: - "Vertex Sets" - Mirrors the vertex sets defined by GMsh

7174:   Level: intermediate

7176: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMHasLabel()`, `DMGetLabelByNum()`, `DMAddLabel()`, `DMCreateLabel()`, `DMPlexGetDepthLabel()`, `DMPlexGetCellType()`
7177: @*/
7178: PetscErrorCode DMGetLabel(DM dm, const char name[], DMLabel *label)
7179: {
7180:   DMLabelLink next = dm->labels;
7181:   PetscBool   hasLabel;
7182:   const char *lname;

7184:   PetscFunctionBegin;
7186:   PetscAssertPointer(name, 2);
7187:   PetscAssertPointer(label, 3);
7188:   *label = NULL;
7189:   while (next) {
7190:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7191:     PetscCall(PetscStrcmp(name, lname, &hasLabel));
7192:     if (hasLabel) {
7193:       *label = next->label;
7194:       break;
7195:     }
7196:     next = next->next;
7197:   }
7198:   PetscFunctionReturn(PETSC_SUCCESS);
7199: }

7201: /*@
7202:   DMGetLabelByNum - Return the nth label on a `DM`

7204:   Not Collective

7206:   Input Parameters:
7207: + dm - The `DM` object
7208: - n  - the label number

7210:   Output Parameter:
7211: . label - the label

7213:   Level: intermediate

7215: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMAddLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7216: @*/
7217: PetscErrorCode DMGetLabelByNum(DM dm, PetscInt n, DMLabel *label)
7218: {
7219:   DMLabelLink next = dm->labels;
7220:   PetscInt    l    = 0;

7222:   PetscFunctionBegin;
7224:   PetscAssertPointer(label, 3);
7225:   while (next) {
7226:     if (l == n) {
7227:       *label = next->label;
7228:       PetscFunctionReturn(PETSC_SUCCESS);
7229:     }
7230:     ++l;
7231:     next = next->next;
7232:   }
7233:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %" PetscInt_FMT " does not exist in this DM", n);
7234: }

7236: /*@
7237:   DMAddLabel - Add the label to this `DM`

7239:   Not Collective

7241:   Input Parameters:
7242: + dm    - The `DM` object
7243: - label - The `DMLabel`

7245:   Level: developer

7247: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7248: @*/
7249: PetscErrorCode DMAddLabel(DM dm, DMLabel label)
7250: {
7251:   DMLabelLink l, *p, tmpLabel;
7252:   PetscBool   hasLabel;
7253:   const char *lname;
7254:   PetscBool   flg;

7256:   PetscFunctionBegin;
7258:   PetscCall(PetscObjectGetName((PetscObject)label, &lname));
7259:   PetscCall(DMHasLabel(dm, lname, &hasLabel));
7260:   PetscCheck(!hasLabel, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %s already exists in this DM", lname);
7261:   PetscCall(PetscCalloc1(1, &tmpLabel));
7262:   tmpLabel->label  = label;
7263:   tmpLabel->output = PETSC_TRUE;
7264:   for (p = &dm->labels; (l = *p); p = &l->next) { }
7265:   *p = tmpLabel;
7266:   PetscCall(PetscObjectReference((PetscObject)label));
7267:   PetscCall(PetscStrcmp(lname, "depth", &flg));
7268:   if (flg) dm->depthLabel = label;
7269:   PetscCall(PetscStrcmp(lname, "celltype", &flg));
7270:   if (flg) dm->celltypeLabel = label;
7271:   PetscFunctionReturn(PETSC_SUCCESS);
7272: }

7274: // PetscClangLinter pragma ignore: -fdoc-section-header-unknown
7275: /*@
7276:   DMSetLabel - Replaces the label of a given name, or ignores it if the name is not present

7278:   Not Collective

7280:   Input Parameters:
7281: + dm    - The `DM` object
7282: - label - The `DMLabel`, having the same name, to substitute

7284:   Default labels in a `DMPLEX`:
7285: + "depth"       - Holds the depth (co-dimension) of each mesh point
7286: . "celltype"    - Holds the topological type of each cell
7287: . "ghost"       - If the DM is distributed with overlap, this marks the cells and faces in the overlap
7288: . "Cell Sets"   - Mirrors the cell sets defined by GMsh and ExodusII
7289: . "Face Sets"   - Mirrors the face sets defined by GMsh and ExodusII
7290: - "Vertex Sets" - Mirrors the vertex sets defined by GMsh

7292:   Level: intermediate

7294: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMPlexGetDepthLabel()`, `DMPlexGetCellType()`
7295: @*/
7296: PetscErrorCode DMSetLabel(DM dm, DMLabel label)
7297: {
7298:   DMLabelLink next = dm->labels;
7299:   PetscBool   hasLabel, flg;
7300:   const char *name, *lname;

7302:   PetscFunctionBegin;
7305:   PetscCall(PetscObjectGetName((PetscObject)label, &name));
7306:   while (next) {
7307:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7308:     PetscCall(PetscStrcmp(name, lname, &hasLabel));
7309:     if (hasLabel) {
7310:       PetscCall(PetscObjectReference((PetscObject)label));
7311:       PetscCall(PetscStrcmp(lname, "depth", &flg));
7312:       if (flg) dm->depthLabel = label;
7313:       PetscCall(PetscStrcmp(lname, "celltype", &flg));
7314:       if (flg) dm->celltypeLabel = label;
7315:       PetscCall(DMLabelDestroy(&next->label));
7316:       next->label = label;
7317:       break;
7318:     }
7319:     next = next->next;
7320:   }
7321:   PetscFunctionReturn(PETSC_SUCCESS);
7322: }

7324: /*@
7325:   DMRemoveLabel - Remove the label given by name from this `DM`

7327:   Not Collective

7329:   Input Parameters:
7330: + dm   - The `DM` object
7331: - name - The label name

7333:   Output Parameter:
7334: . label - The `DMLabel`, or `NULL` if the label is absent. Pass in `NULL` to call `DMLabelDestroy()` on the label, otherwise the
7335:           caller is responsible for calling `DMLabelDestroy()`.

7337:   Level: developer

7339: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMLabelDestroy()`, `DMRemoveLabelBySelf()`
7340: @*/
7341: PetscErrorCode DMRemoveLabel(DM dm, const char name[], DMLabel *label)
7342: {
7343:   DMLabelLink link, *pnext;
7344:   PetscBool   hasLabel;
7345:   const char *lname;

7347:   PetscFunctionBegin;
7349:   PetscAssertPointer(name, 2);
7350:   if (label) {
7351:     PetscAssertPointer(label, 3);
7352:     *label = NULL;
7353:   }
7354:   for (pnext = &dm->labels; (link = *pnext); pnext = &link->next) {
7355:     PetscCall(PetscObjectGetName((PetscObject)link->label, &lname));
7356:     PetscCall(PetscStrcmp(name, lname, &hasLabel));
7357:     if (hasLabel) {
7358:       *pnext = link->next; /* Remove from list */
7359:       PetscCall(PetscStrcmp(name, "depth", &hasLabel));
7360:       if (hasLabel) dm->depthLabel = NULL;
7361:       PetscCall(PetscStrcmp(name, "celltype", &hasLabel));
7362:       if (hasLabel) dm->celltypeLabel = NULL;
7363:       if (label) *label = link->label;
7364:       else PetscCall(DMLabelDestroy(&link->label));
7365:       PetscCall(PetscFree(link));
7366:       break;
7367:     }
7368:   }
7369:   PetscFunctionReturn(PETSC_SUCCESS);
7370: }

7372: /*@
7373:   DMRemoveLabelBySelf - Remove the label from this `DM`

7375:   Not Collective

7377:   Input Parameters:
7378: + dm           - The `DM` object
7379: . label        - The `DMLabel` to be removed from the `DM`
7380: - failNotFound - Should it fail if the label is not found in the `DM`?

7382:   Level: developer

7384:   Note:
7385:   Only exactly the same instance is removed if found, name match is ignored.
7386:   If the `DM` has an exclusive reference to the label, the label gets destroyed and
7387:   *label nullified.

7389: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabel()` `DMGetLabelValue()`, `DMSetLabelValue()`, `DMLabelDestroy()`, `DMRemoveLabel()`
7390: @*/
7391: PetscErrorCode DMRemoveLabelBySelf(DM dm, DMLabel *label, PetscBool failNotFound)
7392: {
7393:   DMLabelLink link, *pnext;
7394:   PetscBool   hasLabel = PETSC_FALSE;

7396:   PetscFunctionBegin;
7398:   PetscAssertPointer(label, 2);
7399:   if (!*label && !failNotFound) PetscFunctionReturn(PETSC_SUCCESS);
7402:   for (pnext = &dm->labels; (link = *pnext); pnext = &link->next) {
7403:     if (*label == link->label) {
7404:       hasLabel = PETSC_TRUE;
7405:       *pnext   = link->next; /* Remove from list */
7406:       if (*label == dm->depthLabel) dm->depthLabel = NULL;
7407:       if (*label == dm->celltypeLabel) dm->celltypeLabel = NULL;
7408:       if (((PetscObject)link->label)->refct < 2) *label = NULL; /* nullify if exclusive reference */
7409:       PetscCall(DMLabelDestroy(&link->label));
7410:       PetscCall(PetscFree(link));
7411:       break;
7412:     }
7413:   }
7414:   PetscCheck(hasLabel || !failNotFound, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Given label not found in DM");
7415:   PetscFunctionReturn(PETSC_SUCCESS);
7416: }

7418: /*@
7419:   DMGetLabelOutput - Get the output flag for a given label

7421:   Not Collective

7423:   Input Parameters:
7424: + dm   - The `DM` object
7425: - name - The label name

7427:   Output Parameter:
7428: . output - The flag for output

7430:   Level: developer

7432: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMSetLabelOutput()`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7433: @*/
7434: PetscErrorCode DMGetLabelOutput(DM dm, const char name[], PetscBool *output)
7435: {
7436:   DMLabelLink next = dm->labels;
7437:   const char *lname;

7439:   PetscFunctionBegin;
7441:   PetscAssertPointer(name, 2);
7442:   PetscAssertPointer(output, 3);
7443:   while (next) {
7444:     PetscBool flg;

7446:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7447:     PetscCall(PetscStrcmp(name, lname, &flg));
7448:     if (flg) {
7449:       *output = next->output;
7450:       PetscFunctionReturn(PETSC_SUCCESS);
7451:     }
7452:     next = next->next;
7453:   }
7454:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No label named %s was present in this dm", name);
7455: }

7457: /*@
7458:   DMSetLabelOutput - Set if a given label should be saved to a `PetscViewer` in calls to `DMView()`

7460:   Not Collective

7462:   Input Parameters:
7463: + dm     - The `DM` object
7464: . name   - The label name
7465: - output - `PETSC_TRUE` to save the label to the viewer

7467:   Level: developer

7469: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetOutputFlag()`, `DMGetLabelOutput()`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7470: @*/
7471: PetscErrorCode DMSetLabelOutput(DM dm, const char name[], PetscBool output)
7472: {
7473:   DMLabelLink next = dm->labels;
7474:   const char *lname;

7476:   PetscFunctionBegin;
7478:   PetscAssertPointer(name, 2);
7479:   while (next) {
7480:     PetscBool flg;

7482:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7483:     PetscCall(PetscStrcmp(name, lname, &flg));
7484:     if (flg) {
7485:       next->output = output;
7486:       PetscFunctionReturn(PETSC_SUCCESS);
7487:     }
7488:     next = next->next;
7489:   }
7490:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No label named %s was present in this dm", name);
7491: }

7493: /*@
7494:   DMCopyLabels - Copy labels from one `DM` mesh to another `DM` with a superset of the points

7496:   Collective

7498:   Input Parameters:
7499: + dmA   - The `DM` object with initial labels
7500: . dmB   - The `DM` object to which labels are copied
7501: . mode  - Copy labels by pointers (`PETSC_OWN_POINTER`) or duplicate them (`PETSC_COPY_VALUES`)
7502: . all   - Copy all labels including "depth", "dim", and "celltype" (`PETSC_TRUE`) which are otherwise ignored (`PETSC_FALSE`)
7503: - emode - How to behave when a `DMLabel` in the source and destination `DM`s with the same name is encountered (see `DMCopyLabelsMode`)

7505:   Level: intermediate

7507:   Note:
7508:   This is typically used when interpolating or otherwise adding to a mesh, or testing.

7510: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMAddLabel()`, `DMCopyLabelsMode`
7511: @*/
7512: PetscErrorCode DMCopyLabels(DM dmA, DM dmB, PetscCopyMode mode, PetscBool all, DMCopyLabelsMode emode)
7513: {
7514:   DMLabel     label, labelNew, labelOld;
7515:   const char *name;
7516:   PetscBool   flg;
7517:   DMLabelLink link;

7519:   PetscFunctionBegin;
7524:   PetscCheck(mode != PETSC_USE_POINTER, PetscObjectComm((PetscObject)dmA), PETSC_ERR_SUP, "PETSC_USE_POINTER not supported for objects");
7525:   if (dmA == dmB) PetscFunctionReturn(PETSC_SUCCESS);
7526:   for (link = dmA->labels; link; link = link->next) {
7527:     label = link->label;
7528:     PetscCall(PetscObjectGetName((PetscObject)label, &name));
7529:     if (!all) {
7530:       PetscCall(PetscStrcmp(name, "depth", &flg));
7531:       if (flg) continue;
7532:       PetscCall(PetscStrcmp(name, "dim", &flg));
7533:       if (flg) continue;
7534:       PetscCall(PetscStrcmp(name, "celltype", &flg));
7535:       if (flg) continue;
7536:     }
7537:     PetscCall(DMGetLabel(dmB, name, &labelOld));
7538:     if (labelOld) {
7539:       switch (emode) {
7540:       case DM_COPY_LABELS_KEEP:
7541:         continue;
7542:       case DM_COPY_LABELS_REPLACE:
7543:         PetscCall(DMRemoveLabelBySelf(dmB, &labelOld, PETSC_TRUE));
7544:         break;
7545:       case DM_COPY_LABELS_FAIL:
7546:         SETERRQ(PetscObjectComm((PetscObject)dmA), PETSC_ERR_ARG_OUTOFRANGE, "Label %s already exists in destination DM", name);
7547:       default:
7548:         SETERRQ(PetscObjectComm((PetscObject)dmA), PETSC_ERR_ARG_OUTOFRANGE, "Unhandled DMCopyLabelsMode %d", (int)emode);
7549:       }
7550:     }
7551:     if (mode == PETSC_COPY_VALUES) {
7552:       PetscCall(DMLabelDuplicate(label, &labelNew));
7553:     } else {
7554:       labelNew = label;
7555:     }
7556:     PetscCall(DMAddLabel(dmB, labelNew));
7557:     if (mode == PETSC_COPY_VALUES) PetscCall(DMLabelDestroy(&labelNew));
7558:   }
7559:   PetscFunctionReturn(PETSC_SUCCESS);
7560: }

7562: /*@C
7563:   DMCompareLabels - Compare labels between two `DM` objects

7565:   Collective; No Fortran Support

7567:   Input Parameters:
7568: + dm0 - First `DM` object
7569: - dm1 - Second `DM` object

7571:   Output Parameters:
7572: + equal   - (Optional) Flag whether labels of dm0 and dm1 are the same
7573: - message - (Optional) Message describing the difference, or `NULL` if there is no difference

7575:   Level: intermediate

7577:   Notes:
7578:   The output flag equal will be the same on all processes.

7580:   If equal is passed as `NULL` and difference is found, an error is thrown on all processes.

7582:   Make sure to pass equal is `NULL` on all processes or none of them.

7584:   The output message is set independently on each rank.

7586:   message must be freed with `PetscFree()`

7588:   If message is passed as `NULL` and a difference is found, the difference description is printed to stderr in synchronized manner.

7590:   Make sure to pass message as `NULL` on all processes or no processes.

7592:   Labels are matched by name. If the number of labels and their names are equal,
7593:   `DMLabelCompare()` is used to compare each pair of labels with the same name.

7595:   Developer Note:
7596:   Can automatically generate the Fortran stub because `message` must be freed with `PetscFree()`

7598: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMAddLabel()`, `DMCopyLabelsMode`, `DMLabelCompare()`
7599: @*/
7600: PetscErrorCode DMCompareLabels(DM dm0, DM dm1, PetscBool *equal, char **message)
7601: {
7602:   PetscInt    n, i;
7603:   char        msg[PETSC_MAX_PATH_LEN] = "";
7604:   PetscBool   eq;
7605:   MPI_Comm    comm;
7606:   PetscMPIInt rank;

7608:   PetscFunctionBegin;
7611:   PetscCheckSameComm(dm0, 1, dm1, 2);
7612:   if (equal) PetscAssertPointer(equal, 3);
7613:   if (message) PetscAssertPointer(message, 4);
7614:   PetscCall(PetscObjectGetComm((PetscObject)dm0, &comm));
7615:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7616:   {
7617:     PetscInt n1;

7619:     PetscCall(DMGetNumLabels(dm0, &n));
7620:     PetscCall(DMGetNumLabels(dm1, &n1));
7621:     eq = (PetscBool)(n == n1);
7622:     if (!eq) PetscCall(PetscSNPrintf(msg, sizeof(msg), "Number of labels in dm0 = %" PetscInt_FMT " != %" PetscInt_FMT " = Number of labels in dm1", n, n1));
7623:     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPIU_BOOL, MPI_LAND, comm));
7624:     if (!eq) goto finish;
7625:   }
7626:   for (i = 0; i < n; i++) {
7627:     DMLabel     l0, l1;
7628:     const char *name;
7629:     char       *msgInner;

7631:     /* Ignore label order */
7632:     PetscCall(DMGetLabelByNum(dm0, i, &l0));
7633:     PetscCall(PetscObjectGetName((PetscObject)l0, &name));
7634:     PetscCall(DMGetLabel(dm1, name, &l1));
7635:     if (!l1) {
7636:       PetscCall(PetscSNPrintf(msg, sizeof(msg), "Label \"%s\" (#%" PetscInt_FMT " in dm0) not found in dm1", name, i));
7637:       eq = PETSC_FALSE;
7638:       break;
7639:     }
7640:     PetscCall(DMLabelCompare(comm, l0, l1, &eq, &msgInner));
7641:     PetscCall(PetscStrncpy(msg, msgInner, sizeof(msg)));
7642:     PetscCall(PetscFree(msgInner));
7643:     if (!eq) break;
7644:   }
7645:   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPIU_BOOL, MPI_LAND, comm));
7646: finish:
7647:   /* If message output arg not set, print to stderr */
7648:   if (message) {
7649:     *message = NULL;
7650:     if (msg[0]) PetscCall(PetscStrallocpy(msg, message));
7651:   } else {
7652:     if (msg[0]) PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] %s\n", rank, msg));
7653:     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR));
7654:   }
7655:   /* If same output arg not ser and labels are not equal, throw error */
7656:   if (equal) *equal = eq;
7657:   else PetscCheck(eq, comm, PETSC_ERR_ARG_INCOMP, "DMLabels are not the same in dm0 and dm1");
7658:   PetscFunctionReturn(PETSC_SUCCESS);
7659: }

7661: PetscErrorCode DMSetLabelValue_Fast(DM dm, DMLabel *label, const char name[], PetscInt point, PetscInt value)
7662: {
7663:   PetscFunctionBegin;
7664:   PetscAssertPointer(label, 2);
7665:   if (!*label) {
7666:     PetscCall(DMCreateLabel(dm, name));
7667:     PetscCall(DMGetLabel(dm, name, label));
7668:   }
7669:   PetscCall(DMLabelSetValue(*label, point, value));
7670:   PetscFunctionReturn(PETSC_SUCCESS);
7671: }

7673: /*
7674:   Many mesh programs, such as Triangle and TetGen, allow only a single label for each mesh point. Therefore, we would
7675:   like to encode all label IDs using a single, universal label. We can do this by assigning an integer to every
7676:   (label, id) pair in the DM.

7678:   However, a mesh point can have multiple labels, so we must separate all these values. We will assign a bit range to
7679:   each label.
7680: */
7681: PetscErrorCode DMUniversalLabelCreate(DM dm, DMUniversalLabel *universal)
7682: {
7683:   DMUniversalLabel ul;
7684:   PetscBool       *active;
7685:   PetscInt         pStart, pEnd, p, Nl, l, m;

7687:   PetscFunctionBegin;
7688:   PetscCall(PetscMalloc1(1, &ul));
7689:   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "universal", &ul->label));
7690:   PetscCall(DMGetNumLabels(dm, &Nl));
7691:   PetscCall(PetscCalloc1(Nl, &active));
7692:   ul->Nl = 0;
7693:   for (l = 0; l < Nl; ++l) {
7694:     PetscBool   isdepth, iscelltype;
7695:     const char *name;

7697:     PetscCall(DMGetLabelName(dm, l, &name));
7698:     PetscCall(PetscStrncmp(name, "depth", 6, &isdepth));
7699:     PetscCall(PetscStrncmp(name, "celltype", 9, &iscelltype));
7700:     active[l] = !(isdepth || iscelltype) ? PETSC_TRUE : PETSC_FALSE;
7701:     if (active[l]) ++ul->Nl;
7702:   }
7703:   PetscCall(PetscCalloc5(ul->Nl, &ul->names, ul->Nl, &ul->indices, ul->Nl + 1, &ul->offsets, ul->Nl + 1, &ul->bits, ul->Nl, &ul->masks));
7704:   ul->Nv = 0;
7705:   for (l = 0, m = 0; l < Nl; ++l) {
7706:     DMLabel     label;
7707:     PetscInt    nv;
7708:     const char *name;

7710:     if (!active[l]) continue;
7711:     PetscCall(DMGetLabelName(dm, l, &name));
7712:     PetscCall(DMGetLabelByNum(dm, l, &label));
7713:     PetscCall(DMLabelGetNumValues(label, &nv));
7714:     PetscCall(PetscStrallocpy(name, &ul->names[m]));
7715:     ul->indices[m] = l;
7716:     ul->Nv += nv;
7717:     ul->offsets[m + 1] = nv;
7718:     ul->bits[m + 1]    = PetscCeilReal(PetscLog2Real(nv + 1));
7719:     ++m;
7720:   }
7721:   for (l = 1; l <= ul->Nl; ++l) {
7722:     ul->offsets[l] = ul->offsets[l - 1] + ul->offsets[l];
7723:     ul->bits[l]    = ul->bits[l - 1] + ul->bits[l];
7724:   }
7725:   for (l = 0; l < ul->Nl; ++l) {
7726:     PetscInt b;

7728:     ul->masks[l] = 0;
7729:     for (b = ul->bits[l]; b < ul->bits[l + 1]; ++b) ul->masks[l] |= 1 << b;
7730:   }
7731:   PetscCall(PetscMalloc1(ul->Nv, &ul->values));
7732:   for (l = 0, m = 0; l < Nl; ++l) {
7733:     DMLabel         label;
7734:     IS              valueIS;
7735:     const PetscInt *varr;
7736:     PetscInt        nv, v;

7738:     if (!active[l]) continue;
7739:     PetscCall(DMGetLabelByNum(dm, l, &label));
7740:     PetscCall(DMLabelGetNumValues(label, &nv));
7741:     PetscCall(DMLabelGetValueIS(label, &valueIS));
7742:     PetscCall(ISGetIndices(valueIS, &varr));
7743:     for (v = 0; v < nv; ++v) ul->values[ul->offsets[m] + v] = varr[v];
7744:     PetscCall(ISRestoreIndices(valueIS, &varr));
7745:     PetscCall(ISDestroy(&valueIS));
7746:     PetscCall(PetscSortInt(nv, &ul->values[ul->offsets[m]]));
7747:     ++m;
7748:   }
7749:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
7750:   for (p = pStart; p < pEnd; ++p) {
7751:     PetscInt  uval   = 0;
7752:     PetscBool marked = PETSC_FALSE;

7754:     for (l = 0, m = 0; l < Nl; ++l) {
7755:       DMLabel  label;
7756:       PetscInt val, defval, loc, nv;

7758:       if (!active[l]) continue;
7759:       PetscCall(DMGetLabelByNum(dm, l, &label));
7760:       PetscCall(DMLabelGetValue(label, p, &val));
7761:       PetscCall(DMLabelGetDefaultValue(label, &defval));
7762:       if (val == defval) {
7763:         ++m;
7764:         continue;
7765:       }
7766:       nv     = ul->offsets[m + 1] - ul->offsets[m];
7767:       marked = PETSC_TRUE;
7768:       PetscCall(PetscFindInt(val, nv, &ul->values[ul->offsets[m]], &loc));
7769:       PetscCheck(loc >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Label value %" PetscInt_FMT " not found in compression array", val);
7770:       uval += (loc + 1) << ul->bits[m];
7771:       ++m;
7772:     }
7773:     if (marked) PetscCall(DMLabelSetValue(ul->label, p, uval));
7774:   }
7775:   PetscCall(PetscFree(active));
7776:   *universal = ul;
7777:   PetscFunctionReturn(PETSC_SUCCESS);
7778: }

7780: PetscErrorCode DMUniversalLabelDestroy(DMUniversalLabel *universal)
7781: {
7782:   PetscInt l;

7784:   PetscFunctionBegin;
7785:   for (l = 0; l < (*universal)->Nl; ++l) PetscCall(PetscFree((*universal)->names[l]));
7786:   PetscCall(DMLabelDestroy(&(*universal)->label));
7787:   PetscCall(PetscFree5((*universal)->names, (*universal)->indices, (*universal)->offsets, (*universal)->bits, (*universal)->masks));
7788:   PetscCall(PetscFree((*universal)->values));
7789:   PetscCall(PetscFree(*universal));
7790:   *universal = NULL;
7791:   PetscFunctionReturn(PETSC_SUCCESS);
7792: }

7794: PetscErrorCode DMUniversalLabelGetLabel(DMUniversalLabel ul, DMLabel *ulabel)
7795: {
7796:   PetscFunctionBegin;
7797:   PetscAssertPointer(ulabel, 2);
7798:   *ulabel = ul->label;
7799:   PetscFunctionReturn(PETSC_SUCCESS);
7800: }

7802: PetscErrorCode DMUniversalLabelCreateLabels(DMUniversalLabel ul, PetscBool preserveOrder, DM dm)
7803: {
7804:   PetscInt Nl = ul->Nl, l;

7806:   PetscFunctionBegin;
7808:   for (l = 0; l < Nl; ++l) {
7809:     if (preserveOrder) PetscCall(DMCreateLabelAtIndex(dm, ul->indices[l], ul->names[l]));
7810:     else PetscCall(DMCreateLabel(dm, ul->names[l]));
7811:   }
7812:   if (preserveOrder) {
7813:     for (l = 0; l < ul->Nl; ++l) {
7814:       const char *name;
7815:       PetscBool   match;

7817:       PetscCall(DMGetLabelName(dm, ul->indices[l], &name));
7818:       PetscCall(PetscStrcmp(name, ul->names[l], &match));
7819:       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]);
7820:     }
7821:   }
7822:   PetscFunctionReturn(PETSC_SUCCESS);
7823: }

7825: PetscErrorCode DMUniversalLabelSetLabelValue(DMUniversalLabel ul, DM dm, PetscBool useIndex, PetscInt p, PetscInt value)
7826: {
7827:   PetscInt l;

7829:   PetscFunctionBegin;
7830:   for (l = 0; l < ul->Nl; ++l) {
7831:     DMLabel  label;
7832:     PetscInt lval = (value & ul->masks[l]) >> ul->bits[l];

7834:     if (lval) {
7835:       if (useIndex) PetscCall(DMGetLabelByNum(dm, ul->indices[l], &label));
7836:       else PetscCall(DMGetLabel(dm, ul->names[l], &label));
7837:       PetscCall(DMLabelSetValue(label, p, ul->values[ul->offsets[l] + lval - 1]));
7838:     }
7839:   }
7840:   PetscFunctionReturn(PETSC_SUCCESS);
7841: }

7843: /*@
7844:   DMGetCoarseDM - Get the coarse `DM`from which this `DM` was obtained by refinement

7846:   Not Collective

7848:   Input Parameter:
7849: . dm - The `DM` object

7851:   Output Parameter:
7852: . cdm - The coarse `DM`

7854:   Level: intermediate

7856: .seealso: [](ch_dmbase), `DM`, `DMSetCoarseDM()`, `DMCoarsen()`
7857: @*/
7858: PetscErrorCode DMGetCoarseDM(DM dm, DM *cdm)
7859: {
7860:   PetscFunctionBegin;
7862:   PetscAssertPointer(cdm, 2);
7863:   *cdm = dm->coarseMesh;
7864:   PetscFunctionReturn(PETSC_SUCCESS);
7865: }

7867: /*@
7868:   DMSetCoarseDM - Set the coarse `DM` from which this `DM` was obtained by refinement

7870:   Input Parameters:
7871: + dm  - The `DM` object
7872: - cdm - The coarse `DM`

7874:   Level: intermediate

7876:   Note:
7877:   Normally this is set automatically by `DMRefine()`

7879: .seealso: [](ch_dmbase), `DM`, `DMGetCoarseDM()`, `DMCoarsen()`, `DMSetRefine()`, `DMSetFineDM()`
7880: @*/
7881: PetscErrorCode DMSetCoarseDM(DM dm, DM cdm)
7882: {
7883:   PetscFunctionBegin;
7886:   if (dm == cdm) cdm = NULL;
7887:   PetscCall(PetscObjectReference((PetscObject)cdm));
7888:   PetscCall(DMDestroy(&dm->coarseMesh));
7889:   dm->coarseMesh = cdm;
7890:   PetscFunctionReturn(PETSC_SUCCESS);
7891: }

7893: /*@
7894:   DMGetFineDM - Get the fine mesh from which this `DM` was obtained by coarsening

7896:   Input Parameter:
7897: . dm - The `DM` object

7899:   Output Parameter:
7900: . fdm - The fine `DM`

7902:   Level: intermediate

7904: .seealso: [](ch_dmbase), `DM`, `DMSetFineDM()`, `DMCoarsen()`, `DMRefine()`
7905: @*/
7906: PetscErrorCode DMGetFineDM(DM dm, DM *fdm)
7907: {
7908:   PetscFunctionBegin;
7910:   PetscAssertPointer(fdm, 2);
7911:   *fdm = dm->fineMesh;
7912:   PetscFunctionReturn(PETSC_SUCCESS);
7913: }

7915: /*@
7916:   DMSetFineDM - Set the fine mesh from which this was obtained by coarsening

7918:   Input Parameters:
7919: + dm  - The `DM` object
7920: - fdm - The fine `DM`

7922:   Level: developer

7924:   Note:
7925:   Normally this is set automatically by `DMCoarsen()`

7927: .seealso: [](ch_dmbase), `DM`, `DMGetFineDM()`, `DMCoarsen()`, `DMRefine()`
7928: @*/
7929: PetscErrorCode DMSetFineDM(DM dm, DM fdm)
7930: {
7931:   PetscFunctionBegin;
7934:   if (dm == fdm) fdm = NULL;
7935:   PetscCall(PetscObjectReference((PetscObject)fdm));
7936:   PetscCall(DMDestroy(&dm->fineMesh));
7937:   dm->fineMesh = fdm;
7938:   PetscFunctionReturn(PETSC_SUCCESS);
7939: }

7941: /*@C
7942:   DMAddBoundary - Add a boundary condition to a model represented by a `DM`

7944:   Collective

7946:   Input Parameters:
7947: + dm       - The `DM`, with a `PetscDS` that matches the problem being constrained
7948: . type     - The type of condition, e.g. `DM_BC_ESSENTIAL_ANALYTIC`, `DM_BC_ESSENTIAL_FIELD` (Dirichlet), or `DM_BC_NATURAL` (Neumann)
7949: . name     - The BC name
7950: . label    - The label defining constrained points
7951: . Nv       - The number of `DMLabel` values for constrained points
7952: . values   - An array of values for constrained points
7953: . field    - The field to constrain
7954: . Nc       - The number of constrained field components (0 will constrain all fields)
7955: . comps    - An array of constrained component numbers
7956: . bcFunc   - A pointwise function giving boundary values
7957: . bcFunc_t - A pointwise function giving the time deriative of the boundary values, or NULL
7958: - ctx      - An optional user context for bcFunc

7960:   Output Parameter:
7961: . bd - (Optional) Boundary number

7963:   Options Database Keys:
7964: + -bc_<boundary name> <num>      - Overrides the boundary ids
7965: - -bc_<boundary name>_comp <num> - Overrides the boundary components

7967:   Level: intermediate

7969:   Notes:
7970:   Both bcFunc and bcFunc_t will depend on the boundary condition type. If the type if `DM_BC_ESSENTIAL`, then the calling sequence is\:
7971: .vb
7972:  void bcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[])
7973: .ve

7975:   If the type is `DM_BC_ESSENTIAL_FIELD` or other _FIELD value, then the calling sequence is\:

7977: .vb
7978:   void bcFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
7979:               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
7980:               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
7981:               PetscReal time, const PetscReal x[], PetscScalar bcval[])
7982: .ve
7983: + dim - the spatial dimension
7984: . Nf - the number of fields
7985: . uOff - the offset into u[] and u_t[] for each field
7986: . uOff_x - the offset into u_x[] for each field
7987: . u - each field evaluated at the current point
7988: . u_t - the time derivative of each field evaluated at the current point
7989: . u_x - the gradient of each field evaluated at the current point
7990: . aOff - the offset into a[] and a_t[] for each auxiliary field
7991: . aOff_x - the offset into a_x[] for each auxiliary field
7992: . a - each auxiliary field evaluated at the current point
7993: . a_t - the time derivative of each auxiliary field evaluated at the current point
7994: . a_x - the gradient of auxiliary each field evaluated at the current point
7995: . t - current time
7996: . x - coordinates of the current point
7997: . numConstants - number of constant parameters
7998: . constants - constant parameters
7999: - bcval - output values at the current point

8001: .seealso: [](ch_dmbase), `DM`, `DSGetBoundary()`, `PetscDSAddBoundary()`
8002: @*/
8003: 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)
8004: {
8005:   PetscDS ds;

8007:   PetscFunctionBegin;
8014:   PetscCheck(!dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot add boundary to DM after creating local section");
8015:   PetscCall(DMGetDS(dm, &ds));
8016:   /* Complete label */
8017:   if (label) {
8018:     PetscObject  obj;
8019:     PetscClassId id;

8021:     PetscCall(DMGetField(dm, field, NULL, &obj));
8022:     PetscCall(PetscObjectGetClassId(obj, &id));
8023:     if (id == PETSCFE_CLASSID) {
8024:       DM plex;

8026:       PetscCall(DMConvert(dm, DMPLEX, &plex));
8027:       if (plex) PetscCall(DMPlexLabelComplete(plex, label));
8028:       PetscCall(DMDestroy(&plex));
8029:     }
8030:   }
8031:   PetscCall(PetscDSAddBoundary(ds, type, name, label, Nv, values, field, Nc, comps, bcFunc, bcFunc_t, ctx, bd));
8032:   PetscFunctionReturn(PETSC_SUCCESS);
8033: }

8035: /* TODO Remove this since now the structures are the same */
8036: static PetscErrorCode DMPopulateBoundary(DM dm)
8037: {
8038:   PetscDS     ds;
8039:   DMBoundary *lastnext;
8040:   DSBoundary  dsbound;

8042:   PetscFunctionBegin;
8043:   PetscCall(DMGetDS(dm, &ds));
8044:   dsbound = ds->boundary;
8045:   if (dm->boundary) {
8046:     DMBoundary next = dm->boundary;

8048:     /* quick check to see if the PetscDS has changed */
8049:     if (next->dsboundary == dsbound) PetscFunctionReturn(PETSC_SUCCESS);
8050:     /* the PetscDS has changed: tear down and rebuild */
8051:     while (next) {
8052:       DMBoundary b = next;

8054:       next = b->next;
8055:       PetscCall(PetscFree(b));
8056:     }
8057:     dm->boundary = NULL;
8058:   }

8060:   lastnext = &dm->boundary;
8061:   while (dsbound) {
8062:     DMBoundary dmbound;

8064:     PetscCall(PetscNew(&dmbound));
8065:     dmbound->dsboundary = dsbound;
8066:     dmbound->label      = dsbound->label;
8067:     /* push on the back instead of the front so that it is in the same order as in the PetscDS */
8068:     *lastnext = dmbound;
8069:     lastnext  = &dmbound->next;
8070:     dsbound   = dsbound->next;
8071:   }
8072:   PetscFunctionReturn(PETSC_SUCCESS);
8073: }

8075: /* TODO: missing manual page */
8076: PetscErrorCode DMIsBoundaryPoint(DM dm, PetscInt point, PetscBool *isBd)
8077: {
8078:   DMBoundary b;

8080:   PetscFunctionBegin;
8082:   PetscAssertPointer(isBd, 3);
8083:   *isBd = PETSC_FALSE;
8084:   PetscCall(DMPopulateBoundary(dm));
8085:   b = dm->boundary;
8086:   while (b && !(*isBd)) {
8087:     DMLabel    label = b->label;
8088:     DSBoundary dsb   = b->dsboundary;
8089:     PetscInt   i;

8091:     if (label) {
8092:       for (i = 0; i < dsb->Nv && !(*isBd); ++i) PetscCall(DMLabelStratumHasPoint(label, dsb->values[i], point, isBd));
8093:     }
8094:     b = b->next;
8095:   }
8096:   PetscFunctionReturn(PETSC_SUCCESS);
8097: }

8099: /*@C
8100:   DMProjectFunction - This projects the given function into the function space provided by a `DM`, putting the coefficients in a global vector.

8102:   Collective

8104:   Input Parameters:
8105: + dm    - The `DM`
8106: . time  - The time
8107: . funcs - The coordinate functions to evaluate, one per field
8108: . ctxs  - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
8109: - mode  - The insertion mode for values

8111:   Output Parameter:
8112: . X - vector

8114:   Calling sequence of `funcs`:
8115: + dim  - The spatial dimension
8116: . time - The time at which to sample
8117: . x    - The coordinates
8118: . Nc   - The number of components
8119: . u    - The output field values
8120: - ctx  - optional user-defined function context

8122:   Level: developer

8124:   Developer Notes:
8125:   This API is specific to only particular usage of `DM`

8127:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8129: .seealso: [](ch_dmbase), `DM`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
8130: @*/
8131: 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)
8132: {
8133:   Vec localX;

8135:   PetscFunctionBegin;
8137:   PetscCall(PetscLogEventBegin(DM_ProjectFunction, dm, X, 0, 0));
8138:   PetscCall(DMGetLocalVector(dm, &localX));
8139:   PetscCall(VecSet(localX, 0.));
8140:   PetscCall(DMProjectFunctionLocal(dm, time, funcs, ctxs, mode, localX));
8141:   PetscCall(DMLocalToGlobalBegin(dm, localX, mode, X));
8142:   PetscCall(DMLocalToGlobalEnd(dm, localX, mode, X));
8143:   PetscCall(DMRestoreLocalVector(dm, &localX));
8144:   PetscCall(PetscLogEventEnd(DM_ProjectFunction, dm, X, 0, 0));
8145:   PetscFunctionReturn(PETSC_SUCCESS);
8146: }

8148: /*@C
8149:   DMProjectFunctionLocal - This projects the given function into the function space provided by a `DM`, putting the coefficients in a local vector.

8151:   Not Collective

8153:   Input Parameters:
8154: + dm    - The `DM`
8155: . time  - The time
8156: . funcs - The coordinate functions to evaluate, one per field
8157: . ctxs  - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
8158: - mode  - The insertion mode for values

8160:   Output Parameter:
8161: . localX - vector

8163:   Calling sequence of `funcs`:
8164: + dim  - The spatial dimension
8165: . time - The current timestep
8166: . x    - The coordinates
8167: . Nc   - The number of components
8168: . u    - The output field values
8169: - ctx  - optional user-defined function context

8171:   Level: developer

8173:   Developer Notes:
8174:   This API is specific to only particular usage of `DM`

8176:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8178: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
8179: @*/
8180: 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)
8181: {
8182:   PetscFunctionBegin;
8185:   PetscUseTypeMethod(dm, projectfunctionlocal, time, funcs, ctxs, mode, localX);
8186:   PetscFunctionReturn(PETSC_SUCCESS);
8187: }

8189: /*@C
8190:   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.

8192:   Collective

8194:   Input Parameters:
8195: + dm     - The `DM`
8196: . time   - The time
8197: . numIds - The number of ids
8198: . ids    - The ids
8199: . Nc     - The number of components
8200: . comps  - The components
8201: . label  - The `DMLabel` selecting the portion of the mesh for projection
8202: . funcs  - The coordinate functions to evaluate, one per field
8203: . ctxs   - Optional array of contexts to pass to each coordinate function.  ctxs may be null.
8204: - mode   - The insertion mode for values

8206:   Output Parameter:
8207: . X - vector

8209:   Calling sequence of `funcs`:
8210: + dim  - The spatial dimension
8211: . time - The current timestep
8212: . x    - The coordinates
8213: . Nc   - The number of components
8214: . u    - The output field values
8215: - ctx  - optional user-defined function context

8217:   Level: developer

8219:   Developer Notes:
8220:   This API is specific to only particular usage of `DM`

8222:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8224: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabelLocal()`, `DMComputeL2Diff()`
8225: @*/
8226: 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)
8227: {
8228:   Vec localX;

8230:   PetscFunctionBegin;
8232:   PetscCall(DMGetLocalVector(dm, &localX));
8233:   PetscCall(VecSet(localX, 0.));
8234:   PetscCall(DMProjectFunctionLabelLocal(dm, time, label, numIds, ids, Nc, comps, funcs, ctxs, mode, localX));
8235:   PetscCall(DMLocalToGlobalBegin(dm, localX, mode, X));
8236:   PetscCall(DMLocalToGlobalEnd(dm, localX, mode, X));
8237:   PetscCall(DMRestoreLocalVector(dm, &localX));
8238:   PetscFunctionReturn(PETSC_SUCCESS);
8239: }

8241: /*@C
8242:   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.

8244:   Not Collective

8246:   Input Parameters:
8247: + dm     - The `DM`
8248: . time   - The time
8249: . label  - The `DMLabel` selecting the portion of the mesh for projection
8250: . numIds - The number of ids
8251: . ids    - The ids
8252: . Nc     - The number of components
8253: . comps  - The components
8254: . funcs  - The coordinate functions to evaluate, one per field
8255: . ctxs   - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
8256: - mode   - The insertion mode for values

8258:   Output Parameter:
8259: . localX - vector

8261:   Calling sequence of `funcs`:
8262: + dim  - The spatial dimension
8263: . time - The current time
8264: . x    - The coordinates
8265: . Nc   - The number of components
8266: . u    - The output field values
8267: - ctx  - optional user-defined function context

8269:   Level: developer

8271:   Developer Notes:
8272:   This API is specific to only particular usage of `DM`

8274:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8276: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
8277: @*/
8278: 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)
8279: {
8280:   PetscFunctionBegin;
8283:   PetscUseTypeMethod(dm, projectfunctionlabellocal, time, label, numIds, ids, Nc, comps, funcs, ctxs, mode, localX);
8284:   PetscFunctionReturn(PETSC_SUCCESS);
8285: }

8287: /*@C
8288:   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.

8290:   Not Collective

8292:   Input Parameters:
8293: + dm     - The `DM`
8294: . time   - The time
8295: . localU - The input field vector; may be `NULL` if projection is defined purely by coordinates
8296: . funcs  - The functions to evaluate, one per field
8297: - mode   - The insertion mode for values

8299:   Output Parameter:
8300: . localX - The output vector

8302:   Calling sequence of `funcs`:
8303: + dim          - The spatial dimension
8304: . Nf           - The number of input fields
8305: . NfAux        - The number of input auxiliary fields
8306: . uOff         - The offset of each field in u[]
8307: . uOff_x       - The offset of each field in u_x[]
8308: . u            - The field values at this point in space
8309: . u_t          - The field time derivative at this point in space (or NULL)
8310: . u_x          - The field derivatives at this point in space
8311: . aOff         - The offset of each auxiliary field in u[]
8312: . aOff_x       - The offset of each auxiliary field in u_x[]
8313: . a            - The auxiliary field values at this point in space
8314: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8315: . a_x          - The auxiliary field derivatives at this point in space
8316: . t            - The current time
8317: . x            - The coordinates of this point
8318: . numConstants - The number of constants
8319: . constants    - The value of each constant
8320: - f            - The value of the function at this point in space

8322:   Level: intermediate

8324:   Note:
8325:   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.
8326:   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
8327:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8328:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8330:   Developer Notes:
8331:   This API is specific to only particular usage of `DM`

8333:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8335: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabelLocal()`,
8336: `DMProjectFunction()`, `DMComputeL2Diff()`
8337: @*/
8338: 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)
8339: {
8340:   PetscFunctionBegin;
8344:   PetscUseTypeMethod(dm, projectfieldlocal, time, localU, funcs, mode, localX);
8345:   PetscFunctionReturn(PETSC_SUCCESS);
8346: }

8348: /*@C
8349:   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.

8351:   Not Collective

8353:   Input Parameters:
8354: + dm     - The `DM`
8355: . time   - The time
8356: . label  - The `DMLabel` marking the portion of the domain to output
8357: . numIds - The number of label ids to use
8358: . ids    - The label ids to use for marking
8359: . Nc     - The number of components to set in the output, or `PETSC_DETERMINE` for all components
8360: . comps  - The components to set in the output, or `NULL` for all components
8361: . localU - The input field vector
8362: . funcs  - The functions to evaluate, one per field
8363: - mode   - The insertion mode for values

8365:   Output Parameter:
8366: . localX - The output vector

8368:   Calling sequence of `funcs`:
8369: + dim          - The spatial dimension
8370: . Nf           - The number of input fields
8371: . NfAux        - The number of input auxiliary fields
8372: . uOff         - The offset of each field in u[]
8373: . uOff_x       - The offset of each field in u_x[]
8374: . u            - The field values at this point in space
8375: . u_t          - The field time derivative at this point in space (or NULL)
8376: . u_x          - The field derivatives at this point in space
8377: . aOff         - The offset of each auxiliary field in u[]
8378: . aOff_x       - The offset of each auxiliary field in u_x[]
8379: . a            - The auxiliary field values at this point in space
8380: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8381: . a_x          - The auxiliary field derivatives at this point in space
8382: . t            - The current time
8383: . x            - The coordinates of this point
8384: . numConstants - The number of constants
8385: . constants    - The value of each constant
8386: - f            - The value of the function at this point in space

8388:   Level: intermediate

8390:   Note:
8391:   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.
8392:   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
8393:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8394:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8396:   Developer Notes:
8397:   This API is specific to only particular usage of `DM`

8399:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8401: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabel()`, `DMProjectFunction()`, `DMComputeL2Diff()`
8402: @*/
8403: 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)
8404: {
8405:   PetscFunctionBegin;
8409:   PetscUseTypeMethod(dm, projectfieldlabellocal, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX);
8410:   PetscFunctionReturn(PETSC_SUCCESS);
8411: }

8413: /*@C
8414:   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.

8416:   Not Collective

8418:   Input Parameters:
8419: + dm     - The `DM`
8420: . time   - The time
8421: . label  - The `DMLabel` marking the portion of the domain to output
8422: . numIds - The number of label ids to use
8423: . ids    - The label ids to use for marking
8424: . Nc     - The number of components to set in the output, or `PETSC_DETERMINE` for all components
8425: . comps  - The components to set in the output, or `NULL` for all components
8426: . U      - The input field vector
8427: . funcs  - The functions to evaluate, one per field
8428: - mode   - The insertion mode for values

8430:   Output Parameter:
8431: . X - The output vector

8433:   Calling sequence of `funcs`:
8434: + dim          - The spatial dimension
8435: . Nf           - The number of input fields
8436: . NfAux        - The number of input auxiliary fields
8437: . uOff         - The offset of each field in u[]
8438: . uOff_x       - The offset of each field in u_x[]
8439: . u            - The field values at this point in space
8440: . u_t          - The field time derivative at this point in space (or NULL)
8441: . u_x          - The field derivatives at this point in space
8442: . aOff         - The offset of each auxiliary field in u[]
8443: . aOff_x       - The offset of each auxiliary field in u_x[]
8444: . a            - The auxiliary field values at this point in space
8445: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8446: . a_x          - The auxiliary field derivatives at this point in space
8447: . t            - The current time
8448: . x            - The coordinates of this point
8449: . numConstants - The number of constants
8450: . constants    - The value of each constant
8451: - f            - The value of the function at this point in space

8453:   Level: intermediate

8455:   Note:
8456:   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.
8457:   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
8458:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8459:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8461:   Developer Notes:
8462:   This API is specific to only particular usage of `DM`

8464:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8466: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabelLocal()`, `DMProjectFunction()`, `DMComputeL2Diff()`
8467: @*/
8468: 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)
8469: {
8470:   DM  dmIn;
8471:   Vec localU, localX;

8473:   PetscFunctionBegin;
8475:   PetscCall(VecGetDM(U, &dmIn));
8476:   PetscCall(DMGetLocalVector(dmIn, &localU));
8477:   PetscCall(DMGetLocalVector(dm, &localX));
8478:   PetscCall(VecSet(localX, 0.));
8479:   PetscCall(DMGlobalToLocalBegin(dmIn, U, mode, localU));
8480:   PetscCall(DMGlobalToLocalEnd(dmIn, U, mode, localU));
8481:   PetscCall(DMProjectFieldLabelLocal(dm, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX));
8482:   PetscCall(DMLocalToGlobalBegin(dm, localX, mode, X));
8483:   PetscCall(DMLocalToGlobalEnd(dm, localX, mode, X));
8484:   PetscCall(DMRestoreLocalVector(dm, &localX));
8485:   PetscCall(DMRestoreLocalVector(dmIn, &localU));
8486:   PetscFunctionReturn(PETSC_SUCCESS);
8487: }

8489: /*@C
8490:   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.

8492:   Not Collective

8494:   Input Parameters:
8495: + dm     - The `DM`
8496: . time   - The time
8497: . label  - The `DMLabel` marking the portion of the domain boundary to output
8498: . numIds - The number of label ids to use
8499: . ids    - The label ids to use for marking
8500: . Nc     - The number of components to set in the output, or `PETSC_DETERMINE` for all components
8501: . comps  - The components to set in the output, or `NULL` for all components
8502: . localU - The input field vector
8503: . funcs  - The functions to evaluate, one per field
8504: - mode   - The insertion mode for values

8506:   Output Parameter:
8507: . localX - The output vector

8509:   Calling sequence of `funcs`:
8510: + dim          - The spatial dimension
8511: . Nf           - The number of input fields
8512: . NfAux        - The number of input auxiliary fields
8513: . uOff         - The offset of each field in u[]
8514: . uOff_x       - The offset of each field in u_x[]
8515: . u            - The field values at this point in space
8516: . u_t          - The field time derivative at this point in space (or NULL)
8517: . u_x          - The field derivatives at this point in space
8518: . aOff         - The offset of each auxiliary field in u[]
8519: . aOff_x       - The offset of each auxiliary field in u_x[]
8520: . a            - The auxiliary field values at this point in space
8521: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8522: . a_x          - The auxiliary field derivatives at this point in space
8523: . t            - The current time
8524: . x            - The coordinates of this point
8525: . n            - The face normal
8526: . numConstants - The number of constants
8527: . constants    - The value of each constant
8528: - f            - The value of the function at this point in space

8530:   Level: intermediate

8532:   Note:
8533:   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.
8534:   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
8535:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8536:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8538:   Developer Notes:
8539:   This API is specific to only particular usage of `DM`

8541:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8543: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabelLocal()`, `DMProjectFunction()`, `DMComputeL2Diff()`
8544: @*/
8545: 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)
8546: {
8547:   PetscFunctionBegin;
8551:   PetscUseTypeMethod(dm, projectbdfieldlabellocal, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX);
8552:   PetscFunctionReturn(PETSC_SUCCESS);
8553: }

8555: /*@C
8556:   DMComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.

8558:   Collective

8560:   Input Parameters:
8561: + dm    - The `DM`
8562: . time  - The time
8563: . funcs - The functions to evaluate for each field component
8564: . ctxs  - Optional array of contexts to pass to each function, or NULL.
8565: - X     - The coefficient vector u_h, a global vector

8567:   Output Parameter:
8568: . diff - The diff ||u - u_h||_2

8570:   Level: developer

8572:   Developer Notes:
8573:   This API is specific to only particular usage of `DM`

8575:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8577: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMComputeL2FieldDiff()`, `DMComputeL2GradientDiff()`
8578: @*/
8579: PetscErrorCode DMComputeL2Diff(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
8580: {
8581:   PetscFunctionBegin;
8584:   PetscUseTypeMethod(dm, computel2diff, time, funcs, ctxs, X, diff);
8585:   PetscFunctionReturn(PETSC_SUCCESS);
8586: }

8588: /*@C
8589:   DMComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h.

8591:   Collective

8593:   Input Parameters:
8594: + dm    - The `DM`
8595: . time  - The time
8596: . funcs - The gradient functions to evaluate for each field component
8597: . ctxs  - Optional array of contexts to pass to each function, or NULL.
8598: . X     - The coefficient vector u_h, a global vector
8599: - n     - The vector to project along

8601:   Output Parameter:
8602: . diff - The diff ||(grad u - grad u_h) . n||_2

8604:   Level: developer

8606:   Developer Notes:
8607:   This API is specific to only particular usage of `DM`

8609:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8611: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMComputeL2Diff()`, `DMComputeL2FieldDiff()`
8612: @*/
8613: 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)
8614: {
8615:   PetscFunctionBegin;
8618:   PetscUseTypeMethod(dm, computel2gradientdiff, time, funcs, ctxs, X, n, diff);
8619:   PetscFunctionReturn(PETSC_SUCCESS);
8620: }

8622: /*@C
8623:   DMComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components.

8625:   Collective

8627:   Input Parameters:
8628: + dm    - The `DM`
8629: . time  - The time
8630: . funcs - The functions to evaluate for each field component
8631: . ctxs  - Optional array of contexts to pass to each function, or NULL.
8632: - X     - The coefficient vector u_h, a global vector

8634:   Output Parameter:
8635: . diff - The array of differences, ||u^f - u^f_h||_2

8637:   Level: developer

8639:   Developer Notes:
8640:   This API is specific to only particular usage of `DM`

8642:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8644: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMComputeL2GradientDiff()`
8645: @*/
8646: PetscErrorCode DMComputeL2FieldDiff(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
8647: {
8648:   PetscFunctionBegin;
8651:   PetscUseTypeMethod(dm, computel2fielddiff, time, funcs, ctxs, X, diff);
8652:   PetscFunctionReturn(PETSC_SUCCESS);
8653: }

8655: /*@C
8656:   DMGetNeighbors - Gets an array containing the MPI ranks of all the processes neighbors

8658:   Not Collective

8660:   Input Parameter:
8661: . dm - The `DM`

8663:   Output Parameters:
8664: + nranks - the number of neighbours
8665: - ranks  - the neighbors ranks

8667:   Level: beginner

8669:   Note:
8670:   Do not free the array, it is freed when the `DM` is destroyed.

8672: .seealso: [](ch_dmbase), `DM`, `DMDAGetNeighbors()`, `PetscSFGetRootRanks()`
8673: @*/
8674: PetscErrorCode DMGetNeighbors(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
8675: {
8676:   PetscFunctionBegin;
8678:   PetscUseTypeMethod(dm, getneighbors, nranks, ranks);
8679:   PetscFunctionReturn(PETSC_SUCCESS);
8680: }

8682: #include <petsc/private/matimpl.h>

8684: /*
8685:     Converts the input vector to a ghosted vector and then calls the standard coloring code.
8686:     This must be a different function because it requires DM which is not defined in the Mat library
8687: */
8688: static PetscErrorCode MatFDColoringApply_AIJDM(Mat J, MatFDColoring coloring, Vec x1, void *sctx)
8689: {
8690:   PetscFunctionBegin;
8691:   if (coloring->ctype == IS_COLORING_LOCAL) {
8692:     Vec x1local;
8693:     DM  dm;
8694:     PetscCall(MatGetDM(J, &dm));
8695:     PetscCheck(dm, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_INCOMP, "IS_COLORING_LOCAL requires a DM");
8696:     PetscCall(DMGetLocalVector(dm, &x1local));
8697:     PetscCall(DMGlobalToLocalBegin(dm, x1, INSERT_VALUES, x1local));
8698:     PetscCall(DMGlobalToLocalEnd(dm, x1, INSERT_VALUES, x1local));
8699:     x1 = x1local;
8700:   }
8701:   PetscCall(MatFDColoringApply_AIJ(J, coloring, x1, sctx));
8702:   if (coloring->ctype == IS_COLORING_LOCAL) {
8703:     DM dm;
8704:     PetscCall(MatGetDM(J, &dm));
8705:     PetscCall(DMRestoreLocalVector(dm, &x1));
8706:   }
8707:   PetscFunctionReturn(PETSC_SUCCESS);
8708: }

8710: /*@
8711:   MatFDColoringUseDM - allows a `MatFDColoring` object to use the `DM` associated with the matrix to compute a `IS_COLORING_LOCAL` coloring

8713:   Input Parameters:
8714: + coloring   - The matrix to get the `DM` from
8715: - fdcoloring - the `MatFDColoring` object

8717:   Level: advanced

8719:   Developer Note:
8720:   This routine exists because the PETSc `Mat` library does not know about the `DM` objects

8722: .seealso: [](ch_dmbase), `DM`, `MatFDColoring`, `MatFDColoringCreate()`, `ISColoringType`
8723: @*/
8724: PetscErrorCode MatFDColoringUseDM(Mat coloring, MatFDColoring fdcoloring)
8725: {
8726:   PetscFunctionBegin;
8727:   coloring->ops->fdcoloringapply = MatFDColoringApply_AIJDM;
8728:   PetscFunctionReturn(PETSC_SUCCESS);
8729: }

8731: /*@
8732:   DMGetCompatibility - determine if two `DM`s are compatible

8734:   Collective

8736:   Input Parameters:
8737: + dm1 - the first `DM`
8738: - dm2 - the second `DM`

8740:   Output Parameters:
8741: + compatible - whether or not the two `DM`s are compatible
8742: - set        - whether or not the compatible value was actually determined and set

8744:   Level: advanced

8746:   Notes:
8747:   Two `DM`s are deemed compatible if they represent the same parallel decomposition
8748:   of the same topology. This implies that the section (field data) on one
8749:   "makes sense" with respect to the topology and parallel decomposition of the other.
8750:   Loosely speaking, compatible `DM`s represent the same domain and parallel
8751:   decomposition, but hold different data.

8753:   Typically, one would confirm compatibility if intending to simultaneously iterate
8754:   over a pair of vectors obtained from different `DM`s.

8756:   For example, two `DMDA` objects are compatible if they have the same local
8757:   and global sizes and the same stencil width. They can have different numbers
8758:   of degrees of freedom per node. Thus, one could use the node numbering from
8759:   either `DM` in bounds for a loop over vectors derived from either `DM`.

8761:   Consider the operation of summing data living on a 2-dof `DMDA` to data living
8762:   on a 1-dof `DMDA`, which should be compatible, as in the following snippet.
8763: .vb
8764:   ...
8765:   PetscCall(DMGetCompatibility(da1,da2,&compatible,&set));
8766:   if (set && compatible)  {
8767:     PetscCall(DMDAVecGetArrayDOF(da1,vec1,&arr1));
8768:     PetscCall(DMDAVecGetArrayDOF(da2,vec2,&arr2));
8769:     PetscCall(DMDAGetCorners(da1,&x,&y,NULL,&m,&n,NULL));
8770:     for (j=y; j<y+n; ++j) {
8771:       for (i=x; i<x+m, ++i) {
8772:         arr1[j][i][0] = arr2[j][i][0] + arr2[j][i][1];
8773:       }
8774:     }
8775:     PetscCall(DMDAVecRestoreArrayDOF(da1,vec1,&arr1));
8776:     PetscCall(DMDAVecRestoreArrayDOF(da2,vec2,&arr2));
8777:   } else {
8778:     SETERRQ(PetscObjectComm((PetscObject)da1,PETSC_ERR_ARG_INCOMP,"DMDA objects incompatible");
8779:   }
8780:   ...
8781: .ve

8783:   Checking compatibility might be expensive for a given implementation of `DM`,
8784:   or might be impossible to unambiguously confirm or deny. For this reason,
8785:   this function may decline to determine compatibility, and hence users should
8786:   always check the "set" output parameter.

8788:   A `DM` is always compatible with itself.

8790:   In the current implementation, `DM`s which live on "unequal" communicators
8791:   (MPI_UNEQUAL in the terminology of MPI_Comm_compare()) are always deemed
8792:   incompatible.

8794:   This function is labeled "Collective," as information about all subdomains
8795:   is required on each rank. However, in `DM` implementations which store all this
8796:   information locally, this function may be merely "Logically Collective".

8798:   Developer Note:
8799:   Compatibility is assumed to be a symmetric concept; `DM` A is compatible with `DM` B
8800:   iff B is compatible with A. Thus, this function checks the implementations
8801:   of both dm and dmc (if they are of different types), attempting to determine
8802:   compatibility. It is left to `DM` implementers to ensure that symmetry is
8803:   preserved. The simplest way to do this is, when implementing type-specific
8804:   logic for this function, is to check for existing logic in the implementation
8805:   of other `DM` types and let *set = PETSC_FALSE if found.

8807: .seealso: [](ch_dmbase), `DM`, `DMDACreateCompatibleDMDA()`, `DMStagCreateCompatibleDMStag()`
8808: @*/
8809: PetscErrorCode DMGetCompatibility(DM dm1, DM dm2, PetscBool *compatible, PetscBool *set)
8810: {
8811:   PetscMPIInt compareResult;
8812:   DMType      type, type2;
8813:   PetscBool   sameType;

8815:   PetscFunctionBegin;

8819:   /* Declare a DM compatible with itself */
8820:   if (dm1 == dm2) {
8821:     *set        = PETSC_TRUE;
8822:     *compatible = PETSC_TRUE;
8823:     PetscFunctionReturn(PETSC_SUCCESS);
8824:   }

8826:   /* Declare a DM incompatible with a DM that lives on an "unequal"
8827:      communicator. Note that this does not preclude compatibility with
8828:      DMs living on "congruent" or "similar" communicators, but this must be
8829:      determined by the implementation-specific logic */
8830:   PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dm1), PetscObjectComm((PetscObject)dm2), &compareResult));
8831:   if (compareResult == MPI_UNEQUAL) {
8832:     *set        = PETSC_TRUE;
8833:     *compatible = PETSC_FALSE;
8834:     PetscFunctionReturn(PETSC_SUCCESS);
8835:   }

8837:   /* Pass to the implementation-specific routine, if one exists. */
8838:   if (dm1->ops->getcompatibility) {
8839:     PetscUseTypeMethod(dm1, getcompatibility, dm2, compatible, set);
8840:     if (*set) PetscFunctionReturn(PETSC_SUCCESS);
8841:   }

8843:   /* If dm1 and dm2 are of different types, then attempt to check compatibility
8844:      with an implementation of this function from dm2 */
8845:   PetscCall(DMGetType(dm1, &type));
8846:   PetscCall(DMGetType(dm2, &type2));
8847:   PetscCall(PetscStrcmp(type, type2, &sameType));
8848:   if (!sameType && dm2->ops->getcompatibility) {
8849:     PetscUseTypeMethod(dm2, getcompatibility, dm1, compatible, set); /* Note argument order */
8850:   } else {
8851:     *set = PETSC_FALSE;
8852:   }
8853:   PetscFunctionReturn(PETSC_SUCCESS);
8854: }

8856: /*@C
8857:   DMMonitorSet - Sets an additional monitor function that is to be used after a solve to monitor discretization performance.

8859:   Logically Collective

8861:   Input Parameters:
8862: + dm             - the `DM`
8863: . f              - the monitor function
8864: . mctx           - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired)
8865: - monitordestroy - [optional] routine that frees monitor context (may be `NULL`)

8867:   Options Database Key:
8868: . -dm_monitor_cancel - cancels all monitors that have been hardwired into a code by calls to `DMMonitorSet()`, but
8869:                             does not cancel those set via the options database.

8871:   Level: intermediate

8873:   Note:
8874:   Several different monitoring routines may be set by calling
8875:   `DMMonitorSet()` multiple times or with `DMMonitorSetFromOptions()`; all will be called in the
8876:   order in which they were set.

8878:   Fortran Note:
8879:   Only a single monitor function can be set for each `DM` object

8881:   Developer Note:
8882:   This API has a generic name but seems specific to a very particular aspect of the use of `DM`

8884: .seealso: [](ch_dmbase), `DM`, `DMMonitorCancel()`, `DMMonitorSetFromOptions()`, `DMMonitor()`
8885: @*/
8886: PetscErrorCode DMMonitorSet(DM dm, PetscErrorCode (*f)(DM, void *), void *mctx, PetscErrorCode (*monitordestroy)(void **))
8887: {
8888:   PetscInt m;

8890:   PetscFunctionBegin;
8892:   for (m = 0; m < dm->numbermonitors; ++m) {
8893:     PetscBool identical;

8895:     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, mctx, monitordestroy, (PetscErrorCode(*)(void))dm->monitor[m], dm->monitorcontext[m], dm->monitordestroy[m], &identical));
8896:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
8897:   }
8898:   PetscCheck(dm->numbermonitors < MAXDMMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
8899:   dm->monitor[dm->numbermonitors]          = f;
8900:   dm->monitordestroy[dm->numbermonitors]   = monitordestroy;
8901:   dm->monitorcontext[dm->numbermonitors++] = (void *)mctx;
8902:   PetscFunctionReturn(PETSC_SUCCESS);
8903: }

8905: /*@
8906:   DMMonitorCancel - Clears all the monitor functions for a `DM` object.

8908:   Logically Collective

8910:   Input Parameter:
8911: . dm - the DM

8913:   Options Database Key:
8914: . -dm_monitor_cancel - cancels all monitors that have been hardwired
8915:   into a code by calls to `DMonitorSet()`, but does not cancel those
8916:   set via the options database

8918:   Level: intermediate

8920:   Note:
8921:   There is no way to clear one specific monitor from a `DM` object.

8923: .seealso: [](ch_dmbase), `DM`, `DMMonitorSet()`, `DMMonitorSetFromOptions()`, `DMMonitor()`
8924: @*/
8925: PetscErrorCode DMMonitorCancel(DM dm)
8926: {
8927:   PetscInt m;

8929:   PetscFunctionBegin;
8931:   for (m = 0; m < dm->numbermonitors; ++m) {
8932:     if (dm->monitordestroy[m]) PetscCall((*dm->monitordestroy[m])(&dm->monitorcontext[m]));
8933:   }
8934:   dm->numbermonitors = 0;
8935:   PetscFunctionReturn(PETSC_SUCCESS);
8936: }

8938: /*@C
8939:   DMMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user

8941:   Collective

8943:   Input Parameters:
8944: + dm           - `DM` object you wish to monitor
8945: . name         - the monitor type one is seeking
8946: . help         - message indicating what monitoring is done
8947: . manual       - manual page for the monitor
8948: . monitor      - the monitor function
8949: - 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

8951:   Output Parameter:
8952: . flg - Flag set if the monitor was created

8954:   Level: developer

8956: .seealso: [](ch_dmbase), `DM`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
8957:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
8958:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
8959:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
8960:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
8961:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
8962:           `PetscOptionsFList()`, `PetscOptionsEList()`, `DMMonitor()`, `DMMonitorSet()`
8963: @*/
8964: PetscErrorCode DMMonitorSetFromOptions(DM dm, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(DM, void *), PetscErrorCode (*monitorsetup)(DM, PetscViewerAndFormat *), PetscBool *flg)
8965: {
8966:   PetscViewer       viewer;
8967:   PetscViewerFormat format;

8969:   PetscFunctionBegin;
8971:   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)dm), ((PetscObject)dm)->options, ((PetscObject)dm)->prefix, name, &viewer, &format, flg));
8972:   if (*flg) {
8973:     PetscViewerAndFormat *vf;

8975:     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
8976:     PetscCall(PetscViewerDestroy(&viewer));
8977:     if (monitorsetup) PetscCall((*monitorsetup)(dm, vf));
8978:     PetscCall(DMMonitorSet(dm, (PetscErrorCode(*)(DM, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
8979:   }
8980:   PetscFunctionReturn(PETSC_SUCCESS);
8981: }

8983: /*@
8984:   DMMonitor - runs the user provided monitor routines, if they exist

8986:   Collective

8988:   Input Parameter:
8989: . dm - The `DM`

8991:   Level: developer

8993:   Developer Note:
8994:   Note should indicate when during the life of the `DM` the monitor is run. It appears to be
8995:   related to the discretization process seems rather specialized since some `DM` have no
8996:   concept of discretization.

8998: .seealso: [](ch_dmbase), `DM`, `DMMonitorSet()`, `DMMonitorSetFromOptions()`
8999: @*/
9000: PetscErrorCode DMMonitor(DM dm)
9001: {
9002:   PetscInt m;

9004:   PetscFunctionBegin;
9005:   if (!dm) PetscFunctionReturn(PETSC_SUCCESS);
9007:   for (m = 0; m < dm->numbermonitors; ++m) PetscCall((*dm->monitor[m])(dm, dm->monitorcontext[m]));
9008:   PetscFunctionReturn(PETSC_SUCCESS);
9009: }

9011: /*@
9012:   DMComputeError - Computes the error assuming the user has provided the exact solution functions

9014:   Collective

9016:   Input Parameters:
9017: + dm  - The `DM`
9018: - sol - The solution vector

9020:   Input/Output Parameter:
9021: . errors - An array of length Nf, the number of fields, or `NULL` for no output; on output
9022:            contains the error in each field

9024:   Output Parameter:
9025: . errorVec - A vector to hold the cellwise error (may be `NULL`)

9027:   Level: developer

9029:   Note:
9030:   The exact solutions come from the `PetscDS` object, and the time comes from `DMGetOutputSequenceNumber()`.

9032: .seealso: [](ch_dmbase), `DM`, `DMMonitorSet()`, `DMGetRegionNumDS()`, `PetscDSGetExactSolution()`, `DMGetOutputSequenceNumber()`
9033: @*/
9034: PetscErrorCode DMComputeError(DM dm, Vec sol, PetscReal errors[], Vec *errorVec)
9035: {
9036:   PetscErrorCode (**exactSol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
9037:   void    **ctxs;
9038:   PetscReal time;
9039:   PetscInt  Nf, f, Nds, s;

9041:   PetscFunctionBegin;
9042:   PetscCall(DMGetNumFields(dm, &Nf));
9043:   PetscCall(PetscCalloc2(Nf, &exactSol, Nf, &ctxs));
9044:   PetscCall(DMGetNumDS(dm, &Nds));
9045:   for (s = 0; s < Nds; ++s) {
9046:     PetscDS         ds;
9047:     DMLabel         label;
9048:     IS              fieldIS;
9049:     const PetscInt *fields;
9050:     PetscInt        dsNf;

9052:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL));
9053:     PetscCall(PetscDSGetNumFields(ds, &dsNf));
9054:     if (fieldIS) PetscCall(ISGetIndices(fieldIS, &fields));
9055:     for (f = 0; f < dsNf; ++f) {
9056:       const PetscInt field = fields[f];
9057:       PetscCall(PetscDSGetExactSolution(ds, field, &exactSol[field], &ctxs[field]));
9058:     }
9059:     if (fieldIS) PetscCall(ISRestoreIndices(fieldIS, &fields));
9060:   }
9061:   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);
9062:   PetscCall(DMGetOutputSequenceNumber(dm, NULL, &time));
9063:   if (errors) PetscCall(DMComputeL2FieldDiff(dm, time, exactSol, ctxs, sol, errors));
9064:   if (errorVec) {
9065:     DM             edm;
9066:     DMPolytopeType ct;
9067:     PetscBool      simplex;
9068:     PetscInt       dim, cStart, Nf;

9070:     PetscCall(DMClone(dm, &edm));
9071:     PetscCall(DMGetDimension(edm, &dim));
9072:     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
9073:     PetscCall(DMPlexGetCellType(dm, cStart, &ct));
9074:     simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
9075:     PetscCall(DMGetNumFields(dm, &Nf));
9076:     for (f = 0; f < Nf; ++f) {
9077:       PetscFE         fe, efe;
9078:       PetscQuadrature q;
9079:       const char     *name;

9081:       PetscCall(DMGetField(dm, f, NULL, (PetscObject *)&fe));
9082:       PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, Nf, simplex, 0, PETSC_DETERMINE, &efe));
9083:       PetscCall(PetscObjectGetName((PetscObject)fe, &name));
9084:       PetscCall(PetscObjectSetName((PetscObject)efe, name));
9085:       PetscCall(PetscFEGetQuadrature(fe, &q));
9086:       PetscCall(PetscFESetQuadrature(efe, q));
9087:       PetscCall(DMSetField(edm, f, NULL, (PetscObject)efe));
9088:       PetscCall(PetscFEDestroy(&efe));
9089:     }
9090:     PetscCall(DMCreateDS(edm));

9092:     PetscCall(DMCreateGlobalVector(edm, errorVec));
9093:     PetscCall(PetscObjectSetName((PetscObject)*errorVec, "Error"));
9094:     PetscCall(DMPlexComputeL2DiffVec(dm, time, exactSol, ctxs, sol, *errorVec));
9095:     PetscCall(DMDestroy(&edm));
9096:   }
9097:   PetscCall(PetscFree2(exactSol, ctxs));
9098:   PetscFunctionReturn(PETSC_SUCCESS);
9099: }

9101: /*@
9102:   DMGetNumAuxiliaryVec - Get the number of auxiliary vectors associated with this `DM`

9104:   Not Collective

9106:   Input Parameter:
9107: . dm - The `DM`

9109:   Output Parameter:
9110: . numAux - The number of auxiliary data vectors

9112:   Level: advanced

9114: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMSetAuxiliaryVec()`, `DMGetAuxiliaryLabels()`, `DMGetAuxiliaryVec()`
9115: @*/
9116: PetscErrorCode DMGetNumAuxiliaryVec(DM dm, PetscInt *numAux)
9117: {
9118:   PetscFunctionBegin;
9120:   PetscCall(PetscHMapAuxGetSize(dm->auxData, numAux));
9121:   PetscFunctionReturn(PETSC_SUCCESS);
9122: }

9124: /*@
9125:   DMGetAuxiliaryVec - Get the auxiliary vector for region specified by the given label and value, and equation part

9127:   Not Collective

9129:   Input Parameters:
9130: + dm    - The `DM`
9131: . label - The `DMLabel`
9132: . value - The label value indicating the region
9133: - part  - The equation part, or 0 if unused

9135:   Output Parameter:
9136: . aux - The `Vec` holding auxiliary field data

9138:   Level: advanced

9140:   Note:
9141:   If no auxiliary vector is found for this (label, value), (NULL, 0, 0) is checked as well.

9143: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMSetAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryLabels()`
9144: @*/
9145: PetscErrorCode DMGetAuxiliaryVec(DM dm, DMLabel label, PetscInt value, PetscInt part, Vec *aux)
9146: {
9147:   PetscHashAuxKey key, wild = {NULL, 0, 0};
9148:   PetscBool       has;

9150:   PetscFunctionBegin;
9153:   key.label = label;
9154:   key.value = value;
9155:   key.part  = part;
9156:   PetscCall(PetscHMapAuxHas(dm->auxData, key, &has));
9157:   if (has) PetscCall(PetscHMapAuxGet(dm->auxData, key, aux));
9158:   else PetscCall(PetscHMapAuxGet(dm->auxData, wild, aux));
9159:   PetscFunctionReturn(PETSC_SUCCESS);
9160: }

9162: /*@
9163:   DMSetAuxiliaryVec - Set an auxiliary vector for region specified by the given label and value, and equation part

9165:   Not Collective because auxiliary vectors are not parallel

9167:   Input Parameters:
9168: + dm    - The `DM`
9169: . label - The `DMLabel`
9170: . value - The label value indicating the region
9171: . part  - The equation part, or 0 if unused
9172: - aux   - The `Vec` holding auxiliary field data

9174:   Level: advanced

9176: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMGetAuxiliaryLabels()`, `DMCopyAuxiliaryVec()`
9177: @*/
9178: PetscErrorCode DMSetAuxiliaryVec(DM dm, DMLabel label, PetscInt value, PetscInt part, Vec aux)
9179: {
9180:   Vec             old;
9181:   PetscHashAuxKey key;

9183:   PetscFunctionBegin;
9186:   key.label = label;
9187:   key.value = value;
9188:   key.part  = part;
9189:   PetscCall(PetscHMapAuxGet(dm->auxData, key, &old));
9190:   PetscCall(PetscObjectReference((PetscObject)aux));
9191:   if (!aux) PetscCall(PetscHMapAuxDel(dm->auxData, key));
9192:   else PetscCall(PetscHMapAuxSet(dm->auxData, key, aux));
9193:   PetscCall(VecDestroy(&old));
9194:   PetscFunctionReturn(PETSC_SUCCESS);
9195: }

9197: /*@
9198:   DMGetAuxiliaryLabels - Get the labels, values, and parts for all auxiliary vectors in this `DM`

9200:   Not Collective

9202:   Input Parameter:
9203: . dm - The `DM`

9205:   Output Parameters:
9206: + labels - The `DMLabel`s for each `Vec`
9207: . values - The label values for each `Vec`
9208: - parts  - The equation parts for each `Vec`

9210:   Level: advanced

9212:   Note:
9213:   The arrays passed in must be at least as large as `DMGetNumAuxiliaryVec()`.

9215: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMSetAuxiliaryVec()`, `DMCopyAuxiliaryVec()`
9216: @*/
9217: PetscErrorCode DMGetAuxiliaryLabels(DM dm, DMLabel labels[], PetscInt values[], PetscInt parts[])
9218: {
9219:   PetscHashAuxKey *keys;
9220:   PetscInt         n, i, off = 0;

9222:   PetscFunctionBegin;
9224:   PetscAssertPointer(labels, 2);
9225:   PetscAssertPointer(values, 3);
9226:   PetscAssertPointer(parts, 4);
9227:   PetscCall(DMGetNumAuxiliaryVec(dm, &n));
9228:   PetscCall(PetscMalloc1(n, &keys));
9229:   PetscCall(PetscHMapAuxGetKeys(dm->auxData, &off, keys));
9230:   for (i = 0; i < n; ++i) {
9231:     labels[i] = keys[i].label;
9232:     values[i] = keys[i].value;
9233:     parts[i]  = keys[i].part;
9234:   }
9235:   PetscCall(PetscFree(keys));
9236:   PetscFunctionReturn(PETSC_SUCCESS);
9237: }

9239: /*@
9240:   DMCopyAuxiliaryVec - Copy the auxiliary vector data on a `DM` to a new `DM`

9242:   Not Collective

9244:   Input Parameter:
9245: . dm - The `DM`

9247:   Output Parameter:
9248: . dmNew - The new `DM`, now with the same auxiliary data

9250:   Level: advanced

9252:   Note:
9253:   This is a shallow copy of the auxiliary vectors

9255: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMSetAuxiliaryVec()`
9256: @*/
9257: PetscErrorCode DMCopyAuxiliaryVec(DM dm, DM dmNew)
9258: {
9259:   PetscFunctionBegin;
9262:   if (dm == dmNew) PetscFunctionReturn(PETSC_SUCCESS);
9263:   PetscCall(DMClearAuxiliaryVec(dmNew));

9265:   PetscCall(PetscHMapAuxDestroy(&dmNew->auxData));
9266:   PetscCall(PetscHMapAuxDuplicate(dm->auxData, &dmNew->auxData));
9267:   {
9268:     Vec     *auxData;
9269:     PetscInt n, i, off = 0;

9271:     PetscCall(PetscHMapAuxGetSize(dmNew->auxData, &n));
9272:     PetscCall(PetscMalloc1(n, &auxData));
9273:     PetscCall(PetscHMapAuxGetVals(dmNew->auxData, &off, auxData));
9274:     for (i = 0; i < n; ++i) PetscCall(PetscObjectReference((PetscObject)auxData[i]));
9275:     PetscCall(PetscFree(auxData));
9276:   }
9277:   PetscFunctionReturn(PETSC_SUCCESS);
9278: }

9280: /*@
9281:   DMClearAuxiliaryVec - Destroys the auxiliary vector information and creates a new empty one

9283:   Not Collective

9285:   Input Parameter:
9286: . dm - The `DM`

9288:   Level: advanced

9290: .seealso: [](ch_dmbase), `DM`, `DMCopyAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMSetAuxiliaryVec()`
9291: @*/
9292: PetscErrorCode DMClearAuxiliaryVec(DM dm)
9293: {
9294:   Vec     *auxData;
9295:   PetscInt n, i, off = 0;

9297:   PetscFunctionBegin;
9298:   PetscCall(PetscHMapAuxGetSize(dm->auxData, &n));
9299:   PetscCall(PetscMalloc1(n, &auxData));
9300:   PetscCall(PetscHMapAuxGetVals(dm->auxData, &off, auxData));
9301:   for (i = 0; i < n; ++i) PetscCall(VecDestroy(&auxData[i]));
9302:   PetscCall(PetscFree(auxData));
9303:   PetscCall(PetscHMapAuxDestroy(&dm->auxData));
9304:   PetscCall(PetscHMapAuxCreate(&dm->auxData));
9305:   PetscFunctionReturn(PETSC_SUCCESS);
9306: }

9308: /*@
9309:   DMPolytopeMatchOrientation - Determine an orientation (transformation) that takes the source face arrangement to the target face arrangement

9311:   Not Collective

9313:   Input Parameters:
9314: + ct         - The `DMPolytopeType`
9315: . sourceCone - The source arrangement of faces
9316: - targetCone - The target arrangement of faces

9318:   Output Parameters:
9319: + ornt  - The orientation (transformation) which will take the source arrangement to the target arrangement
9320: - found - Flag indicating that a suitable orientation was found

9322:   Level: advanced

9324:   Note:
9325:   An arrangement is a face order combined with an orientation for each face

9327:   Each orientation (transformation) is labeled with an integer from negative `DMPolytopeTypeGetNumArrangements(ct)`/2 to `DMPolytopeTypeGetNumArrangements(ct)`/2
9328:   that labels each arrangement (face ordering plus orientation for each face).

9330:   See `DMPolytopeMatchVertexOrientation()` to find a new vertex orientation that takes the source vertex arrangement to the target vertex arrangement

9332: .seealso: [](ch_dmbase), `DM`, `DMPolytopeGetOrientation()`, `DMPolytopeMatchVertexOrientation()`, `DMPolytopeGetVertexOrientation()`
9333: @*/
9334: PetscErrorCode DMPolytopeMatchOrientation(DMPolytopeType ct, const PetscInt sourceCone[], const PetscInt targetCone[], PetscInt *ornt, PetscBool *found)
9335: {
9336:   const PetscInt cS = DMPolytopeTypeGetConeSize(ct);
9337:   const PetscInt nO = DMPolytopeTypeGetNumArrangements(ct) / 2;
9338:   PetscInt       o, c;

9340:   PetscFunctionBegin;
9341:   if (!nO) {
9342:     *ornt  = 0;
9343:     *found = PETSC_TRUE;
9344:     PetscFunctionReturn(PETSC_SUCCESS);
9345:   }
9346:   for (o = -nO; o < nO; ++o) {
9347:     const PetscInt *arr = DMPolytopeTypeGetArrangement(ct, o);

9349:     for (c = 0; c < cS; ++c)
9350:       if (sourceCone[arr[c * 2]] != targetCone[c]) break;
9351:     if (c == cS) {
9352:       *ornt = o;
9353:       break;
9354:     }
9355:   }
9356:   *found = o == nO ? PETSC_FALSE : PETSC_TRUE;
9357:   PetscFunctionReturn(PETSC_SUCCESS);
9358: }

9360: /*@
9361:   DMPolytopeGetOrientation - Determine an orientation (transformation) that takes the source face arrangement to the target face arrangement

9363:   Not Collective

9365:   Input Parameters:
9366: + ct         - The `DMPolytopeType`
9367: . sourceCone - The source arrangement of faces
9368: - targetCone - The target arrangement of faces

9370:   Output Parameter:
9371: . ornt - The orientation (transformation) which will take the source arrangement to the target arrangement

9373:   Level: advanced

9375:   Note:
9376:   This function is the same as `DMPolytopeMatchOrientation()` except it will generate an error if no suitable orientation can be found.

9378:   Developer Note:
9379:   It is unclear why this function needs to exist since one can simply call `DMPolytopeMatchOrientation()` and error if none is found

9381: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMPolytopeMatchOrientation()`, `DMPolytopeGetVertexOrientation()`, `DMPolytopeMatchVertexOrientation()`
9382: @*/
9383: PetscErrorCode DMPolytopeGetOrientation(DMPolytopeType ct, const PetscInt sourceCone[], const PetscInt targetCone[], PetscInt *ornt)
9384: {
9385:   PetscBool found;

9387:   PetscFunctionBegin;
9388:   PetscCall(DMPolytopeMatchOrientation(ct, sourceCone, targetCone, ornt, &found));
9389:   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find orientation for %s", DMPolytopeTypes[ct]);
9390:   PetscFunctionReturn(PETSC_SUCCESS);
9391: }

9393: /*@
9394:   DMPolytopeMatchVertexOrientation - Determine an orientation (transformation) that takes the source vertex arrangement to the target vertex arrangement

9396:   Not Collective

9398:   Input Parameters:
9399: + ct         - The `DMPolytopeType`
9400: . sourceVert - The source arrangement of vertices
9401: - targetVert - The target arrangement of vertices

9403:   Output Parameters:
9404: + ornt  - The orientation (transformation) which will take the source arrangement to the target arrangement
9405: - found - Flag indicating that a suitable orientation was found

9407:   Level: advanced

9409:   Notes:
9410:   An arrangement is a vertex order

9412:   Each orientation (transformation) is labeled with an integer from negative `DMPolytopeTypeGetNumArrangements(ct)`/2 to `DMPolytopeTypeGetNumArrangements(ct)`/2
9413:   that labels each arrangement (vertex ordering).

9415:   See `DMPolytopeMatchOrientation()` to find a new face orientation that takes the source face arrangement to the target face arrangement

9417: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMPolytopeGetOrientation()`, `DMPolytopeMatchOrientation()`, `DMPolytopeTypeGetNumVertices()`, `DMPolytopeTypeGetVertexArrangement()`
9418: @*/
9419: PetscErrorCode DMPolytopeMatchVertexOrientation(DMPolytopeType ct, const PetscInt sourceVert[], const PetscInt targetVert[], PetscInt *ornt, PetscBool *found)
9420: {
9421:   const PetscInt cS = DMPolytopeTypeGetNumVertices(ct);
9422:   const PetscInt nO = DMPolytopeTypeGetNumArrangements(ct) / 2;
9423:   PetscInt       o, c;

9425:   PetscFunctionBegin;
9426:   if (!nO) {
9427:     *ornt  = 0;
9428:     *found = PETSC_TRUE;
9429:     PetscFunctionReturn(PETSC_SUCCESS);
9430:   }
9431:   for (o = -nO; o < nO; ++o) {
9432:     const PetscInt *arr = DMPolytopeTypeGetVertexArrangement(ct, o);

9434:     for (c = 0; c < cS; ++c)
9435:       if (sourceVert[arr[c]] != targetVert[c]) break;
9436:     if (c == cS) {
9437:       *ornt = o;
9438:       break;
9439:     }
9440:   }
9441:   *found = o == nO ? PETSC_FALSE : PETSC_TRUE;
9442:   PetscFunctionReturn(PETSC_SUCCESS);
9443: }

9445: /*@
9446:   DMPolytopeGetVertexOrientation - Determine an orientation (transformation) that takes the source vertex arrangement to the target vertex arrangement

9448:   Not Collective

9450:   Input Parameters:
9451: + ct         - The `DMPolytopeType`
9452: . sourceCone - The source arrangement of vertices
9453: - targetCone - The target arrangement of vertices

9455:   Output Parameter:
9456: . ornt - The orientation (transformation) which will take the source arrangement to the target arrangement

9458:   Level: advanced

9460:   Note:
9461:   This function is the same as `DMPolytopeMatchVertexOrientation()` except it errors if not orientation is possible.

9463:   Developer Note:
9464:   It is unclear why this function needs to exist since one can simply call `DMPolytopeMatchVertexOrientation()` and error if none is found

9466: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMPolytopeMatchVertexOrientation()`, `DMPolytopeGetOrientation()`
9467: @*/
9468: PetscErrorCode DMPolytopeGetVertexOrientation(DMPolytopeType ct, const PetscInt sourceCone[], const PetscInt targetCone[], PetscInt *ornt)
9469: {
9470:   PetscBool found;

9472:   PetscFunctionBegin;
9473:   PetscCall(DMPolytopeMatchVertexOrientation(ct, sourceCone, targetCone, ornt, &found));
9474:   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find orientation for %s", DMPolytopeTypes[ct]);
9475:   PetscFunctionReturn(PETSC_SUCCESS);
9476: }

9478: /*@
9479:   DMPolytopeInCellTest - Check whether a point lies inside the reference cell of given type

9481:   Not Collective

9483:   Input Parameters:
9484: + ct    - The `DMPolytopeType`
9485: - point - Coordinates of the point

9487:   Output Parameter:
9488: . inside - Flag indicating whether the point is inside the reference cell of given type

9490:   Level: advanced

9492: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMLocatePoints()`
9493: @*/
9494: PetscErrorCode DMPolytopeInCellTest(DMPolytopeType ct, const PetscReal point[], PetscBool *inside)
9495: {
9496:   PetscReal sum = 0.0;
9497:   PetscInt  d;

9499:   PetscFunctionBegin;
9500:   *inside = PETSC_TRUE;
9501:   switch (ct) {
9502:   case DM_POLYTOPE_TRIANGLE:
9503:   case DM_POLYTOPE_TETRAHEDRON:
9504:     for (d = 0; d < DMPolytopeTypeGetDim(ct); ++d) {
9505:       if (point[d] < -1.0) {
9506:         *inside = PETSC_FALSE;
9507:         break;
9508:       }
9509:       sum += point[d];
9510:     }
9511:     if (sum > PETSC_SMALL) {
9512:       *inside = PETSC_FALSE;
9513:       break;
9514:     }
9515:     break;
9516:   case DM_POLYTOPE_QUADRILATERAL:
9517:   case DM_POLYTOPE_HEXAHEDRON:
9518:     for (d = 0; d < DMPolytopeTypeGetDim(ct); ++d)
9519:       if (PetscAbsReal(point[d]) > 1. + PETSC_SMALL) {
9520:         *inside = PETSC_FALSE;
9521:         break;
9522:       }
9523:     break;
9524:   default:
9525:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported polytope type %s", DMPolytopeTypes[ct]);
9526:   }
9527:   PetscFunctionReturn(PETSC_SUCCESS);
9528: }

9530: /*@
9531:   DMReorderSectionSetDefault - Set flag indicating whether the local section should be reordered by default

9533:   Logically collective

9535:   Input Parameters:
9536: + dm      - The DM
9537: - reorder - Flag for reordering

9539:   Level: intermediate

9541: .seealso: `DMReorderSectionGetDefault()`
9542: @*/
9543: PetscErrorCode DMReorderSectionSetDefault(DM dm, DMReorderDefaultFlag reorder)
9544: {
9545:   PetscFunctionBegin;
9547:   PetscTryMethod(dm, "DMReorderSectionSetDefault_C", (DM, DMReorderDefaultFlag), (dm, reorder));
9548:   PetscFunctionReturn(PETSC_SUCCESS);
9549: }

9551: /*@
9552:   DMReorderSectionGetDefault - Get flag indicating whether the local section should be reordered by default

9554:   Not collective

9556:   Input Parameter:
9557: . dm - The DM

9559:   Output Parameter:
9560: . reorder - Flag for reordering

9562:   Level: intermediate

9564: .seealso: `DMReorderSetDefault()`
9565: @*/
9566: PetscErrorCode DMReorderSectionGetDefault(DM dm, DMReorderDefaultFlag *reorder)
9567: {
9568:   PetscFunctionBegin;
9570:   PetscAssertPointer(reorder, 2);
9571:   *reorder = DM_REORDER_DEFAULT_NOTSET;
9572:   PetscTryMethod(dm, "DMReorderSectionGetDefault_C", (DM, DMReorderDefaultFlag *), (dm, reorder));
9573:   PetscFunctionReturn(PETSC_SUCCESS);
9574: }

9576: /*@
9577:   DMReorderSectionSetType - Set the type of local section reordering

9579:   Logically collective

9581:   Input Parameters:
9582: + dm      - The DM
9583: - reorder - The reordering method

9585:   Level: intermediate

9587: .seealso: `DMReorderSectionGetType()`, `DMReorderSectionSetDefault()`
9588: @*/
9589: PetscErrorCode DMReorderSectionSetType(DM dm, MatOrderingType reorder)
9590: {
9591:   PetscFunctionBegin;
9593:   PetscTryMethod(dm, "DMReorderSectionSetType_C", (DM, MatOrderingType), (dm, reorder));
9594:   PetscFunctionReturn(PETSC_SUCCESS);
9595: }

9597: /*@
9598:   DMReorderSectionGetType - Get the reordering type for the local section

9600:   Not collective

9602:   Input Parameter:
9603: . dm - The DM

9605:   Output Parameter:
9606: . reorder - The reordering method

9608:   Level: intermediate

9610: .seealso: `DMReorderSetDefault()`, `DMReorderSectionGetDefault()`
9611: @*/
9612: PetscErrorCode DMReorderSectionGetType(DM dm, MatOrderingType *reorder)
9613: {
9614:   PetscFunctionBegin;
9616:   PetscAssertPointer(reorder, 2);
9617:   *reorder = NULL;
9618:   PetscTryMethod(dm, "DMReorderSectionGetType_C", (DM, MatOrderingType *), (dm, reorder));
9619:   PetscFunctionReturn(PETSC_SUCCESS);
9620: }