xref: /petsc/src/dm/impls/plex/plexgeometry.c (revision 97f1a2182ca1afe151ff8b33ee7644d78f3891d9)
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;
30a1e44745SMatthew G. Knepley   PetscScalar    *coords = NULL;
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);
3969d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(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;
88fb150da6SMatthew G. Knepley   const PetscInt faces[24] = {0, 3, 2, 1,  5, 4, 7, 6,  3, 0, 4, 5,
89fb150da6SMatthew G. Knepley                               1, 2, 6, 7,  3, 5, 6, 2,  0, 1, 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);
9669d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(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;
157785e854fSJed Brown   ierr       = PetscMalloc1(numPoints, &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;
18535953a02SMatthew G. Knepley         case 6:
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]);
2138b49ba18SBarry Smith   const PetscReal r = PetscSqrtReal(x*x + y*y), c = x/r, s = y/r;
21417fe8556SMatthew G. Knepley 
21517fe8556SMatthew G. Knepley   PetscFunctionBegin;
2161c99cf0cSGeoffrey Irving   R[0] = c; R[1] = -s;
2171c99cf0cSGeoffrey Irving   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 */
22899dec3a6SMatthew G. Knepley static PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscInt coordSize, PetscScalar coords[], PetscReal R[])
229ccd2543fSMatthew G Knepley {
2301ee9d5ecSMatthew G. Knepley   PetscReal      x1[3],  x2[3], n[3], norm;
23199dec3a6SMatthew G. Knepley   PetscReal      x1p[3], x2p[3], xnp[3];
2324a217a95SMatthew G. Knepley   PetscReal      sqrtz, alpha;
233ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
23499dec3a6SMatthew G. Knepley   PetscInt       d, e, p;
235ccd2543fSMatthew G Knepley 
236ccd2543fSMatthew G Knepley   PetscFunctionBegin;
237ccd2543fSMatthew G Knepley   /* 0) Calculate normal vector */
238ccd2543fSMatthew G Knepley   for (d = 0; d < dim; ++d) {
2391ee9d5ecSMatthew G. Knepley     x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]);
2401ee9d5ecSMatthew G. Knepley     x2[d] = PetscRealPart(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];
2458b49ba18SBarry Smith   norm = PetscSqrtReal(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   */
2598b49ba18SBarry Smith   sqrtz = PetscSqrtReal(1.0 - n[2]*n[2]);
26073868372SMatthew G. Knepley   /* Check for n = z */
26173868372SMatthew G. Knepley   if (sqrtz < 1.0e-10) {
2621ee9d5ecSMatthew G. Knepley     if (n[2] < 0.0) {
26399dec3a6SMatthew G. Knepley       if (coordSize > 9) {
26499dec3a6SMatthew G. Knepley         coords[2] = PetscRealPart(coords[3*dim+0] - coords[0*dim+0]);
26599dec3a6SMatthew G. Knepley         coords[3] = PetscRealPart(coords[3*dim+0] - coords[0*dim+0]);
26699dec3a6SMatthew G. Knepley         coords[4] = x2[0];
26799dec3a6SMatthew G. Knepley         coords[5] = x2[1];
26899dec3a6SMatthew G. Knepley         coords[6] = x1[0];
26999dec3a6SMatthew G. Knepley         coords[7] = x1[1];
27099dec3a6SMatthew G. Knepley       } else {
27173868372SMatthew G. Knepley         coords[2] = x2[0];
27273868372SMatthew G. Knepley         coords[3] = x2[1];
27373868372SMatthew G. Knepley         coords[4] = x1[0];
27473868372SMatthew G. Knepley         coords[5] = x1[1];
27599dec3a6SMatthew G. Knepley       }
276b7ad821dSMatthew G. Knepley       R[0] = 1.0; R[1] = 0.0; R[2] = 0.0;
277b7ad821dSMatthew G. Knepley       R[3] = 0.0; R[4] = 1.0; R[5] = 0.0;
278b7ad821dSMatthew G. Knepley       R[6] = 0.0; R[7] = 0.0; R[8] = -1.0;
27973868372SMatthew G. Knepley     } else {
28099dec3a6SMatthew G. Knepley       for (p = 3; p < coordSize/3; ++p) {
28199dec3a6SMatthew G. Knepley         coords[p*2+0] = PetscRealPart(coords[p*dim+0] - coords[0*dim+0]);
28299dec3a6SMatthew G. Knepley         coords[p*2+1] = PetscRealPart(coords[p*dim+1] - coords[0*dim+1]);
28399dec3a6SMatthew G. Knepley       }
28473868372SMatthew G. Knepley       coords[2] = x1[0];
28573868372SMatthew G. Knepley       coords[3] = x1[1];
28673868372SMatthew G. Knepley       coords[4] = x2[0];
28773868372SMatthew G. Knepley       coords[5] = x2[1];
288b7ad821dSMatthew G. Knepley       R[0] = 1.0; R[1] = 0.0; R[2] = 0.0;
289b7ad821dSMatthew G. Knepley       R[3] = 0.0; R[4] = 1.0; R[5] = 0.0;
290b7ad821dSMatthew G. Knepley       R[6] = 0.0; R[7] = 0.0; R[8] = 1.0;
29173868372SMatthew G. Knepley     }
29299dec3a6SMatthew G. Knepley     coords[0] = 0.0;
29399dec3a6SMatthew G. Knepley     coords[1] = 0.0;
29473868372SMatthew G. Knepley     PetscFunctionReturn(0);
29573868372SMatthew G. Knepley   }
296da18b5e6SMatthew G Knepley   alpha = 1.0/sqrtz;
297ccd2543fSMatthew G Knepley   R[0] =  alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz;
298ccd2543fSMatthew G Knepley   R[3] = -alpha*n[1];      R[4] = alpha*n[0];      R[5] = 0.0;
299ccd2543fSMatthew G Knepley   R[6] =  n[0];            R[7] = n[1];            R[8] = n[2];
300ccd2543fSMatthew G Knepley   for (d = 0; d < dim; ++d) {
301ccd2543fSMatthew G Knepley     x1p[d] = 0.0;
302ccd2543fSMatthew G Knepley     x2p[d] = 0.0;
303ccd2543fSMatthew G Knepley     for (e = 0; e < dim; ++e) {
304ccd2543fSMatthew G Knepley       x1p[d] += R[d*dim+e]*x1[e];
305ccd2543fSMatthew G Knepley       x2p[d] += R[d*dim+e]*x2[e];
306ccd2543fSMatthew G Knepley     }
307ccd2543fSMatthew G Knepley   }
3088763be8eSMatthew G. Knepley   if (PetscAbsReal(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
3098763be8eSMatthew G. Knepley   if (PetscAbsReal(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
310ccd2543fSMatthew G Knepley   /* 2) Project to (x, y) */
31199dec3a6SMatthew G. Knepley   for (p = 3; p < coordSize/3; ++p) {
31299dec3a6SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
31399dec3a6SMatthew G. Knepley       xnp[d] = 0.0;
31499dec3a6SMatthew G. Knepley       for (e = 0; e < dim; ++e) {
31599dec3a6SMatthew G. Knepley         xnp[d] += R[d*dim+e]*PetscRealPart(coords[p*dim+e] - coords[0*dim+e]);
31699dec3a6SMatthew G. Knepley       }
31799dec3a6SMatthew G. Knepley       if (d < dim-1) coords[p*2+d] = xnp[d];
31899dec3a6SMatthew G. Knepley     }
31999dec3a6SMatthew G. Knepley   }
320ccd2543fSMatthew G Knepley   coords[0] = 0.0;
321ccd2543fSMatthew G Knepley   coords[1] = 0.0;
322ccd2543fSMatthew G Knepley   coords[2] = x1p[0];
323ccd2543fSMatthew G Knepley   coords[3] = x1p[1];
324ccd2543fSMatthew G Knepley   coords[4] = x2p[0];
325ccd2543fSMatthew G Knepley   coords[5] = x2p[1];
3267f07f362SMatthew G. Knepley   /* Output R^T which rotates \hat z to the input normal */
3277f07f362SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
3287f07f362SMatthew G. Knepley     for (e = d+1; e < dim; ++e) {
3297f07f362SMatthew G. Knepley       PetscReal tmp;
3307f07f362SMatthew G. Knepley 
3317f07f362SMatthew G. Knepley       tmp        = R[d*dim+e];
3327f07f362SMatthew G. Knepley       R[d*dim+e] = R[e*dim+d];
3337f07f362SMatthew G. Knepley       R[e*dim+d] = tmp;
3347f07f362SMatthew G. Knepley     }
3357f07f362SMatthew G. Knepley   }
336ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
337ccd2543fSMatthew G Knepley }
338ccd2543fSMatthew G Knepley 
339ccd2543fSMatthew G Knepley #undef __FUNCT__
340834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Triangle_Internal"
3416322fe33SJed Brown PETSC_UNUSED
342834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[])
343834e62ceSMatthew G. Knepley {
344834e62ceSMatthew G. Knepley   /* Signed volume is 1/2 the determinant
345834e62ceSMatthew G. Knepley 
346834e62ceSMatthew G. Knepley    |  1  1  1 |
347834e62ceSMatthew G. Knepley    | x0 x1 x2 |
348834e62ceSMatthew G. Knepley    | y0 y1 y2 |
349834e62ceSMatthew G. Knepley 
350834e62ceSMatthew G. Knepley      but if x0,y0 is the origin, we have
351834e62ceSMatthew G. Knepley 
352834e62ceSMatthew G. Knepley    | x1 x2 |
353834e62ceSMatthew G. Knepley    | y1 y2 |
354834e62ceSMatthew G. Knepley   */
355834e62ceSMatthew G. Knepley   const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1];
356834e62ceSMatthew G. Knepley   const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1];
357834e62ceSMatthew G. Knepley   PetscReal       M[4], detM;
358834e62ceSMatthew G. Knepley   M[0] = x1; M[1] = x2;
35986623015SMatthew G. Knepley   M[2] = y1; M[3] = y2;
360923591dfSMatthew G. Knepley   DMPlex_Det2D_Internal(&detM, M);
361834e62ceSMatthew G. Knepley   *vol = 0.5*detM;
362834e62ceSMatthew G. Knepley   PetscLogFlops(5.0);
363834e62ceSMatthew G. Knepley }
364834e62ceSMatthew G. Knepley 
365834e62ceSMatthew G. Knepley #undef __FUNCT__
366834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Triangle_Origin_Internal"
367834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[])
368834e62ceSMatthew G. Knepley {
369923591dfSMatthew G. Knepley   DMPlex_Det2D_Internal(vol, coords);
370834e62ceSMatthew G. Knepley   *vol *= 0.5;
371834e62ceSMatthew G. Knepley }
372834e62ceSMatthew G. Knepley 
373834e62ceSMatthew G. Knepley #undef __FUNCT__
374834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Tetrahedron_Internal"
3756322fe33SJed Brown PETSC_UNUSED
376834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[])
377834e62ceSMatthew G. Knepley {
378834e62ceSMatthew G. Knepley   /* Signed volume is 1/6th of the determinant
379834e62ceSMatthew G. Knepley 
380834e62ceSMatthew G. Knepley    |  1  1  1  1 |
381834e62ceSMatthew G. Knepley    | x0 x1 x2 x3 |
382834e62ceSMatthew G. Knepley    | y0 y1 y2 y3 |
383834e62ceSMatthew G. Knepley    | z0 z1 z2 z3 |
384834e62ceSMatthew G. Knepley 
385834e62ceSMatthew G. Knepley      but if x0,y0,z0 is the origin, we have
386834e62ceSMatthew G. Knepley 
387834e62ceSMatthew G. Knepley    | x1 x2 x3 |
388834e62ceSMatthew G. Knepley    | y1 y2 y3 |
389834e62ceSMatthew G. Knepley    | z1 z2 z3 |
390834e62ceSMatthew G. Knepley   */
391834e62ceSMatthew G. Knepley   const PetscReal x1 = coords[3] - coords[0], y1 = coords[4]  - coords[1], z1 = coords[5]  - coords[2];
392834e62ceSMatthew G. Knepley   const PetscReal x2 = coords[6] - coords[0], y2 = coords[7]  - coords[1], z2 = coords[8]  - coords[2];
393834e62ceSMatthew G. Knepley   const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2];
394834e62ceSMatthew G. Knepley   PetscReal       M[9], detM;
395834e62ceSMatthew G. Knepley   M[0] = x1; M[1] = x2; M[2] = x3;
396834e62ceSMatthew G. Knepley   M[3] = y1; M[4] = y2; M[5] = y3;
397834e62ceSMatthew G. Knepley   M[6] = z1; M[7] = z2; M[8] = z3;
398923591dfSMatthew G. Knepley   DMPlex_Det3D_Internal(&detM, M);
399b7ad821dSMatthew G. Knepley   *vol = -0.16666666666666666666666*detM;
400834e62ceSMatthew G. Knepley   PetscLogFlops(10.0);
401834e62ceSMatthew G. Knepley }
402834e62ceSMatthew G. Knepley 
403834e62ceSMatthew G. Knepley #undef __FUNCT__
4040ec8681fSMatthew G. Knepley #define __FUNCT__ "Volume_Tetrahedron_Origin_Internal"
4050ec8681fSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[])
4060ec8681fSMatthew G. Knepley {
407923591dfSMatthew G. Knepley   DMPlex_Det3D_Internal(vol, coords);
408b7ad821dSMatthew G. Knepley   *vol *= -0.16666666666666666666666;
4090ec8681fSMatthew G. Knepley }
4100ec8681fSMatthew G. Knepley 
4110ec8681fSMatthew G. Knepley #undef __FUNCT__
41217fe8556SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeLineGeometry_Internal"
41317fe8556SMatthew G. Knepley static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
41417fe8556SMatthew G. Knepley {
41517fe8556SMatthew G. Knepley   PetscSection   coordSection;
41617fe8556SMatthew G. Knepley   Vec            coordinates;
417a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
4187f07f362SMatthew G. Knepley   PetscInt       numCoords, d;
41917fe8556SMatthew G. Knepley   PetscErrorCode ierr;
42017fe8556SMatthew G. Knepley 
42117fe8556SMatthew G. Knepley   PetscFunctionBegin;
42217fe8556SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
42369d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
42417fe8556SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
4257f07f362SMatthew G. Knepley   *detJ = 0.0;
42617fe8556SMatthew G. Knepley   if (numCoords == 4) {
4277f07f362SMatthew G. Knepley     const PetscInt dim = 2;
4287f07f362SMatthew G. Knepley     PetscReal      R[4], J0;
4297f07f362SMatthew G. Knepley 
4307f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
4317f07f362SMatthew G. Knepley     ierr = DMPlexComputeProjection2Dto1D_Internal(coords, R);CHKERRQ(ierr);
43217fe8556SMatthew G. Knepley     if (J)    {
4337f07f362SMatthew G. Knepley       J0   = 0.5*PetscRealPart(coords[1]);
4347f07f362SMatthew G. Knepley       J[0] = R[0]*J0; J[1] = R[1];
4357f07f362SMatthew G. Knepley       J[2] = R[2]*J0; J[3] = R[3];
436923591dfSMatthew G. Knepley       DMPlex_Det2D_Internal(detJ, J);
43717fe8556SMatthew G. Knepley     }
438923591dfSMatthew G. Knepley     if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
4397f07f362SMatthew G. Knepley   } else if (numCoords == 2) {
4407f07f362SMatthew G. Knepley     const PetscInt dim = 1;
4417f07f362SMatthew G. Knepley 
4427f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
4437f07f362SMatthew G. Knepley     if (J)    {
4447f07f362SMatthew G. Knepley       J[0]  = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0]));
44517fe8556SMatthew G. Knepley       *detJ = J[0];
44617fe8556SMatthew G. Knepley       PetscLogFlops(2.0);
44717fe8556SMatthew G. Knepley     }
4487f07f362SMatthew G. Knepley     if (invJ) {invJ[0] = 1.0/J[0]; PetscLogFlops(1.0);}
4497f07f362SMatthew G. Knepley   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %d != 2", numCoords);
45017fe8556SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
45117fe8556SMatthew G. Knepley   PetscFunctionReturn(0);
45217fe8556SMatthew G. Knepley }
45317fe8556SMatthew G. Knepley 
45417fe8556SMatthew G. Knepley #undef __FUNCT__
455ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal"
456ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
457ccd2543fSMatthew G Knepley {
458ccd2543fSMatthew G Knepley   PetscSection   coordSection;
459ccd2543fSMatthew G Knepley   Vec            coordinates;
460a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
4617f07f362SMatthew G. Knepley   PetscInt       numCoords, d, f, g;
462ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
463ccd2543fSMatthew G Knepley 
464ccd2543fSMatthew G Knepley   PetscFunctionBegin;
465ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
46669d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
467ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
4687f07f362SMatthew G. Knepley   *detJ = 0.0;
469ccd2543fSMatthew G Knepley   if (numCoords == 9) {
4707f07f362SMatthew G. Knepley     const PetscInt dim = 3;
4717f07f362SMatthew 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};
4727f07f362SMatthew G. Knepley 
4737f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
47499dec3a6SMatthew G. Knepley     ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr);
4757f07f362SMatthew G. Knepley     if (J)    {
476b7ad821dSMatthew G. Knepley       const PetscInt pdim = 2;
477b7ad821dSMatthew G. Knepley 
478b7ad821dSMatthew G. Knepley       for (d = 0; d < pdim; d++) {
479b7ad821dSMatthew G. Knepley         for (f = 0; f < pdim; f++) {
480b7ad821dSMatthew G. Knepley           J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
481ccd2543fSMatthew G Knepley         }
4827f07f362SMatthew G. Knepley       }
4837f07f362SMatthew G. Knepley       PetscLogFlops(8.0);
484923591dfSMatthew G. Knepley       DMPlex_Det3D_Internal(detJ, J0);
4857f07f362SMatthew G. Knepley       for (d = 0; d < dim; d++) {
4867f07f362SMatthew G. Knepley         for (f = 0; f < dim; f++) {
4877f07f362SMatthew G. Knepley           J[d*dim+f] = 0.0;
4887f07f362SMatthew G. Knepley           for (g = 0; g < dim; g++) {
4897f07f362SMatthew G. Knepley             J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
4907f07f362SMatthew G. Knepley           }
4917f07f362SMatthew G. Knepley         }
4927f07f362SMatthew G. Knepley       }
4937f07f362SMatthew G. Knepley       PetscLogFlops(18.0);
4947f07f362SMatthew G. Knepley     }
495923591dfSMatthew G. Knepley     if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
4967f07f362SMatthew G. Knepley   } else if (numCoords == 6) {
4977f07f362SMatthew G. Knepley     const PetscInt dim = 2;
4987f07f362SMatthew G. Knepley 
4997f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
500ccd2543fSMatthew G Knepley     if (J)    {
501ccd2543fSMatthew G Knepley       for (d = 0; d < dim; d++) {
502ccd2543fSMatthew G Knepley         for (f = 0; f < dim; f++) {
503ccd2543fSMatthew G Knepley           J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
504ccd2543fSMatthew G Knepley         }
505ccd2543fSMatthew G Knepley       }
5067f07f362SMatthew G. Knepley       PetscLogFlops(8.0);
507923591dfSMatthew G. Knepley       DMPlex_Det2D_Internal(detJ, J);
508ccd2543fSMatthew G Knepley     }
509923591dfSMatthew G. Knepley     if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
5107f07f362SMatthew G. Knepley   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %d != 6", numCoords);
511ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
512ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
513ccd2543fSMatthew G Knepley }
514ccd2543fSMatthew G Knepley 
515ccd2543fSMatthew G Knepley #undef __FUNCT__
516ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal"
517ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
518ccd2543fSMatthew G Knepley {
519ccd2543fSMatthew G Knepley   PetscSection   coordSection;
520ccd2543fSMatthew G Knepley   Vec            coordinates;
521a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
52299dec3a6SMatthew G. Knepley   PetscInt       numCoords, d, f, g;
523ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
524ccd2543fSMatthew G Knepley 
525ccd2543fSMatthew G Knepley   PetscFunctionBegin;
526ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
52769d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
52899dec3a6SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
5297f07f362SMatthew G. Knepley   *detJ = 0.0;
53099dec3a6SMatthew G. Knepley   if (numCoords == 12) {
53199dec3a6SMatthew G. Knepley     const PetscInt dim = 3;
53299dec3a6SMatthew 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};
53399dec3a6SMatthew G. Knepley 
53499dec3a6SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
53599dec3a6SMatthew G. Knepley     ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr);
53699dec3a6SMatthew G. Knepley     if (J)    {
53799dec3a6SMatthew G. Knepley       const PetscInt pdim = 2;
53899dec3a6SMatthew G. Knepley 
53999dec3a6SMatthew G. Knepley       for (d = 0; d < pdim; d++) {
54099dec3a6SMatthew G. Knepley         J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
54199dec3a6SMatthew G. Knepley         J0[d*dim+1] = 0.5*(PetscRealPart(coords[3*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
54299dec3a6SMatthew G. Knepley       }
54399dec3a6SMatthew G. Knepley       PetscLogFlops(8.0);
544923591dfSMatthew G. Knepley       DMPlex_Det3D_Internal(detJ, J0);
54599dec3a6SMatthew G. Knepley       for (d = 0; d < dim; d++) {
54699dec3a6SMatthew G. Knepley         for (f = 0; f < dim; f++) {
54799dec3a6SMatthew G. Knepley           J[d*dim+f] = 0.0;
54899dec3a6SMatthew G. Knepley           for (g = 0; g < dim; g++) {
54999dec3a6SMatthew G. Knepley             J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
55099dec3a6SMatthew G. Knepley           }
55199dec3a6SMatthew G. Knepley         }
55299dec3a6SMatthew G. Knepley       }
55399dec3a6SMatthew G. Knepley       PetscLogFlops(18.0);
55499dec3a6SMatthew G. Knepley     }
555923591dfSMatthew G. Knepley     if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
556*97f1a218SMatthew G. Knepley   } else if ((numCoords == 8) || (numCoords == 16)) {
55799dec3a6SMatthew G. Knepley     const PetscInt dim = 2;
55899dec3a6SMatthew G. Knepley 
5597f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
560ccd2543fSMatthew G Knepley     if (J)    {
561ccd2543fSMatthew G Knepley       for (d = 0; d < dim; d++) {
56299dec3a6SMatthew G. Knepley         J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
56399dec3a6SMatthew G. Knepley         J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
564ccd2543fSMatthew G Knepley       }
5657f07f362SMatthew G. Knepley       PetscLogFlops(8.0);
566923591dfSMatthew G. Knepley       DMPlex_Det2D_Internal(detJ, J);
567ccd2543fSMatthew G Knepley     }
568923591dfSMatthew G. Knepley     if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
56999dec3a6SMatthew G. Knepley   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %d != 6", numCoords);
57099dec3a6SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
571ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
572ccd2543fSMatthew G Knepley }
573ccd2543fSMatthew G Knepley 
574ccd2543fSMatthew G Knepley #undef __FUNCT__
575ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal"
576ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
577ccd2543fSMatthew G Knepley {
578ccd2543fSMatthew G Knepley   PetscSection   coordSection;
579ccd2543fSMatthew G Knepley   Vec            coordinates;
580a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
581ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
58299dec3a6SMatthew G. Knepley   PetscInt       d;
583ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
584ccd2543fSMatthew G Knepley 
585ccd2543fSMatthew G Knepley   PetscFunctionBegin;
586ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
58769d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
588ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
5897f07f362SMatthew G. Knepley   *detJ = 0.0;
5907f07f362SMatthew G. Knepley   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
591ccd2543fSMatthew G Knepley   if (J)    {
592ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
593f0df753eSMatthew G. Knepley       /* I orient with outward face normals */
594f0df753eSMatthew G. Knepley       J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d]));
595f0df753eSMatthew G. Knepley       J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
596f0df753eSMatthew G. Knepley       J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
597ccd2543fSMatthew G Knepley     }
5987f07f362SMatthew G. Knepley     PetscLogFlops(18.0);
599923591dfSMatthew G. Knepley     DMPlex_Det3D_Internal(detJ, J);
600ccd2543fSMatthew G Knepley   }
601923591dfSMatthew G. Knepley   if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
602ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
603ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
604ccd2543fSMatthew G Knepley }
605ccd2543fSMatthew G Knepley 
606ccd2543fSMatthew G Knepley #undef __FUNCT__
607ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal"
608ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
609ccd2543fSMatthew G Knepley {
610ccd2543fSMatthew G Knepley   PetscSection   coordSection;
611ccd2543fSMatthew G Knepley   Vec            coordinates;
612a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
613ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
614ccd2543fSMatthew G Knepley   PetscInt       d;
615ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
616ccd2543fSMatthew G Knepley 
617ccd2543fSMatthew G Knepley   PetscFunctionBegin;
618ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
61969d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
620ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
6217f07f362SMatthew G. Knepley   *detJ = 0.0;
6227f07f362SMatthew G. Knepley   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
623ccd2543fSMatthew G Knepley   if (J)    {
624ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
625f0df753eSMatthew G. Knepley       J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
626f0df753eSMatthew G. Knepley       J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
627f0df753eSMatthew G. Knepley       J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d]));
628ccd2543fSMatthew G Knepley     }
6297f07f362SMatthew G. Knepley     PetscLogFlops(18.0);
630923591dfSMatthew G. Knepley     DMPlex_Det3D_Internal(detJ, J);
631ccd2543fSMatthew G Knepley   }
632923591dfSMatthew G. Knepley   if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
633ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
634ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
635ccd2543fSMatthew G Knepley }
636ccd2543fSMatthew G Knepley 
637ccd2543fSMatthew G Knepley #undef __FUNCT__
638ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeCellGeometry"
639ccd2543fSMatthew G Knepley /*@C
640ccd2543fSMatthew G Knepley   DMPlexComputeCellGeometry - Compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell
641ccd2543fSMatthew G Knepley 
642ccd2543fSMatthew G Knepley   Collective on DM
643ccd2543fSMatthew G Knepley 
644ccd2543fSMatthew G Knepley   Input Arguments:
645ccd2543fSMatthew G Knepley + dm   - the DM
646ccd2543fSMatthew G Knepley - cell - the cell
647ccd2543fSMatthew G Knepley 
648ccd2543fSMatthew G Knepley   Output Arguments:
649ccd2543fSMatthew G Knepley + v0   - the translation part of this affine transform
650ccd2543fSMatthew G Knepley . J    - the Jacobian of the transform from the reference element
651ccd2543fSMatthew G Knepley . invJ - the inverse of the Jacobian
652ccd2543fSMatthew G Knepley - detJ - the Jacobian determinant
653ccd2543fSMatthew G Knepley 
654ccd2543fSMatthew G Knepley   Level: advanced
655ccd2543fSMatthew G Knepley 
656ccd2543fSMatthew G Knepley   Fortran Notes:
657ccd2543fSMatthew G Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
658ccd2543fSMatthew G Knepley   include petsc.h90 in your code.
659ccd2543fSMatthew G Knepley 
66069d8a9ceSMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetCoordinateVec()
661ccd2543fSMatthew G Knepley @*/
662ccd2543fSMatthew G Knepley PetscErrorCode DMPlexComputeCellGeometry(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
663ccd2543fSMatthew G Knepley {
66449dc4407SMatthew G. Knepley   PetscInt       depth, dim, coneSize;
665ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
666ccd2543fSMatthew G Knepley 
667ccd2543fSMatthew G Knepley   PetscFunctionBegin;
668139a35ccSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
669ccd2543fSMatthew G Knepley   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
67049dc4407SMatthew G. Knepley   if (depth == 1) {
6715743f1d7SMatthew G. Knepley     ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
672ccd2543fSMatthew G Knepley     switch (dim) {
67317fe8556SMatthew G. Knepley     case 1:
67417fe8556SMatthew G. Knepley       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
67517fe8556SMatthew G. Knepley       break;
676ccd2543fSMatthew G Knepley     case 2:
677ccd2543fSMatthew G Knepley       switch (coneSize) {
678ccd2543fSMatthew G Knepley       case 3:
679ccd2543fSMatthew G Knepley         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
680ccd2543fSMatthew G Knepley         break;
681ccd2543fSMatthew G Knepley       case 4:
682ccd2543fSMatthew G Knepley         ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
683ccd2543fSMatthew G Knepley         break;
684ccd2543fSMatthew G Knepley       default:
685ccd2543fSMatthew G Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
686ccd2543fSMatthew G Knepley       }
687ccd2543fSMatthew G Knepley       break;
688ccd2543fSMatthew G Knepley     case 3:
689ccd2543fSMatthew G Knepley       switch (coneSize) {
690ccd2543fSMatthew G Knepley       case 4:
691ccd2543fSMatthew G Knepley         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
692ccd2543fSMatthew G Knepley         break;
693ccd2543fSMatthew G Knepley       case 8:
694ccd2543fSMatthew G Knepley         ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
695ccd2543fSMatthew G Knepley         break;
696ccd2543fSMatthew G Knepley       default:
697ccd2543fSMatthew G Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
698ccd2543fSMatthew G Knepley       }
699ccd2543fSMatthew G Knepley       break;
700ccd2543fSMatthew G Knepley     default:
701ccd2543fSMatthew G Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
702ccd2543fSMatthew G Knepley     }
70349dc4407SMatthew G. Knepley   } else {
7045743f1d7SMatthew G. Knepley     /* We need to keep a pointer to the depth label */
7055743f1d7SMatthew G. Knepley     ierr = DMPlexGetLabelValue(dm, "depth", cell, &dim);CHKERRQ(ierr);
70649dc4407SMatthew G. Knepley     /* Cone size is now the number of faces */
70749dc4407SMatthew G. Knepley     switch (dim) {
70817fe8556SMatthew G. Knepley     case 1:
70917fe8556SMatthew G. Knepley       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
71017fe8556SMatthew G. Knepley       break;
71149dc4407SMatthew G. Knepley     case 2:
71249dc4407SMatthew G. Knepley       switch (coneSize) {
71349dc4407SMatthew G. Knepley       case 3:
71449dc4407SMatthew G. Knepley         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
71549dc4407SMatthew G. Knepley         break;
71649dc4407SMatthew G. Knepley       case 4:
71749dc4407SMatthew G. Knepley         ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
71849dc4407SMatthew G. Knepley         break;
71949dc4407SMatthew G. Knepley       default:
72049dc4407SMatthew G. Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
72149dc4407SMatthew G. Knepley       }
72249dc4407SMatthew G. Knepley       break;
72349dc4407SMatthew G. Knepley     case 3:
72449dc4407SMatthew G. Knepley       switch (coneSize) {
72549dc4407SMatthew G. Knepley       case 4:
72649dc4407SMatthew G. Knepley         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
72749dc4407SMatthew G. Knepley         break;
72849dc4407SMatthew G. Knepley       case 6:
72949dc4407SMatthew G. Knepley         ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
73049dc4407SMatthew G. Knepley         break;
73149dc4407SMatthew G. Knepley       default:
73249dc4407SMatthew G. Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
73349dc4407SMatthew G. Knepley       }
73449dc4407SMatthew G. Knepley       break;
73549dc4407SMatthew G. Knepley     default:
73649dc4407SMatthew G. Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
73749dc4407SMatthew G. Knepley     }
73849dc4407SMatthew G. Knepley   }
739ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
740ccd2543fSMatthew G Knepley }
741834e62ceSMatthew G. Knepley 
742834e62ceSMatthew G. Knepley #undef __FUNCT__
743cc08537eSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_1D_Internal"
744011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
745cc08537eSMatthew G. Knepley {
746cc08537eSMatthew G. Knepley   PetscSection   coordSection;
747cc08537eSMatthew G. Knepley   Vec            coordinates;
748a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
749cc08537eSMatthew G. Knepley   PetscInt       coordSize;
750cc08537eSMatthew G. Knepley   PetscErrorCode ierr;
751cc08537eSMatthew G. Knepley 
752cc08537eSMatthew G. Knepley   PetscFunctionBegin;
753cc08537eSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
75469d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
755cc08537eSMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
756011ea5d8SMatthew G. Knepley   if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now");
757cc08537eSMatthew G. Knepley   if (centroid) {
7581ee9d5ecSMatthew G. Knepley     centroid[0] = 0.5*PetscRealPart(coords[0] + coords[dim+0]);
7591ee9d5ecSMatthew G. Knepley     centroid[1] = 0.5*PetscRealPart(coords[1] + coords[dim+1]);
760cc08537eSMatthew G. Knepley   }
761cc08537eSMatthew G. Knepley   if (normal) {
762a60a936bSMatthew G. Knepley     PetscReal norm;
763a60a936bSMatthew G. Knepley 
7640194f449SMatthew G. Knepley     normal[0] = -PetscRealPart(coords[1] - coords[dim+1]);
7650194f449SMatthew G. Knepley     normal[1] =  PetscRealPart(coords[0] - coords[dim+0]);
766a60a936bSMatthew G. Knepley     norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1]);
767a60a936bSMatthew G. Knepley     normal[0] /= norm;
768a60a936bSMatthew G. Knepley     normal[1] /= norm;
769cc08537eSMatthew G. Knepley   }
770cc08537eSMatthew G. Knepley   if (vol) {
7718b49ba18SBarry Smith     *vol = PetscSqrtReal(PetscSqr(PetscRealPart(coords[0] - coords[dim+0])) + PetscSqr(PetscRealPart(coords[1] - coords[dim+1])));
772cc08537eSMatthew G. Knepley   }
773cc08537eSMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
774cc08537eSMatthew G. Knepley   PetscFunctionReturn(0);
775cc08537eSMatthew G. Knepley }
776cc08537eSMatthew G. Knepley 
777cc08537eSMatthew G. Knepley #undef __FUNCT__
778cc08537eSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_2D_Internal"
779cc08537eSMatthew G. Knepley /* Centroid_i = (\sum_n A_n Cn_i ) / A */
780011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
781cc08537eSMatthew G. Knepley {
782cc08537eSMatthew G. Knepley   PetscSection   coordSection;
783cc08537eSMatthew G. Knepley   Vec            coordinates;
784cc08537eSMatthew G. Knepley   PetscScalar   *coords = NULL;
7850a1d6728SMatthew G. Knepley   PetscReal      vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9];
7860a1d6728SMatthew G. Knepley   PetscInt       tdim = 2, coordSize, numCorners, p, d, e;
787cc08537eSMatthew G. Knepley   PetscErrorCode ierr;
788cc08537eSMatthew G. Knepley 
789cc08537eSMatthew G. Knepley   PetscFunctionBegin;
790cc08537eSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
7910a1d6728SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr);
79269d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
793cc08537eSMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
7940a1d6728SMatthew G. Knepley   dim  = coordSize/numCorners;
795011ea5d8SMatthew G. Knepley   if (normal) {
796011ea5d8SMatthew G. Knepley     if (dim > 2) {
7971ee9d5ecSMatthew G. Knepley       const PetscReal x0 = PetscRealPart(coords[dim+0] - coords[0]), x1 = PetscRealPart(coords[dim*2+0] - coords[0]);
7981ee9d5ecSMatthew G. Knepley       const PetscReal y0 = PetscRealPart(coords[dim+1] - coords[1]), y1 = PetscRealPart(coords[dim*2+1] - coords[1]);
7991ee9d5ecSMatthew G. Knepley       const PetscReal z0 = PetscRealPart(coords[dim+2] - coords[2]), z1 = PetscRealPart(coords[dim*2+2] - coords[2]);
8000a1d6728SMatthew G. Knepley       PetscReal       norm;
8010a1d6728SMatthew G. Knepley 
8021ee9d5ecSMatthew G. Knepley       v0[0]     = PetscRealPart(coords[0]);
8031ee9d5ecSMatthew G. Knepley       v0[1]     = PetscRealPart(coords[1]);
8041ee9d5ecSMatthew G. Knepley       v0[2]     = PetscRealPart(coords[2]);
8050a1d6728SMatthew G. Knepley       normal[0] = y0*z1 - z0*y1;
8060a1d6728SMatthew G. Knepley       normal[1] = z0*x1 - x0*z1;
8070a1d6728SMatthew G. Knepley       normal[2] = x0*y1 - y0*x1;
8088b49ba18SBarry Smith       norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]);
8090a1d6728SMatthew G. Knepley       normal[0] /= norm;
8100a1d6728SMatthew G. Knepley       normal[1] /= norm;
8110a1d6728SMatthew G. Knepley       normal[2] /= norm;
812011ea5d8SMatthew G. Knepley     } else {
813011ea5d8SMatthew G. Knepley       for (d = 0; d < dim; ++d) normal[d] = 0.0;
814011ea5d8SMatthew G. Knepley     }
815011ea5d8SMatthew G. Knepley   }
81699dec3a6SMatthew G. Knepley   if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D_Internal(coordSize, coords, R);CHKERRQ(ierr);}
8170a1d6728SMatthew G. Knepley   for (p = 0; p < numCorners; ++p) {
8180a1d6728SMatthew G. Knepley     /* Need to do this copy to get types right */
8190a1d6728SMatthew G. Knepley     for (d = 0; d < tdim; ++d) {
8201ee9d5ecSMatthew G. Knepley       ctmp[d]      = PetscRealPart(coords[p*tdim+d]);
8211ee9d5ecSMatthew G. Knepley       ctmp[tdim+d] = PetscRealPart(coords[((p+1)%numCorners)*tdim+d]);
8220a1d6728SMatthew G. Knepley     }
8230a1d6728SMatthew G. Knepley     Volume_Triangle_Origin_Internal(&vtmp, ctmp);
8240a1d6728SMatthew G. Knepley     vsum += vtmp;
8250a1d6728SMatthew G. Knepley     for (d = 0; d < tdim; ++d) {
8260a1d6728SMatthew G. Knepley       csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp;
8270a1d6728SMatthew G. Knepley     }
8280a1d6728SMatthew G. Knepley   }
8290a1d6728SMatthew G. Knepley   for (d = 0; d < tdim; ++d) {
8300a1d6728SMatthew G. Knepley     csum[d] /= (tdim+1)*vsum;
8310a1d6728SMatthew G. Knepley   }
8320a1d6728SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
833ee6bbdb2SSatish Balay   if (vol) *vol = PetscAbsReal(vsum);
8340a1d6728SMatthew G. Knepley   if (centroid) {
8350a1d6728SMatthew G. Knepley     if (dim > 2) {
8360a1d6728SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
8370a1d6728SMatthew G. Knepley         centroid[d] = v0[d];
8380a1d6728SMatthew G. Knepley         for (e = 0; e < dim; ++e) {
8390a1d6728SMatthew G. Knepley           centroid[d] += R[d*dim+e]*csum[e];
8400a1d6728SMatthew G. Knepley         }
8410a1d6728SMatthew G. Knepley       }
8420a1d6728SMatthew G. Knepley     } else for (d = 0; d < dim; ++d) centroid[d] = csum[d];
8430a1d6728SMatthew G. Knepley   }
844cc08537eSMatthew G. Knepley   PetscFunctionReturn(0);
845cc08537eSMatthew G. Knepley }
846cc08537eSMatthew G. Knepley 
847cc08537eSMatthew G. Knepley #undef __FUNCT__
8480ec8681fSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_3D_Internal"
8490ec8681fSMatthew G. Knepley /* Centroid_i = (\sum_n V_n Cn_i ) / V */
850011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
8510ec8681fSMatthew G. Knepley {
8520ec8681fSMatthew G. Knepley   PetscSection    coordSection;
8530ec8681fSMatthew G. Knepley   Vec             coordinates;
8540ec8681fSMatthew G. Knepley   PetscScalar    *coords = NULL;
85586623015SMatthew G. Knepley   PetscReal       vsum = 0.0, vtmp, coordsTmp[3*3];
856a7df9edeSMatthew G. Knepley   const PetscInt *faces, *facesO;
8570ec8681fSMatthew G. Knepley   PetscInt        numFaces, f, coordSize, numCorners, p, d;
8580ec8681fSMatthew G. Knepley   PetscErrorCode  ierr;
8590ec8681fSMatthew G. Knepley 
8600ec8681fSMatthew G. Knepley   PetscFunctionBegin;
8610ec8681fSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
86269d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
8630ec8681fSMatthew G. Knepley 
864d9a81ebdSMatthew G. Knepley   if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0;
8650ec8681fSMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr);
8660ec8681fSMatthew G. Knepley   ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr);
867a7df9edeSMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr);
8680ec8681fSMatthew G. Knepley   for (f = 0; f < numFaces; ++f) {
869011ea5d8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr);
8700ec8681fSMatthew G. Knepley     numCorners = coordSize/dim;
8710ec8681fSMatthew G. Knepley     switch (numCorners) {
8720ec8681fSMatthew G. Knepley     case 3:
8730ec8681fSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
8741ee9d5ecSMatthew G. Knepley         coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]);
8751ee9d5ecSMatthew G. Knepley         coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]);
8761ee9d5ecSMatthew G. Knepley         coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]);
8770ec8681fSMatthew G. Knepley       }
8780ec8681fSMatthew G. Knepley       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
879a7df9edeSMatthew G. Knepley       if (facesO[f] < 0) vtmp = -vtmp;
8800ec8681fSMatthew G. Knepley       vsum += vtmp;
8814f25033aSJed Brown       if (centroid) {           /* Centroid of OABC = (a+b+c)/4 */
8820ec8681fSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
8831ee9d5ecSMatthew G. Knepley           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
8840ec8681fSMatthew G. Knepley         }
8850ec8681fSMatthew G. Knepley       }
8860ec8681fSMatthew G. Knepley       break;
8870ec8681fSMatthew G. Knepley     case 4:
8880ec8681fSMatthew G. Knepley       /* DO FOR PYRAMID */
8890ec8681fSMatthew G. Knepley       /* First tet */
8900ec8681fSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
8911ee9d5ecSMatthew G. Knepley         coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]);
8921ee9d5ecSMatthew G. Knepley         coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]);
8931ee9d5ecSMatthew G. Knepley         coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]);
8940ec8681fSMatthew G. Knepley       }
8950ec8681fSMatthew G. Knepley       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
896a7df9edeSMatthew G. Knepley       if (facesO[f] < 0) vtmp = -vtmp;
8970ec8681fSMatthew G. Knepley       vsum += vtmp;
8980ec8681fSMatthew G. Knepley       if (centroid) {
8990ec8681fSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
9000ec8681fSMatthew G. Knepley           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
9010ec8681fSMatthew G. Knepley         }
9020ec8681fSMatthew G. Knepley       }
9030ec8681fSMatthew G. Knepley       /* Second tet */
9040ec8681fSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
9051ee9d5ecSMatthew G. Knepley         coordsTmp[0*dim+d] = PetscRealPart(coords[1*dim+d]);
9061ee9d5ecSMatthew G. Knepley         coordsTmp[1*dim+d] = PetscRealPart(coords[2*dim+d]);
9071ee9d5ecSMatthew G. Knepley         coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]);
9080ec8681fSMatthew G. Knepley       }
9090ec8681fSMatthew G. Knepley       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
910a7df9edeSMatthew G. Knepley       if (facesO[f] < 0) vtmp = -vtmp;
9110ec8681fSMatthew G. Knepley       vsum += vtmp;
9120ec8681fSMatthew G. Knepley       if (centroid) {
9130ec8681fSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
9140ec8681fSMatthew G. Knepley           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
9150ec8681fSMatthew G. Knepley         }
9160ec8681fSMatthew G. Knepley       }
9170ec8681fSMatthew G. Knepley       break;
9180ec8681fSMatthew G. Knepley     default:
9190ec8681fSMatthew G. Knepley       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %d vertices", numCorners);
9200ec8681fSMatthew G. Knepley     }
9214f25033aSJed Brown     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr);
9220ec8681fSMatthew G. Knepley   }
9238763be8eSMatthew G. Knepley   if (vol)     *vol = PetscAbsReal(vsum);
9240ec8681fSMatthew G. Knepley   if (normal)   for (d = 0; d < dim; ++d) normal[d]    = 0.0;
925d9a81ebdSMatthew G. Knepley   if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4);
9260ec8681fSMatthew G. Knepley   PetscFunctionReturn(0);
9270ec8681fSMatthew G. Knepley }
9280ec8681fSMatthew G. Knepley 
9290ec8681fSMatthew G. Knepley #undef __FUNCT__
930834e62ceSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeCellGeometryFVM"
931834e62ceSMatthew G. Knepley /*@C
932834e62ceSMatthew G. Knepley   DMPlexComputeCellGeometryFVM - Compute the volume for a given cell
933834e62ceSMatthew G. Knepley 
934834e62ceSMatthew G. Knepley   Collective on DM
935834e62ceSMatthew G. Knepley 
936834e62ceSMatthew G. Knepley   Input Arguments:
937834e62ceSMatthew G. Knepley + dm   - the DM
938834e62ceSMatthew G. Knepley - cell - the cell
939834e62ceSMatthew G. Knepley 
940834e62ceSMatthew G. Knepley   Output Arguments:
941834e62ceSMatthew G. Knepley + volume   - the cell volume
942cc08537eSMatthew G. Knepley . centroid - the cell centroid
943cc08537eSMatthew G. Knepley - normal - the cell normal, if appropriate
944834e62ceSMatthew G. Knepley 
945834e62ceSMatthew G. Knepley   Level: advanced
946834e62ceSMatthew G. Knepley 
947834e62ceSMatthew G. Knepley   Fortran Notes:
948834e62ceSMatthew G. Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
949834e62ceSMatthew G. Knepley   include petsc.h90 in your code.
950834e62ceSMatthew G. Knepley 
95169d8a9ceSMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetCoordinateVec()
952834e62ceSMatthew G. Knepley @*/
953cc08537eSMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
954834e62ceSMatthew G. Knepley {
9550ec8681fSMatthew G. Knepley   PetscInt       depth, dim;
956834e62ceSMatthew G. Knepley   PetscErrorCode ierr;
957834e62ceSMatthew G. Knepley 
958834e62ceSMatthew G. Knepley   PetscFunctionBegin;
959834e62ceSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
960834e62ceSMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
961834e62ceSMatthew G. Knepley   if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated");
962834e62ceSMatthew G. Knepley   /* We need to keep a pointer to the depth label */
963011ea5d8SMatthew G. Knepley   ierr = DMPlexGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr);
964834e62ceSMatthew G. Knepley   /* Cone size is now the number of faces */
965011ea5d8SMatthew G. Knepley   switch (depth) {
966cc08537eSMatthew G. Knepley   case 1:
967011ea5d8SMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
968cc08537eSMatthew G. Knepley     break;
969834e62ceSMatthew G. Knepley   case 2:
970011ea5d8SMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
971834e62ceSMatthew G. Knepley     break;
972834e62ceSMatthew G. Knepley   case 3:
973011ea5d8SMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
974834e62ceSMatthew G. Knepley     break;
975834e62ceSMatthew G. Knepley   default:
976834e62ceSMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
977834e62ceSMatthew G. Knepley   }
978834e62ceSMatthew G. Knepley   PetscFunctionReturn(0);
979834e62ceSMatthew G. Knepley }
980