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