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