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, §ion);
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, §ion);
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: // If the size of the vector is 0, it will not get the right block size
562: PetscCall(VecGetBlockSize(dm->coordinates[0].x, &bs));
563: PetscCall(VecSetBlockSize(dm->coordinates[0].xl, bs));
564: PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].xl, "coordinates"));
565: PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
566: PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
567: }
568: PetscFunctionReturn(PETSC_SUCCESS);
569: }
571: /*@
572: DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the `DM`.
574: Collective the first time it is called
576: Input Parameter:
577: . dm - the `DM`
579: Output Parameter:
580: . c - coordinate vector
582: Level: intermediate
584: Notes:
585: This is a borrowed reference, so the user should NOT destroy `c`
587: Each process has the local and ghost coordinates
589: For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
590: and (x_0,y_0,z_0,x_1,y_1,z_1...)
592: .seealso: `DM`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`, `DMGetCoordinatesLocalNoncollective()`
593: @*/
594: PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c)
595: {
596: PetscFunctionBegin;
598: PetscAssertPointer(c, 2);
599: PetscCall(DMGetCoordinatesLocalSetUp(dm));
600: *c = dm->coordinates[0].xl;
601: PetscFunctionReturn(PETSC_SUCCESS);
602: }
604: /*@
605: DMGetCoordinatesLocalNoncollective - Non-collective version of `DMGetCoordinatesLocal()`. Fails if global coordinates have been set and `DMGetCoordinatesLocalSetUp()` not called.
607: Not Collective
609: Input Parameter:
610: . dm - the `DM`
612: Output Parameter:
613: . c - coordinate vector
615: Level: advanced
617: Note:
618: A previous call to `DMGetCoordinatesLocal()` or `DMGetCoordinatesLocalSetUp()` ensures that a call to this function will not error.
620: .seealso: `DM`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinatesLocal()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
621: @*/
622: PetscErrorCode DMGetCoordinatesLocalNoncollective(DM dm, Vec *c)
623: {
624: PetscFunctionBegin;
626: PetscAssertPointer(c, 2);
627: PetscCheck(dm->coordinates[0].xl || !dm->coordinates[0].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called");
628: *c = dm->coordinates[0].xl;
629: PetscFunctionReturn(PETSC_SUCCESS);
630: }
632: /*@
633: DMGetCoordinatesLocalTuple - Gets a local vector with the coordinates of specified points and the section describing its layout.
635: Not Collective
637: Input Parameters:
638: + dm - the `DM`
639: - p - the `IS` of points whose coordinates will be returned
641: Output Parameters:
642: + pCoordSection - the `PetscSection` describing the layout of pCoord, i.e. each point corresponds to one point in `p`, and DOFs correspond to coordinates
643: - pCoord - the `Vec` with coordinates of points in `p`
645: Level: advanced
647: Notes:
648: `DMGetCoordinatesLocalSetUp()` must be called first. This function employs `DMGetCoordinatesLocalNoncollective()` so it is not collective.
650: This creates a new vector, so the user SHOULD destroy this vector
652: Each process has the local and ghost coordinates
654: For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
655: and (x_0,y_0,z_0,x_1,y_1,z_1...)
657: .seealso: `DM`, `DMDA`, `DMSetCoordinatesLocal()`, `DMGetCoordinatesLocal()`, `DMGetCoordinatesLocalNoncollective()`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
658: @*/
659: PetscErrorCode DMGetCoordinatesLocalTuple(DM dm, IS p, PetscSection *pCoordSection, Vec *pCoord)
660: {
661: DM cdm;
662: PetscSection cs, newcs;
663: Vec coords;
664: const PetscScalar *arr;
665: PetscScalar *newarr = NULL;
666: PetscInt n;
668: PetscFunctionBegin;
671: if (pCoordSection) PetscAssertPointer(pCoordSection, 3);
672: if (pCoord) PetscAssertPointer(pCoord, 4);
673: PetscCall(DMGetCoordinateDM(dm, &cdm));
674: PetscCall(DMGetLocalSection(cdm, &cs));
675: PetscCall(DMGetCoordinatesLocal(dm, &coords));
676: PetscCheck(coords, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called or coordinates not set");
677: PetscCheck(cdm && cs, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM not supported");
678: PetscCall(VecGetArrayRead(coords, &arr));
679: PetscCall(PetscSectionExtractDofsFromArray(cs, MPIU_SCALAR, arr, p, &newcs, pCoord ? ((void **)&newarr) : NULL));
680: PetscCall(VecRestoreArrayRead(coords, &arr));
681: if (pCoord) {
682: PetscCall(PetscSectionGetStorageSize(newcs, &n));
683: /* set array in two steps to mimic PETSC_OWN_POINTER */
684: PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)p), 1, n, NULL, pCoord));
685: PetscCall(VecReplaceArray(*pCoord, newarr));
686: } else {
687: PetscCall(PetscFree(newarr));
688: }
689: if (pCoordSection) {
690: *pCoordSection = newcs;
691: } else PetscCall(PetscSectionDestroy(&newcs));
692: PetscFunctionReturn(PETSC_SUCCESS);
693: }
695: /*@
696: DMSetCoordinatesLocal - Sets into the `DM` a local vector, including ghost points, that holds the coordinates
698: Not Collective
700: Input Parameters:
701: + dm - the `DM`
702: - c - coordinate vector
704: Level: intermediate
706: Notes:
707: The coordinates of ghost points can be set using `DMSetCoordinates()`
708: followed by `DMGetCoordinatesLocal()`. This is intended to enable the
709: setting of ghost coordinates outside of the domain.
711: The vector `c` should be destroyed by the caller.
713: .seealso: `DM`, `DMGetCoordinatesLocal()`, `DMSetCoordinates()`, `DMGetCoordinates()`, `DMGetCoordinateDM()`
714: @*/
715: PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c)
716: {
717: PetscFunctionBegin;
720: PetscCall(PetscObjectReference((PetscObject)c));
721: PetscCall(VecDestroy(&dm->coordinates[0].xl));
722: dm->coordinates[0].xl = c;
723: PetscCall(VecDestroy(&dm->coordinates[0].x));
724: PetscFunctionReturn(PETSC_SUCCESS);
725: }
727: /*@
728: DMGetCellCoordinatesLocalSetUp - Prepares a local vector of cellwise coordinates, so that `DMGetCellCoordinatesLocalNoncollective()` can be used as non-collective afterwards.
730: Collective
732: Input Parameter:
733: . dm - the `DM`
735: Level: advanced
737: .seealso: `DM`, `DMGetCellCoordinatesLocalNoncollective()`
738: @*/
739: PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM dm)
740: {
741: PetscFunctionBegin;
743: if (!dm->coordinates[1].xl && dm->coordinates[1].x) {
744: DM cdm = NULL;
746: PetscCall(DMGetCellCoordinateDM(dm, &cdm));
747: PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[1].xl));
748: PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].xl, "DG coordinates"));
749: PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
750: PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
751: }
752: PetscFunctionReturn(PETSC_SUCCESS);
753: }
755: /*@
756: DMGetCellCoordinatesLocal - Gets a local vector with the cellwise coordinates associated with the `DM`.
758: Collective
760: Input Parameter:
761: . dm - the `DM`
763: Output Parameter:
764: . c - coordinate vector
766: Level: intermediate
768: Notes:
769: This is a borrowed reference, so the user should NOT destroy this vector
771: Each process has the local and ghost coordinates
773: .seealso: `DM`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`, `DMGetCellCoordinatesLocalNoncollective()`
774: @*/
775: PetscErrorCode DMGetCellCoordinatesLocal(DM dm, Vec *c)
776: {
777: PetscFunctionBegin;
779: PetscAssertPointer(c, 2);
780: PetscCall(DMGetCellCoordinatesLocalSetUp(dm));
781: *c = dm->coordinates[1].xl;
782: PetscFunctionReturn(PETSC_SUCCESS);
783: }
785: /*@
786: DMGetCellCoordinatesLocalNoncollective - Non-collective version of `DMGetCellCoordinatesLocal()`. Fails if global cellwise coordinates have been set and `DMGetCellCoordinatesLocalSetUp()` not called.
788: Not Collective
790: Input Parameter:
791: . dm - the `DM`
793: Output Parameter:
794: . c - cellwise coordinate vector
796: Level: advanced
798: .seealso: `DM`, `DMGetCellCoordinatesLocalSetUp()`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`
799: @*/
800: PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM dm, Vec *c)
801: {
802: PetscFunctionBegin;
804: PetscAssertPointer(c, 2);
805: PetscCheck(dm->coordinates[1].xl || !dm->coordinates[1].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCellCoordinatesLocalSetUp() has not been called");
806: *c = dm->coordinates[1].xl;
807: PetscFunctionReturn(PETSC_SUCCESS);
808: }
810: /*@
811: DMSetCellCoordinatesLocal - Sets into the `DM` a local vector including ghost points that holds the cellwise coordinates
813: Not Collective
815: Input Parameters:
816: + dm - the `DM`
817: - c - cellwise coordinate vector
819: Level: intermediate
821: Notes:
822: The coordinates of ghost points can be set using `DMSetCoordinates()`
823: followed by `DMGetCoordinatesLocal()`. This is intended to enable the
824: setting of ghost coordinates outside of the domain.
826: The vector `c` should be destroyed by the caller.
828: .seealso: `DM`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinates()`, `DMGetCellCoordinates()`, `DMGetCellCoordinateDM()`
829: @*/
830: PetscErrorCode DMSetCellCoordinatesLocal(DM dm, Vec c)
831: {
832: PetscFunctionBegin;
835: PetscCall(PetscObjectReference((PetscObject)c));
836: PetscCall(VecDestroy(&dm->coordinates[1].xl));
837: dm->coordinates[1].xl = c;
838: PetscCall(VecDestroy(&dm->coordinates[1].x));
839: PetscFunctionReturn(PETSC_SUCCESS);
840: }
842: PetscErrorCode DMGetCoordinateField(DM dm, DMField *field)
843: {
844: PetscFunctionBegin;
846: PetscAssertPointer(field, 2);
847: if (!dm->coordinates[0].field) PetscTryTypeMethod(dm, createcoordinatefield, &dm->coordinates[0].field);
848: *field = dm->coordinates[0].field;
849: PetscFunctionReturn(PETSC_SUCCESS);
850: }
852: PetscErrorCode DMSetCoordinateField(DM dm, DMField field)
853: {
854: PetscFunctionBegin;
857: PetscCall(PetscObjectReference((PetscObject)field));
858: PetscCall(DMFieldDestroy(&dm->coordinates[0].field));
859: dm->coordinates[0].field = field;
860: PetscFunctionReturn(PETSC_SUCCESS);
861: }
863: PetscErrorCode DMGetLocalBoundingBox_Coordinates(DM dm, PetscReal lmin[], PetscReal lmax[], PetscInt cs[], PetscInt ce[])
864: {
865: Vec coords = NULL;
866: PetscReal min[3] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL};
867: PetscReal max[3] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
868: PetscInt cdim, i, j;
869: PetscMPIInt size;
871: PetscFunctionBegin;
872: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
873: PetscCall(DMGetCoordinateDim(dm, &cdim));
874: if (size == 1) {
875: const PetscReal *L, *Lstart;
877: PetscCall(DMGetPeriodicity(dm, NULL, &Lstart, &L));
878: if (L) {
879: for (PetscInt d = 0; d < cdim; ++d)
880: if (L[d] > 0.0) {
881: min[d] = Lstart[d];
882: max[d] = Lstart[d] + L[d];
883: }
884: }
885: }
886: PetscCall(DMGetCoordinatesLocal(dm, &coords));
887: if (coords) {
888: const PetscScalar *local_coords;
889: PetscInt N, Ni;
891: for (j = cdim; j < 3; ++j) {
892: min[j] = 0;
893: max[j] = 0;
894: }
895: PetscCall(VecGetArrayRead(coords, &local_coords));
896: PetscCall(VecGetLocalSize(coords, &N));
897: Ni = N / cdim;
898: for (i = 0; i < Ni; ++i) {
899: for (j = 0; j < cdim; ++j) {
900: min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j]));
901: max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j]));
902: }
903: }
904: PetscCall(VecRestoreArrayRead(coords, &local_coords));
905: PetscCall(DMGetCellCoordinatesLocal(dm, &coords));
906: if (coords) {
907: PetscCall(VecGetArrayRead(coords, &local_coords));
908: PetscCall(VecGetLocalSize(coords, &N));
909: Ni = N / cdim;
910: for (i = 0; i < Ni; ++i) {
911: for (j = 0; j < cdim; ++j) {
912: min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j]));
913: max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j]));
914: }
915: }
916: PetscCall(VecRestoreArrayRead(coords, &local_coords));
917: }
918: if (lmin) PetscCall(PetscArraycpy(lmin, min, cdim));
919: if (lmax) PetscCall(PetscArraycpy(lmax, max, cdim));
920: }
921: PetscFunctionReturn(PETSC_SUCCESS);
922: }
924: /*@
925: DMGetLocalBoundingBox - Returns the bounding box for the piece of the `DM` on this process.
927: Not Collective
929: Input Parameter:
930: . dm - the `DM`
932: Output Parameters:
933: + lmin - local minimum coordinates (length coord dim, optional)
934: - lmax - local maximum coordinates (length coord dim, optional)
936: Level: beginner
938: Note:
939: If the `DM` is a `DMDA` and has no coordinates, the index bounds are returned instead.
941: .seealso: `DM`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetBoundingBox()`
942: @*/
943: PetscErrorCode DMGetLocalBoundingBox(DM dm, PetscReal lmin[], PetscReal lmax[])
944: {
945: PetscFunctionBegin;
947: PetscUseTypeMethod(dm, getlocalboundingbox, lmin, lmax, NULL, NULL);
948: PetscFunctionReturn(PETSC_SUCCESS);
949: }
951: /*@
952: DMGetBoundingBox - Returns the global bounding box for the `DM`.
954: Collective
956: Input Parameter:
957: . dm - the `DM`
959: Output Parameters:
960: + gmin - global minimum coordinates (length coord dim, optional)
961: - gmax - global maximum coordinates (length coord dim, optional)
963: Level: beginner
965: .seealso: `DM`, `DMGetLocalBoundingBox()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`
966: @*/
967: PetscErrorCode DMGetBoundingBox(DM dm, PetscReal gmin[], PetscReal gmax[])
968: {
969: PetscReal lmin[3], lmax[3];
970: const PetscReal *L, *Lstart;
971: PetscInt cdim;
972: PetscMPIInt count;
974: PetscFunctionBegin;
976: PetscCall(DMGetCoordinateDim(dm, &cdim));
977: PetscCall(PetscMPIIntCast(cdim, &count));
978: PetscCall(DMGetLocalBoundingBox(dm, lmin, lmax));
979: if (gmin) PetscCall(MPIU_Allreduce(lmin, gmin, count, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm)));
980: if (gmax) PetscCall(MPIU_Allreduce(lmax, gmax, count, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)dm)));
981: PetscCall(DMGetPeriodicity(dm, NULL, &Lstart, &L));
982: if (L) {
983: for (PetscInt d = 0; d < cdim; ++d)
984: if (L[d] > 0.0) {
985: gmin[d] = Lstart[d];
986: gmax[d] = Lstart[d] + L[d];
987: }
988: }
989: PetscFunctionReturn(PETSC_SUCCESS);
990: }
992: static PetscErrorCode DMCreateAffineCoordinates_Internal(DM dm)
993: {
994: DM cdm;
995: PetscFE feLinear;
996: DMPolytopeType ct;
997: PetscInt dim, dE, height, cStart, cEnd, gct;
999: PetscFunctionBegin;
1000: PetscCall(DMGetCoordinateDM(dm, &cdm));
1001: PetscCall(DMGetDimension(dm, &dim));
1002: PetscCall(DMGetCoordinateDim(dm, &dE));
1003: PetscCall(DMPlexGetVTKCellHeight(dm, &height));
1004: PetscCall(DMPlexGetHeightStratum(dm, height, &cStart, &cEnd));
1005: if (cEnd > cStart) PetscCall(DMPlexGetCellType(dm, cStart, &ct));
1006: else ct = DM_POLYTOPE_UNKNOWN;
1007: gct = (PetscInt)ct;
1008: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &gct, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm)));
1009: ct = (DMPolytopeType)gct;
1010: // Work around current bug in PetscDualSpaceSetUp_Lagrange()
1011: // Can be seen in plex_tutorials-ex10_1
1012: if (ct != DM_POLYTOPE_SEG_PRISM_TENSOR && ct != DM_POLYTOPE_TRI_PRISM_TENSOR && ct != DM_POLYTOPE_QUAD_PRISM_TENSOR) {
1013: PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, dE, ct, 1, -1, &feLinear));
1014: PetscCall(DMSetField(cdm, 0, NULL, (PetscObject)feLinear));
1015: PetscCall(PetscFEDestroy(&feLinear));
1016: PetscCall(DMCreateDS(cdm));
1017: }
1018: PetscFunctionReturn(PETSC_SUCCESS);
1019: }
1021: PetscErrorCode DMGetCoordinateDegree_Internal(DM dm, PetscInt *degree)
1022: {
1023: DM cdm;
1024: PetscFE fe;
1025: PetscSpace sp;
1026: PetscClassId id;
1028: PetscFunctionBegin;
1029: *degree = 1;
1030: PetscCall(DMGetCoordinateDM(dm, &cdm));
1031: PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe));
1032: PetscCall(PetscObjectGetClassId((PetscObject)fe, &id));
1033: if (id != PETSCFE_CLASSID) PetscFunctionReturn(PETSC_SUCCESS);
1034: PetscCall(PetscFEGetBasisSpace(fe, &sp));
1035: PetscCall(PetscSpaceGetDegree(sp, degree, NULL));
1036: PetscFunctionReturn(PETSC_SUCCESS);
1037: }
1039: /*@
1040: DMSetCoordinateDisc - Set a coordinate space
1042: Input Parameters:
1043: + dm - The `DM` object
1044: . disc - The new coordinate discretization or `NULL` to ensure a coordinate discretization exists
1045: - project - Project coordinates to new discretization
1047: Level: intermediate
1049: Notes:
1050: 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.
1052: This function takes the current mesh coordinates, which are discretized using some `PetscFE` space, and projects this function into a new `PetscFE` space.
1053: The coordinate projection is done on the continuous coordinates, but the discontinuous coordinates are not updated.
1055: Developer Note:
1056: With more effort, we could directly project the discontinuous coordinates also.
1058: .seealso: `DM`, `PetscFE`, `DMGetCoordinateField()`
1059: @*/
1060: PetscErrorCode DMSetCoordinateDisc(DM dm, PetscFE disc, PetscBool project)
1061: {
1062: DM cdmOld, cdmNew;
1063: PetscFE discOld;
1064: PetscClassId classid;
1065: Vec coordsOld, coordsNew;
1066: PetscBool same_space = PETSC_TRUE;
1067: const char *prefix;
1069: PetscFunctionBegin;
1073: PetscCall(DMGetCoordinateDM(dm, &cdmOld));
1074: /* Check current discretization is compatible */
1075: PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
1076: PetscCall(PetscObjectGetClassId((PetscObject)discOld, &classid));
1077: if (classid != PETSCFE_CLASSID) {
1078: if (classid == PETSC_CONTAINER_CLASSID) {
1079: PetscCall(DMCreateAffineCoordinates_Internal(dm));
1080: PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
1081: } else {
1082: const char *discname;
1084: PetscCall(PetscObjectGetType((PetscObject)discOld, &discname));
1085: SETERRQ(PetscObjectComm((PetscObject)discOld), PETSC_ERR_SUP, "Discretization type %s not supported", discname);
1086: }
1087: }
1088: // Linear space has been created by now
1089: if (!disc) PetscFunctionReturn(PETSC_SUCCESS);
1090: // Check if the new space is the same as the old modulo quadrature
1091: {
1092: PetscDualSpace dsOld, ds;
1093: PetscCall(PetscFEGetDualSpace(discOld, &dsOld));
1094: PetscCall(PetscFEGetDualSpace(disc, &ds));
1095: PetscCall(PetscDualSpaceEqual(dsOld, ds, &same_space));
1096: }
1097: // Make a fresh clone of the coordinate DM
1098: PetscCall(DMClone(cdmOld, &cdmNew));
1099: cdmNew->cloneOpts = PETSC_TRUE;
1100: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)cdmOld, &prefix));
1101: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)cdmNew, prefix));
1102: PetscCall(DMSetField(cdmNew, 0, NULL, (PetscObject)disc));
1103: PetscCall(DMCreateDS(cdmNew));
1104: {
1105: PetscDS ds, nds;
1107: PetscCall(DMGetDS(cdmOld, &ds));
1108: PetscCall(DMGetDS(cdmNew, &nds));
1109: PetscCall(PetscDSCopyConstants(ds, nds));
1110: }
1111: if (cdmOld->periodic.setup) {
1112: cdmNew->periodic.setup = cdmOld->periodic.setup;
1113: PetscCall(cdmNew->periodic.setup(cdmNew));
1114: }
1115: if (dm->setfromoptionscalled) PetscCall(DMSetFromOptions(cdmNew));
1116: if (project) {
1117: PetscCall(DMGetCoordinates(dm, &coordsOld));
1118: PetscCall(DMCreateGlobalVector(cdmNew, &coordsNew));
1119: if (same_space) {
1120: // Need to copy so that the new vector has the right dm
1121: PetscCall(VecCopy(coordsOld, coordsNew));
1122: } else {
1123: Mat In;
1125: PetscCall(DMCreateInterpolation(cdmOld, cdmNew, &In, NULL));
1126: PetscCall(MatMult(In, coordsOld, coordsNew));
1127: PetscCall(MatDestroy(&In));
1128: }
1129: PetscCall(DMSetCoordinates(dm, coordsNew));
1130: PetscCall(VecDestroy(&coordsNew));
1131: }
1132: /* Set new coordinate structures */
1133: PetscCall(DMSetCoordinateField(dm, NULL));
1134: PetscCall(DMSetCoordinateDM(dm, cdmNew));
1135: PetscCall(DMDestroy(&cdmNew));
1136: PetscFunctionReturn(PETSC_SUCCESS);
1137: }
1139: /*@
1140: DMLocatePoints - Locate the points in `v` in the mesh and return a `PetscSF` of the containing cells
1142: Collective
1144: Input Parameters:
1145: + dm - The `DM`
1146: - ltype - The type of point location, e.g. `DM_POINTLOCATION_NONE` or `DM_POINTLOCATION_NEAREST`
1148: Input/Output Parameters:
1149: + v - The `Vec` of points, on output contains the nearest mesh points to the given points if `DM_POINTLOCATION_NEAREST` is used
1150: - cellSF - Points to either `NULL`, or a `PetscSF` with guesses for which cells contain each point;
1151: on output, the `PetscSF` containing the MPI ranks and local indices of the containing points
1153: Level: developer
1155: Notes:
1156: To do a search of the local cells of the mesh, `v` should have `PETSC_COMM_SELF` as its communicator.
1157: To do a search of all the cells in the distributed mesh, `v` should have the same MPI communicator as `dm`.
1159: Points will only be located in owned cells, not overlap cells arising from `DMPlexDistribute()` or other overlapping distributions.
1161: If *cellSF is `NULL` on input, a `PetscSF` will be created.
1162: If *cellSF is not `NULL` on input, it should point to an existing `PetscSF`, whose graph will be used as initial guesses.
1164: An array that maps each point to its containing cell can be obtained with
1165: .vb
1166: const PetscSFNode *cells;
1167: PetscInt nFound;
1168: const PetscInt *found;
1170: PetscSFGetGraph(cellSF,NULL,&nFound,&found,&cells);
1171: .ve
1173: 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
1174: 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.
1176: .seealso: `DM`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMPointLocationType`
1177: @*/
1178: PetscErrorCode DMLocatePoints(DM dm, Vec v, DMPointLocationType ltype, PetscSF *cellSF)
1179: {
1180: PetscFunctionBegin;
1183: PetscAssertPointer(cellSF, 4);
1184: if (*cellSF) {
1185: PetscMPIInt result;
1188: PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)v), PetscObjectComm((PetscObject)*cellSF), &result));
1189: PetscCheck(result == MPI_IDENT || result == MPI_CONGRUENT, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "cellSF must have a communicator congruent to v's");
1190: } else {
1191: PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)v), cellSF));
1192: }
1193: PetscCall(PetscLogEventBegin(DM_LocatePoints, dm, 0, 0, 0));
1194: PetscUseTypeMethod(dm, locatepoints, v, ltype, *cellSF);
1195: PetscCall(PetscLogEventEnd(DM_LocatePoints, dm, 0, 0, 0));
1196: PetscFunctionReturn(PETSC_SUCCESS);
1197: }