xref: /petsc/src/dm/impls/plex/plexgeometry.c (revision bfcb38ea38335faa6e7f8d97f6bc6ce9aa2a1dd1)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 #include <petsc/private/petscfeimpl.h>  /*I      "petscfe.h"       I*/
3 #include <petscblaslapack.h>
4 #include <petsctime.h>
5 
6 static PetscErrorCode DMPlexGetLineIntersection_2D_Internal(const PetscReal segmentA[], const PetscReal segmentB[], PetscReal intersection[], PetscBool *hasIntersection)
7 {
8   const PetscReal p0_x  = segmentA[0*2+0];
9   const PetscReal p0_y  = segmentA[0*2+1];
10   const PetscReal p1_x  = segmentA[1*2+0];
11   const PetscReal p1_y  = segmentA[1*2+1];
12   const PetscReal p2_x  = segmentB[0*2+0];
13   const PetscReal p2_y  = segmentB[0*2+1];
14   const PetscReal p3_x  = segmentB[1*2+0];
15   const PetscReal p3_y  = segmentB[1*2+1];
16   const PetscReal s1_x  = p1_x - p0_x;
17   const PetscReal s1_y  = p1_y - p0_y;
18   const PetscReal s2_x  = p3_x - p2_x;
19   const PetscReal s2_y  = p3_y - p2_y;
20   const PetscReal denom = (-s2_x * s1_y + s1_x * s2_y);
21 
22   PetscFunctionBegin;
23   *hasIntersection = PETSC_FALSE;
24   /* Non-parallel lines */
25   if (denom != 0.0) {
26     const PetscReal s = (-s1_y * (p0_x - p2_x) + s1_x * (p0_y - p2_y)) / denom;
27     const PetscReal t = ( s2_x * (p0_y - p2_y) - s2_y * (p0_x - p2_x)) / denom;
28 
29     if (s >= 0 && s <= 1 && t >= 0 && t <= 1) {
30       *hasIntersection = PETSC_TRUE;
31       if (intersection) {
32         intersection[0] = p0_x + (t * s1_x);
33         intersection[1] = p0_y + (t * s1_y);
34       }
35     }
36   }
37   PetscFunctionReturn(0);
38 }
39 
40 static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
41 {
42   const PetscInt  embedDim = 2;
43   const PetscReal eps      = PETSC_SQRT_MACHINE_EPSILON;
44   PetscReal       x        = PetscRealPart(point[0]);
45   PetscReal       y        = PetscRealPart(point[1]);
46   PetscReal       v0[2], J[4], invJ[4], detJ;
47   PetscReal       xi, eta;
48   PetscErrorCode  ierr;
49 
50   PetscFunctionBegin;
51   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
52   xi  = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]);
53   eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]);
54 
55   if ((xi >= -eps) && (eta >= -eps) && (xi + eta <= 2.0+eps)) *cell = c;
56   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
57   PetscFunctionReturn(0);
58 }
59 
60 static PetscErrorCode DMPlexClosestPoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscReal cpoint[])
61 {
62   const PetscInt  embedDim = 2;
63   PetscReal       x        = PetscRealPart(point[0]);
64   PetscReal       y        = PetscRealPart(point[1]);
65   PetscReal       v0[2], J[4], invJ[4], detJ;
66   PetscReal       xi, eta, r;
67   PetscErrorCode  ierr;
68 
69   PetscFunctionBegin;
70   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
71   xi  = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]);
72   eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]);
73 
74   xi  = PetscMax(xi,  0.0);
75   eta = PetscMax(eta, 0.0);
76   if (xi + eta > 2.0) {
77     r    = (xi + eta)/2.0;
78     xi  /= r;
79     eta /= r;
80   }
81   cpoint[0] = J[0*embedDim+0]*xi + J[0*embedDim+1]*eta + v0[0];
82   cpoint[1] = J[1*embedDim+0]*xi + J[1*embedDim+1]*eta + v0[1];
83   PetscFunctionReturn(0);
84 }
85 
86 static PetscErrorCode DMPlexLocatePoint_General_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
87 {
88   PetscSection       coordSection;
89   Vec             coordsLocal;
90   PetscScalar    *coords = NULL;
91   const PetscInt  faces[8]  = {0, 1, 1, 2, 2, 3, 3, 0};
92   PetscReal       x         = PetscRealPart(point[0]);
93   PetscReal       y         = PetscRealPart(point[1]);
94   PetscInt        crossings = 0, f;
95   PetscErrorCode  ierr;
96 
97   PetscFunctionBegin;
98   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
99   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
100   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
101   for (f = 0; f < 4; ++f) {
102     PetscReal x_i   = PetscRealPart(coords[faces[2*f+0]*2+0]);
103     PetscReal y_i   = PetscRealPart(coords[faces[2*f+0]*2+1]);
104     PetscReal x_j   = PetscRealPart(coords[faces[2*f+1]*2+0]);
105     PetscReal y_j   = PetscRealPart(coords[faces[2*f+1]*2+1]);
106     PetscReal slope = (y_j - y_i) / (x_j - x_i);
107     PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE;
108     PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE;
109     PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE;
110     if ((cond1 || cond2)  && above) ++crossings;
111   }
112   if (crossings % 2) *cell = c;
113   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
114   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
115   PetscFunctionReturn(0);
116 }
117 
118 static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
119 {
120   const PetscInt embedDim = 3;
121   PetscReal      v0[3], J[9], invJ[9], detJ;
122   PetscReal      x = PetscRealPart(point[0]);
123   PetscReal      y = PetscRealPart(point[1]);
124   PetscReal      z = PetscRealPart(point[2]);
125   PetscReal      xi, eta, zeta;
126   PetscErrorCode ierr;
127 
128   PetscFunctionBegin;
129   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
130   xi   = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]);
131   eta  = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]);
132   zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]);
133 
134   if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) *cell = c;
135   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
136   PetscFunctionReturn(0);
137 }
138 
139 static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
140 {
141   PetscSection   coordSection;
142   Vec            coordsLocal;
143   PetscScalar   *coords = NULL;
144   const PetscInt faces[24] = {0, 3, 2, 1,  5, 4, 7, 6,  3, 0, 4, 5,
145                               1, 2, 6, 7,  3, 5, 6, 2,  0, 1, 7, 4};
146   PetscBool      found = PETSC_TRUE;
147   PetscInt       f;
148   PetscErrorCode ierr;
149 
150   PetscFunctionBegin;
151   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
152   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
153   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
154   for (f = 0; f < 6; ++f) {
155     /* Check the point is under plane */
156     /*   Get face normal */
157     PetscReal v_i[3];
158     PetscReal v_j[3];
159     PetscReal normal[3];
160     PetscReal pp[3];
161     PetscReal dot;
162 
163     v_i[0]    = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]);
164     v_i[1]    = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]);
165     v_i[2]    = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]);
166     v_j[0]    = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]);
167     v_j[1]    = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]);
168     v_j[2]    = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]);
169     normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1];
170     normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2];
171     normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0];
172     pp[0]     = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]);
173     pp[1]     = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]);
174     pp[2]     = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]);
175     dot       = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2];
176 
177     /* Check that projected point is in face (2D location problem) */
178     if (dot < 0.0) {
179       found = PETSC_FALSE;
180       break;
181     }
182   }
183   if (found) *cell = c;
184   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
185   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
186   PetscFunctionReturn(0);
187 }
188 
189 static PetscErrorCode PetscGridHashInitialize_Internal(PetscGridHash box, PetscInt dim, const PetscScalar point[])
190 {
191   PetscInt d;
192 
193   PetscFunctionBegin;
194   box->dim = dim;
195   for (d = 0; d < dim; ++d) box->lower[d] = box->upper[d] = PetscRealPart(point[d]);
196   PetscFunctionReturn(0);
197 }
198 
199 PetscErrorCode PetscGridHashCreate(MPI_Comm comm, PetscInt dim, const PetscScalar point[], PetscGridHash *box)
200 {
201   PetscErrorCode ierr;
202 
203   PetscFunctionBegin;
204   ierr = PetscMalloc1(1, box);CHKERRQ(ierr);
205   ierr = PetscGridHashInitialize_Internal(*box, dim, point);CHKERRQ(ierr);
206   PetscFunctionReturn(0);
207 }
208 
209 PetscErrorCode PetscGridHashEnlarge(PetscGridHash box, const PetscScalar point[])
210 {
211   PetscInt d;
212 
213   PetscFunctionBegin;
214   for (d = 0; d < box->dim; ++d) {
215     box->lower[d] = PetscMin(box->lower[d], PetscRealPart(point[d]));
216     box->upper[d] = PetscMax(box->upper[d], PetscRealPart(point[d]));
217   }
218   PetscFunctionReturn(0);
219 }
220 
221 /*
222   PetscGridHashSetGrid - Divide the grid into boxes
223 
224   Not collective
225 
226   Input Parameters:
227 + box - The grid hash object
228 . n   - The number of boxes in each dimension, or PETSC_DETERMINE
229 - h   - The box size in each dimension, only used if n[d] == PETSC_DETERMINE
230 
231   Level: developer
232 
233 .seealso: PetscGridHashCreate()
234 */
235 PetscErrorCode PetscGridHashSetGrid(PetscGridHash box, const PetscInt n[], const PetscReal h[])
236 {
237   PetscInt d;
238 
239   PetscFunctionBegin;
240   for (d = 0; d < box->dim; ++d) {
241     box->extent[d] = box->upper[d] - box->lower[d];
242     if (n[d] == PETSC_DETERMINE) {
243       box->h[d] = h[d];
244       box->n[d] = PetscCeilReal(box->extent[d]/h[d]);
245     } else {
246       box->n[d] = n[d];
247       box->h[d] = box->extent[d]/n[d];
248     }
249   }
250   PetscFunctionReturn(0);
251 }
252 
253 /*
254   PetscGridHashGetEnclosingBox - Find the grid boxes containing each input point
255 
256   Not collective
257 
258   Input Parameters:
259 + box       - The grid hash object
260 . numPoints - The number of input points
261 - points    - The input point coordinates
262 
263   Output Parameters:
264 + dboxes    - An array of numPoints*dim integers expressing the enclosing box as (i_0, i_1, ..., i_dim)
265 - boxes     - An array of numPoints integers expressing the enclosing box as single number, or NULL
266 
267   Level: developer
268 
269 .seealso: PetscGridHashCreate()
270 */
271 PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[])
272 {
273   const PetscReal *lower = box->lower;
274   const PetscReal *upper = box->upper;
275   const PetscReal *h     = box->h;
276   const PetscInt  *n     = box->n;
277   const PetscInt   dim   = box->dim;
278   PetscInt         d, p;
279 
280   PetscFunctionBegin;
281   for (p = 0; p < numPoints; ++p) {
282     for (d = 0; d < dim; ++d) {
283       PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]);
284 
285       if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1;
286       if (dbox == -1   && PetscAbsReal(PetscRealPart(points[p*dim+d]) - lower[d]) < 1.0e-9) dbox = 0;
287       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",
288                                              p, PetscRealPart(points[p*dim+0]), dim > 1 ? PetscRealPart(points[p*dim+1]) : 0.0, dim > 2 ? PetscRealPart(points[p*dim+2]) : 0.0);
289       dboxes[p*dim+d] = dbox;
290     }
291     if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1];
292   }
293   PetscFunctionReturn(0);
294 }
295 
296 /*
297  PetscGridHashGetEnclosingBoxQuery - Find the grid boxes containing each input point
298 
299  Not collective
300 
301   Input Parameters:
302 + box       - The grid hash object
303 . numPoints - The number of input points
304 - points    - The input point coordinates
305 
306   Output Parameters:
307 + dboxes    - An array of numPoints*dim integers expressing the enclosing box as (i_0, i_1, ..., i_dim)
308 . boxes     - An array of numPoints integers expressing the enclosing box as single number, or NULL
309 - found     - Flag indicating if point was located within a box
310 
311   Level: developer
312 
313 .seealso: PetscGridHashGetEnclosingBox()
314 */
315 PetscErrorCode PetscGridHashGetEnclosingBoxQuery(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[],PetscBool *found)
316 {
317   const PetscReal *lower = box->lower;
318   const PetscReal *upper = box->upper;
319   const PetscReal *h     = box->h;
320   const PetscInt  *n     = box->n;
321   const PetscInt   dim   = box->dim;
322   PetscInt         d, p;
323 
324   PetscFunctionBegin;
325   *found = PETSC_FALSE;
326   for (p = 0; p < numPoints; ++p) {
327     for (d = 0; d < dim; ++d) {
328       PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]);
329 
330       if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1;
331       if (dbox < 0 || dbox >= n[d]) {
332         PetscFunctionReturn(0);
333       }
334       dboxes[p*dim+d] = dbox;
335     }
336     if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1];
337   }
338   *found = PETSC_TRUE;
339   PetscFunctionReturn(0);
340 }
341 
342 PetscErrorCode PetscGridHashDestroy(PetscGridHash *box)
343 {
344   PetscErrorCode ierr;
345 
346   PetscFunctionBegin;
347   if (*box) {
348     ierr = PetscSectionDestroy(&(*box)->cellSection);CHKERRQ(ierr);
349     ierr = ISDestroy(&(*box)->cells);CHKERRQ(ierr);
350     ierr = DMLabelDestroy(&(*box)->cellsSparse);CHKERRQ(ierr);
351   }
352   ierr = PetscFree(*box);CHKERRQ(ierr);
353   PetscFunctionReturn(0);
354 }
355 
356 PetscErrorCode DMPlexLocatePoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cellStart, PetscInt *cell)
357 {
358   PetscInt       coneSize;
359   PetscErrorCode ierr;
360 
361   PetscFunctionBegin;
362   switch (dim) {
363   case 2:
364     ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr);
365     switch (coneSize) {
366     case 3:
367       ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
368       break;
369     case 4:
370       ierr = DMPlexLocatePoint_General_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
371       break;
372     default:
373       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize);
374     }
375     break;
376   case 3:
377     ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr);
378     switch (coneSize) {
379     case 4:
380       ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
381       break;
382     case 6:
383       ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
384       break;
385     default:
386       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize);
387     }
388     break;
389   default:
390     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %D", dim);
391   }
392   PetscFunctionReturn(0);
393 }
394 
395 /*
396   DMPlexClosestPoint_Internal - Returns the closest point in the cell to the given point
397 */
398 PetscErrorCode DMPlexClosestPoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cell, PetscReal cpoint[])
399 {
400   PetscInt       coneSize;
401   PetscErrorCode ierr;
402 
403   PetscFunctionBegin;
404   switch (dim) {
405   case 2:
406     ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
407     switch (coneSize) {
408     case 3:
409       ierr = DMPlexClosestPoint_Simplex_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
410       break;
411 #if 0
412     case 4:
413       ierr = DMPlexClosestPoint_General_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
414       break;
415 #endif
416     default:
417       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize);
418     }
419     break;
420 #if 0
421   case 3:
422     ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
423     switch (coneSize) {
424     case 4:
425       ierr = DMPlexClosestPoint_Simplex_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
426       break;
427     case 6:
428       ierr = DMPlexClosestPoint_General_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
429       break;
430     default:
431       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize);
432     }
433     break;
434 #endif
435   default:
436     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for mesh dimension %D", dim);
437   }
438   PetscFunctionReturn(0);
439 }
440 
441 /*
442   DMPlexComputeGridHash_Internal - Create a grid hash structure covering the Plex
443 
444   Collective on DM
445 
446   Input Parameter:
447 . dm - The Plex
448 
449   Output Parameter:
450 . localBox - The grid hash object
451 
452   Level: developer
453 
454 .seealso: PetscGridHashCreate(), PetscGridHashGetEnclosingBox()
455 */
456 PetscErrorCode DMPlexComputeGridHash_Internal(DM dm, PetscGridHash *localBox)
457 {
458   MPI_Comm           comm;
459   PetscGridHash      lbox;
460   Vec                coordinates;
461   PetscSection       coordSection;
462   Vec                coordsLocal;
463   const PetscScalar *coords;
464   PetscInt          *dboxes, *boxes;
465   PetscInt           n[3] = {10, 10, 10};
466   PetscInt           dim, N, cStart, cEnd, cMax, c, i;
467   PetscErrorCode     ierr;
468 
469   PetscFunctionBegin;
470   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
471   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
472   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
473   if (dim != 2) SETERRQ(comm, PETSC_ERR_SUP, "I have only coded this for 2D");
474   ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr);
475   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
476   ierr = PetscGridHashCreate(comm, dim, coords, &lbox);CHKERRQ(ierr);
477   for (i = 0; i < N; i += dim) {ierr = PetscGridHashEnlarge(lbox, &coords[i]);CHKERRQ(ierr);}
478   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
479   ierr = PetscOptionsGetInt(NULL,NULL,"-dm_plex_hash_box_nijk",&n[0],NULL);CHKERRQ(ierr);
480   n[1] = n[0];
481   n[2] = n[0];
482   ierr = PetscGridHashSetGrid(lbox, n, NULL);CHKERRQ(ierr);
483 #if 0
484   /* Could define a custom reduction to merge these */
485   ierr = MPIU_Allreduce(lbox->lower, gbox->lower, 3, MPIU_REAL, MPI_MIN, comm);CHKERRQ(ierr);
486   ierr = MPIU_Allreduce(lbox->upper, gbox->upper, 3, MPIU_REAL, MPI_MAX, comm);CHKERRQ(ierr);
487 #endif
488   /* Is there a reason to snap the local bounding box to a division of the global box? */
489   /* Should we compute all overlaps of local boxes? We could do this with a rendevouz scheme partitioning the global box */
490   /* Create label */
491   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
492   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
493   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
494   ierr = DMLabelCreate(PETSC_COMM_SELF, "cells", &lbox->cellsSparse);CHKERRQ(ierr);
495   ierr = DMLabelCreateIndex(lbox->cellsSparse, cStart, cEnd);CHKERRQ(ierr);
496   /* Compute boxes which overlap each cell: https://stackoverflow.com/questions/13790208/triangle-square-intersection-test-in-2d */
497   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
498   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
499   ierr = PetscCalloc2(16 * dim, &dboxes, 16, &boxes);CHKERRQ(ierr);
500   for (c = cStart; c < cEnd; ++c) {
501     const PetscReal *h       = lbox->h;
502     PetscScalar     *ccoords = NULL;
503     PetscInt         csize   = 0;
504     PetscScalar      point[3];
505     PetscInt         dlim[6], d, e, i, j, k;
506 
507     /* Find boxes enclosing each vertex */
508     ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, &csize, &ccoords);CHKERRQ(ierr);
509     ierr = PetscGridHashGetEnclosingBox(lbox, csize/dim, ccoords, dboxes, boxes);CHKERRQ(ierr);
510     /* Mark cells containing the vertices */
511     for (e = 0; e < csize/dim; ++e) {ierr = DMLabelSetValue(lbox->cellsSparse, c, boxes[e]);CHKERRQ(ierr);}
512     /* Get grid of boxes containing these */
513     for (d = 0;   d < dim; ++d) {dlim[d*2+0] = dlim[d*2+1] = dboxes[d];}
514     for (d = dim; d < 3;   ++d) {dlim[d*2+0] = dlim[d*2+1] = 0;}
515     for (e = 1; e < dim+1; ++e) {
516       for (d = 0; d < dim; ++d) {
517         dlim[d*2+0] = PetscMin(dlim[d*2+0], dboxes[e*dim+d]);
518         dlim[d*2+1] = PetscMax(dlim[d*2+1], dboxes[e*dim+d]);
519       }
520     }
521     /* Check for intersection of box with cell */
522     for (k = dlim[2*2+0], point[2] = lbox->lower[2] + k*h[2]; k <= dlim[2*2+1]; ++k, point[2] += h[2]) {
523       for (j = dlim[1*2+0], point[1] = lbox->lower[1] + j*h[1]; j <= dlim[1*2+1]; ++j, point[1] += h[1]) {
524         for (i = dlim[0*2+0], point[0] = lbox->lower[0] + i*h[0]; i <= dlim[0*2+1]; ++i, point[0] += h[0]) {
525           const PetscInt box = (k*lbox->n[1] + j)*lbox->n[0] + i;
526           PetscScalar    cpoint[3];
527           PetscInt       cell, edge, ii, jj, kk;
528 
529           /* Check whether cell contains any vertex of these subboxes TODO vectorize this */
530           for (kk = 0, cpoint[2] = point[2]; kk < (dim > 2 ? 2 : 1); ++kk, cpoint[2] += h[2]) {
531             for (jj = 0, cpoint[1] = point[1]; jj < (dim > 1 ? 2 : 1); ++jj, cpoint[1] += h[1]) {
532               for (ii = 0, cpoint[0] = point[0]; ii < 2; ++ii, cpoint[0] += h[0]) {
533 
534                 ierr = DMPlexLocatePoint_Internal(dm, dim, cpoint, c, &cell);CHKERRQ(ierr);
535                 if (cell >= 0) { ierr = DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); ii = jj = kk = 2;}
536               }
537             }
538           }
539           /* Check whether cell edge intersects any edge of these subboxes TODO vectorize this */
540           for (edge = 0; edge < dim+1; ++edge) {
541             PetscReal segA[6], segB[6];
542 
543             if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected dim %d > 3",dim);
544             for (d = 0; d < dim; ++d) {segA[d] = PetscRealPart(ccoords[edge*dim+d]); segA[dim+d] = PetscRealPart(ccoords[((edge+1)%(dim+1))*dim+d]);}
545             for (kk = 0; kk < (dim > 2 ? 2 : 1); ++kk) {
546               if (dim > 2) {segB[2]     = PetscRealPart(point[2]);
547                             segB[dim+2] = PetscRealPart(point[2]) + kk*h[2];}
548               for (jj = 0; jj < (dim > 1 ? 2 : 1); ++jj) {
549                 if (dim > 1) {segB[1]     = PetscRealPart(point[1]);
550                               segB[dim+1] = PetscRealPart(point[1]) + jj*h[1];}
551                 for (ii = 0; ii < 2; ++ii) {
552                   PetscBool intersects;
553 
554                   segB[0]     = PetscRealPart(point[0]);
555                   segB[dim+0] = PetscRealPart(point[0]) + ii*h[0];
556                   ierr = DMPlexGetLineIntersection_2D_Internal(segA, segB, NULL, &intersects);CHKERRQ(ierr);
557                   if (intersects) { ierr = DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); edge = ii = jj = kk = dim+1;}
558                 }
559               }
560             }
561           }
562         }
563       }
564     }
565     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &ccoords);CHKERRQ(ierr);
566   }
567   ierr = PetscFree2(dboxes, boxes);CHKERRQ(ierr);
568   ierr = DMLabelConvertToSection(lbox->cellsSparse, &lbox->cellSection, &lbox->cells);CHKERRQ(ierr);
569   ierr = DMLabelDestroy(&lbox->cellsSparse);CHKERRQ(ierr);
570   *localBox = lbox;
571   PetscFunctionReturn(0);
572 }
573 
574 PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, DMPointLocationType ltype, PetscSF cellSF)
575 {
576   DM_Plex        *mesh = (DM_Plex *) dm->data;
577   PetscBool       hash = mesh->useHashLocation, reuse = PETSC_FALSE;
578   PetscInt        bs, numPoints, p, numFound, *found = NULL;
579   PetscInt        dim, cStart, cEnd, cMax, numCells, c, d;
580   const PetscInt *boxCells;
581   PetscSFNode    *cells;
582   PetscScalar    *a;
583   PetscMPIInt     result;
584   PetscLogDouble  t0,t1;
585   PetscReal       gmin[3],gmax[3];
586   PetscInt        terminating_query_type[] = { 0, 0, 0 };
587   PetscErrorCode  ierr;
588 
589   PetscFunctionBegin;
590   ierr = PetscTime(&t0);CHKERRQ(ierr);
591   if (ltype == DM_POINTLOCATION_NEAREST && !hash) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Nearest point location only supported with grid hashing. Use -dm_plex_hash_location to enable it.");
592   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
593   ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr);
594   ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)cellSF),PETSC_COMM_SELF,&result);CHKERRQ(ierr);
595   if (result != MPI_IDENT && result != MPI_CONGRUENT) SETERRQ(PetscObjectComm((PetscObject)cellSF),PETSC_ERR_SUP, "Trying parallel point location: only local point location supported");
596   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);
597   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
598   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
599   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
600   ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr);
601   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
602   numPoints /= bs;
603   {
604     const PetscSFNode *sf_cells;
605 
606     ierr = PetscSFGetGraph(cellSF,NULL,NULL,NULL,&sf_cells);CHKERRQ(ierr);
607     if (sf_cells) {
608       ierr = PetscInfo(dm,"[DMLocatePoints_Plex] Re-using existing StarForest node list\n");CHKERRQ(ierr);
609       cells = (PetscSFNode*)sf_cells;
610       reuse = PETSC_TRUE;
611     } else {
612       ierr = PetscInfo(dm,"[DMLocatePoints_Plex] Creating and initializing new StarForest node list\n");CHKERRQ(ierr);
613       ierr = PetscMalloc1(numPoints, &cells);CHKERRQ(ierr);
614       /* initialize cells if created */
615       for (p=0; p<numPoints; p++) {
616         cells[p].rank  = 0;
617         cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
618       }
619     }
620   }
621   /* define domain bounding box */
622   {
623     Vec coorglobal;
624 
625     ierr = DMGetCoordinates(dm,&coorglobal);CHKERRQ(ierr);
626     ierr = VecStrideMaxAll(coorglobal,NULL,gmax);CHKERRQ(ierr);
627     ierr = VecStrideMinAll(coorglobal,NULL,gmin);CHKERRQ(ierr);
628   }
629   if (hash) {
630     if (!mesh->lbox) {ierr = PetscInfo(dm, "Initializing grid hashing");CHKERRQ(ierr);ierr = DMPlexComputeGridHash_Internal(dm, &mesh->lbox);CHKERRQ(ierr);}
631     /* Designate the local box for each point */
632     /* Send points to correct process */
633     /* Search cells that lie in each subbox */
634     /*   Should we bin points before doing search? */
635     ierr = ISGetIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr);
636   }
637   for (p = 0, numFound = 0; p < numPoints; ++p) {
638     const PetscScalar *point = &a[p*bs];
639     PetscInt           dbin[3] = {-1,-1,-1}, bin, cell = -1, cellOffset;
640     PetscBool          point_outside_domain = PETSC_FALSE;
641 
642     /* check bounding box of domain */
643     for (d=0; d<dim; d++) {
644       if (PetscRealPart(point[d]) < gmin[d]) { point_outside_domain = PETSC_TRUE; break; }
645       if (PetscRealPart(point[d]) > gmax[d]) { point_outside_domain = PETSC_TRUE; break; }
646     }
647     if (point_outside_domain) {
648       cells[p].rank = 0;
649       cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
650       terminating_query_type[0]++;
651       continue;
652     }
653 
654     /* check initial values in cells[].index - abort early if found */
655     if (cells[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
656       c = cells[p].index;
657       cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
658       ierr = DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);CHKERRQ(ierr);
659       if (cell >= 0) {
660         cells[p].rank = 0;
661         cells[p].index = cell;
662         numFound++;
663       }
664     }
665     if (cells[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
666       terminating_query_type[1]++;
667       continue;
668     }
669 
670     if (hash) {
671       PetscBool found_box;
672 
673       /* allow for case that point is outside box - abort early */
674       ierr = PetscGridHashGetEnclosingBoxQuery(mesh->lbox, 1, point, dbin, &bin,&found_box);CHKERRQ(ierr);
675       if (found_box) {
676         /* TODO Lay an interface over this so we can switch between Section (dense) and Label (sparse) */
677         ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr);
678         ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr);
679         for (c = cellOffset; c < cellOffset + numCells; ++c) {
680           ierr = DMPlexLocatePoint_Internal(dm, dim, point, boxCells[c], &cell);CHKERRQ(ierr);
681           if (cell >= 0) {
682             cells[p].rank = 0;
683             cells[p].index = cell;
684             numFound++;
685             terminating_query_type[2]++;
686             break;
687           }
688         }
689       }
690     } else {
691       for (c = cStart; c < cEnd; ++c) {
692         ierr = DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);CHKERRQ(ierr);
693         if (cell >= 0) {
694           cells[p].rank = 0;
695           cells[p].index = cell;
696           numFound++;
697           terminating_query_type[2]++;
698           break;
699         }
700       }
701     }
702   }
703   if (hash) {ierr = ISRestoreIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr);}
704   if (ltype == DM_POINTLOCATION_NEAREST && hash && numFound < numPoints) {
705     for (p = 0; p < numPoints; p++) {
706       const PetscScalar *point = &a[p*bs];
707       PetscReal          cpoint[3], diff[3], dist, distMax = PETSC_MAX_REAL;
708       PetscInt           dbin[3] = {-1,-1,-1}, bin, cellOffset, d;
709 
710       if (cells[p].index < 0) {
711         ++numFound;
712         ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr);
713         ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr);
714         ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr);
715         for (c = cellOffset; c < cellOffset + numCells; ++c) {
716           ierr = DMPlexClosestPoint_Internal(dm, dim, point, boxCells[c], cpoint);CHKERRQ(ierr);
717           for (d = 0; d < dim; ++d) diff[d] = cpoint[d] - PetscRealPart(point[d]);
718           dist = DMPlex_NormD_Internal(dim, diff);
719           if (dist < distMax) {
720             for (d = 0; d < dim; ++d) a[p*bs+d] = cpoint[d];
721             cells[p].rank  = 0;
722             cells[p].index = boxCells[c];
723             distMax = dist;
724             break;
725           }
726         }
727       }
728     }
729   }
730   /* This code is only be relevant when interfaced to parallel point location */
731   /* Check for highest numbered proc that claims a point (do we care?) */
732   if (ltype == DM_POINTLOCATION_REMOVE && numFound < numPoints) {
733     ierr = PetscMalloc1(numFound,&found);CHKERRQ(ierr);
734     for (p = 0, numFound = 0; p < numPoints; p++) {
735       if (cells[p].rank >= 0 && cells[p].index >= 0) {
736         if (numFound < p) {
737           cells[numFound] = cells[p];
738         }
739         found[numFound++] = p;
740       }
741     }
742   }
743   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
744   if (!reuse) {
745     ierr = PetscSFSetGraph(cellSF, cEnd - cStart, numFound, found, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER);CHKERRQ(ierr);
746   }
747   ierr = PetscTime(&t1);CHKERRQ(ierr);
748   if (hash) {
749     ierr = PetscInfo3(dm,"[DMLocatePoints_Plex] terminating_query_type : %D [outside domain] : %D [inside intial cell] : %D [hash]\n",terminating_query_type[0],terminating_query_type[1],terminating_query_type[2]);CHKERRQ(ierr);
750   } else {
751     ierr = PetscInfo3(dm,"[DMLocatePoints_Plex] terminating_query_type : %D [outside domain] : %D [inside intial cell] : %D [brute-force]\n",terminating_query_type[0],terminating_query_type[1],terminating_query_type[2]);CHKERRQ(ierr);
752   }
753   ierr = PetscInfo3(dm,"[DMLocatePoints_Plex] npoints %D : time(rank0) %1.2e (sec): points/sec %1.4e\n",numPoints,t1-t0,(double)((double)numPoints/(t1-t0)));CHKERRQ(ierr);
754   PetscFunctionReturn(0);
755 }
756 
757 /*@C
758   DMPlexComputeProjection2Dto1D - Rewrite coordinates to be the 1D projection of the 2D coordinates
759 
760   Not collective
761 
762   Input Parameter:
763 . coords - The coordinates of a segment
764 
765   Output Parameters:
766 + coords - The new y-coordinate, and 0 for x
767 - R - The rotation which accomplishes the projection
768 
769   Level: developer
770 
771 .seealso: DMPlexComputeProjection3Dto1D(), DMPlexComputeProjection3Dto2D()
772 @*/
773 PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar coords[], PetscReal R[])
774 {
775   const PetscReal x = PetscRealPart(coords[2] - coords[0]);
776   const PetscReal y = PetscRealPart(coords[3] - coords[1]);
777   const PetscReal r = PetscSqrtReal(x*x + y*y), c = x/r, s = y/r;
778 
779   PetscFunctionBegin;
780   R[0] = c; R[1] = -s;
781   R[2] = s; R[3] =  c;
782   coords[0] = 0.0;
783   coords[1] = r;
784   PetscFunctionReturn(0);
785 }
786 
787 /*@C
788   DMPlexComputeProjection3Dto1D - Rewrite coordinates to be the 1D projection of the 3D coordinates
789 
790   Not collective
791 
792   Input Parameter:
793 . coords - The coordinates of a segment
794 
795   Output Parameters:
796 + coords - The new y-coordinate, and 0 for x and z
797 - R - The rotation which accomplishes the projection
798 
799   Note: This uses the basis completion described by Frisvad in http://www.imm.dtu.dk/~jerf/papers/abstracts/onb.html, DOI:10.1080/2165347X.2012.689606
800 
801   Level: developer
802 
803 .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto2D()
804 @*/
805 PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar coords[], PetscReal R[])
806 {
807   PetscReal      x    = PetscRealPart(coords[3] - coords[0]);
808   PetscReal      y    = PetscRealPart(coords[4] - coords[1]);
809   PetscReal      z    = PetscRealPart(coords[5] - coords[2]);
810   PetscReal      r    = PetscSqrtReal(x*x + y*y + z*z);
811   PetscReal      rinv = 1. / r;
812   PetscFunctionBegin;
813 
814   x *= rinv; y *= rinv; z *= rinv;
815   if (x > 0.) {
816     PetscReal inv1pX   = 1./ (1. + x);
817 
818     R[0] = x; R[1] = -y;              R[2] = -z;
819     R[3] = y; R[4] = 1. - y*y*inv1pX; R[5] =     -y*z*inv1pX;
820     R[6] = z; R[7] =     -y*z*inv1pX; R[8] = 1. - z*z*inv1pX;
821   }
822   else {
823     PetscReal inv1mX   = 1./ (1. - x);
824 
825     R[0] = x; R[1] = z;               R[2] = y;
826     R[3] = y; R[4] =     -y*z*inv1mX; R[5] = 1. - y*y*inv1mX;
827     R[6] = z; R[7] = 1. - z*z*inv1mX; R[8] =     -y*z*inv1mX;
828   }
829   coords[0] = 0.0;
830   coords[1] = r;
831   PetscFunctionReturn(0);
832 }
833 
834 /*@
835   DMPlexComputeProjection3Dto2D - Rewrite coordinates to be the 2D projection of the 3D coordinates
836 
837   Not collective
838 
839   Input Parameter:
840 . coords - The coordinates of a segment
841 
842   Output Parameters:
843 + coords - The new y- and z-coordinates, and 0 for x
844 - R - The rotation which accomplishes the projection
845 
846   Level: developer
847 
848 .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto1D()
849 @*/
850 PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt coordSize, PetscScalar coords[], PetscReal R[])
851 {
852   PetscReal      x1[3],  x2[3], n[3], norm;
853   PetscReal      x1p[3], x2p[3], xnp[3];
854   PetscReal      sqrtz, alpha;
855   const PetscInt dim = 3;
856   PetscInt       d, e, p;
857 
858   PetscFunctionBegin;
859   /* 0) Calculate normal vector */
860   for (d = 0; d < dim; ++d) {
861     x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]);
862     x2[d] = PetscRealPart(coords[2*dim+d] - coords[0*dim+d]);
863   }
864   n[0] = x1[1]*x2[2] - x1[2]*x2[1];
865   n[1] = x1[2]*x2[0] - x1[0]*x2[2];
866   n[2] = x1[0]*x2[1] - x1[1]*x2[0];
867   norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]);
868   n[0] /= norm;
869   n[1] /= norm;
870   n[2] /= norm;
871   /* 1) Take the normal vector and rotate until it is \hat z
872 
873     Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then
874 
875     R = /  alpha nx nz  alpha ny nz -1/alpha \
876         | -alpha ny     alpha nx        0    |
877         \     nx            ny         nz    /
878 
879     will rotate the normal vector to \hat z
880   */
881   sqrtz = PetscSqrtReal(1.0 - n[2]*n[2]);
882   /* Check for n = z */
883   if (sqrtz < 1.0e-10) {
884     const PetscInt s = PetscSign(n[2]);
885     /* If nz < 0, rotate 180 degrees around x-axis */
886     for (p = 3; p < coordSize/3; ++p) {
887       coords[p*2+0] = PetscRealPart(coords[p*dim+0] - coords[0*dim+0]);
888       coords[p*2+1] = (PetscRealPart(coords[p*dim+1] - coords[0*dim+1])) * s;
889     }
890     coords[0] = 0.0;
891     coords[1] = 0.0;
892     coords[2] = x1[0];
893     coords[3] = x1[1] * s;
894     coords[4] = x2[0];
895     coords[5] = x2[1] * s;
896     R[0] = 1.0;     R[1] = 0.0;     R[2] = 0.0;
897     R[3] = 0.0;     R[4] = 1.0 * s; R[5] = 0.0;
898     R[6] = 0.0;     R[7] = 0.0;     R[8] = 1.0 * s;
899     PetscFunctionReturn(0);
900   }
901   alpha = 1.0/sqrtz;
902   R[0] =  alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz;
903   R[3] = -alpha*n[1];      R[4] = alpha*n[0];      R[5] = 0.0;
904   R[6] =  n[0];            R[7] = n[1];            R[8] = n[2];
905   for (d = 0; d < dim; ++d) {
906     x1p[d] = 0.0;
907     x2p[d] = 0.0;
908     for (e = 0; e < dim; ++e) {
909       x1p[d] += R[d*dim+e]*x1[e];
910       x2p[d] += R[d*dim+e]*x2[e];
911     }
912   }
913   if (PetscAbsReal(x1p[2]) > 10. * PETSC_SMALL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
914   if (PetscAbsReal(x2p[2]) > 10. * PETSC_SMALL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
915   /* 2) Project to (x, y) */
916   for (p = 3; p < coordSize/3; ++p) {
917     for (d = 0; d < dim; ++d) {
918       xnp[d] = 0.0;
919       for (e = 0; e < dim; ++e) {
920         xnp[d] += R[d*dim+e]*PetscRealPart(coords[p*dim+e] - coords[0*dim+e]);
921       }
922       if (d < dim-1) coords[p*2+d] = xnp[d];
923     }
924   }
925   coords[0] = 0.0;
926   coords[1] = 0.0;
927   coords[2] = x1p[0];
928   coords[3] = x1p[1];
929   coords[4] = x2p[0];
930   coords[5] = x2p[1];
931   /* Output R^T which rotates \hat z to the input normal */
932   for (d = 0; d < dim; ++d) {
933     for (e = d+1; e < dim; ++e) {
934       PetscReal tmp;
935 
936       tmp        = R[d*dim+e];
937       R[d*dim+e] = R[e*dim+d];
938       R[e*dim+d] = tmp;
939     }
940   }
941   PetscFunctionReturn(0);
942 }
943 
944 PETSC_UNUSED
945 PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[])
946 {
947   /* Signed volume is 1/2 the determinant
948 
949    |  1  1  1 |
950    | x0 x1 x2 |
951    | y0 y1 y2 |
952 
953      but if x0,y0 is the origin, we have
954 
955    | x1 x2 |
956    | y1 y2 |
957   */
958   const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1];
959   const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1];
960   PetscReal       M[4], detM;
961   M[0] = x1; M[1] = x2;
962   M[2] = y1; M[3] = y2;
963   DMPlex_Det2D_Internal(&detM, M);
964   *vol = 0.5*detM;
965   (void)PetscLogFlops(5.0);
966 }
967 
968 PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[])
969 {
970   DMPlex_Det2D_Internal(vol, coords);
971   *vol *= 0.5;
972 }
973 
974 PETSC_UNUSED
975 PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[])
976 {
977   /* Signed volume is 1/6th of the determinant
978 
979    |  1  1  1  1 |
980    | x0 x1 x2 x3 |
981    | y0 y1 y2 y3 |
982    | z0 z1 z2 z3 |
983 
984      but if x0,y0,z0 is the origin, we have
985 
986    | x1 x2 x3 |
987    | y1 y2 y3 |
988    | z1 z2 z3 |
989   */
990   const PetscReal x1 = coords[3] - coords[0], y1 = coords[4]  - coords[1], z1 = coords[5]  - coords[2];
991   const PetscReal x2 = coords[6] - coords[0], y2 = coords[7]  - coords[1], z2 = coords[8]  - coords[2];
992   const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2];
993   const PetscReal onesixth = ((PetscReal)1./(PetscReal)6.);
994   PetscReal       M[9], detM;
995   M[0] = x1; M[1] = x2; M[2] = x3;
996   M[3] = y1; M[4] = y2; M[5] = y3;
997   M[6] = z1; M[7] = z2; M[8] = z3;
998   DMPlex_Det3D_Internal(&detM, M);
999   *vol = -onesixth*detM;
1000   (void)PetscLogFlops(10.0);
1001 }
1002 
1003 PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[])
1004 {
1005   const PetscReal onesixth = ((PetscReal)1./(PetscReal)6.);
1006   DMPlex_Det3D_Internal(vol, coords);
1007   *vol *= -onesixth;
1008 }
1009 
1010 static PetscErrorCode DMPlexComputePointGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1011 {
1012   PetscSection   coordSection;
1013   Vec            coordinates;
1014   const PetscScalar *coords;
1015   PetscInt       dim, d, off;
1016   PetscErrorCode ierr;
1017 
1018   PetscFunctionBegin;
1019   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1020   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1021   ierr = PetscSectionGetDof(coordSection,e,&dim);CHKERRQ(ierr);
1022   if (!dim) PetscFunctionReturn(0);
1023   ierr = PetscSectionGetOffset(coordSection,e,&off);CHKERRQ(ierr);
1024   ierr = VecGetArrayRead(coordinates,&coords);CHKERRQ(ierr);
1025   if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[off + d]);}
1026   ierr = VecRestoreArrayRead(coordinates,&coords);CHKERRQ(ierr);
1027   *detJ = 1.;
1028   if (J) {
1029     for (d = 0; d < dim * dim; d++) J[d] = 0.;
1030     for (d = 0; d < dim; d++) J[d * dim + d] = 1.;
1031     if (invJ) {
1032       for (d = 0; d < dim * dim; d++) invJ[d] = 0.;
1033       for (d = 0; d < dim; d++) invJ[d * dim + d] = 1.;
1034     }
1035   }
1036   PetscFunctionReturn(0);
1037 }
1038 
1039 static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1040 {
1041   PetscSection   coordSection;
1042   Vec            coordinates;
1043   PetscScalar   *coords = NULL;
1044   PetscInt       numCoords, d, pStart, pEnd, numSelfCoords = 0;
1045   PetscErrorCode ierr;
1046 
1047   PetscFunctionBegin;
1048   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1049   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1050   ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr);
1051   if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);}
1052   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
1053   numCoords = numSelfCoords ? numSelfCoords : numCoords;
1054   if (invJ && !J) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL");
1055   *detJ = 0.0;
1056   if (numCoords == 6) {
1057     const PetscInt dim = 3;
1058     PetscReal      R[9], J0;
1059 
1060     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1061     ierr = DMPlexComputeProjection3Dto1D(coords, R);CHKERRQ(ierr);
1062     if (J)    {
1063       J0   = 0.5*PetscRealPart(coords[1]);
1064       J[0] = R[0]*J0; J[1] = R[1]; J[2] = R[2];
1065       J[3] = R[3]*J0; J[4] = R[4]; J[5] = R[5];
1066       J[6] = R[6]*J0; J[7] = R[7]; J[8] = R[8];
1067       DMPlex_Det3D_Internal(detJ, J);
1068       if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1069     }
1070   } else if (numCoords == 4) {
1071     const PetscInt dim = 2;
1072     PetscReal      R[4], J0;
1073 
1074     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1075     ierr = DMPlexComputeProjection2Dto1D(coords, R);CHKERRQ(ierr);
1076     if (J)    {
1077       J0   = 0.5*PetscRealPart(coords[1]);
1078       J[0] = R[0]*J0; J[1] = R[1];
1079       J[2] = R[2]*J0; J[3] = R[3];
1080       DMPlex_Det2D_Internal(detJ, J);
1081       if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1082     }
1083   } else if (numCoords == 2) {
1084     const PetscInt dim = 1;
1085 
1086     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1087     if (J)    {
1088       J[0]  = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0]));
1089       *detJ = J[0];
1090       ierr = PetscLogFlops(2.0);CHKERRQ(ierr);
1091       if (invJ) {invJ[0] = 1.0/J[0]; ierr = PetscLogFlops(1.0);CHKERRQ(ierr);}
1092     }
1093   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %D != 2", numCoords);
1094   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
1095   PetscFunctionReturn(0);
1096 }
1097 
1098 static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1099 {
1100   PetscSection   coordSection;
1101   Vec            coordinates;
1102   PetscScalar   *coords = NULL;
1103   PetscInt       numCoords, numSelfCoords = 0, d, f, g, pStart, pEnd;
1104   PetscErrorCode ierr;
1105 
1106   PetscFunctionBegin;
1107   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1108   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1109   ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr);
1110   if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);}
1111   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
1112   numCoords = numSelfCoords ? numSelfCoords : numCoords;
1113   *detJ = 0.0;
1114   if (numCoords == 9) {
1115     const PetscInt dim = 3;
1116     PetscReal      R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0};
1117 
1118     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1119     ierr = DMPlexComputeProjection3Dto2D(numCoords, coords, R);CHKERRQ(ierr);
1120     if (J)    {
1121       const PetscInt pdim = 2;
1122 
1123       for (d = 0; d < pdim; d++) {
1124         for (f = 0; f < pdim; f++) {
1125           J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
1126         }
1127       }
1128       ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1129       DMPlex_Det3D_Internal(detJ, J0);
1130       for (d = 0; d < dim; d++) {
1131         for (f = 0; f < dim; f++) {
1132           J[d*dim+f] = 0.0;
1133           for (g = 0; g < dim; g++) {
1134             J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
1135           }
1136         }
1137       }
1138       ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
1139     }
1140     if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1141   } else if (numCoords == 6) {
1142     const PetscInt dim = 2;
1143 
1144     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1145     if (J)    {
1146       for (d = 0; d < dim; d++) {
1147         for (f = 0; f < dim; f++) {
1148           J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
1149         }
1150       }
1151       ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1152       DMPlex_Det2D_Internal(detJ, J);
1153     }
1154     if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1155   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %D != 6 or 9", numCoords);
1156   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
1157   PetscFunctionReturn(0);
1158 }
1159 
1160 static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1161 {
1162   PetscSection   coordSection;
1163   Vec            coordinates;
1164   PetscScalar   *coords = NULL;
1165   PetscInt       numCoords, numSelfCoords = 0, d, f, g, pStart, pEnd;
1166   PetscErrorCode ierr;
1167 
1168   PetscFunctionBegin;
1169   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1170   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1171   ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr);
1172   if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);}
1173   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
1174   numCoords = numSelfCoords ? numSelfCoords : numCoords;
1175   if (!Nq) {
1176     *detJ = 0.0;
1177     if (numCoords == 12) {
1178       const PetscInt dim = 3;
1179       PetscReal      R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0};
1180 
1181       if (v)   {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);}
1182       ierr = DMPlexComputeProjection3Dto2D(numCoords, coords, R);CHKERRQ(ierr);
1183       if (J)    {
1184         const PetscInt pdim = 2;
1185 
1186         for (d = 0; d < pdim; d++) {
1187           J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
1188           J0[d*dim+1] = 0.5*(PetscRealPart(coords[2*pdim+d]) - PetscRealPart(coords[1*pdim+d]));
1189         }
1190         ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1191         DMPlex_Det3D_Internal(detJ, J0);
1192         for (d = 0; d < dim; d++) {
1193           for (f = 0; f < dim; f++) {
1194             J[d*dim+f] = 0.0;
1195             for (g = 0; g < dim; g++) {
1196               J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
1197             }
1198           }
1199         }
1200         ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
1201       }
1202       if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1203     } else if (numCoords == 8) {
1204       const PetscInt dim = 2;
1205 
1206       if (v)   {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);}
1207       if (J)    {
1208         for (d = 0; d < dim; d++) {
1209           J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
1210           J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1211         }
1212         ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1213         DMPlex_Det2D_Internal(detJ, J);
1214       }
1215       if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1216     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords);
1217   } else {
1218     const PetscInt Nv = 4;
1219     const PetscInt dimR = 2;
1220     const PetscInt zToPlex[4] = {0, 1, 3, 2};
1221     PetscReal zOrder[12];
1222     PetscReal zCoeff[12];
1223     PetscInt  i, j, k, l, dim;
1224 
1225     if (numCoords == 12) {
1226       dim = 3;
1227     } else if (numCoords == 8) {
1228       dim = 2;
1229     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords);
1230     for (i = 0; i < Nv; i++) {
1231       PetscInt zi = zToPlex[i];
1232 
1233       for (j = 0; j < dim; j++) {
1234         zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]);
1235       }
1236     }
1237     for (j = 0; j < dim; j++) {
1238       zCoeff[dim * 0 + j] = 0.25 * (  zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1239       zCoeff[dim * 1 + j] = 0.25 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1240       zCoeff[dim * 2 + j] = 0.25 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1241       zCoeff[dim * 3 + j] = 0.25 * (  zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1242     }
1243     for (i = 0; i < Nq; i++) {
1244       PetscReal xi = points[dimR * i], eta = points[dimR * i + 1];
1245 
1246       if (v) {
1247         PetscReal extPoint[4];
1248 
1249         extPoint[0] = 1.;
1250         extPoint[1] = xi;
1251         extPoint[2] = eta;
1252         extPoint[3] = xi * eta;
1253         for (j = 0; j < dim; j++) {
1254           PetscReal val = 0.;
1255 
1256           for (k = 0; k < Nv; k++) {
1257             val += extPoint[k] * zCoeff[dim * k + j];
1258           }
1259           v[i * dim + j] = val;
1260         }
1261       }
1262       if (J) {
1263         PetscReal extJ[8];
1264 
1265         extJ[0] = 0.;
1266         extJ[1] = 0.;
1267         extJ[2] = 1.;
1268         extJ[3] = 0.;
1269         extJ[4] = 0.;
1270         extJ[5] = 1.;
1271         extJ[6] = eta;
1272         extJ[7] = xi;
1273         for (j = 0; j < dim; j++) {
1274           for (k = 0; k < dimR; k++) {
1275             PetscReal val = 0.;
1276 
1277             for (l = 0; l < Nv; l++) {
1278               val += zCoeff[dim * l + j] * extJ[dimR * l + k];
1279             }
1280             J[i * dim * dim + dim * j + k] = val;
1281           }
1282         }
1283         if (dim == 3) { /* put the cross product in the third component of the Jacobian */
1284           PetscReal x, y, z;
1285           PetscReal *iJ = &J[i * dim * dim];
1286           PetscReal norm;
1287 
1288           x = iJ[1 * dim + 0] * iJ[2 * dim + 1] - iJ[1 * dim + 1] * iJ[2 * dim + 0];
1289           y = iJ[0 * dim + 1] * iJ[2 * dim + 0] - iJ[0 * dim + 0] * iJ[2 * dim + 1];
1290           z = iJ[0 * dim + 0] * iJ[1 * dim + 1] - iJ[0 * dim + 1] * iJ[1 * dim + 0];
1291           norm = PetscSqrtReal(x * x + y * y + z * z);
1292           iJ[2] = x / norm;
1293           iJ[5] = y / norm;
1294           iJ[8] = z / norm;
1295           DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]);
1296           if (invJ) {DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);}
1297         } else {
1298           DMPlex_Det2D_Internal(&detJ[i], &J[i * dim * dim]);
1299           if (invJ) {DMPlex_Invert2D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);}
1300         }
1301       }
1302     }
1303   }
1304   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
1305   PetscFunctionReturn(0);
1306 }
1307 
1308 static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1309 {
1310   PetscSection   coordSection;
1311   Vec            coordinates;
1312   PetscScalar   *coords = NULL;
1313   const PetscInt dim = 3;
1314   PetscInt       d;
1315   PetscErrorCode ierr;
1316 
1317   PetscFunctionBegin;
1318   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1319   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1320   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1321   *detJ = 0.0;
1322   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1323   if (J)    {
1324     for (d = 0; d < dim; d++) {
1325       /* I orient with outward face normals */
1326       J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d]));
1327       J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
1328       J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1329     }
1330     ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
1331     DMPlex_Det3D_Internal(detJ, J);
1332   }
1333   if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1334   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1335   PetscFunctionReturn(0);
1336 }
1337 
1338 static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1339 {
1340   PetscSection   coordSection;
1341   Vec            coordinates;
1342   PetscScalar   *coords = NULL;
1343   const PetscInt dim = 3;
1344   PetscInt       d;
1345   PetscErrorCode ierr;
1346 
1347   PetscFunctionBegin;
1348   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1349   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1350   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1351   if (!Nq) {
1352     *detJ = 0.0;
1353     if (v)   {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);}
1354     if (J)    {
1355       for (d = 0; d < dim; d++) {
1356         J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1357         J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
1358         J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d]));
1359       }
1360       ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
1361       DMPlex_Det3D_Internal(detJ, J);
1362     }
1363     if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1364   } else {
1365     const PetscInt Nv = 8;
1366     const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6};
1367     const PetscInt dim = 3;
1368     const PetscInt dimR = 3;
1369     PetscReal zOrder[24];
1370     PetscReal zCoeff[24];
1371     PetscInt  i, j, k, l;
1372 
1373     for (i = 0; i < Nv; i++) {
1374       PetscInt zi = zToPlex[i];
1375 
1376       for (j = 0; j < dim; j++) {
1377         zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]);
1378       }
1379     }
1380     for (j = 0; j < dim; j++) {
1381       zCoeff[dim * 0 + j] = 0.125 * (  zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j] + zOrder[dim * 4 + j] + zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1382       zCoeff[dim * 1 + j] = 0.125 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j] - zOrder[dim * 4 + j] + zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1383       zCoeff[dim * 2 + j] = 0.125 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j] - zOrder[dim * 4 + j] - zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1384       zCoeff[dim * 3 + j] = 0.125 * (  zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j] + zOrder[dim * 4 + j] - zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1385       zCoeff[dim * 4 + j] = 0.125 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] - zOrder[dim * 3 + j] + zOrder[dim * 4 + j] + zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1386       zCoeff[dim * 5 + j] = 0.125 * (+ zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] - zOrder[dim * 3 + j] - zOrder[dim * 4 + j] + zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1387       zCoeff[dim * 6 + j] = 0.125 * (+ zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] - zOrder[dim * 3 + j] - zOrder[dim * 4 + j] - zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1388       zCoeff[dim * 7 + j] = 0.125 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] - zOrder[dim * 3 + j] + zOrder[dim * 4 + j] - zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1389     }
1390     for (i = 0; i < Nq; i++) {
1391       PetscReal xi = points[dimR * i], eta = points[dimR * i + 1], theta = points[dimR * i + 2];
1392 
1393       if (v) {
1394         PetscReal extPoint[8];
1395 
1396         extPoint[0] = 1.;
1397         extPoint[1] = xi;
1398         extPoint[2] = eta;
1399         extPoint[3] = xi * eta;
1400         extPoint[4] = theta;
1401         extPoint[5] = theta * xi;
1402         extPoint[6] = theta * eta;
1403         extPoint[7] = theta * eta * xi;
1404         for (j = 0; j < dim; j++) {
1405           PetscReal val = 0.;
1406 
1407           for (k = 0; k < Nv; k++) {
1408             val += extPoint[k] * zCoeff[dim * k + j];
1409           }
1410           v[i * dim + j] = val;
1411         }
1412       }
1413       if (J) {
1414         PetscReal extJ[24];
1415 
1416         extJ[0]  = 0.         ; extJ[1]  = 0.        ; extJ[2]  = 0.      ;
1417         extJ[3]  = 1.         ; extJ[4]  = 0.        ; extJ[5]  = 0.      ;
1418         extJ[6]  = 0.         ; extJ[7]  = 1.        ; extJ[8]  = 0.      ;
1419         extJ[9]  = eta        ; extJ[10] = xi        ; extJ[11] = 0.      ;
1420         extJ[12] = 0.         ; extJ[13] = 0.        ; extJ[14] = 1.      ;
1421         extJ[15] = theta      ; extJ[16] = 0.        ; extJ[17] = xi      ;
1422         extJ[18] = 0.         ; extJ[19] = theta     ; extJ[20] = eta     ;
1423         extJ[21] = theta * eta; extJ[22] = theta * xi; extJ[23] = eta * xi;
1424 
1425         for (j = 0; j < dim; j++) {
1426           for (k = 0; k < dimR; k++) {
1427             PetscReal val = 0.;
1428 
1429             for (l = 0; l < Nv; l++) {
1430               val += zCoeff[dim * l + j] * extJ[dimR * l + k];
1431             }
1432             J[i * dim * dim + dim * j + k] = val;
1433           }
1434         }
1435         DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]);
1436         if (invJ) {DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);}
1437       }
1438     }
1439   }
1440   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1441   PetscFunctionReturn(0);
1442 }
1443 
1444 static PetscErrorCode DMPlexComputeCellGeometryFEM_Implicit(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
1445 {
1446   PetscInt        depth, dim, coordDim, coneSize, i;
1447   PetscInt        Nq = 0;
1448   const PetscReal *points = NULL;
1449   DMLabel         depthLabel;
1450   PetscReal       xi0[3] = {-1.,-1.,-1.}, v0[3], J0[9], detJ0;
1451   PetscBool       isAffine = PETSC_TRUE;
1452   PetscErrorCode  ierr;
1453 
1454   PetscFunctionBegin;
1455   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1456   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
1457   ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr);
1458   ierr = DMLabelGetValue(depthLabel, cell, &dim);CHKERRQ(ierr);
1459   if (depth == 1 && dim == 1) {
1460     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1461   }
1462   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
1463   if (coordDim > 3) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported coordinate dimension %D > 3", coordDim);
1464   if (quad) {ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, &points, NULL);CHKERRQ(ierr);}
1465   switch (dim) {
1466   case 0:
1467     ierr = DMPlexComputePointGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr);
1468     isAffine = PETSC_FALSE;
1469     break;
1470   case 1:
1471     if (Nq) {
1472       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr);
1473     } else {
1474       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr);
1475     }
1476     break;
1477   case 2:
1478     switch (coneSize) {
1479     case 3:
1480       if (Nq) {
1481         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr);
1482       } else {
1483         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr);
1484       }
1485       break;
1486     case 4:
1487       ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, Nq, points, v, J, invJ, detJ);CHKERRQ(ierr);
1488       isAffine = PETSC_FALSE;
1489       break;
1490     default:
1491       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell);
1492     }
1493     break;
1494   case 3:
1495     switch (coneSize) {
1496     case 4:
1497       if (Nq) {
1498         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr);
1499       } else {
1500         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr);
1501       }
1502       break;
1503     case 6: /* Faces */
1504     case 8: /* Vertices */
1505       ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, Nq, points, v, J, invJ, detJ);CHKERRQ(ierr);
1506       isAffine = PETSC_FALSE;
1507       break;
1508     default:
1509       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell);
1510     }
1511     break;
1512   default:
1513     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
1514   }
1515   if (isAffine && Nq) {
1516     if (v) {
1517       for (i = 0; i < Nq; i++) {
1518         CoordinatesRefToReal(coordDim, dim, xi0, v0, J0, &points[dim * i], &v[coordDim * i]);
1519       }
1520     }
1521     if (detJ) {
1522       for (i = 0; i < Nq; i++) {
1523         detJ[i] = detJ0;
1524       }
1525     }
1526     if (J) {
1527       PetscInt k;
1528 
1529       for (i = 0, k = 0; i < Nq; i++) {
1530         PetscInt j;
1531 
1532         for (j = 0; j < coordDim * coordDim; j++, k++) {
1533           J[k] = J0[j];
1534         }
1535       }
1536     }
1537     if (invJ) {
1538       PetscInt k;
1539       switch (coordDim) {
1540       case 0:
1541         break;
1542       case 1:
1543         invJ[0] = 1./J0[0];
1544         break;
1545       case 2:
1546         DMPlex_Invert2D_Internal(invJ, J0, detJ0);
1547         break;
1548       case 3:
1549         DMPlex_Invert3D_Internal(invJ, J0, detJ0);
1550         break;
1551       }
1552       for (i = 1, k = coordDim * coordDim; i < Nq; i++) {
1553         PetscInt j;
1554 
1555         for (j = 0; j < coordDim * coordDim; j++, k++) {
1556           invJ[k] = invJ[j];
1557         }
1558       }
1559     }
1560   }
1561   PetscFunctionReturn(0);
1562 }
1563 
1564 /*@C
1565   DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell
1566 
1567   Collective on DM
1568 
1569   Input Arguments:
1570 + dm   - the DM
1571 - cell - the cell
1572 
1573   Output Arguments:
1574 + v0   - the translation part of this affine transform
1575 . J    - the Jacobian of the transform from the reference element
1576 . invJ - the inverse of the Jacobian
1577 - detJ - the Jacobian determinant
1578 
1579   Level: advanced
1580 
1581   Fortran Notes:
1582   Since it returns arrays, this routine is only available in Fortran 90, and you must
1583   include petsc.h90 in your code.
1584 
1585 .seealso: DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinates()
1586 @*/
1587 PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
1588 {
1589   PetscErrorCode ierr;
1590 
1591   PetscFunctionBegin;
1592   ierr = DMPlexComputeCellGeometryFEM_Implicit(dm,cell,NULL,v0,J,invJ,detJ);CHKERRQ(ierr);
1593   PetscFunctionReturn(0);
1594 }
1595 
1596 static PetscErrorCode DMPlexComputeCellGeometryFEM_FE(DM dm, PetscFE fe, PetscInt point, PetscQuadrature quad, PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1597 {
1598   PetscQuadrature  feQuad;
1599   PetscSection     coordSection;
1600   Vec              coordinates;
1601   PetscScalar     *coords = NULL;
1602   const PetscReal *quadPoints;
1603   PetscReal       *basisDer, *basis, detJt;
1604   PetscInt         dim, cdim, pdim, qdim, Nq, numCoords, q;
1605   PetscErrorCode   ierr;
1606 
1607   PetscFunctionBegin;
1608   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1609   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1610   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr);
1611   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1612   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
1613   if (!quad) { /* use the first point of the first functional of the dual space */
1614     PetscDualSpace dsp;
1615 
1616     ierr = PetscFEGetDualSpace(fe, &dsp);CHKERRQ(ierr);
1617     ierr = PetscDualSpaceGetFunctional(dsp, 0, &quad);CHKERRQ(ierr);
1618     ierr = PetscQuadratureGetData(quad, &qdim, NULL, &Nq, &quadPoints, NULL);CHKERRQ(ierr);
1619     Nq = 1;
1620   } else {
1621     ierr = PetscQuadratureGetData(quad, &qdim, NULL, &Nq, &quadPoints, NULL);CHKERRQ(ierr);
1622   }
1623   ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr);
1624   ierr = PetscFEGetQuadrature(fe, &feQuad);CHKERRQ(ierr);
1625   if (feQuad == quad) {
1626     ierr = PetscFEGetDefaultTabulation(fe, &basis, J ? &basisDer : NULL, NULL);CHKERRQ(ierr);
1627     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);
1628   } else {
1629     ierr = PetscFEGetTabulation(fe, Nq, quadPoints, &basis, J ? &basisDer : NULL, NULL);CHKERRQ(ierr);
1630   }
1631   if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %d != quadrature dimension %d", dim, qdim);
1632   if (v) {
1633     ierr = PetscMemzero(v, Nq*cdim*sizeof(PetscReal));CHKERRQ(ierr);
1634     for (q = 0; q < Nq; ++q) {
1635       PetscInt i, k;
1636 
1637       for (k = 0; k < pdim; ++k)
1638         for (i = 0; i < cdim; ++i)
1639           v[q*cdim + i] += basis[q*pdim + k] * PetscRealPart(coords[k*cdim + i]);
1640       ierr = PetscLogFlops(2.0*pdim*cdim);CHKERRQ(ierr);
1641     }
1642   }
1643   if (J) {
1644     ierr = PetscMemzero(J, Nq*cdim*cdim*sizeof(PetscReal));CHKERRQ(ierr);
1645     for (q = 0; q < Nq; ++q) {
1646       PetscInt i, j, k, c, r;
1647 
1648       /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */
1649       for (k = 0; k < pdim; ++k)
1650         for (j = 0; j < dim; ++j)
1651           for (i = 0; i < cdim; ++i)
1652             J[(q*cdim + i)*cdim + j] += basisDer[(q*pdim + k)*dim + j] * PetscRealPart(coords[k*cdim + i]);
1653       ierr = PetscLogFlops(2.0*pdim*dim*cdim);CHKERRQ(ierr);
1654       if (cdim > dim) {
1655         for (c = dim; c < cdim; ++c)
1656           for (r = 0; r < cdim; ++r)
1657             J[r*cdim+c] = r == c ? 1.0 : 0.0;
1658       }
1659       if (!detJ && !invJ) continue;
1660       detJt = 0.;
1661       switch (cdim) {
1662       case 3:
1663         DMPlex_Det3D_Internal(&detJt, &J[q*cdim*dim]);
1664         if (invJ) {DMPlex_Invert3D_Internal(&invJ[q*cdim*dim], &J[q*cdim*dim], detJt);}
1665         break;
1666       case 2:
1667         DMPlex_Det2D_Internal(&detJt, &J[q*cdim*dim]);
1668         if (invJ) {DMPlex_Invert2D_Internal(&invJ[q*cdim*dim], &J[q*cdim*dim], detJt);}
1669         break;
1670       case 1:
1671         detJt = J[q*cdim*dim];
1672         if (invJ) invJ[q*cdim*dim] = 1.0/detJt;
1673       }
1674       if (detJ) detJ[q] = detJt;
1675     }
1676   }
1677   else if (detJ || invJ) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Need J to compute invJ or detJ");
1678   if (feQuad != quad) {
1679     ierr = PetscFERestoreTabulation(fe, Nq, quadPoints, &basis, J ? &basisDer : NULL, NULL);CHKERRQ(ierr);
1680   }
1681   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr);
1682   PetscFunctionReturn(0);
1683 }
1684 
1685 /*@C
1686   DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell
1687 
1688   Collective on DM
1689 
1690   Input Arguments:
1691 + dm   - the DM
1692 . cell - the cell
1693 - quad - the quadrature containing the points in the reference element where the geometry will be evaluated.  If quad == NULL, geometry will be
1694          evaluated at the first vertex of the reference element
1695 
1696   Output Arguments:
1697 + v    - the image of the transformed quadrature points, otherwise the image of the first vertex in the closure of the reference element
1698 . J    - the Jacobian of the transform from the reference element at each quadrature point
1699 . invJ - the inverse of the Jacobian at each quadrature point
1700 - detJ - the Jacobian determinant at each quadrature point
1701 
1702   Level: advanced
1703 
1704   Fortran Notes:
1705   Since it returns arrays, this routine is only available in Fortran 90, and you must
1706   include petsc.h90 in your code.
1707 
1708 .seealso: DMGetCoordinateSection(), DMGetCoordinates()
1709 @*/
1710 PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
1711 {
1712   DM             cdm;
1713   PetscFE        fe = NULL;
1714   PetscErrorCode ierr;
1715 
1716   PetscFunctionBegin;
1717   PetscValidPointer(detJ, 7);
1718   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
1719   if (cdm) {
1720     PetscClassId id;
1721     PetscInt     numFields;
1722     PetscDS      prob;
1723     PetscObject  disc;
1724 
1725     ierr = DMGetNumFields(cdm, &numFields);CHKERRQ(ierr);
1726     if (numFields) {
1727       ierr = DMGetDS(cdm, &prob);CHKERRQ(ierr);
1728       ierr = PetscDSGetDiscretization(prob,0,&disc);CHKERRQ(ierr);
1729       ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
1730       if (id == PETSCFE_CLASSID) {
1731         fe = (PetscFE) disc;
1732       }
1733     }
1734   }
1735   if (!fe) {ierr = DMPlexComputeCellGeometryFEM_Implicit(dm, cell, quad, v, J, invJ, detJ);CHKERRQ(ierr);}
1736   else     {ierr = DMPlexComputeCellGeometryFEM_FE(dm, fe, cell, quad, v, J, invJ, detJ);CHKERRQ(ierr);}
1737   PetscFunctionReturn(0);
1738 }
1739 
1740 static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1741 {
1742   PetscSection   coordSection;
1743   Vec            coordinates;
1744   PetscScalar   *coords = NULL;
1745   PetscScalar    tmp[2];
1746   PetscInt       coordSize;
1747   PetscErrorCode ierr;
1748 
1749   PetscFunctionBegin;
1750   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1751   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1752   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1753   if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now");
1754   ierr = DMLocalizeCoordinate_Internal(dm, dim, coords, &coords[dim], tmp);CHKERRQ(ierr);
1755   if (centroid) {
1756     centroid[0] = 0.5*PetscRealPart(coords[0] + tmp[0]);
1757     centroid[1] = 0.5*PetscRealPart(coords[1] + tmp[1]);
1758   }
1759   if (normal) {
1760     PetscReal norm;
1761 
1762     normal[0]  = -PetscRealPart(coords[1] - tmp[1]);
1763     normal[1]  =  PetscRealPart(coords[0] - tmp[0]);
1764     norm       = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1]);
1765     normal[0] /= norm;
1766     normal[1] /= norm;
1767   }
1768   if (vol) {
1769     *vol = PetscSqrtReal(PetscSqr(PetscRealPart(coords[0] - tmp[0])) + PetscSqr(PetscRealPart(coords[1] - tmp[1])));
1770   }
1771   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1772   PetscFunctionReturn(0);
1773 }
1774 
1775 /* Centroid_i = (\sum_n A_n Cn_i ) / A */
1776 static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1777 {
1778   DMLabel        depth;
1779   PetscSection   coordSection;
1780   Vec            coordinates;
1781   PetscScalar   *coords = NULL;
1782   PetscReal      vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9];
1783   PetscBool      isHybrid = PETSC_FALSE;
1784   PetscInt       fv[4] = {0, 1, 2, 3};
1785   PetscInt       cdepth, pEndInterior[4], tdim = 2, coordSize, numCorners, p, d, e;
1786   PetscErrorCode ierr;
1787 
1788   PetscFunctionBegin;
1789   /* Must check for hybrid cells because prisms have a different orientation scheme */
1790   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
1791   ierr = DMLabelGetValue(depth, cell, &cdepth);CHKERRQ(ierr);
1792   ierr = DMPlexGetHybridBounds(dm, &pEndInterior[3], &pEndInterior[2], &pEndInterior[1], &pEndInterior[0]);CHKERRQ(ierr);
1793   if ((pEndInterior[cdepth]) >=0 && (cell >= pEndInterior[cdepth])) isHybrid = PETSC_TRUE;
1794 
1795   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1796   ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr);
1797   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1798   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1799   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
1800   /* Side faces for hybrid cells are are stored as tensor products */
1801   if (isHybrid && numCorners == 4) {fv[2] = 3; fv[3] = 2;}
1802 
1803   if (dim > 2 && centroid) {
1804     v0[0] = PetscRealPart(coords[0]);
1805     v0[1] = PetscRealPart(coords[1]);
1806     v0[2] = PetscRealPart(coords[2]);
1807   }
1808   if (normal) {
1809     if (dim > 2) {
1810       const PetscReal x0 = PetscRealPart(coords[dim*fv[1]+0] - coords[0]), x1 = PetscRealPart(coords[dim*fv[2]+0] - coords[0]);
1811       const PetscReal y0 = PetscRealPart(coords[dim*fv[1]+1] - coords[1]), y1 = PetscRealPart(coords[dim*fv[2]+1] - coords[1]);
1812       const PetscReal z0 = PetscRealPart(coords[dim*fv[1]+2] - coords[2]), z1 = PetscRealPart(coords[dim*fv[2]+2] - coords[2]);
1813       PetscReal       norm;
1814 
1815       normal[0] = y0*z1 - z0*y1;
1816       normal[1] = z0*x1 - x0*z1;
1817       normal[2] = x0*y1 - y0*x1;
1818       norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]);
1819       normal[0] /= norm;
1820       normal[1] /= norm;
1821       normal[2] /= norm;
1822     } else {
1823       for (d = 0; d < dim; ++d) normal[d] = 0.0;
1824     }
1825   }
1826   if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D(coordSize, coords, R);CHKERRQ(ierr);}
1827   for (p = 0; p < numCorners; ++p) {
1828     const PetscInt pi  = p < 4 ? fv[p] : p;
1829     const PetscInt pin = p < 3 ? fv[(p+1)%numCorners] : (p+1)%numCorners;
1830     /* Need to do this copy to get types right */
1831     for (d = 0; d < tdim; ++d) {
1832       ctmp[d]      = PetscRealPart(coords[pi*tdim+d]);
1833       ctmp[tdim+d] = PetscRealPart(coords[pin*tdim+d]);
1834     }
1835     Volume_Triangle_Origin_Internal(&vtmp, ctmp);
1836     vsum += vtmp;
1837     for (d = 0; d < tdim; ++d) {
1838       csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp;
1839     }
1840   }
1841   for (d = 0; d < tdim; ++d) {
1842     csum[d] /= (tdim+1)*vsum;
1843   }
1844   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1845   if (vol) *vol = PetscAbsReal(vsum);
1846   if (centroid) {
1847     if (dim > 2) {
1848       for (d = 0; d < dim; ++d) {
1849         centroid[d] = v0[d];
1850         for (e = 0; e < dim; ++e) {
1851           centroid[d] += R[d*dim+e]*csum[e];
1852         }
1853       }
1854     } else for (d = 0; d < dim; ++d) centroid[d] = csum[d];
1855   }
1856   PetscFunctionReturn(0);
1857 }
1858 
1859 /* Centroid_i = (\sum_n V_n Cn_i ) / V */
1860 static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1861 {
1862   DMLabel         depth;
1863   PetscSection    coordSection;
1864   Vec             coordinates;
1865   PetscScalar    *coords = NULL;
1866   PetscReal       vsum = 0.0, vtmp, coordsTmp[3*3];
1867   const PetscInt *faces, *facesO;
1868   PetscBool       isHybrid = PETSC_FALSE;
1869   PetscInt        pEndInterior[4], cdepth, numFaces, f, coordSize, numCorners, p, d;
1870   PetscErrorCode  ierr;
1871 
1872   PetscFunctionBegin;
1873   if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No support for dim %D > 3",dim);
1874   /* Must check for hybrid cells because prisms have a different orientation scheme */
1875   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
1876   ierr = DMLabelGetValue(depth, cell, &cdepth);CHKERRQ(ierr);
1877   ierr = DMPlexGetHybridBounds(dm, &pEndInterior[3], &pEndInterior[2], &pEndInterior[1], &pEndInterior[0]);CHKERRQ(ierr);
1878   if ((pEndInterior[cdepth]) >=0 && (cell >= pEndInterior[cdepth])) isHybrid = PETSC_TRUE;
1879 
1880   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1881   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1882 
1883   if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0;
1884   ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr);
1885   ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr);
1886   ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr);
1887   for (f = 0; f < numFaces; ++f) {
1888     PetscBool flip = isHybrid && f == 0 ? PETSC_TRUE : PETSC_FALSE; /* The first hybrid face is reversed */
1889 
1890     ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr);
1891     numCorners = coordSize/dim;
1892     switch (numCorners) {
1893     case 3:
1894       for (d = 0; d < dim; ++d) {
1895         coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]);
1896         coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]);
1897         coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]);
1898       }
1899       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
1900       if (facesO[f] < 0 || flip) vtmp = -vtmp;
1901       vsum += vtmp;
1902       if (centroid) {           /* Centroid of OABC = (a+b+c)/4 */
1903         for (d = 0; d < dim; ++d) {
1904           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
1905         }
1906       }
1907       break;
1908     case 4:
1909     {
1910       PetscInt fv[4] = {0, 1, 2, 3};
1911 
1912       /* Side faces for hybrid cells are are stored as tensor products */
1913       if (isHybrid && f > 1) {fv[2] = 3; fv[3] = 2;}
1914       /* DO FOR PYRAMID */
1915       /* First tet */
1916       for (d = 0; d < dim; ++d) {
1917         coordsTmp[0*dim+d] = PetscRealPart(coords[fv[0]*dim+d]);
1918         coordsTmp[1*dim+d] = PetscRealPart(coords[fv[1]*dim+d]);
1919         coordsTmp[2*dim+d] = PetscRealPart(coords[fv[3]*dim+d]);
1920       }
1921       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
1922       if (facesO[f] < 0 || flip) vtmp = -vtmp;
1923       vsum += vtmp;
1924       if (centroid) {
1925         for (d = 0; d < dim; ++d) {
1926           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
1927         }
1928       }
1929       /* Second tet */
1930       for (d = 0; d < dim; ++d) {
1931         coordsTmp[0*dim+d] = PetscRealPart(coords[fv[1]*dim+d]);
1932         coordsTmp[1*dim+d] = PetscRealPart(coords[fv[2]*dim+d]);
1933         coordsTmp[2*dim+d] = PetscRealPart(coords[fv[3]*dim+d]);
1934       }
1935       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
1936       if (facesO[f] < 0 || flip) vtmp = -vtmp;
1937       vsum += vtmp;
1938       if (centroid) {
1939         for (d = 0; d < dim; ++d) {
1940           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
1941         }
1942       }
1943       break;
1944     }
1945     default:
1946       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %D vertices", numCorners);
1947     }
1948     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr);
1949   }
1950   if (vol)     *vol = PetscAbsReal(vsum);
1951   if (normal)   for (d = 0; d < dim; ++d) normal[d]    = 0.0;
1952   if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4);
1953   PetscFunctionReturn(0);
1954 }
1955 
1956 /*@C
1957   DMPlexComputeCellGeometryFVM - Compute the volume for a given cell
1958 
1959   Collective on DM
1960 
1961   Input Arguments:
1962 + dm   - the DM
1963 - cell - the cell
1964 
1965   Output Arguments:
1966 + volume   - the cell volume
1967 . centroid - the cell centroid
1968 - normal - the cell normal, if appropriate
1969 
1970   Level: advanced
1971 
1972   Fortran Notes:
1973   Since it returns arrays, this routine is only available in Fortran 90, and you must
1974   include petsc.h90 in your code.
1975 
1976 .seealso: DMGetCoordinateSection(), DMGetCoordinates()
1977 @*/
1978 PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1979 {
1980   PetscInt       depth, dim;
1981   PetscErrorCode ierr;
1982 
1983   PetscFunctionBegin;
1984   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1985   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1986   if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated");
1987   /* We need to keep a pointer to the depth label */
1988   ierr = DMGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr);
1989   /* Cone size is now the number of faces */
1990   switch (depth) {
1991   case 1:
1992     ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
1993     break;
1994   case 2:
1995     ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
1996     break;
1997   case 3:
1998     ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
1999     break;
2000   default:
2001     SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D (depth %D) for element geometry computation", dim, depth);
2002   }
2003   PetscFunctionReturn(0);
2004 }
2005 
2006 /*@
2007   DMPlexComputeGeometryFEM - Precompute cell geometry for the entire mesh
2008 
2009   Collective on dm
2010 
2011   Input Parameter:
2012 . dm - The DMPlex
2013 
2014   Output Parameter:
2015 . cellgeom - A vector with the cell geometry data for each cell
2016 
2017   Level: beginner
2018 
2019 @*/
2020 PetscErrorCode DMPlexComputeGeometryFEM(DM dm, Vec *cellgeom)
2021 {
2022   DM             dmCell;
2023   Vec            coordinates;
2024   PetscSection   coordSection, sectionCell;
2025   PetscScalar   *cgeom;
2026   PetscInt       cStart, cEnd, cMax, c;
2027   PetscErrorCode ierr;
2028 
2029   PetscFunctionBegin;
2030   ierr = DMClone(dm, &dmCell);CHKERRQ(ierr);
2031   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2032   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2033   ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr);
2034   ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr);
2035   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionCell);CHKERRQ(ierr);
2036   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2037   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
2038   cEnd = cMax < 0 ? cEnd : cMax;
2039   ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr);
2040   /* TODO This needs to be multiplied by Nq for non-affine */
2041   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFEGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);}
2042   ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr);
2043   ierr = DMSetSection(dmCell, sectionCell);CHKERRQ(ierr);
2044   ierr = PetscSectionDestroy(&sectionCell);CHKERRQ(ierr);
2045   ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr);
2046   ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr);
2047   for (c = cStart; c < cEnd; ++c) {
2048     PetscFEGeom *cg;
2049 
2050     ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
2051     ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr);
2052     ierr = DMPlexComputeCellGeometryFEM(dmCell, c, NULL, cg->v, cg->J, cg->invJ, cg->detJ);CHKERRQ(ierr);
2053     if (*cg->detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cg->detJ, c);
2054   }
2055   PetscFunctionReturn(0);
2056 }
2057 
2058 /*@
2059   DMPlexComputeGeometryFVM - Computes the cell and face geometry for a finite volume method
2060 
2061   Input Parameter:
2062 . dm - The DM
2063 
2064   Output Parameters:
2065 + cellgeom - A Vec of PetscFVCellGeom data
2066 . facegeom - A Vec of PetscFVFaceGeom data
2067 
2068   Level: developer
2069 
2070 .seealso: PetscFVFaceGeom, PetscFVCellGeom, DMPlexComputeGeometryFEM()
2071 @*/
2072 PetscErrorCode DMPlexComputeGeometryFVM(DM dm, Vec *cellgeom, Vec *facegeom)
2073 {
2074   DM             dmFace, dmCell;
2075   DMLabel        ghostLabel;
2076   PetscSection   sectionFace, sectionCell;
2077   PetscSection   coordSection;
2078   Vec            coordinates;
2079   PetscScalar   *fgeom, *cgeom;
2080   PetscReal      minradius, gminradius;
2081   PetscInt       dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f;
2082   PetscErrorCode ierr;
2083 
2084   PetscFunctionBegin;
2085   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2086   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2087   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2088   /* Make cell centroids and volumes */
2089   ierr = DMClone(dm, &dmCell);CHKERRQ(ierr);
2090   ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr);
2091   ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr);
2092   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionCell);CHKERRQ(ierr);
2093   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2094   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
2095   ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr);
2096   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVCellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);}
2097   ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr);
2098   ierr = DMSetSection(dmCell, sectionCell);CHKERRQ(ierr);
2099   ierr = PetscSectionDestroy(&sectionCell);CHKERRQ(ierr);
2100   ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr);
2101   if (cEndInterior < 0) {
2102     cEndInterior = cEnd;
2103   }
2104   ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr);
2105   for (c = cStart; c < cEndInterior; ++c) {
2106     PetscFVCellGeom *cg;
2107 
2108     ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
2109     ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr);
2110     ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr);
2111   }
2112   /* Compute face normals and minimum cell radius */
2113   ierr = DMClone(dm, &dmFace);CHKERRQ(ierr);
2114   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionFace);CHKERRQ(ierr);
2115   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
2116   ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr);
2117   for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVFaceGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);}
2118   ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr);
2119   ierr = DMSetSection(dmFace, sectionFace);CHKERRQ(ierr);
2120   ierr = PetscSectionDestroy(&sectionFace);CHKERRQ(ierr);
2121   ierr = DMCreateLocalVector(dmFace, facegeom);CHKERRQ(ierr);
2122   ierr = VecGetArray(*facegeom, &fgeom);CHKERRQ(ierr);
2123   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
2124   minradius = PETSC_MAX_REAL;
2125   for (f = fStart; f < fEnd; ++f) {
2126     PetscFVFaceGeom *fg;
2127     PetscReal        area;
2128     PetscInt         ghost = -1, d, numChildren;
2129 
2130     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);}
2131     ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr);
2132     if (ghost >= 0 || numChildren) continue;
2133     ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr);
2134     ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr);
2135     for (d = 0; d < dim; ++d) fg->normal[d] *= area;
2136     /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */
2137     {
2138       PetscFVCellGeom *cL, *cR;
2139       PetscInt         ncells;
2140       const PetscInt  *cells;
2141       PetscReal       *lcentroid, *rcentroid;
2142       PetscReal        l[3], r[3], v[3];
2143 
2144       ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr);
2145       ierr = DMPlexGetSupportSize(dm, f, &ncells);CHKERRQ(ierr);
2146       ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr);
2147       lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid;
2148       if (ncells > 1) {
2149         ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr);
2150         rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid;
2151       }
2152       else {
2153         rcentroid = fg->centroid;
2154       }
2155       ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, lcentroid, l);CHKERRQ(ierr);
2156       ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, rcentroid, r);CHKERRQ(ierr);
2157       DMPlex_WaxpyD_Internal(dim, -1, l, r, v);
2158       if (DMPlex_DotRealD_Internal(dim, fg->normal, v) < 0) {
2159         for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d];
2160       }
2161       if (DMPlex_DotRealD_Internal(dim, fg->normal, v) <= 0) {
2162         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]);
2163         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]);
2164         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f);
2165       }
2166       if (cells[0] < cEndInterior) {
2167         DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cL->centroid, v);
2168         minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v));
2169       }
2170       if (ncells > 1 && cells[1] < cEndInterior) {
2171         DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cR->centroid, v);
2172         minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v));
2173       }
2174     }
2175   }
2176   ierr = MPIU_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
2177   ierr = DMPlexSetMinRadius(dm, gminradius);CHKERRQ(ierr);
2178   /* Compute centroids of ghost cells */
2179   for (c = cEndInterior; c < cEnd; ++c) {
2180     PetscFVFaceGeom *fg;
2181     const PetscInt  *cone,    *support;
2182     PetscInt         coneSize, supportSize, s;
2183 
2184     ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr);
2185     if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize);
2186     ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr);
2187     ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr);
2188     if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 2", cone[0], supportSize);
2189     ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr);
2190     ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr);
2191     for (s = 0; s < 2; ++s) {
2192       /* Reflect ghost centroid across plane of face */
2193       if (support[s] == c) {
2194         PetscFVCellGeom       *ci;
2195         PetscFVCellGeom       *cg;
2196         PetscReal              c2f[3], a;
2197 
2198         ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr);
2199         DMPlex_WaxpyD_Internal(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */
2200         a    = DMPlex_DotRealD_Internal(dim, c2f, fg->normal)/DMPlex_DotRealD_Internal(dim, fg->normal, fg->normal);
2201         ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr);
2202         DMPlex_WaxpyD_Internal(dim, 2*a, fg->normal, ci->centroid, cg->centroid);
2203         cg->volume = ci->volume;
2204       }
2205     }
2206   }
2207   ierr = VecRestoreArray(*facegeom, &fgeom);CHKERRQ(ierr);
2208   ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr);
2209   ierr = DMDestroy(&dmCell);CHKERRQ(ierr);
2210   ierr = DMDestroy(&dmFace);CHKERRQ(ierr);
2211   PetscFunctionReturn(0);
2212 }
2213 
2214 /*@C
2215   DMPlexGetMinRadius - Returns the minimum distance from any cell centroid to a face
2216 
2217   Not collective
2218 
2219   Input Argument:
2220 . dm - the DM
2221 
2222   Output Argument:
2223 . minradius - the minium cell radius
2224 
2225   Level: developer
2226 
2227 .seealso: DMGetCoordinates()
2228 @*/
2229 PetscErrorCode DMPlexGetMinRadius(DM dm, PetscReal *minradius)
2230 {
2231   PetscFunctionBegin;
2232   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2233   PetscValidPointer(minradius,2);
2234   *minradius = ((DM_Plex*) dm->data)->minradius;
2235   PetscFunctionReturn(0);
2236 }
2237 
2238 /*@C
2239   DMPlexSetMinRadius - Sets the minimum distance from the cell centroid to a face
2240 
2241   Logically collective
2242 
2243   Input Arguments:
2244 + dm - the DM
2245 - minradius - the minium cell radius
2246 
2247   Level: developer
2248 
2249 .seealso: DMSetCoordinates()
2250 @*/
2251 PetscErrorCode DMPlexSetMinRadius(DM dm, PetscReal minradius)
2252 {
2253   PetscFunctionBegin;
2254   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2255   ((DM_Plex*) dm->data)->minradius = minradius;
2256   PetscFunctionReturn(0);
2257 }
2258 
2259 static PetscErrorCode BuildGradientReconstruction_Internal(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
2260 {
2261   DMLabel        ghostLabel;
2262   PetscScalar   *dx, *grad, **gref;
2263   PetscInt       dim, cStart, cEnd, c, cEndInterior, maxNumFaces;
2264   PetscErrorCode ierr;
2265 
2266   PetscFunctionBegin;
2267   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2268   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2269   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
2270   ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr);
2271   ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr);
2272   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
2273   ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr);
2274   for (c = cStart; c < cEndInterior; c++) {
2275     const PetscInt        *faces;
2276     PetscInt               numFaces, usedFaces, f, d;
2277     PetscFVCellGeom        *cg;
2278     PetscBool              boundary;
2279     PetscInt               ghost;
2280 
2281     ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
2282     ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr);
2283     ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr);
2284     if (numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces);
2285     for (f = 0, usedFaces = 0; f < numFaces; ++f) {
2286       PetscFVCellGeom       *cg1;
2287       PetscFVFaceGeom       *fg;
2288       const PetscInt        *fcells;
2289       PetscInt               ncell, side;
2290 
2291       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
2292       ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
2293       if ((ghost >= 0) || boundary) continue;
2294       ierr  = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr);
2295       side  = (c != fcells[0]); /* c is on left=0 or right=1 of face */
2296       ncell = fcells[!side];    /* the neighbor */
2297       ierr  = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr);
2298       ierr  = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr);
2299       for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d];
2300       gref[usedFaces++] = fg->grad[side];  /* Gradient reconstruction term will go here */
2301     }
2302     if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?");
2303     ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr);
2304     for (f = 0, usedFaces = 0; f < numFaces; ++f) {
2305       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
2306       ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
2307       if ((ghost >= 0) || boundary) continue;
2308       for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d];
2309       ++usedFaces;
2310     }
2311   }
2312   ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr);
2313   PetscFunctionReturn(0);
2314 }
2315 
2316 static PetscErrorCode BuildGradientReconstruction_Internal_Tree(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
2317 {
2318   DMLabel        ghostLabel;
2319   PetscScalar   *dx, *grad, **gref;
2320   PetscInt       dim, cStart, cEnd, c, cEndInterior, fStart, fEnd, f, nStart, nEnd, maxNumFaces = 0;
2321   PetscSection   neighSec;
2322   PetscInt     (*neighbors)[2];
2323   PetscInt      *counter;
2324   PetscErrorCode ierr;
2325 
2326   PetscFunctionBegin;
2327   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2328   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2329   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
2330   if (cEndInterior < 0) {
2331     cEndInterior = cEnd;
2332   }
2333   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&neighSec);CHKERRQ(ierr);
2334   ierr = PetscSectionSetChart(neighSec,cStart,cEndInterior);CHKERRQ(ierr);
2335   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
2336   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
2337   for (f = fStart; f < fEnd; f++) {
2338     const PetscInt        *fcells;
2339     PetscBool              boundary;
2340     PetscInt               ghost = -1;
2341     PetscInt               numChildren, numCells, c;
2342 
2343     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);}
2344     ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr);
2345     ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr);
2346     if ((ghost >= 0) || boundary || numChildren) continue;
2347     ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr);
2348     if (numCells == 2) {
2349       ierr = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr);
2350       for (c = 0; c < 2; c++) {
2351         PetscInt cell = fcells[c];
2352 
2353         if (cell >= cStart && cell < cEndInterior) {
2354           ierr = PetscSectionAddDof(neighSec,cell,1);CHKERRQ(ierr);
2355         }
2356       }
2357     }
2358   }
2359   ierr = PetscSectionSetUp(neighSec);CHKERRQ(ierr);
2360   ierr = PetscSectionGetMaxDof(neighSec,&maxNumFaces);CHKERRQ(ierr);
2361   ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr);
2362   nStart = 0;
2363   ierr = PetscSectionGetStorageSize(neighSec,&nEnd);CHKERRQ(ierr);
2364   ierr = PetscMalloc1((nEnd-nStart),&neighbors);CHKERRQ(ierr);
2365   ierr = PetscCalloc1((cEndInterior-cStart),&counter);CHKERRQ(ierr);
2366   for (f = fStart; f < fEnd; f++) {
2367     const PetscInt        *fcells;
2368     PetscBool              boundary;
2369     PetscInt               ghost = -1;
2370     PetscInt               numChildren, numCells, c;
2371 
2372     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);}
2373     ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr);
2374     ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr);
2375     if ((ghost >= 0) || boundary || numChildren) continue;
2376     ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr);
2377     if (numCells == 2) {
2378       ierr  = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr);
2379       for (c = 0; c < 2; c++) {
2380         PetscInt cell = fcells[c], off;
2381 
2382         if (cell >= cStart && cell < cEndInterior) {
2383           ierr = PetscSectionGetOffset(neighSec,cell,&off);CHKERRQ(ierr);
2384           off += counter[cell - cStart]++;
2385           neighbors[off][0] = f;
2386           neighbors[off][1] = fcells[1 - c];
2387         }
2388       }
2389     }
2390   }
2391   ierr = PetscFree(counter);CHKERRQ(ierr);
2392   ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr);
2393   for (c = cStart; c < cEndInterior; c++) {
2394     PetscInt               numFaces, f, d, off, ghost = -1;
2395     PetscFVCellGeom        *cg;
2396 
2397     ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
2398     ierr = PetscSectionGetDof(neighSec, c, &numFaces);CHKERRQ(ierr);
2399     ierr = PetscSectionGetOffset(neighSec, c, &off);CHKERRQ(ierr);
2400     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, c, &ghost);CHKERRQ(ierr);}
2401     if (ghost < 0 && numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces);
2402     for (f = 0; f < numFaces; ++f) {
2403       PetscFVCellGeom       *cg1;
2404       PetscFVFaceGeom       *fg;
2405       const PetscInt        *fcells;
2406       PetscInt               ncell, side, nface;
2407 
2408       nface = neighbors[off + f][0];
2409       ncell = neighbors[off + f][1];
2410       ierr  = DMPlexGetSupport(dm,nface,&fcells);CHKERRQ(ierr);
2411       side  = (c != fcells[0]);
2412       ierr  = DMPlexPointLocalRef(dmFace, nface, fgeom, &fg);CHKERRQ(ierr);
2413       ierr  = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr);
2414       for (d = 0; d < dim; ++d) dx[f*dim+d] = cg1->centroid[d] - cg->centroid[d];
2415       gref[f] = fg->grad[side];  /* Gradient reconstruction term will go here */
2416     }
2417     ierr = PetscFVComputeGradient(fvm, numFaces, dx, grad);CHKERRQ(ierr);
2418     for (f = 0; f < numFaces; ++f) {
2419       for (d = 0; d < dim; ++d) gref[f][d] = grad[f*dim+d];
2420     }
2421   }
2422   ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr);
2423   ierr = PetscSectionDestroy(&neighSec);CHKERRQ(ierr);
2424   ierr = PetscFree(neighbors);CHKERRQ(ierr);
2425   PetscFunctionReturn(0);
2426 }
2427 
2428 /*@
2429   DMPlexComputeGradientFVM - Compute geometric factors for gradient reconstruction, which are stored in the geometry data, and compute layout for gradient data
2430 
2431   Collective on DM
2432 
2433   Input Arguments:
2434 + dm  - The DM
2435 . fvm - The PetscFV
2436 . faceGeometry - The face geometry from DMPlexComputeFaceGeometryFVM()
2437 - cellGeometry - The face geometry from DMPlexComputeCellGeometryFVM()
2438 
2439   Output Parameters:
2440 + faceGeometry - The geometric factors for gradient calculation are inserted
2441 - dmGrad - The DM describing the layout of gradient data
2442 
2443   Level: developer
2444 
2445 .seealso: DMPlexGetFaceGeometryFVM(), DMPlexGetCellGeometryFVM()
2446 @*/
2447 PetscErrorCode DMPlexComputeGradientFVM(DM dm, PetscFV fvm, Vec faceGeometry, Vec cellGeometry, DM *dmGrad)
2448 {
2449   DM             dmFace, dmCell;
2450   PetscScalar   *fgeom, *cgeom;
2451   PetscSection   sectionGrad, parentSection;
2452   PetscInt       dim, pdim, cStart, cEnd, cEndInterior, c;
2453   PetscErrorCode ierr;
2454 
2455   PetscFunctionBegin;
2456   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2457   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
2458   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2459   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
2460   /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */
2461   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
2462   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
2463   ierr = VecGetArray(faceGeometry, &fgeom);CHKERRQ(ierr);
2464   ierr = VecGetArray(cellGeometry, &cgeom);CHKERRQ(ierr);
2465   ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
2466   if (!parentSection) {
2467     ierr = BuildGradientReconstruction_Internal(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr);
2468   } else {
2469     ierr = BuildGradientReconstruction_Internal_Tree(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr);
2470   }
2471   ierr = VecRestoreArray(faceGeometry, &fgeom);CHKERRQ(ierr);
2472   ierr = VecRestoreArray(cellGeometry, &cgeom);CHKERRQ(ierr);
2473   /* Create storage for gradients */
2474   ierr = DMClone(dm, dmGrad);CHKERRQ(ierr);
2475   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionGrad);CHKERRQ(ierr);
2476   ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr);
2477   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);}
2478   ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr);
2479   ierr = DMSetSection(*dmGrad, sectionGrad);CHKERRQ(ierr);
2480   ierr = PetscSectionDestroy(&sectionGrad);CHKERRQ(ierr);
2481   PetscFunctionReturn(0);
2482 }
2483 
2484 /*@
2485   DMPlexGetDataFVM - Retrieve precomputed cell geometry
2486 
2487   Collective on DM
2488 
2489   Input Arguments:
2490 + dm  - The DM
2491 - fvm - The PetscFV
2492 
2493   Output Parameters:
2494 + cellGeometry - The cell geometry
2495 . faceGeometry - The face geometry
2496 - dmGrad       - The gradient matrices
2497 
2498   Level: developer
2499 
2500 .seealso: DMPlexComputeGeometryFVM()
2501 @*/
2502 PetscErrorCode DMPlexGetDataFVM(DM dm, PetscFV fv, Vec *cellgeom, Vec *facegeom, DM *gradDM)
2503 {
2504   PetscObject    cellgeomobj, facegeomobj;
2505   PetscErrorCode ierr;
2506 
2507   PetscFunctionBegin;
2508   ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);CHKERRQ(ierr);
2509   if (!cellgeomobj) {
2510     Vec cellgeomInt, facegeomInt;
2511 
2512     ierr = DMPlexComputeGeometryFVM(dm, &cellgeomInt, &facegeomInt);CHKERRQ(ierr);
2513     ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_cellgeom_fvm",(PetscObject)cellgeomInt);CHKERRQ(ierr);
2514     ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_facegeom_fvm",(PetscObject)facegeomInt);CHKERRQ(ierr);
2515     ierr = VecDestroy(&cellgeomInt);CHKERRQ(ierr);
2516     ierr = VecDestroy(&facegeomInt);CHKERRQ(ierr);
2517     ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);CHKERRQ(ierr);
2518   }
2519   ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_facegeom_fvm", &facegeomobj);CHKERRQ(ierr);
2520   if (cellgeom) *cellgeom = (Vec) cellgeomobj;
2521   if (facegeom) *facegeom = (Vec) facegeomobj;
2522   if (gradDM) {
2523     PetscObject gradobj;
2524     PetscBool   computeGradients;
2525 
2526     ierr = PetscFVGetComputeGradients(fv,&computeGradients);CHKERRQ(ierr);
2527     if (!computeGradients) {
2528       *gradDM = NULL;
2529       PetscFunctionReturn(0);
2530     }
2531     ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);CHKERRQ(ierr);
2532     if (!gradobj) {
2533       DM dmGradInt;
2534 
2535       ierr = DMPlexComputeGradientFVM(dm,fv,(Vec) facegeomobj,(Vec) cellgeomobj,&dmGradInt);CHKERRQ(ierr);
2536       ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_dmgrad_fvm", (PetscObject)dmGradInt);CHKERRQ(ierr);
2537       ierr = DMDestroy(&dmGradInt);CHKERRQ(ierr);
2538       ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);CHKERRQ(ierr);
2539     }
2540     *gradDM = (DM) gradobj;
2541   }
2542   PetscFunctionReturn(0);
2543 }
2544 
2545 static PetscErrorCode DMPlexCoordinatesToReference_NewtonUpdate(PetscInt dimC, PetscInt dimR, PetscScalar *J, PetscScalar *invJ, PetscScalar *work,  PetscReal *resNeg, PetscReal *guess)
2546 {
2547   PetscInt l, m;
2548 
2549   PetscFunctionBeginHot;
2550   if (dimC == dimR && dimR <= 3) {
2551     /* invert Jacobian, multiply */
2552     PetscScalar det, idet;
2553 
2554     switch (dimR) {
2555     case 1:
2556       invJ[0] = 1./ J[0];
2557       break;
2558     case 2:
2559       det = J[0] * J[3] - J[1] * J[2];
2560       idet = 1./det;
2561       invJ[0] =  J[3] * idet;
2562       invJ[1] = -J[1] * idet;
2563       invJ[2] = -J[2] * idet;
2564       invJ[3] =  J[0] * idet;
2565       break;
2566     case 3:
2567       {
2568         invJ[0] = J[4] * J[8] - J[5] * J[7];
2569         invJ[1] = J[2] * J[7] - J[1] * J[8];
2570         invJ[2] = J[1] * J[5] - J[2] * J[4];
2571         det = invJ[0] * J[0] + invJ[1] * J[3] + invJ[2] * J[6];
2572         idet = 1./det;
2573         invJ[0] *= idet;
2574         invJ[1] *= idet;
2575         invJ[2] *= idet;
2576         invJ[3]  = idet * (J[5] * J[6] - J[3] * J[8]);
2577         invJ[4]  = idet * (J[0] * J[8] - J[2] * J[6]);
2578         invJ[5]  = idet * (J[2] * J[3] - J[0] * J[5]);
2579         invJ[6]  = idet * (J[3] * J[7] - J[4] * J[6]);
2580         invJ[7]  = idet * (J[1] * J[6] - J[0] * J[7]);
2581         invJ[8]  = idet * (J[0] * J[4] - J[1] * J[3]);
2582       }
2583       break;
2584     }
2585     for (l = 0; l < dimR; l++) {
2586       for (m = 0; m < dimC; m++) {
2587         guess[l] += PetscRealPart(invJ[l * dimC + m]) * resNeg[m];
2588       }
2589     }
2590   } else {
2591 #if defined(PETSC_USE_COMPLEX)
2592     char transpose = 'C';
2593 #else
2594     char transpose = 'T';
2595 #endif
2596     PetscBLASInt m = dimR;
2597     PetscBLASInt n = dimC;
2598     PetscBLASInt one = 1;
2599     PetscBLASInt worksize = dimR * dimC, info;
2600 
2601     for (l = 0; l < dimC; l++) {invJ[l] = resNeg[l];}
2602 
2603     PetscStackCallBLAS("LAPACKgels",LAPACKgels_(&transpose,&m,&n,&one,J,&m,invJ,&n,work,&worksize, &info));
2604     if (info != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELS");
2605 
2606     for (l = 0; l < dimR; l++) {guess[l] += PetscRealPart(invJ[l]);}
2607   }
2608   PetscFunctionReturn(0);
2609 }
2610 
2611 static PetscErrorCode DMPlexCoordinatesToReference_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[], Vec coords, PetscInt dimC, PetscInt dimR)
2612 {
2613   PetscInt       coordSize, i, j, k, l, m, maxIts = 7, numV = (1 << dimR);
2614   PetscScalar    *coordsScalar = NULL;
2615   PetscReal      *cellData, *cellCoords, *cellCoeffs, *extJ, *resNeg;
2616   PetscScalar    *J, *invJ, *work;
2617   PetscErrorCode ierr;
2618 
2619   PetscFunctionBegin;
2620   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2621   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
2622   if (coordSize < dimC * numV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expecting at least %D coordinates, got %D",dimC * (1 << dimR), coordSize);
2623   ierr = DMGetWorkArray(dm, 2 * coordSize + dimR + dimC, MPIU_REAL, &cellData);CHKERRQ(ierr);
2624   ierr = DMGetWorkArray(dm, 3 * dimR * dimC, MPIU_SCALAR, &J);CHKERRQ(ierr);
2625   cellCoords = &cellData[0];
2626   cellCoeffs = &cellData[coordSize];
2627   extJ       = &cellData[2 * coordSize];
2628   resNeg     = &cellData[2 * coordSize + dimR];
2629   invJ       = &J[dimR * dimC];
2630   work       = &J[2 * dimR * dimC];
2631   if (dimR == 2) {
2632     const PetscInt zToPlex[4] = {0, 1, 3, 2};
2633 
2634     for (i = 0; i < 4; i++) {
2635       PetscInt plexI = zToPlex[i];
2636 
2637       for (j = 0; j < dimC; j++) {
2638         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
2639       }
2640     }
2641   } else if (dimR == 3) {
2642     const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6};
2643 
2644     for (i = 0; i < 8; i++) {
2645       PetscInt plexI = zToPlex[i];
2646 
2647       for (j = 0; j < dimC; j++) {
2648         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
2649       }
2650     }
2651   } else {
2652     for (i = 0; i < coordSize; i++) {cellCoords[i] = PetscRealPart(coordsScalar[i]);}
2653   }
2654   /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */
2655   for (i = 0; i < dimR; i++) {
2656     PetscReal *swap;
2657 
2658     for (j = 0; j < (numV / 2); j++) {
2659       for (k = 0; k < dimC; k++) {
2660         cellCoeffs[dimC * j + k]                = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]);
2661         cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]);
2662       }
2663     }
2664 
2665     if (i < dimR - 1) {
2666       swap = cellCoeffs;
2667       cellCoeffs = cellCoords;
2668       cellCoords = swap;
2669     }
2670   }
2671   ierr = PetscMemzero(refCoords,numPoints * dimR * sizeof (PetscReal));CHKERRQ(ierr);
2672   for (j = 0; j < numPoints; j++) {
2673     for (i = 0; i < maxIts; i++) {
2674       PetscReal *guess = &refCoords[dimR * j];
2675 
2676       /* compute -residual and Jacobian */
2677       for (k = 0; k < dimC; k++) {resNeg[k] = realCoords[dimC * j + k];}
2678       for (k = 0; k < dimC * dimR; k++) {J[k] = 0.;}
2679       for (k = 0; k < numV; k++) {
2680         PetscReal extCoord = 1.;
2681         for (l = 0; l < dimR; l++) {
2682           PetscReal coord = guess[l];
2683           PetscInt  dep   = (k & (1 << l)) >> l;
2684 
2685           extCoord *= dep * coord + !dep;
2686           extJ[l] = dep;
2687 
2688           for (m = 0; m < dimR; m++) {
2689             PetscReal coord = guess[m];
2690             PetscInt  dep   = ((k & (1 << m)) >> m) && (m != l);
2691             PetscReal mult  = dep * coord + !dep;
2692 
2693             extJ[l] *= mult;
2694           }
2695         }
2696         for (l = 0; l < dimC; l++) {
2697           PetscReal coeff = cellCoeffs[dimC * k + l];
2698 
2699           resNeg[l] -= coeff * extCoord;
2700           for (m = 0; m < dimR; m++) {
2701             J[dimR * l + m] += coeff * extJ[m];
2702           }
2703         }
2704       }
2705 #if 0 && defined(PETSC_USE_DEBUG)
2706       {
2707         PetscReal maxAbs = 0.;
2708 
2709         for (l = 0; l < dimC; l++) {
2710           maxAbs = PetscMax(maxAbs,PetscAbsReal(resNeg[l]));
2711         }
2712         ierr = PetscInfo4(dm,"cell %D, point %D, iter %D: res %g\n",cell,j,i,maxAbs);CHKERRQ(ierr);
2713       }
2714 #endif
2715 
2716       ierr = DMPlexCoordinatesToReference_NewtonUpdate(dimC,dimR,J,invJ,work,resNeg,guess);CHKERRQ(ierr);
2717     }
2718   }
2719   ierr = DMRestoreWorkArray(dm, 3 * dimR * dimC, MPIU_SCALAR, &J);CHKERRQ(ierr);
2720   ierr = DMRestoreWorkArray(dm, 2 * coordSize + dimR + dimC, MPIU_REAL, &cellData);CHKERRQ(ierr);
2721   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
2722   PetscFunctionReturn(0);
2723 }
2724 
2725 static PetscErrorCode DMPlexReferenceToCoordinates_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[], Vec coords, PetscInt dimC, PetscInt dimR)
2726 {
2727   PetscInt       coordSize, i, j, k, l, numV = (1 << dimR);
2728   PetscScalar    *coordsScalar = NULL;
2729   PetscReal      *cellData, *cellCoords, *cellCoeffs;
2730   PetscErrorCode ierr;
2731 
2732   PetscFunctionBegin;
2733   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2734   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
2735   if (coordSize < dimC * numV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expecting at least %D coordinates, got %D",dimC * (1 << dimR), coordSize);
2736   ierr = DMGetWorkArray(dm, 2 * coordSize, MPIU_REAL, &cellData);CHKERRQ(ierr);
2737   cellCoords = &cellData[0];
2738   cellCoeffs = &cellData[coordSize];
2739   if (dimR == 2) {
2740     const PetscInt zToPlex[4] = {0, 1, 3, 2};
2741 
2742     for (i = 0; i < 4; i++) {
2743       PetscInt plexI = zToPlex[i];
2744 
2745       for (j = 0; j < dimC; j++) {
2746         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
2747       }
2748     }
2749   } else if (dimR == 3) {
2750     const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6};
2751 
2752     for (i = 0; i < 8; i++) {
2753       PetscInt plexI = zToPlex[i];
2754 
2755       for (j = 0; j < dimC; j++) {
2756         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
2757       }
2758     }
2759   } else {
2760     for (i = 0; i < coordSize; i++) {cellCoords[i] = PetscRealPart(coordsScalar[i]);}
2761   }
2762   /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */
2763   for (i = 0; i < dimR; i++) {
2764     PetscReal *swap;
2765 
2766     for (j = 0; j < (numV / 2); j++) {
2767       for (k = 0; k < dimC; k++) {
2768         cellCoeffs[dimC * j + k]                = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]);
2769         cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]);
2770       }
2771     }
2772 
2773     if (i < dimR - 1) {
2774       swap = cellCoeffs;
2775       cellCoeffs = cellCoords;
2776       cellCoords = swap;
2777     }
2778   }
2779   ierr = PetscMemzero(realCoords,numPoints * dimC * sizeof (PetscReal));CHKERRQ(ierr);
2780   for (j = 0; j < numPoints; j++) {
2781     const PetscReal *guess  = &refCoords[dimR * j];
2782     PetscReal       *mapped = &realCoords[dimC * j];
2783 
2784     for (k = 0; k < numV; k++) {
2785       PetscReal extCoord = 1.;
2786       for (l = 0; l < dimR; l++) {
2787         PetscReal coord = guess[l];
2788         PetscInt  dep   = (k & (1 << l)) >> l;
2789 
2790         extCoord *= dep * coord + !dep;
2791       }
2792       for (l = 0; l < dimC; l++) {
2793         PetscReal coeff = cellCoeffs[dimC * k + l];
2794 
2795         mapped[l] += coeff * extCoord;
2796       }
2797     }
2798   }
2799   ierr = DMRestoreWorkArray(dm, 2 * coordSize, MPIU_REAL, &cellData);CHKERRQ(ierr);
2800   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
2801   PetscFunctionReturn(0);
2802 }
2803 
2804 /* TODO: TOBY please fix this for Nc > 1 */
2805 static PetscErrorCode DMPlexCoordinatesToReference_FE(DM dm, PetscFE fe, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[], Vec coords, PetscInt Nc, PetscInt dimR)
2806 {
2807   PetscInt       numComp, pdim, i, j, k, l, m, maxIter = 7, coordSize;
2808   PetscScalar    *nodes = NULL;
2809   PetscReal      *invV, *modes;
2810   PetscReal      *B, *D, *resNeg;
2811   PetscScalar    *J, *invJ, *work;
2812   PetscErrorCode ierr;
2813 
2814   PetscFunctionBegin;
2815   ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr);
2816   ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
2817   if (numComp != Nc) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"coordinate discretization must have as many components (%D) as embedding dimension (!= %D)",numComp,Nc);
2818   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
2819   /* convert nodes to values in the stable evaluation basis */
2820   ierr = DMGetWorkArray(dm,pdim,MPIU_REAL,&modes);CHKERRQ(ierr);
2821   invV = fe->invV;
2822   for (i = 0; i < pdim; ++i) {
2823     modes[i] = 0.;
2824     for (j = 0; j < pdim; ++j) {
2825       modes[i] += invV[i * pdim + j] * PetscRealPart(nodes[j]);
2826     }
2827   }
2828   ierr   = DMGetWorkArray(dm,pdim * Nc + pdim * Nc * dimR + Nc,MPIU_REAL,&B);CHKERRQ(ierr);
2829   D      = &B[pdim*Nc];
2830   resNeg = &D[pdim*Nc * dimR];
2831   ierr = DMGetWorkArray(dm,3 * Nc * dimR,MPIU_SCALAR,&J);CHKERRQ(ierr);
2832   invJ = &J[Nc * dimR];
2833   work = &invJ[Nc * dimR];
2834   for (i = 0; i < numPoints * dimR; i++) {refCoords[i] = 0.;}
2835   for (j = 0; j < numPoints; j++) {
2836       for (i = 0; i < maxIter; i++) { /* we could batch this so that we're not making big B and D arrays all the time */
2837       PetscReal *guess = &refCoords[j * dimR];
2838       ierr = PetscSpaceEvaluate(fe->basisSpace, 1, guess, B, D, NULL);CHKERRQ(ierr);
2839       for (k = 0; k < Nc; k++) {resNeg[k] = realCoords[j * Nc + k];}
2840       for (k = 0; k < Nc * dimR; k++) {J[k] = 0.;}
2841       for (k = 0; k < pdim; k++) {
2842         for (l = 0; l < Nc; l++) {
2843           resNeg[l] -= modes[k] * B[k * Nc + l];
2844           for (m = 0; m < dimR; m++) {
2845             J[l * dimR + m] += modes[k] * D[(k * Nc + l) * dimR + m];
2846           }
2847         }
2848       }
2849 #if 0 && defined(PETSC_USE_DEBUG)
2850       {
2851         PetscReal maxAbs = 0.;
2852 
2853         for (l = 0; l < Nc; l++) {
2854           maxAbs = PetscMax(maxAbs,PetscAbsReal(resNeg[l]));
2855         }
2856         ierr = PetscInfo4(dm,"cell %D, point %D, iter %D: res %g\n",cell,j,i,maxAbs);CHKERRQ(ierr);
2857       }
2858 #endif
2859       ierr = DMPlexCoordinatesToReference_NewtonUpdate(Nc,dimR,J,invJ,work,resNeg,guess);CHKERRQ(ierr);
2860     }
2861   }
2862   ierr = DMRestoreWorkArray(dm,3 * Nc * dimR,MPIU_SCALAR,&J);CHKERRQ(ierr);
2863   ierr = DMRestoreWorkArray(dm,pdim * Nc + pdim * Nc * dimR + Nc,MPIU_REAL,&B);CHKERRQ(ierr);
2864   ierr = DMRestoreWorkArray(dm,pdim,MPIU_REAL,&modes);CHKERRQ(ierr);
2865   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
2866   PetscFunctionReturn(0);
2867 }
2868 
2869 /* TODO: TOBY please fix this for Nc > 1 */
2870 static PetscErrorCode DMPlexReferenceToCoordinates_FE(DM dm, PetscFE fe, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[], Vec coords, PetscInt Nc, PetscInt dimR)
2871 {
2872   PetscInt       numComp, pdim, i, j, k, l, coordSize;
2873   PetscScalar    *nodes = NULL;
2874   PetscReal      *invV, *modes;
2875   PetscReal      *B;
2876   PetscErrorCode ierr;
2877 
2878   PetscFunctionBegin;
2879   ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr);
2880   ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
2881   if (numComp != Nc) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"coordinate discretization must have as many components (%D) as embedding dimension (!= %D)",numComp,Nc);
2882   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
2883   /* convert nodes to values in the stable evaluation basis */
2884   ierr = DMGetWorkArray(dm,pdim,MPIU_REAL,&modes);CHKERRQ(ierr);
2885   invV = fe->invV;
2886   for (i = 0; i < pdim; ++i) {
2887     modes[i] = 0.;
2888     for (j = 0; j < pdim; ++j) {
2889       modes[i] += invV[i * pdim + j] * PetscRealPart(nodes[j]);
2890     }
2891   }
2892   ierr = DMGetWorkArray(dm,numPoints * pdim * Nc,MPIU_REAL,&B);CHKERRQ(ierr);
2893   ierr = PetscSpaceEvaluate(fe->basisSpace, numPoints, refCoords, B, NULL, NULL);CHKERRQ(ierr);
2894   for (i = 0; i < numPoints * Nc; i++) {realCoords[i] = 0.;}
2895   for (j = 0; j < numPoints; j++) {
2896     PetscReal *mapped = &realCoords[j * Nc];
2897 
2898     for (k = 0; k < pdim; k++) {
2899       for (l = 0; l < Nc; l++) {
2900         mapped[l] += modes[k] * B[(j * pdim + k) * Nc + l];
2901       }
2902     }
2903   }
2904   ierr = DMRestoreWorkArray(dm,numPoints * pdim * Nc,MPIU_REAL,&B);CHKERRQ(ierr);
2905   ierr = DMRestoreWorkArray(dm,pdim,MPIU_REAL,&modes);CHKERRQ(ierr);
2906   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
2907   PetscFunctionReturn(0);
2908 }
2909 
2910 /*@
2911   DMPlexCoordinatesToReference - Pull coordinates back from the mesh to the reference element using a single element
2912   map.  This inversion will be accurate inside the reference element, but may be inaccurate for mappings that do not
2913   extend uniquely outside the reference cell (e.g, most non-affine maps)
2914 
2915   Not collective
2916 
2917   Input Parameters:
2918 + dm         - The mesh, with coordinate maps defined either by a PetscDS for the coordinate DM (see DMGetCoordinateDM()) or
2919                implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or
2920                as a multilinear map for tensor-product elements
2921 . cell       - the cell whose map is used.
2922 . numPoints  - the number of points to locate
2923 - realCoords - (numPoints x coordinate dimension) array of coordinates (see DMGetCoordinateDim())
2924 
2925   Output Parameters:
2926 . refCoords  - (numPoints x dimension) array of reference coordinates (see DMGetDimension())
2927 
2928   Level: intermediate
2929 
2930 .seealso: DMPlexReferenceToCoordinates()
2931 @*/
2932 PetscErrorCode DMPlexCoordinatesToReference(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[])
2933 {
2934   PetscInt       dimC, dimR, depth, cStart, cEnd, cEndInterior, i;
2935   DM             coordDM = NULL;
2936   Vec            coords;
2937   PetscFE        fe = NULL;
2938   PetscErrorCode ierr;
2939 
2940   PetscFunctionBegin;
2941   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2942   ierr = DMGetDimension(dm,&dimR);CHKERRQ(ierr);
2943   ierr = DMGetCoordinateDim(dm,&dimC);CHKERRQ(ierr);
2944   if (dimR <= 0 || dimC <= 0 || numPoints <= 0) PetscFunctionReturn(0);
2945   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
2946   ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
2947   ierr = DMGetCoordinateDM(dm,&coordDM);CHKERRQ(ierr);
2948   if (coordDM) {
2949     PetscInt coordFields;
2950 
2951     ierr = DMGetNumFields(coordDM,&coordFields);CHKERRQ(ierr);
2952     if (coordFields) {
2953       PetscClassId id;
2954       PetscObject  disc;
2955 
2956       ierr = DMGetField(coordDM,0,NULL,&disc);CHKERRQ(ierr);
2957       ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
2958       if (id == PETSCFE_CLASSID) {
2959         fe = (PetscFE) disc;
2960       }
2961     }
2962   }
2963   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
2964   ierr = DMPlexGetHybridBounds(dm,&cEndInterior,NULL,NULL,NULL);CHKERRQ(ierr);
2965   cEnd = cEndInterior > 0 ? cEndInterior : cEnd;
2966   if (cell < cStart || cell >= cEnd) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D not in cell range [%D,%D)",cell,cStart,cEnd);
2967   if (!fe) { /* implicit discretization: affine or multilinear */
2968     PetscInt  coneSize;
2969     PetscBool isSimplex, isTensor;
2970 
2971     ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr);
2972     isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE;
2973     isTensor  = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE;
2974     if (isSimplex) {
2975       PetscReal detJ, *v0, *J, *invJ;
2976 
2977       ierr = DMGetWorkArray(dm,dimC + 2 * dimC * dimC, MPIU_REAL, &v0);CHKERRQ(ierr);
2978       J    = &v0[dimC];
2979       invJ = &J[dimC * dimC];
2980       ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, &detJ);CHKERRQ(ierr);
2981       for (i = 0; i < numPoints; i++) { /* Apply the inverse affine transformation for each point */
2982         const PetscReal x0[3] = {-1.,-1.,-1.};
2983 
2984         CoordinatesRealToRef(dimC, dimR, x0, v0, invJ, &realCoords[dimC * i], &refCoords[dimR * i]);
2985       }
2986       ierr = DMRestoreWorkArray(dm,dimC + 2 * dimC * dimC, MPIU_REAL, &v0);CHKERRQ(ierr);
2987     } else if (isTensor) {
2988       ierr = DMPlexCoordinatesToReference_Tensor(coordDM, cell, numPoints, realCoords, refCoords, coords, dimC, dimR);CHKERRQ(ierr);
2989     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unrecognized cone size %D",coneSize);
2990   } else {
2991     ierr = DMPlexCoordinatesToReference_FE(coordDM, fe, cell, numPoints, realCoords, refCoords, coords, dimC, dimR);CHKERRQ(ierr);
2992   }
2993   PetscFunctionReturn(0);
2994 }
2995 
2996 /*@
2997   DMPlexReferenceToCoordinates - Map references coordinates to coordinates in the the mesh for a single element map.
2998 
2999   Not collective
3000 
3001   Input Parameters:
3002 + dm         - The mesh, with coordinate maps defined either by a PetscDS for the coordinate DM (see DMGetCoordinateDM()) or
3003                implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or
3004                as a multilinear map for tensor-product elements
3005 . cell       - the cell whose map is used.
3006 . numPoints  - the number of points to locate
3007 + refCoords  - (numPoints x dimension) array of reference coordinates (see DMGetDimension())
3008 
3009   Output Parameters:
3010 . realCoords - (numPoints x coordinate dimension) array of coordinates (see DMGetCoordinateDim())
3011 
3012    Level: intermediate
3013 
3014 .seealso: DMPlexCoordinatesToReference()
3015 @*/
3016 PetscErrorCode DMPlexReferenceToCoordinates(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[])
3017 {
3018   PetscInt       dimC, dimR, depth, cStart, cEnd, cEndInterior, i;
3019   DM             coordDM = NULL;
3020   Vec            coords;
3021   PetscFE        fe = NULL;
3022   PetscErrorCode ierr;
3023 
3024   PetscFunctionBegin;
3025   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3026   ierr = DMGetDimension(dm,&dimR);CHKERRQ(ierr);
3027   ierr = DMGetCoordinateDim(dm,&dimC);CHKERRQ(ierr);
3028   if (dimR <= 0 || dimC <= 0 || numPoints <= 0) PetscFunctionReturn(0);
3029   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
3030   ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
3031   ierr = DMGetCoordinateDM(dm,&coordDM);CHKERRQ(ierr);
3032   if (coordDM) {
3033     PetscInt coordFields;
3034 
3035     ierr = DMGetNumFields(coordDM,&coordFields);CHKERRQ(ierr);
3036     if (coordFields) {
3037       PetscClassId id;
3038       PetscObject  disc;
3039 
3040       ierr = DMGetField(coordDM,0,NULL,&disc);CHKERRQ(ierr);
3041       ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
3042       if (id == PETSCFE_CLASSID) {
3043         fe = (PetscFE) disc;
3044       }
3045     }
3046   }
3047   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
3048   ierr = DMPlexGetHybridBounds(dm,&cEndInterior,NULL,NULL,NULL);CHKERRQ(ierr);
3049   cEnd = cEndInterior > 0 ? cEndInterior : cEnd;
3050   if (cell < cStart || cell >= cEnd) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D not in cell range [%D,%D)",cell,cStart,cEnd);
3051   if (!fe) { /* implicit discretization: affine or multilinear */
3052     PetscInt  coneSize;
3053     PetscBool isSimplex, isTensor;
3054 
3055     ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr);
3056     isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE;
3057     isTensor  = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE;
3058     if (isSimplex) {
3059       PetscReal detJ, *v0, *J;
3060 
3061       ierr = DMGetWorkArray(dm,dimC + 2 * dimC * dimC, MPIU_REAL, &v0);CHKERRQ(ierr);
3062       J    = &v0[dimC];
3063       ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, NULL, &detJ);CHKERRQ(ierr);
3064       for (i = 0; i < numPoints; i++) { /* Apply the affine transformation for each point */
3065         const PetscReal xi0[3] = {-1.,-1.,-1.};
3066 
3067         CoordinatesRefToReal(dimC, dimR, xi0, v0, J, &refCoords[dimR * i], &realCoords[dimC * i]);
3068       }
3069       ierr = DMRestoreWorkArray(dm,dimC + 2 * dimC * dimC, MPIU_REAL, &v0);CHKERRQ(ierr);
3070     } else if (isTensor) {
3071       ierr = DMPlexReferenceToCoordinates_Tensor(coordDM, cell, numPoints, refCoords, realCoords, coords, dimC, dimR);CHKERRQ(ierr);
3072     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unrecognized cone size %D",coneSize);
3073   } else {
3074     ierr = DMPlexReferenceToCoordinates_FE(coordDM, fe, cell, numPoints, refCoords, realCoords, coords, dimC, dimR);CHKERRQ(ierr);
3075   }
3076   PetscFunctionReturn(0);
3077 }
3078