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