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