Actual source code: dm.c

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

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

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

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

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

 31:   Collective

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

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

 39:   Level: beginner

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

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

 47:   `DM` is an orphan initialism or orphan acronym, the letters have no meaning and never did.

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

 56:   PetscFunctionBegin;
 57:   PetscAssertPointer(dm, 2);

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

102:   *dm = v;
103:   PetscFunctionReturn(PETSC_SUCCESS);
104: }

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

109:   Collective

111:   Input Parameter:
112: . dm - The original `DM` object

114:   Output Parameter:
115: . newdm - The new `DM` object

117:   Level: beginner

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

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

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

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

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

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

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

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

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

221:   Logically Collective

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

227:   Options Database Key:
228: . -dm_vec_type ctype - the type of vector to create

230:   Level: intermediate

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

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

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

251:   Logically Collective

253:   Input Parameter:
254: . da - initial distributed array

256:   Output Parameter:
257: . ctype - the vector type

259:   Level: intermediate

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

271: /*@
272:   VecGetDM - Gets the `DM` defining the data layout of the vector

274:   Not Collective

276:   Input Parameter:
277: . v - The `Vec`

279:   Output Parameter:
280: . dm - The `DM`

282:   Level: intermediate

284:   Note:
285:   A `Vec` may not have a `DM` associated with it.

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

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

301:   Not Collective

303:   Input Parameters:
304: + v  - The `Vec`
305: - dm - The `DM`

307:   Level: developer

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

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

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

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

328:   Logically Collective

330:   Input Parameters:
331: + dm    - the `DM` context
332: - ctype - the matrix type

334:   Options Database Key:
335: . -dm_is_coloring_type - global or local

337:   Level: intermediate

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

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

353:   Logically Collective

355:   Input Parameter:
356: . dm - the `DM` context

358:   Output Parameter:
359: . ctype - the matrix type

361:   Options Database Key:
362: . -dm_is_coloring_type - global or local

364:   Level: intermediate

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

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

380:   Logically Collective

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

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

389:   Level: intermediate

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

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

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

409:   Logically Collective

411:   Input Parameter:
412: . dm - the `DM` context

414:   Output Parameter:
415: . ctype - the matrix type

417:   Level: intermediate

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

429: /*@
430:   MatGetDM - Gets the `DM` defining the data layout of the matrix

432:   Not Collective

434:   Input Parameter:
435: . A - The `Mat`

437:   Output Parameter:
438: . dm - The `DM`

440:   Level: intermediate

442:   Note:
443:   A matrix may not have a `DM` associated with it

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

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

459: /*@
460:   MatSetDM - Sets the `DM` defining the data layout of the matrix

462:   Not Collective

464:   Input Parameters:
465: + A  - The `Mat`
466: - dm - The `DM`

468:   Level: developer

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

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

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

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

491:   Logically Collective

493:   Input Parameters:
494: + dm     - the `DM` context
495: - prefix - the prefix to prepend

497:   Level: advanced

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

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

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

519:   Logically Collective

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

525:   Level: advanced

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

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

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

546:   Not Collective

548:   Input Parameter:
549: . dm - the `DM` context

551:   Output Parameter:
552: . prefix - pointer to the prefix string used is returned

554:   Level: advanced

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

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

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

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

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

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

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

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

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

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

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

636: /*@
637:   DMDestroy - Destroys a `DM`.

639:   Collective

641:   Input Parameter:
642: . dm - the `DM` object to destroy

644:   Level: developer

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

652:   PetscFunctionBegin;
653:   if (!*dm) PetscFunctionReturn(PETSC_SUCCESS);

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

666:   PetscCall(DMClearGlobalVectors(*dm));
667:   PetscCall(DMClearLocalVectors(*dm));
668:   PetscCall(DMClearNamedGlobalVectors(*dm));
669:   PetscCall(DMClearNamedLocalVectors(*dm));

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

733:       next = b->next;
734:       PetscCall(PetscFree(b));
735:     }
736:   }

738:   PetscCall(PetscObjectDestroy(&(*dm)->dmksp));
739:   PetscCall(PetscObjectDestroy(&(*dm)->dmsnes));
740:   PetscCall(PetscObjectDestroy(&(*dm)->dmts));

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

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

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

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

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

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

799:   Collective

801:   Input Parameter:
802: . dm - the `DM` object to setup

804:   Level: intermediate

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

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

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

824:   Collective

826:   Input Parameter:
827: . dm - the `DM` object to set options for

829:   Options Database Keys:
830: + -dm_preallocate_only                               - Only preallocate the matrix for `DMCreateMatrix()` and `DMCreateMassMatrix()`, but do not fill it with zeros
831: . -dm_vec_type <type>                                - type of vector to create inside `DM`
832: . -dm_mat_type <type>                                - type of matrix to create inside `DM`
833: . -dm_is_coloring_type                               - <global or local>
834: . -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`
835: . -dm_plex_option_phases <ph0_, ph1_, ...>           - List of prefixes for option processing phases
836: . -dm_plex_filename <str>                            - File containing a mesh
837: . -dm_plex_boundary_filename <str>                   - File containing a mesh boundary
838: . -dm_plex_name <str>                                - Name of the mesh in the file
839: . -dm_plex_shape <shape>                             - The domain shape, such as `BOX`, `SPHERE`, etc.
840: . -dm_plex_cell <ct>                                 - Cell shape
841: . -dm_plex_reference_cell_domain <bool>              - Use a reference cell domain
842: . -dm_plex_dim <dim>                                 - Set the topological dimension
843: . -dm_plex_simplex <bool>                            - `PETSC_TRUE` for simplex elements, `PETSC_FALSE` for tensor elements
844: . -dm_plex_interpolate <bool>                        - `PETSC_TRUE` turns on topological interpolation (creating edges and faces)
845: . -dm_plex_orient <bool>                             - `PETSC_TRUE` turns on topological orientation (flipping 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:   Notes:

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

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

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

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

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

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

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

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

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

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

1017: /*@
1018:   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,
1019:   that is it has no ghost locations.

1021:   Collective

1023:   Input Parameter:
1024: . dm - the `DM` object

1026:   Output Parameter:
1027: . vec - the global vector

1029:   Level: beginner

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

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

1049: /*@
1050:   DMCreateLocalVector - Creates a local vector from a `DM` object.

1052:   Not Collective

1054:   Input Parameter:
1055: . dm - the `DM` object

1057:   Output Parameter:
1058: . vec - the local vector

1060:   Level: beginner

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

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

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

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

1086:   Collective

1088:   Input Parameter:
1089: . dm - the `DM` that provides the mapping

1091:   Output Parameter:
1092: . ltog - the mapping

1094:   Level: advanced

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

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

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

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

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

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

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

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

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

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

1176:   Not Collective

1178:   Input Parameter:
1179: . dm - the `DM` with block structure

1181:   Output Parameter:
1182: . bs - the block size, 1 implies no exploitable block structure

1184:   Level: intermediate

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

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

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

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

1208:   Collective

1210:   Input Parameters:
1211: + dmc - the `DM` object
1212: - dmf - the second, finer `DM` object

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

1218:   Level: developer

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

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

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

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

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

1250:   Output Parameter:
1251: . scale - the scaled vector

1253:   Level: advanced

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

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

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

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

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

1298:   Collective

1300:   Input Parameters:
1301: + dmc - the `DM` object
1302: - dmf - the second, finer `DM` object

1304:   Output Parameter:
1305: . mat - the restriction

1307:   Level: developer

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

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

1327: /*@
1328:   DMCreateInjection - Gets injection matrix between two `DM` objects.

1330:   Collective

1332:   Input Parameters:
1333: + dac - the `DM` object
1334: - daf - the second, finer `DM` object

1336:   Output Parameter:
1337: . mat - the injection

1339:   Level: developer

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

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

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

1367: /*@
1368:   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
1369:   a Galerkin finite element model on the `DM`

1371:   Collective

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

1377:   Output Parameter:
1378: . mat - the mass matrix

1380:   Level: developer

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

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

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

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

1405:   Collective

1407:   Input Parameter:
1408: . dm - the `DM` object

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

1414:   Level: developer

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

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

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

1435:   Collective

1437:   Input Parameters:
1438: + dm    - the `DM` object
1439: - ctype - `IS_COLORING_LOCAL` or `IS_COLORING_GLOBAL`

1441:   Output Parameter:
1442: . coloring - the coloring

1444:   Level: developer

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

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

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

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

1468:   Collective

1470:   Input Parameter:
1471: . dm - the `DM` object

1473:   Output Parameter:
1474: . mat - the empty Jacobian

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

1479:   Level: beginner

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

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

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

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

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

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

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

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

1541:   Logically Collective

1543:   Input Parameters:
1544: + dm   - the `DM`
1545: - skip - `PETSC_TRUE` to skip preallocation

1547:   Level: developer

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

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

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

1567:   Logically Collective

1569:   Input Parameters:
1570: + dm   - the `DM`
1571: - only - `PETSC_TRUE` if only want preallocation

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

1576:   Level: developer

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

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

1592:   Logically Collective

1594:   Input Parameters:
1595: + dm   - the `DM`
1596: - only - `PETSC_TRUE` if you only want matrix nonzero structure

1598:   Level: developer

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

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

1613:   Logically Collective

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

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

1622:   Level: advanced

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

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

1637:   Not Collective

1639:   Input Parameter:
1640: . dm - the `DM`

1642:   Output Parameter:
1643: . btype - block by topological point or field node

1645:   Level: advanced

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

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

1661:   Not Collective

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

1668:   Output Parameter:
1669: . mem - the work array

1671:   Level: developer

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

1676:   The array may contain nonzero values

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

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

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

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

1724:   Not Collective

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

1731:   Output Parameter:
1732: . mem - the work array

1734:   Level: developer

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

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

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

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

1767:   Logically Collective; No Fortran Support

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

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

1780:   Level: intermediate

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

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

1796:   Not Collective; No Fortran Support

1798:   Input Parameters:
1799: + dm    - The `DM`
1800: - field - The field number for the nullspace

1802:   Output Parameter:
1803: . nullsp - A callback to create the nullspace

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

1811:   Level: intermediate

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

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

1828:   Logically Collective; No Fortran Support

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

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

1841:   Level: intermediate

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

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

1858:   Not Collective; No Fortran Support

1860:   Input Parameters:
1861: + dm    - The `DM`
1862: - field - The field number for the nullspace

1864:   Output Parameter:
1865: . nullsp - A callback to create the near-nullspace

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

1873:   Level: intermediate

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

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

1891:   Not Collective; No Fortran Support

1893:   Input Parameter:
1894: . dm - the `DM` object

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

1901:   Level: intermediate

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

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

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

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

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

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

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

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

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

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

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

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

2016:   Not Collective; No Fortran Support

2018:   Input Parameter:
2019: . dm - the `DM` object

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

2027:   Level: intermediate

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

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

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

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

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

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

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

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

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

2106:   Not collective

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

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

2117:   Level: intermediate

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

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

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

2138:   Not collective

2140:   Input Parameters:
2141: + dms - The `DM` objects
2142: - n   - The number of `DM`s

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

2148:   Level: intermediate

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

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

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

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

2177:   Not Collective

2179:   Input Parameter:
2180: . dm - the `DM` object

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

2189:   Level: intermediate

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

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

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

2202:   Developer Notes:
2203:   The `dmlist` is for the inner subdomains or the outer subdomains or all subdomains?

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

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

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

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

2263:   Not Collective

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

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

2275:   Level: developer

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

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

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

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

2301:   Collective

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

2307:   Output Parameter:
2308: . dmf - the refined `DM`, or `NULL`

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

2313:   Level: developer

2315:   Note:
2316:   If no refinement was done, the return value is `NULL`

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

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

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

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

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

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

2350:   Logically Collective; No Fortran Support

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

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

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

2369:   Level: advanced

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

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

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

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

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

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

2403:   Logically Collective; No Fortran Support

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

2411:   Level: advanced

2413:   Note:
2414:   This function does nothing if the hook is not in the list.

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

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

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

2438:   Collective if any hooks are

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

2445:   Level: developer

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

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

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

2464: /*@
2465:   DMInterpolateSolution - Interpolates a solution from a coarse mesh to a fine mesh.

2467:   Collective

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

2477:   Output Parameter:
2478: . fineSol - the interpolation of coarseSol to the fine mesh

2480:   Level: developer

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

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

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

2497:   PetscFunctionBegin;

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

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

2515:   Not Collective

2517:   Input Parameter:
2518: . dm - the `DM` object

2520:   Output Parameter:
2521: . level - number of refinements

2523:   Level: developer

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

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

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

2541:   Not Collective

2543:   Input Parameters:
2544: + dm    - the `DM` object
2545: - level - number of refinements

2547:   Level: advanced

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

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

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

2564: /*@
2565:   DMExtrude - Extrude a `DM` object from a surface

2567:   Collective

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

2573:   Output Parameter:
2574: . dme - the extruded `DM`, or `NULL`

2576:   Level: developer

2578:   Note:
2579:   If no extrusion was done, the return value is `NULL`

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

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

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

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

2618:   Input Parameter:
2619: . dm - The `DM`

2621:   Output Parameter:
2622: . flg - `PETSC_TRUE` if a basis transformation should be done

2624:   Level: developer

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

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

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

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

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

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

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

2704:   Logically Collective

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

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

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

2726:   Level: advanced

2728:   Note:
2729:   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.

2731: .seealso: [](ch_dmbase), `DM`, `DMGlobalToLocal()`, `DMRefineHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`
2732: @*/
2733: 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)
2734: {
2735:   DMGlobalToLocalHookLink link, *p;

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

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

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

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

2784: /*@
2785:   DMGlobalToLocal - update local vectors from global vector

2787:   Neighbor-wise Collective

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

2795:   Level: beginner

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

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

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

2815: /*@
2816:   DMGlobalToLocalBegin - Begins updating local vectors from global vector

2818:   Neighbor-wise Collective

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

2826:   Level: intermediate

2828:   Notes:
2829:   The operation is completed with `DMGlobalToLocalEnd()`

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

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

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

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

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

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

2867: /*@
2868:   DMGlobalToLocalEnd - Ends updating local vectors from global vector

2870:   Neighbor-wise Collective

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

2878:   Level: intermediate

2880:   Note:
2881:   See `DMGlobalToLocalBegin()` for details.

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

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

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

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

2920:   Logically Collective

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

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

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

2942:   Level: advanced

2944: .seealso: [](ch_dmbase), `DM`, `DMLocalToGlobal()`, `DMRefineHookAdd()`, `DMGlobalToLocalHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`
2945: @*/
2946: 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)
2947: {
2948:   DMLocalToGlobalHookLink link, *p;

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

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

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

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

3003:   Neighbor-wise Collective

3005:   Input Parameters:
3006: + dm   - the `DM` object
3007: . l    - the local vector
3008: . 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.
3009: - g    - the global vector

3011:   Level: beginner

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

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

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

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

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

3033: /*@
3034:   DMLocalToGlobalBegin - begins updating global vectors from local vectors

3036:   Neighbor-wise Collective

3038:   Input Parameters:
3039: + dm   - the `DM` object
3040: . l    - the local vector
3041: . 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.
3042: - g    - the global vector

3044:   Level: intermediate

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

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

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

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

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

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

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

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

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

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

3169: /*@
3170:   DMLocalToGlobalEnd - updates global vectors from local vectors

3172:   Neighbor-wise Collective

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

3180:   Level: intermediate

3182:   Note:
3183:   See `DMLocalToGlobalBegin()` for full details

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

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

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

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

3246:   Neighbor-wise Collective

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

3253:   Output Parameter:
3254: . l - the local vector with correct ghost values

3256:   Level: intermediate

3258:   Note:
3259:   Must be followed by `DMLocalToLocalEnd()`.

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

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

3277:   Neighbor-wise Collective

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

3284:   Output Parameter:
3285: . l - the local vector with correct ghost values

3287:   Level: intermediate

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

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

3304:   Collective

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

3310:   Output Parameter:
3311: . dmc - the coarsened `DM`

3313:   Level: developer

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

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

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

3347:   Logically Collective; No Fortran Support

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

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

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

3368:   Level: advanced

3370:   Notes:
3371:   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`.

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

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

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

3380: .seealso: [](ch_dmbase), `DM`, `DMCoarsenHookRemove()`, `DMRefineHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`
3381: @*/
3382: 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)
3383: {
3384:   DMCoarsenHookLink link, *p;

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

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

3403:   Logically Collective; No Fortran Support

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

3411:   Level: advanced

3413:   Notes:
3414:   This function does nothing if the `coarsenhook` is not in the list.

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

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

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

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

3440:   Collective if any hooks are

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

3449:   Level: developer

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

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

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

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

3470:   Logically Collective; No Fortran Support

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

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

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

3490:   Level: advanced

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

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

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

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

3503: .seealso: [](ch_dmbase), `DM`, `DMSubDomainHookRemove()`, `DMRefineHookAdd()`, `SNESFASGetInterpolation()`, `SNESFASGetInjection()`, `PetscObjectCompose()`, `PetscContainerCreate()`, `DMCreateDomainDecomposition()`
3504: @*/
3505: 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)
3506: {
3507:   DMSubDomainHookLink link, *p;

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

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

3526:   Logically Collective; No Fortran Support

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

3534:   Level: advanced

3536:   Note:
3537:   See `DMSubDomainHookAdd()` for the calling sequences of `ddhook` and `restricthook`

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

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

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

3562:   Collective if any hooks are

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

3570:   Level: developer

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

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

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

3588:   Not Collective

3590:   Input Parameter:
3591: . dm - the `DM` object

3593:   Output Parameter:
3594: . level - number of coarsenings

3596:   Level: developer

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

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

3612:   Collective

3614:   Input Parameters:
3615: + dm    - the `DM` object
3616: - level - number of coarsenings

3618:   Level: developer

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

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

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

3636:   Collective

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

3642:   Output Parameter:
3643: . dmf - the refined `DM` hierarchy

3645:   Level: developer

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

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

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

3668:   Collective

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

3674:   Output Parameter:
3675: . dmc - the coarsened `DM` hierarchy

3677:   Level: developer

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

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

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

3700:   Logically Collective if the function is collective

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

3706:   Level: intermediate

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

3719: /*@
3720:   DMSetApplicationContext - Set a user context into a `DM` object

3722:   Not Collective

3724:   Input Parameters:
3725: + dm  - the `DM` object
3726: - ctx - the user context

3728:   Level: intermediate

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

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

3745: /*@
3746:   DMGetApplicationContext - Gets a user context from a `DM` object

3748:   Not Collective

3750:   Input Parameter:
3751: . dm - the `DM` object

3753:   Output Parameter:
3754: . ctx - the user context

3756:   Level: intermediate

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

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

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

3774:   Logically Collective

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

3780:   Level: intermediate

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

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

3796:   Not Collective

3798:   Input Parameter:
3799: . dm - the `DM` object to destroy

3801:   Output Parameter:
3802: . flg - `PETSC_TRUE` if the variable bounds function exists

3804:   Level: developer

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

3817: /*@
3818:   DMComputeVariableBounds - compute variable bounds used by `SNESVI`.

3820:   Logically Collective

3822:   Input Parameter:
3823: . dm - the `DM` object

3825:   Output Parameters:
3826: + xl - lower bound
3827: - xu - upper bound

3829:   Level: advanced

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

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

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

3849:   Not Collective

3851:   Input Parameter:
3852: . dm - the DM object

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

3857:   Level: developer

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

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

3873:   Not Collective

3875:   Input Parameter:
3876: . dm - the `DM` object

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

3881:   Level: developer

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

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

3897:   Not Collective

3899:   Input Parameter:
3900: . dm - the `DM` object

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

3905:   Level: developer

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

3919: PetscFunctionList DMList              = NULL;
3920: PetscBool         DMRegisterAllCalled = PETSC_FALSE;

3922: /*@
3923:   DMSetType - Builds a `DM`, for a particular `DM` implementation.

3925:   Collective

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

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

3934:   Level: intermediate

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

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

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

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

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

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

3965:   Not Collective

3967:   Input Parameter:
3968: . dm - The `DM`

3970:   Output Parameter:
3971: . type - The `DMType` name

3973:   Level: intermediate

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

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

3990:   Collective

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

3996:   Output Parameter:
3997: . M - pointer to new `DM`

3999:   Level: intermediate

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

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

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

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

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

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

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

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

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

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

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

4096: /*--------------------------------------------------------------------------------------------------------------------*/

4098: /*@C
4099:   DMRegister -  Adds a new `DM` type implementation

4101:   Not Collective, No Fortran Support

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

4107:   Level: advanced

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

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

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

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

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

4140:   Collective

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

4148:   Level: intermediate

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

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

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

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

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

4186: /* FEM Support */

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

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

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

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

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

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

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

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

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

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

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

4265:   Input Parameter:
4266: . dm - The `DM`

4268:   Output Parameter:
4269: . section - The `PetscSection`

4271:   Options Database Key:
4272: . -dm_petscsection_view - View the section created by the `DM`

4274:   Level: intermediate

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

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

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

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

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

4317:   Input Parameters:
4318: + dm      - The `DM`
4319: - section - The `PetscSection`

4321:   Level: intermediate

4323:   Note:
4324:   Any existing Section will be destroyed

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

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

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

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

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

4368:   Input Parameter:
4369: . dm - The `DM`

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

4375:   Level: developer

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

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

4391:   not Collective

4393:   Input Parameter:
4394: . dm - The `DM`

4396:   Output Parameters:
4397: + 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.
4398: . 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.
4399: - bias    - Vector containing bias to be added to constrained dofs

4401:   Level: advanced

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

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

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

4422:   Collective

4424:   Input Parameters:
4425: + dm      - The `DM`
4426: . 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).
4427: . 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).
4428: - 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).

4430:   Level: advanced

4432:   Notes:
4433:   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()`.

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

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

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

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

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

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

4483:   Level: intermediate

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

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

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

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

4552: PetscErrorCode DMGetIsoperiodicPointSF_Internal(DM dm, PetscSF *sf)
4553: {
4554:   PetscErrorCode (*f)(DM, PetscSF *);

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

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

4568:   Collective

4570:   Input Parameter:
4571: . dm - The `DM`

4573:   Output Parameter:
4574: . section - The `PetscSection`

4576:   Level: intermediate

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

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

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

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

4608:   Input Parameters:
4609: + dm      - The `DM`
4610: - section - The PetscSection, or `NULL`

4612:   Level: intermediate

4614:   Note:
4615:   Any existing `PetscSection` will be destroyed

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

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

4642:   Input Parameter:
4643: . dm - The `DM`

4645:   Output Parameter:
4646: . sf - The `PetscSF`

4648:   Level: intermediate

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

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

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

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

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

4683:   Input Parameters:
4684: + dm - The `DM`
4685: - sf - The `PetscSF`

4687:   Level: intermediate

4689:   Note:
4690:   Any previous `PetscSF` is destroyed

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

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

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

4714:   Level: developer

4716:   Note:
4717:   One usually uses `DMGetSectionSF()` to obtain the `PetscSF`

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

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

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

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

4740:   Input Parameter:
4741: . dm - The `DM`

4743:   Output Parameter:
4744: . sf - The `PetscSF`

4746:   Level: intermediate

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

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

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

4765:   Collective

4767:   Input Parameters:
4768: + dm - The `DM`
4769: - sf - The `PetscSF`

4771:   Level: intermediate

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

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

4789:   Input Parameter:
4790: . dm - The `DM`

4792:   Output Parameter:
4793: . sf - The `PetscSF`

4795:   Level: intermediate

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

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

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

4814:   Input Parameters:
4815: + dm - The DM
4816: - sf - The PetscSF

4818:   Level: intermediate

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

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

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

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

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

4869: /*@
4870:   DMClearFields - Remove all fields from the `DM`

4872:   Logically Collective

4874:   Input Parameter:
4875: . dm - The `DM`

4877:   Level: intermediate

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

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

4897: /*@
4898:   DMGetNumFields - Get the number of fields in the `DM`

4900:   Not Collective

4902:   Input Parameter:
4903: . dm - The `DM`

4905:   Output Parameter:
4906: . numFields - The number of fields

4908:   Level: intermediate

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

4921: /*@
4922:   DMSetNumFields - Set the number of fields in the `DM`

4924:   Logically Collective

4926:   Input Parameters:
4927: + dm        - The `DM`
4928: - numFields - The number of fields

4930:   Level: intermediate

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

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

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

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

4954:   Not Collective

4956:   Input Parameters:
4957: + dm - The `DM`
4958: - f  - The field number

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

4964:   Level: intermediate

4966: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMSetField()`
4967: @*/
4968: PetscErrorCode DMGetField(DM dm, PetscInt f, DMLabel *label, PetscObject *disc)
4969: {
4970:   PetscFunctionBegin;
4972:   PetscAssertPointer(disc, 4);
4973:   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);
4974:   if (label) *label = dm->fields[f].label;
4975:   if (disc) *disc = dm->fields[f].disc;
4976:   PetscFunctionReturn(PETSC_SUCCESS);
4977: }

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

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

4997:   Logically Collective

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

5005:   Level: intermediate

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

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

5026:   Logically Collective

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

5033:   Level: intermediate

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

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

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

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

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

5065:   Logically Collective

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

5072:   Level: intermediate

5074: .seealso: [](ch_dmbase), `DM`, `DMGetFieldAvoidTensor()`, `DMSetField()`, `DMGetField()`
5075: @*/
5076: PetscErrorCode DMSetFieldAvoidTensor(DM dm, PetscInt f, PetscBool avoidTensor)
5077: {
5078:   PetscFunctionBegin;
5079:   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);
5080:   dm->fields[f].avoidTensor = avoidTensor;
5081:   PetscFunctionReturn(PETSC_SUCCESS);
5082: }

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

5087:   Not Collective

5089:   Input Parameters:
5090: + dm - The `DM`
5091: - f  - The field index

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

5096:   Level: intermediate

5098: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMSetField()`, `DMGetField()`, `DMSetFieldAvoidTensor()`
5099: @*/
5100: PetscErrorCode DMGetFieldAvoidTensor(DM dm, PetscInt f, PetscBool *avoidTensor)
5101: {
5102:   PetscFunctionBegin;
5103:   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);
5104:   *avoidTensor = dm->fields[f].avoidTensor;
5105:   PetscFunctionReturn(PETSC_SUCCESS);
5106: }

5108: /*@
5109:   DMCopyFields - Copy the discretizations for the `DM` into another `DM`

5111:   Collective

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

5118:   Output Parameter:
5119: . newdm - The `DM`

5121:   Level: advanced

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

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

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

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

5156: /*@
5157:   DMGetAdjacency - Returns the flags for determining variable influence

5159:   Not Collective

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

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

5169:   Level: developer

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

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

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

5201: /*@
5202:   DMSetAdjacency - Set the flags for determining variable influence

5204:   Not Collective

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

5212:   Level: developer

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

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

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

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

5245:   Not collective

5247:   Input Parameter:
5248: . dm - The `DM` object

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

5254:   Level: developer

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

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

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

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

5285:   Not Collective

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

5292:   Level: developer

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

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

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

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

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

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

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

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

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

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

5426: /*@
5427:   DMGetNumDS - Get the number of discrete systems in the `DM`

5429:   Not Collective

5431:   Input Parameter:
5432: . dm - The `DM`

5434:   Output Parameter:
5435: . Nds - The number of `PetscDS` objects

5437:   Level: intermediate

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

5450: /*@
5451:   DMClearDS - Remove all discrete systems from the `DM`

5453:   Logically Collective

5455:   Input Parameter:
5456: . dm - The `DM`

5458:   Level: intermediate

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

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

5480: /*@
5481:   DMGetDS - Get the default `PetscDS`

5483:   Not Collective

5485:   Input Parameter:
5486: . dm - The `DM`

5488:   Output Parameter:
5489: . ds - The default `PetscDS`

5491:   Level: intermediate

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

5505: /*@
5506:   DMGetCellDS - Get the `PetscDS` defined on a given cell

5508:   Not Collective

5510:   Input Parameters:
5511: + dm    - The `DM`
5512: - point - Cell for the `PetscDS`

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

5518:   Level: developer

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

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

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

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

5555:   Not Collective

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

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

5566:   Level: advanced

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

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

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

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

5608:   Collective

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

5617:   Level: advanced

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

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

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

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

5664:   Not Collective

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

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

5676:   Level: advanced

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

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

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

5710:   Not Collective

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

5720:   Level: advanced

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

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

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

5760:   Not Collective

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

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

5769:   Level: advanced

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

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

5789: /*@
5790:   DMCreateFEDefault - Create a `PetscFE` based on the celltype for the mesh

5792:   Not Collective

5794:   Input Parameters:
5795: + dm     - The `DM`
5796: . Nc     - The number of components for the field
5797: . prefix - The options prefix for the output `PetscFE`, or `NULL`
5798: - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree

5800:   Output Parameter:
5801: . fem - The `PetscFE`

5803:   Level: intermediate

5805:   Note:
5806:   This is a convenience method that just calls `PetscFECreateByCell()` underneath.

5808: .seealso: [](ch_dmbase), `DM`, `PetscFECreateByCell()`, `DMAddField()`, `DMCreateDS()`, `DMGetCellDS()`, `DMGetRegionDS()`
5809: @*/
5810: PetscErrorCode DMCreateFEDefault(DM dm, PetscInt Nc, const char prefix[], PetscInt qorder, PetscFE *fem)
5811: {
5812:   DMPolytopeType ct;
5813:   PetscInt       dim, cStart;

5815:   PetscFunctionBegin;
5818:   if (prefix) PetscAssertPointer(prefix, 3);
5820:   PetscAssertPointer(fem, 5);
5821:   PetscCall(DMGetDimension(dm, &dim));
5822:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
5823:   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
5824:   PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, Nc, ct, prefix, qorder, fem));
5825:   PetscFunctionReturn(PETSC_SUCCESS);
5826: }

5828: /*@
5829:   DMCreateDS - Create the discrete systems for the `DM` based upon the fields added to the `DM`

5831:   Collective

5833:   Input Parameter:
5834: . dm - The `DM`

5836:   Options Database Key:
5837: . -dm_petscds_view - View all the `PetscDS` objects in this `DM`

5839:   Level: intermediate

5841: .seealso: [](ch_dmbase), `DM`, `DMSetField`, `DMAddField()`, `DMGetDS()`, `DMGetCellDS()`, `DMGetRegionDS()`, `DMSetRegionDS()`
5842: @*/
5843: PetscErrorCode DMCreateDS(DM dm)
5844: {
5845:   MPI_Comm  comm;
5846:   PetscDS   dsDef;
5847:   DMLabel  *labelSet;
5848:   PetscInt  dE, Nf = dm->Nf, f, s, Nl, l, Ndef, k;
5849:   PetscBool doSetup = PETSC_TRUE, flg;

5851:   PetscFunctionBegin;
5853:   if (!dm->fields) PetscFunctionReturn(PETSC_SUCCESS);
5854:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
5855:   PetscCall(DMGetCoordinateDim(dm, &dE));
5856:   /* Determine how many regions we have */
5857:   PetscCall(PetscMalloc1(Nf, &labelSet));
5858:   Nl   = 0;
5859:   Ndef = 0;
5860:   for (f = 0; f < Nf; ++f) {
5861:     DMLabel  label = dm->fields[f].label;
5862:     PetscInt l;

5864: #ifdef PETSC_HAVE_LIBCEED
5865:     /* Move CEED context to discretizations */
5866:     {
5867:       PetscClassId id;

5869:       PetscCall(PetscObjectGetClassId(dm->fields[f].disc, &id));
5870:       if (id == PETSCFE_CLASSID) {
5871:         Ceed ceed;

5873:         PetscCall(DMGetCeed(dm, &ceed));
5874:         PetscCall(PetscFESetCeed((PetscFE)dm->fields[f].disc, ceed));
5875:       }
5876:     }
5877: #endif
5878:     if (!label) {
5879:       ++Ndef;
5880:       continue;
5881:     }
5882:     for (l = 0; l < Nl; ++l)
5883:       if (label == labelSet[l]) break;
5884:     if (l < Nl) continue;
5885:     labelSet[Nl++] = label;
5886:   }
5887:   /* Create default DS if there are no labels to intersect with */
5888:   PetscCall(DMGetRegionDS(dm, NULL, NULL, &dsDef, NULL));
5889:   if (!dsDef && Ndef && !Nl) {
5890:     IS        fields;
5891:     PetscInt *fld, nf;

5893:     for (f = 0, nf = 0; f < Nf; ++f)
5894:       if (!dm->fields[f].label) ++nf;
5895:     PetscCheck(nf, comm, PETSC_ERR_PLIB, "All fields have labels, but we are trying to create a default DS");
5896:     PetscCall(PetscMalloc1(nf, &fld));
5897:     for (f = 0, nf = 0; f < Nf; ++f)
5898:       if (!dm->fields[f].label) fld[nf++] = f;
5899:     PetscCall(ISCreate(PETSC_COMM_SELF, &fields));
5900:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fields, "dm_fields_"));
5901:     PetscCall(ISSetType(fields, ISGENERAL));
5902:     PetscCall(ISGeneralSetIndices(fields, nf, fld, PETSC_OWN_POINTER));

5904:     PetscCall(PetscDSCreate(PETSC_COMM_SELF, &dsDef));
5905:     PetscCall(DMSetRegionDS(dm, NULL, fields, dsDef, NULL));
5906:     PetscCall(PetscDSDestroy(&dsDef));
5907:     PetscCall(ISDestroy(&fields));
5908:   }
5909:   PetscCall(DMGetRegionDS(dm, NULL, NULL, &dsDef, NULL));
5910:   if (dsDef) PetscCall(PetscDSSetCoordinateDimension(dsDef, dE));
5911:   /* Intersect labels with default fields */
5912:   if (Ndef && Nl) {
5913:     DM              plex;
5914:     DMLabel         cellLabel;
5915:     IS              fieldIS, allcellIS, defcellIS = NULL;
5916:     PetscInt       *fields;
5917:     const PetscInt *cells;
5918:     PetscInt        depth, nf = 0, n, c;

5920:     PetscCall(DMConvert(dm, DMPLEX, &plex));
5921:     PetscCall(DMPlexGetDepth(plex, &depth));
5922:     PetscCall(DMGetStratumIS(plex, "dim", depth, &allcellIS));
5923:     if (!allcellIS) PetscCall(DMGetStratumIS(plex, "depth", depth, &allcellIS));
5924:     /* TODO This looks like it only works for one label */
5925:     for (l = 0; l < Nl; ++l) {
5926:       DMLabel label = labelSet[l];
5927:       IS      pointIS;

5929:       PetscCall(ISDestroy(&defcellIS));
5930:       PetscCall(DMLabelGetStratumIS(label, 1, &pointIS));
5931:       PetscCall(ISDifference(allcellIS, pointIS, &defcellIS));
5932:       PetscCall(ISDestroy(&pointIS));
5933:     }
5934:     PetscCall(ISDestroy(&allcellIS));

5936:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "defaultCells", &cellLabel));
5937:     PetscCall(ISGetLocalSize(defcellIS, &n));
5938:     PetscCall(ISGetIndices(defcellIS, &cells));
5939:     for (c = 0; c < n; ++c) PetscCall(DMLabelSetValue(cellLabel, cells[c], 1));
5940:     PetscCall(ISRestoreIndices(defcellIS, &cells));
5941:     PetscCall(ISDestroy(&defcellIS));
5942:     PetscCall(DMPlexLabelComplete(plex, cellLabel));

5944:     PetscCall(PetscMalloc1(Ndef, &fields));
5945:     for (f = 0; f < Nf; ++f)
5946:       if (!dm->fields[f].label) fields[nf++] = f;
5947:     PetscCall(ISCreate(PETSC_COMM_SELF, &fieldIS));
5948:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fieldIS, "dm_fields_"));
5949:     PetscCall(ISSetType(fieldIS, ISGENERAL));
5950:     PetscCall(ISGeneralSetIndices(fieldIS, nf, fields, PETSC_OWN_POINTER));

5952:     PetscCall(PetscDSCreate(PETSC_COMM_SELF, &dsDef));
5953:     PetscCall(DMSetRegionDS(dm, cellLabel, fieldIS, dsDef, NULL));
5954:     PetscCall(PetscDSSetCoordinateDimension(dsDef, dE));
5955:     PetscCall(DMLabelDestroy(&cellLabel));
5956:     PetscCall(PetscDSDestroy(&dsDef));
5957:     PetscCall(ISDestroy(&fieldIS));
5958:     PetscCall(DMDestroy(&plex));
5959:   }
5960:   /* Create label DSes
5961:      - WE ONLY SUPPORT IDENTICAL OR DISJOINT LABELS
5962:   */
5963:   /* TODO Should check that labels are disjoint */
5964:   for (l = 0; l < Nl; ++l) {
5965:     DMLabel   label = labelSet[l];
5966:     PetscDS   ds, dsIn = NULL;
5967:     IS        fields;
5968:     PetscInt *fld, nf;

5970:     PetscCall(PetscDSCreate(PETSC_COMM_SELF, &ds));
5971:     for (f = 0, nf = 0; f < Nf; ++f)
5972:       if (label == dm->fields[f].label || !dm->fields[f].label) ++nf;
5973:     PetscCall(PetscMalloc1(nf, &fld));
5974:     for (f = 0, nf = 0; f < Nf; ++f)
5975:       if (label == dm->fields[f].label || !dm->fields[f].label) fld[nf++] = f;
5976:     PetscCall(ISCreate(PETSC_COMM_SELF, &fields));
5977:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fields, "dm_fields_"));
5978:     PetscCall(ISSetType(fields, ISGENERAL));
5979:     PetscCall(ISGeneralSetIndices(fields, nf, fld, PETSC_OWN_POINTER));
5980:     PetscCall(PetscDSSetCoordinateDimension(ds, dE));
5981:     {
5982:       DMPolytopeType ct;
5983:       PetscInt       lStart, lEnd;
5984:       PetscBool      isCohesiveLocal = PETSC_FALSE, isCohesive;

5986:       PetscCall(DMLabelGetBounds(label, &lStart, &lEnd));
5987:       if (lStart >= 0) {
5988:         PetscCall(DMPlexGetCellType(dm, lStart, &ct));
5989:         switch (ct) {
5990:         case DM_POLYTOPE_POINT_PRISM_TENSOR:
5991:         case DM_POLYTOPE_SEG_PRISM_TENSOR:
5992:         case DM_POLYTOPE_TRI_PRISM_TENSOR:
5993:         case DM_POLYTOPE_QUAD_PRISM_TENSOR:
5994:           isCohesiveLocal = PETSC_TRUE;
5995:           break;
5996:         default:
5997:           break;
5998:         }
5999:       }
6000:       PetscCallMPI(MPIU_Allreduce(&isCohesiveLocal, &isCohesive, 1, MPIU_BOOL, MPI_LOR, comm));
6001:       if (isCohesive) {
6002:         PetscCall(PetscDSCreate(PETSC_COMM_SELF, &dsIn));
6003:         PetscCall(PetscDSSetCoordinateDimension(dsIn, dE));
6004:       }
6005:       for (f = 0, nf = 0; f < Nf; ++f) {
6006:         if (label == dm->fields[f].label || !dm->fields[f].label) {
6007:           if (label == dm->fields[f].label) {
6008:             PetscCall(PetscDSSetDiscretization(ds, nf, NULL));
6009:             PetscCall(PetscDSSetCohesive(ds, nf, isCohesive));
6010:             if (dsIn) {
6011:               PetscCall(PetscDSSetDiscretization(dsIn, nf, NULL));
6012:               PetscCall(PetscDSSetCohesive(dsIn, nf, isCohesive));
6013:             }
6014:           }
6015:           ++nf;
6016:         }
6017:       }
6018:     }
6019:     PetscCall(DMSetRegionDS(dm, label, fields, ds, dsIn));
6020:     PetscCall(ISDestroy(&fields));
6021:     PetscCall(PetscDSDestroy(&ds));
6022:     PetscCall(PetscDSDestroy(&dsIn));
6023:   }
6024:   PetscCall(PetscFree(labelSet));
6025:   /* Set fields in DSes */
6026:   for (s = 0; s < dm->Nds; ++s) {
6027:     PetscDS         ds     = dm->probs[s].ds;
6028:     PetscDS         dsIn   = dm->probs[s].dsIn;
6029:     IS              fields = dm->probs[s].fields;
6030:     const PetscInt *fld;
6031:     PetscInt        nf, dsnf;
6032:     PetscBool       isCohesive;

6034:     PetscCall(PetscDSGetNumFields(ds, &dsnf));
6035:     PetscCall(PetscDSIsCohesive(ds, &isCohesive));
6036:     PetscCall(ISGetLocalSize(fields, &nf));
6037:     PetscCall(ISGetIndices(fields, &fld));
6038:     for (f = 0; f < nf; ++f) {
6039:       PetscObject  disc = dm->fields[fld[f]].disc;
6040:       PetscBool    isCohesiveField;
6041:       PetscClassId id;

6043:       /* Handle DS with no fields */
6044:       if (dsnf) PetscCall(PetscDSGetCohesive(ds, f, &isCohesiveField));
6045:       /* If this is a cohesive cell, then regular fields need the lower dimensional discretization */
6046:       if (isCohesive) {
6047:         if (!isCohesiveField) {
6048:           PetscObject bdDisc;

6050:           PetscCall(PetscFEGetHeightSubspace((PetscFE)disc, 1, (PetscFE *)&bdDisc));
6051:           PetscCall(PetscDSSetDiscretization(ds, f, bdDisc));
6052:           PetscCall(PetscDSSetDiscretization(dsIn, f, disc));
6053:         } else {
6054:           PetscCall(PetscDSSetDiscretization(ds, f, disc));
6055:           PetscCall(PetscDSSetDiscretization(dsIn, f, disc));
6056:         }
6057:       } else {
6058:         PetscCall(PetscDSSetDiscretization(ds, f, disc));
6059:       }
6060:       /* We allow people to have placeholder fields and construct the Section by hand */
6061:       PetscCall(PetscObjectGetClassId(disc, &id));
6062:       if ((id != PETSCFE_CLASSID) && (id != PETSCFV_CLASSID)) doSetup = PETSC_FALSE;
6063:     }
6064:     PetscCall(ISRestoreIndices(fields, &fld));
6065:   }
6066:   /* Allow k-jet tabulation */
6067:   PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)dm)->prefix, "-dm_ds_jet_degree", &k, &flg));
6068:   if (flg) {
6069:     for (s = 0; s < dm->Nds; ++s) {
6070:       PetscDS  ds   = dm->probs[s].ds;
6071:       PetscDS  dsIn = dm->probs[s].dsIn;
6072:       PetscInt Nf, f;

6074:       PetscCall(PetscDSGetNumFields(ds, &Nf));
6075:       for (f = 0; f < Nf; ++f) {
6076:         PetscCall(PetscDSSetJetDegree(ds, f, k));
6077:         if (dsIn) PetscCall(PetscDSSetJetDegree(dsIn, f, k));
6078:       }
6079:     }
6080:   }
6081:   /* Setup DSes */
6082:   if (doSetup) {
6083:     for (s = 0; s < dm->Nds; ++s) {
6084:       if (dm->setfromoptionscalled) {
6085:         PetscCall(PetscDSSetFromOptions(dm->probs[s].ds));
6086:         if (dm->probs[s].dsIn) PetscCall(PetscDSSetFromOptions(dm->probs[s].dsIn));
6087:       }
6088:       PetscCall(PetscDSSetUp(dm->probs[s].ds));
6089:       if (dm->probs[s].dsIn) PetscCall(PetscDSSetUp(dm->probs[s].dsIn));
6090:     }
6091:   }
6092:   PetscFunctionReturn(PETSC_SUCCESS);
6093: }

6095: /*@
6096:   DMUseTensorOrder - Use a tensor product closure ordering for the default section

6098:   Input Parameters:
6099: + dm     - The DM
6100: - tensor - Flag for tensor order

6102:   Level: developer

6104: .seealso: `DMPlexSetClosurePermutationTensor()`, `PetscSectionResetClosurePermutation()`
6105: @*/
6106: PetscErrorCode DMUseTensorOrder(DM dm, PetscBool tensor)
6107: {
6108:   PetscInt  Nf;
6109:   PetscBool reorder = PETSC_TRUE, isPlex;

6111:   PetscFunctionBegin;
6112:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex));
6113:   PetscCall(DMGetNumFields(dm, &Nf));
6114:   for (PetscInt f = 0; f < Nf; ++f) {
6115:     PetscObject  obj;
6116:     PetscClassId id;

6118:     PetscCall(DMGetField(dm, f, NULL, &obj));
6119:     PetscCall(PetscObjectGetClassId(obj, &id));
6120:     if (id == PETSCFE_CLASSID) {
6121:       PetscSpace sp;
6122:       PetscBool  tensor;

6124:       PetscCall(PetscFEGetBasisSpace((PetscFE)obj, &sp));
6125:       PetscCall(PetscSpacePolynomialGetTensor(sp, &tensor));
6126:       reorder = reorder && tensor ? PETSC_TRUE : PETSC_FALSE;
6127:     } else reorder = PETSC_FALSE;
6128:   }
6129:   if (tensor) {
6130:     if (reorder && isPlex) PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
6131:   } else {
6132:     PetscSection s;

6134:     PetscCall(DMGetLocalSection(dm, &s));
6135:     if (s) PetscCall(PetscSectionResetClosurePermutation(s));
6136:   }
6137:   PetscFunctionReturn(PETSC_SUCCESS);
6138: }

6140: /*@
6141:   DMComputeExactSolution - Compute the exact solution for a given `DM`, using the `PetscDS` information.

6143:   Collective

6145:   Input Parameters:
6146: + dm   - The `DM`
6147: - time - The time

6149:   Output Parameters:
6150: + u   - The vector will be filled with exact solution values, or `NULL`
6151: - u_t - The vector will be filled with the time derivative of exact solution values, or `NULL`

6153:   Level: developer

6155:   Note:
6156:   The user must call `PetscDSSetExactSolution()` before using this routine

6158: .seealso: [](ch_dmbase), `DM`, `PetscDSSetExactSolution()`
6159: @*/
6160: PetscErrorCode DMComputeExactSolution(DM dm, PetscReal time, Vec u, Vec u_t)
6161: {
6162:   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
6163:   void   **ectxs;
6164:   Vec      locu, locu_t;
6165:   PetscInt Nf, Nds, s;

6167:   PetscFunctionBegin;
6169:   if (u) {
6171:     PetscCall(DMGetLocalVector(dm, &locu));
6172:     PetscCall(VecSet(locu, 0.));
6173:   }
6174:   if (u_t) {
6176:     PetscCall(DMGetLocalVector(dm, &locu_t));
6177:     PetscCall(VecSet(locu_t, 0.));
6178:   }
6179:   PetscCall(DMGetNumFields(dm, &Nf));
6180:   PetscCall(PetscMalloc2(Nf, &exacts, Nf, &ectxs));
6181:   PetscCall(DMGetNumDS(dm, &Nds));
6182:   for (s = 0; s < Nds; ++s) {
6183:     PetscDS         ds;
6184:     DMLabel         label;
6185:     IS              fieldIS;
6186:     const PetscInt *fields, id = 1;
6187:     PetscInt        dsNf, f;

6189:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL));
6190:     PetscCall(PetscDSGetNumFields(ds, &dsNf));
6191:     PetscCall(ISGetIndices(fieldIS, &fields));
6192:     PetscCall(PetscArrayzero(exacts, Nf));
6193:     PetscCall(PetscArrayzero(ectxs, Nf));
6194:     if (u) {
6195:       for (f = 0; f < dsNf; ++f) PetscCall(PetscDSGetExactSolution(ds, fields[f], &exacts[fields[f]], &ectxs[fields[f]]));
6196:       if (label) PetscCall(DMProjectFunctionLabelLocal(dm, time, label, 1, &id, 0, NULL, exacts, ectxs, INSERT_ALL_VALUES, locu));
6197:       else PetscCall(DMProjectFunctionLocal(dm, time, exacts, ectxs, INSERT_ALL_VALUES, locu));
6198:     }
6199:     if (u_t) {
6200:       PetscCall(PetscArrayzero(exacts, Nf));
6201:       PetscCall(PetscArrayzero(ectxs, Nf));
6202:       for (f = 0; f < dsNf; ++f) PetscCall(PetscDSGetExactSolutionTimeDerivative(ds, fields[f], &exacts[fields[f]], &ectxs[fields[f]]));
6203:       if (label) PetscCall(DMProjectFunctionLabelLocal(dm, time, label, 1, &id, 0, NULL, exacts, ectxs, INSERT_ALL_VALUES, locu_t));
6204:       else PetscCall(DMProjectFunctionLocal(dm, time, exacts, ectxs, INSERT_ALL_VALUES, locu_t));
6205:     }
6206:     PetscCall(ISRestoreIndices(fieldIS, &fields));
6207:   }
6208:   if (u) {
6209:     PetscCall(PetscObjectSetName((PetscObject)u, "Exact Solution"));
6210:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)u, "exact_"));
6211:   }
6212:   if (u_t) {
6213:     PetscCall(PetscObjectSetName((PetscObject)u, "Exact Solution Time Derivative"));
6214:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)u_t, "exact_t_"));
6215:   }
6216:   PetscCall(PetscFree2(exacts, ectxs));
6217:   if (u) {
6218:     PetscCall(DMLocalToGlobalBegin(dm, locu, INSERT_ALL_VALUES, u));
6219:     PetscCall(DMLocalToGlobalEnd(dm, locu, INSERT_ALL_VALUES, u));
6220:     PetscCall(DMRestoreLocalVector(dm, &locu));
6221:   }
6222:   if (u_t) {
6223:     PetscCall(DMLocalToGlobalBegin(dm, locu_t, INSERT_ALL_VALUES, u_t));
6224:     PetscCall(DMLocalToGlobalEnd(dm, locu_t, INSERT_ALL_VALUES, u_t));
6225:     PetscCall(DMRestoreLocalVector(dm, &locu_t));
6226:   }
6227:   PetscFunctionReturn(PETSC_SUCCESS);
6228: }

6230: static PetscErrorCode DMTransferDS_Internal(DM dm, DMLabel label, IS fields, PetscInt minDegree, PetscInt maxDegree, PetscDS ds, PetscDS dsIn)
6231: {
6232:   PetscDS dsNew, dsInNew = NULL;

6234:   PetscFunctionBegin;
6235:   PetscCall(PetscDSCreate(PetscObjectComm((PetscObject)ds), &dsNew));
6236:   PetscCall(PetscDSCopy(ds, minDegree, maxDegree, dm, dsNew));
6237:   if (dsIn) {
6238:     PetscCall(PetscDSCreate(PetscObjectComm((PetscObject)dsIn), &dsInNew));
6239:     PetscCall(PetscDSCopy(dsIn, minDegree, maxDegree, dm, dsInNew));
6240:   }
6241:   PetscCall(DMSetRegionDS(dm, label, fields, dsNew, dsInNew));
6242:   PetscCall(PetscDSDestroy(&dsNew));
6243:   PetscCall(PetscDSDestroy(&dsInNew));
6244:   PetscFunctionReturn(PETSC_SUCCESS);
6245: }

6247: /*@
6248:   DMCopyDS - Copy the discrete systems for the `DM` into another `DM`

6250:   Collective

6252:   Input Parameters:
6253: + dm        - The `DM`
6254: . minDegree - Minimum degree for a discretization, or `PETSC_DETERMINE` for no limit
6255: - maxDegree - Maximum degree for a discretization, or `PETSC_DETERMINE` for no limit

6257:   Output Parameter:
6258: . newdm - The `DM`

6260:   Level: advanced

6262: .seealso: [](ch_dmbase), `DM`, `DMCopyFields()`, `DMAddField()`, `DMGetDS()`, `DMGetCellDS()`, `DMGetRegionDS()`, `DMSetRegionDS()`
6263: @*/
6264: PetscErrorCode DMCopyDS(DM dm, PetscInt minDegree, PetscInt maxDegree, DM newdm)
6265: {
6266:   PetscInt Nds, s;

6268:   PetscFunctionBegin;
6269:   if (dm == newdm) PetscFunctionReturn(PETSC_SUCCESS);
6270:   PetscCall(DMGetNumDS(dm, &Nds));
6271:   PetscCall(DMClearDS(newdm));
6272:   for (s = 0; s < Nds; ++s) {
6273:     DMLabel  label;
6274:     IS       fields;
6275:     PetscDS  ds, dsIn, newds;
6276:     PetscInt Nbd, bd;

6278:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fields, &ds, &dsIn));
6279:     /* TODO: We need to change all keys from labels in the old DM to labels in the new DM */
6280:     PetscCall(DMTransferDS_Internal(newdm, label, fields, minDegree, maxDegree, ds, dsIn));
6281:     /* Complete new labels in the new DS */
6282:     PetscCall(DMGetRegionDS(newdm, label, NULL, &newds, NULL));
6283:     PetscCall(PetscDSGetNumBoundary(newds, &Nbd));
6284:     for (bd = 0; bd < Nbd; ++bd) {
6285:       PetscWeakForm wf;
6286:       DMLabel       label;
6287:       PetscInt      field;

6289:       PetscCall(PetscDSGetBoundary(newds, bd, &wf, NULL, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL));
6290:       PetscCall(PetscWeakFormReplaceLabel(wf, label));
6291:     }
6292:   }
6293:   PetscCall(DMCompleteBCLabels_Internal(newdm));
6294:   PetscFunctionReturn(PETSC_SUCCESS);
6295: }

6297: /*@
6298:   DMCopyDisc - Copy the fields and discrete systems for the `DM` into another `DM`

6300:   Collective

6302:   Input Parameter:
6303: . dm - The `DM`

6305:   Output Parameter:
6306: . newdm - The `DM`

6308:   Level: advanced

6310:   Developer Note:
6311:   Really ugly name, nothing in PETSc is called a `Disc` plus it is an ugly abbreviation

6313: .seealso: [](ch_dmbase), `DM`, `DMCopyFields()`, `DMCopyDS()`
6314: @*/
6315: PetscErrorCode DMCopyDisc(DM dm, DM newdm)
6316: {
6317:   PetscFunctionBegin;
6318:   PetscCall(DMCopyFields(dm, PETSC_DETERMINE, PETSC_DETERMINE, newdm));
6319:   PetscCall(DMCopyDS(dm, PETSC_DETERMINE, PETSC_DETERMINE, newdm));
6320:   PetscFunctionReturn(PETSC_SUCCESS);
6321: }

6323: /*@
6324:   DMGetDimension - Return the topological dimension of the `DM`

6326:   Not Collective

6328:   Input Parameter:
6329: . dm - The `DM`

6331:   Output Parameter:
6332: . dim - The topological dimension

6334:   Level: beginner

6336: .seealso: [](ch_dmbase), `DM`, `DMSetDimension()`, `DMCreate()`
6337: @*/
6338: PetscErrorCode DMGetDimension(DM dm, PetscInt *dim)
6339: {
6340:   PetscFunctionBegin;
6342:   PetscAssertPointer(dim, 2);
6343:   *dim = dm->dim;
6344:   PetscFunctionReturn(PETSC_SUCCESS);
6345: }

6347: /*@
6348:   DMSetDimension - Set the topological dimension of the `DM`

6350:   Collective

6352:   Input Parameters:
6353: + dm  - The `DM`
6354: - dim - The topological dimension

6356:   Level: beginner

6358: .seealso: [](ch_dmbase), `DM`, `DMGetDimension()`, `DMCreate()`
6359: @*/
6360: PetscErrorCode DMSetDimension(DM dm, PetscInt dim)
6361: {
6362:   PetscDS  ds;
6363:   PetscInt Nds, n;

6365:   PetscFunctionBegin;
6368:   dm->dim = dim;
6369:   if (dm->dim >= 0) {
6370:     PetscCall(DMGetNumDS(dm, &Nds));
6371:     for (n = 0; n < Nds; ++n) {
6372:       PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds, NULL));
6373:       if (ds->dimEmbed < 0) PetscCall(PetscDSSetCoordinateDimension(ds, dim));
6374:     }
6375:   }
6376:   PetscFunctionReturn(PETSC_SUCCESS);
6377: }

6379: /*@
6380:   DMGetDimPoints - Get the half-open interval for all points of a given dimension

6382:   Collective

6384:   Input Parameters:
6385: + dm  - the `DM`
6386: - dim - the dimension

6388:   Output Parameters:
6389: + pStart - The first point of the given dimension
6390: - pEnd   - The first point following points of the given dimension

6392:   Level: intermediate

6394:   Note:
6395:   The points are vertices in the Hasse diagram encoding the topology. This is explained in
6396:   https://arxiv.org/abs/0908.4427. If no points exist of this dimension in the storage scheme,
6397:   then the interval is empty.

6399: .seealso: [](ch_dmbase), `DM`, `DMPLEX`, `DMPlexGetDepthStratum()`, `DMPlexGetHeightStratum()`
6400: @*/
6401: PetscErrorCode DMGetDimPoints(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
6402: {
6403:   PetscInt d;

6405:   PetscFunctionBegin;
6407:   PetscCall(DMGetDimension(dm, &d));
6408:   PetscCheck((dim >= 0) && (dim <= d), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim);
6409:   PetscUseTypeMethod(dm, getdimpoints, dim, pStart, pEnd);
6410:   PetscFunctionReturn(PETSC_SUCCESS);
6411: }

6413: /*@
6414:   DMGetOutputDM - Retrieve the `DM` associated with the layout for output

6416:   Collective

6418:   Input Parameter:
6419: . dm - The original `DM`

6421:   Output Parameter:
6422: . odm - The `DM` which provides the layout for output

6424:   Level: intermediate

6426:   Note:
6427:   In some situations the vector obtained with `DMCreateGlobalVector()` excludes points for degrees of freedom that are associated with fixed (Dirichelet) boundary
6428:   conditions since the algebraic solver does not solve for those variables. The output `DM` includes these excluded points and its global vector contains the
6429:   locations for those dof so that they can be output to a file or other viewer along with the unconstrained dof.

6431: .seealso: [](ch_dmbase), `DM`, `VecView()`, `DMGetGlobalSection()`, `DMCreateGlobalVector()`, `PetscSectionHasConstraints()`, `DMSetGlobalSection()`
6432: @*/
6433: PetscErrorCode DMGetOutputDM(DM dm, DM *odm)
6434: {
6435:   PetscSection section;
6436:   IS           perm;
6437:   PetscBool    hasConstraints, newDM, gnewDM;
6438:   PetscInt     num_face_sfs = 0;

6440:   PetscFunctionBegin;
6442:   PetscAssertPointer(odm, 2);
6443:   PetscCall(DMGetLocalSection(dm, &section));
6444:   PetscCall(PetscSectionHasConstraints(section, &hasConstraints));
6445:   PetscCall(PetscSectionGetPermutation(section, &perm));
6446:   PetscCall(DMPlexGetIsoperiodicFaceSF(dm, &num_face_sfs, NULL));
6447:   newDM = hasConstraints || perm || (num_face_sfs > 0) ? PETSC_TRUE : PETSC_FALSE;
6448:   PetscCallMPI(MPIU_Allreduce(&newDM, &gnewDM, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)dm)));
6449:   if (!gnewDM) {
6450:     *odm = dm;
6451:     PetscFunctionReturn(PETSC_SUCCESS);
6452:   }
6453:   if (!dm->dmBC) {
6454:     PetscSection newSection, gsection;
6455:     PetscSF      sf, sfNatural;
6456:     PetscBool    usePerm = dm->ignorePermOutput ? PETSC_FALSE : PETSC_TRUE;

6458:     PetscCall(DMClone(dm, &dm->dmBC));
6459:     PetscCall(DMCopyDisc(dm, dm->dmBC));
6460:     PetscCall(PetscSectionClone(section, &newSection));
6461:     PetscCall(DMSetLocalSection(dm->dmBC, newSection));
6462:     PetscCall(PetscSectionDestroy(&newSection));
6463:     PetscCall(DMGetNaturalSF(dm, &sfNatural));
6464:     PetscCall(DMSetNaturalSF(dm->dmBC, sfNatural));
6465:     PetscCall(DMGetPointSF(dm->dmBC, &sf));
6466:     PetscCall(PetscSectionCreateGlobalSection(section, sf, usePerm, PETSC_TRUE, PETSC_FALSE, &gsection));
6467:     PetscCall(DMSetGlobalSection(dm->dmBC, gsection));
6468:     PetscCall(PetscSectionDestroy(&gsection));
6469:   }
6470:   *odm = dm->dmBC;
6471:   PetscFunctionReturn(PETSC_SUCCESS);
6472: }

6474: /*@
6475:   DMGetOutputSequenceNumber - Retrieve the sequence number/value for output

6477:   Input Parameter:
6478: . dm - The original `DM`

6480:   Output Parameters:
6481: + num - The output sequence number
6482: - val - The output sequence value

6484:   Level: intermediate

6486:   Note:
6487:   This is intended for output that should appear in sequence, for instance
6488:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6490:   Developer Note:
6491:   The `DM` serves as a convenient place to store the current iteration value. The iteration is not
6492:   not directly related to the `DM`.

6494: .seealso: [](ch_dmbase), `DM`, `VecView()`
6495: @*/
6496: PetscErrorCode DMGetOutputSequenceNumber(DM dm, PetscInt *num, PetscReal *val)
6497: {
6498:   PetscFunctionBegin;
6500:   if (num) {
6501:     PetscAssertPointer(num, 2);
6502:     *num = dm->outputSequenceNum;
6503:   }
6504:   if (val) {
6505:     PetscAssertPointer(val, 3);
6506:     *val = dm->outputSequenceVal;
6507:   }
6508:   PetscFunctionReturn(PETSC_SUCCESS);
6509: }

6511: /*@
6512:   DMSetOutputSequenceNumber - Set the sequence number/value for output

6514:   Input Parameters:
6515: + dm  - The original `DM`
6516: . num - The output sequence number
6517: - val - The output sequence value

6519:   Level: intermediate

6521:   Note:
6522:   This is intended for output that should appear in sequence, for instance
6523:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6525: .seealso: [](ch_dmbase), `DM`, `VecView()`
6526: @*/
6527: PetscErrorCode DMSetOutputSequenceNumber(DM dm, PetscInt num, PetscReal val)
6528: {
6529:   PetscFunctionBegin;
6531:   dm->outputSequenceNum = num;
6532:   dm->outputSequenceVal = val;
6533:   PetscFunctionReturn(PETSC_SUCCESS);
6534: }

6536: /*@
6537:   DMOutputSequenceLoad - Retrieve the sequence value from a `PetscViewer`

6539:   Input Parameters:
6540: + dm     - The original `DM`
6541: . viewer - The `PetscViewer` to get it from
6542: . name   - The sequence name
6543: - num    - The output sequence number

6545:   Output Parameter:
6546: . val - The output sequence value

6548:   Level: intermediate

6550:   Note:
6551:   This is intended for output that should appear in sequence, for instance
6552:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6554:   Developer Note:
6555:   It is unclear at the user API level why a `DM` is needed as input

6557: .seealso: [](ch_dmbase), `DM`, `DMGetOutputSequenceNumber()`, `DMSetOutputSequenceNumber()`, `VecView()`
6558: @*/
6559: PetscErrorCode DMOutputSequenceLoad(DM dm, PetscViewer viewer, const char name[], PetscInt num, PetscReal *val)
6560: {
6561:   PetscBool ishdf5;

6563:   PetscFunctionBegin;
6566:   PetscAssertPointer(name, 3);
6567:   PetscAssertPointer(val, 5);
6568:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
6569:   if (ishdf5) {
6570: #if defined(PETSC_HAVE_HDF5)
6571:     PetscScalar value;

6573:     PetscCall(DMSequenceLoad_HDF5_Internal(dm, name, num, &value, viewer));
6574:     *val = PetscRealPart(value);
6575: #endif
6576:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerHDF5Open()");
6577:   PetscFunctionReturn(PETSC_SUCCESS);
6578: }

6580: /*@
6581:   DMGetOutputSequenceLength - Retrieve the number of sequence values from a `PetscViewer`

6583:   Input Parameters:
6584: + dm     - The original `DM`
6585: . viewer - The `PetscViewer` to get it from
6586: - name   - The sequence name

6588:   Output Parameter:
6589: . len - The length of the output sequence

6591:   Level: intermediate

6593:   Note:
6594:   This is intended for output that should appear in sequence, for instance
6595:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6597:   Developer Note:
6598:   It is unclear at the user API level why a `DM` is needed as input

6600: .seealso: [](ch_dmbase), `DM`, `DMGetOutputSequenceNumber()`, `DMSetOutputSequenceNumber()`, `VecView()`
6601: @*/
6602: PetscErrorCode DMGetOutputSequenceLength(DM dm, PetscViewer viewer, const char name[], PetscInt *len)
6603: {
6604:   PetscBool ishdf5;

6606:   PetscFunctionBegin;
6609:   PetscAssertPointer(name, 3);
6610:   PetscAssertPointer(len, 4);
6611:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
6612:   if (ishdf5) {
6613: #if defined(PETSC_HAVE_HDF5)
6614:     PetscCall(DMSequenceGetLength_HDF5_Internal(dm, name, len, viewer));
6615: #endif
6616:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerHDF5Open()");
6617:   PetscFunctionReturn(PETSC_SUCCESS);
6618: }

6620: /*@
6621:   DMGetUseNatural - Get the flag for creating a mapping to the natural order when a `DM` is (re)distributed in parallel

6623:   Not Collective

6625:   Input Parameter:
6626: . dm - The `DM`

6628:   Output Parameter:
6629: . useNatural - `PETSC_TRUE` to build the mapping to a natural order during distribution

6631:   Level: beginner

6633: .seealso: [](ch_dmbase), `DM`, `DMSetUseNatural()`, `DMCreate()`
6634: @*/
6635: PetscErrorCode DMGetUseNatural(DM dm, PetscBool *useNatural)
6636: {
6637:   PetscFunctionBegin;
6639:   PetscAssertPointer(useNatural, 2);
6640:   *useNatural = dm->useNatural;
6641:   PetscFunctionReturn(PETSC_SUCCESS);
6642: }

6644: /*@
6645:   DMSetUseNatural - Set the flag for creating a mapping to the natural order when a `DM` is (re)distributed in parallel

6647:   Collective

6649:   Input Parameters:
6650: + dm         - The `DM`
6651: - useNatural - `PETSC_TRUE` to build the mapping to a natural order during distribution

6653:   Level: beginner

6655:   Note:
6656:   This also causes the map to be build after `DMCreateSubDM()` and `DMCreateSuperDM()`

6658: .seealso: [](ch_dmbase), `DM`, `DMGetUseNatural()`, `DMCreate()`, `DMPlexDistribute()`, `DMCreateSubDM()`, `DMCreateSuperDM()`
6659: @*/
6660: PetscErrorCode DMSetUseNatural(DM dm, PetscBool useNatural)
6661: {
6662:   PetscFunctionBegin;
6665:   dm->useNatural = useNatural;
6666:   PetscFunctionReturn(PETSC_SUCCESS);
6667: }

6669: /*@
6670:   DMCreateLabel - Create a label of the given name if it does not already exist in the `DM`

6672:   Not Collective

6674:   Input Parameters:
6675: + dm   - The `DM` object
6676: - name - The label name

6678:   Level: intermediate

6680: .seealso: [](ch_dmbase), `DM`, `DMLabelCreate()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6681: @*/
6682: PetscErrorCode DMCreateLabel(DM dm, const char name[])
6683: {
6684:   PetscBool flg;
6685:   DMLabel   label;

6687:   PetscFunctionBegin;
6689:   PetscAssertPointer(name, 2);
6690:   PetscCall(DMHasLabel(dm, name, &flg));
6691:   if (!flg) {
6692:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, &label));
6693:     PetscCall(DMAddLabel(dm, label));
6694:     PetscCall(DMLabelDestroy(&label));
6695:   }
6696:   PetscFunctionReturn(PETSC_SUCCESS);
6697: }

6699: /*@
6700:   DMCreateLabelAtIndex - Create a label of the given name at the given index. If it already exists in the `DM`, move it to this index.

6702:   Not Collective

6704:   Input Parameters:
6705: + dm   - The `DM` object
6706: . l    - The index for the label
6707: - name - The label name

6709:   Level: intermediate

6711: .seealso: [](ch_dmbase), `DM`, `DMCreateLabel()`, `DMLabelCreate()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6712: @*/
6713: PetscErrorCode DMCreateLabelAtIndex(DM dm, PetscInt l, const char name[])
6714: {
6715:   DMLabelLink orig, prev = NULL;
6716:   DMLabel     label;
6717:   PetscInt    Nl, m;
6718:   PetscBool   flg, match;
6719:   const char *lname;

6721:   PetscFunctionBegin;
6723:   PetscAssertPointer(name, 3);
6724:   PetscCall(DMHasLabel(dm, name, &flg));
6725:   if (!flg) {
6726:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, &label));
6727:     PetscCall(DMAddLabel(dm, label));
6728:     PetscCall(DMLabelDestroy(&label));
6729:   }
6730:   PetscCall(DMGetNumLabels(dm, &Nl));
6731:   PetscCheck(l < Nl, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label index %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", l, Nl);
6732:   for (m = 0, orig = dm->labels; m < Nl; ++m, prev = orig, orig = orig->next) {
6733:     PetscCall(PetscObjectGetName((PetscObject)orig->label, &lname));
6734:     PetscCall(PetscStrcmp(name, lname, &match));
6735:     if (match) break;
6736:   }
6737:   if (m == l) PetscFunctionReturn(PETSC_SUCCESS);
6738:   if (!m) dm->labels = orig->next;
6739:   else prev->next = orig->next;
6740:   if (!l) {
6741:     orig->next = dm->labels;
6742:     dm->labels = orig;
6743:   } else {
6744:     for (m = 0, prev = dm->labels; m < l - 1; ++m, prev = prev->next);
6745:     orig->next = prev->next;
6746:     prev->next = orig;
6747:   }
6748:   PetscFunctionReturn(PETSC_SUCCESS);
6749: }

6751: /*@
6752:   DMGetLabelValue - Get the value in a `DMLabel` for the given point, with -1 as the default

6754:   Not Collective

6756:   Input Parameters:
6757: + dm    - The `DM` object
6758: . name  - The label name
6759: - point - The mesh point

6761:   Output Parameter:
6762: . value - The label value for this point, or -1 if the point is not in the label

6764:   Level: beginner

6766: .seealso: [](ch_dmbase), `DM`, `DMLabelGetValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6767: @*/
6768: PetscErrorCode DMGetLabelValue(DM dm, const char name[], PetscInt point, PetscInt *value)
6769: {
6770:   DMLabel label;

6772:   PetscFunctionBegin;
6774:   PetscAssertPointer(name, 2);
6775:   PetscCall(DMGetLabel(dm, name, &label));
6776:   PetscCheck(label, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No label named %s was found", name);
6777:   PetscCall(DMLabelGetValue(label, point, value));
6778:   PetscFunctionReturn(PETSC_SUCCESS);
6779: }

6781: /*@
6782:   DMSetLabelValue - Add a point to a `DMLabel` with given value

6784:   Not Collective

6786:   Input Parameters:
6787: + dm    - The `DM` object
6788: . name  - The label name
6789: . point - The mesh point
6790: - value - The label value for this point

6792:   Output Parameter:

6794:   Level: beginner

6796: .seealso: [](ch_dmbase), `DM`, `DMLabelSetValue()`, `DMGetStratumIS()`, `DMClearLabelValue()`
6797: @*/
6798: PetscErrorCode DMSetLabelValue(DM dm, const char name[], PetscInt point, PetscInt value)
6799: {
6800:   DMLabel label;

6802:   PetscFunctionBegin;
6804:   PetscAssertPointer(name, 2);
6805:   PetscCall(DMGetLabel(dm, name, &label));
6806:   if (!label) {
6807:     PetscCall(DMCreateLabel(dm, name));
6808:     PetscCall(DMGetLabel(dm, name, &label));
6809:   }
6810:   PetscCall(DMLabelSetValue(label, point, value));
6811:   PetscFunctionReturn(PETSC_SUCCESS);
6812: }

6814: /*@
6815:   DMClearLabelValue - Remove a point from a `DMLabel` with given value

6817:   Not Collective

6819:   Input Parameters:
6820: + dm    - The `DM` object
6821: . name  - The label name
6822: . point - The mesh point
6823: - value - The label value for this point

6825:   Level: beginner

6827: .seealso: [](ch_dmbase), `DM`, `DMLabelClearValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6828: @*/
6829: PetscErrorCode DMClearLabelValue(DM dm, const char name[], PetscInt point, PetscInt value)
6830: {
6831:   DMLabel label;

6833:   PetscFunctionBegin;
6835:   PetscAssertPointer(name, 2);
6836:   PetscCall(DMGetLabel(dm, name, &label));
6837:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6838:   PetscCall(DMLabelClearValue(label, point, value));
6839:   PetscFunctionReturn(PETSC_SUCCESS);
6840: }

6842: /*@
6843:   DMGetLabelSize - Get the value of `DMLabelGetNumValues()` of a `DMLabel` in the `DM`

6845:   Not Collective

6847:   Input Parameters:
6848: + dm   - The `DM` object
6849: - name - The label name

6851:   Output Parameter:
6852: . size - The number of different integer ids, or 0 if the label does not exist

6854:   Level: beginner

6856:   Developer Note:
6857:   This should be renamed to something like `DMGetLabelNumValues()` or removed.

6859: .seealso: [](ch_dmbase), `DM`, `DMLabelGetNumValues()`, `DMSetLabelValue()`, `DMGetLabel()`
6860: @*/
6861: PetscErrorCode DMGetLabelSize(DM dm, const char name[], PetscInt *size)
6862: {
6863:   DMLabel label;

6865:   PetscFunctionBegin;
6867:   PetscAssertPointer(name, 2);
6868:   PetscAssertPointer(size, 3);
6869:   PetscCall(DMGetLabel(dm, name, &label));
6870:   *size = 0;
6871:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6872:   PetscCall(DMLabelGetNumValues(label, size));
6873:   PetscFunctionReturn(PETSC_SUCCESS);
6874: }

6876: /*@
6877:   DMGetLabelIdIS - Get the `DMLabelGetValueIS()` from a `DMLabel` in the `DM`

6879:   Not Collective

6881:   Input Parameters:
6882: + dm   - The `DM` object
6883: - name - The label name

6885:   Output Parameter:
6886: . ids - The integer ids, or `NULL` if the label does not exist

6888:   Level: beginner

6890: .seealso: [](ch_dmbase), `DM`, `DMLabelGetValueIS()`, `DMGetLabelSize()`
6891: @*/
6892: PetscErrorCode DMGetLabelIdIS(DM dm, const char name[], IS *ids)
6893: {
6894:   DMLabel label;

6896:   PetscFunctionBegin;
6898:   PetscAssertPointer(name, 2);
6899:   PetscAssertPointer(ids, 3);
6900:   PetscCall(DMGetLabel(dm, name, &label));
6901:   *ids = NULL;
6902:   if (label) {
6903:     PetscCall(DMLabelGetValueIS(label, ids));
6904:   } else {
6905:     /* returning an empty IS */
6906:     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_USE_POINTER, ids));
6907:   }
6908:   PetscFunctionReturn(PETSC_SUCCESS);
6909: }

6911: /*@
6912:   DMGetStratumSize - Get the number of points in a label stratum

6914:   Not Collective

6916:   Input Parameters:
6917: + dm    - The `DM` object
6918: . name  - The label name of the stratum
6919: - value - The stratum value

6921:   Output Parameter:
6922: . size - The number of points, also called the stratum size

6924:   Level: beginner

6926: .seealso: [](ch_dmbase), `DM`, `DMLabelGetStratumSize()`, `DMGetLabelSize()`, `DMGetLabelIds()`
6927: @*/
6928: PetscErrorCode DMGetStratumSize(DM dm, const char name[], PetscInt value, PetscInt *size)
6929: {
6930:   DMLabel label;

6932:   PetscFunctionBegin;
6934:   PetscAssertPointer(name, 2);
6935:   PetscAssertPointer(size, 4);
6936:   PetscCall(DMGetLabel(dm, name, &label));
6937:   *size = 0;
6938:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6939:   PetscCall(DMLabelGetStratumSize(label, value, size));
6940:   PetscFunctionReturn(PETSC_SUCCESS);
6941: }

6943: /*@
6944:   DMGetStratumIS - Get the points in a label stratum

6946:   Not Collective

6948:   Input Parameters:
6949: + dm    - The `DM` object
6950: . name  - The label name
6951: - value - The stratum value

6953:   Output Parameter:
6954: . points - The stratum points, or `NULL` if the label does not exist or does not have that value

6956:   Level: beginner

6958: .seealso: [](ch_dmbase), `DM`, `DMLabelGetStratumIS()`, `DMGetStratumSize()`
6959: @*/
6960: PetscErrorCode DMGetStratumIS(DM dm, const char name[], PetscInt value, IS *points)
6961: {
6962:   DMLabel label;

6964:   PetscFunctionBegin;
6966:   PetscAssertPointer(name, 2);
6967:   PetscAssertPointer(points, 4);
6968:   PetscCall(DMGetLabel(dm, name, &label));
6969:   *points = NULL;
6970:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6971:   PetscCall(DMLabelGetStratumIS(label, value, points));
6972:   PetscFunctionReturn(PETSC_SUCCESS);
6973: }

6975: /*@
6976:   DMSetStratumIS - Set the points in a label stratum

6978:   Not Collective

6980:   Input Parameters:
6981: + dm     - The `DM` object
6982: . name   - The label name
6983: . value  - The stratum value
6984: - points - The stratum points

6986:   Level: beginner

6988: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMClearLabelStratum()`, `DMLabelClearStratum()`, `DMLabelSetStratumIS()`, `DMGetStratumSize()`
6989: @*/
6990: PetscErrorCode DMSetStratumIS(DM dm, const char name[], PetscInt value, IS points)
6991: {
6992:   DMLabel label;

6994:   PetscFunctionBegin;
6996:   PetscAssertPointer(name, 2);
6998:   PetscCall(DMGetLabel(dm, name, &label));
6999:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
7000:   PetscCall(DMLabelSetStratumIS(label, value, points));
7001:   PetscFunctionReturn(PETSC_SUCCESS);
7002: }

7004: /*@
7005:   DMClearLabelStratum - Remove all points from a stratum from a `DMLabel`

7007:   Not Collective

7009:   Input Parameters:
7010: + dm    - The `DM` object
7011: . name  - The label name
7012: - value - The label value for this point

7014:   Output Parameter:

7016:   Level: beginner

7018: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMLabelClearStratum()`, `DMSetLabelValue()`, `DMGetStratumIS()`, `DMClearLabelValue()`
7019: @*/
7020: PetscErrorCode DMClearLabelStratum(DM dm, const char name[], PetscInt value)
7021: {
7022:   DMLabel label;

7024:   PetscFunctionBegin;
7026:   PetscAssertPointer(name, 2);
7027:   PetscCall(DMGetLabel(dm, name, &label));
7028:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
7029:   PetscCall(DMLabelClearStratum(label, value));
7030:   PetscFunctionReturn(PETSC_SUCCESS);
7031: }

7033: /*@
7034:   DMGetNumLabels - Return the number of labels defined by on the `DM`

7036:   Not Collective

7038:   Input Parameter:
7039: . dm - The `DM` object

7041:   Output Parameter:
7042: . numLabels - the number of Labels

7044:   Level: intermediate

7046: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetLabelByNum()`, `DMGetLabelName()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7047: @*/
7048: PetscErrorCode DMGetNumLabels(DM dm, PetscInt *numLabels)
7049: {
7050:   DMLabelLink next = dm->labels;
7051:   PetscInt    n    = 0;

7053:   PetscFunctionBegin;
7055:   PetscAssertPointer(numLabels, 2);
7056:   while (next) {
7057:     ++n;
7058:     next = next->next;
7059:   }
7060:   *numLabels = n;
7061:   PetscFunctionReturn(PETSC_SUCCESS);
7062: }

7064: /*@
7065:   DMGetLabelName - Return the name of nth label

7067:   Not Collective

7069:   Input Parameters:
7070: + dm - The `DM` object
7071: - n  - the label number

7073:   Output Parameter:
7074: . name - the label name

7076:   Level: intermediate

7078:   Developer Note:
7079:   Some of the functions that appropriate on labels using their number have the suffix ByNum, others do not.

7081: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetLabelByNum()`, `DMGetLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7082: @*/
7083: PetscErrorCode DMGetLabelName(DM dm, PetscInt n, const char *name[])
7084: {
7085:   DMLabelLink next = dm->labels;
7086:   PetscInt    l    = 0;

7088:   PetscFunctionBegin;
7090:   PetscAssertPointer(name, 3);
7091:   while (next) {
7092:     if (l == n) {
7093:       PetscCall(PetscObjectGetName((PetscObject)next->label, name));
7094:       PetscFunctionReturn(PETSC_SUCCESS);
7095:     }
7096:     ++l;
7097:     next = next->next;
7098:   }
7099:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %" PetscInt_FMT " does not exist in this DM", n);
7100: }

7102: /*@
7103:   DMHasLabel - Determine whether the `DM` has a label of a given name

7105:   Not Collective

7107:   Input Parameters:
7108: + dm   - The `DM` object
7109: - name - The label name

7111:   Output Parameter:
7112: . hasLabel - `PETSC_TRUE` if the label is present

7114:   Level: intermediate

7116: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetLabel()`, `DMGetLabelByNum()`, `DMCreateLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7117: @*/
7118: PetscErrorCode DMHasLabel(DM dm, const char name[], PetscBool *hasLabel)
7119: {
7120:   DMLabelLink next = dm->labels;
7121:   const char *lname;

7123:   PetscFunctionBegin;
7125:   PetscAssertPointer(name, 2);
7126:   PetscAssertPointer(hasLabel, 3);
7127:   *hasLabel = PETSC_FALSE;
7128:   while (next) {
7129:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7130:     PetscCall(PetscStrcmp(name, lname, hasLabel));
7131:     if (*hasLabel) break;
7132:     next = next->next;
7133:   }
7134:   PetscFunctionReturn(PETSC_SUCCESS);
7135: }

7137: // PetscClangLinter pragma ignore: -fdoc-section-header-unknown
7138: /*@
7139:   DMGetLabel - Return the label of a given name, or `NULL`, from a `DM`

7141:   Not Collective

7143:   Input Parameters:
7144: + dm   - The `DM` object
7145: - name - The label name

7147:   Output Parameter:
7148: . label - The `DMLabel`, or `NULL` if the label is absent

7150:   Default labels in a `DMPLEX`:
7151: + "depth"       - Holds the depth (co-dimension) of each mesh point
7152: . "celltype"    - Holds the topological type of each cell
7153: . "ghost"       - If the DM is distributed with overlap, this marks the cells and faces in the overlap
7154: . "Cell Sets"   - Mirrors the cell sets defined by GMsh and ExodusII
7155: . "Face Sets"   - Mirrors the face sets defined by GMsh and ExodusII
7156: - "Vertex Sets" - Mirrors the vertex sets defined by GMsh

7158:   Level: intermediate

7160: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMHasLabel()`, `DMGetLabelByNum()`, `DMAddLabel()`, `DMCreateLabel()`, `DMPlexGetDepthLabel()`, `DMPlexGetCellType()`
7161: @*/
7162: PetscErrorCode DMGetLabel(DM dm, const char name[], DMLabel *label)
7163: {
7164:   DMLabelLink next = dm->labels;
7165:   PetscBool   hasLabel;
7166:   const char *lname;

7168:   PetscFunctionBegin;
7170:   PetscAssertPointer(name, 2);
7171:   PetscAssertPointer(label, 3);
7172:   *label = NULL;
7173:   while (next) {
7174:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7175:     PetscCall(PetscStrcmp(name, lname, &hasLabel));
7176:     if (hasLabel) {
7177:       *label = next->label;
7178:       break;
7179:     }
7180:     next = next->next;
7181:   }
7182:   PetscFunctionReturn(PETSC_SUCCESS);
7183: }

7185: /*@
7186:   DMGetLabelByNum - Return the nth label on a `DM`

7188:   Not Collective

7190:   Input Parameters:
7191: + dm - The `DM` object
7192: - n  - the label number

7194:   Output Parameter:
7195: . label - the label

7197:   Level: intermediate

7199: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMAddLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7200: @*/
7201: PetscErrorCode DMGetLabelByNum(DM dm, PetscInt n, DMLabel *label)
7202: {
7203:   DMLabelLink next = dm->labels;
7204:   PetscInt    l    = 0;

7206:   PetscFunctionBegin;
7208:   PetscAssertPointer(label, 3);
7209:   while (next) {
7210:     if (l == n) {
7211:       *label = next->label;
7212:       PetscFunctionReturn(PETSC_SUCCESS);
7213:     }
7214:     ++l;
7215:     next = next->next;
7216:   }
7217:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %" PetscInt_FMT " does not exist in this DM", n);
7218: }

7220: /*@
7221:   DMAddLabel - Add the label to this `DM`

7223:   Not Collective

7225:   Input Parameters:
7226: + dm    - The `DM` object
7227: - label - The `DMLabel`

7229:   Level: developer

7231: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7232: @*/
7233: PetscErrorCode DMAddLabel(DM dm, DMLabel label)
7234: {
7235:   DMLabelLink l, *p, tmpLabel;
7236:   PetscBool   hasLabel;
7237:   const char *lname;
7238:   PetscBool   flg;

7240:   PetscFunctionBegin;
7242:   PetscCall(PetscObjectGetName((PetscObject)label, &lname));
7243:   PetscCall(DMHasLabel(dm, lname, &hasLabel));
7244:   PetscCheck(!hasLabel, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %s already exists in this DM", lname);
7245:   PetscCall(PetscCalloc1(1, &tmpLabel));
7246:   tmpLabel->label  = label;
7247:   tmpLabel->output = PETSC_TRUE;
7248:   for (p = &dm->labels; (l = *p); p = &l->next) { }
7249:   *p = tmpLabel;
7250:   PetscCall(PetscObjectReference((PetscObject)label));
7251:   PetscCall(PetscStrcmp(lname, "depth", &flg));
7252:   if (flg) dm->depthLabel = label;
7253:   PetscCall(PetscStrcmp(lname, "celltype", &flg));
7254:   if (flg) dm->celltypeLabel = label;
7255:   PetscFunctionReturn(PETSC_SUCCESS);
7256: }

7258: // PetscClangLinter pragma ignore: -fdoc-section-header-unknown
7259: /*@
7260:   DMSetLabel - Replaces the label of a given name, or ignores it if the name is not present

7262:   Not Collective

7264:   Input Parameters:
7265: + dm    - The `DM` object
7266: - label - The `DMLabel`, having the same name, to substitute

7268:   Default labels in a `DMPLEX`:
7269: + "depth"       - Holds the depth (co-dimension) of each mesh point
7270: . "celltype"    - Holds the topological type of each cell
7271: . "ghost"       - If the DM is distributed with overlap, this marks the cells and faces in the overlap
7272: . "Cell Sets"   - Mirrors the cell sets defined by GMsh and ExodusII
7273: . "Face Sets"   - Mirrors the face sets defined by GMsh and ExodusII
7274: - "Vertex Sets" - Mirrors the vertex sets defined by GMsh

7276:   Level: intermediate

7278: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMPlexGetDepthLabel()`, `DMPlexGetCellType()`
7279: @*/
7280: PetscErrorCode DMSetLabel(DM dm, DMLabel label)
7281: {
7282:   DMLabelLink next = dm->labels;
7283:   PetscBool   hasLabel, flg;
7284:   const char *name, *lname;

7286:   PetscFunctionBegin;
7289:   PetscCall(PetscObjectGetName((PetscObject)label, &name));
7290:   while (next) {
7291:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7292:     PetscCall(PetscStrcmp(name, lname, &hasLabel));
7293:     if (hasLabel) {
7294:       PetscCall(PetscObjectReference((PetscObject)label));
7295:       PetscCall(PetscStrcmp(lname, "depth", &flg));
7296:       if (flg) dm->depthLabel = label;
7297:       PetscCall(PetscStrcmp(lname, "celltype", &flg));
7298:       if (flg) dm->celltypeLabel = label;
7299:       PetscCall(DMLabelDestroy(&next->label));
7300:       next->label = label;
7301:       break;
7302:     }
7303:     next = next->next;
7304:   }
7305:   PetscFunctionReturn(PETSC_SUCCESS);
7306: }

7308: /*@
7309:   DMRemoveLabel - Remove the label given by name from this `DM`

7311:   Not Collective

7313:   Input Parameters:
7314: + dm   - The `DM` object
7315: - name - The label name

7317:   Output Parameter:
7318: . label - The `DMLabel`, or `NULL` if the label is absent. Pass in `NULL` to call `DMLabelDestroy()` on the label, otherwise the
7319:           caller is responsible for calling `DMLabelDestroy()`.

7321:   Level: developer

7323: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMLabelDestroy()`, `DMRemoveLabelBySelf()`
7324: @*/
7325: PetscErrorCode DMRemoveLabel(DM dm, const char name[], DMLabel *label)
7326: {
7327:   DMLabelLink link, *pnext;
7328:   PetscBool   hasLabel;
7329:   const char *lname;

7331:   PetscFunctionBegin;
7333:   PetscAssertPointer(name, 2);
7334:   if (label) {
7335:     PetscAssertPointer(label, 3);
7336:     *label = NULL;
7337:   }
7338:   for (pnext = &dm->labels; (link = *pnext); pnext = &link->next) {
7339:     PetscCall(PetscObjectGetName((PetscObject)link->label, &lname));
7340:     PetscCall(PetscStrcmp(name, lname, &hasLabel));
7341:     if (hasLabel) {
7342:       *pnext = link->next; /* Remove from list */
7343:       PetscCall(PetscStrcmp(name, "depth", &hasLabel));
7344:       if (hasLabel) dm->depthLabel = NULL;
7345:       PetscCall(PetscStrcmp(name, "celltype", &hasLabel));
7346:       if (hasLabel) dm->celltypeLabel = NULL;
7347:       if (label) *label = link->label;
7348:       else PetscCall(DMLabelDestroy(&link->label));
7349:       PetscCall(PetscFree(link));
7350:       break;
7351:     }
7352:   }
7353:   PetscFunctionReturn(PETSC_SUCCESS);
7354: }

7356: /*@
7357:   DMRemoveLabelBySelf - Remove the label from this `DM`

7359:   Not Collective

7361:   Input Parameters:
7362: + dm           - The `DM` object
7363: . label        - The `DMLabel` to be removed from the `DM`
7364: - failNotFound - Should it fail if the label is not found in the `DM`?

7366:   Level: developer

7368:   Note:
7369:   Only exactly the same instance is removed if found, name match is ignored.
7370:   If the `DM` has an exclusive reference to the label, the label gets destroyed and
7371:   *label nullified.

7373: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabel()` `DMGetLabelValue()`, `DMSetLabelValue()`, `DMLabelDestroy()`, `DMRemoveLabel()`
7374: @*/
7375: PetscErrorCode DMRemoveLabelBySelf(DM dm, DMLabel *label, PetscBool failNotFound)
7376: {
7377:   DMLabelLink link, *pnext;
7378:   PetscBool   hasLabel = PETSC_FALSE;

7380:   PetscFunctionBegin;
7382:   PetscAssertPointer(label, 2);
7383:   if (!*label && !failNotFound) PetscFunctionReturn(PETSC_SUCCESS);
7386:   for (pnext = &dm->labels; (link = *pnext); pnext = &link->next) {
7387:     if (*label == link->label) {
7388:       hasLabel = PETSC_TRUE;
7389:       *pnext   = link->next; /* Remove from list */
7390:       if (*label == dm->depthLabel) dm->depthLabel = NULL;
7391:       if (*label == dm->celltypeLabel) dm->celltypeLabel = NULL;
7392:       if (((PetscObject)link->label)->refct < 2) *label = NULL; /* nullify if exclusive reference */
7393:       PetscCall(DMLabelDestroy(&link->label));
7394:       PetscCall(PetscFree(link));
7395:       break;
7396:     }
7397:   }
7398:   PetscCheck(hasLabel || !failNotFound, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Given label not found in DM");
7399:   PetscFunctionReturn(PETSC_SUCCESS);
7400: }

7402: /*@
7403:   DMGetLabelOutput - Get the output flag for a given label

7405:   Not Collective

7407:   Input Parameters:
7408: + dm   - The `DM` object
7409: - name - The label name

7411:   Output Parameter:
7412: . output - The flag for output

7414:   Level: developer

7416: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMSetLabelOutput()`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7417: @*/
7418: PetscErrorCode DMGetLabelOutput(DM dm, const char name[], PetscBool *output)
7419: {
7420:   DMLabelLink next = dm->labels;
7421:   const char *lname;

7423:   PetscFunctionBegin;
7425:   PetscAssertPointer(name, 2);
7426:   PetscAssertPointer(output, 3);
7427:   while (next) {
7428:     PetscBool flg;

7430:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7431:     PetscCall(PetscStrcmp(name, lname, &flg));
7432:     if (flg) {
7433:       *output = next->output;
7434:       PetscFunctionReturn(PETSC_SUCCESS);
7435:     }
7436:     next = next->next;
7437:   }
7438:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No label named %s was present in this dm", name);
7439: }

7441: /*@
7442:   DMSetLabelOutput - Set if a given label should be saved to a `PetscViewer` in calls to `DMView()`

7444:   Not Collective

7446:   Input Parameters:
7447: + dm     - The `DM` object
7448: . name   - The label name
7449: - output - `PETSC_TRUE` to save the label to the viewer

7451:   Level: developer

7453: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetOutputFlag()`, `DMGetLabelOutput()`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7454: @*/
7455: PetscErrorCode DMSetLabelOutput(DM dm, const char name[], PetscBool output)
7456: {
7457:   DMLabelLink next = dm->labels;
7458:   const char *lname;

7460:   PetscFunctionBegin;
7462:   PetscAssertPointer(name, 2);
7463:   while (next) {
7464:     PetscBool flg;

7466:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7467:     PetscCall(PetscStrcmp(name, lname, &flg));
7468:     if (flg) {
7469:       next->output = output;
7470:       PetscFunctionReturn(PETSC_SUCCESS);
7471:     }
7472:     next = next->next;
7473:   }
7474:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No label named %s was present in this dm", name);
7475: }

7477: /*@
7478:   DMCopyLabels - Copy labels from one `DM` mesh to another `DM` with a superset of the points

7480:   Collective

7482:   Input Parameters:
7483: + dmA   - The `DM` object with initial labels
7484: . dmB   - The `DM` object to which labels are copied
7485: . mode  - Copy labels by pointers (`PETSC_OWN_POINTER`) or duplicate them (`PETSC_COPY_VALUES`)
7486: . all   - Copy all labels including "depth", "dim", and "celltype" (`PETSC_TRUE`) which are otherwise ignored (`PETSC_FALSE`)
7487: - emode - How to behave when a `DMLabel` in the source and destination `DM`s with the same name is encountered (see `DMCopyLabelsMode`)

7489:   Level: intermediate

7491:   Note:
7492:   This is typically used when interpolating or otherwise adding to a mesh, or testing.

7494: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMAddLabel()`, `DMCopyLabelsMode`
7495: @*/
7496: PetscErrorCode DMCopyLabels(DM dmA, DM dmB, PetscCopyMode mode, PetscBool all, DMCopyLabelsMode emode)
7497: {
7498:   DMLabel     label, labelNew, labelOld;
7499:   const char *name;
7500:   PetscBool   flg;
7501:   DMLabelLink link;

7503:   PetscFunctionBegin;
7508:   PetscCheck(mode != PETSC_USE_POINTER, PetscObjectComm((PetscObject)dmA), PETSC_ERR_SUP, "PETSC_USE_POINTER not supported for objects");
7509:   if (dmA == dmB) PetscFunctionReturn(PETSC_SUCCESS);
7510:   for (link = dmA->labels; link; link = link->next) {
7511:     label = link->label;
7512:     PetscCall(PetscObjectGetName((PetscObject)label, &name));
7513:     if (!all) {
7514:       PetscCall(PetscStrcmp(name, "depth", &flg));
7515:       if (flg) continue;
7516:       PetscCall(PetscStrcmp(name, "dim", &flg));
7517:       if (flg) continue;
7518:       PetscCall(PetscStrcmp(name, "celltype", &flg));
7519:       if (flg) continue;
7520:     }
7521:     PetscCall(DMGetLabel(dmB, name, &labelOld));
7522:     if (labelOld) {
7523:       switch (emode) {
7524:       case DM_COPY_LABELS_KEEP:
7525:         continue;
7526:       case DM_COPY_LABELS_REPLACE:
7527:         PetscCall(DMRemoveLabelBySelf(dmB, &labelOld, PETSC_TRUE));
7528:         break;
7529:       case DM_COPY_LABELS_FAIL:
7530:         SETERRQ(PetscObjectComm((PetscObject)dmA), PETSC_ERR_ARG_OUTOFRANGE, "Label %s already exists in destination DM", name);
7531:       default:
7532:         SETERRQ(PetscObjectComm((PetscObject)dmA), PETSC_ERR_ARG_OUTOFRANGE, "Unhandled DMCopyLabelsMode %d", (int)emode);
7533:       }
7534:     }
7535:     if (mode == PETSC_COPY_VALUES) {
7536:       PetscCall(DMLabelDuplicate(label, &labelNew));
7537:     } else {
7538:       labelNew = label;
7539:     }
7540:     PetscCall(DMAddLabel(dmB, labelNew));
7541:     if (mode == PETSC_COPY_VALUES) PetscCall(DMLabelDestroy(&labelNew));
7542:   }
7543:   PetscFunctionReturn(PETSC_SUCCESS);
7544: }

7546: /*@C
7547:   DMCompareLabels - Compare labels between two `DM` objects

7549:   Collective; No Fortran Support

7551:   Input Parameters:
7552: + dm0 - First `DM` object
7553: - dm1 - Second `DM` object

7555:   Output Parameters:
7556: + equal   - (Optional) Flag whether labels of dm0 and dm1 are the same
7557: - message - (Optional) Message describing the difference, or `NULL` if there is no difference

7559:   Level: intermediate

7561:   Notes:
7562:   The output flag equal will be the same on all processes.

7564:   If equal is passed as `NULL` and difference is found, an error is thrown on all processes.

7566:   Make sure to pass equal is `NULL` on all processes or none of them.

7568:   The output message is set independently on each rank.

7570:   message must be freed with `PetscFree()`

7572:   If message is passed as `NULL` and a difference is found, the difference description is printed to stderr in synchronized manner.

7574:   Make sure to pass message as `NULL` on all processes or no processes.

7576:   Labels are matched by name. If the number of labels and their names are equal,
7577:   `DMLabelCompare()` is used to compare each pair of labels with the same name.

7579:   Developer Note:
7580:   Can automatically generate the Fortran stub because `message` must be freed with `PetscFree()`

7582: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMAddLabel()`, `DMCopyLabelsMode`, `DMLabelCompare()`
7583: @*/
7584: PetscErrorCode DMCompareLabels(DM dm0, DM dm1, PetscBool *equal, char **message)
7585: {
7586:   PetscInt    n, i;
7587:   char        msg[PETSC_MAX_PATH_LEN] = "";
7588:   PetscBool   eq;
7589:   MPI_Comm    comm;
7590:   PetscMPIInt rank;

7592:   PetscFunctionBegin;
7595:   PetscCheckSameComm(dm0, 1, dm1, 2);
7596:   if (equal) PetscAssertPointer(equal, 3);
7597:   if (message) PetscAssertPointer(message, 4);
7598:   PetscCall(PetscObjectGetComm((PetscObject)dm0, &comm));
7599:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7600:   {
7601:     PetscInt n1;

7603:     PetscCall(DMGetNumLabels(dm0, &n));
7604:     PetscCall(DMGetNumLabels(dm1, &n1));
7605:     eq = (PetscBool)(n == n1);
7606:     if (!eq) PetscCall(PetscSNPrintf(msg, sizeof(msg), "Number of labels in dm0 = %" PetscInt_FMT " != %" PetscInt_FMT " = Number of labels in dm1", n, n1));
7607:     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPIU_BOOL, MPI_LAND, comm));
7608:     if (!eq) goto finish;
7609:   }
7610:   for (i = 0; i < n; i++) {
7611:     DMLabel     l0, l1;
7612:     const char *name;
7613:     char       *msgInner;

7615:     /* Ignore label order */
7616:     PetscCall(DMGetLabelByNum(dm0, i, &l0));
7617:     PetscCall(PetscObjectGetName((PetscObject)l0, &name));
7618:     PetscCall(DMGetLabel(dm1, name, &l1));
7619:     if (!l1) {
7620:       PetscCall(PetscSNPrintf(msg, sizeof(msg), "Label \"%s\" (#%" PetscInt_FMT " in dm0) not found in dm1", name, i));
7621:       eq = PETSC_FALSE;
7622:       break;
7623:     }
7624:     PetscCall(DMLabelCompare(comm, l0, l1, &eq, &msgInner));
7625:     PetscCall(PetscStrncpy(msg, msgInner, sizeof(msg)));
7626:     PetscCall(PetscFree(msgInner));
7627:     if (!eq) break;
7628:   }
7629:   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPIU_BOOL, MPI_LAND, comm));
7630: finish:
7631:   /* If message output arg not set, print to stderr */
7632:   if (message) {
7633:     *message = NULL;
7634:     if (msg[0]) PetscCall(PetscStrallocpy(msg, message));
7635:   } else {
7636:     if (msg[0]) PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] %s\n", rank, msg));
7637:     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR));
7638:   }
7639:   /* If same output arg not ser and labels are not equal, throw error */
7640:   if (equal) *equal = eq;
7641:   else PetscCheck(eq, comm, PETSC_ERR_ARG_INCOMP, "DMLabels are not the same in dm0 and dm1");
7642:   PetscFunctionReturn(PETSC_SUCCESS);
7643: }

7645: PetscErrorCode DMSetLabelValue_Fast(DM dm, DMLabel *label, const char name[], PetscInt point, PetscInt value)
7646: {
7647:   PetscFunctionBegin;
7648:   PetscAssertPointer(label, 2);
7649:   if (!*label) {
7650:     PetscCall(DMCreateLabel(dm, name));
7651:     PetscCall(DMGetLabel(dm, name, label));
7652:   }
7653:   PetscCall(DMLabelSetValue(*label, point, value));
7654:   PetscFunctionReturn(PETSC_SUCCESS);
7655: }

7657: /*
7658:   Many mesh programs, such as Triangle and TetGen, allow only a single label for each mesh point. Therefore, we would
7659:   like to encode all label IDs using a single, universal label. We can do this by assigning an integer to every
7660:   (label, id) pair in the DM.

7662:   However, a mesh point can have multiple labels, so we must separate all these values. We will assign a bit range to
7663:   each label.
7664: */
7665: PetscErrorCode DMUniversalLabelCreate(DM dm, DMUniversalLabel *universal)
7666: {
7667:   DMUniversalLabel ul;
7668:   PetscBool       *active;
7669:   PetscInt         pStart, pEnd, p, Nl, l, m;

7671:   PetscFunctionBegin;
7672:   PetscCall(PetscMalloc1(1, &ul));
7673:   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "universal", &ul->label));
7674:   PetscCall(DMGetNumLabels(dm, &Nl));
7675:   PetscCall(PetscCalloc1(Nl, &active));
7676:   ul->Nl = 0;
7677:   for (l = 0; l < Nl; ++l) {
7678:     PetscBool   isdepth, iscelltype;
7679:     const char *name;

7681:     PetscCall(DMGetLabelName(dm, l, &name));
7682:     PetscCall(PetscStrncmp(name, "depth", 6, &isdepth));
7683:     PetscCall(PetscStrncmp(name, "celltype", 9, &iscelltype));
7684:     active[l] = !(isdepth || iscelltype) ? PETSC_TRUE : PETSC_FALSE;
7685:     if (active[l]) ++ul->Nl;
7686:   }
7687:   PetscCall(PetscCalloc5(ul->Nl, &ul->names, ul->Nl, &ul->indices, ul->Nl + 1, &ul->offsets, ul->Nl + 1, &ul->bits, ul->Nl, &ul->masks));
7688:   ul->Nv = 0;
7689:   for (l = 0, m = 0; l < Nl; ++l) {
7690:     DMLabel     label;
7691:     PetscInt    nv;
7692:     const char *name;

7694:     if (!active[l]) continue;
7695:     PetscCall(DMGetLabelName(dm, l, &name));
7696:     PetscCall(DMGetLabelByNum(dm, l, &label));
7697:     PetscCall(DMLabelGetNumValues(label, &nv));
7698:     PetscCall(PetscStrallocpy(name, &ul->names[m]));
7699:     ul->indices[m] = l;
7700:     ul->Nv += nv;
7701:     ul->offsets[m + 1] = nv;
7702:     ul->bits[m + 1]    = PetscCeilReal(PetscLog2Real(nv + 1));
7703:     ++m;
7704:   }
7705:   for (l = 1; l <= ul->Nl; ++l) {
7706:     ul->offsets[l] = ul->offsets[l - 1] + ul->offsets[l];
7707:     ul->bits[l]    = ul->bits[l - 1] + ul->bits[l];
7708:   }
7709:   for (l = 0; l < ul->Nl; ++l) {
7710:     PetscInt b;

7712:     ul->masks[l] = 0;
7713:     for (b = ul->bits[l]; b < ul->bits[l + 1]; ++b) ul->masks[l] |= 1 << b;
7714:   }
7715:   PetscCall(PetscMalloc1(ul->Nv, &ul->values));
7716:   for (l = 0, m = 0; l < Nl; ++l) {
7717:     DMLabel         label;
7718:     IS              valueIS;
7719:     const PetscInt *varr;
7720:     PetscInt        nv, v;

7722:     if (!active[l]) continue;
7723:     PetscCall(DMGetLabelByNum(dm, l, &label));
7724:     PetscCall(DMLabelGetNumValues(label, &nv));
7725:     PetscCall(DMLabelGetValueIS(label, &valueIS));
7726:     PetscCall(ISGetIndices(valueIS, &varr));
7727:     for (v = 0; v < nv; ++v) ul->values[ul->offsets[m] + v] = varr[v];
7728:     PetscCall(ISRestoreIndices(valueIS, &varr));
7729:     PetscCall(ISDestroy(&valueIS));
7730:     PetscCall(PetscSortInt(nv, &ul->values[ul->offsets[m]]));
7731:     ++m;
7732:   }
7733:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
7734:   for (p = pStart; p < pEnd; ++p) {
7735:     PetscInt  uval   = 0;
7736:     PetscBool marked = PETSC_FALSE;

7738:     for (l = 0, m = 0; l < Nl; ++l) {
7739:       DMLabel  label;
7740:       PetscInt val, defval, loc, nv;

7742:       if (!active[l]) continue;
7743:       PetscCall(DMGetLabelByNum(dm, l, &label));
7744:       PetscCall(DMLabelGetValue(label, p, &val));
7745:       PetscCall(DMLabelGetDefaultValue(label, &defval));
7746:       if (val == defval) {
7747:         ++m;
7748:         continue;
7749:       }
7750:       nv     = ul->offsets[m + 1] - ul->offsets[m];
7751:       marked = PETSC_TRUE;
7752:       PetscCall(PetscFindInt(val, nv, &ul->values[ul->offsets[m]], &loc));
7753:       PetscCheck(loc >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Label value %" PetscInt_FMT " not found in compression array", val);
7754:       uval += (loc + 1) << ul->bits[m];
7755:       ++m;
7756:     }
7757:     if (marked) PetscCall(DMLabelSetValue(ul->label, p, uval));
7758:   }
7759:   PetscCall(PetscFree(active));
7760:   *universal = ul;
7761:   PetscFunctionReturn(PETSC_SUCCESS);
7762: }

7764: PetscErrorCode DMUniversalLabelDestroy(DMUniversalLabel *universal)
7765: {
7766:   PetscInt l;

7768:   PetscFunctionBegin;
7769:   for (l = 0; l < (*universal)->Nl; ++l) PetscCall(PetscFree((*universal)->names[l]));
7770:   PetscCall(DMLabelDestroy(&(*universal)->label));
7771:   PetscCall(PetscFree5((*universal)->names, (*universal)->indices, (*universal)->offsets, (*universal)->bits, (*universal)->masks));
7772:   PetscCall(PetscFree((*universal)->values));
7773:   PetscCall(PetscFree(*universal));
7774:   *universal = NULL;
7775:   PetscFunctionReturn(PETSC_SUCCESS);
7776: }

7778: PetscErrorCode DMUniversalLabelGetLabel(DMUniversalLabel ul, DMLabel *ulabel)
7779: {
7780:   PetscFunctionBegin;
7781:   PetscAssertPointer(ulabel, 2);
7782:   *ulabel = ul->label;
7783:   PetscFunctionReturn(PETSC_SUCCESS);
7784: }

7786: PetscErrorCode DMUniversalLabelCreateLabels(DMUniversalLabel ul, PetscBool preserveOrder, DM dm)
7787: {
7788:   PetscInt Nl = ul->Nl, l;

7790:   PetscFunctionBegin;
7792:   for (l = 0; l < Nl; ++l) {
7793:     if (preserveOrder) PetscCall(DMCreateLabelAtIndex(dm, ul->indices[l], ul->names[l]));
7794:     else PetscCall(DMCreateLabel(dm, ul->names[l]));
7795:   }
7796:   if (preserveOrder) {
7797:     for (l = 0; l < ul->Nl; ++l) {
7798:       const char *name;
7799:       PetscBool   match;

7801:       PetscCall(DMGetLabelName(dm, ul->indices[l], &name));
7802:       PetscCall(PetscStrcmp(name, ul->names[l], &match));
7803:       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]);
7804:     }
7805:   }
7806:   PetscFunctionReturn(PETSC_SUCCESS);
7807: }

7809: PetscErrorCode DMUniversalLabelSetLabelValue(DMUniversalLabel ul, DM dm, PetscBool useIndex, PetscInt p, PetscInt value)
7810: {
7811:   PetscInt l;

7813:   PetscFunctionBegin;
7814:   for (l = 0; l < ul->Nl; ++l) {
7815:     DMLabel  label;
7816:     PetscInt lval = (value & ul->masks[l]) >> ul->bits[l];

7818:     if (lval) {
7819:       if (useIndex) PetscCall(DMGetLabelByNum(dm, ul->indices[l], &label));
7820:       else PetscCall(DMGetLabel(dm, ul->names[l], &label));
7821:       PetscCall(DMLabelSetValue(label, p, ul->values[ul->offsets[l] + lval - 1]));
7822:     }
7823:   }
7824:   PetscFunctionReturn(PETSC_SUCCESS);
7825: }

7827: /*@
7828:   DMGetCoarseDM - Get the coarse `DM`from which this `DM` was obtained by refinement

7830:   Not Collective

7832:   Input Parameter:
7833: . dm - The `DM` object

7835:   Output Parameter:
7836: . cdm - The coarse `DM`

7838:   Level: intermediate

7840: .seealso: [](ch_dmbase), `DM`, `DMSetCoarseDM()`, `DMCoarsen()`
7841: @*/
7842: PetscErrorCode DMGetCoarseDM(DM dm, DM *cdm)
7843: {
7844:   PetscFunctionBegin;
7846:   PetscAssertPointer(cdm, 2);
7847:   *cdm = dm->coarseMesh;
7848:   PetscFunctionReturn(PETSC_SUCCESS);
7849: }

7851: /*@
7852:   DMSetCoarseDM - Set the coarse `DM` from which this `DM` was obtained by refinement

7854:   Input Parameters:
7855: + dm  - The `DM` object
7856: - cdm - The coarse `DM`

7858:   Level: intermediate

7860:   Note:
7861:   Normally this is set automatically by `DMRefine()`

7863: .seealso: [](ch_dmbase), `DM`, `DMGetCoarseDM()`, `DMCoarsen()`, `DMSetRefine()`, `DMSetFineDM()`
7864: @*/
7865: PetscErrorCode DMSetCoarseDM(DM dm, DM cdm)
7866: {
7867:   PetscFunctionBegin;
7870:   if (dm == cdm) cdm = NULL;
7871:   PetscCall(PetscObjectReference((PetscObject)cdm));
7872:   PetscCall(DMDestroy(&dm->coarseMesh));
7873:   dm->coarseMesh = cdm;
7874:   PetscFunctionReturn(PETSC_SUCCESS);
7875: }

7877: /*@
7878:   DMGetFineDM - Get the fine mesh from which this `DM` was obtained by coarsening

7880:   Input Parameter:
7881: . dm - The `DM` object

7883:   Output Parameter:
7884: . fdm - The fine `DM`

7886:   Level: intermediate

7888: .seealso: [](ch_dmbase), `DM`, `DMSetFineDM()`, `DMCoarsen()`, `DMRefine()`
7889: @*/
7890: PetscErrorCode DMGetFineDM(DM dm, DM *fdm)
7891: {
7892:   PetscFunctionBegin;
7894:   PetscAssertPointer(fdm, 2);
7895:   *fdm = dm->fineMesh;
7896:   PetscFunctionReturn(PETSC_SUCCESS);
7897: }

7899: /*@
7900:   DMSetFineDM - Set the fine mesh from which this was obtained by coarsening

7902:   Input Parameters:
7903: + dm  - The `DM` object
7904: - fdm - The fine `DM`

7906:   Level: developer

7908:   Note:
7909:   Normally this is set automatically by `DMCoarsen()`

7911: .seealso: [](ch_dmbase), `DM`, `DMGetFineDM()`, `DMCoarsen()`, `DMRefine()`
7912: @*/
7913: PetscErrorCode DMSetFineDM(DM dm, DM fdm)
7914: {
7915:   PetscFunctionBegin;
7918:   if (dm == fdm) fdm = NULL;
7919:   PetscCall(PetscObjectReference((PetscObject)fdm));
7920:   PetscCall(DMDestroy(&dm->fineMesh));
7921:   dm->fineMesh = fdm;
7922:   PetscFunctionReturn(PETSC_SUCCESS);
7923: }

7925: /*@C
7926:   DMAddBoundary - Add a boundary condition to a model represented by a `DM`

7928:   Collective

7930:   Input Parameters:
7931: + dm       - The `DM`, with a `PetscDS` that matches the problem being constrained
7932: . type     - The type of condition, e.g. `DM_BC_ESSENTIAL_ANALYTIC`, `DM_BC_ESSENTIAL_FIELD` (Dirichlet), or `DM_BC_NATURAL` (Neumann)
7933: . name     - The BC name
7934: . label    - The label defining constrained points
7935: . Nv       - The number of `DMLabel` values for constrained points
7936: . values   - An array of values for constrained points
7937: . field    - The field to constrain
7938: . Nc       - The number of constrained field components (0 will constrain all fields)
7939: . comps    - An array of constrained component numbers
7940: . bcFunc   - A pointwise function giving boundary values
7941: . bcFunc_t - A pointwise function giving the time deriative of the boundary values, or NULL
7942: - ctx      - An optional user context for bcFunc

7944:   Output Parameter:
7945: . bd - (Optional) Boundary number

7947:   Options Database Keys:
7948: + -bc_<boundary name> <num>      - Overrides the boundary ids
7949: - -bc_<boundary name>_comp <num> - Overrides the boundary components

7951:   Level: intermediate

7953:   Notes:
7954:   Both bcFunc and bcFunc_t will depend on the boundary condition type. If the type if `DM_BC_ESSENTIAL`, then the calling sequence is\:
7955: .vb
7956:  void bcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[])
7957: .ve

7959:   If the type is `DM_BC_ESSENTIAL_FIELD` or other _FIELD value, then the calling sequence is\:

7961: .vb
7962:   void bcFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
7963:               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
7964:               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
7965:               PetscReal time, const PetscReal x[], PetscScalar bcval[])
7966: .ve
7967: + dim - the spatial dimension
7968: . Nf - the number of fields
7969: . uOff - the offset into u[] and u_t[] for each field
7970: . uOff_x - the offset into u_x[] for each field
7971: . u - each field evaluated at the current point
7972: . u_t - the time derivative of each field evaluated at the current point
7973: . u_x - the gradient of each field evaluated at the current point
7974: . aOff - the offset into a[] and a_t[] for each auxiliary field
7975: . aOff_x - the offset into a_x[] for each auxiliary field
7976: . a - each auxiliary field evaluated at the current point
7977: . a_t - the time derivative of each auxiliary field evaluated at the current point
7978: . a_x - the gradient of auxiliary each field evaluated at the current point
7979: . t - current time
7980: . x - coordinates of the current point
7981: . numConstants - number of constant parameters
7982: . constants - constant parameters
7983: - bcval - output values at the current point

7985: .seealso: [](ch_dmbase), `DM`, `DSGetBoundary()`, `PetscDSAddBoundary()`
7986: @*/
7987: 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)
7988: {
7989:   PetscDS ds;

7991:   PetscFunctionBegin;
7998:   PetscCheck(!dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot add boundary to DM after creating local section");
7999:   PetscCall(DMGetDS(dm, &ds));
8000:   /* Complete label */
8001:   if (label) {
8002:     PetscObject  obj;
8003:     PetscClassId id;

8005:     PetscCall(DMGetField(dm, field, NULL, &obj));
8006:     PetscCall(PetscObjectGetClassId(obj, &id));
8007:     if (id == PETSCFE_CLASSID) {
8008:       DM plex;

8010:       PetscCall(DMConvert(dm, DMPLEX, &plex));
8011:       if (plex) PetscCall(DMPlexLabelComplete(plex, label));
8012:       PetscCall(DMDestroy(&plex));
8013:     }
8014:   }
8015:   PetscCall(PetscDSAddBoundary(ds, type, name, label, Nv, values, field, Nc, comps, bcFunc, bcFunc_t, ctx, bd));
8016:   PetscFunctionReturn(PETSC_SUCCESS);
8017: }

8019: /* TODO Remove this since now the structures are the same */
8020: static PetscErrorCode DMPopulateBoundary(DM dm)
8021: {
8022:   PetscDS     ds;
8023:   DMBoundary *lastnext;
8024:   DSBoundary  dsbound;

8026:   PetscFunctionBegin;
8027:   PetscCall(DMGetDS(dm, &ds));
8028:   dsbound = ds->boundary;
8029:   if (dm->boundary) {
8030:     DMBoundary next = dm->boundary;

8032:     /* quick check to see if the PetscDS has changed */
8033:     if (next->dsboundary == dsbound) PetscFunctionReturn(PETSC_SUCCESS);
8034:     /* the PetscDS has changed: tear down and rebuild */
8035:     while (next) {
8036:       DMBoundary b = next;

8038:       next = b->next;
8039:       PetscCall(PetscFree(b));
8040:     }
8041:     dm->boundary = NULL;
8042:   }

8044:   lastnext = &dm->boundary;
8045:   while (dsbound) {
8046:     DMBoundary dmbound;

8048:     PetscCall(PetscNew(&dmbound));
8049:     dmbound->dsboundary = dsbound;
8050:     dmbound->label      = dsbound->label;
8051:     /* push on the back instead of the front so that it is in the same order as in the PetscDS */
8052:     *lastnext = dmbound;
8053:     lastnext  = &dmbound->next;
8054:     dsbound   = dsbound->next;
8055:   }
8056:   PetscFunctionReturn(PETSC_SUCCESS);
8057: }

8059: /* TODO: missing manual page */
8060: PetscErrorCode DMIsBoundaryPoint(DM dm, PetscInt point, PetscBool *isBd)
8061: {
8062:   DMBoundary b;

8064:   PetscFunctionBegin;
8066:   PetscAssertPointer(isBd, 3);
8067:   *isBd = PETSC_FALSE;
8068:   PetscCall(DMPopulateBoundary(dm));
8069:   b = dm->boundary;
8070:   while (b && !*isBd) {
8071:     DMLabel    label = b->label;
8072:     DSBoundary dsb   = b->dsboundary;
8073:     PetscInt   i;

8075:     if (label) {
8076:       for (i = 0; i < dsb->Nv && !*isBd; ++i) PetscCall(DMLabelStratumHasPoint(label, dsb->values[i], point, isBd));
8077:     }
8078:     b = b->next;
8079:   }
8080:   PetscFunctionReturn(PETSC_SUCCESS);
8081: }

8083: /*@C
8084:   DMProjectFunction - This projects the given function into the function space provided by a `DM`, putting the coefficients in a global vector.

8086:   Collective

8088:   Input Parameters:
8089: + dm    - The `DM`
8090: . time  - The time
8091: . funcs - The coordinate functions to evaluate, one per field
8092: . ctxs  - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
8093: - mode  - The insertion mode for values

8095:   Output Parameter:
8096: . X - vector

8098:   Calling sequence of `funcs`:
8099: + dim  - The spatial dimension
8100: . time - The time at which to sample
8101: . x    - The coordinates
8102: . Nc   - The number of components
8103: . u    - The output field values
8104: - ctx  - optional user-defined function context

8106:   Level: developer

8108:   Developer Notes:
8109:   This API is specific to only particular usage of `DM`

8111:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8113: .seealso: [](ch_dmbase), `DM`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
8114: @*/
8115: 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)
8116: {
8117:   Vec localX;

8119:   PetscFunctionBegin;
8121:   PetscCall(PetscLogEventBegin(DM_ProjectFunction, dm, X, 0, 0));
8122:   PetscCall(DMGetLocalVector(dm, &localX));
8123:   PetscCall(VecSet(localX, 0.));
8124:   PetscCall(DMProjectFunctionLocal(dm, time, funcs, ctxs, mode, localX));
8125:   PetscCall(DMLocalToGlobalBegin(dm, localX, mode, X));
8126:   PetscCall(DMLocalToGlobalEnd(dm, localX, mode, X));
8127:   PetscCall(DMRestoreLocalVector(dm, &localX));
8128:   PetscCall(PetscLogEventEnd(DM_ProjectFunction, dm, X, 0, 0));
8129:   PetscFunctionReturn(PETSC_SUCCESS);
8130: }

8132: /*@C
8133:   DMProjectFunctionLocal - This projects the given function into the function space provided by a `DM`, putting the coefficients in a local vector.

8135:   Not Collective

8137:   Input Parameters:
8138: + dm    - The `DM`
8139: . time  - The time
8140: . funcs - The coordinate functions to evaluate, one per field
8141: . ctxs  - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
8142: - mode  - The insertion mode for values

8144:   Output Parameter:
8145: . localX - vector

8147:   Calling sequence of `funcs`:
8148: + dim  - The spatial dimension
8149: . time - The current timestep
8150: . x    - The coordinates
8151: . Nc   - The number of components
8152: . u    - The output field values
8153: - ctx  - optional user-defined function context

8155:   Level: developer

8157:   Developer Notes:
8158:   This API is specific to only particular usage of `DM`

8160:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8162: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
8163: @*/
8164: 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)
8165: {
8166:   PetscFunctionBegin;
8169:   PetscUseTypeMethod(dm, projectfunctionlocal, time, funcs, ctxs, mode, localX);
8170:   PetscFunctionReturn(PETSC_SUCCESS);
8171: }

8173: /*@C
8174:   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.

8176:   Collective

8178:   Input Parameters:
8179: + dm     - The `DM`
8180: . time   - The time
8181: . numIds - The number of ids
8182: . ids    - The ids
8183: . Nc     - The number of components
8184: . comps  - The components
8185: . label  - The `DMLabel` selecting the portion of the mesh for projection
8186: . funcs  - The coordinate functions to evaluate, one per field
8187: . ctxs   - Optional array of contexts to pass to each coordinate function.  ctxs may be null.
8188: - mode   - The insertion mode for values

8190:   Output Parameter:
8191: . X - vector

8193:   Calling sequence of `funcs`:
8194: + dim  - The spatial dimension
8195: . time - The current timestep
8196: . x    - The coordinates
8197: . Nc   - The number of components
8198: . u    - The output field values
8199: - ctx  - optional user-defined function context

8201:   Level: developer

8203:   Developer Notes:
8204:   This API is specific to only particular usage of `DM`

8206:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8208: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabelLocal()`, `DMComputeL2Diff()`
8209: @*/
8210: 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)
8211: {
8212:   Vec localX;

8214:   PetscFunctionBegin;
8216:   PetscCall(DMGetLocalVector(dm, &localX));
8217:   PetscCall(VecSet(localX, 0.));
8218:   PetscCall(DMProjectFunctionLabelLocal(dm, time, label, numIds, ids, Nc, comps, funcs, ctxs, mode, localX));
8219:   PetscCall(DMLocalToGlobalBegin(dm, localX, mode, X));
8220:   PetscCall(DMLocalToGlobalEnd(dm, localX, mode, X));
8221:   PetscCall(DMRestoreLocalVector(dm, &localX));
8222:   PetscFunctionReturn(PETSC_SUCCESS);
8223: }

8225: /*@C
8226:   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.

8228:   Not Collective

8230:   Input Parameters:
8231: + dm     - The `DM`
8232: . time   - The time
8233: . label  - The `DMLabel` selecting the portion of the mesh for projection
8234: . numIds - The number of ids
8235: . ids    - The ids
8236: . Nc     - The number of components
8237: . comps  - The components
8238: . funcs  - The coordinate functions to evaluate, one per field
8239: . ctxs   - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
8240: - mode   - The insertion mode for values

8242:   Output Parameter:
8243: . localX - vector

8245:   Calling sequence of `funcs`:
8246: + dim  - The spatial dimension
8247: . time - The current time
8248: . x    - The coordinates
8249: . Nc   - The number of components
8250: . u    - The output field values
8251: - ctx  - optional user-defined function context

8253:   Level: developer

8255:   Developer Notes:
8256:   This API is specific to only particular usage of `DM`

8258:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8260: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
8261: @*/
8262: 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)
8263: {
8264:   PetscFunctionBegin;
8267:   PetscUseTypeMethod(dm, projectfunctionlabellocal, time, label, numIds, ids, Nc, comps, funcs, ctxs, mode, localX);
8268:   PetscFunctionReturn(PETSC_SUCCESS);
8269: }

8271: /*@C
8272:   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.

8274:   Not Collective

8276:   Input Parameters:
8277: + dm     - The `DM`
8278: . time   - The time
8279: . localU - The input field vector; may be `NULL` if projection is defined purely by coordinates
8280: . funcs  - The functions to evaluate, one per field
8281: - mode   - The insertion mode for values

8283:   Output Parameter:
8284: . localX - The output vector

8286:   Calling sequence of `funcs`:
8287: + dim          - The spatial dimension
8288: . Nf           - The number of input fields
8289: . NfAux        - The number of input auxiliary fields
8290: . uOff         - The offset of each field in u[]
8291: . uOff_x       - The offset of each field in u_x[]
8292: . u            - The field values at this point in space
8293: . u_t          - The field time derivative at this point in space (or NULL)
8294: . u_x          - The field derivatives at this point in space
8295: . aOff         - The offset of each auxiliary field in u[]
8296: . aOff_x       - The offset of each auxiliary field in u_x[]
8297: . a            - The auxiliary field values at this point in space
8298: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8299: . a_x          - The auxiliary field derivatives at this point in space
8300: . t            - The current time
8301: . x            - The coordinates of this point
8302: . numConstants - The number of constants
8303: . constants    - The value of each constant
8304: - f            - The value of the function at this point in space

8306:   Level: intermediate

8308:   Note:
8309:   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.
8310:   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
8311:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8312:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8314:   Developer Notes:
8315:   This API is specific to only particular usage of `DM`

8317:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8319: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabelLocal()`,
8320: `DMProjectFunction()`, `DMComputeL2Diff()`
8321: @*/
8322: 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)
8323: {
8324:   PetscFunctionBegin;
8328:   PetscUseTypeMethod(dm, projectfieldlocal, time, localU, funcs, mode, localX);
8329:   PetscFunctionReturn(PETSC_SUCCESS);
8330: }

8332: /*@C
8333:   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.

8335:   Not Collective

8337:   Input Parameters:
8338: + dm     - The `DM`
8339: . time   - The time
8340: . label  - The `DMLabel` marking the portion of the domain to output
8341: . numIds - The number of label ids to use
8342: . ids    - The label ids to use for marking
8343: . Nc     - The number of components to set in the output, or `PETSC_DETERMINE` for all components
8344: . comps  - The components to set in the output, or `NULL` for all components
8345: . localU - The input field vector
8346: . funcs  - The functions to evaluate, one per field
8347: - mode   - The insertion mode for values

8349:   Output Parameter:
8350: . localX - The output vector

8352:   Calling sequence of `funcs`:
8353: + dim          - The spatial dimension
8354: . Nf           - The number of input fields
8355: . NfAux        - The number of input auxiliary fields
8356: . uOff         - The offset of each field in u[]
8357: . uOff_x       - The offset of each field in u_x[]
8358: . u            - The field values at this point in space
8359: . u_t          - The field time derivative at this point in space (or NULL)
8360: . u_x          - The field derivatives at this point in space
8361: . aOff         - The offset of each auxiliary field in u[]
8362: . aOff_x       - The offset of each auxiliary field in u_x[]
8363: . a            - The auxiliary field values at this point in space
8364: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8365: . a_x          - The auxiliary field derivatives at this point in space
8366: . t            - The current time
8367: . x            - The coordinates of this point
8368: . numConstants - The number of constants
8369: . constants    - The value of each constant
8370: - f            - The value of the function at this point in space

8372:   Level: intermediate

8374:   Note:
8375:   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.
8376:   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
8377:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8378:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8380:   Developer Notes:
8381:   This API is specific to only particular usage of `DM`

8383:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8385: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabel()`, `DMProjectFunction()`, `DMComputeL2Diff()`
8386: @*/
8387: 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)
8388: {
8389:   PetscFunctionBegin;
8393:   PetscUseTypeMethod(dm, projectfieldlabellocal, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX);
8394:   PetscFunctionReturn(PETSC_SUCCESS);
8395: }

8397: /*@C
8398:   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.

8400:   Not Collective

8402:   Input Parameters:
8403: + dm     - The `DM`
8404: . time   - The time
8405: . label  - The `DMLabel` marking the portion of the domain to output
8406: . numIds - The number of label ids to use
8407: . ids    - The label ids to use for marking
8408: . Nc     - The number of components to set in the output, or `PETSC_DETERMINE` for all components
8409: . comps  - The components to set in the output, or `NULL` for all components
8410: . U      - The input field vector
8411: . funcs  - The functions to evaluate, one per field
8412: - mode   - The insertion mode for values

8414:   Output Parameter:
8415: . X - The output vector

8417:   Calling sequence of `funcs`:
8418: + dim          - The spatial dimension
8419: . Nf           - The number of input fields
8420: . NfAux        - The number of input auxiliary fields
8421: . uOff         - The offset of each field in u[]
8422: . uOff_x       - The offset of each field in u_x[]
8423: . u            - The field values at this point in space
8424: . u_t          - The field time derivative at this point in space (or NULL)
8425: . u_x          - The field derivatives at this point in space
8426: . aOff         - The offset of each auxiliary field in u[]
8427: . aOff_x       - The offset of each auxiliary field in u_x[]
8428: . a            - The auxiliary field values at this point in space
8429: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8430: . a_x          - The auxiliary field derivatives at this point in space
8431: . t            - The current time
8432: . x            - The coordinates of this point
8433: . numConstants - The number of constants
8434: . constants    - The value of each constant
8435: - f            - The value of the function at this point in space

8437:   Level: intermediate

8439:   Note:
8440:   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.
8441:   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
8442:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8443:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8445:   Developer Notes:
8446:   This API is specific to only particular usage of `DM`

8448:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8450: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabelLocal()`, `DMProjectFunction()`, `DMComputeL2Diff()`
8451: @*/
8452: 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)
8453: {
8454:   DM  dmIn;
8455:   Vec localU, localX;

8457:   PetscFunctionBegin;
8459:   PetscCall(VecGetDM(U, &dmIn));
8460:   PetscCall(DMGetLocalVector(dmIn, &localU));
8461:   PetscCall(DMGetLocalVector(dm, &localX));
8462:   PetscCall(VecSet(localX, 0.));
8463:   PetscCall(DMGlobalToLocalBegin(dmIn, U, mode, localU));
8464:   PetscCall(DMGlobalToLocalEnd(dmIn, U, mode, localU));
8465:   PetscCall(DMProjectFieldLabelLocal(dm, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX));
8466:   PetscCall(DMLocalToGlobalBegin(dm, localX, mode, X));
8467:   PetscCall(DMLocalToGlobalEnd(dm, localX, mode, X));
8468:   PetscCall(DMRestoreLocalVector(dm, &localX));
8469:   PetscCall(DMRestoreLocalVector(dmIn, &localU));
8470:   PetscFunctionReturn(PETSC_SUCCESS);
8471: }

8473: /*@C
8474:   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.

8476:   Not Collective

8478:   Input Parameters:
8479: + dm     - The `DM`
8480: . time   - The time
8481: . label  - The `DMLabel` marking the portion of the domain boundary to output
8482: . numIds - The number of label ids to use
8483: . ids    - The label ids to use for marking
8484: . Nc     - The number of components to set in the output, or `PETSC_DETERMINE` for all components
8485: . comps  - The components to set in the output, or `NULL` for all components
8486: . localU - The input field vector
8487: . funcs  - The functions to evaluate, one per field
8488: - mode   - The insertion mode for values

8490:   Output Parameter:
8491: . localX - The output vector

8493:   Calling sequence of `funcs`:
8494: + dim          - The spatial dimension
8495: . Nf           - The number of input fields
8496: . NfAux        - The number of input auxiliary fields
8497: . uOff         - The offset of each field in u[]
8498: . uOff_x       - The offset of each field in u_x[]
8499: . u            - The field values at this point in space
8500: . u_t          - The field time derivative at this point in space (or NULL)
8501: . u_x          - The field derivatives at this point in space
8502: . aOff         - The offset of each auxiliary field in u[]
8503: . aOff_x       - The offset of each auxiliary field in u_x[]
8504: . a            - The auxiliary field values at this point in space
8505: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8506: . a_x          - The auxiliary field derivatives at this point in space
8507: . t            - The current time
8508: . x            - The coordinates of this point
8509: . n            - The face normal
8510: . numConstants - The number of constants
8511: . constants    - The value of each constant
8512: - f            - The value of the function at this point in space

8514:   Level: intermediate

8516:   Note:
8517:   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.
8518:   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
8519:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8520:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8522:   Developer Notes:
8523:   This API is specific to only particular usage of `DM`

8525:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8527: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabelLocal()`, `DMProjectFunction()`, `DMComputeL2Diff()`
8528: @*/
8529: 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)
8530: {
8531:   PetscFunctionBegin;
8535:   PetscUseTypeMethod(dm, projectbdfieldlabellocal, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX);
8536:   PetscFunctionReturn(PETSC_SUCCESS);
8537: }

8539: /*@C
8540:   DMComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.

8542:   Collective

8544:   Input Parameters:
8545: + dm    - The `DM`
8546: . time  - The time
8547: . funcs - The functions to evaluate for each field component
8548: . ctxs  - Optional array of contexts to pass to each function, or NULL.
8549: - X     - The coefficient vector u_h, a global vector

8551:   Output Parameter:
8552: . diff - The diff ||u - u_h||_2

8554:   Level: developer

8556:   Developer Notes:
8557:   This API is specific to only particular usage of `DM`

8559:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8561: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMComputeL2FieldDiff()`, `DMComputeL2GradientDiff()`
8562: @*/
8563: PetscErrorCode DMComputeL2Diff(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
8564: {
8565:   PetscFunctionBegin;
8568:   PetscUseTypeMethod(dm, computel2diff, time, funcs, ctxs, X, diff);
8569:   PetscFunctionReturn(PETSC_SUCCESS);
8570: }

8572: /*@C
8573:   DMComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h.

8575:   Collective

8577:   Input Parameters:
8578: + dm    - The `DM`
8579: . time  - The time
8580: . funcs - The gradient functions to evaluate for each field component
8581: . ctxs  - Optional array of contexts to pass to each function, or NULL.
8582: . X     - The coefficient vector u_h, a global vector
8583: - n     - The vector to project along

8585:   Output Parameter:
8586: . diff - The diff ||(grad u - grad u_h) . n||_2

8588:   Level: developer

8590:   Developer Notes:
8591:   This API is specific to only particular usage of `DM`

8593:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8595: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMComputeL2Diff()`, `DMComputeL2FieldDiff()`
8596: @*/
8597: 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)
8598: {
8599:   PetscFunctionBegin;
8602:   PetscUseTypeMethod(dm, computel2gradientdiff, time, funcs, ctxs, X, n, diff);
8603:   PetscFunctionReturn(PETSC_SUCCESS);
8604: }

8606: /*@C
8607:   DMComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components.

8609:   Collective

8611:   Input Parameters:
8612: + dm    - The `DM`
8613: . time  - The time
8614: . funcs - The functions to evaluate for each field component
8615: . ctxs  - Optional array of contexts to pass to each function, or NULL.
8616: - X     - The coefficient vector u_h, a global vector

8618:   Output Parameter:
8619: . diff - The array of differences, ||u^f - u^f_h||_2

8621:   Level: developer

8623:   Developer Notes:
8624:   This API is specific to only particular usage of `DM`

8626:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8628: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMComputeL2GradientDiff()`
8629: @*/
8630: PetscErrorCode DMComputeL2FieldDiff(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
8631: {
8632:   PetscFunctionBegin;
8635:   PetscUseTypeMethod(dm, computel2fielddiff, time, funcs, ctxs, X, diff);
8636:   PetscFunctionReturn(PETSC_SUCCESS);
8637: }

8639: /*@C
8640:   DMGetNeighbors - Gets an array containing the MPI ranks of all the processes neighbors

8642:   Not Collective

8644:   Input Parameter:
8645: . dm - The `DM`

8647:   Output Parameters:
8648: + nranks - the number of neighbours
8649: - ranks  - the neighbors ranks

8651:   Level: beginner

8653:   Note:
8654:   Do not free the array, it is freed when the `DM` is destroyed.

8656: .seealso: [](ch_dmbase), `DM`, `DMDAGetNeighbors()`, `PetscSFGetRootRanks()`
8657: @*/
8658: PetscErrorCode DMGetNeighbors(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
8659: {
8660:   PetscFunctionBegin;
8662:   PetscUseTypeMethod(dm, getneighbors, nranks, ranks);
8663:   PetscFunctionReturn(PETSC_SUCCESS);
8664: }

8666: #include <petsc/private/matimpl.h>

8668: /*
8669:     Converts the input vector to a ghosted vector and then calls the standard coloring code.
8670:     This must be a different function because it requires DM which is not defined in the Mat library
8671: */
8672: static PetscErrorCode MatFDColoringApply_AIJDM(Mat J, MatFDColoring coloring, Vec x1, void *sctx)
8673: {
8674:   PetscFunctionBegin;
8675:   if (coloring->ctype == IS_COLORING_LOCAL) {
8676:     Vec x1local;
8677:     DM  dm;
8678:     PetscCall(MatGetDM(J, &dm));
8679:     PetscCheck(dm, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_INCOMP, "IS_COLORING_LOCAL requires a DM");
8680:     PetscCall(DMGetLocalVector(dm, &x1local));
8681:     PetscCall(DMGlobalToLocalBegin(dm, x1, INSERT_VALUES, x1local));
8682:     PetscCall(DMGlobalToLocalEnd(dm, x1, INSERT_VALUES, x1local));
8683:     x1 = x1local;
8684:   }
8685:   PetscCall(MatFDColoringApply_AIJ(J, coloring, x1, sctx));
8686:   if (coloring->ctype == IS_COLORING_LOCAL) {
8687:     DM dm;
8688:     PetscCall(MatGetDM(J, &dm));
8689:     PetscCall(DMRestoreLocalVector(dm, &x1));
8690:   }
8691:   PetscFunctionReturn(PETSC_SUCCESS);
8692: }

8694: /*@
8695:   MatFDColoringUseDM - allows a `MatFDColoring` object to use the `DM` associated with the matrix to compute a `IS_COLORING_LOCAL` coloring

8697:   Input Parameters:
8698: + coloring   - The matrix to get the `DM` from
8699: - fdcoloring - the `MatFDColoring` object

8701:   Level: advanced

8703:   Developer Note:
8704:   This routine exists because the PETSc `Mat` library does not know about the `DM` objects

8706: .seealso: [](ch_dmbase), `DM`, `MatFDColoring`, `MatFDColoringCreate()`, `ISColoringType`
8707: @*/
8708: PetscErrorCode MatFDColoringUseDM(Mat coloring, MatFDColoring fdcoloring)
8709: {
8710:   PetscFunctionBegin;
8711:   coloring->ops->fdcoloringapply = MatFDColoringApply_AIJDM;
8712:   PetscFunctionReturn(PETSC_SUCCESS);
8713: }

8715: /*@
8716:   DMGetCompatibility - determine if two `DM`s are compatible

8718:   Collective

8720:   Input Parameters:
8721: + dm1 - the first `DM`
8722: - dm2 - the second `DM`

8724:   Output Parameters:
8725: + compatible - whether or not the two `DM`s are compatible
8726: - set        - whether or not the compatible value was actually determined and set

8728:   Level: advanced

8730:   Notes:
8731:   Two `DM`s are deemed compatible if they represent the same parallel decomposition
8732:   of the same topology. This implies that the section (field data) on one
8733:   "makes sense" with respect to the topology and parallel decomposition of the other.
8734:   Loosely speaking, compatible `DM`s represent the same domain and parallel
8735:   decomposition, but hold different data.

8737:   Typically, one would confirm compatibility if intending to simultaneously iterate
8738:   over a pair of vectors obtained from different `DM`s.

8740:   For example, two `DMDA` objects are compatible if they have the same local
8741:   and global sizes and the same stencil width. They can have different numbers
8742:   of degrees of freedom per node. Thus, one could use the node numbering from
8743:   either `DM` in bounds for a loop over vectors derived from either `DM`.

8745:   Consider the operation of summing data living on a 2-dof `DMDA` to data living
8746:   on a 1-dof `DMDA`, which should be compatible, as in the following snippet.
8747: .vb
8748:   ...
8749:   PetscCall(DMGetCompatibility(da1,da2,&compatible,&set));
8750:   if (set && compatible)  {
8751:     PetscCall(DMDAVecGetArrayDOF(da1,vec1,&arr1));
8752:     PetscCall(DMDAVecGetArrayDOF(da2,vec2,&arr2));
8753:     PetscCall(DMDAGetCorners(da1,&x,&y,NULL,&m,&n,NULL));
8754:     for (j=y; j<y+n; ++j) {
8755:       for (i=x; i<x+m, ++i) {
8756:         arr1[j][i][0] = arr2[j][i][0] + arr2[j][i][1];
8757:       }
8758:     }
8759:     PetscCall(DMDAVecRestoreArrayDOF(da1,vec1,&arr1));
8760:     PetscCall(DMDAVecRestoreArrayDOF(da2,vec2,&arr2));
8761:   } else {
8762:     SETERRQ(PetscObjectComm((PetscObject)da1,PETSC_ERR_ARG_INCOMP,"DMDA objects incompatible");
8763:   }
8764:   ...
8765: .ve

8767:   Checking compatibility might be expensive for a given implementation of `DM`,
8768:   or might be impossible to unambiguously confirm or deny. For this reason,
8769:   this function may decline to determine compatibility, and hence users should
8770:   always check the "set" output parameter.

8772:   A `DM` is always compatible with itself.

8774:   In the current implementation, `DM`s which live on "unequal" communicators
8775:   (MPI_UNEQUAL in the terminology of MPI_Comm_compare()) are always deemed
8776:   incompatible.

8778:   This function is labeled "Collective," as information about all subdomains
8779:   is required on each rank. However, in `DM` implementations which store all this
8780:   information locally, this function may be merely "Logically Collective".

8782:   Developer Note:
8783:   Compatibility is assumed to be a symmetric concept; `DM` A is compatible with `DM` B
8784:   iff B is compatible with A. Thus, this function checks the implementations
8785:   of both dm and dmc (if they are of different types), attempting to determine
8786:   compatibility. It is left to `DM` implementers to ensure that symmetry is
8787:   preserved. The simplest way to do this is, when implementing type-specific
8788:   logic for this function, is to check for existing logic in the implementation
8789:   of other `DM` types and let *set = PETSC_FALSE if found.

8791: .seealso: [](ch_dmbase), `DM`, `DMDACreateCompatibleDMDA()`, `DMStagCreateCompatibleDMStag()`
8792: @*/
8793: PetscErrorCode DMGetCompatibility(DM dm1, DM dm2, PetscBool *compatible, PetscBool *set)
8794: {
8795:   PetscMPIInt compareResult;
8796:   DMType      type, type2;
8797:   PetscBool   sameType;

8799:   PetscFunctionBegin;

8803:   /* Declare a DM compatible with itself */
8804:   if (dm1 == dm2) {
8805:     *set        = PETSC_TRUE;
8806:     *compatible = PETSC_TRUE;
8807:     PetscFunctionReturn(PETSC_SUCCESS);
8808:   }

8810:   /* Declare a DM incompatible with a DM that lives on an "unequal"
8811:      communicator. Note that this does not preclude compatibility with
8812:      DMs living on "congruent" or "similar" communicators, but this must be
8813:      determined by the implementation-specific logic */
8814:   PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dm1), PetscObjectComm((PetscObject)dm2), &compareResult));
8815:   if (compareResult == MPI_UNEQUAL) {
8816:     *set        = PETSC_TRUE;
8817:     *compatible = PETSC_FALSE;
8818:     PetscFunctionReturn(PETSC_SUCCESS);
8819:   }

8821:   /* Pass to the implementation-specific routine, if one exists. */
8822:   if (dm1->ops->getcompatibility) {
8823:     PetscUseTypeMethod(dm1, getcompatibility, dm2, compatible, set);
8824:     if (*set) PetscFunctionReturn(PETSC_SUCCESS);
8825:   }

8827:   /* If dm1 and dm2 are of different types, then attempt to check compatibility
8828:      with an implementation of this function from dm2 */
8829:   PetscCall(DMGetType(dm1, &type));
8830:   PetscCall(DMGetType(dm2, &type2));
8831:   PetscCall(PetscStrcmp(type, type2, &sameType));
8832:   if (!sameType && dm2->ops->getcompatibility) {
8833:     PetscUseTypeMethod(dm2, getcompatibility, dm1, compatible, set); /* Note argument order */
8834:   } else {
8835:     *set = PETSC_FALSE;
8836:   }
8837:   PetscFunctionReturn(PETSC_SUCCESS);
8838: }

8840: /*@C
8841:   DMMonitorSet - Sets an additional monitor function that is to be used after a solve to monitor discretization performance.

8843:   Logically Collective

8845:   Input Parameters:
8846: + dm             - the `DM`
8847: . f              - the monitor function
8848: . mctx           - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired)
8849: - monitordestroy - [optional] routine that frees monitor context (may be `NULL`), see `PetscCtxDestroyFn` for the calling sequence

8851:   Options Database Key:
8852: . -dm_monitor_cancel - cancels all monitors that have been hardwired into a code by calls to `DMMonitorSet()`, but
8853:                        does not cancel those set via the options database.

8855:   Level: intermediate

8857:   Note:
8858:   Several different monitoring routines may be set by calling
8859:   `DMMonitorSet()` multiple times or with `DMMonitorSetFromOptions()`; all will be called in the
8860:   order in which they were set.

8862:   Fortran Note:
8863:   Only a single monitor function can be set for each `DM` object

8865:   Developer Note:
8866:   This API has a generic name but seems specific to a very particular aspect of the use of `DM`

8868: .seealso: [](ch_dmbase), `DM`, `DMMonitorCancel()`, `DMMonitorSetFromOptions()`, `DMMonitor()`, `PetscCtxDestroyFn`
8869: @*/
8870: PetscErrorCode DMMonitorSet(DM dm, PetscErrorCode (*f)(DM, void *), void *mctx, PetscCtxDestroyFn *monitordestroy)
8871: {
8872:   PetscInt m;

8874:   PetscFunctionBegin;
8876:   for (m = 0; m < dm->numbermonitors; ++m) {
8877:     PetscBool identical;

8879:     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))f, mctx, monitordestroy, (PetscErrorCode (*)(void))dm->monitor[m], dm->monitorcontext[m], dm->monitordestroy[m], &identical));
8880:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
8881:   }
8882:   PetscCheck(dm->numbermonitors < MAXDMMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
8883:   dm->monitor[dm->numbermonitors]          = f;
8884:   dm->monitordestroy[dm->numbermonitors]   = monitordestroy;
8885:   dm->monitorcontext[dm->numbermonitors++] = mctx;
8886:   PetscFunctionReturn(PETSC_SUCCESS);
8887: }

8889: /*@
8890:   DMMonitorCancel - Clears all the monitor functions for a `DM` object.

8892:   Logically Collective

8894:   Input Parameter:
8895: . dm - the DM

8897:   Options Database Key:
8898: . -dm_monitor_cancel - cancels all monitors that have been hardwired
8899:   into a code by calls to `DMonitorSet()`, but does not cancel those
8900:   set via the options database

8902:   Level: intermediate

8904:   Note:
8905:   There is no way to clear one specific monitor from a `DM` object.

8907: .seealso: [](ch_dmbase), `DM`, `DMMonitorSet()`, `DMMonitorSetFromOptions()`, `DMMonitor()`
8908: @*/
8909: PetscErrorCode DMMonitorCancel(DM dm)
8910: {
8911:   PetscInt m;

8913:   PetscFunctionBegin;
8915:   for (m = 0; m < dm->numbermonitors; ++m) {
8916:     if (dm->monitordestroy[m]) PetscCall((*dm->monitordestroy[m])(&dm->monitorcontext[m]));
8917:   }
8918:   dm->numbermonitors = 0;
8919:   PetscFunctionReturn(PETSC_SUCCESS);
8920: }

8922: /*@C
8923:   DMMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user

8925:   Collective

8927:   Input Parameters:
8928: + dm           - `DM` object you wish to monitor
8929: . name         - the monitor type one is seeking
8930: . help         - message indicating what monitoring is done
8931: . manual       - manual page for the monitor
8932: . monitor      - the monitor function, this must use a `PetscViewerFormat` as its context
8933: - 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

8935:   Output Parameter:
8936: . flg - Flag set if the monitor was created

8938:   Level: developer

8940: .seealso: [](ch_dmbase), `DM`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
8941:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
8942:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
8943:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
8944:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
8945:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
8946:           `PetscOptionsFList()`, `PetscOptionsEList()`, `DMMonitor()`, `DMMonitorSet()`
8947: @*/
8948: PetscErrorCode DMMonitorSetFromOptions(DM dm, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(DM, void *), PetscErrorCode (*monitorsetup)(DM, PetscViewerAndFormat *), PetscBool *flg)
8949: {
8950:   PetscViewer       viewer;
8951:   PetscViewerFormat format;

8953:   PetscFunctionBegin;
8955:   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)dm), ((PetscObject)dm)->options, ((PetscObject)dm)->prefix, name, &viewer, &format, flg));
8956:   if (*flg) {
8957:     PetscViewerAndFormat *vf;

8959:     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
8960:     PetscCall(PetscViewerDestroy(&viewer));
8961:     if (monitorsetup) PetscCall((*monitorsetup)(dm, vf));
8962:     PetscCall(DMMonitorSet(dm, monitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy));
8963:   }
8964:   PetscFunctionReturn(PETSC_SUCCESS);
8965: }

8967: /*@
8968:   DMMonitor - runs the user provided monitor routines, if they exist

8970:   Collective

8972:   Input Parameter:
8973: . dm - The `DM`

8975:   Level: developer

8977:   Developer Note:
8978:   Note should indicate when during the life of the `DM` the monitor is run. It appears to be
8979:   related to the discretization process seems rather specialized since some `DM` have no
8980:   concept of discretization.

8982: .seealso: [](ch_dmbase), `DM`, `DMMonitorSet()`, `DMMonitorSetFromOptions()`
8983: @*/
8984: PetscErrorCode DMMonitor(DM dm)
8985: {
8986:   PetscInt m;

8988:   PetscFunctionBegin;
8989:   if (!dm) PetscFunctionReturn(PETSC_SUCCESS);
8991:   for (m = 0; m < dm->numbermonitors; ++m) PetscCall((*dm->monitor[m])(dm, dm->monitorcontext[m]));
8992:   PetscFunctionReturn(PETSC_SUCCESS);
8993: }

8995: /*@
8996:   DMComputeError - Computes the error assuming the user has provided the exact solution functions

8998:   Collective

9000:   Input Parameters:
9001: + dm  - The `DM`
9002: - sol - The solution vector

9004:   Input/Output Parameter:
9005: . errors - An array of length Nf, the number of fields, or `NULL` for no output; on output
9006:            contains the error in each field

9008:   Output Parameter:
9009: . errorVec - A vector to hold the cellwise error (may be `NULL`)

9011:   Level: developer

9013:   Note:
9014:   The exact solutions come from the `PetscDS` object, and the time comes from `DMGetOutputSequenceNumber()`.

9016: .seealso: [](ch_dmbase), `DM`, `DMMonitorSet()`, `DMGetRegionNumDS()`, `PetscDSGetExactSolution()`, `DMGetOutputSequenceNumber()`
9017: @*/
9018: PetscErrorCode DMComputeError(DM dm, Vec sol, PetscReal errors[], Vec *errorVec)
9019: {
9020:   PetscErrorCode (**exactSol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
9021:   void    **ctxs;
9022:   PetscReal time;
9023:   PetscInt  Nf, f, Nds, s;

9025:   PetscFunctionBegin;
9026:   PetscCall(DMGetNumFields(dm, &Nf));
9027:   PetscCall(PetscCalloc2(Nf, &exactSol, Nf, &ctxs));
9028:   PetscCall(DMGetNumDS(dm, &Nds));
9029:   for (s = 0; s < Nds; ++s) {
9030:     PetscDS         ds;
9031:     DMLabel         label;
9032:     IS              fieldIS;
9033:     const PetscInt *fields;
9034:     PetscInt        dsNf;

9036:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL));
9037:     PetscCall(PetscDSGetNumFields(ds, &dsNf));
9038:     if (fieldIS) PetscCall(ISGetIndices(fieldIS, &fields));
9039:     for (f = 0; f < dsNf; ++f) {
9040:       const PetscInt field = fields[f];
9041:       PetscCall(PetscDSGetExactSolution(ds, field, &exactSol[field], &ctxs[field]));
9042:     }
9043:     if (fieldIS) PetscCall(ISRestoreIndices(fieldIS, &fields));
9044:   }
9045:   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);
9046:   PetscCall(DMGetOutputSequenceNumber(dm, NULL, &time));
9047:   if (errors) PetscCall(DMComputeL2FieldDiff(dm, time, exactSol, ctxs, sol, errors));
9048:   if (errorVec) {
9049:     DM             edm;
9050:     DMPolytopeType ct;
9051:     PetscBool      simplex;
9052:     PetscInt       dim, cStart, Nf;

9054:     PetscCall(DMClone(dm, &edm));
9055:     PetscCall(DMGetDimension(edm, &dim));
9056:     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
9057:     PetscCall(DMPlexGetCellType(dm, cStart, &ct));
9058:     simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
9059:     PetscCall(DMGetNumFields(dm, &Nf));
9060:     for (f = 0; f < Nf; ++f) {
9061:       PetscFE         fe, efe;
9062:       PetscQuadrature q;
9063:       const char     *name;

9065:       PetscCall(DMGetField(dm, f, NULL, (PetscObject *)&fe));
9066:       PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, Nf, simplex, 0, PETSC_DETERMINE, &efe));
9067:       PetscCall(PetscObjectGetName((PetscObject)fe, &name));
9068:       PetscCall(PetscObjectSetName((PetscObject)efe, name));
9069:       PetscCall(PetscFEGetQuadrature(fe, &q));
9070:       PetscCall(PetscFESetQuadrature(efe, q));
9071:       PetscCall(DMSetField(edm, f, NULL, (PetscObject)efe));
9072:       PetscCall(PetscFEDestroy(&efe));
9073:     }
9074:     PetscCall(DMCreateDS(edm));

9076:     PetscCall(DMCreateGlobalVector(edm, errorVec));
9077:     PetscCall(PetscObjectSetName((PetscObject)*errorVec, "Error"));
9078:     PetscCall(DMPlexComputeL2DiffVec(dm, time, exactSol, ctxs, sol, *errorVec));
9079:     PetscCall(DMDestroy(&edm));
9080:   }
9081:   PetscCall(PetscFree2(exactSol, ctxs));
9082:   PetscFunctionReturn(PETSC_SUCCESS);
9083: }

9085: /*@
9086:   DMGetNumAuxiliaryVec - Get the number of auxiliary vectors associated with this `DM`

9088:   Not Collective

9090:   Input Parameter:
9091: . dm - The `DM`

9093:   Output Parameter:
9094: . numAux - The number of auxiliary data vectors

9096:   Level: advanced

9098: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMSetAuxiliaryVec()`, `DMGetAuxiliaryLabels()`, `DMGetAuxiliaryVec()`
9099: @*/
9100: PetscErrorCode DMGetNumAuxiliaryVec(DM dm, PetscInt *numAux)
9101: {
9102:   PetscFunctionBegin;
9104:   PetscCall(PetscHMapAuxGetSize(dm->auxData, numAux));
9105:   PetscFunctionReturn(PETSC_SUCCESS);
9106: }

9108: /*@
9109:   DMGetAuxiliaryVec - Get the auxiliary vector for region specified by the given label and value, and equation part

9111:   Not Collective

9113:   Input Parameters:
9114: + dm    - The `DM`
9115: . label - The `DMLabel`
9116: . value - The label value indicating the region
9117: - part  - The equation part, or 0 if unused

9119:   Output Parameter:
9120: . aux - The `Vec` holding auxiliary field data

9122:   Level: advanced

9124:   Note:
9125:   If no auxiliary vector is found for this (label, value), (NULL, 0, 0) is checked as well.

9127: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMSetAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryLabels()`
9128: @*/
9129: PetscErrorCode DMGetAuxiliaryVec(DM dm, DMLabel label, PetscInt value, PetscInt part, Vec *aux)
9130: {
9131:   PetscHashAuxKey key, wild = {NULL, 0, 0};
9132:   PetscBool       has;

9134:   PetscFunctionBegin;
9137:   key.label = label;
9138:   key.value = value;
9139:   key.part  = part;
9140:   PetscCall(PetscHMapAuxHas(dm->auxData, key, &has));
9141:   if (has) PetscCall(PetscHMapAuxGet(dm->auxData, key, aux));
9142:   else PetscCall(PetscHMapAuxGet(dm->auxData, wild, aux));
9143:   PetscFunctionReturn(PETSC_SUCCESS);
9144: }

9146: /*@
9147:   DMSetAuxiliaryVec - Set an auxiliary vector for region specified by the given label and value, and equation part

9149:   Not Collective because auxiliary vectors are not parallel

9151:   Input Parameters:
9152: + dm    - The `DM`
9153: . label - The `DMLabel`
9154: . value - The label value indicating the region
9155: . part  - The equation part, or 0 if unused
9156: - aux   - The `Vec` holding auxiliary field data

9158:   Level: advanced

9160: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMGetAuxiliaryLabels()`, `DMCopyAuxiliaryVec()`
9161: @*/
9162: PetscErrorCode DMSetAuxiliaryVec(DM dm, DMLabel label, PetscInt value, PetscInt part, Vec aux)
9163: {
9164:   Vec             old;
9165:   PetscHashAuxKey key;

9167:   PetscFunctionBegin;
9170:   key.label = label;
9171:   key.value = value;
9172:   key.part  = part;
9173:   PetscCall(PetscHMapAuxGet(dm->auxData, key, &old));
9174:   PetscCall(PetscObjectReference((PetscObject)aux));
9175:   if (!aux) PetscCall(PetscHMapAuxDel(dm->auxData, key));
9176:   else PetscCall(PetscHMapAuxSet(dm->auxData, key, aux));
9177:   PetscCall(VecDestroy(&old));
9178:   PetscFunctionReturn(PETSC_SUCCESS);
9179: }

9181: /*@
9182:   DMGetAuxiliaryLabels - Get the labels, values, and parts for all auxiliary vectors in this `DM`

9184:   Not Collective

9186:   Input Parameter:
9187: . dm - The `DM`

9189:   Output Parameters:
9190: + labels - The `DMLabel`s for each `Vec`
9191: . values - The label values for each `Vec`
9192: - parts  - The equation parts for each `Vec`

9194:   Level: advanced

9196:   Note:
9197:   The arrays passed in must be at least as large as `DMGetNumAuxiliaryVec()`.

9199: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMSetAuxiliaryVec()`, `DMCopyAuxiliaryVec()`
9200: @*/
9201: PetscErrorCode DMGetAuxiliaryLabels(DM dm, DMLabel labels[], PetscInt values[], PetscInt parts[])
9202: {
9203:   PetscHashAuxKey *keys;
9204:   PetscInt         n, i, off = 0;

9206:   PetscFunctionBegin;
9208:   PetscAssertPointer(labels, 2);
9209:   PetscAssertPointer(values, 3);
9210:   PetscAssertPointer(parts, 4);
9211:   PetscCall(DMGetNumAuxiliaryVec(dm, &n));
9212:   PetscCall(PetscMalloc1(n, &keys));
9213:   PetscCall(PetscHMapAuxGetKeys(dm->auxData, &off, keys));
9214:   for (i = 0; i < n; ++i) {
9215:     labels[i] = keys[i].label;
9216:     values[i] = keys[i].value;
9217:     parts[i]  = keys[i].part;
9218:   }
9219:   PetscCall(PetscFree(keys));
9220:   PetscFunctionReturn(PETSC_SUCCESS);
9221: }

9223: /*@
9224:   DMCopyAuxiliaryVec - Copy the auxiliary vector data on a `DM` to a new `DM`

9226:   Not Collective

9228:   Input Parameter:
9229: . dm - The `DM`

9231:   Output Parameter:
9232: . dmNew - The new `DM`, now with the same auxiliary data

9234:   Level: advanced

9236:   Note:
9237:   This is a shallow copy of the auxiliary vectors

9239: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMSetAuxiliaryVec()`
9240: @*/
9241: PetscErrorCode DMCopyAuxiliaryVec(DM dm, DM dmNew)
9242: {
9243:   PetscFunctionBegin;
9246:   if (dm == dmNew) PetscFunctionReturn(PETSC_SUCCESS);
9247:   PetscCall(DMClearAuxiliaryVec(dmNew));

9249:   PetscCall(PetscHMapAuxDestroy(&dmNew->auxData));
9250:   PetscCall(PetscHMapAuxDuplicate(dm->auxData, &dmNew->auxData));
9251:   {
9252:     Vec     *auxData;
9253:     PetscInt n, i, off = 0;

9255:     PetscCall(PetscHMapAuxGetSize(dmNew->auxData, &n));
9256:     PetscCall(PetscMalloc1(n, &auxData));
9257:     PetscCall(PetscHMapAuxGetVals(dmNew->auxData, &off, auxData));
9258:     for (i = 0; i < n; ++i) PetscCall(PetscObjectReference((PetscObject)auxData[i]));
9259:     PetscCall(PetscFree(auxData));
9260:   }
9261:   PetscFunctionReturn(PETSC_SUCCESS);
9262: }

9264: /*@
9265:   DMClearAuxiliaryVec - Destroys the auxiliary vector information and creates a new empty one

9267:   Not Collective

9269:   Input Parameter:
9270: . dm - The `DM`

9272:   Level: advanced

9274: .seealso: [](ch_dmbase), `DM`, `DMCopyAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMSetAuxiliaryVec()`
9275: @*/
9276: PetscErrorCode DMClearAuxiliaryVec(DM dm)
9277: {
9278:   Vec     *auxData;
9279:   PetscInt n, i, off = 0;

9281:   PetscFunctionBegin;
9282:   PetscCall(PetscHMapAuxGetSize(dm->auxData, &n));
9283:   PetscCall(PetscMalloc1(n, &auxData));
9284:   PetscCall(PetscHMapAuxGetVals(dm->auxData, &off, auxData));
9285:   for (i = 0; i < n; ++i) PetscCall(VecDestroy(&auxData[i]));
9286:   PetscCall(PetscFree(auxData));
9287:   PetscCall(PetscHMapAuxDestroy(&dm->auxData));
9288:   PetscCall(PetscHMapAuxCreate(&dm->auxData));
9289:   PetscFunctionReturn(PETSC_SUCCESS);
9290: }

9292: /*@
9293:   DMPolytopeMatchOrientation - Determine an orientation (transformation) that takes the source face arrangement to the target face arrangement

9295:   Not Collective

9297:   Input Parameters:
9298: + ct         - The `DMPolytopeType`
9299: . sourceCone - The source arrangement of faces
9300: - targetCone - The target arrangement of faces

9302:   Output Parameters:
9303: + ornt  - The orientation (transformation) which will take the source arrangement to the target arrangement
9304: - found - Flag indicating that a suitable orientation was found

9306:   Level: advanced

9308:   Note:
9309:   An arrangement is a face order combined with an orientation for each face

9311:   Each orientation (transformation) is labeled with an integer from negative `DMPolytopeTypeGetNumArrangements(ct)`/2 to `DMPolytopeTypeGetNumArrangements(ct)`/2
9312:   that labels each arrangement (face ordering plus orientation for each face).

9314:   See `DMPolytopeMatchVertexOrientation()` to find a new vertex orientation that takes the source vertex arrangement to the target vertex arrangement

9316: .seealso: [](ch_dmbase), `DM`, `DMPolytopeGetOrientation()`, `DMPolytopeMatchVertexOrientation()`, `DMPolytopeGetVertexOrientation()`
9317: @*/
9318: PetscErrorCode DMPolytopeMatchOrientation(DMPolytopeType ct, const PetscInt sourceCone[], const PetscInt targetCone[], PetscInt *ornt, PetscBool *found)
9319: {
9320:   const PetscInt cS = DMPolytopeTypeGetConeSize(ct);
9321:   const PetscInt nO = DMPolytopeTypeGetNumArrangements(ct) / 2;
9322:   PetscInt       o, c;

9324:   PetscFunctionBegin;
9325:   if (!nO) {
9326:     *ornt  = 0;
9327:     *found = PETSC_TRUE;
9328:     PetscFunctionReturn(PETSC_SUCCESS);
9329:   }
9330:   for (o = -nO; o < nO; ++o) {
9331:     const PetscInt *arr = DMPolytopeTypeGetArrangement(ct, o);

9333:     for (c = 0; c < cS; ++c)
9334:       if (sourceCone[arr[c * 2]] != targetCone[c]) break;
9335:     if (c == cS) {
9336:       *ornt = o;
9337:       break;
9338:     }
9339:   }
9340:   *found = o == nO ? PETSC_FALSE : PETSC_TRUE;
9341:   PetscFunctionReturn(PETSC_SUCCESS);
9342: }

9344: /*@
9345:   DMPolytopeGetOrientation - Determine an orientation (transformation) that takes the source face arrangement to the target face arrangement

9347:   Not Collective

9349:   Input Parameters:
9350: + ct         - The `DMPolytopeType`
9351: . sourceCone - The source arrangement of faces
9352: - targetCone - The target arrangement of faces

9354:   Output Parameter:
9355: . ornt - The orientation (transformation) which will take the source arrangement to the target arrangement

9357:   Level: advanced

9359:   Note:
9360:   This function is the same as `DMPolytopeMatchOrientation()` except it will generate an error if no suitable orientation can be found.

9362:   Developer Note:
9363:   It is unclear why this function needs to exist since one can simply call `DMPolytopeMatchOrientation()` and error if none is found

9365: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMPolytopeMatchOrientation()`, `DMPolytopeGetVertexOrientation()`, `DMPolytopeMatchVertexOrientation()`
9366: @*/
9367: PetscErrorCode DMPolytopeGetOrientation(DMPolytopeType ct, const PetscInt sourceCone[], const PetscInt targetCone[], PetscInt *ornt)
9368: {
9369:   PetscBool found;

9371:   PetscFunctionBegin;
9372:   PetscCall(DMPolytopeMatchOrientation(ct, sourceCone, targetCone, ornt, &found));
9373:   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find orientation for %s", DMPolytopeTypes[ct]);
9374:   PetscFunctionReturn(PETSC_SUCCESS);
9375: }

9377: /*@
9378:   DMPolytopeMatchVertexOrientation - Determine an orientation (transformation) that takes the source vertex arrangement to the target vertex arrangement

9380:   Not Collective

9382:   Input Parameters:
9383: + ct         - The `DMPolytopeType`
9384: . sourceVert - The source arrangement of vertices
9385: - targetVert - The target arrangement of vertices

9387:   Output Parameters:
9388: + ornt  - The orientation (transformation) which will take the source arrangement to the target arrangement
9389: - found - Flag indicating that a suitable orientation was found

9391:   Level: advanced

9393:   Notes:
9394:   An arrangement is a vertex order

9396:   Each orientation (transformation) is labeled with an integer from negative `DMPolytopeTypeGetNumArrangements(ct)`/2 to `DMPolytopeTypeGetNumArrangements(ct)`/2
9397:   that labels each arrangement (vertex ordering).

9399:   See `DMPolytopeMatchOrientation()` to find a new face orientation that takes the source face arrangement to the target face arrangement

9401: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMPolytopeGetOrientation()`, `DMPolytopeMatchOrientation()`, `DMPolytopeTypeGetNumVertices()`, `DMPolytopeTypeGetVertexArrangement()`
9402: @*/
9403: PetscErrorCode DMPolytopeMatchVertexOrientation(DMPolytopeType ct, const PetscInt sourceVert[], const PetscInt targetVert[], PetscInt *ornt, PetscBool *found)
9404: {
9405:   const PetscInt cS = DMPolytopeTypeGetNumVertices(ct);
9406:   const PetscInt nO = DMPolytopeTypeGetNumArrangements(ct) / 2;
9407:   PetscInt       o, c;

9409:   PetscFunctionBegin;
9410:   if (!nO) {
9411:     *ornt  = 0;
9412:     *found = PETSC_TRUE;
9413:     PetscFunctionReturn(PETSC_SUCCESS);
9414:   }
9415:   for (o = -nO; o < nO; ++o) {
9416:     const PetscInt *arr = DMPolytopeTypeGetVertexArrangement(ct, o);

9418:     for (c = 0; c < cS; ++c)
9419:       if (sourceVert[arr[c]] != targetVert[c]) break;
9420:     if (c == cS) {
9421:       *ornt = o;
9422:       break;
9423:     }
9424:   }
9425:   *found = o == nO ? PETSC_FALSE : PETSC_TRUE;
9426:   PetscFunctionReturn(PETSC_SUCCESS);
9427: }

9429: /*@
9430:   DMPolytopeGetVertexOrientation - Determine an orientation (transformation) that takes the source vertex arrangement to the target vertex arrangement

9432:   Not Collective

9434:   Input Parameters:
9435: + ct         - The `DMPolytopeType`
9436: . sourceCone - The source arrangement of vertices
9437: - targetCone - The target arrangement of vertices

9439:   Output Parameter:
9440: . ornt - The orientation (transformation) which will take the source arrangement to the target arrangement

9442:   Level: advanced

9444:   Note:
9445:   This function is the same as `DMPolytopeMatchVertexOrientation()` except it errors if not orientation is possible.

9447:   Developer Note:
9448:   It is unclear why this function needs to exist since one can simply call `DMPolytopeMatchVertexOrientation()` and error if none is found

9450: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMPolytopeMatchVertexOrientation()`, `DMPolytopeGetOrientation()`
9451: @*/
9452: PetscErrorCode DMPolytopeGetVertexOrientation(DMPolytopeType ct, const PetscInt sourceCone[], const PetscInt targetCone[], PetscInt *ornt)
9453: {
9454:   PetscBool found;

9456:   PetscFunctionBegin;
9457:   PetscCall(DMPolytopeMatchVertexOrientation(ct, sourceCone, targetCone, ornt, &found));
9458:   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find orientation for %s", DMPolytopeTypes[ct]);
9459:   PetscFunctionReturn(PETSC_SUCCESS);
9460: }

9462: /*@
9463:   DMPolytopeInCellTest - Check whether a point lies inside the reference cell of given type

9465:   Not Collective

9467:   Input Parameters:
9468: + ct    - The `DMPolytopeType`
9469: - point - Coordinates of the point

9471:   Output Parameter:
9472: . inside - Flag indicating whether the point is inside the reference cell of given type

9474:   Level: advanced

9476: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMLocatePoints()`
9477: @*/
9478: PetscErrorCode DMPolytopeInCellTest(DMPolytopeType ct, const PetscReal point[], PetscBool *inside)
9479: {
9480:   PetscReal sum = 0.0;
9481:   PetscInt  d;

9483:   PetscFunctionBegin;
9484:   *inside = PETSC_TRUE;
9485:   switch (ct) {
9486:   case DM_POLYTOPE_TRIANGLE:
9487:   case DM_POLYTOPE_TETRAHEDRON:
9488:     for (d = 0; d < DMPolytopeTypeGetDim(ct); ++d) {
9489:       if (point[d] < -1.0) {
9490:         *inside = PETSC_FALSE;
9491:         break;
9492:       }
9493:       sum += point[d];
9494:     }
9495:     if (sum > PETSC_SMALL) {
9496:       *inside = PETSC_FALSE;
9497:       break;
9498:     }
9499:     break;
9500:   case DM_POLYTOPE_QUADRILATERAL:
9501:   case DM_POLYTOPE_HEXAHEDRON:
9502:     for (d = 0; d < DMPolytopeTypeGetDim(ct); ++d)
9503:       if (PetscAbsReal(point[d]) > 1. + PETSC_SMALL) {
9504:         *inside = PETSC_FALSE;
9505:         break;
9506:       }
9507:     break;
9508:   default:
9509:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported polytope type %s", DMPolytopeTypes[ct]);
9510:   }
9511:   PetscFunctionReturn(PETSC_SUCCESS);
9512: }

9514: /*@
9515:   DMReorderSectionSetDefault - Set flag indicating whether the local section should be reordered by default

9517:   Logically collective

9519:   Input Parameters:
9520: + dm      - The DM
9521: - reorder - Flag for reordering

9523:   Level: intermediate

9525: .seealso: `DMReorderSectionGetDefault()`
9526: @*/
9527: PetscErrorCode DMReorderSectionSetDefault(DM dm, DMReorderDefaultFlag reorder)
9528: {
9529:   PetscFunctionBegin;
9531:   PetscTryMethod(dm, "DMReorderSectionSetDefault_C", (DM, DMReorderDefaultFlag), (dm, reorder));
9532:   PetscFunctionReturn(PETSC_SUCCESS);
9533: }

9535: /*@
9536:   DMReorderSectionGetDefault - Get flag indicating whether the local section should be reordered by default

9538:   Not collective

9540:   Input Parameter:
9541: . dm - The DM

9543:   Output Parameter:
9544: . reorder - Flag for reordering

9546:   Level: intermediate

9548: .seealso: `DMReorderSetDefault()`
9549: @*/
9550: PetscErrorCode DMReorderSectionGetDefault(DM dm, DMReorderDefaultFlag *reorder)
9551: {
9552:   PetscFunctionBegin;
9554:   PetscAssertPointer(reorder, 2);
9555:   *reorder = DM_REORDER_DEFAULT_NOTSET;
9556:   PetscTryMethod(dm, "DMReorderSectionGetDefault_C", (DM, DMReorderDefaultFlag *), (dm, reorder));
9557:   PetscFunctionReturn(PETSC_SUCCESS);
9558: }

9560: /*@
9561:   DMReorderSectionSetType - Set the type of local section reordering

9563:   Logically collective

9565:   Input Parameters:
9566: + dm      - The DM
9567: - reorder - The reordering method

9569:   Level: intermediate

9571: .seealso: `DMReorderSectionGetType()`, `DMReorderSectionSetDefault()`
9572: @*/
9573: PetscErrorCode DMReorderSectionSetType(DM dm, MatOrderingType reorder)
9574: {
9575:   PetscFunctionBegin;
9577:   PetscTryMethod(dm, "DMReorderSectionSetType_C", (DM, MatOrderingType), (dm, reorder));
9578:   PetscFunctionReturn(PETSC_SUCCESS);
9579: }

9581: /*@
9582:   DMReorderSectionGetType - Get the reordering type for the local section

9584:   Not collective

9586:   Input Parameter:
9587: . dm - The DM

9589:   Output Parameter:
9590: . reorder - The reordering method

9592:   Level: intermediate

9594: .seealso: `DMReorderSetDefault()`, `DMReorderSectionGetDefault()`
9595: @*/
9596: PetscErrorCode DMReorderSectionGetType(DM dm, MatOrderingType *reorder)
9597: {
9598:   PetscFunctionBegin;
9600:   PetscAssertPointer(reorder, 2);
9601:   *reorder = NULL;
9602:   PetscTryMethod(dm, "DMReorderSectionGetType_C", (DM, MatOrderingType *), (dm, reorder));
9603:   PetscFunctionReturn(PETSC_SUCCESS);
9604: }