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 if the global vector with coordinates has not been set yet but the local vector with coordinates has been set
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), see `DMGetCoordinatesLocal()`.
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: Does not work for `DMSTAG`
420: .seealso: `DM`, `DMDA`, `DMSetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()`
421: @*/
422: PetscErrorCode DMGetCoordinates(DM dm, Vec *c)
423: {
424: PetscFunctionBegin;
426: PetscAssertPointer(c, 2);
427: if (!dm->coordinates[0].x && dm->coordinates[0].xl) {
428: DM cdm = NULL;
430: PetscCall(DMGetCoordinateDM(dm, &cdm));
431: PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[0].x));
432: PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].x, "coordinates"));
433: PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x));
434: PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x));
435: }
436: *c = dm->coordinates[0].x;
437: PetscFunctionReturn(PETSC_SUCCESS);
438: }
440: /*@
441: DMSetCoordinates - Sets into the `DM` a global vector that holds the coordinates
443: Logically Collective
445: Input Parameters:
446: + dm - the `DM`
447: - c - coordinate vector
449: Level: intermediate
451: Notes:
452: The coordinates do not include those for ghost points, which are in the local vector.
454: The vector `c` can be destroyed after the call
456: .seealso: `DM`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()`
457: @*/
458: PetscErrorCode DMSetCoordinates(DM dm, Vec c)
459: {
460: PetscFunctionBegin;
463: PetscCall(PetscObjectReference((PetscObject)c));
464: PetscCall(VecDestroy(&dm->coordinates[0].x));
465: dm->coordinates[0].x = c;
466: PetscCall(VecDestroy(&dm->coordinates[0].xl));
467: PetscCall(DMCoarsenHookAdd(dm, DMRestrictHook_Coordinates, NULL, NULL));
468: PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_Coordinates, NULL, NULL));
469: PetscFunctionReturn(PETSC_SUCCESS);
470: }
472: /*@
473: DMGetCellCoordinates - Gets a global vector with the cellwise coordinates associated with the `DM`.
475: Collective
477: Input Parameter:
478: . dm - the `DM`
480: Output Parameter:
481: . c - global coordinate vector
483: Level: intermediate
485: Notes:
486: This is a borrowed reference, so the user should NOT destroy this vector. When the `DM` is
487: destroyed `c` will no longer be valid.
489: Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates).
491: .seealso: `DM`, `DMGetCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()`
492: @*/
493: PetscErrorCode DMGetCellCoordinates(DM dm, Vec *c)
494: {
495: PetscFunctionBegin;
497: PetscAssertPointer(c, 2);
498: if (!dm->coordinates[1].x && dm->coordinates[1].xl) {
499: DM cdm = NULL;
501: PetscCall(DMGetCellCoordinateDM(dm, &cdm));
502: PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[1].x));
503: PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].x, "DG coordinates"));
504: PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x));
505: PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x));
506: }
507: *c = dm->coordinates[1].x;
508: PetscFunctionReturn(PETSC_SUCCESS);
509: }
511: /*@
512: DMSetCellCoordinates - Sets into the `DM` a global vector that holds the cellwise coordinates
514: Collective
516: Input Parameters:
517: + dm - the `DM`
518: - c - cellwise coordinate vector
520: Level: intermediate
522: Notes:
523: The coordinates do not include those for ghost points, which are in the local vector.
525: The vector `c` should be destroyed by the caller.
527: .seealso: `DM`, `DMGetCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()`
528: @*/
529: PetscErrorCode DMSetCellCoordinates(DM dm, Vec c)
530: {
531: PetscFunctionBegin;
534: PetscCall(PetscObjectReference((PetscObject)c));
535: PetscCall(VecDestroy(&dm->coordinates[1].x));
536: dm->coordinates[1].x = c;
537: PetscCall(VecDestroy(&dm->coordinates[1].xl));
538: PetscFunctionReturn(PETSC_SUCCESS);
539: }
541: /*@
542: DMGetCoordinatesLocalSetUp - Prepares a local vector of coordinates, so that `DMGetCoordinatesLocalNoncollective()` can be used as non-collective afterwards.
544: Collective
546: Input Parameter:
547: . dm - the `DM`
549: Level: advanced
551: .seealso: `DM`, `DMSetCoordinates()`, `DMGetCoordinatesLocalNoncollective()`
552: @*/
553: PetscErrorCode DMGetCoordinatesLocalSetUp(DM dm)
554: {
555: PetscFunctionBegin;
557: if (!dm->coordinates[0].xl && dm->coordinates[0].x) {
558: DM cdm = NULL;
559: PetscInt bs;
561: PetscCall(DMGetCoordinateDM(dm, &cdm));
562: PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[0].xl));
563: PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].xl, "Local Coordinates"));
564: // If the size of the vector is 0, it will not get the right block size
565: PetscCall(VecGetBlockSize(dm->coordinates[0].x, &bs));
566: PetscCall(VecSetBlockSize(dm->coordinates[0].xl, bs));
567: PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].xl, "coordinates"));
568: PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
569: PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
570: }
571: PetscFunctionReturn(PETSC_SUCCESS);
572: }
574: /*@
575: DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the `DM`.
577: Collective the first time it is called
579: Input Parameter:
580: . dm - the `DM`
582: Output Parameter:
583: . c - coordinate vector
585: Level: intermediate
587: Notes:
588: This is a borrowed reference, so the user should NOT destroy `c`
590: Each process has the local and ghost coordinates
592: For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
593: and (x_0,y_0,z_0,x_1,y_1,z_1...)
595: .seealso: `DM`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`, `DMGetCoordinatesLocalNoncollective()`
596: @*/
597: PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c)
598: {
599: PetscFunctionBegin;
601: PetscAssertPointer(c, 2);
602: PetscCall(DMGetCoordinatesLocalSetUp(dm));
603: *c = dm->coordinates[0].xl;
604: PetscFunctionReturn(PETSC_SUCCESS);
605: }
607: /*@
608: DMGetCoordinatesLocalNoncollective - Non-collective version of `DMGetCoordinatesLocal()`. Fails if global coordinates have been set and `DMGetCoordinatesLocalSetUp()` not called.
610: Not Collective
612: Input Parameter:
613: . dm - the `DM`
615: Output Parameter:
616: . c - coordinate vector
618: Level: advanced
620: Note:
621: A previous call to `DMGetCoordinatesLocal()` or `DMGetCoordinatesLocalSetUp()` ensures that a call to this function will not error.
623: .seealso: `DM`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinatesLocal()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
624: @*/
625: PetscErrorCode DMGetCoordinatesLocalNoncollective(DM dm, Vec *c)
626: {
627: PetscFunctionBegin;
629: PetscAssertPointer(c, 2);
630: PetscCheck(dm->coordinates[0].xl || !dm->coordinates[0].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called");
631: *c = dm->coordinates[0].xl;
632: PetscFunctionReturn(PETSC_SUCCESS);
633: }
635: /*@
636: DMGetCoordinatesLocalTuple - Gets a local vector with the coordinates of specified points and the section describing its layout.
638: Not Collective
640: Input Parameters:
641: + dm - the `DM`
642: - p - the `IS` of points whose coordinates will be returned
644: Output Parameters:
645: + pCoordSection - the `PetscSection` describing the layout of pCoord, i.e. each point corresponds to one point in `p`, and DOFs correspond to coordinates
646: - pCoord - the `Vec` with coordinates of points in `p`
648: Level: advanced
650: Notes:
651: `DMGetCoordinatesLocalSetUp()` must be called first. This function employs `DMGetCoordinatesLocalNoncollective()` so it is not collective.
653: This creates a new vector, so the user SHOULD destroy this vector
655: Each process has the local and ghost coordinates
657: For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
658: and (x_0,y_0,z_0,x_1,y_1,z_1...)
660: .seealso: `DM`, `DMDA`, `DMSetCoordinatesLocal()`, `DMGetCoordinatesLocal()`, `DMGetCoordinatesLocalNoncollective()`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
661: @*/
662: PetscErrorCode DMGetCoordinatesLocalTuple(DM dm, IS p, PetscSection *pCoordSection, Vec *pCoord)
663: {
664: DM cdm;
665: PetscSection cs, newcs;
666: Vec coords;
667: const PetscScalar *arr;
668: PetscScalar *newarr = NULL;
669: PetscInt n;
671: PetscFunctionBegin;
674: if (pCoordSection) PetscAssertPointer(pCoordSection, 3);
675: if (pCoord) PetscAssertPointer(pCoord, 4);
676: PetscCall(DMGetCoordinateDM(dm, &cdm));
677: PetscCall(DMGetLocalSection(cdm, &cs));
678: PetscCall(DMGetCoordinatesLocal(dm, &coords));
679: PetscCheck(coords, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called or coordinates not set");
680: PetscCheck(cdm && cs, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM not supported");
681: PetscCall(VecGetArrayRead(coords, &arr));
682: PetscCall(PetscSectionExtractDofsFromArray(cs, MPIU_SCALAR, arr, p, &newcs, pCoord ? ((void **)&newarr) : NULL));
683: PetscCall(VecRestoreArrayRead(coords, &arr));
684: if (pCoord) {
685: PetscCall(PetscSectionGetStorageSize(newcs, &n));
686: /* set array in two steps to mimic PETSC_OWN_POINTER */
687: PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)p), 1, n, NULL, pCoord));
688: PetscCall(VecReplaceArray(*pCoord, newarr));
689: } else {
690: PetscCall(PetscFree(newarr));
691: }
692: if (pCoordSection) {
693: *pCoordSection = newcs;
694: } else PetscCall(PetscSectionDestroy(&newcs));
695: PetscFunctionReturn(PETSC_SUCCESS);
696: }
698: /*@
699: DMSetCoordinatesLocal - Sets into the `DM` a local vector, including ghost points, that holds the coordinates
701: Not Collective
703: Input Parameters:
704: + dm - the `DM`
705: - c - coordinate vector
707: Level: intermediate
709: Notes:
710: The coordinates of ghost points can be set using `DMSetCoordinates()`
711: followed by `DMGetCoordinatesLocal()`. This is intended to enable the
712: setting of ghost coordinates outside of the domain.
714: The vector `c` should be destroyed by the caller.
716: .seealso: `DM`, `DMGetCoordinatesLocal()`, `DMSetCoordinates()`, `DMGetCoordinates()`, `DMGetCoordinateDM()`
717: @*/
718: PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c)
719: {
720: PetscFunctionBegin;
723: PetscCall(PetscObjectReference((PetscObject)c));
724: PetscCall(VecDestroy(&dm->coordinates[0].xl));
725: dm->coordinates[0].xl = c;
726: PetscCall(VecDestroy(&dm->coordinates[0].x));
727: PetscFunctionReturn(PETSC_SUCCESS);
728: }
730: /*@
731: DMGetCellCoordinatesLocalSetUp - Prepares a local vector of cellwise coordinates, so that `DMGetCellCoordinatesLocalNoncollective()` can be used as non-collective afterwards.
733: Collective
735: Input Parameter:
736: . dm - the `DM`
738: Level: advanced
740: .seealso: `DM`, `DMGetCellCoordinatesLocalNoncollective()`
741: @*/
742: PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM dm)
743: {
744: PetscFunctionBegin;
746: if (!dm->coordinates[1].xl && dm->coordinates[1].x) {
747: DM cdm = NULL;
749: PetscCall(DMGetCellCoordinateDM(dm, &cdm));
750: PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[1].xl));
751: PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].xl, "DG coordinates"));
752: PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
753: PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
754: }
755: PetscFunctionReturn(PETSC_SUCCESS);
756: }
758: /*@
759: DMGetCellCoordinatesLocal - Gets a local vector with the cellwise coordinates associated with the `DM`.
761: Collective
763: Input Parameter:
764: . dm - the `DM`
766: Output Parameter:
767: . c - coordinate vector
769: Level: intermediate
771: Notes:
772: This is a borrowed reference, so the user should NOT destroy this vector
774: Each process has the local and ghost coordinates
776: .seealso: `DM`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`, `DMGetCellCoordinatesLocalNoncollective()`
777: @*/
778: PetscErrorCode DMGetCellCoordinatesLocal(DM dm, Vec *c)
779: {
780: PetscFunctionBegin;
782: PetscAssertPointer(c, 2);
783: PetscCall(DMGetCellCoordinatesLocalSetUp(dm));
784: *c = dm->coordinates[1].xl;
785: PetscFunctionReturn(PETSC_SUCCESS);
786: }
788: /*@
789: DMGetCellCoordinatesLocalNoncollective - Non-collective version of `DMGetCellCoordinatesLocal()`. Fails if global cellwise coordinates have been set and `DMGetCellCoordinatesLocalSetUp()` not called.
791: Not Collective
793: Input Parameter:
794: . dm - the `DM`
796: Output Parameter:
797: . c - cellwise coordinate vector
799: Level: advanced
801: .seealso: `DM`, `DMGetCellCoordinatesLocalSetUp()`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`
802: @*/
803: PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM dm, Vec *c)
804: {
805: PetscFunctionBegin;
807: PetscAssertPointer(c, 2);
808: PetscCheck(dm->coordinates[1].xl || !dm->coordinates[1].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCellCoordinatesLocalSetUp() has not been called");
809: *c = dm->coordinates[1].xl;
810: PetscFunctionReturn(PETSC_SUCCESS);
811: }
813: /*@
814: DMSetCellCoordinatesLocal - Sets into the `DM` a local vector including ghost points that holds the cellwise coordinates
816: Not Collective
818: Input Parameters:
819: + dm - the `DM`
820: - c - cellwise coordinate vector
822: Level: intermediate
824: Notes:
825: The coordinates of ghost points can be set using `DMSetCoordinates()`
826: followed by `DMGetCoordinatesLocal()`. This is intended to enable the
827: setting of ghost coordinates outside of the domain.
829: The vector `c` should be destroyed by the caller.
831: .seealso: `DM`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinates()`, `DMGetCellCoordinates()`, `DMGetCellCoordinateDM()`
832: @*/
833: PetscErrorCode DMSetCellCoordinatesLocal(DM dm, Vec c)
834: {
835: PetscFunctionBegin;
838: PetscCall(PetscObjectReference((PetscObject)c));
839: PetscCall(VecDestroy(&dm->coordinates[1].xl));
840: dm->coordinates[1].xl = c;
841: PetscCall(VecDestroy(&dm->coordinates[1].x));
842: PetscFunctionReturn(PETSC_SUCCESS);
843: }
845: PetscErrorCode DMGetCoordinateField(DM dm, DMField *field)
846: {
847: PetscFunctionBegin;
849: PetscAssertPointer(field, 2);
850: if (!dm->coordinates[0].field) PetscTryTypeMethod(dm, createcoordinatefield, &dm->coordinates[0].field);
851: *field = dm->coordinates[0].field;
852: PetscFunctionReturn(PETSC_SUCCESS);
853: }
855: PetscErrorCode DMSetCoordinateField(DM dm, DMField field)
856: {
857: PetscFunctionBegin;
860: PetscCall(PetscObjectReference((PetscObject)field));
861: PetscCall(DMFieldDestroy(&dm->coordinates[0].field));
862: dm->coordinates[0].field = field;
863: PetscFunctionReturn(PETSC_SUCCESS);
864: }
866: PetscErrorCode DMGetLocalBoundingBox_Coordinates(DM dm, PetscReal lmin[], PetscReal lmax[], PetscInt cs[], PetscInt ce[])
867: {
868: Vec coords = NULL;
869: PetscReal min[3] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL};
870: PetscReal max[3] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
871: PetscInt cdim, i, j;
872: PetscMPIInt size;
874: PetscFunctionBegin;
875: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
876: PetscCall(DMGetCoordinateDim(dm, &cdim));
877: if (size == 1) {
878: const PetscReal *L, *Lstart;
880: PetscCall(DMGetPeriodicity(dm, NULL, &Lstart, &L));
881: if (L) {
882: for (PetscInt d = 0; d < cdim; ++d)
883: if (L[d] > 0.0) {
884: min[d] = Lstart[d];
885: max[d] = Lstart[d] + L[d];
886: }
887: }
888: }
889: PetscCall(DMGetCoordinatesLocal(dm, &coords));
890: if (coords) {
891: const PetscScalar *local_coords;
892: PetscInt N, Ni;
894: for (j = cdim; j < 3; ++j) {
895: min[j] = 0;
896: max[j] = 0;
897: }
898: PetscCall(VecGetArrayRead(coords, &local_coords));
899: PetscCall(VecGetLocalSize(coords, &N));
900: Ni = N / cdim;
901: for (i = 0; i < Ni; ++i) {
902: for (j = 0; j < cdim; ++j) {
903: min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j]));
904: max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j]));
905: }
906: }
907: PetscCall(VecRestoreArrayRead(coords, &local_coords));
908: PetscCall(DMGetCellCoordinatesLocal(dm, &coords));
909: if (coords) {
910: PetscCall(VecGetArrayRead(coords, &local_coords));
911: PetscCall(VecGetLocalSize(coords, &N));
912: Ni = N / cdim;
913: for (i = 0; i < Ni; ++i) {
914: for (j = 0; j < cdim; ++j) {
915: min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j]));
916: max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j]));
917: }
918: }
919: PetscCall(VecRestoreArrayRead(coords, &local_coords));
920: }
921: if (lmin) PetscCall(PetscArraycpy(lmin, min, cdim));
922: if (lmax) PetscCall(PetscArraycpy(lmax, max, cdim));
923: }
924: PetscFunctionReturn(PETSC_SUCCESS);
925: }
927: /*@
928: DMGetLocalBoundingBox - Returns the bounding box for the piece of the `DM` on this process.
930: Not Collective
932: Input Parameter:
933: . dm - the `DM`
935: Output Parameters:
936: + lmin - local minimum coordinates (length coord dim, optional)
937: - lmax - local maximum coordinates (length coord dim, optional)
939: Level: beginner
941: Note:
942: If the `DM` is a `DMDA` and has no coordinates, the index bounds are returned instead.
944: .seealso: `DM`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetBoundingBox()`
945: @*/
946: PetscErrorCode DMGetLocalBoundingBox(DM dm, PetscReal lmin[], PetscReal lmax[])
947: {
948: PetscFunctionBegin;
950: PetscUseTypeMethod(dm, getlocalboundingbox, lmin, lmax, NULL, NULL);
951: PetscFunctionReturn(PETSC_SUCCESS);
952: }
954: /*@
955: DMGetBoundingBox - Returns the global bounding box for the `DM`.
957: Collective
959: Input Parameter:
960: . dm - the `DM`
962: Output Parameters:
963: + gmin - global minimum coordinates (length coord dim, optional)
964: - gmax - global maximum coordinates (length coord dim, optional)
966: Level: beginner
968: .seealso: `DM`, `DMGetLocalBoundingBox()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`
969: @*/
970: PetscErrorCode DMGetBoundingBox(DM dm, PetscReal gmin[], PetscReal gmax[])
971: {
972: PetscReal lmin[3], lmax[3];
973: const PetscReal *L, *Lstart;
974: PetscInt cdim;
976: PetscFunctionBegin;
978: PetscCall(DMGetCoordinateDim(dm, &cdim));
979: PetscCall(DMGetLocalBoundingBox(dm, lmin, lmax));
980: if (gmin) PetscCallMPI(MPIU_Allreduce(lmin, gmin, cdim, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm)));
981: if (gmax) PetscCallMPI(MPIU_Allreduce(lmax, gmax, cdim, 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: static void evaluate_coordinates(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 xnew[])
1041: {
1042: for (PetscInt i = 0; i < dim; i++) xnew[i] = x[i];
1043: }
1045: /*@
1046: DMSetCoordinateDisc - Set a coordinate space
1048: Input Parameters:
1049: + dm - The `DM` object
1050: . disc - The new coordinate discretization or `NULL` to ensure a coordinate discretization exists
1051: - project - Project coordinates to new discretization
1053: Level: intermediate
1055: Notes:
1056: 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.
1058: This function takes the current mesh coordinates, which are discretized using some `PetscFE` space, and projects this function into a new `PetscFE` space.
1059: The coordinate projection is done on the continuous coordinates, but the discontinuous coordinates are not updated.
1061: Developer Note:
1062: With more effort, we could directly project the discontinuous coordinates also.
1064: .seealso: `DM`, `PetscFE`, `DMGetCoordinateField()`
1065: @*/
1066: PetscErrorCode DMSetCoordinateDisc(DM dm, PetscFE disc, PetscBool project)
1067: {
1068: DM cdmOld, cdmNew;
1069: PetscFE discOld;
1070: PetscClassId classid;
1071: PetscBool same_space = PETSC_TRUE;
1072: const char *prefix;
1074: PetscFunctionBegin;
1078: PetscCall(DMGetCoordinateDM(dm, &cdmOld));
1079: /* Check current discretization is compatible */
1080: PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
1081: PetscCall(PetscObjectGetClassId((PetscObject)discOld, &classid));
1082: if (classid != PETSCFE_CLASSID) {
1083: if (classid == PETSC_CONTAINER_CLASSID) {
1084: PetscCall(DMCreateAffineCoordinates_Internal(dm));
1085: PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
1086: } else {
1087: const char *discname;
1089: PetscCall(PetscObjectGetType((PetscObject)discOld, &discname));
1090: SETERRQ(PetscObjectComm((PetscObject)discOld), PETSC_ERR_SUP, "Discretization type %s not supported", discname);
1091: }
1092: }
1093: // Linear space has been created by now
1094: if (!disc) PetscFunctionReturn(PETSC_SUCCESS);
1095: // Check if the new space is the same as the old modulo quadrature
1096: {
1097: PetscDualSpace dsOld, ds;
1098: PetscCall(PetscFEGetDualSpace(discOld, &dsOld));
1099: PetscCall(PetscFEGetDualSpace(disc, &ds));
1100: PetscCall(PetscDualSpaceEqual(dsOld, ds, &same_space));
1101: }
1102: // Make a fresh clone of the coordinate DM
1103: PetscCall(DMClone(cdmOld, &cdmNew));
1104: cdmNew->cloneOpts = PETSC_TRUE;
1105: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)cdmOld, &prefix));
1106: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)cdmNew, prefix));
1107: PetscCall(DMSetField(cdmNew, 0, NULL, (PetscObject)disc));
1108: PetscCall(DMCreateDS(cdmNew));
1109: {
1110: PetscDS ds, nds;
1112: PetscCall(DMGetDS(cdmOld, &ds));
1113: PetscCall(DMGetDS(cdmNew, &nds));
1114: PetscCall(PetscDSCopyConstants(ds, nds));
1115: }
1116: if (cdmOld->periodic.setup) {
1117: PetscSF dummy;
1118: // Force IsoperiodicPointSF to be built, required for periodic coordinate setup
1119: PetscCall(DMGetIsoperiodicPointSF_Internal(dm, &dummy));
1120: cdmNew->periodic.setup = cdmOld->periodic.setup;
1121: PetscCall(cdmNew->periodic.setup(cdmNew));
1122: }
1123: if (dm->setfromoptionscalled) PetscCall(DMSetFromOptions(cdmNew));
1124: if (project) {
1125: Vec coordsOld, coordsNew;
1126: PetscInt num_face_sfs = 0;
1128: PetscCall(DMPlexGetIsoperiodicFaceSF(dm, &num_face_sfs, NULL));
1129: if (num_face_sfs) { // Isoperiodicity requires projecting the local coordinates
1130: PetscCall(DMGetCoordinatesLocal(dm, &coordsOld));
1131: PetscCall(DMCreateLocalVector(cdmNew, &coordsNew));
1132: PetscCall(PetscObjectSetName((PetscObject)coordsNew, "coordinates"));
1133: if (same_space) {
1134: // Need to copy so that the new vector has the right dm
1135: PetscCall(VecCopy(coordsOld, coordsNew));
1136: } else {
1137: void (*funcs[])(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]) = {evaluate_coordinates};
1139: // We can't call DMProjectField directly because it depends on KSP for DMGlobalToLocalSolve(), but we can use the core strategy
1140: PetscCall(DMSetCoordinateDM(cdmNew, cdmOld));
1141: // See DMPlexRemapGeometry() for a similar pattern handling the coordinate field
1142: DMField cf;
1143: PetscCall(DMGetCoordinateField(dm, &cf));
1144: cdmNew->coordinates[0].field = cf;
1145: PetscCall(DMProjectFieldLocal(cdmNew, 0.0, NULL, funcs, INSERT_VALUES, coordsNew));
1146: cdmNew->coordinates[0].field = NULL;
1147: PetscCall(DMSetCoordinateDM(cdmNew, NULL));
1148: }
1149: PetscCall(DMSetCoordinatesLocal(dm, coordsNew));
1150: PetscCall(VecDestroy(&coordsNew));
1151: } else {
1152: PetscCall(DMGetCoordinates(dm, &coordsOld));
1153: PetscCall(DMCreateGlobalVector(cdmNew, &coordsNew));
1154: if (same_space) {
1155: // Need to copy so that the new vector has the right dm
1156: PetscCall(VecCopy(coordsOld, coordsNew));
1157: } else {
1158: Mat In;
1160: PetscCall(DMCreateInterpolation(cdmOld, cdmNew, &In, NULL));
1161: PetscCall(MatMult(In, coordsOld, coordsNew));
1162: PetscCall(MatDestroy(&In));
1163: }
1164: PetscCall(DMSetCoordinates(dm, coordsNew));
1165: PetscCall(VecDestroy(&coordsNew));
1166: }
1167: }
1168: /* Set new coordinate structures */
1169: PetscCall(DMSetCoordinateField(dm, NULL));
1170: PetscCall(DMSetCoordinateDM(dm, cdmNew));
1171: PetscCall(DMDestroy(&cdmNew));
1172: PetscFunctionReturn(PETSC_SUCCESS);
1173: }
1175: /*@
1176: DMLocatePoints - Locate the points in `v` in the mesh and return a `PetscSF` of the containing cells
1178: Collective
1180: Input Parameters:
1181: + dm - The `DM`
1182: - ltype - The type of point location, e.g. `DM_POINTLOCATION_NONE` or `DM_POINTLOCATION_NEAREST`
1184: Input/Output Parameters:
1185: + v - The `Vec` of points, on output contains the nearest mesh points to the given points if `DM_POINTLOCATION_NEAREST` is used
1186: - cellSF - Points to either `NULL`, or a `PetscSF` with guesses for which cells contain each point;
1187: on output, the `PetscSF` containing the MPI ranks and local indices of the containing points
1189: Level: developer
1191: Notes:
1192: To do a search of the local cells of the mesh, `v` should have `PETSC_COMM_SELF` as its communicator.
1193: To do a search of all the cells in the distributed mesh, `v` should have the same MPI communicator as `dm`.
1195: Points will only be located in owned cells, not overlap cells arising from `DMPlexDistribute()` or other overlapping distributions.
1197: If *cellSF is `NULL` on input, a `PetscSF` will be created.
1198: If *cellSF is not `NULL` on input, it should point to an existing `PetscSF`, whose graph will be used as initial guesses.
1200: An array that maps each point to its containing cell can be obtained with
1201: .vb
1202: const PetscSFNode *cells;
1203: PetscInt nFound;
1204: const PetscInt *found;
1206: PetscSFGetGraph(cellSF,NULL,&nFound,&found,&cells);
1207: .ve
1209: 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
1210: 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.
1212: .seealso: `DM`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMPointLocationType`
1213: @*/
1214: PetscErrorCode DMLocatePoints(DM dm, Vec v, DMPointLocationType ltype, PetscSF *cellSF)
1215: {
1216: PetscFunctionBegin;
1219: PetscAssertPointer(cellSF, 4);
1220: if (*cellSF) {
1221: PetscMPIInt result;
1224: PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)v), PetscObjectComm((PetscObject)*cellSF), &result));
1225: PetscCheck(result == MPI_IDENT || result == MPI_CONGRUENT, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "cellSF must have a communicator congruent to v's");
1226: } else {
1227: PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)v), cellSF));
1228: }
1229: PetscCall(PetscLogEventBegin(DM_LocatePoints, dm, 0, 0, 0));
1230: PetscUseTypeMethod(dm, locatepoints, v, ltype, *cellSF);
1231: PetscCall(PetscLogEventEnd(DM_LocatePoints, dm, 0, 0, 0));
1232: PetscFunctionReturn(PETSC_SUCCESS);
1233: }