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