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