xref: /petsc/src/dm/impls/plex/plexgeometry.c (revision 17fe8556b0f3b6e25aea499ab158df50dce0bf76)
1ccd2543fSMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2ccd2543fSMatthew G Knepley 
3ccd2543fSMatthew G Knepley #undef __FUNCT__
4ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_Simplex_2D_Internal"
5ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
6ccd2543fSMatthew G Knepley {
7ccd2543fSMatthew G Knepley   const PetscInt embedDim = 2;
8ccd2543fSMatthew G Knepley   PetscReal      x        = PetscRealPart(point[0]);
9ccd2543fSMatthew G Knepley   PetscReal      y        = PetscRealPart(point[1]);
10ccd2543fSMatthew G Knepley   PetscReal      v0[2], J[4], invJ[4], detJ;
11ccd2543fSMatthew G Knepley   PetscReal      xi, eta;
12ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
13ccd2543fSMatthew G Knepley 
14ccd2543fSMatthew G Knepley   PetscFunctionBegin;
15ccd2543fSMatthew G Knepley   ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
16ccd2543fSMatthew G Knepley   xi  = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]);
17ccd2543fSMatthew G Knepley   eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]);
18ccd2543fSMatthew G Knepley 
19ccd2543fSMatthew G Knepley   if ((xi >= 0.0) && (eta >= 0.0) && (xi + eta <= 2.0)) *cell = c;
20ccd2543fSMatthew G Knepley   else *cell = -1;
21ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
22ccd2543fSMatthew G Knepley }
23ccd2543fSMatthew G Knepley 
24ccd2543fSMatthew G Knepley #undef __FUNCT__
25ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_General_2D_Internal"
26ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_General_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
27ccd2543fSMatthew G Knepley {
28ccd2543fSMatthew G Knepley   PetscSection       coordSection;
29ccd2543fSMatthew G Knepley   Vec             coordsLocal;
307c1f9639SMatthew G Knepley   PetscScalar    *coords;
31ccd2543fSMatthew G Knepley   const PetscInt  faces[8]  = {0, 1, 1, 2, 2, 3, 3, 0};
32ccd2543fSMatthew G Knepley   PetscReal       x         = PetscRealPart(point[0]);
33ccd2543fSMatthew G Knepley   PetscReal       y         = PetscRealPart(point[1]);
34ccd2543fSMatthew G Knepley   PetscInt        crossings = 0, f;
35ccd2543fSMatthew G Knepley   PetscErrorCode  ierr;
36ccd2543fSMatthew G Knepley 
37ccd2543fSMatthew G Knepley   PetscFunctionBegin;
38ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
39ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
40ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
41ccd2543fSMatthew G Knepley   for (f = 0; f < 4; ++f) {
42ccd2543fSMatthew G Knepley     PetscReal x_i   = PetscRealPart(coords[faces[2*f+0]*2+0]);
43ccd2543fSMatthew G Knepley     PetscReal y_i   = PetscRealPart(coords[faces[2*f+0]*2+1]);
44ccd2543fSMatthew G Knepley     PetscReal x_j   = PetscRealPart(coords[faces[2*f+1]*2+0]);
45ccd2543fSMatthew G Knepley     PetscReal y_j   = PetscRealPart(coords[faces[2*f+1]*2+1]);
46ccd2543fSMatthew G Knepley     PetscReal slope = (y_j - y_i) / (x_j - x_i);
47ccd2543fSMatthew G Knepley     PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE;
48ccd2543fSMatthew G Knepley     PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE;
49ccd2543fSMatthew G Knepley     PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE;
50ccd2543fSMatthew G Knepley     if ((cond1 || cond2)  && above) ++crossings;
51ccd2543fSMatthew G Knepley   }
52ccd2543fSMatthew G Knepley   if (crossings % 2) *cell = c;
53ccd2543fSMatthew G Knepley   else *cell = -1;
54ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
55ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
56ccd2543fSMatthew G Knepley }
57ccd2543fSMatthew G Knepley 
58ccd2543fSMatthew G Knepley #undef __FUNCT__
59ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_Simplex_3D_Internal"
60ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
61ccd2543fSMatthew G Knepley {
62ccd2543fSMatthew G Knepley   const PetscInt embedDim = 3;
63ccd2543fSMatthew G Knepley   PetscReal      v0[3], J[9], invJ[9], detJ;
64ccd2543fSMatthew G Knepley   PetscReal      x = PetscRealPart(point[0]);
65ccd2543fSMatthew G Knepley   PetscReal      y = PetscRealPart(point[1]);
66ccd2543fSMatthew G Knepley   PetscReal      z = PetscRealPart(point[2]);
67ccd2543fSMatthew G Knepley   PetscReal      xi, eta, zeta;
68ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
69ccd2543fSMatthew G Knepley 
70ccd2543fSMatthew G Knepley   PetscFunctionBegin;
71ccd2543fSMatthew G Knepley   ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
72ccd2543fSMatthew G Knepley   xi   = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]);
73ccd2543fSMatthew G Knepley   eta  = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]);
74ccd2543fSMatthew G Knepley   zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]);
75ccd2543fSMatthew G Knepley 
76ccd2543fSMatthew G Knepley   if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) *cell = c;
77ccd2543fSMatthew G Knepley   else *cell = -1;
78ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
79ccd2543fSMatthew G Knepley }
80ccd2543fSMatthew G Knepley 
81ccd2543fSMatthew G Knepley #undef __FUNCT__
82ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_General_3D_Internal"
83ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
84ccd2543fSMatthew G Knepley {
85ccd2543fSMatthew G Knepley   PetscSection       coordSection;
86ccd2543fSMatthew G Knepley   Vec            coordsLocal;
877c1f9639SMatthew G Knepley   PetscScalar   *coords;
88ccd2543fSMatthew G Knepley   const PetscInt faces[24] = {0, 1, 2, 3,  5, 4, 7, 6,  1, 0, 4, 5,
89ccd2543fSMatthew G Knepley                               3, 2, 6, 7,  1, 5, 6, 2,  0, 3, 7, 4};
90ccd2543fSMatthew G Knepley   PetscBool      found = PETSC_TRUE;
91ccd2543fSMatthew G Knepley   PetscInt       f;
92ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
93ccd2543fSMatthew G Knepley 
94ccd2543fSMatthew G Knepley   PetscFunctionBegin;
95ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
96ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
97ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
98ccd2543fSMatthew G Knepley   for (f = 0; f < 6; ++f) {
99ccd2543fSMatthew G Knepley     /* Check the point is under plane */
100ccd2543fSMatthew G Knepley     /*   Get face normal */
101ccd2543fSMatthew G Knepley     PetscReal v_i[3];
102ccd2543fSMatthew G Knepley     PetscReal v_j[3];
103ccd2543fSMatthew G Knepley     PetscReal normal[3];
104ccd2543fSMatthew G Knepley     PetscReal pp[3];
105ccd2543fSMatthew G Knepley     PetscReal dot;
106ccd2543fSMatthew G Knepley 
107ccd2543fSMatthew G Knepley     v_i[0]    = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]);
108ccd2543fSMatthew G Knepley     v_i[1]    = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]);
109ccd2543fSMatthew G Knepley     v_i[2]    = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]);
110ccd2543fSMatthew G Knepley     v_j[0]    = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]);
111ccd2543fSMatthew G Knepley     v_j[1]    = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]);
112ccd2543fSMatthew G Knepley     v_j[2]    = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]);
113ccd2543fSMatthew G Knepley     normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1];
114ccd2543fSMatthew G Knepley     normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2];
115ccd2543fSMatthew G Knepley     normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0];
116ccd2543fSMatthew G Knepley     pp[0]     = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]);
117ccd2543fSMatthew G Knepley     pp[1]     = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]);
118ccd2543fSMatthew G Knepley     pp[2]     = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]);
119ccd2543fSMatthew G Knepley     dot       = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2];
120ccd2543fSMatthew G Knepley 
121ccd2543fSMatthew G Knepley     /* Check that projected point is in face (2D location problem) */
122ccd2543fSMatthew G Knepley     if (dot < 0.0) {
123ccd2543fSMatthew G Knepley       found = PETSC_FALSE;
124ccd2543fSMatthew G Knepley       break;
125ccd2543fSMatthew G Knepley     }
126ccd2543fSMatthew G Knepley   }
127ccd2543fSMatthew G Knepley   if (found) *cell = c;
128ccd2543fSMatthew G Knepley   else *cell = -1;
129ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
130ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
131ccd2543fSMatthew G Knepley }
132ccd2543fSMatthew G Knepley 
133ccd2543fSMatthew G Knepley #undef __FUNCT__
134ccd2543fSMatthew G Knepley #define __FUNCT__ "DMLocatePoints_Plex"
135ccd2543fSMatthew G Knepley /*
136ccd2543fSMatthew G Knepley  Need to implement using the guess
137ccd2543fSMatthew G Knepley */
138ccd2543fSMatthew G Knepley PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS)
139ccd2543fSMatthew G Knepley {
140ccd2543fSMatthew G Knepley   PetscInt       cell = -1 /*, guess = -1*/;
141ccd2543fSMatthew G Knepley   PetscInt       bs, numPoints, p;
142ccd2543fSMatthew G Knepley   PetscInt       dim, cStart, cEnd, cMax, c, coneSize;
143ccd2543fSMatthew G Knepley   PetscInt      *cells;
144ccd2543fSMatthew G Knepley   PetscScalar   *a;
145ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
146ccd2543fSMatthew G Knepley 
147ccd2543fSMatthew G Knepley   PetscFunctionBegin;
148ccd2543fSMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
149ccd2543fSMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
150ccd2543fSMatthew G Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
151ccd2543fSMatthew G Knepley   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
152ccd2543fSMatthew G Knepley   ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr);
153ccd2543fSMatthew G Knepley   ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr);
154ccd2543fSMatthew G Knepley   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
155ccd2543fSMatthew 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);
156ccd2543fSMatthew G Knepley   numPoints /= bs;
157ccd2543fSMatthew G Knepley   ierr       = PetscMalloc(numPoints * sizeof(PetscInt), &cells);CHKERRQ(ierr);
158ccd2543fSMatthew G Knepley   for (p = 0; p < numPoints; ++p) {
159ccd2543fSMatthew G Knepley     const PetscScalar *point = &a[p*bs];
160ccd2543fSMatthew G Knepley 
161ccd2543fSMatthew G Knepley     switch (dim) {
162ccd2543fSMatthew G Knepley     case 2:
163ccd2543fSMatthew G Knepley       for (c = cStart; c < cEnd; ++c) {
164ccd2543fSMatthew G Knepley         ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
165ccd2543fSMatthew G Knepley         switch (coneSize) {
166ccd2543fSMatthew G Knepley         case 3:
167ccd2543fSMatthew G Knepley           ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, c, &cell);CHKERRQ(ierr);
168ccd2543fSMatthew G Knepley           break;
169ccd2543fSMatthew G Knepley         case 4:
170ccd2543fSMatthew G Knepley           ierr = DMPlexLocatePoint_General_2D_Internal(dm, point, c, &cell);CHKERRQ(ierr);
171ccd2543fSMatthew G Knepley           break;
172ccd2543fSMatthew G Knepley         default:
173ccd2543fSMatthew G Knepley           SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %d", coneSize);
174ccd2543fSMatthew G Knepley         }
175ccd2543fSMatthew G Knepley         if (cell >= 0) break;
176ccd2543fSMatthew G Knepley       }
177ccd2543fSMatthew G Knepley       break;
178ccd2543fSMatthew G Knepley     case 3:
179ccd2543fSMatthew G Knepley       for (c = cStart; c < cEnd; ++c) {
180ccd2543fSMatthew G Knepley         ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
181ccd2543fSMatthew G Knepley         switch (coneSize) {
182ccd2543fSMatthew G Knepley         case 4:
183ccd2543fSMatthew G Knepley           ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, c, &cell);CHKERRQ(ierr);
184ccd2543fSMatthew G Knepley           break;
185ccd2543fSMatthew G Knepley         case 8:
186ccd2543fSMatthew G Knepley           ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, c, &cell);CHKERRQ(ierr);
187ccd2543fSMatthew G Knepley           break;
188ccd2543fSMatthew G Knepley         default:
189ccd2543fSMatthew G Knepley           SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %d", coneSize);
190ccd2543fSMatthew G Knepley         }
191ccd2543fSMatthew G Knepley         if (cell >= 0) break;
192ccd2543fSMatthew G Knepley       }
193ccd2543fSMatthew G Knepley       break;
194ccd2543fSMatthew G Knepley     default:
195ccd2543fSMatthew G Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %d", dim);
196ccd2543fSMatthew G Knepley     }
197ccd2543fSMatthew G Knepley     cells[p] = cell;
198ccd2543fSMatthew G Knepley   }
199ccd2543fSMatthew G Knepley   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
200ccd2543fSMatthew G Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints, cells, PETSC_OWN_POINTER, cellIS);CHKERRQ(ierr);
201ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
202ccd2543fSMatthew G Knepley }
203ccd2543fSMatthew G Knepley 
204ccd2543fSMatthew G Knepley #undef __FUNCT__
205*17fe8556SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeProjection2Dto1D_Internal"
206*17fe8556SMatthew G. Knepley /*
207*17fe8556SMatthew G. Knepley   DMPlexComputeProjection2Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 2D
208*17fe8556SMatthew G. Knepley */
209*17fe8556SMatthew G. Knepley static PetscErrorCode DMPlexComputeProjection2Dto1D_Internal(PetscScalar coords[])
210*17fe8556SMatthew G. Knepley {
211*17fe8556SMatthew G. Knepley   const PetscReal x = PetscRealPart(coords[2] - coords[0]);
212*17fe8556SMatthew G. Knepley   const PetscReal y = PetscRealPart(coords[3] - coords[1]);
213*17fe8556SMatthew G. Knepley 
214*17fe8556SMatthew G. Knepley   PetscFunctionBegin;
215*17fe8556SMatthew G. Knepley   coords[0] = 0.0;
216*17fe8556SMatthew G. Knepley   coords[1] = sqrt(x*x + y*y);
217*17fe8556SMatthew G. Knepley   PetscFunctionReturn(0);
218*17fe8556SMatthew G. Knepley }
219*17fe8556SMatthew G. Knepley 
220*17fe8556SMatthew G. Knepley #undef __FUNCT__
221ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeProjection3Dto2D_Internal"
222ccd2543fSMatthew G Knepley /*
223ccd2543fSMatthew G Knepley   DMPlexComputeProjection3Dto2D_Internal - Rewrite coordinates to be the 2D projection of the 3D
224ccd2543fSMatthew G Knepley */
225ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscScalar coords[])
226ccd2543fSMatthew G Knepley {
227ccd2543fSMatthew G Knepley   PetscScalar    x1[3], x2[3], n[3], norm;
228da18b5e6SMatthew G Knepley   PetscScalar    R[9], x1p[3], x2p[3];
2294a217a95SMatthew G. Knepley   PetscReal      sqrtz, alpha;
230ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
231ccd2543fSMatthew G Knepley   PetscInt       d, e;
232ccd2543fSMatthew G Knepley 
233ccd2543fSMatthew G Knepley   PetscFunctionBegin;
234ccd2543fSMatthew G Knepley   /* 0) Calculate normal vector */
235ccd2543fSMatthew G Knepley   for (d = 0; d < dim; ++d) {
236ccd2543fSMatthew G Knepley     x1[d] = coords[1*dim+d] - coords[0*dim+d];
237ccd2543fSMatthew G Knepley     x2[d] = coords[2*dim+d] - coords[0*dim+d];
238ccd2543fSMatthew G Knepley   }
239ccd2543fSMatthew G Knepley   n[0] = x1[1]*x2[2] - x1[2]*x2[1];
240ccd2543fSMatthew G Knepley   n[1] = x1[2]*x2[0] - x1[0]*x2[2];
241ccd2543fSMatthew G Knepley   n[2] = x1[0]*x2[1] - x1[1]*x2[0];
242ccd2543fSMatthew G Knepley   norm = sqrt(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]);
243ccd2543fSMatthew G Knepley   n[0] /= norm;
244ccd2543fSMatthew G Knepley   n[1] /= norm;
245ccd2543fSMatthew G Knepley   n[2] /= norm;
246ccd2543fSMatthew G Knepley   /* 1) Take the normal vector and rotate until it is \hat z
247ccd2543fSMatthew G Knepley 
248ccd2543fSMatthew G Knepley     Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then
249ccd2543fSMatthew G Knepley 
250ccd2543fSMatthew G Knepley     R = /  alpha nx nz  alpha ny nz -1/alpha \
251ccd2543fSMatthew G Knepley         | -alpha ny     alpha nx        0    |
252ccd2543fSMatthew G Knepley         \     nx            ny         nz    /
253ccd2543fSMatthew G Knepley 
254ccd2543fSMatthew G Knepley     will rotate the normal vector to \hat z
255ccd2543fSMatthew G Knepley   */
2564a217a95SMatthew G. Knepley   sqrtz = sqrt(1.0 - PetscAbsScalar(n[2]*n[2]));
25773868372SMatthew G. Knepley   /* Check for n = z */
25873868372SMatthew G. Knepley   if (sqrtz < 1.0e-10) {
25973868372SMatthew G. Knepley     coords[0] = 0.0;
26073868372SMatthew G. Knepley     coords[1] = 0.0;
2614a217a95SMatthew G. Knepley     if (PetscRealPart(n[2]) < 0.0) {
26273868372SMatthew G. Knepley       coords[2] = x2[0];
26373868372SMatthew G. Knepley       coords[3] = x2[1];
26473868372SMatthew G. Knepley       coords[4] = x1[0];
26573868372SMatthew G. Knepley       coords[5] = x1[1];
26673868372SMatthew G. Knepley     } else {
26773868372SMatthew G. Knepley       coords[2] = x1[0];
26873868372SMatthew G. Knepley       coords[3] = x1[1];
26973868372SMatthew G. Knepley       coords[4] = x2[0];
27073868372SMatthew G. Knepley       coords[5] = x2[1];
27173868372SMatthew G. Knepley     }
27273868372SMatthew G. Knepley     PetscFunctionReturn(0);
27373868372SMatthew G. Knepley   }
274da18b5e6SMatthew G Knepley   alpha = 1.0/sqrtz;
275ccd2543fSMatthew G Knepley   R[0] =  alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz;
276ccd2543fSMatthew G Knepley   R[3] = -alpha*n[1];      R[4] = alpha*n[0];      R[5] = 0.0;
277ccd2543fSMatthew G Knepley   R[6] =  n[0];            R[7] = n[1];            R[8] = n[2];
278ccd2543fSMatthew G Knepley   for (d = 0; d < dim; ++d) {
279ccd2543fSMatthew G Knepley     x1p[d] = 0.0;
280ccd2543fSMatthew G Knepley     x2p[d] = 0.0;
281ccd2543fSMatthew G Knepley     for (e = 0; e < dim; ++e) {
282ccd2543fSMatthew G Knepley       x1p[d] += R[d*dim+e]*x1[e];
283ccd2543fSMatthew G Knepley       x2p[d] += R[d*dim+e]*x2[e];
284ccd2543fSMatthew G Knepley     }
285ccd2543fSMatthew G Knepley   }
28688790e04SMatthew G Knepley   if (PetscAbsScalar(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
28788790e04SMatthew G Knepley   if (PetscAbsScalar(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
288ccd2543fSMatthew G Knepley   /* 2) Project to (x, y) */
289ccd2543fSMatthew G Knepley   coords[0] = 0.0;
290ccd2543fSMatthew G Knepley   coords[1] = 0.0;
291ccd2543fSMatthew G Knepley   coords[2] = x1p[0];
292ccd2543fSMatthew G Knepley   coords[3] = x1p[1];
293ccd2543fSMatthew G Knepley   coords[4] = x2p[0];
294ccd2543fSMatthew G Knepley   coords[5] = x2p[1];
295ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
296ccd2543fSMatthew G Knepley }
297ccd2543fSMatthew G Knepley 
298ccd2543fSMatthew G Knepley #undef __FUNCT__
299*17fe8556SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeLineGeometry_Internal"
300*17fe8556SMatthew G. Knepley static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
301*17fe8556SMatthew G. Knepley {
302*17fe8556SMatthew G. Knepley   PetscSection   coordSection;
303*17fe8556SMatthew G. Knepley   Vec            coordinates;
304*17fe8556SMatthew G. Knepley   PetscScalar   *coords;
305*17fe8556SMatthew G. Knepley   const PetscInt dim = 1;
306*17fe8556SMatthew G. Knepley   PetscInt       numCoords, d, f;
307*17fe8556SMatthew G. Knepley   PetscErrorCode ierr;
308*17fe8556SMatthew G. Knepley 
309*17fe8556SMatthew G. Knepley   PetscFunctionBegin;
310*17fe8556SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
311*17fe8556SMatthew G. Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
312*17fe8556SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
313*17fe8556SMatthew G. Knepley   if (numCoords == 4) {
314*17fe8556SMatthew G. Knepley     ierr = DMPlexComputeProjection2Dto1D_Internal(coords);CHKERRQ(ierr);
315*17fe8556SMatthew G. Knepley   } else if (numCoords != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %d != 2", numCoords);
316*17fe8556SMatthew G. Knepley   if (v0) {
317*17fe8556SMatthew G. Knepley     for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);
318*17fe8556SMatthew G. Knepley   }
319*17fe8556SMatthew G. Knepley   if (J) {
320*17fe8556SMatthew G. Knepley     for (d = 0; d < dim; d++) {
321*17fe8556SMatthew G. Knepley       for (f = 0; f < dim; f++) {
322*17fe8556SMatthew G. Knepley         J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
323*17fe8556SMatthew G. Knepley       }
324*17fe8556SMatthew G. Knepley     }
325*17fe8556SMatthew G. Knepley     *detJ = J[0];
326*17fe8556SMatthew G. Knepley     PetscLogFlops(2.0);
327*17fe8556SMatthew G. Knepley   } else {
328*17fe8556SMatthew G. Knepley     *detJ = 0.0;
329*17fe8556SMatthew G. Knepley   }
330*17fe8556SMatthew G. Knepley   if (invJ) {
331*17fe8556SMatthew G. Knepley     invJ[0] =  1.0/J[0];
332*17fe8556SMatthew G. Knepley     PetscLogFlops(1.0);
333*17fe8556SMatthew G. Knepley   }
334*17fe8556SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
335*17fe8556SMatthew G. Knepley   PetscFunctionReturn(0);
336*17fe8556SMatthew G. Knepley }
337*17fe8556SMatthew G. Knepley 
338*17fe8556SMatthew G. Knepley #undef __FUNCT__
339ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal"
340ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
341ccd2543fSMatthew G Knepley {
342ccd2543fSMatthew G Knepley   PetscSection   coordSection;
343ccd2543fSMatthew G Knepley   Vec            coordinates;
3447c1f9639SMatthew G Knepley   PetscScalar   *coords;
345ccd2543fSMatthew G Knepley   const PetscInt dim = 2;
346ccd2543fSMatthew G Knepley   PetscInt       numCoords, d, f;
347ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
348ccd2543fSMatthew G Knepley 
349ccd2543fSMatthew G Knepley   PetscFunctionBegin;
350ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
351ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
352ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
353ccd2543fSMatthew G Knepley   if (numCoords == 9) {
354ccd2543fSMatthew G Knepley     ierr = DMPlexComputeProjection3Dto2D_Internal(coords);CHKERRQ(ierr);
355ccd2543fSMatthew 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);
356ccd2543fSMatthew G Knepley   if (v0) {
357ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);
358ccd2543fSMatthew G Knepley   }
359ccd2543fSMatthew G Knepley   if (J) {
360ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
361ccd2543fSMatthew G Knepley       for (f = 0; f < dim; f++) {
362ccd2543fSMatthew G Knepley         J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
363ccd2543fSMatthew G Knepley       }
364ccd2543fSMatthew G Knepley     }
365ccd2543fSMatthew G Knepley     *detJ = J[0]*J[3] - J[1]*J[2];
366ccd2543fSMatthew G Knepley #if 0
367ccd2543fSMatthew G Knepley     if (detJ < 0.0) {
368ccd2543fSMatthew G Knepley       const PetscReal xLength = mesh->periodicity[0];
369ccd2543fSMatthew G Knepley 
370ccd2543fSMatthew G Knepley       if (xLength != 0.0) {
371ccd2543fSMatthew G Knepley         PetscReal v0x = coords[0*dim+0];
372ccd2543fSMatthew G Knepley 
373ccd2543fSMatthew G Knepley         if (v0x == 0.0) v0x = v0[0] = xLength;
374ccd2543fSMatthew G Knepley         for (f = 0; f < dim; f++) {
375ccd2543fSMatthew G Knepley           const PetscReal px = coords[(f+1)*dim+0] == 0.0 ? xLength : coords[(f+1)*dim+0];
376ccd2543fSMatthew G Knepley 
377ccd2543fSMatthew G Knepley           J[0*dim+f] = 0.5*(px - v0x);
378ccd2543fSMatthew G Knepley         }
379ccd2543fSMatthew G Knepley       }
380ccd2543fSMatthew G Knepley       detJ = J[0]*J[3] - J[1]*J[2];
381ccd2543fSMatthew G Knepley     }
382ccd2543fSMatthew G Knepley #endif
383ccd2543fSMatthew G Knepley     PetscLogFlops(8.0 + 3.0);
384923cbbb6SMatthew G. Knepley   } else {
385923cbbb6SMatthew G. Knepley     *detJ = 0.0;
386ccd2543fSMatthew G Knepley   }
387ccd2543fSMatthew G Knepley   if (invJ) {
388ccd2543fSMatthew G Knepley     const PetscReal invDet = 1.0/(*detJ);
389ccd2543fSMatthew G Knepley 
390ccd2543fSMatthew G Knepley     invJ[0] =  invDet*J[3];
391ccd2543fSMatthew G Knepley     invJ[1] = -invDet*J[1];
392ccd2543fSMatthew G Knepley     invJ[2] = -invDet*J[2];
393ccd2543fSMatthew G Knepley     invJ[3] =  invDet*J[0];
394ccd2543fSMatthew G Knepley     PetscLogFlops(5.0);
395ccd2543fSMatthew G Knepley   }
396ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
397ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
398ccd2543fSMatthew G Knepley }
399ccd2543fSMatthew G Knepley 
400ccd2543fSMatthew G Knepley #undef __FUNCT__
401ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal"
402ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
403ccd2543fSMatthew G Knepley {
404ccd2543fSMatthew G Knepley   PetscSection   coordSection;
405ccd2543fSMatthew G Knepley   Vec            coordinates;
4067c1f9639SMatthew G Knepley   PetscScalar   *coords;
407ccd2543fSMatthew G Knepley   const PetscInt dim = 2;
408ccd2543fSMatthew G Knepley   PetscInt       d, f;
409ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
410ccd2543fSMatthew G Knepley 
411ccd2543fSMatthew G Knepley   PetscFunctionBegin;
412ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
413ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
414ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
415ccd2543fSMatthew G Knepley   if (v0) {
416ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);
417ccd2543fSMatthew G Knepley   }
418ccd2543fSMatthew G Knepley   if (J) {
419ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
420ccd2543fSMatthew G Knepley       for (f = 0; f < dim; f++) {
421ccd2543fSMatthew G Knepley         J[d*dim+f] = 0.5*(PetscRealPart(coords[(f*2+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
422ccd2543fSMatthew G Knepley       }
423ccd2543fSMatthew G Knepley     }
424ccd2543fSMatthew G Knepley     *detJ = J[0]*J[3] - J[1]*J[2];
425ccd2543fSMatthew G Knepley     PetscLogFlops(8.0 + 3.0);
426923cbbb6SMatthew G. Knepley   } else {
427923cbbb6SMatthew G. Knepley     *detJ = 0.0;
428ccd2543fSMatthew G Knepley   }
429ccd2543fSMatthew G Knepley   if (invJ) {
430ccd2543fSMatthew G Knepley     const PetscReal invDet = 1.0/(*detJ);
431ccd2543fSMatthew G Knepley 
432ccd2543fSMatthew G Knepley     invJ[0] =  invDet*J[3];
433ccd2543fSMatthew G Knepley     invJ[1] = -invDet*J[1];
434ccd2543fSMatthew G Knepley     invJ[2] = -invDet*J[2];
435ccd2543fSMatthew G Knepley     invJ[3] =  invDet*J[0];
436ccd2543fSMatthew G Knepley     PetscLogFlops(5.0);
437ccd2543fSMatthew G Knepley   }
438ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
439ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
440ccd2543fSMatthew G Knepley }
441ccd2543fSMatthew G Knepley 
442ccd2543fSMatthew G Knepley #undef __FUNCT__
443ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal"
444ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
445ccd2543fSMatthew G Knepley {
446ccd2543fSMatthew G Knepley   PetscSection   coordSection;
447ccd2543fSMatthew G Knepley   Vec            coordinates;
4487c1f9639SMatthew G Knepley   PetscScalar   *coords;
449ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
450ccd2543fSMatthew G Knepley   PetscInt       d, f;
451ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
452ccd2543fSMatthew G Knepley 
453ccd2543fSMatthew G Knepley   PetscFunctionBegin;
454ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
455ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
456ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
457ccd2543fSMatthew G Knepley   if (v0) {
458ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);
459ccd2543fSMatthew G Knepley   }
460ccd2543fSMatthew G Knepley   if (J) {
461ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
462ccd2543fSMatthew G Knepley       for (f = 0; f < dim; f++) {
463ccd2543fSMatthew G Knepley         J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
464ccd2543fSMatthew G Knepley       }
465ccd2543fSMatthew G Knepley     }
4662e1b13c2SMatthew G. Knepley     *detJ = (J[0*3+0]*(J[1*3+2]*J[2*3+1] - J[1*3+1]*J[2*3+2]) +
4672e1b13c2SMatthew G. Knepley              J[0*3+1]*(J[1*3+0]*J[2*3+2] - J[1*3+2]*J[2*3+0]) +
4682e1b13c2SMatthew G. Knepley              J[0*3+2]*(J[1*3+1]*J[2*3+0] - J[1*3+0]*J[2*3+1]));
469ccd2543fSMatthew G Knepley     PetscLogFlops(18.0 + 12.0);
470ccd2543fSMatthew G Knepley   }
471ccd2543fSMatthew G Knepley   if (invJ) {
472ccd2543fSMatthew G Knepley     const PetscReal invDet = 1.0/(*detJ);
473ccd2543fSMatthew G Knepley 
474ccd2543fSMatthew G Knepley     invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
475ccd2543fSMatthew G Knepley     invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
476ccd2543fSMatthew G Knepley     invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
477ccd2543fSMatthew G Knepley     invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
478ccd2543fSMatthew G Knepley     invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
479ccd2543fSMatthew G Knepley     invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
480ccd2543fSMatthew G Knepley     invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
481ccd2543fSMatthew G Knepley     invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
482ccd2543fSMatthew G Knepley     invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
483ccd2543fSMatthew G Knepley     PetscLogFlops(37.0);
484ccd2543fSMatthew G Knepley   }
485ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
486ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
487ccd2543fSMatthew G Knepley }
488ccd2543fSMatthew G Knepley 
489ccd2543fSMatthew G Knepley #undef __FUNCT__
490ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal"
491ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
492ccd2543fSMatthew G Knepley {
493ccd2543fSMatthew G Knepley   PetscSection   coordSection;
494ccd2543fSMatthew G Knepley   Vec            coordinates;
4957c1f9639SMatthew G Knepley   PetscScalar   *coords;
496ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
497ccd2543fSMatthew G Knepley   PetscInt       d;
498ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
499ccd2543fSMatthew G Knepley 
500ccd2543fSMatthew G Knepley   PetscFunctionBegin;
501ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
502ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
503ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
504ccd2543fSMatthew G Knepley   if (v0) {
505ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);
506ccd2543fSMatthew G Knepley   }
507ccd2543fSMatthew G Knepley   if (J) {
508ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
509ccd2543fSMatthew G Knepley       J[d*dim+0] = 0.5*(PetscRealPart(coords[(0+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
510ccd2543fSMatthew G Knepley       J[d*dim+1] = 0.5*(PetscRealPart(coords[(1+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
511ccd2543fSMatthew G Knepley       J[d*dim+2] = 0.5*(PetscRealPart(coords[(3+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
512ccd2543fSMatthew G Knepley     }
5132e1b13c2SMatthew G. Knepley     *detJ = (J[0*3+0]*(J[1*3+2]*J[2*3+1] - J[1*3+1]*J[2*3+2]) +
5142e1b13c2SMatthew G. Knepley              J[0*3+1]*(J[1*3+0]*J[2*3+2] - J[1*3+2]*J[2*3+0]) +
5152e1b13c2SMatthew G. Knepley              J[0*3+2]*(J[1*3+1]*J[2*3+0] - J[1*3+0]*J[2*3+1]));
516ccd2543fSMatthew G Knepley     PetscLogFlops(18.0 + 12.0);
517923cbbb6SMatthew G. Knepley   } else {
518923cbbb6SMatthew G. Knepley     *detJ = 0.0;
519ccd2543fSMatthew G Knepley   }
520ccd2543fSMatthew G Knepley   if (invJ) {
521ccd2543fSMatthew G Knepley     const PetscReal invDet = -1.0/(*detJ);
522ccd2543fSMatthew G Knepley 
523ccd2543fSMatthew G Knepley     invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
524ccd2543fSMatthew G Knepley     invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
525ccd2543fSMatthew G Knepley     invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
526ccd2543fSMatthew G Knepley     invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
527ccd2543fSMatthew G Knepley     invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
528ccd2543fSMatthew G Knepley     invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
529ccd2543fSMatthew G Knepley     invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
530ccd2543fSMatthew G Knepley     invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
531ccd2543fSMatthew G Knepley     invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
532ccd2543fSMatthew G Knepley     PetscLogFlops(37.0);
533ccd2543fSMatthew G Knepley   }
534ccd2543fSMatthew G Knepley   *detJ *= 8.0;
535ccd2543fSMatthew G Knepley   ierr   = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
536ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
537ccd2543fSMatthew G Knepley }
538ccd2543fSMatthew G Knepley 
539ccd2543fSMatthew G Knepley #undef __FUNCT__
540ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeCellGeometry"
541ccd2543fSMatthew G Knepley /*@C
542ccd2543fSMatthew G Knepley   DMPlexComputeCellGeometry - Compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell
543ccd2543fSMatthew G Knepley 
544ccd2543fSMatthew G Knepley   Collective on DM
545ccd2543fSMatthew G Knepley 
546ccd2543fSMatthew G Knepley   Input Arguments:
547ccd2543fSMatthew G Knepley + dm   - the DM
548ccd2543fSMatthew G Knepley - cell - the cell
549ccd2543fSMatthew G Knepley 
550ccd2543fSMatthew G Knepley   Output Arguments:
551ccd2543fSMatthew G Knepley + v0   - the translation part of this affine transform
552ccd2543fSMatthew G Knepley . J    - the Jacobian of the transform from the reference element
553ccd2543fSMatthew G Knepley . invJ - the inverse of the Jacobian
554ccd2543fSMatthew G Knepley - detJ - the Jacobian determinant
555ccd2543fSMatthew G Knepley 
556ccd2543fSMatthew G Knepley   Level: advanced
557ccd2543fSMatthew G Knepley 
558ccd2543fSMatthew G Knepley   Fortran Notes:
559ccd2543fSMatthew G Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
560ccd2543fSMatthew G Knepley   include petsc.h90 in your code.
561ccd2543fSMatthew G Knepley 
562ccd2543fSMatthew G Knepley .seealso: DMPlexGetCoordinateSection(), DMPlexGetCoordinateVec()
563ccd2543fSMatthew G Knepley @*/
564ccd2543fSMatthew G Knepley PetscErrorCode DMPlexComputeCellGeometry(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
565ccd2543fSMatthew G Knepley {
56649dc4407SMatthew G. Knepley   PetscInt       depth, dim, coneSize;
567ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
568ccd2543fSMatthew G Knepley 
569ccd2543fSMatthew G Knepley   PetscFunctionBegin;
570139a35ccSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
571ccd2543fSMatthew G Knepley   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
57249dc4407SMatthew G. Knepley   if (depth == 1) {
5735743f1d7SMatthew G. Knepley     ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
574ccd2543fSMatthew G Knepley     switch (dim) {
575*17fe8556SMatthew G. Knepley     case 1:
576*17fe8556SMatthew G. Knepley       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
577*17fe8556SMatthew G. Knepley       break;
578ccd2543fSMatthew G Knepley     case 2:
579ccd2543fSMatthew G Knepley       switch (coneSize) {
580ccd2543fSMatthew G Knepley       case 3:
581ccd2543fSMatthew G Knepley         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
582ccd2543fSMatthew G Knepley         break;
583ccd2543fSMatthew G Knepley       case 4:
584ccd2543fSMatthew G Knepley         ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
585ccd2543fSMatthew G Knepley         break;
586ccd2543fSMatthew G Knepley       default:
587ccd2543fSMatthew G Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
588ccd2543fSMatthew G Knepley       }
589ccd2543fSMatthew G Knepley       break;
590ccd2543fSMatthew G Knepley     case 3:
591ccd2543fSMatthew G Knepley       switch (coneSize) {
592ccd2543fSMatthew G Knepley       case 4:
593ccd2543fSMatthew G Knepley         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
594ccd2543fSMatthew G Knepley         break;
595ccd2543fSMatthew G Knepley       case 8:
596ccd2543fSMatthew G Knepley         ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
597ccd2543fSMatthew G Knepley         break;
598ccd2543fSMatthew G Knepley       default:
599ccd2543fSMatthew G Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
600ccd2543fSMatthew G Knepley       }
601ccd2543fSMatthew G Knepley       break;
602ccd2543fSMatthew G Knepley     default:
603ccd2543fSMatthew G Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
604ccd2543fSMatthew G Knepley     }
60549dc4407SMatthew G. Knepley   } else {
6065743f1d7SMatthew G. Knepley     /* We need to keep a pointer to the depth label */
6075743f1d7SMatthew G. Knepley     ierr = DMPlexGetLabelValue(dm, "depth", cell, &dim);CHKERRQ(ierr);
60849dc4407SMatthew G. Knepley     /* Cone size is now the number of faces */
60949dc4407SMatthew G. Knepley     switch (dim) {
610*17fe8556SMatthew G. Knepley     case 1:
611*17fe8556SMatthew G. Knepley       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
612*17fe8556SMatthew G. Knepley       break;
61349dc4407SMatthew G. Knepley     case 2:
61449dc4407SMatthew G. Knepley       switch (coneSize) {
61549dc4407SMatthew G. Knepley       case 3:
61649dc4407SMatthew G. Knepley         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
61749dc4407SMatthew G. Knepley         break;
61849dc4407SMatthew G. Knepley       case 4:
61949dc4407SMatthew G. Knepley         ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
62049dc4407SMatthew G. Knepley         break;
62149dc4407SMatthew G. Knepley       default:
62249dc4407SMatthew G. Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
62349dc4407SMatthew G. Knepley       }
62449dc4407SMatthew G. Knepley       break;
62549dc4407SMatthew G. Knepley     case 3:
62649dc4407SMatthew G. Knepley       switch (coneSize) {
62749dc4407SMatthew G. Knepley       case 4:
62849dc4407SMatthew G. Knepley         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
62949dc4407SMatthew G. Knepley         break;
63049dc4407SMatthew G. Knepley       case 6:
63149dc4407SMatthew G. Knepley         ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
63249dc4407SMatthew G. Knepley         break;
63349dc4407SMatthew G. Knepley       default:
63449dc4407SMatthew G. Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
63549dc4407SMatthew G. Knepley       }
63649dc4407SMatthew G. Knepley       break;
63749dc4407SMatthew G. Knepley     default:
63849dc4407SMatthew G. Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
63949dc4407SMatthew G. Knepley     }
64049dc4407SMatthew G. Knepley   }
641ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
642ccd2543fSMatthew G Knepley }
643