xref: /petsc/src/dm/interface/dmcoordinates.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1 #include <petsc/private/dmimpl.h> /*I      "petscdm.h"          I*/
2 
3 #include <petscdmplex.h> /* For DMProjectCoordinates() */
4 #include <petscsf.h>     /* For DMLocatePoints() */
5 
6 PetscErrorCode DMRestrictHook_Coordinates(DM dm, DM dmc, void *ctx) {
7   DM  dm_coord, dmc_coord;
8   Vec coords, ccoords;
9   Mat inject;
10   PetscFunctionBegin;
11   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
12   PetscCall(DMGetCoordinateDM(dmc, &dmc_coord));
13   PetscCall(DMGetCoordinates(dm, &coords));
14   PetscCall(DMGetCoordinates(dmc, &ccoords));
15   if (coords && !ccoords) {
16     PetscCall(DMCreateGlobalVector(dmc_coord, &ccoords));
17     PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates"));
18     PetscCall(DMCreateInjection(dmc_coord, dm_coord, &inject));
19     PetscCall(MatRestrict(inject, coords, ccoords));
20     PetscCall(MatDestroy(&inject));
21     PetscCall(DMSetCoordinates(dmc, ccoords));
22     PetscCall(VecDestroy(&ccoords));
23   }
24   PetscFunctionReturn(0);
25 }
26 
27 static PetscErrorCode DMSubDomainHook_Coordinates(DM dm, DM subdm, void *ctx) {
28   DM          dm_coord, subdm_coord;
29   Vec         coords, ccoords, clcoords;
30   VecScatter *scat_i, *scat_g;
31   PetscFunctionBegin;
32   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
33   PetscCall(DMGetCoordinateDM(subdm, &subdm_coord));
34   PetscCall(DMGetCoordinates(dm, &coords));
35   PetscCall(DMGetCoordinates(subdm, &ccoords));
36   if (coords && !ccoords) {
37     PetscCall(DMCreateGlobalVector(subdm_coord, &ccoords));
38     PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates"));
39     PetscCall(DMCreateLocalVector(subdm_coord, &clcoords));
40     PetscCall(PetscObjectSetName((PetscObject)clcoords, "coordinates"));
41     PetscCall(DMCreateDomainDecompositionScatters(dm_coord, 1, &subdm_coord, NULL, &scat_i, &scat_g));
42     PetscCall(VecScatterBegin(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD));
43     PetscCall(VecScatterEnd(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD));
44     PetscCall(VecScatterBegin(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD));
45     PetscCall(VecScatterEnd(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD));
46     PetscCall(DMSetCoordinates(subdm, ccoords));
47     PetscCall(DMSetCoordinatesLocal(subdm, clcoords));
48     PetscCall(VecScatterDestroy(&scat_i[0]));
49     PetscCall(VecScatterDestroy(&scat_g[0]));
50     PetscCall(VecDestroy(&ccoords));
51     PetscCall(VecDestroy(&clcoords));
52     PetscCall(PetscFree(scat_i));
53     PetscCall(PetscFree(scat_g));
54   }
55   PetscFunctionReturn(0);
56 }
57 
58 /*@
59   DMGetCoordinateDM - Gets the DM that prescribes coordinate layout and scatters between global and local coordinates
60 
61   Collective on dm
62 
63   Input Parameter:
64 . dm - the DM
65 
66   Output Parameter:
67 . cdm - coordinate DM
68 
69   Level: intermediate
70 
71 .seealso: `DMSetCoordinateDM()`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`
72 @*/
73 PetscErrorCode DMGetCoordinateDM(DM dm, DM *cdm) {
74   PetscFunctionBegin;
75   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
76   PetscValidPointer(cdm, 2);
77   if (!dm->coordinates[0].dm) {
78     DM cdm;
79 
80     PetscUseTypeMethod(dm, createcoordinatedm, &cdm);
81     PetscCall(PetscObjectSetName((PetscObject)cdm, "coordinateDM"));
82     /* Just in case the DM sets the coordinate DM when creating it (DMP4est can do this, because it may not setup
83      * until the call to CreateCoordinateDM) */
84     PetscCall(DMDestroy(&dm->coordinates[0].dm));
85     dm->coordinates[0].dm = cdm;
86   }
87   *cdm = dm->coordinates[0].dm;
88   PetscFunctionReturn(0);
89 }
90 
91 /*@
92   DMSetCoordinateDM - Sets the DM that prescribes coordinate layout and scatters between global and local coordinates
93 
94   Logically Collective on dm
95 
96   Input Parameters:
97 + dm - the DM
98 - cdm - coordinate DM
99 
100   Level: intermediate
101 
102 .seealso: `DMGetCoordinateDM()`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`
103 @*/
104 PetscErrorCode DMSetCoordinateDM(DM dm, DM cdm) {
105   PetscFunctionBegin;
106   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
107   PetscValidHeaderSpecific(cdm, DM_CLASSID, 2);
108   PetscCall(PetscObjectReference((PetscObject)cdm));
109   PetscCall(DMDestroy(&dm->coordinates[0].dm));
110   dm->coordinates[0].dm = cdm;
111   PetscFunctionReturn(0);
112 }
113 
114 /*@
115   DMGetCellCoordinateDM - Gets the DM that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates
116 
117   Collective on dm
118 
119   Input Parameter:
120 . dm - the DM
121 
122   Output Parameter:
123 . cdm - cellwise coordinate DM, or NULL if they are not defined
124 
125   Note:
126   Call DMLocalizeCoordinates() to automatically create cellwise coordinates for periodic geometries.
127 
128   Level: intermediate
129 
130 .seealso: `DMSetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMLocalizeCoordinates()`
131 @*/
132 PetscErrorCode DMGetCellCoordinateDM(DM dm, DM *cdm) {
133   PetscFunctionBegin;
134   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
135   PetscValidPointer(cdm, 2);
136   *cdm = dm->coordinates[1].dm;
137   PetscFunctionReturn(0);
138 }
139 
140 /*@
141   DMSetCellCoordinateDM - Sets the DM that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates
142 
143   Logically Collective on dm
144 
145   Input Parameters:
146 + dm - the DM
147 - cdm - cellwise coordinate DM
148 
149   Level: intermediate
150 
151 .seealso: `DMGetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`
152 @*/
153 PetscErrorCode DMSetCellCoordinateDM(DM dm, DM cdm) {
154   PetscInt dim;
155 
156   PetscFunctionBegin;
157   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
158   if (cdm) {
159     PetscValidHeaderSpecific(cdm, DM_CLASSID, 2);
160     PetscCall(DMGetCoordinateDim(dm, &dim));
161     dm->coordinates[1].dim = dim;
162   }
163   PetscCall(PetscObjectReference((PetscObject)cdm));
164   PetscCall(DMDestroy(&dm->coordinates[1].dm));
165   dm->coordinates[1].dm = cdm;
166   PetscFunctionReturn(0);
167 }
168 
169 /*@
170   DMGetCoordinateDim - Retrieve the dimension of embedding space for coordinate values.
171 
172   Not Collective
173 
174   Input Parameter:
175 . dm - The DM object
176 
177   Output Parameter:
178 . dim - The embedding dimension
179 
180   Level: intermediate
181 
182 .seealso: `DMSetCoordinateDim()`, `DMGetCoordinateSection()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
183 @*/
184 PetscErrorCode DMGetCoordinateDim(DM dm, PetscInt *dim) {
185   PetscFunctionBegin;
186   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
187   PetscValidIntPointer(dim, 2);
188   if (dm->coordinates[0].dim == PETSC_DEFAULT) dm->coordinates[0].dim = dm->dim;
189   *dim = dm->coordinates[0].dim;
190   PetscFunctionReturn(0);
191 }
192 
193 /*@
194   DMSetCoordinateDim - Set the dimension of the embedding space for coordinate values.
195 
196   Not Collective
197 
198   Input Parameters:
199 + dm  - The DM object
200 - dim - The embedding dimension
201 
202   Level: intermediate
203 
204 .seealso: `DMGetCoordinateDim()`, `DMSetCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()`
205 @*/
206 PetscErrorCode DMSetCoordinateDim(DM dm, PetscInt dim) {
207   PetscDS  ds;
208   PetscInt Nds, n;
209 
210   PetscFunctionBegin;
211   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
212   dm->coordinates[0].dim = dim;
213   if (dm->dim >= 0) {
214     PetscCall(DMGetNumDS(dm, &Nds));
215     for (n = 0; n < Nds; ++n) {
216       PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds));
217       PetscCall(PetscDSSetCoordinateDimension(ds, dim));
218     }
219   }
220   PetscFunctionReturn(0);
221 }
222 
223 /*@
224   DMGetCoordinateSection - Retrieve the layout of coordinate values over the mesh.
225 
226   Collective on dm
227 
228   Input Parameter:
229 . dm - The DM object
230 
231   Output Parameter:
232 . section - The PetscSection object
233 
234   Level: intermediate
235 
236 .seealso: `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
237 @*/
238 PetscErrorCode DMGetCoordinateSection(DM dm, PetscSection *section) {
239   DM cdm;
240 
241   PetscFunctionBegin;
242   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
243   PetscValidPointer(section, 2);
244   PetscCall(DMGetCoordinateDM(dm, &cdm));
245   PetscCall(DMGetLocalSection(cdm, section));
246   PetscFunctionReturn(0);
247 }
248 
249 /*@
250   DMSetCoordinateSection - Set the layout of coordinate values over the mesh.
251 
252   Not Collective
253 
254   Input Parameters:
255 + dm      - The DM object
256 . dim     - The embedding dimension, or PETSC_DETERMINE
257 - section - The PetscSection object
258 
259   Level: intermediate
260 
261 .seealso: `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()`
262 @*/
263 PetscErrorCode DMSetCoordinateSection(DM dm, PetscInt dim, PetscSection section) {
264   DM cdm;
265 
266   PetscFunctionBegin;
267   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
268   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3);
269   PetscCall(DMGetCoordinateDM(dm, &cdm));
270   PetscCall(DMSetLocalSection(cdm, section));
271   if (dim == PETSC_DETERMINE) {
272     PetscInt d = PETSC_DEFAULT;
273     PetscInt pStart, pEnd, vStart, vEnd, v, dd;
274 
275     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
276     PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd));
277     pStart = PetscMax(vStart, pStart);
278     pEnd   = PetscMin(vEnd, pEnd);
279     for (v = pStart; v < pEnd; ++v) {
280       PetscCall(PetscSectionGetDof(section, v, &dd));
281       if (dd) {
282         d = dd;
283         break;
284       }
285     }
286     if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d));
287   }
288   PetscFunctionReturn(0);
289 }
290 
291 /*@
292   DMGetCellCoordinateSection - Retrieve the layout of cellwise coordinate values over the mesh.
293 
294   Collective on dm
295 
296   Input Parameter:
297 . dm - The DM object
298 
299   Output Parameter:
300 . section - The PetscSection object, or NULL if no cellwise coordinates are defined
301 
302   Level: intermediate
303 
304 .seealso: `DMGetCoordinateSection()`, `DMSetCellCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
305 @*/
306 PetscErrorCode DMGetCellCoordinateSection(DM dm, PetscSection *section) {
307   DM cdm;
308 
309   PetscFunctionBegin;
310   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
311   PetscValidPointer(section, 2);
312   *section = NULL;
313   PetscCall(DMGetCellCoordinateDM(dm, &cdm));
314   if (cdm) PetscCall(DMGetLocalSection(cdm, section));
315   PetscFunctionReturn(0);
316 }
317 
318 /*@
319   DMSetCellCoordinateSection - Set the layout of cellwise coordinate values over the mesh.
320 
321   Not Collective
322 
323   Input Parameters:
324 + dm      - The DM object
325 . dim     - The embedding dimension, or PETSC_DETERMINE
326 - section - The PetscSection object for a cellwise layout
327 
328   Level: intermediate
329 
330 .seealso: `DMSetCoordinateSection()`, `DMGetCellCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
331 @*/
332 PetscErrorCode DMSetCellCoordinateSection(DM dm, PetscInt dim, PetscSection section) {
333   DM cdm;
334 
335   PetscFunctionBegin;
336   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
337   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3);
338   PetscCall(DMGetCellCoordinateDM(dm, &cdm));
339   PetscCheck(cdm, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "No DM defined for cellwise coordinates");
340   PetscCall(DMSetLocalSection(cdm, section));
341   if (dim == PETSC_DETERMINE) {
342     PetscInt d = PETSC_DEFAULT;
343     PetscInt pStart, pEnd, vStart, vEnd, v, dd;
344 
345     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
346     PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd));
347     pStart = PetscMax(vStart, pStart);
348     pEnd   = PetscMin(vEnd, pEnd);
349     for (v = pStart; v < pEnd; ++v) {
350       PetscCall(PetscSectionGetDof(section, v, &dd));
351       if (dd) {
352         d = dd;
353         break;
354       }
355     }
356     if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d));
357   }
358   PetscFunctionReturn(0);
359 }
360 
361 /*@
362   DMGetCoordinates - Gets a global vector with the coordinates associated with the DM.
363 
364   Collective on dm
365 
366   Input Parameter:
367 . dm - the DM
368 
369   Output Parameter:
370 . c - global coordinate vector
371 
372   Note:
373   This is a borrowed reference, so the user should NOT destroy this vector. When the DM is
374   destroyed the array will no longer be valid.
375 
376   Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates).
377 
378   For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
379   and (x_0,y_0,z_0,x_1,y_1,z_1...)
380 
381   Level: intermediate
382 
383 .seealso: `DMSetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()`
384 @*/
385 PetscErrorCode DMGetCoordinates(DM dm, Vec *c) {
386   PetscFunctionBegin;
387   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
388   PetscValidPointer(c, 2);
389   if (!dm->coordinates[0].x && dm->coordinates[0].xl) {
390     DM cdm = NULL;
391 
392     PetscCall(DMGetCoordinateDM(dm, &cdm));
393     PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[0].x));
394     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].x, "coordinates"));
395     PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x));
396     PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x));
397   }
398   *c = dm->coordinates[0].x;
399   PetscFunctionReturn(0);
400 }
401 
402 /*@
403   DMSetCoordinates - Sets into the DM a global vector that holds the coordinates
404 
405   Collective on dm
406 
407   Input Parameters:
408 + dm - the DM
409 - c - coordinate vector
410 
411   Notes:
412   The coordinates do include those for ghost points, which are in the local vector.
413 
414   The vector c should be destroyed by the caller.
415 
416   Level: intermediate
417 
418 .seealso: `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()`
419 @*/
420 PetscErrorCode DMSetCoordinates(DM dm, Vec c) {
421   PetscFunctionBegin;
422   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
423   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
424   PetscCall(PetscObjectReference((PetscObject)c));
425   PetscCall(VecDestroy(&dm->coordinates[0].x));
426   dm->coordinates[0].x = c;
427   PetscCall(VecDestroy(&dm->coordinates[0].xl));
428   PetscCall(DMCoarsenHookAdd(dm, DMRestrictHook_Coordinates, NULL, NULL));
429   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_Coordinates, NULL, NULL));
430   PetscFunctionReturn(0);
431 }
432 
433 /*@
434   DMGetCellCoordinates - Gets a global vector with the cellwise coordinates associated with the DM.
435 
436   Collective on dm
437 
438   Input Parameter:
439 . dm - the DM
440 
441   Output Parameter:
442 . c - global coordinate vector
443 
444   Note:
445   This is a borrowed reference, so the user should NOT destroy this vector. When the DM is
446   destroyed the array will no longer be valid.
447 
448   Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates).
449 
450   Level: intermediate
451 
452 .seealso: `DMSetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()`
453 @*/
454 PetscErrorCode DMGetCellCoordinates(DM dm, Vec *c) {
455   PetscFunctionBegin;
456   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
457   PetscValidPointer(c, 2);
458   if (!dm->coordinates[1].x && dm->coordinates[1].xl) {
459     DM cdm = NULL;
460 
461     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
462     PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[1].x));
463     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].x, "DG coordinates"));
464     PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x));
465     PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x));
466   }
467   *c = dm->coordinates[1].x;
468   PetscFunctionReturn(0);
469 }
470 
471 /*@
472   DMSetCellCoordinates - Sets into the DM a global vector that holds the cellwise coordinates
473 
474   Collective on dm
475 
476   Input Parameters:
477 + dm - the DM
478 - c - cellwise coordinate vector
479 
480   Notes:
481   The coordinates do include those for ghost points, which are in the local vector.
482 
483   The vector c should be destroyed by the caller.
484 
485   Level: intermediate
486 
487 .seealso: `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()`
488 @*/
489 PetscErrorCode DMSetCellCoordinates(DM dm, Vec c) {
490   PetscFunctionBegin;
491   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
492   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
493   PetscCall(PetscObjectReference((PetscObject)c));
494   PetscCall(VecDestroy(&dm->coordinates[1].x));
495   dm->coordinates[1].x = c;
496   PetscCall(VecDestroy(&dm->coordinates[1].xl));
497   PetscFunctionReturn(0);
498 }
499 
500 /*@
501   DMGetCoordinatesLocalSetUp - Prepares a local vector of coordinates, so that DMGetCoordinatesLocalNoncollective() can be used as non-collective afterwards.
502 
503   Collective on dm
504 
505   Input Parameter:
506 . dm - the DM
507 
508   Level: advanced
509 
510 .seealso: `DMGetCoordinatesLocalNoncollective()`
511 @*/
512 PetscErrorCode DMGetCoordinatesLocalSetUp(DM dm) {
513   PetscFunctionBegin;
514   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
515   if (!dm->coordinates[0].xl && dm->coordinates[0].x) {
516     DM cdm = NULL;
517 
518     PetscCall(DMGetCoordinateDM(dm, &cdm));
519     PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[0].xl));
520     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].xl, "coordinates"));
521     PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
522     PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
523   }
524   PetscFunctionReturn(0);
525 }
526 
527 /*@
528   DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the DM.
529 
530   Collective on dm
531 
532   Input Parameter:
533 . dm - the DM
534 
535   Output Parameter:
536 . c - coordinate vector
537 
538   Note:
539   This is a borrowed reference, so the user should NOT destroy this vector
540 
541   Each process has the local and ghost coordinates
542 
543   For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
544   and (x_0,y_0,z_0,x_1,y_1,z_1...)
545 
546   Level: intermediate
547 
548 .seealso: `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`, `DMGetCoordinatesLocalNoncollective()`
549 @*/
550 PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c) {
551   PetscFunctionBegin;
552   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
553   PetscValidPointer(c, 2);
554   PetscCall(DMGetCoordinatesLocalSetUp(dm));
555   *c = dm->coordinates[0].xl;
556   PetscFunctionReturn(0);
557 }
558 
559 /*@
560   DMGetCoordinatesLocalNoncollective - Non-collective version of DMGetCoordinatesLocal(). Fails if global coordinates have been set and DMGetCoordinatesLocalSetUp() not called.
561 
562   Not collective
563 
564   Input Parameter:
565 . dm - the DM
566 
567   Output Parameter:
568 . c - coordinate vector
569 
570   Level: advanced
571 
572 .seealso: `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinatesLocal()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
573 @*/
574 PetscErrorCode DMGetCoordinatesLocalNoncollective(DM dm, Vec *c) {
575   PetscFunctionBegin;
576   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
577   PetscValidPointer(c, 2);
578   PetscCheck(dm->coordinates[0].xl || !dm->coordinates[0].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called");
579   *c = dm->coordinates[0].xl;
580   PetscFunctionReturn(0);
581 }
582 
583 /*@
584   DMGetCoordinatesLocalTuple - Gets a local vector with the coordinates of specified points and section describing its layout.
585 
586   Not collective
587 
588   Input Parameters:
589 + dm - the DM
590 - p - the IS of points whose coordinates will be returned
591 
592   Output Parameters:
593 + pCoordSection - the PetscSection describing the layout of pCoord, i.e. each point corresponds to one point in p, and DOFs correspond to coordinates
594 - pCoord - the Vec with coordinates of points in p
595 
596   Note:
597   DMGetCoordinatesLocalSetUp() must be called first. This function employs DMGetCoordinatesLocalNoncollective() so it is not collective.
598 
599   This creates a new vector, so the user SHOULD destroy this vector
600 
601   Each process has the local and ghost coordinates
602 
603   For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
604   and (x_0,y_0,z_0,x_1,y_1,z_1...)
605 
606   Level: advanced
607 
608 .seealso: `DMSetCoordinatesLocal()`, `DMGetCoordinatesLocal()`, `DMGetCoordinatesLocalNoncollective()`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
609 @*/
610 PetscErrorCode DMGetCoordinatesLocalTuple(DM dm, IS p, PetscSection *pCoordSection, Vec *pCoord) {
611   DM                 cdm;
612   PetscSection       cs, newcs;
613   Vec                coords;
614   const PetscScalar *arr;
615   PetscScalar       *newarr = NULL;
616   PetscInt           n;
617 
618   PetscFunctionBegin;
619   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
620   PetscValidHeaderSpecific(p, IS_CLASSID, 2);
621   if (pCoordSection) PetscValidPointer(pCoordSection, 3);
622   if (pCoord) PetscValidPointer(pCoord, 4);
623   PetscCall(DMGetCoordinateDM(dm, &cdm));
624   PetscCall(DMGetLocalSection(cdm, &cs));
625   PetscCall(DMGetCoordinatesLocal(dm, &coords));
626   PetscCheck(coords, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called or coordinates not set");
627   PetscCheck(cdm && cs, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM not supported");
628   PetscCall(VecGetArrayRead(coords, &arr));
629   PetscCall(PetscSectionExtractDofsFromArray(cs, MPIU_SCALAR, arr, p, &newcs, pCoord ? ((void **)&newarr) : NULL));
630   PetscCall(VecRestoreArrayRead(coords, &arr));
631   if (pCoord) {
632     PetscCall(PetscSectionGetStorageSize(newcs, &n));
633     /* set array in two steps to mimic PETSC_OWN_POINTER */
634     PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)p), 1, n, NULL, pCoord));
635     PetscCall(VecReplaceArray(*pCoord, newarr));
636   } else {
637     PetscCall(PetscFree(newarr));
638   }
639   if (pCoordSection) {
640     *pCoordSection = newcs;
641   } else PetscCall(PetscSectionDestroy(&newcs));
642   PetscFunctionReturn(0);
643 }
644 
645 /*@
646   DMSetCoordinatesLocal - Sets into the DM a local vector that holds the coordinates
647 
648   Not collective
649 
650    Input Parameters:
651 +  dm - the DM
652 -  c - coordinate vector
653 
654   Notes:
655   The coordinates of ghost points can be set using DMSetCoordinates()
656   followed by DMGetCoordinatesLocal(). This is intended to enable the
657   setting of ghost coordinates outside of the domain.
658 
659   The vector c should be destroyed by the caller.
660 
661   Level: intermediate
662 
663 .seealso: `DMGetCoordinatesLocal()`, `DMSetCoordinates()`, `DMGetCoordinates()`, `DMGetCoordinateDM()`
664 @*/
665 PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c) {
666   PetscFunctionBegin;
667   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
668   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
669   PetscCall(PetscObjectReference((PetscObject)c));
670   PetscCall(VecDestroy(&dm->coordinates[0].xl));
671   dm->coordinates[0].xl = c;
672   PetscCall(VecDestroy(&dm->coordinates[0].x));
673   PetscFunctionReturn(0);
674 }
675 
676 /*@
677   DMGetCellCoordinatesLocalSetUp - Prepares a local vector of cellwise coordinates, so that DMGetCellCoordinatesLocalNoncollective() can be used as non-collective afterwards.
678 
679   Collective on dm
680 
681   Input Parameter:
682 . dm - the DM
683 
684   Level: advanced
685 
686 .seealso: `DMGetCellCoordinatesLocalNoncollective()`
687 @*/
688 PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM dm) {
689   PetscFunctionBegin;
690   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
691   if (!dm->coordinates[1].xl && dm->coordinates[1].x) {
692     DM cdm = NULL;
693 
694     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
695     PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[1].xl));
696     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].xl, "DG coordinates"));
697     PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
698     PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
699   }
700   PetscFunctionReturn(0);
701 }
702 
703 /*@
704   DMGetCellCoordinatesLocal - Gets a local vector with the cellwise coordinates associated with the DM.
705 
706   Collective on dm
707 
708   Input Parameter:
709 . dm - the DM
710 
711   Output Parameter:
712 . c - coordinate vector
713 
714   Note:
715   This is a borrowed reference, so the user should NOT destroy this vector
716 
717   Each process has the local and ghost coordinates
718 
719   Level: intermediate
720 
721 .seealso: `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`, `DMGetCellCoordinatesLocalNoncollective()`
722 @*/
723 PetscErrorCode DMGetCellCoordinatesLocal(DM dm, Vec *c) {
724   PetscFunctionBegin;
725   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
726   PetscValidPointer(c, 2);
727   PetscCall(DMGetCellCoordinatesLocalSetUp(dm));
728   *c = dm->coordinates[1].xl;
729   PetscFunctionReturn(0);
730 }
731 
732 /*@
733   DMGetCellCoordinatesLocalNoncollective - Non-collective version of DMGetCellCoordinatesLocal(). Fails if global cellwise coordinates have been set and DMGetCellCoordinatesLocalSetUp() not called.
734 
735   Not collective
736 
737   Input Parameter:
738 . dm - the DM
739 
740   Output Parameter:
741 . c - cellwise coordinate vector
742 
743   Level: advanced
744 
745 .seealso: `DMGetCellCoordinatesLocalSetUp()`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`
746 @*/
747 PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM dm, Vec *c) {
748   PetscFunctionBegin;
749   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
750   PetscValidPointer(c, 2);
751   PetscCheck(dm->coordinates[1].xl || !dm->coordinates[1].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCellCoordinatesLocalSetUp() has not been called");
752   *c = dm->coordinates[1].xl;
753   PetscFunctionReturn(0);
754 }
755 
756 /*@
757   DMSetCellCoordinatesLocal - Sets into the DM a local vector that holds the cellwise coordinates
758 
759   Not collective
760 
761    Input Parameters:
762 +  dm - the DM
763 -  c - cellwise coordinate vector
764 
765   Notes:
766   The coordinates of ghost points can be set using DMSetCoordinates()
767   followed by DMGetCoordinatesLocal(). This is intended to enable the
768   setting of ghost coordinates outside of the domain.
769 
770   The vector c should be destroyed by the caller.
771 
772   Level: intermediate
773 
774 .seealso: `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinates()`, `DMGetCellCoordinates()`, `DMGetCellCoordinateDM()`
775 @*/
776 PetscErrorCode DMSetCellCoordinatesLocal(DM dm, Vec c) {
777   PetscFunctionBegin;
778   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
779   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
780   PetscCall(PetscObjectReference((PetscObject)c));
781   PetscCall(VecDestroy(&dm->coordinates[1].xl));
782   dm->coordinates[1].xl = c;
783   PetscCall(VecDestroy(&dm->coordinates[1].x));
784   PetscFunctionReturn(0);
785 }
786 
787 PetscErrorCode DMGetCoordinateField(DM dm, DMField *field) {
788   PetscFunctionBegin;
789   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
790   PetscValidPointer(field, 2);
791   if (!dm->coordinates[0].field) {
792     if (dm->ops->createcoordinatefield) PetscCall((*dm->ops->createcoordinatefield)(dm, &dm->coordinates[0].field));
793   }
794   *field = dm->coordinates[0].field;
795   PetscFunctionReturn(0);
796 }
797 
798 PetscErrorCode DMSetCoordinateField(DM dm, DMField field) {
799   PetscFunctionBegin;
800   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
801   if (field) PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 2);
802   PetscCall(PetscObjectReference((PetscObject)field));
803   PetscCall(DMFieldDestroy(&dm->coordinates[0].field));
804   dm->coordinates[0].field = field;
805   PetscFunctionReturn(0);
806 }
807 
808 /*@
809   DMGetLocalBoundingBox - Returns the bounding box for the piece of the DM on this process.
810 
811   Not collective
812 
813   Input Parameter:
814 . dm - the DM
815 
816   Output Parameters:
817 + lmin - local minimum coordinates (length coord dim, optional)
818 - lmax - local maximim coordinates (length coord dim, optional)
819 
820   Level: beginner
821 
822   Note: If the DM is a DMDA and has no coordinates, the index bounds are returned instead.
823 
824 .seealso: `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetBoundingBox()`
825 @*/
826 PetscErrorCode DMGetLocalBoundingBox(DM dm, PetscReal lmin[], PetscReal lmax[]) {
827   Vec                coords = NULL;
828   PetscReal          min[3] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL};
829   PetscReal          max[3] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
830   const PetscScalar *local_coords;
831   PetscInt           N, Ni;
832   PetscInt           cdim, i, j;
833 
834   PetscFunctionBegin;
835   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
836   PetscCall(DMGetCoordinateDim(dm, &cdim));
837   PetscCall(DMGetCoordinates(dm, &coords));
838   if (coords) {
839     PetscCall(VecGetArrayRead(coords, &local_coords));
840     PetscCall(VecGetLocalSize(coords, &N));
841     Ni = N / cdim;
842     for (i = 0; i < Ni; ++i) {
843       for (j = 0; j < 3; ++j) {
844         min[j] = j < cdim ? PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j])) : 0;
845         max[j] = j < cdim ? PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j])) : 0;
846       }
847     }
848     PetscCall(VecRestoreArrayRead(coords, &local_coords));
849   } else {
850     PetscBool isda;
851 
852     PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda));
853     if (isda) PetscCall(DMGetLocalBoundingIndices_DMDA(dm, min, max));
854   }
855   if (lmin) PetscCall(PetscArraycpy(lmin, min, cdim));
856   if (lmax) PetscCall(PetscArraycpy(lmax, max, cdim));
857   PetscFunctionReturn(0);
858 }
859 
860 /*@
861   DMGetBoundingBox - Returns the global bounding box for the DM.
862 
863   Collective
864 
865   Input Parameter:
866 . dm - the DM
867 
868   Output Parameters:
869 + gmin - global minimum coordinates (length coord dim, optional)
870 - gmax - global maximim coordinates (length coord dim, optional)
871 
872   Level: beginner
873 
874 .seealso: `DMGetLocalBoundingBox()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`
875 @*/
876 PetscErrorCode DMGetBoundingBox(DM dm, PetscReal gmin[], PetscReal gmax[]) {
877   PetscReal   lmin[3], lmax[3];
878   PetscInt    cdim;
879   PetscMPIInt count;
880 
881   PetscFunctionBegin;
882   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
883   PetscCall(DMGetCoordinateDim(dm, &cdim));
884   PetscCall(PetscMPIIntCast(cdim, &count));
885   PetscCall(DMGetLocalBoundingBox(dm, lmin, lmax));
886   if (gmin) PetscCall(MPIU_Allreduce(lmin, gmin, count, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm)));
887   if (gmax) PetscCall(MPIU_Allreduce(lmax, gmax, count, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)dm)));
888   PetscFunctionReturn(0);
889 }
890 
891 /*@
892   DMProjectCoordinates - Project coordinates to a different space
893 
894   Input Parameters:
895 + dm      - The DM object
896 - disc    - The new coordinate discretization or NULL to ensure a coordinate discretization exists
897 
898   Level: intermediate
899 
900 .seealso: `DMGetCoordinateField()`
901 @*/
902 PetscErrorCode DMProjectCoordinates(DM dm, PetscFE disc) {
903   PetscFE      discOld;
904   PetscClassId classid;
905   DM           cdmOld, cdmNew;
906   Vec          coordsOld, coordsNew;
907   Mat          matInterp;
908   PetscBool    same_space = PETSC_TRUE;
909 
910   PetscFunctionBegin;
911   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
912   if (disc) PetscValidHeaderSpecific(disc, PETSCFE_CLASSID, 2);
913 
914   PetscCall(DMGetCoordinateDM(dm, &cdmOld));
915   /* Check current discretization is compatible */
916   PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
917   PetscCall(PetscObjectGetClassId((PetscObject)discOld, &classid));
918   if (classid != PETSCFE_CLASSID) {
919     if (classid == PETSC_CONTAINER_CLASSID) {
920       PetscFE        feLinear;
921       DMPolytopeType ct;
922       PetscInt       dim, dE, cStart, cEnd;
923       PetscBool      simplex;
924 
925       /* Assume linear vertex coordinates */
926       PetscCall(DMGetDimension(dm, &dim));
927       PetscCall(DMGetCoordinateDim(dm, &dE));
928       PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
929       if (cEnd > cStart) {
930         PetscCall(DMPlexGetCellType(dm, cStart, &ct));
931         switch (ct) {
932         case DM_POLYTOPE_TRI_PRISM:
933         case DM_POLYTOPE_TRI_PRISM_TENSOR: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot autoamtically create coordinate space for prisms");
934         default: break;
935         }
936       }
937       PetscCall(DMPlexIsSimplex(dm, &simplex));
938       PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, dE, simplex, 1, -1, &feLinear));
939       PetscCall(DMSetField(cdmOld, 0, NULL, (PetscObject)feLinear));
940       PetscCall(PetscFEDestroy(&feLinear));
941       PetscCall(DMCreateDS(cdmOld));
942       PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
943     } else {
944       const char *discname;
945 
946       PetscCall(PetscObjectGetType((PetscObject)discOld, &discname));
947       SETERRQ(PetscObjectComm((PetscObject)discOld), PETSC_ERR_SUP, "Discretization type %s not supported", discname);
948     }
949   }
950   if (!disc) PetscFunctionReturn(0);
951   { // Check if the new space is the same as the old modulo quadrature
952     PetscDualSpace dsOld, ds;
953     PetscCall(PetscFEGetDualSpace(discOld, &dsOld));
954     PetscCall(PetscFEGetDualSpace(disc, &ds));
955     PetscCall(PetscDualSpaceEqual(dsOld, ds, &same_space));
956   }
957   /* Make a fresh clone of the coordinate DM */
958   PetscCall(DMClone(cdmOld, &cdmNew));
959   PetscCall(DMSetField(cdmNew, 0, NULL, (PetscObject)disc));
960   PetscCall(DMCreateDS(cdmNew));
961   PetscCall(DMGetCoordinates(dm, &coordsOld));
962   if (same_space) {
963     PetscCall(PetscObjectReference((PetscObject)coordsOld));
964     coordsNew = coordsOld;
965   } else { // Project the coordinate vector from old to new space
966     PetscCall(DMCreateGlobalVector(cdmNew, &coordsNew));
967     PetscCall(DMCreateInterpolation(cdmOld, cdmNew, &matInterp, NULL));
968     PetscCall(MatInterpolate(matInterp, coordsOld, coordsNew));
969     PetscCall(MatDestroy(&matInterp));
970   }
971   /* Set new coordinate structures */
972   PetscCall(DMSetCoordinateField(dm, NULL));
973   PetscCall(DMSetCoordinateDM(dm, cdmNew));
974   PetscCall(DMSetCoordinates(dm, coordsNew));
975   PetscCall(VecDestroy(&coordsNew));
976   PetscCall(DMDestroy(&cdmNew));
977   PetscFunctionReturn(0);
978 }
979 
980 /*@
981   DMLocatePoints - Locate the points in v in the mesh and return a PetscSF of the containing cells
982 
983   Collective on v (see explanation below)
984 
985   Input Parameters:
986 + dm - The DM
987 - ltype - The type of point location, e.g. DM_POINTLOCATION_NONE or DM_POINTLOCATION_NEAREST
988 
989   Input/Output Parameters:
990 + v - The Vec of points, on output contains the nearest mesh points to the given points if DM_POINTLOCATION_NEAREST is used
991 - cellSF - Points to either NULL, or a PetscSF with guesses for which cells contain each point;
992            on output, the PetscSF containing the ranks and local indices of the containing points
993 
994   Level: developer
995 
996   Notes:
997   To do a search of the local cells of the mesh, v should have PETSC_COMM_SELF as its communicator.
998   To do a search of all the cells in the distributed mesh, v should have the same communicator as dm.
999 
1000   If *cellSF is NULL on input, a PetscSF will be created.
1001   If *cellSF is not NULL on input, it should point to an existing PetscSF, whose graph will be used as initial guesses.
1002 
1003   An array that maps each point to its containing cell can be obtained with
1004 
1005 $    const PetscSFNode *cells;
1006 $    PetscInt           nFound;
1007 $    const PetscInt    *found;
1008 $
1009 $    PetscSFGetGraph(cellSF,NULL,&nFound,&found,&cells);
1010 
1011   Where cells[i].rank is the rank of the cell containing point found[i] (or i if found == NULL), and cells[i].index is
1012   the index of the cell in its rank's local numbering.
1013 
1014 .seealso: `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMPointLocationType`
1015 @*/
1016 PetscErrorCode DMLocatePoints(DM dm, Vec v, DMPointLocationType ltype, PetscSF *cellSF) {
1017   PetscFunctionBegin;
1018   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1019   PetscValidHeaderSpecific(v, VEC_CLASSID, 2);
1020   PetscValidPointer(cellSF, 4);
1021   if (*cellSF) {
1022     PetscMPIInt result;
1023 
1024     PetscValidHeaderSpecific(*cellSF, PETSCSF_CLASSID, 4);
1025     PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)v), PetscObjectComm((PetscObject)*cellSF), &result));
1026     PetscCheck(result == MPI_IDENT || result == MPI_CONGRUENT, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "cellSF must have a communicator congruent to v's");
1027   } else {
1028     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)v), cellSF));
1029   }
1030   PetscCall(PetscLogEventBegin(DM_LocatePoints, dm, 0, 0, 0));
1031   PetscUseTypeMethod(dm, locatepoints, v, ltype, *cellSF);
1032   PetscCall(PetscLogEventEnd(DM_LocatePoints, dm, 0, 0, 0));
1033   PetscFunctionReturn(0);
1034 }
1035