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