xref: /petsc/src/dm/impls/plex/plexgeometry.c (revision 834e62cec43ffa44d1bbf564ad25bbc4bf1153ae)
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__
20517fe8556SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeProjection2Dto1D_Internal"
20617fe8556SMatthew G. Knepley /*
20717fe8556SMatthew G. Knepley   DMPlexComputeProjection2Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 2D
20817fe8556SMatthew G. Knepley */
2097f07f362SMatthew G. Knepley static PetscErrorCode DMPlexComputeProjection2Dto1D_Internal(PetscScalar coords[], PetscReal R[])
21017fe8556SMatthew G. Knepley {
21117fe8556SMatthew G. Knepley   const PetscReal x = PetscRealPart(coords[2] - coords[0]);
21217fe8556SMatthew G. Knepley   const PetscReal y = PetscRealPart(coords[3] - coords[1]);
2137f07f362SMatthew G. Knepley   const PetscReal r = sqrt(x*x + y*y), c = x/r, s = y/r;
21417fe8556SMatthew G. Knepley 
21517fe8556SMatthew G. Knepley   PetscFunctionBegin;
2167f07f362SMatthew G. Knepley   R[0] =  c; R[1] = s;
2177f07f362SMatthew G. Knepley   R[2] = -s; R[3] = c;
21817fe8556SMatthew G. Knepley   coords[0] = 0.0;
2197f07f362SMatthew G. Knepley   coords[1] = r;
22017fe8556SMatthew G. Knepley   PetscFunctionReturn(0);
22117fe8556SMatthew G. Knepley }
22217fe8556SMatthew G. Knepley 
22317fe8556SMatthew G. Knepley #undef __FUNCT__
224ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeProjection3Dto2D_Internal"
225ccd2543fSMatthew G Knepley /*
226ccd2543fSMatthew G Knepley   DMPlexComputeProjection3Dto2D_Internal - Rewrite coordinates to be the 2D projection of the 3D
227ccd2543fSMatthew G Knepley */
2287f07f362SMatthew G. Knepley static PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscScalar coords[], PetscReal R[])
229ccd2543fSMatthew G Knepley {
230ccd2543fSMatthew G Knepley   PetscScalar    x1[3],  x2[3], n[3], norm;
2317f07f362SMatthew G. Knepley   PetscScalar    x1p[3], x2p[3];
2324a217a95SMatthew G. Knepley   PetscReal      sqrtz, alpha;
233ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
234ccd2543fSMatthew G Knepley   PetscInt       d, e;
235ccd2543fSMatthew G Knepley 
236ccd2543fSMatthew G Knepley   PetscFunctionBegin;
237ccd2543fSMatthew G Knepley   /* 0) Calculate normal vector */
238ccd2543fSMatthew G Knepley   for (d = 0; d < dim; ++d) {
239ccd2543fSMatthew G Knepley     x1[d] = coords[1*dim+d] - coords[0*dim+d];
240ccd2543fSMatthew G Knepley     x2[d] = coords[2*dim+d] - coords[0*dim+d];
241ccd2543fSMatthew G Knepley   }
242ccd2543fSMatthew G Knepley   n[0] = x1[1]*x2[2] - x1[2]*x2[1];
243ccd2543fSMatthew G Knepley   n[1] = x1[2]*x2[0] - x1[0]*x2[2];
244ccd2543fSMatthew G Knepley   n[2] = x1[0]*x2[1] - x1[1]*x2[0];
245ccd2543fSMatthew G Knepley   norm = sqrt(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]);
246ccd2543fSMatthew G Knepley   n[0] /= norm;
247ccd2543fSMatthew G Knepley   n[1] /= norm;
248ccd2543fSMatthew G Knepley   n[2] /= norm;
249ccd2543fSMatthew G Knepley   /* 1) Take the normal vector and rotate until it is \hat z
250ccd2543fSMatthew G Knepley 
251ccd2543fSMatthew G Knepley     Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then
252ccd2543fSMatthew G Knepley 
253ccd2543fSMatthew G Knepley     R = /  alpha nx nz  alpha ny nz -1/alpha \
254ccd2543fSMatthew G Knepley         | -alpha ny     alpha nx        0    |
255ccd2543fSMatthew G Knepley         \     nx            ny         nz    /
256ccd2543fSMatthew G Knepley 
257ccd2543fSMatthew G Knepley     will rotate the normal vector to \hat z
258ccd2543fSMatthew G Knepley   */
2594a217a95SMatthew G. Knepley   sqrtz = sqrt(1.0 - PetscAbsScalar(n[2]*n[2]));
26073868372SMatthew G. Knepley   /* Check for n = z */
26173868372SMatthew G. Knepley   if (sqrtz < 1.0e-10) {
26273868372SMatthew G. Knepley     coords[0] = 0.0;
26373868372SMatthew G. Knepley     coords[1] = 0.0;
2644a217a95SMatthew G. Knepley     if (PetscRealPart(n[2]) < 0.0) {
26573868372SMatthew G. Knepley       coords[2] = x2[0];
26673868372SMatthew G. Knepley       coords[3] = x2[1];
26773868372SMatthew G. Knepley       coords[4] = x1[0];
26873868372SMatthew G. Knepley       coords[5] = x1[1];
26973868372SMatthew G. Knepley     } else {
27073868372SMatthew G. Knepley       coords[2] = x1[0];
27173868372SMatthew G. Knepley       coords[3] = x1[1];
27273868372SMatthew G. Knepley       coords[4] = x2[0];
27373868372SMatthew G. Knepley       coords[5] = x2[1];
27473868372SMatthew G. Knepley     }
27573868372SMatthew G. Knepley     PetscFunctionReturn(0);
27673868372SMatthew G. Knepley   }
277da18b5e6SMatthew G Knepley   alpha = 1.0/sqrtz;
278ccd2543fSMatthew G Knepley   R[0] =  alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz;
279ccd2543fSMatthew G Knepley   R[3] = -alpha*n[1];      R[4] = alpha*n[0];      R[5] = 0.0;
280ccd2543fSMatthew G Knepley   R[6] =  n[0];            R[7] = n[1];            R[8] = n[2];
281ccd2543fSMatthew G Knepley   for (d = 0; d < dim; ++d) {
282ccd2543fSMatthew G Knepley     x1p[d] = 0.0;
283ccd2543fSMatthew G Knepley     x2p[d] = 0.0;
284ccd2543fSMatthew G Knepley     for (e = 0; e < dim; ++e) {
285ccd2543fSMatthew G Knepley       x1p[d] += R[d*dim+e]*x1[e];
286ccd2543fSMatthew G Knepley       x2p[d] += R[d*dim+e]*x2[e];
287ccd2543fSMatthew G Knepley     }
288ccd2543fSMatthew G Knepley   }
28988790e04SMatthew G Knepley   if (PetscAbsScalar(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
29088790e04SMatthew G Knepley   if (PetscAbsScalar(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
291ccd2543fSMatthew G Knepley   /* 2) Project to (x, y) */
292ccd2543fSMatthew G Knepley   coords[0] = 0.0;
293ccd2543fSMatthew G Knepley   coords[1] = 0.0;
294ccd2543fSMatthew G Knepley   coords[2] = x1p[0];
295ccd2543fSMatthew G Knepley   coords[3] = x1p[1];
296ccd2543fSMatthew G Knepley   coords[4] = x2p[0];
297ccd2543fSMatthew G Knepley   coords[5] = x2p[1];
2987f07f362SMatthew G. Knepley   /* Output R^T which rotates \hat z to the input normal */
2997f07f362SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
3007f07f362SMatthew G. Knepley     for (e = d+1; e < dim; ++e) {
3017f07f362SMatthew G. Knepley       PetscReal tmp;
3027f07f362SMatthew G. Knepley 
3037f07f362SMatthew G. Knepley       tmp        = R[d*dim+e];
3047f07f362SMatthew G. Knepley       R[d*dim+e] = R[e*dim+d];
3057f07f362SMatthew G. Knepley       R[e*dim+d] = tmp;
3067f07f362SMatthew G. Knepley     }
3077f07f362SMatthew G. Knepley   }
308ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
309ccd2543fSMatthew G Knepley }
310ccd2543fSMatthew G Knepley 
311ccd2543fSMatthew G Knepley #undef __FUNCT__
3127f07f362SMatthew G. Knepley #define __FUNCT__ "Invert2D_Internal"
3137f07f362SMatthew G. Knepley PETSC_STATIC_INLINE void Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
3147f07f362SMatthew G. Knepley {
3157f07f362SMatthew G. Knepley   const PetscReal invDet = 1.0/detJ;
3167f07f362SMatthew G. Knepley 
3177f07f362SMatthew G. Knepley   invJ[0] =  invDet*J[3];
3187f07f362SMatthew G. Knepley   invJ[1] = -invDet*J[1];
3197f07f362SMatthew G. Knepley   invJ[2] = -invDet*J[2];
3207f07f362SMatthew G. Knepley   invJ[3] =  invDet*J[0];
3217f07f362SMatthew G. Knepley   PetscLogFlops(5.0);
3227f07f362SMatthew G. Knepley }
3237f07f362SMatthew G. Knepley 
3247f07f362SMatthew G. Knepley #undef __FUNCT__
3257f07f362SMatthew G. Knepley #define __FUNCT__ "Invert3D_Internal"
3267f07f362SMatthew G. Knepley PETSC_STATIC_INLINE void Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
3277f07f362SMatthew G. Knepley {
3287f07f362SMatthew G. Knepley   const PetscReal invDet = 1.0/detJ;
3297f07f362SMatthew G. Knepley 
3307f07f362SMatthew G. Knepley   invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
3317f07f362SMatthew G. Knepley   invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
3327f07f362SMatthew G. Knepley   invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
3337f07f362SMatthew G. Knepley   invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
3347f07f362SMatthew G. Knepley   invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
3357f07f362SMatthew G. Knepley   invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
3367f07f362SMatthew G. Knepley   invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
3377f07f362SMatthew G. Knepley   invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
3387f07f362SMatthew G. Knepley   invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
3397f07f362SMatthew G. Knepley   PetscLogFlops(37.0);
3407f07f362SMatthew G. Knepley }
3417f07f362SMatthew G. Knepley 
3427f07f362SMatthew G. Knepley #undef __FUNCT__
3437f07f362SMatthew G. Knepley #define __FUNCT__ "Det2D_Internal"
3447f07f362SMatthew G. Knepley PETSC_STATIC_INLINE void Det2D_Internal(PetscReal *detJ, PetscReal J[])
3457f07f362SMatthew G. Knepley {
3467f07f362SMatthew G. Knepley   *detJ = J[0]*J[3] - J[1]*J[2];
3477f07f362SMatthew G. Knepley   PetscLogFlops(3.0);
3487f07f362SMatthew G. Knepley }
3497f07f362SMatthew G. Knepley 
3507f07f362SMatthew G. Knepley #undef __FUNCT__
3517f07f362SMatthew G. Knepley #define __FUNCT__ "Det3D_Internal"
3527f07f362SMatthew G. Knepley PETSC_STATIC_INLINE void Det3D_Internal(PetscReal *detJ, PetscReal J[])
3537f07f362SMatthew G. Knepley {
3547f07f362SMatthew G. Knepley   *detJ = (J[0*3+0]*(J[1*3+2]*J[2*3+1] - J[1*3+1]*J[2*3+2]) +
3557f07f362SMatthew G. Knepley            J[0*3+1]*(J[1*3+0]*J[2*3+2] - J[1*3+2]*J[2*3+0]) +
3567f07f362SMatthew G. Knepley            J[0*3+2]*(J[1*3+1]*J[2*3+0] - J[1*3+0]*J[2*3+1]));
3577f07f362SMatthew G. Knepley   PetscLogFlops(12.0);
3587f07f362SMatthew G. Knepley }
3597f07f362SMatthew G. Knepley 
3607f07f362SMatthew G. Knepley #undef __FUNCT__
361*834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Triangle_Internal"
362*834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[])
363*834e62ceSMatthew G. Knepley {
364*834e62ceSMatthew G. Knepley   /* Signed volume is 1/2 the determinant
365*834e62ceSMatthew G. Knepley 
366*834e62ceSMatthew G. Knepley    |  1  1  1 |
367*834e62ceSMatthew G. Knepley    | x0 x1 x2 |
368*834e62ceSMatthew G. Knepley    | y0 y1 y2 |
369*834e62ceSMatthew G. Knepley 
370*834e62ceSMatthew G. Knepley      but if x0,y0 is the origin, we have
371*834e62ceSMatthew G. Knepley 
372*834e62ceSMatthew G. Knepley    | x1 x2 |
373*834e62ceSMatthew G. Knepley    | y1 y2 |
374*834e62ceSMatthew G. Knepley   */
375*834e62ceSMatthew G. Knepley   const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1];
376*834e62ceSMatthew G. Knepley   const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1];
377*834e62ceSMatthew G. Knepley   PetscReal       M[4], detM;
378*834e62ceSMatthew G. Knepley   M[0] = x1; M[1] = x2;
379*834e62ceSMatthew G. Knepley   M[3] = y1; M[4] = y2;
380*834e62ceSMatthew G. Knepley   Det2D_Internal(&detM, M);
381*834e62ceSMatthew G. Knepley   *vol = 0.5*detM;
382*834e62ceSMatthew G. Knepley   PetscLogFlops(5.0);
383*834e62ceSMatthew G. Knepley }
384*834e62ceSMatthew G. Knepley 
385*834e62ceSMatthew G. Knepley #undef __FUNCT__
386*834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Triangle_Origin_Internal"
387*834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[])
388*834e62ceSMatthew G. Knepley {
389*834e62ceSMatthew G. Knepley   Det2D_Internal(vol, coords);
390*834e62ceSMatthew G. Knepley   *vol *= 0.5;
391*834e62ceSMatthew G. Knepley }
392*834e62ceSMatthew G. Knepley 
393*834e62ceSMatthew G. Knepley #undef __FUNCT__
394*834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Tetrahedron_Internal"
395*834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[])
396*834e62ceSMatthew G. Knepley {
397*834e62ceSMatthew G. Knepley   /* Signed volume is 1/6th of the determinant
398*834e62ceSMatthew G. Knepley 
399*834e62ceSMatthew G. Knepley    |  1  1  1  1 |
400*834e62ceSMatthew G. Knepley    | x0 x1 x2 x3 |
401*834e62ceSMatthew G. Knepley    | y0 y1 y2 y3 |
402*834e62ceSMatthew G. Knepley    | z0 z1 z2 z3 |
403*834e62ceSMatthew G. Knepley 
404*834e62ceSMatthew G. Knepley      but if x0,y0,z0 is the origin, we have
405*834e62ceSMatthew G. Knepley 
406*834e62ceSMatthew G. Knepley    | x1 x2 x3 |
407*834e62ceSMatthew G. Knepley    | y1 y2 y3 |
408*834e62ceSMatthew G. Knepley    | z1 z2 z3 |
409*834e62ceSMatthew G. Knepley   */
410*834e62ceSMatthew G. Knepley   const PetscReal x1 = coords[3] - coords[0], y1 = coords[4]  - coords[1], z1 = coords[5]  - coords[2];
411*834e62ceSMatthew G. Knepley   const PetscReal x2 = coords[6] - coords[0], y2 = coords[7]  - coords[1], z2 = coords[8]  - coords[2];
412*834e62ceSMatthew G. Knepley   const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2];
413*834e62ceSMatthew G. Knepley   PetscReal       M[9], detM;
414*834e62ceSMatthew G. Knepley   M[0] = x1; M[1] = x2; M[2] = x3;
415*834e62ceSMatthew G. Knepley   M[3] = y1; M[4] = y2; M[5] = y3;
416*834e62ceSMatthew G. Knepley   M[6] = z1; M[7] = z2; M[8] = z3;
417*834e62ceSMatthew G. Knepley   Det3D_Internal(&detM, M);
418*834e62ceSMatthew G. Knepley   *vol = 0.16666666666666666666666*detM;
419*834e62ceSMatthew G. Knepley   PetscLogFlops(10.0);
420*834e62ceSMatthew G. Knepley }
421*834e62ceSMatthew G. Knepley 
422*834e62ceSMatthew G. Knepley #undef __FUNCT__
42317fe8556SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeLineGeometry_Internal"
42417fe8556SMatthew G. Knepley static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
42517fe8556SMatthew G. Knepley {
42617fe8556SMatthew G. Knepley   PetscSection   coordSection;
42717fe8556SMatthew G. Knepley   Vec            coordinates;
42817fe8556SMatthew G. Knepley   PetscScalar   *coords;
4297f07f362SMatthew G. Knepley   PetscInt       numCoords, d;
43017fe8556SMatthew G. Knepley   PetscErrorCode ierr;
43117fe8556SMatthew G. Knepley 
43217fe8556SMatthew G. Knepley   PetscFunctionBegin;
43317fe8556SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
43417fe8556SMatthew G. Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
43517fe8556SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
4367f07f362SMatthew G. Knepley   *detJ = 0.0;
43717fe8556SMatthew G. Knepley   if (numCoords == 4) {
4387f07f362SMatthew G. Knepley     const PetscInt dim = 2;
4397f07f362SMatthew G. Knepley     PetscReal      R[4], J0;
4407f07f362SMatthew G. Knepley 
4417f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
4427f07f362SMatthew G. Knepley     ierr = DMPlexComputeProjection2Dto1D_Internal(coords, R);CHKERRQ(ierr);
44317fe8556SMatthew G. Knepley     if (J)    {
4447f07f362SMatthew G. Knepley       J0   = 0.5*PetscRealPart(coords[1]);
4457f07f362SMatthew G. Knepley       J[0] = R[0]*J0; J[1] = R[1];
4467f07f362SMatthew G. Knepley       J[2] = R[2]*J0; J[3] = R[3];
4477f07f362SMatthew G. Knepley       Det2D_Internal(detJ, J);
44817fe8556SMatthew G. Knepley     }
4497f07f362SMatthew G. Knepley     if (invJ) {Invert2D_Internal(invJ, J, *detJ);}
4507f07f362SMatthew G. Knepley   } else if (numCoords == 2) {
4517f07f362SMatthew G. Knepley     const PetscInt dim = 1;
4527f07f362SMatthew G. Knepley 
4537f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
4547f07f362SMatthew G. Knepley     if (J)    {
4557f07f362SMatthew G. Knepley       J[0]  = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0]));
45617fe8556SMatthew G. Knepley       *detJ = J[0];
45717fe8556SMatthew G. Knepley       PetscLogFlops(2.0);
45817fe8556SMatthew G. Knepley     }
4597f07f362SMatthew G. Knepley     if (invJ) {invJ[0] = 1.0/J[0]; PetscLogFlops(1.0);}
4607f07f362SMatthew G. Knepley   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %d != 2", numCoords);
46117fe8556SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
46217fe8556SMatthew G. Knepley   PetscFunctionReturn(0);
46317fe8556SMatthew G. Knepley }
46417fe8556SMatthew G. Knepley 
46517fe8556SMatthew G. Knepley #undef __FUNCT__
466ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal"
467ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
468ccd2543fSMatthew G Knepley {
469ccd2543fSMatthew G Knepley   PetscSection   coordSection;
470ccd2543fSMatthew G Knepley   Vec            coordinates;
4717c1f9639SMatthew G Knepley   PetscScalar   *coords;
4727f07f362SMatthew G. Knepley   PetscInt       numCoords, d, f, g;
473ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
474ccd2543fSMatthew G Knepley 
475ccd2543fSMatthew G Knepley   PetscFunctionBegin;
476ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
477ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
478ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
4797f07f362SMatthew G. Knepley   *detJ = 0.0;
480ccd2543fSMatthew G Knepley   if (numCoords == 9) {
4817f07f362SMatthew G. Knepley     const PetscInt dim = 3;
4827f07f362SMatthew G. Knepley     PetscReal      R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0};
4837f07f362SMatthew G. Knepley 
4847f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
4857f07f362SMatthew G. Knepley     ierr = DMPlexComputeProjection3Dto2D_Internal(coords, R);CHKERRQ(ierr);
4867f07f362SMatthew G. Knepley     if (J)    {
4877f07f362SMatthew G. Knepley       for (d = 0; d < dim-1; d++) {
4887f07f362SMatthew G. Knepley         for (f = 0; f < dim-1; f++) {
4897f07f362SMatthew G. Knepley           J0[d*(dim+1)+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
490ccd2543fSMatthew G Knepley         }
4917f07f362SMatthew G. Knepley       }
4927f07f362SMatthew G. Knepley       PetscLogFlops(8.0);
4937f07f362SMatthew G. Knepley       for (d = 0; d < dim; d++) {
4947f07f362SMatthew G. Knepley         for (f = 0; f < dim; f++) {
4957f07f362SMatthew G. Knepley           J[d*dim+f] = 0.0;
4967f07f362SMatthew G. Knepley           for (g = 0; g < dim; g++) {
4977f07f362SMatthew G. Knepley             J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
4987f07f362SMatthew G. Knepley           }
4997f07f362SMatthew G. Knepley         }
5007f07f362SMatthew G. Knepley       }
5017f07f362SMatthew G. Knepley       PetscLogFlops(18.0);
5027f07f362SMatthew G. Knepley       Det3D_Internal(detJ, J);
5037f07f362SMatthew G. Knepley     }
5047f07f362SMatthew G. Knepley     if (invJ) {Invert3D_Internal(invJ, J, *detJ);}
5057f07f362SMatthew G. Knepley   } else if (numCoords == 6) {
5067f07f362SMatthew G. Knepley     const PetscInt dim = 2;
5077f07f362SMatthew G. Knepley 
5087f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
509ccd2543fSMatthew G Knepley     if (J)    {
510ccd2543fSMatthew G Knepley       for (d = 0; d < dim; d++) {
511ccd2543fSMatthew G Knepley         for (f = 0; f < dim; f++) {
512ccd2543fSMatthew G Knepley           J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
513ccd2543fSMatthew G Knepley         }
514ccd2543fSMatthew G Knepley       }
5157f07f362SMatthew G. Knepley       PetscLogFlops(8.0);
5167f07f362SMatthew G. Knepley       Det2D_Internal(detJ, J);
517ccd2543fSMatthew G Knepley     }
5187f07f362SMatthew G. Knepley     if (invJ) {Invert2D_Internal(invJ, J, *detJ);}
5197f07f362SMatthew G. Knepley   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %d != 6", numCoords);
520ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
521ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
522ccd2543fSMatthew G Knepley }
523ccd2543fSMatthew G Knepley 
524ccd2543fSMatthew G Knepley #undef __FUNCT__
525ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal"
526ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
527ccd2543fSMatthew G Knepley {
528ccd2543fSMatthew G Knepley   PetscSection   coordSection;
529ccd2543fSMatthew G Knepley   Vec            coordinates;
5307c1f9639SMatthew G Knepley   PetscScalar   *coords;
531ccd2543fSMatthew G Knepley   const PetscInt dim = 2;
532ccd2543fSMatthew G Knepley   PetscInt       d, f;
533ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
534ccd2543fSMatthew G Knepley 
535ccd2543fSMatthew G Knepley   PetscFunctionBegin;
536ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
537ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
538ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
5397f07f362SMatthew G. Knepley   *detJ = 0.0;
5407f07f362SMatthew G. Knepley   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
541ccd2543fSMatthew G Knepley   if (J)    {
542ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
543ccd2543fSMatthew G Knepley       for (f = 0; f < dim; f++) {
544ccd2543fSMatthew G Knepley         J[d*dim+f] = 0.5*(PetscRealPart(coords[(f*2+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
545ccd2543fSMatthew G Knepley       }
546ccd2543fSMatthew G Knepley     }
5477f07f362SMatthew G. Knepley     PetscLogFlops(8.0);
5487f07f362SMatthew G. Knepley     Det2D_Internal(detJ, J);
549ccd2543fSMatthew G Knepley   }
5507f07f362SMatthew G. Knepley   if (invJ) {Invert2D_Internal(invJ, J, *detJ);}
551ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
552ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
553ccd2543fSMatthew G Knepley }
554ccd2543fSMatthew G Knepley 
555ccd2543fSMatthew G Knepley #undef __FUNCT__
556ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal"
557ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
558ccd2543fSMatthew G Knepley {
559ccd2543fSMatthew G Knepley   PetscSection   coordSection;
560ccd2543fSMatthew G Knepley   Vec            coordinates;
5617c1f9639SMatthew G Knepley   PetscScalar   *coords;
562ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
563ccd2543fSMatthew G Knepley   PetscInt       d, f;
564ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
565ccd2543fSMatthew G Knepley 
566ccd2543fSMatthew G Knepley   PetscFunctionBegin;
567ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
568ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
569ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
5707f07f362SMatthew G. Knepley   *detJ = 0.0;
5717f07f362SMatthew G. Knepley   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
572ccd2543fSMatthew G Knepley   if (J)    {
573ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
574ccd2543fSMatthew G Knepley       for (f = 0; f < dim; f++) {
575ccd2543fSMatthew G Knepley         J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
576ccd2543fSMatthew G Knepley       }
577ccd2543fSMatthew G Knepley     }
5787f07f362SMatthew G. Knepley     PetscLogFlops(18.0);
5797f07f362SMatthew G. Knepley     Det3D_Internal(detJ, J);
580ccd2543fSMatthew G Knepley   }
5817f07f362SMatthew G. Knepley   if (invJ) {Invert3D_Internal(invJ, J, *detJ);}
582ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
583ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
584ccd2543fSMatthew G Knepley }
585ccd2543fSMatthew G Knepley 
586ccd2543fSMatthew G Knepley #undef __FUNCT__
587ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal"
588ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
589ccd2543fSMatthew G Knepley {
590ccd2543fSMatthew G Knepley   PetscSection   coordSection;
591ccd2543fSMatthew G Knepley   Vec            coordinates;
5927c1f9639SMatthew G Knepley   PetscScalar   *coords;
593ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
594ccd2543fSMatthew G Knepley   PetscInt       d;
595ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
596ccd2543fSMatthew G Knepley 
597ccd2543fSMatthew G Knepley   PetscFunctionBegin;
598ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
599ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
600ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
6017f07f362SMatthew G. Knepley   *detJ = 0.0;
6027f07f362SMatthew G. Knepley   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
603ccd2543fSMatthew G Knepley   if (J)    {
604ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
605de36ddddSMatthew G. Knepley       J[d*dim+0] = 0.5*(PetscRealPart(coords[(2+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
606ccd2543fSMatthew G Knepley       J[d*dim+1] = 0.5*(PetscRealPart(coords[(1+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
607ccd2543fSMatthew G Knepley       J[d*dim+2] = 0.5*(PetscRealPart(coords[(3+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
608ccd2543fSMatthew G Knepley     }
6097f07f362SMatthew G. Knepley     PetscLogFlops(18.0);
6107f07f362SMatthew G. Knepley     Det3D_Internal(detJ, J);
611ccd2543fSMatthew G Knepley   }
6127f07f362SMatthew G. Knepley   if (invJ) {Invert3D_Internal(invJ, J, *detJ);}
6137f07f362SMatthew G. Knepley   *detJ *= -8.0;
614ccd2543fSMatthew G Knepley   ierr   = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
615ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
616ccd2543fSMatthew G Knepley }
617ccd2543fSMatthew G Knepley 
618ccd2543fSMatthew G Knepley #undef __FUNCT__
619ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeCellGeometry"
620ccd2543fSMatthew G Knepley /*@C
621ccd2543fSMatthew G Knepley   DMPlexComputeCellGeometry - Compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell
622ccd2543fSMatthew G Knepley 
623ccd2543fSMatthew G Knepley   Collective on DM
624ccd2543fSMatthew G Knepley 
625ccd2543fSMatthew G Knepley   Input Arguments:
626ccd2543fSMatthew G Knepley + dm   - the DM
627ccd2543fSMatthew G Knepley - cell - the cell
628ccd2543fSMatthew G Knepley 
629ccd2543fSMatthew G Knepley   Output Arguments:
630ccd2543fSMatthew G Knepley + v0   - the translation part of this affine transform
631ccd2543fSMatthew G Knepley . J    - the Jacobian of the transform from the reference element
632ccd2543fSMatthew G Knepley . invJ - the inverse of the Jacobian
633ccd2543fSMatthew G Knepley - detJ - the Jacobian determinant
634ccd2543fSMatthew G Knepley 
635ccd2543fSMatthew G Knepley   Level: advanced
636ccd2543fSMatthew G Knepley 
637ccd2543fSMatthew G Knepley   Fortran Notes:
638ccd2543fSMatthew G Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
639ccd2543fSMatthew G Knepley   include petsc.h90 in your code.
640ccd2543fSMatthew G Knepley 
641ccd2543fSMatthew G Knepley .seealso: DMPlexGetCoordinateSection(), DMPlexGetCoordinateVec()
642ccd2543fSMatthew G Knepley @*/
643ccd2543fSMatthew G Knepley PetscErrorCode DMPlexComputeCellGeometry(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
644ccd2543fSMatthew G Knepley {
64549dc4407SMatthew G. Knepley   PetscInt       depth, dim, coneSize;
646ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
647ccd2543fSMatthew G Knepley 
648ccd2543fSMatthew G Knepley   PetscFunctionBegin;
649139a35ccSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
650ccd2543fSMatthew G Knepley   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
65149dc4407SMatthew G. Knepley   if (depth == 1) {
6525743f1d7SMatthew G. Knepley     ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
653ccd2543fSMatthew G Knepley     switch (dim) {
65417fe8556SMatthew G. Knepley     case 1:
65517fe8556SMatthew G. Knepley       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
65617fe8556SMatthew G. Knepley       break;
657ccd2543fSMatthew G Knepley     case 2:
658ccd2543fSMatthew G Knepley       switch (coneSize) {
659ccd2543fSMatthew G Knepley       case 3:
660ccd2543fSMatthew G Knepley         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
661ccd2543fSMatthew G Knepley         break;
662ccd2543fSMatthew G Knepley       case 4:
663ccd2543fSMatthew G Knepley         ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
664ccd2543fSMatthew G Knepley         break;
665ccd2543fSMatthew G Knepley       default:
666ccd2543fSMatthew G Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
667ccd2543fSMatthew G Knepley       }
668ccd2543fSMatthew G Knepley       break;
669ccd2543fSMatthew G Knepley     case 3:
670ccd2543fSMatthew G Knepley       switch (coneSize) {
671ccd2543fSMatthew G Knepley       case 4:
672ccd2543fSMatthew G Knepley         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
673ccd2543fSMatthew G Knepley         break;
674ccd2543fSMatthew G Knepley       case 8:
675ccd2543fSMatthew G Knepley         ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
676ccd2543fSMatthew G Knepley         break;
677ccd2543fSMatthew G Knepley       default:
678ccd2543fSMatthew G Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
679ccd2543fSMatthew G Knepley       }
680ccd2543fSMatthew G Knepley       break;
681ccd2543fSMatthew G Knepley     default:
682ccd2543fSMatthew G Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
683ccd2543fSMatthew G Knepley     }
68449dc4407SMatthew G. Knepley   } else {
6855743f1d7SMatthew G. Knepley     /* We need to keep a pointer to the depth label */
6865743f1d7SMatthew G. Knepley     ierr = DMPlexGetLabelValue(dm, "depth", cell, &dim);CHKERRQ(ierr);
68749dc4407SMatthew G. Knepley     /* Cone size is now the number of faces */
68849dc4407SMatthew G. Knepley     switch (dim) {
68917fe8556SMatthew G. Knepley     case 1:
69017fe8556SMatthew G. Knepley       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
69117fe8556SMatthew G. Knepley       break;
69249dc4407SMatthew G. Knepley     case 2:
69349dc4407SMatthew G. Knepley       switch (coneSize) {
69449dc4407SMatthew G. Knepley       case 3:
69549dc4407SMatthew G. Knepley         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
69649dc4407SMatthew G. Knepley         break;
69749dc4407SMatthew G. Knepley       case 4:
69849dc4407SMatthew G. Knepley         ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
69949dc4407SMatthew G. Knepley         break;
70049dc4407SMatthew G. Knepley       default:
70149dc4407SMatthew G. Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
70249dc4407SMatthew G. Knepley       }
70349dc4407SMatthew G. Knepley       break;
70449dc4407SMatthew G. Knepley     case 3:
70549dc4407SMatthew G. Knepley       switch (coneSize) {
70649dc4407SMatthew G. Knepley       case 4:
70749dc4407SMatthew G. Knepley         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
70849dc4407SMatthew G. Knepley         break;
70949dc4407SMatthew G. Knepley       case 6:
71049dc4407SMatthew G. Knepley         ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
71149dc4407SMatthew G. Knepley         break;
71249dc4407SMatthew G. Knepley       default:
71349dc4407SMatthew G. Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
71449dc4407SMatthew G. Knepley       }
71549dc4407SMatthew G. Knepley       break;
71649dc4407SMatthew G. Knepley     default:
71749dc4407SMatthew G. Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
71849dc4407SMatthew G. Knepley     }
71949dc4407SMatthew G. Knepley   }
720ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
721ccd2543fSMatthew G Knepley }
722*834e62ceSMatthew G. Knepley 
723*834e62ceSMatthew G. Knepley #undef __FUNCT__
724*834e62ceSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeCellGeometryFVM"
725*834e62ceSMatthew G. Knepley /*@C
726*834e62ceSMatthew G. Knepley   DMPlexComputeCellGeometryFVM - Compute the volume for a given cell
727*834e62ceSMatthew G. Knepley 
728*834e62ceSMatthew G. Knepley   Collective on DM
729*834e62ceSMatthew G. Knepley 
730*834e62ceSMatthew G. Knepley   Input Arguments:
731*834e62ceSMatthew G. Knepley + dm   - the DM
732*834e62ceSMatthew G. Knepley - cell - the cell
733*834e62ceSMatthew G. Knepley 
734*834e62ceSMatthew G. Knepley   Output Arguments:
735*834e62ceSMatthew G. Knepley + volume   - the cell volume
736*834e62ceSMatthew G. Knepley - centroid - the cell centroid
737*834e62ceSMatthew G. Knepley 
738*834e62ceSMatthew G. Knepley   Level: advanced
739*834e62ceSMatthew G. Knepley 
740*834e62ceSMatthew G. Knepley   Fortran Notes:
741*834e62ceSMatthew G. Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
742*834e62ceSMatthew G. Knepley   include petsc.h90 in your code.
743*834e62ceSMatthew G. Knepley 
744*834e62ceSMatthew G. Knepley .seealso: DMPlexGetCoordinateSection(), DMPlexGetCoordinateVec()
745*834e62ceSMatthew G. Knepley @*/
746*834e62ceSMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[])
747*834e62ceSMatthew G. Knepley {
748*834e62ceSMatthew G. Knepley   PetscInt       depth, dim, coneSize;
749*834e62ceSMatthew G. Knepley   PetscErrorCode ierr;
750*834e62ceSMatthew G. Knepley 
751*834e62ceSMatthew G. Knepley   PetscFunctionBegin;
752*834e62ceSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
753*834e62ceSMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
754*834e62ceSMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
755*834e62ceSMatthew G. Knepley   if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated");
756*834e62ceSMatthew G. Knepley   /* We need to keep a pointer to the depth label */
757*834e62ceSMatthew G. Knepley   ierr = DMPlexGetLabelValue(dm, "depth", cell, &dim);CHKERRQ(ierr);
758*834e62ceSMatthew G. Knepley   /* Cone size is now the number of faces */
759*834e62ceSMatthew G. Knepley   switch (dim) {
760*834e62ceSMatthew G. Knepley   case 2:
761*834e62ceSMatthew G. Knepley   {
762*834e62ceSMatthew G. Knepley     /*
763*834e62ceSMatthew G. Knepley      Centroid_i = (\sum_n A_n Cn_i ) / A
764*834e62ceSMatthew G. Knepley     */
765*834e62ceSMatthew G. Knepley     PetscSection coordSection;
766*834e62ceSMatthew G. Knepley     Vec          coordinates;
767*834e62ceSMatthew G. Knepley     PetscScalar *coords = NULL;
768*834e62ceSMatthew G. Knepley     PetscInt     coordSize, numCorners, p, d;
769*834e62ceSMatthew G. Knepley     PetscReal    vsum, vtmp, ctmp[4];
770*834e62ceSMatthew G. Knepley 
771*834e62ceSMatthew G. Knepley     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
772*834e62ceSMatthew G. Knepley     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
773*834e62ceSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
774*834e62ceSMatthew G. Knepley     numCorners = coordSize/dim;
775*834e62ceSMatthew G. Knepley     for (p = 0; p < numCorners-1; ++p) {
776*834e62ceSMatthew G. Knepley       Volume_Triangle_Origin_Internal(&vtmp, &coords[p*dim]);
777*834e62ceSMatthew G. Knepley       vsum += vtmp;
778*834e62ceSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
779*834e62ceSMatthew G. Knepley         centroid[d] += (coords[p*dim+d] + coords[(p+1)*dim+d])*vtmp;
780*834e62ceSMatthew G. Knepley       }
781*834e62ceSMatthew G. Knepley     }
782*834e62ceSMatthew G. Knepley     for (d = 0; d < dim; ++d) {
783*834e62ceSMatthew G. Knepley       ctmp[d]     = coords[(numCorners-1)*dim+d];
784*834e62ceSMatthew G. Knepley       ctmp[dim+d] = coords[d];
785*834e62ceSMatthew G. Knepley     }
786*834e62ceSMatthew G. Knepley     Volume_Triangle_Origin_Internal(&vtmp, ctmp);
787*834e62ceSMatthew G. Knepley     vsum += vtmp;
788*834e62ceSMatthew G. Knepley     for (d = 0; d < dim; ++d) {
789*834e62ceSMatthew G. Knepley       centroid[d] += (coords[(numCorners-1)*dim+d] + coords[0*dim+d])*vtmp;
790*834e62ceSMatthew G. Knepley       centroid[d] /= (dim+1)*vsum;
791*834e62ceSMatthew G. Knepley     }
792*834e62ceSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
793*834e62ceSMatthew G. Knepley     *vol = PetscAbsScalar(vsum);
794*834e62ceSMatthew G. Knepley   }
795*834e62ceSMatthew G. Knepley   break;
796*834e62ceSMatthew G. Knepley   case 3:
797*834e62ceSMatthew G. Knepley     switch (coneSize) {
798*834e62ceSMatthew G. Knepley     default:
799*834e62ceSMatthew G. Knepley       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
800*834e62ceSMatthew G. Knepley     }
801*834e62ceSMatthew G. Knepley     break;
802*834e62ceSMatthew G. Knepley   default:
803*834e62ceSMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
804*834e62ceSMatthew G. Knepley   }
805*834e62ceSMatthew G. Knepley   PetscFunctionReturn(0);
806*834e62ceSMatthew G. Knepley }
807