xref: /petsc/src/dm/impls/plex/plexgeometry.c (revision cc08537e2ac1abf2e280922c2ca6efaf01b264e6)
1 #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 
3 #undef __FUNCT__
4 #define __FUNCT__ "DMPlexLocatePoint_Simplex_2D_Internal"
5 static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
6 {
7   const PetscInt embedDim = 2;
8   PetscReal      x        = PetscRealPart(point[0]);
9   PetscReal      y        = PetscRealPart(point[1]);
10   PetscReal      v0[2], J[4], invJ[4], detJ;
11   PetscReal      xi, eta;
12   PetscErrorCode ierr;
13 
14   PetscFunctionBegin;
15   ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
16   xi  = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]);
17   eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]);
18 
19   if ((xi >= 0.0) && (eta >= 0.0) && (xi + eta <= 2.0)) *cell = c;
20   else *cell = -1;
21   PetscFunctionReturn(0);
22 }
23 
24 #undef __FUNCT__
25 #define __FUNCT__ "DMPlexLocatePoint_General_2D_Internal"
26 static PetscErrorCode DMPlexLocatePoint_General_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
27 {
28   PetscSection       coordSection;
29   Vec             coordsLocal;
30   PetscScalar    *coords;
31   const PetscInt  faces[8]  = {0, 1, 1, 2, 2, 3, 3, 0};
32   PetscReal       x         = PetscRealPart(point[0]);
33   PetscReal       y         = PetscRealPart(point[1]);
34   PetscInt        crossings = 0, f;
35   PetscErrorCode  ierr;
36 
37   PetscFunctionBegin;
38   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
39   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
40   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
41   for (f = 0; f < 4; ++f) {
42     PetscReal x_i   = PetscRealPart(coords[faces[2*f+0]*2+0]);
43     PetscReal y_i   = PetscRealPart(coords[faces[2*f+0]*2+1]);
44     PetscReal x_j   = PetscRealPart(coords[faces[2*f+1]*2+0]);
45     PetscReal y_j   = PetscRealPart(coords[faces[2*f+1]*2+1]);
46     PetscReal slope = (y_j - y_i) / (x_j - x_i);
47     PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE;
48     PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE;
49     PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE;
50     if ((cond1 || cond2)  && above) ++crossings;
51   }
52   if (crossings % 2) *cell = c;
53   else *cell = -1;
54   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
55   PetscFunctionReturn(0);
56 }
57 
58 #undef __FUNCT__
59 #define __FUNCT__ "DMPlexLocatePoint_Simplex_3D_Internal"
60 static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
61 {
62   const PetscInt embedDim = 3;
63   PetscReal      v0[3], J[9], invJ[9], detJ;
64   PetscReal      x = PetscRealPart(point[0]);
65   PetscReal      y = PetscRealPart(point[1]);
66   PetscReal      z = PetscRealPart(point[2]);
67   PetscReal      xi, eta, zeta;
68   PetscErrorCode ierr;
69 
70   PetscFunctionBegin;
71   ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
72   xi   = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]);
73   eta  = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]);
74   zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]);
75 
76   if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) *cell = c;
77   else *cell = -1;
78   PetscFunctionReturn(0);
79 }
80 
81 #undef __FUNCT__
82 #define __FUNCT__ "DMPlexLocatePoint_General_3D_Internal"
83 static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
84 {
85   PetscSection       coordSection;
86   Vec            coordsLocal;
87   PetscScalar   *coords;
88   const PetscInt faces[24] = {0, 1, 2, 3,  5, 4, 7, 6,  1, 0, 4, 5,
89                               3, 2, 6, 7,  1, 5, 6, 2,  0, 3, 7, 4};
90   PetscBool      found = PETSC_TRUE;
91   PetscInt       f;
92   PetscErrorCode ierr;
93 
94   PetscFunctionBegin;
95   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
96   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
97   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
98   for (f = 0; f < 6; ++f) {
99     /* Check the point is under plane */
100     /*   Get face normal */
101     PetscReal v_i[3];
102     PetscReal v_j[3];
103     PetscReal normal[3];
104     PetscReal pp[3];
105     PetscReal dot;
106 
107     v_i[0]    = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]);
108     v_i[1]    = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]);
109     v_i[2]    = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]);
110     v_j[0]    = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]);
111     v_j[1]    = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]);
112     v_j[2]    = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]);
113     normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1];
114     normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2];
115     normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0];
116     pp[0]     = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]);
117     pp[1]     = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]);
118     pp[2]     = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]);
119     dot       = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2];
120 
121     /* Check that projected point is in face (2D location problem) */
122     if (dot < 0.0) {
123       found = PETSC_FALSE;
124       break;
125     }
126   }
127   if (found) *cell = c;
128   else *cell = -1;
129   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
130   PetscFunctionReturn(0);
131 }
132 
133 #undef __FUNCT__
134 #define __FUNCT__ "DMLocatePoints_Plex"
135 /*
136  Need to implement using the guess
137 */
138 PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS)
139 {
140   PetscInt       cell = -1 /*, guess = -1*/;
141   PetscInt       bs, numPoints, p;
142   PetscInt       dim, cStart, cEnd, cMax, c, coneSize;
143   PetscInt      *cells;
144   PetscScalar   *a;
145   PetscErrorCode ierr;
146 
147   PetscFunctionBegin;
148   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
149   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
150   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
151   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
152   ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr);
153   ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr);
154   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
155   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);
156   numPoints /= bs;
157   ierr       = PetscMalloc(numPoints * sizeof(PetscInt), &cells);CHKERRQ(ierr);
158   for (p = 0; p < numPoints; ++p) {
159     const PetscScalar *point = &a[p*bs];
160 
161     switch (dim) {
162     case 2:
163       for (c = cStart; c < cEnd; ++c) {
164         ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
165         switch (coneSize) {
166         case 3:
167           ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, c, &cell);CHKERRQ(ierr);
168           break;
169         case 4:
170           ierr = DMPlexLocatePoint_General_2D_Internal(dm, point, c, &cell);CHKERRQ(ierr);
171           break;
172         default:
173           SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %d", coneSize);
174         }
175         if (cell >= 0) break;
176       }
177       break;
178     case 3:
179       for (c = cStart; c < cEnd; ++c) {
180         ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
181         switch (coneSize) {
182         case 4:
183           ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, c, &cell);CHKERRQ(ierr);
184           break;
185         case 8:
186           ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, c, &cell);CHKERRQ(ierr);
187           break;
188         default:
189           SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %d", coneSize);
190         }
191         if (cell >= 0) break;
192       }
193       break;
194     default:
195       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %d", dim);
196     }
197     cells[p] = cell;
198   }
199   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
200   ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints, cells, PETSC_OWN_POINTER, cellIS);CHKERRQ(ierr);
201   PetscFunctionReturn(0);
202 }
203 
204 #undef __FUNCT__
205 #define __FUNCT__ "DMPlexComputeProjection2Dto1D_Internal"
206 /*
207   DMPlexComputeProjection2Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 2D
208 */
209 static PetscErrorCode DMPlexComputeProjection2Dto1D_Internal(PetscScalar coords[], PetscReal R[])
210 {
211   const PetscReal x = PetscRealPart(coords[2] - coords[0]);
212   const PetscReal y = PetscRealPart(coords[3] - coords[1]);
213   const PetscReal r = sqrt(x*x + y*y), c = x/r, s = y/r;
214 
215   PetscFunctionBegin;
216   R[0] =  c; R[1] = s;
217   R[2] = -s; R[3] = c;
218   coords[0] = 0.0;
219   coords[1] = r;
220   PetscFunctionReturn(0);
221 }
222 
223 #undef __FUNCT__
224 #define __FUNCT__ "DMPlexComputeProjection3Dto2D_Internal"
225 /*
226   DMPlexComputeProjection3Dto2D_Internal - Rewrite coordinates to be the 2D projection of the 3D
227 */
228 static PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscScalar coords[], PetscReal R[])
229 {
230   PetscScalar    x1[3],  x2[3], n[3], norm;
231   PetscScalar    x1p[3], x2p[3];
232   PetscReal      sqrtz, alpha;
233   const PetscInt dim = 3;
234   PetscInt       d, e;
235 
236   PetscFunctionBegin;
237   /* 0) Calculate normal vector */
238   for (d = 0; d < dim; ++d) {
239     x1[d] = coords[1*dim+d] - coords[0*dim+d];
240     x2[d] = coords[2*dim+d] - coords[0*dim+d];
241   }
242   n[0] = x1[1]*x2[2] - x1[2]*x2[1];
243   n[1] = x1[2]*x2[0] - x1[0]*x2[2];
244   n[2] = x1[0]*x2[1] - x1[1]*x2[0];
245   norm = sqrt(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]);
246   n[0] /= norm;
247   n[1] /= norm;
248   n[2] /= norm;
249   /* 1) Take the normal vector and rotate until it is \hat z
250 
251     Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then
252 
253     R = /  alpha nx nz  alpha ny nz -1/alpha \
254         | -alpha ny     alpha nx        0    |
255         \     nx            ny         nz    /
256 
257     will rotate the normal vector to \hat z
258   */
259   sqrtz = sqrt(1.0 - PetscAbsScalar(n[2]*n[2]));
260   /* Check for n = z */
261   if (sqrtz < 1.0e-10) {
262     coords[0] = 0.0;
263     coords[1] = 0.0;
264     if (PetscRealPart(n[2]) < 0.0) {
265       coords[2] = x2[0];
266       coords[3] = x2[1];
267       coords[4] = x1[0];
268       coords[5] = x1[1];
269     } else {
270       coords[2] = x1[0];
271       coords[3] = x1[1];
272       coords[4] = x2[0];
273       coords[5] = x2[1];
274     }
275     PetscFunctionReturn(0);
276   }
277   alpha = 1.0/sqrtz;
278   R[0] =  alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz;
279   R[3] = -alpha*n[1];      R[4] = alpha*n[0];      R[5] = 0.0;
280   R[6] =  n[0];            R[7] = n[1];            R[8] = n[2];
281   for (d = 0; d < dim; ++d) {
282     x1p[d] = 0.0;
283     x2p[d] = 0.0;
284     for (e = 0; e < dim; ++e) {
285       x1p[d] += R[d*dim+e]*x1[e];
286       x2p[d] += R[d*dim+e]*x2[e];
287     }
288   }
289   if (PetscAbsScalar(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
290   if (PetscAbsScalar(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
291   /* 2) Project to (x, y) */
292   coords[0] = 0.0;
293   coords[1] = 0.0;
294   coords[2] = x1p[0];
295   coords[3] = x1p[1];
296   coords[4] = x2p[0];
297   coords[5] = x2p[1];
298   /* Output R^T which rotates \hat z to the input normal */
299   for (d = 0; d < dim; ++d) {
300     for (e = d+1; e < dim; ++e) {
301       PetscReal tmp;
302 
303       tmp        = R[d*dim+e];
304       R[d*dim+e] = R[e*dim+d];
305       R[e*dim+d] = tmp;
306     }
307   }
308   PetscFunctionReturn(0);
309 }
310 
311 #undef __FUNCT__
312 #define __FUNCT__ "Invert2D_Internal"
313 PETSC_STATIC_INLINE void Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
314 {
315   const PetscReal invDet = 1.0/detJ;
316 
317   invJ[0] =  invDet*J[3];
318   invJ[1] = -invDet*J[1];
319   invJ[2] = -invDet*J[2];
320   invJ[3] =  invDet*J[0];
321   PetscLogFlops(5.0);
322 }
323 
324 #undef __FUNCT__
325 #define __FUNCT__ "Invert3D_Internal"
326 PETSC_STATIC_INLINE void Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
327 {
328   const PetscReal invDet = 1.0/detJ;
329 
330   invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
331   invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
332   invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
333   invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
334   invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
335   invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
336   invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
337   invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
338   invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
339   PetscLogFlops(37.0);
340 }
341 
342 #undef __FUNCT__
343 #define __FUNCT__ "Det2D_Internal"
344 PETSC_STATIC_INLINE void Det2D_Internal(PetscReal *detJ, PetscReal J[])
345 {
346   *detJ = J[0]*J[3] - J[1]*J[2];
347   PetscLogFlops(3.0);
348 }
349 
350 #undef __FUNCT__
351 #define __FUNCT__ "Det3D_Internal"
352 PETSC_STATIC_INLINE void Det3D_Internal(PetscReal *detJ, PetscReal J[])
353 {
354   *detJ = (J[0*3+0]*(J[1*3+2]*J[2*3+1] - J[1*3+1]*J[2*3+2]) +
355            J[0*3+1]*(J[1*3+0]*J[2*3+2] - J[1*3+2]*J[2*3+0]) +
356            J[0*3+2]*(J[1*3+1]*J[2*3+0] - J[1*3+0]*J[2*3+1]));
357   PetscLogFlops(12.0);
358 }
359 
360 #undef __FUNCT__
361 #define __FUNCT__ "Volume_Triangle_Internal"
362 PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[])
363 {
364   /* Signed volume is 1/2 the determinant
365 
366    |  1  1  1 |
367    | x0 x1 x2 |
368    | y0 y1 y2 |
369 
370      but if x0,y0 is the origin, we have
371 
372    | x1 x2 |
373    | y1 y2 |
374   */
375   const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1];
376   const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1];
377   PetscReal       M[4], detM;
378   M[0] = x1; M[1] = x2;
379   M[3] = y1; M[4] = y2;
380   Det2D_Internal(&detM, M);
381   *vol = 0.5*detM;
382   PetscLogFlops(5.0);
383 }
384 
385 #undef __FUNCT__
386 #define __FUNCT__ "Volume_Triangle_Origin_Internal"
387 PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[])
388 {
389   Det2D_Internal(vol, coords);
390   *vol *= 0.5;
391 }
392 
393 #undef __FUNCT__
394 #define __FUNCT__ "Volume_Tetrahedron_Internal"
395 PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[])
396 {
397   /* Signed volume is 1/6th of the determinant
398 
399    |  1  1  1  1 |
400    | x0 x1 x2 x3 |
401    | y0 y1 y2 y3 |
402    | z0 z1 z2 z3 |
403 
404      but if x0,y0,z0 is the origin, we have
405 
406    | x1 x2 x3 |
407    | y1 y2 y3 |
408    | z1 z2 z3 |
409   */
410   const PetscReal x1 = coords[3] - coords[0], y1 = coords[4]  - coords[1], z1 = coords[5]  - coords[2];
411   const PetscReal x2 = coords[6] - coords[0], y2 = coords[7]  - coords[1], z2 = coords[8]  - coords[2];
412   const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2];
413   PetscReal       M[9], detM;
414   M[0] = x1; M[1] = x2; M[2] = x3;
415   M[3] = y1; M[4] = y2; M[5] = y3;
416   M[6] = z1; M[7] = z2; M[8] = z3;
417   Det3D_Internal(&detM, M);
418   *vol = 0.16666666666666666666666*detM;
419   PetscLogFlops(10.0);
420 }
421 
422 #undef __FUNCT__
423 #define __FUNCT__ "DMPlexComputeLineGeometry_Internal"
424 static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
425 {
426   PetscSection   coordSection;
427   Vec            coordinates;
428   PetscScalar   *coords;
429   PetscInt       numCoords, d;
430   PetscErrorCode ierr;
431 
432   PetscFunctionBegin;
433   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
434   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
435   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
436   *detJ = 0.0;
437   if (numCoords == 4) {
438     const PetscInt dim = 2;
439     PetscReal      R[4], J0;
440 
441     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
442     ierr = DMPlexComputeProjection2Dto1D_Internal(coords, R);CHKERRQ(ierr);
443     if (J)    {
444       J0   = 0.5*PetscRealPart(coords[1]);
445       J[0] = R[0]*J0; J[1] = R[1];
446       J[2] = R[2]*J0; J[3] = R[3];
447       Det2D_Internal(detJ, J);
448     }
449     if (invJ) {Invert2D_Internal(invJ, J, *detJ);}
450   } else if (numCoords == 2) {
451     const PetscInt dim = 1;
452 
453     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
454     if (J)    {
455       J[0]  = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0]));
456       *detJ = J[0];
457       PetscLogFlops(2.0);
458     }
459     if (invJ) {invJ[0] = 1.0/J[0]; PetscLogFlops(1.0);}
460   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %d != 2", numCoords);
461   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
462   PetscFunctionReturn(0);
463 }
464 
465 #undef __FUNCT__
466 #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal"
467 static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
468 {
469   PetscSection   coordSection;
470   Vec            coordinates;
471   PetscScalar   *coords;
472   PetscInt       numCoords, d, f, g;
473   PetscErrorCode ierr;
474 
475   PetscFunctionBegin;
476   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
477   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
478   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
479   *detJ = 0.0;
480   if (numCoords == 9) {
481     const PetscInt dim = 3;
482     PetscReal      R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0};
483 
484     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
485     ierr = DMPlexComputeProjection3Dto2D_Internal(coords, R);CHKERRQ(ierr);
486     if (J)    {
487       for (d = 0; d < dim-1; d++) {
488         for (f = 0; f < dim-1; f++) {
489           J0[d*(dim+1)+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
490         }
491       }
492       PetscLogFlops(8.0);
493       for (d = 0; d < dim; d++) {
494         for (f = 0; f < dim; f++) {
495           J[d*dim+f] = 0.0;
496           for (g = 0; g < dim; g++) {
497             J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
498           }
499         }
500       }
501       PetscLogFlops(18.0);
502       Det3D_Internal(detJ, J);
503     }
504     if (invJ) {Invert3D_Internal(invJ, J, *detJ);}
505   } else if (numCoords == 6) {
506     const PetscInt dim = 2;
507 
508     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
509     if (J)    {
510       for (d = 0; d < dim; d++) {
511         for (f = 0; f < dim; f++) {
512           J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
513         }
514       }
515       PetscLogFlops(8.0);
516       Det2D_Internal(detJ, J);
517     }
518     if (invJ) {Invert2D_Internal(invJ, J, *detJ);}
519   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %d != 6", numCoords);
520   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
521   PetscFunctionReturn(0);
522 }
523 
524 #undef __FUNCT__
525 #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal"
526 static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
527 {
528   PetscSection   coordSection;
529   Vec            coordinates;
530   PetscScalar   *coords;
531   const PetscInt dim = 2;
532   PetscInt       d, f;
533   PetscErrorCode ierr;
534 
535   PetscFunctionBegin;
536   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
537   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
538   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
539   *detJ = 0.0;
540   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
541   if (J)    {
542     for (d = 0; d < dim; d++) {
543       for (f = 0; f < dim; f++) {
544         J[d*dim+f] = 0.5*(PetscRealPart(coords[(f*2+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
545       }
546     }
547     PetscLogFlops(8.0);
548     Det2D_Internal(detJ, J);
549   }
550   if (invJ) {Invert2D_Internal(invJ, J, *detJ);}
551   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
552   PetscFunctionReturn(0);
553 }
554 
555 #undef __FUNCT__
556 #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal"
557 static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
558 {
559   PetscSection   coordSection;
560   Vec            coordinates;
561   PetscScalar   *coords;
562   const PetscInt dim = 3;
563   PetscInt       d, f;
564   PetscErrorCode ierr;
565 
566   PetscFunctionBegin;
567   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
568   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
569   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
570   *detJ = 0.0;
571   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
572   if (J)    {
573     for (d = 0; d < dim; d++) {
574       for (f = 0; f < dim; f++) {
575         J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
576       }
577     }
578     PetscLogFlops(18.0);
579     Det3D_Internal(detJ, J);
580   }
581   if (invJ) {Invert3D_Internal(invJ, J, *detJ);}
582   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
583   PetscFunctionReturn(0);
584 }
585 
586 #undef __FUNCT__
587 #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal"
588 static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
589 {
590   PetscSection   coordSection;
591   Vec            coordinates;
592   PetscScalar   *coords;
593   const PetscInt dim = 3;
594   PetscInt       d;
595   PetscErrorCode ierr;
596 
597   PetscFunctionBegin;
598   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
599   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
600   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
601   *detJ = 0.0;
602   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
603   if (J)    {
604     for (d = 0; d < dim; d++) {
605       J[d*dim+0] = 0.5*(PetscRealPart(coords[(2+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
606       J[d*dim+1] = 0.5*(PetscRealPart(coords[(1+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
607       J[d*dim+2] = 0.5*(PetscRealPart(coords[(3+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
608     }
609     PetscLogFlops(18.0);
610     Det3D_Internal(detJ, J);
611   }
612   if (invJ) {Invert3D_Internal(invJ, J, *detJ);}
613   *detJ *= -8.0;
614   ierr   = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
615   PetscFunctionReturn(0);
616 }
617 
618 #undef __FUNCT__
619 #define __FUNCT__ "DMPlexComputeCellGeometry"
620 /*@C
621   DMPlexComputeCellGeometry - Compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell
622 
623   Collective on DM
624 
625   Input Arguments:
626 + dm   - the DM
627 - cell - the cell
628 
629   Output Arguments:
630 + v0   - the translation part of this affine transform
631 . J    - the Jacobian of the transform from the reference element
632 . invJ - the inverse of the Jacobian
633 - detJ - the Jacobian determinant
634 
635   Level: advanced
636 
637   Fortran Notes:
638   Since it returns arrays, this routine is only available in Fortran 90, and you must
639   include petsc.h90 in your code.
640 
641 .seealso: DMPlexGetCoordinateSection(), DMPlexGetCoordinateVec()
642 @*/
643 PetscErrorCode DMPlexComputeCellGeometry(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
644 {
645   PetscInt       depth, dim, coneSize;
646   PetscErrorCode ierr;
647 
648   PetscFunctionBegin;
649   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
650   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
651   if (depth == 1) {
652     ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
653     switch (dim) {
654     case 1:
655       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
656       break;
657     case 2:
658       switch (coneSize) {
659       case 3:
660         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
661         break;
662       case 4:
663         ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
664         break;
665       default:
666         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
667       }
668       break;
669     case 3:
670       switch (coneSize) {
671       case 4:
672         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
673         break;
674       case 8:
675         ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
676         break;
677       default:
678         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
679       }
680       break;
681     default:
682       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
683     }
684   } else {
685     /* We need to keep a pointer to the depth label */
686     ierr = DMPlexGetLabelValue(dm, "depth", cell, &dim);CHKERRQ(ierr);
687     /* Cone size is now the number of faces */
688     switch (dim) {
689     case 1:
690       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
691       break;
692     case 2:
693       switch (coneSize) {
694       case 3:
695         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
696         break;
697       case 4:
698         ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
699         break;
700       default:
701         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
702       }
703       break;
704     case 3:
705       switch (coneSize) {
706       case 4:
707         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
708         break;
709       case 6:
710         ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
711         break;
712       default:
713         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
714       }
715       break;
716     default:
717       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
718     }
719   }
720   PetscFunctionReturn(0);
721 }
722 
723 #undef __FUNCT__
724 #define __FUNCT__ "DMPlexComputeGeometryFVM_1D_Internal"
725 static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
726 {
727   PetscSection   coordSection;
728   Vec            coordinates;
729   PetscScalar   *coords;
730   const PetscInt dim = 2;
731   PetscInt       coordSize;
732   PetscErrorCode ierr;
733 
734   PetscFunctionBegin;
735   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
736   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
737   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
738   if (coordSize != dim*2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now");
739   if (centroid) {
740     centroid[0] = 0.5*(coords[0] + coords[dim+0]);
741     centroid[1] = 0.5*(coords[1] + coords[dim+1]);
742   }
743   if (normal) {
744     normal[0] =  (coords[1] - coords[dim+1]);
745     normal[1] = -(coords[0] - coords[dim+0]);
746   }
747   if (vol) {
748     *vol = sqrt(PetscSqr(coords[0] - coords[dim+0]) + PetscSqr(coords[1] - coords[dim+1]));
749   }
750   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
751   PetscFunctionReturn(0);
752 }
753 
754 #undef __FUNCT__
755 #define __FUNCT__ "DMPlexComputeGeometryFVM_2D_Internal"
756 /* Centroid_i = (\sum_n A_n Cn_i ) / A */
757 static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
758 {
759   PetscSection   coordSection;
760   Vec            coordinates;
761   PetscScalar   *coords = NULL;
762   PetscReal      vsum, csum[3], vtmp, ctmp[6];
763   const PetscInt dim = 2;
764   PetscInt       coordSize, numCorners, p, d;
765   PetscErrorCode ierr;
766 
767   PetscFunctionBegin;
768   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
769   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
770   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
771   numCorners = coordSize/dim;
772   for (p = 0; p < numCorners-1; ++p) {
773     Volume_Triangle_Origin_Internal(&vtmp, &coords[p*dim]);
774     vsum += vtmp;
775     for (d = 0; d < dim; ++d) {
776       csum[d] += (coords[p*dim+d] + coords[(p+1)*dim+d])*vtmp;
777     }
778   }
779   for (d = 0; d < dim; ++d) {
780     ctmp[d]     = coords[(numCorners-1)*dim+d];
781     ctmp[dim+d] = coords[d];
782   }
783   Volume_Triangle_Origin_Internal(&vtmp, ctmp);
784   vsum += vtmp;
785   for (d = 0; d < dim; ++d) {
786     csum[d] += (coords[(numCorners-1)*dim+d] + coords[0*dim+d])*vtmp;
787     csum[d] /= (dim+1)*vsum;
788   }
789   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
790   if (vol) *vol = PetscAbsScalar(vsum);
791   if (centroid) for (d = 0; d < dim; ++d) centroid[d] = csum[d];
792   if (normal)   for (d = 0; d < dim; ++d) normal[d]   = 0.0;
793   PetscFunctionReturn(0);
794 }
795 
796 #undef __FUNCT__
797 #define __FUNCT__ "DMPlexComputeCellGeometryFVM"
798 /*@C
799   DMPlexComputeCellGeometryFVM - Compute the volume for a given cell
800 
801   Collective on DM
802 
803   Input Arguments:
804 + dm   - the DM
805 - cell - the cell
806 
807   Output Arguments:
808 + volume   - the cell volume
809 . centroid - the cell centroid
810 - normal - the cell normal, if appropriate
811 
812   Level: advanced
813 
814   Fortran Notes:
815   Since it returns arrays, this routine is only available in Fortran 90, and you must
816   include petsc.h90 in your code.
817 
818 .seealso: DMPlexGetCoordinateSection(), DMPlexGetCoordinateVec()
819 @*/
820 PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
821 {
822   PetscInt       depth, dim, coneSize;
823   PetscErrorCode ierr;
824 
825   PetscFunctionBegin;
826   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
827   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
828   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
829   if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated");
830   /* We need to keep a pointer to the depth label */
831   ierr = DMPlexGetLabelValue(dm, "depth", cell, &dim);CHKERRQ(ierr);
832   /* Cone size is now the number of faces */
833   switch (dim) {
834   case 1:
835     ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, cell, vol, centroid, normal);CHKERRQ(ierr);
836     break;
837   case 2:
838     ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, cell, vol, centroid, normal);CHKERRQ(ierr);
839     break;
840   case 3:
841     switch (coneSize) {
842     default:
843       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
844     }
845     break;
846   default:
847     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
848   }
849   PetscFunctionReturn(0);
850 }
851