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