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 (i) PetscCall(DMSetCellCoordinateDM(*newdm, ncdm));
183:         else PetscCall(DMSetCoordinateDM(*newdm, ncdm));
184:         PetscCall(DMDestroy(&ncdm));
185:       }
186:     }
187:   }
188:   PetscCall(DMGetCoordinateDim(dm, &cdim));
189:   PetscCall(DMSetCoordinateDim(*newdm, cdim));
190:   PetscCall(DMGetCoordinatesLocal(dm, &coords));
191:   if (coords) {
192:     PetscCall(DMSetCoordinatesLocal(*newdm, coords));
193:   } else {
194:     PetscCall(DMGetCoordinates(dm, &coords));
195:     if (coords) PetscCall(DMSetCoordinates(*newdm, coords));
196:   }
197:   PetscCall(DMGetCellCoordinatesLocal(dm, &coords));
198:   if (coords) {
199:     PetscCall(DMSetCellCoordinatesLocal(*newdm, coords));
200:   } else {
201:     PetscCall(DMGetCellCoordinates(dm, &coords));
202:     if (coords) PetscCall(DMSetCellCoordinates(*newdm, coords));
203:   }
204:   {
205:     const PetscReal *maxCell, *Lstart, *L;

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

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

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

222:   Logically Collective

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

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

231:   Level: intermediate

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

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

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

252:   Logically Collective

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

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

260:   Level: intermediate

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

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

275:   Not Collective

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

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

283:   Level: intermediate

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

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

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

302:   Not Collective

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

308:   Level: developer

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

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

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

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

329:   Logically Collective

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

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

338:   Level: intermediate

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

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

354:   Logically Collective

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

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

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

365:   Level: intermediate

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

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

381:   Logically Collective

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

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

390:   Level: intermediate

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

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

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

410:   Logically Collective

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

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

418:   Level: intermediate

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

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

433:   Not Collective

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

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

441:   Level: intermediate

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

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

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

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

463:   Not Collective

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

469:   Level: developer

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

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

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

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

492:   Logically Collective

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

498:   Level: advanced

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

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

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

520:   Logically Collective

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

526:   Level: advanced

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

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

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

547:   Not Collective

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

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

555:   Level: advanced

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

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

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

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

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

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

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

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

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

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

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

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

640:   Collective

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

645:   Level: developer

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

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

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

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

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

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

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

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

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

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

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

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

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

800:   Collective

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

805:   Level: intermediate

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

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

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

825:   Collective

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

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

883:   Level: intermediate

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

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

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

920:   Collective

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

927:   Level: intermediate

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

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

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

946:   Collective

948:   Input Parameters:
949: + dm - the `DM` object to view
950: - v  - the viewer

952:   Level: beginner

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

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

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

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

990:     PetscCall(PetscViewerBinaryWrite(v, &classid, 1, PETSC_INT));
991:     PetscCall(PetscStrncpy(type, ((PetscObject)dm)->type_name, sizeof(type)));
992:     PetscCall(PetscViewerBinaryWrite(v, type, 256, PETSC_CHAR));
993:   }
994:   PetscTryTypeMethod(dm, view, v);
995:   PetscFunctionReturn(PETSC_SUCCESS);
996: }

998: /*@
999:   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,
1000:   that is it has no ghost locations.

1002:   Collective

1004:   Input Parameter:
1005: . dm - the `DM` object

1007:   Output Parameter:
1008: . vec - the global vector

1010:   Level: beginner

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

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

1030: /*@
1031:   DMCreateLocalVector - Creates a local vector from a `DM` object.

1033:   Not Collective

1035:   Input Parameter:
1036: . dm - the `DM` object

1038:   Output Parameter:
1039: . vec - the local vector

1041:   Level: beginner

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

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

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

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

1067:   Collective

1069:   Input Parameter:
1070: . dm - the `DM` that provides the mapping

1072:   Output Parameter:
1073: . ltog - the mapping

1075:   Level: advanced

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

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

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

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

1092:   PetscFunctionBegin;
1094:   PetscAssertPointer(ltog, 2);
1095:   if (!dm->ltogmap) {
1096:     PetscSection section, sectionGlobal;

1098:     PetscCall(DMGetLocalSection(dm, &section));
1099:     if (section) {
1100:       const PetscInt *cdofs;
1101:       PetscInt       *ltog;
1102:       PetscInt        pStart, pEnd, n, p, k, l;

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

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

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

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

1157:   Not Collective

1159:   Input Parameter:
1160: . dm - the `DM` with block structure

1162:   Output Parameter:
1163: . bs - the block size, 1 implies no exploitable block structure

1165:   Level: intermediate

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

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

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

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

1189:   Collective

1191:   Input Parameters:
1192: + dmc - the `DM` object
1193: - dmf - the second, finer `DM` object

1195:   Output Parameters:
1196: + mat - the interpolation
1197: - vec - the scaling (optional), see `DMCreateInterpolationScale()`

1199:   Level: developer

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

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

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

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

1226:   Input Parameters:
1227: + dac - `DM` that defines a coarse mesh
1228: . daf - `DM` that defines a fine mesh
1229: - mat - the restriction (or interpolation operator) from fine to coarse

1231:   Output Parameter:
1232: . scale - the scaled vector

1234:   Level: advanced

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

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

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

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

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

1279:   Collective

1281:   Input Parameters:
1282: + dmc - the `DM` object
1283: - dmf - the second, finer `DM` object

1285:   Output Parameter:
1286: . mat - the restriction

1288:   Level: developer

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

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

1308: /*@
1309:   DMCreateInjection - Gets injection matrix between two `DM` objects.

1311:   Collective

1313:   Input Parameters:
1314: + dac - the `DM` object
1315: - daf - the second, finer `DM` object

1317:   Output Parameter:
1318: . mat - the injection

1320:   Level: developer

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

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

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

1348: /*@
1349:   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
1350:   a Galerkin finite element model on the `DM`

1352:   Collective

1354:   Input Parameters:
1355: + dmc - the target `DM` object
1356: - dmf - the source `DM` object

1358:   Output Parameter:
1359: . mat - the mass matrix

1361:   Level: developer

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

1366:   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()`

1368: .seealso: [](ch_dmbase), `DM`, `DMCreateMassMatrixLumped()`, `DMCreateMatrix()`, `DMRefine()`, `DMCoarsen()`, `DMCreateRestriction()`, `DMCreateInterpolation()`, `DMCreateInjection()`
1369: @*/
1370: PetscErrorCode DMCreateMassMatrix(DM dmc, DM dmf, Mat *mat)
1371: {
1372:   PetscFunctionBegin;
1375:   PetscAssertPointer(mat, 3);
1376:   PetscCall(PetscLogEventBegin(DM_CreateMassMatrix, 0, 0, 0, 0));
1377:   PetscUseTypeMethod(dmc, createmassmatrix, dmf, mat);
1378:   PetscCall(PetscLogEventEnd(DM_CreateMassMatrix, 0, 0, 0, 0));
1379:   PetscFunctionReturn(PETSC_SUCCESS);
1380: }

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

1385:   Collective

1387:   Input Parameter:
1388: . dm - the `DM` object

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

1393:   Level: developer

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

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

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

1413:   Collective

1415:   Input Parameters:
1416: + dm    - the `DM` object
1417: - ctype - `IS_COLORING_LOCAL` or `IS_COLORING_GLOBAL`

1419:   Output Parameter:
1420: . coloring - the coloring

1422:   Level: developer

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

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

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

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

1446:   Collective

1448:   Input Parameter:
1449: . dm - the `DM` object

1451:   Output Parameter:
1452: . mat - the empty Jacobian

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

1457:   Level: beginner

1459:   Notes:
1460:   This properly preallocates the number of nonzeros in the sparse matrix so you
1461:   do not need to do it yourself.

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

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

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

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

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

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

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

1519:   Logically Collective

1521:   Input Parameters:
1522: + dm   - the `DM`
1523: - skip - `PETSC_TRUE` to skip preallocation

1525:   Level: developer

1527:   Note:
1528:   This is most useful to reduce initialization costs when `MatSetPreallocationCOO()` and
1529:   `MatSetValuesCOO()` will be used.

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

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

1545:   Logically Collective

1547:   Input Parameters:
1548: + dm   - the `DM`
1549: - only - `PETSC_TRUE` if only want preallocation

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

1554:   Level: developer

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

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

1570:   Logically Collective

1572:   Input Parameters:
1573: + dm   - the `DM`
1574: - only - `PETSC_TRUE` if you only want matrix structure

1576:   Level: developer

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

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

1591:   Logically Collective

1593:   Input Parameters:
1594: + dm    - the `DM`
1595: - btype - block by topological point or field node

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

1600:   Level: advanced

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

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

1615:   Not Collective

1617:   Input Parameter:
1618: . dm - the `DM`

1620:   Output Parameter:
1621: . btype - block by topological point or field node

1623:   Level: advanced

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

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

1639:   Not Collective

1641:   Input Parameters:
1642: + dm    - the `DM` object
1643: . count - The minimum size
1644: - dtype - MPI data type, often `MPIU_REAL`, `MPIU_SCALAR`, or `MPIU_INT`)

1646:   Output Parameter:
1647: . mem - the work array

1649:   Level: developer

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

1654:   The array may contain nonzero values

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

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

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

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

1706:   Not Collective

1708:   Input Parameters:
1709: + dm    - the `DM` object
1710: . count - The minimum size
1711: - dtype - MPI data type, often `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_INT`

1713:   Output Parameter:
1714: . mem - the work array

1716:   Level: developer

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

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

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

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

1749:   Logically Collective; No Fortran Support

1751:   Input Parameters:
1752: + dm     - The `DM`
1753: . field  - The field number for the nullspace
1754: - nullsp - A callback to create the nullspace

1756:   Calling sequence of `nullsp`:
1757: + dm        - The present `DM`
1758: . origField - The field number given above, in the original `DM`
1759: . field     - The field number in dm
1760: - nullSpace - The nullspace for the given field

1762:   Level: intermediate

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

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

1778:   Not Collective; No Fortran Support

1780:   Input Parameters:
1781: + dm    - The `DM`
1782: - field - The field number for the nullspace

1784:   Output Parameter:
1785: . nullsp - A callback to create the nullspace

1787:   Calling sequence of `nullsp`:
1788: + dm        - The present DM
1789: . origField - The field number given above, in the original DM
1790: . field     - The field number in dm
1791: - nullSpace - The nullspace for the given field

1793:   Level: intermediate

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

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

1810:   Logically Collective; No Fortran Support

1812:   Input Parameters:
1813: + dm     - The `DM`
1814: . field  - The field number for the nullspace
1815: - nullsp - A callback to create the near-nullspace

1817:   Calling sequence of `nullsp`:
1818: + dm        - The present `DM`
1819: . origField - The field number given above, in the original `DM`
1820: . field     - The field number in dm
1821: - nullSpace - The nullspace for the given field

1823:   Level: intermediate

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

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

1840:   Not Collective; No Fortran Support

1842:   Input Parameters:
1843: + dm    - The `DM`
1844: - field - The field number for the nullspace

1846:   Output Parameter:
1847: . nullsp - A callback to create the near-nullspace

1849:   Calling sequence of `nullsp`:
1850: + dm        - The present `DM`
1851: . origField - The field number given above, in the original `DM`
1852: . field     - The field number in dm
1853: - nullSpace - The nullspace for the given field

1855:   Level: intermediate

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

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

1873:   Not Collective; No Fortran Support

1875:   Input Parameter:
1876: . dm - the `DM` object

1878:   Output Parameters:
1879: + numFields  - The number of fields (or `NULL` if not requested)
1880: . fieldNames - The number of each field (or `NULL` if not requested)
1881: - fields     - The global indices for each field (or `NULL` if not requested)

1883:   Level: intermediate

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

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

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

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

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

1931:       PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof));
1932:       if (gdof > 0) {
1933:         for (f = 0; f < nF; ++f) {
1934:           PetscInt fdof, fcdof, fpdof;

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

1954:       PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof));
1955:       if (gdof > 0) {
1956:         PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
1957:         for (f = 0; f < nF; ++f) {
1958:           PetscInt fdof, fcdof, fc;

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

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

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

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

1998:   Not Collective; No Fortran Support

2000:   Input Parameter:
2001: . dm - the `DM` object

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

2009:   Level: intermediate

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

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

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

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

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

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

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

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

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

2088:   Not collective

2090:   Input Parameters:
2091: + dm        - The `DM` object
2092: . numFields - The number of fields to select
2093: - fields    - The field numbers of the selected fields

2095:   Output Parameters:
2096: + is    - The global indices for all the degrees of freedom in the new sub `DM`
2097: - subdm - The `DM` for the subproblem

2099:   Level: intermediate

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

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

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

2120:   Not collective

2122:   Input Parameters:
2123: + dms - The `DM` objects
2124: - n   - The number of `DM`s

2126:   Output Parameters:
2127: + is      - The global indices for each of subproblem within the super `DM`, or NULL
2128: - superdm - The `DM` for the superproblem

2130:   Level: intermediate

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

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

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

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

2159:   Not Collective

2161:   Input Parameter:
2162: . dm - the `DM` object

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

2171:   Level: intermediate

2173:   Notes:
2174:   Each `IS` contains the global indices of the dofs of the corresponding subdomains with in the
2175:   dofs of the original `DM`. The inner subdomains conceptually define a nonoverlapping
2176:   covering, while outer subdomains can overlap.

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

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

2184:   Developer Notes:
2185:   The `dmlist` is for the inner subdomains or the outer subdomains or all subdomains?

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

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

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

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

2245:   Not Collective

2247:   Input Parameters:
2248: + dm     - the `DM` object
2249: . n      - the number of subdomains
2250: - subdms - the local subdomains

2252:   Output Parameters:
2253: + iscat - scatter from global vector to nonoverlapping global vector entries on subdomain
2254: . oscat - scatter from global vector to overlapping global vector entries on subdomain
2255: - gscat - scatter from global vector to local vector on subdomain (fills in ghosts)

2257:   Level: developer

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

2265:   Developer Note:
2266:   Can the subdms input be anything or are they exactly the `DM` obtained from
2267:   `DMCreateDomainDecomposition()`?

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

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

2283:   Collective

2285:   Input Parameters:
2286: + dm   - the `DM` object
2287: - comm - the communicator to contain the new `DM` object (or `MPI_COMM_NULL`)

2289:   Output Parameter:
2290: . dmf - the refined `DM`, or `NULL`

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

2295:   Level: developer

2297:   Note:
2298:   If no refinement was done, the return value is `NULL`

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

2307:   PetscFunctionBegin;
2309:   PetscCall(PetscLogEventBegin(DM_Refine, dm, 0, 0, 0));
2310:   PetscUseTypeMethod(dm, refine, comm, dmf);
2311:   if (*dmf) {
2312:     (*dmf)->ops->creatematrix = dm->ops->creatematrix;

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

2316:     (*dmf)->ctx       = dm->ctx;
2317:     (*dmf)->leveldown = dm->leveldown;
2318:     (*dmf)->levelup   = dm->levelup + 1;

2320:     PetscCall(DMSetMatType(*dmf, dm->mattype));
2321:     for (link = dm->refinehook; link; link = link->next) {
2322:       if (link->refinehook) PetscCall((*link->refinehook)(dm, *dmf, link->ctx));
2323:     }
2324:   }
2325:   PetscCall(PetscLogEventEnd(DM_Refine, dm, 0, 0, 0));
2326:   PetscFunctionReturn(PETSC_SUCCESS);
2327: }

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

2332:   Logically Collective; No Fortran Support

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

2340:   Calling sequence of `refinehook`:
2341: + coarse - coarse level `DM`
2342: . fine   - fine level `DM` to interpolate problem to
2343: - ctx    - optional user-defined function context

2345:   Calling sequence of `interphook`:
2346: + coarse - coarse level `DM`
2347: . interp - matrix interpolating a coarse-level solution to the finer grid
2348: . fine   - fine level `DM` to update
2349: - ctx    - optional user-defined function context

2351:   Level: advanced

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

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

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

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

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

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

2385:   Logically Collective; No Fortran Support

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

2393:   Level: advanced

2395:   Note:
2396:   This function does nothing if the hook is not in the list.

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

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

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

2420:   Collective if any hooks are

2422:   Input Parameters:
2423: + coarse - coarser `DM` to use as a base
2424: . interp - interpolation matrix, apply using `MatInterpolate()`
2425: - fine   - finer `DM` to update

2427:   Level: developer

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

2433: .seealso: [](ch_dmbase), `DM`, `DMRefineHookAdd()`, `MatInterpolate()`
2434: @*/
2435: PetscErrorCode DMInterpolate(DM coarse, Mat interp, DM fine)
2436: {
2437:   DMRefineHookLink link;

2439:   PetscFunctionBegin;
2440:   for (link = fine->refinehook; link; link = link->next) {
2441:     if (link->interphook) PetscCall((*link->interphook)(coarse, interp, fine, link->ctx));
2442:   }
2443:   PetscFunctionReturn(PETSC_SUCCESS);
2444: }

2446: /*@
2447:   DMInterpolateSolution - Interpolates a solution from a coarse mesh to a fine mesh.

2449:   Collective

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

2459:   Output Parameter:
2460: . fineSol - the interpolation of coarseSol to the fine mesh

2462:   Level: developer

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

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

2473: .seealso: [](ch_dmbase), `DM`, `DMInterpolate()`, `DMCreateInterpolation()`
2474: @*/
2475: PetscErrorCode DMInterpolateSolution(DM coarse, DM fine, Mat interp, Vec coarseSol, Vec fineSol)
2476: {
2477:   PetscErrorCode (*interpsol)(DM, DM, Mat, Vec, Vec) = NULL;

2479:   PetscFunctionBegin;

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

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

2497:   Not Collective

2499:   Input Parameter:
2500: . dm - the `DM` object

2502:   Output Parameter:
2503: . level - number of refinements

2505:   Level: developer

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

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

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

2523:   Not Collective

2525:   Input Parameters:
2526: + dm    - the `DM` object
2527: - level - number of refinements

2529:   Level: advanced

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

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

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

2546: /*@
2547:   DMExtrude - Extrude a `DM` object from a surface

2549:   Collective

2551:   Input Parameters:
2552: + dm     - the `DM` object
2553: - layers - the number of extruded cell layers

2555:   Output Parameter:
2556: . dme - the extruded `DM`, or `NULL`

2558:   Level: developer

2560:   Note:
2561:   If no extrusion was done, the return value is `NULL`

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

2579: PetscErrorCode DMGetBasisTransformDM_Internal(DM dm, DM *tdm)
2580: {
2581:   PetscFunctionBegin;
2583:   PetscAssertPointer(tdm, 2);
2584:   *tdm = dm->transformDM;
2585:   PetscFunctionReturn(PETSC_SUCCESS);
2586: }

2588: PetscErrorCode DMGetBasisTransformVec_Internal(DM dm, Vec *tv)
2589: {
2590:   PetscFunctionBegin;
2592:   PetscAssertPointer(tv, 2);
2593:   *tv = dm->transform;
2594:   PetscFunctionReturn(PETSC_SUCCESS);
2595: }

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

2600:   Input Parameter:
2601: . dm - The `DM`

2603:   Output Parameter:
2604: . flg - `PETSC_TRUE` if a basis transformation should be done

2606:   Level: developer

2608: .seealso: [](ch_dmbase), `DM`, `DMPlexGlobalToLocalBasis()`, `DMPlexLocalToGlobalBasis()`, `DMPlexCreateBasisRotation()`
2609: @*/
2610: PetscErrorCode DMHasBasisTransform(DM dm, PetscBool *flg)
2611: {
2612:   Vec tv;

2614:   PetscFunctionBegin;
2616:   PetscAssertPointer(flg, 2);
2617:   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
2618:   *flg = tv ? PETSC_TRUE : PETSC_FALSE;
2619:   PetscFunctionReturn(PETSC_SUCCESS);
2620: }

2622: PetscErrorCode DMConstructBasisTransform_Internal(DM dm)
2623: {
2624:   PetscSection s, ts;
2625:   PetscScalar *ta;
2626:   PetscInt     cdim, pStart, pEnd, p, Nf, f, Nc, dof;

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

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

2670: PetscErrorCode DMCopyTransform(DM dm, DM newdm)
2671: {
2672:   PetscFunctionBegin;
2675:   newdm->transformCtx       = dm->transformCtx;
2676:   newdm->transformSetUp     = dm->transformSetUp;
2677:   newdm->transformDestroy   = NULL;
2678:   newdm->transformGetMatrix = dm->transformGetMatrix;
2679:   if (newdm->transformSetUp) PetscCall(DMConstructBasisTransform_Internal(newdm));
2680:   PetscFunctionReturn(PETSC_SUCCESS);
2681: }

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

2686:   Logically Collective

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

2694:   Calling sequence of `beginhook`:
2695: + dm   - global `DM`
2696: . g    - global vector
2697: . mode - mode
2698: . l    - local vector
2699: - ctx  - optional user-defined function context

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

2708:   Level: advanced

2710:   Note:
2711:   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.

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

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

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

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

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

2766: /*@
2767:   DMGlobalToLocal - update local vectors from global vector

2769:   Neighbor-wise Collective

2771:   Input Parameters:
2772: + dm   - the `DM` object
2773: . g    - the global vector
2774: . mode - `INSERT_VALUES` or `ADD_VALUES`
2775: - l    - the local vector

2777:   Level: beginner

2779:   Notes:
2780:   The communication involved in this update can be overlapped with computation by instead using
2781:   `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()`.

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

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

2797: /*@
2798:   DMGlobalToLocalBegin - Begins updating local vectors from global vector

2800:   Neighbor-wise Collective

2802:   Input Parameters:
2803: + dm   - the `DM` object
2804: . g    - the global vector
2805: . mode - `INSERT_VALUES` or `ADD_VALUES`
2806: - l    - the local vector

2808:   Level: intermediate

2810:   Notes:
2811:   The operation is completed with `DMGlobalToLocalEnd()`

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

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

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

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

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

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

2849: /*@
2850:   DMGlobalToLocalEnd - Ends updating local vectors from global vector

2852:   Neighbor-wise Collective

2854:   Input Parameters:
2855: + dm   - the `DM` object
2856: . g    - the global vector
2857: . mode - `INSERT_VALUES` or `ADD_VALUES`
2858: - l    - the local vector

2860:   Level: intermediate

2862:   Note:
2863:   See `DMGlobalToLocalBegin()` for details.

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

2876:   PetscFunctionBegin;
2878:   PetscCall(DMGetSectionSF(dm, &sf));
2879:   PetscCall(DMHasBasisTransform(dm, &transform));
2880:   if (sf) {
2881:     PetscCheck(mode != ADD_VALUES, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %d", (int)mode);

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

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

2902:   Logically Collective

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

2910:   Calling sequence of `beginhook`:
2911: + global - global `DM`
2912: . l      - local vector
2913: . mode   - mode
2914: . g      - global vector
2915: - ctx    - optional user-defined function context

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

2924:   Level: advanced

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

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

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

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

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

2985:   Neighbor-wise Collective

2987:   Input Parameters:
2988: + dm   - the `DM` object
2989: . l    - the local vector
2990: . 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.
2991: - g    - the global vector

2993:   Level: beginner

2995:   Notes:
2996:   The communication involved in this update can be overlapped with computation by using
2997:   `DMLocalToGlobalBegin()` and `DMLocalToGlobalEnd()`.

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

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

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

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

3015: /*@
3016:   DMLocalToGlobalBegin - begins updating global vectors from local vectors

3018:   Neighbor-wise Collective

3020:   Input Parameters:
3021: + dm   - the `DM` object
3022: . l    - the local vector
3023: . 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.
3024: - g    - the global vector

3026:   Level: intermediate

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

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

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

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

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

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

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

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

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

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

3151: /*@
3152:   DMLocalToGlobalEnd - updates global vectors from local vectors

3154:   Neighbor-wise Collective

3156:   Input Parameters:
3157: + dm   - the `DM` object
3158: . l    - the local vector
3159: . mode - `INSERT_VALUES` or `ADD_VALUES`
3160: - g    - the global vector

3162:   Level: intermediate

3164:   Note:
3165:   See `DMLocalToGlobalBegin()` for full details

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

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

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

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

3228:   Neighbor-wise Collective

3230:   Input Parameters:
3231: + dm   - the `DM` object
3232: . g    - the original local vector
3233: - mode - one of `INSERT_VALUES` or `ADD_VALUES`

3235:   Output Parameter:
3236: . l - the local vector with correct ghost values

3238:   Level: intermediate

3240:   Note:
3241:   Must be followed by `DMLocalToLocalEnd()`.

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

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

3259:   Neighbor-wise Collective

3261:   Input Parameters:
3262: + dm   - the `DM` object
3263: . g    - the original local vector
3264: - mode - one of `INSERT_VALUES` or `ADD_VALUES`

3266:   Output Parameter:
3267: . l - the local vector with correct ghost values

3269:   Level: intermediate

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

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

3286:   Collective

3288:   Input Parameters:
3289: + dm   - the `DM` object
3290: - comm - the communicator to contain the new `DM` object (or `MPI_COMM_NULL`)

3292:   Output Parameter:
3293: . dmc - the coarsened `DM`

3295:   Level: developer

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

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

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

3329:   Logically Collective; No Fortran Support

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

3337:   Calling sequence of `coarsenhook`:
3338: + fine   - fine level `DM`
3339: . coarse - coarse level `DM` to restrict problem to
3340: - ctx    - optional user-defined function context

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

3350:   Level: advanced

3352:   Notes:
3353:   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`.

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

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

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

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

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

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

3385:   Logically Collective; No Fortran Support

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

3393:   Level: advanced

3395:   Notes:
3396:   This function does nothing if the `coarsenhook` is not in the list.

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

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

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

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

3422:   Collective if any hooks are

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

3431:   Level: developer

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

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

3442:   PetscFunctionBegin;
3443:   for (link = fine->coarsenhook; link; link = link->next) {
3444:     if (link->restricthook) PetscCall((*link->restricthook)(fine, restrct, rscale, inject, coarse, link->ctx));
3445:   }
3446:   PetscFunctionReturn(PETSC_SUCCESS);
3447: }

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

3452:   Logically Collective; No Fortran Support

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

3460:   Calling sequence of `ddhook`:
3461: + global - global `DM`
3462: . block  - subdomain `DM`
3463: - ctx    - optional user-defined function context

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

3472:   Level: advanced

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

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

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

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

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

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

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

3508:   Logically Collective; No Fortran Support

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

3516:   Level: advanced

3518:   Note:
3519:   See `DMSubDomainHookAdd()` for the calling sequences of `ddhook` and `restricthook`

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

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

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

3544:   Collective if any hooks are

3546:   Input Parameters:
3547: + global   - The global `DM` to use as a base
3548: . oscatter - The scatter from domain global vector filling subdomain global vector with overlap
3549: . gscatter - The scatter from domain global vector filling subdomain local vector with ghosts
3550: - subdm    - The subdomain `DM` to update

3552:   Level: developer

3554: .seealso: [](ch_dmbase), `DM`, `DMCoarsenHookAdd()`, `MatRestrict()`, `DMCreateDomainDecomposition()`
3555: @*/
3556: PetscErrorCode DMSubDomainRestrict(DM global, VecScatter oscatter, VecScatter gscatter, DM subdm)
3557: {
3558:   DMSubDomainHookLink link;

3560:   PetscFunctionBegin;
3561:   for (link = global->subdomainhook; link; link = link->next) {
3562:     if (link->restricthook) PetscCall((*link->restricthook)(global, oscatter, gscatter, subdm, link->ctx));
3563:   }
3564:   PetscFunctionReturn(PETSC_SUCCESS);
3565: }

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

3570:   Not Collective

3572:   Input Parameter:
3573: . dm - the `DM` object

3575:   Output Parameter:
3576: . level - number of coarsenings

3578:   Level: developer

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

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

3594:   Collective

3596:   Input Parameters:
3597: + dm    - the `DM` object
3598: - level - number of coarsenings

3600:   Level: developer

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

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

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

3618:   Collective

3620:   Input Parameters:
3621: + dm      - the `DM` object
3622: - nlevels - the number of levels of refinement

3624:   Output Parameter:
3625: . dmf - the refined `DM` hierarchy

3627:   Level: developer

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

3641:     PetscCall(DMRefine(dm, PetscObjectComm((PetscObject)dm), &dmf[0]));
3642:     for (i = 1; i < nlevels; i++) PetscCall(DMRefine(dmf[i - 1], PetscObjectComm((PetscObject)dm), &dmf[i]));
3643:   } else PetscUseTypeMethod(dm, refinehierarchy, nlevels, dmf);
3644:   PetscFunctionReturn(PETSC_SUCCESS);
3645: }

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

3650:   Collective

3652:   Input Parameters:
3653: + dm      - the `DM` object
3654: - nlevels - the number of levels of coarsening

3656:   Output Parameter:
3657: . dmc - the coarsened `DM` hierarchy

3659:   Level: developer

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

3673:     PetscCall(DMCoarsen(dm, PetscObjectComm((PetscObject)dm), &dmc[0]));
3674:     for (i = 1; i < nlevels; i++) PetscCall(DMCoarsen(dmc[i - 1], PetscObjectComm((PetscObject)dm), &dmc[i]));
3675:   } else PetscUseTypeMethod(dm, coarsenhierarchy, nlevels, dmc);
3676:   PetscFunctionReturn(PETSC_SUCCESS);
3677: }

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

3682:   Logically Collective if the function is collective

3684:   Input Parameters:
3685: + dm      - the `DM` object
3686: - destroy - the destroy function

3688:   Level: intermediate

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

3700: /*@
3701:   DMSetApplicationContext - Set a user context into a `DM` object

3703:   Not Collective

3705:   Input Parameters:
3706: + dm  - the `DM` object
3707: - ctx - the user context

3709:   Level: intermediate

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

3714: .seealso: [](ch_dmbase), `DM`, `DMGetApplicationContext()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`
3715: @*/
3716: PetscErrorCode DMSetApplicationContext(DM dm, void *ctx)
3717: {
3718:   PetscFunctionBegin;
3720:   dm->ctx = ctx;
3721:   PetscFunctionReturn(PETSC_SUCCESS);
3722: }

3724: /*@
3725:   DMGetApplicationContext - Gets a user context from a `DM` object

3727:   Not Collective

3729:   Input Parameter:
3730: . dm - the `DM` object

3732:   Output Parameter:
3733: . ctx - the user context

3735:   Level: intermediate

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

3740: .seealso: [](ch_dmbase), `DM`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`
3741: @*/
3742: PetscErrorCode DMGetApplicationContext(DM dm, void *ctx)
3743: {
3744:   PetscFunctionBegin;
3746:   *(void **)ctx = dm->ctx;
3747:   PetscFunctionReturn(PETSC_SUCCESS);
3748: }

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

3753:   Logically Collective

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

3759:   Level: intermediate

3761: .seealso: [](ch_dmbase), `DM`, `DMComputeVariableBounds()`, `DMHasVariableBounds()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMGetApplicationContext()`,
3762:          `DMSetJacobian()`
3763: @*/
3764: PetscErrorCode DMSetVariableBounds(DM dm, PetscErrorCode (*f)(DM, Vec, Vec))
3765: {
3766:   PetscFunctionBegin;
3768:   dm->ops->computevariablebounds = f;
3769:   PetscFunctionReturn(PETSC_SUCCESS);
3770: }

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

3775:   Not Collective

3777:   Input Parameter:
3778: . dm - the `DM` object to destroy

3780:   Output Parameter:
3781: . flg - `PETSC_TRUE` if the variable bounds function exists

3783:   Level: developer

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

3796: /*@C
3797:   DMComputeVariableBounds - compute variable bounds used by `SNESVI`.

3799:   Logically Collective

3801:   Input Parameter:
3802: . dm - the `DM` object

3804:   Output Parameters:
3805: + xl - lower bound
3806: - xu - upper bound

3808:   Level: advanced

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

3813: .seealso: [](ch_dmbase), `DM`, `DMHasVariableBounds()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateInterpolation()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMGetApplicationContext()`
3814: @*/
3815: PetscErrorCode DMComputeVariableBounds(DM dm, Vec xl, Vec xu)
3816: {
3817:   PetscFunctionBegin;
3821:   PetscUseTypeMethod(dm, computevariablebounds, xl, xu);
3822:   PetscFunctionReturn(PETSC_SUCCESS);
3823: }

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

3828:   Not Collective

3830:   Input Parameter:
3831: . dm - the DM object

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

3836:   Level: developer

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

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

3852:   Not Collective

3854:   Input Parameter:
3855: . dm - the `DM` object

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

3860:   Level: developer

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

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

3876:   Not Collective

3878:   Input Parameter:
3879: . dm - the `DM` object

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

3884:   Level: developer

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

3898: PetscFunctionList DMList              = NULL;
3899: PetscBool         DMRegisterAllCalled = PETSC_FALSE;

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

3904:   Collective

3906:   Input Parameters:
3907: + dm     - The `DM` object
3908: - method - The name of the `DMType`, for example `DMDA`, `DMPLEX`

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

3913:   Level: intermediate

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

3918: .seealso: [](ch_dmbase), `DM`, `DMType`, `DMDA`, `DMPLEX`, `DMGetType()`, `DMCreate()`, `DMDACreate2d()`
3919: @*/
3920: PetscErrorCode DMSetType(DM dm, DMType method)
3921: {
3922:   PetscErrorCode (*r)(DM);
3923:   PetscBool match;

3925:   PetscFunctionBegin;
3927:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, method, &match));
3928:   if (match) PetscFunctionReturn(PETSC_SUCCESS);

3930:   PetscCall(DMRegisterAll());
3931:   PetscCall(PetscFunctionListFind(DMList, method, &r));
3932:   PetscCheck(r, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);

3934:   PetscTryTypeMethod(dm, destroy);
3935:   PetscCall(PetscMemzero(dm->ops, sizeof(*dm->ops)));
3936:   PetscCall(PetscObjectChangeTypeName((PetscObject)dm, method));
3937:   PetscCall((*r)(dm));
3938:   PetscFunctionReturn(PETSC_SUCCESS);
3939: }

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

3944:   Not Collective

3946:   Input Parameter:
3947: . dm - The `DM`

3949:   Output Parameter:
3950: . type - The `DMType` name

3952:   Level: intermediate

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

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

3969:   Collective

3971:   Input Parameters:
3972: + dm      - the `DM`
3973: - newtype - new `DM` type (use "same" for the same type)

3975:   Output Parameter:
3976: . M - pointer to new `DM`

3978:   Level: intermediate

3980:   Note:
3981:   Cannot be used to convert a sequential `DM` to a parallel or a parallel to sequential,
3982:   the MPI communicator of the generated `DM` is always the same as the communicator
3983:   of the input `DM`.

3985: .seealso: [](ch_dmbase), `DM`, `DMSetType()`, `DMCreate()`, `DMClone()`
3986: @*/
3987: PetscErrorCode DMConvert(DM dm, DMType newtype, DM *M)
3988: {
3989:   DM        B;
3990:   char      convname[256];
3991:   PetscBool sametype /*, issame */;

3993:   PetscFunctionBegin;
3996:   PetscAssertPointer(M, 3);
3997:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, newtype, &sametype));
3998:   /* PetscCall(PetscStrcmp(newtype, "same", &issame)); */
3999:   if (sametype) {
4000:     *M = dm;
4001:     PetscCall(PetscObjectReference((PetscObject)dm));
4002:     PetscFunctionReturn(PETSC_SUCCESS);
4003:   } else {
4004:     PetscErrorCode (*conv)(DM, DMType, DM *) = NULL;

4006:     /*
4007:        Order of precedence:
4008:        1) See if a specialized converter is known to the current DM.
4009:        2) See if a specialized converter is known to the desired DM class.
4010:        3) See if a good general converter is registered for the desired class
4011:        4) See if a good general converter is known for the current matrix.
4012:        5) Use a really basic converter.
4013:     */

4015:     /* 1) See if a specialized converter is known to the current DM and the desired class */
4016:     PetscCall(PetscStrncpy(convname, "DMConvert_", sizeof(convname)));
4017:     PetscCall(PetscStrlcat(convname, ((PetscObject)dm)->type_name, sizeof(convname)));
4018:     PetscCall(PetscStrlcat(convname, "_", sizeof(convname)));
4019:     PetscCall(PetscStrlcat(convname, newtype, sizeof(convname)));
4020:     PetscCall(PetscStrlcat(convname, "_C", sizeof(convname)));
4021:     PetscCall(PetscObjectQueryFunction((PetscObject)dm, convname, &conv));
4022:     if (conv) goto foundconv;

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

4038: #if 0
4039:     /* 3) See if a good general converter is registered for the desired class */
4040:     conv = B->ops->convertfrom;
4041:     PetscCall(DMDestroy(&B));
4042:     if (conv) goto foundconv;

4044:     /* 4) See if a good general converter is known for the current matrix */
4045:     if (dm->ops->convert) {
4046:       conv = dm->ops->convert;
4047:     }
4048:     if (conv) goto foundconv;
4049: #endif

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

4054:   foundconv:
4055:     PetscCall(PetscLogEventBegin(DM_Convert, dm, 0, 0, 0));
4056:     PetscCall((*conv)(dm, newtype, M));
4057:     /* Things that are independent of DM type: We should consult DMClone() here */
4058:     {
4059:       const PetscReal *maxCell, *Lstart, *L;

4061:       PetscCall(DMGetPeriodicity(dm, &maxCell, &Lstart, &L));
4062:       PetscCall(DMSetPeriodicity(*M, maxCell, Lstart, L));
4063:       (*M)->prealloc_only = dm->prealloc_only;
4064:       PetscCall(PetscFree((*M)->vectype));
4065:       PetscCall(PetscStrallocpy(dm->vectype, (char **)&(*M)->vectype));
4066:       PetscCall(PetscFree((*M)->mattype));
4067:       PetscCall(PetscStrallocpy(dm->mattype, (char **)&(*M)->mattype));
4068:     }
4069:     PetscCall(PetscLogEventEnd(DM_Convert, dm, 0, 0, 0));
4070:   }
4071:   PetscCall(PetscObjectStateIncrease((PetscObject)*M));
4072:   PetscFunctionReturn(PETSC_SUCCESS);
4073: }

4075: /*--------------------------------------------------------------------------------------------------------------------*/

4077: /*@C
4078:   DMRegister -  Adds a new `DM` type implementation

4080:   Not Collective

4082:   Input Parameters:
4083: + sname    - The name of a new user-defined creation routine
4084: - function - The creation routine itself

4086:   Level: advanced

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

4091:   Example Usage:
4092: .vb
4093:     DMRegister("my_da", MyDMCreate);
4094: .ve

4096:   Then, your `DM` type can be chosen with the procedural interface via
4097: .vb
4098:     DMCreate(MPI_Comm, DM *);
4099:     DMSetType(DM,"my_da");
4100: .ve
4101:   or at runtime via the option
4102: .vb
4103:     -da_type my_da
4104: .ve

4106: .seealso: [](ch_dmbase), `DM`, `DMType`, `DMSetType()`, `DMRegisterAll()`, `DMRegisterDestroy()`
4107: @*/
4108: PetscErrorCode DMRegister(const char sname[], PetscErrorCode (*function)(DM))
4109: {
4110:   PetscFunctionBegin;
4111:   PetscCall(DMInitializePackage());
4112:   PetscCall(PetscFunctionListAdd(&DMList, sname, function));
4113:   PetscFunctionReturn(PETSC_SUCCESS);
4114: }

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

4119:   Collective

4121:   Input Parameters:
4122: + newdm  - the newly loaded `DM`, this needs to have been created with `DMCreate()` or
4123:            some related function before a call to `DMLoad()`.
4124: - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()` or
4125:            `PETSCVIEWERHDF5` file viewer, obtained from `PetscViewerHDF5Open()`

4127:   Level: intermediate

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

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

4136: .seealso: [](ch_dmbase), `DM`, `PetscViewerBinaryOpen()`, `DMView()`, `MatLoad()`, `VecLoad()`
4137: @*/
4138: PetscErrorCode DMLoad(DM newdm, PetscViewer viewer)
4139: {
4140:   PetscBool isbinary, ishdf5;

4142:   PetscFunctionBegin;
4145:   PetscCall(PetscViewerCheckReadable(viewer));
4146:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
4147:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
4148:   PetscCall(PetscLogEventBegin(DM_Load, viewer, 0, 0, 0));
4149:   if (isbinary) {
4150:     PetscInt classid;
4151:     char     type[256];

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

4165: /******************************** FEM Support **********************************/

4167: PetscErrorCode DMPrintCellIndices(PetscInt c, const char name[], PetscInt len, const PetscInt x[])
4168: {
4169:   PetscInt f;

4171:   PetscFunctionBegin;
4172:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " Element %s\n", c, name));
4173:   for (f = 0; f < len; ++f) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  | %" PetscInt_FMT " |\n", x[f]));
4174:   PetscFunctionReturn(PETSC_SUCCESS);
4175: }

4177: PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[])
4178: {
4179:   PetscInt f;

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

4187: PetscErrorCode DMPrintCellVectorReal(PetscInt c, const char name[], PetscInt len, const PetscReal x[])
4188: {
4189:   PetscInt f;

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

4197: PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[])
4198: {
4199:   PetscInt f, g;

4201:   PetscFunctionBegin;
4202:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " Element %s\n", c, name));
4203:   for (f = 0; f < rows; ++f) {
4204:     PetscCall(PetscPrintf(PETSC_COMM_SELF, "  |"));
4205:     for (g = 0; g < cols; ++g) PetscCall(PetscPrintf(PETSC_COMM_SELF, " % 9.5g", (double)PetscRealPart(A[f * cols + g])));
4206:     PetscCall(PetscPrintf(PETSC_COMM_SELF, " |\n"));
4207:   }
4208:   PetscFunctionReturn(PETSC_SUCCESS);
4209: }

4211: PetscErrorCode DMPrintLocalVec(DM dm, const char name[], PetscReal tol, Vec X)
4212: {
4213:   PetscInt           localSize, bs;
4214:   PetscMPIInt        size;
4215:   Vec                x, xglob;
4216:   const PetscScalar *xarray;

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

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

4244:   Input Parameter:
4245: . dm - The `DM`

4247:   Output Parameter:
4248: . section - The `PetscSection`

4250:   Options Database Key:
4251: . -dm_petscsection_view - View the `PetscSection` created by the `DM`

4253:   Level: advanced

4255:   Notes:
4256:   Use `DMGetLocalSection()` in new code.

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

4260: .seealso: [](ch_dmbase), `DM`, `DMGetLocalSection()`, `DMSetLocalSection()`, `DMGetGlobalSection()`
4261: @*/
4262: PetscErrorCode DMGetSection(DM dm, PetscSection *section)
4263: {
4264:   PetscFunctionBegin;
4265:   PetscCall(DMGetLocalSection(dm, section));
4266:   PetscFunctionReturn(PETSC_SUCCESS);
4267: }

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

4272:   Input Parameter:
4273: . dm - The `DM`

4275:   Output Parameter:
4276: . section - The `PetscSection`

4278:   Options Database Key:
4279: . -dm_petscsection_view - View the section created by the `DM`

4281:   Level: intermediate

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

4286: .seealso: [](ch_dmbase), `DM`, `DMSetLocalSection()`, `DMGetGlobalSection()`
4287: @*/
4288: PetscErrorCode DMGetLocalSection(DM dm, PetscSection *section)
4289: {
4290:   PetscFunctionBegin;
4292:   PetscAssertPointer(section, 2);
4293:   if (!dm->localSection && dm->ops->createlocalsection) {
4294:     PetscInt d;

4296:     if (dm->setfromoptionscalled) {
4297:       PetscObject       obj = (PetscObject)dm;
4298:       PetscViewer       viewer;
4299:       PetscViewerFormat format;
4300:       PetscBool         flg;

4302:       PetscCall(PetscOptionsGetViewer(PetscObjectComm(obj), obj->options, obj->prefix, "-dm_petscds_view", &viewer, &format, &flg));
4303:       if (flg) PetscCall(PetscViewerPushFormat(viewer, format));
4304:       for (d = 0; d < dm->Nds; ++d) {
4305:         PetscCall(PetscDSSetFromOptions(dm->probs[d].ds));
4306:         if (flg) PetscCall(PetscDSView(dm->probs[d].ds, viewer));
4307:       }
4308:       if (flg) {
4309:         PetscCall(PetscViewerFlush(viewer));
4310:         PetscCall(PetscViewerPopFormat(viewer));
4311:         PetscCall(PetscOptionsRestoreViewer(&viewer));
4312:       }
4313:     }
4314:     PetscUseTypeMethod(dm, createlocalsection);
4315:     if (dm->localSection) PetscCall(PetscObjectViewFromOptions((PetscObject)dm->localSection, NULL, "-dm_petscsection_view"));
4316:   }
4317:   *section = dm->localSection;
4318:   PetscFunctionReturn(PETSC_SUCCESS);
4319: }

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

4324:   Input Parameters:
4325: + dm      - The `DM`
4326: - section - The `PetscSection`

4328:   Level: advanced

4330:   Notes:
4331:   Use `DMSetLocalSection()` in new code.

4333:   Any existing `PetscSection` will be destroyed

4335: .seealso: [](ch_dmbase), `DM`, `DMSetLocalSection()`, `DMGetLocalSection()`, `DMSetGlobalSection()`
4336: @*/
4337: PetscErrorCode DMSetSection(DM dm, PetscSection section)
4338: {
4339:   PetscFunctionBegin;
4340:   PetscCall(DMSetLocalSection(dm, section));
4341:   PetscFunctionReturn(PETSC_SUCCESS);
4342: }

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

4347:   Input Parameters:
4348: + dm      - The `DM`
4349: - section - The `PetscSection`

4351:   Level: intermediate

4353:   Note:
4354:   Any existing Section will be destroyed

4356: .seealso: [](ch_dmbase), `DM`, `PetscSection`, `DMGetLocalSection()`, `DMSetGlobalSection()`
4357: @*/
4358: PetscErrorCode DMSetLocalSection(DM dm, PetscSection section)
4359: {
4360:   PetscInt numFields = 0;
4361:   PetscInt f;

4363:   PetscFunctionBegin;
4366:   PetscCall(PetscObjectReference((PetscObject)section));
4367:   PetscCall(PetscSectionDestroy(&dm->localSection));
4368:   dm->localSection = section;
4369:   if (section) PetscCall(PetscSectionGetNumFields(dm->localSection, &numFields));
4370:   if (numFields) {
4371:     PetscCall(DMSetNumFields(dm, numFields));
4372:     for (f = 0; f < numFields; ++f) {
4373:       PetscObject disc;
4374:       const char *name;

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

4387:   /* Clear scratch vectors */
4388:   PetscCall(DMClearGlobalVectors(dm));
4389:   PetscCall(DMClearLocalVectors(dm));
4390:   PetscCall(DMClearNamedGlobalVectors(dm));
4391:   PetscCall(DMClearNamedLocalVectors(dm));
4392:   PetscFunctionReturn(PETSC_SUCCESS);
4393: }

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

4398:   Input Parameter:
4399: . dm - The `DM`

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

4405:   Level: developer

4407: .seealso: [](ch_dmbase), `DM`, `PetscSection`, `DMGetLocalSection()`, `DMGetGlobalSection()`
4408: @*/
4409: PetscErrorCode DMCreateSectionPermutation(DM dm, IS *perm, PetscBT *blockStarts)
4410: {
4411:   PetscFunctionBegin;
4412:   *perm        = NULL;
4413:   *blockStarts = NULL;
4414:   PetscTryTypeMethod(dm, createsectionpermutation, perm, blockStarts);
4415:   PetscFunctionReturn(PETSC_SUCCESS);
4416: }

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

4421:   not Collective

4423:   Input Parameter:
4424: . dm - The `DM`

4426:   Output Parameters:
4427: + 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.
4428: . 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.
4429: - bias    - Vector containing bias to be added to constrained dofs

4431:   Level: advanced

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

4436: .seealso: [](ch_dmbase), `DM`, `DMSetDefaultConstraints()`
4437: @*/
4438: PetscErrorCode DMGetDefaultConstraints(DM dm, PetscSection *section, Mat *mat, Vec *bias)
4439: {
4440:   PetscFunctionBegin;
4442:   if (!dm->defaultConstraint.section && !dm->defaultConstraint.mat && dm->ops->createdefaultconstraints) PetscUseTypeMethod(dm, createdefaultconstraints);
4443:   if (section) *section = dm->defaultConstraint.section;
4444:   if (mat) *mat = dm->defaultConstraint.mat;
4445:   if (bias) *bias = dm->defaultConstraint.bias;
4446:   PetscFunctionReturn(PETSC_SUCCESS);
4447: }

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

4452:   Collective

4454:   Input Parameters:
4455: + dm      - The `DM`
4456: . 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).
4457: . 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).
4458: - 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).

4460:   Level: advanced

4462:   Notes:
4463:   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()`.

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

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

4469: .seealso: [](ch_dmbase), `DM`, `DMGetDefaultConstraints()`
4470: @*/
4471: PetscErrorCode DMSetDefaultConstraints(DM dm, PetscSection section, Mat mat, Vec bias)
4472: {
4473:   PetscMPIInt result;

4475:   PetscFunctionBegin;
4477:   if (section) {
4479:     PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, PetscObjectComm((PetscObject)section), &result));
4480:     PetscCheck(result == MPI_CONGRUENT || result == MPI_IDENT, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMECOMM, "constraint section must have local communicator");
4481:   }
4482:   if (mat) {
4484:     PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, PetscObjectComm((PetscObject)mat), &result));
4485:     PetscCheck(result == MPI_CONGRUENT || result == MPI_IDENT, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMECOMM, "constraint matrix must have local communicator");
4486:   }
4487:   if (bias) {
4489:     PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, PetscObjectComm((PetscObject)bias), &result));
4490:     PetscCheck(result == MPI_CONGRUENT || result == MPI_IDENT, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMECOMM, "constraint bias must have local communicator");
4491:   }
4492:   PetscCall(PetscObjectReference((PetscObject)section));
4493:   PetscCall(PetscSectionDestroy(&dm->defaultConstraint.section));
4494:   dm->defaultConstraint.section = section;
4495:   PetscCall(PetscObjectReference((PetscObject)mat));
4496:   PetscCall(MatDestroy(&dm->defaultConstraint.mat));
4497:   dm->defaultConstraint.mat = mat;
4498:   PetscCall(PetscObjectReference((PetscObject)bias));
4499:   PetscCall(VecDestroy(&dm->defaultConstraint.bias));
4500:   dm->defaultConstraint.bias = bias;
4501:   PetscFunctionReturn(PETSC_SUCCESS);
4502: }

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

4508:   Input Parameters:
4509: + dm - The `DM`
4510: . localSection - `PetscSection` describing the local data layout
4511: - globalSection - `PetscSection` describing the global data layout

4513:   Level: intermediate

4515: .seealso: [](ch_dmbase), `DM`, `DMGetSectionSF()`, `DMSetSectionSF()`
4516: */
4517: static PetscErrorCode DMDefaultSectionCheckConsistency_Internal(DM dm, PetscSection localSection, PetscSection globalSection)
4518: {
4519:   MPI_Comm        comm;
4520:   PetscLayout     layout;
4521:   const PetscInt *ranges;
4522:   PetscInt        pStart, pEnd, p, nroots;
4523:   PetscMPIInt     size, rank;
4524:   PetscBool       valid = PETSC_TRUE, gvalid;

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

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

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

4582: static PetscErrorCode DMGetIsoperiodicPointSF_Internal(DM dm, PetscSF *sf)
4583: {
4584:   PetscErrorCode (*f)(DM, PetscSF *);

4586:   PetscFunctionBegin;
4588:   PetscAssertPointer(sf, 2);
4589:   PetscCall(PetscObjectQueryFunction((PetscObject)dm, "DMGetIsoperiodicPointSF_C", &f));
4590:   if (f) PetscCall(f(dm, sf));
4591:   else *sf = dm->sf;
4592:   PetscFunctionReturn(PETSC_SUCCESS);
4593: }

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

4598:   Collective

4600:   Input Parameter:
4601: . dm - The `DM`

4603:   Output Parameter:
4604: . section - The `PetscSection`

4606:   Level: intermediate

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

4611: .seealso: [](ch_dmbase), `DM`, `DMSetLocalSection()`, `DMGetLocalSection()`
4612: @*/
4613: PetscErrorCode DMGetGlobalSection(DM dm, PetscSection *section)
4614: {
4615:   PetscFunctionBegin;
4617:   PetscAssertPointer(section, 2);
4618:   if (!dm->globalSection) {
4619:     PetscSection s;
4620:     PetscSF      sf;

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

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

4638:   Input Parameters:
4639: + dm      - The `DM`
4640: - section - The PetscSection, or `NULL`

4642:   Level: intermediate

4644:   Note:
4645:   Any existing `PetscSection` will be destroyed

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

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

4672:   Input Parameter:
4673: . dm - The `DM`

4675:   Output Parameter:
4676: . sf - The `PetscSF`

4678:   Level: intermediate

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

4683: .seealso: [](ch_dmbase), `DM`, `DMSetSectionSF()`, `DMCreateSectionSF()`
4684: @*/
4685: PetscErrorCode DMGetSectionSF(DM dm, PetscSF *sf)
4686: {
4687:   PetscInt nroots;

4689:   PetscFunctionBegin;
4691:   PetscAssertPointer(sf, 2);
4692:   if (!dm->sectionSF) PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), &dm->sectionSF));
4693:   PetscCall(PetscSFGetGraph(dm->sectionSF, &nroots, NULL, NULL, NULL));
4694:   if (nroots < 0) {
4695:     PetscSection section, gSection;

4697:     PetscCall(DMGetLocalSection(dm, &section));
4698:     if (section) {
4699:       PetscCall(DMGetGlobalSection(dm, &gSection));
4700:       PetscCall(DMCreateSectionSF(dm, section, gSection));
4701:     } else {
4702:       *sf = NULL;
4703:       PetscFunctionReturn(PETSC_SUCCESS);
4704:     }
4705:   }
4706:   *sf = dm->sectionSF;
4707:   PetscFunctionReturn(PETSC_SUCCESS);
4708: }

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

4713:   Input Parameters:
4714: + dm - The `DM`
4715: - sf - The `PetscSF`

4717:   Level: intermediate

4719:   Note:
4720:   Any previous `PetscSF` is destroyed

4722: .seealso: [](ch_dmbase), `DM`, `DMGetSectionSF()`, `DMCreateSectionSF()`
4723: @*/
4724: PetscErrorCode DMSetSectionSF(DM dm, PetscSF sf)
4725: {
4726:   PetscFunctionBegin;
4729:   PetscCall(PetscObjectReference((PetscObject)sf));
4730:   PetscCall(PetscSFDestroy(&dm->sectionSF));
4731:   dm->sectionSF = sf;
4732:   PetscFunctionReturn(PETSC_SUCCESS);
4733: }

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

4739:   Input Parameters:
4740: + dm            - The `DM`
4741: . localSection  - `PetscSection` describing the local data layout
4742: - globalSection - `PetscSection` describing the global data layout

4744:   Level: developer

4746:   Note:
4747:   One usually uses `DMGetSectionSF()` to obtain the `PetscSF`

4749:   Developer Note:
4750:   Since this routine has for arguments the two sections from the `DM` and puts the resulting `PetscSF`
4751:   directly into the `DM`, perhaps this function should not take the local and global sections as
4752:   input and should just obtain them from the `DM`?

4754: .seealso: [](ch_dmbase), `DM`, `DMGetSectionSF()`, `DMSetSectionSF()`, `DMGetLocalSection()`, `DMGetGlobalSection()`
4755: @*/
4756: PetscErrorCode DMCreateSectionSF(DM dm, PetscSection localSection, PetscSection globalSection)
4757: {
4758:   PetscFunctionBegin;
4760:   PetscCall(PetscSFSetGraphSection(dm->sectionSF, localSection, globalSection));
4761:   PetscFunctionReturn(PETSC_SUCCESS);
4762: }

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

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

4769:   Input Parameter:
4770: . dm - The `DM`

4772:   Output Parameter:
4773: . sf - The `PetscSF`

4775:   Level: intermediate

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

4780: .seealso: [](ch_dmbase), `DM`, `DMSetPointSF()`, `DMGetSectionSF()`, `DMSetSectionSF()`, `DMCreateSectionSF()`
4781: @*/
4782: PetscErrorCode DMGetPointSF(DM dm, PetscSF *sf)
4783: {
4784:   PetscFunctionBegin;
4786:   PetscAssertPointer(sf, 2);
4787:   *sf = dm->sf;
4788:   PetscFunctionReturn(PETSC_SUCCESS);
4789: }

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

4794:   Collective

4796:   Input Parameters:
4797: + dm - The `DM`
4798: - sf - The `PetscSF`

4800:   Level: intermediate

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

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

4818:   Input Parameter:
4819: . dm - The `DM`

4821:   Output Parameter:
4822: . sf - The `PetscSF`

4824:   Level: intermediate

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

4829: .seealso: [](ch_dmbase), `DM`, `DMSetNaturalSF()`, `DMSetUseNatural()`, `DMGetUseNatural()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexDistribute()`
4830: @*/
4831: PetscErrorCode DMGetNaturalSF(DM dm, PetscSF *sf)
4832: {
4833:   PetscFunctionBegin;
4835:   PetscAssertPointer(sf, 2);
4836:   *sf = dm->sfNatural;
4837:   PetscFunctionReturn(PETSC_SUCCESS);
4838: }

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

4843:   Input Parameters:
4844: + dm - The DM
4845: - sf - The PetscSF

4847:   Level: intermediate

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

4862: static PetscErrorCode DMSetDefaultAdjacency_Private(DM dm, PetscInt f, PetscObject disc)
4863: {
4864:   PetscClassId id;

4866:   PetscFunctionBegin;
4867:   PetscCall(PetscObjectGetClassId(disc, &id));
4868:   if (id == PETSCFE_CLASSID) {
4869:     PetscCall(DMSetAdjacency(dm, f, PETSC_FALSE, PETSC_TRUE));
4870:   } else if (id == PETSCFV_CLASSID) {
4871:     PetscCall(DMSetAdjacency(dm, f, PETSC_TRUE, PETSC_FALSE));
4872:   } else {
4873:     PetscCall(DMSetAdjacency(dm, f, PETSC_FALSE, PETSC_TRUE));
4874:   }
4875:   PetscFunctionReturn(PETSC_SUCCESS);
4876: }

4878: static PetscErrorCode DMFieldEnlarge_Static(DM dm, PetscInt NfNew)
4879: {
4880:   RegionField *tmpr;
4881:   PetscInt     Nf = dm->Nf, f;

4883:   PetscFunctionBegin;
4884:   if (Nf >= NfNew) PetscFunctionReturn(PETSC_SUCCESS);
4885:   PetscCall(PetscMalloc1(NfNew, &tmpr));
4886:   for (f = 0; f < Nf; ++f) tmpr[f] = dm->fields[f];
4887:   for (f = Nf; f < NfNew; ++f) {
4888:     tmpr[f].disc        = NULL;
4889:     tmpr[f].label       = NULL;
4890:     tmpr[f].avoidTensor = PETSC_FALSE;
4891:   }
4892:   PetscCall(PetscFree(dm->fields));
4893:   dm->Nf     = NfNew;
4894:   dm->fields = tmpr;
4895:   PetscFunctionReturn(PETSC_SUCCESS);
4896: }

4898: /*@
4899:   DMClearFields - Remove all fields from the `DM`

4901:   Logically Collective

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

4906:   Level: intermediate

4908: .seealso: [](ch_dmbase), `DM`, `DMGetNumFields()`, `DMSetNumFields()`, `DMSetField()`
4909: @*/
4910: PetscErrorCode DMClearFields(DM dm)
4911: {
4912:   PetscInt f;

4914:   PetscFunctionBegin;
4916:   for (f = 0; f < dm->Nf; ++f) {
4917:     PetscCall(PetscObjectDestroy(&dm->fields[f].disc));
4918:     PetscCall(DMLabelDestroy(&dm->fields[f].label));
4919:   }
4920:   PetscCall(PetscFree(dm->fields));
4921:   dm->fields = NULL;
4922:   dm->Nf     = 0;
4923:   PetscFunctionReturn(PETSC_SUCCESS);
4924: }

4926: /*@
4927:   DMGetNumFields - Get the number of fields in the `DM`

4929:   Not Collective

4931:   Input Parameter:
4932: . dm - The `DM`

4934:   Output Parameter:
4935: . numFields - The number of fields

4937:   Level: intermediate

4939: .seealso: [](ch_dmbase), `DM`, `DMSetNumFields()`, `DMSetField()`
4940: @*/
4941: PetscErrorCode DMGetNumFields(DM dm, PetscInt *numFields)
4942: {
4943:   PetscFunctionBegin;
4945:   PetscAssertPointer(numFields, 2);
4946:   *numFields = dm->Nf;
4947:   PetscFunctionReturn(PETSC_SUCCESS);
4948: }

4950: /*@
4951:   DMSetNumFields - Set the number of fields in the `DM`

4953:   Logically Collective

4955:   Input Parameters:
4956: + dm        - The `DM`
4957: - numFields - The number of fields

4959:   Level: intermediate

4961: .seealso: [](ch_dmbase), `DM`, `DMGetNumFields()`, `DMSetField()`
4962: @*/
4963: PetscErrorCode DMSetNumFields(DM dm, PetscInt numFields)
4964: {
4965:   PetscInt Nf, f;

4967:   PetscFunctionBegin;
4969:   PetscCall(DMGetNumFields(dm, &Nf));
4970:   for (f = Nf; f < numFields; ++f) {
4971:     PetscContainer obj;

4973:     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)dm), &obj));
4974:     PetscCall(DMAddField(dm, NULL, (PetscObject)obj));
4975:     PetscCall(PetscContainerDestroy(&obj));
4976:   }
4977:   PetscFunctionReturn(PETSC_SUCCESS);
4978: }

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

4983:   Not Collective

4985:   Input Parameters:
4986: + dm - The `DM`
4987: - f  - The field number

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

4993:   Level: intermediate

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

5008: /* Does not clear the DS */
5009: PetscErrorCode DMSetField_Internal(DM dm, PetscInt f, DMLabel label, PetscObject disc)
5010: {
5011:   PetscFunctionBegin;
5012:   PetscCall(DMFieldEnlarge_Static(dm, f + 1));
5013:   PetscCall(DMLabelDestroy(&dm->fields[f].label));
5014:   PetscCall(PetscObjectDestroy(&dm->fields[f].disc));
5015:   dm->fields[f].label = label;
5016:   dm->fields[f].disc  = disc;
5017:   PetscCall(PetscObjectReference((PetscObject)label));
5018:   PetscCall(PetscObjectReference((PetscObject)disc));
5019:   PetscFunctionReturn(PETSC_SUCCESS);
5020: }

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

5026:   Logically Collective

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

5034:   Level: intermediate

5036: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMGetField()`
5037: @*/
5038: PetscErrorCode DMSetField(DM dm, PetscInt f, DMLabel label, PetscObject disc)
5039: {
5040:   PetscFunctionBegin;
5044:   PetscCheck(f >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be non-negative", f);
5045:   PetscCall(DMSetField_Internal(dm, f, label, disc));
5046:   PetscCall(DMSetDefaultAdjacency_Private(dm, f, disc));
5047:   PetscCall(DMClearDS(dm));
5048:   PetscFunctionReturn(PETSC_SUCCESS);
5049: }

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

5055:   Logically Collective

5057:   Input Parameters:
5058: + dm    - The `DM`
5059: . label - The label indicating the support of the field, or `NULL` for the entire mesh
5060: - disc  - The discretization object

5062:   Level: intermediate

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

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

5071: .seealso: [](ch_dmbase), `DM`, `DMSetLabel()`, `DMSetField()`, `DMGetField()`, `PetscFE`
5072: @*/
5073: PetscErrorCode DMAddField(DM dm, DMLabel label, PetscObject disc)
5074: {
5075:   PetscInt Nf = dm->Nf;

5077:   PetscFunctionBegin;
5081:   PetscCall(DMFieldEnlarge_Static(dm, Nf + 1));
5082:   dm->fields[Nf].label = label;
5083:   dm->fields[Nf].disc  = disc;
5084:   PetscCall(PetscObjectReference((PetscObject)label));
5085:   PetscCall(PetscObjectReference((PetscObject)disc));
5086:   PetscCall(DMSetDefaultAdjacency_Private(dm, Nf, disc));
5087:   PetscCall(DMClearDS(dm));
5088:   PetscFunctionReturn(PETSC_SUCCESS);
5089: }

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

5094:   Logically Collective

5096:   Input Parameters:
5097: + dm          - The `DM`
5098: . f           - The field index
5099: - avoidTensor - `PETSC_TRUE` to skip defining the field on tensor cells

5101:   Level: intermediate

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

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

5116:   Not Collective

5118:   Input Parameters:
5119: + dm - The `DM`
5120: - f  - The field index

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

5125:   Level: intermediate

5127: .seealso: [](ch_dmbase), `DM`, `DMAddField()`, `DMSetField()`, `DMGetField()`, `DMSetFieldAvoidTensor()`
5128: @*/
5129: PetscErrorCode DMGetFieldAvoidTensor(DM dm, PetscInt f, PetscBool *avoidTensor)
5130: {
5131:   PetscFunctionBegin;
5132:   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);
5133:   *avoidTensor = dm->fields[f].avoidTensor;
5134:   PetscFunctionReturn(PETSC_SUCCESS);
5135: }

5137: /*@
5138:   DMCopyFields - Copy the discretizations for the `DM` into another `DM`

5140:   Collective

5142:   Input Parameter:
5143: . dm - The `DM`

5145:   Output Parameter:
5146: . newdm - The `DM`

5148:   Level: advanced

5150: .seealso: [](ch_dmbase), `DM`, `DMGetField()`, `DMSetField()`, `DMAddField()`, `DMCopyDS()`, `DMGetDS()`, `DMGetCellDS()`
5151: @*/
5152: PetscErrorCode DMCopyFields(DM dm, DM newdm)
5153: {
5154:   PetscInt Nf, f;

5156:   PetscFunctionBegin;
5157:   if (dm == newdm) PetscFunctionReturn(PETSC_SUCCESS);
5158:   PetscCall(DMGetNumFields(dm, &Nf));
5159:   PetscCall(DMClearFields(newdm));
5160:   for (f = 0; f < Nf; ++f) {
5161:     DMLabel     label;
5162:     PetscObject field;
5163:     PetscBool   useCone, useClosure;

5165:     PetscCall(DMGetField(dm, f, &label, &field));
5166:     PetscCall(DMSetField(newdm, f, label, field));
5167:     PetscCall(DMGetAdjacency(dm, f, &useCone, &useClosure));
5168:     PetscCall(DMSetAdjacency(newdm, f, useCone, useClosure));
5169:   }
5170:   PetscFunctionReturn(PETSC_SUCCESS);
5171: }

5173: /*@
5174:   DMGetAdjacency - Returns the flags for determining variable influence

5176:   Not Collective

5178:   Input Parameters:
5179: + dm - The `DM` object
5180: - f  - The field number, or `PETSC_DEFAULT` for the default adjacency

5182:   Output Parameters:
5183: + useCone    - Flag for variable influence starting with the cone operation
5184: - useClosure - Flag for variable influence using transitive closure

5186:   Level: developer

5188:   Notes:
5189: .vb
5190:      FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
5191:      FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
5192:      FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
5193: .ve
5194:   Further explanation can be found in the User's Manual Section on the Influence of Variables on One Another.

5196: .seealso: [](ch_dmbase), `DM`, `DMSetAdjacency()`, `DMGetField()`, `DMSetField()`
5197: @*/
5198: PetscErrorCode DMGetAdjacency(DM dm, PetscInt f, PetscBool *useCone, PetscBool *useClosure)
5199: {
5200:   PetscFunctionBegin;
5202:   if (useCone) PetscAssertPointer(useCone, 3);
5203:   if (useClosure) PetscAssertPointer(useClosure, 4);
5204:   if (f < 0) {
5205:     if (useCone) *useCone = dm->adjacency[0];
5206:     if (useClosure) *useClosure = dm->adjacency[1];
5207:   } else {
5208:     PetscInt Nf;

5210:     PetscCall(DMGetNumFields(dm, &Nf));
5211:     PetscCheck(f < Nf, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, Nf);
5212:     if (useCone) *useCone = dm->fields[f].adjacency[0];
5213:     if (useClosure) *useClosure = dm->fields[f].adjacency[1];
5214:   }
5215:   PetscFunctionReturn(PETSC_SUCCESS);
5216: }

5218: /*@
5219:   DMSetAdjacency - Set the flags for determining variable influence

5221:   Not Collective

5223:   Input Parameters:
5224: + dm         - The `DM` object
5225: . f          - The field number
5226: . useCone    - Flag for variable influence starting with the cone operation
5227: - useClosure - Flag for variable influence using transitive closure

5229:   Level: developer

5231:   Notes:
5232: .vb
5233:      FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
5234:      FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
5235:      FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
5236: .ve
5237:   Further explanation can be found in the User's Manual Section on the Influence of Variables on One Another.

5239: .seealso: [](ch_dmbase), `DM`, `DMGetAdjacency()`, `DMGetField()`, `DMSetField()`
5240: @*/
5241: PetscErrorCode DMSetAdjacency(DM dm, PetscInt f, PetscBool useCone, PetscBool useClosure)
5242: {
5243:   PetscFunctionBegin;
5245:   if (f < 0) {
5246:     dm->adjacency[0] = useCone;
5247:     dm->adjacency[1] = useClosure;
5248:   } else {
5249:     PetscInt Nf;

5251:     PetscCall(DMGetNumFields(dm, &Nf));
5252:     PetscCheck(f < Nf, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", f, Nf);
5253:     dm->fields[f].adjacency[0] = useCone;
5254:     dm->fields[f].adjacency[1] = useClosure;
5255:   }
5256:   PetscFunctionReturn(PETSC_SUCCESS);
5257: }

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

5262:   Not collective

5264:   Input Parameter:
5265: . dm - The `DM` object

5267:   Output Parameters:
5268: + useCone    - Flag for variable influence starting with the cone operation
5269: - useClosure - Flag for variable influence using transitive closure

5271:   Level: developer

5273:   Notes:
5274: .vb
5275:      FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
5276:      FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
5277:      FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
5278: .ve

5280: .seealso: [](ch_dmbase), `DM`, `DMSetBasicAdjacency()`, `DMGetField()`, `DMSetField()`
5281: @*/
5282: PetscErrorCode DMGetBasicAdjacency(DM dm, PetscBool *useCone, PetscBool *useClosure)
5283: {
5284:   PetscInt Nf;

5286:   PetscFunctionBegin;
5288:   if (useCone) PetscAssertPointer(useCone, 2);
5289:   if (useClosure) PetscAssertPointer(useClosure, 3);
5290:   PetscCall(DMGetNumFields(dm, &Nf));
5291:   if (!Nf) {
5292:     PetscCall(DMGetAdjacency(dm, PETSC_DEFAULT, useCone, useClosure));
5293:   } else {
5294:     PetscCall(DMGetAdjacency(dm, 0, useCone, useClosure));
5295:   }
5296:   PetscFunctionReturn(PETSC_SUCCESS);
5297: }

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

5302:   Not Collective

5304:   Input Parameters:
5305: + dm         - The `DM` object
5306: . useCone    - Flag for variable influence starting with the cone operation
5307: - useClosure - Flag for variable influence using transitive closure

5309:   Level: developer

5311:   Notes:
5312: .vb
5313:      FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
5314:      FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
5315:      FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
5316: .ve

5318: .seealso: [](ch_dmbase), `DM`, `DMGetBasicAdjacency()`, `DMGetField()`, `DMSetField()`
5319: @*/
5320: PetscErrorCode DMSetBasicAdjacency(DM dm, PetscBool useCone, PetscBool useClosure)
5321: {
5322:   PetscInt Nf;

5324:   PetscFunctionBegin;
5326:   PetscCall(DMGetNumFields(dm, &Nf));
5327:   if (!Nf) {
5328:     PetscCall(DMSetAdjacency(dm, PETSC_DEFAULT, useCone, useClosure));
5329:   } else {
5330:     PetscCall(DMSetAdjacency(dm, 0, useCone, useClosure));
5331:   }
5332:   PetscFunctionReturn(PETSC_SUCCESS);
5333: }

5335: PetscErrorCode DMCompleteBCLabels_Internal(DM dm)
5336: {
5337:   DM           plex;
5338:   DMLabel     *labels, *glabels;
5339:   const char **names;
5340:   char        *sendNames, *recvNames;
5341:   PetscInt     Nds, s, maxLabels = 0, maxLen = 0, gmaxLen, Nl = 0, gNl, l, gl, m;
5342:   size_t       len;
5343:   MPI_Comm     comm;
5344:   PetscMPIInt  rank, size, p, *counts, *displs;

5346:   PetscFunctionBegin;
5347:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
5348:   PetscCallMPI(MPI_Comm_size(comm, &size));
5349:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5350:   PetscCall(DMGetNumDS(dm, &Nds));
5351:   for (s = 0; s < Nds; ++s) {
5352:     PetscDS  dsBC;
5353:     PetscInt numBd;

5355:     PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC, NULL));
5356:     PetscCall(PetscDSGetNumBoundary(dsBC, &numBd));
5357:     maxLabels += numBd;
5358:   }
5359:   PetscCall(PetscCalloc1(maxLabels, &labels));
5360:   /* Get list of labels to be completed */
5361:   for (s = 0; s < Nds; ++s) {
5362:     PetscDS  dsBC;
5363:     PetscInt numBd, bd;

5365:     PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC, NULL));
5366:     PetscCall(PetscDSGetNumBoundary(dsBC, &numBd));
5367:     for (bd = 0; bd < numBd; ++bd) {
5368:       DMLabel      label;
5369:       PetscInt     field;
5370:       PetscObject  obj;
5371:       PetscClassId id;

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

5421: static PetscErrorCode DMDSEnlarge_Static(DM dm, PetscInt NdsNew)
5422: {
5423:   DMSpace *tmpd;
5424:   PetscInt Nds = dm->Nds, s;

5426:   PetscFunctionBegin;
5427:   if (Nds >= NdsNew) PetscFunctionReturn(PETSC_SUCCESS);
5428:   PetscCall(PetscMalloc1(NdsNew, &tmpd));
5429:   for (s = 0; s < Nds; ++s) tmpd[s] = dm->probs[s];
5430:   for (s = Nds; s < NdsNew; ++s) {
5431:     tmpd[s].ds     = NULL;
5432:     tmpd[s].label  = NULL;
5433:     tmpd[s].fields = NULL;
5434:   }
5435:   PetscCall(PetscFree(dm->probs));
5436:   dm->Nds   = NdsNew;
5437:   dm->probs = tmpd;
5438:   PetscFunctionReturn(PETSC_SUCCESS);
5439: }

5441: /*@
5442:   DMGetNumDS - Get the number of discrete systems in the `DM`

5444:   Not Collective

5446:   Input Parameter:
5447: . dm - The `DM`

5449:   Output Parameter:
5450: . Nds - The number of `PetscDS` objects

5452:   Level: intermediate

5454: .seealso: [](ch_dmbase), `DM`, `DMGetDS()`, `DMGetCellDS()`
5455: @*/
5456: PetscErrorCode DMGetNumDS(DM dm, PetscInt *Nds)
5457: {
5458:   PetscFunctionBegin;
5460:   PetscAssertPointer(Nds, 2);
5461:   *Nds = dm->Nds;
5462:   PetscFunctionReturn(PETSC_SUCCESS);
5463: }

5465: /*@
5466:   DMClearDS - Remove all discrete systems from the `DM`

5468:   Logically Collective

5470:   Input Parameter:
5471: . dm - The `DM`

5473:   Level: intermediate

5475: .seealso: [](ch_dmbase), `DM`, `DMGetNumDS()`, `DMGetDS()`, `DMSetField()`
5476: @*/
5477: PetscErrorCode DMClearDS(DM dm)
5478: {
5479:   PetscInt s;

5481:   PetscFunctionBegin;
5483:   for (s = 0; s < dm->Nds; ++s) {
5484:     PetscCall(PetscDSDestroy(&dm->probs[s].ds));
5485:     PetscCall(PetscDSDestroy(&dm->probs[s].dsIn));
5486:     PetscCall(DMLabelDestroy(&dm->probs[s].label));
5487:     PetscCall(ISDestroy(&dm->probs[s].fields));
5488:   }
5489:   PetscCall(PetscFree(dm->probs));
5490:   dm->probs = NULL;
5491:   dm->Nds   = 0;
5492:   PetscFunctionReturn(PETSC_SUCCESS);
5493: }

5495: /*@
5496:   DMGetDS - Get the default `PetscDS`

5498:   Not Collective

5500:   Input Parameter:
5501: . dm - The `DM`

5503:   Output Parameter:
5504: . ds - The default `PetscDS`

5506:   Level: intermediate

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

5520: /*@
5521:   DMGetCellDS - Get the `PetscDS` defined on a given cell

5523:   Not Collective

5525:   Input Parameters:
5526: + dm    - The `DM`
5527: - point - Cell for the `PetscDS`

5529:   Output Parameters:
5530: + ds   - The `PetscDS` defined on the given cell
5531: - dsIn - The `PetscDS` for input on the given cell, or NULL if the same ds

5533:   Level: developer

5535: .seealso: [](ch_dmbase), `DM`, `DMGetDS()`, `DMSetRegionDS()`
5536: @*/
5537: PetscErrorCode DMGetCellDS(DM dm, PetscInt point, PetscDS *ds, PetscDS *dsIn)
5538: {
5539:   PetscDS  dsDef = NULL;
5540:   PetscInt s;

5542:   PetscFunctionBeginHot;
5544:   if (ds) PetscAssertPointer(ds, 3);
5545:   if (dsIn) PetscAssertPointer(dsIn, 4);
5546:   PetscCheck(point >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Mesh point cannot be negative: %" PetscInt_FMT, point);
5547:   if (ds) *ds = NULL;
5548:   if (dsIn) *dsIn = NULL;
5549:   for (s = 0; s < dm->Nds; ++s) {
5550:     PetscInt val;

5552:     if (!dm->probs[s].label) {
5553:       dsDef = dm->probs[s].ds;
5554:     } else {
5555:       PetscCall(DMLabelGetValue(dm->probs[s].label, point, &val));
5556:       if (val >= 0) {
5557:         if (ds) *ds = dm->probs[s].ds;
5558:         if (dsIn) *dsIn = dm->probs[s].dsIn;
5559:         break;
5560:       }
5561:     }
5562:   }
5563:   if (ds && !*ds) *ds = dsDef;
5564:   PetscFunctionReturn(PETSC_SUCCESS);
5565: }

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

5570:   Not Collective

5572:   Input Parameters:
5573: + dm    - The `DM`
5574: - label - The `DMLabel` defining the mesh region, or `NULL` for the entire mesh

5576:   Output Parameters:
5577: + fields - The `IS` containing the `DM` field numbers for the fields in this `PetscDS`, or `NULL`
5578: . ds     - The `PetscDS` defined on the given region, or `NULL`
5579: - dsIn   - The `PetscDS` for input in the given region, or `NULL`

5581:   Level: advanced

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

5588: .seealso: [](ch_dmbase), `DM`, `DMGetRegionNumDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5589: @*/
5590: PetscErrorCode DMGetRegionDS(DM dm, DMLabel label, IS *fields, PetscDS *ds, PetscDS *dsIn)
5591: {
5592:   PetscInt Nds = dm->Nds, s;

5594:   PetscFunctionBegin;
5597:   if (fields) {
5598:     PetscAssertPointer(fields, 3);
5599:     *fields = NULL;
5600:   }
5601:   if (ds) {
5602:     PetscAssertPointer(ds, 4);
5603:     *ds = NULL;
5604:   }
5605:   if (dsIn) {
5606:     PetscAssertPointer(dsIn, 5);
5607:     *dsIn = NULL;
5608:   }
5609:   for (s = 0; s < Nds; ++s) {
5610:     if (dm->probs[s].label == label || !dm->probs[s].label) {
5611:       if (fields) *fields = dm->probs[s].fields;
5612:       if (ds) *ds = dm->probs[s].ds;
5613:       if (dsIn) *dsIn = dm->probs[s].dsIn;
5614:       if (dm->probs[s].label) PetscFunctionReturn(PETSC_SUCCESS);
5615:     }
5616:   }
5617:   PetscFunctionReturn(PETSC_SUCCESS);
5618: }

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

5623:   Collective

5625:   Input Parameters:
5626: + dm     - The `DM`
5627: . label  - The `DMLabel` defining the mesh region, or `NULL` for the entire mesh
5628: . fields - The `IS` containing the `DM` field numbers for the fields in this `PetscDS`, or `NULL` for all fields
5629: . ds     - The `PetscDS` defined on the given region
5630: - dsIn   - The `PetscDS` for input on the given cell, or `NULL` if it is the same `PetscDS`

5632:   Level: advanced

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

5638: .seealso: [](ch_dmbase), `DM`, `DMGetRegionDS()`, `DMSetRegionNumDS()`, `DMGetDS()`, `DMGetCellDS()`
5639: @*/
5640: PetscErrorCode DMSetRegionDS(DM dm, DMLabel label, IS fields, PetscDS ds, PetscDS dsIn)
5641: {
5642:   PetscInt Nds = dm->Nds, s;

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

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

5679:   Not Collective

5681:   Input Parameters:
5682: + dm  - The `DM`
5683: - num - The region number, in [0, Nds)

5685:   Output Parameters:
5686: + label  - The region label, or `NULL`
5687: . fields - The `IS` containing the `DM` field numbers for the fields in this `PetscDS`, or `NULL`
5688: . ds     - The `PetscDS` defined on the given region, or `NULL`
5689: - dsIn   - The `PetscDS` for input in the given region, or `NULL`

5691:   Level: advanced

5693: .seealso: [](ch_dmbase), `DM`, `DMGetRegionDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5694: @*/
5695: PetscErrorCode DMGetRegionNumDS(DM dm, PetscInt num, DMLabel *label, IS *fields, PetscDS *ds, PetscDS *dsIn)
5696: {
5697:   PetscInt Nds;

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

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

5725:   Not Collective

5727:   Input Parameters:
5728: + dm     - The `DM`
5729: . num    - The region number, in [0, Nds)
5730: . label  - The region label, or `NULL`
5731: . fields - The `IS` containing the `DM` field numbers for the fields in this `PetscDS`, or `NULL` to prevent setting
5732: . ds     - The `PetscDS` defined on the given region, or `NULL` to prevent setting
5733: - dsIn   - The `PetscDS` for input on the given cell, or `NULL` if it is the same `PetscDS`

5735:   Level: advanced

5737: .seealso: [](ch_dmbase), `DM`, `DMGetRegionDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5738: @*/
5739: PetscErrorCode DMSetRegionNumDS(DM dm, PetscInt num, DMLabel label, IS fields, PetscDS ds, PetscDS dsIn)
5740: {
5741:   PetscInt Nds;

5743:   PetscFunctionBegin;
5746:   PetscCall(DMGetNumDS(dm, &Nds));
5747:   PetscCheck((num >= 0) && (num < Nds), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Region number %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", num, Nds);
5748:   PetscCall(PetscObjectReference((PetscObject)label));
5749:   PetscCall(DMLabelDestroy(&dm->probs[num].label));
5750:   dm->probs[num].label = label;
5751:   if (fields) {
5753:     PetscCall(PetscObjectReference((PetscObject)fields));
5754:     PetscCall(ISDestroy(&dm->probs[num].fields));
5755:     dm->probs[num].fields = fields;
5756:   }
5757:   if (ds) {
5759:     PetscCall(PetscObjectReference((PetscObject)ds));
5760:     PetscCall(PetscDSDestroy(&dm->probs[num].ds));
5761:     dm->probs[num].ds = ds;
5762:   }
5763:   if (dsIn) {
5765:     PetscCall(PetscObjectReference((PetscObject)dsIn));
5766:     PetscCall(PetscDSDestroy(&dm->probs[num].dsIn));
5767:     dm->probs[num].dsIn = dsIn;
5768:   }
5769:   PetscFunctionReturn(PETSC_SUCCESS);
5770: }

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

5775:   Not Collective

5777:   Input Parameters:
5778: + dm - The `DM`
5779: - ds - The `PetscDS` defined on the given region

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

5784:   Level: advanced

5786: .seealso: [](ch_dmbase), `DM`, `DMGetRegionNumDS()`, `DMGetRegionDS()`, `DMSetRegionDS()`, `DMGetDS()`, `DMGetCellDS()`
5787: @*/
5788: PetscErrorCode DMFindRegionNum(DM dm, PetscDS ds, PetscInt *num)
5789: {
5790:   PetscInt Nds, n;

5792:   PetscFunctionBegin;
5795:   PetscAssertPointer(num, 3);
5796:   PetscCall(DMGetNumDS(dm, &Nds));
5797:   for (n = 0; n < Nds; ++n)
5798:     if (ds == dm->probs[n].ds) break;
5799:   if (n >= Nds) *num = -1;
5800:   else *num = n;
5801:   PetscFunctionReturn(PETSC_SUCCESS);
5802: }

5804: /*@C
5805:   DMCreateFEDefault - Create a `PetscFE` based on the celltype for the mesh

5807:   Not Collective

5809:   Input Parameters:
5810: + dm     - The `DM`
5811: . Nc     - The number of components for the field
5812: . prefix - The options prefix for the output `PetscFE`, or `NULL`
5813: - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree

5815:   Output Parameter:
5816: . fem - The `PetscFE`

5818:   Level: intermediate

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

5823: .seealso: [](ch_dmbase), `DM`, `PetscFECreateByCell()`, `DMAddField()`, `DMCreateDS()`, `DMGetCellDS()`, `DMGetRegionDS()`
5824: @*/
5825: PetscErrorCode DMCreateFEDefault(DM dm, PetscInt Nc, const char prefix[], PetscInt qorder, PetscFE *fem)
5826: {
5827:   DMPolytopeType ct;
5828:   PetscInt       dim, cStart;

5830:   PetscFunctionBegin;
5833:   if (prefix) PetscAssertPointer(prefix, 3);
5835:   PetscAssertPointer(fem, 5);
5836:   PetscCall(DMGetDimension(dm, &dim));
5837:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
5838:   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
5839:   PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, Nc, ct, prefix, qorder, fem));
5840:   PetscFunctionReturn(PETSC_SUCCESS);
5841: }

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

5846:   Collective

5848:   Input Parameter:
5849: . dm - The `DM`

5851:   Options Database Key:
5852: . -dm_petscds_view - View all the `PetscDS` objects in this `DM`

5854:   Level: intermediate

5856: .seealso: [](ch_dmbase), `DM`, `DMSetField`, `DMAddField()`, `DMGetDS()`, `DMGetCellDS()`, `DMGetRegionDS()`, `DMSetRegionDS()`
5857: @*/
5858: PetscErrorCode DMCreateDS(DM dm)
5859: {
5860:   MPI_Comm  comm;
5861:   PetscDS   dsDef;
5862:   DMLabel  *labelSet;
5863:   PetscInt  dE, Nf = dm->Nf, f, s, Nl, l, Ndef, k;
5864:   PetscBool doSetup = PETSC_TRUE, flg;

5866:   PetscFunctionBegin;
5868:   if (!dm->fields) PetscFunctionReturn(PETSC_SUCCESS);
5869:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
5870:   PetscCall(DMGetCoordinateDim(dm, &dE));
5871:   /* Determine how many regions we have */
5872:   PetscCall(PetscMalloc1(Nf, &labelSet));
5873:   Nl   = 0;
5874:   Ndef = 0;
5875:   for (f = 0; f < Nf; ++f) {
5876:     DMLabel  label = dm->fields[f].label;
5877:     PetscInt l;

5879: #ifdef PETSC_HAVE_LIBCEED
5880:     /* Move CEED context to discretizations */
5881:     {
5882:       PetscClassId id;

5884:       PetscCall(PetscObjectGetClassId(dm->fields[f].disc, &id));
5885:       if (id == PETSCFE_CLASSID) {
5886:         Ceed ceed;

5888:         PetscCall(DMGetCeed(dm, &ceed));
5889:         PetscCall(PetscFESetCeed((PetscFE)dm->fields[f].disc, ceed));
5890:       }
5891:     }
5892: #endif
5893:     if (!label) {
5894:       ++Ndef;
5895:       continue;
5896:     }
5897:     for (l = 0; l < Nl; ++l)
5898:       if (label == labelSet[l]) break;
5899:     if (l < Nl) continue;
5900:     labelSet[Nl++] = label;
5901:   }
5902:   /* Create default DS if there are no labels to intersect with */
5903:   PetscCall(DMGetRegionDS(dm, NULL, NULL, &dsDef, NULL));
5904:   if (!dsDef && Ndef && !Nl) {
5905:     IS        fields;
5906:     PetscInt *fld, nf;

5908:     for (f = 0, nf = 0; f < Nf; ++f)
5909:       if (!dm->fields[f].label) ++nf;
5910:     PetscCheck(nf, comm, PETSC_ERR_PLIB, "All fields have labels, but we are trying to create a default DS");
5911:     PetscCall(PetscMalloc1(nf, &fld));
5912:     for (f = 0, nf = 0; f < Nf; ++f)
5913:       if (!dm->fields[f].label) fld[nf++] = f;
5914:     PetscCall(ISCreate(PETSC_COMM_SELF, &fields));
5915:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fields, "dm_fields_"));
5916:     PetscCall(ISSetType(fields, ISGENERAL));
5917:     PetscCall(ISGeneralSetIndices(fields, nf, fld, PETSC_OWN_POINTER));

5919:     PetscCall(PetscDSCreate(PETSC_COMM_SELF, &dsDef));
5920:     PetscCall(DMSetRegionDS(dm, NULL, fields, dsDef, NULL));
5921:     PetscCall(PetscDSDestroy(&dsDef));
5922:     PetscCall(ISDestroy(&fields));
5923:   }
5924:   PetscCall(DMGetRegionDS(dm, NULL, NULL, &dsDef, NULL));
5925:   if (dsDef) PetscCall(PetscDSSetCoordinateDimension(dsDef, dE));
5926:   /* Intersect labels with default fields */
5927:   if (Ndef && Nl) {
5928:     DM              plex;
5929:     DMLabel         cellLabel;
5930:     IS              fieldIS, allcellIS, defcellIS = NULL;
5931:     PetscInt       *fields;
5932:     const PetscInt *cells;
5933:     PetscInt        depth, nf = 0, n, c;

5935:     PetscCall(DMConvert(dm, DMPLEX, &plex));
5936:     PetscCall(DMPlexGetDepth(plex, &depth));
5937:     PetscCall(DMGetStratumIS(plex, "dim", depth, &allcellIS));
5938:     if (!allcellIS) PetscCall(DMGetStratumIS(plex, "depth", depth, &allcellIS));
5939:     /* TODO This looks like it only works for one label */
5940:     for (l = 0; l < Nl; ++l) {
5941:       DMLabel label = labelSet[l];
5942:       IS      pointIS;

5944:       PetscCall(ISDestroy(&defcellIS));
5945:       PetscCall(DMLabelGetStratumIS(label, 1, &pointIS));
5946:       PetscCall(ISDifference(allcellIS, pointIS, &defcellIS));
5947:       PetscCall(ISDestroy(&pointIS));
5948:     }
5949:     PetscCall(ISDestroy(&allcellIS));

5951:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "defaultCells", &cellLabel));
5952:     PetscCall(ISGetLocalSize(defcellIS, &n));
5953:     PetscCall(ISGetIndices(defcellIS, &cells));
5954:     for (c = 0; c < n; ++c) PetscCall(DMLabelSetValue(cellLabel, cells[c], 1));
5955:     PetscCall(ISRestoreIndices(defcellIS, &cells));
5956:     PetscCall(ISDestroy(&defcellIS));
5957:     PetscCall(DMPlexLabelComplete(plex, cellLabel));

5959:     PetscCall(PetscMalloc1(Ndef, &fields));
5960:     for (f = 0; f < Nf; ++f)
5961:       if (!dm->fields[f].label) fields[nf++] = f;
5962:     PetscCall(ISCreate(PETSC_COMM_SELF, &fieldIS));
5963:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fieldIS, "dm_fields_"));
5964:     PetscCall(ISSetType(fieldIS, ISGENERAL));
5965:     PetscCall(ISGeneralSetIndices(fieldIS, nf, fields, PETSC_OWN_POINTER));

5967:     PetscCall(PetscDSCreate(PETSC_COMM_SELF, &dsDef));
5968:     PetscCall(DMSetRegionDS(dm, cellLabel, fieldIS, dsDef, NULL));
5969:     PetscCall(PetscDSSetCoordinateDimension(dsDef, dE));
5970:     PetscCall(DMLabelDestroy(&cellLabel));
5971:     PetscCall(PetscDSDestroy(&dsDef));
5972:     PetscCall(ISDestroy(&fieldIS));
5973:     PetscCall(DMDestroy(&plex));
5974:   }
5975:   /* Create label DSes
5976:      - WE ONLY SUPPORT IDENTICAL OR DISJOINT LABELS
5977:   */
5978:   /* TODO Should check that labels are disjoint */
5979:   for (l = 0; l < Nl; ++l) {
5980:     DMLabel   label = labelSet[l];
5981:     PetscDS   ds, dsIn = NULL;
5982:     IS        fields;
5983:     PetscInt *fld, nf;

5985:     PetscCall(PetscDSCreate(PETSC_COMM_SELF, &ds));
5986:     for (f = 0, nf = 0; f < Nf; ++f)
5987:       if (label == dm->fields[f].label || !dm->fields[f].label) ++nf;
5988:     PetscCall(PetscMalloc1(nf, &fld));
5989:     for (f = 0, nf = 0; f < Nf; ++f)
5990:       if (label == dm->fields[f].label || !dm->fields[f].label) fld[nf++] = f;
5991:     PetscCall(ISCreate(PETSC_COMM_SELF, &fields));
5992:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fields, "dm_fields_"));
5993:     PetscCall(ISSetType(fields, ISGENERAL));
5994:     PetscCall(ISGeneralSetIndices(fields, nf, fld, PETSC_OWN_POINTER));
5995:     PetscCall(PetscDSSetCoordinateDimension(ds, dE));
5996:     {
5997:       DMPolytopeType ct;
5998:       PetscInt       lStart, lEnd;
5999:       PetscBool      isCohesiveLocal = PETSC_FALSE, isCohesive;

6001:       PetscCall(DMLabelGetBounds(label, &lStart, &lEnd));
6002:       if (lStart >= 0) {
6003:         PetscCall(DMPlexGetCellType(dm, lStart, &ct));
6004:         switch (ct) {
6005:         case DM_POLYTOPE_POINT_PRISM_TENSOR:
6006:         case DM_POLYTOPE_SEG_PRISM_TENSOR:
6007:         case DM_POLYTOPE_TRI_PRISM_TENSOR:
6008:         case DM_POLYTOPE_QUAD_PRISM_TENSOR:
6009:           isCohesiveLocal = PETSC_TRUE;
6010:           break;
6011:         default:
6012:           break;
6013:         }
6014:       }
6015:       PetscCall(MPIU_Allreduce(&isCohesiveLocal, &isCohesive, 1, MPIU_BOOL, MPI_LOR, comm));
6016:       if (isCohesive) {
6017:         PetscCall(PetscDSCreate(PETSC_COMM_SELF, &dsIn));
6018:         PetscCall(PetscDSSetCoordinateDimension(dsIn, dE));
6019:       }
6020:       for (f = 0, nf = 0; f < Nf; ++f) {
6021:         if (label == dm->fields[f].label || !dm->fields[f].label) {
6022:           if (label == dm->fields[f].label) {
6023:             PetscCall(PetscDSSetDiscretization(ds, nf, NULL));
6024:             PetscCall(PetscDSSetCohesive(ds, nf, isCohesive));
6025:             if (dsIn) {
6026:               PetscCall(PetscDSSetDiscretization(dsIn, nf, NULL));
6027:               PetscCall(PetscDSSetCohesive(dsIn, nf, isCohesive));
6028:             }
6029:           }
6030:           ++nf;
6031:         }
6032:       }
6033:     }
6034:     PetscCall(DMSetRegionDS(dm, label, fields, ds, dsIn));
6035:     PetscCall(ISDestroy(&fields));
6036:     PetscCall(PetscDSDestroy(&ds));
6037:     PetscCall(PetscDSDestroy(&dsIn));
6038:   }
6039:   PetscCall(PetscFree(labelSet));
6040:   /* Set fields in DSes */
6041:   for (s = 0; s < dm->Nds; ++s) {
6042:     PetscDS         ds     = dm->probs[s].ds;
6043:     PetscDS         dsIn   = dm->probs[s].dsIn;
6044:     IS              fields = dm->probs[s].fields;
6045:     const PetscInt *fld;
6046:     PetscInt        nf, dsnf;
6047:     PetscBool       isCohesive;

6049:     PetscCall(PetscDSGetNumFields(ds, &dsnf));
6050:     PetscCall(PetscDSIsCohesive(ds, &isCohesive));
6051:     PetscCall(ISGetLocalSize(fields, &nf));
6052:     PetscCall(ISGetIndices(fields, &fld));
6053:     for (f = 0; f < nf; ++f) {
6054:       PetscObject  disc = dm->fields[fld[f]].disc;
6055:       PetscBool    isCohesiveField;
6056:       PetscClassId id;

6058:       /* Handle DS with no fields */
6059:       if (dsnf) PetscCall(PetscDSGetCohesive(ds, f, &isCohesiveField));
6060:       /* If this is a cohesive cell, then regular fields need the lower dimensional discretization */
6061:       if (isCohesive) {
6062:         if (!isCohesiveField) {
6063:           PetscObject bdDisc;

6065:           PetscCall(PetscFEGetHeightSubspace((PetscFE)disc, 1, (PetscFE *)&bdDisc));
6066:           PetscCall(PetscDSSetDiscretization(ds, f, bdDisc));
6067:           PetscCall(PetscDSSetDiscretization(dsIn, f, disc));
6068:         } else {
6069:           PetscCall(PetscDSSetDiscretization(ds, f, disc));
6070:           PetscCall(PetscDSSetDiscretization(dsIn, f, disc));
6071:         }
6072:       } else {
6073:         PetscCall(PetscDSSetDiscretization(ds, f, disc));
6074:       }
6075:       /* We allow people to have placeholder fields and construct the Section by hand */
6076:       PetscCall(PetscObjectGetClassId(disc, &id));
6077:       if ((id != PETSCFE_CLASSID) && (id != PETSCFV_CLASSID)) doSetup = PETSC_FALSE;
6078:     }
6079:     PetscCall(ISRestoreIndices(fields, &fld));
6080:   }
6081:   /* Allow k-jet tabulation */
6082:   PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)dm)->prefix, "-dm_ds_jet_degree", &k, &flg));
6083:   if (flg) {
6084:     for (s = 0; s < dm->Nds; ++s) {
6085:       PetscDS  ds   = dm->probs[s].ds;
6086:       PetscDS  dsIn = dm->probs[s].dsIn;
6087:       PetscInt Nf, f;

6089:       PetscCall(PetscDSGetNumFields(ds, &Nf));
6090:       for (f = 0; f < Nf; ++f) {
6091:         PetscCall(PetscDSSetJetDegree(ds, f, k));
6092:         if (dsIn) PetscCall(PetscDSSetJetDegree(dsIn, f, k));
6093:       }
6094:     }
6095:   }
6096:   /* Setup DSes */
6097:   if (doSetup) {
6098:     for (s = 0; s < dm->Nds; ++s) {
6099:       if (dm->setfromoptionscalled) {
6100:         PetscCall(PetscDSSetFromOptions(dm->probs[s].ds));
6101:         if (dm->probs[s].dsIn) PetscCall(PetscDSSetFromOptions(dm->probs[s].dsIn));
6102:       }
6103:       PetscCall(PetscDSSetUp(dm->probs[s].ds));
6104:       if (dm->probs[s].dsIn) PetscCall(PetscDSSetUp(dm->probs[s].dsIn));
6105:     }
6106:   }
6107:   PetscFunctionReturn(PETSC_SUCCESS);
6108: }

6110: /*@
6111:   DMUseTensorOrder - Use a tensor product closure ordering for the default section

6113:   Input Parameters:
6114: + dm     - The DM
6115: - tensor - Flag for tensor order

6117:   Level: developer

6119: .seealso: `DMPlexSetClosurePermutationTensor()`, `PetscSectionResetClosurePermutation()`
6120: @*/
6121: PetscErrorCode DMUseTensorOrder(DM dm, PetscBool tensor)
6122: {
6123:   PetscInt  Nf;
6124:   PetscBool reorder = PETSC_TRUE, isPlex;

6126:   PetscFunctionBegin;
6127:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex));
6128:   PetscCall(DMGetNumFields(dm, &Nf));
6129:   for (PetscInt f = 0; f < Nf; ++f) {
6130:     PetscObject  obj;
6131:     PetscClassId id;

6133:     PetscCall(DMGetField(dm, f, NULL, &obj));
6134:     PetscCall(PetscObjectGetClassId(obj, &id));
6135:     if (id == PETSCFE_CLASSID) {
6136:       PetscSpace sp;
6137:       PetscBool  tensor;

6139:       PetscCall(PetscFEGetBasisSpace((PetscFE)obj, &sp));
6140:       PetscCall(PetscSpacePolynomialGetTensor(sp, &tensor));
6141:       reorder = reorder && tensor ? PETSC_TRUE : PETSC_FALSE;
6142:     } else reorder = PETSC_FALSE;
6143:   }
6144:   if (tensor) {
6145:     if (reorder && isPlex) PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
6146:   } else {
6147:     PetscSection s;

6149:     PetscCall(DMGetLocalSection(dm, &s));
6150:     if (s) PetscCall(PetscSectionResetClosurePermutation(s));
6151:   }
6152:   PetscFunctionReturn(PETSC_SUCCESS);
6153: }

6155: /*@
6156:   DMComputeExactSolution - Compute the exact solution for a given `DM`, using the `PetscDS` information.

6158:   Collective

6160:   Input Parameters:
6161: + dm   - The `DM`
6162: - time - The time

6164:   Output Parameters:
6165: + u   - The vector will be filled with exact solution values, or `NULL`
6166: - u_t - The vector will be filled with the time derivative of exact solution values, or `NULL`

6168:   Level: developer

6170:   Note:
6171:   The user must call `PetscDSSetExactSolution()` before using this routine

6173: .seealso: [](ch_dmbase), `DM`, `PetscDSSetExactSolution()`
6174: @*/
6175: PetscErrorCode DMComputeExactSolution(DM dm, PetscReal time, Vec u, Vec u_t)
6176: {
6177:   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
6178:   void   **ectxs;
6179:   Vec      locu, locu_t;
6180:   PetscInt Nf, Nds, s;

6182:   PetscFunctionBegin;
6184:   if (u) {
6186:     PetscCall(DMGetLocalVector(dm, &locu));
6187:     PetscCall(VecSet(locu, 0.));
6188:   }
6189:   if (u_t) {
6191:     PetscCall(DMGetLocalVector(dm, &locu_t));
6192:     PetscCall(VecSet(locu_t, 0.));
6193:   }
6194:   PetscCall(DMGetNumFields(dm, &Nf));
6195:   PetscCall(PetscMalloc2(Nf, &exacts, Nf, &ectxs));
6196:   PetscCall(DMGetNumDS(dm, &Nds));
6197:   for (s = 0; s < Nds; ++s) {
6198:     PetscDS         ds;
6199:     DMLabel         label;
6200:     IS              fieldIS;
6201:     const PetscInt *fields, id = 1;
6202:     PetscInt        dsNf, f;

6204:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL));
6205:     PetscCall(PetscDSGetNumFields(ds, &dsNf));
6206:     PetscCall(ISGetIndices(fieldIS, &fields));
6207:     PetscCall(PetscArrayzero(exacts, Nf));
6208:     PetscCall(PetscArrayzero(ectxs, Nf));
6209:     if (u) {
6210:       for (f = 0; f < dsNf; ++f) PetscCall(PetscDSGetExactSolution(ds, fields[f], &exacts[fields[f]], &ectxs[fields[f]]));
6211:       if (label) PetscCall(DMProjectFunctionLabelLocal(dm, time, label, 1, &id, 0, NULL, exacts, ectxs, INSERT_ALL_VALUES, locu));
6212:       else PetscCall(DMProjectFunctionLocal(dm, time, exacts, ectxs, INSERT_ALL_VALUES, locu));
6213:     }
6214:     if (u_t) {
6215:       PetscCall(PetscArrayzero(exacts, Nf));
6216:       PetscCall(PetscArrayzero(ectxs, Nf));
6217:       for (f = 0; f < dsNf; ++f) PetscCall(PetscDSGetExactSolutionTimeDerivative(ds, fields[f], &exacts[fields[f]], &ectxs[fields[f]]));
6218:       if (label) PetscCall(DMProjectFunctionLabelLocal(dm, time, label, 1, &id, 0, NULL, exacts, ectxs, INSERT_ALL_VALUES, locu_t));
6219:       else PetscCall(DMProjectFunctionLocal(dm, time, exacts, ectxs, INSERT_ALL_VALUES, locu_t));
6220:     }
6221:     PetscCall(ISRestoreIndices(fieldIS, &fields));
6222:   }
6223:   if (u) {
6224:     PetscCall(PetscObjectSetName((PetscObject)u, "Exact Solution"));
6225:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)u, "exact_"));
6226:   }
6227:   if (u_t) {
6228:     PetscCall(PetscObjectSetName((PetscObject)u, "Exact Solution Time Derivative"));
6229:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)u_t, "exact_t_"));
6230:   }
6231:   PetscCall(PetscFree2(exacts, ectxs));
6232:   if (u) {
6233:     PetscCall(DMLocalToGlobalBegin(dm, locu, INSERT_ALL_VALUES, u));
6234:     PetscCall(DMLocalToGlobalEnd(dm, locu, INSERT_ALL_VALUES, u));
6235:     PetscCall(DMRestoreLocalVector(dm, &locu));
6236:   }
6237:   if (u_t) {
6238:     PetscCall(DMLocalToGlobalBegin(dm, locu_t, INSERT_ALL_VALUES, u_t));
6239:     PetscCall(DMLocalToGlobalEnd(dm, locu_t, INSERT_ALL_VALUES, u_t));
6240:     PetscCall(DMRestoreLocalVector(dm, &locu_t));
6241:   }
6242:   PetscFunctionReturn(PETSC_SUCCESS);
6243: }

6245: static PetscErrorCode DMTransferDS_Internal(DM dm, DMLabel label, IS fields, PetscDS ds, PetscDS dsIn)
6246: {
6247:   PetscDS dsNew, dsInNew = NULL;

6249:   PetscFunctionBegin;
6250:   PetscCall(PetscDSCreate(PetscObjectComm((PetscObject)ds), &dsNew));
6251:   PetscCall(PetscDSCopy(ds, dm, dsNew));
6252:   if (dsIn) {
6253:     PetscCall(PetscDSCreate(PetscObjectComm((PetscObject)dsIn), &dsInNew));
6254:     PetscCall(PetscDSCopy(dsIn, dm, dsInNew));
6255:   }
6256:   PetscCall(DMSetRegionDS(dm, label, fields, dsNew, dsInNew));
6257:   PetscCall(PetscDSDestroy(&dsNew));
6258:   PetscCall(PetscDSDestroy(&dsInNew));
6259:   PetscFunctionReturn(PETSC_SUCCESS);
6260: }

6262: /*@
6263:   DMCopyDS - Copy the discrete systems for the `DM` into another `DM`

6265:   Collective

6267:   Input Parameter:
6268: . dm - The `DM`

6270:   Output Parameter:
6271: . newdm - The `DM`

6273:   Level: advanced

6275: .seealso: [](ch_dmbase), `DM`, `DMCopyFields()`, `DMAddField()`, `DMGetDS()`, `DMGetCellDS()`, `DMGetRegionDS()`, `DMSetRegionDS()`
6276: @*/
6277: PetscErrorCode DMCopyDS(DM dm, DM newdm)
6278: {
6279:   PetscInt Nds, s;

6281:   PetscFunctionBegin;
6282:   if (dm == newdm) PetscFunctionReturn(PETSC_SUCCESS);
6283:   PetscCall(DMGetNumDS(dm, &Nds));
6284:   PetscCall(DMClearDS(newdm));
6285:   for (s = 0; s < Nds; ++s) {
6286:     DMLabel  label;
6287:     IS       fields;
6288:     PetscDS  ds, dsIn, newds;
6289:     PetscInt Nbd, bd;

6291:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fields, &ds, &dsIn));
6292:     /* TODO: We need to change all keys from labels in the old DM to labels in the new DM */
6293:     PetscCall(DMTransferDS_Internal(newdm, label, fields, ds, dsIn));
6294:     /* Complete new labels in the new DS */
6295:     PetscCall(DMGetRegionDS(newdm, label, NULL, &newds, NULL));
6296:     PetscCall(PetscDSGetNumBoundary(newds, &Nbd));
6297:     for (bd = 0; bd < Nbd; ++bd) {
6298:       PetscWeakForm wf;
6299:       DMLabel       label;
6300:       PetscInt      field;

6302:       PetscCall(PetscDSGetBoundary(newds, bd, &wf, NULL, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL));
6303:       PetscCall(PetscWeakFormReplaceLabel(wf, label));
6304:     }
6305:   }
6306:   PetscCall(DMCompleteBCLabels_Internal(newdm));
6307:   PetscFunctionReturn(PETSC_SUCCESS);
6308: }

6310: /*@
6311:   DMCopyDisc - Copy the fields and discrete systems for the `DM` into another `DM`

6313:   Collective

6315:   Input Parameter:
6316: . dm - The `DM`

6318:   Output Parameter:
6319: . newdm - The `DM`

6321:   Level: advanced

6323:   Developer Note:
6324:   Really ugly name, nothing in PETSc is called a `Disc` plus it is an ugly abbreviation

6326: .seealso: [](ch_dmbase), `DM`, `DMCopyFields()`, `DMCopyDS()`
6327: @*/
6328: PetscErrorCode DMCopyDisc(DM dm, DM newdm)
6329: {
6330:   PetscFunctionBegin;
6331:   PetscCall(DMCopyFields(dm, newdm));
6332:   PetscCall(DMCopyDS(dm, newdm));
6333:   PetscFunctionReturn(PETSC_SUCCESS);
6334: }

6336: /*@
6337:   DMGetDimension - Return the topological dimension of the `DM`

6339:   Not Collective

6341:   Input Parameter:
6342: . dm - The `DM`

6344:   Output Parameter:
6345: . dim - The topological dimension

6347:   Level: beginner

6349: .seealso: [](ch_dmbase), `DM`, `DMSetDimension()`, `DMCreate()`
6350: @*/
6351: PetscErrorCode DMGetDimension(DM dm, PetscInt *dim)
6352: {
6353:   PetscFunctionBegin;
6355:   PetscAssertPointer(dim, 2);
6356:   *dim = dm->dim;
6357:   PetscFunctionReturn(PETSC_SUCCESS);
6358: }

6360: /*@
6361:   DMSetDimension - Set the topological dimension of the `DM`

6363:   Collective

6365:   Input Parameters:
6366: + dm  - The `DM`
6367: - dim - The topological dimension

6369:   Level: beginner

6371: .seealso: [](ch_dmbase), `DM`, `DMGetDimension()`, `DMCreate()`
6372: @*/
6373: PetscErrorCode DMSetDimension(DM dm, PetscInt dim)
6374: {
6375:   PetscDS  ds;
6376:   PetscInt Nds, n;

6378:   PetscFunctionBegin;
6381:   dm->dim = dim;
6382:   if (dm->dim >= 0) {
6383:     PetscCall(DMGetNumDS(dm, &Nds));
6384:     for (n = 0; n < Nds; ++n) {
6385:       PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds, NULL));
6386:       if (ds->dimEmbed < 0) PetscCall(PetscDSSetCoordinateDimension(ds, dim));
6387:     }
6388:   }
6389:   PetscFunctionReturn(PETSC_SUCCESS);
6390: }

6392: /*@
6393:   DMGetDimPoints - Get the half-open interval for all points of a given dimension

6395:   Collective

6397:   Input Parameters:
6398: + dm  - the `DM`
6399: - dim - the dimension

6401:   Output Parameters:
6402: + pStart - The first point of the given dimension
6403: - pEnd   - The first point following points of the given dimension

6405:   Level: intermediate

6407:   Note:
6408:   The points are vertices in the Hasse diagram encoding the topology. This is explained in
6409:   https://arxiv.org/abs/0908.4427. If no points exist of this dimension in the storage scheme,
6410:   then the interval is empty.

6412: .seealso: [](ch_dmbase), `DM`, `DMPLEX`, `DMPlexGetDepthStratum()`, `DMPlexGetHeightStratum()`
6413: @*/
6414: PetscErrorCode DMGetDimPoints(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
6415: {
6416:   PetscInt d;

6418:   PetscFunctionBegin;
6420:   PetscCall(DMGetDimension(dm, &d));
6421:   PetscCheck((dim >= 0) && (dim <= d), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim);
6422:   PetscUseTypeMethod(dm, getdimpoints, dim, pStart, pEnd);
6423:   PetscFunctionReturn(PETSC_SUCCESS);
6424: }

6426: /*@
6427:   DMGetOutputDM - Retrieve the `DM` associated with the layout for output

6429:   Collective

6431:   Input Parameter:
6432: . dm - The original `DM`

6434:   Output Parameter:
6435: . odm - The `DM` which provides the layout for output

6437:   Level: intermediate

6439:   Note:
6440:   In some situations the vector obtained with `DMCreateGlobalVector()` excludes points for degrees of freedom that are associated with fixed (Dirichelet) boundary
6441:   conditions since the algebraic solver does not solve for those variables. The output `DM` includes these excluded points and its global vector contains the
6442:   locations for those dof so that they can be output to a file or other viewer along with the unconstrained dof.

6444: .seealso: [](ch_dmbase), `DM`, `VecView()`, `DMGetGlobalSection()`, `DMCreateGlobalVector()`, `PetscSectionHasConstraints()`, `DMSetGlobalSection()`
6445: @*/
6446: PetscErrorCode DMGetOutputDM(DM dm, DM *odm)
6447: {
6448:   PetscSection section;
6449:   IS           perm;
6450:   PetscBool    hasConstraints, newDM, gnewDM;

6452:   PetscFunctionBegin;
6454:   PetscAssertPointer(odm, 2);
6455:   PetscCall(DMGetLocalSection(dm, &section));
6456:   PetscCall(PetscSectionHasConstraints(section, &hasConstraints));
6457:   PetscCall(PetscSectionGetPermutation(section, &perm));
6458:   newDM = hasConstraints || perm ? PETSC_TRUE : PETSC_FALSE;
6459:   PetscCall(MPIU_Allreduce(&newDM, &gnewDM, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)dm)));
6460:   if (!gnewDM) {
6461:     *odm = dm;
6462:     PetscFunctionReturn(PETSC_SUCCESS);
6463:   }
6464:   if (!dm->dmBC) {
6465:     PetscSection newSection, gsection;
6466:     PetscSF      sf;
6467:     PetscBool    usePerm = dm->ignorePermOutput ? PETSC_FALSE : PETSC_TRUE;

6469:     PetscCall(DMClone(dm, &dm->dmBC));
6470:     PetscCall(DMCopyDisc(dm, dm->dmBC));
6471:     PetscCall(PetscSectionClone(section, &newSection));
6472:     PetscCall(DMSetLocalSection(dm->dmBC, newSection));
6473:     PetscCall(PetscSectionDestroy(&newSection));
6474:     PetscCall(DMGetPointSF(dm->dmBC, &sf));
6475:     PetscCall(PetscSectionCreateGlobalSection(section, sf, usePerm, PETSC_TRUE, PETSC_FALSE, &gsection));
6476:     PetscCall(DMSetGlobalSection(dm->dmBC, gsection));
6477:     PetscCall(PetscSectionDestroy(&gsection));
6478:   }
6479:   *odm = dm->dmBC;
6480:   PetscFunctionReturn(PETSC_SUCCESS);
6481: }

6483: /*@
6484:   DMGetOutputSequenceNumber - Retrieve the sequence number/value for output

6486:   Input Parameter:
6487: . dm - The original `DM`

6489:   Output Parameters:
6490: + num - The output sequence number
6491: - val - The output sequence value

6493:   Level: intermediate

6495:   Note:
6496:   This is intended for output that should appear in sequence, for instance
6497:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6499:   Developer Note:
6500:   The `DM` serves as a convenient place to store the current iteration value. The iteration is not
6501:   not directly related to the `DM`.

6503: .seealso: [](ch_dmbase), `DM`, `VecView()`
6504: @*/
6505: PetscErrorCode DMGetOutputSequenceNumber(DM dm, PetscInt *num, PetscReal *val)
6506: {
6507:   PetscFunctionBegin;
6509:   if (num) {
6510:     PetscAssertPointer(num, 2);
6511:     *num = dm->outputSequenceNum;
6512:   }
6513:   if (val) {
6514:     PetscAssertPointer(val, 3);
6515:     *val = dm->outputSequenceVal;
6516:   }
6517:   PetscFunctionReturn(PETSC_SUCCESS);
6518: }

6520: /*@
6521:   DMSetOutputSequenceNumber - Set the sequence number/value for output

6523:   Input Parameters:
6524: + dm  - The original `DM`
6525: . num - The output sequence number
6526: - val - The output sequence value

6528:   Level: intermediate

6530:   Note:
6531:   This is intended for output that should appear in sequence, for instance
6532:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6534: .seealso: [](ch_dmbase), `DM`, `VecView()`
6535: @*/
6536: PetscErrorCode DMSetOutputSequenceNumber(DM dm, PetscInt num, PetscReal val)
6537: {
6538:   PetscFunctionBegin;
6540:   dm->outputSequenceNum = num;
6541:   dm->outputSequenceVal = val;
6542:   PetscFunctionReturn(PETSC_SUCCESS);
6543: }

6545: /*@C
6546:   DMOutputSequenceLoad - Retrieve the sequence value from a `PetscViewer`

6548:   Input Parameters:
6549: + dm     - The original `DM`
6550: . viewer - The viewer to get it from
6551: . name   - The sequence name
6552: - num    - The output sequence number

6554:   Output Parameter:
6555: . val - The output sequence value

6557:   Level: intermediate

6559:   Note:
6560:   This is intended for output that should appear in sequence, for instance
6561:   a set of timesteps in an `PETSCVIEWERHDF5` file, or a set of realizations of a stochastic system.

6563:   Developer Note:
6564:   It is unclear at the user API level why a `DM` is needed as input

6566: .seealso: [](ch_dmbase), `DM`, `DMGetOutputSequenceNumber()`, `DMSetOutputSequenceNumber()`, `VecView()`
6567: @*/
6568: PetscErrorCode DMOutputSequenceLoad(DM dm, PetscViewer viewer, const char *name, PetscInt num, PetscReal *val)
6569: {
6570:   PetscBool ishdf5;

6572:   PetscFunctionBegin;
6575:   PetscAssertPointer(val, 5);
6576:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
6577:   if (ishdf5) {
6578: #if defined(PETSC_HAVE_HDF5)
6579:     PetscScalar value;

6581:     PetscCall(DMSequenceLoad_HDF5_Internal(dm, name, num, &value, viewer));
6582:     *val = PetscRealPart(value);
6583: #endif
6584:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerHDF5Open()");
6585:   PetscFunctionReturn(PETSC_SUCCESS);
6586: }

6588: /*@
6589:   DMGetUseNatural - Get the flag for creating a mapping to the natural order when a `DM` is (re)distributed in parallel

6591:   Not Collective

6593:   Input Parameter:
6594: . dm - The `DM`

6596:   Output Parameter:
6597: . useNatural - `PETSC_TRUE` to build the mapping to a natural order during distribution

6599:   Level: beginner

6601: .seealso: [](ch_dmbase), `DM`, `DMSetUseNatural()`, `DMCreate()`
6602: @*/
6603: PetscErrorCode DMGetUseNatural(DM dm, PetscBool *useNatural)
6604: {
6605:   PetscFunctionBegin;
6607:   PetscAssertPointer(useNatural, 2);
6608:   *useNatural = dm->useNatural;
6609:   PetscFunctionReturn(PETSC_SUCCESS);
6610: }

6612: /*@
6613:   DMSetUseNatural - Set the flag for creating a mapping to the natural order when a `DM` is (re)distributed in parallel

6615:   Collective

6617:   Input Parameters:
6618: + dm         - The `DM`
6619: - useNatural - `PETSC_TRUE` to build the mapping to a natural order during distribution

6621:   Level: beginner

6623:   Note:
6624:   This also causes the map to be build after `DMCreateSubDM()` and `DMCreateSuperDM()`

6626: .seealso: [](ch_dmbase), `DM`, `DMGetUseNatural()`, `DMCreate()`, `DMPlexDistribute()`, `DMCreateSubDM()`, `DMCreateSuperDM()`
6627: @*/
6628: PetscErrorCode DMSetUseNatural(DM dm, PetscBool useNatural)
6629: {
6630:   PetscFunctionBegin;
6633:   dm->useNatural = useNatural;
6634:   PetscFunctionReturn(PETSC_SUCCESS);
6635: }

6637: /*@C
6638:   DMCreateLabel - Create a label of the given name if it does not already exist in the `DM`

6640:   Not Collective

6642:   Input Parameters:
6643: + dm   - The `DM` object
6644: - name - The label name

6646:   Level: intermediate

6648: .seealso: [](ch_dmbase), `DM`, `DMLabelCreate()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6649: @*/
6650: PetscErrorCode DMCreateLabel(DM dm, const char name[])
6651: {
6652:   PetscBool flg;
6653:   DMLabel   label;

6655:   PetscFunctionBegin;
6657:   PetscAssertPointer(name, 2);
6658:   PetscCall(DMHasLabel(dm, name, &flg));
6659:   if (!flg) {
6660:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, &label));
6661:     PetscCall(DMAddLabel(dm, label));
6662:     PetscCall(DMLabelDestroy(&label));
6663:   }
6664:   PetscFunctionReturn(PETSC_SUCCESS);
6665: }

6667: /*@C
6668:   DMCreateLabelAtIndex - Create a label of the given name at the given index. If it already exists in the `DM`, move it to this index.

6670:   Not Collective

6672:   Input Parameters:
6673: + dm   - The `DM` object
6674: . l    - The index for the label
6675: - name - The label name

6677:   Level: intermediate

6679: .seealso: [](ch_dmbase), `DM`, `DMCreateLabel()`, `DMLabelCreate()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6680: @*/
6681: PetscErrorCode DMCreateLabelAtIndex(DM dm, PetscInt l, const char name[])
6682: {
6683:   DMLabelLink orig, prev = NULL;
6684:   DMLabel     label;
6685:   PetscInt    Nl, m;
6686:   PetscBool   flg, match;
6687:   const char *lname;

6689:   PetscFunctionBegin;
6691:   PetscAssertPointer(name, 3);
6692:   PetscCall(DMHasLabel(dm, name, &flg));
6693:   if (!flg) {
6694:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, &label));
6695:     PetscCall(DMAddLabel(dm, label));
6696:     PetscCall(DMLabelDestroy(&label));
6697:   }
6698:   PetscCall(DMGetNumLabels(dm, &Nl));
6699:   PetscCheck(l < Nl, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label index %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", l, Nl);
6700:   for (m = 0, orig = dm->labels; m < Nl; ++m, prev = orig, orig = orig->next) {
6701:     PetscCall(PetscObjectGetName((PetscObject)orig->label, &lname));
6702:     PetscCall(PetscStrcmp(name, lname, &match));
6703:     if (match) break;
6704:   }
6705:   if (m == l) PetscFunctionReturn(PETSC_SUCCESS);
6706:   if (!m) dm->labels = orig->next;
6707:   else prev->next = orig->next;
6708:   if (!l) {
6709:     orig->next = dm->labels;
6710:     dm->labels = orig;
6711:   } else {
6712:     for (m = 0, prev = dm->labels; m < l - 1; ++m, prev = prev->next);
6713:     orig->next = prev->next;
6714:     prev->next = orig;
6715:   }
6716:   PetscFunctionReturn(PETSC_SUCCESS);
6717: }

6719: /*@C
6720:   DMGetLabelValue - Get the value in a `DMLabel` for the given point, with -1 as the default

6722:   Not Collective

6724:   Input Parameters:
6725: + dm    - The `DM` object
6726: . name  - The label name
6727: - point - The mesh point

6729:   Output Parameter:
6730: . value - The label value for this point, or -1 if the point is not in the label

6732:   Level: beginner

6734: .seealso: [](ch_dmbase), `DM`, `DMLabelGetValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6735: @*/
6736: PetscErrorCode DMGetLabelValue(DM dm, const char name[], PetscInt point, PetscInt *value)
6737: {
6738:   DMLabel label;

6740:   PetscFunctionBegin;
6742:   PetscAssertPointer(name, 2);
6743:   PetscCall(DMGetLabel(dm, name, &label));
6744:   PetscCheck(label, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No label named %s was found", name);
6745:   PetscCall(DMLabelGetValue(label, point, value));
6746:   PetscFunctionReturn(PETSC_SUCCESS);
6747: }

6749: /*@C
6750:   DMSetLabelValue - Add a point to a `DMLabel` with given value

6752:   Not Collective

6754:   Input Parameters:
6755: + dm    - The `DM` object
6756: . name  - The label name
6757: . point - The mesh point
6758: - value - The label value for this point

6760:   Output Parameter:

6762:   Level: beginner

6764: .seealso: [](ch_dmbase), `DM`, `DMLabelSetValue()`, `DMGetStratumIS()`, `DMClearLabelValue()`
6765: @*/
6766: PetscErrorCode DMSetLabelValue(DM dm, const char name[], PetscInt point, PetscInt value)
6767: {
6768:   DMLabel label;

6770:   PetscFunctionBegin;
6772:   PetscAssertPointer(name, 2);
6773:   PetscCall(DMGetLabel(dm, name, &label));
6774:   if (!label) {
6775:     PetscCall(DMCreateLabel(dm, name));
6776:     PetscCall(DMGetLabel(dm, name, &label));
6777:   }
6778:   PetscCall(DMLabelSetValue(label, point, value));
6779:   PetscFunctionReturn(PETSC_SUCCESS);
6780: }

6782: /*@C
6783:   DMClearLabelValue - Remove a point from a `DMLabel` with given value

6785:   Not Collective

6787:   Input Parameters:
6788: + dm    - The `DM` object
6789: . name  - The label name
6790: . point - The mesh point
6791: - value - The label value for this point

6793:   Level: beginner

6795: .seealso: [](ch_dmbase), `DM`, `DMLabelClearValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
6796: @*/
6797: PetscErrorCode DMClearLabelValue(DM dm, const char name[], PetscInt point, PetscInt value)
6798: {
6799:   DMLabel label;

6801:   PetscFunctionBegin;
6803:   PetscAssertPointer(name, 2);
6804:   PetscCall(DMGetLabel(dm, name, &label));
6805:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6806:   PetscCall(DMLabelClearValue(label, point, value));
6807:   PetscFunctionReturn(PETSC_SUCCESS);
6808: }

6810: /*@C
6811:   DMGetLabelSize - Get the value of `DMLabelGetNumValues()` of a `DMLabel` in the `DM`

6813:   Not Collective

6815:   Input Parameters:
6816: + dm   - The `DM` object
6817: - name - The label name

6819:   Output Parameter:
6820: . size - The number of different integer ids, or 0 if the label does not exist

6822:   Level: beginner

6824:   Developer Note:
6825:   This should be renamed to something like `DMGetLabelNumValues()` or removed.

6827: .seealso: [](ch_dmbase), `DM`, `DMLabelGetNumValues()`, `DMSetLabelValue()`, `DMGetLabel()`
6828: @*/
6829: PetscErrorCode DMGetLabelSize(DM dm, const char name[], PetscInt *size)
6830: {
6831:   DMLabel label;

6833:   PetscFunctionBegin;
6835:   PetscAssertPointer(name, 2);
6836:   PetscAssertPointer(size, 3);
6837:   PetscCall(DMGetLabel(dm, name, &label));
6838:   *size = 0;
6839:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6840:   PetscCall(DMLabelGetNumValues(label, size));
6841:   PetscFunctionReturn(PETSC_SUCCESS);
6842: }

6844: /*@C
6845:   DMGetLabelIdIS - Get the `DMLabelGetValueIS()` from a `DMLabel` in the `DM`

6847:   Not Collective

6849:   Input Parameters:
6850: + dm   - The `DM` object
6851: - name - The label name

6853:   Output Parameter:
6854: . ids - The integer ids, or `NULL` if the label does not exist

6856:   Level: beginner

6858: .seealso: [](ch_dmbase), `DM`, `DMLabelGetValueIS()`, `DMGetLabelSize()`
6859: @*/
6860: PetscErrorCode DMGetLabelIdIS(DM dm, const char name[], IS *ids)
6861: {
6862:   DMLabel label;

6864:   PetscFunctionBegin;
6866:   PetscAssertPointer(name, 2);
6867:   PetscAssertPointer(ids, 3);
6868:   PetscCall(DMGetLabel(dm, name, &label));
6869:   *ids = NULL;
6870:   if (label) {
6871:     PetscCall(DMLabelGetValueIS(label, ids));
6872:   } else {
6873:     /* returning an empty IS */
6874:     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_USE_POINTER, ids));
6875:   }
6876:   PetscFunctionReturn(PETSC_SUCCESS);
6877: }

6879: /*@C
6880:   DMGetStratumSize - Get the number of points in a label stratum

6882:   Not Collective

6884:   Input Parameters:
6885: + dm    - The `DM` object
6886: . name  - The label name
6887: - value - The stratum value

6889:   Output Parameter:
6890: . size - The number of points, also called the stratum size

6892:   Level: beginner

6894: .seealso: [](ch_dmbase), `DM`, `DMLabelGetStratumSize()`, `DMGetLabelSize()`, `DMGetLabelIds()`
6895: @*/
6896: PetscErrorCode DMGetStratumSize(DM dm, const char name[], PetscInt value, PetscInt *size)
6897: {
6898:   DMLabel label;

6900:   PetscFunctionBegin;
6902:   PetscAssertPointer(name, 2);
6903:   PetscAssertPointer(size, 4);
6904:   PetscCall(DMGetLabel(dm, name, &label));
6905:   *size = 0;
6906:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6907:   PetscCall(DMLabelGetStratumSize(label, value, size));
6908:   PetscFunctionReturn(PETSC_SUCCESS);
6909: }

6911: /*@C
6912:   DMGetStratumIS - Get the points in a label stratum

6914:   Not Collective

6916:   Input Parameters:
6917: + dm    - The `DM` object
6918: . name  - The label name
6919: - value - The stratum value

6921:   Output Parameter:
6922: . points - The stratum points, or `NULL` if the label does not exist or does not have that value

6924:   Level: beginner

6926: .seealso: [](ch_dmbase), `DM`, `DMLabelGetStratumIS()`, `DMGetStratumSize()`
6927: @*/
6928: PetscErrorCode DMGetStratumIS(DM dm, const char name[], PetscInt value, IS *points)
6929: {
6930:   DMLabel label;

6932:   PetscFunctionBegin;
6934:   PetscAssertPointer(name, 2);
6935:   PetscAssertPointer(points, 4);
6936:   PetscCall(DMGetLabel(dm, name, &label));
6937:   *points = NULL;
6938:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6939:   PetscCall(DMLabelGetStratumIS(label, value, points));
6940:   PetscFunctionReturn(PETSC_SUCCESS);
6941: }

6943: /*@C
6944:   DMSetStratumIS - Set the points in a label stratum

6946:   Not Collective

6948:   Input Parameters:
6949: + dm     - The `DM` object
6950: . name   - The label name
6951: . value  - The stratum value
6952: - points - The stratum points

6954:   Level: beginner

6956: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMClearLabelStratum()`, `DMLabelClearStratum()`, `DMLabelSetStratumIS()`, `DMGetStratumSize()`
6957: @*/
6958: PetscErrorCode DMSetStratumIS(DM dm, const char name[], PetscInt value, IS points)
6959: {
6960:   DMLabel label;

6962:   PetscFunctionBegin;
6964:   PetscAssertPointer(name, 2);
6966:   PetscCall(DMGetLabel(dm, name, &label));
6967:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6968:   PetscCall(DMLabelSetStratumIS(label, value, points));
6969:   PetscFunctionReturn(PETSC_SUCCESS);
6970: }

6972: /*@C
6973:   DMClearLabelStratum - Remove all points from a stratum from a `DMLabel`

6975:   Not Collective

6977:   Input Parameters:
6978: + dm    - The `DM` object
6979: . name  - The label name
6980: - value - The label value for this point

6982:   Output Parameter:

6984:   Level: beginner

6986: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMLabelClearStratum()`, `DMSetLabelValue()`, `DMGetStratumIS()`, `DMClearLabelValue()`
6987: @*/
6988: PetscErrorCode DMClearLabelStratum(DM dm, const char name[], PetscInt value)
6989: {
6990:   DMLabel label;

6992:   PetscFunctionBegin;
6994:   PetscAssertPointer(name, 2);
6995:   PetscCall(DMGetLabel(dm, name, &label));
6996:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6997:   PetscCall(DMLabelClearStratum(label, value));
6998:   PetscFunctionReturn(PETSC_SUCCESS);
6999: }

7001: /*@
7002:   DMGetNumLabels - Return the number of labels defined by on the `DM`

7004:   Not Collective

7006:   Input Parameter:
7007: . dm - The `DM` object

7009:   Output Parameter:
7010: . numLabels - the number of Labels

7012:   Level: intermediate

7014: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetLabelByNum()`, `DMGetLabelName()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7015: @*/
7016: PetscErrorCode DMGetNumLabels(DM dm, PetscInt *numLabels)
7017: {
7018:   DMLabelLink next = dm->labels;
7019:   PetscInt    n    = 0;

7021:   PetscFunctionBegin;
7023:   PetscAssertPointer(numLabels, 2);
7024:   while (next) {
7025:     ++n;
7026:     next = next->next;
7027:   }
7028:   *numLabels = n;
7029:   PetscFunctionReturn(PETSC_SUCCESS);
7030: }

7032: /*@C
7033:   DMGetLabelName - Return the name of nth label

7035:   Not Collective

7037:   Input Parameters:
7038: + dm - The `DM` object
7039: - n  - the label number

7041:   Output Parameter:
7042: . name - the label name

7044:   Level: intermediate

7046:   Developer Note:
7047:   Some of the functions that appropriate on labels using their number have the suffix ByNum, others do not.

7049: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetLabelByNum()`, `DMGetLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7050: @*/
7051: PetscErrorCode DMGetLabelName(DM dm, PetscInt n, const char **name)
7052: {
7053:   DMLabelLink next = dm->labels;
7054:   PetscInt    l    = 0;

7056:   PetscFunctionBegin;
7058:   PetscAssertPointer(name, 3);
7059:   while (next) {
7060:     if (l == n) {
7061:       PetscCall(PetscObjectGetName((PetscObject)next->label, name));
7062:       PetscFunctionReturn(PETSC_SUCCESS);
7063:     }
7064:     ++l;
7065:     next = next->next;
7066:   }
7067:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %" PetscInt_FMT " does not exist in this DM", n);
7068: }

7070: /*@C
7071:   DMHasLabel - Determine whether the `DM` has a label of a given name

7073:   Not Collective

7075:   Input Parameters:
7076: + dm   - The `DM` object
7077: - name - The label name

7079:   Output Parameter:
7080: . hasLabel - `PETSC_TRUE` if the label is present

7082:   Level: intermediate

7084: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetLabel()`, `DMGetLabelByNum()`, `DMCreateLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7085: @*/
7086: PetscErrorCode DMHasLabel(DM dm, const char name[], PetscBool *hasLabel)
7087: {
7088:   DMLabelLink next = dm->labels;
7089:   const char *lname;

7091:   PetscFunctionBegin;
7093:   PetscAssertPointer(name, 2);
7094:   PetscAssertPointer(hasLabel, 3);
7095:   *hasLabel = PETSC_FALSE;
7096:   while (next) {
7097:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7098:     PetscCall(PetscStrcmp(name, lname, hasLabel));
7099:     if (*hasLabel) break;
7100:     next = next->next;
7101:   }
7102:   PetscFunctionReturn(PETSC_SUCCESS);
7103: }

7105: // PetscClangLinter pragma ignore: -fdoc-section-header-unknown
7106: /*@C
7107:   DMGetLabel - Return the label of a given name, or `NULL`, from a `DM`

7109:   Not Collective

7111:   Input Parameters:
7112: + dm   - The `DM` object
7113: - name - The label name

7115:   Output Parameter:
7116: . label - The `DMLabel`, or `NULL` if the label is absent

7118:   Default labels in a `DMPLEX`:
7119: + "depth"       - Holds the depth (co-dimension) of each mesh point
7120: . "celltype"    - Holds the topological type of each cell
7121: . "ghost"       - If the DM is distributed with overlap, this marks the cells and faces in the overlap
7122: . "Cell Sets"   - Mirrors the cell sets defined by GMsh and ExodusII
7123: . "Face Sets"   - Mirrors the face sets defined by GMsh and ExodusII
7124: - "Vertex Sets" - Mirrors the vertex sets defined by GMsh

7126:   Level: intermediate

7128: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMHasLabel()`, `DMGetLabelByNum()`, `DMAddLabel()`, `DMCreateLabel()`, `DMPlexGetDepthLabel()`, `DMPlexGetCellType()`
7129: @*/
7130: PetscErrorCode DMGetLabel(DM dm, const char name[], DMLabel *label)
7131: {
7132:   DMLabelLink next = dm->labels;
7133:   PetscBool   hasLabel;
7134:   const char *lname;

7136:   PetscFunctionBegin;
7138:   PetscAssertPointer(name, 2);
7139:   PetscAssertPointer(label, 3);
7140:   *label = NULL;
7141:   while (next) {
7142:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7143:     PetscCall(PetscStrcmp(name, lname, &hasLabel));
7144:     if (hasLabel) {
7145:       *label = next->label;
7146:       break;
7147:     }
7148:     next = next->next;
7149:   }
7150:   PetscFunctionReturn(PETSC_SUCCESS);
7151: }

7153: /*@C
7154:   DMGetLabelByNum - Return the nth label on a `DM`

7156:   Not Collective

7158:   Input Parameters:
7159: + dm - The `DM` object
7160: - n  - the label number

7162:   Output Parameter:
7163: . label - the label

7165:   Level: intermediate

7167: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMAddLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7168: @*/
7169: PetscErrorCode DMGetLabelByNum(DM dm, PetscInt n, DMLabel *label)
7170: {
7171:   DMLabelLink next = dm->labels;
7172:   PetscInt    l    = 0;

7174:   PetscFunctionBegin;
7176:   PetscAssertPointer(label, 3);
7177:   while (next) {
7178:     if (l == n) {
7179:       *label = next->label;
7180:       PetscFunctionReturn(PETSC_SUCCESS);
7181:     }
7182:     ++l;
7183:     next = next->next;
7184:   }
7185:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %" PetscInt_FMT " does not exist in this DM", n);
7186: }

7188: /*@C
7189:   DMAddLabel - Add the label to this `DM`

7191:   Not Collective

7193:   Input Parameters:
7194: + dm    - The `DM` object
7195: - label - The `DMLabel`

7197:   Level: developer

7199: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7200: @*/
7201: PetscErrorCode DMAddLabel(DM dm, DMLabel label)
7202: {
7203:   DMLabelLink l, *p, tmpLabel;
7204:   PetscBool   hasLabel;
7205:   const char *lname;
7206:   PetscBool   flg;

7208:   PetscFunctionBegin;
7210:   PetscCall(PetscObjectGetName((PetscObject)label, &lname));
7211:   PetscCall(DMHasLabel(dm, lname, &hasLabel));
7212:   PetscCheck(!hasLabel, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %s already exists in this DM", lname);
7213:   PetscCall(PetscCalloc1(1, &tmpLabel));
7214:   tmpLabel->label  = label;
7215:   tmpLabel->output = PETSC_TRUE;
7216:   for (p = &dm->labels; (l = *p); p = &l->next) { }
7217:   *p = tmpLabel;
7218:   PetscCall(PetscObjectReference((PetscObject)label));
7219:   PetscCall(PetscStrcmp(lname, "depth", &flg));
7220:   if (flg) dm->depthLabel = label;
7221:   PetscCall(PetscStrcmp(lname, "celltype", &flg));
7222:   if (flg) dm->celltypeLabel = label;
7223:   PetscFunctionReturn(PETSC_SUCCESS);
7224: }

7226: // PetscClangLinter pragma ignore: -fdoc-section-header-unknown
7227: /*@C
7228:   DMSetLabel - Replaces the label of a given name, or ignores it if the name is not present

7230:   Not Collective

7232:   Input Parameters:
7233: + dm    - The `DM` object
7234: - label - The `DMLabel`, having the same name, to substitute

7236:   Default labels in a `DMPLEX`:
7237: + "depth"       - Holds the depth (co-dimension) of each mesh point
7238: . "celltype"    - Holds the topological type of each cell
7239: . "ghost"       - If the DM is distributed with overlap, this marks the cells and faces in the overlap
7240: . "Cell Sets"   - Mirrors the cell sets defined by GMsh and ExodusII
7241: . "Face Sets"   - Mirrors the face sets defined by GMsh and ExodusII
7242: - "Vertex Sets" - Mirrors the vertex sets defined by GMsh

7244:   Level: intermediate

7246: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMPlexGetDepthLabel()`, `DMPlexGetCellType()`
7247: @*/
7248: PetscErrorCode DMSetLabel(DM dm, DMLabel label)
7249: {
7250:   DMLabelLink next = dm->labels;
7251:   PetscBool   hasLabel, flg;
7252:   const char *name, *lname;

7254:   PetscFunctionBegin;
7257:   PetscCall(PetscObjectGetName((PetscObject)label, &name));
7258:   while (next) {
7259:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7260:     PetscCall(PetscStrcmp(name, lname, &hasLabel));
7261:     if (hasLabel) {
7262:       PetscCall(PetscObjectReference((PetscObject)label));
7263:       PetscCall(PetscStrcmp(lname, "depth", &flg));
7264:       if (flg) dm->depthLabel = label;
7265:       PetscCall(PetscStrcmp(lname, "celltype", &flg));
7266:       if (flg) dm->celltypeLabel = label;
7267:       PetscCall(DMLabelDestroy(&next->label));
7268:       next->label = label;
7269:       break;
7270:     }
7271:     next = next->next;
7272:   }
7273:   PetscFunctionReturn(PETSC_SUCCESS);
7274: }

7276: /*@C
7277:   DMRemoveLabel - Remove the label given by name from this `DM`

7279:   Not Collective

7281:   Input Parameters:
7282: + dm   - The `DM` object
7283: - name - The label name

7285:   Output Parameter:
7286: . label - The `DMLabel`, or `NULL` if the label is absent. Pass in `NULL` to call `DMLabelDestroy()` on the label, otherwise the
7287:           caller is responsible for calling `DMLabelDestroy()`.

7289:   Level: developer

7291: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMLabelDestroy()`, `DMRemoveLabelBySelf()`
7292: @*/
7293: PetscErrorCode DMRemoveLabel(DM dm, const char name[], DMLabel *label)
7294: {
7295:   DMLabelLink link, *pnext;
7296:   PetscBool   hasLabel;
7297:   const char *lname;

7299:   PetscFunctionBegin;
7301:   PetscAssertPointer(name, 2);
7302:   if (label) {
7303:     PetscAssertPointer(label, 3);
7304:     *label = NULL;
7305:   }
7306:   for (pnext = &dm->labels; (link = *pnext); pnext = &link->next) {
7307:     PetscCall(PetscObjectGetName((PetscObject)link->label, &lname));
7308:     PetscCall(PetscStrcmp(name, lname, &hasLabel));
7309:     if (hasLabel) {
7310:       *pnext = link->next; /* Remove from list */
7311:       PetscCall(PetscStrcmp(name, "depth", &hasLabel));
7312:       if (hasLabel) dm->depthLabel = NULL;
7313:       PetscCall(PetscStrcmp(name, "celltype", &hasLabel));
7314:       if (hasLabel) dm->celltypeLabel = NULL;
7315:       if (label) *label = link->label;
7316:       else PetscCall(DMLabelDestroy(&link->label));
7317:       PetscCall(PetscFree(link));
7318:       break;
7319:     }
7320:   }
7321:   PetscFunctionReturn(PETSC_SUCCESS);
7322: }

7324: /*@
7325:   DMRemoveLabelBySelf - Remove the label from this `DM`

7327:   Not Collective

7329:   Input Parameters:
7330: + dm           - The `DM` object
7331: . label        - The `DMLabel` to be removed from the `DM`
7332: - failNotFound - Should it fail if the label is not found in the `DM`?

7334:   Level: developer

7336:   Note:
7337:   Only exactly the same instance is removed if found, name match is ignored.
7338:   If the `DM` has an exclusive reference to the label, the label gets destroyed and
7339:   *label nullified.

7341: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabel()` `DMGetLabelValue()`, `DMSetLabelValue()`, `DMLabelDestroy()`, `DMRemoveLabel()`
7342: @*/
7343: PetscErrorCode DMRemoveLabelBySelf(DM dm, DMLabel *label, PetscBool failNotFound)
7344: {
7345:   DMLabelLink link, *pnext;
7346:   PetscBool   hasLabel = PETSC_FALSE;

7348:   PetscFunctionBegin;
7350:   PetscAssertPointer(label, 2);
7351:   if (!*label && !failNotFound) PetscFunctionReturn(PETSC_SUCCESS);
7354:   for (pnext = &dm->labels; (link = *pnext); pnext = &link->next) {
7355:     if (*label == link->label) {
7356:       hasLabel = PETSC_TRUE;
7357:       *pnext   = link->next; /* Remove from list */
7358:       if (*label == dm->depthLabel) dm->depthLabel = NULL;
7359:       if (*label == dm->celltypeLabel) dm->celltypeLabel = NULL;
7360:       if (((PetscObject)link->label)->refct < 2) *label = NULL; /* nullify if exclusive reference */
7361:       PetscCall(DMLabelDestroy(&link->label));
7362:       PetscCall(PetscFree(link));
7363:       break;
7364:     }
7365:   }
7366:   PetscCheck(hasLabel || !failNotFound, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Given label not found in DM");
7367:   PetscFunctionReturn(PETSC_SUCCESS);
7368: }

7370: /*@C
7371:   DMGetLabelOutput - Get the output flag for a given label

7373:   Not Collective

7375:   Input Parameters:
7376: + dm   - The `DM` object
7377: - name - The label name

7379:   Output Parameter:
7380: . output - The flag for output

7382:   Level: developer

7384: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMSetLabelOutput()`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7385: @*/
7386: PetscErrorCode DMGetLabelOutput(DM dm, const char name[], PetscBool *output)
7387: {
7388:   DMLabelLink next = dm->labels;
7389:   const char *lname;

7391:   PetscFunctionBegin;
7393:   PetscAssertPointer(name, 2);
7394:   PetscAssertPointer(output, 3);
7395:   while (next) {
7396:     PetscBool flg;

7398:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7399:     PetscCall(PetscStrcmp(name, lname, &flg));
7400:     if (flg) {
7401:       *output = next->output;
7402:       PetscFunctionReturn(PETSC_SUCCESS);
7403:     }
7404:     next = next->next;
7405:   }
7406:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No label named %s was present in this dm", name);
7407: }

7409: /*@C
7410:   DMSetLabelOutput - Set if a given label should be saved to a `PetscViewer` in calls to `DMView()`

7412:   Not Collective

7414:   Input Parameters:
7415: + dm     - The `DM` object
7416: . name   - The label name
7417: - output - `PETSC_TRUE` to save the label to the viewer

7419:   Level: developer

7421: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMGetOutputFlag()`, `DMGetLabelOutput()`, `DMCreateLabel()`, `DMHasLabel()`, `DMGetLabelValue()`, `DMSetLabelValue()`, `DMGetStratumIS()`
7422: @*/
7423: PetscErrorCode DMSetLabelOutput(DM dm, const char name[], PetscBool output)
7424: {
7425:   DMLabelLink next = dm->labels;
7426:   const char *lname;

7428:   PetscFunctionBegin;
7430:   PetscAssertPointer(name, 2);
7431:   while (next) {
7432:     PetscBool flg;

7434:     PetscCall(PetscObjectGetName((PetscObject)next->label, &lname));
7435:     PetscCall(PetscStrcmp(name, lname, &flg));
7436:     if (flg) {
7437:       next->output = output;
7438:       PetscFunctionReturn(PETSC_SUCCESS);
7439:     }
7440:     next = next->next;
7441:   }
7442:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No label named %s was present in this dm", name);
7443: }

7445: /*@
7446:   DMCopyLabels - Copy labels from one `DM` mesh to another `DM` with a superset of the points

7448:   Collective

7450:   Input Parameters:
7451: + dmA   - The `DM` object with initial labels
7452: . dmB   - The `DM` object to which labels are copied
7453: . mode  - Copy labels by pointers (`PETSC_OWN_POINTER`) or duplicate them (`PETSC_COPY_VALUES`)
7454: . all   - Copy all labels including "depth", "dim", and "celltype" (`PETSC_TRUE`) which are otherwise ignored (`PETSC_FALSE`)
7455: - emode - How to behave when a `DMLabel` in the source and destination `DM`s with the same name is encountered (see `DMCopyLabelsMode`)

7457:   Level: intermediate

7459:   Note:
7460:   This is typically used when interpolating or otherwise adding to a mesh, or testing.

7462: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMAddLabel()`, `DMCopyLabelsMode`
7463: @*/
7464: PetscErrorCode DMCopyLabels(DM dmA, DM dmB, PetscCopyMode mode, PetscBool all, DMCopyLabelsMode emode)
7465: {
7466:   DMLabel     label, labelNew, labelOld;
7467:   const char *name;
7468:   PetscBool   flg;
7469:   DMLabelLink link;

7471:   PetscFunctionBegin;
7476:   PetscCheck(mode != PETSC_USE_POINTER, PetscObjectComm((PetscObject)dmA), PETSC_ERR_SUP, "PETSC_USE_POINTER not supported for objects");
7477:   if (dmA == dmB) PetscFunctionReturn(PETSC_SUCCESS);
7478:   for (link = dmA->labels; link; link = link->next) {
7479:     label = link->label;
7480:     PetscCall(PetscObjectGetName((PetscObject)label, &name));
7481:     if (!all) {
7482:       PetscCall(PetscStrcmp(name, "depth", &flg));
7483:       if (flg) continue;
7484:       PetscCall(PetscStrcmp(name, "dim", &flg));
7485:       if (flg) continue;
7486:       PetscCall(PetscStrcmp(name, "celltype", &flg));
7487:       if (flg) continue;
7488:     }
7489:     PetscCall(DMGetLabel(dmB, name, &labelOld));
7490:     if (labelOld) {
7491:       switch (emode) {
7492:       case DM_COPY_LABELS_KEEP:
7493:         continue;
7494:       case DM_COPY_LABELS_REPLACE:
7495:         PetscCall(DMRemoveLabelBySelf(dmB, &labelOld, PETSC_TRUE));
7496:         break;
7497:       case DM_COPY_LABELS_FAIL:
7498:         SETERRQ(PetscObjectComm((PetscObject)dmA), PETSC_ERR_ARG_OUTOFRANGE, "Label %s already exists in destination DM", name);
7499:       default:
7500:         SETERRQ(PetscObjectComm((PetscObject)dmA), PETSC_ERR_ARG_OUTOFRANGE, "Unhandled DMCopyLabelsMode %d", (int)emode);
7501:       }
7502:     }
7503:     if (mode == PETSC_COPY_VALUES) {
7504:       PetscCall(DMLabelDuplicate(label, &labelNew));
7505:     } else {
7506:       labelNew = label;
7507:     }
7508:     PetscCall(DMAddLabel(dmB, labelNew));
7509:     if (mode == PETSC_COPY_VALUES) PetscCall(DMLabelDestroy(&labelNew));
7510:   }
7511:   PetscFunctionReturn(PETSC_SUCCESS);
7512: }

7514: /*@C
7515:   DMCompareLabels - Compare labels between two `DM` objects

7517:   Collective; No Fortran Support

7519:   Input Parameters:
7520: + dm0 - First `DM` object
7521: - dm1 - Second `DM` object

7523:   Output Parameters:
7524: + equal   - (Optional) Flag whether labels of dm0 and dm1 are the same
7525: - message - (Optional) Message describing the difference, or `NULL` if there is no difference

7527:   Level: intermediate

7529:   Notes:
7530:   The output flag equal will be the same on all processes.

7532:   If equal is passed as `NULL` and difference is found, an error is thrown on all processes.

7534:   Make sure to pass equal is `NULL` on all processes or none of them.

7536:   The output message is set independently on each rank.

7538:   message must be freed with `PetscFree()`

7540:   If message is passed as `NULL` and a difference is found, the difference description is printed to stderr in synchronized manner.

7542:   Make sure to pass message as `NULL` on all processes or no processes.

7544:   Labels are matched by name. If the number of labels and their names are equal,
7545:   `DMLabelCompare()` is used to compare each pair of labels with the same name.

7547: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMAddLabel()`, `DMCopyLabelsMode`, `DMLabelCompare()`
7548: @*/
7549: PetscErrorCode DMCompareLabels(DM dm0, DM dm1, PetscBool *equal, char **message)
7550: {
7551:   PetscInt    n, i;
7552:   char        msg[PETSC_MAX_PATH_LEN] = "";
7553:   PetscBool   eq;
7554:   MPI_Comm    comm;
7555:   PetscMPIInt rank;

7557:   PetscFunctionBegin;
7560:   PetscCheckSameComm(dm0, 1, dm1, 2);
7561:   if (equal) PetscAssertPointer(equal, 3);
7562:   if (message) PetscAssertPointer(message, 4);
7563:   PetscCall(PetscObjectGetComm((PetscObject)dm0, &comm));
7564:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7565:   {
7566:     PetscInt n1;

7568:     PetscCall(DMGetNumLabels(dm0, &n));
7569:     PetscCall(DMGetNumLabels(dm1, &n1));
7570:     eq = (PetscBool)(n == n1);
7571:     if (!eq) PetscCall(PetscSNPrintf(msg, sizeof(msg), "Number of labels in dm0 = %" PetscInt_FMT " != %" PetscInt_FMT " = Number of labels in dm1", n, n1));
7572:     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPIU_BOOL, MPI_LAND, comm));
7573:     if (!eq) goto finish;
7574:   }
7575:   for (i = 0; i < n; i++) {
7576:     DMLabel     l0, l1;
7577:     const char *name;
7578:     char       *msgInner;

7580:     /* Ignore label order */
7581:     PetscCall(DMGetLabelByNum(dm0, i, &l0));
7582:     PetscCall(PetscObjectGetName((PetscObject)l0, &name));
7583:     PetscCall(DMGetLabel(dm1, name, &l1));
7584:     if (!l1) {
7585:       PetscCall(PetscSNPrintf(msg, sizeof(msg), "Label \"%s\" (#%" PetscInt_FMT " in dm0) not found in dm1", name, i));
7586:       eq = PETSC_FALSE;
7587:       break;
7588:     }
7589:     PetscCall(DMLabelCompare(comm, l0, l1, &eq, &msgInner));
7590:     PetscCall(PetscStrncpy(msg, msgInner, sizeof(msg)));
7591:     PetscCall(PetscFree(msgInner));
7592:     if (!eq) break;
7593:   }
7594:   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPIU_BOOL, MPI_LAND, comm));
7595: finish:
7596:   /* If message output arg not set, print to stderr */
7597:   if (message) {
7598:     *message = NULL;
7599:     if (msg[0]) PetscCall(PetscStrallocpy(msg, message));
7600:   } else {
7601:     if (msg[0]) PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] %s\n", rank, msg));
7602:     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR));
7603:   }
7604:   /* If same output arg not ser and labels are not equal, throw error */
7605:   if (equal) *equal = eq;
7606:   else PetscCheck(eq, comm, PETSC_ERR_ARG_INCOMP, "DMLabels are not the same in dm0 and dm1");
7607:   PetscFunctionReturn(PETSC_SUCCESS);
7608: }

7610: PetscErrorCode DMSetLabelValue_Fast(DM dm, DMLabel *label, const char name[], PetscInt point, PetscInt value)
7611: {
7612:   PetscFunctionBegin;
7613:   PetscAssertPointer(label, 2);
7614:   if (!*label) {
7615:     PetscCall(DMCreateLabel(dm, name));
7616:     PetscCall(DMGetLabel(dm, name, label));
7617:   }
7618:   PetscCall(DMLabelSetValue(*label, point, value));
7619:   PetscFunctionReturn(PETSC_SUCCESS);
7620: }

7622: /*
7623:   Many mesh programs, such as Triangle and TetGen, allow only a single label for each mesh point. Therefore, we would
7624:   like to encode all label IDs using a single, universal label. We can do this by assigning an integer to every
7625:   (label, id) pair in the DM.

7627:   However, a mesh point can have multiple labels, so we must separate all these values. We will assign a bit range to
7628:   each label.
7629: */
7630: PetscErrorCode DMUniversalLabelCreate(DM dm, DMUniversalLabel *universal)
7631: {
7632:   DMUniversalLabel ul;
7633:   PetscBool       *active;
7634:   PetscInt         pStart, pEnd, p, Nl, l, m;

7636:   PetscFunctionBegin;
7637:   PetscCall(PetscMalloc1(1, &ul));
7638:   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "universal", &ul->label));
7639:   PetscCall(DMGetNumLabels(dm, &Nl));
7640:   PetscCall(PetscCalloc1(Nl, &active));
7641:   ul->Nl = 0;
7642:   for (l = 0; l < Nl; ++l) {
7643:     PetscBool   isdepth, iscelltype;
7644:     const char *name;

7646:     PetscCall(DMGetLabelName(dm, l, &name));
7647:     PetscCall(PetscStrncmp(name, "depth", 6, &isdepth));
7648:     PetscCall(PetscStrncmp(name, "celltype", 9, &iscelltype));
7649:     active[l] = !(isdepth || iscelltype) ? PETSC_TRUE : PETSC_FALSE;
7650:     if (active[l]) ++ul->Nl;
7651:   }
7652:   PetscCall(PetscCalloc5(ul->Nl, &ul->names, ul->Nl, &ul->indices, ul->Nl + 1, &ul->offsets, ul->Nl + 1, &ul->bits, ul->Nl, &ul->masks));
7653:   ul->Nv = 0;
7654:   for (l = 0, m = 0; l < Nl; ++l) {
7655:     DMLabel     label;
7656:     PetscInt    nv;
7657:     const char *name;

7659:     if (!active[l]) continue;
7660:     PetscCall(DMGetLabelName(dm, l, &name));
7661:     PetscCall(DMGetLabelByNum(dm, l, &label));
7662:     PetscCall(DMLabelGetNumValues(label, &nv));
7663:     PetscCall(PetscStrallocpy(name, &ul->names[m]));
7664:     ul->indices[m] = l;
7665:     ul->Nv += nv;
7666:     ul->offsets[m + 1] = nv;
7667:     ul->bits[m + 1]    = PetscCeilReal(PetscLog2Real(nv + 1));
7668:     ++m;
7669:   }
7670:   for (l = 1; l <= ul->Nl; ++l) {
7671:     ul->offsets[l] = ul->offsets[l - 1] + ul->offsets[l];
7672:     ul->bits[l]    = ul->bits[l - 1] + ul->bits[l];
7673:   }
7674:   for (l = 0; l < ul->Nl; ++l) {
7675:     PetscInt b;

7677:     ul->masks[l] = 0;
7678:     for (b = ul->bits[l]; b < ul->bits[l + 1]; ++b) ul->masks[l] |= 1 << b;
7679:   }
7680:   PetscCall(PetscMalloc1(ul->Nv, &ul->values));
7681:   for (l = 0, m = 0; l < Nl; ++l) {
7682:     DMLabel         label;
7683:     IS              valueIS;
7684:     const PetscInt *varr;
7685:     PetscInt        nv, v;

7687:     if (!active[l]) continue;
7688:     PetscCall(DMGetLabelByNum(dm, l, &label));
7689:     PetscCall(DMLabelGetNumValues(label, &nv));
7690:     PetscCall(DMLabelGetValueIS(label, &valueIS));
7691:     PetscCall(ISGetIndices(valueIS, &varr));
7692:     for (v = 0; v < nv; ++v) ul->values[ul->offsets[m] + v] = varr[v];
7693:     PetscCall(ISRestoreIndices(valueIS, &varr));
7694:     PetscCall(ISDestroy(&valueIS));
7695:     PetscCall(PetscSortInt(nv, &ul->values[ul->offsets[m]]));
7696:     ++m;
7697:   }
7698:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
7699:   for (p = pStart; p < pEnd; ++p) {
7700:     PetscInt  uval   = 0;
7701:     PetscBool marked = PETSC_FALSE;

7703:     for (l = 0, m = 0; l < Nl; ++l) {
7704:       DMLabel  label;
7705:       PetscInt val, defval, loc, nv;

7707:       if (!active[l]) continue;
7708:       PetscCall(DMGetLabelByNum(dm, l, &label));
7709:       PetscCall(DMLabelGetValue(label, p, &val));
7710:       PetscCall(DMLabelGetDefaultValue(label, &defval));
7711:       if (val == defval) {
7712:         ++m;
7713:         continue;
7714:       }
7715:       nv     = ul->offsets[m + 1] - ul->offsets[m];
7716:       marked = PETSC_TRUE;
7717:       PetscCall(PetscFindInt(val, nv, &ul->values[ul->offsets[m]], &loc));
7718:       PetscCheck(loc >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Label value %" PetscInt_FMT " not found in compression array", val);
7719:       uval += (loc + 1) << ul->bits[m];
7720:       ++m;
7721:     }
7722:     if (marked) PetscCall(DMLabelSetValue(ul->label, p, uval));
7723:   }
7724:   PetscCall(PetscFree(active));
7725:   *universal = ul;
7726:   PetscFunctionReturn(PETSC_SUCCESS);
7727: }

7729: PetscErrorCode DMUniversalLabelDestroy(DMUniversalLabel *universal)
7730: {
7731:   PetscInt l;

7733:   PetscFunctionBegin;
7734:   for (l = 0; l < (*universal)->Nl; ++l) PetscCall(PetscFree((*universal)->names[l]));
7735:   PetscCall(DMLabelDestroy(&(*universal)->label));
7736:   PetscCall(PetscFree5((*universal)->names, (*universal)->indices, (*universal)->offsets, (*universal)->bits, (*universal)->masks));
7737:   PetscCall(PetscFree((*universal)->values));
7738:   PetscCall(PetscFree(*universal));
7739:   *universal = NULL;
7740:   PetscFunctionReturn(PETSC_SUCCESS);
7741: }

7743: PetscErrorCode DMUniversalLabelGetLabel(DMUniversalLabel ul, DMLabel *ulabel)
7744: {
7745:   PetscFunctionBegin;
7746:   PetscAssertPointer(ulabel, 2);
7747:   *ulabel = ul->label;
7748:   PetscFunctionReturn(PETSC_SUCCESS);
7749: }

7751: PetscErrorCode DMUniversalLabelCreateLabels(DMUniversalLabel ul, PetscBool preserveOrder, DM dm)
7752: {
7753:   PetscInt Nl = ul->Nl, l;

7755:   PetscFunctionBegin;
7757:   for (l = 0; l < Nl; ++l) {
7758:     if (preserveOrder) PetscCall(DMCreateLabelAtIndex(dm, ul->indices[l], ul->names[l]));
7759:     else PetscCall(DMCreateLabel(dm, ul->names[l]));
7760:   }
7761:   if (preserveOrder) {
7762:     for (l = 0; l < ul->Nl; ++l) {
7763:       const char *name;
7764:       PetscBool   match;

7766:       PetscCall(DMGetLabelName(dm, ul->indices[l], &name));
7767:       PetscCall(PetscStrcmp(name, ul->names[l], &match));
7768:       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]);
7769:     }
7770:   }
7771:   PetscFunctionReturn(PETSC_SUCCESS);
7772: }

7774: PetscErrorCode DMUniversalLabelSetLabelValue(DMUniversalLabel ul, DM dm, PetscBool useIndex, PetscInt p, PetscInt value)
7775: {
7776:   PetscInt l;

7778:   PetscFunctionBegin;
7779:   for (l = 0; l < ul->Nl; ++l) {
7780:     DMLabel  label;
7781:     PetscInt lval = (value & ul->masks[l]) >> ul->bits[l];

7783:     if (lval) {
7784:       if (useIndex) PetscCall(DMGetLabelByNum(dm, ul->indices[l], &label));
7785:       else PetscCall(DMGetLabel(dm, ul->names[l], &label));
7786:       PetscCall(DMLabelSetValue(label, p, ul->values[ul->offsets[l] + lval - 1]));
7787:     }
7788:   }
7789:   PetscFunctionReturn(PETSC_SUCCESS);
7790: }

7792: /*@
7793:   DMGetCoarseDM - Get the coarse `DM`from which this `DM` was obtained by refinement

7795:   Not Collective

7797:   Input Parameter:
7798: . dm - The `DM` object

7800:   Output Parameter:
7801: . cdm - The coarse `DM`

7803:   Level: intermediate

7805: .seealso: [](ch_dmbase), `DM`, `DMSetCoarseDM()`, `DMCoarsen()`
7806: @*/
7807: PetscErrorCode DMGetCoarseDM(DM dm, DM *cdm)
7808: {
7809:   PetscFunctionBegin;
7811:   PetscAssertPointer(cdm, 2);
7812:   *cdm = dm->coarseMesh;
7813:   PetscFunctionReturn(PETSC_SUCCESS);
7814: }

7816: /*@
7817:   DMSetCoarseDM - Set the coarse `DM` from which this `DM` was obtained by refinement

7819:   Input Parameters:
7820: + dm  - The `DM` object
7821: - cdm - The coarse `DM`

7823:   Level: intermediate

7825:   Note:
7826:   Normally this is set automatically by `DMRefine()`

7828: .seealso: [](ch_dmbase), `DM`, `DMGetCoarseDM()`, `DMCoarsen()`, `DMSetRefine()`, `DMSetFineDM()`
7829: @*/
7830: PetscErrorCode DMSetCoarseDM(DM dm, DM cdm)
7831: {
7832:   PetscFunctionBegin;
7835:   if (dm == cdm) cdm = NULL;
7836:   PetscCall(PetscObjectReference((PetscObject)cdm));
7837:   PetscCall(DMDestroy(&dm->coarseMesh));
7838:   dm->coarseMesh = cdm;
7839:   PetscFunctionReturn(PETSC_SUCCESS);
7840: }

7842: /*@
7843:   DMGetFineDM - Get the fine mesh from which this `DM` was obtained by coarsening

7845:   Input Parameter:
7846: . dm - The `DM` object

7848:   Output Parameter:
7849: . fdm - The fine `DM`

7851:   Level: intermediate

7853: .seealso: [](ch_dmbase), `DM`, `DMSetFineDM()`, `DMCoarsen()`, `DMRefine()`
7854: @*/
7855: PetscErrorCode DMGetFineDM(DM dm, DM *fdm)
7856: {
7857:   PetscFunctionBegin;
7859:   PetscAssertPointer(fdm, 2);
7860:   *fdm = dm->fineMesh;
7861:   PetscFunctionReturn(PETSC_SUCCESS);
7862: }

7864: /*@
7865:   DMSetFineDM - Set the fine mesh from which this was obtained by coarsening

7867:   Input Parameters:
7868: + dm  - The `DM` object
7869: - fdm - The fine `DM`

7871:   Level: developer

7873:   Note:
7874:   Normally this is set automatically by `DMCoarsen()`

7876: .seealso: [](ch_dmbase), `DM`, `DMGetFineDM()`, `DMCoarsen()`, `DMRefine()`
7877: @*/
7878: PetscErrorCode DMSetFineDM(DM dm, DM fdm)
7879: {
7880:   PetscFunctionBegin;
7883:   if (dm == fdm) fdm = NULL;
7884:   PetscCall(PetscObjectReference((PetscObject)fdm));
7885:   PetscCall(DMDestroy(&dm->fineMesh));
7886:   dm->fineMesh = fdm;
7887:   PetscFunctionReturn(PETSC_SUCCESS);
7888: }

7890: /*@C
7891:   DMAddBoundary - Add a boundary condition to a model represented by a `DM`

7893:   Collective

7895:   Input Parameters:
7896: + dm       - The `DM`, with a `PetscDS` that matches the problem being constrained
7897: . type     - The type of condition, e.g. `DM_BC_ESSENTIAL_ANALYTIC`, `DM_BC_ESSENTIAL_FIELD` (Dirichlet), or `DM_BC_NATURAL` (Neumann)
7898: . name     - The BC name
7899: . label    - The label defining constrained points
7900: . Nv       - The number of `DMLabel` values for constrained points
7901: . values   - An array of values for constrained points
7902: . field    - The field to constrain
7903: . Nc       - The number of constrained field components (0 will constrain all fields)
7904: . comps    - An array of constrained component numbers
7905: . bcFunc   - A pointwise function giving boundary values
7906: . bcFunc_t - A pointwise function giving the time deriative of the boundary values, or NULL
7907: - ctx      - An optional user context for bcFunc

7909:   Output Parameter:
7910: . bd - (Optional) Boundary number

7912:   Options Database Keys:
7913: + -bc_<boundary name> <num>      - Overrides the boundary ids
7914: - -bc_<boundary name>_comp <num> - Overrides the boundary components

7916:   Level: intermediate

7918:   Notes:
7919:   Both bcFunc abd bcFunc_t will depend on the boundary condition type. If the type if `DM_BC_ESSENTIAL`, then the calling sequence is\:
7920: .vb
7921:  void bcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[])
7922: .ve

7924:   If the type is `DM_BC_ESSENTIAL_FIELD` or other _FIELD value, then the calling sequence is\:

7926: .vb
7927:   void bcFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
7928:               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
7929:               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
7930:               PetscReal time, const PetscReal x[], PetscScalar bcval[])
7931: .ve
7932: + dim - the spatial dimension
7933: . Nf - the number of fields
7934: . uOff - the offset into u[] and u_t[] for each field
7935: . uOff_x - the offset into u_x[] for each field
7936: . u - each field evaluated at the current point
7937: . u_t - the time derivative of each field evaluated at the current point
7938: . u_x - the gradient of each field evaluated at the current point
7939: . aOff - the offset into a[] and a_t[] for each auxiliary field
7940: . aOff_x - the offset into a_x[] for each auxiliary field
7941: . a - each auxiliary field evaluated at the current point
7942: . a_t - the time derivative of each auxiliary field evaluated at the current point
7943: . a_x - the gradient of auxiliary each field evaluated at the current point
7944: . t - current time
7945: . x - coordinates of the current point
7946: . numConstants - number of constant parameters
7947: . constants - constant parameters
7948: - bcval - output values at the current point

7950: .seealso: [](ch_dmbase), `DM`, `DSGetBoundary()`, `PetscDSAddBoundary()`
7951: @*/
7952: 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)
7953: {
7954:   PetscDS ds;

7956:   PetscFunctionBegin;
7963:   PetscCheck(!dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot add boundary to DM after creating local section");
7964:   PetscCall(DMGetDS(dm, &ds));
7965:   /* Complete label */
7966:   if (label) {
7967:     PetscObject  obj;
7968:     PetscClassId id;

7970:     PetscCall(DMGetField(dm, field, NULL, &obj));
7971:     PetscCall(PetscObjectGetClassId(obj, &id));
7972:     if (id == PETSCFE_CLASSID) {
7973:       DM plex;

7975:       PetscCall(DMConvert(dm, DMPLEX, &plex));
7976:       if (plex) PetscCall(DMPlexLabelComplete(plex, label));
7977:       PetscCall(DMDestroy(&plex));
7978:     }
7979:   }
7980:   PetscCall(PetscDSAddBoundary(ds, type, name, label, Nv, values, field, Nc, comps, bcFunc, bcFunc_t, ctx, bd));
7981:   PetscFunctionReturn(PETSC_SUCCESS);
7982: }

7984: /* TODO Remove this since now the structures are the same */
7985: static PetscErrorCode DMPopulateBoundary(DM dm)
7986: {
7987:   PetscDS     ds;
7988:   DMBoundary *lastnext;
7989:   DSBoundary  dsbound;

7991:   PetscFunctionBegin;
7992:   PetscCall(DMGetDS(dm, &ds));
7993:   dsbound = ds->boundary;
7994:   if (dm->boundary) {
7995:     DMBoundary next = dm->boundary;

7997:     /* quick check to see if the PetscDS has changed */
7998:     if (next->dsboundary == dsbound) PetscFunctionReturn(PETSC_SUCCESS);
7999:     /* the PetscDS has changed: tear down and rebuild */
8000:     while (next) {
8001:       DMBoundary b = next;

8003:       next = b->next;
8004:       PetscCall(PetscFree(b));
8005:     }
8006:     dm->boundary = NULL;
8007:   }

8009:   lastnext = &dm->boundary;
8010:   while (dsbound) {
8011:     DMBoundary dmbound;

8013:     PetscCall(PetscNew(&dmbound));
8014:     dmbound->dsboundary = dsbound;
8015:     dmbound->label      = dsbound->label;
8016:     /* push on the back instead of the front so that it is in the same order as in the PetscDS */
8017:     *lastnext = dmbound;
8018:     lastnext  = &dmbound->next;
8019:     dsbound   = dsbound->next;
8020:   }
8021:   PetscFunctionReturn(PETSC_SUCCESS);
8022: }

8024: /* TODO: missing manual page */
8025: PetscErrorCode DMIsBoundaryPoint(DM dm, PetscInt point, PetscBool *isBd)
8026: {
8027:   DMBoundary b;

8029:   PetscFunctionBegin;
8031:   PetscAssertPointer(isBd, 3);
8032:   *isBd = PETSC_FALSE;
8033:   PetscCall(DMPopulateBoundary(dm));
8034:   b = dm->boundary;
8035:   while (b && !(*isBd)) {
8036:     DMLabel    label = b->label;
8037:     DSBoundary dsb   = b->dsboundary;
8038:     PetscInt   i;

8040:     if (label) {
8041:       for (i = 0; i < dsb->Nv && !(*isBd); ++i) PetscCall(DMLabelStratumHasPoint(label, dsb->values[i], point, isBd));
8042:     }
8043:     b = b->next;
8044:   }
8045:   PetscFunctionReturn(PETSC_SUCCESS);
8046: }

8048: /*@C
8049:   DMProjectFunction - This projects the given function into the function space provided by a `DM`, putting the coefficients in a global vector.

8051:   Collective

8053:   Input Parameters:
8054: + dm    - The `DM`
8055: . time  - The time
8056: . funcs - The coordinate functions to evaluate, one per field
8057: . ctxs  - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
8058: - mode  - The insertion mode for values

8060:   Output Parameter:
8061: . X - vector

8063:   Calling sequence of `funcs`:
8064: + dim  - The spatial dimension
8065: . time - The time at which to sample
8066: . x    - The coordinates
8067: . Nc   - The number of components
8068: . u    - The output field values
8069: - ctx  - optional user-defined function context

8071:   Level: developer

8073:   Developer Notes:
8074:   This API is specific to only particular usage of `DM`

8076:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8078: .seealso: [](ch_dmbase), `DM`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
8079: @*/
8080: 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)
8081: {
8082:   Vec localX;

8084:   PetscFunctionBegin;
8086:   PetscCall(PetscLogEventBegin(DM_ProjectFunction, dm, X, 0, 0));
8087:   PetscCall(DMGetLocalVector(dm, &localX));
8088:   PetscCall(VecSet(localX, 0.));
8089:   PetscCall(DMProjectFunctionLocal(dm, time, funcs, ctxs, mode, localX));
8090:   PetscCall(DMLocalToGlobalBegin(dm, localX, mode, X));
8091:   PetscCall(DMLocalToGlobalEnd(dm, localX, mode, X));
8092:   PetscCall(DMRestoreLocalVector(dm, &localX));
8093:   PetscCall(PetscLogEventEnd(DM_ProjectFunction, dm, X, 0, 0));
8094:   PetscFunctionReturn(PETSC_SUCCESS);
8095: }

8097: /*@C
8098:   DMProjectFunctionLocal - This projects the given function into the function space provided by a `DM`, putting the coefficients in a local vector.

8100:   Not Collective

8102:   Input Parameters:
8103: + dm    - The `DM`
8104: . time  - The time
8105: . funcs - The coordinate functions to evaluate, one per field
8106: . ctxs  - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
8107: - mode  - The insertion mode for values

8109:   Output Parameter:
8110: . localX - vector

8112:   Calling sequence of `funcs`:
8113: + dim  - The spatial dimension
8114: . time - The current timestep
8115: . x    - The coordinates
8116: . Nc   - The number of components
8117: . u    - The output field values
8118: - ctx  - optional user-defined function context

8120:   Level: developer

8122:   Developer Notes:
8123:   This API is specific to only particular usage of `DM`

8125:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8127: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
8128: @*/
8129: 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)
8130: {
8131:   PetscFunctionBegin;
8134:   PetscUseTypeMethod(dm, projectfunctionlocal, time, funcs, ctxs, mode, localX);
8135:   PetscFunctionReturn(PETSC_SUCCESS);
8136: }

8138: /*@C
8139:   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.

8141:   Collective

8143:   Input Parameters:
8144: + dm     - The `DM`
8145: . time   - The time
8146: . numIds - The number of ids
8147: . ids    - The ids
8148: . Nc     - The number of components
8149: . comps  - The components
8150: . label  - The `DMLabel` selecting the portion of the mesh for projection
8151: . funcs  - The coordinate functions to evaluate, one per field
8152: . ctxs   - Optional array of contexts to pass to each coordinate function.  ctxs may be null.
8153: - mode   - The insertion mode for values

8155:   Output Parameter:
8156: . X - vector

8158:   Calling sequence of `funcs`:
8159: + dim  - The spatial dimension
8160: . time - The current timestep
8161: . x    - The coordinates
8162: . Nc   - The number of components
8163: . u    - The output field values
8164: - ctx  - optional user-defined function context

8166:   Level: developer

8168:   Developer Notes:
8169:   This API is specific to only particular usage of `DM`

8171:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8173: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabelLocal()`, `DMComputeL2Diff()`
8174: @*/
8175: 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)
8176: {
8177:   Vec localX;

8179:   PetscFunctionBegin;
8181:   PetscCall(DMGetLocalVector(dm, &localX));
8182:   PetscCall(VecSet(localX, 0.));
8183:   PetscCall(DMProjectFunctionLabelLocal(dm, time, label, numIds, ids, Nc, comps, funcs, ctxs, mode, localX));
8184:   PetscCall(DMLocalToGlobalBegin(dm, localX, mode, X));
8185:   PetscCall(DMLocalToGlobalEnd(dm, localX, mode, X));
8186:   PetscCall(DMRestoreLocalVector(dm, &localX));
8187:   PetscFunctionReturn(PETSC_SUCCESS);
8188: }

8190: /*@C
8191:   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.

8193:   Not Collective

8195:   Input Parameters:
8196: + dm     - The `DM`
8197: . time   - The time
8198: . label  - The `DMLabel` selecting the portion of the mesh for projection
8199: . numIds - The number of ids
8200: . ids    - The ids
8201: . Nc     - The number of components
8202: . comps  - The components
8203: . funcs  - The coordinate functions to evaluate, one per field
8204: . ctxs   - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
8205: - mode   - The insertion mode for values

8207:   Output Parameter:
8208: . localX - vector

8210:   Calling sequence of `funcs`:
8211: + dim  - The spatial dimension
8212: . time - The current time
8213: . x    - The coordinates
8214: . Nc   - The number of components
8215: . u    - The output field values
8216: - ctx  - optional user-defined function context

8218:   Level: developer

8220:   Developer Notes:
8221:   This API is specific to only particular usage of `DM`

8223:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8225: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
8226: @*/
8227: 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)
8228: {
8229:   PetscFunctionBegin;
8232:   PetscUseTypeMethod(dm, projectfunctionlabellocal, time, label, numIds, ids, Nc, comps, funcs, ctxs, mode, localX);
8233:   PetscFunctionReturn(PETSC_SUCCESS);
8234: }

8236: /*@C
8237:   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.

8239:   Not Collective

8241:   Input Parameters:
8242: + dm     - The `DM`
8243: . time   - The time
8244: . localU - The input field vector; may be `NULL` if projection is defined purely by coordinates
8245: . funcs  - The functions to evaluate, one per field
8246: - mode   - The insertion mode for values

8248:   Output Parameter:
8249: . localX - The output vector

8251:   Calling sequence of `funcs`:
8252: + dim          - The spatial dimension
8253: . Nf           - The number of input fields
8254: . NfAux        - The number of input auxiliary fields
8255: . uOff         - The offset of each field in u[]
8256: . uOff_x       - The offset of each field in u_x[]
8257: . u            - The field values at this point in space
8258: . u_t          - The field time derivative at this point in space (or NULL)
8259: . u_x          - The field derivatives at this point in space
8260: . aOff         - The offset of each auxiliary field in u[]
8261: . aOff_x       - The offset of each auxiliary field in u_x[]
8262: . a            - The auxiliary field values at this point in space
8263: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8264: . a_x          - The auxiliary field derivatives at this point in space
8265: . t            - The current time
8266: . x            - The coordinates of this point
8267: . numConstants - The number of constants
8268: . constants    - The value of each constant
8269: - f            - The value of the function at this point in space

8271:   Level: intermediate

8273:   Note:
8274:   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.
8275:   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
8276:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8277:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8279:   Developer Notes:
8280:   This API is specific to only particular usage of `DM`

8282:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8284: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabelLocal()`,
8285: `DMProjectFunction()`, `DMComputeL2Diff()`
8286: @*/
8287: 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)
8288: {
8289:   PetscFunctionBegin;
8293:   PetscUseTypeMethod(dm, projectfieldlocal, time, localU, funcs, mode, localX);
8294:   PetscFunctionReturn(PETSC_SUCCESS);
8295: }

8297: /*@C
8298:   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.

8300:   Not Collective

8302:   Input Parameters:
8303: + dm     - The `DM`
8304: . time   - The time
8305: . label  - The `DMLabel` marking the portion of the domain to output
8306: . numIds - The number of label ids to use
8307: . ids    - The label ids to use for marking
8308: . Nc     - The number of components to set in the output, or `PETSC_DETERMINE` for all components
8309: . comps  - The components to set in the output, or `NULL` for all components
8310: . localU - The input field vector
8311: . funcs  - The functions to evaluate, one per field
8312: - mode   - The insertion mode for values

8314:   Output Parameter:
8315: . localX - The output vector

8317:   Calling sequence of `funcs`:
8318: + dim          - The spatial dimension
8319: . Nf           - The number of input fields
8320: . NfAux        - The number of input auxiliary fields
8321: . uOff         - The offset of each field in u[]
8322: . uOff_x       - The offset of each field in u_x[]
8323: . u            - The field values at this point in space
8324: . u_t          - The field time derivative at this point in space (or NULL)
8325: . u_x          - The field derivatives at this point in space
8326: . aOff         - The offset of each auxiliary field in u[]
8327: . aOff_x       - The offset of each auxiliary field in u_x[]
8328: . a            - The auxiliary field values at this point in space
8329: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8330: . a_x          - The auxiliary field derivatives at this point in space
8331: . t            - The current time
8332: . x            - The coordinates of this point
8333: . numConstants - The number of constants
8334: . constants    - The value of each constant
8335: - f            - The value of the function at this point in space

8337:   Level: intermediate

8339:   Note:
8340:   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.
8341:   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
8342:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8343:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8345:   Developer Notes:
8346:   This API is specific to only particular usage of `DM`

8348:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8350: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabel()`, `DMProjectFunction()`, `DMComputeL2Diff()`
8351: @*/
8352: 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)
8353: {
8354:   PetscFunctionBegin;
8358:   PetscUseTypeMethod(dm, projectfieldlabellocal, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX);
8359:   PetscFunctionReturn(PETSC_SUCCESS);
8360: }

8362: /*@C
8363:   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.

8365:   Not Collective

8367:   Input Parameters:
8368: + dm     - The `DM`
8369: . time   - The time
8370: . label  - The `DMLabel` marking the portion of the domain to output
8371: . numIds - The number of label ids to use
8372: . ids    - The label ids to use for marking
8373: . Nc     - The number of components to set in the output, or `PETSC_DETERMINE` for all components
8374: . comps  - The components to set in the output, or `NULL` for all components
8375: . U      - The input field vector
8376: . funcs  - The functions to evaluate, one per field
8377: - mode   - The insertion mode for values

8379:   Output Parameter:
8380: . X - The output vector

8382:   Calling sequence of `funcs`:
8383: + dim          - The spatial dimension
8384: . Nf           - The number of input fields
8385: . NfAux        - The number of input auxiliary fields
8386: . uOff         - The offset of each field in u[]
8387: . uOff_x       - The offset of each field in u_x[]
8388: . u            - The field values at this point in space
8389: . u_t          - The field time derivative at this point in space (or NULL)
8390: . u_x          - The field derivatives at this point in space
8391: . aOff         - The offset of each auxiliary field in u[]
8392: . aOff_x       - The offset of each auxiliary field in u_x[]
8393: . a            - The auxiliary field values at this point in space
8394: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8395: . a_x          - The auxiliary field derivatives at this point in space
8396: . t            - The current time
8397: . x            - The coordinates of this point
8398: . numConstants - The number of constants
8399: . constants    - The value of each constant
8400: - f            - The value of the function at this point in space

8402:   Level: intermediate

8404:   Note:
8405:   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.
8406:   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
8407:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8408:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8410:   Developer Notes:
8411:   This API is specific to only particular usage of `DM`

8413:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8415: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabelLocal()`, `DMProjectFunction()`, `DMComputeL2Diff()`
8416: @*/
8417: 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)
8418: {
8419:   DM  dmIn;
8420:   Vec localU, localX;

8422:   PetscFunctionBegin;
8424:   PetscCall(VecGetDM(U, &dmIn));
8425:   PetscCall(DMGetLocalVector(dmIn, &localU));
8426:   PetscCall(DMGetLocalVector(dm, &localX));
8427:   PetscCall(VecSet(localX, 0.));
8428:   PetscCall(DMGlobalToLocalBegin(dmIn, U, mode, localU));
8429:   PetscCall(DMGlobalToLocalEnd(dmIn, U, mode, localU));
8430:   PetscCall(DMProjectFieldLabelLocal(dm, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX));
8431:   PetscCall(DMLocalToGlobalBegin(dm, localX, mode, X));
8432:   PetscCall(DMLocalToGlobalEnd(dm, localX, mode, X));
8433:   PetscCall(DMRestoreLocalVector(dm, &localX));
8434:   PetscCall(DMRestoreLocalVector(dmIn, &localU));
8435:   PetscFunctionReturn(PETSC_SUCCESS);
8436: }

8438: /*@C
8439:   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.

8441:   Not Collective

8443:   Input Parameters:
8444: + dm     - The `DM`
8445: . time   - The time
8446: . label  - The `DMLabel` marking the portion of the domain boundary to output
8447: . numIds - The number of label ids to use
8448: . ids    - The label ids to use for marking
8449: . Nc     - The number of components to set in the output, or `PETSC_DETERMINE` for all components
8450: . comps  - The components to set in the output, or `NULL` for all components
8451: . localU - The input field vector
8452: . funcs  - The functions to evaluate, one per field
8453: - mode   - The insertion mode for values

8455:   Output Parameter:
8456: . localX - The output vector

8458:   Calling sequence of `funcs`:
8459: + dim          - The spatial dimension
8460: . Nf           - The number of input fields
8461: . NfAux        - The number of input auxiliary fields
8462: . uOff         - The offset of each field in u[]
8463: . uOff_x       - The offset of each field in u_x[]
8464: . u            - The field values at this point in space
8465: . u_t          - The field time derivative at this point in space (or NULL)
8466: . u_x          - The field derivatives at this point in space
8467: . aOff         - The offset of each auxiliary field in u[]
8468: . aOff_x       - The offset of each auxiliary field in u_x[]
8469: . a            - The auxiliary field values at this point in space
8470: . a_t          - The auxiliary field time derivative at this point in space (or NULL)
8471: . a_x          - The auxiliary field derivatives at this point in space
8472: . t            - The current time
8473: . x            - The coordinates of this point
8474: . n            - The face normal
8475: . numConstants - The number of constants
8476: . constants    - The value of each constant
8477: - f            - The value of the function at this point in space

8479:   Level: intermediate

8481:   Note:
8482:   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.
8483:   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
8484:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary `DM`, attached to the
8485:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

8487:   Developer Notes:
8488:   This API is specific to only particular usage of `DM`

8490:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8492: .seealso: [](ch_dmbase), `DM`, `DMProjectField()`, `DMProjectFieldLabelLocal()`, `DMProjectFunction()`, `DMComputeL2Diff()`
8493: @*/
8494: 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)
8495: {
8496:   PetscFunctionBegin;
8500:   PetscUseTypeMethod(dm, projectbdfieldlabellocal, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX);
8501:   PetscFunctionReturn(PETSC_SUCCESS);
8502: }

8504: /*@C
8505:   DMComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.

8507:   Collective

8509:   Input Parameters:
8510: + dm    - The `DM`
8511: . time  - The time
8512: . funcs - The functions to evaluate for each field component
8513: . ctxs  - Optional array of contexts to pass to each function, or NULL.
8514: - X     - The coefficient vector u_h, a global vector

8516:   Output Parameter:
8517: . diff - The diff ||u - u_h||_2

8519:   Level: developer

8521:   Developer Notes:
8522:   This API is specific to only particular usage of `DM`

8524:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8526: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMComputeL2FieldDiff()`, `DMComputeL2GradientDiff()`
8527: @*/
8528: PetscErrorCode DMComputeL2Diff(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
8529: {
8530:   PetscFunctionBegin;
8533:   PetscUseTypeMethod(dm, computel2diff, time, funcs, ctxs, X, diff);
8534:   PetscFunctionReturn(PETSC_SUCCESS);
8535: }

8537: /*@C
8538:   DMComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h.

8540:   Collective

8542:   Input Parameters:
8543: + dm    - The `DM`
8544: . time  - The time
8545: . funcs - The gradient functions to evaluate for each field component
8546: . ctxs  - Optional array of contexts to pass to each function, or NULL.
8547: . X     - The coefficient vector u_h, a global vector
8548: - n     - The vector to project along

8550:   Output Parameter:
8551: . diff - The diff ||(grad u - grad u_h) . n||_2

8553:   Level: developer

8555:   Developer Notes:
8556:   This API is specific to only particular usage of `DM`

8558:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8560: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMComputeL2Diff()`, `DMComputeL2FieldDiff()`
8561: @*/
8562: 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)
8563: {
8564:   PetscFunctionBegin;
8567:   PetscUseTypeMethod(dm, computel2gradientdiff, time, funcs, ctxs, X, n, diff);
8568:   PetscFunctionReturn(PETSC_SUCCESS);
8569: }

8571: /*@C
8572:   DMComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components.

8574:   Collective

8576:   Input Parameters:
8577: + dm    - The `DM`
8578: . time  - The time
8579: . funcs - The functions to evaluate for each field component
8580: . ctxs  - Optional array of contexts to pass to each function, or NULL.
8581: - X     - The coefficient vector u_h, a global vector

8583:   Output Parameter:
8584: . diff - The array of differences, ||u^f - u^f_h||_2

8586:   Level: developer

8588:   Developer Notes:
8589:   This API is specific to only particular usage of `DM`

8591:   The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.

8593: .seealso: [](ch_dmbase), `DM`, `DMProjectFunction()`, `DMComputeL2GradientDiff()`
8594: @*/
8595: PetscErrorCode DMComputeL2FieldDiff(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
8596: {
8597:   PetscFunctionBegin;
8600:   PetscUseTypeMethod(dm, computel2fielddiff, time, funcs, ctxs, X, diff);
8601:   PetscFunctionReturn(PETSC_SUCCESS);
8602: }

8604: /*@C
8605:   DMGetNeighbors - Gets an array containing the MPI ranks of all the processes neighbors

8607:   Not Collective

8609:   Input Parameter:
8610: . dm - The `DM`

8612:   Output Parameters:
8613: + nranks - the number of neighbours
8614: - ranks  - the neighbors ranks

8616:   Level: beginner

8618:   Note:
8619:   Do not free the array, it is freed when the `DM` is destroyed.

8621: .seealso: [](ch_dmbase), `DM`, `DMDAGetNeighbors()`, `PetscSFGetRootRanks()`
8622: @*/
8623: PetscErrorCode DMGetNeighbors(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
8624: {
8625:   PetscFunctionBegin;
8627:   PetscUseTypeMethod(dm, getneighbors, nranks, ranks);
8628:   PetscFunctionReturn(PETSC_SUCCESS);
8629: }

8631: #include <petsc/private/matimpl.h>

8633: /*
8634:     Converts the input vector to a ghosted vector and then calls the standard coloring code.
8635:     This must be a different function because it requires DM which is not defined in the Mat library
8636: */
8637: static PetscErrorCode MatFDColoringApply_AIJDM(Mat J, MatFDColoring coloring, Vec x1, void *sctx)
8638: {
8639:   PetscFunctionBegin;
8640:   if (coloring->ctype == IS_COLORING_LOCAL) {
8641:     Vec x1local;
8642:     DM  dm;
8643:     PetscCall(MatGetDM(J, &dm));
8644:     PetscCheck(dm, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_INCOMP, "IS_COLORING_LOCAL requires a DM");
8645:     PetscCall(DMGetLocalVector(dm, &x1local));
8646:     PetscCall(DMGlobalToLocalBegin(dm, x1, INSERT_VALUES, x1local));
8647:     PetscCall(DMGlobalToLocalEnd(dm, x1, INSERT_VALUES, x1local));
8648:     x1 = x1local;
8649:   }
8650:   PetscCall(MatFDColoringApply_AIJ(J, coloring, x1, sctx));
8651:   if (coloring->ctype == IS_COLORING_LOCAL) {
8652:     DM dm;
8653:     PetscCall(MatGetDM(J, &dm));
8654:     PetscCall(DMRestoreLocalVector(dm, &x1));
8655:   }
8656:   PetscFunctionReturn(PETSC_SUCCESS);
8657: }

8659: /*@
8660:   MatFDColoringUseDM - allows a `MatFDColoring` object to use the `DM` associated with the matrix to compute a `IS_COLORING_LOCAL` coloring

8662:   Input Parameters:
8663: + coloring   - The matrix to get the `DM` from
8664: - fdcoloring - the `MatFDColoring` object

8666:   Level: advanced

8668:   Developer Note:
8669:   This routine exists because the PETSc `Mat` library does not know about the `DM` objects

8671: .seealso: [](ch_dmbase), `DM`, `MatFDColoring`, `MatFDColoringCreate()`, `ISColoringType`
8672: @*/
8673: PetscErrorCode MatFDColoringUseDM(Mat coloring, MatFDColoring fdcoloring)
8674: {
8675:   PetscFunctionBegin;
8676:   coloring->ops->fdcoloringapply = MatFDColoringApply_AIJDM;
8677:   PetscFunctionReturn(PETSC_SUCCESS);
8678: }

8680: /*@
8681:   DMGetCompatibility - determine if two `DM`s are compatible

8683:   Collective

8685:   Input Parameters:
8686: + dm1 - the first `DM`
8687: - dm2 - the second `DM`

8689:   Output Parameters:
8690: + compatible - whether or not the two `DM`s are compatible
8691: - set        - whether or not the compatible value was actually determined and set

8693:   Level: advanced

8695:   Notes:
8696:   Two `DM`s are deemed compatible if they represent the same parallel decomposition
8697:   of the same topology. This implies that the section (field data) on one
8698:   "makes sense" with respect to the topology and parallel decomposition of the other.
8699:   Loosely speaking, compatible `DM`s represent the same domain and parallel
8700:   decomposition, but hold different data.

8702:   Typically, one would confirm compatibility if intending to simultaneously iterate
8703:   over a pair of vectors obtained from different `DM`s.

8705:   For example, two `DMDA` objects are compatible if they have the same local
8706:   and global sizes and the same stencil width. They can have different numbers
8707:   of degrees of freedom per node. Thus, one could use the node numbering from
8708:   either `DM` in bounds for a loop over vectors derived from either `DM`.

8710:   Consider the operation of summing data living on a 2-dof `DMDA` to data living
8711:   on a 1-dof `DMDA`, which should be compatible, as in the following snippet.
8712: .vb
8713:   ...
8714:   PetscCall(DMGetCompatibility(da1,da2,&compatible,&set));
8715:   if (set && compatible)  {
8716:     PetscCall(DMDAVecGetArrayDOF(da1,vec1,&arr1));
8717:     PetscCall(DMDAVecGetArrayDOF(da2,vec2,&arr2));
8718:     PetscCall(DMDAGetCorners(da1,&x,&y,NULL,&m,&n,NULL));
8719:     for (j=y; j<y+n; ++j) {
8720:       for (i=x; i<x+m, ++i) {
8721:         arr1[j][i][0] = arr2[j][i][0] + arr2[j][i][1];
8722:       }
8723:     }
8724:     PetscCall(DMDAVecRestoreArrayDOF(da1,vec1,&arr1));
8725:     PetscCall(DMDAVecRestoreArrayDOF(da2,vec2,&arr2));
8726:   } else {
8727:     SETERRQ(PetscObjectComm((PetscObject)da1,PETSC_ERR_ARG_INCOMP,"DMDA objects incompatible");
8728:   }
8729:   ...
8730: .ve

8732:   Checking compatibility might be expensive for a given implementation of `DM`,
8733:   or might be impossible to unambiguously confirm or deny. For this reason,
8734:   this function may decline to determine compatibility, and hence users should
8735:   always check the "set" output parameter.

8737:   A `DM` is always compatible with itself.

8739:   In the current implementation, `DM`s which live on "unequal" communicators
8740:   (MPI_UNEQUAL in the terminology of MPI_Comm_compare()) are always deemed
8741:   incompatible.

8743:   This function is labeled "Collective," as information about all subdomains
8744:   is required on each rank. However, in `DM` implementations which store all this
8745:   information locally, this function may be merely "Logically Collective".

8747:   Developer Note:
8748:   Compatibility is assumed to be a symmetric concept; `DM` A is compatible with `DM` B
8749:   iff B is compatible with A. Thus, this function checks the implementations
8750:   of both dm and dmc (if they are of different types), attempting to determine
8751:   compatibility. It is left to `DM` implementers to ensure that symmetry is
8752:   preserved. The simplest way to do this is, when implementing type-specific
8753:   logic for this function, is to check for existing logic in the implementation
8754:   of other `DM` types and let *set = PETSC_FALSE if found.

8756: .seealso: [](ch_dmbase), `DM`, `DMDACreateCompatibleDMDA()`, `DMStagCreateCompatibleDMStag()`
8757: @*/
8758: PetscErrorCode DMGetCompatibility(DM dm1, DM dm2, PetscBool *compatible, PetscBool *set)
8759: {
8760:   PetscMPIInt compareResult;
8761:   DMType      type, type2;
8762:   PetscBool   sameType;

8764:   PetscFunctionBegin;

8768:   /* Declare a DM compatible with itself */
8769:   if (dm1 == dm2) {
8770:     *set        = PETSC_TRUE;
8771:     *compatible = PETSC_TRUE;
8772:     PetscFunctionReturn(PETSC_SUCCESS);
8773:   }

8775:   /* Declare a DM incompatible with a DM that lives on an "unequal"
8776:      communicator. Note that this does not preclude compatibility with
8777:      DMs living on "congruent" or "similar" communicators, but this must be
8778:      determined by the implementation-specific logic */
8779:   PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dm1), PetscObjectComm((PetscObject)dm2), &compareResult));
8780:   if (compareResult == MPI_UNEQUAL) {
8781:     *set        = PETSC_TRUE;
8782:     *compatible = PETSC_FALSE;
8783:     PetscFunctionReturn(PETSC_SUCCESS);
8784:   }

8786:   /* Pass to the implementation-specific routine, if one exists. */
8787:   if (dm1->ops->getcompatibility) {
8788:     PetscUseTypeMethod(dm1, getcompatibility, dm2, compatible, set);
8789:     if (*set) PetscFunctionReturn(PETSC_SUCCESS);
8790:   }

8792:   /* If dm1 and dm2 are of different types, then attempt to check compatibility
8793:      with an implementation of this function from dm2 */
8794:   PetscCall(DMGetType(dm1, &type));
8795:   PetscCall(DMGetType(dm2, &type2));
8796:   PetscCall(PetscStrcmp(type, type2, &sameType));
8797:   if (!sameType && dm2->ops->getcompatibility) {
8798:     PetscUseTypeMethod(dm2, getcompatibility, dm1, compatible, set); /* Note argument order */
8799:   } else {
8800:     *set = PETSC_FALSE;
8801:   }
8802:   PetscFunctionReturn(PETSC_SUCCESS);
8803: }

8805: /*@C
8806:   DMMonitorSet - Sets an additional monitor function that is to be used after a solve to monitor discretization performance.

8808:   Logically Collective

8810:   Input Parameters:
8811: + dm             - the `DM`
8812: . f              - the monitor function
8813: . mctx           - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired)
8814: - monitordestroy - [optional] routine that frees monitor context (may be `NULL`)

8816:   Options Database Key:
8817: . -dm_monitor_cancel - cancels all monitors that have been hardwired into a code by calls to `DMMonitorSet()`, but
8818:                             does not cancel those set via the options database.

8820:   Level: intermediate

8822:   Note:
8823:   Several different monitoring routines may be set by calling
8824:   `DMMonitorSet()` multiple times or with `DMMonitorSetFromOptions()`; all will be called in the
8825:   order in which they were set.

8827:   Fortran Note:
8828:   Only a single monitor function can be set for each `DM` object

8830:   Developer Note:
8831:   This API has a generic name but seems specific to a very particular aspect of the use of `DM`

8833: .seealso: [](ch_dmbase), `DM`, `DMMonitorCancel()`, `DMMonitorSetFromOptions()`, `DMMonitor()`
8834: @*/
8835: PetscErrorCode DMMonitorSet(DM dm, PetscErrorCode (*f)(DM, void *), void *mctx, PetscErrorCode (*monitordestroy)(void **))
8836: {
8837:   PetscInt m;

8839:   PetscFunctionBegin;
8841:   for (m = 0; m < dm->numbermonitors; ++m) {
8842:     PetscBool identical;

8844:     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, mctx, monitordestroy, (PetscErrorCode(*)(void))dm->monitor[m], dm->monitorcontext[m], dm->monitordestroy[m], &identical));
8845:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
8846:   }
8847:   PetscCheck(dm->numbermonitors < MAXDMMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
8848:   dm->monitor[dm->numbermonitors]          = f;
8849:   dm->monitordestroy[dm->numbermonitors]   = monitordestroy;
8850:   dm->monitorcontext[dm->numbermonitors++] = (void *)mctx;
8851:   PetscFunctionReturn(PETSC_SUCCESS);
8852: }

8854: /*@
8855:   DMMonitorCancel - Clears all the monitor functions for a `DM` object.

8857:   Logically Collective

8859:   Input Parameter:
8860: . dm - the DM

8862:   Options Database Key:
8863: . -dm_monitor_cancel - cancels all monitors that have been hardwired
8864:   into a code by calls to `DMonitorSet()`, but does not cancel those
8865:   set via the options database

8867:   Level: intermediate

8869:   Note:
8870:   There is no way to clear one specific monitor from a `DM` object.

8872: .seealso: [](ch_dmbase), `DM`, `DMMonitorSet()`, `DMMonitorSetFromOptions()`, `DMMonitor()`
8873: @*/
8874: PetscErrorCode DMMonitorCancel(DM dm)
8875: {
8876:   PetscInt m;

8878:   PetscFunctionBegin;
8880:   for (m = 0; m < dm->numbermonitors; ++m) {
8881:     if (dm->monitordestroy[m]) PetscCall((*dm->monitordestroy[m])(&dm->monitorcontext[m]));
8882:   }
8883:   dm->numbermonitors = 0;
8884:   PetscFunctionReturn(PETSC_SUCCESS);
8885: }

8887: /*@C
8888:   DMMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user

8890:   Collective

8892:   Input Parameters:
8893: + dm           - `DM` object you wish to monitor
8894: . name         - the monitor type one is seeking
8895: . help         - message indicating what monitoring is done
8896: . manual       - manual page for the monitor
8897: . monitor      - the monitor function
8898: - 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

8900:   Output Parameter:
8901: . flg - Flag set if the monitor was created

8903:   Level: developer

8905: .seealso: [](ch_dmbase), `DM`, `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
8906:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
8907:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
8908:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
8909:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
8910:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
8911:           `PetscOptionsFList()`, `PetscOptionsEList()`, `DMMonitor()`, `DMMonitorSet()`
8912: @*/
8913: PetscErrorCode DMMonitorSetFromOptions(DM dm, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(DM, void *), PetscErrorCode (*monitorsetup)(DM, PetscViewerAndFormat *), PetscBool *flg)
8914: {
8915:   PetscViewer       viewer;
8916:   PetscViewerFormat format;

8918:   PetscFunctionBegin;
8920:   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)dm), ((PetscObject)dm)->options, ((PetscObject)dm)->prefix, name, &viewer, &format, flg));
8921:   if (*flg) {
8922:     PetscViewerAndFormat *vf;

8924:     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
8925:     PetscCall(PetscOptionsRestoreViewer(&viewer));
8926:     if (monitorsetup) PetscCall((*monitorsetup)(dm, vf));
8927:     PetscCall(DMMonitorSet(dm, (PetscErrorCode(*)(DM, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
8928:   }
8929:   PetscFunctionReturn(PETSC_SUCCESS);
8930: }

8932: /*@
8933:   DMMonitor - runs the user provided monitor routines, if they exist

8935:   Collective

8937:   Input Parameter:
8938: . dm - The `DM`

8940:   Level: developer

8942:   Developer Note:
8943:   Note should indicate when during the life of the `DM` the monitor is run. It appears to be
8944:   related to the discretization process seems rather specialized since some `DM` have no
8945:   concept of discretization.

8947: .seealso: [](ch_dmbase), `DM`, `DMMonitorSet()`, `DMMonitorSetFromOptions()`
8948: @*/
8949: PetscErrorCode DMMonitor(DM dm)
8950: {
8951:   PetscInt m;

8953:   PetscFunctionBegin;
8954:   if (!dm) PetscFunctionReturn(PETSC_SUCCESS);
8956:   for (m = 0; m < dm->numbermonitors; ++m) PetscCall((*dm->monitor[m])(dm, dm->monitorcontext[m]));
8957:   PetscFunctionReturn(PETSC_SUCCESS);
8958: }

8960: /*@
8961:   DMComputeError - Computes the error assuming the user has provided the exact solution functions

8963:   Collective

8965:   Input Parameters:
8966: + dm  - The `DM`
8967: - sol - The solution vector

8969:   Input/Output Parameter:
8970: . errors - An array of length Nf, the number of fields, or `NULL` for no output; on output
8971:            contains the error in each field

8973:   Output Parameter:
8974: . errorVec - A vector to hold the cellwise error (may be `NULL`)

8976:   Level: developer

8978:   Note:
8979:   The exact solutions come from the `PetscDS` object, and the time comes from `DMGetOutputSequenceNumber()`.

8981: .seealso: [](ch_dmbase), `DM`, `DMMonitorSet()`, `DMGetRegionNumDS()`, `PetscDSGetExactSolution()`, `DMGetOutputSequenceNumber()`
8982: @*/
8983: PetscErrorCode DMComputeError(DM dm, Vec sol, PetscReal errors[], Vec *errorVec)
8984: {
8985:   PetscErrorCode (**exactSol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
8986:   void    **ctxs;
8987:   PetscReal time;
8988:   PetscInt  Nf, f, Nds, s;

8990:   PetscFunctionBegin;
8991:   PetscCall(DMGetNumFields(dm, &Nf));
8992:   PetscCall(PetscCalloc2(Nf, &exactSol, Nf, &ctxs));
8993:   PetscCall(DMGetNumDS(dm, &Nds));
8994:   for (s = 0; s < Nds; ++s) {
8995:     PetscDS         ds;
8996:     DMLabel         label;
8997:     IS              fieldIS;
8998:     const PetscInt *fields;
8999:     PetscInt        dsNf;

9001:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL));
9002:     PetscCall(PetscDSGetNumFields(ds, &dsNf));
9003:     if (fieldIS) PetscCall(ISGetIndices(fieldIS, &fields));
9004:     for (f = 0; f < dsNf; ++f) {
9005:       const PetscInt field = fields[f];
9006:       PetscCall(PetscDSGetExactSolution(ds, field, &exactSol[field], &ctxs[field]));
9007:     }
9008:     if (fieldIS) PetscCall(ISRestoreIndices(fieldIS, &fields));
9009:   }
9010:   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);
9011:   PetscCall(DMGetOutputSequenceNumber(dm, NULL, &time));
9012:   if (errors) PetscCall(DMComputeL2FieldDiff(dm, time, exactSol, ctxs, sol, errors));
9013:   if (errorVec) {
9014:     DM             edm;
9015:     DMPolytopeType ct;
9016:     PetscBool      simplex;
9017:     PetscInt       dim, cStart, Nf;

9019:     PetscCall(DMClone(dm, &edm));
9020:     PetscCall(DMGetDimension(edm, &dim));
9021:     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
9022:     PetscCall(DMPlexGetCellType(dm, cStart, &ct));
9023:     simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
9024:     PetscCall(DMGetNumFields(dm, &Nf));
9025:     for (f = 0; f < Nf; ++f) {
9026:       PetscFE         fe, efe;
9027:       PetscQuadrature q;
9028:       const char     *name;

9030:       PetscCall(DMGetField(dm, f, NULL, (PetscObject *)&fe));
9031:       PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, Nf, simplex, 0, PETSC_DETERMINE, &efe));
9032:       PetscCall(PetscObjectGetName((PetscObject)fe, &name));
9033:       PetscCall(PetscObjectSetName((PetscObject)efe, name));
9034:       PetscCall(PetscFEGetQuadrature(fe, &q));
9035:       PetscCall(PetscFESetQuadrature(efe, q));
9036:       PetscCall(DMSetField(edm, f, NULL, (PetscObject)efe));
9037:       PetscCall(PetscFEDestroy(&efe));
9038:     }
9039:     PetscCall(DMCreateDS(edm));

9041:     PetscCall(DMCreateGlobalVector(edm, errorVec));
9042:     PetscCall(PetscObjectSetName((PetscObject)*errorVec, "Error"));
9043:     PetscCall(DMPlexComputeL2DiffVec(dm, time, exactSol, ctxs, sol, *errorVec));
9044:     PetscCall(DMDestroy(&edm));
9045:   }
9046:   PetscCall(PetscFree2(exactSol, ctxs));
9047:   PetscFunctionReturn(PETSC_SUCCESS);
9048: }

9050: /*@
9051:   DMGetNumAuxiliaryVec - Get the number of auxiliary vectors associated with this `DM`

9053:   Not Collective

9055:   Input Parameter:
9056: . dm - The `DM`

9058:   Output Parameter:
9059: . numAux - The number of auxiliary data vectors

9061:   Level: advanced

9063: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMSetAuxiliaryVec()`, `DMGetAuxiliaryLabels()`, `DMGetAuxiliaryVec()`
9064: @*/
9065: PetscErrorCode DMGetNumAuxiliaryVec(DM dm, PetscInt *numAux)
9066: {
9067:   PetscFunctionBegin;
9069:   PetscCall(PetscHMapAuxGetSize(dm->auxData, numAux));
9070:   PetscFunctionReturn(PETSC_SUCCESS);
9071: }

9073: /*@
9074:   DMGetAuxiliaryVec - Get the auxiliary vector for region specified by the given label and value, and equation part

9076:   Not Collective

9078:   Input Parameters:
9079: + dm    - The `DM`
9080: . label - The `DMLabel`
9081: . value - The label value indicating the region
9082: - part  - The equation part, or 0 if unused

9084:   Output Parameter:
9085: . aux - The `Vec` holding auxiliary field data

9087:   Level: advanced

9089:   Note:
9090:   If no auxiliary vector is found for this (label, value), (NULL, 0, 0) is checked as well.

9092: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMSetAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryLabels()`
9093: @*/
9094: PetscErrorCode DMGetAuxiliaryVec(DM dm, DMLabel label, PetscInt value, PetscInt part, Vec *aux)
9095: {
9096:   PetscHashAuxKey key, wild = {NULL, 0, 0};
9097:   PetscBool       has;

9099:   PetscFunctionBegin;
9102:   key.label = label;
9103:   key.value = value;
9104:   key.part  = part;
9105:   PetscCall(PetscHMapAuxHas(dm->auxData, key, &has));
9106:   if (has) PetscCall(PetscHMapAuxGet(dm->auxData, key, aux));
9107:   else PetscCall(PetscHMapAuxGet(dm->auxData, wild, aux));
9108:   PetscFunctionReturn(PETSC_SUCCESS);
9109: }

9111: /*@
9112:   DMSetAuxiliaryVec - Set an auxiliary vector for region specified by the given label and value, and equation part

9114:   Not Collective because auxiliary vectors are not parallel

9116:   Input Parameters:
9117: + dm    - The `DM`
9118: . label - The `DMLabel`
9119: . value - The label value indicating the region
9120: . part  - The equation part, or 0 if unused
9121: - aux   - The `Vec` holding auxiliary field data

9123:   Level: advanced

9125: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMGetAuxiliaryLabels()`, `DMCopyAuxiliaryVec()`
9126: @*/
9127: PetscErrorCode DMSetAuxiliaryVec(DM dm, DMLabel label, PetscInt value, PetscInt part, Vec aux)
9128: {
9129:   Vec             old;
9130:   PetscHashAuxKey key;

9132:   PetscFunctionBegin;
9135:   key.label = label;
9136:   key.value = value;
9137:   key.part  = part;
9138:   PetscCall(PetscHMapAuxGet(dm->auxData, key, &old));
9139:   PetscCall(PetscObjectReference((PetscObject)aux));
9140:   if (!aux) PetscCall(PetscHMapAuxDel(dm->auxData, key));
9141:   else PetscCall(PetscHMapAuxSet(dm->auxData, key, aux));
9142:   PetscCall(VecDestroy(&old));
9143:   PetscFunctionReturn(PETSC_SUCCESS);
9144: }

9146: /*@C
9147:   DMGetAuxiliaryLabels - Get the labels, values, and parts for all auxiliary vectors in this `DM`

9149:   Not Collective

9151:   Input Parameter:
9152: . dm - The `DM`

9154:   Output Parameters:
9155: + labels - The `DMLabel`s for each `Vec`
9156: . values - The label values for each `Vec`
9157: - parts  - The equation parts for each `Vec`

9159:   Level: advanced

9161:   Note:
9162:   The arrays passed in must be at least as large as `DMGetNumAuxiliaryVec()`.

9164: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMSetAuxiliaryVec()`, `DMCopyAuxiliaryVec()`
9165: @*/
9166: PetscErrorCode DMGetAuxiliaryLabels(DM dm, DMLabel labels[], PetscInt values[], PetscInt parts[])
9167: {
9168:   PetscHashAuxKey *keys;
9169:   PetscInt         n, i, off = 0;

9171:   PetscFunctionBegin;
9173:   PetscAssertPointer(labels, 2);
9174:   PetscAssertPointer(values, 3);
9175:   PetscAssertPointer(parts, 4);
9176:   PetscCall(DMGetNumAuxiliaryVec(dm, &n));
9177:   PetscCall(PetscMalloc1(n, &keys));
9178:   PetscCall(PetscHMapAuxGetKeys(dm->auxData, &off, keys));
9179:   for (i = 0; i < n; ++i) {
9180:     labels[i] = keys[i].label;
9181:     values[i] = keys[i].value;
9182:     parts[i]  = keys[i].part;
9183:   }
9184:   PetscCall(PetscFree(keys));
9185:   PetscFunctionReturn(PETSC_SUCCESS);
9186: }

9188: /*@
9189:   DMCopyAuxiliaryVec - Copy the auxiliary vector data on a `DM` to a new `DM`

9191:   Not Collective

9193:   Input Parameter:
9194: . dm - The `DM`

9196:   Output Parameter:
9197: . dmNew - The new `DM`, now with the same auxiliary data

9199:   Level: advanced

9201:   Note:
9202:   This is a shallow copy of the auxiliary vectors

9204: .seealso: [](ch_dmbase), `DM`, `DMClearAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMSetAuxiliaryVec()`
9205: @*/
9206: PetscErrorCode DMCopyAuxiliaryVec(DM dm, DM dmNew)
9207: {
9208:   PetscFunctionBegin;
9211:   if (dm == dmNew) PetscFunctionReturn(PETSC_SUCCESS);
9212:   PetscCall(DMClearAuxiliaryVec(dmNew));

9214:   PetscCall(PetscHMapAuxDestroy(&dmNew->auxData));
9215:   PetscCall(PetscHMapAuxDuplicate(dm->auxData, &dmNew->auxData));
9216:   {
9217:     Vec     *auxData;
9218:     PetscInt n, i, off = 0;

9220:     PetscCall(PetscHMapAuxGetSize(dmNew->auxData, &n));
9221:     PetscCall(PetscMalloc1(n, &auxData));
9222:     PetscCall(PetscHMapAuxGetVals(dmNew->auxData, &off, auxData));
9223:     for (i = 0; i < n; ++i) PetscCall(PetscObjectReference((PetscObject)auxData[i]));
9224:     PetscCall(PetscFree(auxData));
9225:   }
9226:   PetscFunctionReturn(PETSC_SUCCESS);
9227: }

9229: /*@
9230:   DMClearAuxiliaryVec - Destroys the auxiliary vector information and creates a new empty one

9232:   Not Collective

9234:   Input Parameter:
9235: . dm - The `DM`

9237:   Level: advanced

9239: .seealso: [](ch_dmbase), `DM`, `DMCopyAuxiliaryVec()`, `DMGetNumAuxiliaryVec()`, `DMGetAuxiliaryVec()`, `DMSetAuxiliaryVec()`
9240: @*/
9241: PetscErrorCode DMClearAuxiliaryVec(DM dm)
9242: {
9243:   Vec     *auxData;
9244:   PetscInt n, i, off = 0;

9246:   PetscFunctionBegin;
9247:   PetscCall(PetscHMapAuxGetSize(dm->auxData, &n));
9248:   PetscCall(PetscMalloc1(n, &auxData));
9249:   PetscCall(PetscHMapAuxGetVals(dm->auxData, &off, auxData));
9250:   for (i = 0; i < n; ++i) PetscCall(VecDestroy(&auxData[i]));
9251:   PetscCall(PetscFree(auxData));
9252:   PetscCall(PetscHMapAuxDestroy(&dm->auxData));
9253:   PetscCall(PetscHMapAuxCreate(&dm->auxData));
9254:   PetscFunctionReturn(PETSC_SUCCESS);
9255: }

9257: /*@C
9258:   DMPolytopeMatchOrientation - Determine an orientation (transformation) that takes the source face arrangement to the target face arrangement

9260:   Not Collective

9262:   Input Parameters:
9263: + ct         - The `DMPolytopeType`
9264: . sourceCone - The source arrangement of faces
9265: - targetCone - The target arrangement of faces

9267:   Output Parameters:
9268: + ornt  - The orientation (transformation) which will take the source arrangement to the target arrangement
9269: - found - Flag indicating that a suitable orientation was found

9271:   Level: advanced

9273:   Note:
9274:   An arrangement is a face order combined with an orientation for each face

9276:   Each orientation (transformation) is labeled with an integer from negative `DMPolytopeTypeGetNumArrangements(ct)`/2 to `DMPolytopeTypeGetNumArrangements(ct)`/2
9277:   that labels each arrangement (face ordering plus orientation for each face).

9279:   See `DMPolytopeMatchVertexOrientation()` to find a new vertex orientation that takes the source vertex arrangement to the target vertex arrangement

9281: .seealso: [](ch_dmbase), `DM`, `DMPolytopeGetOrientation()`, `DMPolytopeMatchVertexOrientation()`, `DMPolytopeGetVertexOrientation()`
9282: @*/
9283: PetscErrorCode DMPolytopeMatchOrientation(DMPolytopeType ct, const PetscInt sourceCone[], const PetscInt targetCone[], PetscInt *ornt, PetscBool *found)
9284: {
9285:   const PetscInt cS = DMPolytopeTypeGetConeSize(ct);
9286:   const PetscInt nO = DMPolytopeTypeGetNumArrangements(ct) / 2;
9287:   PetscInt       o, c;

9289:   PetscFunctionBegin;
9290:   if (!nO) {
9291:     *ornt  = 0;
9292:     *found = PETSC_TRUE;
9293:     PetscFunctionReturn(PETSC_SUCCESS);
9294:   }
9295:   for (o = -nO; o < nO; ++o) {
9296:     const PetscInt *arr = DMPolytopeTypeGetArrangement(ct, o);

9298:     for (c = 0; c < cS; ++c)
9299:       if (sourceCone[arr[c * 2]] != targetCone[c]) break;
9300:     if (c == cS) {
9301:       *ornt = o;
9302:       break;
9303:     }
9304:   }
9305:   *found = o == nO ? PETSC_FALSE : PETSC_TRUE;
9306:   PetscFunctionReturn(PETSC_SUCCESS);
9307: }

9309: /*@C
9310:   DMPolytopeGetOrientation - Determine an orientation (transformation) that takes the source face arrangement to the target face arrangement

9312:   Not Collective

9314:   Input Parameters:
9315: + ct         - The `DMPolytopeType`
9316: . sourceCone - The source arrangement of faces
9317: - targetCone - The target arrangement of faces

9319:   Output Parameter:
9320: . ornt - The orientation (transformation) which will take the source arrangement to the target arrangement

9322:   Level: advanced

9324:   Note:
9325:   This function is the same as `DMPolytopeMatchOrientation()` except it will generate an error if no suitable orientation can be found.

9327:   Developer Note:
9328:   It is unclear why this function needs to exist since one can simply call `DMPolytopeMatchOrientation()` and error if none is found

9330: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMPolytopeMatchOrientation()`, `DMPolytopeGetVertexOrientation()`, `DMPolytopeMatchVertexOrientation()`
9331: @*/
9332: PetscErrorCode DMPolytopeGetOrientation(DMPolytopeType ct, const PetscInt sourceCone[], const PetscInt targetCone[], PetscInt *ornt)
9333: {
9334:   PetscBool found;

9336:   PetscFunctionBegin;
9337:   PetscCall(DMPolytopeMatchOrientation(ct, sourceCone, targetCone, ornt, &found));
9338:   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find orientation for %s", DMPolytopeTypes[ct]);
9339:   PetscFunctionReturn(PETSC_SUCCESS);
9340: }

9342: /*@C
9343:   DMPolytopeMatchVertexOrientation - Determine an orientation (transformation) that takes the source vertex arrangement to the target vertex arrangement

9345:   Not Collective

9347:   Input Parameters:
9348: + ct         - The `DMPolytopeType`
9349: . sourceVert - The source arrangement of vertices
9350: - targetVert - The target arrangement of vertices

9352:   Output Parameters:
9353: + ornt  - The orientation (transformation) which will take the source arrangement to the target arrangement
9354: - found - Flag indicating that a suitable orientation was found

9356:   Level: advanced

9358:   Notes:
9359:   An arrangement is a vertex order

9361:   Each orientation (transformation) is labeled with an integer from negative `DMPolytopeTypeGetNumArrangements(ct)`/2 to `DMPolytopeTypeGetNumArrangements(ct)`/2
9362:   that labels each arrangement (vertex ordering).

9364:   See `DMPolytopeMatchOrientation()` to find a new face orientation that takes the source face arrangement to the target face arrangement

9366: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMPolytopeGetOrientation()`, `DMPolytopeMatchOrientation()`, `DMPolytopeTypeGetNumVertices()`, `DMPolytopeTypeGetVertexArrangement()`
9367: @*/
9368: PetscErrorCode DMPolytopeMatchVertexOrientation(DMPolytopeType ct, const PetscInt sourceVert[], const PetscInt targetVert[], PetscInt *ornt, PetscBool *found)
9369: {
9370:   const PetscInt cS = DMPolytopeTypeGetNumVertices(ct);
9371:   const PetscInt nO = DMPolytopeTypeGetNumArrangements(ct) / 2;
9372:   PetscInt       o, c;

9374:   PetscFunctionBegin;
9375:   if (!nO) {
9376:     *ornt  = 0;
9377:     *found = PETSC_TRUE;
9378:     PetscFunctionReturn(PETSC_SUCCESS);
9379:   }
9380:   for (o = -nO; o < nO; ++o) {
9381:     const PetscInt *arr = DMPolytopeTypeGetVertexArrangement(ct, o);

9383:     for (c = 0; c < cS; ++c)
9384:       if (sourceVert[arr[c]] != targetVert[c]) break;
9385:     if (c == cS) {
9386:       *ornt = o;
9387:       break;
9388:     }
9389:   }
9390:   *found = o == nO ? PETSC_FALSE : PETSC_TRUE;
9391:   PetscFunctionReturn(PETSC_SUCCESS);
9392: }

9394: /*@C
9395:   DMPolytopeGetVertexOrientation - Determine an orientation (transformation) that takes the source vertex arrangement to the target vertex arrangement

9397:   Not Collective

9399:   Input Parameters:
9400: + ct         - The `DMPolytopeType`
9401: . sourceCone - The source arrangement of vertices
9402: - targetCone - The target arrangement of vertices

9404:   Output Parameter:
9405: . ornt - The orientation (transformation) which will take the source arrangement to the target arrangement

9407:   Level: advanced

9409:   Note:
9410:   This function is the same as `DMPolytopeMatchVertexOrientation()` except it errors if not orientation is possible.

9412:   Developer Note:
9413:   It is unclear why this function needs to exist since one can simply call `DMPolytopeMatchVertexOrientation()` and error if none is found

9415: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMPolytopeMatchVertexOrientation()`, `DMPolytopeGetOrientation()`
9416: @*/
9417: PetscErrorCode DMPolytopeGetVertexOrientation(DMPolytopeType ct, const PetscInt sourceCone[], const PetscInt targetCone[], PetscInt *ornt)
9418: {
9419:   PetscBool found;

9421:   PetscFunctionBegin;
9422:   PetscCall(DMPolytopeMatchVertexOrientation(ct, sourceCone, targetCone, ornt, &found));
9423:   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find orientation for %s", DMPolytopeTypes[ct]);
9424:   PetscFunctionReturn(PETSC_SUCCESS);
9425: }

9427: /*@C
9428:   DMPolytopeInCellTest - Check whether a point lies inside the reference cell of given type

9430:   Not Collective

9432:   Input Parameters:
9433: + ct    - The `DMPolytopeType`
9434: - point - Coordinates of the point

9436:   Output Parameter:
9437: . inside - Flag indicating whether the point is inside the reference cell of given type

9439:   Level: advanced

9441: .seealso: [](ch_dmbase), `DM`, `DMPolytopeType`, `DMLocatePoints()`
9442: @*/
9443: PetscErrorCode DMPolytopeInCellTest(DMPolytopeType ct, const PetscReal point[], PetscBool *inside)
9444: {
9445:   PetscReal sum = 0.0;
9446:   PetscInt  d;

9448:   PetscFunctionBegin;
9449:   *inside = PETSC_TRUE;
9450:   switch (ct) {
9451:   case DM_POLYTOPE_TRIANGLE:
9452:   case DM_POLYTOPE_TETRAHEDRON:
9453:     for (d = 0; d < DMPolytopeTypeGetDim(ct); ++d) {
9454:       if (point[d] < -1.0) {
9455:         *inside = PETSC_FALSE;
9456:         break;
9457:       }
9458:       sum += point[d];
9459:     }
9460:     if (sum > PETSC_SMALL) {
9461:       *inside = PETSC_FALSE;
9462:       break;
9463:     }
9464:     break;
9465:   case DM_POLYTOPE_QUADRILATERAL:
9466:   case DM_POLYTOPE_HEXAHEDRON:
9467:     for (d = 0; d < DMPolytopeTypeGetDim(ct); ++d)
9468:       if (PetscAbsReal(point[d]) > 1. + PETSC_SMALL) {
9469:         *inside = PETSC_FALSE;
9470:         break;
9471:       }
9472:     break;
9473:   default:
9474:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported polytope type %s", DMPolytopeTypes[ct]);
9475:   }
9476:   PetscFunctionReturn(PETSC_SUCCESS);
9477: }

9479: /*@
9480:   DMReorderSectionSetDefault - Set flag indicating whether the local section should be reordered by default

9482:   Logically collective

9484:   Input Parameters:
9485: + dm      - The DM
9486: - reorder - Flag for reordering

9488:   Level: intermediate

9490: .seealso: `DMReorderSectionGetDefault()`
9491: @*/
9492: PetscErrorCode DMReorderSectionSetDefault(DM dm, DMReorderDefaultFlag reorder)
9493: {
9494:   PetscFunctionBegin;
9496:   PetscTryMethod(dm, "DMReorderSectionSetDefault_C", (DM, DMReorderDefaultFlag), (dm, reorder));
9497:   PetscFunctionReturn(PETSC_SUCCESS);
9498: }

9500: /*@
9501:   DMReorderSectionGetDefault - Get flag indicating whether the local section should be reordered by default

9503:   Not collective

9505:   Input Parameter:
9506: . dm - The DM

9508:   Output Parameter:
9509: . reorder - Flag for reordering

9511:   Level: intermediate

9513: .seealso: `DMReorderSetDefault()`
9514: @*/
9515: PetscErrorCode DMReorderSectionGetDefault(DM dm, DMReorderDefaultFlag *reorder)
9516: {
9517:   PetscFunctionBegin;
9519:   PetscAssertPointer(reorder, 2);
9520:   *reorder = DM_REORDER_DEFAULT_NOTSET;
9521:   PetscTryMethod(dm, "DMReorderSectionGetDefault_C", (DM, DMReorderDefaultFlag *), (dm, reorder));
9522:   PetscFunctionReturn(PETSC_SUCCESS);
9523: }

9525: /*@C
9526:   DMReorderSectionSetType - Set the type of local section reordering

9528:   Logically collective

9530:   Input Parameters:
9531: + dm      - The DM
9532: - reorder - The reordering method

9534:   Level: intermediate

9536: .seealso: `DMReorderSectionGetType()`, `DMReorderSectionSetDefault()`
9537: @*/
9538: PetscErrorCode DMReorderSectionSetType(DM dm, MatOrderingType reorder)
9539: {
9540:   PetscFunctionBegin;
9542:   PetscTryMethod(dm, "DMReorderSectionSetType_C", (DM, MatOrderingType), (dm, reorder));
9543:   PetscFunctionReturn(PETSC_SUCCESS);
9544: }

9546: /*@C
9547:   DMReorderSectionGetType - Get the reordering type for the local section

9549:   Not collective

9551:   Input Parameter:
9552: . dm - The DM

9554:   Output Parameter:
9555: . reorder - The reordering method

9557:   Level: intermediate

9559: .seealso: `DMReorderSetDefault()`, `DMReorderSectionGetDefault()`
9560: @*/
9561: PetscErrorCode DMReorderSectionGetType(DM dm, MatOrderingType *reorder)
9562: {
9563:   PetscFunctionBegin;
9565:   PetscAssertPointer(reorder, 2);
9566:   *reorder = NULL;
9567:   PetscTryMethod(dm, "DMReorderSectionGetType_C", (DM, MatOrderingType *), (dm, reorder));
9568:   PetscFunctionReturn(PETSC_SUCCESS);
9569: }