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