Actual source code: dm.c

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

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

 15: #if !defined(PETSC_HAVE_WINDOWS_COMPILERS)
 16: #include <petsc/private/valgrind/memcheck.h>
 17: #endif

 19: PetscClassId DM_CLASSID;
 20: PetscClassId DMLABEL_CLASSID;
 21: PetscLogEvent DM_Convert, DM_GlobalToLocal, DM_LocalToGlobal, DM_LocalToLocal, DM_LocatePoints, DM_Coarsen, DM_Refine, DM_CreateInterpolation, DM_CreateRestriction, DM_CreateInjection, DM_CreateMatrix, DM_CreateMassMatrix, DM_Load, DM_AdaptInterpolator, DM_ProjectFunction;

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

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

 35:   Collective

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

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

 43:   Level: beginner

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

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

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

 58:   PetscFunctionBegin;
 59:   PetscAssertPointer(dm, 2);
 60:   *dm = NULL;
 61:   PetscCall(DMInitializePackage());

 63:   PetscCall(PetscHeaderCreate(v, DM_CLASSID, "DM", "Distribution Manager", "DM", comm, DMDestroy, DMView));

 65:   ((PetscObject)v)->non_cyclic_references = &DMCountNonCyclicReferences;

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

107:   *dm = v;
108:   PetscFunctionReturn(PETSC_SUCCESS);
109: }

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

114:   Collective

116:   Input Parameter:
117: . dm - The original `DM` object

119:   Output Parameter:
120: . newdm - The new `DM` object

122:   Level: beginner

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

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

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

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

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

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

211:     PetscCall(DMGetPeriodicity(dm, &maxCell, &Lstart, &L));
212:     PetscCall(DMSetPeriodicity(*newdm, maxCell, Lstart, L));
213:   }
214:   {
215:     PetscBool useCone, useClosure;

217:     PetscCall(DMGetAdjacency(dm, PETSC_DEFAULT, &useCone, &useClosure));
218:     PetscCall(DMSetAdjacency(*newdm, PETSC_DEFAULT, useCone, useClosure));
219:   }
220:   PetscFunctionReturn(PETSC_SUCCESS);
221: }

223: /*@C
224:   DMSetVecType - Sets the type of vector created with `DMCreateLocalVector()` and `DMCreateGlobalVector()`

226:   Logically Collective

228:   Input Parameters:
229: + dm    - initial distributed array
230: - ctype - the vector type, for example `VECSTANDARD`, `VECCUDA`, or `VECVIENNACL`

232:   Options Database Key:
233: . -dm_vec_type ctype - the type of vector to create

235:   Level: intermediate

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

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

253: /*@C
254:   DMGetVecType - Gets the type of vector created with `DMCreateLocalVector()` and `DMCreateGlobalVector()`

256:   Logically Collective

258:   Input Parameter:
259: . da - initial distributed array

261:   Output Parameter:
262: . ctype - the vector type

264:   Level: intermediate

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

276: /*@
277:   VecGetDM - Gets the `DM` defining the data layout of the vector

279:   Not Collective

281:   Input Parameter:
282: . v - The `Vec`

284:   Output Parameter:
285: . dm - The `DM`

287:   Level: intermediate

289:   Note:
290:   A `Vec` may not have a `DM` associated with it.

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

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

306:   Not Collective

308:   Input Parameters:
309: + v  - The `Vec`
310: - dm - The `DM`

312:   Level: developer

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

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

319: .seealso: [](ch_dmbase), `DM`, `VecGetDM()`, `DMGetLocalVector()`, `DMGetGlobalVector()`, `DMSetVecType()`
320: @*/
321: PetscErrorCode VecSetDM(Vec v, DM dm)
322: {
323:   PetscFunctionBegin;
326:   PetscCall(PetscObjectCompose((PetscObject)v, "__PETSc_dm", (PetscObject)dm));
327:   PetscFunctionReturn(PETSC_SUCCESS);
328: }

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

333:   Logically Collective

335:   Input Parameters:
336: + dm    - the `DM` context
337: - ctype - the matrix type

339:   Options Database Key:
340: . -dm_is_coloring_type - global or local

342:   Level: intermediate

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

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

358:   Logically Collective

360:   Input Parameter:
361: . dm - the `DM` context

363:   Output Parameter:
364: . ctype - the matrix type

366:   Options Database Key:
367: . -dm_is_coloring_type - global or local

369:   Level: intermediate

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

382: /*@C
383:   DMSetMatType - Sets the type of matrix created with `DMCreateMatrix()`

385:   Logically Collective

387:   Input Parameters:
388: + dm    - the `DM` context
389: - ctype - the matrix type, for example `MATMPIAIJ`

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

394:   Level: intermediate

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

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

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

414:   Logically Collective

416:   Input Parameter:
417: . dm - the `DM` context

419:   Output Parameter:
420: . ctype - the matrix type

422:   Level: intermediate

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

434: /*@
435:   MatGetDM - Gets the `DM` defining the data layout of the matrix

437:   Not Collective

439:   Input Parameter:
440: . A - The `Mat`

442:   Output Parameter:
443: . dm - The `DM`

445:   Level: intermediate

447:   Note:
448:   A matrix may not have a `DM` associated with it

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

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

464: /*@
465:   MatSetDM - Sets the `DM` defining the data layout of the matrix

467:   Not Collective

469:   Input Parameters:
470: + A  - The `Mat`
471: - dm - The `DM`

473:   Level: developer

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

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

482: .seealso: [](ch_dmbase), `DM`, `MatGetDM()`, `DMCreateMatrix()`, `DMSetMatType()`
483: @*/
484: PetscErrorCode MatSetDM(Mat A, DM dm)
485: {
486:   PetscFunctionBegin;
489:   PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_dm", (PetscObject)dm));
490:   PetscFunctionReturn(PETSC_SUCCESS);
491: }

493: /*@C
494:   DMSetOptionsPrefix - Sets the prefix prepended to all option names when searching through the options database

496:   Logically Collective

498:   Input Parameters:
499: + dm     - the `DM` context
500: - prefix - the prefix to prepend

502:   Level: advanced

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

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

520: /*@C
521:   DMAppendOptionsPrefix - Appends an additional string to an already existing prefix used for searching for
522:   `DM` options in the options database.

524:   Logically Collective

526:   Input Parameters:
527: + dm     - the `DM` context
528: - prefix - the string to append to the current prefix

530:   Level: advanced

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

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

547: /*@C
548:   DMGetOptionsPrefix - Gets the prefix used for searching for all
549:   DM options in the options database.

551:   Not Collective

553:   Input Parameter:
554: . dm - the `DM` context

556:   Output Parameter:
557: . prefix - pointer to the prefix string used is returned

559:   Level: advanced

561:   Fortran Note:
562:   Pass in a string 'prefix' of
563:   sufficient length to hold the prefix.

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

575: static PetscErrorCode DMCountNonCyclicReferences_Internal(DM dm, PetscBool recurseCoarse, PetscBool recurseFine, PetscInt *ncrefct)
576: {
577:   PetscInt refct = ((PetscObject)dm)->refct;

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

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

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

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

611: PetscErrorCode DMDestroyLabelLinkList_Internal(DM dm)
612: {
613:   DMLabelLink next = dm->labels;

615:   PetscFunctionBegin;
616:   /* destroy the labels */
617:   while (next) {
618:     DMLabelLink tmp = next->next;

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

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

641: /*@C
642:   DMDestroy - Destroys a `DM`.

644:   Collective

646:   Input Parameter:
647: . dm - the `DM` object to destroy

649:   Level: developer

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

657:   PetscFunctionBegin;
658:   if (!*dm) PetscFunctionReturn(PETSC_SUCCESS);

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

671:   PetscCall(DMClearGlobalVectors(*dm));
672:   PetscCall(DMClearLocalVectors(*dm));
673:   PetscCall(DMClearNamedGlobalVectors(*dm));
674:   PetscCall(DMClearNamedLocalVectors(*dm));

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

738:       next = b->next;
739:       PetscCall(PetscFree(b));
740:     }
741:   }

743:   PetscCall(PetscObjectDestroy(&(*dm)->dmksp));
744:   PetscCall(PetscObjectDestroy(&(*dm)->dmsnes));
745:   PetscCall(PetscObjectDestroy(&(*dm)->dmts));

747:   if ((*dm)->ctx && (*dm)->ctxdestroy) PetscCall((*(*dm)->ctxdestroy)(&(*dm)->ctx));
748:   PetscCall(MatFDColoringDestroy(&(*dm)->fd));
749:   PetscCall(ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap));
750:   PetscCall(PetscFree((*dm)->vectype));
751:   PetscCall(PetscFree((*dm)->mattype));

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

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

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

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

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

804:   Collective

806:   Input Parameter:
807: . dm - the `DM` object to setup

809:   Level: intermediate

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

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

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

829:   Collective

831:   Input Parameter:
832: . dm - the `DM` object to set options for

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

887:   Level: intermediate

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

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

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

924: /*@C
925:   DMViewFromOptions - View a `DM` in a particular way based on a request in the options database

927:   Collective

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

934:   Level: intermediate

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

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

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

953:   Collective

955:   Input Parameters:
956: + dm - the `DM` object to view
957: - v  - the viewer

959:   Level: beginner

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

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

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

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

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

1005: /*@
1006:   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,
1007:   that is it has no ghost locations.

1009:   Collective

1011:   Input Parameter:
1012: . dm - the `DM` object

1014:   Output Parameter:
1015: . vec - the global vector

1017:   Level: beginner

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

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

1037: /*@
1038:   DMCreateLocalVector - Creates a local vector from a `DM` object.

1040:   Not Collective

1042:   Input Parameter:
1043: . dm - the `DM` object

1045:   Output Parameter:
1046: . vec - the local vector

1048:   Level: beginner

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

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

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

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

1074:   Collective

1076:   Input Parameter:
1077: . dm - the `DM` that provides the mapping

1079:   Output Parameter:
1080: . ltog - the mapping

1082:   Level: advanced

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

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

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

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

1099:   PetscFunctionBegin;
1101:   PetscAssertPointer(ltog, 2);
1102:   if (!dm->ltogmap) {
1103:     PetscSection section, sectionGlobal;

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

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

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

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

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

1164:   Not Collective

1166:   Input Parameter:
1167: . dm - the `DM` with block structure

1169:   Output Parameter:
1170: . bs - the block size, 1 implies no exploitable block structure

1172:   Level: intermediate

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

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

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

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

1196:   Collective

1198:   Input Parameters:
1199: + dmc - the `DM` object
1200: - dmf - the second, finer `DM` object

1202:   Output Parameters:
1203: + mat - the interpolation
1204: - vec - the scaling (optional, pass `NULL` if not needed), see `DMCreateInterpolationScale()`

1206:   Level: developer

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

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

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

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

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

1238:   Output Parameter:
1239: . scale - the scaled vector

1241:   Level: advanced

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

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

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

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

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

1286:   Collective

1288:   Input Parameters:
1289: + dmc - the `DM` object
1290: - dmf - the second, finer `DM` object

1292:   Output Parameter:
1293: . mat - the restriction

1295:   Level: developer

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

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

1315: /*@
1316:   DMCreateInjection - Gets injection matrix between two `DM` objects.

1318:   Collective

1320:   Input Parameters:
1321: + dac - the `DM` object
1322: - daf - the second, finer `DM` object

1324:   Output Parameter:
1325: . mat - the injection

1327:   Level: developer

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

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

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

1355: /*@
1356:   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
1357:   a Galerkin finite element model on the `DM`

1359:   Collective

1361:   Input Parameters:
1362: + dmc - the target `DM` object
1363: - dmf - the source `DM` object

1365:   Output Parameter:
1366: . mat - the mass matrix

1368:   Level: developer

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

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

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

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

1392:   Collective

1394:   Input Parameter:
1395: . dm - the `DM` object

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

1400:   Level: developer

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

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

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

1420:   Collective

1422:   Input Parameters:
1423: + dm    - the `DM` object
1424: - ctype - `IS_COLORING_LOCAL` or `IS_COLORING_GLOBAL`

1426:   Output Parameter:
1427: . coloring - the coloring

1429:   Level: developer

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

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

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

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

1453:   Collective

1455:   Input Parameter:
1456: . dm - the `DM` object

1458:   Output Parameter:
1459: . mat - the empty Jacobian

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

1464:   Level: beginner

1466:   Notes:
1467:   This properly preallocates the number of nonzeros in the sparse matrix so you
1468:   do not need to do it yourself.

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

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

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

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

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

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

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

1526:   Logically Collective

1528:   Input Parameters:
1529: + dm   - the `DM`
1530: - skip - `PETSC_TRUE` to skip preallocation

1532:   Level: developer

1534:   Note:
1535:   This is most useful to reduce initialization costs when `MatSetPreallocationCOO()` and
1536:   `MatSetValuesCOO()` will be used.

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

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

1552:   Logically Collective

1554:   Input Parameters:
1555: + dm   - the `DM`
1556: - only - `PETSC_TRUE` if only want preallocation

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

1561:   Level: developer

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

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

1577:   Logically Collective

1579:   Input Parameters:
1580: + dm   - the `DM`
1581: - only - `PETSC_TRUE` if you only want matrix structure

1583:   Level: developer

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

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

1598:   Logically Collective

1600:   Input Parameters:
1601: + dm    - the `DM`
1602: - btype - block by topological point or field node

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

1607:   Level: advanced

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

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

1622:   Not Collective

1624:   Input Parameter:
1625: . dm - the `DM`

1627:   Output Parameter:
1628: . btype - block by topological point or field node

1630:   Level: advanced

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

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

1646:   Not Collective

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

1653:   Output Parameter:
1654: . mem - the work array

1656:   Level: developer

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

1661:   The array may contain nonzero values

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

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

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

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

1713:   Not Collective

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

1720:   Output Parameter:
1721: . mem - the work array

1723:   Level: developer

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

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

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

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

1756:   Logically Collective; No Fortran Support

1758:   Input Parameters:
1759: + dm     - The `DM`
1760: . field  - The field number for the nullspace
1761: - nullsp - A callback to create the nullspace

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

1769:   Level: intermediate

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

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

1785:   Not Collective; No Fortran Support

1787:   Input Parameters:
1788: + dm    - The `DM`
1789: - field - The field number for the nullspace

1791:   Output Parameter:
1792: . nullsp - A callback to create the nullspace

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

1800:   Level: intermediate

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

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

1817:   Logically Collective; No Fortran Support

1819:   Input Parameters:
1820: + dm     - The `DM`
1821: . field  - The field number for the nullspace
1822: - nullsp - A callback to create the near-nullspace

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

1830:   Level: intermediate

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

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

1847:   Not Collective; No Fortran Support

1849:   Input Parameters:
1850: + dm    - The `DM`
1851: - field - The field number for the nullspace

1853:   Output Parameter:
1854: . nullsp - A callback to create the near-nullspace

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

1862:   Level: intermediate

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

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

1880:   Not Collective; No Fortran Support

1882:   Input Parameter:
1883: . dm - the `DM` object

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

1890:   Level: intermediate

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

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

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

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

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

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

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

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

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

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

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

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

2005:   Not Collective; No Fortran Support

2007:   Input Parameter:
2008: . dm - the `DM` object

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

2016:   Level: intermediate

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

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

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

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

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

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

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

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

2091: /*@C
2092:   DMCreateSubDM - Returns an `IS` and `DM` encapsulating a subproblem defined by the fields passed in.
2093:   The fields are defined by `DMCreateFieldIS()`.

2095:   Not collective

2097:   Input Parameters:
2098: + dm        - The `DM` object
2099: . numFields - The number of fields to select
2100: - fields    - The field numbers of the selected fields

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

2106:   Level: intermediate

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

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

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

2127:   Not collective

2129:   Input Parameters:
2130: + dms - The `DM` objects
2131: - n   - The number of `DM`s

2133:   Output Parameters:
2134: + is      - The global indices for each of subproblem within the super `DM`, or NULL
2135: - superdm - The `DM` for the superproblem

2137:   Level: intermediate

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

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

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

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

2166:   Not Collective

2168:   Input Parameter:
2169: . dm - the `DM` object

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

2178:   Level: intermediate

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

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

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

2191:   Developer Notes:
2192:   The `dmlist` is for the inner subdomains or the outer subdomains or all subdomains?

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

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

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

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

2252:   Not Collective

2254:   Input Parameters:
2255: + dm     - the `DM` object
2256: . n      - the number of subdomains
2257: - subdms - the local subdomains

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

2264:   Level: developer

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

2272:   Developer Note:
2273:   Can the subdms input be anything or are they exactly the `DM` obtained from
2274:   `DMCreateDomainDecomposition()`?

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

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

2290:   Collective

2292:   Input Parameters:
2293: + dm   - the `DM` object
2294: - comm - the communicator to contain the new `DM` object (or `MPI_COMM_NULL`)

2296:   Output Parameter:
2297: . dmf - the refined `DM`, or `NULL`

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

2302:   Level: developer

2304:   Note:
2305:   If no refinement was done, the return value is `NULL`

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

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

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

2323:     (*dmf)->ctx       = dm->ctx;
2324:     (*dmf)->leveldown = dm->leveldown;
2325:     (*dmf)->levelup   = dm->levelup + 1;

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

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

2339:   Logically Collective; No Fortran Support

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

2347:   Calling sequence of `refinehook`:
2348: + coarse - coarse level `DM`
2349: . fine   - fine level `DM` to interpolate problem to
2350: - ctx    - optional user-defined function context

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

2358:   Level: advanced

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

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

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

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

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

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

2392:   Logically Collective; No Fortran Support

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

2400:   Level: advanced

2402:   Note:
2403:   This function does nothing if the hook is not in the list.

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

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

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

2427:   Collective if any hooks are

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

2434:   Level: developer

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

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

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

2453: /*@
2454:   DMInterpolateSolution - Interpolates a solution from a coarse mesh to a fine mesh.

2456:   Collective

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

2466:   Output Parameter:
2467: . fineSol - the interpolation of coarseSol to the fine mesh

2469:   Level: developer

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

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

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

2486:   PetscFunctionBegin;

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

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

2504:   Not Collective

2506:   Input Parameter:
2507: . dm - the `DM` object

2509:   Output Parameter:
2510: . level - number of refinements

2512:   Level: developer

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

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

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

2530:   Not Collective

2532:   Input Parameters:
2533: + dm    - the `DM` object
2534: - level - number of refinements

2536:   Level: advanced

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

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

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

2553: /*@
2554:   DMExtrude - Extrude a `DM` object from a surface

2556:   Collective

2558:   Input Parameters:
2559: + dm     - the `DM` object
2560: - layers - the number of extruded cell layers

2562:   Output Parameter:
2563: . dme - the extruded `DM`, or `NULL`

2565:   Level: developer

2567:   Note:
2568:   If no extrusion was done, the return value is `NULL`

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

2586: PetscErrorCode DMGetBasisTransformDM_Internal(DM dm, DM *tdm)
2587: {
2588:   PetscFunctionBegin;
2590:   PetscAssertPointer(tdm, 2);
2591:   *tdm = dm->transformDM;
2592:   PetscFunctionReturn(PETSC_SUCCESS);
2593: }

2595: PetscErrorCode DMGetBasisTransformVec_Internal(DM dm, Vec *tv)
2596: {
2597:   PetscFunctionBegin;
2599:   PetscAssertPointer(tv, 2);
2600:   *tv = dm->transform;
2601:   PetscFunctionReturn(PETSC_SUCCESS);
2602: }

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

2607:   Input Parameter:
2608: . dm - The `DM`

2610:   Output Parameter:
2611: . flg - `PETSC_TRUE` if a basis transformation should be done

2613:   Level: developer

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

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

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

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

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

2677: PetscErrorCode DMCopyTransform(DM dm, DM newdm)
2678: {
2679:   PetscFunctionBegin;
2682:   newdm->transformCtx       = dm->transformCtx;
2683:   newdm->transformSetUp     = dm->transformSetUp;
2684:   newdm->transformDestroy   = NULL;
2685:   newdm->transformGetMatrix = dm->transformGetMatrix;
2686:   if (newdm->transformSetUp) PetscCall(DMConstructBasisTransform_Internal(newdm));
2687:   PetscFunctionReturn(PETSC_SUCCESS);
2688: }

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

2693:   Logically Collective

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

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

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

2715:   Level: advanced

2717:   Note:
2718:   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.

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

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

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

2745:   PetscFunctionBegin;
2746:   (void)g;
2747:   (void)ctx;
2749:   PetscCall(DMGetDefaultConstraints(dm, &cSec, &cMat, &cBias));
2750:   if (cMat && (mode == INSERT_VALUES || mode == INSERT_ALL_VALUES || mode == INSERT_BC_VALUES)) {
2751:     PetscInt nRows;

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

2773: /*@
2774:   DMGlobalToLocal - update local vectors from global vector

2776:   Neighbor-wise Collective

2778:   Input Parameters:
2779: + dm   - the `DM` object
2780: . g    - the global vector
2781: . mode - `INSERT_VALUES` or `ADD_VALUES`
2782: - l    - the local vector

2784:   Level: beginner

2786:   Notes:
2787:   The communication involved in this update can be overlapped with computation by instead using
2788:   `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()`.

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

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

2804: /*@
2805:   DMGlobalToLocalBegin - Begins updating local vectors from global vector

2807:   Neighbor-wise Collective

2809:   Input Parameters:
2810: + dm   - the `DM` object
2811: . g    - the global vector
2812: . mode - `INSERT_VALUES` or `ADD_VALUES`
2813: - l    - the local vector

2815:   Level: intermediate

2817:   Notes:
2818:   The operation is completed with `DMGlobalToLocalEnd()`

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

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

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

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

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

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

2856: /*@
2857:   DMGlobalToLocalEnd - Ends updating local vectors from global vector

2859:   Neighbor-wise Collective

2861:   Input Parameters:
2862: + dm   - the `DM` object
2863: . g    - the global vector
2864: . mode - `INSERT_VALUES` or `ADD_VALUES`
2865: - l    - the local vector

2867:   Level: intermediate

2869:   Note:
2870:   See `DMGlobalToLocalBegin()` for details.

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

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

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

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

2909:   Logically Collective

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

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

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

2931:   Level: advanced

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

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

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

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

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

2992:   Neighbor-wise Collective

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

3000:   Level: beginner

3002:   Notes:
3003:   The communication involved in this update can be overlapped with computation by using
3004:   `DMLocalToGlobalBegin()` and `DMLocalToGlobalEnd()`.

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

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

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

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

3022: /*@
3023:   DMLocalToGlobalBegin - begins updating global vectors from local vectors

3025:   Neighbor-wise Collective

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

3033:   Level: intermediate

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

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

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

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

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

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

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

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

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

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

3158: /*@
3159:   DMLocalToGlobalEnd - updates global vectors from local vectors

3161:   Neighbor-wise Collective

3163:   Input Parameters:
3164: + dm   - the `DM` object
3165: . l    - the local vector
3166: . mode - `INSERT_VALUES` or `ADD_VALUES`
3167: - g    - the global vector

3169:   Level: intermediate

3171:   Note:
3172:   See `DMLocalToGlobalBegin()` for full details

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

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

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

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

3235:   Neighbor-wise Collective

3237:   Input Parameters:
3238: + dm   - the `DM` object
3239: . g    - the original local vector
3240: - mode - one of `INSERT_VALUES` or `ADD_VALUES`

3242:   Output Parameter:
3243: . l - the local vector with correct ghost values

3245:   Level: intermediate

3247:   Note:
3248:   Must be followed by `DMLocalToLocalEnd()`.

3250: .seealso: [](ch_dmbase), `DM`, `DMLocalToLocalEnd()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateLocalVector()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`
3251: @*/
3252: PetscErrorCode DMLocalToLocalBegin(DM dm, Vec g, InsertMode mode, Vec l)
3253: {
3254:   PetscFunctionBegin;
3258:   PetscUseTypeMethod(dm, localtolocalbegin, g, mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode), l);
3259:   PetscFunctionReturn(PETSC_SUCCESS);
3260: }

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

3266:   Neighbor-wise Collective

3268:   Input Parameters:
3269: + dm   - the `DM` object
3270: . g    - the original local vector
3271: - mode - one of `INSERT_VALUES` or `ADD_VALUES`

3273:   Output Parameter:
3274: . l - the local vector with correct ghost values

3276:   Level: intermediate

3278: .seealso: [](ch_dmbase), `DM`, `DMLocalToLocalBegin()`, `DMCoarsen()`, `DMDestroy()`, `DMView()`, `DMCreateLocalVector()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`
3279: @*/
3280: PetscErrorCode DMLocalToLocalEnd(DM dm, Vec g, InsertMode mode, Vec l)
3281: {
3282:   PetscFunctionBegin;
3286:   PetscUseTypeMethod(dm, localtolocalend, g, mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode), l);
3287:   PetscFunctionReturn(PETSC_SUCCESS);
3288: }

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

3293:   Collective

3295:   Input Parameters:
3296: + dm   - the `DM` object
3297: - comm - the communicator to contain the new `DM` object (or `MPI_COMM_NULL`)

3299:   Output Parameter:
3300: . dmc - the coarsened `DM`

3302:   Level: developer

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

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

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

3336:   Logically Collective; No Fortran Support

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

3344:   Calling sequence of `coarsenhook`:
3345: + fine   - fine level `DM`
3346: . coarse - coarse level `DM` to restrict problem to
3347: - ctx    - optional user-defined function context

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

3357:   Level: advanced

3359:   Notes:
3360:   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`.

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

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

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

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

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

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

3392:   Logically Collective; No Fortran Support

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

3400:   Level: advanced

3402:   Notes:
3403:   This function does nothing if the `coarsenhook` is not in the list.

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

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

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

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

3429:   Collective if any hooks are

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

3438:   Level: developer

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

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

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

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

3459:   Logically Collective; No Fortran Support

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

3467:   Calling sequence of `ddhook`:
3468: + global - global `DM`
3469: . block  - subdomain `DM`
3470: - ctx    - optional user-defined function context

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

3479:   Level: advanced

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

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

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

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

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

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

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

3515:   Logically Collective; No Fortran Support

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

3523:   Level: advanced

3525:   Note:
3526:   See `DMSubDomainHookAdd()` for the calling sequences of `ddhook` and `restricthook`

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

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

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

3551:   Collective if any hooks are

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

3559:   Level: developer

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

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

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

3577:   Not Collective

3579:   Input Parameter:
3580: . dm - the `DM` object

3582:   Output Parameter:
3583: . level - number of coarsenings

3585:   Level: developer

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

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

3601:   Collective

3603:   Input Parameters:
3604: + dm    - the `DM` object
3605: - level - number of coarsenings

3607:   Level: developer

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

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

3622: /*@C
3623:   DMRefineHierarchy - Refines a `DM` object, all levels at once

3625:   Collective

3627:   Input Parameters:
3628: + dm      - the `DM` object
3629: - nlevels - the number of levels of refinement

3631:   Output Parameter:
3632: . dmf - the refined `DM` hierarchy

3634:   Level: developer

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

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

3654: /*@C
3655:   DMCoarsenHierarchy - Coarsens a `DM` object, all levels at once

3657:   Collective

3659:   Input Parameters:
3660: + dm      - the `DM` object
3661: - nlevels - the number of levels of coarsening

3663:   Output Parameter:
3664: . dmc - the coarsened `DM` hierarchy

3666:   Level: developer

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

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

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

3689:   Logically Collective if the function is collective

3691:   Input Parameters:
3692: + dm      - the `DM` object
3693: - destroy - the destroy function

3695:   Level: intermediate

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

3707: /*@
3708:   DMSetApplicationContext - Set a user context into a `DM` object

3710:   Not Collective

3712:   Input Parameters:
3713: + dm  - the `DM` object
3714: - ctx - the user context

3716:   Level: intermediate

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

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

3733: /*@
3734:   DMGetApplicationContext - Gets a user context from a `DM` object

3736:   Not Collective

3738:   Input Parameter:
3739: . dm - the `DM` object

3741:   Output Parameter:
3742: . ctx - the user context

3744:   Level: intermediate

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

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

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

3762:   Logically Collective

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

3768:   Level: intermediate

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

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

3784:   Not Collective

3786:   Input Parameter:
3787: . dm - the `DM` object to destroy

3789:   Output Parameter:
3790: . flg - `PETSC_TRUE` if the variable bounds function exists

3792:   Level: developer

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

3805: /*@C
3806:   DMComputeVariableBounds - compute variable bounds used by `SNESVI`.

3808:   Logically Collective

3810:   Input Parameter:
3811: . dm - the `DM` object

3813:   Output Parameters:
3814: + xl - lower bound
3815: - xu - upper bound

3817:   Level: advanced

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

3822: .seealso: [](ch_dmbase), `DM`, `DMHasVariableBounds()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMGetApplicationContext()`
3823: @*/
3824: PetscErrorCode DMComputeVariableBounds(DM dm, Vec xl, Vec xu)
3825: {
3826:   PetscFunctionBegin;
3830:   PetscUseTypeMethod(dm, computevariablebounds, xl, xu);
3831:   PetscFunctionReturn(PETSC_SUCCESS);
3832: }

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

3837:   Not Collective

3839:   Input Parameter:
3840: . dm - the DM object

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

3845:   Level: developer

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

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

3861:   Not Collective

3863:   Input Parameter:
3864: . dm - the `DM` object

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

3869:   Level: developer

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

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

3885:   Not Collective

3887:   Input Parameter:
3888: . dm - the `DM` object

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

3893:   Level: developer

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

3907: PetscFunctionList DMList              = NULL;
3908: PetscBool         DMRegisterAllCalled = PETSC_FALSE;

3910: /*@C
3911:   DMSetType - Builds a `DM`, for a particular `DM` implementation.

3913:   Collective

3915:   Input Parameters:
3916: + dm     - The `DM` object
3917: - method - The name of the `DMType`, for example `DMDA`, `DMPLEX`

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

3922:   Level: intermediate

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

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

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

3939:   PetscCall(DMRegisterAll());
3940:   PetscCall(PetscFunctionListFind(DMList, method, &r));
3941:   PetscCheck(r, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);

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

3950: /*@C
3951:   DMGetType - Gets the `DM` type name (as a string) from the `DM`.

3953:   Not Collective

3955:   Input Parameter:
3956: . dm - The `DM`

3958:   Output Parameter:
3959: . type - The `DMType` name

3961:   Level: intermediate

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

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

3978:   Collective

3980:   Input Parameters:
3981: + dm      - the `DM`
3982: - newtype - new `DM` type (use "same" for the same type)

3984:   Output Parameter:
3985: . M - pointer to new `DM`

3987:   Level: intermediate

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

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

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

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

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

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

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

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

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

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

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

4084: /*--------------------------------------------------------------------------------------------------------------------*/

4086: /*@C
4087:   DMRegister -  Adds a new `DM` type implementation

4089:   Not Collective

4091:   Input Parameters:
4092: + sname    - The name of a new user-defined creation routine
4093: - function - The creation routine itself

4095:   Level: advanced

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

4100:   Example Usage:
4101: .vb
4102:     DMRegister("my_da", MyDMCreate);
4103: .ve

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

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

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

4128:   Collective

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

4136:   Level: intermediate

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

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

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

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

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

4174: /******************************** FEM Support **********************************/

4176: PetscErrorCode DMPrintCellIndices(PetscInt c, const char name[], PetscInt len, const PetscInt x[])
4177: {
4178:   PetscInt f;

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

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

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

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

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

4206: PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[])
4207: {
4208:   PetscInt f, g;

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

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

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

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

4253:   Input Parameter:
4254: . dm - The `DM`

4256:   Output Parameter:
4257: . section - The `PetscSection`

4259:   Options Database Key:
4260: . -dm_petscsection_view - View the `PetscSection` created by the `DM`

4262:   Level: advanced

4264:   Notes:
4265:   Use `DMGetLocalSection()` in new code.

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

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

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

4281:   Input Parameter:
4282: . dm - The `DM`

4284:   Output Parameter:
4285: . section - The `PetscSection`

4287:   Options Database Key:
4288: . -dm_petscsection_view - View the section created by the `DM`

4290:   Level: intermediate

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

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

4305:     if (dm->setfromoptionscalled) {
4306:       PetscObject       obj = (PetscObject)dm;
4307:       PetscViewer       viewer;
4308:       PetscViewerFormat format;
4309:       PetscBool         flg;

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

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

4333:   Input Parameters:
4334: + dm      - The `DM`
4335: - section - The `PetscSection`

4337:   Level: advanced

4339:   Notes:
4340:   Use `DMSetLocalSection()` in new code.

4342:   Any existing `PetscSection` will be destroyed

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

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

4356:   Input Parameters:
4357: + dm      - The `DM`
4358: - section - The `PetscSection`

4360:   Level: intermediate

4362:   Note:
4363:   Any existing Section will be destroyed

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

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

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

4396:   /* Clear scratch vectors */
4397:   PetscCall(DMClearGlobalVectors(dm));
4398:   PetscCall(DMClearLocalVectors(dm));
4399:   PetscCall(DMClearNamedGlobalVectors(dm));
4400:   PetscCall(DMClearNamedLocalVectors(dm));
4401:   PetscFunctionReturn(PETSC_SUCCESS);
4402: }

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

4407:   Input Parameter:
4408: . dm - The `DM`

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

4414:   Level: developer

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

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

4430:   not Collective

4432:   Input Parameter:
4433: . dm - The `DM`

4435:   Output Parameters:
4436: + 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.
4437: . 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.
4438: - bias    - Vector containing bias to be added to constrained dofs

4440:   Level: advanced

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

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

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

4461:   Collective

4463:   Input Parameters:
4464: + dm      - The `DM`
4465: . 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).
4466: . 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).
4467: - 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).

4469:   Level: advanced

4471:   Notes:
4472:   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()`.

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

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

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

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

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

4517:   Input Parameters:
4518: + dm - The `DM`
4519: . localSection - `PetscSection` describing the local data layout
4520: - globalSection - `PetscSection` describing the global data layout

4522:   Level: intermediate

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

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

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

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

4591: static PetscErrorCode DMGetIsoperiodicPointSF_Internal(DM dm, PetscSF *sf)
4592: {
4593:   PetscErrorCode (*f)(DM, PetscSF *);

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

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

4607:   Collective

4609:   Input Parameter:
4610: . dm - The `DM`

4612:   Output Parameter:
4613: . section - The `PetscSection`

4615:   Level: intermediate

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

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

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

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

4647:   Input Parameters:
4648: + dm      - The `DM`
4649: - section - The PetscSection, or `NULL`

4651:   Level: intermediate

4653:   Note:
4654:   Any existing `PetscSection` will be destroyed

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

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

4681:   Input Parameter:
4682: . dm - The `DM`

4684:   Output Parameter:
4685: . sf - The `PetscSF`

4687:   Level: intermediate

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

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

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

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

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

4722:   Input Parameters:
4723: + dm - The `DM`
4724: - sf - The `PetscSF`

4726:   Level: intermediate

4728:   Note:
4729:   Any previous `PetscSF` is destroyed

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

4744: /*@C
4745:   DMCreateSectionSF - Create the `PetscSF` encoding the parallel dof overlap for the `DM` based upon the `PetscSection`s
4746:   describing the data layout.

4748:   Input Parameters:
4749: + dm            - The `DM`
4750: . localSection  - `PetscSection` describing the local data layout
4751: - globalSection - `PetscSection` describing the global data layout

4753:   Level: developer

4755:   Note:
4756:   One usually uses `DMGetSectionSF()` to obtain the `PetscSF`

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

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

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

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

4779:   Input Parameter:
4780: . dm - The `DM`

4782:   Output Parameter:
4783: . sf - The `PetscSF`

4785:   Level: intermediate

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

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

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

4804:   Collective

4806:   Input Parameters:
4807: + dm - The `DM`
4808: - sf - The `PetscSF`

4810:   Level: intermediate

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

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

4828:   Input Parameter:
4829: . dm - The `DM`

4831:   Output Parameter:
4832: . sf - The `PetscSF`

4834:   Level: intermediate

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

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

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

4853:   Input Parameters:
4854: + dm - The DM
4855: - sf - The PetscSF

4857:   Level: intermediate

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

4872: static PetscErrorCode DMSetDefaultAdjacency_Private(DM dm, PetscInt f, PetscObject disc)
4873: {
4874:   PetscClassId id;

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

4888: static PetscErrorCode DMFieldEnlarge_Static(DM dm, PetscInt NfNew)
4889: {
4890:   RegionField *tmpr;
4891:   PetscInt     Nf = dm->Nf, f;

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

4908: /*@
4909:   DMClearFields - Remove all fields from the `DM`

4911:   Logically Collective

4913:   Input Parameter:
4914: . dm - The `DM`

4916:   Level: intermediate

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

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

4936: /*@
4937:   DMGetNumFields - Get the number of fields in the `DM`

4939:   Not Collective

4941:   Input Parameter:
4942: . dm - The `DM`

4944:   Output Parameter:
4945: . numFields - The number of fields

4947:   Level: intermediate

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

4960: /*@
4961:   DMSetNumFields - Set the number of fields in the `DM`

4963:   Logically Collective

4965:   Input Parameters:
4966: + dm        - The `DM`
4967: - numFields - The number of fields

4969:   Level: intermediate

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

4977:   PetscFunctionBegin;
4979:   PetscCall(DMGetNumFields(dm, &Nf));
4980:   for (f = Nf; f < numFields; ++f) {
4981:     PetscContainer obj;

4983:     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)dm), &obj));
4984:     PetscCall(DMAddField(dm, NULL, (PetscObject)obj));
4985:     PetscCall(PetscContainerDestroy(&obj));
4986:   }
4987:   PetscFunctionReturn(PETSC_SUCCESS);
4988: }

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

4993:   Not Collective

4995:   Input Parameters:
4996: + dm - The `DM`
4997: - f  - The field number

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

5003:   Level: intermediate

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

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

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

5036:   Logically Collective

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

5044:   Level: intermediate

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

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

5065:   Logically Collective

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

5072:   Level: intermediate

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

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

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

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

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

5104:   Logically Collective

5106:   Input Parameters:
5107: + dm          - The `DM`
5108: . f           - The field index
5109: - avoidTensor - `PETSC_TRUE` to skip defining the field on tensor cells

5111:   Level: intermediate

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

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

5126:   Not Collective

5128:   Input Parameters:
5129: + dm - The `DM`
5130: - f  - The field index

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

5135:   Level: intermediate

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

5147: /*@
5148:   DMCopyFields - Copy the discretizations for the `DM` into another `DM`

5150:   Collective

5152:   Input Parameter:
5153: . dm - The `DM`

5155:   Output Parameter:
5156: . newdm - The `DM`

5158:   Level: advanced

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

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

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

5183: /*@
5184:   DMGetAdjacency - Returns the flags for determining variable influence

5186:   Not Collective

5188:   Input Parameters:
5189: + dm - The `DM` object
5190: - f  - The field number, or `PETSC_DEFAULT` for the default adjacency

5192:   Output Parameters:
5193: + useCone    - Flag for variable influence starting with the cone operation
5194: - useClosure - Flag for variable influence using transitive closure

5196:   Level: developer

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

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

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

5228: /*@
5229:   DMSetAdjacency - Set the flags for determining variable influence

5231:   Not Collective

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

5239:   Level: developer

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

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

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

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

5272:   Not collective

5274:   Input Parameter:
5275: . dm - The `DM` object

5277:   Output Parameters:
5278: + useCone    - Flag for variable influence starting with the cone operation
5279: - useClosure - Flag for variable influence using transitive closure

5281:   Level: developer

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

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

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

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

5312:   Not Collective

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

5319:   Level: developer

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

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

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

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

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

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

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

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

5433: static PetscErrorCode DMDSEnlarge_Static(DM dm, PetscInt NdsNew)
5434: {
5435:   DMSpace *tmpd;
5436:   PetscInt Nds = dm->Nds, s;

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

5453: /*@
5454:   DMGetNumDS - Get the number of discrete systems in the `DM`

5456:   Not Collective

5458:   Input Parameter:
5459: . dm - The `DM`

5461:   Output Parameter:
5462: . Nds - The number of `PetscDS` objects

5464:   Level: intermediate

5466: .seealso: [](ch_dmbase), `DM`, `DMGetDS()`, `DMGetCellDS()`
5467: @*/
5468: PetscErrorCode DMGetNumDS(DM dm, PetscInt *Nds)
5469: {
5470:   PetscFunctionBegin;
5472:   PetscAssertPointer(Nds, 2);
5473:   *Nds = dm->Nds;
5474:   PetscFunctionReturn(PETSC_SUCCESS);
5475: }

5477: /*@
5478:   DMClearDS - Remove all discrete systems from the `DM`

5480:   Logically Collective

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

5485:   Level: intermediate

5487: .seealso: [](ch_dmbase), `DM`, `DMGetNumDS()`, `DMGetDS()`, `DMSetField()`
5488: @*/
5489: PetscErrorCode DMClearDS(DM dm)
5490: {
5491:   PetscInt s;

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

5507: /*@
5508:   DMGetDS - Get the default `PetscDS`

5510:   Not Collective

5512:   Input Parameter:
5513: . dm - The `DM`

5515:   Output Parameter:
5516: . ds - The default `PetscDS`

5518:   Level: intermediate

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

5532: /*@
5533:   DMGetCellDS - Get the `PetscDS` defined on a given cell

5535:   Not Collective

5537:   Input Parameters:
5538: + dm    - The `DM`
5539: - point - Cell for the `PetscDS`

5541:   Output Parameters:
5542: + ds   - The `PetscDS` defined on the given cell
5543: - dsIn - The `PetscDS` for input on the given cell, or NULL if the same ds

5545:   Level: developer

5547: .seealso: [](ch_dmbase), `DM`, `DMGetDS()`, `DMSetRegionDS()`
5548: @*/
5549: PetscErrorCode DMGetCellDS(DM dm, PetscInt point, PetscDS *ds, PetscDS *dsIn)
5550: {
5551:   PetscDS  dsDef = NULL;
5552:   PetscInt s;

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

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

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

5582:   Not Collective

5584:   Input Parameters:
5585: + dm    - The `DM`
5586: - label - The `DMLabel` defining the mesh region, or `NULL` for the entire mesh

5588:   Output Parameters:
5589: + fields - The `IS` containing the `DM` field numbers for the fields in this `PetscDS`, or `NULL`
5590: . ds     - The `PetscDS` defined on the given region, or `NULL`
5591: - dsIn   - The `PetscDS` for input in the given region, or `NULL`

5593:   Level: advanced

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

5600: .seealso: [](ch_dmbase), `DM`, `DMGetRegionNumDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5601: @*/
5602: PetscErrorCode DMGetRegionDS(DM dm, DMLabel label, IS *fields, PetscDS *ds, PetscDS *dsIn)
5603: {
5604:   PetscInt Nds = dm->Nds, s;

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

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

5635:   Collective

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

5644:   Level: advanced

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

5650: .seealso: [](ch_dmbase), `DM`, `DMGetRegionDS()`, `DMSetRegionNumDS()`, `DMGetDS()`, `DMGetCellDS()`
5651: @*/
5652: PetscErrorCode DMSetRegionDS(DM dm, DMLabel label, IS fields, PetscDS ds, PetscDS dsIn)
5653: {
5654:   PetscInt Nds = dm->Nds, s;

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

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

5691:   Not Collective

5693:   Input Parameters:
5694: + dm  - The `DM`
5695: - num - The region number, in [0, Nds)

5697:   Output Parameters:
5698: + label  - The region label, or `NULL`
5699: . fields - The `IS` containing the `DM` field numbers for the fields in this `PetscDS`, or `NULL`
5700: . ds     - The `PetscDS` defined on the given region, or `NULL`
5701: - dsIn   - The `PetscDS` for input in the given region, or `NULL`

5703:   Level: advanced

5705: .seealso: [](ch_dmbase), `DM`, `DMGetRegionDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5706: @*/
5707: PetscErrorCode DMGetRegionNumDS(DM dm, PetscInt num, DMLabel *label, IS *fields, PetscDS *ds, PetscDS *dsIn)
5708: {
5709:   PetscInt Nds;

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

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

5737:   Not Collective

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

5747:   Level: advanced

5749: .seealso: [](ch_dmbase), `DM`, `DMGetRegionDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5750: @*/
5751: PetscErrorCode DMSetRegionNumDS(DM dm, PetscInt num, DMLabel label, IS fields, PetscDS ds, PetscDS dsIn)
5752: {
5753:   PetscInt Nds;

5755:   PetscFunctionBegin;
5758:   PetscCall(DMGetNumDS(dm, &Nds));
5759:   PetscCheck((num >= 0) && (num < Nds), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Region number %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", num, Nds);
5760:   PetscCall(PetscObjectReference((PetscObject)label));
5761:   PetscCall(DMLabelDestroy(&dm->probs[num].label));
5762:   dm->probs[num].label = label;
5763:   if (fields) {
5765:     PetscCall(PetscObjectReference((PetscObject)fields));
5766:     PetscCall(ISDestroy(&dm->probs[num].fields));
5767:     dm->probs[num].fields = fields;
5768:   }
5769:   if (ds) {
5771:     PetscCall(PetscObjectReference((PetscObject)ds));
5772:     PetscCall(PetscDSDestroy(&dm->probs[num].ds));
5773:     dm->probs[num].ds = ds;
5774:   }
5775:   if (dsIn) {
5777:     PetscCall(PetscObjectReference((PetscObject)dsIn));
5778:     PetscCall(PetscDSDestroy(&dm->probs[num].dsIn));
5779:     dm->probs[num].dsIn = dsIn;
5780:   }
5781:   PetscFunctionReturn(PETSC_SUCCESS);
5782: }

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

5787:   Not Collective

5789:   Input Parameters:
5790: + dm - The `DM`
5791: - ds - The `PetscDS` defined on the given region

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

5796:   Level: advanced

5798: .seealso: [](ch_dmbase), `DM`, `DMGetRegionNumDS()`, `DMGetRegionDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5799: @*/
5800: PetscErrorCode DMFindRegionNum(DM dm, PetscDS ds, PetscInt *num)
5801: {
5802:   PetscInt Nds, n;

5804:   PetscFunctionBegin;
5807:   PetscAssertPointer(num, 3);
5808:   PetscCall(DMGetNumDS(dm, &Nds));
5809:   for (n = 0; n < Nds; ++n)
5810:     if (ds == dm->probs[n].ds) break;
5811:   if (n >= Nds) *num = -1;
5812:   else *num = n;
5813:   PetscFunctionReturn(PETSC_SUCCESS);
5814: }

5816: /*@C
5817:   DMCreateFEDefault - Create a `PetscFE` based on the celltype for the mesh

5819:   Not Collective

5821:   Input Parameters:
5822: + dm     - The `DM`
5823: . Nc     - The number of components for the field
5824: . prefix - The options prefix for the output `PetscFE`, or `NULL`
5825: - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree

5827:   Output Parameter:
5828: . fem - The `PetscFE`

5830:   Level: intermediate

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

5835: .seealso: [](ch_dmbase), `DM`, `PetscFECreateByCell()`, `DMAddField()`, `DMCreateDS()`, `DMGetCellDS()`, `DMGetRegionDS()`
5836: @*/
5837: PetscErrorCode DMCreateFEDefault(DM dm, PetscInt Nc, const char prefix[], PetscInt qorder, PetscFE *fem)
5838: {
5839:   DMPolytopeType ct;
5840:   PetscInt       dim, cStart;

5842:   PetscFunctionBegin;
5845:   if (prefix) PetscAssertPointer(prefix, 3);
5847:   PetscAssertPointer(fem, 5);
5848:   PetscCall(DMGetDimension(dm, &dim));
5849:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
5850:   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
5851:   PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, Nc, ct, prefix, qorder, fem));
5852:   PetscFunctionReturn(PETSC_SUCCESS);
5853: }

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

5858:   Collective

5860:   Input Parameter:
5861: . dm - The `DM`

5863:   Options Database Key:
5864: . -dm_petscds_view - View all the `PetscDS` objects in this `DM`

5866:   Level: intermediate

5868: .seealso: [](ch_dmbase), `DM`, `DMSetField`, `DMAddField()`, `DMGetDS()`, `DMGetCellDS()`, `DMGetRegionDS()`, `DMSetRegionDS()`
5869: @*/
5870: PetscErrorCode DMCreateDS(DM dm)
5871: {
5872:   MPI_Comm  comm;
5873:   PetscDS   dsDef;
5874:   DMLabel  *labelSet;
5875:   PetscInt  dE, Nf = dm->Nf, f, s, Nl, l, Ndef, k;
5876:   PetscBool doSetup = PETSC_TRUE, flg;

5878:   PetscFunctionBegin;
5880:   if (!dm->fields) PetscFunctionReturn(PETSC_SUCCESS);
5881:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
5882:   PetscCall(DMGetCoordinateDim(dm, &dE));
5883:   /* Determine how many regions we have */
5884:   PetscCall(PetscMalloc1(Nf, &labelSet));
5885:   Nl   = 0;
5886:   Ndef = 0;
5887:   for (f = 0; f < Nf; ++f) {
5888:     DMLabel  label = dm->fields[f].label;
5889:     PetscInt l;

5891: #ifdef PETSC_HAVE_LIBCEED
5892:     /* Move CEED context to discretizations */
5893:     {
5894:       PetscClassId id;

5896:       PetscCall(PetscObjectGetClassId(dm->fields[f].disc, &id));
5897:       if (id == PETSCFE_CLASSID) {
5898:         Ceed ceed;

5900:         PetscCall(DMGetCeed(dm, &ceed));
5901:         PetscCall(PetscFESetCeed((PetscFE)dm->fields[f].disc, ceed));
5902:       }
5903:     }
5904: #endif
5905:     if (!label) {
5906:       ++Ndef;
5907:       continue;
5908:     }
5909:     for (l = 0; l < Nl; ++l)
5910:       if (label == labelSet[l]) break;
5911:     if (l < Nl) continue;
5912:     labelSet[Nl++] = label;
5913:   }
5914:   /* Create default DS if there are no labels to intersect with */
5915:   PetscCall(DMGetRegionDS(dm, NULL, NULL, &dsDef, NULL));
5916:   if (!dsDef && Ndef && !Nl) {
5917:     IS        fields;
5918:     PetscInt *fld, nf;

5920:     for (f = 0, nf = 0; f < Nf; ++f)
5921:       if (!dm->fields[f].label) ++nf;
5922:     PetscCheck(nf, comm, PETSC_ERR_PLIB, "All fields have labels, but we are trying to create a default DS");
5923:     PetscCall(PetscMalloc1(nf, &fld));
5924:     for (f = 0, nf = 0; f < Nf; ++f)
5925:       if (!dm->fields[f].label) fld[nf++] = f;
5926:     PetscCall(ISCreate(PETSC_COMM_SELF, &fields));
5927:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fields, "dm_fields_"));
5928:     PetscCall(ISSetType(fields, ISGENERAL));
5929:     PetscCall(ISGeneralSetIndices(fields, nf, fld, PETSC_OWN_POINTER));

5931:     PetscCall(PetscDSCreate(PETSC_COMM_SELF, &dsDef));
5932:     PetscCall(DMSetRegionDS(dm, NULL, fields, dsDef, NULL));
5933:     PetscCall(PetscDSDestroy(&dsDef));
5934:     PetscCall(ISDestroy(&fields));
5935:   }
5936:   PetscCall(DMGetRegionDS(dm, NULL, NULL, &dsDef, NULL));
5937:   if (dsDef) PetscCall(PetscDSSetCoordinateDimension(dsDef, dE));
5938:   /* Intersect labels with default fields */
5939:   if (Ndef && Nl) {
5940:     DM              plex;
5941:     DMLabel         cellLabel;
5942:     IS              fieldIS, allcellIS, defcellIS = NULL;
5943:     PetscInt       *fields;
5944:     const PetscInt *cells;
5945:     PetscInt        depth, nf = 0, n, c;

5947:     PetscCall(DMConvert(dm, DMPLEX, &plex));
5948:     PetscCall(DMPlexGetDepth(plex, &depth));
5949:     PetscCall(DMGetStratumIS(plex, "dim", depth, &allcellIS));
5950:     if (!allcellIS) PetscCall(DMGetStratumIS(plex, "depth", depth, &allcellIS));
5951:     /* TODO This looks like it only works for one label */
5952:     for (l = 0; l < Nl; ++l) {
5953:       DMLabel label = labelSet[l];
5954:       IS      pointIS;

5956:       PetscCall(ISDestroy(&defcellIS));
5957:       PetscCall(DMLabelGetStratumIS(label, 1, &pointIS));
5958:       PetscCall(ISDifference(allcellIS, pointIS, &defcellIS));
5959:       PetscCall(ISDestroy(&pointIS));
5960:     }
5961:     PetscCall(ISDestroy(&allcellIS));

5963:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "defaultCells", &cellLabel));
5964:     PetscCall(ISGetLocalSize(defcellIS, &n));
5965:     PetscCall(ISGetIndices(defcellIS, &cells));
5966:     for (c = 0; c < n; ++c) PetscCall(DMLabelSetValue(cellLabel, cells[c], 1));
5967:     PetscCall(ISRestoreIndices(defcellIS, &cells));
5968:     PetscCall(ISDestroy(&defcellIS));
5969:     PetscCall(DMPlexLabelComplete(plex, cellLabel));

5971:     PetscCall(PetscMalloc1(Ndef, &fields));
5972:     for (f = 0; f < Nf; ++f)
5973:       if (!dm->fields[f].label) fields[nf++] = f;
5974:     PetscCall(ISCreate(PETSC_COMM_SELF, &fieldIS));
5975:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fieldIS, "dm_fields_"));
5976:     PetscCall(ISSetType(fieldIS, ISGENERAL));
5977:     PetscCall(ISGeneralSetIndices(fieldIS, nf, fields, PETSC_OWN_POINTER));

5979:     PetscCall(PetscDSCreate(PETSC_COMM_SELF, &dsDef));
5980:     PetscCall(DMSetRegionDS(dm, cellLabel, fieldIS, dsDef, NULL));
5981:     PetscCall(PetscDSSetCoordinateDimension(dsDef, dE));
5982:     PetscCall(DMLabelDestroy(&cellLabel));
5983:     PetscCall(PetscDSDestroy(&dsDef));
5984:     PetscCall(ISDestroy(&fieldIS));
5985:     PetscCall(DMDestroy(&plex));
5986:   }
5987:   /* Create label DSes
5988:      - WE ONLY SUPPORT IDENTICAL OR DISJOINT LABELS
5989:   */
5990:   /* TODO Should check that labels are disjoint */
5991:   for (l = 0; l < Nl; ++l) {
5992:     DMLabel   label = labelSet[l];
5993:     PetscDS   ds, dsIn = NULL;
5994:     IS        fields;
5995:     PetscInt *fld, nf;

5997:     PetscCall(PetscDSCreate(PETSC_COMM_SELF, &ds));
5998:     for (f = 0, nf = 0; f < Nf; ++f)
5999:       if (label == dm->fields[f].label || !dm->fields[f].label) ++nf;
6000:     PetscCall(PetscMalloc1(nf, &fld));
6001:     for (f = 0, nf = 0; f < Nf; ++f)
6002:       if (label == dm->fields[f].label || !dm->fields[f].label) fld[nf++] = f;
6003:     PetscCall(ISCreate(PETSC_COMM_SELF, &fields));
6004:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fields, "dm_fields_"));
6005:     PetscCall(ISSetType(fields, ISGENERAL));
6006:     PetscCall(ISGeneralSetIndices(fields, nf, fld, PETSC_OWN_POINTER));
6007:     PetscCall(PetscDSSetCoordinateDimension(ds, dE));
6008:     {
6009:       DMPolytopeType ct;
6010:       PetscInt       lStart, lEnd;
6011:       PetscBool      isCohesiveLocal = PETSC_FALSE, isCohesive;

6013:       PetscCall(DMLabelGetBounds(label, &lStart, &lEnd));
6014:       if (lStart >= 0) {
6015:         PetscCall(DMPlexGetCellType(dm, lStart, &ct));
6016:         switch (ct) {
6017:         case DM_POLYTOPE_POINT_PRISM_TENSOR:
6018:         case DM_POLYTOPE_SEG_PRISM_TENSOR:
6019:         case DM_POLYTOPE_TRI_PRISM_TENSOR:
6020:         case DM_POLYTOPE_QUAD_PRISM_TENSOR:
6021:           isCohesiveLocal = PETSC_TRUE;
6022:           break;
6023:         default:
6024:           break;
6025:         }
6026:       }
6027:       PetscCall(MPIU_Allreduce(&isCohesiveLocal, &isCohesive, 1, MPIU_BOOL, MPI_LOR, comm));
6028:       if (isCohesive) {
6029:         PetscCall(PetscDSCreate(PETSC_COMM_SELF, &dsIn));
6030:         PetscCall(PetscDSSetCoordinateDimension(dsIn, dE));
6031:       }
6032:       for (f = 0, nf = 0; f < Nf; ++f) {
6033:         if (label == dm->fields[f].label || !dm->fields[f].label) {
6034:           if (label == dm->fields[f].label) {
6035:             PetscCall(PetscDSSetDiscretization(ds, nf, NULL));
6036:             PetscCall(PetscDSSetCohesive(ds, nf, isCohesive));
6037:             if (dsIn) {
6038:               PetscCall(PetscDSSetDiscretization(dsIn, nf, NULL));
6039:               PetscCall(PetscDSSetCohesive(dsIn, nf, isCohesive));
6040:             }
6041:           }
6042:           ++nf;
6043:         }
6044:       }
6045:     }
6046:     PetscCall(DMSetRegionDS(dm, label, fields, ds, dsIn));
6047:     PetscCall(ISDestroy(&fields));
6048:     PetscCall(PetscDSDestroy(&ds));
6049:     PetscCall(PetscDSDestroy(&dsIn));
6050:   }
6051:   PetscCall(PetscFree(labelSet));
6052:   /* Set fields in DSes */
6053:   for (s = 0; s < dm->Nds; ++s) {
6054:     PetscDS         ds     = dm->probs[s].ds;
6055:     PetscDS         dsIn   = dm->probs[s].dsIn;
6056:     IS              fields = dm->probs[s].fields;
6057:     const PetscInt *fld;
6058:     PetscInt        nf, dsnf;
6059:     PetscBool       isCohesive;

6061:     PetscCall(PetscDSGetNumFields(ds, &dsnf));
6062:     PetscCall(PetscDSIsCohesive(ds, &isCohesive));
6063:     PetscCall(ISGetLocalSize(fields, &nf));
6064:     PetscCall(ISGetIndices(fields, &fld));
6065:     for (f = 0; f < nf; ++f) {
6066:       PetscObject  disc = dm->fields[fld[f]].disc;
6067:       PetscBool    isCohesiveField;
6068:       PetscClassId id;

6070:       /* Handle DS with no fields */
6071:       if (dsnf) PetscCall(PetscDSGetCohesive(ds, f, &isCohesiveField));
6072:       /* If this is a cohesive cell, then regular fields need the lower dimensional discretization */
6073:       if (isCohesive) {
6074:         if (!isCohesiveField) {
6075:           PetscObject bdDisc;

6077:           PetscCall(PetscFEGetHeightSubspace((PetscFE)disc, 1, (PetscFE *)&bdDisc));
6078:           PetscCall(PetscDSSetDiscretization(ds, f, bdDisc));
6079:           PetscCall(PetscDSSetDiscretization(dsIn, f, disc));
6080:         } else {
6081:           PetscCall(PetscDSSetDiscretization(ds, f, disc));
6082:           PetscCall(PetscDSSetDiscretization(dsIn, f, disc));
6083:         }
6084:       } else {
6085:         PetscCall(PetscDSSetDiscretization(ds, f, disc));
6086:       }
6087:       /* We allow people to have placeholder fields and construct the Section by hand */
6088:       PetscCall(PetscObjectGetClassId(disc, &id));
6089:       if ((id != PETSCFE_CLASSID) && (id != PETSCFV_CLASSID)) doSetup = PETSC_FALSE;
6090:     }
6091:     PetscCall(ISRestoreIndices(fields, &fld));
6092:   }
6093:   /* Allow k-jet tabulation */
6094:   PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)dm)->prefix, "-dm_ds_jet_degree", &k, &flg));
6095:   if (flg) {
6096:     for (s = 0; s < dm->Nds; ++s) {
6097:       PetscDS  ds   = dm->probs[s].ds;
6098:       PetscDS  dsIn = dm->probs[s].dsIn;
6099:       PetscInt Nf, f;

6101:       PetscCall(PetscDSGetNumFields(ds, &Nf));
6102:       for (f = 0; f < Nf; ++f) {
6103:         PetscCall(PetscDSSetJetDegree(ds, f, k));
6104:         if (dsIn) PetscCall(PetscDSSetJetDegree(dsIn, f, k));
6105:       }
6106:     }
6107:   }
6108:   /* Setup DSes */
6109:   if (doSetup) {
6110:     for (s = 0; s < dm->Nds; ++s) {
6111:       if (dm->setfromoptionscalled) {
6112:         PetscCall(PetscDSSetFromOptions(dm->probs[s].ds));
6113:         if (dm->probs[s].dsIn) PetscCall(PetscDSSetFromOptions(dm->probs[s].dsIn));
6114:       }
6115:       PetscCall(PetscDSSetUp(dm->probs[s].ds));
6116:       if (dm->probs[s].dsIn) PetscCall(PetscDSSetUp(dm->probs[s].dsIn));
6117:     }
6118:   }
6119:   PetscFunctionReturn(PETSC_SUCCESS);
6120: }

6122: /*@
6123:   DMUseTensorOrder - Use a tensor product closure ordering for the default section

6125:   Input Parameters:
6126: + dm     - The DM
6127: - tensor - Flag for tensor order

6129:   Level: developer

6131: .seealso: `DMPlexSetClosurePermutationTensor()`, `PetscSectionResetClosurePermutation()`
6132: @*/
6133: PetscErrorCode DMUseTensorOrder(DM dm, PetscBool tensor)
6134: {
6135:   PetscInt  Nf;
6136:   PetscBool reorder = PETSC_TRUE, isPlex;

6138:   PetscFunctionBegin;
6139:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex));
6140:   PetscCall(DMGetNumFields(dm, &Nf));
6141:   for (PetscInt f = 0; f < Nf; ++f) {
6142:     PetscObject  obj;
6143:     PetscClassId id;

6145:     PetscCall(DMGetField(dm, f, NULL, &obj));
6146:     PetscCall(PetscObjectGetClassId(obj, &id));
6147:     if (id == PETSCFE_CLASSID) {
6148:       PetscSpace sp;
6149:       PetscBool  tensor;

6151:       PetscCall(PetscFEGetBasisSpace((PetscFE)obj, &sp));
6152:       PetscCall(PetscSpacePolynomialGetTensor(sp, &tensor));
6153:       reorder = reorder && tensor ? PETSC_TRUE : PETSC_FALSE;
6154:     } else reorder = PETSC_FALSE;
6155:   }
6156:   if (tensor) {
6157:     if (reorder && isPlex) PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
6158:   } else {
6159:     PetscSection s;

6161:     PetscCall(DMGetLocalSection(dm, &s));
6162:     if (s) PetscCall(PetscSectionResetClosurePermutation(s));
6163:   }
6164:   PetscFunctionReturn(PETSC_SUCCESS);
6165: }

6167: /*@
6168:   DMComputeExactSolution - Compute the exact solution for a given `DM`, using the `PetscDS` information.

6170:   Collective

6172:   Input Parameters:
6173: + dm   - The `DM`
6174: - time - The time

6176:   Output Parameters:
6177: + u   - The vector will be filled with exact solution values, or `NULL`
6178: - u_t - The vector will be filled with the time derivative of exact solution values, or `NULL`

6180:   Level: developer

6182:   Note:
6183:   The user must call `PetscDSSetExactSolution()` before using this routine

6185: .seealso: [](ch_dmbase), `DM`, `PetscDSSetExactSolution()`
6186: @*/
6187: PetscErrorCode DMComputeExactSolution(DM dm, PetscReal time, Vec u, Vec u_t)
6188: {
6189:   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
6190:   void   **ectxs;
6191:   Vec      locu, locu_t;
6192:   PetscInt Nf, Nds, s;

6194:   PetscFunctionBegin;
6196:   if (u) {
6198:     PetscCall(DMGetLocalVector(dm, &locu));
6199:     PetscCall(VecSet(locu, 0.));
6200:   }
6201:   if (u_t) {
6203:     PetscCall(DMGetLocalVector(dm, &locu_t));
6204:     PetscCall(VecSet(locu_t, 0.));
6205:   }
6206:   PetscCall(DMGetNumFields(dm, &Nf));
6207:   PetscCall(PetscMalloc2(Nf, &exacts, Nf, &ectxs));
6208:   PetscCall(DMGetNumDS(dm, &Nds));
6209:   for (s = 0; s < Nds; ++s) {
6210:     PetscDS         ds;
6211:     DMLabel         label;
6212:     IS              fieldIS;
6213:     const PetscInt *fields, id = 1;
6214:     PetscInt        dsNf, f;

6216:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL));
6217:     PetscCall(PetscDSGetNumFields(ds, &dsNf));
6218:     PetscCall(ISGetIndices(fieldIS, &fields));
6219:     PetscCall(PetscArrayzero(exacts, Nf));
6220:     PetscCall(PetscArrayzero(ectxs, Nf));
6221:     if (u) {
6222:       for (f = 0; f < dsNf; ++f) PetscCall(PetscDSGetExactSolution(ds, fields[f], &exacts[fields[f]], &ectxs[fields[f]]));
6223:       if (label) PetscCall(DMProjectFunctionLabelLocal(dm, time, label, 1, &id, 0, NULL, exacts, ectxs, INSERT_ALL_VALUES, locu));
6224:       else PetscCall(DMProjectFunctionLocal(dm, time, exacts, ectxs, INSERT_ALL_VALUES, locu));
6225:     }
6226:     if (u_t) {
6227:       PetscCall(PetscArrayzero(exacts, Nf));
6228:       PetscCall(PetscArrayzero(ectxs, Nf));
6229:       for (f = 0; f < dsNf; ++f) PetscCall(PetscDSGetExactSolutionTimeDerivative(ds, fields[f], &exacts[fields[f]], &ectxs[fields[f]]));
6230:       if (label) PetscCall(DMProjectFunctionLabelLocal(dm, time, label, 1, &id, 0, NULL, exacts, ectxs, INSERT_ALL_VALUES, locu_t));
6231:       else PetscCall(DMProjectFunctionLocal(dm, time, exacts, ectxs, INSERT_ALL_VALUES, locu_t));
6232:     }
6233:     PetscCall(ISRestoreIndices(fieldIS, &fields));
6234:   }
6235:   if (u) {
6236:     PetscCall(PetscObjectSetName((PetscObject)u, "Exact Solution"));
6237:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)u, "exact_"));
6238:   }
6239:   if (u_t) {
6240:     PetscCall(PetscObjectSetName((PetscObject)u, "Exact Solution Time Derivative"));
6241:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)u_t, "exact_t_"));
6242:   }
6243:   PetscCall(PetscFree2(exacts, ectxs));
6244:   if (u) {
6245:     PetscCall(DMLocalToGlobalBegin(dm, locu, INSERT_ALL_VALUES, u));
6246:     PetscCall(DMLocalToGlobalEnd(dm, locu, INSERT_ALL_VALUES, u));
6247:     PetscCall(DMRestoreLocalVector(dm, &locu));
6248:   }
6249:   if (u_t) {
6250:     PetscCall(DMLocalToGlobalBegin(dm, locu_t, INSERT_ALL_VALUES, u_t));
6251:     PetscCall(DMLocalToGlobalEnd(dm, locu_t, INSERT_ALL_VALUES, u_t));
6252:     PetscCall(DMRestoreLocalVector(dm, &locu_t));
6253:   }
6254:   PetscFunctionReturn(PETSC_SUCCESS);
6255: }

6257: static PetscErrorCode DMTransferDS_Internal(DM dm, DMLabel label, IS fields, PetscDS ds, PetscDS dsIn)
6258: {
6259:   PetscDS dsNew, dsInNew = NULL;

6261:   PetscFunctionBegin;
6262:   PetscCall(PetscDSCreate(PetscObjectComm((PetscObject)ds), &dsNew));
6263:   PetscCall(PetscDSCopy(ds, dm, dsNew));
6264:   if (dsIn) {
6265:     PetscCall(PetscDSCreate(PetscObjectComm((PetscObject)dsIn), &dsInNew));
6266:     PetscCall(PetscDSCopy(dsIn, dm, dsInNew));
6267:   }
6268:   PetscCall(DMSetRegionDS(dm, label, fields, dsNew, dsInNew));
6269:   PetscCall(PetscDSDestroy(&dsNew));
6270:   PetscCall(PetscDSDestroy(&dsInNew));
6271:   PetscFunctionReturn(PETSC_SUCCESS);
6272: }

6274: /*@
6275:   DMCopyDS - Copy the discrete systems for the `DM` into another `DM`

6277:   Collective

6279:   Input Parameter:
6280: . dm - The `DM`

6282:   Output Parameter:
6283: . newdm - The `DM`

6285:   Level: advanced

6287: .seealso: [](ch_dmbase), `DM`, `DMCopyFields()`, `DMAddField()`, `DMGetDS()`, `DMGetCellDS()`, `DMGetRegionDS()`, `DMSetRegionDS()`
6288: @*/
6289: PetscErrorCode DMCopyDS(DM dm, DM newdm)
6290: {
6291:   PetscInt Nds, s;

6293:   PetscFunctionBegin;
6294:   if (dm == newdm) PetscFunctionReturn(PETSC_SUCCESS);
6295:   PetscCall(DMGetNumDS(dm, &Nds));
6296:   PetscCall(DMClearDS(newdm));
6297:   for (s = 0; s < Nds; ++s) {
6298:     DMLabel  label;
6299:     IS       fields;
6300:     PetscDS  ds, dsIn, newds;
6301:     PetscInt Nbd, bd;

6303:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fields, &ds, &dsIn));
6304:     /* TODO: We need to change all keys from labels in the old DM to labels in the new DM */
6305:     PetscCall(DMTransferDS_Internal(newdm, label, fields, ds, dsIn));
6306:     /* Complete new labels in the new DS */
6307:     PetscCall(DMGetRegionDS(newdm, label, NULL, &newds, NULL));
6308:     PetscCall(PetscDSGetNumBoundary(newds, &Nbd));
6309:     for (bd = 0; bd < Nbd; ++bd) {
6310:       PetscWeakForm wf;
6311:       DMLabel       label;
6312:       PetscInt      field;

6314:       PetscCall(PetscDSGetBoundary(newds, bd, &wf, NULL, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL));
6315:       PetscCall(PetscWeakFormReplaceLabel(wf, label));
6316:     }
6317:   }
6318:   PetscCall(DMCompleteBCLabels_Internal(newdm));
6319:   PetscFunctionReturn(PETSC_SUCCESS);
6320: }

6322: /*@
6323:   DMCopyDisc - Copy the fields and discrete systems for the `DM` into another `DM`

6325:   Collective

6327:   Input Parameter:
6328: . dm - The `DM`

6330:   Output Parameter:
6331: . newdm - The `DM`

6333:   Level: advanced

6335:   Developer Note:
6336:   Really ugly name, nothing in PETSc is called a `Disc` plus it is an ugly abbreviation

6338: .seealso: [](ch_dmbase), `DM`, `DMCopyFields()`, `DMCopyDS()`
6339: @*/
6340: PetscErrorCode DMCopyDisc(DM dm, DM newdm)
6341: {
6342:   PetscFunctionBegin;
6343:   PetscCall(DMCopyFields(dm, newdm));
6344:   PetscCall(DMCopyDS(dm, newdm));
6345:   PetscFunctionReturn(PETSC_SUCCESS);
6346: }

6348: /*@
6349:   DMGetDimension - Return the topological dimension of the `DM`

6351:   Not Collective

6353:   Input Parameter:
6354: . dm - The `DM`

6356:   Output Parameter:
6357: . dim - The topological dimension

6359:   Level: beginner

6361: .seealso: [](ch_dmbase), `DM`, `DMSetDimension()`, `DMCreate()`
6362: @*/
6363: PetscErrorCode DMGetDimension(DM dm, PetscInt *dim)
6364: {
6365:   PetscFunctionBegin;
6367:   PetscAssertPointer(dim, 2);
6368:   *dim = dm->dim;
6369:   PetscFunctionReturn(PETSC_SUCCESS);
6370: }

6372: /*@
6373:   DMSetDimension - Set the topological dimension of the `DM`

6375:   Collective

6377:   Input Parameters:
6378: + dm  - The `DM`
6379: - dim - The topological dimension

6381:   Level: beginner

6383: .seealso: [](ch_dmbase), `DM`, `DMGetDimension()`, `DMCreate()`
6384: @*/
6385: PetscErrorCode DMSetDimension(DM dm, PetscInt dim)
6386: {
6387:   PetscDS  ds;
6388:   PetscInt Nds, n;

6390:   PetscFunctionBegin;
6393:   dm->dim = dim;
6394:   if (dm->dim >= 0) {
6395:     PetscCall(DMGetNumDS(dm, &Nds));
6396:     for (n = 0; n < Nds; ++n) {
6397:       PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds, NULL));
6398:       if (ds->dimEmbed < 0) PetscCall(PetscDSSetCoordinateDimension(ds, dim));
6399:     }
6400:   }
6401:   PetscFunctionReturn(PETSC_SUCCESS);
6402: }

6404: /*@
6405:   DMGetDimPoints - Get the half-open interval for all points of a given dimension

6407:   Collective

6409:   Input Parameters:
6410: + dm  - the `DM`
6411: - dim - the dimension

6413:   Output Parameters:
6414: + pStart - The first point of the given dimension
6415: - pEnd   - The first point following points of the given dimension

6417:   Level: intermediate

6419:   Note:
6420:   The points are vertices in the Hasse diagram encoding the topology. This is explained in
6421:   https://arxiv.org/abs/0908.4427. If no points exist of this dimension in the storage scheme,
6422:   then the interval is empty.

6424: .seealso: [](ch_dmbase), `DM`, `DMPLEX`, `DMPlexGetDepthStratum()`, `DMPlexGetHeightStratum()`
6425: @*/
6426: PetscErrorCode DMGetDimPoints(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
6427: {
6428:   PetscInt d;

6430:   PetscFunctionBegin;
6432:   PetscCall(DMGetDimension(dm, &d));
6433:   PetscCheck((dim >= 0) && (dim <= d), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim);
6434:   PetscUseTypeMethod(dm, getdimpoints, dim, pStart, pEnd);
6435:   PetscFunctionReturn(PETSC_SUCCESS);
6436: }

6438: /*@
6439:   DMGetOutputDM - Retrieve the `DM` associated with the layout for output

6441:   Collective

6443:   Input Parameter:
6444: . dm - The original `DM`

6446:   Output Parameter:
6447: . odm - The `DM` which provides the layout for output

6449:   Level: intermediate

6451:   Note:
6452:   In some situations the vector obtained with `DMCreateGlobalVector()` excludes points for degrees of freedom that are associated with fixed (Dirichelet) boundary
6453:   conditions since the algebraic solver does not solve for those variables. The output `DM` includes these excluded points and its global vector contains the
6454:   locations for those dof so that they can be output to a file or other viewer along with the unconstrained dof.

6456: .seealso: [](ch_dmbase), `DM`, `VecView()`, `DMGetGlobalSection()`, `DMCreateGlobalVector()`, `PetscSectionHasConstraints()`, `DMSetGlobalSection()`
6457: @*/
6458: PetscErrorCode DMGetOutputDM(DM dm, DM *odm)
6459: {
6460:   PetscSection section;
6461:   IS           perm;
6462:   PetscBool    hasConstraints, newDM, gnewDM;

6464:   PetscFunctionBegin;
6466:   PetscAssertPointer(odm, 2);
6467:   PetscCall(DMGetLocalSection(dm, &section));
6468:   PetscCall(PetscSectionHasConstraints(section, &hasConstraints));
6469:   PetscCall(PetscSectionGetPermutation(section, &perm));
6470:   newDM = hasConstraints || perm ? PETSC_TRUE : PETSC_FALSE;
6471:   PetscCall(MPIU_Allreduce(&newDM, &gnewDM, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)dm)));
6472:   if (!gnewDM) {
6473:     *odm = dm;
6474:     PetscFunctionReturn(PETSC_SUCCESS);
6475:   }
6476:   if (!dm->dmBC) {
6477:     PetscSection newSection, gsection;
6478:     PetscSF      sf;
6479:     PetscBool    usePerm = dm->ignorePermOutput ? PETSC_FALSE : PETSC_TRUE;

6481:     PetscCall(DMClone(dm, &dm->dmBC));
6482:     PetscCall(DMCopyDisc(dm, dm->dmBC));
6483:     PetscCall(PetscSectionClone(section, &newSection));
6484:     PetscCall(DMSetLocalSection(dm->dmBC, newSection));
6485:     PetscCall(PetscSectionDestroy(&newSection));
6486:     PetscCall(DMGetPointSF(dm->dmBC, &sf));
6487:     PetscCall(PetscSectionCreateGlobalSection(section, sf, usePerm, PETSC_TRUE, PETSC_FALSE, &gsection));
6488:     PetscCall(DMSetGlobalSection(dm->dmBC, gsection));
6489:     PetscCall(PetscSectionDestroy(&gsection));
6490:   }
6491:   *odm = dm->dmBC;
6492:   PetscFunctionReturn(PETSC_SUCCESS);
6493: }

6495: /*@
6496:   DMGetOutputSequenceNumber - Retrieve the sequence number/value for output

6498:   Input Parameter:
6499: . dm - The original `DM`

6501:   Output Parameters:
6502: + num - The output sequence number
6503: - val - The output sequence value

6505:   Level: intermediate

6507:   Note:
6508:   This is intended for output that should appear in sequence, for instance
6509:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6511:   Developer Note:
6512:   The `DM` serves as a convenient place to store the current iteration value. The iteration is not
6513:   not directly related to the `DM`.

6515: .seealso: [](ch_dmbase), `DM`, `VecView()`
6516: @*/
6517: PetscErrorCode DMGetOutputSequenceNumber(DM dm, PetscInt *num, PetscReal *val)
6518: {
6519:   PetscFunctionBegin;
6521:   if (num) {
6522:     PetscAssertPointer(num, 2);
6523:     *num = dm->outputSequenceNum;
6524:   }
6525:   if (val) {
6526:     PetscAssertPointer(val, 3);
6527:     *val = dm->outputSequenceVal;
6528:   }
6529:   PetscFunctionReturn(PETSC_SUCCESS);
6530: }

6532: /*@
6533:   DMSetOutputSequenceNumber - Set the sequence number/value for output

6535:   Input Parameters:
6536: + dm  - The original `DM`
6537: . num - The output sequence number
6538: - val - The output sequence value

6540:   Level: intermediate

6542:   Note:
6543:   This is intended for output that should appear in sequence, for instance
6544:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6546: .seealso: [](ch_dmbase), `DM`, `VecView()`
6547: @*/
6548: PetscErrorCode DMSetOutputSequenceNumber(DM dm, PetscInt num, PetscReal val)
6549: {
6550:   PetscFunctionBegin;
6552:   dm->outputSequenceNum = num;
6553:   dm->outputSequenceVal = val;
6554:   PetscFunctionReturn(PETSC_SUCCESS);
6555: }

6557: /*@C
6558:   DMOutputSequenceLoad - Retrieve the sequence value from a `PetscViewer`

6560:   Input Parameters:
6561: + dm     - The original `DM`
6562: . viewer - The viewer to get it from
6563: . name   - The sequence name
6564: - num    - The output sequence number

6566:   Output Parameter:
6567: . val - The output sequence value

6569:   Level: intermediate

6571:   Note:
6572:   This is intended for output that should appear in sequence, for instance
6573:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6575:   Developer Note:
6576:   It is unclear at the user API level why a `DM` is needed as input

6578: .seealso: [](ch_dmbase), `DM`, `DMGetOutputSequenceNumber()`, `DMSetOutputSequenceNumber()`, `VecView()`
6579: @*/
6580: PetscErrorCode DMOutputSequenceLoad(DM dm, PetscViewer viewer, const char *name, PetscInt num, PetscReal *val)
6581: {
6582:   PetscBool ishdf5;

6584:   PetscFunctionBegin;
6587:   PetscAssertPointer(val, 5);
6588:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
6589:   if (ishdf5) {
6590: #if defined(PETSC_HAVE_HDF5)
6591:     PetscScalar value;

6593:     PetscCall(DMSequenceLoad_HDF5_Internal(dm, name, num, &value, viewer));
6594:     *val = PetscRealPart(value);
6595: #endif
6596:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerHDF5Open()");
6597:   PetscFunctionReturn(PETSC_SUCCESS);
6598: }

6600: /*@
6601:   DMGetUseNatural - Get the flag for creating a mapping to the natural order when a `DM` is (re)distributed in parallel

6603:   Not Collective

6605:   Input Parameter:
6606: . dm - The `DM`

6608:   Output Parameter:
6609: . useNatural - `PETSC_TRUE` to build the mapping to a natural order during distribution

6611:   Level: beginner

6613: .seealso: [](ch_dmbase), `DM`, `DMSetUseNatural()`, `DMCreate()`
6614: @*/
6615: PetscErrorCode DMGetUseNatural(DM dm, PetscBool *useNatural)
6616: {
6617:   PetscFunctionBegin;
6619:   PetscAssertPointer(useNatural, 2);
6620:   *useNatural = dm->useNatural;
6621:   PetscFunctionReturn(PETSC_SUCCESS);
6622: }

6624: /*@
6625:   DMSetUseNatural - Set the flag for creating a mapping to the natural order when a `DM` is (re)distributed in parallel

6627:   Collective

6629:   Input Parameters:
6630: + dm         - The `DM`
6631: - useNatural - `PETSC_TRUE` to build the mapping to a natural order during distribution

6633:   Level: beginner

6635:   Note:
6636:   This also causes the map to be build after `DMCreateSubDM()` and `DMCreateSuperDM()`

6638: .seealso: [](ch_dmbase), `DM`, `DMGetUseNatural()`, `DMCreate()`, `DMPlexDistribute()`, `DMCreateSubDM()`, `DMCreateSuperDM()`
6639: @*/
6640: PetscErrorCode DMSetUseNatural(DM dm, PetscBool useNatural)
6641: {
6642:   PetscFunctionBegin;
6645:   dm->useNatural = useNatural;
6646:   PetscFunctionReturn(PETSC_SUCCESS);
6647: }

6649: /*@C
6650:   DMCreateLabel - Create a label of the given name if it does not already exist in the `DM`

6652:   Not Collective

6654:   Input Parameters:
6655: + dm   - The `DM` object
6656: - name - The label name

6658:   Level: intermediate

6660: .seealso: [](ch_dmbase), `DM`, `DMLabelCreate()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6661: @*/
6662: PetscErrorCode DMCreateLabel(DM dm, const char name[])
6663: {
6664:   PetscBool flg;
6665:   DMLabel   label;

6667:   PetscFunctionBegin;
6669:   PetscAssertPointer(name, 2);
6670:   PetscCall(DMHasLabel(dm, name, &flg));
6671:   if (!flg) {
6672:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, &label));
6673:     PetscCall(DMAddLabel(dm, label));
6674:     PetscCall(DMLabelDestroy(&label));
6675:   }
6676:   PetscFunctionReturn(PETSC_SUCCESS);
6677: }

6679: /*@C
6680:   DMCreateLabelAtIndex - Create a label of the given name at the given index. If it already exists in the `DM`, move it to this index.

6682:   Not Collective

6684:   Input Parameters:
6685: + dm   - The `DM` object
6686: . l    - The index for the label
6687: - name - The label name

6689:   Level: intermediate

6691: .seealso: [](ch_dmbase), `DM`, `DMCreateLabel()`, `DMLabelCreate()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6692: @*/
6693: PetscErrorCode DMCreateLabelAtIndex(DM dm, PetscInt l, const char name[])
6694: {
6695:   DMLabelLink orig, prev = NULL;
6696:   DMLabel     label;
6697:   PetscInt    Nl, m;
6698:   PetscBool   flg, match;
6699:   const char *lname;

6701:   PetscFunctionBegin;
6703:   PetscAssertPointer(name, 3);
6704:   PetscCall(DMHasLabel(dm, name, &flg));
6705:   if (!flg) {
6706:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, &label));
6707:     PetscCall(DMAddLabel(dm, label));
6708:     PetscCall(DMLabelDestroy(&label));
6709:   }
6710:   PetscCall(DMGetNumLabels(dm, &Nl));
6711:   PetscCheck(l < Nl, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label index %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", l, Nl);
6712:   for (m = 0, orig = dm->labels; m < Nl; ++m, prev = orig, orig = orig->next) {
6713:     PetscCall(PetscObjectGetName((PetscObject)orig->label, &lname));
6714:     PetscCall(PetscStrcmp(name, lname, &match));
6715:     if (match) break;
6716:   }
6717:   if (m == l) PetscFunctionReturn(PETSC_SUCCESS);
6718:   if (!m) dm->labels = orig->next;
6719:   else prev->next = orig->next;
6720:   if (!l) {
6721:     orig->next = dm->labels;
6722:     dm->labels = orig;
6723:   } else {
6724:     for (m = 0, prev = dm->labels; m < l - 1; ++m, prev = prev->next);
6725:     orig->next = prev->next;
6726:     prev->next = orig;
6727:   }
6728:   PetscFunctionReturn(PETSC_SUCCESS);
6729: }

6731: /*@C
6732:   DMGetLabelValue - Get the value in a `DMLabel` for the given point, with -1 as the default

6734:   Not Collective

6736:   Input Parameters:
6737: + dm    - The `DM` object
6738: . name  - The label name
6739: - point - The mesh point

6741:   Output Parameter:
6742: . value - The label value for this point, or -1 if the point is not in the label

6744:   Level: beginner

6746: .seealso: [](ch_dmbase), `DM`, `DMLabelGetValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6747: @*/
6748: PetscErrorCode DMGetLabelValue(DM dm, const char name[], PetscInt point, PetscInt *value)
6749: {
6750:   DMLabel label;

6752:   PetscFunctionBegin;
6754:   PetscAssertPointer(name, 2);
6755:   PetscCall(DMGetLabel(dm, name, &label));
6756:   PetscCheck(label, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No label named %s was found", name);
6757:   PetscCall(DMLabelGetValue(label, point, value));
6758:   PetscFunctionReturn(PETSC_SUCCESS);
6759: }

6761: /*@C
6762:   DMSetLabelValue - Add a point to a `DMLabel` with given value

6764:   Not Collective

6766:   Input Parameters:
6767: + dm    - The `DM` object
6768: . name  - The label name
6769: . point - The mesh point
6770: - value - The label value for this point

6772:   Output Parameter:

6774:   Level: beginner

6776: .seealso: [](ch_dmbase), `DM`, `DMLabelSetValue()`, `DMGetStratumIS()`, `DMClearLabelValue()`
6777: @*/
6778: PetscErrorCode DMSetLabelValue(DM dm, const char name[], PetscInt point, PetscInt value)
6779: {
6780:   DMLabel label;

6782:   PetscFunctionBegin;
6784:   PetscAssertPointer(name, 2);
6785:   PetscCall(DMGetLabel(dm, name, &label));
6786:   if (!label) {
6787:     PetscCall(DMCreateLabel(dm, name));
6788:     PetscCall(DMGetLabel(dm, name, &label));
6789:   }
6790:   PetscCall(DMLabelSetValue(label, point, value));
6791:   PetscFunctionReturn(PETSC_SUCCESS);
6792: }

6794: /*@C
6795:   DMClearLabelValue - Remove a point from a `DMLabel` with given value

6797:   Not Collective

6799:   Input Parameters:
6800: + dm    - The `DM` object
6801: . name  - The label name
6802: . point - The mesh point
6803: - value - The label value for this point

6805:   Level: beginner

6807: .seealso: [](ch_dmbase), `DM`, `DMLabelClearValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6808: @*/
6809: PetscErrorCode DMClearLabelValue(DM dm, const char name[], PetscInt point, PetscInt value)
6810: {
6811:   DMLabel label;

6813:   PetscFunctionBegin;
6815:   PetscAssertPointer(name, 2);
6816:   PetscCall(DMGetLabel(dm, name, &label));
6817:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6818:   PetscCall(DMLabelClearValue(label, point, value));
6819:   PetscFunctionReturn(PETSC_SUCCESS);
6820: }

6822: /*@C
6823:   DMGetLabelSize - Get the value of `DMLabelGetNumValues()` of a `DMLabel` in the `DM`

6825:   Not Collective

6827:   Input Parameters:
6828: + dm   - The `DM` object
6829: - name - The label name

6831:   Output Parameter:
6832: . size - The number of different integer ids, or 0 if the label does not exist

6834:   Level: beginner

6836:   Developer Note:
6837:   This should be renamed to something like `DMGetLabelNumValues()` or removed.

6839: .seealso: [](ch_dmbase), `DM`, `DMLabelGetNumValues()`, `DMSetLabelValue()`, `DMGetLabel()`
6840: @*/
6841: PetscErrorCode DMGetLabelSize(DM dm, const char name[], PetscInt *size)
6842: {
6843:   DMLabel label;

6845:   PetscFunctionBegin;
6847:   PetscAssertPointer(name, 2);
6848:   PetscAssertPointer(size, 3);
6849:   PetscCall(DMGetLabel(dm, name, &label));
6850:   *size = 0;
6851:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6852:   PetscCall(DMLabelGetNumValues(label, size));
6853:   PetscFunctionReturn(PETSC_SUCCESS);
6854: }

6856: /*@C
6857:   DMGetLabelIdIS - Get the `DMLabelGetValueIS()` from a `DMLabel` in the `DM`

6859:   Not Collective

6861:   Input Parameters:
6862: + dm   - The `DM` object
6863: - name - The label name

6865:   Output Parameter:
6866: . ids - The integer ids, or `NULL` if the label does not exist

6868:   Level: beginner

6870: .seealso: [](ch_dmbase), `DM`, `DMLabelGetValueIS()`, `DMGetLabelSize()`
6871: @*/
6872: PetscErrorCode DMGetLabelIdIS(DM dm, const char name[], IS *ids)
6873: {
6874:   DMLabel label;

6876:   PetscFunctionBegin;
6878:   PetscAssertPointer(name, 2);
6879:   PetscAssertPointer(ids, 3);
6880:   PetscCall(DMGetLabel(dm, name, &label));
6881:   *ids = NULL;
6882:   if (label) {
6883:     PetscCall(DMLabelGetValueIS(label, ids));
6884:   } else {
6885:     /* returning an empty IS */
6886:     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_USE_POINTER, ids));
6887:   }
6888:   PetscFunctionReturn(PETSC_SUCCESS);
6889: }

6891: /*@C
6892:   DMGetStratumSize - Get the number of points in a label stratum

6894:   Not Collective

6896:   Input Parameters:
6897: + dm    - The `DM` object
6898: . name  - The label name of the stratum
6899: - value - The stratum value

6901:   Output Parameter:
6902: . size - The number of points, also called the stratum size

6904:   Level: beginner

6906: .seealso: [](ch_dmbase), `DM`, `DMLabelGetStratumSize()`, `DMGetLabelSize()`, `DMGetLabelIds()`
6907: @*/
6908: PetscErrorCode DMGetStratumSize(DM dm, const char name[], PetscInt value, PetscInt *size)
6909: {
6910:   DMLabel label;

6912:   PetscFunctionBegin;
6914:   PetscAssertPointer(name, 2);
6915:   PetscAssertPointer(size, 4);
6916:   PetscCall(DMGetLabel(dm, name, &label));
6917:   *size = 0;
6918:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6919:   PetscCall(DMLabelGetStratumSize(label, value, size));
6920:   PetscFunctionReturn(PETSC_SUCCESS);
6921: }

6923: /*@C
6924:   DMGetStratumIS - Get the points in a label stratum

6926:   Not Collective

6928:   Input Parameters:
6929: + dm    - The `DM` object
6930: . name  - The label name
6931: - value - The stratum value

6933:   Output Parameter:
6934: . points - The stratum points, or `NULL` if the label does not exist or does not have that value

6936:   Level: beginner

6938: .seealso: [](ch_dmbase), `DM`, `DMLabelGetStratumIS()`, `DMGetStratumSize()`
6939: @*/
6940: PetscErrorCode DMGetStratumIS(DM dm, const char name[], PetscInt value, IS *points)
6941: {
6942:   DMLabel label;

6944:   PetscFunctionBegin;
6946:   PetscAssertPointer(name, 2);
6947:   PetscAssertPointer(points, 4);
6948:   PetscCall(DMGetLabel(dm, name, &label));
6949:   *points = NULL;
6950:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6951:   PetscCall(DMLabelGetStratumIS(label, value, points));
6952:   PetscFunctionReturn(PETSC_SUCCESS);
6953: }

6955: /*@C
6956:   DMSetStratumIS - Set the points in a label stratum

6958:   Not Collective

6960:   Input Parameters:
6961: + dm     - The `DM` object
6962: . name   - The label name
6963: . value  - The stratum value
6964: - points - The stratum points

6966:   Level: beginner

6968: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMClearLabelStratum()`, `DMLabelClearStratum()`, `DMLabelSetStratumIS()`, `DMGetStratumSize()`
6969: @*/
6970: PetscErrorCode DMSetStratumIS(DM dm, const char name[], PetscInt value, IS points)
6971: {
6972:   DMLabel label;

6974:   PetscFunctionBegin;
6976:   PetscAssertPointer(name, 2);
6978:   PetscCall(DMGetLabel(dm, name, &label));
6979:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6980:   PetscCall(DMLabelSetStratumIS(label, value, points));
6981:   PetscFunctionReturn(PETSC_SUCCESS);
6982: }

6984: /*@C
6985:   DMClearLabelStratum - Remove all points from a stratum from a `DMLabel`

6987:   Not Collective

6989:   Input Parameters:
6990: + dm    - The `DM` object
6991: . name  - The label name
6992: - value - The label value for this point

6994:   Output Parameter:

6996:   Level: beginner

6998: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMLabelClearStratum()`, `DMSetLabelValue()`, `DMGetStratumIS()`, `DMClearLabelValue()`
6999: @*/
7000: PetscErrorCode DMClearLabelStratum(DM dm, const char name[], PetscInt value)
7001: {
7002:   DMLabel label;

7004:   PetscFunctionBegin;
7006:   PetscAssertPointer(name, 2);
7007:   PetscCall(DMGetLabel(dm, name, &label));
7008:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
7009:   PetscCall(DMLabelClearStratum(label, value));
7010:   PetscFunctionReturn(PETSC_SUCCESS);
7011: }

7013: /*@
7014:   DMGetNumLabels - Return the number of labels defined by on the `DM`

7016:   Not Collective

7018:   Input Parameter:
7019: . dm - The `DM` object

7021:   Output Parameter:
7022: . numLabels - the number of Labels

7024:   Level: intermediate

7026: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetLabelByNum()`, `DMGetLabelName()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7027: @*/
7028: PetscErrorCode DMGetNumLabels(DM dm, PetscInt *numLabels)
7029: {
7030:   DMLabelLink next = dm->labels;
7031:   PetscInt    n    = 0;

7033:   PetscFunctionBegin;
7035:   PetscAssertPointer(numLabels, 2);
7036:   while (next) {
7037:     ++n;
7038:     next = next->next;
7039:   }
7040:   *numLabels = n;
7041:   PetscFunctionReturn(PETSC_SUCCESS);
7042: }

7044: /*@C
7045:   DMGetLabelName - Return the name of nth label

7047:   Not Collective

7049:   Input Parameters:
7050: + dm - The `DM` object
7051: - n  - the label number

7053:   Output Parameter:
7054: . name - the label name

7056:   Level: intermediate

7058:   Developer Note:
7059:   Some of the functions that appropriate on labels using their number have the suffix ByNum, others do not.

7061: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetLabelByNum()`, `DMGetLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7062: @*/
7063: PetscErrorCode DMGetLabelName(DM dm, PetscInt n, const char **name)
7064: {
7065:   DMLabelLink next = dm->labels;
7066:   PetscInt    l    = 0;

7068:   PetscFunctionBegin;
7070:   PetscAssertPointer(name, 3);
7071:   while (next) {
7072:     if (l == n) {
7073:       PetscCall(PetscObjectGetName((PetscObject)next->label, name));
7074:       PetscFunctionReturn(PETSC_SUCCESS);
7075:     }
7076:     ++l;
7077:     next = next->next;
7078:   }
7079:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %" PetscInt_FMT " does not exist in this DM", n);
7080: }

7082: /*@C
7083:   DMHasLabel - Determine whether the `DM` has a label of a given name

7085:   Not Collective

7087:   Input Parameters:
7088: + dm   - The `DM` object
7089: - name - The label name

7091:   Output Parameter:
7092: . hasLabel - `PETSC_TRUE` if the label is present

7094:   Level: intermediate

7096: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetLabel()`, `DMGetLabelByNum()`, `DMCreateLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7097: @*/
7098: PetscErrorCode DMHasLabel(DM dm, const char name[], PetscBool *hasLabel)
7099: {
7100:   DMLabelLink next = dm->labels;
7101:   const char *lname;

7103:   PetscFunctionBegin;
7105:   PetscAssertPointer(name, 2);
7106:   PetscAssertPointer(hasLabel, 3);
7107:   *hasLabel = PETSC_FALSE;
7108:   while (next) {
7109:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7110:     PetscCall(PetscStrcmp(name, lname, hasLabel));
7111:     if (*hasLabel) break;
7112:     next = next->next;
7113:   }
7114:   PetscFunctionReturn(PETSC_SUCCESS);
7115: }

7117: // PetscClangLinter pragma ignore: -fdoc-section-header-unknown
7118: /*@C
7119:   DMGetLabel - Return the label of a given name, or `NULL`, from a `DM`

7121:   Not Collective

7123:   Input Parameters:
7124: + dm   - The `DM` object
7125: - name - The label name

7127:   Output Parameter:
7128: . label - The `DMLabel`, or `NULL` if the label is absent

7130:   Default labels in a `DMPLEX`:
7131: + "depth"       - Holds the depth (co-dimension) of each mesh point
7132: . "celltype"    - Holds the topological type of each cell
7133: . "ghost"       - If the DM is distributed with overlap, this marks the cells and faces in the overlap
7134: . "Cell Sets"   - Mirrors the cell sets defined by GMsh and ExodusII
7135: . "Face Sets"   - Mirrors the face sets defined by GMsh and ExodusII
7136: - "Vertex Sets" - Mirrors the vertex sets defined by GMsh

7138:   Level: intermediate

7140: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMHasLabel()`, `DMGetLabelByNum()`, `DMAddLabel()`, `DMCreateLabel()`, `DMPlexGetDepthLabel()`, `DMPlexGetCellType()`
7141: @*/
7142: PetscErrorCode DMGetLabel(DM dm, const char name[], DMLabel *label)
7143: {
7144:   DMLabelLink next = dm->labels;
7145:   PetscBool   hasLabel;
7146:   const char *lname;

7148:   PetscFunctionBegin;
7150:   PetscAssertPointer(name, 2);
7151:   PetscAssertPointer(label, 3);
7152:   *label = NULL;
7153:   while (next) {
7154:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7155:     PetscCall(PetscStrcmp(name, lname, &hasLabel));
7156:     if (hasLabel) {
7157:       *label = next->label;
7158:       break;
7159:     }
7160:     next = next->next;
7161:   }
7162:   PetscFunctionReturn(PETSC_SUCCESS);
7163: }

7165: /*@C
7166:   DMGetLabelByNum - Return the nth label on a `DM`

7168:   Not Collective

7170:   Input Parameters:
7171: + dm - The `DM` object
7172: - n  - the label number

7174:   Output Parameter:
7175: . label - the label

7177:   Level: intermediate

7179: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMAddLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7180: @*/
7181: PetscErrorCode DMGetLabelByNum(DM dm, PetscInt n, DMLabel *label)
7182: {
7183:   DMLabelLink next = dm->labels;
7184:   PetscInt    l    = 0;

7186:   PetscFunctionBegin;
7188:   PetscAssertPointer(label, 3);
7189:   while (next) {
7190:     if (l == n) {
7191:       *label = next->label;
7192:       PetscFunctionReturn(PETSC_SUCCESS);
7193:     }
7194:     ++l;
7195:     next = next->next;
7196:   }
7197:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %" PetscInt_FMT " does not exist in this DM", n);
7198: }

7200: /*@C
7201:   DMAddLabel - Add the label to this `DM`

7203:   Not Collective

7205:   Input Parameters:
7206: + dm    - The `DM` object
7207: - label - The `DMLabel`

7209:   Level: developer

7211: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7212: @*/
7213: PetscErrorCode DMAddLabel(DM dm, DMLabel label)
7214: {
7215:   DMLabelLink l, *p, tmpLabel;
7216:   PetscBool   hasLabel;
7217:   const char *lname;
7218:   PetscBool   flg;

7220:   PetscFunctionBegin;
7222:   PetscCall(PetscObjectGetName((PetscObject)label, &lname));
7223:   PetscCall(DMHasLabel(dm, lname, &hasLabel));
7224:   PetscCheck(!hasLabel, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %s already exists in this DM", lname);
7225:   PetscCall(PetscCalloc1(1, &tmpLabel));
7226:   tmpLabel->label  = label;
7227:   tmpLabel->output = PETSC_TRUE;
7228:   for (p = &dm->labels; (l = *p); p = &l->next) { }
7229:   *p = tmpLabel;
7230:   PetscCall(PetscObjectReference((PetscObject)label));
7231:   PetscCall(PetscStrcmp(lname, "depth", &flg));
7232:   if (flg) dm->depthLabel = label;
7233:   PetscCall(PetscStrcmp(lname, "celltype", &flg));
7234:   if (flg) dm->celltypeLabel = label;
7235:   PetscFunctionReturn(PETSC_SUCCESS);
7236: }

7238: // PetscClangLinter pragma ignore: -fdoc-section-header-unknown
7239: /*@C
7240:   DMSetLabel - Replaces the label of a given name, or ignores it if the name is not present

7242:   Not Collective

7244:   Input Parameters:
7245: + dm    - The `DM` object
7246: - label - The `DMLabel`, having the same name, to substitute

7248:   Default labels in a `DMPLEX`:
7249: + "depth"       - Holds the depth (co-dimension) of each mesh point
7250: . "celltype"    - Holds the topological type of each cell
7251: . "ghost"       - If the DM is distributed with overlap, this marks the cells and faces in the overlap
7252: . "Cell Sets"   - Mirrors the cell sets defined by GMsh and ExodusII
7253: . "Face Sets"   - Mirrors the face sets defined by GMsh and ExodusII
7254: - "Vertex Sets" - Mirrors the vertex sets defined by GMsh

7256:   Level: intermediate

7258: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMPlexGetDepthLabel()`, `DMPlexGetCellType()`
7259: @*/
7260: PetscErrorCode DMSetLabel(DM dm, DMLabel label)
7261: {
7262:   DMLabelLink next = dm->labels;
7263:   PetscBool   hasLabel, flg;
7264:   const char *name, *lname;

7266:   PetscFunctionBegin;
7269:   PetscCall(PetscObjectGetName((PetscObject)label, &name));
7270:   while (next) {
7271:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7272:     PetscCall(PetscStrcmp(name, lname, &hasLabel));
7273:     if (hasLabel) {
7274:       PetscCall(PetscObjectReference((PetscObject)label));
7275:       PetscCall(PetscStrcmp(lname, "depth", &flg));
7276:       if (flg) dm->depthLabel = label;
7277:       PetscCall(PetscStrcmp(lname, "celltype", &flg));
7278:       if (flg) dm->celltypeLabel = label;
7279:       PetscCall(DMLabelDestroy(&next->label));
7280:       next->label = label;
7281:       break;
7282:     }
7283:     next = next->next;
7284:   }
7285:   PetscFunctionReturn(PETSC_SUCCESS);
7286: }

7288: /*@C
7289:   DMRemoveLabel - Remove the label given by name from this `DM`

7291:   Not Collective

7293:   Input Parameters:
7294: + dm   - The `DM` object
7295: - name - The label name

7297:   Output Parameter:
7298: . label - The `DMLabel`, or `NULL` if the label is absent. Pass in `NULL` to call `DMLabelDestroy()` on the label, otherwise the
7299:           caller is responsible for calling `DMLabelDestroy()`.

7301:   Level: developer

7303: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMLabelDestroy()`, `DMRemoveLabelBySelf()`
7304: @*/
7305: PetscErrorCode DMRemoveLabel(DM dm, const char name[], DMLabel *label)
7306: {
7307:   DMLabelLink link, *pnext;
7308:   PetscBool   hasLabel;
7309:   const char *lname;

7311:   PetscFunctionBegin;
7313:   PetscAssertPointer(name, 2);
7314:   if (label) {
7315:     PetscAssertPointer(label, 3);
7316:     *label = NULL;
7317:   }
7318:   for (pnext = &dm->labels; (link = *pnext); pnext = &link->next) {
7319:     PetscCall(PetscObjectGetName((PetscObject)link->label, &lname));
7320:     PetscCall(PetscStrcmp(name, lname, &hasLabel));
7321:     if (hasLabel) {
7322:       *pnext = link->next; /* Remove from list */
7323:       PetscCall(PetscStrcmp(name, "depth", &hasLabel));
7324:       if (hasLabel) dm->depthLabel = NULL;
7325:       PetscCall(PetscStrcmp(name, "celltype", &hasLabel));
7326:       if (hasLabel) dm->celltypeLabel = NULL;
7327:       if (label) *label = link->label;
7328:       else PetscCall(DMLabelDestroy(&link->label));
7329:       PetscCall(PetscFree(link));
7330:       break;
7331:     }
7332:   }
7333:   PetscFunctionReturn(PETSC_SUCCESS);
7334: }

7336: /*@
7337:   DMRemoveLabelBySelf - Remove the label from this `DM`

7339:   Not Collective

7341:   Input Parameters:
7342: + dm           - The `DM` object
7343: . label        - The `DMLabel` to be removed from the `DM`
7344: - failNotFound - Should it fail if the label is not found in the `DM`?

7346:   Level: developer

7348:   Note:
7349:   Only exactly the same instance is removed if found, name match is ignored.
7350:   If the `DM` has an exclusive reference to the label, the label gets destroyed and
7351:   *label nullified.

7353: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabel()` `DMGetLabelValue()`, `DMSetLabelValue()`, `DMLabelDestroy()`, `DMRemoveLabel()`
7354: @*/
7355: PetscErrorCode DMRemoveLabelBySelf(DM dm, DMLabel *label, PetscBool failNotFound)
7356: {
7357:   DMLabelLink link, *pnext;
7358:   PetscBool   hasLabel = PETSC_FALSE;

7360:   PetscFunctionBegin;
7362:   PetscAssertPointer(label, 2);
7363:   if (!*label && !failNotFound) PetscFunctionReturn(PETSC_SUCCESS);
7366:   for (pnext = &dm->labels; (link = *pnext); pnext = &link->next) {
7367:     if (*label == link->label) {
7368:       hasLabel = PETSC_TRUE;
7369:       *pnext   = link->next; /* Remove from list */
7370:       if (*label == dm->depthLabel) dm->depthLabel = NULL;
7371:       if (*label == dm->celltypeLabel) dm->celltypeLabel = NULL;
7372:       if (((PetscObject)link->label)->refct < 2) *label = NULL; /* nullify if exclusive reference */
7373:       PetscCall(DMLabelDestroy(&link->label));
7374:       PetscCall(PetscFree(link));
7375:       break;
7376:     }
7377:   }
7378:   PetscCheck(hasLabel || !failNotFound, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Given label not found in DM");
7379:   PetscFunctionReturn(PETSC_SUCCESS);
7380: }

7382: /*@C
7383:   DMGetLabelOutput - Get the output flag for a given label

7385:   Not Collective

7387:   Input Parameters:
7388: + dm   - The `DM` object
7389: - name - The label name

7391:   Output Parameter:
7392: . output - The flag for output

7394:   Level: developer

7396: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMSetLabelOutput()`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7397: @*/
7398: PetscErrorCode DMGetLabelOutput(DM dm, const char name[], PetscBool *output)
7399: {
7400:   DMLabelLink next = dm->labels;
7401:   const char *lname;

7403:   PetscFunctionBegin;
7405:   PetscAssertPointer(name, 2);
7406:   PetscAssertPointer(output, 3);
7407:   while (next) {
7408:     PetscBool flg;

7410:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7411:     PetscCall(PetscStrcmp(name, lname, &flg));
7412:     if (flg) {
7413:       *output = next->output;
7414:       PetscFunctionReturn(PETSC_SUCCESS);
7415:     }
7416:     next = next->next;
7417:   }
7418:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No label named %s was present in this dm", name);
7419: }

7421: /*@C
7422:   DMSetLabelOutput - Set if a given label should be saved to a `PetscViewer` in calls to `DMView()`

7424:   Not Collective

7426:   Input Parameters:
7427: + dm     - The `DM` object
7428: . name   - The label name
7429: - output - `PETSC_TRUE` to save the label to the viewer

7431:   Level: developer

7433: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetOutputFlag()`, `DMGetLabelOutput()`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7434: @*/
7435: PetscErrorCode DMSetLabelOutput(DM dm, const char name[], PetscBool output)
7436: {
7437:   DMLabelLink next = dm->labels;
7438:   const char *lname;

7440:   PetscFunctionBegin;
7442:   PetscAssertPointer(name, 2);
7443:   while (next) {
7444:     PetscBool flg;

7446:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7447:     PetscCall(PetscStrcmp(name, lname, &flg));
7448:     if (flg) {
7449:       next->output = output;
7450:       PetscFunctionReturn(PETSC_SUCCESS);
7451:     }
7452:     next = next->next;
7453:   }
7454:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No label named %s was present in this dm", name);
7455: }

7457: /*@
7458:   DMCopyLabels - Copy labels from one `DM` mesh to another `DM` with a superset of the points

7460:   Collective

7462:   Input Parameters:
7463: + dmA   - The `DM` object with initial labels
7464: . dmB   - The `DM` object to which labels are copied
7465: . mode  - Copy labels by pointers (`PETSC_OWN_POINTER`) or duplicate them (`PETSC_COPY_VALUES`)
7466: . all   - Copy all labels including "depth", "dim", and "celltype" (`PETSC_TRUE`) which are otherwise ignored (`PETSC_FALSE`)
7467: - emode - How to behave when a `DMLabel` in the source and destination `DM`s with the same name is encountered (see `DMCopyLabelsMode`)

7469:   Level: intermediate

7471:   Note:
7472:   This is typically used when interpolating or otherwise adding to a mesh, or testing.

7474: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMAddLabel()`, `DMCopyLabelsMode`
7475: @*/
7476: PetscErrorCode DMCopyLabels(DM dmA, DM dmB, PetscCopyMode mode, PetscBool all, DMCopyLabelsMode emode)
7477: {
7478:   DMLabel     label, labelNew, labelOld;
7479:   const char *name;
7480:   PetscBool   flg;
7481:   DMLabelLink link;

7483:   PetscFunctionBegin;
7488:   PetscCheck(mode != PETSC_USE_POINTER, PetscObjectComm((PetscObject)dmA), PETSC_ERR_SUP, "PETSC_USE_POINTER not supported for objects");
7489:   if (dmA == dmB) PetscFunctionReturn(PETSC_SUCCESS);
7490:   for (link = dmA->labels; link; link = link->next) {
7491:     label = link->label;
7492:     PetscCall(PetscObjectGetName((PetscObject)label, &name));
7493:     if (!all) {
7494:       PetscCall(PetscStrcmp(name, "depth", &flg));
7495:       if (flg) continue;
7496:       PetscCall(PetscStrcmp(name, "dim", &flg));
7497:       if (flg) continue;
7498:       PetscCall(PetscStrcmp(name, "celltype", &flg));
7499:       if (flg) continue;
7500:     }
7501:     PetscCall(DMGetLabel(dmB, name, &labelOld));
7502:     if (labelOld) {
7503:       switch (emode) {
7504:       case DM_COPY_LABELS_KEEP:
7505:         continue;
7506:       case DM_COPY_LABELS_REPLACE:
7507:         PetscCall(DMRemoveLabelBySelf(dmB, &labelOld, PETSC_TRUE));
7508:         break;
7509:       case DM_COPY_LABELS_FAIL:
7510:         SETERRQ(PetscObjectComm((PetscObject)dmA), PETSC_ERR_ARG_OUTOFRANGE, "Label %s already exists in destination DM", name);
7511:       default:
7512:         SETERRQ(PetscObjectComm((PetscObject)dmA), PETSC_ERR_ARG_OUTOFRANGE, "Unhandled DMCopyLabelsMode %d", (int)emode);
7513:       }
7514:     }
7515:     if (mode == PETSC_COPY_VALUES) {
7516:       PetscCall(DMLabelDuplicate(label, &labelNew));
7517:     } else {
7518:       labelNew = label;
7519:     }
7520:     PetscCall(DMAddLabel(dmB, labelNew));
7521:     if (mode == PETSC_COPY_VALUES) PetscCall(DMLabelDestroy(&labelNew));
7522:   }
7523:   PetscFunctionReturn(PETSC_SUCCESS);
7524: }

7526: /*@C
7527:   DMCompareLabels - Compare labels between two `DM` objects

7529:   Collective; No Fortran Support

7531:   Input Parameters:
7532: + dm0 - First `DM` object
7533: - dm1 - Second `DM` object

7535:   Output Parameters:
7536: + equal   - (Optional) Flag whether labels of dm0 and dm1 are the same
7537: - message - (Optional) Message describing the difference, or `NULL` if there is no difference

7539:   Level: intermediate

7541:   Notes:
7542:   The output flag equal will be the same on all processes.

7544:   If equal is passed as `NULL` and difference is found, an error is thrown on all processes.

7546:   Make sure to pass equal is `NULL` on all processes or none of them.

7548:   The output message is set independently on each rank.

7550:   message must be freed with `PetscFree()`

7552:   If message is passed as `NULL` and a difference is found, the difference description is printed to stderr in synchronized manner.

7554:   Make sure to pass message as `NULL` on all processes or no processes.

7556:   Labels are matched by name. If the number of labels and their names are equal,
7557:   `DMLabelCompare()` is used to compare each pair of labels with the same name.

7559: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMAddLabel()`, `DMCopyLabelsMode`, `DMLabelCompare()`
7560: @*/
7561: PetscErrorCode DMCompareLabels(DM dm0, DM dm1, PetscBool *equal, char **message)
7562: {
7563:   PetscInt    n, i;
7564:   char        msg[PETSC_MAX_PATH_LEN] = "";
7565:   PetscBool   eq;
7566:   MPI_Comm    comm;
7567:   PetscMPIInt rank;

7569:   PetscFunctionBegin;
7572:   PetscCheckSameComm(dm0, 1, dm1, 2);
7573:   if (equal) PetscAssertPointer(equal, 3);
7574:   if (message) PetscAssertPointer(message, 4);
7575:   PetscCall(PetscObjectGetComm((PetscObject)dm0, &comm));
7576:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7577:   {
7578:     PetscInt n1;

7580:     PetscCall(DMGetNumLabels(dm0, &n));
7581:     PetscCall(DMGetNumLabels(dm1, &n1));
7582:     eq = (PetscBool)(n == n1);
7583:     if (!eq) PetscCall(PetscSNPrintf(msg, sizeof(msg), "Number of labels in dm0 = %" PetscInt_FMT " != %" PetscInt_FMT " = Number of labels in dm1", n, n1));
7584:     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPIU_BOOL, MPI_LAND, comm));
7585:     if (!eq) goto finish;
7586:   }
7587:   for (i = 0; i < n; i++) {
7588:     DMLabel     l0, l1;
7589:     const char *name;
7590:     char       *msgInner;

7592:     /* Ignore label order */
7593:     PetscCall(DMGetLabelByNum(dm0, i, &l0));
7594:     PetscCall(PetscObjectGetName((PetscObject)l0, &name));
7595:     PetscCall(DMGetLabel(dm1, name, &l1));
7596:     if (!l1) {
7597:       PetscCall(PetscSNPrintf(msg, sizeof(msg), "Label \"%s\" (#%" PetscInt_FMT " in dm0) not found in dm1", name, i));
7598:       eq = PETSC_FALSE;
7599:       break;
7600:     }
7601:     PetscCall(DMLabelCompare(comm, l0, l1, &eq, &msgInner));
7602:     PetscCall(PetscStrncpy(msg, msgInner, sizeof(msg)));
7603:     PetscCall(PetscFree(msgInner));
7604:     if (!eq) break;
7605:   }
7606:   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPIU_BOOL, MPI_LAND, comm));
7607: finish:
7608:   /* If message output arg not set, print to stderr */
7609:   if (message) {
7610:     *message = NULL;
7611:     if (msg[0]) PetscCall(PetscStrallocpy(msg, message));
7612:   } else {
7613:     if (msg[0]) PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] %s\n", rank, msg));
7614:     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR));
7615:   }
7616:   /* If same output arg not ser and labels are not equal, throw error */
7617:   if (equal) *equal = eq;
7618:   else PetscCheck(eq, comm, PETSC_ERR_ARG_INCOMP, "DMLabels are not the same in dm0 and dm1");
7619:   PetscFunctionReturn(PETSC_SUCCESS);
7620: }

7622: PetscErrorCode DMSetLabelValue_Fast(DM dm, DMLabel *label, const char name[], PetscInt point, PetscInt value)
7623: {
7624:   PetscFunctionBegin;
7625:   PetscAssertPointer(label, 2);
7626:   if (!*label) {
7627:     PetscCall(DMCreateLabel(dm, name));
7628:     PetscCall(DMGetLabel(dm, name, label));
7629:   }
7630:   PetscCall(DMLabelSetValue(*label, point, value));
7631:   PetscFunctionReturn(PETSC_SUCCESS);
7632: }

7634: /*
7635:   Many mesh programs, such as Triangle and TetGen, allow only a single label for each mesh point. Therefore, we would
7636:   like to encode all label IDs using a single, universal label. We can do this by assigning an integer to every
7637:   (label, id) pair in the DM.

7639:   However, a mesh point can have multiple labels, so we must separate all these values. We will assign a bit range to
7640:   each label.
7641: */
7642: PetscErrorCode DMUniversalLabelCreate(DM dm, DMUniversalLabel *universal)
7643: {
7644:   DMUniversalLabel ul;
7645:   PetscBool       *active;
7646:   PetscInt         pStart, pEnd, p, Nl, l, m;

7648:   PetscFunctionBegin;
7649:   PetscCall(PetscMalloc1(1, &ul));
7650:   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "universal", &ul->label));
7651:   PetscCall(DMGetNumLabels(dm, &Nl));
7652:   PetscCall(PetscCalloc1(Nl, &active));
7653:   ul->Nl = 0;
7654:   for (l = 0; l < Nl; ++l) {
7655:     PetscBool   isdepth, iscelltype;
7656:     const char *name;

7658:     PetscCall(DMGetLabelName(dm, l, &name));
7659:     PetscCall(PetscStrncmp(name, "depth", 6, &isdepth));
7660:     PetscCall(PetscStrncmp(name, "celltype", 9, &iscelltype));
7661:     active[l] = !(isdepth || iscelltype) ? PETSC_TRUE : PETSC_FALSE;
7662:     if (active[l]) ++ul->Nl;
7663:   }
7664:   PetscCall(PetscCalloc5(ul->Nl, &ul->names, ul->Nl, &ul->indices, ul->Nl + 1, &ul->offsets, ul->Nl + 1, &ul->bits, ul->Nl, &ul->masks));
7665:   ul->Nv = 0;
7666:   for (l = 0, m = 0; l < Nl; ++l) {
7667:     DMLabel     label;
7668:     PetscInt    nv;
7669:     const char *name;

7671:     if (!active[l]) continue;
7672:     PetscCall(DMGetLabelName(dm, l, &name));
7673:     PetscCall(DMGetLabelByNum(dm, l, &label));
7674:     PetscCall(DMLabelGetNumValues(label, &nv));
7675:     PetscCall(PetscStrallocpy(name, &ul->names[m]));
7676:     ul->indices[m] = l;
7677:     ul->Nv += nv;
7678:     ul->offsets[m + 1] = nv;
7679:     ul->bits[m + 1]    = PetscCeilReal(PetscLog2Real(nv + 1));
7680:     ++m;
7681:   }
7682:   for (l = 1; l <= ul->Nl; ++l) {
7683:     ul->offsets[l] = ul->offsets[l - 1] + ul->offsets[l];
7684:     ul->bits[l]    = ul->bits[l - 1] + ul->bits[l];
7685:   }
7686:   for (l = 0; l < ul->Nl; ++l) {
7687:     PetscInt b;

7689:     ul->masks[l] = 0;
7690:     for (b = ul->bits[l]; b < ul->bits[l + 1]; ++b) ul->masks[l] |= 1 << b;
7691:   }
7692:   PetscCall(PetscMalloc1(ul->Nv, &ul->values));
7693:   for (l = 0, m = 0; l < Nl; ++l) {
7694:     DMLabel         label;
7695:     IS              valueIS;
7696:     const PetscInt *varr;
7697:     PetscInt        nv, v;

7699:     if (!active[l]) continue;
7700:     PetscCall(DMGetLabelByNum(dm, l, &label));
7701:     PetscCall(DMLabelGetNumValues(label, &nv));
7702:     PetscCall(DMLabelGetValueIS(label, &valueIS));
7703:     PetscCall(ISGetIndices(valueIS, &varr));
7704:     for (v = 0; v < nv; ++v) ul->values[ul->offsets[m] + v] = varr[v];
7705:     PetscCall(ISRestoreIndices(valueIS, &varr));
7706:     PetscCall(ISDestroy(&valueIS));
7707:     PetscCall(PetscSortInt(nv, &ul->values[ul->offsets[m]]));
7708:     ++m;
7709:   }
7710:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
7711:   for (p = pStart; p < pEnd; ++p) {
7712:     PetscInt  uval   = 0;
7713:     PetscBool marked = PETSC_FALSE;

7715:     for (l = 0, m = 0; l < Nl; ++l) {
7716:       DMLabel  label;
7717:       PetscInt val, defval, loc, nv;

7719:       if (!active[l]) continue;
7720:       PetscCall(DMGetLabelByNum(dm, l, &label));
7721:       PetscCall(DMLabelGetValue(label, p, &val));
7722:       PetscCall(DMLabelGetDefaultValue(label, &defval));
7723:       if (val == defval) {
7724:         ++m;
7725:         continue;
7726:       }
7727:       nv     = ul->offsets[m + 1] - ul->offsets[m];
7728:       marked = PETSC_TRUE;
7729:       PetscCall(PetscFindInt(val, nv, &ul->values[ul->offsets[m]], &loc));
7730:       PetscCheck(loc >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Label value %" PetscInt_FMT " not found in compression array", val);
7731:       uval += (loc + 1) << ul->bits[m];
7732:       ++m;
7733:     }
7734:     if (marked) PetscCall(DMLabelSetValue(ul->label, p, uval));
7735:   }
7736:   PetscCall(PetscFree(active));
7737:   *universal = ul;
7738:   PetscFunctionReturn(PETSC_SUCCESS);
7739: }

7741: PetscErrorCode DMUniversalLabelDestroy(DMUniversalLabel *universal)
7742: {
7743:   PetscInt l;

7745:   PetscFunctionBegin;
7746:   for (l = 0; l < (*universal)->Nl; ++l) PetscCall(PetscFree((*universal)->names[l]));
7747:   PetscCall(DMLabelDestroy(&(*universal)->label));
7748:   PetscCall(PetscFree5((*universal)->names, (*universal)->indices, (*universal)->offsets, (*universal)->bits, (*universal)->masks));
7749:   PetscCall(PetscFree((*universal)->values));
7750:   PetscCall(PetscFree(*universal));
7751:   *universal = NULL;
7752:   PetscFunctionReturn(PETSC_SUCCESS);
7753: }

7755: PetscErrorCode DMUniversalLabelGetLabel(DMUniversalLabel ul, DMLabel *ulabel)
7756: {
7757:   PetscFunctionBegin;
7758:   PetscAssertPointer(ulabel, 2);
7759:   *ulabel = ul->label;
7760:   PetscFunctionReturn(PETSC_SUCCESS);
7761: }

7763: PetscErrorCode DMUniversalLabelCreateLabels(DMUniversalLabel ul, PetscBool preserveOrder, DM dm)
7764: {
7765:   PetscInt Nl = ul->Nl, l;

7767:   PetscFunctionBegin;
7769:   for (l = 0; l < Nl; ++l) {
7770:     if (preserveOrder) PetscCall(DMCreateLabelAtIndex(dm, ul->indices[l], ul->names[l]));
7771:     else PetscCall(DMCreateLabel(dm, ul->names[l]));
7772:   }
7773:   if (preserveOrder) {
7774:     for (l = 0; l < ul->Nl; ++l) {
7775:       const char *name;
7776:       PetscBool   match;

7778:       PetscCall(DMGetLabelName(dm, ul->indices[l], &name));
7779:       PetscCall(PetscStrcmp(name, ul->names[l], &match));
7780:       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]);
7781:     }
7782:   }
7783:   PetscFunctionReturn(PETSC_SUCCESS);
7784: }

7786: PetscErrorCode DMUniversalLabelSetLabelValue(DMUniversalLabel ul, DM dm, PetscBool useIndex, PetscInt p, PetscInt value)
7787: {
7788:   PetscInt l;

7790:   PetscFunctionBegin;
7791:   for (l = 0; l < ul->Nl; ++l) {
7792:     DMLabel  label;
7793:     PetscInt lval = (value & ul->masks[l]) >> ul->bits[l];

7795:     if (lval) {
7796:       if (useIndex) PetscCall(DMGetLabelByNum(dm, ul->indices[l], &label));
7797:       else PetscCall(DMGetLabel(dm, ul->names[l], &label));
7798:       PetscCall(DMLabelSetValue(label, p, ul->values[ul->offsets[l] + lval - 1]));
7799:     }
7800:   }
7801:   PetscFunctionReturn(PETSC_SUCCESS);
7802: }

7804: /*@
7805:   DMGetCoarseDM - Get the coarse `DM`from which this `DM` was obtained by refinement

7807:   Not Collective

7809:   Input Parameter:
7810: . dm - The `DM` object

7812:   Output Parameter:
7813: . cdm - The coarse `DM`

7815:   Level: intermediate

7817: .seealso: [](ch_dmbase), `DM`, `DMSetCoarseDM()`, `DMCoarsen()`
7818: @*/
7819: PetscErrorCode DMGetCoarseDM(DM dm, DM *cdm)
7820: {
7821:   PetscFunctionBegin;
7823:   PetscAssertPointer(cdm, 2);
7824:   *cdm = dm->coarseMesh;
7825:   PetscFunctionReturn(PETSC_SUCCESS);
7826: }

7828: /*@
7829:   DMSetCoarseDM - Set the coarse `DM` from which this `DM` was obtained by refinement

7831:   Input Parameters:
7832: + dm  - The `DM` object
7833: - cdm - The coarse `DM`

7835:   Level: intermediate

7837:   Note:
7838:   Normally this is set automatically by `DMRefine()`

7840: .seealso: [](ch_dmbase), `DM`, `DMGetCoarseDM()`, `DMCoarsen()`, `DMSetRefine()`, `DMSetFineDM()`
7841: @*/
7842: PetscErrorCode DMSetCoarseDM(DM dm, DM cdm)
7843: {
7844:   PetscFunctionBegin;
7847:   if (dm == cdm) cdm = NULL;
7848:   PetscCall(PetscObjectReference((PetscObject)cdm));
7849:   PetscCall(DMDestroy(&dm->coarseMesh));
7850:   dm->coarseMesh = cdm;
7851:   PetscFunctionReturn(PETSC_SUCCESS);
7852: }

7854: /*@
7855:   DMGetFineDM - Get the fine mesh from which this `DM` was obtained by coarsening

7857:   Input Parameter:
7858: . dm - The `DM` object

7860:   Output Parameter:
7861: . fdm - The fine `DM`

7863:   Level: intermediate

7865: .seealso: [](ch_dmbase), `DM`, `DMSetFineDM()`, `DMCoarsen()`, `DMRefine()`
7866: @*/
7867: PetscErrorCode DMGetFineDM(DM dm, DM *fdm)
7868: {
7869:   PetscFunctionBegin;
7871:   PetscAssertPointer(fdm, 2);
7872:   *fdm = dm->fineMesh;
7873:   PetscFunctionReturn(PETSC_SUCCESS);
7874: }

7876: /*@
7877:   DMSetFineDM - Set the fine mesh from which this was obtained by coarsening

7879:   Input Parameters:
7880: + dm  - The `DM` object
7881: - fdm - The fine `DM`

7883:   Level: developer

7885:   Note:
7886:   Normally this is set automatically by `DMCoarsen()`

7888: .seealso: [](ch_dmbase), `DM`, `DMGetFineDM()`, `DMCoarsen()`, `DMRefine()`
7889: @*/
7890: PetscErrorCode DMSetFineDM(DM dm, DM fdm)
7891: {
7892:   PetscFunctionBegin;
7895:   if (dm == fdm) fdm = NULL;
7896:   PetscCall(PetscObjectReference((PetscObject)fdm));
7897:   PetscCall(DMDestroy(&dm->fineMesh));
7898:   dm->fineMesh = fdm;
7899:   PetscFunctionReturn(PETSC_SUCCESS);
7900: }

7902: /*@C
7903:   DMAddBoundary - Add a boundary condition to a model represented by a `DM`

7905:   Collective

7907:   Input Parameters:
7908: + dm       - The `DM`, with a `PetscDS` that matches the problem being constrained
7909: . type     - The type of condition, e.g. `DM_BC_ESSENTIAL_ANALYTIC`, `DM_BC_ESSENTIAL_FIELD` (Dirichlet), or `DM_BC_NATURAL` (Neumann)
7910: . name     - The BC name
7911: . label    - The label defining constrained points
7912: . Nv       - The number of `DMLabel` values for constrained points
7913: . values   - An array of values for constrained points
7914: . field    - The field to constrain
7915: . Nc       - The number of constrained field components (0 will constrain all fields)
7916: . comps    - An array of constrained component numbers
7917: . bcFunc   - A pointwise function giving boundary values
7918: . bcFunc_t - A pointwise function giving the time deriative of the boundary values, or NULL
7919: - ctx      - An optional user context for bcFunc

7921:   Output Parameter:
7922: . bd - (Optional) Boundary number

7924:   Options Database Keys:
7925: + -bc_<boundary name> <num>      - Overrides the boundary ids
7926: - -bc_<boundary name>_comp <num> - Overrides the boundary components

7928:   Level: intermediate

7930:   Notes:
7931:   Both bcFunc and bcFunc_t will depend on the boundary condition type. If the type if `DM_BC_ESSENTIAL`, then the calling sequence is\:
7932: .vb
7933:  void bcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[])
7934: .ve

7936:   If the type is `DM_BC_ESSENTIAL_FIELD` or other _FIELD value, then the calling sequence is\:

7938: .vb
7939:   void bcFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
7940:               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
7941:               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
7942:               PetscReal time, const PetscReal x[], PetscScalar bcval[])
7943: .ve
7944: + dim - the spatial dimension
7945: . Nf - the number of fields
7946: . uOff - the offset into u[] and u_t[] for each field
7947: . uOff_x - the offset into u_x[] for each field
7948: . u - each field evaluated at the current point
7949: . u_t - the time derivative of each field evaluated at the current point
7950: . u_x - the gradient of each field evaluated at the current point
7951: . aOff - the offset into a[] and a_t[] for each auxiliary field
7952: . aOff_x - the offset into a_x[] for each auxiliary field
7953: . a - each auxiliary field evaluated at the current point
7954: . a_t - the time derivative of each auxiliary field evaluated at the current point
7955: . a_x - the gradient of auxiliary each field evaluated at the current point
7956: . t - current time
7957: . x - coordinates of the current point
7958: . numConstants - number of constant parameters
7959: . constants - constant parameters
7960: - bcval - output values at the current point

7962: .seealso: [](ch_dmbase), `DM`, `DSGetBoundary()`, `PetscDSAddBoundary()`
7963: @*/
7964: 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)
7965: {
7966:   PetscDS ds;

7968:   PetscFunctionBegin;
7975:   PetscCheck(!dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot add boundary to DM after creating local section");
7976:   PetscCall(DMGetDS(dm, &ds));
7977:   /* Complete label */
7978:   if (label) {
7979:     PetscObject  obj;
7980:     PetscClassId id;

7982:     PetscCall(DMGetField(dm, field, NULL, &obj));
7983:     PetscCall(PetscObjectGetClassId(obj, &id));
7984:     if (id == PETSCFE_CLASSID) {
7985:       DM plex;

7987:       PetscCall(DMConvert(dm, DMPLEX, &plex));
7988:       if (plex) PetscCall(DMPlexLabelComplete(plex, label));
7989:       PetscCall(DMDestroy(&plex));
7990:     }
7991:   }
7992:   PetscCall(PetscDSAddBoundary(ds, type, name, label, Nv, values, field, Nc, comps, bcFunc, bcFunc_t, ctx, bd));
7993:   PetscFunctionReturn(PETSC_SUCCESS);
7994: }

7996: /* TODO Remove this since now the structures are the same */
7997: static PetscErrorCode DMPopulateBoundary(DM dm)
7998: {
7999:   PetscDS     ds;
8000:   DMBoundary *lastnext;
8001:   DSBoundary  dsbound;

8003:   PetscFunctionBegin;
8004:   PetscCall(DMGetDS(dm, &ds));
8005:   dsbound = ds->boundary;
8006:   if (dm->boundary) {
8007:     DMBoundary next = dm->boundary;

8009:     /* quick check to see if the PetscDS has changed */
8010:     if (next->dsboundary == dsbound) PetscFunctionReturn(PETSC_SUCCESS);
8011:     /* the PetscDS has changed: tear down and rebuild */
8012:     while (next) {
8013:       DMBoundary b = next;

8015:       next = b->next;
8016:       PetscCall(PetscFree(b));
8017:     }
8018:     dm->boundary = NULL;
8019:   }

8021:   lastnext = &dm->boundary;
8022:   while (dsbound) {
8023:     DMBoundary dmbound;

8025:     PetscCall(PetscNew(&dmbound));
8026:     dmbound->dsboundary = dsbound;
8027:     dmbound->label      = dsbound->label;
8028:     /* push on the back instead of the front so that it is in the same order as in the PetscDS */
8029:     *lastnext = dmbound;
8030:     lastnext  = &dmbound->next;
8031:     dsbound   = dsbound->next;
8032:   }
8033:   PetscFunctionReturn(PETSC_SUCCESS);
8034: }

8036: /* TODO: missing manual page */
8037: PetscErrorCode DMIsBoundaryPoint(DM dm, PetscInt point, PetscBool *isBd)
8038: {
8039:   DMBoundary b;

8041:   PetscFunctionBegin;
8043:   PetscAssertPointer(isBd, 3);
8044:   *isBd = PETSC_FALSE;
8045:   PetscCall(DMPopulateBoundary(dm));
8046:   b = dm->boundary;
8047:   while (b && !(*isBd)) {
8048:     DMLabel    label = b->label;
8049:     DSBoundary dsb   = b->dsboundary;
8050:     PetscInt   i;

8052:     if (label) {
8053:       for (i = 0; i < dsb->Nv && !(*isBd); ++i) PetscCall(DMLabelStratumHasPoint(label, dsb->values[i], point, isBd));
8054:     }
8055:     b = b->next;
8056:   }
8057:   PetscFunctionReturn(PETSC_SUCCESS);
8058: }

8060: /*@C
8061:   DMProjectFunction - This projects the given function into the function space provided by a `DM`, putting the coefficients in a global vector.

8063:   Collective

8065:   Input Parameters:
8066: + dm    - The `DM`
8067: . time  - The time
8068: . funcs - The coordinate functions to evaluate, one per field
8069: . ctxs  - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
8070: - mode  - The insertion mode for values

8072:   Output Parameter:
8073: . X - vector

8075:   Calling sequence of `funcs`:
8076: + dim  - The spatial dimension
8077: . time - The time at which to sample
8078: . x    - The coordinates
8079: . Nc   - The number of components
8080: . u    - The output field values
8081: - ctx  - optional user-defined function context

8083:   Level: developer

8085:   Developer Notes:
8086:   This API is specific to only particular usage of `DM`

8088:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8090: .seealso: [](ch_dmbase), `DM`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
8091: @*/
8092: 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)
8093: {
8094:   Vec localX;

8096:   PetscFunctionBegin;
8098:   PetscCall(PetscLogEventBegin(DM_ProjectFunction, dm, X, 0, 0));
8099:   PetscCall(DMGetLocalVector(dm, &localX));
8100:   PetscCall(VecSet(localX, 0.));
8101:   PetscCall(DMProjectFunctionLocal(dm, time, funcs, ctxs, mode, localX));
8102:   PetscCall(DMLocalToGlobalBegin(dm, localX, mode, X));
8103:   PetscCall(DMLocalToGlobalEnd(dm, localX, mode, X));
8104:   PetscCall(DMRestoreLocalVector(dm, &localX));
8105:   PetscCall(PetscLogEventEnd(DM_ProjectFunction, dm, X, 0, 0));
8106:   PetscFunctionReturn(PETSC_SUCCESS);
8107: }

8109: /*@C
8110:   DMProjectFunctionLocal - This projects the given function into the function space provided by a `DM`, putting the coefficients in a local vector.

8112:   Not Collective

8114:   Input Parameters:
8115: + dm    - The `DM`
8116: . time  - The time
8117: . funcs - The coordinate functions to evaluate, one per field
8118: . ctxs  - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
8119: - mode  - The insertion mode for values

8121:   Output Parameter:
8122: . localX - vector

8124:   Calling sequence of `funcs`:
8125: + dim  - The spatial dimension
8126: . time - The current timestep
8127: . x    - The coordinates
8128: . Nc   - The number of components
8129: . u    - The output field values
8130: - ctx  - optional user-defined function context

8132:   Level: developer

8134:   Developer Notes:
8135:   This API is specific to only particular usage of `DM`

8137:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8139: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
8140: @*/
8141: 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)
8142: {
8143:   PetscFunctionBegin;
8146:   PetscUseTypeMethod(dm, projectfunctionlocal, time, funcs, ctxs, mode, localX);
8147:   PetscFunctionReturn(PETSC_SUCCESS);
8148: }

8150: /*@C
8151:   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.

8153:   Collective

8155:   Input Parameters:
8156: + dm     - The `DM`
8157: . time   - The time
8158: . numIds - The number of ids
8159: . ids    - The ids
8160: . Nc     - The number of components
8161: . comps  - The components
8162: . label  - The `DMLabel` selecting the portion of the mesh for projection
8163: . funcs  - The coordinate functions to evaluate, one per field
8164: . ctxs   - Optional array of contexts to pass to each coordinate function.  ctxs may be null.
8165: - mode   - The insertion mode for values

8167:   Output Parameter:
8168: . X - vector

8170:   Calling sequence of `funcs`:
8171: + dim  - The spatial dimension
8172: . time - The current timestep
8173: . x    - The coordinates
8174: . Nc   - The number of components
8175: . u    - The output field values
8176: - ctx  - optional user-defined function context

8178:   Level: developer

8180:   Developer Notes:
8181:   This API is specific to only particular usage of `DM`

8183:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8185: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabelLocal()`, `DMComputeL2Diff()`
8186: @*/
8187: 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)
8188: {
8189:   Vec localX;

8191:   PetscFunctionBegin;
8193:   PetscCall(DMGetLocalVector(dm, &localX));
8194:   PetscCall(VecSet(localX, 0.));
8195:   PetscCall(DMProjectFunctionLabelLocal(dm, time, label, numIds, ids, Nc, comps, funcs, ctxs, mode, localX));
8196:   PetscCall(DMLocalToGlobalBegin(dm, localX, mode, X));
8197:   PetscCall(DMLocalToGlobalEnd(dm, localX, mode, X));
8198:   PetscCall(DMRestoreLocalVector(dm, &localX));
8199:   PetscFunctionReturn(PETSC_SUCCESS);
8200: }

8202: /*@C
8203:   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.

8205:   Not Collective

8207:   Input Parameters:
8208: + dm     - The `DM`
8209: . time   - The time
8210: . label  - The `DMLabel` selecting the portion of the mesh for projection
8211: . numIds - The number of ids
8212: . ids    - The ids
8213: . Nc     - The number of components
8214: . comps  - The components
8215: . funcs  - The coordinate functions to evaluate, one per field
8216: . ctxs   - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
8217: - mode   - The insertion mode for values

8219:   Output Parameter:
8220: . localX - vector

8222:   Calling sequence of `funcs`:
8223: + dim  - The spatial dimension
8224: . time - The current time
8225: . x    - The coordinates
8226: . Nc   - The number of components
8227: . u    - The output field values
8228: - ctx  - optional user-defined function context

8230:   Level: developer

8232:   Developer Notes:
8233:   This API is specific to only particular usage of `DM`

8235:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8237: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
8238: @*/
8239: 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)
8240: {
8241:   PetscFunctionBegin;
8244:   PetscUseTypeMethod(dm, projectfunctionlabellocal, time, label, numIds, ids, Nc, comps, funcs, ctxs, mode, localX);
8245:   PetscFunctionReturn(PETSC_SUCCESS);
8246: }

8248: /*@C
8249:   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.

8251:   Not Collective

8253:   Input Parameters:
8254: + dm     - The `DM`
8255: . time   - The time
8256: . localU - The input field vector; may be `NULL` if projection is defined purely by coordinates
8257: . funcs  - The functions to evaluate, one per field
8258: - mode   - The insertion mode for values

8260:   Output Parameter:
8261: . localX - The output vector

8263:   Calling sequence of `funcs`:
8264: + dim          - The spatial dimension
8265: . Nf           - The number of input fields
8266: . NfAux        - The number of input auxiliary fields
8267: . uOff         - The offset of each field in u[]
8268: . uOff_x       - The offset of each field in u_x[]
8269: . u            - The field values at this point in space
8270: . u_t          - The field time derivative at this point in space (or NULL)
8271: . u_x          - The field derivatives at this point in space
8272: . aOff         - The offset of each auxiliary field in u[]
8273: . aOff_x       - The offset of each auxiliary field in u_x[]
8274: . a            - The auxiliary field values at this point in space
8275: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8276: . a_x          - The auxiliary field derivatives at this point in space
8277: . t            - The current time
8278: . x            - The coordinates of this point
8279: . numConstants - The number of constants
8280: . constants    - The value of each constant
8281: - f            - The value of the function at this point in space

8283:   Level: intermediate

8285:   Note:
8286:   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.
8287:   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
8288:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8289:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8291:   Developer Notes:
8292:   This API is specific to only particular usage of `DM`

8294:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8296: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabelLocal()`,
8297: `DMProjectFunction()`, `DMComputeL2Diff()`
8298: @*/
8299: 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)
8300: {
8301:   PetscFunctionBegin;
8305:   PetscUseTypeMethod(dm, projectfieldlocal, time, localU, funcs, mode, localX);
8306:   PetscFunctionReturn(PETSC_SUCCESS);
8307: }

8309: /*@C
8310:   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.

8312:   Not Collective

8314:   Input Parameters:
8315: + dm     - The `DM`
8316: . time   - The time
8317: . label  - The `DMLabel` marking the portion of the domain to output
8318: . numIds - The number of label ids to use
8319: . ids    - The label ids to use for marking
8320: . Nc     - The number of components to set in the output, or `PETSC_DETERMINE` for all components
8321: . comps  - The components to set in the output, or `NULL` for all components
8322: . localU - The input field vector
8323: . funcs  - The functions to evaluate, one per field
8324: - mode   - The insertion mode for values

8326:   Output Parameter:
8327: . localX - The output vector

8329:   Calling sequence of `funcs`:
8330: + dim          - The spatial dimension
8331: . Nf           - The number of input fields
8332: . NfAux        - The number of input auxiliary fields
8333: . uOff         - The offset of each field in u[]
8334: . uOff_x       - The offset of each field in u_x[]
8335: . u            - The field values at this point in space
8336: . u_t          - The field time derivative at this point in space (or NULL)
8337: . u_x          - The field derivatives at this point in space
8338: . aOff         - The offset of each auxiliary field in u[]
8339: . aOff_x       - The offset of each auxiliary field in u_x[]
8340: . a            - The auxiliary field values at this point in space
8341: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8342: . a_x          - The auxiliary field derivatives at this point in space
8343: . t            - The current time
8344: . x            - The coordinates of this point
8345: . numConstants - The number of constants
8346: . constants    - The value of each constant
8347: - f            - The value of the function at this point in space

8349:   Level: intermediate

8351:   Note:
8352:   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.
8353:   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
8354:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8355:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8357:   Developer Notes:
8358:   This API is specific to only particular usage of `DM`

8360:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8362: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabel()`, `DMProjectFunction()`, `DMComputeL2Diff()`
8363: @*/
8364: 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)
8365: {
8366:   PetscFunctionBegin;
8370:   PetscUseTypeMethod(dm, projectfieldlabellocal, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX);
8371:   PetscFunctionReturn(PETSC_SUCCESS);
8372: }

8374: /*@C
8375:   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.

8377:   Not Collective

8379:   Input Parameters:
8380: + dm     - The `DM`
8381: . time   - The time
8382: . label  - The `DMLabel` marking the portion of the domain to output
8383: . numIds - The number of label ids to use
8384: . ids    - The label ids to use for marking
8385: . Nc     - The number of components to set in the output, or `PETSC_DETERMINE` for all components
8386: . comps  - The components to set in the output, or `NULL` for all components
8387: . U      - The input field vector
8388: . funcs  - The functions to evaluate, one per field
8389: - mode   - The insertion mode for values

8391:   Output Parameter:
8392: . X - The output vector

8394:   Calling sequence of `funcs`:
8395: + dim          - The spatial dimension
8396: . Nf           - The number of input fields
8397: . NfAux        - The number of input auxiliary fields
8398: . uOff         - The offset of each field in u[]
8399: . uOff_x       - The offset of each field in u_x[]
8400: . u            - The field values at this point in space
8401: . u_t          - The field time derivative at this point in space (or NULL)
8402: . u_x          - The field derivatives at this point in space
8403: . aOff         - The offset of each auxiliary field in u[]
8404: . aOff_x       - The offset of each auxiliary field in u_x[]
8405: . a            - The auxiliary field values at this point in space
8406: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8407: . a_x          - The auxiliary field derivatives at this point in space
8408: . t            - The current time
8409: . x            - The coordinates of this point
8410: . numConstants - The number of constants
8411: . constants    - The value of each constant
8412: - f            - The value of the function at this point in space

8414:   Level: intermediate

8416:   Note:
8417:   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.
8418:   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
8419:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8420:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8422:   Developer Notes:
8423:   This API is specific to only particular usage of `DM`

8425:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8427: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabelLocal()`, `DMProjectFunction()`, `DMComputeL2Diff()`
8428: @*/
8429: 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)
8430: {
8431:   DM  dmIn;
8432:   Vec localU, localX;

8434:   PetscFunctionBegin;
8436:   PetscCall(VecGetDM(U, &dmIn));
8437:   PetscCall(DMGetLocalVector(dmIn, &localU));
8438:   PetscCall(DMGetLocalVector(dm, &localX));
8439:   PetscCall(VecSet(localX, 0.));
8440:   PetscCall(DMGlobalToLocalBegin(dmIn, U, mode, localU));
8441:   PetscCall(DMGlobalToLocalEnd(dmIn, U, mode, localU));
8442:   PetscCall(DMProjectFieldLabelLocal(dm, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX));
8443:   PetscCall(DMLocalToGlobalBegin(dm, localX, mode, X));
8444:   PetscCall(DMLocalToGlobalEnd(dm, localX, mode, X));
8445:   PetscCall(DMRestoreLocalVector(dm, &localX));
8446:   PetscCall(DMRestoreLocalVector(dmIn, &localU));
8447:   PetscFunctionReturn(PETSC_SUCCESS);
8448: }

8450: /*@C
8451:   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.

8453:   Not Collective

8455:   Input Parameters:
8456: + dm     - The `DM`
8457: . time   - The time
8458: . label  - The `DMLabel` marking the portion of the domain boundary to output
8459: . numIds - The number of label ids to use
8460: . ids    - The label ids to use for marking
8461: . Nc     - The number of components to set in the output, or `PETSC_DETERMINE` for all components
8462: . comps  - The components to set in the output, or `NULL` for all components
8463: . localU - The input field vector
8464: . funcs  - The functions to evaluate, one per field
8465: - mode   - The insertion mode for values

8467:   Output Parameter:
8468: . localX - The output vector

8470:   Calling sequence of `funcs`:
8471: + dim          - The spatial dimension
8472: . Nf           - The number of input fields
8473: . NfAux        - The number of input auxiliary fields
8474: . uOff         - The offset of each field in u[]
8475: . uOff_x       - The offset of each field in u_x[]
8476: . u            - The field values at this point in space
8477: . u_t          - The field time derivative at this point in space (or NULL)
8478: . u_x          - The field derivatives at this point in space
8479: . aOff         - The offset of each auxiliary field in u[]
8480: . aOff_x       - The offset of each auxiliary field in u_x[]
8481: . a            - The auxiliary field values at this point in space
8482: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8483: . a_x          - The auxiliary field derivatives at this point in space
8484: . t            - The current time
8485: . x            - The coordinates of this point
8486: . n            - The face normal
8487: . numConstants - The number of constants
8488: . constants    - The value of each constant
8489: - f            - The value of the function at this point in space

8491:   Level: intermediate

8493:   Note:
8494:   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.
8495:   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
8496:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8497:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8499:   Developer Notes:
8500:   This API is specific to only particular usage of `DM`

8502:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8504: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabelLocal()`, `DMProjectFunction()`, `DMComputeL2Diff()`
8505: @*/
8506: 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)
8507: {
8508:   PetscFunctionBegin;
8512:   PetscUseTypeMethod(dm, projectbdfieldlabellocal, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX);
8513:   PetscFunctionReturn(PETSC_SUCCESS);
8514: }

8516: /*@C
8517:   DMComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.

8519:   Collective

8521:   Input Parameters:
8522: + dm    - The `DM`
8523: . time  - The time
8524: . funcs - The functions to evaluate for each field component
8525: . ctxs  - Optional array of contexts to pass to each function, or NULL.
8526: - X     - The coefficient vector u_h, a global vector

8528:   Output Parameter:
8529: . diff - The diff ||u - u_h||_2

8531:   Level: developer

8533:   Developer Notes:
8534:   This API is specific to only particular usage of `DM`

8536:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8538: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMComputeL2FieldDiff()`, `DMComputeL2GradientDiff()`
8539: @*/
8540: PetscErrorCode DMComputeL2Diff(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
8541: {
8542:   PetscFunctionBegin;
8545:   PetscUseTypeMethod(dm, computel2diff, time, funcs, ctxs, X, diff);
8546:   PetscFunctionReturn(PETSC_SUCCESS);
8547: }

8549: /*@C
8550:   DMComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h.

8552:   Collective

8554:   Input Parameters:
8555: + dm    - The `DM`
8556: . time  - The time
8557: . funcs - The gradient functions to evaluate for each field component
8558: . ctxs  - Optional array of contexts to pass to each function, or NULL.
8559: . X     - The coefficient vector u_h, a global vector
8560: - n     - The vector to project along

8562:   Output Parameter:
8563: . diff - The diff ||(grad u - grad u_h) . n||_2

8565:   Level: developer

8567:   Developer Notes:
8568:   This API is specific to only particular usage of `DM`

8570:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8572: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMComputeL2Diff()`, `DMComputeL2FieldDiff()`
8573: @*/
8574: 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)
8575: {
8576:   PetscFunctionBegin;
8579:   PetscUseTypeMethod(dm, computel2gradientdiff, time, funcs, ctxs, X, n, diff);
8580:   PetscFunctionReturn(PETSC_SUCCESS);
8581: }

8583: /*@C
8584:   DMComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components.

8586:   Collective

8588:   Input Parameters:
8589: + dm    - The `DM`
8590: . time  - The time
8591: . funcs - The functions to evaluate for each field component
8592: . ctxs  - Optional array of contexts to pass to each function, or NULL.
8593: - X     - The coefficient vector u_h, a global vector

8595:   Output Parameter:
8596: . diff - The array of differences, ||u^f - u^f_h||_2

8598:   Level: developer

8600:   Developer Notes:
8601:   This API is specific to only particular usage of `DM`

8603:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8605: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMComputeL2GradientDiff()`
8606: @*/
8607: PetscErrorCode DMComputeL2FieldDiff(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
8608: {
8609:   PetscFunctionBegin;
8612:   PetscUseTypeMethod(dm, computel2fielddiff, time, funcs, ctxs, X, diff);
8613:   PetscFunctionReturn(PETSC_SUCCESS);
8614: }

8616: /*@C
8617:   DMGetNeighbors - Gets an array containing the MPI ranks of all the processes neighbors

8619:   Not Collective

8621:   Input Parameter:
8622: . dm - The `DM`

8624:   Output Parameters:
8625: + nranks - the number of neighbours
8626: - ranks  - the neighbors ranks

8628:   Level: beginner

8630:   Note:
8631:   Do not free the array, it is freed when the `DM` is destroyed.

8633: .seealso: [](ch_dmbase), `DM`, `DMDAGetNeighbors()`, `PetscSFGetRootRanks()`
8634: @*/
8635: PetscErrorCode DMGetNeighbors(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
8636: {
8637:   PetscFunctionBegin;
8639:   PetscUseTypeMethod(dm, getneighbors, nranks, ranks);
8640:   PetscFunctionReturn(PETSC_SUCCESS);
8641: }

8643: #include <petsc/private/matimpl.h>

8645: /*
8646:     Converts the input vector to a ghosted vector and then calls the standard coloring code.
8647:     This must be a different function because it requires DM which is not defined in the Mat library
8648: */
8649: static PetscErrorCode MatFDColoringApply_AIJDM(Mat J, MatFDColoring coloring, Vec x1, void *sctx)
8650: {
8651:   PetscFunctionBegin;
8652:   if (coloring->ctype == IS_COLORING_LOCAL) {
8653:     Vec x1local;
8654:     DM  dm;
8655:     PetscCall(MatGetDM(J, &dm));
8656:     PetscCheck(dm, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_INCOMP, "IS_COLORING_LOCAL requires a DM");
8657:     PetscCall(DMGetLocalVector(dm, &x1local));
8658:     PetscCall(DMGlobalToLocalBegin(dm, x1, INSERT_VALUES, x1local));
8659:     PetscCall(DMGlobalToLocalEnd(dm, x1, INSERT_VALUES, x1local));
8660:     x1 = x1local;
8661:   }
8662:   PetscCall(MatFDColoringApply_AIJ(J, coloring, x1, sctx));
8663:   if (coloring->ctype == IS_COLORING_LOCAL) {
8664:     DM dm;
8665:     PetscCall(MatGetDM(J, &dm));
8666:     PetscCall(DMRestoreLocalVector(dm, &x1));
8667:   }
8668:   PetscFunctionReturn(PETSC_SUCCESS);
8669: }

8671: /*@
8672:   MatFDColoringUseDM - allows a `MatFDColoring` object to use the `DM` associated with the matrix to compute a `IS_COLORING_LOCAL` coloring

8674:   Input Parameters:
8675: + coloring   - The matrix to get the `DM` from
8676: - fdcoloring - the `MatFDColoring` object

8678:   Level: advanced

8680:   Developer Note:
8681:   This routine exists because the PETSc `Mat` library does not know about the `DM` objects

8683: .seealso: [](ch_dmbase), `DM`, `MatFDColoring`, `MatFDColoringCreate()`, `ISColoringType`
8684: @*/
8685: PetscErrorCode MatFDColoringUseDM(Mat coloring, MatFDColoring fdcoloring)
8686: {
8687:   PetscFunctionBegin;
8688:   coloring->ops->fdcoloringapply = MatFDColoringApply_AIJDM;
8689:   PetscFunctionReturn(PETSC_SUCCESS);
8690: }

8692: /*@
8693:   DMGetCompatibility - determine if two `DM`s are compatible

8695:   Collective

8697:   Input Parameters:
8698: + dm1 - the first `DM`
8699: - dm2 - the second `DM`

8701:   Output Parameters:
8702: + compatible - whether or not the two `DM`s are compatible
8703: - set        - whether or not the compatible value was actually determined and set

8705:   Level: advanced

8707:   Notes:
8708:   Two `DM`s are deemed compatible if they represent the same parallel decomposition
8709:   of the same topology. This implies that the section (field data) on one
8710:   "makes sense" with respect to the topology and parallel decomposition of the other.
8711:   Loosely speaking, compatible `DM`s represent the same domain and parallel
8712:   decomposition, but hold different data.

8714:   Typically, one would confirm compatibility if intending to simultaneously iterate
8715:   over a pair of vectors obtained from different `DM`s.

8717:   For example, two `DMDA` objects are compatible if they have the same local
8718:   and global sizes and the same stencil width. They can have different numbers
8719:   of degrees of freedom per node. Thus, one could use the node numbering from
8720:   either `DM` in bounds for a loop over vectors derived from either `DM`.

8722:   Consider the operation of summing data living on a 2-dof `DMDA` to data living
8723:   on a 1-dof `DMDA`, which should be compatible, as in the following snippet.
8724: .vb
8725:   ...
8726:   PetscCall(DMGetCompatibility(da1,da2,&compatible,&set));
8727:   if (set && compatible)  {
8728:     PetscCall(DMDAVecGetArrayDOF(da1,vec1,&arr1));
8729:     PetscCall(DMDAVecGetArrayDOF(da2,vec2,&arr2));
8730:     PetscCall(DMDAGetCorners(da1,&x,&y,NULL,&m,&n,NULL));
8731:     for (j=y; j<y+n; ++j) {
8732:       for (i=x; i<x+m, ++i) {
8733:         arr1[j][i][0] = arr2[j][i][0] + arr2[j][i][1];
8734:       }
8735:     }
8736:     PetscCall(DMDAVecRestoreArrayDOF(da1,vec1,&arr1));
8737:     PetscCall(DMDAVecRestoreArrayDOF(da2,vec2,&arr2));
8738:   } else {
8739:     SETERRQ(PetscObjectComm((PetscObject)da1,PETSC_ERR_ARG_INCOMP,"DMDA objects incompatible");
8740:   }
8741:   ...
8742: .ve

8744:   Checking compatibility might be expensive for a given implementation of `DM`,
8745:   or might be impossible to unambiguously confirm or deny. For this reason,
8746:   this function may decline to determine compatibility, and hence users should
8747:   always check the "set" output parameter.

8749:   A `DM` is always compatible with itself.

8751:   In the current implementation, `DM`s which live on "unequal" communicators
8752:   (MPI_UNEQUAL in the terminology of MPI_Comm_compare()) are always deemed
8753:   incompatible.

8755:   This function is labeled "Collective," as information about all subdomains
8756:   is required on each rank. However, in `DM` implementations which store all this
8757:   information locally, this function may be merely "Logically Collective".

8759:   Developer Note:
8760:   Compatibility is assumed to be a symmetric concept; `DM` A is compatible with `DM` B
8761:   iff B is compatible with A. Thus, this function checks the implementations
8762:   of both dm and dmc (if they are of different types), attempting to determine
8763:   compatibility. It is left to `DM` implementers to ensure that symmetry is
8764:   preserved. The simplest way to do this is, when implementing type-specific
8765:   logic for this function, is to check for existing logic in the implementation
8766:   of other `DM` types and let *set = PETSC_FALSE if found.

8768: .seealso: [](ch_dmbase), `DM`, `DMDACreateCompatibleDMDA()`, `DMStagCreateCompatibleDMStag()`
8769: @*/
8770: PetscErrorCode DMGetCompatibility(DM dm1, DM dm2, PetscBool *compatible, PetscBool *set)
8771: {
8772:   PetscMPIInt compareResult;
8773:   DMType      type, type2;
8774:   PetscBool   sameType;

8776:   PetscFunctionBegin;

8780:   /* Declare a DM compatible with itself */
8781:   if (dm1 == dm2) {
8782:     *set        = PETSC_TRUE;
8783:     *compatible = PETSC_TRUE;
8784:     PetscFunctionReturn(PETSC_SUCCESS);
8785:   }

8787:   /* Declare a DM incompatible with a DM that lives on an "unequal"
8788:      communicator. Note that this does not preclude compatibility with
8789:      DMs living on "congruent" or "similar" communicators, but this must be
8790:      determined by the implementation-specific logic */
8791:   PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dm1), PetscObjectComm((PetscObject)dm2), &compareResult));
8792:   if (compareResult == MPI_UNEQUAL) {
8793:     *set        = PETSC_TRUE;
8794:     *compatible = PETSC_FALSE;
8795:     PetscFunctionReturn(PETSC_SUCCESS);
8796:   }

8798:   /* Pass to the implementation-specific routine, if one exists. */
8799:   if (dm1->ops->getcompatibility) {
8800:     PetscUseTypeMethod(dm1, getcompatibility, dm2, compatible, set);
8801:     if (*set) PetscFunctionReturn(PETSC_SUCCESS);
8802:   }

8804:   /* If dm1 and dm2 are of different types, then attempt to check compatibility
8805:      with an implementation of this function from dm2 */
8806:   PetscCall(DMGetType(dm1, &type));
8807:   PetscCall(DMGetType(dm2, &type2));
8808:   PetscCall(PetscStrcmp(type, type2, &sameType));
8809:   if (!sameType && dm2->ops->getcompatibility) {
8810:     PetscUseTypeMethod(dm2, getcompatibility, dm1, compatible, set); /* Note argument order */
8811:   } else {
8812:     *set = PETSC_FALSE;
8813:   }
8814:   PetscFunctionReturn(PETSC_SUCCESS);
8815: }

8817: /*@C
8818:   DMMonitorSet - Sets an additional monitor function that is to be used after a solve to monitor discretization performance.

8820:   Logically Collective

8822:   Input Parameters:
8823: + dm             - the `DM`
8824: . f              - the monitor function
8825: . mctx           - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired)
8826: - monitordestroy - [optional] routine that frees monitor context (may be `NULL`)

8828:   Options Database Key:
8829: . -dm_monitor_cancel - cancels all monitors that have been hardwired into a code by calls to `DMMonitorSet()`, but
8830:                             does not cancel those set via the options database.

8832:   Level: intermediate

8834:   Note:
8835:   Several different monitoring routines may be set by calling
8836:   `DMMonitorSet()` multiple times or with `DMMonitorSetFromOptions()`; all will be called in the
8837:   order in which they were set.

8839:   Fortran Note:
8840:   Only a single monitor function can be set for each `DM` object

8842:   Developer Note:
8843:   This API has a generic name but seems specific to a very particular aspect of the use of `DM`

8845: .seealso: [](ch_dmbase), `DM`, `DMMonitorCancel()`, `DMMonitorSetFromOptions()`, `DMMonitor()`
8846: @*/
8847: PetscErrorCode DMMonitorSet(DM dm, PetscErrorCode (*f)(DM, void *), void *mctx, PetscErrorCode (*monitordestroy)(void **))
8848: {
8849:   PetscInt m;

8851:   PetscFunctionBegin;
8853:   for (m = 0; m < dm->numbermonitors; ++m) {
8854:     PetscBool identical;

8856:     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, mctx, monitordestroy, (PetscErrorCode(*)(void))dm->monitor[m], dm->monitorcontext[m], dm->monitordestroy[m], &identical));
8857:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
8858:   }
8859:   PetscCheck(dm->numbermonitors < MAXDMMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
8860:   dm->monitor[dm->numbermonitors]          = f;
8861:   dm->monitordestroy[dm->numbermonitors]   = monitordestroy;
8862:   dm->monitorcontext[dm->numbermonitors++] = (void *)mctx;
8863:   PetscFunctionReturn(PETSC_SUCCESS);
8864: }

8866: /*@
8867:   DMMonitorCancel - Clears all the monitor functions for a `DM` object.

8869:   Logically Collective

8871:   Input Parameter:
8872: . dm - the DM

8874:   Options Database Key:
8875: . -dm_monitor_cancel - cancels all monitors that have been hardwired
8876:   into a code by calls to `DMonitorSet()`, but does not cancel those
8877:   set via the options database

8879:   Level: intermediate

8881:   Note:
8882:   There is no way to clear one specific monitor from a `DM` object.

8884: .seealso: [](ch_dmbase), `DM`, `DMMonitorSet()`, `DMMonitorSetFromOptions()`, `DMMonitor()`
8885: @*/
8886: PetscErrorCode DMMonitorCancel(DM dm)
8887: {
8888:   PetscInt m;

8890:   PetscFunctionBegin;
8892:   for (m = 0; m < dm->numbermonitors; ++m) {
8893:     if (dm->monitordestroy[m]) PetscCall((*dm->monitordestroy[m])(&dm->monitorcontext[m]));
8894:   }
8895:   dm->numbermonitors = 0;
8896:   PetscFunctionReturn(PETSC_SUCCESS);
8897: }

8899: /*@C
8900:   DMMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user

8902:   Collective

8904:   Input Parameters:
8905: + dm           - `DM` object you wish to monitor
8906: . name         - the monitor type one is seeking
8907: . help         - message indicating what monitoring is done
8908: . manual       - manual page for the monitor
8909: . monitor      - the monitor function
8910: - 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

8912:   Output Parameter:
8913: . flg - Flag set if the monitor was created

8915:   Level: developer

8917: .seealso: [](ch_dmbase), `DM`, `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
8918:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
8919:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
8920:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
8921:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
8922:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
8923:           `PetscOptionsFList()`, `PetscOptionsEList()`, `DMMonitor()`, `DMMonitorSet()`
8924: @*/
8925: PetscErrorCode DMMonitorSetFromOptions(DM dm, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(DM, void *), PetscErrorCode (*monitorsetup)(DM, PetscViewerAndFormat *), PetscBool *flg)
8926: {
8927:   PetscViewer       viewer;
8928:   PetscViewerFormat format;

8930:   PetscFunctionBegin;
8932:   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)dm), ((PetscObject)dm)->options, ((PetscObject)dm)->prefix, name, &viewer, &format, flg));
8933:   if (*flg) {
8934:     PetscViewerAndFormat *vf;

8936:     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
8937:     PetscCall(PetscOptionsRestoreViewer(&viewer));
8938:     if (monitorsetup) PetscCall((*monitorsetup)(dm, vf));
8939:     PetscCall(DMMonitorSet(dm, (PetscErrorCode(*)(DM, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
8940:   }
8941:   PetscFunctionReturn(PETSC_SUCCESS);
8942: }

8944: /*@
8945:   DMMonitor - runs the user provided monitor routines, if they exist

8947:   Collective

8949:   Input Parameter:
8950: . dm - The `DM`

8952:   Level: developer

8954:   Developer Note:
8955:   Note should indicate when during the life of the `DM` the monitor is run. It appears to be
8956:   related to the discretization process seems rather specialized since some `DM` have no
8957:   concept of discretization.

8959: .seealso: [](ch_dmbase), `DM`, `DMMonitorSet()`, `DMMonitorSetFromOptions()`
8960: @*/
8961: PetscErrorCode DMMonitor(DM dm)
8962: {
8963:   PetscInt m;

8965:   PetscFunctionBegin;
8966:   if (!dm) PetscFunctionReturn(PETSC_SUCCESS);
8968:   for (m = 0; m < dm->numbermonitors; ++m) PetscCall((*dm->monitor[m])(dm, dm->monitorcontext[m]));
8969:   PetscFunctionReturn(PETSC_SUCCESS);
8970: }

8972: /*@
8973:   DMComputeError - Computes the error assuming the user has provided the exact solution functions

8975:   Collective

8977:   Input Parameters:
8978: + dm  - The `DM`
8979: - sol - The solution vector

8981:   Input/Output Parameter:
8982: . errors - An array of length Nf, the number of fields, or `NULL` for no output; on output
8983:            contains the error in each field

8985:   Output Parameter:
8986: . errorVec - A vector to hold the cellwise error (may be `NULL`)

8988:   Level: developer

8990:   Note:
8991:   The exact solutions come from the `PetscDS` object, and the time comes from `DMGetOutputSequenceNumber()`.

8993: .seealso: [](ch_dmbase), `DM`, `DMMonitorSet()`, `DMGetRegionNumDS()`, `PetscDSGetExactSolution()`, `DMGetOutputSequenceNumber()`
8994: @*/
8995: PetscErrorCode DMComputeError(DM dm, Vec sol, PetscReal errors[], Vec *errorVec)
8996: {
8997:   PetscErrorCode (**exactSol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
8998:   void    **ctxs;
8999:   PetscReal time;
9000:   PetscInt  Nf, f, Nds, s;

9002:   PetscFunctionBegin;
9003:   PetscCall(DMGetNumFields(dm, &Nf));
9004:   PetscCall(PetscCalloc2(Nf, &exactSol, Nf, &ctxs));
9005:   PetscCall(DMGetNumDS(dm, &Nds));
9006:   for (s = 0; s < Nds; ++s) {
9007:     PetscDS         ds;
9008:     DMLabel         label;
9009:     IS              fieldIS;
9010:     const PetscInt *fields;
9011:     PetscInt        dsNf;

9013:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL));
9014:     PetscCall(PetscDSGetNumFields(ds, &dsNf));
9015:     if (fieldIS) PetscCall(ISGetIndices(fieldIS, &fields));
9016:     for (f = 0; f < dsNf; ++f) {
9017:       const PetscInt field = fields[f];
9018:       PetscCall(PetscDSGetExactSolution(ds, field, &exactSol[field], &ctxs[field]));
9019:     }
9020:     if (fieldIS) PetscCall(ISRestoreIndices(fieldIS, &fields));
9021:   }
9022:   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);
9023:   PetscCall(DMGetOutputSequenceNumber(dm, NULL, &time));
9024:   if (errors) PetscCall(DMComputeL2FieldDiff(dm, time, exactSol, ctxs, sol, errors));
9025:   if (errorVec) {
9026:     DM             edm;
9027:     DMPolytopeType ct;
9028:     PetscBool      simplex;
9029:     PetscInt       dim, cStart, Nf;

9031:     PetscCall(DMClone(dm, &edm));
9032:     PetscCall(DMGetDimension(edm, &dim));
9033:     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
9034:     PetscCall(DMPlexGetCellType(dm, cStart, &ct));
9035:     simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
9036:     PetscCall(DMGetNumFields(dm, &Nf));
9037:     for (f = 0; f < Nf; ++f) {
9038:       PetscFE         fe, efe;
9039:       PetscQuadrature q;
9040:       const char     *name;

9042:       PetscCall(DMGetField(dm, f, NULL, (PetscObject *)&fe));
9043:       PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, Nf, simplex, 0, PETSC_DETERMINE, &efe));
9044:       PetscCall(PetscObjectGetName((PetscObject)fe, &name));
9045:       PetscCall(PetscObjectSetName((PetscObject)efe, name));
9046:       PetscCall(PetscFEGetQuadrature(fe, &q));
9047:       PetscCall(PetscFESetQuadrature(efe, q));
9048:       PetscCall(DMSetField(edm, f, NULL, (PetscObject)efe));
9049:       PetscCall(PetscFEDestroy(&efe));
9050:     }
9051:     PetscCall(DMCreateDS(edm));

9053:     PetscCall(DMCreateGlobalVector(edm, errorVec));
9054:     PetscCall(PetscObjectSetName((PetscObject)*errorVec, "Error"));
9055:     PetscCall(DMPlexComputeL2DiffVec(dm, time, exactSol, ctxs, sol, *errorVec));
9056:     PetscCall(DMDestroy(&edm));
9057:   }
9058:   PetscCall(PetscFree2(exactSol, ctxs));
9059:   PetscFunctionReturn(PETSC_SUCCESS);
9060: }

9062: /*@
9063:   DMGetNumAuxiliaryVec - Get the number of auxiliary vectors associated with this `DM`

9065:   Not Collective

9067:   Input Parameter:
9068: . dm - The `DM`

9070:   Output Parameter:
9071: . numAux - The number of auxiliary data vectors

9073:   Level: advanced

9075: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMSetAuxiliaryVec()`, `DMGetAuxiliaryLabels()`, `DMGetAuxiliaryVec()`
9076: @*/
9077: PetscErrorCode DMGetNumAuxiliaryVec(DM dm, PetscInt *numAux)
9078: {
9079:   PetscFunctionBegin;
9081:   PetscCall(PetscHMapAuxGetSize(dm->auxData, numAux));
9082:   PetscFunctionReturn(PETSC_SUCCESS);
9083: }

9085: /*@
9086:   DMGetAuxiliaryVec - Get the auxiliary vector for region specified by the given label and value, and equation part

9088:   Not Collective

9090:   Input Parameters:
9091: + dm    - The `DM`
9092: . label - The `DMLabel`
9093: . value - The label value indicating the region
9094: - part  - The equation part, or 0 if unused

9096:   Output Parameter:
9097: . aux - The `Vec` holding auxiliary field data

9099:   Level: advanced

9101:   Note:
9102:   If no auxiliary vector is found for this (label, value), (NULL, 0, 0) is checked as well.

9104: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMSetAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryLabels()`
9105: @*/
9106: PetscErrorCode DMGetAuxiliaryVec(DM dm, DMLabel label, PetscInt value, PetscInt part, Vec *aux)
9107: {
9108:   PetscHashAuxKey key, wild = {NULL, 0, 0};
9109:   PetscBool       has;

9111:   PetscFunctionBegin;
9114:   key.label = label;
9115:   key.value = value;
9116:   key.part  = part;
9117:   PetscCall(PetscHMapAuxHas(dm->auxData, key, &has));
9118:   if (has) PetscCall(PetscHMapAuxGet(dm->auxData, key, aux));
9119:   else PetscCall(PetscHMapAuxGet(dm->auxData, wild, aux));
9120:   PetscFunctionReturn(PETSC_SUCCESS);
9121: }

9123: /*@
9124:   DMSetAuxiliaryVec - Set an auxiliary vector for region specified by the given label and value, and equation part

9126:   Not Collective because auxiliary vectors are not parallel

9128:   Input Parameters:
9129: + dm    - The `DM`
9130: . label - The `DMLabel`
9131: . value - The label value indicating the region
9132: . part  - The equation part, or 0 if unused
9133: - aux   - The `Vec` holding auxiliary field data

9135:   Level: advanced

9137: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMGetAuxiliaryLabels()`, `DMCopyAuxiliaryVec()`
9138: @*/
9139: PetscErrorCode DMSetAuxiliaryVec(DM dm, DMLabel label, PetscInt value, PetscInt part, Vec aux)
9140: {
9141:   Vec             old;
9142:   PetscHashAuxKey key;

9144:   PetscFunctionBegin;
9147:   key.label = label;
9148:   key.value = value;
9149:   key.part  = part;
9150:   PetscCall(PetscHMapAuxGet(dm->auxData, key, &old));
9151:   PetscCall(PetscObjectReference((PetscObject)aux));
9152:   if (!aux) PetscCall(PetscHMapAuxDel(dm->auxData, key));
9153:   else PetscCall(PetscHMapAuxSet(dm->auxData, key, aux));
9154:   PetscCall(VecDestroy(&old));
9155:   PetscFunctionReturn(PETSC_SUCCESS);
9156: }

9158: /*@C
9159:   DMGetAuxiliaryLabels - Get the labels, values, and parts for all auxiliary vectors in this `DM`

9161:   Not Collective

9163:   Input Parameter:
9164: . dm - The `DM`

9166:   Output Parameters:
9167: + labels - The `DMLabel`s for each `Vec`
9168: . values - The label values for each `Vec`
9169: - parts  - The equation parts for each `Vec`

9171:   Level: advanced

9173:   Note:
9174:   The arrays passed in must be at least as large as `DMGetNumAuxiliaryVec()`.

9176: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMSetAuxiliaryVec()`, `DMCopyAuxiliaryVec()`
9177: @*/
9178: PetscErrorCode DMGetAuxiliaryLabels(DM dm, DMLabel labels[], PetscInt values[], PetscInt parts[])
9179: {
9180:   PetscHashAuxKey *keys;
9181:   PetscInt         n, i, off = 0;

9183:   PetscFunctionBegin;
9185:   PetscAssertPointer(labels, 2);
9186:   PetscAssertPointer(values, 3);
9187:   PetscAssertPointer(parts, 4);
9188:   PetscCall(DMGetNumAuxiliaryVec(dm, &n));
9189:   PetscCall(PetscMalloc1(n, &keys));
9190:   PetscCall(PetscHMapAuxGetKeys(dm->auxData, &off, keys));
9191:   for (i = 0; i < n; ++i) {
9192:     labels[i] = keys[i].label;
9193:     values[i] = keys[i].value;
9194:     parts[i]  = keys[i].part;
9195:   }
9196:   PetscCall(PetscFree(keys));
9197:   PetscFunctionReturn(PETSC_SUCCESS);
9198: }

9200: /*@
9201:   DMCopyAuxiliaryVec - Copy the auxiliary vector data on a `DM` to a new `DM`

9203:   Not Collective

9205:   Input Parameter:
9206: . dm - The `DM`

9208:   Output Parameter:
9209: . dmNew - The new `DM`, now with the same auxiliary data

9211:   Level: advanced

9213:   Note:
9214:   This is a shallow copy of the auxiliary vectors

9216: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMSetAuxiliaryVec()`
9217: @*/
9218: PetscErrorCode DMCopyAuxiliaryVec(DM dm, DM dmNew)
9219: {
9220:   PetscFunctionBegin;
9223:   if (dm == dmNew) PetscFunctionReturn(PETSC_SUCCESS);
9224:   PetscCall(DMClearAuxiliaryVec(dmNew));

9226:   PetscCall(PetscHMapAuxDestroy(&dmNew->auxData));
9227:   PetscCall(PetscHMapAuxDuplicate(dm->auxData, &dmNew->auxData));
9228:   {
9229:     Vec     *auxData;
9230:     PetscInt n, i, off = 0;

9232:     PetscCall(PetscHMapAuxGetSize(dmNew->auxData, &n));
9233:     PetscCall(PetscMalloc1(n, &auxData));
9234:     PetscCall(PetscHMapAuxGetVals(dmNew->auxData, &off, auxData));
9235:     for (i = 0; i < n; ++i) PetscCall(PetscObjectReference((PetscObject)auxData[i]));
9236:     PetscCall(PetscFree(auxData));
9237:   }
9238:   PetscFunctionReturn(PETSC_SUCCESS);
9239: }

9241: /*@
9242:   DMClearAuxiliaryVec - Destroys the auxiliary vector information and creates a new empty one

9244:   Not Collective

9246:   Input Parameter:
9247: . dm - The `DM`

9249:   Level: advanced

9251: .seealso: [](ch_dmbase), `DM`, `DMCopyAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMSetAuxiliaryVec()`
9252: @*/
9253: PetscErrorCode DMClearAuxiliaryVec(DM dm)
9254: {
9255:   Vec     *auxData;
9256:   PetscInt n, i, off = 0;

9258:   PetscFunctionBegin;
9259:   PetscCall(PetscHMapAuxGetSize(dm->auxData, &n));
9260:   PetscCall(PetscMalloc1(n, &auxData));
9261:   PetscCall(PetscHMapAuxGetVals(dm->auxData, &off, auxData));
9262:   for (i = 0; i < n; ++i) PetscCall(VecDestroy(&auxData[i]));
9263:   PetscCall(PetscFree(auxData));
9264:   PetscCall(PetscHMapAuxDestroy(&dm->auxData));
9265:   PetscCall(PetscHMapAuxCreate(&dm->auxData));
9266:   PetscFunctionReturn(PETSC_SUCCESS);
9267: }

9269: /*@C
9270:   DMPolytopeMatchOrientation - Determine an orientation (transformation) that takes the source face arrangement to the target face arrangement

9272:   Not Collective

9274:   Input Parameters:
9275: + ct         - The `DMPolytopeType`
9276: . sourceCone - The source arrangement of faces
9277: - targetCone - The target arrangement of faces

9279:   Output Parameters:
9280: + ornt  - The orientation (transformation) which will take the source arrangement to the target arrangement
9281: - found - Flag indicating that a suitable orientation was found

9283:   Level: advanced

9285:   Note:
9286:   An arrangement is a face order combined with an orientation for each face

9288:   Each orientation (transformation) is labeled with an integer from negative `DMPolytopeTypeGetNumArrangements(ct)`/2 to `DMPolytopeTypeGetNumArrangements(ct)`/2
9289:   that labels each arrangement (face ordering plus orientation for each face).

9291:   See `DMPolytopeMatchVertexOrientation()` to find a new vertex orientation that takes the source vertex arrangement to the target vertex arrangement

9293: .seealso: [](ch_dmbase), `DM`, `DMPolytopeGetOrientation()`, `DMPolytopeMatchVertexOrientation()`, `DMPolytopeGetVertexOrientation()`
9294: @*/
9295: PetscErrorCode DMPolytopeMatchOrientation(DMPolytopeType ct, const PetscInt sourceCone[], const PetscInt targetCone[], PetscInt *ornt, PetscBool *found)
9296: {
9297:   const PetscInt cS = DMPolytopeTypeGetConeSize(ct);
9298:   const PetscInt nO = DMPolytopeTypeGetNumArrangements(ct) / 2;
9299:   PetscInt       o, c;

9301:   PetscFunctionBegin;
9302:   if (!nO) {
9303:     *ornt  = 0;
9304:     *found = PETSC_TRUE;
9305:     PetscFunctionReturn(PETSC_SUCCESS);
9306:   }
9307:   for (o = -nO; o < nO; ++o) {
9308:     const PetscInt *arr = DMPolytopeTypeGetArrangement(ct, o);

9310:     for (c = 0; c < cS; ++c)
9311:       if (sourceCone[arr[c * 2]] != targetCone[c]) break;
9312:     if (c == cS) {
9313:       *ornt = o;
9314:       break;
9315:     }
9316:   }
9317:   *found = o == nO ? PETSC_FALSE : PETSC_TRUE;
9318:   PetscFunctionReturn(PETSC_SUCCESS);
9319: }

9321: /*@C
9322:   DMPolytopeGetOrientation - Determine an orientation (transformation) that takes the source face arrangement to the target face arrangement

9324:   Not Collective

9326:   Input Parameters:
9327: + ct         - The `DMPolytopeType`
9328: . sourceCone - The source arrangement of faces
9329: - targetCone - The target arrangement of faces

9331:   Output Parameter:
9332: . ornt - The orientation (transformation) which will take the source arrangement to the target arrangement

9334:   Level: advanced

9336:   Note:
9337:   This function is the same as `DMPolytopeMatchOrientation()` except it will generate an error if no suitable orientation can be found.

9339:   Developer Note:
9340:   It is unclear why this function needs to exist since one can simply call `DMPolytopeMatchOrientation()` and error if none is found

9342: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMPolytopeMatchOrientation()`, `DMPolytopeGetVertexOrientation()`, `DMPolytopeMatchVertexOrientation()`
9343: @*/
9344: PetscErrorCode DMPolytopeGetOrientation(DMPolytopeType ct, const PetscInt sourceCone[], const PetscInt targetCone[], PetscInt *ornt)
9345: {
9346:   PetscBool found;

9348:   PetscFunctionBegin;
9349:   PetscCall(DMPolytopeMatchOrientation(ct, sourceCone, targetCone, ornt, &found));
9350:   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find orientation for %s", DMPolytopeTypes[ct]);
9351:   PetscFunctionReturn(PETSC_SUCCESS);
9352: }

9354: /*@C
9355:   DMPolytopeMatchVertexOrientation - Determine an orientation (transformation) that takes the source vertex arrangement to the target vertex arrangement

9357:   Not Collective

9359:   Input Parameters:
9360: + ct         - The `DMPolytopeType`
9361: . sourceVert - The source arrangement of vertices
9362: - targetVert - The target arrangement of vertices

9364:   Output Parameters:
9365: + ornt  - The orientation (transformation) which will take the source arrangement to the target arrangement
9366: - found - Flag indicating that a suitable orientation was found

9368:   Level: advanced

9370:   Notes:
9371:   An arrangement is a vertex order

9373:   Each orientation (transformation) is labeled with an integer from negative `DMPolytopeTypeGetNumArrangements(ct)`/2 to `DMPolytopeTypeGetNumArrangements(ct)`/2
9374:   that labels each arrangement (vertex ordering).

9376:   See `DMPolytopeMatchOrientation()` to find a new face orientation that takes the source face arrangement to the target face arrangement

9378: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMPolytopeGetOrientation()`, `DMPolytopeMatchOrientation()`, `DMPolytopeTypeGetNumVertices()`, `DMPolytopeTypeGetVertexArrangement()`
9379: @*/
9380: PetscErrorCode DMPolytopeMatchVertexOrientation(DMPolytopeType ct, const PetscInt sourceVert[], const PetscInt targetVert[], PetscInt *ornt, PetscBool *found)
9381: {
9382:   const PetscInt cS = DMPolytopeTypeGetNumVertices(ct);
9383:   const PetscInt nO = DMPolytopeTypeGetNumArrangements(ct) / 2;
9384:   PetscInt       o, c;

9386:   PetscFunctionBegin;
9387:   if (!nO) {
9388:     *ornt  = 0;
9389:     *found = PETSC_TRUE;
9390:     PetscFunctionReturn(PETSC_SUCCESS);
9391:   }
9392:   for (o = -nO; o < nO; ++o) {
9393:     const PetscInt *arr = DMPolytopeTypeGetVertexArrangement(ct, o);

9395:     for (c = 0; c < cS; ++c)
9396:       if (sourceVert[arr[c]] != targetVert[c]) break;
9397:     if (c == cS) {
9398:       *ornt = o;
9399:       break;
9400:     }
9401:   }
9402:   *found = o == nO ? PETSC_FALSE : PETSC_TRUE;
9403:   PetscFunctionReturn(PETSC_SUCCESS);
9404: }

9406: /*@C
9407:   DMPolytopeGetVertexOrientation - Determine an orientation (transformation) that takes the source vertex arrangement to the target vertex arrangement

9409:   Not Collective

9411:   Input Parameters:
9412: + ct         - The `DMPolytopeType`
9413: . sourceCone - The source arrangement of vertices
9414: - targetCone - The target arrangement of vertices

9416:   Output Parameter:
9417: . ornt - The orientation (transformation) which will take the source arrangement to the target arrangement

9419:   Level: advanced

9421:   Note:
9422:   This function is the same as `DMPolytopeMatchVertexOrientation()` except it errors if not orientation is possible.

9424:   Developer Note:
9425:   It is unclear why this function needs to exist since one can simply call `DMPolytopeMatchVertexOrientation()` and error if none is found

9427: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMPolytopeMatchVertexOrientation()`, `DMPolytopeGetOrientation()`
9428: @*/
9429: PetscErrorCode DMPolytopeGetVertexOrientation(DMPolytopeType ct, const PetscInt sourceCone[], const PetscInt targetCone[], PetscInt *ornt)
9430: {
9431:   PetscBool found;

9433:   PetscFunctionBegin;
9434:   PetscCall(DMPolytopeMatchVertexOrientation(ct, sourceCone, targetCone, ornt, &found));
9435:   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find orientation for %s", DMPolytopeTypes[ct]);
9436:   PetscFunctionReturn(PETSC_SUCCESS);
9437: }

9439: /*@C
9440:   DMPolytopeInCellTest - Check whether a point lies inside the reference cell of given type

9442:   Not Collective

9444:   Input Parameters:
9445: + ct    - The `DMPolytopeType`
9446: - point - Coordinates of the point

9448:   Output Parameter:
9449: . inside - Flag indicating whether the point is inside the reference cell of given type

9451:   Level: advanced

9453: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMLocatePoints()`
9454: @*/
9455: PetscErrorCode DMPolytopeInCellTest(DMPolytopeType ct, const PetscReal point[], PetscBool *inside)
9456: {
9457:   PetscReal sum = 0.0;
9458:   PetscInt  d;

9460:   PetscFunctionBegin;
9461:   *inside = PETSC_TRUE;
9462:   switch (ct) {
9463:   case DM_POLYTOPE_TRIANGLE:
9464:   case DM_POLYTOPE_TETRAHEDRON:
9465:     for (d = 0; d < DMPolytopeTypeGetDim(ct); ++d) {
9466:       if (point[d] < -1.0) {
9467:         *inside = PETSC_FALSE;
9468:         break;
9469:       }
9470:       sum += point[d];
9471:     }
9472:     if (sum > PETSC_SMALL) {
9473:       *inside = PETSC_FALSE;
9474:       break;
9475:     }
9476:     break;
9477:   case DM_POLYTOPE_QUADRILATERAL:
9478:   case DM_POLYTOPE_HEXAHEDRON:
9479:     for (d = 0; d < DMPolytopeTypeGetDim(ct); ++d)
9480:       if (PetscAbsReal(point[d]) > 1. + PETSC_SMALL) {
9481:         *inside = PETSC_FALSE;
9482:         break;
9483:       }
9484:     break;
9485:   default:
9486:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported polytope type %s", DMPolytopeTypes[ct]);
9487:   }
9488:   PetscFunctionReturn(PETSC_SUCCESS);
9489: }

9491: /*@
9492:   DMReorderSectionSetDefault - Set flag indicating whether the local section should be reordered by default

9494:   Logically collective

9496:   Input Parameters:
9497: + dm      - The DM
9498: - reorder - Flag for reordering

9500:   Level: intermediate

9502: .seealso: `DMReorderSectionGetDefault()`
9503: @*/
9504: PetscErrorCode DMReorderSectionSetDefault(DM dm, DMReorderDefaultFlag reorder)
9505: {
9506:   PetscFunctionBegin;
9508:   PetscTryMethod(dm, "DMReorderSectionSetDefault_C", (DM, DMReorderDefaultFlag), (dm, reorder));
9509:   PetscFunctionReturn(PETSC_SUCCESS);
9510: }

9512: /*@
9513:   DMReorderSectionGetDefault - Get flag indicating whether the local section should be reordered by default

9515:   Not collective

9517:   Input Parameter:
9518: . dm - The DM

9520:   Output Parameter:
9521: . reorder - Flag for reordering

9523:   Level: intermediate

9525: .seealso: `DMReorderSetDefault()`
9526: @*/
9527: PetscErrorCode DMReorderSectionGetDefault(DM dm, DMReorderDefaultFlag *reorder)
9528: {
9529:   PetscFunctionBegin;
9531:   PetscAssertPointer(reorder, 2);
9532:   *reorder = DM_REORDER_DEFAULT_NOTSET;
9533:   PetscTryMethod(dm, "DMReorderSectionGetDefault_C", (DM, DMReorderDefaultFlag *), (dm, reorder));
9534:   PetscFunctionReturn(PETSC_SUCCESS);
9535: }

9537: /*@C
9538:   DMReorderSectionSetType - Set the type of local section reordering

9540:   Logically collective

9542:   Input Parameters:
9543: + dm      - The DM
9544: - reorder - The reordering method

9546:   Level: intermediate

9548: .seealso: `DMReorderSectionGetType()`, `DMReorderSectionSetDefault()`
9549: @*/
9550: PetscErrorCode DMReorderSectionSetType(DM dm, MatOrderingType reorder)
9551: {
9552:   PetscFunctionBegin;
9554:   PetscTryMethod(dm, "DMReorderSectionSetType_C", (DM, MatOrderingType), (dm, reorder));
9555:   PetscFunctionReturn(PETSC_SUCCESS);
9556: }

9558: /*@C
9559:   DMReorderSectionGetType - Get the reordering type for the local section

9561:   Not collective

9563:   Input Parameter:
9564: . dm - The DM

9566:   Output Parameter:
9567: . reorder - The reordering method

9569:   Level: intermediate

9571: .seealso: `DMReorderSetDefault()`, `DMReorderSectionGetDefault()`
9572: @*/
9573: PetscErrorCode DMReorderSectionGetType(DM dm, MatOrderingType *reorder)
9574: {
9575:   PetscFunctionBegin;
9577:   PetscAssertPointer(reorder, 2);
9578:   *reorder = NULL;
9579:   PetscTryMethod(dm, "DMReorderSectionGetType_C", (DM, MatOrderingType *), (dm, reorder));
9580:   PetscFunctionReturn(PETSC_SUCCESS);
9581: }