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