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