xref: /petsc/src/dm/impls/plex/plexgeometry.c (revision a7df9edef2185f07a4a9e227064efeea538de3b2)
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 */
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];
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   */
2598763be8eSMatthew G. Knepley   sqrtz = sqrt(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__
3407f07f362SMatthew G. Knepley #define __FUNCT__ "Invert2D_Internal"
3417f07f362SMatthew G. Knepley PETSC_STATIC_INLINE void Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
3427f07f362SMatthew G. Knepley {
3437f07f362SMatthew G. Knepley   const PetscReal invDet = 1.0/detJ;
3447f07f362SMatthew G. Knepley 
3457f07f362SMatthew G. Knepley   invJ[0] =  invDet*J[3];
3467f07f362SMatthew G. Knepley   invJ[1] = -invDet*J[1];
3477f07f362SMatthew G. Knepley   invJ[2] = -invDet*J[2];
3487f07f362SMatthew G. Knepley   invJ[3] =  invDet*J[0];
3497f07f362SMatthew G. Knepley   PetscLogFlops(5.0);
3507f07f362SMatthew G. Knepley }
3517f07f362SMatthew G. Knepley 
3527f07f362SMatthew G. Knepley #undef __FUNCT__
3537f07f362SMatthew G. Knepley #define __FUNCT__ "Invert3D_Internal"
3547f07f362SMatthew G. Knepley PETSC_STATIC_INLINE void Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
3557f07f362SMatthew G. Knepley {
3567f07f362SMatthew G. Knepley   const PetscReal invDet = 1.0/detJ;
3577f07f362SMatthew G. Knepley 
3587f07f362SMatthew G. Knepley   invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
3597f07f362SMatthew G. Knepley   invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
3607f07f362SMatthew G. Knepley   invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
3617f07f362SMatthew G. Knepley   invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
3627f07f362SMatthew G. Knepley   invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
3637f07f362SMatthew G. Knepley   invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
3647f07f362SMatthew G. Knepley   invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
3657f07f362SMatthew G. Knepley   invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
3667f07f362SMatthew G. Knepley   invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
3677f07f362SMatthew G. Knepley   PetscLogFlops(37.0);
3687f07f362SMatthew G. Knepley }
3697f07f362SMatthew G. Knepley 
3707f07f362SMatthew G. Knepley #undef __FUNCT__
3717f07f362SMatthew G. Knepley #define __FUNCT__ "Det2D_Internal"
3727f07f362SMatthew G. Knepley PETSC_STATIC_INLINE void Det2D_Internal(PetscReal *detJ, PetscReal J[])
3737f07f362SMatthew G. Knepley {
3747f07f362SMatthew G. Knepley   *detJ = J[0]*J[3] - J[1]*J[2];
3757f07f362SMatthew G. Knepley   PetscLogFlops(3.0);
3767f07f362SMatthew G. Knepley }
3777f07f362SMatthew G. Knepley 
3787f07f362SMatthew G. Knepley #undef __FUNCT__
3797f07f362SMatthew G. Knepley #define __FUNCT__ "Det3D_Internal"
3807f07f362SMatthew G. Knepley PETSC_STATIC_INLINE void Det3D_Internal(PetscReal *detJ, PetscReal J[])
3817f07f362SMatthew G. Knepley {
382b7ad821dSMatthew G. Knepley   *detJ = (J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
383b7ad821dSMatthew G. Knepley            J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
384b7ad821dSMatthew G. Knepley            J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
3857f07f362SMatthew G. Knepley   PetscLogFlops(12.0);
3867f07f362SMatthew G. Knepley }
3877f07f362SMatthew G. Knepley 
3887f07f362SMatthew G. Knepley #undef __FUNCT__
389834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Triangle_Internal"
390834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[])
391834e62ceSMatthew G. Knepley {
392834e62ceSMatthew G. Knepley   /* Signed volume is 1/2 the determinant
393834e62ceSMatthew G. Knepley 
394834e62ceSMatthew G. Knepley    |  1  1  1 |
395834e62ceSMatthew G. Knepley    | x0 x1 x2 |
396834e62ceSMatthew G. Knepley    | y0 y1 y2 |
397834e62ceSMatthew G. Knepley 
398834e62ceSMatthew G. Knepley      but if x0,y0 is the origin, we have
399834e62ceSMatthew G. Knepley 
400834e62ceSMatthew G. Knepley    | x1 x2 |
401834e62ceSMatthew G. Knepley    | y1 y2 |
402834e62ceSMatthew G. Knepley   */
403834e62ceSMatthew G. Knepley   const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1];
404834e62ceSMatthew G. Knepley   const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1];
405834e62ceSMatthew G. Knepley   PetscReal       M[4], detM;
406834e62ceSMatthew G. Knepley   M[0] = x1; M[1] = x2;
40786623015SMatthew G. Knepley   M[2] = y1; M[3] = y2;
408834e62ceSMatthew G. Knepley   Det2D_Internal(&detM, M);
409834e62ceSMatthew G. Knepley   *vol = 0.5*detM;
410834e62ceSMatthew G. Knepley   PetscLogFlops(5.0);
411834e62ceSMatthew G. Knepley }
412834e62ceSMatthew G. Knepley 
413834e62ceSMatthew G. Knepley #undef __FUNCT__
414834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Triangle_Origin_Internal"
415834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[])
416834e62ceSMatthew G. Knepley {
417834e62ceSMatthew G. Knepley   Det2D_Internal(vol, coords);
418834e62ceSMatthew G. Knepley   *vol *= 0.5;
419834e62ceSMatthew G. Knepley }
420834e62ceSMatthew G. Knepley 
421834e62ceSMatthew G. Knepley #undef __FUNCT__
422834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Tetrahedron_Internal"
423834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[])
424834e62ceSMatthew G. Knepley {
425834e62ceSMatthew G. Knepley   /* Signed volume is 1/6th of the determinant
426834e62ceSMatthew G. Knepley 
427834e62ceSMatthew G. Knepley    |  1  1  1  1 |
428834e62ceSMatthew G. Knepley    | x0 x1 x2 x3 |
429834e62ceSMatthew G. Knepley    | y0 y1 y2 y3 |
430834e62ceSMatthew G. Knepley    | z0 z1 z2 z3 |
431834e62ceSMatthew G. Knepley 
432834e62ceSMatthew G. Knepley      but if x0,y0,z0 is the origin, we have
433834e62ceSMatthew G. Knepley 
434834e62ceSMatthew G. Knepley    | x1 x2 x3 |
435834e62ceSMatthew G. Knepley    | y1 y2 y3 |
436834e62ceSMatthew G. Knepley    | z1 z2 z3 |
437834e62ceSMatthew G. Knepley   */
438834e62ceSMatthew G. Knepley   const PetscReal x1 = coords[3] - coords[0], y1 = coords[4]  - coords[1], z1 = coords[5]  - coords[2];
439834e62ceSMatthew G. Knepley   const PetscReal x2 = coords[6] - coords[0], y2 = coords[7]  - coords[1], z2 = coords[8]  - coords[2];
440834e62ceSMatthew G. Knepley   const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2];
441834e62ceSMatthew G. Knepley   PetscReal       M[9], detM;
442834e62ceSMatthew G. Knepley   M[0] = x1; M[1] = x2; M[2] = x3;
443834e62ceSMatthew G. Knepley   M[3] = y1; M[4] = y2; M[5] = y3;
444834e62ceSMatthew G. Knepley   M[6] = z1; M[7] = z2; M[8] = z3;
445834e62ceSMatthew G. Knepley   Det3D_Internal(&detM, M);
446b7ad821dSMatthew G. Knepley   *vol = -0.16666666666666666666666*detM;
447834e62ceSMatthew G. Knepley   PetscLogFlops(10.0);
448834e62ceSMatthew G. Knepley }
449834e62ceSMatthew G. Knepley 
450834e62ceSMatthew G. Knepley #undef __FUNCT__
4510ec8681fSMatthew G. Knepley #define __FUNCT__ "Volume_Tetrahedron_Origin_Internal"
4520ec8681fSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[])
4530ec8681fSMatthew G. Knepley {
4540ec8681fSMatthew G. Knepley   Det3D_Internal(vol, coords);
455b7ad821dSMatthew G. Knepley   *vol *= -0.16666666666666666666666;
4560ec8681fSMatthew G. Knepley }
4570ec8681fSMatthew G. Knepley 
4580ec8681fSMatthew G. Knepley #undef __FUNCT__
45917fe8556SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeLineGeometry_Internal"
46017fe8556SMatthew G. Knepley static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
46117fe8556SMatthew G. Knepley {
46217fe8556SMatthew G. Knepley   PetscSection   coordSection;
46317fe8556SMatthew G. Knepley   Vec            coordinates;
46417fe8556SMatthew G. Knepley   PetscScalar   *coords;
4657f07f362SMatthew G. Knepley   PetscInt       numCoords, d;
46617fe8556SMatthew G. Knepley   PetscErrorCode ierr;
46717fe8556SMatthew G. Knepley 
46817fe8556SMatthew G. Knepley   PetscFunctionBegin;
46917fe8556SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
47017fe8556SMatthew G. Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
47117fe8556SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
4727f07f362SMatthew G. Knepley   *detJ = 0.0;
47317fe8556SMatthew G. Knepley   if (numCoords == 4) {
4747f07f362SMatthew G. Knepley     const PetscInt dim = 2;
4757f07f362SMatthew G. Knepley     PetscReal      R[4], J0;
4767f07f362SMatthew G. Knepley 
4777f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
4787f07f362SMatthew G. Knepley     ierr = DMPlexComputeProjection2Dto1D_Internal(coords, R);CHKERRQ(ierr);
47917fe8556SMatthew G. Knepley     if (J)    {
4807f07f362SMatthew G. Knepley       J0   = 0.5*PetscRealPart(coords[1]);
4817f07f362SMatthew G. Knepley       J[0] = R[0]*J0; J[1] = R[1];
4827f07f362SMatthew G. Knepley       J[2] = R[2]*J0; J[3] = R[3];
4837f07f362SMatthew G. Knepley       Det2D_Internal(detJ, J);
48417fe8556SMatthew G. Knepley     }
4857f07f362SMatthew G. Knepley     if (invJ) {Invert2D_Internal(invJ, J, *detJ);}
4867f07f362SMatthew G. Knepley   } else if (numCoords == 2) {
4877f07f362SMatthew G. Knepley     const PetscInt dim = 1;
4887f07f362SMatthew G. Knepley 
4897f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
4907f07f362SMatthew G. Knepley     if (J)    {
4917f07f362SMatthew G. Knepley       J[0]  = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0]));
49217fe8556SMatthew G. Knepley       *detJ = J[0];
49317fe8556SMatthew G. Knepley       PetscLogFlops(2.0);
49417fe8556SMatthew G. Knepley     }
4957f07f362SMatthew G. Knepley     if (invJ) {invJ[0] = 1.0/J[0]; PetscLogFlops(1.0);}
4967f07f362SMatthew G. Knepley   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %d != 2", numCoords);
49717fe8556SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
49817fe8556SMatthew G. Knepley   PetscFunctionReturn(0);
49917fe8556SMatthew G. Knepley }
50017fe8556SMatthew G. Knepley 
50117fe8556SMatthew G. Knepley #undef __FUNCT__
502ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal"
503ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
504ccd2543fSMatthew G Knepley {
505ccd2543fSMatthew G Knepley   PetscSection   coordSection;
506ccd2543fSMatthew G Knepley   Vec            coordinates;
5077c1f9639SMatthew G Knepley   PetscScalar   *coords;
5087f07f362SMatthew G. Knepley   PetscInt       numCoords, d, f, g;
509ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
510ccd2543fSMatthew G Knepley 
511ccd2543fSMatthew G Knepley   PetscFunctionBegin;
512ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
513ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
514ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
5157f07f362SMatthew G. Knepley   *detJ = 0.0;
516ccd2543fSMatthew G Knepley   if (numCoords == 9) {
5177f07f362SMatthew G. Knepley     const PetscInt dim = 3;
5187f07f362SMatthew 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};
5197f07f362SMatthew G. Knepley 
5207f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
52199dec3a6SMatthew G. Knepley     ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr);
5227f07f362SMatthew G. Knepley     if (J)    {
523b7ad821dSMatthew G. Knepley       const PetscInt pdim = 2;
524b7ad821dSMatthew G. Knepley 
525b7ad821dSMatthew G. Knepley       for (d = 0; d < pdim; d++) {
526b7ad821dSMatthew G. Knepley         for (f = 0; f < pdim; f++) {
527b7ad821dSMatthew G. Knepley           J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
528ccd2543fSMatthew G Knepley         }
5297f07f362SMatthew G. Knepley       }
5307f07f362SMatthew G. Knepley       PetscLogFlops(8.0);
53187d60e71SMatthew G. Knepley       Det3D_Internal(detJ, J0);
5327f07f362SMatthew G. Knepley       for (d = 0; d < dim; d++) {
5337f07f362SMatthew G. Knepley         for (f = 0; f < dim; f++) {
5347f07f362SMatthew G. Knepley           J[d*dim+f] = 0.0;
5357f07f362SMatthew G. Knepley           for (g = 0; g < dim; g++) {
5367f07f362SMatthew G. Knepley             J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
5377f07f362SMatthew G. Knepley           }
5387f07f362SMatthew G. Knepley         }
5397f07f362SMatthew G. Knepley       }
5407f07f362SMatthew G. Knepley       PetscLogFlops(18.0);
5417f07f362SMatthew G. Knepley     }
5427f07f362SMatthew G. Knepley     if (invJ) {Invert3D_Internal(invJ, J, *detJ);}
5437f07f362SMatthew G. Knepley   } else if (numCoords == 6) {
5447f07f362SMatthew G. Knepley     const PetscInt dim = 2;
5457f07f362SMatthew G. Knepley 
5467f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
547ccd2543fSMatthew G Knepley     if (J)    {
548ccd2543fSMatthew G Knepley       for (d = 0; d < dim; d++) {
549ccd2543fSMatthew G Knepley         for (f = 0; f < dim; f++) {
550ccd2543fSMatthew G Knepley           J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
551ccd2543fSMatthew G Knepley         }
552ccd2543fSMatthew G Knepley       }
5537f07f362SMatthew G. Knepley       PetscLogFlops(8.0);
5547f07f362SMatthew G. Knepley       Det2D_Internal(detJ, J);
555ccd2543fSMatthew G Knepley     }
5567f07f362SMatthew G. Knepley     if (invJ) {Invert2D_Internal(invJ, J, *detJ);}
5577f07f362SMatthew G. Knepley   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %d != 6", numCoords);
558ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
559ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
560ccd2543fSMatthew G Knepley }
561ccd2543fSMatthew G Knepley 
562ccd2543fSMatthew G Knepley #undef __FUNCT__
563ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal"
564ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
565ccd2543fSMatthew G Knepley {
566ccd2543fSMatthew G Knepley   PetscSection   coordSection;
567ccd2543fSMatthew G Knepley   Vec            coordinates;
5687c1f9639SMatthew G Knepley   PetscScalar   *coords;
56999dec3a6SMatthew G. Knepley   PetscInt       numCoords, d, f, g;
570ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
571ccd2543fSMatthew G Knepley 
572ccd2543fSMatthew G Knepley   PetscFunctionBegin;
573ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
574ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
57599dec3a6SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
5767f07f362SMatthew G. Knepley   *detJ = 0.0;
57799dec3a6SMatthew G. Knepley   if (numCoords == 12) {
57899dec3a6SMatthew G. Knepley     const PetscInt dim = 3;
57999dec3a6SMatthew 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};
58099dec3a6SMatthew G. Knepley 
58199dec3a6SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
58299dec3a6SMatthew G. Knepley     ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr);
58399dec3a6SMatthew G. Knepley     if (J)    {
58499dec3a6SMatthew G. Knepley       const PetscInt pdim = 2;
58599dec3a6SMatthew G. Knepley 
58699dec3a6SMatthew G. Knepley       for (d = 0; d < pdim; d++) {
58799dec3a6SMatthew G. Knepley         J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
58899dec3a6SMatthew G. Knepley         J0[d*dim+1] = 0.5*(PetscRealPart(coords[3*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
58999dec3a6SMatthew G. Knepley       }
59099dec3a6SMatthew G. Knepley       PetscLogFlops(8.0);
59199dec3a6SMatthew G. Knepley       Det3D_Internal(detJ, J0);
59299dec3a6SMatthew G. Knepley       for (d = 0; d < dim; d++) {
59399dec3a6SMatthew G. Knepley         for (f = 0; f < dim; f++) {
59499dec3a6SMatthew G. Knepley           J[d*dim+f] = 0.0;
59599dec3a6SMatthew G. Knepley           for (g = 0; g < dim; g++) {
59699dec3a6SMatthew G. Knepley             J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
59799dec3a6SMatthew G. Knepley           }
59899dec3a6SMatthew G. Knepley         }
59999dec3a6SMatthew G. Knepley       }
60099dec3a6SMatthew G. Knepley       PetscLogFlops(18.0);
60199dec3a6SMatthew G. Knepley     }
60299dec3a6SMatthew G. Knepley     if (invJ) {Invert3D_Internal(invJ, J, *detJ);}
60399dec3a6SMatthew G. Knepley   } else if (numCoords == 8) {
60499dec3a6SMatthew G. Knepley     const PetscInt dim = 2;
60599dec3a6SMatthew G. Knepley 
6067f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
607ccd2543fSMatthew G Knepley     if (J)    {
608ccd2543fSMatthew G Knepley       for (d = 0; d < dim; d++) {
60999dec3a6SMatthew G. Knepley         J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
61099dec3a6SMatthew G. Knepley         J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
611ccd2543fSMatthew G Knepley       }
6127f07f362SMatthew G. Knepley       PetscLogFlops(8.0);
6137f07f362SMatthew G. Knepley       Det2D_Internal(detJ, J);
614ccd2543fSMatthew G Knepley     }
6157f07f362SMatthew G. Knepley     if (invJ) {Invert2D_Internal(invJ, J, *detJ);}
61699dec3a6SMatthew G. Knepley   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %d != 6", numCoords);
61799dec3a6SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
618ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
619ccd2543fSMatthew G Knepley }
620ccd2543fSMatthew G Knepley 
621ccd2543fSMatthew G Knepley #undef __FUNCT__
622ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal"
623ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
624ccd2543fSMatthew G Knepley {
625ccd2543fSMatthew G Knepley   PetscSection   coordSection;
626ccd2543fSMatthew G Knepley   Vec            coordinates;
6277c1f9639SMatthew G Knepley   PetscScalar   *coords;
628ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
62999dec3a6SMatthew G. Knepley   PetscInt       d;
630ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
631ccd2543fSMatthew G Knepley 
632ccd2543fSMatthew G Knepley   PetscFunctionBegin;
633ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
634ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
635ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
6367f07f362SMatthew G. Knepley   *detJ = 0.0;
6377f07f362SMatthew G. Knepley   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
638ccd2543fSMatthew G Knepley   if (J)    {
639ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
640f0df753eSMatthew G. Knepley       /* I orient with outward face normals */
641f0df753eSMatthew G. Knepley       J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d]));
642f0df753eSMatthew G. Knepley       J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
643f0df753eSMatthew G. Knepley       J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
644ccd2543fSMatthew G Knepley     }
6457f07f362SMatthew G. Knepley     PetscLogFlops(18.0);
6467f07f362SMatthew G. Knepley     Det3D_Internal(detJ, J);
647ccd2543fSMatthew G Knepley   }
648f0df753eSMatthew G. Knepley   if (invJ) {Invert3D_Internal(invJ, J, *detJ);}
649ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
650ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
651ccd2543fSMatthew G Knepley }
652ccd2543fSMatthew G Knepley 
653ccd2543fSMatthew G Knepley #undef __FUNCT__
654ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal"
655ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
656ccd2543fSMatthew G Knepley {
657ccd2543fSMatthew G Knepley   PetscSection   coordSection;
658ccd2543fSMatthew G Knepley   Vec            coordinates;
6597c1f9639SMatthew G Knepley   PetscScalar   *coords;
660ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
661ccd2543fSMatthew G Knepley   PetscInt       d;
662ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
663ccd2543fSMatthew G Knepley 
664ccd2543fSMatthew G Knepley   PetscFunctionBegin;
665ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
666ccd2543fSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
667ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
6687f07f362SMatthew G. Knepley   *detJ = 0.0;
6697f07f362SMatthew G. Knepley   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
670ccd2543fSMatthew G Knepley   if (J)    {
671ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
672f0df753eSMatthew G. Knepley       J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
673f0df753eSMatthew G. Knepley       J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
674f0df753eSMatthew G. Knepley       J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d]));
675ccd2543fSMatthew G Knepley     }
6767f07f362SMatthew G. Knepley     PetscLogFlops(18.0);
6777f07f362SMatthew G. Knepley     Det3D_Internal(detJ, J);
678ccd2543fSMatthew G Knepley   }
6797f07f362SMatthew G. Knepley   if (invJ) {Invert3D_Internal(invJ, J, *detJ);}
680ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
681ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
682ccd2543fSMatthew G Knepley }
683ccd2543fSMatthew G Knepley 
684ccd2543fSMatthew G Knepley #undef __FUNCT__
685ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeCellGeometry"
686ccd2543fSMatthew G Knepley /*@C
687ccd2543fSMatthew G Knepley   DMPlexComputeCellGeometry - Compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell
688ccd2543fSMatthew G Knepley 
689ccd2543fSMatthew G Knepley   Collective on DM
690ccd2543fSMatthew G Knepley 
691ccd2543fSMatthew G Knepley   Input Arguments:
692ccd2543fSMatthew G Knepley + dm   - the DM
693ccd2543fSMatthew G Knepley - cell - the cell
694ccd2543fSMatthew G Knepley 
695ccd2543fSMatthew G Knepley   Output Arguments:
696ccd2543fSMatthew G Knepley + v0   - the translation part of this affine transform
697ccd2543fSMatthew G Knepley . J    - the Jacobian of the transform from the reference element
698ccd2543fSMatthew G Knepley . invJ - the inverse of the Jacobian
699ccd2543fSMatthew G Knepley - detJ - the Jacobian determinant
700ccd2543fSMatthew G Knepley 
701ccd2543fSMatthew G Knepley   Level: advanced
702ccd2543fSMatthew G Knepley 
703ccd2543fSMatthew G Knepley   Fortran Notes:
704ccd2543fSMatthew G Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
705ccd2543fSMatthew G Knepley   include petsc.h90 in your code.
706ccd2543fSMatthew G Knepley 
707ccd2543fSMatthew G Knepley .seealso: DMPlexGetCoordinateSection(), DMPlexGetCoordinateVec()
708ccd2543fSMatthew G Knepley @*/
709ccd2543fSMatthew G Knepley PetscErrorCode DMPlexComputeCellGeometry(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
710ccd2543fSMatthew G Knepley {
71149dc4407SMatthew G. Knepley   PetscInt       depth, dim, coneSize;
712ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
713ccd2543fSMatthew G Knepley 
714ccd2543fSMatthew G Knepley   PetscFunctionBegin;
715139a35ccSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
716ccd2543fSMatthew G Knepley   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
71749dc4407SMatthew G. Knepley   if (depth == 1) {
7185743f1d7SMatthew G. Knepley     ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
719ccd2543fSMatthew G Knepley     switch (dim) {
72017fe8556SMatthew G. Knepley     case 1:
72117fe8556SMatthew G. Knepley       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
72217fe8556SMatthew G. Knepley       break;
723ccd2543fSMatthew G Knepley     case 2:
724ccd2543fSMatthew G Knepley       switch (coneSize) {
725ccd2543fSMatthew G Knepley       case 3:
726ccd2543fSMatthew G Knepley         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
727ccd2543fSMatthew G Knepley         break;
728ccd2543fSMatthew G Knepley       case 4:
729ccd2543fSMatthew G Knepley         ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
730ccd2543fSMatthew G Knepley         break;
731ccd2543fSMatthew G Knepley       default:
732ccd2543fSMatthew G Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
733ccd2543fSMatthew G Knepley       }
734ccd2543fSMatthew G Knepley       break;
735ccd2543fSMatthew G Knepley     case 3:
736ccd2543fSMatthew G Knepley       switch (coneSize) {
737ccd2543fSMatthew G Knepley       case 4:
738ccd2543fSMatthew G Knepley         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
739ccd2543fSMatthew G Knepley         break;
740ccd2543fSMatthew G Knepley       case 8:
741ccd2543fSMatthew G Knepley         ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
742ccd2543fSMatthew G Knepley         break;
743ccd2543fSMatthew G Knepley       default:
744ccd2543fSMatthew G Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
745ccd2543fSMatthew G Knepley       }
746ccd2543fSMatthew G Knepley       break;
747ccd2543fSMatthew G Knepley     default:
748ccd2543fSMatthew G Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
749ccd2543fSMatthew G Knepley     }
75049dc4407SMatthew G. Knepley   } else {
7515743f1d7SMatthew G. Knepley     /* We need to keep a pointer to the depth label */
7525743f1d7SMatthew G. Knepley     ierr = DMPlexGetLabelValue(dm, "depth", cell, &dim);CHKERRQ(ierr);
75349dc4407SMatthew G. Knepley     /* Cone size is now the number of faces */
75449dc4407SMatthew G. Knepley     switch (dim) {
75517fe8556SMatthew G. Knepley     case 1:
75617fe8556SMatthew G. Knepley       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
75717fe8556SMatthew G. Knepley       break;
75849dc4407SMatthew G. Knepley     case 2:
75949dc4407SMatthew G. Knepley       switch (coneSize) {
76049dc4407SMatthew G. Knepley       case 3:
76149dc4407SMatthew G. Knepley         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
76249dc4407SMatthew G. Knepley         break;
76349dc4407SMatthew G. Knepley       case 4:
76449dc4407SMatthew G. Knepley         ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
76549dc4407SMatthew G. Knepley         break;
76649dc4407SMatthew G. Knepley       default:
76749dc4407SMatthew G. Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
76849dc4407SMatthew G. Knepley       }
76949dc4407SMatthew G. Knepley       break;
77049dc4407SMatthew G. Knepley     case 3:
77149dc4407SMatthew G. Knepley       switch (coneSize) {
77249dc4407SMatthew G. Knepley       case 4:
77349dc4407SMatthew G. Knepley         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
77449dc4407SMatthew G. Knepley         break;
77549dc4407SMatthew G. Knepley       case 6:
77649dc4407SMatthew G. Knepley         ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);
77749dc4407SMatthew G. Knepley         break;
77849dc4407SMatthew G. Knepley       default:
77949dc4407SMatthew G. Knepley         SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell);
78049dc4407SMatthew G. Knepley       }
78149dc4407SMatthew G. Knepley       break;
78249dc4407SMatthew G. Knepley     default:
78349dc4407SMatthew G. Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
78449dc4407SMatthew G. Knepley     }
78549dc4407SMatthew G. Knepley   }
786ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
787ccd2543fSMatthew G Knepley }
788834e62ceSMatthew G. Knepley 
789834e62ceSMatthew G. Knepley #undef __FUNCT__
790cc08537eSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_1D_Internal"
791011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
792cc08537eSMatthew G. Knepley {
793cc08537eSMatthew G. Knepley   PetscSection   coordSection;
794cc08537eSMatthew G. Knepley   Vec            coordinates;
795cc08537eSMatthew G. Knepley   PetscScalar   *coords;
796cc08537eSMatthew G. Knepley   PetscInt       coordSize;
797cc08537eSMatthew G. Knepley   PetscErrorCode ierr;
798cc08537eSMatthew G. Knepley 
799cc08537eSMatthew G. Knepley   PetscFunctionBegin;
800cc08537eSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
801cc08537eSMatthew G. Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
802cc08537eSMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
803011ea5d8SMatthew G. Knepley   if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now");
804cc08537eSMatthew G. Knepley   if (centroid) {
8051ee9d5ecSMatthew G. Knepley     centroid[0] = 0.5*PetscRealPart(coords[0] + coords[dim+0]);
8061ee9d5ecSMatthew G. Knepley     centroid[1] = 0.5*PetscRealPart(coords[1] + coords[dim+1]);
807cc08537eSMatthew G. Knepley   }
808cc08537eSMatthew G. Knepley   if (normal) {
8091ee9d5ecSMatthew G. Knepley     normal[0] =  PetscRealPart(coords[1] - coords[dim+1]);
8101ee9d5ecSMatthew G. Knepley     normal[1] = -PetscRealPart(coords[0] - coords[dim+0]);
811cc08537eSMatthew G. Knepley   }
812cc08537eSMatthew G. Knepley   if (vol) {
8131ee9d5ecSMatthew G. Knepley     *vol = sqrt(PetscSqr(PetscRealPart(coords[0] - coords[dim+0])) + PetscSqr(PetscRealPart(coords[1] - coords[dim+1])));
814cc08537eSMatthew G. Knepley   }
815cc08537eSMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
816cc08537eSMatthew G. Knepley   PetscFunctionReturn(0);
817cc08537eSMatthew G. Knepley }
818cc08537eSMatthew G. Knepley 
819cc08537eSMatthew G. Knepley #undef __FUNCT__
820cc08537eSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_2D_Internal"
821cc08537eSMatthew G. Knepley /* Centroid_i = (\sum_n A_n Cn_i ) / A */
822011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
823cc08537eSMatthew G. Knepley {
824cc08537eSMatthew G. Knepley   PetscSection   coordSection;
825cc08537eSMatthew G. Knepley   Vec            coordinates;
826cc08537eSMatthew G. Knepley   PetscScalar   *coords = NULL;
8270a1d6728SMatthew G. Knepley   PetscReal      vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9];
8280a1d6728SMatthew G. Knepley   PetscInt       tdim = 2, coordSize, numCorners, p, d, e;
829cc08537eSMatthew G. Knepley   PetscErrorCode ierr;
830cc08537eSMatthew G. Knepley 
831cc08537eSMatthew G. Knepley   PetscFunctionBegin;
832cc08537eSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
8330a1d6728SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr);
834cc08537eSMatthew G. Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
835cc08537eSMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
8360a1d6728SMatthew G. Knepley   dim  = coordSize/numCorners;
837011ea5d8SMatthew G. Knepley   if (normal) {
838011ea5d8SMatthew G. Knepley     if (dim > 2) {
8391ee9d5ecSMatthew G. Knepley       const PetscReal x0 = PetscRealPart(coords[dim+0] - coords[0]), x1 = PetscRealPart(coords[dim*2+0] - coords[0]);
8401ee9d5ecSMatthew G. Knepley       const PetscReal y0 = PetscRealPart(coords[dim+1] - coords[1]), y1 = PetscRealPart(coords[dim*2+1] - coords[1]);
8411ee9d5ecSMatthew G. Knepley       const PetscReal z0 = PetscRealPart(coords[dim+2] - coords[2]), z1 = PetscRealPart(coords[dim*2+2] - coords[2]);
8420a1d6728SMatthew G. Knepley       PetscReal       norm;
8430a1d6728SMatthew G. Knepley 
8441ee9d5ecSMatthew G. Knepley       v0[0]     = PetscRealPart(coords[0]);
8451ee9d5ecSMatthew G. Knepley       v0[1]     = PetscRealPart(coords[1]);
8461ee9d5ecSMatthew G. Knepley       v0[2]     = PetscRealPart(coords[2]);
8470a1d6728SMatthew G. Knepley       normal[0] = y0*z1 - z0*y1;
8480a1d6728SMatthew G. Knepley       normal[1] = z0*x1 - x0*z1;
8490a1d6728SMatthew G. Knepley       normal[2] = x0*y1 - y0*x1;
8500a1d6728SMatthew G. Knepley       norm = sqrt(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]);
8510a1d6728SMatthew G. Knepley       normal[0] /= norm;
8520a1d6728SMatthew G. Knepley       normal[1] /= norm;
8530a1d6728SMatthew G. Knepley       normal[2] /= norm;
854011ea5d8SMatthew G. Knepley     } else {
855011ea5d8SMatthew G. Knepley       for (d = 0; d < dim; ++d) normal[d] = 0.0;
856011ea5d8SMatthew G. Knepley     }
857011ea5d8SMatthew G. Knepley   }
85899dec3a6SMatthew G. Knepley   if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D_Internal(coordSize, coords, R);CHKERRQ(ierr);}
8590a1d6728SMatthew G. Knepley   for (p = 0; p < numCorners; ++p) {
8600a1d6728SMatthew G. Knepley     /* Need to do this copy to get types right */
8610a1d6728SMatthew G. Knepley     for (d = 0; d < tdim; ++d) {
8621ee9d5ecSMatthew G. Knepley       ctmp[d]      = PetscRealPart(coords[p*tdim+d]);
8631ee9d5ecSMatthew G. Knepley       ctmp[tdim+d] = PetscRealPart(coords[((p+1)%numCorners)*tdim+d]);
8640a1d6728SMatthew G. Knepley     }
8650a1d6728SMatthew G. Knepley     Volume_Triangle_Origin_Internal(&vtmp, ctmp);
8660a1d6728SMatthew G. Knepley     vsum += vtmp;
8670a1d6728SMatthew G. Knepley     for (d = 0; d < tdim; ++d) {
8680a1d6728SMatthew G. Knepley       csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp;
8690a1d6728SMatthew G. Knepley     }
8700a1d6728SMatthew G. Knepley   }
8710a1d6728SMatthew G. Knepley   for (d = 0; d < tdim; ++d) {
8720a1d6728SMatthew G. Knepley     csum[d] /= (tdim+1)*vsum;
8730a1d6728SMatthew G. Knepley   }
8740a1d6728SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
875ee6bbdb2SSatish Balay   if (vol) *vol = PetscAbsReal(vsum);
8760a1d6728SMatthew G. Knepley   if (centroid) {
8770a1d6728SMatthew G. Knepley     if (dim > 2) {
8780a1d6728SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
8790a1d6728SMatthew G. Knepley         centroid[d] = v0[d];
8800a1d6728SMatthew G. Knepley         for (e = 0; e < dim; ++e) {
8810a1d6728SMatthew G. Knepley           centroid[d] += R[d*dim+e]*csum[e];
8820a1d6728SMatthew G. Knepley         }
8830a1d6728SMatthew G. Knepley       }
8840a1d6728SMatthew G. Knepley     } else for (d = 0; d < dim; ++d) centroid[d] = csum[d];
8850a1d6728SMatthew G. Knepley   }
886cc08537eSMatthew G. Knepley   PetscFunctionReturn(0);
887cc08537eSMatthew G. Knepley }
888cc08537eSMatthew G. Knepley 
889cc08537eSMatthew G. Knepley #undef __FUNCT__
8900ec8681fSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_3D_Internal"
8910ec8681fSMatthew G. Knepley /* Centroid_i = (\sum_n V_n Cn_i ) / V */
892011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
8930ec8681fSMatthew G. Knepley {
8940ec8681fSMatthew G. Knepley   PetscSection    coordSection;
8950ec8681fSMatthew G. Knepley   Vec             coordinates;
8960ec8681fSMatthew G. Knepley   PetscScalar    *coords = NULL;
89786623015SMatthew G. Knepley   PetscReal       vsum = 0.0, vtmp, coordsTmp[3*3];
898*a7df9edeSMatthew G. Knepley   const PetscInt *faces, *facesO;
8990ec8681fSMatthew G. Knepley   PetscInt        numFaces, f, coordSize, numCorners, p, d;
9000ec8681fSMatthew G. Knepley   PetscErrorCode  ierr;
9010ec8681fSMatthew G. Knepley 
9020ec8681fSMatthew G. Knepley   PetscFunctionBegin;
9030ec8681fSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
9040ec8681fSMatthew G. Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
9050ec8681fSMatthew G. Knepley 
906d9a81ebdSMatthew G. Knepley   if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0;
9070ec8681fSMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr);
9080ec8681fSMatthew G. Knepley   ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr);
909*a7df9edeSMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr);
9100ec8681fSMatthew G. Knepley   for (f = 0; f < numFaces; ++f) {
911011ea5d8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr);
9120ec8681fSMatthew G. Knepley     numCorners = coordSize/dim;
9130ec8681fSMatthew G. Knepley     switch (numCorners) {
9140ec8681fSMatthew G. Knepley     case 3:
9150ec8681fSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
9161ee9d5ecSMatthew G. Knepley         coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]);
9171ee9d5ecSMatthew G. Knepley         coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]);
9181ee9d5ecSMatthew G. Knepley         coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]);
9190ec8681fSMatthew G. Knepley       }
9200ec8681fSMatthew G. Knepley       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
921*a7df9edeSMatthew G. Knepley       if (facesO[f] < 0) vtmp = -vtmp;
9220ec8681fSMatthew G. Knepley       vsum += vtmp;
9230ec8681fSMatthew G. Knepley       if (centroid) {
9240ec8681fSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
9251ee9d5ecSMatthew G. Knepley           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
9260ec8681fSMatthew G. Knepley         }
9270ec8681fSMatthew G. Knepley       }
9280ec8681fSMatthew G. Knepley       break;
9290ec8681fSMatthew G. Knepley     case 4:
9300ec8681fSMatthew G. Knepley       /* DO FOR PYRAMID */
9310ec8681fSMatthew G. Knepley       /* First tet */
9320ec8681fSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
9331ee9d5ecSMatthew G. Knepley         coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]);
9341ee9d5ecSMatthew G. Knepley         coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]);
9351ee9d5ecSMatthew G. Knepley         coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]);
9360ec8681fSMatthew G. Knepley       }
9370ec8681fSMatthew G. Knepley       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
938*a7df9edeSMatthew G. Knepley       if (facesO[f] < 0) vtmp = -vtmp;
9390ec8681fSMatthew G. Knepley       vsum += vtmp;
9400ec8681fSMatthew G. Knepley       if (centroid) {
9410ec8681fSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
9420ec8681fSMatthew G. Knepley           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
9430ec8681fSMatthew G. Knepley         }
9440ec8681fSMatthew G. Knepley       }
9450ec8681fSMatthew G. Knepley       /* Second tet */
9460ec8681fSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
9471ee9d5ecSMatthew G. Knepley         coordsTmp[0*dim+d] = PetscRealPart(coords[1*dim+d]);
9481ee9d5ecSMatthew G. Knepley         coordsTmp[1*dim+d] = PetscRealPart(coords[2*dim+d]);
9491ee9d5ecSMatthew G. Knepley         coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]);
9500ec8681fSMatthew G. Knepley       }
9510ec8681fSMatthew G. Knepley       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
952*a7df9edeSMatthew G. Knepley       if (facesO[f] < 0) vtmp = -vtmp;
9530ec8681fSMatthew G. Knepley       vsum += vtmp;
9540ec8681fSMatthew G. Knepley       if (centroid) {
9550ec8681fSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
9560ec8681fSMatthew G. Knepley           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
9570ec8681fSMatthew G. Knepley         }
9580ec8681fSMatthew G. Knepley       }
9590ec8681fSMatthew G. Knepley       break;
9600ec8681fSMatthew G. Knepley     default:
9610ec8681fSMatthew G. Knepley       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %d vertices", numCorners);
9620ec8681fSMatthew G. Knepley     }
9630ec8681fSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
9640ec8681fSMatthew G. Knepley   }
9658763be8eSMatthew G. Knepley   if (vol)     *vol = PetscAbsReal(vsum);
9660ec8681fSMatthew G. Knepley   if (normal)   for (d = 0; d < dim; ++d) normal[d]    = 0.0;
967d9a81ebdSMatthew G. Knepley   if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4);
9680ec8681fSMatthew G. Knepley   PetscFunctionReturn(0);
9690ec8681fSMatthew G. Knepley }
9700ec8681fSMatthew G. Knepley 
9710ec8681fSMatthew G. Knepley #undef __FUNCT__
972834e62ceSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeCellGeometryFVM"
973834e62ceSMatthew G. Knepley /*@C
974834e62ceSMatthew G. Knepley   DMPlexComputeCellGeometryFVM - Compute the volume for a given cell
975834e62ceSMatthew G. Knepley 
976834e62ceSMatthew G. Knepley   Collective on DM
977834e62ceSMatthew G. Knepley 
978834e62ceSMatthew G. Knepley   Input Arguments:
979834e62ceSMatthew G. Knepley + dm   - the DM
980834e62ceSMatthew G. Knepley - cell - the cell
981834e62ceSMatthew G. Knepley 
982834e62ceSMatthew G. Knepley   Output Arguments:
983834e62ceSMatthew G. Knepley + volume   - the cell volume
984cc08537eSMatthew G. Knepley . centroid - the cell centroid
985cc08537eSMatthew G. Knepley - normal - the cell normal, if appropriate
986834e62ceSMatthew G. Knepley 
987834e62ceSMatthew G. Knepley   Level: advanced
988834e62ceSMatthew G. Knepley 
989834e62ceSMatthew G. Knepley   Fortran Notes:
990834e62ceSMatthew G. Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
991834e62ceSMatthew G. Knepley   include petsc.h90 in your code.
992834e62ceSMatthew G. Knepley 
993834e62ceSMatthew G. Knepley .seealso: DMPlexGetCoordinateSection(), DMPlexGetCoordinateVec()
994834e62ceSMatthew G. Knepley @*/
995cc08537eSMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
996834e62ceSMatthew G. Knepley {
9970ec8681fSMatthew G. Knepley   PetscInt       depth, dim;
998834e62ceSMatthew G. Knepley   PetscErrorCode ierr;
999834e62ceSMatthew G. Knepley 
1000834e62ceSMatthew G. Knepley   PetscFunctionBegin;
1001834e62ceSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1002834e62ceSMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1003834e62ceSMatthew G. Knepley   if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated");
1004834e62ceSMatthew G. Knepley   /* We need to keep a pointer to the depth label */
1005011ea5d8SMatthew G. Knepley   ierr = DMPlexGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr);
1006834e62ceSMatthew G. Knepley   /* Cone size is now the number of faces */
1007011ea5d8SMatthew G. Knepley   switch (depth) {
1008cc08537eSMatthew G. Knepley   case 1:
1009011ea5d8SMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
1010cc08537eSMatthew G. Knepley     break;
1011834e62ceSMatthew G. Knepley   case 2:
1012011ea5d8SMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
1013834e62ceSMatthew G. Knepley     break;
1014834e62ceSMatthew G. Knepley   case 3:
1015011ea5d8SMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
1016834e62ceSMatthew G. Knepley     break;
1017834e62ceSMatthew G. Knepley   default:
1018834e62ceSMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
1019834e62ceSMatthew G. Knepley   }
1020834e62ceSMatthew G. Knepley   PetscFunctionReturn(0);
1021834e62ceSMatthew G. Knepley }
1022