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