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