xref: /petsc/src/dm/impls/plex/plexgeometry.c (revision e9b685f588e4684eded62e6afe66c807c541da1d)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 
3 #undef __FUNCT__
4 #define __FUNCT__ "DMPlexGetLineIntersection_2D_Internal"
5 static PetscErrorCode DMPlexGetLineIntersection_2D_Internal(const PetscReal segmentA[], const PetscReal segmentB[], PetscReal intersection[], PetscBool *hasIntersection)
6 {
7   const PetscReal p0_x  = segmentA[0*2+0];
8   const PetscReal p0_y  = segmentA[0*2+1];
9   const PetscReal p1_x  = segmentA[1*2+0];
10   const PetscReal p1_y  = segmentA[1*2+1];
11   const PetscReal p2_x  = segmentB[0*2+0];
12   const PetscReal p2_y  = segmentB[0*2+1];
13   const PetscReal p3_x  = segmentB[1*2+0];
14   const PetscReal p3_y  = segmentB[1*2+1];
15   const PetscReal s1_x  = p1_x - p0_x;
16   const PetscReal s1_y  = p1_y - p0_y;
17   const PetscReal s2_x  = p3_x - p2_x;
18   const PetscReal s2_y  = p3_y - p2_y;
19   const PetscReal denom = (-s2_x * s1_y + s1_x * s2_y);
20 
21   PetscFunctionBegin;
22   *hasIntersection = PETSC_FALSE;
23   /* Non-parallel lines */
24   if (denom != 0.0) {
25     const PetscReal s = (-s1_y * (p0_x - p2_x) + s1_x * (p0_y - p2_y)) / denom;
26     const PetscReal t = ( s2_x * (p0_y - p2_y) - s2_y * (p0_x - p2_x)) / denom;
27 
28     if (s >= 0 && s <= 1 && t >= 0 && t <= 1) {
29       *hasIntersection = PETSC_TRUE;
30       if (intersection) {
31         intersection[0] = p0_x + (t * s1_x);
32         intersection[1] = p0_y + (t * s1_y);
33       }
34     }
35   }
36   PetscFunctionReturn(0);
37 }
38 
39 #undef __FUNCT__
40 #define __FUNCT__ "DMPlexLocatePoint_Simplex_2D_Internal"
41 static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
42 {
43   const PetscInt  embedDim = 2;
44   const PetscReal eps      = PETSC_SQRT_MACHINE_EPSILON;
45   PetscReal       x        = PetscRealPart(point[0]);
46   PetscReal       y        = PetscRealPart(point[1]);
47   PetscReal       v0[2], J[4], invJ[4], detJ;
48   PetscReal       xi, eta;
49   PetscErrorCode  ierr;
50 
51   PetscFunctionBegin;
52   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
53   xi  = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]);
54   eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]);
55 
56   if ((xi >= -eps) && (eta >= -eps) && (xi + eta <= 2.0+eps)) *cell = c;
57   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
58   PetscFunctionReturn(0);
59 }
60 
61 #undef __FUNCT__
62 #define __FUNCT__ "DMPlexClosestPoint_Simplex_2D_Internal"
63 static PetscErrorCode DMPlexClosestPoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscReal cpoint[])
64 {
65   const PetscInt  embedDim = 2;
66   PetscReal       x        = PetscRealPart(point[0]);
67   PetscReal       y        = PetscRealPart(point[1]);
68   PetscReal       v0[2], J[4], invJ[4], detJ;
69   PetscReal       xi, eta, r;
70   PetscErrorCode  ierr;
71 
72   PetscFunctionBegin;
73   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
74   xi  = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]);
75   eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]);
76 
77   xi  = PetscMax(xi,  0.0);
78   eta = PetscMax(eta, 0.0);
79   r   = (xi + eta)/2.0;
80   if (xi + eta > 2.0) {
81     r    = (xi + eta)/2.0;
82     xi  /= r;
83     eta /= r;
84   }
85   cpoint[0] = J[0*embedDim+0]*xi + J[0*embedDim+1]*eta + v0[0];
86   cpoint[1] = J[1*embedDim+0]*xi + J[1*embedDim+1]*eta + v0[1];
87   PetscFunctionReturn(0);
88 }
89 
90 #undef __FUNCT__
91 #define __FUNCT__ "DMPlexLocatePoint_General_2D_Internal"
92 static PetscErrorCode DMPlexLocatePoint_General_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
93 {
94   PetscSection       coordSection;
95   Vec             coordsLocal;
96   PetscScalar    *coords = NULL;
97   const PetscInt  faces[8]  = {0, 1, 1, 2, 2, 3, 3, 0};
98   PetscReal       x         = PetscRealPart(point[0]);
99   PetscReal       y         = PetscRealPart(point[1]);
100   PetscInt        crossings = 0, f;
101   PetscErrorCode  ierr;
102 
103   PetscFunctionBegin;
104   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
105   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
106   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
107   for (f = 0; f < 4; ++f) {
108     PetscReal x_i   = PetscRealPart(coords[faces[2*f+0]*2+0]);
109     PetscReal y_i   = PetscRealPart(coords[faces[2*f+0]*2+1]);
110     PetscReal x_j   = PetscRealPart(coords[faces[2*f+1]*2+0]);
111     PetscReal y_j   = PetscRealPart(coords[faces[2*f+1]*2+1]);
112     PetscReal slope = (y_j - y_i) / (x_j - x_i);
113     PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE;
114     PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE;
115     PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE;
116     if ((cond1 || cond2)  && above) ++crossings;
117   }
118   if (crossings % 2) *cell = c;
119   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
120   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
121   PetscFunctionReturn(0);
122 }
123 
124 #undef __FUNCT__
125 #define __FUNCT__ "DMPlexLocatePoint_Simplex_3D_Internal"
126 static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
127 {
128   const PetscInt embedDim = 3;
129   PetscReal      v0[3], J[9], invJ[9], detJ;
130   PetscReal      x = PetscRealPart(point[0]);
131   PetscReal      y = PetscRealPart(point[1]);
132   PetscReal      z = PetscRealPart(point[2]);
133   PetscReal      xi, eta, zeta;
134   PetscErrorCode ierr;
135 
136   PetscFunctionBegin;
137   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
138   xi   = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]);
139   eta  = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]);
140   zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]);
141 
142   if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) *cell = c;
143   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
144   PetscFunctionReturn(0);
145 }
146 
147 #undef __FUNCT__
148 #define __FUNCT__ "DMPlexLocatePoint_General_3D_Internal"
149 static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
150 {
151   PetscSection   coordSection;
152   Vec            coordsLocal;
153   PetscScalar   *coords;
154   const PetscInt faces[24] = {0, 3, 2, 1,  5, 4, 7, 6,  3, 0, 4, 5,
155                               1, 2, 6, 7,  3, 5, 6, 2,  0, 1, 7, 4};
156   PetscBool      found = PETSC_TRUE;
157   PetscInt       f;
158   PetscErrorCode ierr;
159 
160   PetscFunctionBegin;
161   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
162   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
163   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
164   for (f = 0; f < 6; ++f) {
165     /* Check the point is under plane */
166     /*   Get face normal */
167     PetscReal v_i[3];
168     PetscReal v_j[3];
169     PetscReal normal[3];
170     PetscReal pp[3];
171     PetscReal dot;
172 
173     v_i[0]    = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]);
174     v_i[1]    = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]);
175     v_i[2]    = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]);
176     v_j[0]    = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]);
177     v_j[1]    = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]);
178     v_j[2]    = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]);
179     normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1];
180     normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2];
181     normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0];
182     pp[0]     = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]);
183     pp[1]     = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]);
184     pp[2]     = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]);
185     dot       = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2];
186 
187     /* Check that projected point is in face (2D location problem) */
188     if (dot < 0.0) {
189       found = PETSC_FALSE;
190       break;
191     }
192   }
193   if (found) *cell = c;
194   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
195   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
196   PetscFunctionReturn(0);
197 }
198 
199 #undef __FUNCT__
200 #define __FUNCT__ "PetscGridHashInitialize_Internal"
201 static PetscErrorCode PetscGridHashInitialize_Internal(PetscGridHash box, PetscInt dim, const PetscScalar point[])
202 {
203   PetscInt d;
204 
205   PetscFunctionBegin;
206   box->dim = dim;
207   for (d = 0; d < dim; ++d) box->lower[d] = box->upper[d] = PetscRealPart(point[d]);
208   PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNCT__
212 #define __FUNCT__ "PetscGridHashCreate"
213 PetscErrorCode PetscGridHashCreate(MPI_Comm comm, PetscInt dim, const PetscScalar point[], PetscGridHash *box)
214 {
215   PetscErrorCode ierr;
216 
217   PetscFunctionBegin;
218   ierr = PetscMalloc1(1, box);CHKERRQ(ierr);
219   ierr = PetscGridHashInitialize_Internal(*box, dim, point);CHKERRQ(ierr);
220   PetscFunctionReturn(0);
221 }
222 
223 #undef __FUNCT__
224 #define __FUNCT__ "PetscGridHashEnlarge"
225 PetscErrorCode PetscGridHashEnlarge(PetscGridHash box, const PetscScalar point[])
226 {
227   PetscInt d;
228 
229   PetscFunctionBegin;
230   for (d = 0; d < box->dim; ++d) {
231     box->lower[d] = PetscMin(box->lower[d], PetscRealPart(point[d]));
232     box->upper[d] = PetscMax(box->upper[d], PetscRealPart(point[d]));
233   }
234   PetscFunctionReturn(0);
235 }
236 
237 #undef __FUNCT__
238 #define __FUNCT__ "PetscGridHashSetGrid"
239 /*
240   PetscGridHashSetGrid - Divide the grid into boxes
241 
242   Not collective
243 
244   Input Parameters:
245 + box - The grid hash object
246 . n   - The number of boxes in each dimension, or PETSC_DETERMINE
247 - h   - The box size in each dimension, only used if n[d] == PETSC_DETERMINE
248 
249   Level: developer
250 
251 .seealso: PetscGridHashCreate()
252 */
253 PetscErrorCode PetscGridHashSetGrid(PetscGridHash box, const PetscInt n[], const PetscReal h[])
254 {
255   PetscInt d;
256 
257   PetscFunctionBegin;
258   for (d = 0; d < box->dim; ++d) {
259     box->extent[d] = box->upper[d] - box->lower[d];
260     if (n[d] == PETSC_DETERMINE) {
261       box->h[d] = h[d];
262       box->n[d] = PetscCeilReal(box->extent[d]/h[d]);
263     } else {
264       box->n[d] = n[d];
265       box->h[d] = box->extent[d]/n[d];
266     }
267   }
268   PetscFunctionReturn(0);
269 }
270 
271 #undef __FUNCT__
272 #define __FUNCT__ "PetscGridHashGetEnclosingBox"
273 /*
274   PetscGridHashGetEnclosingBox - Find the grid boxes containing each input point
275 
276   Not collective
277 
278   Input Parameters:
279 + box       - The grid hash object
280 . numPoints - The number of input points
281 - points    - The input point coordinates
282 
283   Output Parameters:
284 + dboxes    - An array of numPoints*dim integers expressing the enclosing box as (i_0, i_1, ..., i_dim)
285 - boxes     - An array of numPoints integers expressing the enclosing box as single number, or NULL
286 
287   Level: developer
288 
289 .seealso: PetscGridHashCreate()
290 */
291 PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[])
292 {
293   const PetscReal *lower = box->lower;
294   const PetscReal *upper = box->upper;
295   const PetscReal *h     = box->h;
296   const PetscInt  *n     = box->n;
297   const PetscInt   dim   = box->dim;
298   PetscInt         d, p;
299 
300   PetscFunctionBegin;
301   for (p = 0; p < numPoints; ++p) {
302     for (d = 0; d < dim; ++d) {
303       PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]);
304 
305       if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1;
306       if (dbox < 0 || dbox >= n[d]) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input point %d (%g, %g, %g) is outside of our bounding box",
307                                              p, PetscRealPart(points[p*dim+0]), dim > 1 ? PetscRealPart(points[p*dim+1]) : 0.0, dim > 2 ? PetscRealPart(points[p*dim+2]) : 0.0);
308       dboxes[p*dim+d] = dbox;
309     }
310     if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1];
311   }
312   PetscFunctionReturn(0);
313 }
314 
315 #undef __FUNCT__
316 #define __FUNCT__ "PetscGridHashDestroy"
317 PetscErrorCode PetscGridHashDestroy(PetscGridHash *box)
318 {
319   PetscErrorCode ierr;
320 
321   PetscFunctionBegin;
322   if (*box) {
323     ierr = PetscSectionDestroy(&(*box)->cellSection);CHKERRQ(ierr);
324     ierr = ISDestroy(&(*box)->cells);CHKERRQ(ierr);
325     ierr = DMLabelDestroy(&(*box)->cellsSparse);CHKERRQ(ierr);
326   }
327   ierr = PetscFree(*box);CHKERRQ(ierr);
328   PetscFunctionReturn(0);
329 }
330 
331 #undef __FUNCT__
332 #define __FUNCT__ "DMPlexLocatePoint_Internal"
333 PetscErrorCode DMPlexLocatePoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cellStart, PetscInt *cell)
334 {
335   PetscInt       coneSize;
336   PetscErrorCode ierr;
337 
338   PetscFunctionBegin;
339   switch (dim) {
340   case 2:
341     ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr);
342     switch (coneSize) {
343     case 3:
344       ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
345       break;
346     case 4:
347       ierr = DMPlexLocatePoint_General_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
348       break;
349     default:
350       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize);
351     }
352     break;
353   case 3:
354     ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr);
355     switch (coneSize) {
356     case 4:
357       ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
358       break;
359     case 6:
360       ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
361       break;
362     default:
363       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize);
364     }
365     break;
366   default:
367     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %D", dim);
368   }
369   PetscFunctionReturn(0);
370 }
371 
372 #undef __FUNCT__
373 #define __FUNCT__ "DMPlexClosestPoint_Internal"
374 /*
375   DMPlexClosestPoint_Internal - Returns the closest point in the cell to the given point
376 */
377 PetscErrorCode DMPlexClosestPoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cell, PetscReal cpoint[])
378 {
379   PetscInt       coneSize;
380   PetscErrorCode ierr;
381 
382   PetscFunctionBegin;
383   switch (dim) {
384   case 2:
385     ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
386     switch (coneSize) {
387     case 3:
388       ierr = DMPlexClosestPoint_Simplex_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
389       break;
390 #if 0
391     case 4:
392       ierr = DMPlexClosestPoint_General_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
393       break;
394 #endif
395     default:
396       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize);
397     }
398     break;
399 #if 0
400   case 3:
401     ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
402     switch (coneSize) {
403     case 4:
404       ierr = DMPlexClosestPoint_Simplex_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
405       break;
406     case 6:
407       ierr = DMPlexClosestPoint_General_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
408       break;
409     default:
410       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize);
411     }
412     break;
413 #endif
414   default:
415     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for mesh dimension %D", dim);
416   }
417   PetscFunctionReturn(0);
418 }
419 
420 #undef __FUNCT__
421 #define __FUNCT__ "DMPlexComputeGridHash_Internal"
422 /*
423   DMPlexComputeGridHash_Internal - Create a grid hash structure covering the Plex
424 
425   Collective on DM
426 
427   Input Parameter:
428 . dm - The Plex
429 
430   Output Parameter:
431 . localBox - The grid hash object
432 
433   Level: developer
434 
435 .seealso: PetscGridHashCreate(), PetscGridHashGetEnclosingBox()
436 */
437 PetscErrorCode DMPlexComputeGridHash_Internal(DM dm, PetscGridHash *localBox)
438 {
439   MPI_Comm           comm;
440   PetscGridHash      lbox;
441   Vec                coordinates;
442   PetscSection       coordSection;
443   Vec                coordsLocal;
444   const PetscScalar *coords;
445   PetscInt          *dboxes, *boxes;
446   PetscInt           n[3] = {10, 10, 10};
447   PetscInt           dim, N, cStart, cEnd, cMax, c, i;
448   PetscErrorCode     ierr;
449 
450   PetscFunctionBegin;
451   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
452   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
453   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
454   if (dim != 2) SETERRQ(comm, PETSC_ERR_SUP, "I have only coded this for 2D");
455   ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr);
456   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
457   ierr = PetscGridHashCreate(comm, dim, coords, &lbox);CHKERRQ(ierr);
458   for (i = 0; i < N; i += dim) {ierr = PetscGridHashEnlarge(lbox, &coords[i]);CHKERRQ(ierr);}
459   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
460   ierr = PetscGridHashSetGrid(lbox, n, NULL);CHKERRQ(ierr);
461 #if 0
462   /* Could define a custom reduction to merge these */
463   ierr = MPIU_Allreduce(lbox->lower, gbox->lower, 3, MPIU_REAL, MPI_MIN, comm);CHKERRQ(ierr);
464   ierr = MPIU_Allreduce(lbox->upper, gbox->upper, 3, MPIU_REAL, MPI_MAX, comm);CHKERRQ(ierr);
465 #endif
466   /* Is there a reason to snap the local bounding box to a division of the global box? */
467   /* Should we compute all overlaps of local boxes? We could do this with a rendevouz scheme partitioning the global box */
468   /* Create label */
469   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
470   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
471   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
472   ierr = DMLabelCreate("cells", &lbox->cellsSparse);CHKERRQ(ierr);
473   ierr = DMLabelCreateIndex(lbox->cellsSparse, cStart, cEnd);CHKERRQ(ierr);
474   /* Compute boxes which overlap each cell: http://stackoverflow.com/questions/13790208/triangle-square-intersection-test-in-2d */
475   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
476   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
477   ierr = PetscCalloc2(16 * dim, &dboxes, 16, &boxes);CHKERRQ(ierr);
478   for (c = cStart; c < cEnd; ++c) {
479     const PetscReal *h       = lbox->h;
480     PetscScalar     *ccoords = NULL;
481     PetscInt         csize   = 0;
482     PetscScalar      point[3];
483     PetscInt         dlim[6], d, e, i, j, k;
484 
485     /* Find boxes enclosing each vertex */
486     ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, &csize, &ccoords);CHKERRQ(ierr);
487     ierr = PetscGridHashGetEnclosingBox(lbox, csize/dim, ccoords, dboxes, boxes);CHKERRQ(ierr);
488     /* Mark cells containing the vertices */
489     for (e = 0; e < csize/dim; ++e) {ierr = DMLabelSetValue(lbox->cellsSparse, c, boxes[e]);CHKERRQ(ierr);}
490     /* Get grid of boxes containing these */
491     for (d = 0;   d < dim; ++d) {dlim[d*2+0] = dlim[d*2+1] = dboxes[d];}
492     for (d = dim; d < 3;   ++d) {dlim[d*2+0] = dlim[d*2+1] = 0;}
493     for (e = 1; e < dim+1; ++e) {
494       for (d = 0; d < dim; ++d) {
495         dlim[d*2+0] = PetscMin(dlim[d*2+0], dboxes[e*dim+d]);
496         dlim[d*2+1] = PetscMax(dlim[d*2+1], dboxes[e*dim+d]);
497       }
498     }
499     /* Check for intersection of box with cell */
500     for (k = dlim[2*2+0], point[2] = lbox->lower[2] + k*h[2]; k <= dlim[2*2+1]; ++k, point[2] += h[2]) {
501       for (j = dlim[1*2+0], point[1] = lbox->lower[1] + j*h[1]; j <= dlim[1*2+1]; ++j, point[1] += h[1]) {
502         for (i = dlim[0*2+0], point[0] = lbox->lower[0] + i*h[0]; i <= dlim[0*2+1]; ++i, point[0] += h[0]) {
503           const PetscInt box = (k*lbox->n[1] + j)*lbox->n[0] + i;
504           PetscScalar    cpoint[3];
505           PetscInt       cell, edge, ii, jj, kk;
506 
507           /* Check whether cell contains any vertex of these subboxes TODO vectorize this */
508           for (kk = 0, cpoint[2] = point[2]; kk < (dim > 2 ? 2 : 1); ++kk, cpoint[2] += h[2]) {
509             for (jj = 0, cpoint[1] = point[1]; jj < (dim > 1 ? 2 : 1); ++jj, cpoint[1] += h[1]) {
510               for (ii = 0, cpoint[0] = point[0]; ii < 2; ++ii, cpoint[0] += h[0]) {
511 
512                 ierr = DMPlexLocatePoint_Internal(dm, dim, cpoint, c, &cell);CHKERRQ(ierr);
513                 if (cell >= 0) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); ii = jj = kk = 2;}
514               }
515             }
516           }
517           /* Check whether cell edge intersects any edge of these subboxes TODO vectorize this */
518           for (edge = 0; edge < dim+1; ++edge) {
519             PetscReal segA[6], segB[6];
520 
521             for (d = 0; d < dim; ++d) {segA[d] = PetscRealPart(ccoords[edge*dim+d]); segA[dim+d] = PetscRealPart(ccoords[((edge+1)%(dim+1))*dim+d]);}
522             for (kk = 0; kk < (dim > 2 ? 2 : 1); ++kk) {
523               if (dim > 2) {segB[2]     = PetscRealPart(point[2]);
524                             segB[dim+2] = PetscRealPart(point[2]) + kk*h[2];}
525               for (jj = 0; jj < (dim > 1 ? 2 : 1); ++jj) {
526                 if (dim > 1) {segB[1]     = PetscRealPart(point[1]);
527                               segB[dim+1] = PetscRealPart(point[1]) + jj*h[1];}
528                 for (ii = 0; ii < 2; ++ii) {
529                   PetscBool intersects;
530 
531                   segB[0]     = PetscRealPart(point[0]);
532                   segB[dim+0] = PetscRealPart(point[0]) + ii*h[0];
533                   ierr = DMPlexGetLineIntersection_2D_Internal(segA, segB, NULL, &intersects);CHKERRQ(ierr);
534                   if (intersects) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); edge = ii = jj = kk = dim+1;}
535                 }
536               }
537             }
538           }
539         }
540       }
541     }
542     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &ccoords);CHKERRQ(ierr);
543   }
544   ierr = PetscFree2(dboxes, boxes);CHKERRQ(ierr);
545   ierr = DMLabelConvertToSection(lbox->cellsSparse, &lbox->cellSection, &lbox->cells);CHKERRQ(ierr);
546   ierr = DMLabelDestroy(&lbox->cellsSparse);CHKERRQ(ierr);
547   *localBox = lbox;
548   PetscFunctionReturn(0);
549 }
550 
551 #undef __FUNCT__
552 #define __FUNCT__ "DMLocatePoints_Plex"
553 PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, DMPointLocationType ltype, PetscSF cellSF)
554 {
555   DM_Plex        *mesh = (DM_Plex *) dm->data;
556   PetscBool       hash = mesh->useHashLocation;
557   PetscInt        bs, numPoints, p, numFound, *found = NULL;
558   PetscInt        dim, cStart, cEnd, cMax, numCells, c;
559   const PetscInt *boxCells;
560   PetscSFNode    *cells;
561   PetscScalar    *a;
562   PetscMPIInt     result;
563   PetscErrorCode  ierr;
564 
565   PetscFunctionBegin;
566   if (ltype == DM_POINTLOCATION_NEAREST && !hash) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Nearest point location only supported with grid hashing.");
567   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
568   ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr);
569   ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)cellSF),PETSC_COMM_SELF,&result);CHKERRQ(ierr);
570   if (result != MPI_IDENT && result != MPI_CONGRUENT) SETERRQ(PetscObjectComm((PetscObject)cellSF),PETSC_ERR_SUP, "Trying parallel point location: only local point location supported");
571   if (bs != dim) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Block size for point vector %D must be the mesh coordinate dimension %D", bs, dim);
572   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
573   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
574   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
575   ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr);
576   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
577   numPoints /= bs;
578   ierr = PetscMalloc1(numPoints, &cells);CHKERRQ(ierr);
579   if (hash) {
580     if (!mesh->lbox) {ierr = PetscInfo(dm, "Initializing grid hashing");CHKERRQ(ierr);ierr = DMPlexComputeGridHash_Internal(dm, &mesh->lbox);CHKERRQ(ierr);}
581     /* Designate the local box for each point */
582     /* Send points to correct process */
583     /* Search cells that lie in each subbox */
584     /*   Should we bin points before doing search? */
585     ierr = ISGetIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr);
586   }
587   for (p = 0, numFound = 0; p < numPoints; ++p) {
588     const PetscScalar *point = &a[p*bs];
589     PetscInt           dbin[3], bin, cell = -1, cellOffset;
590 
591     cells[p].rank  = 0;
592     cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
593     if (hash) {
594       ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr);
595       /* TODO Lay an interface over this so we can switch between Section (dense) and Label (sparse) */
596       ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr);
597       ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr);
598       for (c = cellOffset; c < cellOffset + numCells; ++c) {
599         ierr = DMPlexLocatePoint_Internal(dm, dim, point, boxCells[c], &cell);CHKERRQ(ierr);
600         if (cell >= 0) {
601           cells[p].rank = 0;
602           cells[p].index = cell;
603           numFound++;
604           break;
605         }
606       }
607     } else {
608       for (c = cStart; c < cEnd; ++c) {
609         ierr = DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);CHKERRQ(ierr);
610         if (cell >= 0) {
611           cells[p].rank = 0;
612           cells[p].index = cell;
613           numFound++;
614           break;
615         }
616       }
617     }
618   }
619   if (hash) {ierr = ISRestoreIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr);}
620   if (ltype == DM_POINTLOCATION_NEAREST && hash && numFound < numPoints) {
621     for (p = 0; p < numPoints; p++) {
622       const PetscScalar *point = &a[p*bs];
623       PetscReal          cpoint[3], diff[3], dist, distMax = PETSC_MAX_REAL;
624       PetscInt           dbin[3], bin, cellOffset, d;
625 
626       if (cells[p].index < 0) {
627         ++numFound;
628         ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr);
629         ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr);
630         ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr);
631         for (c = cellOffset; c < cellOffset + numCells; ++c) {
632           ierr = DMPlexClosestPoint_Internal(dm, dim, point, boxCells[c], cpoint);CHKERRQ(ierr);
633           for (d = 0; d < dim; ++d) diff[d] = cpoint[d] - PetscRealPart(point[d]);
634           dist = DMPlex_NormD_Internal(dim, diff);
635           if (dist < distMax) {
636             for (d = 0; d < dim; ++d) a[p*bs+d] = cpoint[d];
637             cells[p].rank  = 0;
638             cells[p].index = boxCells[c];
639             distMax = dist;
640             break;
641           }
642         }
643       }
644     }
645   }
646   /* This code is only be relevant when interfaced to parallel point location */
647   /* Check for highest numbered proc that claims a point (do we care?) */
648   if (ltype == DM_POINTLOCATION_REMOVE && numFound < numPoints) {
649     ierr = PetscMalloc1(numFound,&found);CHKERRQ(ierr);
650     for (p = 0, numFound = 0; p < numPoints; p++) {
651       if (cells[p].rank >= 0 && cells[p].index >= 0) {
652         if (numFound < p) {
653           cells[numFound] = cells[p];
654         }
655         found[numFound++] = p;
656       }
657     }
658   }
659   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
660   ierr = PetscSFSetGraph(cellSF, cEnd - cStart, numFound, found, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER);CHKERRQ(ierr);
661   PetscFunctionReturn(0);
662 }
663 
664 #undef __FUNCT__
665 #define __FUNCT__ "DMPlexComputeProjection2Dto1D_Internal"
666 /*
667   DMPlexComputeProjection2Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 2D
668 */
669 PetscErrorCode DMPlexComputeProjection2Dto1D_Internal(PetscScalar coords[], PetscReal R[])
670 {
671   const PetscReal x = PetscRealPart(coords[2] - coords[0]);
672   const PetscReal y = PetscRealPart(coords[3] - coords[1]);
673   const PetscReal r = PetscSqrtReal(x*x + y*y), c = x/r, s = y/r;
674 
675   PetscFunctionBegin;
676   R[0] = c; R[1] = -s;
677   R[2] = s; R[3] =  c;
678   coords[0] = 0.0;
679   coords[1] = r;
680   PetscFunctionReturn(0);
681 }
682 
683 #undef __FUNCT__
684 #define __FUNCT__ "DMPlexComputeProjection3Dto1D_Internal"
685 /*
686   DMPlexComputeProjection3Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 3D
687 
688   This uses the basis completion described by Frisvad,
689 
690   http://www.imm.dtu.dk/~jerf/papers/abstracts/onb.html
691   DOI:10.1080/2165347X.2012.689606
692 */
693 PetscErrorCode DMPlexComputeProjection3Dto1D_Internal(PetscScalar coords[], PetscReal R[])
694 {
695   PetscReal      x    = PetscRealPart(coords[3] - coords[0]);
696   PetscReal      y    = PetscRealPart(coords[4] - coords[1]);
697   PetscReal      z    = PetscRealPart(coords[5] - coords[2]);
698   PetscReal      r    = PetscSqrtReal(x*x + y*y + z*z);
699   PetscReal      rinv = 1. / r;
700   PetscFunctionBegin;
701 
702   x *= rinv; y *= rinv; z *= rinv;
703   if (x > 0.) {
704     PetscReal inv1pX   = 1./ (1. + x);
705 
706     R[0] = x; R[1] = -y;              R[2] = -z;
707     R[3] = y; R[4] = 1. - y*y*inv1pX; R[5] =     -y*z*inv1pX;
708     R[6] = z; R[7] =     -y*z*inv1pX; R[8] = 1. - z*z*inv1pX;
709   }
710   else {
711     PetscReal inv1mX   = 1./ (1. - x);
712 
713     R[0] = x; R[1] = z;               R[2] = y;
714     R[3] = y; R[4] =     -y*z*inv1mX; R[5] = 1. - y*y*inv1mX;
715     R[6] = z; R[7] = 1. - z*z*inv1mX; R[8] =     -y*z*inv1mX;
716   }
717   coords[0] = 0.0;
718   coords[1] = r;
719   PetscFunctionReturn(0);
720 }
721 
722 #undef __FUNCT__
723 #define __FUNCT__ "DMPlexComputeProjection3Dto2D_Internal"
724 /*
725   DMPlexComputeProjection3Dto2D_Internal - Rewrite coordinates to be the 2D projection of the 3D
726 */
727 PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscInt coordSize, PetscScalar coords[], PetscReal R[])
728 {
729   PetscReal      x1[3],  x2[3], n[3], norm;
730   PetscReal      x1p[3], x2p[3], xnp[3];
731   PetscReal      sqrtz, alpha;
732   const PetscInt dim = 3;
733   PetscInt       d, e, p;
734 
735   PetscFunctionBegin;
736   /* 0) Calculate normal vector */
737   for (d = 0; d < dim; ++d) {
738     x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]);
739     x2[d] = PetscRealPart(coords[2*dim+d] - coords[0*dim+d]);
740   }
741   n[0] = x1[1]*x2[2] - x1[2]*x2[1];
742   n[1] = x1[2]*x2[0] - x1[0]*x2[2];
743   n[2] = x1[0]*x2[1] - x1[1]*x2[0];
744   norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]);
745   n[0] /= norm;
746   n[1] /= norm;
747   n[2] /= norm;
748   /* 1) Take the normal vector and rotate until it is \hat z
749 
750     Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then
751 
752     R = /  alpha nx nz  alpha ny nz -1/alpha \
753         | -alpha ny     alpha nx        0    |
754         \     nx            ny         nz    /
755 
756     will rotate the normal vector to \hat z
757   */
758   sqrtz = PetscSqrtReal(1.0 - n[2]*n[2]);
759   /* Check for n = z */
760   if (sqrtz < 1.0e-10) {
761     const PetscInt s = PetscSign(n[2]);
762     /* If nz < 0, rotate 180 degrees around x-axis */
763     for (p = 3; p < coordSize/3; ++p) {
764       coords[p*2+0] = PetscRealPart(coords[p*dim+0] - coords[0*dim+0]);
765       coords[p*2+1] = (PetscRealPart(coords[p*dim+1] - coords[0*dim+1])) * s;
766     }
767     coords[0] = 0.0;
768     coords[1] = 0.0;
769     coords[2] = x1[0];
770     coords[3] = x1[1] * s;
771     coords[4] = x2[0];
772     coords[5] = x2[1] * s;
773     R[0] = 1.0;     R[1] = 0.0;     R[2] = 0.0;
774     R[3] = 0.0;     R[4] = 1.0 * s; R[5] = 0.0;
775     R[6] = 0.0;     R[7] = 0.0;     R[8] = 1.0 * s;
776     PetscFunctionReturn(0);
777   }
778   alpha = 1.0/sqrtz;
779   R[0] =  alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz;
780   R[3] = -alpha*n[1];      R[4] = alpha*n[0];      R[5] = 0.0;
781   R[6] =  n[0];            R[7] = n[1];            R[8] = n[2];
782   for (d = 0; d < dim; ++d) {
783     x1p[d] = 0.0;
784     x2p[d] = 0.0;
785     for (e = 0; e < dim; ++e) {
786       x1p[d] += R[d*dim+e]*x1[e];
787       x2p[d] += R[d*dim+e]*x2[e];
788     }
789   }
790   if (PetscAbsReal(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
791   if (PetscAbsReal(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
792   /* 2) Project to (x, y) */
793   for (p = 3; p < coordSize/3; ++p) {
794     for (d = 0; d < dim; ++d) {
795       xnp[d] = 0.0;
796       for (e = 0; e < dim; ++e) {
797         xnp[d] += R[d*dim+e]*PetscRealPart(coords[p*dim+e] - coords[0*dim+e]);
798       }
799       if (d < dim-1) coords[p*2+d] = xnp[d];
800     }
801   }
802   coords[0] = 0.0;
803   coords[1] = 0.0;
804   coords[2] = x1p[0];
805   coords[3] = x1p[1];
806   coords[4] = x2p[0];
807   coords[5] = x2p[1];
808   /* Output R^T which rotates \hat z to the input normal */
809   for (d = 0; d < dim; ++d) {
810     for (e = d+1; e < dim; ++e) {
811       PetscReal tmp;
812 
813       tmp        = R[d*dim+e];
814       R[d*dim+e] = R[e*dim+d];
815       R[e*dim+d] = tmp;
816     }
817   }
818   PetscFunctionReturn(0);
819 }
820 
821 #undef __FUNCT__
822 #define __FUNCT__ "Volume_Triangle_Internal"
823 PETSC_UNUSED
824 PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[])
825 {
826   /* Signed volume is 1/2 the determinant
827 
828    |  1  1  1 |
829    | x0 x1 x2 |
830    | y0 y1 y2 |
831 
832      but if x0,y0 is the origin, we have
833 
834    | x1 x2 |
835    | y1 y2 |
836   */
837   const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1];
838   const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1];
839   PetscReal       M[4], detM;
840   M[0] = x1; M[1] = x2;
841   M[2] = y1; M[3] = y2;
842   DMPlex_Det2D_Internal(&detM, M);
843   *vol = 0.5*detM;
844   (void)PetscLogFlops(5.0);
845 }
846 
847 #undef __FUNCT__
848 #define __FUNCT__ "Volume_Triangle_Origin_Internal"
849 PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[])
850 {
851   DMPlex_Det2D_Internal(vol, coords);
852   *vol *= 0.5;
853 }
854 
855 #undef __FUNCT__
856 #define __FUNCT__ "Volume_Tetrahedron_Internal"
857 PETSC_UNUSED
858 PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[])
859 {
860   /* Signed volume is 1/6th of the determinant
861 
862    |  1  1  1  1 |
863    | x0 x1 x2 x3 |
864    | y0 y1 y2 y3 |
865    | z0 z1 z2 z3 |
866 
867      but if x0,y0,z0 is the origin, we have
868 
869    | x1 x2 x3 |
870    | y1 y2 y3 |
871    | z1 z2 z3 |
872   */
873   const PetscReal x1 = coords[3] - coords[0], y1 = coords[4]  - coords[1], z1 = coords[5]  - coords[2];
874   const PetscReal x2 = coords[6] - coords[0], y2 = coords[7]  - coords[1], z2 = coords[8]  - coords[2];
875   const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2];
876   PetscReal       M[9], detM;
877   M[0] = x1; M[1] = x2; M[2] = x3;
878   M[3] = y1; M[4] = y2; M[5] = y3;
879   M[6] = z1; M[7] = z2; M[8] = z3;
880   DMPlex_Det3D_Internal(&detM, M);
881   *vol = -0.16666666666666666666666*detM;
882   (void)PetscLogFlops(10.0);
883 }
884 
885 #undef __FUNCT__
886 #define __FUNCT__ "Volume_Tetrahedron_Origin_Internal"
887 PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[])
888 {
889   DMPlex_Det3D_Internal(vol, coords);
890   *vol *= -0.16666666666666666666666;
891 }
892 
893 #undef __FUNCT__
894 #define __FUNCT__ "DMPlexComputeLineGeometry_Internal"
895 static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
896 {
897   PetscSection   coordSection;
898   Vec            coordinates;
899   PetscScalar   *coords = NULL;
900   PetscInt       numCoords, d, pStart, pEnd, numSelfCoords = 0;
901   PetscErrorCode ierr;
902 
903   PetscFunctionBegin;
904   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
905   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
906   ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr);
907   if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);}
908   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
909   numCoords = numSelfCoords ? numSelfCoords : numCoords;
910   if (invJ && !J) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL");
911   *detJ = 0.0;
912   if (numCoords == 6) {
913     const PetscInt dim = 3;
914     PetscReal      R[9], J0;
915 
916     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
917     ierr = DMPlexComputeProjection3Dto1D_Internal(coords, R);CHKERRQ(ierr);
918     if (J)    {
919       J0   = 0.5*PetscRealPart(coords[1]);
920       J[0] = R[0]*J0; J[1] = R[1]; J[2] = R[2];
921       J[3] = R[3]*J0; J[4] = R[4]; J[5] = R[5];
922       J[6] = R[6]*J0; J[7] = R[7]; J[8] = R[8];
923       DMPlex_Det3D_Internal(detJ, J);
924       if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
925     }
926   } else if (numCoords == 4) {
927     const PetscInt dim = 2;
928     PetscReal      R[4], J0;
929 
930     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
931     ierr = DMPlexComputeProjection2Dto1D_Internal(coords, R);CHKERRQ(ierr);
932     if (J)    {
933       J0   = 0.5*PetscRealPart(coords[1]);
934       J[0] = R[0]*J0; J[1] = R[1];
935       J[2] = R[2]*J0; J[3] = R[3];
936       DMPlex_Det2D_Internal(detJ, J);
937       if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
938     }
939   } else if (numCoords == 2) {
940     const PetscInt dim = 1;
941 
942     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
943     if (J)    {
944       J[0]  = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0]));
945       *detJ = J[0];
946       ierr = PetscLogFlops(2.0);CHKERRQ(ierr);
947       if (invJ) {invJ[0] = 1.0/J[0]; ierr = PetscLogFlops(1.0);CHKERRQ(ierr);}
948     }
949   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %D != 2", numCoords);
950   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
951   PetscFunctionReturn(0);
952 }
953 
954 #undef __FUNCT__
955 #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal"
956 static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
957 {
958   PetscSection   coordSection;
959   Vec            coordinates;
960   PetscScalar   *coords = NULL;
961   PetscInt       numCoords, d, f, g;
962   PetscErrorCode ierr;
963 
964   PetscFunctionBegin;
965   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
966   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
967   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
968   *detJ = 0.0;
969   if (numCoords == 9) {
970     const PetscInt dim = 3;
971     PetscReal      R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0};
972 
973     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
974     ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr);
975     if (J)    {
976       const PetscInt pdim = 2;
977 
978       for (d = 0; d < pdim; d++) {
979         for (f = 0; f < pdim; f++) {
980           J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
981         }
982       }
983       ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
984       DMPlex_Det3D_Internal(detJ, J0);
985       for (d = 0; d < dim; d++) {
986         for (f = 0; f < dim; f++) {
987           J[d*dim+f] = 0.0;
988           for (g = 0; g < dim; g++) {
989             J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
990           }
991         }
992       }
993       ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
994     }
995     if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
996   } else if (numCoords == 6) {
997     const PetscInt dim = 2;
998 
999     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1000     if (J)    {
1001       for (d = 0; d < dim; d++) {
1002         for (f = 0; f < dim; f++) {
1003           J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
1004         }
1005       }
1006       ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1007       DMPlex_Det2D_Internal(detJ, J);
1008     }
1009     if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1010   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %D != 6 or 9", numCoords);
1011   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
1012   PetscFunctionReturn(0);
1013 }
1014 
1015 #undef __FUNCT__
1016 #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal"
1017 static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1018 {
1019   PetscSection   coordSection;
1020   Vec            coordinates;
1021   PetscScalar   *coords = NULL;
1022   PetscInt       numCoords, numSelfCoords = 0, d, f, g, pStart, pEnd;
1023   PetscErrorCode ierr;
1024 
1025   PetscFunctionBegin;
1026   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1027   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1028   ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr);
1029   if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);}
1030   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
1031   numCoords = numSelfCoords ? numSelfCoords : numCoords;
1032   *detJ = 0.0;
1033   if (numCoords == 12) {
1034     const PetscInt dim = 3;
1035     PetscReal      R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0};
1036 
1037     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1038     ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr);
1039     if (J)    {
1040       const PetscInt pdim = 2;
1041 
1042       for (d = 0; d < pdim; d++) {
1043         J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
1044         J0[d*dim+1] = 0.5*(PetscRealPart(coords[3*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
1045       }
1046       ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1047       DMPlex_Det3D_Internal(detJ, J0);
1048       for (d = 0; d < dim; d++) {
1049         for (f = 0; f < dim; f++) {
1050           J[d*dim+f] = 0.0;
1051           for (g = 0; g < dim; g++) {
1052             J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
1053           }
1054         }
1055       }
1056       ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
1057     }
1058     if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1059   } else if (numCoords == 8) {
1060     const PetscInt dim = 2;
1061 
1062     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1063     if (J)    {
1064       for (d = 0; d < dim; d++) {
1065         J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
1066         J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1067       }
1068       ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1069       DMPlex_Det2D_Internal(detJ, J);
1070     }
1071     if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1072   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords);
1073   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
1074   PetscFunctionReturn(0);
1075 }
1076 
1077 #undef __FUNCT__
1078 #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal"
1079 static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1080 {
1081   PetscSection   coordSection;
1082   Vec            coordinates;
1083   PetscScalar   *coords = NULL;
1084   const PetscInt dim = 3;
1085   PetscInt       d;
1086   PetscErrorCode ierr;
1087 
1088   PetscFunctionBegin;
1089   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1090   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1091   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1092   *detJ = 0.0;
1093   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1094   if (J)    {
1095     for (d = 0; d < dim; d++) {
1096       /* I orient with outward face normals */
1097       J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d]));
1098       J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
1099       J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1100     }
1101     ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
1102     DMPlex_Det3D_Internal(detJ, J);
1103   }
1104   if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1105   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1106   PetscFunctionReturn(0);
1107 }
1108 
1109 #undef __FUNCT__
1110 #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal"
1111 static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1112 {
1113   PetscSection   coordSection;
1114   Vec            coordinates;
1115   PetscScalar   *coords = NULL;
1116   const PetscInt dim = 3;
1117   PetscInt       d;
1118   PetscErrorCode ierr;
1119 
1120   PetscFunctionBegin;
1121   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1122   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1123   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1124   *detJ = 0.0;
1125   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1126   if (J)    {
1127     for (d = 0; d < dim; d++) {
1128       J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1129       J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
1130       J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d]));
1131     }
1132     ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
1133     DMPlex_Det3D_Internal(detJ, J);
1134   }
1135   if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1136   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1137   PetscFunctionReturn(0);
1138 }
1139 
1140 #undef __FUNCT__
1141 #define __FUNCT__ "DMPlexComputeCellGeometryAffineFEM"
1142 /*@C
1143   DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell
1144 
1145   Collective on DM
1146 
1147   Input Arguments:
1148 + dm   - the DM
1149 - cell - the cell
1150 
1151   Output Arguments:
1152 + v0   - the translation part of this affine transform
1153 . J    - the Jacobian of the transform from the reference element
1154 . invJ - the inverse of the Jacobian
1155 - detJ - the Jacobian determinant
1156 
1157   Level: advanced
1158 
1159   Fortran Notes:
1160   Since it returns arrays, this routine is only available in Fortran 90, and you must
1161   include petsc.h90 in your code.
1162 
1163 .seealso: DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinateVec()
1164 @*/
1165 PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
1166 {
1167   PetscInt       depth, dim, coneSize;
1168   PetscErrorCode ierr;
1169 
1170   PetscFunctionBegin;
1171   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1172   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
1173   if (depth == 1) {
1174     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1175   } else {
1176     DMLabel depth;
1177 
1178     ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
1179     ierr = DMLabelGetValue(depth, cell, &dim);CHKERRQ(ierr);
1180   }
1181   switch (dim) {
1182   case 1:
1183     ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
1184     break;
1185   case 2:
1186     switch (coneSize) {
1187     case 3:
1188       ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
1189       break;
1190     case 4:
1191       ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
1192       break;
1193     default:
1194       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell);
1195     }
1196     break;
1197   case 3:
1198     switch (coneSize) {
1199     case 4:
1200       ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
1201       break;
1202     case 6: /* Faces */
1203     case 8: /* Vertices */
1204       ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
1205       break;
1206     default:
1207         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell);
1208     }
1209       break;
1210   default:
1211     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
1212   }
1213   PetscFunctionReturn(0);
1214 }
1215 
1216 #undef __FUNCT__
1217 #define __FUNCT__ "DMPlexComputeIsoparametricGeometry_Internal"
1218 static PetscErrorCode DMPlexComputeIsoparametricGeometry_Internal(DM dm, PetscFE fe, PetscInt point, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1219 {
1220   PetscQuadrature  quad;
1221   PetscSection     coordSection;
1222   Vec              coordinates;
1223   PetscScalar     *coords = NULL;
1224   const PetscReal *quadPoints;
1225   PetscReal       *basisDer;
1226   PetscInt         dim, cdim, pdim, qdim, Nq, numCoords, d, q;
1227   PetscErrorCode   ierr;
1228 
1229   PetscFunctionBegin;
1230   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1231   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1232   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr);
1233   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1234   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
1235   ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
1236   ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr);
1237   ierr = PetscQuadratureGetData(quad, &qdim, &Nq, &quadPoints, NULL);CHKERRQ(ierr);
1238   ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
1239   *detJ = 0.0;
1240   if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %d != quadrature dimension %d", dim, qdim);
1241   if (numCoords != pdim*cdim) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "There are %d coordinates for point %d != %d*%d", numCoords, point, pdim, cdim);
1242   if (v0) {for (d = 0; d < cdim; d++) v0[d] = PetscRealPart(coords[d]);}
1243   if (J) {
1244     ierr = PetscMemzero(J, Nq*cdim*dim*sizeof(PetscReal));CHKERRQ(ierr);
1245     for (q = 0; q < Nq; ++q) {
1246       PetscInt i, j, k, c, r;
1247 
1248       /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */
1249       for (k = 0; k < pdim; ++k)
1250         for (j = 0; j < dim; ++j)
1251           for (i = 0; i < cdim; ++i)
1252             J[(q*cdim + i)*dim + j] += basisDer[(q*pdim + k)*dim + j] * PetscRealPart(coords[k*cdim + i]);
1253       ierr = PetscLogFlops(2.0*pdim*dim*cdim);CHKERRQ(ierr);
1254       if (cdim > dim) {
1255         for (c = dim; c < cdim; ++c)
1256           for (r = 0; r < cdim; ++r)
1257             J[r*cdim+c] = r == c ? 1.0 : 0.0;
1258       }
1259       switch (cdim) {
1260       case 3:
1261         DMPlex_Det3D_Internal(detJ, J);
1262         if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1263         break;
1264       case 2:
1265         DMPlex_Det2D_Internal(detJ, J);
1266         if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1267         break;
1268       case 1:
1269         *detJ = J[0];
1270         if (invJ) invJ[0] = 1.0/J[0];
1271       }
1272     }
1273   }
1274   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr);
1275   PetscFunctionReturn(0);
1276 }
1277 
1278 #undef __FUNCT__
1279 #define __FUNCT__ "DMPlexComputeCellGeometryFEM"
1280 /*@C
1281   DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell
1282 
1283   Collective on DM
1284 
1285   Input Arguments:
1286 + dm   - the DM
1287 . cell - the cell
1288 - fe   - the finite element containing the quadrature
1289 
1290   Output Arguments:
1291 + v0   - the translation part of this transform
1292 . J    - the Jacobian of the transform from the reference element at each quadrature point
1293 . invJ - the inverse of the Jacobian at each quadrature point
1294 - detJ - the Jacobian determinant at each quadrature point
1295 
1296   Level: advanced
1297 
1298   Fortran Notes:
1299   Since it returns arrays, this routine is only available in Fortran 90, and you must
1300   include petsc.h90 in your code.
1301 
1302 .seealso: DMGetCoordinateSection(), DMGetCoordinateVec()
1303 @*/
1304 PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscFE fe, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
1305 {
1306   PetscErrorCode ierr;
1307 
1308   PetscFunctionBegin;
1309   if (!fe) {ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);}
1310   else     {ierr = DMPlexComputeIsoparametricGeometry_Internal(dm, fe, cell, v0, J, invJ, detJ);CHKERRQ(ierr);}
1311   PetscFunctionReturn(0);
1312 }
1313 
1314 #undef __FUNCT__
1315 #define __FUNCT__ "DMPlexComputeGeometryFVM_1D_Internal"
1316 static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1317 {
1318   PetscSection   coordSection;
1319   Vec            coordinates;
1320   PetscScalar   *coords = NULL;
1321   PetscScalar    tmp[2];
1322   PetscInt       coordSize;
1323   PetscErrorCode ierr;
1324 
1325   PetscFunctionBegin;
1326   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1327   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1328   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1329   if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now");
1330   ierr = DMLocalizeCoordinate_Internal(dm, dim, coords, &coords[dim], tmp);CHKERRQ(ierr);
1331   if (centroid) {
1332     centroid[0] = 0.5*PetscRealPart(coords[0] + tmp[0]);
1333     centroid[1] = 0.5*PetscRealPart(coords[1] + tmp[1]);
1334   }
1335   if (normal) {
1336     PetscReal norm;
1337 
1338     normal[0]  = -PetscRealPart(coords[1] - tmp[1]);
1339     normal[1]  =  PetscRealPart(coords[0] - tmp[0]);
1340     norm       = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1]);
1341     normal[0] /= norm;
1342     normal[1] /= norm;
1343   }
1344   if (vol) {
1345     *vol = PetscSqrtReal(PetscSqr(PetscRealPart(coords[0] - tmp[0])) + PetscSqr(PetscRealPart(coords[1] - tmp[1])));
1346   }
1347   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1348   PetscFunctionReturn(0);
1349 }
1350 
1351 #undef __FUNCT__
1352 #define __FUNCT__ "DMPlexComputeGeometryFVM_2D_Internal"
1353 /* Centroid_i = (\sum_n A_n Cn_i ) / A */
1354 static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1355 {
1356   PetscSection   coordSection;
1357   Vec            coordinates;
1358   PetscScalar   *coords = NULL;
1359   PetscReal      vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9];
1360   PetscInt       tdim = 2, coordSize, numCorners, p, d, e;
1361   PetscErrorCode ierr;
1362 
1363   PetscFunctionBegin;
1364   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1365   ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr);
1366   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1367   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1368   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
1369   if (dim > 2 && centroid) {
1370     v0[0] = PetscRealPart(coords[0]);
1371     v0[1] = PetscRealPart(coords[1]);
1372     v0[2] = PetscRealPart(coords[2]);
1373   }
1374   if (normal) {
1375     if (dim > 2) {
1376       const PetscReal x0 = PetscRealPart(coords[dim+0] - coords[0]), x1 = PetscRealPart(coords[dim*2+0] - coords[0]);
1377       const PetscReal y0 = PetscRealPart(coords[dim+1] - coords[1]), y1 = PetscRealPart(coords[dim*2+1] - coords[1]);
1378       const PetscReal z0 = PetscRealPart(coords[dim+2] - coords[2]), z1 = PetscRealPart(coords[dim*2+2] - coords[2]);
1379       PetscReal       norm;
1380 
1381       normal[0] = y0*z1 - z0*y1;
1382       normal[1] = z0*x1 - x0*z1;
1383       normal[2] = x0*y1 - y0*x1;
1384       norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]);
1385       normal[0] /= norm;
1386       normal[1] /= norm;
1387       normal[2] /= norm;
1388     } else {
1389       for (d = 0; d < dim; ++d) normal[d] = 0.0;
1390     }
1391   }
1392   if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D_Internal(coordSize, coords, R);CHKERRQ(ierr);}
1393   for (p = 0; p < numCorners; ++p) {
1394     /* Need to do this copy to get types right */
1395     for (d = 0; d < tdim; ++d) {
1396       ctmp[d]      = PetscRealPart(coords[p*tdim+d]);
1397       ctmp[tdim+d] = PetscRealPart(coords[((p+1)%numCorners)*tdim+d]);
1398     }
1399     Volume_Triangle_Origin_Internal(&vtmp, ctmp);
1400     vsum += vtmp;
1401     for (d = 0; d < tdim; ++d) {
1402       csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp;
1403     }
1404   }
1405   for (d = 0; d < tdim; ++d) {
1406     csum[d] /= (tdim+1)*vsum;
1407   }
1408   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1409   if (vol) *vol = PetscAbsReal(vsum);
1410   if (centroid) {
1411     if (dim > 2) {
1412       for (d = 0; d < dim; ++d) {
1413         centroid[d] = v0[d];
1414         for (e = 0; e < dim; ++e) {
1415           centroid[d] += R[d*dim+e]*csum[e];
1416         }
1417       }
1418     } else for (d = 0; d < dim; ++d) centroid[d] = csum[d];
1419   }
1420   PetscFunctionReturn(0);
1421 }
1422 
1423 #undef __FUNCT__
1424 #define __FUNCT__ "DMPlexComputeGeometryFVM_3D_Internal"
1425 /* Centroid_i = (\sum_n V_n Cn_i ) / V */
1426 static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1427 {
1428   PetscSection    coordSection;
1429   Vec             coordinates;
1430   PetscScalar    *coords = NULL;
1431   PetscReal       vsum = 0.0, vtmp, coordsTmp[3*3];
1432   const PetscInt *faces, *facesO;
1433   PetscInt        numFaces, f, coordSize, numCorners, p, d;
1434   PetscErrorCode  ierr;
1435 
1436   PetscFunctionBegin;
1437   if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No support for dim %D > 3",dim);
1438   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1439   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1440 
1441   if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0;
1442   ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr);
1443   ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr);
1444   ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr);
1445   for (f = 0; f < numFaces; ++f) {
1446     ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr);
1447     numCorners = coordSize/dim;
1448     switch (numCorners) {
1449     case 3:
1450       for (d = 0; d < dim; ++d) {
1451         coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]);
1452         coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]);
1453         coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]);
1454       }
1455       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
1456       if (facesO[f] < 0) vtmp = -vtmp;
1457       vsum += vtmp;
1458       if (centroid) {           /* Centroid of OABC = (a+b+c)/4 */
1459         for (d = 0; d < dim; ++d) {
1460           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
1461         }
1462       }
1463       break;
1464     case 4:
1465       /* DO FOR PYRAMID */
1466       /* First tet */
1467       for (d = 0; d < dim; ++d) {
1468         coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]);
1469         coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]);
1470         coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]);
1471       }
1472       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
1473       if (facesO[f] < 0) vtmp = -vtmp;
1474       vsum += vtmp;
1475       if (centroid) {
1476         for (d = 0; d < dim; ++d) {
1477           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
1478         }
1479       }
1480       /* Second tet */
1481       for (d = 0; d < dim; ++d) {
1482         coordsTmp[0*dim+d] = PetscRealPart(coords[1*dim+d]);
1483         coordsTmp[1*dim+d] = PetscRealPart(coords[2*dim+d]);
1484         coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]);
1485       }
1486       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
1487       if (facesO[f] < 0) vtmp = -vtmp;
1488       vsum += vtmp;
1489       if (centroid) {
1490         for (d = 0; d < dim; ++d) {
1491           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
1492         }
1493       }
1494       break;
1495     default:
1496       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %D vertices", numCorners);
1497     }
1498     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr);
1499   }
1500   if (vol)     *vol = PetscAbsReal(vsum);
1501   if (normal)   for (d = 0; d < dim; ++d) normal[d]    = 0.0;
1502   if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4);
1503   PetscFunctionReturn(0);
1504 }
1505 
1506 #undef __FUNCT__
1507 #define __FUNCT__ "DMPlexComputeCellGeometryFVM"
1508 /*@C
1509   DMPlexComputeCellGeometryFVM - Compute the volume for a given cell
1510 
1511   Collective on DM
1512 
1513   Input Arguments:
1514 + dm   - the DM
1515 - cell - the cell
1516 
1517   Output Arguments:
1518 + volume   - the cell volume
1519 . centroid - the cell centroid
1520 - normal - the cell normal, if appropriate
1521 
1522   Level: advanced
1523 
1524   Fortran Notes:
1525   Since it returns arrays, this routine is only available in Fortran 90, and you must
1526   include petsc.h90 in your code.
1527 
1528 .seealso: DMGetCoordinateSection(), DMGetCoordinateVec()
1529 @*/
1530 PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1531 {
1532   PetscInt       depth, dim;
1533   PetscErrorCode ierr;
1534 
1535   PetscFunctionBegin;
1536   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1537   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1538   if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated");
1539   /* We need to keep a pointer to the depth label */
1540   ierr = DMGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr);
1541   /* Cone size is now the number of faces */
1542   switch (depth) {
1543   case 1:
1544     ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
1545     break;
1546   case 2:
1547     ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
1548     break;
1549   case 3:
1550     ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
1551     break;
1552   default:
1553     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
1554   }
1555   PetscFunctionReturn(0);
1556 }
1557 
1558 #undef __FUNCT__
1559 #define __FUNCT__ "DMPlexComputeGeometryFEM"
1560 /* This should also take a PetscFE argument I think */
1561 PetscErrorCode DMPlexComputeGeometryFEM(DM dm, Vec *cellgeom)
1562 {
1563   DM             dmCell;
1564   Vec            coordinates;
1565   PetscSection   coordSection, sectionCell;
1566   PetscScalar   *cgeom;
1567   PetscInt       cStart, cEnd, cMax, c;
1568   PetscErrorCode ierr;
1569 
1570   PetscFunctionBegin;
1571   ierr = DMClone(dm, &dmCell);CHKERRQ(ierr);
1572   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1573   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1574   ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr);
1575   ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr);
1576   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionCell);CHKERRQ(ierr);
1577   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1578   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
1579   cEnd = cMax < 0 ? cEnd : cMax;
1580   ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr);
1581   /* TODO This needs to be multiplied by Nq for non-affine */
1582   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFECellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);}
1583   ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr);
1584   ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr);
1585   ierr = PetscSectionDestroy(&sectionCell);CHKERRQ(ierr);
1586   ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr);
1587   ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr);
1588   for (c = cStart; c < cEnd; ++c) {
1589     PetscFECellGeom *cg;
1590 
1591     ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
1592     ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr);
1593     ierr = DMPlexComputeCellGeometryFEM(dmCell, c, NULL, cg->v0, cg->J, cg->invJ, &cg->detJ);CHKERRQ(ierr);
1594     if (cg->detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cg->detJ, c);
1595   }
1596   ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr);
1597   ierr = DMDestroy(&dmCell);CHKERRQ(ierr);
1598   PetscFunctionReturn(0);
1599 }
1600 
1601 #undef __FUNCT__
1602 #define __FUNCT__ "DMPlexComputeGeometryFVM"
1603 /*@
1604   DMPlexComputeGeometryFVM - Computes the cell and face geometry for a finite volume method
1605 
1606   Input Parameter:
1607 . dm - The DM
1608 
1609   Output Parameters:
1610 + cellgeom - A Vec of PetscFVCellGeom data
1611 . facegeom - A Vec of PetscFVFaceGeom data
1612 
1613   Level: developer
1614 
1615 .seealso: PetscFVFaceGeom, PetscFVCellGeom, DMPlexComputeGeometryFEM()
1616 @*/
1617 PetscErrorCode DMPlexComputeGeometryFVM(DM dm, Vec *cellgeom, Vec *facegeom)
1618 {
1619   DM             dmFace, dmCell;
1620   DMLabel        ghostLabel;
1621   PetscSection   sectionFace, sectionCell;
1622   PetscSection   coordSection;
1623   Vec            coordinates;
1624   PetscScalar   *fgeom, *cgeom;
1625   PetscReal      minradius, gminradius;
1626   PetscInt       dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f;
1627   PetscErrorCode ierr;
1628 
1629   PetscFunctionBegin;
1630   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1631   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1632   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1633   /* Make cell centroids and volumes */
1634   ierr = DMClone(dm, &dmCell);CHKERRQ(ierr);
1635   ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr);
1636   ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr);
1637   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionCell);CHKERRQ(ierr);
1638   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1639   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1640   ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr);
1641   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVCellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);}
1642   ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr);
1643   ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr);
1644   ierr = PetscSectionDestroy(&sectionCell);CHKERRQ(ierr);
1645   ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr);
1646   if (cEndInterior < 0) {
1647     cEndInterior = cEnd;
1648   }
1649   ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr);
1650   for (c = cStart; c < cEndInterior; ++c) {
1651     PetscFVCellGeom *cg;
1652 
1653     ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
1654     ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr);
1655     ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr);
1656   }
1657   /* Compute face normals and minimum cell radius */
1658   ierr = DMClone(dm, &dmFace);CHKERRQ(ierr);
1659   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionFace);CHKERRQ(ierr);
1660   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
1661   ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr);
1662   for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVFaceGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);}
1663   ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr);
1664   ierr = DMSetDefaultSection(dmFace, sectionFace);CHKERRQ(ierr);
1665   ierr = PetscSectionDestroy(&sectionFace);CHKERRQ(ierr);
1666   ierr = DMCreateLocalVector(dmFace, facegeom);CHKERRQ(ierr);
1667   ierr = VecGetArray(*facegeom, &fgeom);CHKERRQ(ierr);
1668   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
1669   minradius = PETSC_MAX_REAL;
1670   for (f = fStart; f < fEnd; ++f) {
1671     PetscFVFaceGeom *fg;
1672     PetscReal        area;
1673     PetscInt         ghost = -1, d, numChildren;
1674 
1675     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);}
1676     ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr);
1677     if (ghost >= 0 || numChildren) continue;
1678     ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr);
1679     ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr);
1680     for (d = 0; d < dim; ++d) fg->normal[d] *= area;
1681     /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */
1682     {
1683       PetscFVCellGeom *cL, *cR;
1684       PetscInt         ncells;
1685       const PetscInt  *cells;
1686       PetscReal       *lcentroid, *rcentroid;
1687       PetscReal        l[3], r[3], v[3];
1688 
1689       ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr);
1690       ierr = DMPlexGetSupportSize(dm, f, &ncells);CHKERRQ(ierr);
1691       ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr);
1692       lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid;
1693       if (ncells > 1) {
1694         ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr);
1695         rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid;
1696       }
1697       else {
1698         rcentroid = fg->centroid;
1699       }
1700       ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, lcentroid, l);CHKERRQ(ierr);
1701       ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, rcentroid, r);CHKERRQ(ierr);
1702       DMPlex_WaxpyD_Internal(dim, -1, l, r, v);
1703       if (DMPlex_DotRealD_Internal(dim, fg->normal, v) < 0) {
1704         for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d];
1705       }
1706       if (DMPlex_DotRealD_Internal(dim, fg->normal, v) <= 0) {
1707         if (dim == 2) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed, normal (%g,%g) v (%g,%g)", f, (double) fg->normal[0], (double) fg->normal[1], (double) v[0], (double) v[1]);
1708         if (dim == 3) SETERRQ7(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed, normal (%g,%g,%g) v (%g,%g,%g)", f, (double) fg->normal[0], (double) fg->normal[1], (double) fg->normal[2], (double) v[0], (double) v[1], (double) v[2]);
1709         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f);
1710       }
1711       if (cells[0] < cEndInterior) {
1712         DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cL->centroid, v);
1713         minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v));
1714       }
1715       if (ncells > 1 && cells[1] < cEndInterior) {
1716         DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cR->centroid, v);
1717         minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v));
1718       }
1719     }
1720   }
1721   ierr = MPIU_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
1722   ierr = DMPlexSetMinRadius(dm, gminradius);CHKERRQ(ierr);
1723   /* Compute centroids of ghost cells */
1724   for (c = cEndInterior; c < cEnd; ++c) {
1725     PetscFVFaceGeom *fg;
1726     const PetscInt  *cone,    *support;
1727     PetscInt         coneSize, supportSize, s;
1728 
1729     ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr);
1730     if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize);
1731     ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr);
1732     ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr);
1733     if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 2", cone[0], supportSize);
1734     ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr);
1735     ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr);
1736     for (s = 0; s < 2; ++s) {
1737       /* Reflect ghost centroid across plane of face */
1738       if (support[s] == c) {
1739         PetscFVCellGeom       *ci;
1740         PetscFVCellGeom       *cg;
1741         PetscReal              c2f[3], a;
1742 
1743         ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr);
1744         DMPlex_WaxpyD_Internal(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */
1745         a    = DMPlex_DotRealD_Internal(dim, c2f, fg->normal)/DMPlex_DotRealD_Internal(dim, fg->normal, fg->normal);
1746         ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr);
1747         DMPlex_WaxpyD_Internal(dim, 2*a, fg->normal, ci->centroid, cg->centroid);
1748         cg->volume = ci->volume;
1749       }
1750     }
1751   }
1752   ierr = VecRestoreArray(*facegeom, &fgeom);CHKERRQ(ierr);
1753   ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr);
1754   ierr = DMDestroy(&dmCell);CHKERRQ(ierr);
1755   ierr = DMDestroy(&dmFace);CHKERRQ(ierr);
1756   PetscFunctionReturn(0);
1757 }
1758 
1759 #undef __FUNCT__
1760 #define __FUNCT__ "DMPlexGetMinRadius"
1761 /*@C
1762   DMPlexGetMinRadius - Returns the minimum distance from any cell centroid to a face
1763 
1764   Not collective
1765 
1766   Input Argument:
1767 . dm - the DM
1768 
1769   Output Argument:
1770 . minradius - the minium cell radius
1771 
1772   Level: developer
1773 
1774 .seealso: DMGetCoordinates()
1775 @*/
1776 PetscErrorCode DMPlexGetMinRadius(DM dm, PetscReal *minradius)
1777 {
1778   PetscFunctionBegin;
1779   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1780   PetscValidPointer(minradius,2);
1781   *minradius = ((DM_Plex*) dm->data)->minradius;
1782   PetscFunctionReturn(0);
1783 }
1784 
1785 #undef __FUNCT__
1786 #define __FUNCT__ "DMPlexSetMinRadius"
1787 /*@C
1788   DMPlexSetMinRadius - Sets the minimum distance from the cell centroid to a face
1789 
1790   Logically collective
1791 
1792   Input Arguments:
1793 + dm - the DM
1794 - minradius - the minium cell radius
1795 
1796   Level: developer
1797 
1798 .seealso: DMSetCoordinates()
1799 @*/
1800 PetscErrorCode DMPlexSetMinRadius(DM dm, PetscReal minradius)
1801 {
1802   PetscFunctionBegin;
1803   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1804   ((DM_Plex*) dm->data)->minradius = minradius;
1805   PetscFunctionReturn(0);
1806 }
1807 
1808 #undef __FUNCT__
1809 #define __FUNCT__ "BuildGradientReconstruction_Internal"
1810 static PetscErrorCode BuildGradientReconstruction_Internal(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
1811 {
1812   DMLabel        ghostLabel;
1813   PetscScalar   *dx, *grad, **gref;
1814   PetscInt       dim, cStart, cEnd, c, cEndInterior, maxNumFaces;
1815   PetscErrorCode ierr;
1816 
1817   PetscFunctionBegin;
1818   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1819   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1820   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1821   ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr);
1822   ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr);
1823   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
1824   ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr);
1825   for (c = cStart; c < cEndInterior; c++) {
1826     const PetscInt        *faces;
1827     PetscInt               numFaces, usedFaces, f, d;
1828     PetscFVCellGeom        *cg;
1829     PetscBool              boundary;
1830     PetscInt               ghost;
1831 
1832     ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
1833     ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr);
1834     ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr);
1835     if (numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces);
1836     for (f = 0, usedFaces = 0; f < numFaces; ++f) {
1837       PetscFVCellGeom       *cg1;
1838       PetscFVFaceGeom       *fg;
1839       const PetscInt        *fcells;
1840       PetscInt               ncell, side;
1841 
1842       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
1843       ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
1844       if ((ghost >= 0) || boundary) continue;
1845       ierr  = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr);
1846       side  = (c != fcells[0]); /* c is on left=0 or right=1 of face */
1847       ncell = fcells[!side];    /* the neighbor */
1848       ierr  = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr);
1849       ierr  = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr);
1850       for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d];
1851       gref[usedFaces++] = fg->grad[side];  /* Gradient reconstruction term will go here */
1852     }
1853     if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?");
1854     ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr);
1855     for (f = 0, usedFaces = 0; f < numFaces; ++f) {
1856       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
1857       ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
1858       if ((ghost >= 0) || boundary) continue;
1859       for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d];
1860       ++usedFaces;
1861     }
1862   }
1863   ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr);
1864   PetscFunctionReturn(0);
1865 }
1866 
1867 #undef __FUNCT__
1868 #define __FUNCT__ "BuildGradientReconstruction_Internal_Tree"
1869 static PetscErrorCode BuildGradientReconstruction_Internal_Tree(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
1870 {
1871   DMLabel        ghostLabel;
1872   PetscScalar   *dx, *grad, **gref;
1873   PetscInt       dim, cStart, cEnd, c, cEndInterior, fStart, fEnd, f, nStart, nEnd, maxNumFaces = 0;
1874   PetscSection   neighSec;
1875   PetscInt     (*neighbors)[2];
1876   PetscInt      *counter;
1877   PetscErrorCode ierr;
1878 
1879   PetscFunctionBegin;
1880   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1881   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1882   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1883   if (cEndInterior < 0) {
1884     cEndInterior = cEnd;
1885   }
1886   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&neighSec);CHKERRQ(ierr);
1887   ierr = PetscSectionSetChart(neighSec,cStart,cEndInterior);CHKERRQ(ierr);
1888   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
1889   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
1890   for (f = fStart; f < fEnd; f++) {
1891     const PetscInt        *fcells;
1892     PetscBool              boundary;
1893     PetscInt               ghost = -1;
1894     PetscInt               numChildren, numCells, c;
1895 
1896     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);}
1897     ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr);
1898     ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr);
1899     if ((ghost >= 0) || boundary || numChildren) continue;
1900     ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr);
1901     if (numCells == 2) {
1902       ierr = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr);
1903       for (c = 0; c < 2; c++) {
1904         PetscInt cell = fcells[c];
1905 
1906         if (cell >= cStart && cell < cEndInterior) {
1907           ierr = PetscSectionAddDof(neighSec,cell,1);CHKERRQ(ierr);
1908         }
1909       }
1910     }
1911   }
1912   ierr = PetscSectionSetUp(neighSec);CHKERRQ(ierr);
1913   ierr = PetscSectionGetMaxDof(neighSec,&maxNumFaces);CHKERRQ(ierr);
1914   ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr);
1915   nStart = 0;
1916   ierr = PetscSectionGetStorageSize(neighSec,&nEnd);CHKERRQ(ierr);
1917   ierr = PetscMalloc1((nEnd-nStart),&neighbors);CHKERRQ(ierr);
1918   ierr = PetscCalloc1((cEndInterior-cStart),&counter);CHKERRQ(ierr);
1919   for (f = fStart; f < fEnd; f++) {
1920     const PetscInt        *fcells;
1921     PetscBool              boundary;
1922     PetscInt               ghost = -1;
1923     PetscInt               numChildren, numCells, c;
1924 
1925     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);}
1926     ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr);
1927     ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr);
1928     if ((ghost >= 0) || boundary || numChildren) continue;
1929     ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr);
1930     if (numCells == 2) {
1931       ierr  = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr);
1932       for (c = 0; c < 2; c++) {
1933         PetscInt cell = fcells[c], off;
1934 
1935         if (cell >= cStart && cell < cEndInterior) {
1936           ierr = PetscSectionGetOffset(neighSec,cell,&off);CHKERRQ(ierr);
1937           off += counter[cell - cStart]++;
1938           neighbors[off][0] = f;
1939           neighbors[off][1] = fcells[1 - c];
1940         }
1941       }
1942     }
1943   }
1944   ierr = PetscFree(counter);CHKERRQ(ierr);
1945   ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr);
1946   for (c = cStart; c < cEndInterior; c++) {
1947     PetscInt               numFaces, f, d, off, ghost = -1;
1948     PetscFVCellGeom        *cg;
1949 
1950     ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
1951     ierr = PetscSectionGetDof(neighSec, c, &numFaces);CHKERRQ(ierr);
1952     ierr = PetscSectionGetOffset(neighSec, c, &off);CHKERRQ(ierr);
1953     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, c, &ghost);CHKERRQ(ierr);}
1954     if (ghost < 0 && numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces);
1955     for (f = 0; f < numFaces; ++f) {
1956       PetscFVCellGeom       *cg1;
1957       PetscFVFaceGeom       *fg;
1958       const PetscInt        *fcells;
1959       PetscInt               ncell, side, nface;
1960 
1961       nface = neighbors[off + f][0];
1962       ncell = neighbors[off + f][1];
1963       ierr  = DMPlexGetSupport(dm,nface,&fcells);CHKERRQ(ierr);
1964       side  = (c != fcells[0]);
1965       ierr  = DMPlexPointLocalRef(dmFace, nface, fgeom, &fg);CHKERRQ(ierr);
1966       ierr  = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr);
1967       for (d = 0; d < dim; ++d) dx[f*dim+d] = cg1->centroid[d] - cg->centroid[d];
1968       gref[f] = fg->grad[side];  /* Gradient reconstruction term will go here */
1969     }
1970     ierr = PetscFVComputeGradient(fvm, numFaces, dx, grad);CHKERRQ(ierr);
1971     for (f = 0; f < numFaces; ++f) {
1972       for (d = 0; d < dim; ++d) gref[f][d] = grad[f*dim+d];
1973     }
1974   }
1975   ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr);
1976   ierr = PetscSectionDestroy(&neighSec);CHKERRQ(ierr);
1977   ierr = PetscFree(neighbors);CHKERRQ(ierr);
1978   PetscFunctionReturn(0);
1979 }
1980 
1981 #undef __FUNCT__
1982 #define __FUNCT__ "DMPlexComputeGradientFVM"
1983 /*@
1984   DMPlexComputeGradientFVM - Compute geometric factors for gradient reconstruction, which are stored in the geometry data, and compute layout for gradient data
1985 
1986   Collective on DM
1987 
1988   Input Arguments:
1989 + dm  - The DM
1990 . fvm - The PetscFV
1991 . faceGeometry - The face geometry from DMPlexComputeFaceGeometryFVM()
1992 - cellGeometry - The face geometry from DMPlexComputeCellGeometryFVM()
1993 
1994   Output Parameters:
1995 + faceGeometry - The geometric factors for gradient calculation are inserted
1996 - dmGrad - The DM describing the layout of gradient data
1997 
1998   Level: developer
1999 
2000 .seealso: DMPlexGetFaceGeometryFVM(), DMPlexGetCellGeometryFVM()
2001 @*/
2002 PetscErrorCode DMPlexComputeGradientFVM(DM dm, PetscFV fvm, Vec faceGeometry, Vec cellGeometry, DM *dmGrad)
2003 {
2004   DM             dmFace, dmCell;
2005   PetscScalar   *fgeom, *cgeom;
2006   PetscSection   sectionGrad, parentSection;
2007   PetscInt       dim, pdim, cStart, cEnd, cEndInterior, c;
2008   PetscErrorCode ierr;
2009 
2010   PetscFunctionBegin;
2011   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2012   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
2013   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2014   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
2015   /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */
2016   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
2017   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
2018   ierr = VecGetArray(faceGeometry, &fgeom);CHKERRQ(ierr);
2019   ierr = VecGetArray(cellGeometry, &cgeom);CHKERRQ(ierr);
2020   ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
2021   if (!parentSection) {
2022     ierr = BuildGradientReconstruction_Internal(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr);
2023   } else {
2024     ierr = BuildGradientReconstruction_Internal_Tree(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr);
2025   }
2026   ierr = VecRestoreArray(faceGeometry, &fgeom);CHKERRQ(ierr);
2027   ierr = VecRestoreArray(cellGeometry, &cgeom);CHKERRQ(ierr);
2028   /* Create storage for gradients */
2029   ierr = DMClone(dm, dmGrad);CHKERRQ(ierr);
2030   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionGrad);CHKERRQ(ierr);
2031   ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr);
2032   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);}
2033   ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr);
2034   ierr = DMSetDefaultSection(*dmGrad, sectionGrad);CHKERRQ(ierr);
2035   ierr = PetscSectionDestroy(&sectionGrad);CHKERRQ(ierr);
2036   PetscFunctionReturn(0);
2037 }
2038