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