Actual source code: dmcoordinates.c

  1: #include <petsc/private/dmimpl.h>

  3: #include <petscdmplex.h>
  4: #include <petscsf.h>

  6: PetscErrorCode DMRestrictHook_Coordinates(DM dm, DM dmc, void *ctx)
  7: {
  8:   DM  dm_coord, dmc_coord;
  9:   Vec coords, ccoords;
 10:   Mat inject;

 12:   PetscFunctionBegin;
 13:   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
 14:   PetscCall(DMGetCoordinateDM(dmc, &dmc_coord));
 15:   PetscCall(DMGetCoordinates(dm, &coords));
 16:   PetscCall(DMGetCoordinates(dmc, &ccoords));
 17:   if (coords && !ccoords) {
 18:     PetscCall(DMCreateGlobalVector(dmc_coord, &ccoords));
 19:     PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates"));
 20:     PetscCall(DMCreateInjection(dmc_coord, dm_coord, &inject));
 21:     PetscCall(MatRestrict(inject, coords, ccoords));
 22:     PetscCall(MatDestroy(&inject));
 23:     PetscCall(DMSetCoordinates(dmc, ccoords));
 24:     PetscCall(VecDestroy(&ccoords));
 25:   }
 26:   PetscFunctionReturn(PETSC_SUCCESS);
 27: }

 29: static PetscErrorCode DMSubDomainHook_Coordinates(DM dm, DM subdm, void *ctx)
 30: {
 31:   DM          dm_coord, subdm_coord;
 32:   Vec         coords, ccoords, clcoords;
 33:   VecScatter *scat_i, *scat_g;

 35:   PetscFunctionBegin;
 36:   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
 37:   PetscCall(DMGetCoordinateDM(subdm, &subdm_coord));
 38:   PetscCall(DMGetCoordinates(dm, &coords));
 39:   PetscCall(DMGetCoordinates(subdm, &ccoords));
 40:   if (coords && !ccoords) {
 41:     PetscCall(DMCreateGlobalVector(subdm_coord, &ccoords));
 42:     PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates"));
 43:     PetscCall(DMCreateLocalVector(subdm_coord, &clcoords));
 44:     PetscCall(PetscObjectSetName((PetscObject)clcoords, "coordinates"));
 45:     PetscCall(DMCreateDomainDecompositionScatters(dm_coord, 1, &subdm_coord, NULL, &scat_i, &scat_g));
 46:     PetscCall(VecScatterBegin(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD));
 47:     PetscCall(VecScatterEnd(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD));
 48:     PetscCall(VecScatterBegin(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD));
 49:     PetscCall(VecScatterEnd(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD));
 50:     PetscCall(DMSetCoordinates(subdm, ccoords));
 51:     PetscCall(DMSetCoordinatesLocal(subdm, clcoords));
 52:     PetscCall(VecScatterDestroy(&scat_i[0]));
 53:     PetscCall(VecScatterDestroy(&scat_g[0]));
 54:     PetscCall(VecDestroy(&ccoords));
 55:     PetscCall(VecDestroy(&clcoords));
 56:     PetscCall(PetscFree(scat_i));
 57:     PetscCall(PetscFree(scat_g));
 58:   }
 59:   PetscFunctionReturn(PETSC_SUCCESS);
 60: }

 62: /*@
 63:   DMGetCoordinateDM - Gets the `DM` that prescribes coordinate layout and scatters between global and local coordinates

 65:   Collective

 67:   Input Parameter:
 68: . dm - the `DM`

 70:   Output Parameter:
 71: . cdm - coordinate `DM`

 73:   Level: intermediate

 75: .seealso: `DM`, `DMSetCoordinateDM()`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGSetCellCoordinateDM()`,

 77: @*/
 78: PetscErrorCode DMGetCoordinateDM(DM dm, DM *cdm)
 79: {
 80:   PetscFunctionBegin;
 82:   PetscAssertPointer(cdm, 2);
 83:   if (!dm->coordinates[0].dm) {
 84:     DM cdm;

 86:     PetscUseTypeMethod(dm, createcoordinatedm, &cdm);
 87:     PetscCall(PetscObjectSetName((PetscObject)cdm, "coordinateDM"));
 88:     /* Just in case the DM sets the coordinate DM when creating it (DMP4est can do this, because it may not setup
 89:      * until the call to CreateCoordinateDM) */
 90:     PetscCall(DMDestroy(&dm->coordinates[0].dm));
 91:     dm->coordinates[0].dm = cdm;
 92:   }
 93:   *cdm = dm->coordinates[0].dm;
 94:   PetscFunctionReturn(PETSC_SUCCESS);
 95: }

 97: /*@
 98:   DMSetCoordinateDM - Sets the `DM` that prescribes coordinate layout and scatters between global and local coordinates

100:   Logically Collective

102:   Input Parameters:
103: + dm  - the `DM`
104: - cdm - coordinate `DM`

106:   Level: intermediate

108: .seealso: `DM`, `DMGetCoordinateDM()`, `DMSetCoordinates()`, `DMGetCellCoordinateDM()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`,
109:           `DMGSetCellCoordinateDM()`
110: @*/
111: PetscErrorCode DMSetCoordinateDM(DM dm, DM cdm)
112: {
113:   PetscFunctionBegin;
116:   PetscCall(PetscObjectReference((PetscObject)cdm));
117:   PetscCall(DMDestroy(&dm->coordinates[0].dm));
118:   dm->coordinates[0].dm = cdm;
119:   PetscFunctionReturn(PETSC_SUCCESS);
120: }

122: /*@
123:   DMGetCellCoordinateDM - Gets the `DM` that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates

125:   Collective

127:   Input Parameter:
128: . dm - the `DM`

130:   Output Parameter:
131: . cdm - cellwise coordinate `DM`, or `NULL` if they are not defined

133:   Level: intermediate

135:   Note:
136:   Call `DMLocalizeCoordinates()` to automatically create cellwise coordinates for periodic geometries.

138: .seealso: `DM`, `DMSetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`,
139:           `DMLocalizeCoordinates()`, `DMSetCoordinateDM()`, `DMGetCoordinateDM()`
140: @*/
141: PetscErrorCode DMGetCellCoordinateDM(DM dm, DM *cdm)
142: {
143:   PetscFunctionBegin;
145:   PetscAssertPointer(cdm, 2);
146:   *cdm = dm->coordinates[1].dm;
147:   PetscFunctionReturn(PETSC_SUCCESS);
148: }

150: /*@
151:   DMSetCellCoordinateDM - Sets the `DM` that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates

153:   Logically Collective

155:   Input Parameters:
156: + dm  - the `DM`
157: - cdm - cellwise coordinate `DM`

159:   Level: intermediate

161:   Note:
162:   As opposed to `DMSetCoordinateDM()` these coordinates are useful for discontinuous Galerkin methods since they support coordinate fields that are discontinuous at cell boundaries.

164: .seealso: `DMGetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`,
165:           `DMSetCoordinateDM()`, `DMGetCoordinateDM()`
166: @*/
167: PetscErrorCode DMSetCellCoordinateDM(DM dm, DM cdm)
168: {
169:   PetscInt dim;

171:   PetscFunctionBegin;
173:   if (cdm) {
175:     PetscCall(DMGetCoordinateDim(dm, &dim));
176:     dm->coordinates[1].dim = dim;
177:   }
178:   PetscCall(PetscObjectReference((PetscObject)cdm));
179:   PetscCall(DMDestroy(&dm->coordinates[1].dm));
180:   dm->coordinates[1].dm = cdm;
181:   PetscFunctionReturn(PETSC_SUCCESS);
182: }

184: /*@
185:   DMGetCoordinateDim - Retrieve the dimension of the embedding space for coordinate values. For example a mesh on the surface of a sphere would have a 3 dimensional embedding space

187:   Not Collective

189:   Input Parameter:
190: . dm - The `DM` object

192:   Output Parameter:
193: . dim - The embedding dimension

195:   Level: intermediate

197: .seealso: `DM`, `DMSetCoordinateDim()`, `DMGetCoordinateSection()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
198: @*/
199: PetscErrorCode DMGetCoordinateDim(DM dm, PetscInt *dim)
200: {
201:   PetscFunctionBegin;
203:   PetscAssertPointer(dim, 2);
204:   if (dm->coordinates[0].dim == PETSC_DEFAULT) dm->coordinates[0].dim = dm->dim;
205:   *dim = dm->coordinates[0].dim;
206:   PetscFunctionReturn(PETSC_SUCCESS);
207: }

209: /*@
210:   DMSetCoordinateDim - Set the dimension of the embedding space for coordinate values.

212:   Not Collective

214:   Input Parameters:
215: + dm  - The `DM` object
216: - dim - The embedding dimension

218:   Level: intermediate

220: .seealso: `DM`, `DMGetCoordinateDim()`, `DMSetCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()`
221: @*/
222: PetscErrorCode DMSetCoordinateDim(DM dm, PetscInt dim)
223: {
224:   PetscDS  ds;
225:   PetscInt Nds, n;

227:   PetscFunctionBegin;
229:   dm->coordinates[0].dim = dim;
230:   if (dm->dim >= 0) {
231:     PetscCall(DMGetNumDS(dm, &Nds));
232:     for (n = 0; n < Nds; ++n) {
233:       PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds, NULL));
234:       PetscCall(PetscDSSetCoordinateDimension(ds, dim));
235:     }
236:   }
237:   PetscFunctionReturn(PETSC_SUCCESS);
238: }

240: /*@
241:   DMGetCoordinateSection - Retrieve the `PetscSection` of coordinate values over the mesh.

243:   Collective

245:   Input Parameter:
246: . dm - The `DM` object

248:   Output Parameter:
249: . section - The `PetscSection` object

251:   Level: intermediate

253:   Note:
254:   This just retrieves the local section from the coordinate `DM`. In other words,
255: .vb
256:   DMGetCoordinateDM(dm, &cdm);
257:   DMGetLocalSection(cdm, &section);
258: .ve

260: .seealso: `DM`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
261: @*/
262: PetscErrorCode DMGetCoordinateSection(DM dm, PetscSection *section)
263: {
264:   DM cdm;

266:   PetscFunctionBegin;
268:   PetscAssertPointer(section, 2);
269:   PetscCall(DMGetCoordinateDM(dm, &cdm));
270:   PetscCall(DMGetLocalSection(cdm, section));
271:   PetscFunctionReturn(PETSC_SUCCESS);
272: }

274: /*@
275:   DMSetCoordinateSection - Set the `PetscSection` of coordinate values over the mesh.

277:   Not Collective

279:   Input Parameters:
280: + dm      - The `DM` object
281: . dim     - The embedding dimension, or `PETSC_DETERMINE`
282: - section - The `PetscSection` object

284:   Level: intermediate

286: .seealso: `DM`, `DMGetCoordinateDim()`, `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()`
287: @*/
288: PetscErrorCode DMSetCoordinateSection(DM dm, PetscInt dim, PetscSection section)
289: {
290:   DM cdm;

292:   PetscFunctionBegin;
295:   PetscCall(DMGetCoordinateDM(dm, &cdm));
296:   PetscCall(DMSetLocalSection(cdm, section));
297:   if (dim == PETSC_DETERMINE) {
298:     PetscInt d = PETSC_DEFAULT;
299:     PetscInt pStart, pEnd, vStart, vEnd, v, dd;

301:     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
302:     PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd));
303:     pStart = PetscMax(vStart, pStart);
304:     pEnd   = PetscMin(vEnd, pEnd);
305:     for (v = pStart; v < pEnd; ++v) {
306:       PetscCall(PetscSectionGetDof(section, v, &dd));
307:       if (dd) {
308:         d = dd;
309:         break;
310:       }
311:     }
312:     if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d));
313:   }
314:   PetscFunctionReturn(PETSC_SUCCESS);
315: }

317: /*@
318:   DMGetCellCoordinateSection - Retrieve the `PetscSection` of cellwise coordinate values over the mesh.

320:   Collective

322:   Input Parameter:
323: . dm - The `DM` object

325:   Output Parameter:
326: . section - The `PetscSection` object, or `NULL` if no cellwise coordinates are defined

328:   Level: intermediate

330:   Note:
331:   This just retrieves the local section from the cell coordinate `DM`. In other words,
332: .vb
333:   DMGetCellCoordinateDM(dm, &cdm);
334:   DMGetLocalSection(cdm, &section);
335: .ve

337: .seealso: `DM`, `DMGetCoordinateSection()`, `DMSetCellCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
338: @*/
339: PetscErrorCode DMGetCellCoordinateSection(DM dm, PetscSection *section)
340: {
341:   DM cdm;

343:   PetscFunctionBegin;
345:   PetscAssertPointer(section, 2);
346:   *section = NULL;
347:   PetscCall(DMGetCellCoordinateDM(dm, &cdm));
348:   if (cdm) PetscCall(DMGetLocalSection(cdm, section));
349:   PetscFunctionReturn(PETSC_SUCCESS);
350: }

352: /*@
353:   DMSetCellCoordinateSection - Set the `PetscSection` of cellwise coordinate values over the mesh.

355:   Not Collective

357:   Input Parameters:
358: + dm      - The `DM` object
359: . dim     - The embedding dimension, or `PETSC_DETERMINE`
360: - section - The `PetscSection` object for a cellwise layout

362:   Level: intermediate

364: .seealso: `DM`, `DMGetCoordinateDim()`, `DMSetCoordinateSection()`, `DMGetCellCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
365: @*/
366: PetscErrorCode DMSetCellCoordinateSection(DM dm, PetscInt dim, PetscSection section)
367: {
368:   DM cdm;

370:   PetscFunctionBegin;
373:   PetscCall(DMGetCellCoordinateDM(dm, &cdm));
374:   PetscCheck(cdm, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "No DM defined for cellwise coordinates");
375:   PetscCall(DMSetLocalSection(cdm, section));
376:   if (dim == PETSC_DETERMINE) {
377:     PetscInt d = PETSC_DEFAULT;
378:     PetscInt pStart, pEnd, vStart, vEnd, v, dd;

380:     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
381:     PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd));
382:     pStart = PetscMax(vStart, pStart);
383:     pEnd   = PetscMin(vEnd, pEnd);
384:     for (v = pStart; v < pEnd; ++v) {
385:       PetscCall(PetscSectionGetDof(section, v, &dd));
386:       if (dd) {
387:         d = dd;
388:         break;
389:       }
390:     }
391:     if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d));
392:   }
393:   PetscFunctionReturn(PETSC_SUCCESS);
394: }

396: /*@
397:   DMGetCoordinates - Gets a global vector with the coordinates associated with the `DM`.

399:   Collective

401:   Input Parameter:
402: . dm - the `DM`

404:   Output Parameter:
405: . c - global coordinate vector

407:   Level: intermediate

409:   Notes:
410:   This is a borrowed reference, so the user should NOT destroy this vector. When the `DM` is
411:   destroyed `c` will no longer be valid.

413:   Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates).

415:   For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
416:   and (x_0,y_0,z_0,x_1,y_1,z_1...)

418: .seealso: `DM`, `DMDA`, `DMSetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()`
419: @*/
420: PetscErrorCode DMGetCoordinates(DM dm, Vec *c)
421: {
422:   PetscFunctionBegin;
424:   PetscAssertPointer(c, 2);
425:   if (!dm->coordinates[0].x && dm->coordinates[0].xl) {
426:     DM cdm = NULL;

428:     PetscCall(DMGetCoordinateDM(dm, &cdm));
429:     PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[0].x));
430:     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].x, "coordinates"));
431:     PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x));
432:     PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x));
433:   }
434:   *c = dm->coordinates[0].x;
435:   PetscFunctionReturn(PETSC_SUCCESS);
436: }

438: /*@
439:   DMSetCoordinates - Sets into the `DM` a global vector that holds the coordinates

441:   Collective

443:   Input Parameters:
444: + dm - the `DM`
445: - c  - coordinate vector

447:   Level: intermediate

449:   Notes:
450:   The coordinates do not include those for ghost points, which are in the local vector.

452:   The vector `c` can be destroyed after the call

454: .seealso: `DM`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()`
455: @*/
456: PetscErrorCode DMSetCoordinates(DM dm, Vec c)
457: {
458:   PetscFunctionBegin;
461:   PetscCall(PetscObjectReference((PetscObject)c));
462:   PetscCall(VecDestroy(&dm->coordinates[0].x));
463:   dm->coordinates[0].x = c;
464:   PetscCall(VecDestroy(&dm->coordinates[0].xl));
465:   PetscCall(DMCoarsenHookAdd(dm, DMRestrictHook_Coordinates, NULL, NULL));
466:   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_Coordinates, NULL, NULL));
467:   PetscFunctionReturn(PETSC_SUCCESS);
468: }

470: /*@
471:   DMGetCellCoordinates - Gets a global vector with the cellwise coordinates associated with the `DM`.

473:   Collective

475:   Input Parameter:
476: . dm - the `DM`

478:   Output Parameter:
479: . c - global coordinate vector

481:   Level: intermediate

483:   Notes:
484:   This is a borrowed reference, so the user should NOT destroy this vector. When the `DM` is
485:   destroyed `c` will no longer be valid.

487:   Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates).

489: .seealso: `DM`, `DMGetCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()`
490: @*/
491: PetscErrorCode DMGetCellCoordinates(DM dm, Vec *c)
492: {
493:   PetscFunctionBegin;
495:   PetscAssertPointer(c, 2);
496:   if (!dm->coordinates[1].x && dm->coordinates[1].xl) {
497:     DM cdm = NULL;

499:     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
500:     PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[1].x));
501:     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].x, "DG coordinates"));
502:     PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x));
503:     PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x));
504:   }
505:   *c = dm->coordinates[1].x;
506:   PetscFunctionReturn(PETSC_SUCCESS);
507: }

509: /*@
510:   DMSetCellCoordinates - Sets into the `DM` a global vector that holds the cellwise coordinates

512:   Collective

514:   Input Parameters:
515: + dm - the `DM`
516: - c  - cellwise coordinate vector

518:   Level: intermediate

520:   Notes:
521:   The coordinates do not include those for ghost points, which are in the local vector.

523:   The vector `c` should be destroyed by the caller.

525: .seealso: `DM`, `DMGetCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()`
526: @*/
527: PetscErrorCode DMSetCellCoordinates(DM dm, Vec c)
528: {
529:   PetscFunctionBegin;
532:   PetscCall(PetscObjectReference((PetscObject)c));
533:   PetscCall(VecDestroy(&dm->coordinates[1].x));
534:   dm->coordinates[1].x = c;
535:   PetscCall(VecDestroy(&dm->coordinates[1].xl));
536:   PetscFunctionReturn(PETSC_SUCCESS);
537: }

539: /*@
540:   DMGetCoordinatesLocalSetUp - Prepares a local vector of coordinates, so that `DMGetCoordinatesLocalNoncollective()` can be used as non-collective afterwards.

542:   Collective

544:   Input Parameter:
545: . dm - the `DM`

547:   Level: advanced

549: .seealso: `DM`, `DMSetCoordinates()`, `DMGetCoordinatesLocalNoncollective()`
550: @*/
551: PetscErrorCode DMGetCoordinatesLocalSetUp(DM dm)
552: {
553:   PetscFunctionBegin;
555:   if (!dm->coordinates[0].xl && dm->coordinates[0].x) {
556:     DM       cdm = NULL;
557:     PetscInt bs;

559:     PetscCall(DMGetCoordinateDM(dm, &cdm));
560:     PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[0].xl));
561:     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].xl, "Local Coordinates"));
562:     // If the size of the vector is 0, it will not get the right block size
563:     PetscCall(VecGetBlockSize(dm->coordinates[0].x, &bs));
564:     PetscCall(VecSetBlockSize(dm->coordinates[0].xl, bs));
565:     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].xl, "coordinates"));
566:     PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
567:     PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
568:   }
569:   PetscFunctionReturn(PETSC_SUCCESS);
570: }

572: /*@
573:   DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the `DM`.

575:   Collective the first time it is called

577:   Input Parameter:
578: . dm - the `DM`

580:   Output Parameter:
581: . c - coordinate vector

583:   Level: intermediate

585:   Notes:
586:   This is a borrowed reference, so the user should NOT destroy `c`

588:   Each process has the local and ghost coordinates

590:   For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
591:   and (x_0,y_0,z_0,x_1,y_1,z_1...)

593: .seealso: `DM`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`, `DMGetCoordinatesLocalNoncollective()`
594: @*/
595: PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c)
596: {
597:   PetscFunctionBegin;
599:   PetscAssertPointer(c, 2);
600:   PetscCall(DMGetCoordinatesLocalSetUp(dm));
601:   *c = dm->coordinates[0].xl;
602:   PetscFunctionReturn(PETSC_SUCCESS);
603: }

605: /*@
606:   DMGetCoordinatesLocalNoncollective - Non-collective version of `DMGetCoordinatesLocal()`. Fails if global coordinates have been set and `DMGetCoordinatesLocalSetUp()` not called.

608:   Not Collective

610:   Input Parameter:
611: . dm - the `DM`

613:   Output Parameter:
614: . c - coordinate vector

616:   Level: advanced

618:   Note:
619:   A previous call to  `DMGetCoordinatesLocal()` or `DMGetCoordinatesLocalSetUp()` ensures that a call to this function will not error.

621: .seealso: `DM`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinatesLocal()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
622: @*/
623: PetscErrorCode DMGetCoordinatesLocalNoncollective(DM dm, Vec *c)
624: {
625:   PetscFunctionBegin;
627:   PetscAssertPointer(c, 2);
628:   PetscCheck(dm->coordinates[0].xl || !dm->coordinates[0].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called");
629:   *c = dm->coordinates[0].xl;
630:   PetscFunctionReturn(PETSC_SUCCESS);
631: }

633: /*@
634:   DMGetCoordinatesLocalTuple - Gets a local vector with the coordinates of specified points and the section describing its layout.

636:   Not Collective

638:   Input Parameters:
639: + dm - the `DM`
640: - p  - the `IS` of points whose coordinates will be returned

642:   Output Parameters:
643: + pCoordSection - the `PetscSection` describing the layout of pCoord, i.e. each point corresponds to one point in `p`, and DOFs correspond to coordinates
644: - pCoord        - the `Vec` with coordinates of points in `p`

646:   Level: advanced

648:   Notes:
649:   `DMGetCoordinatesLocalSetUp()` must be called first. This function employs `DMGetCoordinatesLocalNoncollective()` so it is not collective.

651:   This creates a new vector, so the user SHOULD destroy this vector

653:   Each process has the local and ghost coordinates

655:   For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
656:   and (x_0,y_0,z_0,x_1,y_1,z_1...)

658: .seealso: `DM`, `DMDA`, `DMSetCoordinatesLocal()`, `DMGetCoordinatesLocal()`, `DMGetCoordinatesLocalNoncollective()`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
659: @*/
660: PetscErrorCode DMGetCoordinatesLocalTuple(DM dm, IS p, PetscSection *pCoordSection, Vec *pCoord)
661: {
662:   DM                 cdm;
663:   PetscSection       cs, newcs;
664:   Vec                coords;
665:   const PetscScalar *arr;
666:   PetscScalar       *newarr = NULL;
667:   PetscInt           n;

669:   PetscFunctionBegin;
672:   if (pCoordSection) PetscAssertPointer(pCoordSection, 3);
673:   if (pCoord) PetscAssertPointer(pCoord, 4);
674:   PetscCall(DMGetCoordinateDM(dm, &cdm));
675:   PetscCall(DMGetLocalSection(cdm, &cs));
676:   PetscCall(DMGetCoordinatesLocal(dm, &coords));
677:   PetscCheck(coords, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called or coordinates not set");
678:   PetscCheck(cdm && cs, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM not supported");
679:   PetscCall(VecGetArrayRead(coords, &arr));
680:   PetscCall(PetscSectionExtractDofsFromArray(cs, MPIU_SCALAR, arr, p, &newcs, pCoord ? ((void **)&newarr) : NULL));
681:   PetscCall(VecRestoreArrayRead(coords, &arr));
682:   if (pCoord) {
683:     PetscCall(PetscSectionGetStorageSize(newcs, &n));
684:     /* set array in two steps to mimic PETSC_OWN_POINTER */
685:     PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)p), 1, n, NULL, pCoord));
686:     PetscCall(VecReplaceArray(*pCoord, newarr));
687:   } else {
688:     PetscCall(PetscFree(newarr));
689:   }
690:   if (pCoordSection) {
691:     *pCoordSection = newcs;
692:   } else PetscCall(PetscSectionDestroy(&newcs));
693:   PetscFunctionReturn(PETSC_SUCCESS);
694: }

696: /*@
697:   DMSetCoordinatesLocal - Sets into the `DM` a local vector, including ghost points, that holds the coordinates

699:   Not Collective

701:   Input Parameters:
702: + dm - the `DM`
703: - c  - coordinate vector

705:   Level: intermediate

707:   Notes:
708:   The coordinates of ghost points can be set using `DMSetCoordinates()`
709:   followed by `DMGetCoordinatesLocal()`. This is intended to enable the
710:   setting of ghost coordinates outside of the domain.

712:   The vector `c` should be destroyed by the caller.

714: .seealso: `DM`, `DMGetCoordinatesLocal()`, `DMSetCoordinates()`, `DMGetCoordinates()`, `DMGetCoordinateDM()`
715: @*/
716: PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c)
717: {
718:   PetscFunctionBegin;
721:   PetscCall(PetscObjectReference((PetscObject)c));
722:   PetscCall(VecDestroy(&dm->coordinates[0].xl));
723:   dm->coordinates[0].xl = c;
724:   PetscCall(VecDestroy(&dm->coordinates[0].x));
725:   PetscFunctionReturn(PETSC_SUCCESS);
726: }

728: /*@
729:   DMGetCellCoordinatesLocalSetUp - Prepares a local vector of cellwise coordinates, so that `DMGetCellCoordinatesLocalNoncollective()` can be used as non-collective afterwards.

731:   Collective

733:   Input Parameter:
734: . dm - the `DM`

736:   Level: advanced

738: .seealso: `DM`, `DMGetCellCoordinatesLocalNoncollective()`
739: @*/
740: PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM dm)
741: {
742:   PetscFunctionBegin;
744:   if (!dm->coordinates[1].xl && dm->coordinates[1].x) {
745:     DM cdm = NULL;

747:     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
748:     PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[1].xl));
749:     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].xl, "DG coordinates"));
750:     PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
751:     PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
752:   }
753:   PetscFunctionReturn(PETSC_SUCCESS);
754: }

756: /*@
757:   DMGetCellCoordinatesLocal - Gets a local vector with the cellwise coordinates associated with the `DM`.

759:   Collective

761:   Input Parameter:
762: . dm - the `DM`

764:   Output Parameter:
765: . c - coordinate vector

767:   Level: intermediate

769:   Notes:
770:   This is a borrowed reference, so the user should NOT destroy this vector

772:   Each process has the local and ghost coordinates

774: .seealso: `DM`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`, `DMGetCellCoordinatesLocalNoncollective()`
775: @*/
776: PetscErrorCode DMGetCellCoordinatesLocal(DM dm, Vec *c)
777: {
778:   PetscFunctionBegin;
780:   PetscAssertPointer(c, 2);
781:   PetscCall(DMGetCellCoordinatesLocalSetUp(dm));
782:   *c = dm->coordinates[1].xl;
783:   PetscFunctionReturn(PETSC_SUCCESS);
784: }

786: /*@
787:   DMGetCellCoordinatesLocalNoncollective - Non-collective version of `DMGetCellCoordinatesLocal()`. Fails if global cellwise coordinates have been set and `DMGetCellCoordinatesLocalSetUp()` not called.

789:   Not Collective

791:   Input Parameter:
792: . dm - the `DM`

794:   Output Parameter:
795: . c - cellwise coordinate vector

797:   Level: advanced

799: .seealso: `DM`, `DMGetCellCoordinatesLocalSetUp()`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`
800: @*/
801: PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM dm, Vec *c)
802: {
803:   PetscFunctionBegin;
805:   PetscAssertPointer(c, 2);
806:   PetscCheck(dm->coordinates[1].xl || !dm->coordinates[1].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCellCoordinatesLocalSetUp() has not been called");
807:   *c = dm->coordinates[1].xl;
808:   PetscFunctionReturn(PETSC_SUCCESS);
809: }

811: /*@
812:   DMSetCellCoordinatesLocal - Sets into the `DM` a local vector including ghost points that holds the cellwise coordinates

814:   Not Collective

816:   Input Parameters:
817: + dm - the `DM`
818: - c  - cellwise coordinate vector

820:   Level: intermediate

822:   Notes:
823:   The coordinates of ghost points can be set using `DMSetCoordinates()`
824:   followed by `DMGetCoordinatesLocal()`. This is intended to enable the
825:   setting of ghost coordinates outside of the domain.

827:   The vector `c` should be destroyed by the caller.

829: .seealso: `DM`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinates()`, `DMGetCellCoordinates()`, `DMGetCellCoordinateDM()`
830: @*/
831: PetscErrorCode DMSetCellCoordinatesLocal(DM dm, Vec c)
832: {
833:   PetscFunctionBegin;
836:   PetscCall(PetscObjectReference((PetscObject)c));
837:   PetscCall(VecDestroy(&dm->coordinates[1].xl));
838:   dm->coordinates[1].xl = c;
839:   PetscCall(VecDestroy(&dm->coordinates[1].x));
840:   PetscFunctionReturn(PETSC_SUCCESS);
841: }

843: PetscErrorCode DMGetCoordinateField(DM dm, DMField *field)
844: {
845:   PetscFunctionBegin;
847:   PetscAssertPointer(field, 2);
848:   if (!dm->coordinates[0].field) PetscTryTypeMethod(dm, createcoordinatefield, &dm->coordinates[0].field);
849:   *field = dm->coordinates[0].field;
850:   PetscFunctionReturn(PETSC_SUCCESS);
851: }

853: PetscErrorCode DMSetCoordinateField(DM dm, DMField field)
854: {
855:   PetscFunctionBegin;
858:   PetscCall(PetscObjectReference((PetscObject)field));
859:   PetscCall(DMFieldDestroy(&dm->coordinates[0].field));
860:   dm->coordinates[0].field = field;
861:   PetscFunctionReturn(PETSC_SUCCESS);
862: }

864: PetscErrorCode DMGetLocalBoundingBox_Coordinates(DM dm, PetscReal lmin[], PetscReal lmax[], PetscInt cs[], PetscInt ce[])
865: {
866:   Vec         coords = NULL;
867:   PetscReal   min[3] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL};
868:   PetscReal   max[3] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
869:   PetscInt    cdim, i, j;
870:   PetscMPIInt size;

872:   PetscFunctionBegin;
873:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
874:   PetscCall(DMGetCoordinateDim(dm, &cdim));
875:   if (size == 1) {
876:     const PetscReal *L, *Lstart;

878:     PetscCall(DMGetPeriodicity(dm, NULL, &Lstart, &L));
879:     if (L) {
880:       for (PetscInt d = 0; d < cdim; ++d)
881:         if (L[d] > 0.0) {
882:           min[d] = Lstart[d];
883:           max[d] = Lstart[d] + L[d];
884:         }
885:     }
886:   }
887:   PetscCall(DMGetCoordinatesLocal(dm, &coords));
888:   if (coords) {
889:     const PetscScalar *local_coords;
890:     PetscInt           N, Ni;

892:     for (j = cdim; j < 3; ++j) {
893:       min[j] = 0;
894:       max[j] = 0;
895:     }
896:     PetscCall(VecGetArrayRead(coords, &local_coords));
897:     PetscCall(VecGetLocalSize(coords, &N));
898:     Ni = N / cdim;
899:     for (i = 0; i < Ni; ++i) {
900:       for (j = 0; j < cdim; ++j) {
901:         min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j]));
902:         max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j]));
903:       }
904:     }
905:     PetscCall(VecRestoreArrayRead(coords, &local_coords));
906:     PetscCall(DMGetCellCoordinatesLocal(dm, &coords));
907:     if (coords) {
908:       PetscCall(VecGetArrayRead(coords, &local_coords));
909:       PetscCall(VecGetLocalSize(coords, &N));
910:       Ni = N / cdim;
911:       for (i = 0; i < Ni; ++i) {
912:         for (j = 0; j < cdim; ++j) {
913:           min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j]));
914:           max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j]));
915:         }
916:       }
917:       PetscCall(VecRestoreArrayRead(coords, &local_coords));
918:     }
919:     if (lmin) PetscCall(PetscArraycpy(lmin, min, cdim));
920:     if (lmax) PetscCall(PetscArraycpy(lmax, max, cdim));
921:   }
922:   PetscFunctionReturn(PETSC_SUCCESS);
923: }

925: /*@
926:   DMGetLocalBoundingBox - Returns the bounding box for the piece of the `DM` on this process.

928:   Not Collective

930:   Input Parameter:
931: . dm - the `DM`

933:   Output Parameters:
934: + lmin - local minimum coordinates (length coord dim, optional)
935: - lmax - local maximum coordinates (length coord dim, optional)

937:   Level: beginner

939:   Note:
940:   If the `DM` is a `DMDA` and has no coordinates, the index bounds are returned instead.

942: .seealso: `DM`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetBoundingBox()`
943: @*/
944: PetscErrorCode DMGetLocalBoundingBox(DM dm, PetscReal lmin[], PetscReal lmax[])
945: {
946:   PetscFunctionBegin;
948:   PetscUseTypeMethod(dm, getlocalboundingbox, lmin, lmax, NULL, NULL);
949:   PetscFunctionReturn(PETSC_SUCCESS);
950: }

952: /*@
953:   DMGetBoundingBox - Returns the global bounding box for the `DM`.

955:   Collective

957:   Input Parameter:
958: . dm - the `DM`

960:   Output Parameters:
961: + gmin - global minimum coordinates (length coord dim, optional)
962: - gmax - global maximum coordinates (length coord dim, optional)

964:   Level: beginner

966: .seealso: `DM`, `DMGetLocalBoundingBox()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`
967: @*/
968: PetscErrorCode DMGetBoundingBox(DM dm, PetscReal gmin[], PetscReal gmax[])
969: {
970:   PetscReal        lmin[3], lmax[3];
971:   const PetscReal *L, *Lstart;
972:   PetscInt         cdim;
973:   PetscMPIInt      count;

975:   PetscFunctionBegin;
977:   PetscCall(DMGetCoordinateDim(dm, &cdim));
978:   PetscCall(PetscMPIIntCast(cdim, &count));
979:   PetscCall(DMGetLocalBoundingBox(dm, lmin, lmax));
980:   if (gmin) PetscCallMPI(MPIU_Allreduce(lmin, gmin, count, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm)));
981:   if (gmax) PetscCallMPI(MPIU_Allreduce(lmax, gmax, count, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)dm)));
982:   PetscCall(DMGetPeriodicity(dm, NULL, &Lstart, &L));
983:   if (L) {
984:     for (PetscInt d = 0; d < cdim; ++d)
985:       if (L[d] > 0.0) {
986:         gmin[d] = Lstart[d];
987:         gmax[d] = Lstart[d] + L[d];
988:       }
989:   }
990:   PetscFunctionReturn(PETSC_SUCCESS);
991: }

993: static PetscErrorCode DMCreateAffineCoordinates_Internal(DM dm)
994: {
995:   DM             cdm;
996:   PetscFE        feLinear;
997:   DMPolytopeType ct;
998:   PetscInt       dim, dE, height, cStart, cEnd, gct;

1000:   PetscFunctionBegin;
1001:   PetscCall(DMGetCoordinateDM(dm, &cdm));
1002:   PetscCall(DMGetDimension(dm, &dim));
1003:   PetscCall(DMGetCoordinateDim(dm, &dE));
1004:   PetscCall(DMPlexGetVTKCellHeight(dm, &height));
1005:   PetscCall(DMPlexGetHeightStratum(dm, height, &cStart, &cEnd));
1006:   if (cEnd > cStart) PetscCall(DMPlexGetCellType(dm, cStart, &ct));
1007:   else ct = DM_POLYTOPE_UNKNOWN;
1008:   gct = (PetscInt)ct;
1009:   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &gct, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm)));
1010:   ct = (DMPolytopeType)gct;
1011:   // Work around current bug in PetscDualSpaceSetUp_Lagrange()
1012:   //   Can be seen in plex_tutorials-ex10_1
1013:   if (ct != DM_POLYTOPE_SEG_PRISM_TENSOR && ct != DM_POLYTOPE_TRI_PRISM_TENSOR && ct != DM_POLYTOPE_QUAD_PRISM_TENSOR) {
1014:     PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, dE, ct, 1, -1, &feLinear));
1015:     PetscCall(DMSetField(cdm, 0, NULL, (PetscObject)feLinear));
1016:     PetscCall(PetscFEDestroy(&feLinear));
1017:     PetscCall(DMCreateDS(cdm));
1018:   }
1019:   PetscFunctionReturn(PETSC_SUCCESS);
1020: }

1022: PetscErrorCode DMGetCoordinateDegree_Internal(DM dm, PetscInt *degree)
1023: {
1024:   DM           cdm;
1025:   PetscFE      fe;
1026:   PetscSpace   sp;
1027:   PetscClassId id;

1029:   PetscFunctionBegin;
1030:   *degree = 1;
1031:   PetscCall(DMGetCoordinateDM(dm, &cdm));
1032:   PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe));
1033:   PetscCall(PetscObjectGetClassId((PetscObject)fe, &id));
1034:   if (id != PETSCFE_CLASSID) PetscFunctionReturn(PETSC_SUCCESS);
1035:   PetscCall(PetscFEGetBasisSpace(fe, &sp));
1036:   PetscCall(PetscSpaceGetDegree(sp, degree, NULL));
1037:   PetscFunctionReturn(PETSC_SUCCESS);
1038: }

1040: /*@
1041:   DMSetCoordinateDisc - Set a coordinate space

1043:   Input Parameters:
1044: + dm      - The `DM` object
1045: . disc    - The new coordinate discretization or `NULL` to ensure a coordinate discretization exists
1046: - project - Project coordinates to new discretization

1048:   Level: intermediate

1050:   Notes:
1051:   A `PetscFE` defines an approximation space using a `PetscSpace`, which represents the basis functions, and a `PetscDualSpace`, which defines the interpolation operation in the space.

1053:   This function takes the current mesh coordinates, which are discretized using some `PetscFE` space, and projects this function into a new `PetscFE` space.
1054:   The coordinate projection is done on the continuous coordinates, but the discontinuous coordinates are not updated.

1056:   Developer Note:
1057:   With more effort, we could directly project the discontinuous coordinates also.

1059: .seealso: `DM`, `PetscFE`, `DMGetCoordinateField()`
1060: @*/
1061: PetscErrorCode DMSetCoordinateDisc(DM dm, PetscFE disc, PetscBool project)
1062: {
1063:   DM           cdmOld, cdmNew;
1064:   PetscFE      discOld;
1065:   PetscClassId classid;
1066:   Vec          coordsOld, coordsNew;
1067:   PetscBool    same_space = PETSC_TRUE;
1068:   const char  *prefix;

1070:   PetscFunctionBegin;

1074:   PetscCall(DMGetCoordinateDM(dm, &cdmOld));
1075:   /* Check current discretization is compatible */
1076:   PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
1077:   PetscCall(PetscObjectGetClassId((PetscObject)discOld, &classid));
1078:   if (classid != PETSCFE_CLASSID) {
1079:     if (classid == PETSC_CONTAINER_CLASSID) {
1080:       PetscCall(DMCreateAffineCoordinates_Internal(dm));
1081:       PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
1082:     } else {
1083:       const char *discname;

1085:       PetscCall(PetscObjectGetType((PetscObject)discOld, &discname));
1086:       SETERRQ(PetscObjectComm((PetscObject)discOld), PETSC_ERR_SUP, "Discretization type %s not supported", discname);
1087:     }
1088:   }
1089:   // Linear space has been created by now
1090:   if (!disc) PetscFunctionReturn(PETSC_SUCCESS);
1091:   // Check if the new space is the same as the old modulo quadrature
1092:   {
1093:     PetscDualSpace dsOld, ds;
1094:     PetscCall(PetscFEGetDualSpace(discOld, &dsOld));
1095:     PetscCall(PetscFEGetDualSpace(disc, &ds));
1096:     PetscCall(PetscDualSpaceEqual(dsOld, ds, &same_space));
1097:   }
1098:   // Make a fresh clone of the coordinate DM
1099:   PetscCall(DMClone(cdmOld, &cdmNew));
1100:   cdmNew->cloneOpts = PETSC_TRUE;
1101:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)cdmOld, &prefix));
1102:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)cdmNew, prefix));
1103:   PetscCall(DMSetField(cdmNew, 0, NULL, (PetscObject)disc));
1104:   PetscCall(DMCreateDS(cdmNew));
1105:   {
1106:     PetscDS ds, nds;

1108:     PetscCall(DMGetDS(cdmOld, &ds));
1109:     PetscCall(DMGetDS(cdmNew, &nds));
1110:     PetscCall(PetscDSCopyConstants(ds, nds));
1111:   }
1112:   if (cdmOld->periodic.setup) {
1113:     cdmNew->periodic.setup = cdmOld->periodic.setup;
1114:     PetscCall(cdmNew->periodic.setup(cdmNew));
1115:   }
1116:   if (dm->setfromoptionscalled) PetscCall(DMSetFromOptions(cdmNew));
1117:   if (project) {
1118:     PetscCall(DMGetCoordinates(dm, &coordsOld));
1119:     PetscCall(DMCreateGlobalVector(cdmNew, &coordsNew));
1120:     if (same_space) {
1121:       // Need to copy so that the new vector has the right dm
1122:       PetscCall(VecCopy(coordsOld, coordsNew));
1123:     } else {
1124:       Mat In;

1126:       PetscCall(DMCreateInterpolation(cdmOld, cdmNew, &In, NULL));
1127:       PetscCall(MatMult(In, coordsOld, coordsNew));
1128:       PetscCall(MatDestroy(&In));
1129:     }
1130:     PetscCall(DMSetCoordinates(dm, coordsNew));
1131:     PetscCall(VecDestroy(&coordsNew));
1132:   }
1133:   /* Set new coordinate structures */
1134:   PetscCall(DMSetCoordinateField(dm, NULL));
1135:   PetscCall(DMSetCoordinateDM(dm, cdmNew));
1136:   PetscCall(DMDestroy(&cdmNew));
1137:   PetscFunctionReturn(PETSC_SUCCESS);
1138: }

1140: /*@
1141:   DMLocatePoints - Locate the points in `v` in the mesh and return a `PetscSF` of the containing cells

1143:   Collective

1145:   Input Parameters:
1146: + dm    - The `DM`
1147: - ltype - The type of point location, e.g. `DM_POINTLOCATION_NONE` or `DM_POINTLOCATION_NEAREST`

1149:   Input/Output Parameters:
1150: + v      - The `Vec` of points, on output contains the nearest mesh points to the given points if `DM_POINTLOCATION_NEAREST` is used
1151: - cellSF - Points to either `NULL`, or a `PetscSF` with guesses for which cells contain each point;
1152:            on output, the `PetscSF` containing the MPI ranks and local indices of the containing points

1154:   Level: developer

1156:   Notes:
1157:   To do a search of the local cells of the mesh, `v` should have `PETSC_COMM_SELF` as its communicator.
1158:   To do a search of all the cells in the distributed mesh, `v` should have the same MPI communicator as `dm`.

1160:   Points will only be located in owned cells, not overlap cells arising from `DMPlexDistribute()` or other overlapping distributions.

1162:   If *cellSF is `NULL` on input, a `PetscSF` will be created.
1163:   If *cellSF is not `NULL` on input, it should point to an existing `PetscSF`, whose graph will be used as initial guesses.

1165:   An array that maps each point to its containing cell can be obtained with
1166: .vb
1167:     const PetscSFNode *cells;
1168:     PetscInt           nFound;
1169:     const PetscInt    *found;

1171:     PetscSFGetGraph(cellSF,NULL,&nFound,&found,&cells);
1172: .ve

1174:   Where cells[i].rank is the MPI rank of the process owning the cell containing point found[i] (or i if found == NULL), and cells[i].index is
1175:   the index of the cell in its MPI process' local numbering. This rank is in the communicator for `v`, so if `v` is on `PETSC_COMM_SELF` then the rank will always be 0.

1177: .seealso: `DM`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMPointLocationType`
1178: @*/
1179: PetscErrorCode DMLocatePoints(DM dm, Vec v, DMPointLocationType ltype, PetscSF *cellSF)
1180: {
1181:   PetscFunctionBegin;
1184:   PetscAssertPointer(cellSF, 4);
1185:   if (*cellSF) {
1186:     PetscMPIInt result;

1189:     PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)v), PetscObjectComm((PetscObject)*cellSF), &result));
1190:     PetscCheck(result == MPI_IDENT || result == MPI_CONGRUENT, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "cellSF must have a communicator congruent to v's");
1191:   } else {
1192:     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)v), cellSF));
1193:   }
1194:   PetscCall(PetscLogEventBegin(DM_LocatePoints, dm, 0, 0, 0));
1195:   PetscUseTypeMethod(dm, locatepoints, v, ltype, *cellSF);
1196:   PetscCall(PetscLogEventEnd(DM_LocatePoints, dm, 0, 0, 0));
1197:   PetscFunctionReturn(PETSC_SUCCESS);
1198: }