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