xref: /petsc/src/dm/impls/plex/plexgeometry.c (revision ccd2543ff18b66b3097fd04d6e61c25760cf3b9e)
1*ccd2543fSMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2*ccd2543fSMatthew G Knepley 
3*ccd2543fSMatthew G Knepley #undef __FUNCT__
4*ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_Simplex_2D_Internal"
5*ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
6*ccd2543fSMatthew G Knepley {
7*ccd2543fSMatthew G Knepley   const PetscInt embedDim = 2;
8*ccd2543fSMatthew G Knepley   PetscReal      x        = PetscRealPart(point[0]);
9*ccd2543fSMatthew G Knepley   PetscReal      y        = PetscRealPart(point[1]);
10*ccd2543fSMatthew G Knepley   PetscReal      v0[2], J[4], invJ[4], detJ;
11*ccd2543fSMatthew G Knepley   PetscReal      xi, eta;
12*ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
13*ccd2543fSMatthew G Knepley 
14*ccd2543fSMatthew G Knepley   PetscFunctionBegin;
15*ccd2543fSMatthew G Knepley   ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
16*ccd2543fSMatthew G Knepley   xi  = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]);
17*ccd2543fSMatthew G Knepley   eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]);
18*ccd2543fSMatthew G Knepley 
19*ccd2543fSMatthew G Knepley   if ((xi >= 0.0) && (eta >= 0.0) && (xi + eta <= 2.0)) *cell = c;
20*ccd2543fSMatthew G Knepley   else *cell = -1;
21*ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
22*ccd2543fSMatthew G Knepley }
23*ccd2543fSMatthew G Knepley 
24*ccd2543fSMatthew G Knepley #undef __FUNCT__
25*ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_General_2D_Internal"
26*ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_General_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
27*ccd2543fSMatthew G Knepley {
28*ccd2543fSMatthew G Knepley   PetscSection       coordSection;
29*ccd2543fSMatthew G Knepley   Vec                coordsLocal;
30*ccd2543fSMatthew G Knepley   const PetscScalar *coords;
31*ccd2543fSMatthew G Knepley   const PetscInt     faces[8]  = {0, 1, 1, 2, 2, 3, 3, 0};
32*ccd2543fSMatthew G Knepley   PetscReal          x         = PetscRealPart(point[0]);
33*ccd2543fSMatthew G Knepley   PetscReal          y         = PetscRealPart(point[1]);
34*ccd2543fSMatthew G Knepley   PetscInt           crossings = 0, f;
35*ccd2543fSMatthew G Knepley   PetscErrorCode     ierr;
36*ccd2543fSMatthew G Knepley 
37*ccd2543fSMatthew G Knepley   PetscFunctionBegin;
38*ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
39*ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
40*ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
41*ccd2543fSMatthew G Knepley   for (f = 0; f < 4; ++f) {
42*ccd2543fSMatthew G Knepley     PetscReal x_i   = PetscRealPart(coords[faces[2*f+0]*2+0]);
43*ccd2543fSMatthew G Knepley     PetscReal y_i   = PetscRealPart(coords[faces[2*f+0]*2+1]);
44*ccd2543fSMatthew G Knepley     PetscReal x_j   = PetscRealPart(coords[faces[2*f+1]*2+0]);
45*ccd2543fSMatthew G Knepley     PetscReal y_j   = PetscRealPart(coords[faces[2*f+1]*2+1]);
46*ccd2543fSMatthew G Knepley     PetscReal slope = (y_j - y_i) / (x_j - x_i);
47*ccd2543fSMatthew G Knepley     PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE;
48*ccd2543fSMatthew G Knepley     PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE;
49*ccd2543fSMatthew G Knepley     PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE;
50*ccd2543fSMatthew G Knepley     if ((cond1 || cond2)  && above) ++crossings;
51*ccd2543fSMatthew G Knepley   }
52*ccd2543fSMatthew G Knepley   if (crossings % 2) *cell = c;
53*ccd2543fSMatthew G Knepley   else *cell = -1;
54*ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
55*ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
56*ccd2543fSMatthew G Knepley }
57*ccd2543fSMatthew G Knepley 
58*ccd2543fSMatthew G Knepley #undef __FUNCT__
59*ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_Simplex_3D_Internal"
60*ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
61*ccd2543fSMatthew G Knepley {
62*ccd2543fSMatthew G Knepley   const PetscInt embedDim = 3;
63*ccd2543fSMatthew G Knepley   PetscReal      v0[3], J[9], invJ[9], detJ;
64*ccd2543fSMatthew G Knepley   PetscReal      x = PetscRealPart(point[0]);
65*ccd2543fSMatthew G Knepley   PetscReal      y = PetscRealPart(point[1]);
66*ccd2543fSMatthew G Knepley   PetscReal      z = PetscRealPart(point[2]);
67*ccd2543fSMatthew G Knepley   PetscReal      xi, eta, zeta;
68*ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
69*ccd2543fSMatthew G Knepley 
70*ccd2543fSMatthew G Knepley   PetscFunctionBegin;
71*ccd2543fSMatthew G Knepley   ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
72*ccd2543fSMatthew G Knepley   xi   = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]);
73*ccd2543fSMatthew G Knepley   eta  = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]);
74*ccd2543fSMatthew G Knepley   zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]);
75*ccd2543fSMatthew G Knepley 
76*ccd2543fSMatthew G Knepley   if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) *cell = c;
77*ccd2543fSMatthew G Knepley   else *cell = -1;
78*ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
79*ccd2543fSMatthew G Knepley }
80*ccd2543fSMatthew G Knepley 
81*ccd2543fSMatthew G Knepley #undef __FUNCT__
82*ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_General_3D_Internal"
83*ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
84*ccd2543fSMatthew G Knepley {
85*ccd2543fSMatthew G Knepley   PetscSection       coordSection;
86*ccd2543fSMatthew G Knepley   Vec                coordsLocal;
87*ccd2543fSMatthew G Knepley   const PetscScalar *coords;
88*ccd2543fSMatthew G Knepley   const PetscInt     faces[24] = {0, 1, 2, 3,  5, 4, 7, 6,  1, 0, 4, 5,
89*ccd2543fSMatthew G Knepley                                   3, 2, 6, 7,  1, 5, 6, 2,  0, 3, 7, 4};
90*ccd2543fSMatthew G Knepley   PetscBool          found = PETSC_TRUE;
91*ccd2543fSMatthew G Knepley   PetscInt           f;
92*ccd2543fSMatthew G Knepley   PetscErrorCode     ierr;
93*ccd2543fSMatthew G Knepley 
94*ccd2543fSMatthew G Knepley   PetscFunctionBegin;
95*ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
96*ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
97*ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
98*ccd2543fSMatthew G Knepley   for (f = 0; f < 6; ++f) {
99*ccd2543fSMatthew G Knepley     /* Check the point is under plane */
100*ccd2543fSMatthew G Knepley     /*   Get face normal */
101*ccd2543fSMatthew G Knepley     PetscReal v_i[3];
102*ccd2543fSMatthew G Knepley     PetscReal v_j[3];
103*ccd2543fSMatthew G Knepley     PetscReal normal[3];
104*ccd2543fSMatthew G Knepley     PetscReal pp[3];
105*ccd2543fSMatthew G Knepley     PetscReal dot;
106*ccd2543fSMatthew G Knepley 
107*ccd2543fSMatthew G Knepley     v_i[0]    = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]);
108*ccd2543fSMatthew G Knepley     v_i[1]    = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]);
109*ccd2543fSMatthew G Knepley     v_i[2]    = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]);
110*ccd2543fSMatthew G Knepley     v_j[0]    = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]);
111*ccd2543fSMatthew G Knepley     v_j[1]    = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]);
112*ccd2543fSMatthew G Knepley     v_j[2]    = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]);
113*ccd2543fSMatthew G Knepley     normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1];
114*ccd2543fSMatthew G Knepley     normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2];
115*ccd2543fSMatthew G Knepley     normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0];
116*ccd2543fSMatthew G Knepley     pp[0]     = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]);
117*ccd2543fSMatthew G Knepley     pp[1]     = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]);
118*ccd2543fSMatthew G Knepley     pp[2]     = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]);
119*ccd2543fSMatthew G Knepley     dot       = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2];
120*ccd2543fSMatthew G Knepley 
121*ccd2543fSMatthew G Knepley     /* Check that projected point is in face (2D location problem) */
122*ccd2543fSMatthew G Knepley     if (dot < 0.0) {
123*ccd2543fSMatthew G Knepley       found = PETSC_FALSE;
124*ccd2543fSMatthew G Knepley       break;
125*ccd2543fSMatthew G Knepley     }
126*ccd2543fSMatthew G Knepley   }
127*ccd2543fSMatthew G Knepley   if (found) *cell = c;
128*ccd2543fSMatthew G Knepley   else *cell = -1;
129*ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
130*ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
131*ccd2543fSMatthew G Knepley }
132*ccd2543fSMatthew G Knepley 
133*ccd2543fSMatthew G Knepley #undef __FUNCT__
134*ccd2543fSMatthew G Knepley #define __FUNCT__ "DMLocatePoints_Plex"
135*ccd2543fSMatthew G Knepley /*
136*ccd2543fSMatthew G Knepley  Need to implement using the guess
137*ccd2543fSMatthew G Knepley */
138*ccd2543fSMatthew G Knepley PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS)
139*ccd2543fSMatthew G Knepley {
140*ccd2543fSMatthew G Knepley   PetscInt       cell = -1 /*, guess = -1*/;
141*ccd2543fSMatthew G Knepley   PetscInt       bs, numPoints, p;
142*ccd2543fSMatthew G Knepley   PetscInt       dim, cStart, cEnd, cMax, c, coneSize;
143*ccd2543fSMatthew G Knepley   PetscInt      *cells;
144*ccd2543fSMatthew G Knepley   PetscScalar   *a;
145*ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
146*ccd2543fSMatthew G Knepley 
147*ccd2543fSMatthew G Knepley   PetscFunctionBegin;
148*ccd2543fSMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
149*ccd2543fSMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
150*ccd2543fSMatthew G Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
151*ccd2543fSMatthew G Knepley   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
152*ccd2543fSMatthew G Knepley   ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr);
153*ccd2543fSMatthew G Knepley   ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr);
154*ccd2543fSMatthew G Knepley   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
155*ccd2543fSMatthew G Knepley   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*ccd2543fSMatthew G Knepley   numPoints /= bs;
157*ccd2543fSMatthew G Knepley   ierr       = PetscMalloc(numPoints * sizeof(PetscInt), &cells);CHKERRQ(ierr);
158*ccd2543fSMatthew G Knepley   for (p = 0; p < numPoints; ++p) {
159*ccd2543fSMatthew G Knepley     const PetscScalar *point = &a[p*bs];
160*ccd2543fSMatthew G Knepley 
161*ccd2543fSMatthew G Knepley     switch (dim) {
162*ccd2543fSMatthew G Knepley     case 2:
163*ccd2543fSMatthew G Knepley       for (c = cStart; c < cEnd; ++c) {
164*ccd2543fSMatthew G Knepley         ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
165*ccd2543fSMatthew G Knepley         switch (coneSize) {
166*ccd2543fSMatthew G Knepley         case 3:
167*ccd2543fSMatthew G Knepley           ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, c, &cell);CHKERRQ(ierr);
168*ccd2543fSMatthew G Knepley           break;
169*ccd2543fSMatthew G Knepley         case 4:
170*ccd2543fSMatthew G Knepley           ierr = DMPlexLocatePoint_General_2D_Internal(dm, point, c, &cell);CHKERRQ(ierr);
171*ccd2543fSMatthew G Knepley           break;
172*ccd2543fSMatthew G Knepley         default:
173*ccd2543fSMatthew G Knepley           SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %d", coneSize);
174*ccd2543fSMatthew G Knepley         }
175*ccd2543fSMatthew G Knepley         if (cell >= 0) break;
176*ccd2543fSMatthew G Knepley       }
177*ccd2543fSMatthew G Knepley       break;
178*ccd2543fSMatthew G Knepley     case 3:
179*ccd2543fSMatthew G Knepley       for (c = cStart; c < cEnd; ++c) {
180*ccd2543fSMatthew G Knepley         ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
181*ccd2543fSMatthew G Knepley         switch (coneSize) {
182*ccd2543fSMatthew G Knepley         case 4:
183*ccd2543fSMatthew G Knepley           ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, c, &cell);CHKERRQ(ierr);
184*ccd2543fSMatthew G Knepley           break;
185*ccd2543fSMatthew G Knepley         case 8:
186*ccd2543fSMatthew G Knepley           ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, c, &cell);CHKERRQ(ierr);
187*ccd2543fSMatthew G Knepley           break;
188*ccd2543fSMatthew G Knepley         default:
189*ccd2543fSMatthew G Knepley           SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %d", coneSize);
190*ccd2543fSMatthew G Knepley         }
191*ccd2543fSMatthew G Knepley         if (cell >= 0) break;
192*ccd2543fSMatthew G Knepley       }
193*ccd2543fSMatthew G Knepley       break;
194*ccd2543fSMatthew G Knepley     default:
195*ccd2543fSMatthew G Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %d", dim);
196*ccd2543fSMatthew G Knepley     }
197*ccd2543fSMatthew G Knepley     cells[p] = cell;
198*ccd2543fSMatthew G Knepley   }
199*ccd2543fSMatthew G Knepley   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
200*ccd2543fSMatthew G Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints, cells, PETSC_OWN_POINTER, cellIS);CHKERRQ(ierr);
201*ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
202*ccd2543fSMatthew G Knepley }
203*ccd2543fSMatthew G Knepley 
204*ccd2543fSMatthew G Knepley #undef __FUNCT__
205*ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeProjection3Dto2D_Internal"
206*ccd2543fSMatthew G Knepley /*
207*ccd2543fSMatthew G Knepley   DMPlexComputeProjection3Dto2D_Internal - Rewrite coordinates to be the 2D projection of the 3D
208*ccd2543fSMatthew G Knepley */
209*ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscScalar coords[])
210*ccd2543fSMatthew G Knepley {
211*ccd2543fSMatthew G Knepley   PetscScalar    x1[3], x2[3], n[3], norm;
212*ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
213*ccd2543fSMatthew G Knepley   PetscInt       d, e;
214*ccd2543fSMatthew G Knepley 
215*ccd2543fSMatthew G Knepley   PetscFunctionBegin;
216*ccd2543fSMatthew G Knepley   /* 0) Calculate normal vector */
217*ccd2543fSMatthew G Knepley   for (d = 0; d < dim; ++d) {
218*ccd2543fSMatthew G Knepley     x1[d] = coords[1*dim+d] - coords[0*dim+d];
219*ccd2543fSMatthew G Knepley     x2[d] = coords[2*dim+d] - coords[0*dim+d];
220*ccd2543fSMatthew G Knepley   }
221*ccd2543fSMatthew G Knepley   n[0] = x1[1]*x2[2] - x1[2]*x2[1];
222*ccd2543fSMatthew G Knepley   n[1] = x1[2]*x2[0] - x1[0]*x2[2];
223*ccd2543fSMatthew G Knepley   n[2] = x1[0]*x2[1] - x1[1]*x2[0];
224*ccd2543fSMatthew G Knepley   norm = sqrt(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]);
225*ccd2543fSMatthew G Knepley   n[0] /= norm;
226*ccd2543fSMatthew G Knepley   n[1] /= norm;
227*ccd2543fSMatthew G Knepley   n[2] /= norm;
228*ccd2543fSMatthew G Knepley   /* 1) Take the normal vector and rotate until it is \hat z
229*ccd2543fSMatthew G Knepley 
230*ccd2543fSMatthew G Knepley     Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then
231*ccd2543fSMatthew G Knepley 
232*ccd2543fSMatthew G Knepley     R = /  alpha nx nz  alpha ny nz -1/alpha \
233*ccd2543fSMatthew G Knepley         | -alpha ny     alpha nx        0    |
234*ccd2543fSMatthew G Knepley         \     nx            ny         nz    /
235*ccd2543fSMatthew G Knepley 
236*ccd2543fSMatthew G Knepley     will rotate the normal vector to \hat z
237*ccd2543fSMatthew G Knepley   */
238*ccd2543fSMatthew G Knepley   PetscScalar R[9], x1p[3], x2p[3];
239*ccd2543fSMatthew G Knepley   PetscScalar sqrtz = sqrt(1.0 - n[2]*n[2]), alpha = 1.0/sqrtz;
240*ccd2543fSMatthew G Knepley 
241*ccd2543fSMatthew G Knepley   R[0] =  alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz;
242*ccd2543fSMatthew G Knepley   R[3] = -alpha*n[1];      R[4] = alpha*n[0];      R[5] = 0.0;
243*ccd2543fSMatthew G Knepley   R[6] =  n[0];            R[7] = n[1];            R[8] = n[2];
244*ccd2543fSMatthew G Knepley   for (d = 0; d < dim; ++d) {
245*ccd2543fSMatthew G Knepley     x1p[d] = 0.0;
246*ccd2543fSMatthew G Knepley     x2p[d] = 0.0;
247*ccd2543fSMatthew G Knepley     for (e = 0; e < dim; ++e) {
248*ccd2543fSMatthew G Knepley       x1p[d] += R[d*dim+e]*x1[e];
249*ccd2543fSMatthew G Knepley       x2p[d] += R[d*dim+e]*x2[e];
250*ccd2543fSMatthew G Knepley     }
251*ccd2543fSMatthew G Knepley   }
252*ccd2543fSMatthew G Knepley   if (fabs(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
253*ccd2543fSMatthew G Knepley   if (fabs(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
254*ccd2543fSMatthew G Knepley   /* 2) Project to (x, y) */
255*ccd2543fSMatthew G Knepley   coords[0] = 0.0;
256*ccd2543fSMatthew G Knepley   coords[1] = 0.0;
257*ccd2543fSMatthew G Knepley   coords[2] = x1p[0];
258*ccd2543fSMatthew G Knepley   coords[3] = x1p[1];
259*ccd2543fSMatthew G Knepley   coords[4] = x2p[0];
260*ccd2543fSMatthew G Knepley   coords[5] = x2p[1];
261*ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
262*ccd2543fSMatthew G Knepley }
263*ccd2543fSMatthew G Knepley 
264*ccd2543fSMatthew G Knepley #undef __FUNCT__
265*ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal"
266*ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
267*ccd2543fSMatthew G Knepley {
268*ccd2543fSMatthew G Knepley   PetscSection       coordSection;
269*ccd2543fSMatthew G Knepley   Vec                coordinates;
270*ccd2543fSMatthew G Knepley   const PetscScalar *coords;
271*ccd2543fSMatthew G Knepley   const PetscInt     dim = 2;
272*ccd2543fSMatthew G Knepley   PetscInt           numCoords, d, f;
273*ccd2543fSMatthew G Knepley   PetscErrorCode     ierr;
274*ccd2543fSMatthew G Knepley 
275*ccd2543fSMatthew G Knepley   PetscFunctionBegin;
276*ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
277*ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
278*ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
279*ccd2543fSMatthew G Knepley   if (numCoords == 9) {
280*ccd2543fSMatthew G Knepley     ierr = DMPlexComputeProjection3Dto2D_Internal(coords);CHKERRQ(ierr);
281*ccd2543fSMatthew G Knepley   } else if (numCoords != 6) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %d != 6", numCoords);
282*ccd2543fSMatthew G Knepley   if (v0) {
283*ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);
284*ccd2543fSMatthew G Knepley   }
285*ccd2543fSMatthew G Knepley   if (J) {
286*ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
287*ccd2543fSMatthew G Knepley       for (f = 0; f < dim; f++) {
288*ccd2543fSMatthew G Knepley         J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
289*ccd2543fSMatthew G Knepley       }
290*ccd2543fSMatthew G Knepley     }
291*ccd2543fSMatthew G Knepley     *detJ = J[0]*J[3] - J[1]*J[2];
292*ccd2543fSMatthew G Knepley #if 0
293*ccd2543fSMatthew G Knepley     if (detJ < 0.0) {
294*ccd2543fSMatthew G Knepley       const PetscReal xLength = mesh->periodicity[0];
295*ccd2543fSMatthew G Knepley 
296*ccd2543fSMatthew G Knepley       if (xLength != 0.0) {
297*ccd2543fSMatthew G Knepley         PetscReal v0x = coords[0*dim+0];
298*ccd2543fSMatthew G Knepley 
299*ccd2543fSMatthew G Knepley         if (v0x == 0.0) v0x = v0[0] = xLength;
300*ccd2543fSMatthew G Knepley         for (f = 0; f < dim; f++) {
301*ccd2543fSMatthew G Knepley           const PetscReal px = coords[(f+1)*dim+0] == 0.0 ? xLength : coords[(f+1)*dim+0];
302*ccd2543fSMatthew G Knepley 
303*ccd2543fSMatthew G Knepley           J[0*dim+f] = 0.5*(px - v0x);
304*ccd2543fSMatthew G Knepley         }
305*ccd2543fSMatthew G Knepley       }
306*ccd2543fSMatthew G Knepley       detJ = J[0]*J[3] - J[1]*J[2];
307*ccd2543fSMatthew G Knepley     }
308*ccd2543fSMatthew G Knepley #endif
309*ccd2543fSMatthew G Knepley     PetscLogFlops(8.0 + 3.0);
310*ccd2543fSMatthew G Knepley   }
311*ccd2543fSMatthew G Knepley   if (invJ) {
312*ccd2543fSMatthew G Knepley     const PetscReal invDet = 1.0/(*detJ);
313*ccd2543fSMatthew G Knepley 
314*ccd2543fSMatthew G Knepley     invJ[0] =  invDet*J[3];
315*ccd2543fSMatthew G Knepley     invJ[1] = -invDet*J[1];
316*ccd2543fSMatthew G Knepley     invJ[2] = -invDet*J[2];
317*ccd2543fSMatthew G Knepley     invJ[3] =  invDet*J[0];
318*ccd2543fSMatthew G Knepley     PetscLogFlops(5.0);
319*ccd2543fSMatthew G Knepley   }
320*ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
321*ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
322*ccd2543fSMatthew G Knepley }
323*ccd2543fSMatthew G Knepley 
324*ccd2543fSMatthew G Knepley #undef __FUNCT__
325*ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal"
326*ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
327*ccd2543fSMatthew G Knepley {
328*ccd2543fSMatthew G Knepley   PetscSection       coordSection;
329*ccd2543fSMatthew G Knepley   Vec                coordinates;
330*ccd2543fSMatthew G Knepley   const PetscScalar *coords;
331*ccd2543fSMatthew G Knepley   const PetscInt     dim = 2;
332*ccd2543fSMatthew G Knepley   PetscInt           d, f;
333*ccd2543fSMatthew G Knepley   PetscErrorCode     ierr;
334*ccd2543fSMatthew G Knepley 
335*ccd2543fSMatthew G Knepley   PetscFunctionBegin;
336*ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
337*ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
338*ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
339*ccd2543fSMatthew G Knepley   if (v0) {
340*ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);
341*ccd2543fSMatthew G Knepley   }
342*ccd2543fSMatthew G Knepley   if (J) {
343*ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
344*ccd2543fSMatthew G Knepley       for (f = 0; f < dim; f++) {
345*ccd2543fSMatthew G Knepley         J[d*dim+f] = 0.5*(PetscRealPart(coords[(f*2+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
346*ccd2543fSMatthew G Knepley       }
347*ccd2543fSMatthew G Knepley     }
348*ccd2543fSMatthew G Knepley     *detJ = J[0]*J[3] - J[1]*J[2];
349*ccd2543fSMatthew G Knepley     PetscLogFlops(8.0 + 3.0);
350*ccd2543fSMatthew G Knepley   }
351*ccd2543fSMatthew G Knepley   if (invJ) {
352*ccd2543fSMatthew G Knepley     const PetscReal invDet = 1.0/(*detJ);
353*ccd2543fSMatthew G Knepley 
354*ccd2543fSMatthew G Knepley     invJ[0] =  invDet*J[3];
355*ccd2543fSMatthew G Knepley     invJ[1] = -invDet*J[1];
356*ccd2543fSMatthew G Knepley     invJ[2] = -invDet*J[2];
357*ccd2543fSMatthew G Knepley     invJ[3] =  invDet*J[0];
358*ccd2543fSMatthew G Knepley     PetscLogFlops(5.0);
359*ccd2543fSMatthew G Knepley   }
360*ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
361*ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
362*ccd2543fSMatthew G Knepley }
363*ccd2543fSMatthew G Knepley 
364*ccd2543fSMatthew G Knepley #undef __FUNCT__
365*ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal"
366*ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
367*ccd2543fSMatthew G Knepley {
368*ccd2543fSMatthew G Knepley   PetscSection       coordSection;
369*ccd2543fSMatthew G Knepley   Vec                coordinates;
370*ccd2543fSMatthew G Knepley   const PetscScalar *coords;
371*ccd2543fSMatthew G Knepley   const PetscInt     dim = 3;
372*ccd2543fSMatthew G Knepley   PetscInt           d, f;
373*ccd2543fSMatthew G Knepley   PetscErrorCode     ierr;
374*ccd2543fSMatthew G Knepley 
375*ccd2543fSMatthew G Knepley   PetscFunctionBegin;
376*ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
377*ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
378*ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
379*ccd2543fSMatthew G Knepley   if (v0) {
380*ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);
381*ccd2543fSMatthew G Knepley   }
382*ccd2543fSMatthew G Knepley   if (J) {
383*ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
384*ccd2543fSMatthew G Knepley       for (f = 0; f < dim; f++) {
385*ccd2543fSMatthew G Knepley         J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
386*ccd2543fSMatthew G Knepley       }
387*ccd2543fSMatthew G Knepley     }
388*ccd2543fSMatthew G Knepley     /* ??? This does not work with CTetGen: The minus sign is here since I orient the first face to get the outward normal */
389*ccd2543fSMatthew G Knepley     *detJ = (J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
390*ccd2543fSMatthew G Knepley              J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
391*ccd2543fSMatthew G Knepley              J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
392*ccd2543fSMatthew G Knepley     PetscLogFlops(18.0 + 12.0);
393*ccd2543fSMatthew G Knepley   }
394*ccd2543fSMatthew G Knepley   if (invJ) {
395*ccd2543fSMatthew G Knepley     const PetscReal invDet = 1.0/(*detJ);
396*ccd2543fSMatthew G Knepley 
397*ccd2543fSMatthew G Knepley     invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
398*ccd2543fSMatthew G Knepley     invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
399*ccd2543fSMatthew G Knepley     invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
400*ccd2543fSMatthew G Knepley     invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
401*ccd2543fSMatthew G Knepley     invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
402*ccd2543fSMatthew G Knepley     invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
403*ccd2543fSMatthew G Knepley     invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
404*ccd2543fSMatthew G Knepley     invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
405*ccd2543fSMatthew G Knepley     invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
406*ccd2543fSMatthew G Knepley     PetscLogFlops(37.0);
407*ccd2543fSMatthew G Knepley   }
408*ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
409*ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
410*ccd2543fSMatthew G Knepley }
411*ccd2543fSMatthew G Knepley 
412*ccd2543fSMatthew G Knepley #undef __FUNCT__
413*ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal"
414*ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
415*ccd2543fSMatthew G Knepley {
416*ccd2543fSMatthew G Knepley   PetscSection       coordSection;
417*ccd2543fSMatthew G Knepley   Vec                coordinates;
418*ccd2543fSMatthew G Knepley   const PetscScalar *coords;
419*ccd2543fSMatthew G Knepley   const PetscInt     dim = 3;
420*ccd2543fSMatthew G Knepley   PetscInt           d;
421*ccd2543fSMatthew G Knepley   PetscErrorCode     ierr;
422*ccd2543fSMatthew G Knepley 
423*ccd2543fSMatthew G Knepley   PetscFunctionBegin;
424*ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
425*ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
426*ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
427*ccd2543fSMatthew G Knepley   if (v0) {
428*ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);
429*ccd2543fSMatthew G Knepley   }
430*ccd2543fSMatthew G Knepley   if (J) {
431*ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
432*ccd2543fSMatthew G Knepley       J[d*dim+0] = 0.5*(PetscRealPart(coords[(0+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
433*ccd2543fSMatthew G Knepley       J[d*dim+1] = 0.5*(PetscRealPart(coords[(1+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
434*ccd2543fSMatthew G Knepley       J[d*dim+2] = 0.5*(PetscRealPart(coords[(3+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
435*ccd2543fSMatthew G Knepley     }
436*ccd2543fSMatthew G Knepley     *detJ = (J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
437*ccd2543fSMatthew G Knepley              J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
438*ccd2543fSMatthew G Knepley              J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
439*ccd2543fSMatthew G Knepley     PetscLogFlops(18.0 + 12.0);
440*ccd2543fSMatthew G Knepley   }
441*ccd2543fSMatthew G Knepley   if (invJ) {
442*ccd2543fSMatthew G Knepley     const PetscReal invDet = -1.0/(*detJ);
443*ccd2543fSMatthew G Knepley 
444*ccd2543fSMatthew G Knepley     invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
445*ccd2543fSMatthew G Knepley     invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
446*ccd2543fSMatthew G Knepley     invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
447*ccd2543fSMatthew G Knepley     invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
448*ccd2543fSMatthew G Knepley     invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
449*ccd2543fSMatthew G Knepley     invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
450*ccd2543fSMatthew G Knepley     invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
451*ccd2543fSMatthew G Knepley     invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
452*ccd2543fSMatthew G Knepley     invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
453*ccd2543fSMatthew G Knepley     PetscLogFlops(37.0);
454*ccd2543fSMatthew G Knepley   }
455*ccd2543fSMatthew G Knepley   *detJ *= 8.0;
456*ccd2543fSMatthew G Knepley   ierr   = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
457*ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
458*ccd2543fSMatthew G Knepley }
459*ccd2543fSMatthew G Knepley 
460*ccd2543fSMatthew G Knepley #undef __FUNCT__
461*ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeCellGeometry"
462*ccd2543fSMatthew G Knepley /*@C
463*ccd2543fSMatthew G Knepley   DMPlexComputeCellGeometry - Compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell
464*ccd2543fSMatthew G Knepley 
465*ccd2543fSMatthew G Knepley   Collective on DM
466*ccd2543fSMatthew G Knepley 
467*ccd2543fSMatthew G Knepley   Input Arguments:
468*ccd2543fSMatthew G Knepley + dm   - the DM
469*ccd2543fSMatthew G Knepley - cell - the cell
470*ccd2543fSMatthew G Knepley 
471*ccd2543fSMatthew G Knepley   Output Arguments:
472*ccd2543fSMatthew G Knepley + v0   - the translation part of this affine transform
473*ccd2543fSMatthew G Knepley . J    - the Jacobian of the transform from the reference element
474*ccd2543fSMatthew G Knepley . invJ - the inverse of the Jacobian
475*ccd2543fSMatthew G Knepley - detJ - the Jacobian determinant
476*ccd2543fSMatthew G Knepley 
477*ccd2543fSMatthew G Knepley   Level: advanced
478*ccd2543fSMatthew G Knepley 
479*ccd2543fSMatthew G Knepley   Fortran Notes:
480*ccd2543fSMatthew G Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
481*ccd2543fSMatthew G Knepley   include petsc.h90 in your code.
482*ccd2543fSMatthew G Knepley 
483*ccd2543fSMatthew G Knepley .seealso: DMPlexGetCoordinateSection(), DMPlexGetCoordinateVec()
484*ccd2543fSMatthew G Knepley @*/
485*ccd2543fSMatthew G Knepley PetscErrorCode DMPlexComputeCellGeometry(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
486*ccd2543fSMatthew G Knepley {
487*ccd2543fSMatthew G Knepley   PetscInt       dim, coneSize;
488*ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
489*ccd2543fSMatthew G Knepley 
490*ccd2543fSMatthew G Knepley   PetscFunctionBegin;
491*ccd2543fSMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
492*ccd2543fSMatthew G Knepley   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
493*ccd2543fSMatthew G Knepley   switch (dim) {
494*ccd2543fSMatthew G Knepley   case 2:
495*ccd2543fSMatthew G Knepley     switch (coneSize) {
496*ccd2543fSMatthew G Knepley     case 3:
497*ccd2543fSMatthew G Knepley       ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
498*ccd2543fSMatthew G Knepley       break;
499*ccd2543fSMatthew G Knepley     case 4:
500*ccd2543fSMatthew G Knepley       ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
501*ccd2543fSMatthew G Knepley       break;
502*ccd2543fSMatthew G Knepley     default:
503*ccd2543fSMatthew G Knepley       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
504*ccd2543fSMatthew G Knepley     }
505*ccd2543fSMatthew G Knepley     break;
506*ccd2543fSMatthew G Knepley   case 3:
507*ccd2543fSMatthew G Knepley     switch (coneSize) {
508*ccd2543fSMatthew G Knepley     case 4:
509*ccd2543fSMatthew G Knepley       ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
510*ccd2543fSMatthew G Knepley       break;
511*ccd2543fSMatthew G Knepley     case 8:
512*ccd2543fSMatthew G Knepley       ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
513*ccd2543fSMatthew G Knepley       break;
514*ccd2543fSMatthew G Knepley     default:
515*ccd2543fSMatthew G Knepley       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
516*ccd2543fSMatthew G Knepley     }
517*ccd2543fSMatthew G Knepley     break;
518*ccd2543fSMatthew G Knepley   default:
519*ccd2543fSMatthew G Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
520*ccd2543fSMatthew G Knepley   }
521*ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
522*ccd2543fSMatthew G Knepley }
523