xref: /petsc/src/dm/impls/plex/plexgeometry.c (revision c501906f0abd97a72ec10222ec0d348467dbcc70)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
29d150b73SToby Isaac #include <petsc/private/petscfeimpl.h>  /*I      "petscfe.h"       I*/
39d150b73SToby Isaac #include <petscblaslapack.h>
4ccd2543fSMatthew G Knepley 
5fea14342SMatthew G. Knepley static PetscErrorCode DMPlexGetLineIntersection_2D_Internal(const PetscReal segmentA[], const PetscReal segmentB[], PetscReal intersection[], PetscBool *hasIntersection)
6fea14342SMatthew G. Knepley {
7fea14342SMatthew G. Knepley   const PetscReal p0_x  = segmentA[0*2+0];
8fea14342SMatthew G. Knepley   const PetscReal p0_y  = segmentA[0*2+1];
9fea14342SMatthew G. Knepley   const PetscReal p1_x  = segmentA[1*2+0];
10fea14342SMatthew G. Knepley   const PetscReal p1_y  = segmentA[1*2+1];
11fea14342SMatthew G. Knepley   const PetscReal p2_x  = segmentB[0*2+0];
12fea14342SMatthew G. Knepley   const PetscReal p2_y  = segmentB[0*2+1];
13fea14342SMatthew G. Knepley   const PetscReal p3_x  = segmentB[1*2+0];
14fea14342SMatthew G. Knepley   const PetscReal p3_y  = segmentB[1*2+1];
15fea14342SMatthew G. Knepley   const PetscReal s1_x  = p1_x - p0_x;
16fea14342SMatthew G. Knepley   const PetscReal s1_y  = p1_y - p0_y;
17fea14342SMatthew G. Knepley   const PetscReal s2_x  = p3_x - p2_x;
18fea14342SMatthew G. Knepley   const PetscReal s2_y  = p3_y - p2_y;
19fea14342SMatthew G. Knepley   const PetscReal denom = (-s2_x * s1_y + s1_x * s2_y);
20fea14342SMatthew G. Knepley 
21fea14342SMatthew G. Knepley   PetscFunctionBegin;
22fea14342SMatthew G. Knepley   *hasIntersection = PETSC_FALSE;
23fea14342SMatthew G. Knepley   /* Non-parallel lines */
24fea14342SMatthew G. Knepley   if (denom != 0.0) {
25fea14342SMatthew G. Knepley     const PetscReal s = (-s1_y * (p0_x - p2_x) + s1_x * (p0_y - p2_y)) / denom;
26fea14342SMatthew G. Knepley     const PetscReal t = ( s2_x * (p0_y - p2_y) - s2_y * (p0_x - p2_x)) / denom;
27fea14342SMatthew G. Knepley 
28fea14342SMatthew G. Knepley     if (s >= 0 && s <= 1 && t >= 0 && t <= 1) {
29fea14342SMatthew G. Knepley       *hasIntersection = PETSC_TRUE;
30fea14342SMatthew G. Knepley       if (intersection) {
31fea14342SMatthew G. Knepley         intersection[0] = p0_x + (t * s1_x);
32fea14342SMatthew G. Knepley         intersection[1] = p0_y + (t * s1_y);
33fea14342SMatthew G. Knepley       }
34fea14342SMatthew G. Knepley     }
35fea14342SMatthew G. Knepley   }
36fea14342SMatthew G. Knepley   PetscFunctionReturn(0);
37fea14342SMatthew G. Knepley }
38fea14342SMatthew G. Knepley 
39ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
40ccd2543fSMatthew G Knepley {
41ccd2543fSMatthew G Knepley   const PetscInt  embedDim = 2;
42f5ebc837SMatthew G. Knepley   const PetscReal eps      = PETSC_SQRT_MACHINE_EPSILON;
43ccd2543fSMatthew G Knepley   PetscReal       x        = PetscRealPart(point[0]);
44ccd2543fSMatthew G Knepley   PetscReal       y        = PetscRealPart(point[1]);
45ccd2543fSMatthew G Knepley   PetscReal       v0[2], J[4], invJ[4], detJ;
46ccd2543fSMatthew G Knepley   PetscReal       xi, eta;
47ccd2543fSMatthew G Knepley   PetscErrorCode  ierr;
48ccd2543fSMatthew G Knepley 
49ccd2543fSMatthew G Knepley   PetscFunctionBegin;
508e0841e0SMatthew G. Knepley   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
51ccd2543fSMatthew G Knepley   xi  = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]);
52ccd2543fSMatthew G Knepley   eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]);
53ccd2543fSMatthew G Knepley 
54f5ebc837SMatthew G. Knepley   if ((xi >= -eps) && (eta >= -eps) && (xi + eta <= 2.0+eps)) *cell = c;
55c1496c66SMatthew G. Knepley   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
56ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
57ccd2543fSMatthew G Knepley }
58ccd2543fSMatthew G Knepley 
5962a38674SMatthew G. Knepley static PetscErrorCode DMPlexClosestPoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscReal cpoint[])
6062a38674SMatthew G. Knepley {
6162a38674SMatthew G. Knepley   const PetscInt  embedDim = 2;
6262a38674SMatthew G. Knepley   PetscReal       x        = PetscRealPart(point[0]);
6362a38674SMatthew G. Knepley   PetscReal       y        = PetscRealPart(point[1]);
6462a38674SMatthew G. Knepley   PetscReal       v0[2], J[4], invJ[4], detJ;
6562a38674SMatthew G. Knepley   PetscReal       xi, eta, r;
6662a38674SMatthew G. Knepley   PetscErrorCode  ierr;
6762a38674SMatthew G. Knepley 
6862a38674SMatthew G. Knepley   PetscFunctionBegin;
6962a38674SMatthew G. Knepley   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
7062a38674SMatthew G. Knepley   xi  = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]);
7162a38674SMatthew G. Knepley   eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]);
7262a38674SMatthew G. Knepley 
7362a38674SMatthew G. Knepley   xi  = PetscMax(xi,  0.0);
7462a38674SMatthew G. Knepley   eta = PetscMax(eta, 0.0);
7562a38674SMatthew G. Knepley   r   = (xi + eta)/2.0;
7662a38674SMatthew G. Knepley   if (xi + eta > 2.0) {
7762a38674SMatthew G. Knepley     r    = (xi + eta)/2.0;
7862a38674SMatthew G. Knepley     xi  /= r;
7962a38674SMatthew G. Knepley     eta /= r;
8062a38674SMatthew G. Knepley   }
8162a38674SMatthew G. Knepley   cpoint[0] = J[0*embedDim+0]*xi + J[0*embedDim+1]*eta + v0[0];
8262a38674SMatthew G. Knepley   cpoint[1] = J[1*embedDim+0]*xi + J[1*embedDim+1]*eta + v0[1];
8362a38674SMatthew G. Knepley   PetscFunctionReturn(0);
8462a38674SMatthew G. Knepley }
8562a38674SMatthew G. Knepley 
86ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_General_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
87ccd2543fSMatthew G Knepley {
88ccd2543fSMatthew G Knepley   PetscSection       coordSection;
89ccd2543fSMatthew G Knepley   Vec             coordsLocal;
90a1e44745SMatthew G. Knepley   PetscScalar    *coords = NULL;
91ccd2543fSMatthew G Knepley   const PetscInt  faces[8]  = {0, 1, 1, 2, 2, 3, 3, 0};
92ccd2543fSMatthew G Knepley   PetscReal       x         = PetscRealPart(point[0]);
93ccd2543fSMatthew G Knepley   PetscReal       y         = PetscRealPart(point[1]);
94ccd2543fSMatthew G Knepley   PetscInt        crossings = 0, f;
95ccd2543fSMatthew G Knepley   PetscErrorCode  ierr;
96ccd2543fSMatthew G Knepley 
97ccd2543fSMatthew G Knepley   PetscFunctionBegin;
98ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
9969d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
100ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
101ccd2543fSMatthew G Knepley   for (f = 0; f < 4; ++f) {
102ccd2543fSMatthew G Knepley     PetscReal x_i   = PetscRealPart(coords[faces[2*f+0]*2+0]);
103ccd2543fSMatthew G Knepley     PetscReal y_i   = PetscRealPart(coords[faces[2*f+0]*2+1]);
104ccd2543fSMatthew G Knepley     PetscReal x_j   = PetscRealPart(coords[faces[2*f+1]*2+0]);
105ccd2543fSMatthew G Knepley     PetscReal y_j   = PetscRealPart(coords[faces[2*f+1]*2+1]);
106ccd2543fSMatthew G Knepley     PetscReal slope = (y_j - y_i) / (x_j - x_i);
107ccd2543fSMatthew G Knepley     PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE;
108ccd2543fSMatthew G Knepley     PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE;
109ccd2543fSMatthew G Knepley     PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE;
110ccd2543fSMatthew G Knepley     if ((cond1 || cond2)  && above) ++crossings;
111ccd2543fSMatthew G Knepley   }
112ccd2543fSMatthew G Knepley   if (crossings % 2) *cell = c;
113c1496c66SMatthew G. Knepley   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
114ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
115ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
116ccd2543fSMatthew G Knepley }
117ccd2543fSMatthew G Knepley 
118ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
119ccd2543fSMatthew G Knepley {
120ccd2543fSMatthew G Knepley   const PetscInt embedDim = 3;
121ccd2543fSMatthew G Knepley   PetscReal      v0[3], J[9], invJ[9], detJ;
122ccd2543fSMatthew G Knepley   PetscReal      x = PetscRealPart(point[0]);
123ccd2543fSMatthew G Knepley   PetscReal      y = PetscRealPart(point[1]);
124ccd2543fSMatthew G Knepley   PetscReal      z = PetscRealPart(point[2]);
125ccd2543fSMatthew G Knepley   PetscReal      xi, eta, zeta;
126ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
127ccd2543fSMatthew G Knepley 
128ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1298e0841e0SMatthew G. Knepley   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
130ccd2543fSMatthew G Knepley   xi   = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]);
131ccd2543fSMatthew G Knepley   eta  = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]);
132ccd2543fSMatthew G Knepley   zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]);
133ccd2543fSMatthew G Knepley 
134ccd2543fSMatthew G Knepley   if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) *cell = c;
135c1496c66SMatthew G. Knepley   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
136ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
137ccd2543fSMatthew G Knepley }
138ccd2543fSMatthew G Knepley 
139ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
140ccd2543fSMatthew G Knepley {
141ccd2543fSMatthew G Knepley   PetscSection   coordSection;
142ccd2543fSMatthew G Knepley   Vec            coordsLocal;
1437c1f9639SMatthew G Knepley   PetscScalar   *coords;
144fb150da6SMatthew G. Knepley   const PetscInt faces[24] = {0, 3, 2, 1,  5, 4, 7, 6,  3, 0, 4, 5,
145fb150da6SMatthew G. Knepley                               1, 2, 6, 7,  3, 5, 6, 2,  0, 1, 7, 4};
146ccd2543fSMatthew G Knepley   PetscBool      found = PETSC_TRUE;
147ccd2543fSMatthew G Knepley   PetscInt       f;
148ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
149ccd2543fSMatthew G Knepley 
150ccd2543fSMatthew G Knepley   PetscFunctionBegin;
151ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
15269d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
153ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
154ccd2543fSMatthew G Knepley   for (f = 0; f < 6; ++f) {
155ccd2543fSMatthew G Knepley     /* Check the point is under plane */
156ccd2543fSMatthew G Knepley     /*   Get face normal */
157ccd2543fSMatthew G Knepley     PetscReal v_i[3];
158ccd2543fSMatthew G Knepley     PetscReal v_j[3];
159ccd2543fSMatthew G Knepley     PetscReal normal[3];
160ccd2543fSMatthew G Knepley     PetscReal pp[3];
161ccd2543fSMatthew G Knepley     PetscReal dot;
162ccd2543fSMatthew G Knepley 
163ccd2543fSMatthew G Knepley     v_i[0]    = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]);
164ccd2543fSMatthew G Knepley     v_i[1]    = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]);
165ccd2543fSMatthew G Knepley     v_i[2]    = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]);
166ccd2543fSMatthew G Knepley     v_j[0]    = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]);
167ccd2543fSMatthew G Knepley     v_j[1]    = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]);
168ccd2543fSMatthew G Knepley     v_j[2]    = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]);
169ccd2543fSMatthew G Knepley     normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1];
170ccd2543fSMatthew G Knepley     normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2];
171ccd2543fSMatthew G Knepley     normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0];
172ccd2543fSMatthew G Knepley     pp[0]     = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]);
173ccd2543fSMatthew G Knepley     pp[1]     = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]);
174ccd2543fSMatthew G Knepley     pp[2]     = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]);
175ccd2543fSMatthew G Knepley     dot       = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2];
176ccd2543fSMatthew G Knepley 
177ccd2543fSMatthew G Knepley     /* Check that projected point is in face (2D location problem) */
178ccd2543fSMatthew G Knepley     if (dot < 0.0) {
179ccd2543fSMatthew G Knepley       found = PETSC_FALSE;
180ccd2543fSMatthew G Knepley       break;
181ccd2543fSMatthew G Knepley     }
182ccd2543fSMatthew G Knepley   }
183ccd2543fSMatthew G Knepley   if (found) *cell = c;
184c1496c66SMatthew G. Knepley   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
185ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
186ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
187ccd2543fSMatthew G Knepley }
188ccd2543fSMatthew G Knepley 
189c4eade1cSMatthew G. Knepley static PetscErrorCode PetscGridHashInitialize_Internal(PetscGridHash box, PetscInt dim, const PetscScalar point[])
190c4eade1cSMatthew G. Knepley {
191c4eade1cSMatthew G. Knepley   PetscInt d;
192c4eade1cSMatthew G. Knepley 
193c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
194c4eade1cSMatthew G. Knepley   box->dim = dim;
195c4eade1cSMatthew G. Knepley   for (d = 0; d < dim; ++d) box->lower[d] = box->upper[d] = PetscRealPart(point[d]);
196c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
197c4eade1cSMatthew G. Knepley }
198c4eade1cSMatthew G. Knepley 
199c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashCreate(MPI_Comm comm, PetscInt dim, const PetscScalar point[], PetscGridHash *box)
200c4eade1cSMatthew G. Knepley {
201c4eade1cSMatthew G. Knepley   PetscErrorCode ierr;
202c4eade1cSMatthew G. Knepley 
203c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
204c4eade1cSMatthew G. Knepley   ierr = PetscMalloc1(1, box);CHKERRQ(ierr);
205c4eade1cSMatthew G. Knepley   ierr = PetscGridHashInitialize_Internal(*box, dim, point);CHKERRQ(ierr);
206c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
207c4eade1cSMatthew G. Knepley }
208c4eade1cSMatthew G. Knepley 
209c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashEnlarge(PetscGridHash box, const PetscScalar point[])
210c4eade1cSMatthew G. Knepley {
211c4eade1cSMatthew G. Knepley   PetscInt d;
212c4eade1cSMatthew G. Knepley 
213c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
214c4eade1cSMatthew G. Knepley   for (d = 0; d < box->dim; ++d) {
215c4eade1cSMatthew G. Knepley     box->lower[d] = PetscMin(box->lower[d], PetscRealPart(point[d]));
216c4eade1cSMatthew G. Knepley     box->upper[d] = PetscMax(box->upper[d], PetscRealPart(point[d]));
217c4eade1cSMatthew G. Knepley   }
218c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
219c4eade1cSMatthew G. Knepley }
220c4eade1cSMatthew G. Knepley 
22162a38674SMatthew G. Knepley /*
22262a38674SMatthew G. Knepley   PetscGridHashSetGrid - Divide the grid into boxes
22362a38674SMatthew G. Knepley 
22462a38674SMatthew G. Knepley   Not collective
22562a38674SMatthew G. Knepley 
22662a38674SMatthew G. Knepley   Input Parameters:
22762a38674SMatthew G. Knepley + box - The grid hash object
22862a38674SMatthew G. Knepley . n   - The number of boxes in each dimension, or PETSC_DETERMINE
22962a38674SMatthew G. Knepley - h   - The box size in each dimension, only used if n[d] == PETSC_DETERMINE
23062a38674SMatthew G. Knepley 
23162a38674SMatthew G. Knepley   Level: developer
23262a38674SMatthew G. Knepley 
23362a38674SMatthew G. Knepley .seealso: PetscGridHashCreate()
23462a38674SMatthew G. Knepley */
235c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashSetGrid(PetscGridHash box, const PetscInt n[], const PetscReal h[])
236c4eade1cSMatthew G. Knepley {
237c4eade1cSMatthew G. Knepley   PetscInt d;
238c4eade1cSMatthew G. Knepley 
239c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
240c4eade1cSMatthew G. Knepley   for (d = 0; d < box->dim; ++d) {
241c4eade1cSMatthew G. Knepley     box->extent[d] = box->upper[d] - box->lower[d];
242c4eade1cSMatthew G. Knepley     if (n[d] == PETSC_DETERMINE) {
243c4eade1cSMatthew G. Knepley       box->h[d] = h[d];
244c4eade1cSMatthew G. Knepley       box->n[d] = PetscCeilReal(box->extent[d]/h[d]);
245c4eade1cSMatthew G. Knepley     } else {
246c4eade1cSMatthew G. Knepley       box->n[d] = n[d];
247c4eade1cSMatthew G. Knepley       box->h[d] = box->extent[d]/n[d];
248c4eade1cSMatthew G. Knepley     }
249c4eade1cSMatthew G. Knepley   }
250c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
251c4eade1cSMatthew G. Knepley }
252c4eade1cSMatthew G. Knepley 
25362a38674SMatthew G. Knepley /*
25462a38674SMatthew G. Knepley   PetscGridHashGetEnclosingBox - Find the grid boxes containing each input point
25562a38674SMatthew G. Knepley 
25662a38674SMatthew G. Knepley   Not collective
25762a38674SMatthew G. Knepley 
25862a38674SMatthew G. Knepley   Input Parameters:
25962a38674SMatthew G. Knepley + box       - The grid hash object
26062a38674SMatthew G. Knepley . numPoints - The number of input points
26162a38674SMatthew G. Knepley - points    - The input point coordinates
26262a38674SMatthew G. Knepley 
26362a38674SMatthew G. Knepley   Output Parameters:
26462a38674SMatthew G. Knepley + dboxes    - An array of numPoints*dim integers expressing the enclosing box as (i_0, i_1, ..., i_dim)
26562a38674SMatthew G. Knepley - boxes     - An array of numPoints integers expressing the enclosing box as single number, or NULL
26662a38674SMatthew G. Knepley 
26762a38674SMatthew G. Knepley   Level: developer
26862a38674SMatthew G. Knepley 
26962a38674SMatthew G. Knepley .seealso: PetscGridHashCreate()
27062a38674SMatthew G. Knepley */
2711c6dfc3eSMatthew G. Knepley PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[])
272c4eade1cSMatthew G. Knepley {
273c4eade1cSMatthew G. Knepley   const PetscReal *lower = box->lower;
274c4eade1cSMatthew G. Knepley   const PetscReal *upper = box->upper;
275c4eade1cSMatthew G. Knepley   const PetscReal *h     = box->h;
276c4eade1cSMatthew G. Knepley   const PetscInt  *n     = box->n;
277c4eade1cSMatthew G. Knepley   const PetscInt   dim   = box->dim;
278c4eade1cSMatthew G. Knepley   PetscInt         d, p;
279c4eade1cSMatthew G. Knepley 
280c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
281c4eade1cSMatthew G. Knepley   for (p = 0; p < numPoints; ++p) {
282c4eade1cSMatthew G. Knepley     for (d = 0; d < dim; ++d) {
2831c6dfc3eSMatthew G. Knepley       PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]);
284c4eade1cSMatthew G. Knepley 
2851c6dfc3eSMatthew G. Knepley       if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1;
286c4eade1cSMatthew G. Knepley       if (dbox < 0 || dbox >= n[d]) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input point %d (%g, %g, %g) is outside of our bounding box",
2871c6dfc3eSMatthew G. Knepley                                              p, PetscRealPart(points[p*dim+0]), dim > 1 ? PetscRealPart(points[p*dim+1]) : 0.0, dim > 2 ? PetscRealPart(points[p*dim+2]) : 0.0);
288c4eade1cSMatthew G. Knepley       dboxes[p*dim+d] = dbox;
289c4eade1cSMatthew G. Knepley     }
290c4eade1cSMatthew G. Knepley     if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1];
291c4eade1cSMatthew G. Knepley   }
292c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
293c4eade1cSMatthew G. Knepley }
294c4eade1cSMatthew G. Knepley 
295c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashDestroy(PetscGridHash *box)
296c4eade1cSMatthew G. Knepley {
297c4eade1cSMatthew G. Knepley   PetscErrorCode ierr;
298c4eade1cSMatthew G. Knepley 
299c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
300c4eade1cSMatthew G. Knepley   if (*box) {
301c4eade1cSMatthew G. Knepley     ierr = PetscSectionDestroy(&(*box)->cellSection);CHKERRQ(ierr);
302c4eade1cSMatthew G. Knepley     ierr = ISDestroy(&(*box)->cells);CHKERRQ(ierr);
303c4eade1cSMatthew G. Knepley     ierr = DMLabelDestroy(&(*box)->cellsSparse);CHKERRQ(ierr);
304c4eade1cSMatthew G. Knepley   }
305c4eade1cSMatthew G. Knepley   ierr = PetscFree(*box);CHKERRQ(ierr);
306c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
307c4eade1cSMatthew G. Knepley }
308c4eade1cSMatthew G. Knepley 
309cafe43deSMatthew G. Knepley PetscErrorCode DMPlexLocatePoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cellStart, PetscInt *cell)
310cafe43deSMatthew G. Knepley {
311cafe43deSMatthew G. Knepley   PetscInt       coneSize;
312cafe43deSMatthew G. Knepley   PetscErrorCode ierr;
313cafe43deSMatthew G. Knepley 
314cafe43deSMatthew G. Knepley   PetscFunctionBegin;
315cafe43deSMatthew G. Knepley   switch (dim) {
316cafe43deSMatthew G. Knepley   case 2:
317cafe43deSMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr);
318cafe43deSMatthew G. Knepley     switch (coneSize) {
319cafe43deSMatthew G. Knepley     case 3:
320cafe43deSMatthew G. Knepley       ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
321cafe43deSMatthew G. Knepley       break;
322cafe43deSMatthew G. Knepley     case 4:
323cafe43deSMatthew G. Knepley       ierr = DMPlexLocatePoint_General_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
324cafe43deSMatthew G. Knepley       break;
325cafe43deSMatthew G. Knepley     default:
326cafe43deSMatthew G. Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize);
327cafe43deSMatthew G. Knepley     }
328cafe43deSMatthew G. Knepley     break;
329cafe43deSMatthew G. Knepley   case 3:
330cafe43deSMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr);
331cafe43deSMatthew G. Knepley     switch (coneSize) {
332cafe43deSMatthew G. Knepley     case 4:
333cafe43deSMatthew G. Knepley       ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
334cafe43deSMatthew G. Knepley       break;
335cafe43deSMatthew G. Knepley     case 6:
336cafe43deSMatthew G. Knepley       ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
337cafe43deSMatthew G. Knepley       break;
338cafe43deSMatthew G. Knepley     default:
339cafe43deSMatthew G. Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize);
340cafe43deSMatthew G. Knepley     }
341cafe43deSMatthew G. Knepley     break;
342cafe43deSMatthew G. Knepley   default:
343cafe43deSMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %D", dim);
344cafe43deSMatthew G. Knepley   }
345cafe43deSMatthew G. Knepley   PetscFunctionReturn(0);
346cafe43deSMatthew G. Knepley }
347cafe43deSMatthew G. Knepley 
34862a38674SMatthew G. Knepley /*
34962a38674SMatthew G. Knepley   DMPlexClosestPoint_Internal - Returns the closest point in the cell to the given point
35062a38674SMatthew G. Knepley */
35162a38674SMatthew G. Knepley PetscErrorCode DMPlexClosestPoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cell, PetscReal cpoint[])
35262a38674SMatthew G. Knepley {
35362a38674SMatthew G. Knepley   PetscInt       coneSize;
35462a38674SMatthew G. Knepley   PetscErrorCode ierr;
35562a38674SMatthew G. Knepley 
35662a38674SMatthew G. Knepley   PetscFunctionBegin;
35762a38674SMatthew G. Knepley   switch (dim) {
35862a38674SMatthew G. Knepley   case 2:
35962a38674SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
36062a38674SMatthew G. Knepley     switch (coneSize) {
36162a38674SMatthew G. Knepley     case 3:
36262a38674SMatthew G. Knepley       ierr = DMPlexClosestPoint_Simplex_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
36362a38674SMatthew G. Knepley       break;
36462a38674SMatthew G. Knepley #if 0
36562a38674SMatthew G. Knepley     case 4:
36662a38674SMatthew G. Knepley       ierr = DMPlexClosestPoint_General_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
36762a38674SMatthew G. Knepley       break;
36862a38674SMatthew G. Knepley #endif
36962a38674SMatthew G. Knepley     default:
37062a38674SMatthew G. Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize);
37162a38674SMatthew G. Knepley     }
37262a38674SMatthew G. Knepley     break;
37362a38674SMatthew G. Knepley #if 0
37462a38674SMatthew G. Knepley   case 3:
37562a38674SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
37662a38674SMatthew G. Knepley     switch (coneSize) {
37762a38674SMatthew G. Knepley     case 4:
37862a38674SMatthew G. Knepley       ierr = DMPlexClosestPoint_Simplex_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
37962a38674SMatthew G. Knepley       break;
38062a38674SMatthew G. Knepley     case 6:
38162a38674SMatthew G. Knepley       ierr = DMPlexClosestPoint_General_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
38262a38674SMatthew G. Knepley       break;
38362a38674SMatthew G. Knepley     default:
38462a38674SMatthew G. Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize);
38562a38674SMatthew G. Knepley     }
38662a38674SMatthew G. Knepley     break;
38762a38674SMatthew G. Knepley #endif
38862a38674SMatthew G. Knepley   default:
38962a38674SMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for mesh dimension %D", dim);
39062a38674SMatthew G. Knepley   }
39162a38674SMatthew G. Knepley   PetscFunctionReturn(0);
39262a38674SMatthew G. Knepley }
39362a38674SMatthew G. Knepley 
39462a38674SMatthew G. Knepley /*
39562a38674SMatthew G. Knepley   DMPlexComputeGridHash_Internal - Create a grid hash structure covering the Plex
39662a38674SMatthew G. Knepley 
39762a38674SMatthew G. Knepley   Collective on DM
39862a38674SMatthew G. Knepley 
39962a38674SMatthew G. Knepley   Input Parameter:
40062a38674SMatthew G. Knepley . dm - The Plex
40162a38674SMatthew G. Knepley 
40262a38674SMatthew G. Knepley   Output Parameter:
40362a38674SMatthew G. Knepley . localBox - The grid hash object
40462a38674SMatthew G. Knepley 
40562a38674SMatthew G. Knepley   Level: developer
40662a38674SMatthew G. Knepley 
40762a38674SMatthew G. Knepley .seealso: PetscGridHashCreate(), PetscGridHashGetEnclosingBox()
40862a38674SMatthew G. Knepley */
409cafe43deSMatthew G. Knepley PetscErrorCode DMPlexComputeGridHash_Internal(DM dm, PetscGridHash *localBox)
410cafe43deSMatthew G. Knepley {
411cafe43deSMatthew G. Knepley   MPI_Comm           comm;
412cafe43deSMatthew G. Knepley   PetscGridHash      lbox;
413cafe43deSMatthew G. Knepley   Vec                coordinates;
414cafe43deSMatthew G. Knepley   PetscSection       coordSection;
415cafe43deSMatthew G. Knepley   Vec                coordsLocal;
416cafe43deSMatthew G. Knepley   const PetscScalar *coords;
417722d0f5cSMatthew G. Knepley   PetscInt          *dboxes, *boxes;
418cafe43deSMatthew G. Knepley   PetscInt           n[3] = {10, 10, 10};
4191d0c6c94SMatthew G. Knepley   PetscInt           dim, N, cStart, cEnd, cMax, c, i;
420cafe43deSMatthew G. Knepley   PetscErrorCode     ierr;
421cafe43deSMatthew G. Knepley 
422cafe43deSMatthew G. Knepley   PetscFunctionBegin;
423cafe43deSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
424cafe43deSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
425cafe43deSMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
4265b3353d8SMatthew G. Knepley   if (dim != 2) SETERRQ(comm, PETSC_ERR_SUP, "I have only coded this for 2D");
427cafe43deSMatthew G. Knepley   ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr);
428cafe43deSMatthew G. Knepley   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
429cafe43deSMatthew G. Knepley   ierr = PetscGridHashCreate(comm, dim, coords, &lbox);CHKERRQ(ierr);
430cafe43deSMatthew G. Knepley   for (i = 0; i < N; i += dim) {ierr = PetscGridHashEnlarge(lbox, &coords[i]);CHKERRQ(ierr);}
431cafe43deSMatthew G. Knepley   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
432cafe43deSMatthew G. Knepley   ierr = PetscGridHashSetGrid(lbox, n, NULL);CHKERRQ(ierr);
433cafe43deSMatthew G. Knepley #if 0
434cafe43deSMatthew G. Knepley   /* Could define a custom reduction to merge these */
435b2566f29SBarry Smith   ierr = MPIU_Allreduce(lbox->lower, gbox->lower, 3, MPIU_REAL, MPI_MIN, comm);CHKERRQ(ierr);
436b2566f29SBarry Smith   ierr = MPIU_Allreduce(lbox->upper, gbox->upper, 3, MPIU_REAL, MPI_MAX, comm);CHKERRQ(ierr);
437cafe43deSMatthew G. Knepley #endif
438cafe43deSMatthew G. Knepley   /* Is there a reason to snap the local bounding box to a division of the global box? */
439cafe43deSMatthew G. Knepley   /* Should we compute all overlaps of local boxes? We could do this with a rendevouz scheme partitioning the global box */
440cafe43deSMatthew G. Knepley   /* Create label */
441cafe43deSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
4421d0c6c94SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
4431d0c6c94SMatthew G. Knepley   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
444cafe43deSMatthew G. Knepley   ierr = DMLabelCreate("cells", &lbox->cellsSparse);CHKERRQ(ierr);
445cafe43deSMatthew G. Knepley   ierr = DMLabelCreateIndex(lbox->cellsSparse, cStart, cEnd);CHKERRQ(ierr);
446722d0f5cSMatthew G. Knepley   /* Compute boxes which overlap each cell: http://stackoverflow.com/questions/13790208/triangle-square-intersection-test-in-2d */
447cafe43deSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
448cafe43deSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
44938353de4SMatthew G. Knepley   ierr = PetscCalloc2(16 * dim, &dboxes, 16, &boxes);CHKERRQ(ierr);
450cafe43deSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
451cafe43deSMatthew G. Knepley     const PetscReal *h       = lbox->h;
452cafe43deSMatthew G. Knepley     PetscScalar     *ccoords = NULL;
45338353de4SMatthew G. Knepley     PetscInt         csize   = 0;
454cafe43deSMatthew G. Knepley     PetscScalar      point[3];
455cafe43deSMatthew G. Knepley     PetscInt         dlim[6], d, e, i, j, k;
456cafe43deSMatthew G. Knepley 
457cafe43deSMatthew G. Knepley     /* Find boxes enclosing each vertex */
45838353de4SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, &csize, &ccoords);CHKERRQ(ierr);
45938353de4SMatthew G. Knepley     ierr = PetscGridHashGetEnclosingBox(lbox, csize/dim, ccoords, dboxes, boxes);CHKERRQ(ierr);
460722d0f5cSMatthew G. Knepley     /* Mark cells containing the vertices */
46138353de4SMatthew G. Knepley     for (e = 0; e < csize/dim; ++e) {ierr = DMLabelSetValue(lbox->cellsSparse, c, boxes[e]);CHKERRQ(ierr);}
462cafe43deSMatthew G. Knepley     /* Get grid of boxes containing these */
463cafe43deSMatthew G. Knepley     for (d = 0;   d < dim; ++d) {dlim[d*2+0] = dlim[d*2+1] = dboxes[d];}
4642291669eSMatthew G. Knepley     for (d = dim; d < 3;   ++d) {dlim[d*2+0] = dlim[d*2+1] = 0;}
465cafe43deSMatthew G. Knepley     for (e = 1; e < dim+1; ++e) {
466cafe43deSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
467cafe43deSMatthew G. Knepley         dlim[d*2+0] = PetscMin(dlim[d*2+0], dboxes[e*dim+d]);
468cafe43deSMatthew G. Knepley         dlim[d*2+1] = PetscMax(dlim[d*2+1], dboxes[e*dim+d]);
469cafe43deSMatthew G. Knepley       }
470cafe43deSMatthew G. Knepley     }
471fea14342SMatthew G. Knepley     /* Check for intersection of box with cell */
472cafe43deSMatthew G. Knepley     for (k = dlim[2*2+0], point[2] = lbox->lower[2] + k*h[2]; k <= dlim[2*2+1]; ++k, point[2] += h[2]) {
473cafe43deSMatthew G. Knepley       for (j = dlim[1*2+0], point[1] = lbox->lower[1] + j*h[1]; j <= dlim[1*2+1]; ++j, point[1] += h[1]) {
474cafe43deSMatthew G. Knepley         for (i = dlim[0*2+0], point[0] = lbox->lower[0] + i*h[0]; i <= dlim[0*2+1]; ++i, point[0] += h[0]) {
475cafe43deSMatthew G. Knepley           const PetscInt box = (k*lbox->n[1] + j)*lbox->n[0] + i;
476cafe43deSMatthew G. Knepley           PetscScalar    cpoint[3];
477fea14342SMatthew G. Knepley           PetscInt       cell, edge, ii, jj, kk;
478cafe43deSMatthew G. Knepley 
479fea14342SMatthew G. Knepley           /* Check whether cell contains any vertex of these subboxes TODO vectorize this */
480cafe43deSMatthew G. Knepley           for (kk = 0, cpoint[2] = point[2]; kk < (dim > 2 ? 2 : 1); ++kk, cpoint[2] += h[2]) {
481cafe43deSMatthew G. Knepley             for (jj = 0, cpoint[1] = point[1]; jj < (dim > 1 ? 2 : 1); ++jj, cpoint[1] += h[1]) {
482cafe43deSMatthew G. Knepley               for (ii = 0, cpoint[0] = point[0]; ii < 2; ++ii, cpoint[0] += h[0]) {
483cafe43deSMatthew G. Knepley 
484cafe43deSMatthew G. Knepley                 ierr = DMPlexLocatePoint_Internal(dm, dim, cpoint, c, &cell);CHKERRQ(ierr);
485cafe43deSMatthew G. Knepley                 if (cell >= 0) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); ii = jj = kk = 2;}
486cafe43deSMatthew G. Knepley               }
487cafe43deSMatthew G. Knepley             }
488cafe43deSMatthew G. Knepley           }
489fea14342SMatthew G. Knepley           /* Check whether cell edge intersects any edge of these subboxes TODO vectorize this */
490fea14342SMatthew G. Knepley           for (edge = 0; edge < dim+1; ++edge) {
491fea14342SMatthew G. Knepley             PetscReal segA[6], segB[6];
492fea14342SMatthew G. Knepley 
493fea14342SMatthew G. Knepley             for (d = 0; d < dim; ++d) {segA[d] = PetscRealPart(ccoords[edge*dim+d]); segA[dim+d] = PetscRealPart(ccoords[((edge+1)%(dim+1))*dim+d]);}
494fea14342SMatthew G. Knepley             for (kk = 0; kk < (dim > 2 ? 2 : 1); ++kk) {
4959a128ed2SMatthew G. Knepley               if (dim > 2) {segB[2]     = PetscRealPart(point[2]);
4969a128ed2SMatthew G. Knepley                             segB[dim+2] = PetscRealPart(point[2]) + kk*h[2];}
497fea14342SMatthew G. Knepley               for (jj = 0; jj < (dim > 1 ? 2 : 1); ++jj) {
4989a128ed2SMatthew G. Knepley                 if (dim > 1) {segB[1]     = PetscRealPart(point[1]);
4999a128ed2SMatthew G. Knepley                               segB[dim+1] = PetscRealPart(point[1]) + jj*h[1];}
500fea14342SMatthew G. Knepley                 for (ii = 0; ii < 2; ++ii) {
501fea14342SMatthew G. Knepley                   PetscBool intersects;
502fea14342SMatthew G. Knepley 
5039a128ed2SMatthew G. Knepley                   segB[0]     = PetscRealPart(point[0]);
5049a128ed2SMatthew G. Knepley                   segB[dim+0] = PetscRealPart(point[0]) + ii*h[0];
505fea14342SMatthew G. Knepley                   ierr = DMPlexGetLineIntersection_2D_Internal(segA, segB, NULL, &intersects);CHKERRQ(ierr);
506fea14342SMatthew G. Knepley                   if (intersects) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); edge = ii = jj = kk = dim+1;}
507cafe43deSMatthew G. Knepley                 }
508cafe43deSMatthew G. Knepley               }
509cafe43deSMatthew G. Knepley             }
510cafe43deSMatthew G. Knepley           }
511fea14342SMatthew G. Knepley         }
512fea14342SMatthew G. Knepley       }
513fea14342SMatthew G. Knepley     }
514fea14342SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &ccoords);CHKERRQ(ierr);
515fea14342SMatthew G. Knepley   }
516722d0f5cSMatthew G. Knepley   ierr = PetscFree2(dboxes, boxes);CHKERRQ(ierr);
517cafe43deSMatthew G. Knepley   ierr = DMLabelConvertToSection(lbox->cellsSparse, &lbox->cellSection, &lbox->cells);CHKERRQ(ierr);
518cafe43deSMatthew G. Knepley   ierr = DMLabelDestroy(&lbox->cellsSparse);CHKERRQ(ierr);
519cafe43deSMatthew G. Knepley   *localBox = lbox;
520cafe43deSMatthew G. Knepley   PetscFunctionReturn(0);
521cafe43deSMatthew G. Knepley }
522cafe43deSMatthew G. Knepley 
52362a38674SMatthew G. Knepley PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, DMPointLocationType ltype, PetscSF cellSF)
524ccd2543fSMatthew G Knepley {
525cafe43deSMatthew G. Knepley   DM_Plex        *mesh = (DM_Plex *) dm->data;
526953fc75cSMatthew G. Knepley   PetscBool       hash = mesh->useHashLocation;
5273a93e3b7SToby Isaac   PetscInt        bs, numPoints, p, numFound, *found = NULL;
5281318edbeSMatthew G. Knepley   PetscInt        dim, cStart, cEnd, cMax, numCells, c;
529cafe43deSMatthew G. Knepley   const PetscInt *boxCells;
5303a93e3b7SToby Isaac   PetscSFNode    *cells;
531ccd2543fSMatthew G Knepley   PetscScalar    *a;
5323a93e3b7SToby Isaac   PetscMPIInt     result;
533ccd2543fSMatthew G Knepley   PetscErrorCode  ierr;
534ccd2543fSMatthew G Knepley 
535ccd2543fSMatthew G Knepley   PetscFunctionBegin;
536080342d1SMatthew G. Knepley   if (ltype == DM_POINTLOCATION_NEAREST && !hash) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Nearest point location only supported with grid hashing. Use -dm_plex_hash_location to enable it.");
537cafe43deSMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
538cafe43deSMatthew G. Knepley   ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr);
5393a93e3b7SToby Isaac   ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)cellSF),PETSC_COMM_SELF,&result);CHKERRQ(ierr);
5403a93e3b7SToby Isaac   if (result != MPI_IDENT && result != MPI_CONGRUENT) SETERRQ(PetscObjectComm((PetscObject)cellSF),PETSC_ERR_SUP, "Trying parallel point location: only local point location supported");
541cafe43deSMatthew 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);
542ccd2543fSMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
543ccd2543fSMatthew G Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
544ccd2543fSMatthew G Knepley   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
545ccd2543fSMatthew G Knepley   ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr);
546ccd2543fSMatthew G Knepley   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
547ccd2543fSMatthew G Knepley   numPoints /= bs;
548785e854fSJed Brown   ierr = PetscMalloc1(numPoints, &cells);CHKERRQ(ierr);
549953fc75cSMatthew G. Knepley   if (hash) {
550ac6ec2abSMatthew G. Knepley     if (!mesh->lbox) {ierr = PetscInfo(dm, "Initializing grid hashing");CHKERRQ(ierr);ierr = DMPlexComputeGridHash_Internal(dm, &mesh->lbox);CHKERRQ(ierr);}
551cafe43deSMatthew G. Knepley     /* Designate the local box for each point */
552cafe43deSMatthew G. Knepley     /* Send points to correct process */
553cafe43deSMatthew G. Knepley     /* Search cells that lie in each subbox */
554cafe43deSMatthew G. Knepley     /*   Should we bin points before doing search? */
555cafe43deSMatthew G. Knepley     ierr = ISGetIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr);
556953fc75cSMatthew G. Knepley   }
5573a93e3b7SToby Isaac   for (p = 0, numFound = 0; p < numPoints; ++p) {
558ccd2543fSMatthew G Knepley     const PetscScalar *point = &a[p*bs];
559953fc75cSMatthew G. Knepley     PetscInt           dbin[3], bin, cell = -1, cellOffset;
560ccd2543fSMatthew G Knepley 
561e9b685f5SMatthew G. Knepley     cells[p].rank  = 0;
562e9b685f5SMatthew G. Knepley     cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
563953fc75cSMatthew G. Knepley     if (hash) {
564cafe43deSMatthew G. Knepley       ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr);
565cafe43deSMatthew G. Knepley       /* TODO Lay an interface over this so we can switch between Section (dense) and Label (sparse) */
566cafe43deSMatthew G. Knepley       ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr);
567cafe43deSMatthew G. Knepley       ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr);
568cafe43deSMatthew G. Knepley       for (c = cellOffset; c < cellOffset + numCells; ++c) {
569cafe43deSMatthew G. Knepley         ierr = DMPlexLocatePoint_Internal(dm, dim, point, boxCells[c], &cell);CHKERRQ(ierr);
5703a93e3b7SToby Isaac         if (cell >= 0) {
5713a93e3b7SToby Isaac           cells[p].rank = 0;
5723a93e3b7SToby Isaac           cells[p].index = cell;
5733a93e3b7SToby Isaac           numFound++;
5743a93e3b7SToby Isaac           break;
575ccd2543fSMatthew G Knepley         }
5763a93e3b7SToby Isaac       }
577953fc75cSMatthew G. Knepley     } else {
578953fc75cSMatthew G. Knepley       for (c = cStart; c < cEnd; ++c) {
579953fc75cSMatthew G. Knepley         ierr = DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);CHKERRQ(ierr);
5803a93e3b7SToby Isaac         if (cell >= 0) {
5813a93e3b7SToby Isaac           cells[p].rank = 0;
5823a93e3b7SToby Isaac           cells[p].index = cell;
5833a93e3b7SToby Isaac           numFound++;
5843a93e3b7SToby Isaac           break;
585953fc75cSMatthew G. Knepley         }
586953fc75cSMatthew G. Knepley       }
5873a93e3b7SToby Isaac     }
588ccd2543fSMatthew G Knepley   }
589953fc75cSMatthew G. Knepley   if (hash) {ierr = ISRestoreIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr);}
59062a38674SMatthew G. Knepley   if (ltype == DM_POINTLOCATION_NEAREST && hash && numFound < numPoints) {
59162a38674SMatthew G. Knepley     for (p = 0; p < numPoints; p++) {
59262a38674SMatthew G. Knepley       const PetscScalar *point = &a[p*bs];
59362a38674SMatthew G. Knepley       PetscReal          cpoint[3], diff[3], dist, distMax = PETSC_MAX_REAL;
594b716b415SMatthew G. Knepley       PetscInt           dbin[3], bin, cellOffset, d;
59562a38674SMatthew G. Knepley 
596e9b685f5SMatthew G. Knepley       if (cells[p].index < 0) {
59762a38674SMatthew G. Knepley         ++numFound;
59862a38674SMatthew G. Knepley         ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr);
59962a38674SMatthew G. Knepley         ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr);
60062a38674SMatthew G. Knepley         ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr);
60162a38674SMatthew G. Knepley         for (c = cellOffset; c < cellOffset + numCells; ++c) {
60262a38674SMatthew G. Knepley           ierr = DMPlexClosestPoint_Internal(dm, dim, point, boxCells[c], cpoint);CHKERRQ(ierr);
603b716b415SMatthew G. Knepley           for (d = 0; d < dim; ++d) diff[d] = cpoint[d] - PetscRealPart(point[d]);
60462a38674SMatthew G. Knepley           dist = DMPlex_NormD_Internal(dim, diff);
60562a38674SMatthew G. Knepley           if (dist < distMax) {
60662a38674SMatthew G. Knepley             for (d = 0; d < dim; ++d) a[p*bs+d] = cpoint[d];
60762a38674SMatthew G. Knepley             cells[p].rank  = 0;
60862a38674SMatthew G. Knepley             cells[p].index = boxCells[c];
60962a38674SMatthew G. Knepley             distMax = dist;
61062a38674SMatthew G. Knepley             break;
61162a38674SMatthew G. Knepley           }
61262a38674SMatthew G. Knepley         }
61362a38674SMatthew G. Knepley       }
61462a38674SMatthew G. Knepley     }
61562a38674SMatthew G. Knepley   }
61662a38674SMatthew G. Knepley   /* This code is only be relevant when interfaced to parallel point location */
617cafe43deSMatthew G. Knepley   /* Check for highest numbered proc that claims a point (do we care?) */
6182d1fa6caSMatthew G. Knepley   if (ltype == DM_POINTLOCATION_REMOVE && numFound < numPoints) {
6193a93e3b7SToby Isaac     ierr = PetscMalloc1(numFound,&found);CHKERRQ(ierr);
6203a93e3b7SToby Isaac     for (p = 0, numFound = 0; p < numPoints; p++) {
6213a93e3b7SToby Isaac       if (cells[p].rank >= 0 && cells[p].index >= 0) {
6223a93e3b7SToby Isaac         if (numFound < p) {
6233a93e3b7SToby Isaac           cells[numFound] = cells[p];
6243a93e3b7SToby Isaac         }
6253a93e3b7SToby Isaac         found[numFound++] = p;
6263a93e3b7SToby Isaac       }
6273a93e3b7SToby Isaac     }
6283a93e3b7SToby Isaac   }
62962a38674SMatthew G. Knepley   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
6303a93e3b7SToby Isaac   ierr = PetscSFSetGraph(cellSF, cEnd - cStart, numFound, found, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER);CHKERRQ(ierr);
631ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
632ccd2543fSMatthew G Knepley }
633ccd2543fSMatthew G Knepley 
634741bfc07SMatthew G. Knepley /*@C
635741bfc07SMatthew G. Knepley   DMPlexComputeProjection2Dto1D - Rewrite coordinates to be the 1D projection of the 2D coordinates
636741bfc07SMatthew G. Knepley 
637741bfc07SMatthew G. Knepley   Not collective
638741bfc07SMatthew G. Knepley 
639741bfc07SMatthew G. Knepley   Input Parameter:
640741bfc07SMatthew G. Knepley . coords - The coordinates of a segment
641741bfc07SMatthew G. Knepley 
642741bfc07SMatthew G. Knepley   Output Parameters:
643741bfc07SMatthew G. Knepley + coords - The new y-coordinate, and 0 for x
644741bfc07SMatthew G. Knepley - R - The rotation which accomplishes the projection
645741bfc07SMatthew G. Knepley 
646741bfc07SMatthew G. Knepley   Level: developer
647741bfc07SMatthew G. Knepley 
648741bfc07SMatthew G. Knepley .seealso: DMPlexComputeProjection3Dto1D(), DMPlexComputeProjection3Dto2D()
649741bfc07SMatthew G. Knepley @*/
650741bfc07SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar coords[], PetscReal R[])
65117fe8556SMatthew G. Knepley {
65217fe8556SMatthew G. Knepley   const PetscReal x = PetscRealPart(coords[2] - coords[0]);
65317fe8556SMatthew G. Knepley   const PetscReal y = PetscRealPart(coords[3] - coords[1]);
6548b49ba18SBarry Smith   const PetscReal r = PetscSqrtReal(x*x + y*y), c = x/r, s = y/r;
65517fe8556SMatthew G. Knepley 
65617fe8556SMatthew G. Knepley   PetscFunctionBegin;
6571c99cf0cSGeoffrey Irving   R[0] = c; R[1] = -s;
6581c99cf0cSGeoffrey Irving   R[2] = s; R[3] =  c;
65917fe8556SMatthew G. Knepley   coords[0] = 0.0;
6607f07f362SMatthew G. Knepley   coords[1] = r;
66117fe8556SMatthew G. Knepley   PetscFunctionReturn(0);
66217fe8556SMatthew G. Knepley }
66317fe8556SMatthew G. Knepley 
664741bfc07SMatthew G. Knepley /*@C
665741bfc07SMatthew G. Knepley   DMPlexComputeProjection3Dto1D - Rewrite coordinates to be the 1D projection of the 3D coordinates
66628dbe442SToby Isaac 
667741bfc07SMatthew G. Knepley   Not collective
66828dbe442SToby Isaac 
669741bfc07SMatthew G. Knepley   Input Parameter:
670741bfc07SMatthew G. Knepley . coords - The coordinates of a segment
671741bfc07SMatthew G. Knepley 
672741bfc07SMatthew G. Knepley   Output Parameters:
673741bfc07SMatthew G. Knepley + coords - The new y-coordinate, and 0 for x and z
674741bfc07SMatthew G. Knepley - R - The rotation which accomplishes the projection
675741bfc07SMatthew G. Knepley 
676741bfc07SMatthew G. Knepley   Note: This uses the basis completion described by Frisvad in http://www.imm.dtu.dk/~jerf/papers/abstracts/onb.html, DOI:10.1080/2165347X.2012.689606
677741bfc07SMatthew G. Knepley 
678741bfc07SMatthew G. Knepley   Level: developer
679741bfc07SMatthew G. Knepley 
680741bfc07SMatthew G. Knepley .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto2D()
681741bfc07SMatthew G. Knepley @*/
682741bfc07SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar coords[], PetscReal R[])
68328dbe442SToby Isaac {
68428dbe442SToby Isaac   PetscReal      x    = PetscRealPart(coords[3] - coords[0]);
68528dbe442SToby Isaac   PetscReal      y    = PetscRealPart(coords[4] - coords[1]);
68628dbe442SToby Isaac   PetscReal      z    = PetscRealPart(coords[5] - coords[2]);
68728dbe442SToby Isaac   PetscReal      r    = PetscSqrtReal(x*x + y*y + z*z);
68828dbe442SToby Isaac   PetscReal      rinv = 1. / r;
68928dbe442SToby Isaac   PetscFunctionBegin;
69028dbe442SToby Isaac 
69128dbe442SToby Isaac   x *= rinv; y *= rinv; z *= rinv;
69228dbe442SToby Isaac   if (x > 0.) {
69328dbe442SToby Isaac     PetscReal inv1pX   = 1./ (1. + x);
69428dbe442SToby Isaac 
69528dbe442SToby Isaac     R[0] = x; R[1] = -y;              R[2] = -z;
69628dbe442SToby Isaac     R[3] = y; R[4] = 1. - y*y*inv1pX; R[5] =     -y*z*inv1pX;
69728dbe442SToby Isaac     R[6] = z; R[7] =     -y*z*inv1pX; R[8] = 1. - z*z*inv1pX;
69828dbe442SToby Isaac   }
69928dbe442SToby Isaac   else {
70028dbe442SToby Isaac     PetscReal inv1mX   = 1./ (1. - x);
70128dbe442SToby Isaac 
70228dbe442SToby Isaac     R[0] = x; R[1] = z;               R[2] = y;
70328dbe442SToby Isaac     R[3] = y; R[4] =     -y*z*inv1mX; R[5] = 1. - y*y*inv1mX;
70428dbe442SToby Isaac     R[6] = z; R[7] = 1. - z*z*inv1mX; R[8] =     -y*z*inv1mX;
70528dbe442SToby Isaac   }
70628dbe442SToby Isaac   coords[0] = 0.0;
70728dbe442SToby Isaac   coords[1] = r;
70828dbe442SToby Isaac   PetscFunctionReturn(0);
70928dbe442SToby Isaac }
71028dbe442SToby Isaac 
711741bfc07SMatthew G. Knepley /*@
712741bfc07SMatthew G. Knepley   DMPlexComputeProjection3Dto2D - Rewrite coordinates to be the 2D projection of the 3D coordinates
713741bfc07SMatthew G. Knepley 
714741bfc07SMatthew G. Knepley   Not collective
715741bfc07SMatthew G. Knepley 
716741bfc07SMatthew G. Knepley   Input Parameter:
717741bfc07SMatthew G. Knepley . coords - The coordinates of a segment
718741bfc07SMatthew G. Knepley 
719741bfc07SMatthew G. Knepley   Output Parameters:
720741bfc07SMatthew G. Knepley + coords - The new y- and z-coordinates, and 0 for x
721741bfc07SMatthew G. Knepley - R - The rotation which accomplishes the projection
722741bfc07SMatthew G. Knepley 
723741bfc07SMatthew G. Knepley   Level: developer
724741bfc07SMatthew G. Knepley 
725741bfc07SMatthew G. Knepley .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto1D()
726741bfc07SMatthew G. Knepley @*/
727741bfc07SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt coordSize, PetscScalar coords[], PetscReal R[])
728ccd2543fSMatthew G Knepley {
7291ee9d5ecSMatthew G. Knepley   PetscReal      x1[3],  x2[3], n[3], norm;
73099dec3a6SMatthew G. Knepley   PetscReal      x1p[3], x2p[3], xnp[3];
7314a217a95SMatthew G. Knepley   PetscReal      sqrtz, alpha;
732ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
73399dec3a6SMatthew G. Knepley   PetscInt       d, e, p;
734ccd2543fSMatthew G Knepley 
735ccd2543fSMatthew G Knepley   PetscFunctionBegin;
736ccd2543fSMatthew G Knepley   /* 0) Calculate normal vector */
737ccd2543fSMatthew G Knepley   for (d = 0; d < dim; ++d) {
7381ee9d5ecSMatthew G. Knepley     x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]);
7391ee9d5ecSMatthew G. Knepley     x2[d] = PetscRealPart(coords[2*dim+d] - coords[0*dim+d]);
740ccd2543fSMatthew G Knepley   }
741ccd2543fSMatthew G Knepley   n[0] = x1[1]*x2[2] - x1[2]*x2[1];
742ccd2543fSMatthew G Knepley   n[1] = x1[2]*x2[0] - x1[0]*x2[2];
743ccd2543fSMatthew G Knepley   n[2] = x1[0]*x2[1] - x1[1]*x2[0];
7448b49ba18SBarry Smith   norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]);
745ccd2543fSMatthew G Knepley   n[0] /= norm;
746ccd2543fSMatthew G Knepley   n[1] /= norm;
747ccd2543fSMatthew G Knepley   n[2] /= norm;
748ccd2543fSMatthew G Knepley   /* 1) Take the normal vector and rotate until it is \hat z
749ccd2543fSMatthew G Knepley 
750ccd2543fSMatthew G Knepley     Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then
751ccd2543fSMatthew G Knepley 
752ccd2543fSMatthew G Knepley     R = /  alpha nx nz  alpha ny nz -1/alpha \
753ccd2543fSMatthew G Knepley         | -alpha ny     alpha nx        0    |
754ccd2543fSMatthew G Knepley         \     nx            ny         nz    /
755ccd2543fSMatthew G Knepley 
756ccd2543fSMatthew G Knepley     will rotate the normal vector to \hat z
757ccd2543fSMatthew G Knepley   */
7588b49ba18SBarry Smith   sqrtz = PetscSqrtReal(1.0 - n[2]*n[2]);
75973868372SMatthew G. Knepley   /* Check for n = z */
76073868372SMatthew G. Knepley   if (sqrtz < 1.0e-10) {
7617df32b8bSSanderA     const PetscInt s = PetscSign(n[2]);
7627df32b8bSSanderA     /* If nz < 0, rotate 180 degrees around x-axis */
76399dec3a6SMatthew G. Knepley     for (p = 3; p < coordSize/3; ++p) {
76499dec3a6SMatthew G. Knepley       coords[p*2+0] = PetscRealPart(coords[p*dim+0] - coords[0*dim+0]);
7657df32b8bSSanderA       coords[p*2+1] = (PetscRealPart(coords[p*dim+1] - coords[0*dim+1])) * s;
76673868372SMatthew G. Knepley     }
76799dec3a6SMatthew G. Knepley     coords[0] = 0.0;
76899dec3a6SMatthew G. Knepley     coords[1] = 0.0;
7697df32b8bSSanderA     coords[2] = x1[0];
7707df32b8bSSanderA     coords[3] = x1[1] * s;
7717df32b8bSSanderA     coords[4] = x2[0];
7727df32b8bSSanderA     coords[5] = x2[1] * s;
7737df32b8bSSanderA     R[0] = 1.0;     R[1] = 0.0;     R[2] = 0.0;
7747df32b8bSSanderA     R[3] = 0.0;     R[4] = 1.0 * s; R[5] = 0.0;
7757df32b8bSSanderA     R[6] = 0.0;     R[7] = 0.0;     R[8] = 1.0 * s;
77673868372SMatthew G. Knepley     PetscFunctionReturn(0);
77773868372SMatthew G. Knepley   }
778da18b5e6SMatthew G Knepley   alpha = 1.0/sqrtz;
779ccd2543fSMatthew G Knepley   R[0] =  alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz;
780ccd2543fSMatthew G Knepley   R[3] = -alpha*n[1];      R[4] = alpha*n[0];      R[5] = 0.0;
781ccd2543fSMatthew G Knepley   R[6] =  n[0];            R[7] = n[1];            R[8] = n[2];
782ccd2543fSMatthew G Knepley   for (d = 0; d < dim; ++d) {
783ccd2543fSMatthew G Knepley     x1p[d] = 0.0;
784ccd2543fSMatthew G Knepley     x2p[d] = 0.0;
785ccd2543fSMatthew G Knepley     for (e = 0; e < dim; ++e) {
786ccd2543fSMatthew G Knepley       x1p[d] += R[d*dim+e]*x1[e];
787ccd2543fSMatthew G Knepley       x2p[d] += R[d*dim+e]*x2[e];
788ccd2543fSMatthew G Knepley     }
789ccd2543fSMatthew G Knepley   }
79048120919SToby Isaac   if (PetscAbsReal(x1p[2]) > 10. * PETSC_SMALL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
79148120919SToby Isaac   if (PetscAbsReal(x2p[2]) > 10. * PETSC_SMALL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
792ccd2543fSMatthew G Knepley   /* 2) Project to (x, y) */
79399dec3a6SMatthew G. Knepley   for (p = 3; p < coordSize/3; ++p) {
79499dec3a6SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
79599dec3a6SMatthew G. Knepley       xnp[d] = 0.0;
79699dec3a6SMatthew G. Knepley       for (e = 0; e < dim; ++e) {
79799dec3a6SMatthew G. Knepley         xnp[d] += R[d*dim+e]*PetscRealPart(coords[p*dim+e] - coords[0*dim+e]);
79899dec3a6SMatthew G. Knepley       }
79999dec3a6SMatthew G. Knepley       if (d < dim-1) coords[p*2+d] = xnp[d];
80099dec3a6SMatthew G. Knepley     }
80199dec3a6SMatthew G. Knepley   }
802ccd2543fSMatthew G Knepley   coords[0] = 0.0;
803ccd2543fSMatthew G Knepley   coords[1] = 0.0;
804ccd2543fSMatthew G Knepley   coords[2] = x1p[0];
805ccd2543fSMatthew G Knepley   coords[3] = x1p[1];
806ccd2543fSMatthew G Knepley   coords[4] = x2p[0];
807ccd2543fSMatthew G Knepley   coords[5] = x2p[1];
8087f07f362SMatthew G. Knepley   /* Output R^T which rotates \hat z to the input normal */
8097f07f362SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
8107f07f362SMatthew G. Knepley     for (e = d+1; e < dim; ++e) {
8117f07f362SMatthew G. Knepley       PetscReal tmp;
8127f07f362SMatthew G. Knepley 
8137f07f362SMatthew G. Knepley       tmp        = R[d*dim+e];
8147f07f362SMatthew G. Knepley       R[d*dim+e] = R[e*dim+d];
8157f07f362SMatthew G. Knepley       R[e*dim+d] = tmp;
8167f07f362SMatthew G. Knepley     }
8177f07f362SMatthew G. Knepley   }
818ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
819ccd2543fSMatthew G Knepley }
820ccd2543fSMatthew G Knepley 
8216322fe33SJed Brown PETSC_UNUSED
822834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[])
823834e62ceSMatthew G. Knepley {
824834e62ceSMatthew G. Knepley   /* Signed volume is 1/2 the determinant
825834e62ceSMatthew G. Knepley 
826834e62ceSMatthew G. Knepley    |  1  1  1 |
827834e62ceSMatthew G. Knepley    | x0 x1 x2 |
828834e62ceSMatthew G. Knepley    | y0 y1 y2 |
829834e62ceSMatthew G. Knepley 
830834e62ceSMatthew G. Knepley      but if x0,y0 is the origin, we have
831834e62ceSMatthew G. Knepley 
832834e62ceSMatthew G. Knepley    | x1 x2 |
833834e62ceSMatthew G. Knepley    | y1 y2 |
834834e62ceSMatthew G. Knepley   */
835834e62ceSMatthew G. Knepley   const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1];
836834e62ceSMatthew G. Knepley   const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1];
837834e62ceSMatthew G. Knepley   PetscReal       M[4], detM;
838834e62ceSMatthew G. Knepley   M[0] = x1; M[1] = x2;
83986623015SMatthew G. Knepley   M[2] = y1; M[3] = y2;
840923591dfSMatthew G. Knepley   DMPlex_Det2D_Internal(&detM, M);
841834e62ceSMatthew G. Knepley   *vol = 0.5*detM;
8423bc0b13bSBarry Smith   (void)PetscLogFlops(5.0);
843834e62ceSMatthew G. Knepley }
844834e62ceSMatthew G. Knepley 
845834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[])
846834e62ceSMatthew G. Knepley {
847923591dfSMatthew G. Knepley   DMPlex_Det2D_Internal(vol, coords);
848834e62ceSMatthew G. Knepley   *vol *= 0.5;
849834e62ceSMatthew G. Knepley }
850834e62ceSMatthew G. Knepley 
8516322fe33SJed Brown PETSC_UNUSED
852834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[])
853834e62ceSMatthew G. Knepley {
854834e62ceSMatthew G. Knepley   /* Signed volume is 1/6th of the determinant
855834e62ceSMatthew G. Knepley 
856834e62ceSMatthew G. Knepley    |  1  1  1  1 |
857834e62ceSMatthew G. Knepley    | x0 x1 x2 x3 |
858834e62ceSMatthew G. Knepley    | y0 y1 y2 y3 |
859834e62ceSMatthew G. Knepley    | z0 z1 z2 z3 |
860834e62ceSMatthew G. Knepley 
861834e62ceSMatthew G. Knepley      but if x0,y0,z0 is the origin, we have
862834e62ceSMatthew G. Knepley 
863834e62ceSMatthew G. Knepley    | x1 x2 x3 |
864834e62ceSMatthew G. Knepley    | y1 y2 y3 |
865834e62ceSMatthew G. Knepley    | z1 z2 z3 |
866834e62ceSMatthew G. Knepley   */
867834e62ceSMatthew G. Knepley   const PetscReal x1 = coords[3] - coords[0], y1 = coords[4]  - coords[1], z1 = coords[5]  - coords[2];
868834e62ceSMatthew G. Knepley   const PetscReal x2 = coords[6] - coords[0], y2 = coords[7]  - coords[1], z2 = coords[8]  - coords[2];
869834e62ceSMatthew G. Knepley   const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2];
870834e62ceSMatthew G. Knepley   PetscReal       M[9], detM;
871834e62ceSMatthew G. Knepley   M[0] = x1; M[1] = x2; M[2] = x3;
872834e62ceSMatthew G. Knepley   M[3] = y1; M[4] = y2; M[5] = y3;
873834e62ceSMatthew G. Knepley   M[6] = z1; M[7] = z2; M[8] = z3;
874923591dfSMatthew G. Knepley   DMPlex_Det3D_Internal(&detM, M);
875b7ad821dSMatthew G. Knepley   *vol = -0.16666666666666666666666*detM;
8763bc0b13bSBarry Smith   (void)PetscLogFlops(10.0);
877834e62ceSMatthew G. Knepley }
878834e62ceSMatthew G. Knepley 
8790ec8681fSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[])
8800ec8681fSMatthew G. Knepley {
881923591dfSMatthew G. Knepley   DMPlex_Det3D_Internal(vol, coords);
882b7ad821dSMatthew G. Knepley   *vol *= -0.16666666666666666666666;
8830ec8681fSMatthew G. Knepley }
8840ec8681fSMatthew G. Knepley 
885cb92db44SToby Isaac static PetscErrorCode DMPlexComputePointGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
886cb92db44SToby Isaac {
887cb92db44SToby Isaac   PetscSection   coordSection;
888cb92db44SToby Isaac   Vec            coordinates;
889cb92db44SToby Isaac   const PetscScalar *coords;
890cb92db44SToby Isaac   PetscInt       dim, d, off;
891cb92db44SToby Isaac   PetscErrorCode ierr;
892cb92db44SToby Isaac 
893cb92db44SToby Isaac   PetscFunctionBegin;
894cb92db44SToby Isaac   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
895cb92db44SToby Isaac   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
896cb92db44SToby Isaac   ierr = PetscSectionGetDof(coordSection,e,&dim);CHKERRQ(ierr);
897cb92db44SToby Isaac   if (!dim) PetscFunctionReturn(0);
898cb92db44SToby Isaac   ierr = PetscSectionGetOffset(coordSection,e,&off);CHKERRQ(ierr);
899cb92db44SToby Isaac   ierr = VecGetArrayRead(coordinates,&coords);CHKERRQ(ierr);
900cb92db44SToby Isaac   if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[off + d]);}
901cb92db44SToby Isaac   ierr = VecRestoreArrayRead(coordinates,&coords);CHKERRQ(ierr);
902cb92db44SToby Isaac   *detJ = 1.;
903cb92db44SToby Isaac   if (J) {
904cb92db44SToby Isaac     for (d = 0; d < dim * dim; d++) J[d] = 0.;
905cb92db44SToby Isaac     for (d = 0; d < dim; d++) J[d * dim + d] = 1.;
906cb92db44SToby Isaac     if (invJ) {
907cb92db44SToby Isaac       for (d = 0; d < dim * dim; d++) invJ[d] = 0.;
908cb92db44SToby Isaac       for (d = 0; d < dim; d++) invJ[d * dim + d] = 1.;
909cb92db44SToby Isaac     }
910cb92db44SToby Isaac   }
911cb92db44SToby Isaac   PetscFunctionReturn(0);
912cb92db44SToby Isaac }
913cb92db44SToby Isaac 
91417fe8556SMatthew G. Knepley static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
91517fe8556SMatthew G. Knepley {
91617fe8556SMatthew G. Knepley   PetscSection   coordSection;
91717fe8556SMatthew G. Knepley   Vec            coordinates;
918a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
9198bf5c034SToby Isaac   PetscInt       numCoords, d, pStart, pEnd, numSelfCoords = 0;
92017fe8556SMatthew G. Knepley   PetscErrorCode ierr;
92117fe8556SMatthew G. Knepley 
92217fe8556SMatthew G. Knepley   PetscFunctionBegin;
92317fe8556SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
92469d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
9258bf5c034SToby Isaac   ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr);
9268bf5c034SToby Isaac   if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);}
92717fe8556SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
9288bf5c034SToby Isaac   numCoords = numSelfCoords ? numSelfCoords : numCoords;
929adac9986SMatthew G. Knepley   if (invJ && !J) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL");
9307f07f362SMatthew G. Knepley   *detJ = 0.0;
93128dbe442SToby Isaac   if (numCoords == 6) {
93228dbe442SToby Isaac     const PetscInt dim = 3;
93328dbe442SToby Isaac     PetscReal      R[9], J0;
93428dbe442SToby Isaac 
93528dbe442SToby Isaac     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
936741bfc07SMatthew G. Knepley     ierr = DMPlexComputeProjection3Dto1D(coords, R);CHKERRQ(ierr);
93728dbe442SToby Isaac     if (J)    {
93828dbe442SToby Isaac       J0   = 0.5*PetscRealPart(coords[1]);
93928dbe442SToby Isaac       J[0] = R[0]*J0; J[1] = R[1]; J[2] = R[2];
94028dbe442SToby Isaac       J[3] = R[3]*J0; J[4] = R[4]; J[5] = R[5];
94128dbe442SToby Isaac       J[6] = R[6]*J0; J[7] = R[7]; J[8] = R[8];
94228dbe442SToby Isaac       DMPlex_Det3D_Internal(detJ, J);
94328dbe442SToby Isaac       if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
944adac9986SMatthew G. Knepley     }
94528dbe442SToby Isaac   } else if (numCoords == 4) {
9467f07f362SMatthew G. Knepley     const PetscInt dim = 2;
9477f07f362SMatthew G. Knepley     PetscReal      R[4], J0;
9487f07f362SMatthew G. Knepley 
9497f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
950741bfc07SMatthew G. Knepley     ierr = DMPlexComputeProjection2Dto1D(coords, R);CHKERRQ(ierr);
95117fe8556SMatthew G. Knepley     if (J)    {
9527f07f362SMatthew G. Knepley       J0   = 0.5*PetscRealPart(coords[1]);
9537f07f362SMatthew G. Knepley       J[0] = R[0]*J0; J[1] = R[1];
9547f07f362SMatthew G. Knepley       J[2] = R[2]*J0; J[3] = R[3];
955923591dfSMatthew G. Knepley       DMPlex_Det2D_Internal(detJ, J);
956923591dfSMatthew G. Knepley       if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
957adac9986SMatthew G. Knepley     }
9587f07f362SMatthew G. Knepley   } else if (numCoords == 2) {
9597f07f362SMatthew G. Knepley     const PetscInt dim = 1;
9607f07f362SMatthew G. Knepley 
9617f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
9627f07f362SMatthew G. Knepley     if (J)    {
9637f07f362SMatthew G. Knepley       J[0]  = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0]));
96417fe8556SMatthew G. Knepley       *detJ = J[0];
9653bc0b13bSBarry Smith       ierr = PetscLogFlops(2.0);CHKERRQ(ierr);
9663bc0b13bSBarry Smith       if (invJ) {invJ[0] = 1.0/J[0]; ierr = PetscLogFlops(1.0);CHKERRQ(ierr);}
967adac9986SMatthew G. Knepley     }
968796f034aSJed Brown   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %D != 2", numCoords);
96917fe8556SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
97017fe8556SMatthew G. Knepley   PetscFunctionReturn(0);
97117fe8556SMatthew G. Knepley }
97217fe8556SMatthew G. Knepley 
973ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
974ccd2543fSMatthew G Knepley {
975ccd2543fSMatthew G Knepley   PetscSection   coordSection;
976ccd2543fSMatthew G Knepley   Vec            coordinates;
977a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
9787f07f362SMatthew G. Knepley   PetscInt       numCoords, d, f, g;
979ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
980ccd2543fSMatthew G Knepley 
981ccd2543fSMatthew G Knepley   PetscFunctionBegin;
982ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
98369d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
984ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
9857f07f362SMatthew G. Knepley   *detJ = 0.0;
986ccd2543fSMatthew G Knepley   if (numCoords == 9) {
9877f07f362SMatthew G. Knepley     const PetscInt dim = 3;
9887f07f362SMatthew 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};
9897f07f362SMatthew G. Knepley 
9907f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
991741bfc07SMatthew G. Knepley     ierr = DMPlexComputeProjection3Dto2D(numCoords, coords, R);CHKERRQ(ierr);
9927f07f362SMatthew G. Knepley     if (J)    {
993b7ad821dSMatthew G. Knepley       const PetscInt pdim = 2;
994b7ad821dSMatthew G. Knepley 
995b7ad821dSMatthew G. Knepley       for (d = 0; d < pdim; d++) {
996b7ad821dSMatthew G. Knepley         for (f = 0; f < pdim; f++) {
997b7ad821dSMatthew G. Knepley           J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
998ccd2543fSMatthew G Knepley         }
9997f07f362SMatthew G. Knepley       }
10003bc0b13bSBarry Smith       ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1001923591dfSMatthew G. Knepley       DMPlex_Det3D_Internal(detJ, J0);
10027f07f362SMatthew G. Knepley       for (d = 0; d < dim; d++) {
10037f07f362SMatthew G. Knepley         for (f = 0; f < dim; f++) {
10047f07f362SMatthew G. Knepley           J[d*dim+f] = 0.0;
10057f07f362SMatthew G. Knepley           for (g = 0; g < dim; g++) {
10067f07f362SMatthew G. Knepley             J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
10077f07f362SMatthew G. Knepley           }
10087f07f362SMatthew G. Knepley         }
10097f07f362SMatthew G. Knepley       }
10103bc0b13bSBarry Smith       ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
10117f07f362SMatthew G. Knepley     }
1012923591dfSMatthew G. Knepley     if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
10137f07f362SMatthew G. Knepley   } else if (numCoords == 6) {
10147f07f362SMatthew G. Knepley     const PetscInt dim = 2;
10157f07f362SMatthew G. Knepley 
10167f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1017ccd2543fSMatthew G Knepley     if (J)    {
1018ccd2543fSMatthew G Knepley       for (d = 0; d < dim; d++) {
1019ccd2543fSMatthew G Knepley         for (f = 0; f < dim; f++) {
1020ccd2543fSMatthew G Knepley           J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
1021ccd2543fSMatthew G Knepley         }
1022ccd2543fSMatthew G Knepley       }
10233bc0b13bSBarry Smith       ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1024923591dfSMatthew G. Knepley       DMPlex_Det2D_Internal(detJ, J);
1025ccd2543fSMatthew G Knepley     }
1026923591dfSMatthew G. Knepley     if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1027796f034aSJed Brown   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %D != 6 or 9", numCoords);
1028ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
1029ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1030ccd2543fSMatthew G Knepley }
1031ccd2543fSMatthew G Knepley 
1032dfccc68fSToby Isaac static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1033ccd2543fSMatthew G Knepley {
1034ccd2543fSMatthew G Knepley   PetscSection   coordSection;
1035ccd2543fSMatthew G Knepley   Vec            coordinates;
1036a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
10370d29256aSToby Isaac   PetscInt       numCoords, numSelfCoords = 0, d, f, g, pStart, pEnd;
1038ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
1039ccd2543fSMatthew G Knepley 
1040ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1041ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
104269d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
10430d29256aSToby Isaac   ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr);
10440d29256aSToby Isaac   if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);}
104599dec3a6SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
104671f58de1SToby Isaac   numCoords = numSelfCoords ? numSelfCoords : numCoords;
1047dfccc68fSToby Isaac   if (!Nq) {
10487f07f362SMatthew G. Knepley     *detJ = 0.0;
104999dec3a6SMatthew G. Knepley     if (numCoords == 12) {
105099dec3a6SMatthew G. Knepley       const PetscInt dim = 3;
105199dec3a6SMatthew 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};
105299dec3a6SMatthew G. Knepley 
1053dfccc68fSToby Isaac       if (v)   {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);}
1054741bfc07SMatthew G. Knepley       ierr = DMPlexComputeProjection3Dto2D(numCoords, coords, R);CHKERRQ(ierr);
105599dec3a6SMatthew G. Knepley       if (J)    {
105699dec3a6SMatthew G. Knepley         const PetscInt pdim = 2;
105799dec3a6SMatthew G. Knepley 
105899dec3a6SMatthew G. Knepley         for (d = 0; d < pdim; d++) {
105999dec3a6SMatthew G. Knepley           J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
106099dec3a6SMatthew G. Knepley           J0[d*dim+1] = 0.5*(PetscRealPart(coords[3*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
106199dec3a6SMatthew G. Knepley         }
10623bc0b13bSBarry Smith         ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1063923591dfSMatthew G. Knepley         DMPlex_Det3D_Internal(detJ, J0);
106499dec3a6SMatthew G. Knepley         for (d = 0; d < dim; d++) {
106599dec3a6SMatthew G. Knepley           for (f = 0; f < dim; f++) {
106699dec3a6SMatthew G. Knepley             J[d*dim+f] = 0.0;
106799dec3a6SMatthew G. Knepley             for (g = 0; g < dim; g++) {
106899dec3a6SMatthew G. Knepley               J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
106999dec3a6SMatthew G. Knepley             }
107099dec3a6SMatthew G. Knepley           }
107199dec3a6SMatthew G. Knepley         }
10723bc0b13bSBarry Smith         ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
107399dec3a6SMatthew G. Knepley       }
1074923591dfSMatthew G. Knepley       if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
107571f58de1SToby Isaac     } else if (numCoords == 8) {
107699dec3a6SMatthew G. Knepley       const PetscInt dim = 2;
107799dec3a6SMatthew G. Knepley 
1078dfccc68fSToby Isaac       if (v)   {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);}
1079ccd2543fSMatthew G Knepley       if (J)    {
1080ccd2543fSMatthew G Knepley         for (d = 0; d < dim; d++) {
108199dec3a6SMatthew G. Knepley           J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
108299dec3a6SMatthew G. Knepley           J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1083ccd2543fSMatthew G Knepley         }
10843bc0b13bSBarry Smith         ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1085923591dfSMatthew G. Knepley         DMPlex_Det2D_Internal(detJ, J);
1086ccd2543fSMatthew G Knepley       }
1087923591dfSMatthew G. Knepley       if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1088796f034aSJed Brown     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords);
1089dfccc68fSToby Isaac   } else {
1090dfccc68fSToby Isaac     const PetscInt Nv = 4;
1091dfccc68fSToby Isaac     const PetscInt dimR = 2;
1092dfccc68fSToby Isaac     const PetscInt zToPlex[4] = {0, 1, 3, 2};
1093dfccc68fSToby Isaac     PetscReal zOrder[12];
1094dfccc68fSToby Isaac     PetscReal zCoeff[12];
1095dfccc68fSToby Isaac     PetscInt  i, j, k, l, dim;
1096dfccc68fSToby Isaac 
1097dfccc68fSToby Isaac     if (numCoords == 12) {
1098dfccc68fSToby Isaac       dim = 3;
1099dfccc68fSToby Isaac     } else if (numCoords == 8) {
1100dfccc68fSToby Isaac       dim = 2;
1101dfccc68fSToby Isaac     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords);
1102dfccc68fSToby Isaac     for (i = 0; i < Nv; i++) {
1103dfccc68fSToby Isaac       PetscInt zi = zToPlex[i];
1104dfccc68fSToby Isaac 
1105dfccc68fSToby Isaac       for (j = 0; j < dim; j++) {
1106dfccc68fSToby Isaac         zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]);
1107dfccc68fSToby Isaac       }
1108dfccc68fSToby Isaac     }
1109dfccc68fSToby Isaac     for (j = 0; j < dim; j++) {
1110dfccc68fSToby Isaac       zCoeff[dim * 0 + j] = 0.25 * (  zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1111dfccc68fSToby Isaac       zCoeff[dim * 1 + j] = 0.25 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1112dfccc68fSToby Isaac       zCoeff[dim * 2 + j] = 0.25 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1113dfccc68fSToby Isaac       zCoeff[dim * 3 + j] = 0.25 * (  zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1114dfccc68fSToby Isaac     }
1115dfccc68fSToby Isaac     for (i = 0; i < Nq; i++) {
1116dfccc68fSToby Isaac       PetscReal xi = points[dimR * i], eta = points[dimR * i + 1];
1117dfccc68fSToby Isaac 
1118dfccc68fSToby Isaac       if (v) {
1119dfccc68fSToby Isaac         PetscReal extPoint[4];
1120dfccc68fSToby Isaac 
1121dfccc68fSToby Isaac         extPoint[0] = 1.;
1122dfccc68fSToby Isaac         extPoint[1] = xi;
1123dfccc68fSToby Isaac         extPoint[2] = eta;
1124dfccc68fSToby Isaac         extPoint[3] = xi * eta;
1125dfccc68fSToby Isaac         for (j = 0; j < dim; j++) {
1126dfccc68fSToby Isaac           PetscReal val = 0.;
1127dfccc68fSToby Isaac 
1128dfccc68fSToby Isaac           for (k = 0; k < Nv; k++) {
1129dfccc68fSToby Isaac             val += extPoint[k] * zCoeff[dim * k + j];
1130dfccc68fSToby Isaac           }
1131dfccc68fSToby Isaac           v[i * dim + j] = val;
1132dfccc68fSToby Isaac         }
1133dfccc68fSToby Isaac       }
1134dfccc68fSToby Isaac       if (J) {
1135dfccc68fSToby Isaac         PetscReal extJ[8];
1136dfccc68fSToby Isaac 
1137dfccc68fSToby Isaac         extJ[0] = 0.;
1138dfccc68fSToby Isaac         extJ[1] = 0.;
1139dfccc68fSToby Isaac         extJ[2] = 1.;
1140dfccc68fSToby Isaac         extJ[3] = 0.;
1141dfccc68fSToby Isaac         extJ[4] = 0.;
1142dfccc68fSToby Isaac         extJ[5] = 1.;
1143dfccc68fSToby Isaac         extJ[6] = eta;
1144dfccc68fSToby Isaac         extJ[7] = xi;
1145dfccc68fSToby Isaac         for (j = 0; j < dim; j++) {
1146dfccc68fSToby Isaac           for (k = 0; k < dimR; k++) {
1147dfccc68fSToby Isaac             PetscReal val = 0.;
1148dfccc68fSToby Isaac 
1149dfccc68fSToby Isaac             for (l = 0; l < Nv; l++) {
1150dfccc68fSToby Isaac               val += zCoeff[dim * l + j] * extJ[dimR * l + k];
1151dfccc68fSToby Isaac             }
1152dfccc68fSToby Isaac             J[i * dim * dim + dim * j + k] = val;
1153dfccc68fSToby Isaac           }
1154dfccc68fSToby Isaac         }
1155dfccc68fSToby Isaac         if (dim == 3) { /* put the cross product in the third component of the Jacobian */
1156dfccc68fSToby Isaac           PetscReal x, y, z;
1157dfccc68fSToby Isaac           PetscReal *iJ = &J[i * dim * dim];
1158dfccc68fSToby Isaac           PetscReal norm;
1159dfccc68fSToby Isaac 
1160dfccc68fSToby Isaac           x = iJ[1 * dim + 0] * iJ[2 * dim + 1] - iJ[1 * dim + 1] * iJ[2 * dim + 0];
1161dfccc68fSToby Isaac           y = iJ[0 * dim + 1] * iJ[2 * dim + 0] - iJ[0 * dim + 0] * iJ[2 * dim + 1];
1162dfccc68fSToby Isaac           z = iJ[0 * dim + 0] * iJ[1 * dim + 1] - iJ[0 * dim + 1] * iJ[1 * dim + 0];
1163dfccc68fSToby Isaac           norm = PetscSqrtReal(x * x + y * y + z * z);
1164dfccc68fSToby Isaac           iJ[2] = x / norm;
1165dfccc68fSToby Isaac           iJ[5] = y / norm;
1166dfccc68fSToby Isaac           iJ[8] = z / norm;
1167dfccc68fSToby Isaac           DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]);
1168dfccc68fSToby Isaac           if (invJ) {DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);}
1169dfccc68fSToby Isaac         } else {
1170dfccc68fSToby Isaac           DMPlex_Det2D_Internal(&detJ[i], &J[i * dim * dim]);
1171dfccc68fSToby Isaac           if (invJ) {DMPlex_Invert2D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);}
1172dfccc68fSToby Isaac         }
1173dfccc68fSToby Isaac       }
1174dfccc68fSToby Isaac     }
1175dfccc68fSToby Isaac   }
117699dec3a6SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
1177ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1178ccd2543fSMatthew G Knepley }
1179ccd2543fSMatthew G Knepley 
1180ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1181ccd2543fSMatthew G Knepley {
1182ccd2543fSMatthew G Knepley   PetscSection   coordSection;
1183ccd2543fSMatthew G Knepley   Vec            coordinates;
1184a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
1185ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
118699dec3a6SMatthew G. Knepley   PetscInt       d;
1187ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
1188ccd2543fSMatthew G Knepley 
1189ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1190ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
119169d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1192ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
11937f07f362SMatthew G. Knepley   *detJ = 0.0;
11947f07f362SMatthew G. Knepley   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1195ccd2543fSMatthew G Knepley   if (J)    {
1196ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
1197f0df753eSMatthew G. Knepley       /* I orient with outward face normals */
1198f0df753eSMatthew G. Knepley       J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d]));
1199f0df753eSMatthew G. Knepley       J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
1200f0df753eSMatthew G. Knepley       J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1201ccd2543fSMatthew G Knepley     }
12023bc0b13bSBarry Smith     ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
1203923591dfSMatthew G. Knepley     DMPlex_Det3D_Internal(detJ, J);
1204ccd2543fSMatthew G Knepley   }
1205923591dfSMatthew G. Knepley   if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1206ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1207ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1208ccd2543fSMatthew G Knepley }
1209ccd2543fSMatthew G Knepley 
1210dfccc68fSToby Isaac static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1211ccd2543fSMatthew G Knepley {
1212ccd2543fSMatthew G Knepley   PetscSection   coordSection;
1213ccd2543fSMatthew G Knepley   Vec            coordinates;
1214a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
1215ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
1216ccd2543fSMatthew G Knepley   PetscInt       d;
1217ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
1218ccd2543fSMatthew G Knepley 
1219ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1220ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
122169d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1222ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1223dfccc68fSToby Isaac   if (!Nq) {
12247f07f362SMatthew G. Knepley     *detJ = 0.0;
1225dfccc68fSToby Isaac     if (v)   {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);}
1226ccd2543fSMatthew G Knepley     if (J)    {
1227ccd2543fSMatthew G Knepley       for (d = 0; d < dim; d++) {
1228f0df753eSMatthew G. Knepley         J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1229f0df753eSMatthew G. Knepley         J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
1230f0df753eSMatthew G. Knepley         J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d]));
1231ccd2543fSMatthew G Knepley       }
12323bc0b13bSBarry Smith       ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
1233923591dfSMatthew G. Knepley       DMPlex_Det3D_Internal(detJ, J);
1234ccd2543fSMatthew G Knepley     }
1235923591dfSMatthew G. Knepley     if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1236dfccc68fSToby Isaac   } else {
1237dfccc68fSToby Isaac     const PetscInt Nv = 8;
1238dfccc68fSToby Isaac     const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6};
1239dfccc68fSToby Isaac     const PetscInt dim = 3;
1240dfccc68fSToby Isaac     const PetscInt dimR = 3;
1241dfccc68fSToby Isaac     PetscReal zOrder[24];
1242dfccc68fSToby Isaac     PetscReal zCoeff[24];
1243dfccc68fSToby Isaac     PetscInt  i, j, k, l;
1244dfccc68fSToby Isaac 
1245dfccc68fSToby Isaac     for (i = 0; i < Nv; i++) {
1246dfccc68fSToby Isaac       PetscInt zi = zToPlex[i];
1247dfccc68fSToby Isaac 
1248dfccc68fSToby Isaac       for (j = 0; j < dim; j++) {
1249dfccc68fSToby Isaac         zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]);
1250dfccc68fSToby Isaac       }
1251dfccc68fSToby Isaac     }
1252dfccc68fSToby Isaac     for (j = 0; j < dim; j++) {
1253dfccc68fSToby Isaac       zCoeff[dim * 0 + j] = 0.125 * (  zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j] + zOrder[dim * 4 + j] + zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1254dfccc68fSToby Isaac       zCoeff[dim * 1 + j] = 0.125 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j] - zOrder[dim * 4 + j] + zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1255dfccc68fSToby Isaac       zCoeff[dim * 2 + j] = 0.125 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j] - zOrder[dim * 4 + j] - zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1256dfccc68fSToby Isaac       zCoeff[dim * 3 + j] = 0.125 * (  zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j] + zOrder[dim * 4 + j] - zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1257dfccc68fSToby Isaac       zCoeff[dim * 4 + j] = 0.125 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] - zOrder[dim * 3 + j] + zOrder[dim * 4 + j] + zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1258dfccc68fSToby Isaac       zCoeff[dim * 5 + j] = 0.125 * (+ zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] - zOrder[dim * 3 + j] - zOrder[dim * 4 + j] + zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1259dfccc68fSToby Isaac       zCoeff[dim * 6 + j] = 0.125 * (+ zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] - zOrder[dim * 3 + j] - zOrder[dim * 4 + j] - zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1260dfccc68fSToby Isaac       zCoeff[dim * 7 + j] = 0.125 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] - zOrder[dim * 3 + j] + zOrder[dim * 4 + j] - zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1261dfccc68fSToby Isaac     }
1262dfccc68fSToby Isaac     for (i = 0; i < Nq; i++) {
1263dfccc68fSToby Isaac       PetscReal xi = points[dimR * i], eta = points[dimR * i + 1], theta = points[dimR * i + 2];
1264dfccc68fSToby Isaac 
1265dfccc68fSToby Isaac       if (v) {
126691d2b7ceSToby Isaac         PetscReal extPoint[8];
1267dfccc68fSToby Isaac 
1268dfccc68fSToby Isaac         extPoint[0] = 1.;
1269dfccc68fSToby Isaac         extPoint[1] = xi;
1270dfccc68fSToby Isaac         extPoint[2] = eta;
1271dfccc68fSToby Isaac         extPoint[3] = xi * eta;
1272dfccc68fSToby Isaac         extPoint[4] = theta;
1273dfccc68fSToby Isaac         extPoint[5] = theta * xi;
1274dfccc68fSToby Isaac         extPoint[6] = theta * eta;
1275dfccc68fSToby Isaac         extPoint[7] = theta * eta * xi;
1276dfccc68fSToby Isaac         for (j = 0; j < dim; j++) {
1277dfccc68fSToby Isaac           PetscReal val = 0.;
1278dfccc68fSToby Isaac 
1279dfccc68fSToby Isaac           for (k = 0; k < Nv; k++) {
1280dfccc68fSToby Isaac             val += extPoint[k] * zCoeff[dim * k + j];
1281dfccc68fSToby Isaac           }
1282dfccc68fSToby Isaac           v[i * dim + j] = val;
1283dfccc68fSToby Isaac         }
1284dfccc68fSToby Isaac       }
1285dfccc68fSToby Isaac       if (J) {
1286dfccc68fSToby Isaac         PetscReal extJ[24];
1287dfccc68fSToby Isaac 
1288dfccc68fSToby Isaac         extJ[0]  = 0.         ; extJ[1]  = 0.        ; extJ[2]  = 0.      ;
1289dfccc68fSToby Isaac         extJ[3]  = 1.         ; extJ[4]  = 0.        ; extJ[5]  = 0.      ;
1290dfccc68fSToby Isaac         extJ[6]  = 0.         ; extJ[7]  = 1.        ; extJ[8]  = 0.      ;
1291dfccc68fSToby Isaac         extJ[9]  = eta        ; extJ[10] = xi        ; extJ[11] = 0.      ;
1292dfccc68fSToby Isaac         extJ[12] = 0.         ; extJ[13] = 0.        ; extJ[14] = 1.      ;
1293dfccc68fSToby Isaac         extJ[15] = theta      ; extJ[16] = 0.        ; extJ[17] = xi      ;
1294dfccc68fSToby Isaac         extJ[18] = 0.         ; extJ[19] = theta     ; extJ[20] = eta     ;
1295dfccc68fSToby Isaac         extJ[21] = theta * eta; extJ[22] = theta * xi; extJ[23] = eta * xi;
1296dfccc68fSToby Isaac 
1297dfccc68fSToby Isaac         for (j = 0; j < dim; j++) {
1298dfccc68fSToby Isaac           for (k = 0; k < dimR; k++) {
1299dfccc68fSToby Isaac             PetscReal val = 0.;
1300dfccc68fSToby Isaac 
1301dfccc68fSToby Isaac             for (l = 0; l < Nv; l++) {
1302dfccc68fSToby Isaac               val += zCoeff[dim * l + j] * extJ[dimR * l + k];
1303dfccc68fSToby Isaac             }
1304dfccc68fSToby Isaac             J[i * dim * dim + dim * j + k] = val;
1305dfccc68fSToby Isaac           }
1306dfccc68fSToby Isaac         }
1307dfccc68fSToby Isaac         DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]);
1308dfccc68fSToby Isaac         if (invJ) {DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);}
1309dfccc68fSToby Isaac       }
1310dfccc68fSToby Isaac     }
1311dfccc68fSToby Isaac   }
1312ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1313ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1314ccd2543fSMatthew G Knepley }
1315ccd2543fSMatthew G Knepley 
1316dfccc68fSToby Isaac static PetscErrorCode DMPlexComputeCellGeometryFEM_Implicit(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
1317dfccc68fSToby Isaac {
1318dfccc68fSToby Isaac   PetscInt        depth, dim, coordDim, coneSize, i;
1319dfccc68fSToby Isaac   PetscInt        Nq = 0;
1320dfccc68fSToby Isaac   const PetscReal *points = NULL;
1321dfccc68fSToby Isaac   DMLabel         depthLabel;
13227318780aSToby Isaac   PetscReal       v0[3] = {-1.}, J0[9] = {-1.}, detJ0 = -1.;
1323dfccc68fSToby Isaac   PetscBool       isAffine = PETSC_TRUE;
1324dfccc68fSToby Isaac   PetscErrorCode  ierr;
1325dfccc68fSToby Isaac 
1326dfccc68fSToby Isaac   PetscFunctionBegin;
1327dfccc68fSToby Isaac   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1328dfccc68fSToby Isaac   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
1329dfccc68fSToby Isaac   ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr);
1330dfccc68fSToby Isaac   ierr = DMLabelGetValue(depthLabel, cell, &dim);CHKERRQ(ierr);
1331dfccc68fSToby Isaac   if (depth == 1 && dim == 1) {
1332dfccc68fSToby Isaac     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1333dfccc68fSToby Isaac   }
1334dfccc68fSToby Isaac   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
1335dfccc68fSToby Isaac   if (coordDim > 3) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported coordinate dimension %D > 3", coordDim);
1336dfccc68fSToby Isaac   if (quad) {ierr = PetscQuadratureGetData(quad, NULL, &Nq, &points, NULL);CHKERRQ(ierr);}
1337dfccc68fSToby Isaac   switch (dim) {
1338dfccc68fSToby Isaac   case 0:
1339dfccc68fSToby Isaac     ierr = DMPlexComputePointGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr);
1340dfccc68fSToby Isaac     isAffine = PETSC_FALSE;
1341dfccc68fSToby Isaac     break;
1342dfccc68fSToby Isaac   case 1:
13437318780aSToby Isaac     if (Nq) {
13447318780aSToby Isaac       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr);
13457318780aSToby Isaac     } else {
13467318780aSToby Isaac       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr);
13477318780aSToby Isaac     }
1348dfccc68fSToby Isaac     break;
1349dfccc68fSToby Isaac   case 2:
1350dfccc68fSToby Isaac     switch (coneSize) {
1351dfccc68fSToby Isaac     case 3:
13527318780aSToby Isaac       if (Nq) {
13537318780aSToby Isaac         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr);
13547318780aSToby Isaac       } else {
13557318780aSToby Isaac         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr);
13567318780aSToby Isaac       }
1357dfccc68fSToby Isaac       break;
1358dfccc68fSToby Isaac     case 4:
1359dfccc68fSToby Isaac       ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, Nq, points, v, J, invJ, detJ);CHKERRQ(ierr);
1360dfccc68fSToby Isaac       isAffine = PETSC_FALSE;
1361dfccc68fSToby Isaac       break;
1362dfccc68fSToby Isaac     default:
1363dfccc68fSToby Isaac       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell);
1364dfccc68fSToby Isaac     }
1365dfccc68fSToby Isaac     break;
1366dfccc68fSToby Isaac   case 3:
1367dfccc68fSToby Isaac     switch (coneSize) {
1368dfccc68fSToby Isaac     case 4:
13697318780aSToby Isaac       if (Nq) {
13707318780aSToby Isaac         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr);
13717318780aSToby Isaac       } else {
13727318780aSToby Isaac         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr);
13737318780aSToby Isaac       }
1374dfccc68fSToby Isaac       break;
1375dfccc68fSToby Isaac     case 6: /* Faces */
1376dfccc68fSToby Isaac     case 8: /* Vertices */
1377dfccc68fSToby Isaac       ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, Nq, points, v, J, invJ, detJ);CHKERRQ(ierr);
1378dfccc68fSToby Isaac       isAffine = PETSC_FALSE;
1379dfccc68fSToby Isaac       break;
1380dfccc68fSToby Isaac     default:
1381dfccc68fSToby Isaac       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell);
1382dfccc68fSToby Isaac     }
1383dfccc68fSToby Isaac     break;
1384dfccc68fSToby Isaac   default:
1385dfccc68fSToby Isaac     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
1386dfccc68fSToby Isaac   }
13877318780aSToby Isaac   if (isAffine && Nq) {
1388dfccc68fSToby Isaac     if (v) {
1389dfccc68fSToby Isaac       for (i = 0; i < Nq; i++) {
13907318780aSToby Isaac         CoordinatesRefToReal(coordDim,dim,v0, J0, &points[dim * i], &v[coordDim * i]);
1391dfccc68fSToby Isaac       }
1392dfccc68fSToby Isaac     }
13937318780aSToby Isaac     if (detJ) {
13947318780aSToby Isaac       for (i = 0; i < Nq; i++) {
13957318780aSToby Isaac         detJ[i] = detJ0;
1396dfccc68fSToby Isaac       }
13977318780aSToby Isaac     }
13987318780aSToby Isaac     if (J) {
13997318780aSToby Isaac       PetscInt k;
14007318780aSToby Isaac 
14017318780aSToby Isaac       for (i = 0, k = 0; i < Nq; i++) {
1402dfccc68fSToby Isaac         PetscInt j;
1403dfccc68fSToby Isaac 
14047318780aSToby Isaac         for (j = 0; j < coordDim * coordDim; j++, k++) {
14057318780aSToby Isaac           J[k] = J0[j];
14067318780aSToby Isaac         }
14077318780aSToby Isaac       }
14087318780aSToby Isaac     }
14097318780aSToby Isaac     if (invJ) {
14107318780aSToby Isaac       PetscInt k;
14117318780aSToby Isaac       switch (coordDim) {
14127318780aSToby Isaac       case 0:
14137318780aSToby Isaac         break;
14147318780aSToby Isaac       case 1:
14157318780aSToby Isaac         invJ[0] = 1./J0[0];
14167318780aSToby Isaac         break;
14177318780aSToby Isaac       case 2:
14187318780aSToby Isaac         DMPlex_Invert2D_Internal(invJ, J0, detJ0);
14197318780aSToby Isaac         break;
14207318780aSToby Isaac       case 3:
14217318780aSToby Isaac         DMPlex_Invert3D_Internal(invJ, J0, detJ0);
14227318780aSToby Isaac         break;
14237318780aSToby Isaac       }
14247318780aSToby Isaac       for (i = 1, k = coordDim * coordDim; i < Nq; i++) {
14257318780aSToby Isaac         PetscInt j;
14267318780aSToby Isaac 
14277318780aSToby Isaac         for (j = 0; j < coordDim * coordDim; j++, k++) {
14287318780aSToby Isaac           invJ[k] = invJ[j];
14297318780aSToby Isaac         }
1430dfccc68fSToby Isaac       }
1431dfccc68fSToby Isaac     }
1432dfccc68fSToby Isaac   }
1433dfccc68fSToby Isaac   PetscFunctionReturn(0);
1434dfccc68fSToby Isaac }
1435dfccc68fSToby Isaac 
1436ccd2543fSMatthew G Knepley /*@C
14378e0841e0SMatthew G. Knepley   DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell
1438ccd2543fSMatthew G Knepley 
1439ccd2543fSMatthew G Knepley   Collective on DM
1440ccd2543fSMatthew G Knepley 
1441ccd2543fSMatthew G Knepley   Input Arguments:
1442ccd2543fSMatthew G Knepley + dm   - the DM
1443ccd2543fSMatthew G Knepley - cell - the cell
1444ccd2543fSMatthew G Knepley 
1445ccd2543fSMatthew G Knepley   Output Arguments:
1446ccd2543fSMatthew G Knepley + v0   - the translation part of this affine transform
1447ccd2543fSMatthew G Knepley . J    - the Jacobian of the transform from the reference element
1448ccd2543fSMatthew G Knepley . invJ - the inverse of the Jacobian
1449ccd2543fSMatthew G Knepley - detJ - the Jacobian determinant
1450ccd2543fSMatthew G Knepley 
1451ccd2543fSMatthew G Knepley   Level: advanced
1452ccd2543fSMatthew G Knepley 
1453ccd2543fSMatthew G Knepley   Fortran Notes:
1454ccd2543fSMatthew G Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
1455ccd2543fSMatthew G Knepley   include petsc.h90 in your code.
1456ccd2543fSMatthew G Knepley 
14578e0841e0SMatthew G. Knepley .seealso: DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinateVec()
1458ccd2543fSMatthew G Knepley @*/
14598e0841e0SMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
1460ccd2543fSMatthew G Knepley {
1461ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
1462ccd2543fSMatthew G Knepley 
1463ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1464dfccc68fSToby Isaac   ierr = DMPlexComputeCellGeometryFEM_Implicit(dm,cell,NULL,v0,J,invJ,detJ);CHKERRQ(ierr);
14658e0841e0SMatthew G. Knepley   PetscFunctionReturn(0);
14668e0841e0SMatthew G. Knepley }
14678e0841e0SMatthew G. Knepley 
14688e0841e0SMatthew G. Knepley #undef __FUNCT__
1469dfccc68fSToby Isaac #define __FUNCT__ "DMPlexComputeCellGeometryFEM_FE"
1470dfccc68fSToby Isaac static PetscErrorCode DMPlexComputeCellGeometryFEM_FE(DM dm, PetscFE fe, PetscInt point, PetscQuadrature quad, PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
14718e0841e0SMatthew G. Knepley {
1472dfccc68fSToby Isaac   PetscQuadrature  feQuad;
14738e0841e0SMatthew G. Knepley   PetscSection     coordSection;
14748e0841e0SMatthew G. Knepley   Vec              coordinates;
14758e0841e0SMatthew G. Knepley   PetscScalar     *coords = NULL;
14768e0841e0SMatthew G. Knepley   const PetscReal *quadPoints;
1477f960e424SToby Isaac   PetscReal       *basisDer, *basis, detJt;
1478f960e424SToby Isaac   PetscInt         dim, cdim, pdim, qdim, Nq, numCoords, q;
14798e0841e0SMatthew G. Knepley   PetscErrorCode   ierr;
14808e0841e0SMatthew G. Knepley 
14818e0841e0SMatthew G. Knepley   PetscFunctionBegin;
14828e0841e0SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
14838e0841e0SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
14848e0841e0SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr);
14858e0841e0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
14868e0841e0SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
1487dfccc68fSToby Isaac   if (!quad) { /* use the first point of the first functional of the dual space */
1488dfccc68fSToby Isaac     PetscDualSpace dsp;
1489dfccc68fSToby Isaac 
1490dfccc68fSToby Isaac     ierr = PetscFEGetDualSpace(fe, &dsp);CHKERRQ(ierr);
1491dfccc68fSToby Isaac     ierr = PetscDualSpaceGetFunctional(dsp, 0, &quad);CHKERRQ(ierr);
14928e0841e0SMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, &qdim, &Nq, &quadPoints, NULL);CHKERRQ(ierr);
1493dfccc68fSToby Isaac     Nq = 1;
1494dfccc68fSToby Isaac   } else {
1495dfccc68fSToby Isaac     ierr = PetscQuadratureGetData(quad, &qdim, &Nq, &quadPoints, NULL);CHKERRQ(ierr);
1496dfccc68fSToby Isaac   }
149791d2b7ceSToby Isaac   ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr);
1498dfccc68fSToby Isaac   ierr = PetscFEGetQuadrature(fe, &feQuad);CHKERRQ(ierr);
1499dfccc68fSToby Isaac   if (feQuad == quad) {
1500f960e424SToby Isaac     ierr = PetscFEGetDefaultTabulation(fe, &basis, &basisDer, NULL);CHKERRQ(ierr);
15018e0841e0SMatthew G. Knepley     if (numCoords != pdim*cdim) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "There are %d coordinates for point %d != %d*%d", numCoords, point, pdim, cdim);
1502dfccc68fSToby Isaac   } else {
1503dfccc68fSToby Isaac     ierr = PetscFEGetTabulation(fe, Nq, quadPoints, &basis, &basisDer, NULL);CHKERRQ(ierr);
1504dfccc68fSToby Isaac   }
1505dfccc68fSToby Isaac   if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %d != quadrature dimension %d", dim, qdim);
1506dfccc68fSToby Isaac   if (v) {
1507dfccc68fSToby Isaac     ierr = PetscMemzero(v, Nq*cdim*sizeof(PetscReal));CHKERRQ(ierr);
1508f960e424SToby Isaac     for (q = 0; q < Nq; ++q) {
1509f960e424SToby Isaac       PetscInt i, k;
1510f960e424SToby Isaac 
1511f960e424SToby Isaac       for (k = 0; k < pdim; ++k)
1512f960e424SToby Isaac         for (i = 0; i < cdim; ++i)
1513dfccc68fSToby Isaac           v[q*cdim + i] += basis[q*pdim + k] * PetscRealPart(coords[k*cdim + i]);
1514f960e424SToby Isaac       ierr = PetscLogFlops(2.0*pdim*cdim);CHKERRQ(ierr);
1515f960e424SToby Isaac     }
1516f960e424SToby Isaac   }
15178e0841e0SMatthew G. Knepley   if (J) {
1518dfccc68fSToby Isaac     ierr = PetscMemzero(J, Nq*cdim*cdim*sizeof(PetscReal));CHKERRQ(ierr);
15198e0841e0SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
15208e0841e0SMatthew G. Knepley       PetscInt i, j, k, c, r;
15218e0841e0SMatthew G. Knepley 
15228e0841e0SMatthew G. Knepley       /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */
15238e0841e0SMatthew G. Knepley       for (k = 0; k < pdim; ++k)
15248e0841e0SMatthew G. Knepley         for (j = 0; j < dim; ++j)
15258e0841e0SMatthew G. Knepley           for (i = 0; i < cdim; ++i)
1526dfccc68fSToby Isaac             J[(q*cdim + i)*cdim + j] += basisDer[(q*pdim + k)*dim + j] * PetscRealPart(coords[k*cdim + i]);
15273bc0b13bSBarry Smith       ierr = PetscLogFlops(2.0*pdim*dim*cdim);CHKERRQ(ierr);
15288e0841e0SMatthew G. Knepley       if (cdim > dim) {
15298e0841e0SMatthew G. Knepley         for (c = dim; c < cdim; ++c)
15308e0841e0SMatthew G. Knepley           for (r = 0; r < cdim; ++r)
15318e0841e0SMatthew G. Knepley             J[r*cdim+c] = r == c ? 1.0 : 0.0;
15328e0841e0SMatthew G. Knepley       }
1533f960e424SToby Isaac       if (!detJ && !invJ) continue;
1534a63b72c6SToby Isaac       detJt = 0.;
15358e0841e0SMatthew G. Knepley       switch (cdim) {
15368e0841e0SMatthew G. Knepley       case 3:
1537037dc194SToby Isaac         DMPlex_Det3D_Internal(&detJt, &J[q*cdim*dim]);
1538037dc194SToby Isaac         if (invJ) {DMPlex_Invert3D_Internal(&invJ[q*cdim*dim], &J[q*cdim*dim], detJt);}
153917fe8556SMatthew G. Knepley         break;
154049dc4407SMatthew G. Knepley       case 2:
15419f328543SToby Isaac         DMPlex_Det2D_Internal(&detJt, &J[q*cdim*dim]);
1542037dc194SToby Isaac         if (invJ) {DMPlex_Invert2D_Internal(&invJ[q*cdim*dim], &J[q*cdim*dim], detJt);}
154349dc4407SMatthew G. Knepley         break;
15448e0841e0SMatthew G. Knepley       case 1:
1545037dc194SToby Isaac         detJt = J[q*cdim*dim];
1546037dc194SToby Isaac         if (invJ) invJ[q*cdim*dim] = 1.0/detJt;
154749dc4407SMatthew G. Knepley       }
1548f960e424SToby Isaac       if (detJ) detJ[q] = detJt;
154949dc4407SMatthew G. Knepley     }
155049dc4407SMatthew G. Knepley   }
1551037dc194SToby Isaac   else if (detJ || invJ) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Need J to compute invJ or detJ");
1552dfccc68fSToby Isaac   if (feQuad != quad) {
1553dfccc68fSToby Isaac     ierr = PetscFERestoreTabulation(fe, Nq, quadPoints, &basis, &basisDer, NULL);CHKERRQ(ierr);
1554dfccc68fSToby Isaac   }
15558e0841e0SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr);
15568e0841e0SMatthew G. Knepley   PetscFunctionReturn(0);
15578e0841e0SMatthew G. Knepley }
15588e0841e0SMatthew G. Knepley 
15598e0841e0SMatthew G. Knepley /*@C
15608e0841e0SMatthew G. Knepley   DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell
15618e0841e0SMatthew G. Knepley 
15628e0841e0SMatthew G. Knepley   Collective on DM
15638e0841e0SMatthew G. Knepley 
15648e0841e0SMatthew G. Knepley   Input Arguments:
15658e0841e0SMatthew G. Knepley + dm   - the DM
15668e0841e0SMatthew G. Knepley . cell - the cell
1567dfccc68fSToby Isaac - quad - the quadrature containing the points in the reference element where the geometry will be evaluated.  If quad == NULL, geometry will be
1568dfccc68fSToby Isaac          evaluated at the first vertex of the reference element
15698e0841e0SMatthew G. Knepley 
15708e0841e0SMatthew G. Knepley   Output Arguments:
1571dfccc68fSToby Isaac + v    - the image of the transformed quadrature points, otherwise the image of the first vertex in the closure of the reference element
15728e0841e0SMatthew G. Knepley . J    - the Jacobian of the transform from the reference element at each quadrature point
15738e0841e0SMatthew G. Knepley . invJ - the inverse of the Jacobian at each quadrature point
15748e0841e0SMatthew G. Knepley - detJ - the Jacobian determinant at each quadrature point
15758e0841e0SMatthew G. Knepley 
15768e0841e0SMatthew G. Knepley   Level: advanced
15778e0841e0SMatthew G. Knepley 
15788e0841e0SMatthew G. Knepley   Fortran Notes:
15798e0841e0SMatthew G. Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
15808e0841e0SMatthew G. Knepley   include petsc.h90 in your code.
15818e0841e0SMatthew G. Knepley 
15828e0841e0SMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetCoordinateVec()
15838e0841e0SMatthew G. Knepley @*/
1584dfccc68fSToby Isaac PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
15858e0841e0SMatthew G. Knepley {
1586dfccc68fSToby Isaac   PetscFE        fe = NULL;
15878e0841e0SMatthew G. Knepley   PetscErrorCode ierr;
15888e0841e0SMatthew G. Knepley 
15898e0841e0SMatthew G. Knepley   PetscFunctionBegin;
1590dfccc68fSToby Isaac   if (dm->coordinateDM) {
1591dfccc68fSToby Isaac     PetscClassId id;
1592dfccc68fSToby Isaac     PetscInt     numFields;
1593dfccc68fSToby Isaac     PetscDS      prob = dm->coordinateDM->prob;
1594dfccc68fSToby Isaac     PetscObject  disc;
1595dfccc68fSToby Isaac 
1596dfccc68fSToby Isaac     ierr = PetscDSGetNumFields(prob, &numFields);CHKERRQ(ierr);
1597dfccc68fSToby Isaac     if (numFields) {
1598dfccc68fSToby Isaac       ierr = PetscDSGetDiscretization(prob,0,&disc);CHKERRQ(ierr);
1599dfccc68fSToby Isaac       ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
1600dfccc68fSToby Isaac       if (id == PETSCFE_CLASSID) {
1601dfccc68fSToby Isaac         fe = (PetscFE) disc;
1602dfccc68fSToby Isaac       }
1603dfccc68fSToby Isaac     }
1604dfccc68fSToby Isaac   }
1605dfccc68fSToby Isaac   if (!fe) {ierr = DMPlexComputeCellGeometryFEM_Implicit(dm, cell, quad, v, J, invJ, detJ);CHKERRQ(ierr);}
1606dfccc68fSToby Isaac   else     {ierr = DMPlexComputeCellGeometryFEM_FE(dm, fe, cell, quad, v, J, invJ, detJ);CHKERRQ(ierr);}
1607ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1608ccd2543fSMatthew G Knepley }
1609834e62ceSMatthew G. Knepley 
1610011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1611cc08537eSMatthew G. Knepley {
1612cc08537eSMatthew G. Knepley   PetscSection   coordSection;
1613cc08537eSMatthew G. Knepley   Vec            coordinates;
1614a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
161506e2781eSMatthew G. Knepley   PetscScalar    tmp[2];
1616cc08537eSMatthew G. Knepley   PetscInt       coordSize;
1617cc08537eSMatthew G. Knepley   PetscErrorCode ierr;
1618cc08537eSMatthew G. Knepley 
1619cc08537eSMatthew G. Knepley   PetscFunctionBegin;
1620cc08537eSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
162169d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1622cc08537eSMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1623011ea5d8SMatthew G. Knepley   if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now");
16242e17dfb7SMatthew G. Knepley   ierr = DMLocalizeCoordinate_Internal(dm, dim, coords, &coords[dim], tmp);CHKERRQ(ierr);
1625cc08537eSMatthew G. Knepley   if (centroid) {
162606e2781eSMatthew G. Knepley     centroid[0] = 0.5*PetscRealPart(coords[0] + tmp[0]);
162706e2781eSMatthew G. Knepley     centroid[1] = 0.5*PetscRealPart(coords[1] + tmp[1]);
1628cc08537eSMatthew G. Knepley   }
1629cc08537eSMatthew G. Knepley   if (normal) {
1630a60a936bSMatthew G. Knepley     PetscReal norm;
1631a60a936bSMatthew G. Knepley 
163206e2781eSMatthew G. Knepley     normal[0]  = -PetscRealPart(coords[1] - tmp[1]);
163306e2781eSMatthew G. Knepley     normal[1]  =  PetscRealPart(coords[0] - tmp[0]);
1634a60a936bSMatthew G. Knepley     norm       = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1]);
1635a60a936bSMatthew G. Knepley     normal[0] /= norm;
1636a60a936bSMatthew G. Knepley     normal[1] /= norm;
1637cc08537eSMatthew G. Knepley   }
1638cc08537eSMatthew G. Knepley   if (vol) {
163906e2781eSMatthew G. Knepley     *vol = PetscSqrtReal(PetscSqr(PetscRealPart(coords[0] - tmp[0])) + PetscSqr(PetscRealPart(coords[1] - tmp[1])));
1640cc08537eSMatthew G. Knepley   }
1641cc08537eSMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1642cc08537eSMatthew G. Knepley   PetscFunctionReturn(0);
1643cc08537eSMatthew G. Knepley }
1644cc08537eSMatthew G. Knepley 
1645cc08537eSMatthew G. Knepley /* Centroid_i = (\sum_n A_n Cn_i ) / A */
1646011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1647cc08537eSMatthew G. Knepley {
1648cc08537eSMatthew G. Knepley   PetscSection   coordSection;
1649cc08537eSMatthew G. Knepley   Vec            coordinates;
1650cc08537eSMatthew G. Knepley   PetscScalar   *coords = NULL;
16510a1d6728SMatthew G. Knepley   PetscReal      vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9];
16520a1d6728SMatthew G. Knepley   PetscInt       tdim = 2, coordSize, numCorners, p, d, e;
1653cc08537eSMatthew G. Knepley   PetscErrorCode ierr;
1654cc08537eSMatthew G. Knepley 
1655cc08537eSMatthew G. Knepley   PetscFunctionBegin;
1656cc08537eSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
16570a1d6728SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr);
165869d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1659cc08537eSMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
16600bce18caSMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
1661ceee4971SMatthew G. Knepley   if (dim > 2 && centroid) {
1662ceee4971SMatthew G. Knepley     v0[0] = PetscRealPart(coords[0]);
1663ceee4971SMatthew G. Knepley     v0[1] = PetscRealPart(coords[1]);
1664ceee4971SMatthew G. Knepley     v0[2] = PetscRealPart(coords[2]);
1665ceee4971SMatthew G. Knepley   }
1666011ea5d8SMatthew G. Knepley   if (normal) {
1667011ea5d8SMatthew G. Knepley     if (dim > 2) {
16681ee9d5ecSMatthew G. Knepley       const PetscReal x0 = PetscRealPart(coords[dim+0] - coords[0]), x1 = PetscRealPart(coords[dim*2+0] - coords[0]);
16691ee9d5ecSMatthew G. Knepley       const PetscReal y0 = PetscRealPart(coords[dim+1] - coords[1]), y1 = PetscRealPart(coords[dim*2+1] - coords[1]);
16701ee9d5ecSMatthew G. Knepley       const PetscReal z0 = PetscRealPart(coords[dim+2] - coords[2]), z1 = PetscRealPart(coords[dim*2+2] - coords[2]);
16710a1d6728SMatthew G. Knepley       PetscReal       norm;
16720a1d6728SMatthew G. Knepley 
16730a1d6728SMatthew G. Knepley       normal[0] = y0*z1 - z0*y1;
16740a1d6728SMatthew G. Knepley       normal[1] = z0*x1 - x0*z1;
16750a1d6728SMatthew G. Knepley       normal[2] = x0*y1 - y0*x1;
16768b49ba18SBarry Smith       norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]);
16770a1d6728SMatthew G. Knepley       normal[0] /= norm;
16780a1d6728SMatthew G. Knepley       normal[1] /= norm;
16790a1d6728SMatthew G. Knepley       normal[2] /= norm;
1680011ea5d8SMatthew G. Knepley     } else {
1681011ea5d8SMatthew G. Knepley       for (d = 0; d < dim; ++d) normal[d] = 0.0;
1682011ea5d8SMatthew G. Knepley     }
1683011ea5d8SMatthew G. Knepley   }
1684741bfc07SMatthew G. Knepley   if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D(coordSize, coords, R);CHKERRQ(ierr);}
16850a1d6728SMatthew G. Knepley   for (p = 0; p < numCorners; ++p) {
16860a1d6728SMatthew G. Knepley     /* Need to do this copy to get types right */
16870a1d6728SMatthew G. Knepley     for (d = 0; d < tdim; ++d) {
16881ee9d5ecSMatthew G. Knepley       ctmp[d]      = PetscRealPart(coords[p*tdim+d]);
16891ee9d5ecSMatthew G. Knepley       ctmp[tdim+d] = PetscRealPart(coords[((p+1)%numCorners)*tdim+d]);
16900a1d6728SMatthew G. Knepley     }
16910a1d6728SMatthew G. Knepley     Volume_Triangle_Origin_Internal(&vtmp, ctmp);
16920a1d6728SMatthew G. Knepley     vsum += vtmp;
16930a1d6728SMatthew G. Knepley     for (d = 0; d < tdim; ++d) {
16940a1d6728SMatthew G. Knepley       csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp;
16950a1d6728SMatthew G. Knepley     }
16960a1d6728SMatthew G. Knepley   }
16970a1d6728SMatthew G. Knepley   for (d = 0; d < tdim; ++d) {
16980a1d6728SMatthew G. Knepley     csum[d] /= (tdim+1)*vsum;
16990a1d6728SMatthew G. Knepley   }
17000a1d6728SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1701ee6bbdb2SSatish Balay   if (vol) *vol = PetscAbsReal(vsum);
17020a1d6728SMatthew G. Knepley   if (centroid) {
17030a1d6728SMatthew G. Knepley     if (dim > 2) {
17040a1d6728SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
17050a1d6728SMatthew G. Knepley         centroid[d] = v0[d];
17060a1d6728SMatthew G. Knepley         for (e = 0; e < dim; ++e) {
17070a1d6728SMatthew G. Knepley           centroid[d] += R[d*dim+e]*csum[e];
17080a1d6728SMatthew G. Knepley         }
17090a1d6728SMatthew G. Knepley       }
17100a1d6728SMatthew G. Knepley     } else for (d = 0; d < dim; ++d) centroid[d] = csum[d];
17110a1d6728SMatthew G. Knepley   }
1712cc08537eSMatthew G. Knepley   PetscFunctionReturn(0);
1713cc08537eSMatthew G. Knepley }
1714cc08537eSMatthew G. Knepley 
17150ec8681fSMatthew G. Knepley /* Centroid_i = (\sum_n V_n Cn_i ) / V */
1716011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
17170ec8681fSMatthew G. Knepley {
17180ec8681fSMatthew G. Knepley   PetscSection    coordSection;
17190ec8681fSMatthew G. Knepley   Vec             coordinates;
17200ec8681fSMatthew G. Knepley   PetscScalar    *coords = NULL;
172186623015SMatthew G. Knepley   PetscReal       vsum = 0.0, vtmp, coordsTmp[3*3];
1722a7df9edeSMatthew G. Knepley   const PetscInt *faces, *facesO;
17230ec8681fSMatthew G. Knepley   PetscInt        numFaces, f, coordSize, numCorners, p, d;
17240ec8681fSMatthew G. Knepley   PetscErrorCode  ierr;
17250ec8681fSMatthew G. Knepley 
17260ec8681fSMatthew G. Knepley   PetscFunctionBegin;
1727f6dae198SJed Brown   if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No support for dim %D > 3",dim);
17280ec8681fSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
172969d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
17300ec8681fSMatthew G. Knepley 
1731d9a81ebdSMatthew G. Knepley   if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0;
17320ec8681fSMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr);
17330ec8681fSMatthew G. Knepley   ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr);
1734a7df9edeSMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr);
17350ec8681fSMatthew G. Knepley   for (f = 0; f < numFaces; ++f) {
1736011ea5d8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr);
17370ec8681fSMatthew G. Knepley     numCorners = coordSize/dim;
17380ec8681fSMatthew G. Knepley     switch (numCorners) {
17390ec8681fSMatthew G. Knepley     case 3:
17400ec8681fSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
17411ee9d5ecSMatthew G. Knepley         coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]);
17421ee9d5ecSMatthew G. Knepley         coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]);
17431ee9d5ecSMatthew G. Knepley         coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]);
17440ec8681fSMatthew G. Knepley       }
17450ec8681fSMatthew G. Knepley       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
1746a7df9edeSMatthew G. Knepley       if (facesO[f] < 0) vtmp = -vtmp;
17470ec8681fSMatthew G. Knepley       vsum += vtmp;
17484f25033aSJed Brown       if (centroid) {           /* Centroid of OABC = (a+b+c)/4 */
17490ec8681fSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
17501ee9d5ecSMatthew G. Knepley           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
17510ec8681fSMatthew G. Knepley         }
17520ec8681fSMatthew G. Knepley       }
17530ec8681fSMatthew G. Knepley       break;
17540ec8681fSMatthew G. Knepley     case 4:
17550ec8681fSMatthew G. Knepley       /* DO FOR PYRAMID */
17560ec8681fSMatthew G. Knepley       /* First tet */
17570ec8681fSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
17581ee9d5ecSMatthew G. Knepley         coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]);
17591ee9d5ecSMatthew G. Knepley         coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]);
17601ee9d5ecSMatthew G. Knepley         coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]);
17610ec8681fSMatthew G. Knepley       }
17620ec8681fSMatthew G. Knepley       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
1763a7df9edeSMatthew G. Knepley       if (facesO[f] < 0) vtmp = -vtmp;
17640ec8681fSMatthew G. Knepley       vsum += vtmp;
17650ec8681fSMatthew G. Knepley       if (centroid) {
17660ec8681fSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
17670ec8681fSMatthew G. Knepley           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
17680ec8681fSMatthew G. Knepley         }
17690ec8681fSMatthew G. Knepley       }
17700ec8681fSMatthew G. Knepley       /* Second tet */
17710ec8681fSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
17721ee9d5ecSMatthew G. Knepley         coordsTmp[0*dim+d] = PetscRealPart(coords[1*dim+d]);
17731ee9d5ecSMatthew G. Knepley         coordsTmp[1*dim+d] = PetscRealPart(coords[2*dim+d]);
17741ee9d5ecSMatthew G. Knepley         coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]);
17750ec8681fSMatthew G. Knepley       }
17760ec8681fSMatthew G. Knepley       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
1777a7df9edeSMatthew G. Knepley       if (facesO[f] < 0) vtmp = -vtmp;
17780ec8681fSMatthew G. Knepley       vsum += vtmp;
17790ec8681fSMatthew G. Knepley       if (centroid) {
17800ec8681fSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
17810ec8681fSMatthew G. Knepley           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
17820ec8681fSMatthew G. Knepley         }
17830ec8681fSMatthew G. Knepley       }
17840ec8681fSMatthew G. Knepley       break;
17850ec8681fSMatthew G. Knepley     default:
1786796f034aSJed Brown       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %D vertices", numCorners);
17870ec8681fSMatthew G. Knepley     }
17884f25033aSJed Brown     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr);
17890ec8681fSMatthew G. Knepley   }
17908763be8eSMatthew G. Knepley   if (vol)     *vol = PetscAbsReal(vsum);
17910ec8681fSMatthew G. Knepley   if (normal)   for (d = 0; d < dim; ++d) normal[d]    = 0.0;
1792d9a81ebdSMatthew G. Knepley   if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4);
17930ec8681fSMatthew G. Knepley   PetscFunctionReturn(0);
17940ec8681fSMatthew G. Knepley }
17950ec8681fSMatthew G. Knepley 
1796834e62ceSMatthew G. Knepley /*@C
1797834e62ceSMatthew G. Knepley   DMPlexComputeCellGeometryFVM - Compute the volume for a given cell
1798834e62ceSMatthew G. Knepley 
1799834e62ceSMatthew G. Knepley   Collective on DM
1800834e62ceSMatthew G. Knepley 
1801834e62ceSMatthew G. Knepley   Input Arguments:
1802834e62ceSMatthew G. Knepley + dm   - the DM
1803834e62ceSMatthew G. Knepley - cell - the cell
1804834e62ceSMatthew G. Knepley 
1805834e62ceSMatthew G. Knepley   Output Arguments:
1806834e62ceSMatthew G. Knepley + volume   - the cell volume
1807cc08537eSMatthew G. Knepley . centroid - the cell centroid
1808cc08537eSMatthew G. Knepley - normal - the cell normal, if appropriate
1809834e62ceSMatthew G. Knepley 
1810834e62ceSMatthew G. Knepley   Level: advanced
1811834e62ceSMatthew G. Knepley 
1812834e62ceSMatthew G. Knepley   Fortran Notes:
1813834e62ceSMatthew G. Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
1814834e62ceSMatthew G. Knepley   include petsc.h90 in your code.
1815834e62ceSMatthew G. Knepley 
181669d8a9ceSMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetCoordinateVec()
1817834e62ceSMatthew G. Knepley @*/
1818cc08537eSMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1819834e62ceSMatthew G. Knepley {
18200ec8681fSMatthew G. Knepley   PetscInt       depth, dim;
1821834e62ceSMatthew G. Knepley   PetscErrorCode ierr;
1822834e62ceSMatthew G. Knepley 
1823834e62ceSMatthew G. Knepley   PetscFunctionBegin;
1824834e62ceSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1825c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1826834e62ceSMatthew G. Knepley   if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated");
1827834e62ceSMatthew G. Knepley   /* We need to keep a pointer to the depth label */
1828c58f1c22SToby Isaac   ierr = DMGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr);
1829834e62ceSMatthew G. Knepley   /* Cone size is now the number of faces */
1830011ea5d8SMatthew G. Knepley   switch (depth) {
1831cc08537eSMatthew G. Knepley   case 1:
1832011ea5d8SMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
1833cc08537eSMatthew G. Knepley     break;
1834834e62ceSMatthew G. Knepley   case 2:
1835011ea5d8SMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
1836834e62ceSMatthew G. Knepley     break;
1837834e62ceSMatthew G. Knepley   case 3:
1838011ea5d8SMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
1839834e62ceSMatthew G. Knepley     break;
1840834e62ceSMatthew G. Knepley   default:
1841834e62ceSMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
1842834e62ceSMatthew G. Knepley   }
1843834e62ceSMatthew G. Knepley   PetscFunctionReturn(0);
1844834e62ceSMatthew G. Knepley }
1845113c68e6SMatthew G. Knepley 
1846*c501906fSMatthew G. Knepley /*@
1847*c501906fSMatthew G. Knepley   DMPlexComputeGeometryFEM - Precompute cell geometry for the entire mesh
1848*c501906fSMatthew G. Knepley 
1849*c501906fSMatthew G. Knepley   Collective on dm
1850*c501906fSMatthew G. Knepley 
1851*c501906fSMatthew G. Knepley   Input Parameter:
1852*c501906fSMatthew G. Knepley . dm - The DMPlex
1853*c501906fSMatthew G. Knepley 
1854*c501906fSMatthew G. Knepley   Output Parameter:
1855*c501906fSMatthew G. Knepley . cellgeom - A vector with the cell geometry data for each cell
1856*c501906fSMatthew G. Knepley 
1857*c501906fSMatthew G. Knepley   Level: beginner
1858*c501906fSMatthew G. Knepley 
1859*c501906fSMatthew G. Knepley .keywords: DMPlexComputeCellGeometryFEM()
1860*c501906fSMatthew G. Knepley @*/
1861c0d900a5SMatthew G. Knepley PetscErrorCode DMPlexComputeGeometryFEM(DM dm, Vec *cellgeom)
1862c0d900a5SMatthew G. Knepley {
1863c0d900a5SMatthew G. Knepley   DM             dmCell;
1864c0d900a5SMatthew G. Knepley   Vec            coordinates;
1865c0d900a5SMatthew G. Knepley   PetscSection   coordSection, sectionCell;
1866c0d900a5SMatthew G. Knepley   PetscScalar   *cgeom;
1867c0d900a5SMatthew G. Knepley   PetscInt       cStart, cEnd, cMax, c;
1868c0d900a5SMatthew G. Knepley   PetscErrorCode ierr;
1869c0d900a5SMatthew G. Knepley 
1870c0d900a5SMatthew G. Knepley   PetscFunctionBegin;
1871c0d900a5SMatthew G. Knepley   ierr = DMClone(dm, &dmCell);CHKERRQ(ierr);
1872c0d900a5SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1873c0d900a5SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1874c0d900a5SMatthew G. Knepley   ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr);
1875c0d900a5SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr);
1876c0d900a5SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionCell);CHKERRQ(ierr);
1877c0d900a5SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1878c0d900a5SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
1879c0d900a5SMatthew G. Knepley   cEnd = cMax < 0 ? cEnd : cMax;
1880c0d900a5SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr);
1881c0d900a5SMatthew G. Knepley   /* TODO This needs to be multiplied by Nq for non-affine */
18829e5edeeeSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFECellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);}
1883c0d900a5SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr);
1884c0d900a5SMatthew G. Knepley   ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr);
1885c0d900a5SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionCell);CHKERRQ(ierr);
1886c0d900a5SMatthew G. Knepley   ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr);
1887c0d900a5SMatthew G. Knepley   ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr);
1888c0d900a5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1889c0d900a5SMatthew G. Knepley     PetscFECellGeom *cg;
1890c0d900a5SMatthew G. Knepley 
1891c0d900a5SMatthew G. Knepley     ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
1892c0d900a5SMatthew G. Knepley     ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr);
1893c0d900a5SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dmCell, c, NULL, cg->v0, cg->J, cg->invJ, &cg->detJ);CHKERRQ(ierr);
1894c0d900a5SMatthew G. Knepley     if (cg->detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cg->detJ, c);
1895c0d900a5SMatthew G. Knepley   }
1896c0d900a5SMatthew G. Knepley   ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr);
1897c0d900a5SMatthew G. Knepley   ierr = DMDestroy(&dmCell);CHKERRQ(ierr);
1898c0d900a5SMatthew G. Knepley   PetscFunctionReturn(0);
1899c0d900a5SMatthew G. Knepley }
1900c0d900a5SMatthew G. Knepley 
1901891a9168SMatthew G. Knepley /*@
1902891a9168SMatthew G. Knepley   DMPlexComputeGeometryFVM - Computes the cell and face geometry for a finite volume method
1903891a9168SMatthew G. Knepley 
1904891a9168SMatthew G. Knepley   Input Parameter:
1905891a9168SMatthew G. Knepley . dm - The DM
1906891a9168SMatthew G. Knepley 
1907891a9168SMatthew G. Knepley   Output Parameters:
1908891a9168SMatthew G. Knepley + cellgeom - A Vec of PetscFVCellGeom data
1909891a9168SMatthew G. Knepley . facegeom - A Vec of PetscFVFaceGeom data
1910891a9168SMatthew G. Knepley 
1911891a9168SMatthew G. Knepley   Level: developer
1912891a9168SMatthew G. Knepley 
1913891a9168SMatthew G. Knepley .seealso: PetscFVFaceGeom, PetscFVCellGeom, DMPlexComputeGeometryFEM()
1914891a9168SMatthew G. Knepley @*/
1915113c68e6SMatthew G. Knepley PetscErrorCode DMPlexComputeGeometryFVM(DM dm, Vec *cellgeom, Vec *facegeom)
1916113c68e6SMatthew G. Knepley {
1917113c68e6SMatthew G. Knepley   DM             dmFace, dmCell;
1918113c68e6SMatthew G. Knepley   DMLabel        ghostLabel;
1919113c68e6SMatthew G. Knepley   PetscSection   sectionFace, sectionCell;
1920113c68e6SMatthew G. Knepley   PetscSection   coordSection;
1921113c68e6SMatthew G. Knepley   Vec            coordinates;
1922113c68e6SMatthew G. Knepley   PetscScalar   *fgeom, *cgeom;
1923113c68e6SMatthew G. Knepley   PetscReal      minradius, gminradius;
1924113c68e6SMatthew G. Knepley   PetscInt       dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f;
1925113c68e6SMatthew G. Knepley   PetscErrorCode ierr;
1926113c68e6SMatthew G. Knepley 
1927113c68e6SMatthew G. Knepley   PetscFunctionBegin;
1928113c68e6SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1929113c68e6SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1930113c68e6SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1931113c68e6SMatthew G. Knepley   /* Make cell centroids and volumes */
1932113c68e6SMatthew G. Knepley   ierr = DMClone(dm, &dmCell);CHKERRQ(ierr);
1933113c68e6SMatthew G. Knepley   ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr);
1934113c68e6SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr);
1935113c68e6SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionCell);CHKERRQ(ierr);
1936113c68e6SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1937113c68e6SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1938113c68e6SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr);
19399e5edeeeSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVCellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);}
1940113c68e6SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr);
1941113c68e6SMatthew G. Knepley   ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr);
1942113c68e6SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionCell);CHKERRQ(ierr);
1943113c68e6SMatthew G. Knepley   ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr);
194406348e87SToby Isaac   if (cEndInterior < 0) {
194506348e87SToby Isaac     cEndInterior = cEnd;
194606348e87SToby Isaac   }
1947113c68e6SMatthew G. Knepley   ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr);
1948113c68e6SMatthew G. Knepley   for (c = cStart; c < cEndInterior; ++c) {
1949113c68e6SMatthew G. Knepley     PetscFVCellGeom *cg;
1950113c68e6SMatthew G. Knepley 
1951113c68e6SMatthew G. Knepley     ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
1952113c68e6SMatthew G. Knepley     ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr);
1953113c68e6SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr);
1954113c68e6SMatthew G. Knepley   }
1955113c68e6SMatthew G. Knepley   /* Compute face normals and minimum cell radius */
1956113c68e6SMatthew G. Knepley   ierr = DMClone(dm, &dmFace);CHKERRQ(ierr);
1957113c68e6SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionFace);CHKERRQ(ierr);
1958113c68e6SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
1959113c68e6SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr);
19609e5edeeeSMatthew G. Knepley   for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVFaceGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);}
1961113c68e6SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr);
1962113c68e6SMatthew G. Knepley   ierr = DMSetDefaultSection(dmFace, sectionFace);CHKERRQ(ierr);
1963113c68e6SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionFace);CHKERRQ(ierr);
1964113c68e6SMatthew G. Knepley   ierr = DMCreateLocalVector(dmFace, facegeom);CHKERRQ(ierr);
1965113c68e6SMatthew G. Knepley   ierr = VecGetArray(*facegeom, &fgeom);CHKERRQ(ierr);
1966c58f1c22SToby Isaac   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
1967113c68e6SMatthew G. Knepley   minradius = PETSC_MAX_REAL;
1968113c68e6SMatthew G. Knepley   for (f = fStart; f < fEnd; ++f) {
1969113c68e6SMatthew G. Knepley     PetscFVFaceGeom *fg;
1970113c68e6SMatthew G. Knepley     PetscReal        area;
197150d63984SToby Isaac     PetscInt         ghost = -1, d, numChildren;
1972113c68e6SMatthew G. Knepley 
19739ac3fadcSMatthew G. Knepley     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);}
197450d63984SToby Isaac     ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr);
197550d63984SToby Isaac     if (ghost >= 0 || numChildren) continue;
1976113c68e6SMatthew G. Knepley     ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr);
1977113c68e6SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr);
1978113c68e6SMatthew G. Knepley     for (d = 0; d < dim; ++d) fg->normal[d] *= area;
1979113c68e6SMatthew G. Knepley     /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */
1980113c68e6SMatthew G. Knepley     {
1981113c68e6SMatthew G. Knepley       PetscFVCellGeom *cL, *cR;
198206348e87SToby Isaac       PetscInt         ncells;
1983113c68e6SMatthew G. Knepley       const PetscInt  *cells;
1984113c68e6SMatthew G. Knepley       PetscReal       *lcentroid, *rcentroid;
19850453c0cdSMatthew G. Knepley       PetscReal        l[3], r[3], v[3];
1986113c68e6SMatthew G. Knepley 
1987113c68e6SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr);
198806348e87SToby Isaac       ierr = DMPlexGetSupportSize(dm, f, &ncells);CHKERRQ(ierr);
1989113c68e6SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr);
1990113c68e6SMatthew G. Knepley       lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid;
199106348e87SToby Isaac       if (ncells > 1) {
199206348e87SToby Isaac         ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr);
1993113c68e6SMatthew G. Knepley         rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid;
199406348e87SToby Isaac       }
199506348e87SToby Isaac       else {
199606348e87SToby Isaac         rcentroid = fg->centroid;
199706348e87SToby Isaac       }
19982e17dfb7SMatthew G. Knepley       ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, lcentroid, l);CHKERRQ(ierr);
19992e17dfb7SMatthew G. Knepley       ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, rcentroid, r);CHKERRQ(ierr);
20000453c0cdSMatthew G. Knepley       DMPlex_WaxpyD_Internal(dim, -1, l, r, v);
2001113c68e6SMatthew G. Knepley       if (DMPlex_DotRealD_Internal(dim, fg->normal, v) < 0) {
2002113c68e6SMatthew G. Knepley         for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d];
2003113c68e6SMatthew G. Knepley       }
2004113c68e6SMatthew G. Knepley       if (DMPlex_DotRealD_Internal(dim, fg->normal, v) <= 0) {
2005113c68e6SMatthew G. Knepley         if (dim == 2) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed, normal (%g,%g) v (%g,%g)", f, (double) fg->normal[0], (double) fg->normal[1], (double) v[0], (double) v[1]);
2006113c68e6SMatthew G. Knepley         if (dim == 3) SETERRQ7(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed, normal (%g,%g,%g) v (%g,%g,%g)", f, (double) fg->normal[0], (double) fg->normal[1], (double) fg->normal[2], (double) v[0], (double) v[1], (double) v[2]);
2007113c68e6SMatthew G. Knepley         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f);
2008113c68e6SMatthew G. Knepley       }
2009113c68e6SMatthew G. Knepley       if (cells[0] < cEndInterior) {
2010113c68e6SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cL->centroid, v);
2011113c68e6SMatthew G. Knepley         minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v));
2012113c68e6SMatthew G. Knepley       }
201306348e87SToby Isaac       if (ncells > 1 && cells[1] < cEndInterior) {
2014113c68e6SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cR->centroid, v);
2015113c68e6SMatthew G. Knepley         minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v));
2016113c68e6SMatthew G. Knepley       }
2017113c68e6SMatthew G. Knepley     }
2018113c68e6SMatthew G. Knepley   }
2019b2566f29SBarry Smith   ierr = MPIU_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
2020113c68e6SMatthew G. Knepley   ierr = DMPlexSetMinRadius(dm, gminradius);CHKERRQ(ierr);
2021113c68e6SMatthew G. Knepley   /* Compute centroids of ghost cells */
2022113c68e6SMatthew G. Knepley   for (c = cEndInterior; c < cEnd; ++c) {
2023113c68e6SMatthew G. Knepley     PetscFVFaceGeom *fg;
2024113c68e6SMatthew G. Knepley     const PetscInt  *cone,    *support;
2025113c68e6SMatthew G. Knepley     PetscInt         coneSize, supportSize, s;
2026113c68e6SMatthew G. Knepley 
2027113c68e6SMatthew G. Knepley     ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr);
2028113c68e6SMatthew G. Knepley     if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize);
2029113c68e6SMatthew G. Knepley     ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr);
2030113c68e6SMatthew G. Knepley     ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr);
203150d63984SToby Isaac     if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 2", cone[0], supportSize);
2032113c68e6SMatthew G. Knepley     ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr);
2033113c68e6SMatthew G. Knepley     ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr);
2034113c68e6SMatthew G. Knepley     for (s = 0; s < 2; ++s) {
2035113c68e6SMatthew G. Knepley       /* Reflect ghost centroid across plane of face */
2036113c68e6SMatthew G. Knepley       if (support[s] == c) {
2037640bce14SSatish Balay         PetscFVCellGeom       *ci;
2038113c68e6SMatthew G. Knepley         PetscFVCellGeom       *cg;
2039113c68e6SMatthew G. Knepley         PetscReal              c2f[3], a;
2040113c68e6SMatthew G. Knepley 
2041113c68e6SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr);
2042113c68e6SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */
2043113c68e6SMatthew G. Knepley         a    = DMPlex_DotRealD_Internal(dim, c2f, fg->normal)/DMPlex_DotRealD_Internal(dim, fg->normal, fg->normal);
2044113c68e6SMatthew G. Knepley         ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr);
2045113c68e6SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, 2*a, fg->normal, ci->centroid, cg->centroid);
2046113c68e6SMatthew G. Knepley         cg->volume = ci->volume;
2047113c68e6SMatthew G. Knepley       }
2048113c68e6SMatthew G. Knepley     }
2049113c68e6SMatthew G. Knepley   }
2050113c68e6SMatthew G. Knepley   ierr = VecRestoreArray(*facegeom, &fgeom);CHKERRQ(ierr);
2051113c68e6SMatthew G. Knepley   ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr);
2052113c68e6SMatthew G. Knepley   ierr = DMDestroy(&dmCell);CHKERRQ(ierr);
2053113c68e6SMatthew G. Knepley   ierr = DMDestroy(&dmFace);CHKERRQ(ierr);
2054113c68e6SMatthew G. Knepley   PetscFunctionReturn(0);
2055113c68e6SMatthew G. Knepley }
2056113c68e6SMatthew G. Knepley 
2057113c68e6SMatthew G. Knepley /*@C
2058113c68e6SMatthew G. Knepley   DMPlexGetMinRadius - Returns the minimum distance from any cell centroid to a face
2059113c68e6SMatthew G. Knepley 
2060113c68e6SMatthew G. Knepley   Not collective
2061113c68e6SMatthew G. Knepley 
2062113c68e6SMatthew G. Knepley   Input Argument:
2063113c68e6SMatthew G. Knepley . dm - the DM
2064113c68e6SMatthew G. Knepley 
2065113c68e6SMatthew G. Knepley   Output Argument:
2066113c68e6SMatthew G. Knepley . minradius - the minium cell radius
2067113c68e6SMatthew G. Knepley 
2068113c68e6SMatthew G. Knepley   Level: developer
2069113c68e6SMatthew G. Knepley 
2070113c68e6SMatthew G. Knepley .seealso: DMGetCoordinates()
2071113c68e6SMatthew G. Knepley @*/
2072113c68e6SMatthew G. Knepley PetscErrorCode DMPlexGetMinRadius(DM dm, PetscReal *minradius)
2073113c68e6SMatthew G. Knepley {
2074113c68e6SMatthew G. Knepley   PetscFunctionBegin;
2075113c68e6SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2076113c68e6SMatthew G. Knepley   PetscValidPointer(minradius,2);
2077113c68e6SMatthew G. Knepley   *minradius = ((DM_Plex*) dm->data)->minradius;
2078113c68e6SMatthew G. Knepley   PetscFunctionReturn(0);
2079113c68e6SMatthew G. Knepley }
2080113c68e6SMatthew G. Knepley 
2081113c68e6SMatthew G. Knepley /*@C
2082113c68e6SMatthew G. Knepley   DMPlexSetMinRadius - Sets the minimum distance from the cell centroid to a face
2083113c68e6SMatthew G. Knepley 
2084113c68e6SMatthew G. Knepley   Logically collective
2085113c68e6SMatthew G. Knepley 
2086113c68e6SMatthew G. Knepley   Input Arguments:
2087113c68e6SMatthew G. Knepley + dm - the DM
2088113c68e6SMatthew G. Knepley - minradius - the minium cell radius
2089113c68e6SMatthew G. Knepley 
2090113c68e6SMatthew G. Knepley   Level: developer
2091113c68e6SMatthew G. Knepley 
2092113c68e6SMatthew G. Knepley .seealso: DMSetCoordinates()
2093113c68e6SMatthew G. Knepley @*/
2094113c68e6SMatthew G. Knepley PetscErrorCode DMPlexSetMinRadius(DM dm, PetscReal minradius)
2095113c68e6SMatthew G. Knepley {
2096113c68e6SMatthew G. Knepley   PetscFunctionBegin;
2097113c68e6SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2098113c68e6SMatthew G. Knepley   ((DM_Plex*) dm->data)->minradius = minradius;
2099113c68e6SMatthew G. Knepley   PetscFunctionReturn(0);
2100113c68e6SMatthew G. Knepley }
2101856ac710SMatthew G. Knepley 
2102856ac710SMatthew G. Knepley static PetscErrorCode BuildGradientReconstruction_Internal(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
2103856ac710SMatthew G. Knepley {
2104856ac710SMatthew G. Knepley   DMLabel        ghostLabel;
2105856ac710SMatthew G. Knepley   PetscScalar   *dx, *grad, **gref;
2106856ac710SMatthew G. Knepley   PetscInt       dim, cStart, cEnd, c, cEndInterior, maxNumFaces;
2107856ac710SMatthew G. Knepley   PetscErrorCode ierr;
2108856ac710SMatthew G. Knepley 
2109856ac710SMatthew G. Knepley   PetscFunctionBegin;
2110856ac710SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2111856ac710SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2112856ac710SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
2113856ac710SMatthew G. Knepley   ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr);
2114856ac710SMatthew G. Knepley   ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr);
2115c58f1c22SToby Isaac   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
2116856ac710SMatthew G. Knepley   ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr);
2117856ac710SMatthew G. Knepley   for (c = cStart; c < cEndInterior; c++) {
2118856ac710SMatthew G. Knepley     const PetscInt        *faces;
2119856ac710SMatthew G. Knepley     PetscInt               numFaces, usedFaces, f, d;
2120640bce14SSatish Balay     PetscFVCellGeom        *cg;
2121856ac710SMatthew G. Knepley     PetscBool              boundary;
2122856ac710SMatthew G. Knepley     PetscInt               ghost;
2123856ac710SMatthew G. Knepley 
2124856ac710SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
2125856ac710SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr);
2126856ac710SMatthew G. Knepley     ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr);
2127856ac710SMatthew G. Knepley     if (numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces);
2128856ac710SMatthew G. Knepley     for (f = 0, usedFaces = 0; f < numFaces; ++f) {
2129640bce14SSatish Balay       PetscFVCellGeom       *cg1;
2130856ac710SMatthew G. Knepley       PetscFVFaceGeom       *fg;
2131856ac710SMatthew G. Knepley       const PetscInt        *fcells;
2132856ac710SMatthew G. Knepley       PetscInt               ncell, side;
2133856ac710SMatthew G. Knepley 
2134856ac710SMatthew G. Knepley       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
2135a6ba4734SToby Isaac       ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
2136856ac710SMatthew G. Knepley       if ((ghost >= 0) || boundary) continue;
2137856ac710SMatthew G. Knepley       ierr  = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr);
2138856ac710SMatthew G. Knepley       side  = (c != fcells[0]); /* c is on left=0 or right=1 of face */
2139856ac710SMatthew G. Knepley       ncell = fcells[!side];    /* the neighbor */
2140856ac710SMatthew G. Knepley       ierr  = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr);
2141856ac710SMatthew G. Knepley       ierr  = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr);
2142856ac710SMatthew G. Knepley       for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d];
2143856ac710SMatthew G. Knepley       gref[usedFaces++] = fg->grad[side];  /* Gradient reconstruction term will go here */
2144856ac710SMatthew G. Knepley     }
2145856ac710SMatthew G. Knepley     if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?");
2146856ac710SMatthew G. Knepley     ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr);
2147856ac710SMatthew G. Knepley     for (f = 0, usedFaces = 0; f < numFaces; ++f) {
2148856ac710SMatthew G. Knepley       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
2149a6ba4734SToby Isaac       ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
2150856ac710SMatthew G. Knepley       if ((ghost >= 0) || boundary) continue;
2151856ac710SMatthew G. Knepley       for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d];
2152856ac710SMatthew G. Knepley       ++usedFaces;
2153856ac710SMatthew G. Knepley     }
2154856ac710SMatthew G. Knepley   }
2155856ac710SMatthew G. Knepley   ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr);
2156856ac710SMatthew G. Knepley   PetscFunctionReturn(0);
2157856ac710SMatthew G. Knepley }
2158856ac710SMatthew G. Knepley 
2159b81db932SToby Isaac static PetscErrorCode BuildGradientReconstruction_Internal_Tree(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
2160b81db932SToby Isaac {
2161b81db932SToby Isaac   DMLabel        ghostLabel;
2162b81db932SToby Isaac   PetscScalar   *dx, *grad, **gref;
2163b81db932SToby Isaac   PetscInt       dim, cStart, cEnd, c, cEndInterior, fStart, fEnd, f, nStart, nEnd, maxNumFaces = 0;
2164b81db932SToby Isaac   PetscSection   neighSec;
2165b81db932SToby Isaac   PetscInt     (*neighbors)[2];
2166b81db932SToby Isaac   PetscInt      *counter;
2167b81db932SToby Isaac   PetscErrorCode ierr;
2168b81db932SToby Isaac 
2169b81db932SToby Isaac   PetscFunctionBegin;
2170b81db932SToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2171b81db932SToby Isaac   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2172b81db932SToby Isaac   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
21735bc680faSToby Isaac   if (cEndInterior < 0) {
21745bc680faSToby Isaac     cEndInterior = cEnd;
21755bc680faSToby Isaac   }
2176b81db932SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&neighSec);CHKERRQ(ierr);
2177b81db932SToby Isaac   ierr = PetscSectionSetChart(neighSec,cStart,cEndInterior);CHKERRQ(ierr);
2178b81db932SToby Isaac   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
2179c58f1c22SToby Isaac   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
2180b81db932SToby Isaac   for (f = fStart; f < fEnd; f++) {
2181b81db932SToby Isaac     const PetscInt        *fcells;
2182b81db932SToby Isaac     PetscBool              boundary;
21835bc680faSToby Isaac     PetscInt               ghost = -1;
2184b81db932SToby Isaac     PetscInt               numChildren, numCells, c;
2185b81db932SToby Isaac 
218606348e87SToby Isaac     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);}
2187a6ba4734SToby Isaac     ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr);
2188b81db932SToby Isaac     ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr);
2189b81db932SToby Isaac     if ((ghost >= 0) || boundary || numChildren) continue;
2190b81db932SToby Isaac     ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr);
219106348e87SToby Isaac     if (numCells == 2) {
2192b81db932SToby Isaac       ierr = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr);
2193b81db932SToby Isaac       for (c = 0; c < 2; c++) {
2194b81db932SToby Isaac         PetscInt cell = fcells[c];
2195b81db932SToby Isaac 
2196e6885bbbSToby Isaac         if (cell >= cStart && cell < cEndInterior) {
2197b81db932SToby Isaac           ierr = PetscSectionAddDof(neighSec,cell,1);CHKERRQ(ierr);
2198b81db932SToby Isaac         }
2199b81db932SToby Isaac       }
2200b81db932SToby Isaac     }
220106348e87SToby Isaac   }
2202b81db932SToby Isaac   ierr = PetscSectionSetUp(neighSec);CHKERRQ(ierr);
2203b81db932SToby Isaac   ierr = PetscSectionGetMaxDof(neighSec,&maxNumFaces);CHKERRQ(ierr);
2204b81db932SToby Isaac   ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr);
2205b81db932SToby Isaac   nStart = 0;
2206b81db932SToby Isaac   ierr = PetscSectionGetStorageSize(neighSec,&nEnd);CHKERRQ(ierr);
2207b81db932SToby Isaac   ierr = PetscMalloc1((nEnd-nStart),&neighbors);CHKERRQ(ierr);
2208b81db932SToby Isaac   ierr = PetscCalloc1((cEndInterior-cStart),&counter);CHKERRQ(ierr);
2209b81db932SToby Isaac   for (f = fStart; f < fEnd; f++) {
2210b81db932SToby Isaac     const PetscInt        *fcells;
2211b81db932SToby Isaac     PetscBool              boundary;
22125bc680faSToby Isaac     PetscInt               ghost = -1;
2213b81db932SToby Isaac     PetscInt               numChildren, numCells, c;
2214b81db932SToby Isaac 
221506348e87SToby Isaac     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);}
2216a6ba4734SToby Isaac     ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr);
2217b81db932SToby Isaac     ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr);
2218b81db932SToby Isaac     if ((ghost >= 0) || boundary || numChildren) continue;
2219b81db932SToby Isaac     ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr);
222006348e87SToby Isaac     if (numCells == 2) {
2221b81db932SToby Isaac       ierr  = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr);
2222b81db932SToby Isaac       for (c = 0; c < 2; c++) {
2223b81db932SToby Isaac         PetscInt cell = fcells[c], off;
2224b81db932SToby Isaac 
2225e6885bbbSToby Isaac         if (cell >= cStart && cell < cEndInterior) {
2226b81db932SToby Isaac           ierr = PetscSectionGetOffset(neighSec,cell,&off);CHKERRQ(ierr);
2227b81db932SToby Isaac           off += counter[cell - cStart]++;
2228b81db932SToby Isaac           neighbors[off][0] = f;
2229b81db932SToby Isaac           neighbors[off][1] = fcells[1 - c];
2230b81db932SToby Isaac         }
2231b81db932SToby Isaac       }
2232b81db932SToby Isaac     }
223306348e87SToby Isaac   }
2234b81db932SToby Isaac   ierr = PetscFree(counter);CHKERRQ(ierr);
2235b81db932SToby Isaac   ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr);
2236b81db932SToby Isaac   for (c = cStart; c < cEndInterior; c++) {
2237317218b9SToby Isaac     PetscInt               numFaces, f, d, off, ghost = -1;
2238640bce14SSatish Balay     PetscFVCellGeom        *cg;
2239b81db932SToby Isaac 
2240b81db932SToby Isaac     ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
2241b81db932SToby Isaac     ierr = PetscSectionGetDof(neighSec, c, &numFaces);CHKERRQ(ierr);
2242b81db932SToby Isaac     ierr = PetscSectionGetOffset(neighSec, c, &off);CHKERRQ(ierr);
2243317218b9SToby Isaac     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, c, &ghost);CHKERRQ(ierr);}
2244317218b9SToby Isaac     if (ghost < 0 && numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces);
2245b81db932SToby Isaac     for (f = 0; f < numFaces; ++f) {
2246640bce14SSatish Balay       PetscFVCellGeom       *cg1;
2247b81db932SToby Isaac       PetscFVFaceGeom       *fg;
2248b81db932SToby Isaac       const PetscInt        *fcells;
2249b81db932SToby Isaac       PetscInt               ncell, side, nface;
2250b81db932SToby Isaac 
2251b81db932SToby Isaac       nface = neighbors[off + f][0];
2252b81db932SToby Isaac       ncell = neighbors[off + f][1];
2253b81db932SToby Isaac       ierr  = DMPlexGetSupport(dm,nface,&fcells);CHKERRQ(ierr);
2254b81db932SToby Isaac       side  = (c != fcells[0]);
2255b81db932SToby Isaac       ierr  = DMPlexPointLocalRef(dmFace, nface, fgeom, &fg);CHKERRQ(ierr);
2256b81db932SToby Isaac       ierr  = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr);
2257b81db932SToby Isaac       for (d = 0; d < dim; ++d) dx[f*dim+d] = cg1->centroid[d] - cg->centroid[d];
2258b81db932SToby Isaac       gref[f] = fg->grad[side];  /* Gradient reconstruction term will go here */
2259b81db932SToby Isaac     }
2260b81db932SToby Isaac     ierr = PetscFVComputeGradient(fvm, numFaces, dx, grad);CHKERRQ(ierr);
2261b81db932SToby Isaac     for (f = 0; f < numFaces; ++f) {
2262b81db932SToby Isaac       for (d = 0; d < dim; ++d) gref[f][d] = grad[f*dim+d];
2263b81db932SToby Isaac     }
2264b81db932SToby Isaac   }
2265b81db932SToby Isaac   ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr);
22665fe94518SToby Isaac   ierr = PetscSectionDestroy(&neighSec);CHKERRQ(ierr);
2267b81db932SToby Isaac   ierr = PetscFree(neighbors);CHKERRQ(ierr);
2268b81db932SToby Isaac   PetscFunctionReturn(0);
2269b81db932SToby Isaac }
2270b81db932SToby Isaac 
2271856ac710SMatthew G. Knepley /*@
2272856ac710SMatthew G. Knepley   DMPlexComputeGradientFVM - Compute geometric factors for gradient reconstruction, which are stored in the geometry data, and compute layout for gradient data
2273856ac710SMatthew G. Knepley 
2274856ac710SMatthew G. Knepley   Collective on DM
2275856ac710SMatthew G. Knepley 
2276856ac710SMatthew G. Knepley   Input Arguments:
2277856ac710SMatthew G. Knepley + dm  - The DM
2278856ac710SMatthew G. Knepley . fvm - The PetscFV
22798f9f38e3SMatthew G. Knepley . faceGeometry - The face geometry from DMPlexComputeFaceGeometryFVM()
22808f9f38e3SMatthew G. Knepley - cellGeometry - The face geometry from DMPlexComputeCellGeometryFVM()
2281856ac710SMatthew G. Knepley 
2282856ac710SMatthew G. Knepley   Output Parameters:
2283856ac710SMatthew G. Knepley + faceGeometry - The geometric factors for gradient calculation are inserted
2284856ac710SMatthew G. Knepley - dmGrad - The DM describing the layout of gradient data
2285856ac710SMatthew G. Knepley 
2286856ac710SMatthew G. Knepley   Level: developer
2287856ac710SMatthew G. Knepley 
2288856ac710SMatthew G. Knepley .seealso: DMPlexGetFaceGeometryFVM(), DMPlexGetCellGeometryFVM()
2289856ac710SMatthew G. Knepley @*/
2290856ac710SMatthew G. Knepley PetscErrorCode DMPlexComputeGradientFVM(DM dm, PetscFV fvm, Vec faceGeometry, Vec cellGeometry, DM *dmGrad)
2291856ac710SMatthew G. Knepley {
2292856ac710SMatthew G. Knepley   DM             dmFace, dmCell;
2293856ac710SMatthew G. Knepley   PetscScalar   *fgeom, *cgeom;
2294b81db932SToby Isaac   PetscSection   sectionGrad, parentSection;
2295856ac710SMatthew G. Knepley   PetscInt       dim, pdim, cStart, cEnd, cEndInterior, c;
2296856ac710SMatthew G. Knepley   PetscErrorCode ierr;
2297856ac710SMatthew G. Knepley 
2298856ac710SMatthew G. Knepley   PetscFunctionBegin;
2299856ac710SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2300856ac710SMatthew G. Knepley   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
2301856ac710SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2302856ac710SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
2303856ac710SMatthew G. Knepley   /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */
2304856ac710SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
2305856ac710SMatthew G. Knepley   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
2306856ac710SMatthew G. Knepley   ierr = VecGetArray(faceGeometry, &fgeom);CHKERRQ(ierr);
2307856ac710SMatthew G. Knepley   ierr = VecGetArray(cellGeometry, &cgeom);CHKERRQ(ierr);
2308b81db932SToby Isaac   ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
2309b81db932SToby Isaac   if (!parentSection) {
2310856ac710SMatthew G. Knepley     ierr = BuildGradientReconstruction_Internal(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr);
2311b5a3613cSMatthew G. Knepley   } else {
2312b81db932SToby Isaac     ierr = BuildGradientReconstruction_Internal_Tree(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr);
2313b81db932SToby Isaac   }
2314856ac710SMatthew G. Knepley   ierr = VecRestoreArray(faceGeometry, &fgeom);CHKERRQ(ierr);
2315856ac710SMatthew G. Knepley   ierr = VecRestoreArray(cellGeometry, &cgeom);CHKERRQ(ierr);
2316856ac710SMatthew G. Knepley   /* Create storage for gradients */
2317856ac710SMatthew G. Knepley   ierr = DMClone(dm, dmGrad);CHKERRQ(ierr);
2318856ac710SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionGrad);CHKERRQ(ierr);
2319856ac710SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr);
2320856ac710SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);}
2321856ac710SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr);
2322856ac710SMatthew G. Knepley   ierr = DMSetDefaultSection(*dmGrad, sectionGrad);CHKERRQ(ierr);
2323856ac710SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionGrad);CHKERRQ(ierr);
2324856ac710SMatthew G. Knepley   PetscFunctionReturn(0);
2325856ac710SMatthew G. Knepley }
2326b27d5b9eSToby Isaac 
2327*c501906fSMatthew G. Knepley /*@
2328*c501906fSMatthew G. Knepley   DMPlexGetDataFVM - Retrieve precomputed cell geometry
2329*c501906fSMatthew G. Knepley 
2330*c501906fSMatthew G. Knepley   Collective on DM
2331*c501906fSMatthew G. Knepley 
2332*c501906fSMatthew G. Knepley   Input Arguments:
2333*c501906fSMatthew G. Knepley + dm  - The DM
2334*c501906fSMatthew G. Knepley - fvm - The PetscFV
2335*c501906fSMatthew G. Knepley 
2336*c501906fSMatthew G. Knepley   Output Parameters:
2337*c501906fSMatthew G. Knepley + cellGeometry - The cell geometry
2338*c501906fSMatthew G. Knepley . faceGeometry - The face geometry
2339*c501906fSMatthew G. Knepley - dmGrad       - The gradient matrices
2340*c501906fSMatthew G. Knepley 
2341*c501906fSMatthew G. Knepley   Level: developer
2342*c501906fSMatthew G. Knepley 
2343*c501906fSMatthew G. Knepley .seealso: DMPlexComputeGeometryFVM()
2344*c501906fSMatthew G. Knepley @*/
2345b27d5b9eSToby Isaac PetscErrorCode DMPlexGetDataFVM(DM dm, PetscFV fv, Vec *cellgeom, Vec *facegeom, DM *gradDM)
2346b27d5b9eSToby Isaac {
2347b27d5b9eSToby Isaac   PetscObject    cellgeomobj, facegeomobj;
2348b27d5b9eSToby Isaac   PetscErrorCode ierr;
2349b27d5b9eSToby Isaac 
2350b27d5b9eSToby Isaac   PetscFunctionBegin;
2351b27d5b9eSToby Isaac   ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);CHKERRQ(ierr);
2352b27d5b9eSToby Isaac   if (!cellgeomobj) {
2353b27d5b9eSToby Isaac     Vec cellgeomInt, facegeomInt;
2354b27d5b9eSToby Isaac 
2355b27d5b9eSToby Isaac     ierr = DMPlexComputeGeometryFVM(dm, &cellgeomInt, &facegeomInt);CHKERRQ(ierr);
2356b27d5b9eSToby Isaac     ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_cellgeom_fvm",(PetscObject)cellgeomInt);CHKERRQ(ierr);
2357b27d5b9eSToby Isaac     ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_facegeom_fvm",(PetscObject)facegeomInt);CHKERRQ(ierr);
2358b27d5b9eSToby Isaac     ierr = VecDestroy(&cellgeomInt);CHKERRQ(ierr);
2359b27d5b9eSToby Isaac     ierr = VecDestroy(&facegeomInt);CHKERRQ(ierr);
2360b27d5b9eSToby Isaac     ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);CHKERRQ(ierr);
2361b27d5b9eSToby Isaac   }
2362b27d5b9eSToby Isaac   ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_facegeom_fvm", &facegeomobj);CHKERRQ(ierr);
2363b27d5b9eSToby Isaac   if (cellgeom) *cellgeom = (Vec) cellgeomobj;
2364b27d5b9eSToby Isaac   if (facegeom) *facegeom = (Vec) facegeomobj;
2365b27d5b9eSToby Isaac   if (gradDM) {
2366b27d5b9eSToby Isaac     PetscObject gradobj;
2367b27d5b9eSToby Isaac     PetscBool   computeGradients;
2368b27d5b9eSToby Isaac 
2369b27d5b9eSToby Isaac     ierr = PetscFVGetComputeGradients(fv,&computeGradients);CHKERRQ(ierr);
2370b27d5b9eSToby Isaac     if (!computeGradients) {
2371b27d5b9eSToby Isaac       *gradDM = NULL;
2372b27d5b9eSToby Isaac       PetscFunctionReturn(0);
2373b27d5b9eSToby Isaac     }
2374b27d5b9eSToby Isaac     ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);CHKERRQ(ierr);
2375b27d5b9eSToby Isaac     if (!gradobj) {
2376b27d5b9eSToby Isaac       DM dmGradInt;
2377b27d5b9eSToby Isaac 
2378b27d5b9eSToby Isaac       ierr = DMPlexComputeGradientFVM(dm,fv,(Vec) facegeomobj,(Vec) cellgeomobj,&dmGradInt);CHKERRQ(ierr);
2379b27d5b9eSToby Isaac       ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_dmgrad_fvm", (PetscObject)dmGradInt);CHKERRQ(ierr);
2380b27d5b9eSToby Isaac       ierr = DMDestroy(&dmGradInt);CHKERRQ(ierr);
2381b27d5b9eSToby Isaac       ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);CHKERRQ(ierr);
2382b27d5b9eSToby Isaac     }
2383b27d5b9eSToby Isaac     *gradDM = (DM) gradobj;
2384b27d5b9eSToby Isaac   }
2385b27d5b9eSToby Isaac   PetscFunctionReturn(0);
2386b27d5b9eSToby Isaac }
2387d6143a4eSToby Isaac 
23889d150b73SToby Isaac static PetscErrorCode DMPlexCoordinatesToReference_NewtonUpdate(PetscInt dimC, PetscInt dimR, PetscScalar *J, PetscScalar *invJ, PetscScalar *work,  PetscReal *resNeg, PetscReal *guess)
23899d150b73SToby Isaac {
23909d150b73SToby Isaac   PetscInt l, m;
23919d150b73SToby Isaac 
2392cd345991SToby Isaac   PetscFunctionBeginHot;
23939d150b73SToby Isaac   if (dimC == dimR && dimR <= 3) {
23949d150b73SToby Isaac     /* invert Jacobian, multiply */
23959d150b73SToby Isaac     PetscScalar det, idet;
23969d150b73SToby Isaac 
23979d150b73SToby Isaac     switch (dimR) {
23989d150b73SToby Isaac     case 1:
23999d150b73SToby Isaac       invJ[0] = 1./ J[0];
24009d150b73SToby Isaac       break;
24019d150b73SToby Isaac     case 2:
24029d150b73SToby Isaac       det = J[0] * J[3] - J[1] * J[2];
24039d150b73SToby Isaac       idet = 1./det;
24049d150b73SToby Isaac       invJ[0] =  J[3] * idet;
24059d150b73SToby Isaac       invJ[1] = -J[1] * idet;
24069d150b73SToby Isaac       invJ[2] = -J[2] * idet;
24079d150b73SToby Isaac       invJ[3] =  J[0] * idet;
24089d150b73SToby Isaac       break;
24099d150b73SToby Isaac     case 3:
24109d150b73SToby Isaac       {
24119d150b73SToby Isaac         invJ[0] = J[4] * J[8] - J[5] * J[7];
24129d150b73SToby Isaac         invJ[1] = J[2] * J[7] - J[1] * J[8];
24139d150b73SToby Isaac         invJ[2] = J[1] * J[5] - J[2] * J[4];
24149d150b73SToby Isaac         det = invJ[0] * J[0] + invJ[1] * J[3] + invJ[2] * J[6];
24159d150b73SToby Isaac         idet = 1./det;
24169d150b73SToby Isaac         invJ[0] *= idet;
24179d150b73SToby Isaac         invJ[1] *= idet;
24189d150b73SToby Isaac         invJ[2] *= idet;
24199d150b73SToby Isaac         invJ[3]  = idet * (J[5] * J[6] - J[3] * J[8]);
24209d150b73SToby Isaac         invJ[4]  = idet * (J[0] * J[8] - J[2] * J[6]);
24219d150b73SToby Isaac         invJ[5]  = idet * (J[2] * J[3] - J[0] * J[5]);
24229d150b73SToby Isaac         invJ[6]  = idet * (J[3] * J[7] - J[4] * J[6]);
24239d150b73SToby Isaac         invJ[7]  = idet * (J[1] * J[6] - J[0] * J[7]);
24249d150b73SToby Isaac         invJ[8]  = idet * (J[0] * J[4] - J[1] * J[3]);
24259d150b73SToby Isaac       }
24269d150b73SToby Isaac       break;
24279d150b73SToby Isaac     }
24289d150b73SToby Isaac     for (l = 0; l < dimR; l++) {
24299d150b73SToby Isaac       for (m = 0; m < dimC; m++) {
2430c6e120d1SToby Isaac         guess[l] += PetscRealPart(invJ[l * dimC + m]) * resNeg[m];
24319d150b73SToby Isaac       }
24329d150b73SToby Isaac     }
24339d150b73SToby Isaac   } else {
24349d150b73SToby Isaac #if defined(PETSC_USE_COMPLEX)
24359d150b73SToby Isaac     char transpose = 'C';
24369d150b73SToby Isaac #else
24379d150b73SToby Isaac     char transpose = 'T';
24389d150b73SToby Isaac #endif
24399d150b73SToby Isaac     PetscBLASInt m = dimR;
24409d150b73SToby Isaac     PetscBLASInt n = dimC;
24419d150b73SToby Isaac     PetscBLASInt one = 1;
24429d150b73SToby Isaac     PetscBLASInt worksize = dimR * dimC, info;
24439d150b73SToby Isaac 
24449d150b73SToby Isaac     for (l = 0; l < dimC; l++) {invJ[l] = resNeg[l];}
24459d150b73SToby Isaac 
24469d150b73SToby Isaac     PetscStackCallBLAS("LAPACKgels",LAPACKgels_(&transpose,&m,&n,&one,J,&m,invJ,&n,work,&worksize, &info));
24479d150b73SToby Isaac     if (info != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELS");
24489d150b73SToby Isaac 
2449c6e120d1SToby Isaac     for (l = 0; l < dimR; l++) {guess[l] += PetscRealPart(invJ[l]);}
24509d150b73SToby Isaac   }
24519d150b73SToby Isaac   PetscFunctionReturn(0);
24529d150b73SToby Isaac }
24539d150b73SToby Isaac 
24549d150b73SToby Isaac static PetscErrorCode DMPlexCoordinatesToReference_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[], Vec coords, PetscInt dimC, PetscInt dimR)
24559d150b73SToby Isaac {
2456c0cbe899SToby Isaac   PetscInt       coordSize, i, j, k, l, m, maxIts = 7, numV = (1 << dimR);
24579d150b73SToby Isaac   PetscScalar    *coordsScalar = NULL;
24589d150b73SToby Isaac   PetscReal      *cellData, *cellCoords, *cellCoeffs, *extJ, *resNeg;
24599d150b73SToby Isaac   PetscScalar    *J, *invJ, *work;
24609d150b73SToby Isaac   PetscErrorCode ierr;
24619d150b73SToby Isaac 
24629d150b73SToby Isaac   PetscFunctionBegin;
24639d150b73SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
24649d150b73SToby Isaac   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
24659d150b73SToby Isaac   if (coordSize < dimC * numV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expecting at least %D coordinates, got %D",dimC * (1 << dimR), coordSize);CHKERRQ(ierr);
24669d150b73SToby Isaac   ierr = DMGetWorkArray(dm, 2 * coordSize + dimR + dimC, PETSC_REAL, &cellData);CHKERRQ(ierr);
24679d150b73SToby Isaac   ierr = DMGetWorkArray(dm, 3 * dimR * dimC, PETSC_SCALAR, &J);CHKERRQ(ierr);
24689d150b73SToby Isaac   cellCoords = &cellData[0];
24699d150b73SToby Isaac   cellCoeffs = &cellData[coordSize];
24709d150b73SToby Isaac   extJ       = &cellData[2 * coordSize];
24719d150b73SToby Isaac   resNeg     = &cellData[2 * coordSize + dimR];
24729d150b73SToby Isaac   invJ       = &J[dimR * dimC];
24739d150b73SToby Isaac   work       = &J[2 * dimR * dimC];
24749d150b73SToby Isaac   if (dimR == 2) {
24759d150b73SToby Isaac     const PetscInt zToPlex[4] = {0, 1, 3, 2};
24769d150b73SToby Isaac 
24779d150b73SToby Isaac     for (i = 0; i < 4; i++) {
24789d150b73SToby Isaac       PetscInt plexI = zToPlex[i];
24799d150b73SToby Isaac 
24809d150b73SToby Isaac       for (j = 0; j < dimC; j++) {
24819d150b73SToby Isaac         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
24829d150b73SToby Isaac       }
24839d150b73SToby Isaac     }
24849d150b73SToby Isaac   } else if (dimR == 3) {
24859d150b73SToby Isaac     const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6};
24869d150b73SToby Isaac 
24879d150b73SToby Isaac     for (i = 0; i < 8; i++) {
24889d150b73SToby Isaac       PetscInt plexI = zToPlex[i];
24899d150b73SToby Isaac 
24909d150b73SToby Isaac       for (j = 0; j < dimC; j++) {
24919d150b73SToby Isaac         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
24929d150b73SToby Isaac       }
24939d150b73SToby Isaac     }
24949d150b73SToby Isaac   } else {
24959d150b73SToby Isaac     for (i = 0; i < coordSize; i++) {cellCoords[i] = PetscRealPart(coordsScalar[i]);}
24969d150b73SToby Isaac   }
24979d150b73SToby Isaac   /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */
24989d150b73SToby Isaac   for (i = 0; i < dimR; i++) {
24999d150b73SToby Isaac     PetscReal *swap;
25009d150b73SToby Isaac 
25019d150b73SToby Isaac     for (j = 0; j < (numV / 2); j++) {
25029d150b73SToby Isaac       for (k = 0; k < dimC; k++) {
25039d150b73SToby Isaac         cellCoeffs[dimC * j + k]                = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]);
25049d150b73SToby Isaac         cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]);
25059d150b73SToby Isaac       }
25069d150b73SToby Isaac     }
25079d150b73SToby Isaac 
25089d150b73SToby Isaac     if (i < dimR - 1) {
25099d150b73SToby Isaac       swap = cellCoeffs;
25109d150b73SToby Isaac       cellCoeffs = cellCoords;
25119d150b73SToby Isaac       cellCoords = swap;
25129d150b73SToby Isaac     }
25139d150b73SToby Isaac   }
25149d150b73SToby Isaac   ierr = PetscMemzero(refCoords,numPoints * dimR * sizeof (PetscReal));CHKERRQ(ierr);
25159d150b73SToby Isaac   for (j = 0; j < numPoints; j++) {
25169d150b73SToby Isaac     for (i = 0; i < maxIts; i++) {
25179d150b73SToby Isaac       PetscReal *guess = &refCoords[dimR * j];
25189d150b73SToby Isaac 
25199d150b73SToby Isaac       /* compute -residual and Jacobian */
25209d150b73SToby Isaac       for (k = 0; k < dimC; k++) {resNeg[k] = realCoords[dimC * j + k];}
25219d150b73SToby Isaac       for (k = 0; k < dimC * dimR; k++) {J[k] = 0.;}
25229d150b73SToby Isaac       for (k = 0; k < numV; k++) {
25239d150b73SToby Isaac         PetscReal extCoord = 1.;
25249d150b73SToby Isaac         for (l = 0; l < dimR; l++) {
25259d150b73SToby Isaac           PetscReal coord = guess[l];
25269d150b73SToby Isaac           PetscInt  dep   = (k & (1 << l)) >> l;
25279d150b73SToby Isaac 
25289d150b73SToby Isaac           extCoord *= dep * coord + !dep;
25299d150b73SToby Isaac           extJ[l] = dep;
25309d150b73SToby Isaac 
25319d150b73SToby Isaac           for (m = 0; m < dimR; m++) {
25329d150b73SToby Isaac             PetscReal coord = guess[m];
25339d150b73SToby Isaac             PetscInt  dep   = ((k & (1 << m)) >> m) && (m != l);
25349d150b73SToby Isaac             PetscReal mult  = dep * coord + !dep;
25359d150b73SToby Isaac 
25369d150b73SToby Isaac             extJ[l] *= mult;
25379d150b73SToby Isaac           }
25389d150b73SToby Isaac         }
25399d150b73SToby Isaac         for (l = 0; l < dimC; l++) {
25409d150b73SToby Isaac           PetscReal coeff = cellCoeffs[dimC * k + l];
25419d150b73SToby Isaac 
25429d150b73SToby Isaac           resNeg[l] -= coeff * extCoord;
25439d150b73SToby Isaac           for (m = 0; m < dimR; m++) {
25449d150b73SToby Isaac             J[dimR * l + m] += coeff * extJ[m];
25459d150b73SToby Isaac           }
25469d150b73SToby Isaac         }
25479d150b73SToby Isaac       }
25480611203eSToby Isaac #if 0 && defined(PETSC_USE_DEBUG)
25490611203eSToby Isaac       {
25500611203eSToby Isaac         PetscReal maxAbs = 0.;
25510611203eSToby Isaac 
25520611203eSToby Isaac         for (l = 0; l < dimC; l++) {
25530611203eSToby Isaac           maxAbs = PetscMax(maxAbs,PetscAbsReal(resNeg[l]));
25540611203eSToby Isaac         }
25550611203eSToby Isaac         ierr = PetscInfo4(dm,"cell %D, point %D, iter %D: res %g\n",cell,j,i,maxAbs);CHKERRQ(ierr);
25560611203eSToby Isaac       }
25570611203eSToby Isaac #endif
25589d150b73SToby Isaac 
25599d150b73SToby Isaac       ierr = DMPlexCoordinatesToReference_NewtonUpdate(dimC,dimR,J,invJ,work,resNeg,guess);CHKERRQ(ierr);
25609d150b73SToby Isaac     }
25619d150b73SToby Isaac   }
25629d150b73SToby Isaac   ierr = DMRestoreWorkArray(dm, 3 * dimR * dimC, PETSC_SCALAR, &J);CHKERRQ(ierr);
25639d150b73SToby Isaac   ierr = DMRestoreWorkArray(dm, 2 * coordSize + dimR + dimC, PETSC_REAL, &cellData);CHKERRQ(ierr);
25649d150b73SToby Isaac   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
25659d150b73SToby Isaac   PetscFunctionReturn(0);
25669d150b73SToby Isaac }
25679d150b73SToby Isaac 
25689d150b73SToby Isaac static PetscErrorCode DMPlexReferenceToCoordinates_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[], Vec coords, PetscInt dimC, PetscInt dimR)
25699d150b73SToby Isaac {
25709d150b73SToby Isaac   PetscInt       coordSize, i, j, k, l, numV = (1 << dimR);
25719d150b73SToby Isaac   PetscScalar    *coordsScalar = NULL;
25729d150b73SToby Isaac   PetscReal      *cellData, *cellCoords, *cellCoeffs;
25739d150b73SToby Isaac   PetscErrorCode ierr;
25749d150b73SToby Isaac 
25759d150b73SToby Isaac   PetscFunctionBegin;
25769d150b73SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
25779d150b73SToby Isaac   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
25789d150b73SToby Isaac   if (coordSize < dimC * numV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expecting at least %D coordinates, got %D",dimC * (1 << dimR), coordSize);CHKERRQ(ierr);
25799d150b73SToby Isaac   ierr = DMGetWorkArray(dm, 2 * coordSize, PETSC_REAL, &cellData);CHKERRQ(ierr);
25809d150b73SToby Isaac   cellCoords = &cellData[0];
25819d150b73SToby Isaac   cellCoeffs = &cellData[coordSize];
25829d150b73SToby Isaac   if (dimR == 2) {
25839d150b73SToby Isaac     const PetscInt zToPlex[4] = {0, 1, 3, 2};
25849d150b73SToby Isaac 
25859d150b73SToby Isaac     for (i = 0; i < 4; i++) {
25869d150b73SToby Isaac       PetscInt plexI = zToPlex[i];
25879d150b73SToby Isaac 
25889d150b73SToby Isaac       for (j = 0; j < dimC; j++) {
25899d150b73SToby Isaac         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
25909d150b73SToby Isaac       }
25919d150b73SToby Isaac     }
25929d150b73SToby Isaac   } else if (dimR == 3) {
25939d150b73SToby Isaac     const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6};
25949d150b73SToby Isaac 
25959d150b73SToby Isaac     for (i = 0; i < 8; i++) {
25969d150b73SToby Isaac       PetscInt plexI = zToPlex[i];
25979d150b73SToby Isaac 
25989d150b73SToby Isaac       for (j = 0; j < dimC; j++) {
25999d150b73SToby Isaac         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
26009d150b73SToby Isaac       }
26019d150b73SToby Isaac     }
26029d150b73SToby Isaac   } else {
26039d150b73SToby Isaac     for (i = 0; i < coordSize; i++) {cellCoords[i] = PetscRealPart(coordsScalar[i]);}
26049d150b73SToby Isaac   }
26059d150b73SToby Isaac   /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */
26069d150b73SToby Isaac   for (i = 0; i < dimR; i++) {
26079d150b73SToby Isaac     PetscReal *swap;
26089d150b73SToby Isaac 
26099d150b73SToby Isaac     for (j = 0; j < (numV / 2); j++) {
26109d150b73SToby Isaac       for (k = 0; k < dimC; k++) {
26119d150b73SToby Isaac         cellCoeffs[dimC * j + k]                = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]);
26129d150b73SToby Isaac         cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]);
26139d150b73SToby Isaac       }
26149d150b73SToby Isaac     }
26159d150b73SToby Isaac 
26169d150b73SToby Isaac     if (i < dimR - 1) {
26179d150b73SToby Isaac       swap = cellCoeffs;
26189d150b73SToby Isaac       cellCoeffs = cellCoords;
26199d150b73SToby Isaac       cellCoords = swap;
26209d150b73SToby Isaac     }
26219d150b73SToby Isaac   }
26229d150b73SToby Isaac   ierr = PetscMemzero(realCoords,numPoints * dimC * sizeof (PetscReal));CHKERRQ(ierr);
26239d150b73SToby Isaac   for (j = 0; j < numPoints; j++) {
26249d150b73SToby Isaac     const PetscReal *guess  = &refCoords[dimR * j];
26259d150b73SToby Isaac     PetscReal       *mapped = &realCoords[dimC * j];
26269d150b73SToby Isaac 
26279d150b73SToby Isaac     for (k = 0; k < numV; k++) {
26289d150b73SToby Isaac       PetscReal extCoord = 1.;
26299d150b73SToby Isaac       for (l = 0; l < dimR; l++) {
26309d150b73SToby Isaac         PetscReal coord = guess[l];
26319d150b73SToby Isaac         PetscInt  dep   = (k & (1 << l)) >> l;
26329d150b73SToby Isaac 
26339d150b73SToby Isaac         extCoord *= dep * coord + !dep;
26349d150b73SToby Isaac       }
26359d150b73SToby Isaac       for (l = 0; l < dimC; l++) {
26369d150b73SToby Isaac         PetscReal coeff = cellCoeffs[dimC * k + l];
26379d150b73SToby Isaac 
26389d150b73SToby Isaac         mapped[l] += coeff * extCoord;
26399d150b73SToby Isaac       }
26409d150b73SToby Isaac     }
26419d150b73SToby Isaac   }
26429d150b73SToby Isaac   ierr = DMRestoreWorkArray(dm, 2 * coordSize, PETSC_REAL, &cellData);CHKERRQ(ierr);
26439d150b73SToby Isaac   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
26449d150b73SToby Isaac   PetscFunctionReturn(0);
26459d150b73SToby Isaac }
26469d150b73SToby Isaac 
26479d150b73SToby Isaac static PetscErrorCode DMPlexCoordinatesToReference_FE(DM dm, PetscFE fe, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[], Vec coords, PetscInt dimC, PetscInt dimR)
26489d150b73SToby Isaac {
2649c0cbe899SToby Isaac   PetscInt       numComp, numDof, i, j, k, l, m, maxIter = 7, coordSize;
2650c6e120d1SToby Isaac   PetscScalar    *nodes = NULL;
2651c6e120d1SToby Isaac   PetscReal      *invV, *modes;
2652c6e120d1SToby Isaac   PetscReal      *B, *D, *resNeg;
2653c6e120d1SToby Isaac   PetscScalar    *J, *invJ, *work;
26549d150b73SToby Isaac   PetscErrorCode ierr;
26559d150b73SToby Isaac 
26569d150b73SToby Isaac   PetscFunctionBegin;
26579d150b73SToby Isaac   ierr = PetscFEGetDimension(fe, &numDof);CHKERRQ(ierr);
26589d150b73SToby Isaac   ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
26599d150b73SToby Isaac   if (numComp != dimC) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"coordinate discretization must have as many components (%D) as embedding dimension (!= %D)",numComp,dimC);
26609d150b73SToby Isaac   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
26619d150b73SToby Isaac   /* convert nodes to values in the stable evaluation basis */
2662c6e120d1SToby Isaac   ierr = DMGetWorkArray(dm,dimC * numDof,PETSC_REAL,&modes);CHKERRQ(ierr);
26639d150b73SToby Isaac   invV = fe->invV;
26649d150b73SToby Isaac   for (i = 0; i < numDof; i++) {
26659d150b73SToby Isaac     for (j = 0; j < dimC; j++) {
26669d150b73SToby Isaac       modes[i * dimC + j] = 0.;
26679d150b73SToby Isaac       for (k = 0; k < numDof; k++) {
2668c6e120d1SToby Isaac         modes[i * dimC + j] += invV[i * numDof + k] * PetscRealPart(nodes[k * dimC + j]);
26699d150b73SToby Isaac       }
26709d150b73SToby Isaac     }
26719d150b73SToby Isaac   }
2672c6e120d1SToby Isaac   ierr   = DMGetWorkArray(dm,numDof + numDof * dimR + dimC,PETSC_REAL,&B);CHKERRQ(ierr);
2673c6e120d1SToby Isaac   D      = &B[numDof];
2674c6e120d1SToby Isaac   resNeg = &D[numDof * dimR];
2675c6e120d1SToby Isaac   ierr = DMGetWorkArray(dm,3 * dimC * dimR,PETSC_SCALAR,&J);CHKERRQ(ierr);
26769d150b73SToby Isaac   invJ = &J[dimC * dimR];
26779d150b73SToby Isaac   work = &invJ[dimC * dimR];
26789d150b73SToby Isaac   for (i = 0; i < numPoints * dimR; i++) {refCoords[i] = 0.;}
26799d150b73SToby Isaac   for (j = 0; j < numPoints; j++) {
26809b1f03cbSToby Isaac       for (i = 0; i < maxIter; i++) { /* we could batch this so that we're not making big B and D arrays all the time */
26819d150b73SToby Isaac       PetscReal *guess = &refCoords[j * dimR];
26829d150b73SToby Isaac       ierr = PetscSpaceEvaluate(fe->basisSpace, 1, guess, B, D, NULL);CHKERRQ(ierr);
26839d150b73SToby Isaac       for (k = 0; k < dimC; k++) {resNeg[k] = realCoords[j * dimC + k];}
26849d150b73SToby Isaac       for (k = 0; k < dimC * dimR; k++) {J[k] = 0.;}
26859d150b73SToby Isaac       for (k = 0; k < numDof; k++) {
26869d150b73SToby Isaac         for (l = 0; l < dimC; l++) {
26879d150b73SToby Isaac           resNeg[l] -= modes[k * dimC + l] * B[k];
26889d150b73SToby Isaac           for (m = 0; m < dimR; m++) {
26899d150b73SToby Isaac             J[l * dimR + m] += modes[k * dimC + l] * D[k * dimR + m];
26909d150b73SToby Isaac           }
26919d150b73SToby Isaac         }
26929d150b73SToby Isaac       }
26930611203eSToby Isaac #if 0 && defined(PETSC_USE_DEBUG)
26940611203eSToby Isaac       {
26950611203eSToby Isaac         PetscReal maxAbs = 0.;
26960611203eSToby Isaac 
26970611203eSToby Isaac         for (l = 0; l < dimC; l++) {
26980611203eSToby Isaac           maxAbs = PetscMax(maxAbs,PetscAbsReal(resNeg[l]));
26990611203eSToby Isaac         }
27000611203eSToby Isaac         ierr = PetscInfo4(dm,"cell %D, point %D, iter %D: res %g\n",cell,j,i,maxAbs);CHKERRQ(ierr);
27010611203eSToby Isaac       }
27020611203eSToby Isaac #endif
27039d150b73SToby Isaac       ierr = DMPlexCoordinatesToReference_NewtonUpdate(dimC,dimR,J,invJ,work,resNeg,guess);CHKERRQ(ierr);
27049d150b73SToby Isaac     }
27059d150b73SToby Isaac   }
2706c6e120d1SToby Isaac   ierr = DMRestoreWorkArray(dm,3 * dimC * dimR,PETSC_SCALAR,&J);CHKERRQ(ierr);
2707c6e120d1SToby Isaac   ierr = DMRestoreWorkArray(dm,numDof + numDof * dimR + dimC,PETSC_REAL,&B);CHKERRQ(ierr);
2708c6e120d1SToby Isaac   ierr = DMRestoreWorkArray(dm,dimC * numDof,PETSC_REAL,&modes);CHKERRQ(ierr);
27099d150b73SToby Isaac   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
27109d150b73SToby Isaac   PetscFunctionReturn(0);
27119d150b73SToby Isaac }
27129d150b73SToby Isaac 
27139d150b73SToby Isaac static PetscErrorCode DMPlexReferenceToCoordinates_FE(DM dm, PetscFE fe, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[], Vec coords, PetscInt dimC, PetscInt dimR)
27149d150b73SToby Isaac {
27159d150b73SToby Isaac   PetscInt       numComp, numDof, i, j, k, l, coordSize;
2716c6e120d1SToby Isaac   PetscScalar    *nodes = NULL;
2717c6e120d1SToby Isaac   PetscReal      *invV, *modes;
27189d150b73SToby Isaac   PetscReal      *B;
27199d150b73SToby Isaac   PetscErrorCode ierr;
27209d150b73SToby Isaac 
27219d150b73SToby Isaac   PetscFunctionBegin;
27229d150b73SToby Isaac   ierr = PetscFEGetDimension(fe, &numDof);CHKERRQ(ierr);
27239d150b73SToby Isaac   ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
27249d150b73SToby Isaac   if (numComp != dimC) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"coordinate discretization must have as many components (%D) as embedding dimension (!= %D)",numComp,dimC);
27259d150b73SToby Isaac   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
27269d150b73SToby Isaac   /* convert nodes to values in the stable evaluation basis */
2727c6e120d1SToby Isaac   ierr = DMGetWorkArray(dm,dimC * numDof,PETSC_REAL,&modes);CHKERRQ(ierr);
27289d150b73SToby Isaac   invV = fe->invV;
27299d150b73SToby Isaac   for (i = 0; i < numDof; i++) {
27309d150b73SToby Isaac     for (j = 0; j < dimC; j++) {
27319d150b73SToby Isaac       modes[i * dimC + j] = 0.;
27329d150b73SToby Isaac       for (k = 0; k < numDof; k++) {
2733c6e120d1SToby Isaac         modes[i * dimC + j] += invV[i * numDof + k] * PetscRealPart(nodes[k * dimC + j]);
27349d150b73SToby Isaac       }
27359d150b73SToby Isaac     }
27369d150b73SToby Isaac   }
27379b1f03cbSToby Isaac   ierr = DMGetWorkArray(dm,numDof,PETSC_REAL,&B);CHKERRQ(ierr);
27389d150b73SToby Isaac   for (i = 0; i < numPoints * dimC; i++) {realCoords[i] = 0.;}
27399d150b73SToby Isaac   for (j = 0; j < numPoints; j++) {
27409d150b73SToby Isaac     const PetscReal *guess  = &refCoords[j * dimR];
27419d150b73SToby Isaac     PetscReal       *mapped = &realCoords[j * dimC];
27429d150b73SToby Isaac 
27439d150b73SToby Isaac     ierr = PetscSpaceEvaluate(fe->basisSpace, 1, guess, B, NULL, NULL);CHKERRQ(ierr);
27449d150b73SToby Isaac     for (k = 0; k < numDof; k++) {
27459d150b73SToby Isaac       for (l = 0; l < dimC; l++) {
27469d150b73SToby Isaac         mapped[l] += modes[k * dimC + l] * B[k];
27479d150b73SToby Isaac       }
27489d150b73SToby Isaac     }
27499d150b73SToby Isaac   }
27509b1f03cbSToby Isaac   ierr = DMRestoreWorkArray(dm,numDof,PETSC_REAL,&B);CHKERRQ(ierr);
2751c6e120d1SToby Isaac   ierr = DMRestoreWorkArray(dm,dimC * numDof,PETSC_REAL,&modes);CHKERRQ(ierr);
27529d150b73SToby Isaac   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
27539d150b73SToby Isaac   PetscFunctionReturn(0);
27549d150b73SToby Isaac }
27559d150b73SToby Isaac 
2756d6143a4eSToby Isaac /*@
2757d6143a4eSToby Isaac   DMPlexCoordinatesToReference - Pull coordinates back from the mesh to the reference element using a single element
2758d6143a4eSToby Isaac   map.  This inversion will be accurate inside the reference element, but may be inaccurate for mappings that do not
2759d6143a4eSToby Isaac   extend uniquely outside the reference cell (e.g, most non-affine maps)
2760d6143a4eSToby Isaac 
2761d6143a4eSToby Isaac   Not collective
2762d6143a4eSToby Isaac 
2763d6143a4eSToby Isaac   Input Parameters:
2764d6143a4eSToby Isaac + dm         - The mesh, with coordinate maps defined either by a PetscDS for the coordinate DM (see DMGetCoordinateDM()) or
2765d6143a4eSToby Isaac                implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or
2766d6143a4eSToby Isaac                as a multilinear map for tensor-product elements
2767d6143a4eSToby Isaac . cell       - the cell whose map is used.
2768d6143a4eSToby Isaac . numPoints  - the number of points to locate
27691b266c99SBarry Smith - realCoords - (numPoints x coordinate dimension) array of coordinates (see DMGetCoordinateDim())
2770d6143a4eSToby Isaac 
2771d6143a4eSToby Isaac   Output Parameters:
2772d6143a4eSToby Isaac . refCoords  - (numPoints x dimension) array of reference coordinates (see DMGetDimension())
27731b266c99SBarry Smith 
27741b266c99SBarry Smith   Level: intermediate
2775d6143a4eSToby Isaac @*/
2776d6143a4eSToby Isaac PetscErrorCode DMPlexCoordinatesToReference(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[])
2777d6143a4eSToby Isaac {
27789d150b73SToby Isaac   PetscInt       dimC, dimR, depth, cStart, cEnd, cEndInterior, i;
27799d150b73SToby Isaac   DM             coordDM = NULL;
27809d150b73SToby Isaac   Vec            coords;
27819d150b73SToby Isaac   PetscFE        fe = NULL;
27829d150b73SToby Isaac   PetscErrorCode ierr;
27839d150b73SToby Isaac 
2784d6143a4eSToby Isaac   PetscFunctionBegin;
27859d150b73SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
27869d150b73SToby Isaac   ierr = DMGetDimension(dm,&dimR);CHKERRQ(ierr);
27879d150b73SToby Isaac   ierr = DMGetCoordinateDim(dm,&dimC);CHKERRQ(ierr);
27889d150b73SToby Isaac   if (dimR <= 0 || dimC <= 0 || numPoints <= 0) PetscFunctionReturn(0);
27899d150b73SToby Isaac   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
27909d150b73SToby Isaac   ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
27919d150b73SToby Isaac   ierr = DMGetCoordinateDM(dm,&coordDM);CHKERRQ(ierr);
27929d150b73SToby Isaac   if (coordDM) {
27939d150b73SToby Isaac     PetscInt coordFields;
27949d150b73SToby Isaac 
27959d150b73SToby Isaac     ierr = DMGetNumFields(coordDM,&coordFields);CHKERRQ(ierr);
27969d150b73SToby Isaac     if (coordFields) {
27979d150b73SToby Isaac       PetscClassId id;
27989d150b73SToby Isaac       PetscObject  disc;
27999d150b73SToby Isaac 
28009d150b73SToby Isaac       ierr = DMGetField(coordDM,0,&disc);CHKERRQ(ierr);
28019d150b73SToby Isaac       ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
28029d150b73SToby Isaac       if (id == PETSCFE_CLASSID) {
28039d150b73SToby Isaac         fe = (PetscFE) disc;
28049d150b73SToby Isaac       }
28059d150b73SToby Isaac     }
28069d150b73SToby Isaac   }
28079d150b73SToby Isaac   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
28089d150b73SToby Isaac   ierr = DMPlexGetHybridBounds(dm,&cEndInterior,NULL,NULL,NULL);CHKERRQ(ierr);
28099d150b73SToby Isaac   cEnd = cEndInterior > 0 ? cEndInterior : cEnd;
28109d150b73SToby Isaac   if (cell < cStart || cell >= cEnd) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D not in cell range [%D,%D)",cell,cStart,cEnd);CHKERRQ(ierr);
28119d150b73SToby Isaac   if (!fe) { /* implicit discretization: affine or multilinear */
28129d150b73SToby Isaac     PetscInt  coneSize;
28139d150b73SToby Isaac     PetscBool isSimplex, isTensor;
28149d150b73SToby Isaac 
28159d150b73SToby Isaac     ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr);
28169d150b73SToby Isaac     isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE;
28179d150b73SToby Isaac     isTensor  = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE;
28189d150b73SToby Isaac     if (isSimplex) {
28199d150b73SToby Isaac       PetscReal detJ, *v0, *J, *invJ;
28209d150b73SToby Isaac 
28219d150b73SToby Isaac       ierr = DMGetWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);CHKERRQ(ierr);
28229d150b73SToby Isaac       J    = &v0[dimC];
28239d150b73SToby Isaac       invJ = &J[dimC * dimC];
28249d150b73SToby Isaac       ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, &detJ);CHKERRQ(ierr);
28259d150b73SToby Isaac       for (i = 0; i < numPoints; i++) { /* Apply the inverse affine transformation for each point */
28269d150b73SToby Isaac         CoordinatesRealToRef(dimC, dimR, v0, invJ, &realCoords[dimC * i], &refCoords[dimR * i]);
28279d150b73SToby Isaac       }
28289d150b73SToby Isaac       ierr = DMRestoreWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);CHKERRQ(ierr);
28299d150b73SToby Isaac     } else if (isTensor) {
28309d150b73SToby Isaac       ierr = DMPlexCoordinatesToReference_Tensor(coordDM, cell, numPoints, realCoords, refCoords, coords, dimC, dimR);CHKERRQ(ierr);
28319d150b73SToby Isaac     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unrecognized cone size %D",coneSize);
28329d150b73SToby Isaac   } else {
28339d150b73SToby Isaac     ierr = DMPlexCoordinatesToReference_FE(coordDM, fe, cell, numPoints, realCoords, refCoords, coords, dimC, dimR);CHKERRQ(ierr);
28349d150b73SToby Isaac   }
28359d150b73SToby Isaac   PetscFunctionReturn(0);
28369d150b73SToby Isaac }
28379d150b73SToby Isaac 
28389d150b73SToby Isaac /*@
28399d150b73SToby Isaac   DMPlexReferenceToCoordinates - Map references coordinates to coordinates in the the mesh for a single element map.
28409d150b73SToby Isaac 
28419d150b73SToby Isaac   Not collective
28429d150b73SToby Isaac 
28439d150b73SToby Isaac   Input Parameters:
28449d150b73SToby Isaac + dm         - The mesh, with coordinate maps defined either by a PetscDS for the coordinate DM (see DMGetCoordinateDM()) or
28459d150b73SToby Isaac                implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or
28469d150b73SToby Isaac                as a multilinear map for tensor-product elements
28479d150b73SToby Isaac . cell       - the cell whose map is used.
28489d150b73SToby Isaac . numPoints  - the number of points to locate
28499d150b73SToby Isaac + refCoords  - (numPoints x dimension) array of reference coordinates (see DMGetDimension())
28509d150b73SToby Isaac 
28519d150b73SToby Isaac   Output Parameters:
28529d150b73SToby Isaac . realCoords - (numPoints x coordinate dimension) array of coordinates (see DMGetCoordinateDim())
28531b266c99SBarry Smith 
28541b266c99SBarry Smith    Level: intermediate
28559d150b73SToby Isaac @*/
28569d150b73SToby Isaac PetscErrorCode DMPlexReferenceToCoordinates(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[])
28579d150b73SToby Isaac {
28589d150b73SToby Isaac   PetscInt       dimC, dimR, depth, cStart, cEnd, cEndInterior, i;
28599d150b73SToby Isaac   DM             coordDM = NULL;
28609d150b73SToby Isaac   Vec            coords;
28619d150b73SToby Isaac   PetscFE        fe = NULL;
28629d150b73SToby Isaac   PetscErrorCode ierr;
28639d150b73SToby Isaac 
28649d150b73SToby Isaac   PetscFunctionBegin;
28659d150b73SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
28669d150b73SToby Isaac   ierr = DMGetDimension(dm,&dimR);CHKERRQ(ierr);
28679d150b73SToby Isaac   ierr = DMGetCoordinateDim(dm,&dimC);CHKERRQ(ierr);
28689d150b73SToby Isaac   if (dimR <= 0 || dimC <= 0 || numPoints <= 0) PetscFunctionReturn(0);
28699d150b73SToby Isaac   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
28709d150b73SToby Isaac   ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
28719d150b73SToby Isaac   ierr = DMGetCoordinateDM(dm,&coordDM);CHKERRQ(ierr);
28729d150b73SToby Isaac   if (coordDM) {
28739d150b73SToby Isaac     PetscInt coordFields;
28749d150b73SToby Isaac 
28759d150b73SToby Isaac     ierr = DMGetNumFields(coordDM,&coordFields);CHKERRQ(ierr);
28769d150b73SToby Isaac     if (coordFields) {
28779d150b73SToby Isaac       PetscClassId id;
28789d150b73SToby Isaac       PetscObject  disc;
28799d150b73SToby Isaac 
28809d150b73SToby Isaac       ierr = DMGetField(coordDM,0,&disc);CHKERRQ(ierr);
28819d150b73SToby Isaac       ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
28829d150b73SToby Isaac       if (id == PETSCFE_CLASSID) {
28839d150b73SToby Isaac         fe = (PetscFE) disc;
28849d150b73SToby Isaac       }
28859d150b73SToby Isaac     }
28869d150b73SToby Isaac   }
28879d150b73SToby Isaac   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
28889d150b73SToby Isaac   ierr = DMPlexGetHybridBounds(dm,&cEndInterior,NULL,NULL,NULL);CHKERRQ(ierr);
28899d150b73SToby Isaac   cEnd = cEndInterior > 0 ? cEndInterior : cEnd;
28909d150b73SToby Isaac   if (cell < cStart || cell >= cEnd) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D not in cell range [%D,%D)",cell,cStart,cEnd);CHKERRQ(ierr);
28919d150b73SToby Isaac   if (!fe) { /* implicit discretization: affine or multilinear */
28929d150b73SToby Isaac     PetscInt  coneSize;
28939d150b73SToby Isaac     PetscBool isSimplex, isTensor;
28949d150b73SToby Isaac 
28959d150b73SToby Isaac     ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr);
28969d150b73SToby Isaac     isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE;
28979d150b73SToby Isaac     isTensor  = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE;
28989d150b73SToby Isaac     if (isSimplex) {
28999d150b73SToby Isaac       PetscReal detJ, *v0, *J;
29009d150b73SToby Isaac 
29019d150b73SToby Isaac       ierr = DMGetWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);CHKERRQ(ierr);
29029d150b73SToby Isaac       J    = &v0[dimC];
29039d150b73SToby Isaac       ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, NULL, &detJ);CHKERRQ(ierr);
29049d150b73SToby Isaac       for (i = 0; i < numPoints; i++) { /* Apply the inverse affine transformation for each point */
29059d150b73SToby Isaac         CoordinatesRefToReal(dimC, dimR, v0, J, &refCoords[dimR * i], &realCoords[dimC * i]);
29069d150b73SToby Isaac       }
29079d150b73SToby Isaac       ierr = DMRestoreWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);CHKERRQ(ierr);
29089d150b73SToby Isaac     } else if (isTensor) {
29099d150b73SToby Isaac       ierr = DMPlexReferenceToCoordinates_Tensor(coordDM, cell, numPoints, refCoords, realCoords, coords, dimC, dimR);CHKERRQ(ierr);
29109d150b73SToby Isaac     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unrecognized cone size %D",coneSize);
29119d150b73SToby Isaac   } else {
29129d150b73SToby Isaac     ierr = DMPlexReferenceToCoordinates_FE(coordDM, fe, cell, numPoints, refCoords, realCoords, coords, dimC, dimR);CHKERRQ(ierr);
29139d150b73SToby Isaac   }
2914d6143a4eSToby Isaac   PetscFunctionReturn(0);
2915d6143a4eSToby Isaac }
2916