xref: /petsc/src/dm/impls/plex/plexgeometry.c (revision eca9f51847c038c68d1253765f4fe9f2802a3d24)
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>
4af74b616SDave May #include <petsctime.h>
5ccd2543fSMatthew G Knepley 
63985bb02SVaclav Hapla /*@
73985bb02SVaclav Hapla   DMPlexFindVertices - Try to find DAG points based on their coordinates.
83985bb02SVaclav Hapla 
93985bb02SVaclav Hapla   Not Collective (provided DMGetCoordinatesLocalSetUp() has been called already)
103985bb02SVaclav Hapla 
113985bb02SVaclav Hapla   Input Parameters:
123985bb02SVaclav Hapla + dm - The DMPlex object
133985bb02SVaclav Hapla . npoints - The number of sought points
143985bb02SVaclav Hapla . coords - The array of coordinates of the sought points
153985bb02SVaclav Hapla - eps - The tolerance or PETSC_DEFAULT
163985bb02SVaclav Hapla 
173985bb02SVaclav Hapla   Output Parameters:
183985bb02SVaclav Hapla . dagPoints - The array of found DAG points, or -1 if not found
193985bb02SVaclav Hapla 
203985bb02SVaclav Hapla   Level: intermediate
213985bb02SVaclav Hapla 
223985bb02SVaclav Hapla   Notes:
233985bb02SVaclav Hapla   The length of the array coords must be npoints * dim where dim is the spatial dimension returned by DMGetDimension().
243985bb02SVaclav Hapla 
253985bb02SVaclav Hapla   The output array dagPoints is NOT newly allocated; the user must pass an array of length npoints.
263985bb02SVaclav Hapla 
273985bb02SVaclav Hapla   Each rank does the search independently; a nonnegative value is returned only if this rank's local DMPlex portion contains the point.
283985bb02SVaclav Hapla 
293985bb02SVaclav Hapla   The tolerance is interpreted as the maximum Euclidean (L2) distance of the sought point from the specified coordinates.
303985bb02SVaclav Hapla 
31335ef845SVaclav Hapla   Complexity of this function is currently O(mn) with m number of vertices to find and n number of vertices in the local mesh. This could probably be improved.
32335ef845SVaclav Hapla 
333985bb02SVaclav Hapla .seealso: DMPlexCreate(), DMGetCoordinates()
343985bb02SVaclav Hapla @*/
353985bb02SVaclav Hapla PetscErrorCode DMPlexFindVertices(DM dm, PetscInt npoints, const PetscReal coord[], PetscReal eps, PetscInt dagPoints[])
363985bb02SVaclav Hapla {
37335ef845SVaclav Hapla   PetscInt          c, dim, i, j, o, p, vStart, vEnd;
383985bb02SVaclav Hapla   Vec               allCoordsVec;
393985bb02SVaclav Hapla   const PetscScalar *allCoords;
403985bb02SVaclav Hapla   PetscReal         norm;
413985bb02SVaclav Hapla   PetscErrorCode    ierr;
423985bb02SVaclav Hapla 
433985bb02SVaclav Hapla   PetscFunctionBegin;
443985bb02SVaclav Hapla   if (eps < 0) eps = PETSC_SQRT_MACHINE_EPSILON;
453985bb02SVaclav Hapla   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
463985bb02SVaclav Hapla   ierr = DMGetCoordinatesLocal(dm, &allCoordsVec);CHKERRQ(ierr);
473985bb02SVaclav Hapla   ierr = VecGetArrayRead(allCoordsVec, &allCoords);CHKERRQ(ierr);
483985bb02SVaclav Hapla   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
49335ef845SVaclav Hapla #if defined(PETSC_USE_DEBUG)
50335ef845SVaclav Hapla   /* check coordinate section is consistent with DM dimension */
51335ef845SVaclav Hapla   {
52335ef845SVaclav Hapla     PetscSection      cs;
53335ef845SVaclav Hapla     PetscInt          ndof;
54335ef845SVaclav Hapla 
55335ef845SVaclav Hapla     ierr = DMGetCoordinateSection(dm, &cs);CHKERRQ(ierr);
563985bb02SVaclav Hapla     for (p = vStart; p < vEnd; p++) {
573985bb02SVaclav Hapla       ierr = PetscSectionGetDof(cs, p, &ndof);CHKERRQ(ierr);
583985bb02SVaclav Hapla       if (PetscUnlikely(ndof != dim)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "point %D: ndof = %D != %D = dim", p, ndof, dim);
59335ef845SVaclav Hapla     }
60335ef845SVaclav Hapla   }
61335ef845SVaclav Hapla #endif
62*eca9f518SVaclav Hapla   if (eps == 0.0) {
63*eca9f518SVaclav Hapla     for (i=0,j=0; i < npoints; i++,j+=dim) {
64*eca9f518SVaclav Hapla       dagPoints[i] = -1;
65*eca9f518SVaclav Hapla       for (p = vStart,o=0; p < vEnd; p++,o+=dim) {
66*eca9f518SVaclav Hapla         for (c = 0; c < dim; c++) {
67*eca9f518SVaclav Hapla           if (coord[j+c] != PetscRealPart(allCoords[o+c])) break;
68*eca9f518SVaclav Hapla         }
69*eca9f518SVaclav Hapla         if (c == dim) {
70*eca9f518SVaclav Hapla           dagPoints[i] = p;
71*eca9f518SVaclav Hapla           break;
72*eca9f518SVaclav Hapla         }
73*eca9f518SVaclav Hapla       }
74*eca9f518SVaclav Hapla     }
75*eca9f518SVaclav Hapla     ierr = VecRestoreArrayRead(allCoordsVec, &allCoords);CHKERRQ(ierr);
76*eca9f518SVaclav Hapla     PetscFunctionReturn(0);
77*eca9f518SVaclav Hapla   }
78335ef845SVaclav Hapla   for (i=0,j=0; i < npoints; i++,j+=dim) {
79335ef845SVaclav Hapla     dagPoints[i] = -1;
80335ef845SVaclav Hapla     for (p = vStart,o=0; p < vEnd; p++,o+=dim) {
813985bb02SVaclav Hapla       norm = 0.0;
823985bb02SVaclav Hapla       for (c = 0; c < dim; c++) {
833985bb02SVaclav Hapla         norm += PetscSqr(coord[j+c] - PetscRealPart(allCoords[o+c]));
843985bb02SVaclav Hapla       }
853985bb02SVaclav Hapla       norm = PetscSqrtReal(norm);
863985bb02SVaclav Hapla       if (norm <= eps) {
873985bb02SVaclav Hapla         dagPoints[i] = p;
883985bb02SVaclav Hapla         break;
893985bb02SVaclav Hapla       }
903985bb02SVaclav Hapla     }
913985bb02SVaclav Hapla   }
923985bb02SVaclav Hapla   ierr = VecRestoreArrayRead(allCoordsVec, &allCoords);CHKERRQ(ierr);
933985bb02SVaclav Hapla   PetscFunctionReturn(0);
943985bb02SVaclav Hapla }
953985bb02SVaclav Hapla 
96fea14342SMatthew G. Knepley static PetscErrorCode DMPlexGetLineIntersection_2D_Internal(const PetscReal segmentA[], const PetscReal segmentB[], PetscReal intersection[], PetscBool *hasIntersection)
97fea14342SMatthew G. Knepley {
98fea14342SMatthew G. Knepley   const PetscReal p0_x  = segmentA[0*2+0];
99fea14342SMatthew G. Knepley   const PetscReal p0_y  = segmentA[0*2+1];
100fea14342SMatthew G. Knepley   const PetscReal p1_x  = segmentA[1*2+0];
101fea14342SMatthew G. Knepley   const PetscReal p1_y  = segmentA[1*2+1];
102fea14342SMatthew G. Knepley   const PetscReal p2_x  = segmentB[0*2+0];
103fea14342SMatthew G. Knepley   const PetscReal p2_y  = segmentB[0*2+1];
104fea14342SMatthew G. Knepley   const PetscReal p3_x  = segmentB[1*2+0];
105fea14342SMatthew G. Knepley   const PetscReal p3_y  = segmentB[1*2+1];
106fea14342SMatthew G. Knepley   const PetscReal s1_x  = p1_x - p0_x;
107fea14342SMatthew G. Knepley   const PetscReal s1_y  = p1_y - p0_y;
108fea14342SMatthew G. Knepley   const PetscReal s2_x  = p3_x - p2_x;
109fea14342SMatthew G. Knepley   const PetscReal s2_y  = p3_y - p2_y;
110fea14342SMatthew G. Knepley   const PetscReal denom = (-s2_x * s1_y + s1_x * s2_y);
111fea14342SMatthew G. Knepley 
112fea14342SMatthew G. Knepley   PetscFunctionBegin;
113fea14342SMatthew G. Knepley   *hasIntersection = PETSC_FALSE;
114fea14342SMatthew G. Knepley   /* Non-parallel lines */
115fea14342SMatthew G. Knepley   if (denom != 0.0) {
116fea14342SMatthew G. Knepley     const PetscReal s = (-s1_y * (p0_x - p2_x) + s1_x * (p0_y - p2_y)) / denom;
117fea14342SMatthew G. Knepley     const PetscReal t = ( s2_x * (p0_y - p2_y) - s2_y * (p0_x - p2_x)) / denom;
118fea14342SMatthew G. Knepley 
119fea14342SMatthew G. Knepley     if (s >= 0 && s <= 1 && t >= 0 && t <= 1) {
120fea14342SMatthew G. Knepley       *hasIntersection = PETSC_TRUE;
121fea14342SMatthew G. Knepley       if (intersection) {
122fea14342SMatthew G. Knepley         intersection[0] = p0_x + (t * s1_x);
123fea14342SMatthew G. Knepley         intersection[1] = p0_y + (t * s1_y);
124fea14342SMatthew G. Knepley       }
125fea14342SMatthew G. Knepley     }
126fea14342SMatthew G. Knepley   }
127fea14342SMatthew G. Knepley   PetscFunctionReturn(0);
128fea14342SMatthew G. Knepley }
129fea14342SMatthew G. Knepley 
130ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
131ccd2543fSMatthew G Knepley {
132ccd2543fSMatthew G Knepley   const PetscInt  embedDim = 2;
133f5ebc837SMatthew G. Knepley   const PetscReal eps      = PETSC_SQRT_MACHINE_EPSILON;
134ccd2543fSMatthew G Knepley   PetscReal       x        = PetscRealPart(point[0]);
135ccd2543fSMatthew G Knepley   PetscReal       y        = PetscRealPart(point[1]);
136ccd2543fSMatthew G Knepley   PetscReal       v0[2], J[4], invJ[4], detJ;
137ccd2543fSMatthew G Knepley   PetscReal       xi, eta;
138ccd2543fSMatthew G Knepley   PetscErrorCode  ierr;
139ccd2543fSMatthew G Knepley 
140ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1418e0841e0SMatthew G. Knepley   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
142ccd2543fSMatthew G Knepley   xi  = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]);
143ccd2543fSMatthew G Knepley   eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]);
144ccd2543fSMatthew G Knepley 
145f5ebc837SMatthew G. Knepley   if ((xi >= -eps) && (eta >= -eps) && (xi + eta <= 2.0+eps)) *cell = c;
146c1496c66SMatthew G. Knepley   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
147ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
148ccd2543fSMatthew G Knepley }
149ccd2543fSMatthew G Knepley 
15062a38674SMatthew G. Knepley static PetscErrorCode DMPlexClosestPoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscReal cpoint[])
15162a38674SMatthew G. Knepley {
15262a38674SMatthew G. Knepley   const PetscInt  embedDim = 2;
15362a38674SMatthew G. Knepley   PetscReal       x        = PetscRealPart(point[0]);
15462a38674SMatthew G. Knepley   PetscReal       y        = PetscRealPart(point[1]);
15562a38674SMatthew G. Knepley   PetscReal       v0[2], J[4], invJ[4], detJ;
15662a38674SMatthew G. Knepley   PetscReal       xi, eta, r;
15762a38674SMatthew G. Knepley   PetscErrorCode  ierr;
15862a38674SMatthew G. Knepley 
15962a38674SMatthew G. Knepley   PetscFunctionBegin;
16062a38674SMatthew G. Knepley   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
16162a38674SMatthew G. Knepley   xi  = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]);
16262a38674SMatthew G. Knepley   eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]);
16362a38674SMatthew G. Knepley 
16462a38674SMatthew G. Knepley   xi  = PetscMax(xi,  0.0);
16562a38674SMatthew G. Knepley   eta = PetscMax(eta, 0.0);
16662a38674SMatthew G. Knepley   if (xi + eta > 2.0) {
16762a38674SMatthew G. Knepley     r    = (xi + eta)/2.0;
16862a38674SMatthew G. Knepley     xi  /= r;
16962a38674SMatthew G. Knepley     eta /= r;
17062a38674SMatthew G. Knepley   }
17162a38674SMatthew G. Knepley   cpoint[0] = J[0*embedDim+0]*xi + J[0*embedDim+1]*eta + v0[0];
17262a38674SMatthew G. Knepley   cpoint[1] = J[1*embedDim+0]*xi + J[1*embedDim+1]*eta + v0[1];
17362a38674SMatthew G. Knepley   PetscFunctionReturn(0);
17462a38674SMatthew G. Knepley }
17562a38674SMatthew G. Knepley 
176ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_General_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
177ccd2543fSMatthew G Knepley {
178ccd2543fSMatthew G Knepley   PetscSection       coordSection;
179ccd2543fSMatthew G Knepley   Vec             coordsLocal;
180a1e44745SMatthew G. Knepley   PetscScalar    *coords = NULL;
181ccd2543fSMatthew G Knepley   const PetscInt  faces[8]  = {0, 1, 1, 2, 2, 3, 3, 0};
182ccd2543fSMatthew G Knepley   PetscReal       x         = PetscRealPart(point[0]);
183ccd2543fSMatthew G Knepley   PetscReal       y         = PetscRealPart(point[1]);
184ccd2543fSMatthew G Knepley   PetscInt        crossings = 0, f;
185ccd2543fSMatthew G Knepley   PetscErrorCode  ierr;
186ccd2543fSMatthew G Knepley 
187ccd2543fSMatthew G Knepley   PetscFunctionBegin;
188ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
18969d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
190ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
191ccd2543fSMatthew G Knepley   for (f = 0; f < 4; ++f) {
192ccd2543fSMatthew G Knepley     PetscReal x_i   = PetscRealPart(coords[faces[2*f+0]*2+0]);
193ccd2543fSMatthew G Knepley     PetscReal y_i   = PetscRealPart(coords[faces[2*f+0]*2+1]);
194ccd2543fSMatthew G Knepley     PetscReal x_j   = PetscRealPart(coords[faces[2*f+1]*2+0]);
195ccd2543fSMatthew G Knepley     PetscReal y_j   = PetscRealPart(coords[faces[2*f+1]*2+1]);
196ccd2543fSMatthew G Knepley     PetscReal slope = (y_j - y_i) / (x_j - x_i);
197ccd2543fSMatthew G Knepley     PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE;
198ccd2543fSMatthew G Knepley     PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE;
199ccd2543fSMatthew G Knepley     PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE;
200ccd2543fSMatthew G Knepley     if ((cond1 || cond2)  && above) ++crossings;
201ccd2543fSMatthew G Knepley   }
202ccd2543fSMatthew G Knepley   if (crossings % 2) *cell = c;
203c1496c66SMatthew G. Knepley   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
204ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
205ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
206ccd2543fSMatthew G Knepley }
207ccd2543fSMatthew G Knepley 
208ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
209ccd2543fSMatthew G Knepley {
210ccd2543fSMatthew G Knepley   const PetscInt embedDim = 3;
211ccd2543fSMatthew G Knepley   PetscReal      v0[3], J[9], invJ[9], detJ;
212ccd2543fSMatthew G Knepley   PetscReal      x = PetscRealPart(point[0]);
213ccd2543fSMatthew G Knepley   PetscReal      y = PetscRealPart(point[1]);
214ccd2543fSMatthew G Knepley   PetscReal      z = PetscRealPart(point[2]);
215ccd2543fSMatthew G Knepley   PetscReal      xi, eta, zeta;
216ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
217ccd2543fSMatthew G Knepley 
218ccd2543fSMatthew G Knepley   PetscFunctionBegin;
2198e0841e0SMatthew G. Knepley   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
220ccd2543fSMatthew G Knepley   xi   = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]);
221ccd2543fSMatthew G Knepley   eta  = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]);
222ccd2543fSMatthew G Knepley   zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]);
223ccd2543fSMatthew G Knepley 
224ccd2543fSMatthew G Knepley   if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) *cell = c;
225c1496c66SMatthew G. Knepley   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
226ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
227ccd2543fSMatthew G Knepley }
228ccd2543fSMatthew G Knepley 
229ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
230ccd2543fSMatthew G Knepley {
231ccd2543fSMatthew G Knepley   PetscSection   coordSection;
232ccd2543fSMatthew G Knepley   Vec            coordsLocal;
233872a9804SMatthew G. Knepley   PetscScalar   *coords = NULL;
234fb150da6SMatthew G. Knepley   const PetscInt faces[24] = {0, 3, 2, 1,  5, 4, 7, 6,  3, 0, 4, 5,
235fb150da6SMatthew G. Knepley                               1, 2, 6, 7,  3, 5, 6, 2,  0, 1, 7, 4};
236ccd2543fSMatthew G Knepley   PetscBool      found = PETSC_TRUE;
237ccd2543fSMatthew G Knepley   PetscInt       f;
238ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
239ccd2543fSMatthew G Knepley 
240ccd2543fSMatthew G Knepley   PetscFunctionBegin;
241ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
24269d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
243ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
244ccd2543fSMatthew G Knepley   for (f = 0; f < 6; ++f) {
245ccd2543fSMatthew G Knepley     /* Check the point is under plane */
246ccd2543fSMatthew G Knepley     /*   Get face normal */
247ccd2543fSMatthew G Knepley     PetscReal v_i[3];
248ccd2543fSMatthew G Knepley     PetscReal v_j[3];
249ccd2543fSMatthew G Knepley     PetscReal normal[3];
250ccd2543fSMatthew G Knepley     PetscReal pp[3];
251ccd2543fSMatthew G Knepley     PetscReal dot;
252ccd2543fSMatthew G Knepley 
253ccd2543fSMatthew G Knepley     v_i[0]    = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]);
254ccd2543fSMatthew G Knepley     v_i[1]    = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]);
255ccd2543fSMatthew G Knepley     v_i[2]    = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]);
256ccd2543fSMatthew G Knepley     v_j[0]    = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]);
257ccd2543fSMatthew G Knepley     v_j[1]    = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]);
258ccd2543fSMatthew G Knepley     v_j[2]    = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]);
259ccd2543fSMatthew G Knepley     normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1];
260ccd2543fSMatthew G Knepley     normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2];
261ccd2543fSMatthew G Knepley     normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0];
262ccd2543fSMatthew G Knepley     pp[0]     = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]);
263ccd2543fSMatthew G Knepley     pp[1]     = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]);
264ccd2543fSMatthew G Knepley     pp[2]     = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]);
265ccd2543fSMatthew G Knepley     dot       = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2];
266ccd2543fSMatthew G Knepley 
267ccd2543fSMatthew G Knepley     /* Check that projected point is in face (2D location problem) */
268ccd2543fSMatthew G Knepley     if (dot < 0.0) {
269ccd2543fSMatthew G Knepley       found = PETSC_FALSE;
270ccd2543fSMatthew G Knepley       break;
271ccd2543fSMatthew G Knepley     }
272ccd2543fSMatthew G Knepley   }
273ccd2543fSMatthew G Knepley   if (found) *cell = c;
274c1496c66SMatthew G. Knepley   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
275ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
276ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
277ccd2543fSMatthew G Knepley }
278ccd2543fSMatthew G Knepley 
279c4eade1cSMatthew G. Knepley static PetscErrorCode PetscGridHashInitialize_Internal(PetscGridHash box, PetscInt dim, const PetscScalar point[])
280c4eade1cSMatthew G. Knepley {
281c4eade1cSMatthew G. Knepley   PetscInt d;
282c4eade1cSMatthew G. Knepley 
283c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
284c4eade1cSMatthew G. Knepley   box->dim = dim;
285c4eade1cSMatthew G. Knepley   for (d = 0; d < dim; ++d) box->lower[d] = box->upper[d] = PetscRealPart(point[d]);
286c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
287c4eade1cSMatthew G. Knepley }
288c4eade1cSMatthew G. Knepley 
289c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashCreate(MPI_Comm comm, PetscInt dim, const PetscScalar point[], PetscGridHash *box)
290c4eade1cSMatthew G. Knepley {
291c4eade1cSMatthew G. Knepley   PetscErrorCode ierr;
292c4eade1cSMatthew G. Knepley 
293c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
294c4eade1cSMatthew G. Knepley   ierr = PetscMalloc1(1, box);CHKERRQ(ierr);
295c4eade1cSMatthew G. Knepley   ierr = PetscGridHashInitialize_Internal(*box, dim, point);CHKERRQ(ierr);
296c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
297c4eade1cSMatthew G. Knepley }
298c4eade1cSMatthew G. Knepley 
299c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashEnlarge(PetscGridHash box, const PetscScalar point[])
300c4eade1cSMatthew G. Knepley {
301c4eade1cSMatthew G. Knepley   PetscInt d;
302c4eade1cSMatthew G. Knepley 
303c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
304c4eade1cSMatthew G. Knepley   for (d = 0; d < box->dim; ++d) {
305c4eade1cSMatthew G. Knepley     box->lower[d] = PetscMin(box->lower[d], PetscRealPart(point[d]));
306c4eade1cSMatthew G. Knepley     box->upper[d] = PetscMax(box->upper[d], PetscRealPart(point[d]));
307c4eade1cSMatthew G. Knepley   }
308c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
309c4eade1cSMatthew G. Knepley }
310c4eade1cSMatthew G. Knepley 
31162a38674SMatthew G. Knepley /*
31262a38674SMatthew G. Knepley   PetscGridHashSetGrid - Divide the grid into boxes
31362a38674SMatthew G. Knepley 
31462a38674SMatthew G. Knepley   Not collective
31562a38674SMatthew G. Knepley 
31662a38674SMatthew G. Knepley   Input Parameters:
31762a38674SMatthew G. Knepley + box - The grid hash object
31862a38674SMatthew G. Knepley . n   - The number of boxes in each dimension, or PETSC_DETERMINE
31962a38674SMatthew G. Knepley - h   - The box size in each dimension, only used if n[d] == PETSC_DETERMINE
32062a38674SMatthew G. Knepley 
32162a38674SMatthew G. Knepley   Level: developer
32262a38674SMatthew G. Knepley 
32362a38674SMatthew G. Knepley .seealso: PetscGridHashCreate()
32462a38674SMatthew G. Knepley */
325c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashSetGrid(PetscGridHash box, const PetscInt n[], const PetscReal h[])
326c4eade1cSMatthew G. Knepley {
327c4eade1cSMatthew G. Knepley   PetscInt d;
328c4eade1cSMatthew G. Knepley 
329c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
330c4eade1cSMatthew G. Knepley   for (d = 0; d < box->dim; ++d) {
331c4eade1cSMatthew G. Knepley     box->extent[d] = box->upper[d] - box->lower[d];
332c4eade1cSMatthew G. Knepley     if (n[d] == PETSC_DETERMINE) {
333c4eade1cSMatthew G. Knepley       box->h[d] = h[d];
334c4eade1cSMatthew G. Knepley       box->n[d] = PetscCeilReal(box->extent[d]/h[d]);
335c4eade1cSMatthew G. Knepley     } else {
336c4eade1cSMatthew G. Knepley       box->n[d] = n[d];
337c4eade1cSMatthew G. Knepley       box->h[d] = box->extent[d]/n[d];
338c4eade1cSMatthew G. Knepley     }
339c4eade1cSMatthew G. Knepley   }
340c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
341c4eade1cSMatthew G. Knepley }
342c4eade1cSMatthew G. Knepley 
34362a38674SMatthew G. Knepley /*
34462a38674SMatthew G. Knepley   PetscGridHashGetEnclosingBox - Find the grid boxes containing each input point
34562a38674SMatthew G. Knepley 
34662a38674SMatthew G. Knepley   Not collective
34762a38674SMatthew G. Knepley 
34862a38674SMatthew G. Knepley   Input Parameters:
34962a38674SMatthew G. Knepley + box       - The grid hash object
35062a38674SMatthew G. Knepley . numPoints - The number of input points
35162a38674SMatthew G. Knepley - points    - The input point coordinates
35262a38674SMatthew G. Knepley 
35362a38674SMatthew G. Knepley   Output Parameters:
35462a38674SMatthew G. Knepley + dboxes    - An array of numPoints*dim integers expressing the enclosing box as (i_0, i_1, ..., i_dim)
35562a38674SMatthew G. Knepley - boxes     - An array of numPoints integers expressing the enclosing box as single number, or NULL
35662a38674SMatthew G. Knepley 
35762a38674SMatthew G. Knepley   Level: developer
35862a38674SMatthew G. Knepley 
35962a38674SMatthew G. Knepley .seealso: PetscGridHashCreate()
36062a38674SMatthew G. Knepley */
3611c6dfc3eSMatthew G. Knepley PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[])
362c4eade1cSMatthew G. Knepley {
363c4eade1cSMatthew G. Knepley   const PetscReal *lower = box->lower;
364c4eade1cSMatthew G. Knepley   const PetscReal *upper = box->upper;
365c4eade1cSMatthew G. Knepley   const PetscReal *h     = box->h;
366c4eade1cSMatthew G. Knepley   const PetscInt  *n     = box->n;
367c4eade1cSMatthew G. Knepley   const PetscInt   dim   = box->dim;
368c4eade1cSMatthew G. Knepley   PetscInt         d, p;
369c4eade1cSMatthew G. Knepley 
370c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
371c4eade1cSMatthew G. Knepley   for (p = 0; p < numPoints; ++p) {
372c4eade1cSMatthew G. Knepley     for (d = 0; d < dim; ++d) {
3731c6dfc3eSMatthew G. Knepley       PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]);
374c4eade1cSMatthew G. Knepley 
3751c6dfc3eSMatthew G. Knepley       if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1;
3762a705cacSMatthew G. Knepley       if (dbox == -1   && PetscAbsReal(PetscRealPart(points[p*dim+d]) - lower[d]) < 1.0e-9) dbox = 0;
377c4eade1cSMatthew 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",
378087ef6b2SMatthew G. Knepley                                              p, (double) PetscRealPart(points[p*dim+0]), dim > 1 ? (double) PetscRealPart(points[p*dim+1]) : 0.0, dim > 2 ? (double) PetscRealPart(points[p*dim+2]) : 0.0);
379c4eade1cSMatthew G. Knepley       dboxes[p*dim+d] = dbox;
380c4eade1cSMatthew G. Knepley     }
381c4eade1cSMatthew G. Knepley     if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1];
382c4eade1cSMatthew G. Knepley   }
383c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
384c4eade1cSMatthew G. Knepley }
385c4eade1cSMatthew G. Knepley 
386af74b616SDave May /*
387af74b616SDave May  PetscGridHashGetEnclosingBoxQuery - Find the grid boxes containing each input point
388af74b616SDave May 
389af74b616SDave May  Not collective
390af74b616SDave May 
391af74b616SDave May   Input Parameters:
392af74b616SDave May + box       - The grid hash object
393af74b616SDave May . numPoints - The number of input points
394af74b616SDave May - points    - The input point coordinates
395af74b616SDave May 
396af74b616SDave May   Output Parameters:
397af74b616SDave May + dboxes    - An array of numPoints*dim integers expressing the enclosing box as (i_0, i_1, ..., i_dim)
398af74b616SDave May . boxes     - An array of numPoints integers expressing the enclosing box as single number, or NULL
399af74b616SDave May - found     - Flag indicating if point was located within a box
400af74b616SDave May 
401af74b616SDave May   Level: developer
402af74b616SDave May 
403af74b616SDave May .seealso: PetscGridHashGetEnclosingBox()
404af74b616SDave May */
405af74b616SDave May PetscErrorCode PetscGridHashGetEnclosingBoxQuery(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[],PetscBool *found)
406af74b616SDave May {
407af74b616SDave May   const PetscReal *lower = box->lower;
408af74b616SDave May   const PetscReal *upper = box->upper;
409af74b616SDave May   const PetscReal *h     = box->h;
410af74b616SDave May   const PetscInt  *n     = box->n;
411af74b616SDave May   const PetscInt   dim   = box->dim;
412af74b616SDave May   PetscInt         d, p;
413af74b616SDave May 
414af74b616SDave May   PetscFunctionBegin;
415af74b616SDave May   *found = PETSC_FALSE;
416af74b616SDave May   for (p = 0; p < numPoints; ++p) {
417af74b616SDave May     for (d = 0; d < dim; ++d) {
418af74b616SDave May       PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]);
419af74b616SDave May 
420af74b616SDave May       if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1;
421af74b616SDave May       if (dbox < 0 || dbox >= n[d]) {
422af74b616SDave May         PetscFunctionReturn(0);
423af74b616SDave May       }
424af74b616SDave May       dboxes[p*dim+d] = dbox;
425af74b616SDave May     }
426af74b616SDave May     if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1];
427af74b616SDave May   }
428af74b616SDave May   *found = PETSC_TRUE;
429af74b616SDave May   PetscFunctionReturn(0);
430af74b616SDave May }
431af74b616SDave May 
432c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashDestroy(PetscGridHash *box)
433c4eade1cSMatthew G. Knepley {
434c4eade1cSMatthew G. Knepley   PetscErrorCode ierr;
435c4eade1cSMatthew G. Knepley 
436c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
437c4eade1cSMatthew G. Knepley   if (*box) {
438c4eade1cSMatthew G. Knepley     ierr = PetscSectionDestroy(&(*box)->cellSection);CHKERRQ(ierr);
439c4eade1cSMatthew G. Knepley     ierr = ISDestroy(&(*box)->cells);CHKERRQ(ierr);
440c4eade1cSMatthew G. Knepley     ierr = DMLabelDestroy(&(*box)->cellsSparse);CHKERRQ(ierr);
441c4eade1cSMatthew G. Knepley   }
442c4eade1cSMatthew G. Knepley   ierr = PetscFree(*box);CHKERRQ(ierr);
443c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
444c4eade1cSMatthew G. Knepley }
445c4eade1cSMatthew G. Knepley 
446cafe43deSMatthew G. Knepley PetscErrorCode DMPlexLocatePoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cellStart, PetscInt *cell)
447cafe43deSMatthew G. Knepley {
448cafe43deSMatthew G. Knepley   PetscInt       coneSize;
449cafe43deSMatthew G. Knepley   PetscErrorCode ierr;
450cafe43deSMatthew G. Knepley 
451cafe43deSMatthew G. Knepley   PetscFunctionBegin;
452cafe43deSMatthew G. Knepley   switch (dim) {
453cafe43deSMatthew G. Knepley   case 2:
454cafe43deSMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr);
455cafe43deSMatthew G. Knepley     switch (coneSize) {
456cafe43deSMatthew G. Knepley     case 3:
457cafe43deSMatthew G. Knepley       ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
458cafe43deSMatthew G. Knepley       break;
459cafe43deSMatthew G. Knepley     case 4:
460cafe43deSMatthew G. Knepley       ierr = DMPlexLocatePoint_General_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
461cafe43deSMatthew G. Knepley       break;
462cafe43deSMatthew G. Knepley     default:
463cafe43deSMatthew G. Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize);
464cafe43deSMatthew G. Knepley     }
465cafe43deSMatthew G. Knepley     break;
466cafe43deSMatthew G. Knepley   case 3:
467cafe43deSMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr);
468cafe43deSMatthew G. Knepley     switch (coneSize) {
469cafe43deSMatthew G. Knepley     case 4:
470cafe43deSMatthew G. Knepley       ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
471cafe43deSMatthew G. Knepley       break;
472cafe43deSMatthew G. Knepley     case 6:
473cafe43deSMatthew G. Knepley       ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);
474cafe43deSMatthew G. Knepley       break;
475cafe43deSMatthew G. Knepley     default:
476cafe43deSMatthew G. Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize);
477cafe43deSMatthew G. Knepley     }
478cafe43deSMatthew G. Knepley     break;
479cafe43deSMatthew G. Knepley   default:
480cafe43deSMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %D", dim);
481cafe43deSMatthew G. Knepley   }
482cafe43deSMatthew G. Knepley   PetscFunctionReturn(0);
483cafe43deSMatthew G. Knepley }
484cafe43deSMatthew G. Knepley 
48562a38674SMatthew G. Knepley /*
48662a38674SMatthew G. Knepley   DMPlexClosestPoint_Internal - Returns the closest point in the cell to the given point
48762a38674SMatthew G. Knepley */
48862a38674SMatthew G. Knepley PetscErrorCode DMPlexClosestPoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cell, PetscReal cpoint[])
48962a38674SMatthew G. Knepley {
49062a38674SMatthew G. Knepley   PetscInt       coneSize;
49162a38674SMatthew G. Knepley   PetscErrorCode ierr;
49262a38674SMatthew G. Knepley 
49362a38674SMatthew G. Knepley   PetscFunctionBegin;
49462a38674SMatthew G. Knepley   switch (dim) {
49562a38674SMatthew G. Knepley   case 2:
49662a38674SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
49762a38674SMatthew G. Knepley     switch (coneSize) {
49862a38674SMatthew G. Knepley     case 3:
49962a38674SMatthew G. Knepley       ierr = DMPlexClosestPoint_Simplex_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
50062a38674SMatthew G. Knepley       break;
50162a38674SMatthew G. Knepley #if 0
50262a38674SMatthew G. Knepley     case 4:
50362a38674SMatthew G. Knepley       ierr = DMPlexClosestPoint_General_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
50462a38674SMatthew G. Knepley       break;
50562a38674SMatthew G. Knepley #endif
50662a38674SMatthew G. Knepley     default:
50762a38674SMatthew G. Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize);
50862a38674SMatthew G. Knepley     }
50962a38674SMatthew G. Knepley     break;
51062a38674SMatthew G. Knepley #if 0
51162a38674SMatthew G. Knepley   case 3:
51262a38674SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
51362a38674SMatthew G. Knepley     switch (coneSize) {
51462a38674SMatthew G. Knepley     case 4:
51562a38674SMatthew G. Knepley       ierr = DMPlexClosestPoint_Simplex_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
51662a38674SMatthew G. Knepley       break;
51762a38674SMatthew G. Knepley     case 6:
51862a38674SMatthew G. Knepley       ierr = DMPlexClosestPoint_General_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);
51962a38674SMatthew G. Knepley       break;
52062a38674SMatthew G. Knepley     default:
52162a38674SMatthew G. Knepley       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize);
52262a38674SMatthew G. Knepley     }
52362a38674SMatthew G. Knepley     break;
52462a38674SMatthew G. Knepley #endif
52562a38674SMatthew G. Knepley   default:
52662a38674SMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for mesh dimension %D", dim);
52762a38674SMatthew G. Knepley   }
52862a38674SMatthew G. Knepley   PetscFunctionReturn(0);
52962a38674SMatthew G. Knepley }
53062a38674SMatthew G. Knepley 
53162a38674SMatthew G. Knepley /*
53262a38674SMatthew G. Knepley   DMPlexComputeGridHash_Internal - Create a grid hash structure covering the Plex
53362a38674SMatthew G. Knepley 
534d083f849SBarry Smith   Collective on dm
53562a38674SMatthew G. Knepley 
53662a38674SMatthew G. Knepley   Input Parameter:
53762a38674SMatthew G. Knepley . dm - The Plex
53862a38674SMatthew G. Knepley 
53962a38674SMatthew G. Knepley   Output Parameter:
54062a38674SMatthew G. Knepley . localBox - The grid hash object
54162a38674SMatthew G. Knepley 
54262a38674SMatthew G. Knepley   Level: developer
54362a38674SMatthew G. Knepley 
54462a38674SMatthew G. Knepley .seealso: PetscGridHashCreate(), PetscGridHashGetEnclosingBox()
54562a38674SMatthew G. Knepley */
546cafe43deSMatthew G. Knepley PetscErrorCode DMPlexComputeGridHash_Internal(DM dm, PetscGridHash *localBox)
547cafe43deSMatthew G. Knepley {
548cafe43deSMatthew G. Knepley   MPI_Comm           comm;
549cafe43deSMatthew G. Knepley   PetscGridHash      lbox;
550cafe43deSMatthew G. Knepley   Vec                coordinates;
551cafe43deSMatthew G. Knepley   PetscSection       coordSection;
552cafe43deSMatthew G. Knepley   Vec                coordsLocal;
553cafe43deSMatthew G. Knepley   const PetscScalar *coords;
554722d0f5cSMatthew G. Knepley   PetscInt          *dboxes, *boxes;
555cafe43deSMatthew G. Knepley   PetscInt           n[3] = {10, 10, 10};
5561d0c6c94SMatthew G. Knepley   PetscInt           dim, N, cStart, cEnd, cMax, c, i;
557cafe43deSMatthew G. Knepley   PetscErrorCode     ierr;
558cafe43deSMatthew G. Knepley 
559cafe43deSMatthew G. Knepley   PetscFunctionBegin;
560cafe43deSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
561cafe43deSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
562cafe43deSMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
5635b3353d8SMatthew G. Knepley   if (dim != 2) SETERRQ(comm, PETSC_ERR_SUP, "I have only coded this for 2D");
564cafe43deSMatthew G. Knepley   ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr);
565cafe43deSMatthew G. Knepley   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
566cafe43deSMatthew G. Knepley   ierr = PetscGridHashCreate(comm, dim, coords, &lbox);CHKERRQ(ierr);
567cafe43deSMatthew G. Knepley   for (i = 0; i < N; i += dim) {ierr = PetscGridHashEnlarge(lbox, &coords[i]);CHKERRQ(ierr);}
568cafe43deSMatthew G. Knepley   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
5699cb35068SDave May   ierr = PetscOptionsGetInt(NULL,NULL,"-dm_plex_hash_box_nijk",&n[0],NULL);CHKERRQ(ierr);
5709cb35068SDave May   n[1] = n[0];
5719cb35068SDave May   n[2] = n[0];
572cafe43deSMatthew G. Knepley   ierr = PetscGridHashSetGrid(lbox, n, NULL);CHKERRQ(ierr);
573cafe43deSMatthew G. Knepley #if 0
574cafe43deSMatthew G. Knepley   /* Could define a custom reduction to merge these */
575b2566f29SBarry Smith   ierr = MPIU_Allreduce(lbox->lower, gbox->lower, 3, MPIU_REAL, MPI_MIN, comm);CHKERRQ(ierr);
576b2566f29SBarry Smith   ierr = MPIU_Allreduce(lbox->upper, gbox->upper, 3, MPIU_REAL, MPI_MAX, comm);CHKERRQ(ierr);
577cafe43deSMatthew G. Knepley #endif
578cafe43deSMatthew G. Knepley   /* Is there a reason to snap the local bounding box to a division of the global box? */
579cafe43deSMatthew G. Knepley   /* Should we compute all overlaps of local boxes? We could do this with a rendevouz scheme partitioning the global box */
580cafe43deSMatthew G. Knepley   /* Create label */
581cafe43deSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
5821d0c6c94SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
5831d0c6c94SMatthew G. Knepley   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
584d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, "cells", &lbox->cellsSparse);CHKERRQ(ierr);
585cafe43deSMatthew G. Knepley   ierr = DMLabelCreateIndex(lbox->cellsSparse, cStart, cEnd);CHKERRQ(ierr);
586a8d69d7bSBarry Smith   /* Compute boxes which overlap each cell: https://stackoverflow.com/questions/13790208/triangle-square-intersection-test-in-2d */
587cafe43deSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
588cafe43deSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
58938353de4SMatthew G. Knepley   ierr = PetscCalloc2(16 * dim, &dboxes, 16, &boxes);CHKERRQ(ierr);
590cafe43deSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
591cafe43deSMatthew G. Knepley     const PetscReal *h       = lbox->h;
592cafe43deSMatthew G. Knepley     PetscScalar     *ccoords = NULL;
59338353de4SMatthew G. Knepley     PetscInt         csize   = 0;
594cafe43deSMatthew G. Knepley     PetscScalar      point[3];
595cafe43deSMatthew G. Knepley     PetscInt         dlim[6], d, e, i, j, k;
596cafe43deSMatthew G. Knepley 
597cafe43deSMatthew G. Knepley     /* Find boxes enclosing each vertex */
59838353de4SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, &csize, &ccoords);CHKERRQ(ierr);
59938353de4SMatthew G. Knepley     ierr = PetscGridHashGetEnclosingBox(lbox, csize/dim, ccoords, dboxes, boxes);CHKERRQ(ierr);
600722d0f5cSMatthew G. Knepley     /* Mark cells containing the vertices */
60138353de4SMatthew G. Knepley     for (e = 0; e < csize/dim; ++e) {ierr = DMLabelSetValue(lbox->cellsSparse, c, boxes[e]);CHKERRQ(ierr);}
602cafe43deSMatthew G. Knepley     /* Get grid of boxes containing these */
603cafe43deSMatthew G. Knepley     for (d = 0;   d < dim; ++d) {dlim[d*2+0] = dlim[d*2+1] = dboxes[d];}
6042291669eSMatthew G. Knepley     for (d = dim; d < 3;   ++d) {dlim[d*2+0] = dlim[d*2+1] = 0;}
605cafe43deSMatthew G. Knepley     for (e = 1; e < dim+1; ++e) {
606cafe43deSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
607cafe43deSMatthew G. Knepley         dlim[d*2+0] = PetscMin(dlim[d*2+0], dboxes[e*dim+d]);
608cafe43deSMatthew G. Knepley         dlim[d*2+1] = PetscMax(dlim[d*2+1], dboxes[e*dim+d]);
609cafe43deSMatthew G. Knepley       }
610cafe43deSMatthew G. Knepley     }
611fea14342SMatthew G. Knepley     /* Check for intersection of box with cell */
612cafe43deSMatthew 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]) {
613cafe43deSMatthew 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]) {
614cafe43deSMatthew 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]) {
615cafe43deSMatthew G. Knepley           const PetscInt box = (k*lbox->n[1] + j)*lbox->n[0] + i;
616cafe43deSMatthew G. Knepley           PetscScalar    cpoint[3];
617fea14342SMatthew G. Knepley           PetscInt       cell, edge, ii, jj, kk;
618cafe43deSMatthew G. Knepley 
619fea14342SMatthew G. Knepley           /* Check whether cell contains any vertex of these subboxes TODO vectorize this */
620cafe43deSMatthew G. Knepley           for (kk = 0, cpoint[2] = point[2]; kk < (dim > 2 ? 2 : 1); ++kk, cpoint[2] += h[2]) {
621cafe43deSMatthew G. Knepley             for (jj = 0, cpoint[1] = point[1]; jj < (dim > 1 ? 2 : 1); ++jj, cpoint[1] += h[1]) {
622cafe43deSMatthew G. Knepley               for (ii = 0, cpoint[0] = point[0]; ii < 2; ++ii, cpoint[0] += h[0]) {
623cafe43deSMatthew G. Knepley 
624cafe43deSMatthew G. Knepley                 ierr = DMPlexLocatePoint_Internal(dm, dim, cpoint, c, &cell);CHKERRQ(ierr);
625367003a6SStefano Zampini                 if (cell >= 0) { ierr = DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); ii = jj = kk = 2;}
626cafe43deSMatthew G. Knepley               }
627cafe43deSMatthew G. Knepley             }
628cafe43deSMatthew G. Knepley           }
629fea14342SMatthew G. Knepley           /* Check whether cell edge intersects any edge of these subboxes TODO vectorize this */
630fea14342SMatthew G. Knepley           for (edge = 0; edge < dim+1; ++edge) {
631fea14342SMatthew G. Knepley             PetscReal segA[6], segB[6];
632fea14342SMatthew G. Knepley 
633aab5bcd8SJed Brown             if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected dim %d > 3",dim);
634fea14342SMatthew 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]);}
635fea14342SMatthew G. Knepley             for (kk = 0; kk < (dim > 2 ? 2 : 1); ++kk) {
6369a128ed2SMatthew G. Knepley               if (dim > 2) {segB[2]     = PetscRealPart(point[2]);
6379a128ed2SMatthew G. Knepley                             segB[dim+2] = PetscRealPart(point[2]) + kk*h[2];}
638fea14342SMatthew G. Knepley               for (jj = 0; jj < (dim > 1 ? 2 : 1); ++jj) {
6399a128ed2SMatthew G. Knepley                 if (dim > 1) {segB[1]     = PetscRealPart(point[1]);
6409a128ed2SMatthew G. Knepley                               segB[dim+1] = PetscRealPart(point[1]) + jj*h[1];}
641fea14342SMatthew G. Knepley                 for (ii = 0; ii < 2; ++ii) {
642fea14342SMatthew G. Knepley                   PetscBool intersects;
643fea14342SMatthew G. Knepley 
6449a128ed2SMatthew G. Knepley                   segB[0]     = PetscRealPart(point[0]);
6459a128ed2SMatthew G. Knepley                   segB[dim+0] = PetscRealPart(point[0]) + ii*h[0];
646fea14342SMatthew G. Knepley                   ierr = DMPlexGetLineIntersection_2D_Internal(segA, segB, NULL, &intersects);CHKERRQ(ierr);
647367003a6SStefano Zampini                   if (intersects) { ierr = DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); edge = ii = jj = kk = dim+1;}
648cafe43deSMatthew G. Knepley                 }
649cafe43deSMatthew G. Knepley               }
650cafe43deSMatthew G. Knepley             }
651cafe43deSMatthew G. Knepley           }
652fea14342SMatthew G. Knepley         }
653fea14342SMatthew G. Knepley       }
654fea14342SMatthew G. Knepley     }
655fea14342SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &ccoords);CHKERRQ(ierr);
656fea14342SMatthew G. Knepley   }
657722d0f5cSMatthew G. Knepley   ierr = PetscFree2(dboxes, boxes);CHKERRQ(ierr);
658cafe43deSMatthew G. Knepley   ierr = DMLabelConvertToSection(lbox->cellsSparse, &lbox->cellSection, &lbox->cells);CHKERRQ(ierr);
659cafe43deSMatthew G. Knepley   ierr = DMLabelDestroy(&lbox->cellsSparse);CHKERRQ(ierr);
660cafe43deSMatthew G. Knepley   *localBox = lbox;
661cafe43deSMatthew G. Knepley   PetscFunctionReturn(0);
662cafe43deSMatthew G. Knepley }
663cafe43deSMatthew G. Knepley 
66462a38674SMatthew G. Knepley PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, DMPointLocationType ltype, PetscSF cellSF)
665ccd2543fSMatthew G Knepley {
666cafe43deSMatthew G. Knepley   DM_Plex        *mesh = (DM_Plex *) dm->data;
667af74b616SDave May   PetscBool       hash = mesh->useHashLocation, reuse = PETSC_FALSE;
6683a93e3b7SToby Isaac   PetscInt        bs, numPoints, p, numFound, *found = NULL;
6699cb35068SDave May   PetscInt        dim, cStart, cEnd, cMax, numCells, c, d;
670cafe43deSMatthew G. Knepley   const PetscInt *boxCells;
6713a93e3b7SToby Isaac   PetscSFNode    *cells;
672ccd2543fSMatthew G Knepley   PetscScalar    *a;
6733a93e3b7SToby Isaac   PetscMPIInt     result;
674af74b616SDave May   PetscLogDouble  t0,t1;
6759cb35068SDave May   PetscReal       gmin[3],gmax[3];
6769cb35068SDave May   PetscInt        terminating_query_type[] = { 0, 0, 0 };
677ccd2543fSMatthew G Knepley   PetscErrorCode  ierr;
678ccd2543fSMatthew G Knepley 
679ccd2543fSMatthew G Knepley   PetscFunctionBegin;
680af74b616SDave May   ierr = PetscTime(&t0);CHKERRQ(ierr);
681080342d1SMatthew 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.");
682cafe43deSMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
683cafe43deSMatthew G. Knepley   ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr);
6843a93e3b7SToby Isaac   ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)cellSF),PETSC_COMM_SELF,&result);CHKERRQ(ierr);
6853a93e3b7SToby Isaac   if (result != MPI_IDENT && result != MPI_CONGRUENT) SETERRQ(PetscObjectComm((PetscObject)cellSF),PETSC_ERR_SUP, "Trying parallel point location: only local point location supported");
686cafe43deSMatthew 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);
687ccd2543fSMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
688ccd2543fSMatthew G Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
689ccd2543fSMatthew G Knepley   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
690ccd2543fSMatthew G Knepley   ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr);
691ccd2543fSMatthew G Knepley   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
692ccd2543fSMatthew G Knepley   numPoints /= bs;
693af74b616SDave May   {
694af74b616SDave May     const PetscSFNode *sf_cells;
695af74b616SDave May 
696af74b616SDave May     ierr = PetscSFGetGraph(cellSF,NULL,NULL,NULL,&sf_cells);CHKERRQ(ierr);
697af74b616SDave May     if (sf_cells) {
698af74b616SDave May       ierr = PetscInfo(dm,"[DMLocatePoints_Plex] Re-using existing StarForest node list\n");CHKERRQ(ierr);
699af74b616SDave May       cells = (PetscSFNode*)sf_cells;
700af74b616SDave May       reuse = PETSC_TRUE;
701af74b616SDave May     } else {
702af74b616SDave May       ierr = PetscInfo(dm,"[DMLocatePoints_Plex] Creating and initializing new StarForest node list\n");CHKERRQ(ierr);
703785e854fSJed Brown       ierr = PetscMalloc1(numPoints, &cells);CHKERRQ(ierr);
704af74b616SDave May       /* initialize cells if created */
705af74b616SDave May       for (p=0; p<numPoints; p++) {
706af74b616SDave May         cells[p].rank  = 0;
707af74b616SDave May         cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
708af74b616SDave May       }
709af74b616SDave May     }
710af74b616SDave May   }
7119cb35068SDave May   /* define domain bounding box */
7129cb35068SDave May   {
7139cb35068SDave May     Vec coorglobal;
7149cb35068SDave May 
7159cb35068SDave May     ierr = DMGetCoordinates(dm,&coorglobal);CHKERRQ(ierr);
7169cb35068SDave May     ierr = VecStrideMaxAll(coorglobal,NULL,gmax);CHKERRQ(ierr);
7179cb35068SDave May     ierr = VecStrideMinAll(coorglobal,NULL,gmin);CHKERRQ(ierr);
7189cb35068SDave May   }
719953fc75cSMatthew G. Knepley   if (hash) {
720ac6ec2abSMatthew G. Knepley     if (!mesh->lbox) {ierr = PetscInfo(dm, "Initializing grid hashing");CHKERRQ(ierr);ierr = DMPlexComputeGridHash_Internal(dm, &mesh->lbox);CHKERRQ(ierr);}
721cafe43deSMatthew G. Knepley     /* Designate the local box for each point */
722cafe43deSMatthew G. Knepley     /* Send points to correct process */
723cafe43deSMatthew G. Knepley     /* Search cells that lie in each subbox */
724cafe43deSMatthew G. Knepley     /*   Should we bin points before doing search? */
725cafe43deSMatthew G. Knepley     ierr = ISGetIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr);
726953fc75cSMatthew G. Knepley   }
7273a93e3b7SToby Isaac   for (p = 0, numFound = 0; p < numPoints; ++p) {
728ccd2543fSMatthew G Knepley     const PetscScalar *point = &a[p*bs];
729e56f9228SJed Brown     PetscInt           dbin[3] = {-1,-1,-1}, bin, cell = -1, cellOffset;
7309cb35068SDave May     PetscBool          point_outside_domain = PETSC_FALSE;
731ccd2543fSMatthew G Knepley 
7329cb35068SDave May     /* check bounding box of domain */
7339cb35068SDave May     for (d=0; d<dim; d++) {
734a5f152d1SDave May       if (PetscRealPart(point[d]) < gmin[d]) { point_outside_domain = PETSC_TRUE; break; }
735a5f152d1SDave May       if (PetscRealPart(point[d]) > gmax[d]) { point_outside_domain = PETSC_TRUE; break; }
7369cb35068SDave May     }
7379cb35068SDave May     if (point_outside_domain) {
738e9b685f5SMatthew G. Knepley       cells[p].rank = 0;
739e9b685f5SMatthew G. Knepley       cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
7409cb35068SDave May       terminating_query_type[0]++;
7419cb35068SDave May       continue;
7429cb35068SDave May     }
743ccd2543fSMatthew G Knepley 
744af74b616SDave May     /* check initial values in cells[].index - abort early if found */
745af74b616SDave May     if (cells[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
746af74b616SDave May       c = cells[p].index;
7473a93e3b7SToby Isaac       cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
748af74b616SDave May       ierr = DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);CHKERRQ(ierr);
749af74b616SDave May       if (cell >= 0) {
750af74b616SDave May         cells[p].rank = 0;
751af74b616SDave May         cells[p].index = cell;
752af74b616SDave May         numFound++;
753af74b616SDave May       }
754af74b616SDave May     }
7559cb35068SDave May     if (cells[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
7569cb35068SDave May       terminating_query_type[1]++;
7579cb35068SDave May       continue;
7589cb35068SDave May     }
759af74b616SDave May 
760953fc75cSMatthew G. Knepley     if (hash) {
761af74b616SDave May       PetscBool found_box;
762af74b616SDave May 
763af74b616SDave May       /* allow for case that point is outside box - abort early */
764af74b616SDave May       ierr = PetscGridHashGetEnclosingBoxQuery(mesh->lbox, 1, point, dbin, &bin,&found_box);CHKERRQ(ierr);
765af74b616SDave May       if (found_box) {
766cafe43deSMatthew G. Knepley         /* TODO Lay an interface over this so we can switch between Section (dense) and Label (sparse) */
767cafe43deSMatthew G. Knepley         ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr);
768cafe43deSMatthew G. Knepley         ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr);
769cafe43deSMatthew G. Knepley         for (c = cellOffset; c < cellOffset + numCells; ++c) {
770cafe43deSMatthew G. Knepley           ierr = DMPlexLocatePoint_Internal(dm, dim, point, boxCells[c], &cell);CHKERRQ(ierr);
7713a93e3b7SToby Isaac           if (cell >= 0) {
7723a93e3b7SToby Isaac             cells[p].rank = 0;
7733a93e3b7SToby Isaac             cells[p].index = cell;
7743a93e3b7SToby Isaac             numFound++;
7759cb35068SDave May             terminating_query_type[2]++;
7763a93e3b7SToby Isaac             break;
777ccd2543fSMatthew G Knepley           }
7783a93e3b7SToby Isaac         }
779af74b616SDave May       }
780953fc75cSMatthew G. Knepley     } else {
781953fc75cSMatthew G. Knepley       for (c = cStart; c < cEnd; ++c) {
782953fc75cSMatthew G. Knepley         ierr = DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);CHKERRQ(ierr);
7833a93e3b7SToby Isaac         if (cell >= 0) {
7843a93e3b7SToby Isaac           cells[p].rank = 0;
7853a93e3b7SToby Isaac           cells[p].index = cell;
7863a93e3b7SToby Isaac           numFound++;
7879cb35068SDave May           terminating_query_type[2]++;
7883a93e3b7SToby Isaac           break;
789953fc75cSMatthew G. Knepley         }
790953fc75cSMatthew G. Knepley       }
7913a93e3b7SToby Isaac     }
792ccd2543fSMatthew G Knepley   }
793953fc75cSMatthew G. Knepley   if (hash) {ierr = ISRestoreIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr);}
79462a38674SMatthew G. Knepley   if (ltype == DM_POINTLOCATION_NEAREST && hash && numFound < numPoints) {
79562a38674SMatthew G. Knepley     for (p = 0; p < numPoints; p++) {
79662a38674SMatthew G. Knepley       const PetscScalar *point = &a[p*bs];
79762a38674SMatthew G. Knepley       PetscReal          cpoint[3], diff[3], dist, distMax = PETSC_MAX_REAL;
798e56f9228SJed Brown       PetscInt           dbin[3] = {-1,-1,-1}, bin, cellOffset, d;
79962a38674SMatthew G. Knepley 
800e9b685f5SMatthew G. Knepley       if (cells[p].index < 0) {
80162a38674SMatthew G. Knepley         ++numFound;
80262a38674SMatthew G. Knepley         ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr);
80362a38674SMatthew G. Knepley         ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr);
80462a38674SMatthew G. Knepley         ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr);
80562a38674SMatthew G. Knepley         for (c = cellOffset; c < cellOffset + numCells; ++c) {
80662a38674SMatthew G. Knepley           ierr = DMPlexClosestPoint_Internal(dm, dim, point, boxCells[c], cpoint);CHKERRQ(ierr);
807b716b415SMatthew G. Knepley           for (d = 0; d < dim; ++d) diff[d] = cpoint[d] - PetscRealPart(point[d]);
80862a38674SMatthew G. Knepley           dist = DMPlex_NormD_Internal(dim, diff);
80962a38674SMatthew G. Knepley           if (dist < distMax) {
81062a38674SMatthew G. Knepley             for (d = 0; d < dim; ++d) a[p*bs+d] = cpoint[d];
81162a38674SMatthew G. Knepley             cells[p].rank  = 0;
81262a38674SMatthew G. Knepley             cells[p].index = boxCells[c];
81362a38674SMatthew G. Knepley             distMax = dist;
81462a38674SMatthew G. Knepley             break;
81562a38674SMatthew G. Knepley           }
81662a38674SMatthew G. Knepley         }
81762a38674SMatthew G. Knepley       }
81862a38674SMatthew G. Knepley     }
81962a38674SMatthew G. Knepley   }
82062a38674SMatthew G. Knepley   /* This code is only be relevant when interfaced to parallel point location */
821cafe43deSMatthew G. Knepley   /* Check for highest numbered proc that claims a point (do we care?) */
8222d1fa6caSMatthew G. Knepley   if (ltype == DM_POINTLOCATION_REMOVE && numFound < numPoints) {
8233a93e3b7SToby Isaac     ierr = PetscMalloc1(numFound,&found);CHKERRQ(ierr);
8243a93e3b7SToby Isaac     for (p = 0, numFound = 0; p < numPoints; p++) {
8253a93e3b7SToby Isaac       if (cells[p].rank >= 0 && cells[p].index >= 0) {
8263a93e3b7SToby Isaac         if (numFound < p) {
8273a93e3b7SToby Isaac           cells[numFound] = cells[p];
8283a93e3b7SToby Isaac         }
8293a93e3b7SToby Isaac         found[numFound++] = p;
8303a93e3b7SToby Isaac       }
8313a93e3b7SToby Isaac     }
8323a93e3b7SToby Isaac   }
83362a38674SMatthew G. Knepley   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
834af74b616SDave May   if (!reuse) {
8353a93e3b7SToby Isaac     ierr = PetscSFSetGraph(cellSF, cEnd - cStart, numFound, found, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER);CHKERRQ(ierr);
836af74b616SDave May   }
837af74b616SDave May   ierr = PetscTime(&t1);CHKERRQ(ierr);
8389cb35068SDave May   if (hash) {
839d6e0b8ebSSatish Balay     ierr = PetscInfo3(dm,"[DMLocatePoints_Plex] terminating_query_type : %D [outside domain] : %D [inside intial cell] : %D [hash]\n",terminating_query_type[0],terminating_query_type[1],terminating_query_type[2]);CHKERRQ(ierr);
8409cb35068SDave May   } else {
841d6e0b8ebSSatish Balay     ierr = PetscInfo3(dm,"[DMLocatePoints_Plex] terminating_query_type : %D [outside domain] : %D [inside intial cell] : %D [brute-force]\n",terminating_query_type[0],terminating_query_type[1],terminating_query_type[2]);CHKERRQ(ierr);
8429cb35068SDave May   }
843af74b616SDave May   ierr = PetscInfo3(dm,"[DMLocatePoints_Plex] npoints %D : time(rank0) %1.2e (sec): points/sec %1.4e\n",numPoints,t1-t0,(double)((double)numPoints/(t1-t0)));CHKERRQ(ierr);
844ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
845ccd2543fSMatthew G Knepley }
846ccd2543fSMatthew G Knepley 
847741bfc07SMatthew G. Knepley /*@C
848741bfc07SMatthew G. Knepley   DMPlexComputeProjection2Dto1D - Rewrite coordinates to be the 1D projection of the 2D coordinates
849741bfc07SMatthew G. Knepley 
850741bfc07SMatthew G. Knepley   Not collective
851741bfc07SMatthew G. Knepley 
852741bfc07SMatthew G. Knepley   Input Parameter:
853741bfc07SMatthew G. Knepley . coords - The coordinates of a segment
854741bfc07SMatthew G. Knepley 
855741bfc07SMatthew G. Knepley   Output Parameters:
856741bfc07SMatthew G. Knepley + coords - The new y-coordinate, and 0 for x
857741bfc07SMatthew G. Knepley - R - The rotation which accomplishes the projection
858741bfc07SMatthew G. Knepley 
859741bfc07SMatthew G. Knepley   Level: developer
860741bfc07SMatthew G. Knepley 
861741bfc07SMatthew G. Knepley .seealso: DMPlexComputeProjection3Dto1D(), DMPlexComputeProjection3Dto2D()
862741bfc07SMatthew G. Knepley @*/
863741bfc07SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar coords[], PetscReal R[])
86417fe8556SMatthew G. Knepley {
86517fe8556SMatthew G. Knepley   const PetscReal x = PetscRealPart(coords[2] - coords[0]);
86617fe8556SMatthew G. Knepley   const PetscReal y = PetscRealPart(coords[3] - coords[1]);
8678b49ba18SBarry Smith   const PetscReal r = PetscSqrtReal(x*x + y*y), c = x/r, s = y/r;
86817fe8556SMatthew G. Knepley 
86917fe8556SMatthew G. Knepley   PetscFunctionBegin;
8701c99cf0cSGeoffrey Irving   R[0] = c; R[1] = -s;
8711c99cf0cSGeoffrey Irving   R[2] = s; R[3] =  c;
87217fe8556SMatthew G. Knepley   coords[0] = 0.0;
8737f07f362SMatthew G. Knepley   coords[1] = r;
87417fe8556SMatthew G. Knepley   PetscFunctionReturn(0);
87517fe8556SMatthew G. Knepley }
87617fe8556SMatthew G. Knepley 
877741bfc07SMatthew G. Knepley /*@C
878741bfc07SMatthew G. Knepley   DMPlexComputeProjection3Dto1D - Rewrite coordinates to be the 1D projection of the 3D coordinates
87928dbe442SToby Isaac 
880741bfc07SMatthew G. Knepley   Not collective
88128dbe442SToby Isaac 
882741bfc07SMatthew G. Knepley   Input Parameter:
883741bfc07SMatthew G. Knepley . coords - The coordinates of a segment
884741bfc07SMatthew G. Knepley 
885741bfc07SMatthew G. Knepley   Output Parameters:
886741bfc07SMatthew G. Knepley + coords - The new y-coordinate, and 0 for x and z
887741bfc07SMatthew G. Knepley - R - The rotation which accomplishes the projection
888741bfc07SMatthew G. Knepley 
889741bfc07SMatthew 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
890741bfc07SMatthew G. Knepley 
891741bfc07SMatthew G. Knepley   Level: developer
892741bfc07SMatthew G. Knepley 
893741bfc07SMatthew G. Knepley .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto2D()
894741bfc07SMatthew G. Knepley @*/
895741bfc07SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar coords[], PetscReal R[])
89628dbe442SToby Isaac {
89728dbe442SToby Isaac   PetscReal      x    = PetscRealPart(coords[3] - coords[0]);
89828dbe442SToby Isaac   PetscReal      y    = PetscRealPart(coords[4] - coords[1]);
89928dbe442SToby Isaac   PetscReal      z    = PetscRealPart(coords[5] - coords[2]);
90028dbe442SToby Isaac   PetscReal      r    = PetscSqrtReal(x*x + y*y + z*z);
90128dbe442SToby Isaac   PetscReal      rinv = 1. / r;
90228dbe442SToby Isaac   PetscFunctionBegin;
90328dbe442SToby Isaac 
90428dbe442SToby Isaac   x *= rinv; y *= rinv; z *= rinv;
90528dbe442SToby Isaac   if (x > 0.) {
90628dbe442SToby Isaac     PetscReal inv1pX   = 1./ (1. + x);
90728dbe442SToby Isaac 
90828dbe442SToby Isaac     R[0] = x; R[1] = -y;              R[2] = -z;
90928dbe442SToby Isaac     R[3] = y; R[4] = 1. - y*y*inv1pX; R[5] =     -y*z*inv1pX;
91028dbe442SToby Isaac     R[6] = z; R[7] =     -y*z*inv1pX; R[8] = 1. - z*z*inv1pX;
91128dbe442SToby Isaac   }
91228dbe442SToby Isaac   else {
91328dbe442SToby Isaac     PetscReal inv1mX   = 1./ (1. - x);
91428dbe442SToby Isaac 
91528dbe442SToby Isaac     R[0] = x; R[1] = z;               R[2] = y;
91628dbe442SToby Isaac     R[3] = y; R[4] =     -y*z*inv1mX; R[5] = 1. - y*y*inv1mX;
91728dbe442SToby Isaac     R[6] = z; R[7] = 1. - z*z*inv1mX; R[8] =     -y*z*inv1mX;
91828dbe442SToby Isaac   }
91928dbe442SToby Isaac   coords[0] = 0.0;
92028dbe442SToby Isaac   coords[1] = r;
92128dbe442SToby Isaac   PetscFunctionReturn(0);
92228dbe442SToby Isaac }
92328dbe442SToby Isaac 
924741bfc07SMatthew G. Knepley /*@
925741bfc07SMatthew G. Knepley   DMPlexComputeProjection3Dto2D - Rewrite coordinates to be the 2D projection of the 3D coordinates
926741bfc07SMatthew G. Knepley 
927741bfc07SMatthew G. Knepley   Not collective
928741bfc07SMatthew G. Knepley 
929741bfc07SMatthew G. Knepley   Input Parameter:
930741bfc07SMatthew G. Knepley . coords - The coordinates of a segment
931741bfc07SMatthew G. Knepley 
932741bfc07SMatthew G. Knepley   Output Parameters:
933741bfc07SMatthew G. Knepley + coords - The new y- and z-coordinates, and 0 for x
934741bfc07SMatthew G. Knepley - R - The rotation which accomplishes the projection
935741bfc07SMatthew G. Knepley 
936741bfc07SMatthew G. Knepley   Level: developer
937741bfc07SMatthew G. Knepley 
938741bfc07SMatthew G. Knepley .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto1D()
939741bfc07SMatthew G. Knepley @*/
940741bfc07SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt coordSize, PetscScalar coords[], PetscReal R[])
941ccd2543fSMatthew G Knepley {
9421ee9d5ecSMatthew G. Knepley   PetscReal      x1[3],  x2[3], n[3], norm;
94399dec3a6SMatthew G. Knepley   PetscReal      x1p[3], x2p[3], xnp[3];
9444a217a95SMatthew G. Knepley   PetscReal      sqrtz, alpha;
945ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
94699dec3a6SMatthew G. Knepley   PetscInt       d, e, p;
947ccd2543fSMatthew G Knepley 
948ccd2543fSMatthew G Knepley   PetscFunctionBegin;
949ccd2543fSMatthew G Knepley   /* 0) Calculate normal vector */
950ccd2543fSMatthew G Knepley   for (d = 0; d < dim; ++d) {
9511ee9d5ecSMatthew G. Knepley     x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]);
9521ee9d5ecSMatthew G. Knepley     x2[d] = PetscRealPart(coords[2*dim+d] - coords[0*dim+d]);
953ccd2543fSMatthew G Knepley   }
954ccd2543fSMatthew G Knepley   n[0] = x1[1]*x2[2] - x1[2]*x2[1];
955ccd2543fSMatthew G Knepley   n[1] = x1[2]*x2[0] - x1[0]*x2[2];
956ccd2543fSMatthew G Knepley   n[2] = x1[0]*x2[1] - x1[1]*x2[0];
9578b49ba18SBarry Smith   norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]);
958ccd2543fSMatthew G Knepley   n[0] /= norm;
959ccd2543fSMatthew G Knepley   n[1] /= norm;
960ccd2543fSMatthew G Knepley   n[2] /= norm;
961ccd2543fSMatthew G Knepley   /* 1) Take the normal vector and rotate until it is \hat z
962ccd2543fSMatthew G Knepley 
963ccd2543fSMatthew G Knepley     Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then
964ccd2543fSMatthew G Knepley 
965ccd2543fSMatthew G Knepley     R = /  alpha nx nz  alpha ny nz -1/alpha \
966ccd2543fSMatthew G Knepley         | -alpha ny     alpha nx        0    |
967ccd2543fSMatthew G Knepley         \     nx            ny         nz    /
968ccd2543fSMatthew G Knepley 
969ccd2543fSMatthew G Knepley     will rotate the normal vector to \hat z
970ccd2543fSMatthew G Knepley   */
9718b49ba18SBarry Smith   sqrtz = PetscSqrtReal(1.0 - n[2]*n[2]);
97273868372SMatthew G. Knepley   /* Check for n = z */
97373868372SMatthew G. Knepley   if (sqrtz < 1.0e-10) {
9747df32b8bSSanderA     const PetscInt s = PetscSign(n[2]);
9757df32b8bSSanderA     /* If nz < 0, rotate 180 degrees around x-axis */
97699dec3a6SMatthew G. Knepley     for (p = 3; p < coordSize/3; ++p) {
97799dec3a6SMatthew G. Knepley       coords[p*2+0] = PetscRealPart(coords[p*dim+0] - coords[0*dim+0]);
9787df32b8bSSanderA       coords[p*2+1] = (PetscRealPart(coords[p*dim+1] - coords[0*dim+1])) * s;
97973868372SMatthew G. Knepley     }
98099dec3a6SMatthew G. Knepley     coords[0] = 0.0;
98199dec3a6SMatthew G. Knepley     coords[1] = 0.0;
9827df32b8bSSanderA     coords[2] = x1[0];
9837df32b8bSSanderA     coords[3] = x1[1] * s;
9847df32b8bSSanderA     coords[4] = x2[0];
9857df32b8bSSanderA     coords[5] = x2[1] * s;
9867df32b8bSSanderA     R[0] = 1.0;     R[1] = 0.0;     R[2] = 0.0;
9877df32b8bSSanderA     R[3] = 0.0;     R[4] = 1.0 * s; R[5] = 0.0;
9887df32b8bSSanderA     R[6] = 0.0;     R[7] = 0.0;     R[8] = 1.0 * s;
98973868372SMatthew G. Knepley     PetscFunctionReturn(0);
99073868372SMatthew G. Knepley   }
991da18b5e6SMatthew G Knepley   alpha = 1.0/sqrtz;
992ccd2543fSMatthew G Knepley   R[0] =  alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz;
993ccd2543fSMatthew G Knepley   R[3] = -alpha*n[1];      R[4] = alpha*n[0];      R[5] = 0.0;
994ccd2543fSMatthew G Knepley   R[6] =  n[0];            R[7] = n[1];            R[8] = n[2];
995ccd2543fSMatthew G Knepley   for (d = 0; d < dim; ++d) {
996ccd2543fSMatthew G Knepley     x1p[d] = 0.0;
997ccd2543fSMatthew G Knepley     x2p[d] = 0.0;
998ccd2543fSMatthew G Knepley     for (e = 0; e < dim; ++e) {
999ccd2543fSMatthew G Knepley       x1p[d] += R[d*dim+e]*x1[e];
1000ccd2543fSMatthew G Knepley       x2p[d] += R[d*dim+e]*x2[e];
1001ccd2543fSMatthew G Knepley     }
1002ccd2543fSMatthew G Knepley   }
100348120919SToby Isaac   if (PetscAbsReal(x1p[2]) > 10. * PETSC_SMALL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
100448120919SToby Isaac   if (PetscAbsReal(x2p[2]) > 10. * PETSC_SMALL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
1005ccd2543fSMatthew G Knepley   /* 2) Project to (x, y) */
100699dec3a6SMatthew G. Knepley   for (p = 3; p < coordSize/3; ++p) {
100799dec3a6SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
100899dec3a6SMatthew G. Knepley       xnp[d] = 0.0;
100999dec3a6SMatthew G. Knepley       for (e = 0; e < dim; ++e) {
101099dec3a6SMatthew G. Knepley         xnp[d] += R[d*dim+e]*PetscRealPart(coords[p*dim+e] - coords[0*dim+e]);
101199dec3a6SMatthew G. Knepley       }
101299dec3a6SMatthew G. Knepley       if (d < dim-1) coords[p*2+d] = xnp[d];
101399dec3a6SMatthew G. Knepley     }
101499dec3a6SMatthew G. Knepley   }
1015ccd2543fSMatthew G Knepley   coords[0] = 0.0;
1016ccd2543fSMatthew G Knepley   coords[1] = 0.0;
1017ccd2543fSMatthew G Knepley   coords[2] = x1p[0];
1018ccd2543fSMatthew G Knepley   coords[3] = x1p[1];
1019ccd2543fSMatthew G Knepley   coords[4] = x2p[0];
1020ccd2543fSMatthew G Knepley   coords[5] = x2p[1];
10217f07f362SMatthew G. Knepley   /* Output R^T which rotates \hat z to the input normal */
10227f07f362SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
10237f07f362SMatthew G. Knepley     for (e = d+1; e < dim; ++e) {
10247f07f362SMatthew G. Knepley       PetscReal tmp;
10257f07f362SMatthew G. Knepley 
10267f07f362SMatthew G. Knepley       tmp        = R[d*dim+e];
10277f07f362SMatthew G. Knepley       R[d*dim+e] = R[e*dim+d];
10287f07f362SMatthew G. Knepley       R[e*dim+d] = tmp;
10297f07f362SMatthew G. Knepley     }
10307f07f362SMatthew G. Knepley   }
1031ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1032ccd2543fSMatthew G Knepley }
1033ccd2543fSMatthew G Knepley 
10346322fe33SJed Brown PETSC_UNUSED
1035834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[])
1036834e62ceSMatthew G. Knepley {
1037834e62ceSMatthew G. Knepley   /* Signed volume is 1/2 the determinant
1038834e62ceSMatthew G. Knepley 
1039834e62ceSMatthew G. Knepley    |  1  1  1 |
1040834e62ceSMatthew G. Knepley    | x0 x1 x2 |
1041834e62ceSMatthew G. Knepley    | y0 y1 y2 |
1042834e62ceSMatthew G. Knepley 
1043834e62ceSMatthew G. Knepley      but if x0,y0 is the origin, we have
1044834e62ceSMatthew G. Knepley 
1045834e62ceSMatthew G. Knepley    | x1 x2 |
1046834e62ceSMatthew G. Knepley    | y1 y2 |
1047834e62ceSMatthew G. Knepley   */
1048834e62ceSMatthew G. Knepley   const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1];
1049834e62ceSMatthew G. Knepley   const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1];
1050834e62ceSMatthew G. Knepley   PetscReal       M[4], detM;
1051834e62ceSMatthew G. Knepley   M[0] = x1; M[1] = x2;
105286623015SMatthew G. Knepley   M[2] = y1; M[3] = y2;
1053923591dfSMatthew G. Knepley   DMPlex_Det2D_Internal(&detM, M);
1054834e62ceSMatthew G. Knepley   *vol = 0.5*detM;
10553bc0b13bSBarry Smith   (void)PetscLogFlops(5.0);
1056834e62ceSMatthew G. Knepley }
1057834e62ceSMatthew G. Knepley 
1058834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[])
1059834e62ceSMatthew G. Knepley {
1060923591dfSMatthew G. Knepley   DMPlex_Det2D_Internal(vol, coords);
1061834e62ceSMatthew G. Knepley   *vol *= 0.5;
1062834e62ceSMatthew G. Knepley }
1063834e62ceSMatthew G. Knepley 
10646322fe33SJed Brown PETSC_UNUSED
1065834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[])
1066834e62ceSMatthew G. Knepley {
1067834e62ceSMatthew G. Knepley   /* Signed volume is 1/6th of the determinant
1068834e62ceSMatthew G. Knepley 
1069834e62ceSMatthew G. Knepley    |  1  1  1  1 |
1070834e62ceSMatthew G. Knepley    | x0 x1 x2 x3 |
1071834e62ceSMatthew G. Knepley    | y0 y1 y2 y3 |
1072834e62ceSMatthew G. Knepley    | z0 z1 z2 z3 |
1073834e62ceSMatthew G. Knepley 
1074834e62ceSMatthew G. Knepley      but if x0,y0,z0 is the origin, we have
1075834e62ceSMatthew G. Knepley 
1076834e62ceSMatthew G. Knepley    | x1 x2 x3 |
1077834e62ceSMatthew G. Knepley    | y1 y2 y3 |
1078834e62ceSMatthew G. Knepley    | z1 z2 z3 |
1079834e62ceSMatthew G. Knepley   */
1080834e62ceSMatthew G. Knepley   const PetscReal x1 = coords[3] - coords[0], y1 = coords[4]  - coords[1], z1 = coords[5]  - coords[2];
1081834e62ceSMatthew G. Knepley   const PetscReal x2 = coords[6] - coords[0], y2 = coords[7]  - coords[1], z2 = coords[8]  - coords[2];
1082834e62ceSMatthew G. Knepley   const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2];
10830a3da2c2SToby Isaac   const PetscReal onesixth = ((PetscReal)1./(PetscReal)6.);
1084834e62ceSMatthew G. Knepley   PetscReal       M[9], detM;
1085834e62ceSMatthew G. Knepley   M[0] = x1; M[1] = x2; M[2] = x3;
1086834e62ceSMatthew G. Knepley   M[3] = y1; M[4] = y2; M[5] = y3;
1087834e62ceSMatthew G. Knepley   M[6] = z1; M[7] = z2; M[8] = z3;
1088923591dfSMatthew G. Knepley   DMPlex_Det3D_Internal(&detM, M);
10890a3da2c2SToby Isaac   *vol = -onesixth*detM;
10903bc0b13bSBarry Smith   (void)PetscLogFlops(10.0);
1091834e62ceSMatthew G. Knepley }
1092834e62ceSMatthew G. Knepley 
10930ec8681fSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[])
10940ec8681fSMatthew G. Knepley {
10950a3da2c2SToby Isaac   const PetscReal onesixth = ((PetscReal)1./(PetscReal)6.);
1096923591dfSMatthew G. Knepley   DMPlex_Det3D_Internal(vol, coords);
10970a3da2c2SToby Isaac   *vol *= -onesixth;
10980ec8681fSMatthew G. Knepley }
10990ec8681fSMatthew G. Knepley 
1100cb92db44SToby Isaac static PetscErrorCode DMPlexComputePointGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1101cb92db44SToby Isaac {
1102cb92db44SToby Isaac   PetscSection   coordSection;
1103cb92db44SToby Isaac   Vec            coordinates;
1104cb92db44SToby Isaac   const PetscScalar *coords;
1105cb92db44SToby Isaac   PetscInt       dim, d, off;
1106cb92db44SToby Isaac   PetscErrorCode ierr;
1107cb92db44SToby Isaac 
1108cb92db44SToby Isaac   PetscFunctionBegin;
1109cb92db44SToby Isaac   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1110cb92db44SToby Isaac   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1111cb92db44SToby Isaac   ierr = PetscSectionGetDof(coordSection,e,&dim);CHKERRQ(ierr);
1112cb92db44SToby Isaac   if (!dim) PetscFunctionReturn(0);
1113cb92db44SToby Isaac   ierr = PetscSectionGetOffset(coordSection,e,&off);CHKERRQ(ierr);
1114cb92db44SToby Isaac   ierr = VecGetArrayRead(coordinates,&coords);CHKERRQ(ierr);
1115cb92db44SToby Isaac   if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[off + d]);}
1116cb92db44SToby Isaac   ierr = VecRestoreArrayRead(coordinates,&coords);CHKERRQ(ierr);
1117cb92db44SToby Isaac   *detJ = 1.;
1118cb92db44SToby Isaac   if (J) {
1119cb92db44SToby Isaac     for (d = 0; d < dim * dim; d++) J[d] = 0.;
1120cb92db44SToby Isaac     for (d = 0; d < dim; d++) J[d * dim + d] = 1.;
1121cb92db44SToby Isaac     if (invJ) {
1122cb92db44SToby Isaac       for (d = 0; d < dim * dim; d++) invJ[d] = 0.;
1123cb92db44SToby Isaac       for (d = 0; d < dim; d++) invJ[d * dim + d] = 1.;
1124cb92db44SToby Isaac     }
1125cb92db44SToby Isaac   }
1126cb92db44SToby Isaac   PetscFunctionReturn(0);
1127cb92db44SToby Isaac }
1128cb92db44SToby Isaac 
112917fe8556SMatthew G. Knepley static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
113017fe8556SMatthew G. Knepley {
113117fe8556SMatthew G. Knepley   PetscSection   coordSection;
113217fe8556SMatthew G. Knepley   Vec            coordinates;
1133a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
11348bf5c034SToby Isaac   PetscInt       numCoords, d, pStart, pEnd, numSelfCoords = 0;
113517fe8556SMatthew G. Knepley   PetscErrorCode ierr;
113617fe8556SMatthew G. Knepley 
113717fe8556SMatthew G. Knepley   PetscFunctionBegin;
113817fe8556SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
113969d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
11408bf5c034SToby Isaac   ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr);
11418bf5c034SToby Isaac   if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);}
114217fe8556SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
11438bf5c034SToby Isaac   numCoords = numSelfCoords ? numSelfCoords : numCoords;
1144adac9986SMatthew G. Knepley   if (invJ && !J) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL");
11457f07f362SMatthew G. Knepley   *detJ = 0.0;
114628dbe442SToby Isaac   if (numCoords == 6) {
114728dbe442SToby Isaac     const PetscInt dim = 3;
114828dbe442SToby Isaac     PetscReal      R[9], J0;
114928dbe442SToby Isaac 
115028dbe442SToby Isaac     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1151741bfc07SMatthew G. Knepley     ierr = DMPlexComputeProjection3Dto1D(coords, R);CHKERRQ(ierr);
115228dbe442SToby Isaac     if (J)    {
115328dbe442SToby Isaac       J0   = 0.5*PetscRealPart(coords[1]);
115428dbe442SToby Isaac       J[0] = R[0]*J0; J[1] = R[1]; J[2] = R[2];
115528dbe442SToby Isaac       J[3] = R[3]*J0; J[4] = R[4]; J[5] = R[5];
115628dbe442SToby Isaac       J[6] = R[6]*J0; J[7] = R[7]; J[8] = R[8];
115728dbe442SToby Isaac       DMPlex_Det3D_Internal(detJ, J);
115828dbe442SToby Isaac       if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1159adac9986SMatthew G. Knepley     }
116028dbe442SToby Isaac   } else if (numCoords == 4) {
11617f07f362SMatthew G. Knepley     const PetscInt dim = 2;
11627f07f362SMatthew G. Knepley     PetscReal      R[4], J0;
11637f07f362SMatthew G. Knepley 
11647f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1165741bfc07SMatthew G. Knepley     ierr = DMPlexComputeProjection2Dto1D(coords, R);CHKERRQ(ierr);
116617fe8556SMatthew G. Knepley     if (J)    {
11677f07f362SMatthew G. Knepley       J0   = 0.5*PetscRealPart(coords[1]);
11687f07f362SMatthew G. Knepley       J[0] = R[0]*J0; J[1] = R[1];
11697f07f362SMatthew G. Knepley       J[2] = R[2]*J0; J[3] = R[3];
1170923591dfSMatthew G. Knepley       DMPlex_Det2D_Internal(detJ, J);
1171923591dfSMatthew G. Knepley       if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1172adac9986SMatthew G. Knepley     }
11737f07f362SMatthew G. Knepley   } else if (numCoords == 2) {
11747f07f362SMatthew G. Knepley     const PetscInt dim = 1;
11757f07f362SMatthew G. Knepley 
11767f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
11777f07f362SMatthew G. Knepley     if (J)    {
11787f07f362SMatthew G. Knepley       J[0]  = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0]));
117917fe8556SMatthew G. Knepley       *detJ = J[0];
11803bc0b13bSBarry Smith       ierr = PetscLogFlops(2.0);CHKERRQ(ierr);
11813bc0b13bSBarry Smith       if (invJ) {invJ[0] = 1.0/J[0]; ierr = PetscLogFlops(1.0);CHKERRQ(ierr);}
1182adac9986SMatthew G. Knepley     }
1183796f034aSJed Brown   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %D != 2", numCoords);
118417fe8556SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
118517fe8556SMatthew G. Knepley   PetscFunctionReturn(0);
118617fe8556SMatthew G. Knepley }
118717fe8556SMatthew G. Knepley 
1188ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1189ccd2543fSMatthew G Knepley {
1190ccd2543fSMatthew G Knepley   PetscSection   coordSection;
1191ccd2543fSMatthew G Knepley   Vec            coordinates;
1192a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
119303d7ed2eSStefano Zampini   PetscInt       numCoords, numSelfCoords = 0, d, f, g, pStart, pEnd;
1194ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
1195ccd2543fSMatthew G Knepley 
1196ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1197ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
119869d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
119903d7ed2eSStefano Zampini   ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr);
120003d7ed2eSStefano Zampini   if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);}
1201ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
120203d7ed2eSStefano Zampini   numCoords = numSelfCoords ? numSelfCoords : numCoords;
12037f07f362SMatthew G. Knepley   *detJ = 0.0;
1204ccd2543fSMatthew G Knepley   if (numCoords == 9) {
12057f07f362SMatthew G. Knepley     const PetscInt dim = 3;
12067f07f362SMatthew 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};
12077f07f362SMatthew G. Knepley 
12087f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1209741bfc07SMatthew G. Knepley     ierr = DMPlexComputeProjection3Dto2D(numCoords, coords, R);CHKERRQ(ierr);
12107f07f362SMatthew G. Knepley     if (J)    {
1211b7ad821dSMatthew G. Knepley       const PetscInt pdim = 2;
1212b7ad821dSMatthew G. Knepley 
1213b7ad821dSMatthew G. Knepley       for (d = 0; d < pdim; d++) {
1214b7ad821dSMatthew G. Knepley         for (f = 0; f < pdim; f++) {
1215b7ad821dSMatthew G. Knepley           J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
1216ccd2543fSMatthew G Knepley         }
12177f07f362SMatthew G. Knepley       }
12183bc0b13bSBarry Smith       ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1219923591dfSMatthew G. Knepley       DMPlex_Det3D_Internal(detJ, J0);
12207f07f362SMatthew G. Knepley       for (d = 0; d < dim; d++) {
12217f07f362SMatthew G. Knepley         for (f = 0; f < dim; f++) {
12227f07f362SMatthew G. Knepley           J[d*dim+f] = 0.0;
12237f07f362SMatthew G. Knepley           for (g = 0; g < dim; g++) {
12247f07f362SMatthew G. Knepley             J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
12257f07f362SMatthew G. Knepley           }
12267f07f362SMatthew G. Knepley         }
12277f07f362SMatthew G. Knepley       }
12283bc0b13bSBarry Smith       ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
12297f07f362SMatthew G. Knepley     }
1230923591dfSMatthew G. Knepley     if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
12317f07f362SMatthew G. Knepley   } else if (numCoords == 6) {
12327f07f362SMatthew G. Knepley     const PetscInt dim = 2;
12337f07f362SMatthew G. Knepley 
12347f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1235ccd2543fSMatthew G Knepley     if (J)    {
1236ccd2543fSMatthew G Knepley       for (d = 0; d < dim; d++) {
1237ccd2543fSMatthew G Knepley         for (f = 0; f < dim; f++) {
1238ccd2543fSMatthew G Knepley           J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
1239ccd2543fSMatthew G Knepley         }
1240ccd2543fSMatthew G Knepley       }
12413bc0b13bSBarry Smith       ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1242923591dfSMatthew G. Knepley       DMPlex_Det2D_Internal(detJ, J);
1243ccd2543fSMatthew G Knepley     }
1244923591dfSMatthew G. Knepley     if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1245796f034aSJed Brown   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %D != 6 or 9", numCoords);
1246ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
1247ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1248ccd2543fSMatthew G Knepley }
1249ccd2543fSMatthew G Knepley 
1250dfccc68fSToby Isaac static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1251ccd2543fSMatthew G Knepley {
1252ccd2543fSMatthew G Knepley   PetscSection   coordSection;
1253ccd2543fSMatthew G Knepley   Vec            coordinates;
1254a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
12550d29256aSToby Isaac   PetscInt       numCoords, numSelfCoords = 0, d, f, g, pStart, pEnd;
1256ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
1257ccd2543fSMatthew G Knepley 
1258ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1259ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
126069d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
12610d29256aSToby Isaac   ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr);
12620d29256aSToby Isaac   if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);}
126399dec3a6SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
126471f58de1SToby Isaac   numCoords = numSelfCoords ? numSelfCoords : numCoords;
1265dfccc68fSToby Isaac   if (!Nq) {
12667f07f362SMatthew G. Knepley     *detJ = 0.0;
126799dec3a6SMatthew G. Knepley     if (numCoords == 12) {
126899dec3a6SMatthew G. Knepley       const PetscInt dim = 3;
126999dec3a6SMatthew 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};
127099dec3a6SMatthew G. Knepley 
1271dfccc68fSToby Isaac       if (v)   {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);}
1272741bfc07SMatthew G. Knepley       ierr = DMPlexComputeProjection3Dto2D(numCoords, coords, R);CHKERRQ(ierr);
127399dec3a6SMatthew G. Knepley       if (J)    {
127499dec3a6SMatthew G. Knepley         const PetscInt pdim = 2;
127599dec3a6SMatthew G. Knepley 
127699dec3a6SMatthew G. Knepley         for (d = 0; d < pdim; d++) {
127799dec3a6SMatthew G. Knepley           J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
12780892ce5fSStefano Zampini           J0[d*dim+1] = 0.5*(PetscRealPart(coords[2*pdim+d]) - PetscRealPart(coords[1*pdim+d]));
127999dec3a6SMatthew G. Knepley         }
12803bc0b13bSBarry Smith         ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1281923591dfSMatthew G. Knepley         DMPlex_Det3D_Internal(detJ, J0);
128299dec3a6SMatthew G. Knepley         for (d = 0; d < dim; d++) {
128399dec3a6SMatthew G. Knepley           for (f = 0; f < dim; f++) {
128499dec3a6SMatthew G. Knepley             J[d*dim+f] = 0.0;
128599dec3a6SMatthew G. Knepley             for (g = 0; g < dim; g++) {
128699dec3a6SMatthew G. Knepley               J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
128799dec3a6SMatthew G. Knepley             }
128899dec3a6SMatthew G. Knepley           }
128999dec3a6SMatthew G. Knepley         }
12903bc0b13bSBarry Smith         ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
129199dec3a6SMatthew G. Knepley       }
1292923591dfSMatthew G. Knepley       if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
129371f58de1SToby Isaac     } else if (numCoords == 8) {
129499dec3a6SMatthew G. Knepley       const PetscInt dim = 2;
129599dec3a6SMatthew G. Knepley 
1296dfccc68fSToby Isaac       if (v)   {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);}
1297ccd2543fSMatthew G Knepley       if (J)    {
1298ccd2543fSMatthew G Knepley         for (d = 0; d < dim; d++) {
129999dec3a6SMatthew G. Knepley           J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
130099dec3a6SMatthew G. Knepley           J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1301ccd2543fSMatthew G Knepley         }
13023bc0b13bSBarry Smith         ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1303923591dfSMatthew G. Knepley         DMPlex_Det2D_Internal(detJ, J);
1304ccd2543fSMatthew G Knepley       }
1305923591dfSMatthew G. Knepley       if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1306796f034aSJed Brown     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords);
1307dfccc68fSToby Isaac   } else {
1308dfccc68fSToby Isaac     const PetscInt Nv = 4;
1309dfccc68fSToby Isaac     const PetscInt dimR = 2;
1310dfccc68fSToby Isaac     const PetscInt zToPlex[4] = {0, 1, 3, 2};
1311dfccc68fSToby Isaac     PetscReal zOrder[12];
1312dfccc68fSToby Isaac     PetscReal zCoeff[12];
1313dfccc68fSToby Isaac     PetscInt  i, j, k, l, dim;
1314dfccc68fSToby Isaac 
1315dfccc68fSToby Isaac     if (numCoords == 12) {
1316dfccc68fSToby Isaac       dim = 3;
1317dfccc68fSToby Isaac     } else if (numCoords == 8) {
1318dfccc68fSToby Isaac       dim = 2;
1319dfccc68fSToby Isaac     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords);
1320dfccc68fSToby Isaac     for (i = 0; i < Nv; i++) {
1321dfccc68fSToby Isaac       PetscInt zi = zToPlex[i];
1322dfccc68fSToby Isaac 
1323dfccc68fSToby Isaac       for (j = 0; j < dim; j++) {
1324dfccc68fSToby Isaac         zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]);
1325dfccc68fSToby Isaac       }
1326dfccc68fSToby Isaac     }
1327dfccc68fSToby Isaac     for (j = 0; j < dim; j++) {
1328dfccc68fSToby Isaac       zCoeff[dim * 0 + j] = 0.25 * (  zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1329dfccc68fSToby Isaac       zCoeff[dim * 1 + j] = 0.25 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1330dfccc68fSToby Isaac       zCoeff[dim * 2 + j] = 0.25 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1331dfccc68fSToby Isaac       zCoeff[dim * 3 + j] = 0.25 * (  zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1332dfccc68fSToby Isaac     }
1333dfccc68fSToby Isaac     for (i = 0; i < Nq; i++) {
1334dfccc68fSToby Isaac       PetscReal xi = points[dimR * i], eta = points[dimR * i + 1];
1335dfccc68fSToby Isaac 
1336dfccc68fSToby Isaac       if (v) {
1337dfccc68fSToby Isaac         PetscReal extPoint[4];
1338dfccc68fSToby Isaac 
1339dfccc68fSToby Isaac         extPoint[0] = 1.;
1340dfccc68fSToby Isaac         extPoint[1] = xi;
1341dfccc68fSToby Isaac         extPoint[2] = eta;
1342dfccc68fSToby Isaac         extPoint[3] = xi * eta;
1343dfccc68fSToby Isaac         for (j = 0; j < dim; j++) {
1344dfccc68fSToby Isaac           PetscReal val = 0.;
1345dfccc68fSToby Isaac 
1346dfccc68fSToby Isaac           for (k = 0; k < Nv; k++) {
1347dfccc68fSToby Isaac             val += extPoint[k] * zCoeff[dim * k + j];
1348dfccc68fSToby Isaac           }
1349dfccc68fSToby Isaac           v[i * dim + j] = val;
1350dfccc68fSToby Isaac         }
1351dfccc68fSToby Isaac       }
1352dfccc68fSToby Isaac       if (J) {
1353dfccc68fSToby Isaac         PetscReal extJ[8];
1354dfccc68fSToby Isaac 
1355dfccc68fSToby Isaac         extJ[0] = 0.;
1356dfccc68fSToby Isaac         extJ[1] = 0.;
1357dfccc68fSToby Isaac         extJ[2] = 1.;
1358dfccc68fSToby Isaac         extJ[3] = 0.;
1359dfccc68fSToby Isaac         extJ[4] = 0.;
1360dfccc68fSToby Isaac         extJ[5] = 1.;
1361dfccc68fSToby Isaac         extJ[6] = eta;
1362dfccc68fSToby Isaac         extJ[7] = xi;
1363dfccc68fSToby Isaac         for (j = 0; j < dim; j++) {
1364dfccc68fSToby Isaac           for (k = 0; k < dimR; k++) {
1365dfccc68fSToby Isaac             PetscReal val = 0.;
1366dfccc68fSToby Isaac 
1367dfccc68fSToby Isaac             for (l = 0; l < Nv; l++) {
1368dfccc68fSToby Isaac               val += zCoeff[dim * l + j] * extJ[dimR * l + k];
1369dfccc68fSToby Isaac             }
1370dfccc68fSToby Isaac             J[i * dim * dim + dim * j + k] = val;
1371dfccc68fSToby Isaac           }
1372dfccc68fSToby Isaac         }
1373dfccc68fSToby Isaac         if (dim == 3) { /* put the cross product in the third component of the Jacobian */
1374dfccc68fSToby Isaac           PetscReal x, y, z;
1375dfccc68fSToby Isaac           PetscReal *iJ = &J[i * dim * dim];
1376dfccc68fSToby Isaac           PetscReal norm;
1377dfccc68fSToby Isaac 
1378dfccc68fSToby Isaac           x = iJ[1 * dim + 0] * iJ[2 * dim + 1] - iJ[1 * dim + 1] * iJ[2 * dim + 0];
1379dfccc68fSToby Isaac           y = iJ[0 * dim + 1] * iJ[2 * dim + 0] - iJ[0 * dim + 0] * iJ[2 * dim + 1];
1380dfccc68fSToby Isaac           z = iJ[0 * dim + 0] * iJ[1 * dim + 1] - iJ[0 * dim + 1] * iJ[1 * dim + 0];
1381dfccc68fSToby Isaac           norm = PetscSqrtReal(x * x + y * y + z * z);
1382dfccc68fSToby Isaac           iJ[2] = x / norm;
1383dfccc68fSToby Isaac           iJ[5] = y / norm;
1384dfccc68fSToby Isaac           iJ[8] = z / norm;
1385dfccc68fSToby Isaac           DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]);
1386dfccc68fSToby Isaac           if (invJ) {DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);}
1387dfccc68fSToby Isaac         } else {
1388dfccc68fSToby Isaac           DMPlex_Det2D_Internal(&detJ[i], &J[i * dim * dim]);
1389dfccc68fSToby Isaac           if (invJ) {DMPlex_Invert2D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);}
1390dfccc68fSToby Isaac         }
1391dfccc68fSToby Isaac       }
1392dfccc68fSToby Isaac     }
1393dfccc68fSToby Isaac   }
139499dec3a6SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
1395ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1396ccd2543fSMatthew G Knepley }
1397ccd2543fSMatthew G Knepley 
1398ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1399ccd2543fSMatthew G Knepley {
1400ccd2543fSMatthew G Knepley   PetscSection   coordSection;
1401ccd2543fSMatthew G Knepley   Vec            coordinates;
1402a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
1403ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
140499dec3a6SMatthew G. Knepley   PetscInt       d;
1405ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
1406ccd2543fSMatthew G Knepley 
1407ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1408ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
140969d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1410ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
14117f07f362SMatthew G. Knepley   *detJ = 0.0;
14127f07f362SMatthew G. Knepley   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1413ccd2543fSMatthew G Knepley   if (J)    {
1414ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
1415f0df753eSMatthew G. Knepley       /* I orient with outward face normals */
1416f0df753eSMatthew G. Knepley       J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d]));
1417f0df753eSMatthew G. Knepley       J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
1418f0df753eSMatthew G. Knepley       J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1419ccd2543fSMatthew G Knepley     }
14203bc0b13bSBarry Smith     ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
1421923591dfSMatthew G. Knepley     DMPlex_Det3D_Internal(detJ, J);
1422ccd2543fSMatthew G Knepley   }
1423923591dfSMatthew G. Knepley   if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1424ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1425ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1426ccd2543fSMatthew G Knepley }
1427ccd2543fSMatthew G Knepley 
1428dfccc68fSToby Isaac static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1429ccd2543fSMatthew G Knepley {
1430ccd2543fSMatthew G Knepley   PetscSection   coordSection;
1431ccd2543fSMatthew G Knepley   Vec            coordinates;
1432a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
1433ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
1434ccd2543fSMatthew G Knepley   PetscInt       d;
1435ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
1436ccd2543fSMatthew G Knepley 
1437ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1438ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
143969d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1440ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1441dfccc68fSToby Isaac   if (!Nq) {
14427f07f362SMatthew G. Knepley     *detJ = 0.0;
1443dfccc68fSToby Isaac     if (v)   {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);}
1444ccd2543fSMatthew G Knepley     if (J)    {
1445ccd2543fSMatthew G Knepley       for (d = 0; d < dim; d++) {
1446f0df753eSMatthew G. Knepley         J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1447f0df753eSMatthew G. Knepley         J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
1448f0df753eSMatthew G. Knepley         J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d]));
1449ccd2543fSMatthew G Knepley       }
14503bc0b13bSBarry Smith       ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
1451923591dfSMatthew G. Knepley       DMPlex_Det3D_Internal(detJ, J);
1452ccd2543fSMatthew G Knepley     }
1453923591dfSMatthew G. Knepley     if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1454dfccc68fSToby Isaac   } else {
1455dfccc68fSToby Isaac     const PetscInt Nv = 8;
1456dfccc68fSToby Isaac     const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6};
1457dfccc68fSToby Isaac     const PetscInt dim = 3;
1458dfccc68fSToby Isaac     const PetscInt dimR = 3;
1459dfccc68fSToby Isaac     PetscReal zOrder[24];
1460dfccc68fSToby Isaac     PetscReal zCoeff[24];
1461dfccc68fSToby Isaac     PetscInt  i, j, k, l;
1462dfccc68fSToby Isaac 
1463dfccc68fSToby Isaac     for (i = 0; i < Nv; i++) {
1464dfccc68fSToby Isaac       PetscInt zi = zToPlex[i];
1465dfccc68fSToby Isaac 
1466dfccc68fSToby Isaac       for (j = 0; j < dim; j++) {
1467dfccc68fSToby Isaac         zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]);
1468dfccc68fSToby Isaac       }
1469dfccc68fSToby Isaac     }
1470dfccc68fSToby Isaac     for (j = 0; j < dim; j++) {
1471dfccc68fSToby 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]);
1472dfccc68fSToby 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]);
1473dfccc68fSToby 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]);
1474dfccc68fSToby 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]);
1475dfccc68fSToby 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]);
1476dfccc68fSToby 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]);
1477dfccc68fSToby 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]);
1478dfccc68fSToby 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]);
1479dfccc68fSToby Isaac     }
1480dfccc68fSToby Isaac     for (i = 0; i < Nq; i++) {
1481dfccc68fSToby Isaac       PetscReal xi = points[dimR * i], eta = points[dimR * i + 1], theta = points[dimR * i + 2];
1482dfccc68fSToby Isaac 
1483dfccc68fSToby Isaac       if (v) {
148491d2b7ceSToby Isaac         PetscReal extPoint[8];
1485dfccc68fSToby Isaac 
1486dfccc68fSToby Isaac         extPoint[0] = 1.;
1487dfccc68fSToby Isaac         extPoint[1] = xi;
1488dfccc68fSToby Isaac         extPoint[2] = eta;
1489dfccc68fSToby Isaac         extPoint[3] = xi * eta;
1490dfccc68fSToby Isaac         extPoint[4] = theta;
1491dfccc68fSToby Isaac         extPoint[5] = theta * xi;
1492dfccc68fSToby Isaac         extPoint[6] = theta * eta;
1493dfccc68fSToby Isaac         extPoint[7] = theta * eta * xi;
1494dfccc68fSToby Isaac         for (j = 0; j < dim; j++) {
1495dfccc68fSToby Isaac           PetscReal val = 0.;
1496dfccc68fSToby Isaac 
1497dfccc68fSToby Isaac           for (k = 0; k < Nv; k++) {
1498dfccc68fSToby Isaac             val += extPoint[k] * zCoeff[dim * k + j];
1499dfccc68fSToby Isaac           }
1500dfccc68fSToby Isaac           v[i * dim + j] = val;
1501dfccc68fSToby Isaac         }
1502dfccc68fSToby Isaac       }
1503dfccc68fSToby Isaac       if (J) {
1504dfccc68fSToby Isaac         PetscReal extJ[24];
1505dfccc68fSToby Isaac 
1506dfccc68fSToby Isaac         extJ[0]  = 0.         ; extJ[1]  = 0.        ; extJ[2]  = 0.      ;
1507dfccc68fSToby Isaac         extJ[3]  = 1.         ; extJ[4]  = 0.        ; extJ[5]  = 0.      ;
1508dfccc68fSToby Isaac         extJ[6]  = 0.         ; extJ[7]  = 1.        ; extJ[8]  = 0.      ;
1509dfccc68fSToby Isaac         extJ[9]  = eta        ; extJ[10] = xi        ; extJ[11] = 0.      ;
1510dfccc68fSToby Isaac         extJ[12] = 0.         ; extJ[13] = 0.        ; extJ[14] = 1.      ;
1511dfccc68fSToby Isaac         extJ[15] = theta      ; extJ[16] = 0.        ; extJ[17] = xi      ;
1512dfccc68fSToby Isaac         extJ[18] = 0.         ; extJ[19] = theta     ; extJ[20] = eta     ;
1513dfccc68fSToby Isaac         extJ[21] = theta * eta; extJ[22] = theta * xi; extJ[23] = eta * xi;
1514dfccc68fSToby Isaac 
1515dfccc68fSToby Isaac         for (j = 0; j < dim; j++) {
1516dfccc68fSToby Isaac           for (k = 0; k < dimR; k++) {
1517dfccc68fSToby Isaac             PetscReal val = 0.;
1518dfccc68fSToby Isaac 
1519dfccc68fSToby Isaac             for (l = 0; l < Nv; l++) {
1520dfccc68fSToby Isaac               val += zCoeff[dim * l + j] * extJ[dimR * l + k];
1521dfccc68fSToby Isaac             }
1522dfccc68fSToby Isaac             J[i * dim * dim + dim * j + k] = val;
1523dfccc68fSToby Isaac           }
1524dfccc68fSToby Isaac         }
1525dfccc68fSToby Isaac         DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]);
1526dfccc68fSToby Isaac         if (invJ) {DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);}
1527dfccc68fSToby Isaac       }
1528dfccc68fSToby Isaac     }
1529dfccc68fSToby Isaac   }
1530ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1531ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1532ccd2543fSMatthew G Knepley }
1533ccd2543fSMatthew G Knepley 
1534dfccc68fSToby Isaac static PetscErrorCode DMPlexComputeCellGeometryFEM_Implicit(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
1535dfccc68fSToby Isaac {
1536dfccc68fSToby Isaac   PetscInt        depth, dim, coordDim, coneSize, i;
1537dfccc68fSToby Isaac   PetscInt        Nq = 0;
1538dfccc68fSToby Isaac   const PetscReal *points = NULL;
1539dfccc68fSToby Isaac   DMLabel         depthLabel;
1540c330f8ffSToby Isaac   PetscReal       xi0[3] = {-1.,-1.,-1.}, v0[3], J0[9], detJ0;
1541dfccc68fSToby Isaac   PetscBool       isAffine = PETSC_TRUE;
1542dfccc68fSToby Isaac   PetscErrorCode  ierr;
1543dfccc68fSToby Isaac 
1544dfccc68fSToby Isaac   PetscFunctionBegin;
1545dfccc68fSToby Isaac   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1546dfccc68fSToby Isaac   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
1547dfccc68fSToby Isaac   ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr);
1548dfccc68fSToby Isaac   ierr = DMLabelGetValue(depthLabel, cell, &dim);CHKERRQ(ierr);
1549dfccc68fSToby Isaac   if (depth == 1 && dim == 1) {
1550dfccc68fSToby Isaac     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1551dfccc68fSToby Isaac   }
1552dfccc68fSToby Isaac   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
1553dfccc68fSToby Isaac   if (coordDim > 3) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported coordinate dimension %D > 3", coordDim);
15549c3cf19fSMatthew G. Knepley   if (quad) {ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, &points, NULL);CHKERRQ(ierr);}
1555dfccc68fSToby Isaac   switch (dim) {
1556dfccc68fSToby Isaac   case 0:
1557dfccc68fSToby Isaac     ierr = DMPlexComputePointGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr);
1558dfccc68fSToby Isaac     isAffine = PETSC_FALSE;
1559dfccc68fSToby Isaac     break;
1560dfccc68fSToby Isaac   case 1:
15617318780aSToby Isaac     if (Nq) {
15627318780aSToby Isaac       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr);
15637318780aSToby Isaac     } else {
15647318780aSToby Isaac       ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr);
15657318780aSToby Isaac     }
1566dfccc68fSToby Isaac     break;
1567dfccc68fSToby Isaac   case 2:
1568dfccc68fSToby Isaac     switch (coneSize) {
1569dfccc68fSToby Isaac     case 3:
15707318780aSToby Isaac       if (Nq) {
15717318780aSToby Isaac         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr);
15727318780aSToby Isaac       } else {
15737318780aSToby Isaac         ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr);
15747318780aSToby Isaac       }
1575dfccc68fSToby Isaac       break;
1576dfccc68fSToby Isaac     case 4:
1577dfccc68fSToby Isaac       ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, Nq, points, v, J, invJ, detJ);CHKERRQ(ierr);
1578dfccc68fSToby Isaac       isAffine = PETSC_FALSE;
1579dfccc68fSToby Isaac       break;
1580dfccc68fSToby Isaac     default:
1581dfccc68fSToby Isaac       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell);
1582dfccc68fSToby Isaac     }
1583dfccc68fSToby Isaac     break;
1584dfccc68fSToby Isaac   case 3:
1585dfccc68fSToby Isaac     switch (coneSize) {
1586dfccc68fSToby Isaac     case 4:
15877318780aSToby Isaac       if (Nq) {
15887318780aSToby Isaac         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr);
15897318780aSToby Isaac       } else {
15907318780aSToby Isaac         ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr);
15917318780aSToby Isaac       }
1592dfccc68fSToby Isaac       break;
1593dfccc68fSToby Isaac     case 6: /* Faces */
1594dfccc68fSToby Isaac     case 8: /* Vertices */
1595dfccc68fSToby Isaac       ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, Nq, points, v, J, invJ, detJ);CHKERRQ(ierr);
1596dfccc68fSToby Isaac       isAffine = PETSC_FALSE;
1597dfccc68fSToby Isaac       break;
1598dfccc68fSToby Isaac     default:
1599dfccc68fSToby Isaac       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell);
1600dfccc68fSToby Isaac     }
1601dfccc68fSToby Isaac     break;
1602dfccc68fSToby Isaac   default:
1603dfccc68fSToby Isaac     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
1604dfccc68fSToby Isaac   }
16057318780aSToby Isaac   if (isAffine && Nq) {
1606dfccc68fSToby Isaac     if (v) {
1607dfccc68fSToby Isaac       for (i = 0; i < Nq; i++) {
1608c330f8ffSToby Isaac         CoordinatesRefToReal(coordDim, dim, xi0, v0, J0, &points[dim * i], &v[coordDim * i]);
1609dfccc68fSToby Isaac       }
1610dfccc68fSToby Isaac     }
16117318780aSToby Isaac     if (detJ) {
16127318780aSToby Isaac       for (i = 0; i < Nq; i++) {
16137318780aSToby Isaac         detJ[i] = detJ0;
1614dfccc68fSToby Isaac       }
16157318780aSToby Isaac     }
16167318780aSToby Isaac     if (J) {
16177318780aSToby Isaac       PetscInt k;
16187318780aSToby Isaac 
16197318780aSToby Isaac       for (i = 0, k = 0; i < Nq; i++) {
1620dfccc68fSToby Isaac         PetscInt j;
1621dfccc68fSToby Isaac 
16227318780aSToby Isaac         for (j = 0; j < coordDim * coordDim; j++, k++) {
16237318780aSToby Isaac           J[k] = J0[j];
16247318780aSToby Isaac         }
16257318780aSToby Isaac       }
16267318780aSToby Isaac     }
16277318780aSToby Isaac     if (invJ) {
16287318780aSToby Isaac       PetscInt k;
16297318780aSToby Isaac       switch (coordDim) {
16307318780aSToby Isaac       case 0:
16317318780aSToby Isaac         break;
16327318780aSToby Isaac       case 1:
16337318780aSToby Isaac         invJ[0] = 1./J0[0];
16347318780aSToby Isaac         break;
16357318780aSToby Isaac       case 2:
16367318780aSToby Isaac         DMPlex_Invert2D_Internal(invJ, J0, detJ0);
16377318780aSToby Isaac         break;
16387318780aSToby Isaac       case 3:
16397318780aSToby Isaac         DMPlex_Invert3D_Internal(invJ, J0, detJ0);
16407318780aSToby Isaac         break;
16417318780aSToby Isaac       }
16427318780aSToby Isaac       for (i = 1, k = coordDim * coordDim; i < Nq; i++) {
16437318780aSToby Isaac         PetscInt j;
16447318780aSToby Isaac 
16457318780aSToby Isaac         for (j = 0; j < coordDim * coordDim; j++, k++) {
16467318780aSToby Isaac           invJ[k] = invJ[j];
16477318780aSToby Isaac         }
1648dfccc68fSToby Isaac       }
1649dfccc68fSToby Isaac     }
1650dfccc68fSToby Isaac   }
1651dfccc68fSToby Isaac   PetscFunctionReturn(0);
1652dfccc68fSToby Isaac }
1653dfccc68fSToby Isaac 
1654ccd2543fSMatthew G Knepley /*@C
16558e0841e0SMatthew G. Knepley   DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell
1656ccd2543fSMatthew G Knepley 
1657d083f849SBarry Smith   Collective on dm
1658ccd2543fSMatthew G Knepley 
1659ccd2543fSMatthew G Knepley   Input Arguments:
1660ccd2543fSMatthew G Knepley + dm   - the DM
1661ccd2543fSMatthew G Knepley - cell - the cell
1662ccd2543fSMatthew G Knepley 
1663ccd2543fSMatthew G Knepley   Output Arguments:
1664ccd2543fSMatthew G Knepley + v0   - the translation part of this affine transform
1665ccd2543fSMatthew G Knepley . J    - the Jacobian of the transform from the reference element
1666ccd2543fSMatthew G Knepley . invJ - the inverse of the Jacobian
1667ccd2543fSMatthew G Knepley - detJ - the Jacobian determinant
1668ccd2543fSMatthew G Knepley 
1669ccd2543fSMatthew G Knepley   Level: advanced
1670ccd2543fSMatthew G Knepley 
1671ccd2543fSMatthew G Knepley   Fortran Notes:
1672ccd2543fSMatthew G Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
1673ccd2543fSMatthew G Knepley   include petsc.h90 in your code.
1674ccd2543fSMatthew G Knepley 
1675e8964c0aSStefano Zampini .seealso: DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinates()
1676ccd2543fSMatthew G Knepley @*/
16778e0841e0SMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
1678ccd2543fSMatthew G Knepley {
1679ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
1680ccd2543fSMatthew G Knepley 
1681ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1682dfccc68fSToby Isaac   ierr = DMPlexComputeCellGeometryFEM_Implicit(dm,cell,NULL,v0,J,invJ,detJ);CHKERRQ(ierr);
16838e0841e0SMatthew G. Knepley   PetscFunctionReturn(0);
16848e0841e0SMatthew G. Knepley }
16858e0841e0SMatthew G. Knepley 
1686dfccc68fSToby Isaac static PetscErrorCode DMPlexComputeCellGeometryFEM_FE(DM dm, PetscFE fe, PetscInt point, PetscQuadrature quad, PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
16878e0841e0SMatthew G. Knepley {
1688dfccc68fSToby Isaac   PetscQuadrature  feQuad;
16898e0841e0SMatthew G. Knepley   PetscSection     coordSection;
16908e0841e0SMatthew G. Knepley   Vec              coordinates;
16918e0841e0SMatthew G. Knepley   PetscScalar     *coords = NULL;
16928e0841e0SMatthew G. Knepley   const PetscReal *quadPoints;
1693f960e424SToby Isaac   PetscReal       *basisDer, *basis, detJt;
1694f960e424SToby Isaac   PetscInt         dim, cdim, pdim, qdim, Nq, numCoords, q;
16958e0841e0SMatthew G. Knepley   PetscErrorCode   ierr;
16968e0841e0SMatthew G. Knepley 
16978e0841e0SMatthew G. Knepley   PetscFunctionBegin;
16988e0841e0SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
16998e0841e0SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
17008e0841e0SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr);
17018e0841e0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
17028e0841e0SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
1703dfccc68fSToby Isaac   if (!quad) { /* use the first point of the first functional of the dual space */
1704dfccc68fSToby Isaac     PetscDualSpace dsp;
1705dfccc68fSToby Isaac 
1706dfccc68fSToby Isaac     ierr = PetscFEGetDualSpace(fe, &dsp);CHKERRQ(ierr);
1707dfccc68fSToby Isaac     ierr = PetscDualSpaceGetFunctional(dsp, 0, &quad);CHKERRQ(ierr);
17089c3cf19fSMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, &qdim, NULL, &Nq, &quadPoints, NULL);CHKERRQ(ierr);
1709dfccc68fSToby Isaac     Nq = 1;
1710dfccc68fSToby Isaac   } else {
17119c3cf19fSMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, &qdim, NULL, &Nq, &quadPoints, NULL);CHKERRQ(ierr);
1712dfccc68fSToby Isaac   }
171391d2b7ceSToby Isaac   ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr);
1714dfccc68fSToby Isaac   ierr = PetscFEGetQuadrature(fe, &feQuad);CHKERRQ(ierr);
1715dfccc68fSToby Isaac   if (feQuad == quad) {
17168135c375SStefano Zampini     ierr = PetscFEGetDefaultTabulation(fe, &basis, J ? &basisDer : NULL, NULL);CHKERRQ(ierr);
17178e0841e0SMatthew 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);
1718dfccc68fSToby Isaac   } else {
17198135c375SStefano Zampini     ierr = PetscFEGetTabulation(fe, Nq, quadPoints, &basis, J ? &basisDer : NULL, NULL);CHKERRQ(ierr);
1720dfccc68fSToby Isaac   }
1721dfccc68fSToby Isaac   if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %d != quadrature dimension %d", dim, qdim);
1722dfccc68fSToby Isaac   if (v) {
1723580bdb30SBarry Smith     ierr = PetscArrayzero(v, Nq*cdim);CHKERRQ(ierr);
1724f960e424SToby Isaac     for (q = 0; q < Nq; ++q) {
1725f960e424SToby Isaac       PetscInt i, k;
1726f960e424SToby Isaac 
1727f960e424SToby Isaac       for (k = 0; k < pdim; ++k)
1728f960e424SToby Isaac         for (i = 0; i < cdim; ++i)
1729dfccc68fSToby Isaac           v[q*cdim + i] += basis[q*pdim + k] * PetscRealPart(coords[k*cdim + i]);
1730f960e424SToby Isaac       ierr = PetscLogFlops(2.0*pdim*cdim);CHKERRQ(ierr);
1731f960e424SToby Isaac     }
1732f960e424SToby Isaac   }
17338e0841e0SMatthew G. Knepley   if (J) {
1734580bdb30SBarry Smith     ierr = PetscArrayzero(J, Nq*cdim*cdim);CHKERRQ(ierr);
17358e0841e0SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
17368e0841e0SMatthew G. Knepley       PetscInt i, j, k, c, r;
17378e0841e0SMatthew G. Knepley 
17388e0841e0SMatthew G. Knepley       /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */
17398e0841e0SMatthew G. Knepley       for (k = 0; k < pdim; ++k)
17408e0841e0SMatthew G. Knepley         for (j = 0; j < dim; ++j)
17418e0841e0SMatthew G. Knepley           for (i = 0; i < cdim; ++i)
1742dfccc68fSToby Isaac             J[(q*cdim + i)*cdim + j] += basisDer[(q*pdim + k)*dim + j] * PetscRealPart(coords[k*cdim + i]);
17433bc0b13bSBarry Smith       ierr = PetscLogFlops(2.0*pdim*dim*cdim);CHKERRQ(ierr);
17448e0841e0SMatthew G. Knepley       if (cdim > dim) {
17458e0841e0SMatthew G. Knepley         for (c = dim; c < cdim; ++c)
17468e0841e0SMatthew G. Knepley           for (r = 0; r < cdim; ++r)
17478e0841e0SMatthew G. Knepley             J[r*cdim+c] = r == c ? 1.0 : 0.0;
17488e0841e0SMatthew G. Knepley       }
1749f960e424SToby Isaac       if (!detJ && !invJ) continue;
1750a63b72c6SToby Isaac       detJt = 0.;
17518e0841e0SMatthew G. Knepley       switch (cdim) {
17528e0841e0SMatthew G. Knepley       case 3:
1753037dc194SToby Isaac         DMPlex_Det3D_Internal(&detJt, &J[q*cdim*dim]);
1754037dc194SToby Isaac         if (invJ) {DMPlex_Invert3D_Internal(&invJ[q*cdim*dim], &J[q*cdim*dim], detJt);}
175517fe8556SMatthew G. Knepley         break;
175649dc4407SMatthew G. Knepley       case 2:
17579f328543SToby Isaac         DMPlex_Det2D_Internal(&detJt, &J[q*cdim*dim]);
1758037dc194SToby Isaac         if (invJ) {DMPlex_Invert2D_Internal(&invJ[q*cdim*dim], &J[q*cdim*dim], detJt);}
175949dc4407SMatthew G. Knepley         break;
17608e0841e0SMatthew G. Knepley       case 1:
1761037dc194SToby Isaac         detJt = J[q*cdim*dim];
1762037dc194SToby Isaac         if (invJ) invJ[q*cdim*dim] = 1.0/detJt;
176349dc4407SMatthew G. Knepley       }
1764f960e424SToby Isaac       if (detJ) detJ[q] = detJt;
176549dc4407SMatthew G. Knepley     }
176649dc4407SMatthew G. Knepley   }
1767037dc194SToby Isaac   else if (detJ || invJ) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Need J to compute invJ or detJ");
1768dfccc68fSToby Isaac   if (feQuad != quad) {
17698135c375SStefano Zampini     ierr = PetscFERestoreTabulation(fe, Nq, quadPoints, &basis, J ? &basisDer : NULL, NULL);CHKERRQ(ierr);
1770dfccc68fSToby Isaac   }
17718e0841e0SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr);
17728e0841e0SMatthew G. Knepley   PetscFunctionReturn(0);
17738e0841e0SMatthew G. Knepley }
17748e0841e0SMatthew G. Knepley 
17758e0841e0SMatthew G. Knepley /*@C
17768e0841e0SMatthew G. Knepley   DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell
17778e0841e0SMatthew G. Knepley 
1778d083f849SBarry Smith   Collective on dm
17798e0841e0SMatthew G. Knepley 
17808e0841e0SMatthew G. Knepley   Input Arguments:
17818e0841e0SMatthew G. Knepley + dm   - the DM
17828e0841e0SMatthew G. Knepley . cell - the cell
1783dfccc68fSToby Isaac - quad - the quadrature containing the points in the reference element where the geometry will be evaluated.  If quad == NULL, geometry will be
1784dfccc68fSToby Isaac          evaluated at the first vertex of the reference element
17858e0841e0SMatthew G. Knepley 
17868e0841e0SMatthew G. Knepley   Output Arguments:
1787dfccc68fSToby Isaac + v    - the image of the transformed quadrature points, otherwise the image of the first vertex in the closure of the reference element
17888e0841e0SMatthew G. Knepley . J    - the Jacobian of the transform from the reference element at each quadrature point
17898e0841e0SMatthew G. Knepley . invJ - the inverse of the Jacobian at each quadrature point
17908e0841e0SMatthew G. Knepley - detJ - the Jacobian determinant at each quadrature point
17918e0841e0SMatthew G. Knepley 
17928e0841e0SMatthew G. Knepley   Level: advanced
17938e0841e0SMatthew G. Knepley 
17948e0841e0SMatthew G. Knepley   Fortran Notes:
17958e0841e0SMatthew G. Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
17968e0841e0SMatthew G. Knepley   include petsc.h90 in your code.
17978e0841e0SMatthew G. Knepley 
1798e8964c0aSStefano Zampini .seealso: DMGetCoordinateSection(), DMGetCoordinates()
17998e0841e0SMatthew G. Knepley @*/
1800dfccc68fSToby Isaac PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
18018e0841e0SMatthew G. Knepley {
1802bb4a5db5SMatthew G. Knepley   DM             cdm;
1803dfccc68fSToby Isaac   PetscFE        fe = NULL;
18048e0841e0SMatthew G. Knepley   PetscErrorCode ierr;
18058e0841e0SMatthew G. Knepley 
18068e0841e0SMatthew G. Knepley   PetscFunctionBegin;
18072d89661fSMatthew G. Knepley   PetscValidPointer(detJ, 7);
1808bb4a5db5SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
1809bb4a5db5SMatthew G. Knepley   if (cdm) {
1810dfccc68fSToby Isaac     PetscClassId id;
1811dfccc68fSToby Isaac     PetscInt     numFields;
1812e5e52638SMatthew G. Knepley     PetscDS      prob;
1813dfccc68fSToby Isaac     PetscObject  disc;
1814dfccc68fSToby Isaac 
1815bb4a5db5SMatthew G. Knepley     ierr = DMGetNumFields(cdm, &numFields);CHKERRQ(ierr);
1816dfccc68fSToby Isaac     if (numFields) {
1817bb4a5db5SMatthew G. Knepley       ierr = DMGetDS(cdm, &prob);CHKERRQ(ierr);
1818dfccc68fSToby Isaac       ierr = PetscDSGetDiscretization(prob,0,&disc);CHKERRQ(ierr);
1819dfccc68fSToby Isaac       ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
1820dfccc68fSToby Isaac       if (id == PETSCFE_CLASSID) {
1821dfccc68fSToby Isaac         fe = (PetscFE) disc;
1822dfccc68fSToby Isaac       }
1823dfccc68fSToby Isaac     }
1824dfccc68fSToby Isaac   }
1825dfccc68fSToby Isaac   if (!fe) {ierr = DMPlexComputeCellGeometryFEM_Implicit(dm, cell, quad, v, J, invJ, detJ);CHKERRQ(ierr);}
1826dfccc68fSToby Isaac   else     {ierr = DMPlexComputeCellGeometryFEM_FE(dm, fe, cell, quad, v, J, invJ, detJ);CHKERRQ(ierr);}
1827ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1828ccd2543fSMatthew G Knepley }
1829834e62ceSMatthew G. Knepley 
1830011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1831cc08537eSMatthew G. Knepley {
1832cc08537eSMatthew G. Knepley   PetscSection   coordSection;
1833cc08537eSMatthew G. Knepley   Vec            coordinates;
1834a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
183506e2781eSMatthew G. Knepley   PetscScalar    tmp[2];
1836714b99b6SMatthew G. Knepley   PetscInt       coordSize, d;
1837cc08537eSMatthew G. Knepley   PetscErrorCode ierr;
1838cc08537eSMatthew G. Knepley 
1839cc08537eSMatthew G. Knepley   PetscFunctionBegin;
1840cc08537eSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
184169d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1842cc08537eSMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
18432e17dfb7SMatthew G. Knepley   ierr = DMLocalizeCoordinate_Internal(dm, dim, coords, &coords[dim], tmp);CHKERRQ(ierr);
1844cc08537eSMatthew G. Knepley   if (centroid) {
1845714b99b6SMatthew G. Knepley     for (d = 0; d < dim; ++d) centroid[d] = 0.5*PetscRealPart(coords[d] + tmp[d]);
1846cc08537eSMatthew G. Knepley   }
1847cc08537eSMatthew G. Knepley   if (normal) {
1848a60a936bSMatthew G. Knepley     PetscReal norm;
1849a60a936bSMatthew G. Knepley 
1850714b99b6SMatthew G. Knepley     if (dim != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "We only support 2D edges right now");
185106e2781eSMatthew G. Knepley     normal[0]  = -PetscRealPart(coords[1] - tmp[1]);
185206e2781eSMatthew G. Knepley     normal[1]  =  PetscRealPart(coords[0] - tmp[0]);
1853714b99b6SMatthew G. Knepley     norm       = DMPlex_NormD_Internal(dim, normal);
1854714b99b6SMatthew G. Knepley     for (d = 0; d < dim; ++d) normal[d] /= norm;
1855cc08537eSMatthew G. Knepley   }
1856cc08537eSMatthew G. Knepley   if (vol) {
1857714b99b6SMatthew G. Knepley     *vol = 0.0;
1858714b99b6SMatthew G. Knepley     for (d = 0; d < dim; ++d) *vol += PetscSqr(PetscRealPart(coords[d] - tmp[d]));
1859714b99b6SMatthew G. Knepley     *vol = PetscSqrtReal(*vol);
1860cc08537eSMatthew G. Knepley   }
1861cc08537eSMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1862cc08537eSMatthew G. Knepley   PetscFunctionReturn(0);
1863cc08537eSMatthew G. Knepley }
1864cc08537eSMatthew G. Knepley 
1865cc08537eSMatthew G. Knepley /* Centroid_i = (\sum_n A_n Cn_i ) / A */
1866011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1867cc08537eSMatthew G. Knepley {
1868793a2a13SMatthew G. Knepley   DMLabel        depth;
1869cc08537eSMatthew G. Knepley   PetscSection   coordSection;
1870cc08537eSMatthew G. Knepley   Vec            coordinates;
1871cc08537eSMatthew G. Knepley   PetscScalar   *coords = NULL;
18720a1d6728SMatthew G. Knepley   PetscReal      vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9];
1873793a2a13SMatthew G. Knepley   PetscBool      isHybrid = PETSC_FALSE;
1874793a2a13SMatthew G. Knepley   PetscInt       fv[4] = {0, 1, 2, 3};
1875d80ece95SMatthew G. Knepley   PetscInt       pEndInterior[4], cdepth, tdim = 2, coordSize, numCorners, p, d, e;
1876cc08537eSMatthew G. Knepley   PetscErrorCode ierr;
1877cc08537eSMatthew G. Knepley 
1878cc08537eSMatthew G. Knepley   PetscFunctionBegin;
1879793a2a13SMatthew G. Knepley   /* Must check for hybrid cells because prisms have a different orientation scheme */
1880793a2a13SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
1881793a2a13SMatthew G. Knepley   ierr = DMLabelGetValue(depth, cell, &cdepth);CHKERRQ(ierr);
1882d80ece95SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &pEndInterior[3], &pEndInterior[2], &pEndInterior[1], &pEndInterior[0]);CHKERRQ(ierr);
1883d80ece95SMatthew G. Knepley   if ((pEndInterior[cdepth]) >=0 && (cell >= pEndInterior[cdepth])) isHybrid = PETSC_TRUE;
1884cc08537eSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
18850a1d6728SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr);
188669d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1887cc08537eSMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
18880bce18caSMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
1889793a2a13SMatthew G. Knepley   /* Side faces for hybrid cells are are stored as tensor products */
1890793a2a13SMatthew G. Knepley   if (isHybrid && numCorners == 4) {fv[2] = 3; fv[3] = 2;}
1891793a2a13SMatthew G. Knepley 
1892ceee4971SMatthew G. Knepley   if (dim > 2 && centroid) {
1893ceee4971SMatthew G. Knepley     v0[0] = PetscRealPart(coords[0]);
1894ceee4971SMatthew G. Knepley     v0[1] = PetscRealPart(coords[1]);
1895ceee4971SMatthew G. Knepley     v0[2] = PetscRealPart(coords[2]);
1896ceee4971SMatthew G. Knepley   }
1897011ea5d8SMatthew G. Knepley   if (normal) {
1898011ea5d8SMatthew G. Knepley     if (dim > 2) {
1899793a2a13SMatthew G. Knepley       const PetscReal x0 = PetscRealPart(coords[dim*fv[1]+0] - coords[0]), x1 = PetscRealPart(coords[dim*fv[2]+0] - coords[0]);
1900793a2a13SMatthew G. Knepley       const PetscReal y0 = PetscRealPart(coords[dim*fv[1]+1] - coords[1]), y1 = PetscRealPart(coords[dim*fv[2]+1] - coords[1]);
1901793a2a13SMatthew G. Knepley       const PetscReal z0 = PetscRealPart(coords[dim*fv[1]+2] - coords[2]), z1 = PetscRealPart(coords[dim*fv[2]+2] - coords[2]);
19020a1d6728SMatthew G. Knepley       PetscReal       norm;
19030a1d6728SMatthew G. Knepley 
19040a1d6728SMatthew G. Knepley       normal[0] = y0*z1 - z0*y1;
19050a1d6728SMatthew G. Knepley       normal[1] = z0*x1 - x0*z1;
19060a1d6728SMatthew G. Knepley       normal[2] = x0*y1 - y0*x1;
19078b49ba18SBarry Smith       norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]);
19080a1d6728SMatthew G. Knepley       normal[0] /= norm;
19090a1d6728SMatthew G. Knepley       normal[1] /= norm;
19100a1d6728SMatthew G. Knepley       normal[2] /= norm;
1911011ea5d8SMatthew G. Knepley     } else {
1912011ea5d8SMatthew G. Knepley       for (d = 0; d < dim; ++d) normal[d] = 0.0;
1913011ea5d8SMatthew G. Knepley     }
1914011ea5d8SMatthew G. Knepley   }
1915741bfc07SMatthew G. Knepley   if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D(coordSize, coords, R);CHKERRQ(ierr);}
19160a1d6728SMatthew G. Knepley   for (p = 0; p < numCorners; ++p) {
1917793a2a13SMatthew G. Knepley     const PetscInt pi  = p < 4 ? fv[p] : p;
1918793a2a13SMatthew G. Knepley     const PetscInt pin = p < 3 ? fv[(p+1)%numCorners] : (p+1)%numCorners;
19190a1d6728SMatthew G. Knepley     /* Need to do this copy to get types right */
19200a1d6728SMatthew G. Knepley     for (d = 0; d < tdim; ++d) {
1921793a2a13SMatthew G. Knepley       ctmp[d]      = PetscRealPart(coords[pi*tdim+d]);
1922793a2a13SMatthew G. Knepley       ctmp[tdim+d] = PetscRealPart(coords[pin*tdim+d]);
19230a1d6728SMatthew G. Knepley     }
19240a1d6728SMatthew G. Knepley     Volume_Triangle_Origin_Internal(&vtmp, ctmp);
19250a1d6728SMatthew G. Knepley     vsum += vtmp;
19260a1d6728SMatthew G. Knepley     for (d = 0; d < tdim; ++d) {
19270a1d6728SMatthew G. Knepley       csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp;
19280a1d6728SMatthew G. Knepley     }
19290a1d6728SMatthew G. Knepley   }
19300a1d6728SMatthew G. Knepley   for (d = 0; d < tdim; ++d) {
19310a1d6728SMatthew G. Knepley     csum[d] /= (tdim+1)*vsum;
19320a1d6728SMatthew G. Knepley   }
19330a1d6728SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1934ee6bbdb2SSatish Balay   if (vol) *vol = PetscAbsReal(vsum);
19350a1d6728SMatthew G. Knepley   if (centroid) {
19360a1d6728SMatthew G. Knepley     if (dim > 2) {
19370a1d6728SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
19380a1d6728SMatthew G. Knepley         centroid[d] = v0[d];
19390a1d6728SMatthew G. Knepley         for (e = 0; e < dim; ++e) {
19400a1d6728SMatthew G. Knepley           centroid[d] += R[d*dim+e]*csum[e];
19410a1d6728SMatthew G. Knepley         }
19420a1d6728SMatthew G. Knepley       }
19430a1d6728SMatthew G. Knepley     } else for (d = 0; d < dim; ++d) centroid[d] = csum[d];
19440a1d6728SMatthew G. Knepley   }
1945cc08537eSMatthew G. Knepley   PetscFunctionReturn(0);
1946cc08537eSMatthew G. Knepley }
1947cc08537eSMatthew G. Knepley 
19480ec8681fSMatthew G. Knepley /* Centroid_i = (\sum_n V_n Cn_i ) / V */
1949011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
19500ec8681fSMatthew G. Knepley {
1951793a2a13SMatthew G. Knepley   DMLabel         depth;
19520ec8681fSMatthew G. Knepley   PetscSection    coordSection;
19530ec8681fSMatthew G. Knepley   Vec             coordinates;
19540ec8681fSMatthew G. Knepley   PetscScalar    *coords = NULL;
195586623015SMatthew G. Knepley   PetscReal       vsum = 0.0, vtmp, coordsTmp[3*3];
1956a7df9edeSMatthew G. Knepley   const PetscInt *faces, *facesO;
1957793a2a13SMatthew G. Knepley   PetscBool       isHybrid = PETSC_FALSE;
1958793a2a13SMatthew G. Knepley   PetscInt        pEndInterior[4], cdepth, numFaces, f, coordSize, numCorners, p, d;
19590ec8681fSMatthew G. Knepley   PetscErrorCode  ierr;
19600ec8681fSMatthew G. Knepley 
19610ec8681fSMatthew G. Knepley   PetscFunctionBegin;
1962f6dae198SJed Brown   if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No support for dim %D > 3",dim);
1963793a2a13SMatthew G. Knepley   /* Must check for hybrid cells because prisms have a different orientation scheme */
1964793a2a13SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
1965793a2a13SMatthew G. Knepley   ierr = DMLabelGetValue(depth, cell, &cdepth);CHKERRQ(ierr);
19662b2dc32bSSatish Balay   ierr = DMPlexGetHybridBounds(dm, &pEndInterior[3], &pEndInterior[2], &pEndInterior[1], &pEndInterior[0]);CHKERRQ(ierr);
1967793a2a13SMatthew G. Knepley   if ((pEndInterior[cdepth]) >=0 && (cell >= pEndInterior[cdepth])) isHybrid = PETSC_TRUE;
1968793a2a13SMatthew G. Knepley 
19690ec8681fSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
197069d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
19710ec8681fSMatthew G. Knepley 
1972d9a81ebdSMatthew G. Knepley   if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0;
19730ec8681fSMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr);
19740ec8681fSMatthew G. Knepley   ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr);
1975a7df9edeSMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr);
19760ec8681fSMatthew G. Knepley   for (f = 0; f < numFaces; ++f) {
1977793a2a13SMatthew G. Knepley     PetscBool flip = isHybrid && f == 0 ? PETSC_TRUE : PETSC_FALSE; /* The first hybrid face is reversed */
1978793a2a13SMatthew G. Knepley 
1979011ea5d8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr);
19800ec8681fSMatthew G. Knepley     numCorners = coordSize/dim;
19810ec8681fSMatthew G. Knepley     switch (numCorners) {
19820ec8681fSMatthew G. Knepley     case 3:
19830ec8681fSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
19841ee9d5ecSMatthew G. Knepley         coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]);
19851ee9d5ecSMatthew G. Knepley         coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]);
19861ee9d5ecSMatthew G. Knepley         coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]);
19870ec8681fSMatthew G. Knepley       }
19880ec8681fSMatthew G. Knepley       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
1989793a2a13SMatthew G. Knepley       if (facesO[f] < 0 || flip) vtmp = -vtmp;
19900ec8681fSMatthew G. Knepley       vsum += vtmp;
19914f25033aSJed Brown       if (centroid) {           /* Centroid of OABC = (a+b+c)/4 */
19920ec8681fSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
19931ee9d5ecSMatthew G. Knepley           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
19940ec8681fSMatthew G. Knepley         }
19950ec8681fSMatthew G. Knepley       }
19960ec8681fSMatthew G. Knepley       break;
19970ec8681fSMatthew G. Knepley     case 4:
1998793a2a13SMatthew G. Knepley     {
1999793a2a13SMatthew G. Knepley       PetscInt fv[4] = {0, 1, 2, 3};
2000793a2a13SMatthew G. Knepley 
2001793a2a13SMatthew G. Knepley       /* Side faces for hybrid cells are are stored as tensor products */
2002793a2a13SMatthew G. Knepley       if (isHybrid && f > 1) {fv[2] = 3; fv[3] = 2;}
20030ec8681fSMatthew G. Knepley       /* DO FOR PYRAMID */
20040ec8681fSMatthew G. Knepley       /* First tet */
20050ec8681fSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
2006793a2a13SMatthew G. Knepley         coordsTmp[0*dim+d] = PetscRealPart(coords[fv[0]*dim+d]);
2007793a2a13SMatthew G. Knepley         coordsTmp[1*dim+d] = PetscRealPart(coords[fv[1]*dim+d]);
2008793a2a13SMatthew G. Knepley         coordsTmp[2*dim+d] = PetscRealPart(coords[fv[3]*dim+d]);
20090ec8681fSMatthew G. Knepley       }
20100ec8681fSMatthew G. Knepley       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
2011793a2a13SMatthew G. Knepley       if (facesO[f] < 0 || flip) vtmp = -vtmp;
20120ec8681fSMatthew G. Knepley       vsum += vtmp;
20130ec8681fSMatthew G. Knepley       if (centroid) {
20140ec8681fSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
20150ec8681fSMatthew G. Knepley           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
20160ec8681fSMatthew G. Knepley         }
20170ec8681fSMatthew G. Knepley       }
20180ec8681fSMatthew G. Knepley       /* Second tet */
20190ec8681fSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
2020793a2a13SMatthew G. Knepley         coordsTmp[0*dim+d] = PetscRealPart(coords[fv[1]*dim+d]);
2021793a2a13SMatthew G. Knepley         coordsTmp[1*dim+d] = PetscRealPart(coords[fv[2]*dim+d]);
2022793a2a13SMatthew G. Knepley         coordsTmp[2*dim+d] = PetscRealPart(coords[fv[3]*dim+d]);
20230ec8681fSMatthew G. Knepley       }
20240ec8681fSMatthew G. Knepley       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
2025793a2a13SMatthew G. Knepley       if (facesO[f] < 0 || flip) vtmp = -vtmp;
20260ec8681fSMatthew G. Knepley       vsum += vtmp;
20270ec8681fSMatthew G. Knepley       if (centroid) {
20280ec8681fSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
20290ec8681fSMatthew G. Knepley           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
20300ec8681fSMatthew G. Knepley         }
20310ec8681fSMatthew G. Knepley       }
20320ec8681fSMatthew G. Knepley       break;
2033793a2a13SMatthew G. Knepley     }
20340ec8681fSMatthew G. Knepley     default:
2035796f034aSJed Brown       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %D vertices", numCorners);
20360ec8681fSMatthew G. Knepley     }
20374f25033aSJed Brown     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr);
20380ec8681fSMatthew G. Knepley   }
20398763be8eSMatthew G. Knepley   if (vol)     *vol = PetscAbsReal(vsum);
20400ec8681fSMatthew G. Knepley   if (normal)   for (d = 0; d < dim; ++d) normal[d]    = 0.0;
2041d9a81ebdSMatthew G. Knepley   if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4);
20420ec8681fSMatthew G. Knepley   PetscFunctionReturn(0);
20430ec8681fSMatthew G. Knepley }
20440ec8681fSMatthew G. Knepley 
2045834e62ceSMatthew G. Knepley /*@C
2046834e62ceSMatthew G. Knepley   DMPlexComputeCellGeometryFVM - Compute the volume for a given cell
2047834e62ceSMatthew G. Knepley 
2048d083f849SBarry Smith   Collective on dm
2049834e62ceSMatthew G. Knepley 
2050834e62ceSMatthew G. Knepley   Input Arguments:
2051834e62ceSMatthew G. Knepley + dm   - the DM
2052834e62ceSMatthew G. Knepley - cell - the cell
2053834e62ceSMatthew G. Knepley 
2054834e62ceSMatthew G. Knepley   Output Arguments:
2055834e62ceSMatthew G. Knepley + volume   - the cell volume
2056cc08537eSMatthew G. Knepley . centroid - the cell centroid
2057cc08537eSMatthew G. Knepley - normal - the cell normal, if appropriate
2058834e62ceSMatthew G. Knepley 
2059834e62ceSMatthew G. Knepley   Level: advanced
2060834e62ceSMatthew G. Knepley 
2061834e62ceSMatthew G. Knepley   Fortran Notes:
2062834e62ceSMatthew G. Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
2063834e62ceSMatthew G. Knepley   include petsc.h90 in your code.
2064834e62ceSMatthew G. Knepley 
2065e8964c0aSStefano Zampini .seealso: DMGetCoordinateSection(), DMGetCoordinates()
2066834e62ceSMatthew G. Knepley @*/
2067cc08537eSMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
2068834e62ceSMatthew G. Knepley {
20690ec8681fSMatthew G. Knepley   PetscInt       depth, dim;
2070834e62ceSMatthew G. Knepley   PetscErrorCode ierr;
2071834e62ceSMatthew G. Knepley 
2072834e62ceSMatthew G. Knepley   PetscFunctionBegin;
2073834e62ceSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2074c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2075834e62ceSMatthew G. Knepley   if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated");
2076834e62ceSMatthew G. Knepley   /* We need to keep a pointer to the depth label */
2077c58f1c22SToby Isaac   ierr = DMGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr);
2078834e62ceSMatthew G. Knepley   /* Cone size is now the number of faces */
2079011ea5d8SMatthew G. Knepley   switch (depth) {
2080cc08537eSMatthew G. Knepley   case 1:
2081011ea5d8SMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
2082cc08537eSMatthew G. Knepley     break;
2083834e62ceSMatthew G. Knepley   case 2:
2084011ea5d8SMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
2085834e62ceSMatthew G. Knepley     break;
2086834e62ceSMatthew G. Knepley   case 3:
2087011ea5d8SMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
2088834e62ceSMatthew G. Knepley     break;
2089834e62ceSMatthew G. Knepley   default:
2090b81cf158SMatthew G. Knepley     SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D (depth %D) for element geometry computation", dim, depth);
2091834e62ceSMatthew G. Knepley   }
2092834e62ceSMatthew G. Knepley   PetscFunctionReturn(0);
2093834e62ceSMatthew G. Knepley }
2094113c68e6SMatthew G. Knepley 
2095c501906fSMatthew G. Knepley /*@
2096c501906fSMatthew G. Knepley   DMPlexComputeGeometryFEM - Precompute cell geometry for the entire mesh
2097c501906fSMatthew G. Knepley 
2098c501906fSMatthew G. Knepley   Collective on dm
2099c501906fSMatthew G. Knepley 
2100c501906fSMatthew G. Knepley   Input Parameter:
2101c501906fSMatthew G. Knepley . dm - The DMPlex
2102c501906fSMatthew G. Knepley 
2103c501906fSMatthew G. Knepley   Output Parameter:
2104c501906fSMatthew G. Knepley . cellgeom - A vector with the cell geometry data for each cell
2105c501906fSMatthew G. Knepley 
2106c501906fSMatthew G. Knepley   Level: beginner
2107c501906fSMatthew G. Knepley 
2108c501906fSMatthew G. Knepley @*/
2109c0d900a5SMatthew G. Knepley PetscErrorCode DMPlexComputeGeometryFEM(DM dm, Vec *cellgeom)
2110c0d900a5SMatthew G. Knepley {
2111c0d900a5SMatthew G. Knepley   DM             dmCell;
2112c0d900a5SMatthew G. Knepley   Vec            coordinates;
2113c0d900a5SMatthew G. Knepley   PetscSection   coordSection, sectionCell;
2114c0d900a5SMatthew G. Knepley   PetscScalar   *cgeom;
2115c0d900a5SMatthew G. Knepley   PetscInt       cStart, cEnd, cMax, c;
2116c0d900a5SMatthew G. Knepley   PetscErrorCode ierr;
2117c0d900a5SMatthew G. Knepley 
2118c0d900a5SMatthew G. Knepley   PetscFunctionBegin;
2119c0d900a5SMatthew G. Knepley   ierr = DMClone(dm, &dmCell);CHKERRQ(ierr);
2120c0d900a5SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2121c0d900a5SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2122c0d900a5SMatthew G. Knepley   ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr);
2123c0d900a5SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr);
2124c0d900a5SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionCell);CHKERRQ(ierr);
2125c0d900a5SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2126c0d900a5SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
2127c0d900a5SMatthew G. Knepley   cEnd = cMax < 0 ? cEnd : cMax;
2128c0d900a5SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr);
2129c0d900a5SMatthew G. Knepley   /* TODO This needs to be multiplied by Nq for non-affine */
2130cf0b7c11SKarl Rupp   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFEGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);}
2131c0d900a5SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr);
213292fd8e1eSJed Brown   ierr = DMSetLocalSection(dmCell, sectionCell);CHKERRQ(ierr);
2133c0d900a5SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionCell);CHKERRQ(ierr);
2134c0d900a5SMatthew G. Knepley   ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr);
2135c0d900a5SMatthew G. Knepley   ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr);
2136c0d900a5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
2137cf0b7c11SKarl Rupp     PetscFEGeom *cg;
2138c0d900a5SMatthew G. Knepley 
2139c0d900a5SMatthew G. Knepley     ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
2140580bdb30SBarry Smith     ierr = PetscArrayzero(cg, 1);CHKERRQ(ierr);
2141cf0b7c11SKarl Rupp     ierr = DMPlexComputeCellGeometryFEM(dmCell, c, NULL, cg->v, cg->J, cg->invJ, cg->detJ);CHKERRQ(ierr);
2142087ef6b2SMatthew G. Knepley     if (*cg->detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D", (double) *cg->detJ, c);
2143c0d900a5SMatthew G. Knepley   }
2144c0d900a5SMatthew G. Knepley   PetscFunctionReturn(0);
2145c0d900a5SMatthew G. Knepley }
2146c0d900a5SMatthew G. Knepley 
2147891a9168SMatthew G. Knepley /*@
2148891a9168SMatthew G. Knepley   DMPlexComputeGeometryFVM - Computes the cell and face geometry for a finite volume method
2149891a9168SMatthew G. Knepley 
2150891a9168SMatthew G. Knepley   Input Parameter:
2151891a9168SMatthew G. Knepley . dm - The DM
2152891a9168SMatthew G. Knepley 
2153891a9168SMatthew G. Knepley   Output Parameters:
2154891a9168SMatthew G. Knepley + cellgeom - A Vec of PetscFVCellGeom data
2155a2b725a8SWilliam Gropp - facegeom - A Vec of PetscFVFaceGeom data
2156891a9168SMatthew G. Knepley 
2157891a9168SMatthew G. Knepley   Level: developer
2158891a9168SMatthew G. Knepley 
2159891a9168SMatthew G. Knepley .seealso: PetscFVFaceGeom, PetscFVCellGeom, DMPlexComputeGeometryFEM()
2160891a9168SMatthew G. Knepley @*/
2161113c68e6SMatthew G. Knepley PetscErrorCode DMPlexComputeGeometryFVM(DM dm, Vec *cellgeom, Vec *facegeom)
2162113c68e6SMatthew G. Knepley {
2163113c68e6SMatthew G. Knepley   DM             dmFace, dmCell;
2164113c68e6SMatthew G. Knepley   DMLabel        ghostLabel;
2165113c68e6SMatthew G. Knepley   PetscSection   sectionFace, sectionCell;
2166113c68e6SMatthew G. Knepley   PetscSection   coordSection;
2167113c68e6SMatthew G. Knepley   Vec            coordinates;
2168113c68e6SMatthew G. Knepley   PetscScalar   *fgeom, *cgeom;
2169113c68e6SMatthew G. Knepley   PetscReal      minradius, gminradius;
2170113c68e6SMatthew G. Knepley   PetscInt       dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f;
2171113c68e6SMatthew G. Knepley   PetscErrorCode ierr;
2172113c68e6SMatthew G. Knepley 
2173113c68e6SMatthew G. Knepley   PetscFunctionBegin;
2174113c68e6SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2175113c68e6SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2176113c68e6SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2177113c68e6SMatthew G. Knepley   /* Make cell centroids and volumes */
2178113c68e6SMatthew G. Knepley   ierr = DMClone(dm, &dmCell);CHKERRQ(ierr);
2179113c68e6SMatthew G. Knepley   ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr);
2180113c68e6SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr);
2181113c68e6SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionCell);CHKERRQ(ierr);
2182113c68e6SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2183485ad865SMatthew G. Knepley   ierr = DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);CHKERRQ(ierr);
2184113c68e6SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr);
21859e5edeeeSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVCellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);}
2186113c68e6SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr);
218792fd8e1eSJed Brown   ierr = DMSetLocalSection(dmCell, sectionCell);CHKERRQ(ierr);
2188113c68e6SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionCell);CHKERRQ(ierr);
2189113c68e6SMatthew G. Knepley   ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr);
2190485ad865SMatthew G. Knepley   if (cEndInterior < 0) cEndInterior = cEnd;
2191113c68e6SMatthew G. Knepley   ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr);
2192113c68e6SMatthew G. Knepley   for (c = cStart; c < cEndInterior; ++c) {
2193113c68e6SMatthew G. Knepley     PetscFVCellGeom *cg;
2194113c68e6SMatthew G. Knepley 
2195113c68e6SMatthew G. Knepley     ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
2196580bdb30SBarry Smith     ierr = PetscArrayzero(cg, 1);CHKERRQ(ierr);
2197113c68e6SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr);
2198113c68e6SMatthew G. Knepley   }
2199113c68e6SMatthew G. Knepley   /* Compute face normals and minimum cell radius */
2200113c68e6SMatthew G. Knepley   ierr = DMClone(dm, &dmFace);CHKERRQ(ierr);
2201113c68e6SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionFace);CHKERRQ(ierr);
2202113c68e6SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
2203113c68e6SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr);
22049e5edeeeSMatthew G. Knepley   for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVFaceGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);}
2205113c68e6SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr);
220692fd8e1eSJed Brown   ierr = DMSetLocalSection(dmFace, sectionFace);CHKERRQ(ierr);
2207113c68e6SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionFace);CHKERRQ(ierr);
2208113c68e6SMatthew G. Knepley   ierr = DMCreateLocalVector(dmFace, facegeom);CHKERRQ(ierr);
2209113c68e6SMatthew G. Knepley   ierr = VecGetArray(*facegeom, &fgeom);CHKERRQ(ierr);
2210c58f1c22SToby Isaac   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
2211113c68e6SMatthew G. Knepley   minradius = PETSC_MAX_REAL;
2212113c68e6SMatthew G. Knepley   for (f = fStart; f < fEnd; ++f) {
2213113c68e6SMatthew G. Knepley     PetscFVFaceGeom *fg;
2214113c68e6SMatthew G. Knepley     PetscReal        area;
221550d63984SToby Isaac     PetscInt         ghost = -1, d, numChildren;
2216113c68e6SMatthew G. Knepley 
22179ac3fadcSMatthew G. Knepley     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);}
221850d63984SToby Isaac     ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr);
221950d63984SToby Isaac     if (ghost >= 0 || numChildren) continue;
2220113c68e6SMatthew G. Knepley     ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr);
2221113c68e6SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr);
2222113c68e6SMatthew G. Knepley     for (d = 0; d < dim; ++d) fg->normal[d] *= area;
2223113c68e6SMatthew G. Knepley     /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */
2224113c68e6SMatthew G. Knepley     {
2225113c68e6SMatthew G. Knepley       PetscFVCellGeom *cL, *cR;
222606348e87SToby Isaac       PetscInt         ncells;
2227113c68e6SMatthew G. Knepley       const PetscInt  *cells;
2228113c68e6SMatthew G. Knepley       PetscReal       *lcentroid, *rcentroid;
22290453c0cdSMatthew G. Knepley       PetscReal        l[3], r[3], v[3];
2230113c68e6SMatthew G. Knepley 
2231113c68e6SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr);
223206348e87SToby Isaac       ierr = DMPlexGetSupportSize(dm, f, &ncells);CHKERRQ(ierr);
2233113c68e6SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr);
2234113c68e6SMatthew G. Knepley       lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid;
223506348e87SToby Isaac       if (ncells > 1) {
223606348e87SToby Isaac         ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr);
2237113c68e6SMatthew G. Knepley         rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid;
223806348e87SToby Isaac       }
223906348e87SToby Isaac       else {
224006348e87SToby Isaac         rcentroid = fg->centroid;
224106348e87SToby Isaac       }
22422e17dfb7SMatthew G. Knepley       ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, lcentroid, l);CHKERRQ(ierr);
22432e17dfb7SMatthew G. Knepley       ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, rcentroid, r);CHKERRQ(ierr);
22440453c0cdSMatthew G. Knepley       DMPlex_WaxpyD_Internal(dim, -1, l, r, v);
2245113c68e6SMatthew G. Knepley       if (DMPlex_DotRealD_Internal(dim, fg->normal, v) < 0) {
2246113c68e6SMatthew G. Knepley         for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d];
2247113c68e6SMatthew G. Knepley       }
2248113c68e6SMatthew G. Knepley       if (DMPlex_DotRealD_Internal(dim, fg->normal, v) <= 0) {
2249113c68e6SMatthew 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]);
2250113c68e6SMatthew 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]);
2251113c68e6SMatthew G. Knepley         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f);
2252113c68e6SMatthew G. Knepley       }
2253113c68e6SMatthew G. Knepley       if (cells[0] < cEndInterior) {
2254113c68e6SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cL->centroid, v);
2255113c68e6SMatthew G. Knepley         minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v));
2256113c68e6SMatthew G. Knepley       }
225706348e87SToby Isaac       if (ncells > 1 && cells[1] < cEndInterior) {
2258113c68e6SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cR->centroid, v);
2259113c68e6SMatthew G. Knepley         minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v));
2260113c68e6SMatthew G. Knepley       }
2261113c68e6SMatthew G. Knepley     }
2262113c68e6SMatthew G. Knepley   }
2263b2566f29SBarry Smith   ierr = MPIU_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
2264113c68e6SMatthew G. Knepley   ierr = DMPlexSetMinRadius(dm, gminradius);CHKERRQ(ierr);
2265113c68e6SMatthew G. Knepley   /* Compute centroids of ghost cells */
2266113c68e6SMatthew G. Knepley   for (c = cEndInterior; c < cEnd; ++c) {
2267113c68e6SMatthew G. Knepley     PetscFVFaceGeom *fg;
2268113c68e6SMatthew G. Knepley     const PetscInt  *cone,    *support;
2269113c68e6SMatthew G. Knepley     PetscInt         coneSize, supportSize, s;
2270113c68e6SMatthew G. Knepley 
2271113c68e6SMatthew G. Knepley     ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr);
2272113c68e6SMatthew G. Knepley     if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize);
2273113c68e6SMatthew G. Knepley     ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr);
2274113c68e6SMatthew G. Knepley     ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr);
227550d63984SToby Isaac     if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 2", cone[0], supportSize);
2276113c68e6SMatthew G. Knepley     ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr);
2277113c68e6SMatthew G. Knepley     ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr);
2278113c68e6SMatthew G. Knepley     for (s = 0; s < 2; ++s) {
2279113c68e6SMatthew G. Knepley       /* Reflect ghost centroid across plane of face */
2280113c68e6SMatthew G. Knepley       if (support[s] == c) {
2281640bce14SSatish Balay         PetscFVCellGeom       *ci;
2282113c68e6SMatthew G. Knepley         PetscFVCellGeom       *cg;
2283113c68e6SMatthew G. Knepley         PetscReal              c2f[3], a;
2284113c68e6SMatthew G. Knepley 
2285113c68e6SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr);
2286113c68e6SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */
2287113c68e6SMatthew G. Knepley         a    = DMPlex_DotRealD_Internal(dim, c2f, fg->normal)/DMPlex_DotRealD_Internal(dim, fg->normal, fg->normal);
2288113c68e6SMatthew G. Knepley         ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr);
2289113c68e6SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, 2*a, fg->normal, ci->centroid, cg->centroid);
2290113c68e6SMatthew G. Knepley         cg->volume = ci->volume;
2291113c68e6SMatthew G. Knepley       }
2292113c68e6SMatthew G. Knepley     }
2293113c68e6SMatthew G. Knepley   }
2294113c68e6SMatthew G. Knepley   ierr = VecRestoreArray(*facegeom, &fgeom);CHKERRQ(ierr);
2295113c68e6SMatthew G. Knepley   ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr);
2296113c68e6SMatthew G. Knepley   ierr = DMDestroy(&dmCell);CHKERRQ(ierr);
2297113c68e6SMatthew G. Knepley   ierr = DMDestroy(&dmFace);CHKERRQ(ierr);
2298113c68e6SMatthew G. Knepley   PetscFunctionReturn(0);
2299113c68e6SMatthew G. Knepley }
2300113c68e6SMatthew G. Knepley 
2301113c68e6SMatthew G. Knepley /*@C
2302113c68e6SMatthew G. Knepley   DMPlexGetMinRadius - Returns the minimum distance from any cell centroid to a face
2303113c68e6SMatthew G. Knepley 
2304113c68e6SMatthew G. Knepley   Not collective
2305113c68e6SMatthew G. Knepley 
2306113c68e6SMatthew G. Knepley   Input Argument:
2307113c68e6SMatthew G. Knepley . dm - the DM
2308113c68e6SMatthew G. Knepley 
2309113c68e6SMatthew G. Knepley   Output Argument:
2310113c68e6SMatthew G. Knepley . minradius - the minium cell radius
2311113c68e6SMatthew G. Knepley 
2312113c68e6SMatthew G. Knepley   Level: developer
2313113c68e6SMatthew G. Knepley 
2314113c68e6SMatthew G. Knepley .seealso: DMGetCoordinates()
2315113c68e6SMatthew G. Knepley @*/
2316113c68e6SMatthew G. Knepley PetscErrorCode DMPlexGetMinRadius(DM dm, PetscReal *minradius)
2317113c68e6SMatthew G. Knepley {
2318113c68e6SMatthew G. Knepley   PetscFunctionBegin;
2319113c68e6SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2320113c68e6SMatthew G. Knepley   PetscValidPointer(minradius,2);
2321113c68e6SMatthew G. Knepley   *minradius = ((DM_Plex*) dm->data)->minradius;
2322113c68e6SMatthew G. Knepley   PetscFunctionReturn(0);
2323113c68e6SMatthew G. Knepley }
2324113c68e6SMatthew G. Knepley 
2325113c68e6SMatthew G. Knepley /*@C
2326113c68e6SMatthew G. Knepley   DMPlexSetMinRadius - Sets the minimum distance from the cell centroid to a face
2327113c68e6SMatthew G. Knepley 
2328113c68e6SMatthew G. Knepley   Logically collective
2329113c68e6SMatthew G. Knepley 
2330113c68e6SMatthew G. Knepley   Input Arguments:
2331113c68e6SMatthew G. Knepley + dm - the DM
2332113c68e6SMatthew G. Knepley - minradius - the minium cell radius
2333113c68e6SMatthew G. Knepley 
2334113c68e6SMatthew G. Knepley   Level: developer
2335113c68e6SMatthew G. Knepley 
2336113c68e6SMatthew G. Knepley .seealso: DMSetCoordinates()
2337113c68e6SMatthew G. Knepley @*/
2338113c68e6SMatthew G. Knepley PetscErrorCode DMPlexSetMinRadius(DM dm, PetscReal minradius)
2339113c68e6SMatthew G. Knepley {
2340113c68e6SMatthew G. Knepley   PetscFunctionBegin;
2341113c68e6SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2342113c68e6SMatthew G. Knepley   ((DM_Plex*) dm->data)->minradius = minradius;
2343113c68e6SMatthew G. Knepley   PetscFunctionReturn(0);
2344113c68e6SMatthew G. Knepley }
2345856ac710SMatthew G. Knepley 
2346856ac710SMatthew G. Knepley static PetscErrorCode BuildGradientReconstruction_Internal(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
2347856ac710SMatthew G. Knepley {
2348856ac710SMatthew G. Knepley   DMLabel        ghostLabel;
2349856ac710SMatthew G. Knepley   PetscScalar   *dx, *grad, **gref;
2350856ac710SMatthew G. Knepley   PetscInt       dim, cStart, cEnd, c, cEndInterior, maxNumFaces;
2351856ac710SMatthew G. Knepley   PetscErrorCode ierr;
2352856ac710SMatthew G. Knepley 
2353856ac710SMatthew G. Knepley   PetscFunctionBegin;
2354856ac710SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2355856ac710SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2356485ad865SMatthew G. Knepley   ierr = DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);CHKERRQ(ierr);
2357856ac710SMatthew G. Knepley   ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr);
2358856ac710SMatthew G. Knepley   ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr);
2359c58f1c22SToby Isaac   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
2360856ac710SMatthew G. Knepley   ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr);
2361856ac710SMatthew G. Knepley   for (c = cStart; c < cEndInterior; c++) {
2362856ac710SMatthew G. Knepley     const PetscInt        *faces;
2363856ac710SMatthew G. Knepley     PetscInt               numFaces, usedFaces, f, d;
2364640bce14SSatish Balay     PetscFVCellGeom        *cg;
2365856ac710SMatthew G. Knepley     PetscBool              boundary;
2366856ac710SMatthew G. Knepley     PetscInt               ghost;
2367856ac710SMatthew G. Knepley 
2368856ac710SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
2369856ac710SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr);
2370856ac710SMatthew G. Knepley     ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr);
2371856ac710SMatthew 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);
2372856ac710SMatthew G. Knepley     for (f = 0, usedFaces = 0; f < numFaces; ++f) {
2373640bce14SSatish Balay       PetscFVCellGeom       *cg1;
2374856ac710SMatthew G. Knepley       PetscFVFaceGeom       *fg;
2375856ac710SMatthew G. Knepley       const PetscInt        *fcells;
2376856ac710SMatthew G. Knepley       PetscInt               ncell, side;
2377856ac710SMatthew G. Knepley 
2378856ac710SMatthew G. Knepley       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
2379a6ba4734SToby Isaac       ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
2380856ac710SMatthew G. Knepley       if ((ghost >= 0) || boundary) continue;
2381856ac710SMatthew G. Knepley       ierr  = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr);
2382856ac710SMatthew G. Knepley       side  = (c != fcells[0]); /* c is on left=0 or right=1 of face */
2383856ac710SMatthew G. Knepley       ncell = fcells[!side];    /* the neighbor */
2384856ac710SMatthew G. Knepley       ierr  = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr);
2385856ac710SMatthew G. Knepley       ierr  = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr);
2386856ac710SMatthew G. Knepley       for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d];
2387856ac710SMatthew G. Knepley       gref[usedFaces++] = fg->grad[side];  /* Gradient reconstruction term will go here */
2388856ac710SMatthew G. Knepley     }
2389856ac710SMatthew G. Knepley     if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?");
2390856ac710SMatthew G. Knepley     ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr);
2391856ac710SMatthew G. Knepley     for (f = 0, usedFaces = 0; f < numFaces; ++f) {
2392856ac710SMatthew G. Knepley       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
2393a6ba4734SToby Isaac       ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
2394856ac710SMatthew G. Knepley       if ((ghost >= 0) || boundary) continue;
2395856ac710SMatthew G. Knepley       for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d];
2396856ac710SMatthew G. Knepley       ++usedFaces;
2397856ac710SMatthew G. Knepley     }
2398856ac710SMatthew G. Knepley   }
2399856ac710SMatthew G. Knepley   ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr);
2400856ac710SMatthew G. Knepley   PetscFunctionReturn(0);
2401856ac710SMatthew G. Knepley }
2402856ac710SMatthew G. Knepley 
2403b81db932SToby Isaac static PetscErrorCode BuildGradientReconstruction_Internal_Tree(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
2404b81db932SToby Isaac {
2405b81db932SToby Isaac   DMLabel        ghostLabel;
2406b81db932SToby Isaac   PetscScalar   *dx, *grad, **gref;
2407b81db932SToby Isaac   PetscInt       dim, cStart, cEnd, c, cEndInterior, fStart, fEnd, f, nStart, nEnd, maxNumFaces = 0;
2408b81db932SToby Isaac   PetscSection   neighSec;
2409b81db932SToby Isaac   PetscInt     (*neighbors)[2];
2410b81db932SToby Isaac   PetscInt      *counter;
2411b81db932SToby Isaac   PetscErrorCode ierr;
2412b81db932SToby Isaac 
2413b81db932SToby Isaac   PetscFunctionBegin;
2414b81db932SToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2415b81db932SToby Isaac   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2416485ad865SMatthew G. Knepley   ierr = DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);CHKERRQ(ierr);
2417485ad865SMatthew G. Knepley   if (cEndInterior < 0) cEndInterior = cEnd;
2418b81db932SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&neighSec);CHKERRQ(ierr);
2419b81db932SToby Isaac   ierr = PetscSectionSetChart(neighSec,cStart,cEndInterior);CHKERRQ(ierr);
2420b81db932SToby Isaac   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
2421c58f1c22SToby Isaac   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
2422b81db932SToby Isaac   for (f = fStart; f < fEnd; f++) {
2423b81db932SToby Isaac     const PetscInt        *fcells;
2424b81db932SToby Isaac     PetscBool              boundary;
24255bc680faSToby Isaac     PetscInt               ghost = -1;
2426b81db932SToby Isaac     PetscInt               numChildren, numCells, c;
2427b81db932SToby Isaac 
242806348e87SToby Isaac     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);}
2429a6ba4734SToby Isaac     ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr);
2430b81db932SToby Isaac     ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr);
2431b81db932SToby Isaac     if ((ghost >= 0) || boundary || numChildren) continue;
2432b81db932SToby Isaac     ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr);
243306348e87SToby Isaac     if (numCells == 2) {
2434b81db932SToby Isaac       ierr = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr);
2435b81db932SToby Isaac       for (c = 0; c < 2; c++) {
2436b81db932SToby Isaac         PetscInt cell = fcells[c];
2437b81db932SToby Isaac 
2438e6885bbbSToby Isaac         if (cell >= cStart && cell < cEndInterior) {
2439b81db932SToby Isaac           ierr = PetscSectionAddDof(neighSec,cell,1);CHKERRQ(ierr);
2440b81db932SToby Isaac         }
2441b81db932SToby Isaac       }
2442b81db932SToby Isaac     }
244306348e87SToby Isaac   }
2444b81db932SToby Isaac   ierr = PetscSectionSetUp(neighSec);CHKERRQ(ierr);
2445b81db932SToby Isaac   ierr = PetscSectionGetMaxDof(neighSec,&maxNumFaces);CHKERRQ(ierr);
2446b81db932SToby Isaac   ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr);
2447b81db932SToby Isaac   nStart = 0;
2448b81db932SToby Isaac   ierr = PetscSectionGetStorageSize(neighSec,&nEnd);CHKERRQ(ierr);
2449b81db932SToby Isaac   ierr = PetscMalloc1((nEnd-nStart),&neighbors);CHKERRQ(ierr);
2450b81db932SToby Isaac   ierr = PetscCalloc1((cEndInterior-cStart),&counter);CHKERRQ(ierr);
2451b81db932SToby Isaac   for (f = fStart; f < fEnd; f++) {
2452b81db932SToby Isaac     const PetscInt        *fcells;
2453b81db932SToby Isaac     PetscBool              boundary;
24545bc680faSToby Isaac     PetscInt               ghost = -1;
2455b81db932SToby Isaac     PetscInt               numChildren, numCells, c;
2456b81db932SToby Isaac 
245706348e87SToby Isaac     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);}
2458a6ba4734SToby Isaac     ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr);
2459b81db932SToby Isaac     ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr);
2460b81db932SToby Isaac     if ((ghost >= 0) || boundary || numChildren) continue;
2461b81db932SToby Isaac     ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr);
246206348e87SToby Isaac     if (numCells == 2) {
2463b81db932SToby Isaac       ierr  = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr);
2464b81db932SToby Isaac       for (c = 0; c < 2; c++) {
2465b81db932SToby Isaac         PetscInt cell = fcells[c], off;
2466b81db932SToby Isaac 
2467e6885bbbSToby Isaac         if (cell >= cStart && cell < cEndInterior) {
2468b81db932SToby Isaac           ierr = PetscSectionGetOffset(neighSec,cell,&off);CHKERRQ(ierr);
2469b81db932SToby Isaac           off += counter[cell - cStart]++;
2470b81db932SToby Isaac           neighbors[off][0] = f;
2471b81db932SToby Isaac           neighbors[off][1] = fcells[1 - c];
2472b81db932SToby Isaac         }
2473b81db932SToby Isaac       }
2474b81db932SToby Isaac     }
247506348e87SToby Isaac   }
2476b81db932SToby Isaac   ierr = PetscFree(counter);CHKERRQ(ierr);
2477b81db932SToby Isaac   ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr);
2478b81db932SToby Isaac   for (c = cStart; c < cEndInterior; c++) {
2479317218b9SToby Isaac     PetscInt               numFaces, f, d, off, ghost = -1;
2480640bce14SSatish Balay     PetscFVCellGeom        *cg;
2481b81db932SToby Isaac 
2482b81db932SToby Isaac     ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
2483b81db932SToby Isaac     ierr = PetscSectionGetDof(neighSec, c, &numFaces);CHKERRQ(ierr);
2484b81db932SToby Isaac     ierr = PetscSectionGetOffset(neighSec, c, &off);CHKERRQ(ierr);
2485317218b9SToby Isaac     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, c, &ghost);CHKERRQ(ierr);}
2486317218b9SToby 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);
2487b81db932SToby Isaac     for (f = 0; f < numFaces; ++f) {
2488640bce14SSatish Balay       PetscFVCellGeom       *cg1;
2489b81db932SToby Isaac       PetscFVFaceGeom       *fg;
2490b81db932SToby Isaac       const PetscInt        *fcells;
2491b81db932SToby Isaac       PetscInt               ncell, side, nface;
2492b81db932SToby Isaac 
2493b81db932SToby Isaac       nface = neighbors[off + f][0];
2494b81db932SToby Isaac       ncell = neighbors[off + f][1];
2495b81db932SToby Isaac       ierr  = DMPlexGetSupport(dm,nface,&fcells);CHKERRQ(ierr);
2496b81db932SToby Isaac       side  = (c != fcells[0]);
2497b81db932SToby Isaac       ierr  = DMPlexPointLocalRef(dmFace, nface, fgeom, &fg);CHKERRQ(ierr);
2498b81db932SToby Isaac       ierr  = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr);
2499b81db932SToby Isaac       for (d = 0; d < dim; ++d) dx[f*dim+d] = cg1->centroid[d] - cg->centroid[d];
2500b81db932SToby Isaac       gref[f] = fg->grad[side];  /* Gradient reconstruction term will go here */
2501b81db932SToby Isaac     }
2502b81db932SToby Isaac     ierr = PetscFVComputeGradient(fvm, numFaces, dx, grad);CHKERRQ(ierr);
2503b81db932SToby Isaac     for (f = 0; f < numFaces; ++f) {
2504b81db932SToby Isaac       for (d = 0; d < dim; ++d) gref[f][d] = grad[f*dim+d];
2505b81db932SToby Isaac     }
2506b81db932SToby Isaac   }
2507b81db932SToby Isaac   ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr);
25085fe94518SToby Isaac   ierr = PetscSectionDestroy(&neighSec);CHKERRQ(ierr);
2509b81db932SToby Isaac   ierr = PetscFree(neighbors);CHKERRQ(ierr);
2510b81db932SToby Isaac   PetscFunctionReturn(0);
2511b81db932SToby Isaac }
2512b81db932SToby Isaac 
2513856ac710SMatthew G. Knepley /*@
2514856ac710SMatthew G. Knepley   DMPlexComputeGradientFVM - Compute geometric factors for gradient reconstruction, which are stored in the geometry data, and compute layout for gradient data
2515856ac710SMatthew G. Knepley 
2516d083f849SBarry Smith   Collective on dm
2517856ac710SMatthew G. Knepley 
2518856ac710SMatthew G. Knepley   Input Arguments:
2519856ac710SMatthew G. Knepley + dm  - The DM
2520856ac710SMatthew G. Knepley . fvm - The PetscFV
25218f9f38e3SMatthew G. Knepley . faceGeometry - The face geometry from DMPlexComputeFaceGeometryFVM()
25228f9f38e3SMatthew G. Knepley - cellGeometry - The face geometry from DMPlexComputeCellGeometryFVM()
2523856ac710SMatthew G. Knepley 
2524856ac710SMatthew G. Knepley   Output Parameters:
2525856ac710SMatthew G. Knepley + faceGeometry - The geometric factors for gradient calculation are inserted
2526856ac710SMatthew G. Knepley - dmGrad - The DM describing the layout of gradient data
2527856ac710SMatthew G. Knepley 
2528856ac710SMatthew G. Knepley   Level: developer
2529856ac710SMatthew G. Knepley 
2530856ac710SMatthew G. Knepley .seealso: DMPlexGetFaceGeometryFVM(), DMPlexGetCellGeometryFVM()
2531856ac710SMatthew G. Knepley @*/
2532856ac710SMatthew G. Knepley PetscErrorCode DMPlexComputeGradientFVM(DM dm, PetscFV fvm, Vec faceGeometry, Vec cellGeometry, DM *dmGrad)
2533856ac710SMatthew G. Knepley {
2534856ac710SMatthew G. Knepley   DM             dmFace, dmCell;
2535856ac710SMatthew G. Knepley   PetscScalar   *fgeom, *cgeom;
2536b81db932SToby Isaac   PetscSection   sectionGrad, parentSection;
2537856ac710SMatthew G. Knepley   PetscInt       dim, pdim, cStart, cEnd, cEndInterior, c;
2538856ac710SMatthew G. Knepley   PetscErrorCode ierr;
2539856ac710SMatthew G. Knepley 
2540856ac710SMatthew G. Knepley   PetscFunctionBegin;
2541856ac710SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2542856ac710SMatthew G. Knepley   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
2543856ac710SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2544485ad865SMatthew G. Knepley   ierr = DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);CHKERRQ(ierr);
2545856ac710SMatthew G. Knepley   /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */
2546856ac710SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
2547856ac710SMatthew G. Knepley   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
2548856ac710SMatthew G. Knepley   ierr = VecGetArray(faceGeometry, &fgeom);CHKERRQ(ierr);
2549856ac710SMatthew G. Knepley   ierr = VecGetArray(cellGeometry, &cgeom);CHKERRQ(ierr);
2550b81db932SToby Isaac   ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
2551b81db932SToby Isaac   if (!parentSection) {
2552856ac710SMatthew G. Knepley     ierr = BuildGradientReconstruction_Internal(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr);
2553b5a3613cSMatthew G. Knepley   } else {
2554b81db932SToby Isaac     ierr = BuildGradientReconstruction_Internal_Tree(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr);
2555b81db932SToby Isaac   }
2556856ac710SMatthew G. Knepley   ierr = VecRestoreArray(faceGeometry, &fgeom);CHKERRQ(ierr);
2557856ac710SMatthew G. Knepley   ierr = VecRestoreArray(cellGeometry, &cgeom);CHKERRQ(ierr);
2558856ac710SMatthew G. Knepley   /* Create storage for gradients */
2559856ac710SMatthew G. Knepley   ierr = DMClone(dm, dmGrad);CHKERRQ(ierr);
2560856ac710SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionGrad);CHKERRQ(ierr);
2561856ac710SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr);
2562856ac710SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);}
2563856ac710SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr);
256492fd8e1eSJed Brown   ierr = DMSetLocalSection(*dmGrad, sectionGrad);CHKERRQ(ierr);
2565856ac710SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionGrad);CHKERRQ(ierr);
2566856ac710SMatthew G. Knepley   PetscFunctionReturn(0);
2567856ac710SMatthew G. Knepley }
2568b27d5b9eSToby Isaac 
2569c501906fSMatthew G. Knepley /*@
2570c501906fSMatthew G. Knepley   DMPlexGetDataFVM - Retrieve precomputed cell geometry
2571c501906fSMatthew G. Knepley 
2572d083f849SBarry Smith   Collective on dm
2573c501906fSMatthew G. Knepley 
2574c501906fSMatthew G. Knepley   Input Arguments:
2575c501906fSMatthew G. Knepley + dm  - The DM
2576c501906fSMatthew G. Knepley - fvm - The PetscFV
2577c501906fSMatthew G. Knepley 
2578c501906fSMatthew G. Knepley   Output Parameters:
2579c501906fSMatthew G. Knepley + cellGeometry - The cell geometry
2580c501906fSMatthew G. Knepley . faceGeometry - The face geometry
2581c501906fSMatthew G. Knepley - dmGrad       - The gradient matrices
2582c501906fSMatthew G. Knepley 
2583c501906fSMatthew G. Knepley   Level: developer
2584c501906fSMatthew G. Knepley 
2585c501906fSMatthew G. Knepley .seealso: DMPlexComputeGeometryFVM()
2586c501906fSMatthew G. Knepley @*/
2587b27d5b9eSToby Isaac PetscErrorCode DMPlexGetDataFVM(DM dm, PetscFV fv, Vec *cellgeom, Vec *facegeom, DM *gradDM)
2588b27d5b9eSToby Isaac {
2589b27d5b9eSToby Isaac   PetscObject    cellgeomobj, facegeomobj;
2590b27d5b9eSToby Isaac   PetscErrorCode ierr;
2591b27d5b9eSToby Isaac 
2592b27d5b9eSToby Isaac   PetscFunctionBegin;
2593b27d5b9eSToby Isaac   ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);CHKERRQ(ierr);
2594b27d5b9eSToby Isaac   if (!cellgeomobj) {
2595b27d5b9eSToby Isaac     Vec cellgeomInt, facegeomInt;
2596b27d5b9eSToby Isaac 
2597b27d5b9eSToby Isaac     ierr = DMPlexComputeGeometryFVM(dm, &cellgeomInt, &facegeomInt);CHKERRQ(ierr);
2598b27d5b9eSToby Isaac     ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_cellgeom_fvm",(PetscObject)cellgeomInt);CHKERRQ(ierr);
2599b27d5b9eSToby Isaac     ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_facegeom_fvm",(PetscObject)facegeomInt);CHKERRQ(ierr);
2600b27d5b9eSToby Isaac     ierr = VecDestroy(&cellgeomInt);CHKERRQ(ierr);
2601b27d5b9eSToby Isaac     ierr = VecDestroy(&facegeomInt);CHKERRQ(ierr);
2602b27d5b9eSToby Isaac     ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);CHKERRQ(ierr);
2603b27d5b9eSToby Isaac   }
2604b27d5b9eSToby Isaac   ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_facegeom_fvm", &facegeomobj);CHKERRQ(ierr);
2605b27d5b9eSToby Isaac   if (cellgeom) *cellgeom = (Vec) cellgeomobj;
2606b27d5b9eSToby Isaac   if (facegeom) *facegeom = (Vec) facegeomobj;
2607b27d5b9eSToby Isaac   if (gradDM) {
2608b27d5b9eSToby Isaac     PetscObject gradobj;
2609b27d5b9eSToby Isaac     PetscBool   computeGradients;
2610b27d5b9eSToby Isaac 
2611b27d5b9eSToby Isaac     ierr = PetscFVGetComputeGradients(fv,&computeGradients);CHKERRQ(ierr);
2612b27d5b9eSToby Isaac     if (!computeGradients) {
2613b27d5b9eSToby Isaac       *gradDM = NULL;
2614b27d5b9eSToby Isaac       PetscFunctionReturn(0);
2615b27d5b9eSToby Isaac     }
2616b27d5b9eSToby Isaac     ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);CHKERRQ(ierr);
2617b27d5b9eSToby Isaac     if (!gradobj) {
2618b27d5b9eSToby Isaac       DM dmGradInt;
2619b27d5b9eSToby Isaac 
2620b27d5b9eSToby Isaac       ierr = DMPlexComputeGradientFVM(dm,fv,(Vec) facegeomobj,(Vec) cellgeomobj,&dmGradInt);CHKERRQ(ierr);
2621b27d5b9eSToby Isaac       ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_dmgrad_fvm", (PetscObject)dmGradInt);CHKERRQ(ierr);
2622b27d5b9eSToby Isaac       ierr = DMDestroy(&dmGradInt);CHKERRQ(ierr);
2623b27d5b9eSToby Isaac       ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);CHKERRQ(ierr);
2624b27d5b9eSToby Isaac     }
2625b27d5b9eSToby Isaac     *gradDM = (DM) gradobj;
2626b27d5b9eSToby Isaac   }
2627b27d5b9eSToby Isaac   PetscFunctionReturn(0);
2628b27d5b9eSToby Isaac }
2629d6143a4eSToby Isaac 
26309d150b73SToby Isaac static PetscErrorCode DMPlexCoordinatesToReference_NewtonUpdate(PetscInt dimC, PetscInt dimR, PetscScalar *J, PetscScalar *invJ, PetscScalar *work,  PetscReal *resNeg, PetscReal *guess)
26319d150b73SToby Isaac {
26329d150b73SToby Isaac   PetscInt l, m;
26339d150b73SToby Isaac 
2634cd345991SToby Isaac   PetscFunctionBeginHot;
26359d150b73SToby Isaac   if (dimC == dimR && dimR <= 3) {
26369d150b73SToby Isaac     /* invert Jacobian, multiply */
26379d150b73SToby Isaac     PetscScalar det, idet;
26389d150b73SToby Isaac 
26399d150b73SToby Isaac     switch (dimR) {
26409d150b73SToby Isaac     case 1:
26419d150b73SToby Isaac       invJ[0] = 1./ J[0];
26429d150b73SToby Isaac       break;
26439d150b73SToby Isaac     case 2:
26449d150b73SToby Isaac       det = J[0] * J[3] - J[1] * J[2];
26459d150b73SToby Isaac       idet = 1./det;
26469d150b73SToby Isaac       invJ[0] =  J[3] * idet;
26479d150b73SToby Isaac       invJ[1] = -J[1] * idet;
26489d150b73SToby Isaac       invJ[2] = -J[2] * idet;
26499d150b73SToby Isaac       invJ[3] =  J[0] * idet;
26509d150b73SToby Isaac       break;
26519d150b73SToby Isaac     case 3:
26529d150b73SToby Isaac       {
26539d150b73SToby Isaac         invJ[0] = J[4] * J[8] - J[5] * J[7];
26549d150b73SToby Isaac         invJ[1] = J[2] * J[7] - J[1] * J[8];
26559d150b73SToby Isaac         invJ[2] = J[1] * J[5] - J[2] * J[4];
26569d150b73SToby Isaac         det = invJ[0] * J[0] + invJ[1] * J[3] + invJ[2] * J[6];
26579d150b73SToby Isaac         idet = 1./det;
26589d150b73SToby Isaac         invJ[0] *= idet;
26599d150b73SToby Isaac         invJ[1] *= idet;
26609d150b73SToby Isaac         invJ[2] *= idet;
26619d150b73SToby Isaac         invJ[3]  = idet * (J[5] * J[6] - J[3] * J[8]);
26629d150b73SToby Isaac         invJ[4]  = idet * (J[0] * J[8] - J[2] * J[6]);
26639d150b73SToby Isaac         invJ[5]  = idet * (J[2] * J[3] - J[0] * J[5]);
26649d150b73SToby Isaac         invJ[6]  = idet * (J[3] * J[7] - J[4] * J[6]);
26659d150b73SToby Isaac         invJ[7]  = idet * (J[1] * J[6] - J[0] * J[7]);
26669d150b73SToby Isaac         invJ[8]  = idet * (J[0] * J[4] - J[1] * J[3]);
26679d150b73SToby Isaac       }
26689d150b73SToby Isaac       break;
26699d150b73SToby Isaac     }
26709d150b73SToby Isaac     for (l = 0; l < dimR; l++) {
26719d150b73SToby Isaac       for (m = 0; m < dimC; m++) {
2672c6e120d1SToby Isaac         guess[l] += PetscRealPart(invJ[l * dimC + m]) * resNeg[m];
26739d150b73SToby Isaac       }
26749d150b73SToby Isaac     }
26759d150b73SToby Isaac   } else {
26769d150b73SToby Isaac #if defined(PETSC_USE_COMPLEX)
26779d150b73SToby Isaac     char transpose = 'C';
26789d150b73SToby Isaac #else
26799d150b73SToby Isaac     char transpose = 'T';
26809d150b73SToby Isaac #endif
26819d150b73SToby Isaac     PetscBLASInt m = dimR;
26829d150b73SToby Isaac     PetscBLASInt n = dimC;
26839d150b73SToby Isaac     PetscBLASInt one = 1;
26849d150b73SToby Isaac     PetscBLASInt worksize = dimR * dimC, info;
26859d150b73SToby Isaac 
26869d150b73SToby Isaac     for (l = 0; l < dimC; l++) {invJ[l] = resNeg[l];}
26879d150b73SToby Isaac 
26889d150b73SToby Isaac     PetscStackCallBLAS("LAPACKgels",LAPACKgels_(&transpose,&m,&n,&one,J,&m,invJ,&n,work,&worksize, &info));
26899d150b73SToby Isaac     if (info != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELS");
26909d150b73SToby Isaac 
2691c6e120d1SToby Isaac     for (l = 0; l < dimR; l++) {guess[l] += PetscRealPart(invJ[l]);}
26929d150b73SToby Isaac   }
26939d150b73SToby Isaac   PetscFunctionReturn(0);
26949d150b73SToby Isaac }
26959d150b73SToby Isaac 
26969d150b73SToby Isaac static PetscErrorCode DMPlexCoordinatesToReference_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[], Vec coords, PetscInt dimC, PetscInt dimR)
26979d150b73SToby Isaac {
2698c0cbe899SToby Isaac   PetscInt       coordSize, i, j, k, l, m, maxIts = 7, numV = (1 << dimR);
26999d150b73SToby Isaac   PetscScalar    *coordsScalar = NULL;
27009d150b73SToby Isaac   PetscReal      *cellData, *cellCoords, *cellCoeffs, *extJ, *resNeg;
27019d150b73SToby Isaac   PetscScalar    *J, *invJ, *work;
27029d150b73SToby Isaac   PetscErrorCode ierr;
27039d150b73SToby Isaac 
27049d150b73SToby Isaac   PetscFunctionBegin;
27059d150b73SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
27069d150b73SToby Isaac   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
270713903a91SSatish Balay   if (coordSize < dimC * numV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expecting at least %D coordinates, got %D",dimC * (1 << dimR), coordSize);
270869291d52SBarry Smith   ierr = DMGetWorkArray(dm, 2 * coordSize + dimR + dimC, MPIU_REAL, &cellData);CHKERRQ(ierr);
270969291d52SBarry Smith   ierr = DMGetWorkArray(dm, 3 * dimR * dimC, MPIU_SCALAR, &J);CHKERRQ(ierr);
27109d150b73SToby Isaac   cellCoords = &cellData[0];
27119d150b73SToby Isaac   cellCoeffs = &cellData[coordSize];
27129d150b73SToby Isaac   extJ       = &cellData[2 * coordSize];
27139d150b73SToby Isaac   resNeg     = &cellData[2 * coordSize + dimR];
27149d150b73SToby Isaac   invJ       = &J[dimR * dimC];
27159d150b73SToby Isaac   work       = &J[2 * dimR * dimC];
27169d150b73SToby Isaac   if (dimR == 2) {
27179d150b73SToby Isaac     const PetscInt zToPlex[4] = {0, 1, 3, 2};
27189d150b73SToby Isaac 
27199d150b73SToby Isaac     for (i = 0; i < 4; i++) {
27209d150b73SToby Isaac       PetscInt plexI = zToPlex[i];
27219d150b73SToby Isaac 
27229d150b73SToby Isaac       for (j = 0; j < dimC; j++) {
27239d150b73SToby Isaac         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
27249d150b73SToby Isaac       }
27259d150b73SToby Isaac     }
27269d150b73SToby Isaac   } else if (dimR == 3) {
27279d150b73SToby Isaac     const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6};
27289d150b73SToby Isaac 
27299d150b73SToby Isaac     for (i = 0; i < 8; i++) {
27309d150b73SToby Isaac       PetscInt plexI = zToPlex[i];
27319d150b73SToby Isaac 
27329d150b73SToby Isaac       for (j = 0; j < dimC; j++) {
27339d150b73SToby Isaac         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
27349d150b73SToby Isaac       }
27359d150b73SToby Isaac     }
27369d150b73SToby Isaac   } else {
27379d150b73SToby Isaac     for (i = 0; i < coordSize; i++) {cellCoords[i] = PetscRealPart(coordsScalar[i]);}
27389d150b73SToby Isaac   }
27399d150b73SToby Isaac   /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */
27409d150b73SToby Isaac   for (i = 0; i < dimR; i++) {
27419d150b73SToby Isaac     PetscReal *swap;
27429d150b73SToby Isaac 
27439d150b73SToby Isaac     for (j = 0; j < (numV / 2); j++) {
27449d150b73SToby Isaac       for (k = 0; k < dimC; k++) {
27459d150b73SToby Isaac         cellCoeffs[dimC * j + k]                = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]);
27469d150b73SToby Isaac         cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]);
27479d150b73SToby Isaac       }
27489d150b73SToby Isaac     }
27499d150b73SToby Isaac 
27509d150b73SToby Isaac     if (i < dimR - 1) {
27519d150b73SToby Isaac       swap = cellCoeffs;
27529d150b73SToby Isaac       cellCoeffs = cellCoords;
27539d150b73SToby Isaac       cellCoords = swap;
27549d150b73SToby Isaac     }
27559d150b73SToby Isaac   }
2756580bdb30SBarry Smith   ierr = PetscArrayzero(refCoords,numPoints * dimR);CHKERRQ(ierr);
27579d150b73SToby Isaac   for (j = 0; j < numPoints; j++) {
27589d150b73SToby Isaac     for (i = 0; i < maxIts; i++) {
27599d150b73SToby Isaac       PetscReal *guess = &refCoords[dimR * j];
27609d150b73SToby Isaac 
27619d150b73SToby Isaac       /* compute -residual and Jacobian */
27629d150b73SToby Isaac       for (k = 0; k < dimC; k++) {resNeg[k] = realCoords[dimC * j + k];}
27639d150b73SToby Isaac       for (k = 0; k < dimC * dimR; k++) {J[k] = 0.;}
27649d150b73SToby Isaac       for (k = 0; k < numV; k++) {
27659d150b73SToby Isaac         PetscReal extCoord = 1.;
27669d150b73SToby Isaac         for (l = 0; l < dimR; l++) {
27679d150b73SToby Isaac           PetscReal coord = guess[l];
27689d150b73SToby Isaac           PetscInt  dep   = (k & (1 << l)) >> l;
27699d150b73SToby Isaac 
27709d150b73SToby Isaac           extCoord *= dep * coord + !dep;
27719d150b73SToby Isaac           extJ[l] = dep;
27729d150b73SToby Isaac 
27739d150b73SToby Isaac           for (m = 0; m < dimR; m++) {
27749d150b73SToby Isaac             PetscReal coord = guess[m];
27759d150b73SToby Isaac             PetscInt  dep   = ((k & (1 << m)) >> m) && (m != l);
27769d150b73SToby Isaac             PetscReal mult  = dep * coord + !dep;
27779d150b73SToby Isaac 
27789d150b73SToby Isaac             extJ[l] *= mult;
27799d150b73SToby Isaac           }
27809d150b73SToby Isaac         }
27819d150b73SToby Isaac         for (l = 0; l < dimC; l++) {
27829d150b73SToby Isaac           PetscReal coeff = cellCoeffs[dimC * k + l];
27839d150b73SToby Isaac 
27849d150b73SToby Isaac           resNeg[l] -= coeff * extCoord;
27859d150b73SToby Isaac           for (m = 0; m < dimR; m++) {
27869d150b73SToby Isaac             J[dimR * l + m] += coeff * extJ[m];
27879d150b73SToby Isaac           }
27889d150b73SToby Isaac         }
27899d150b73SToby Isaac       }
27900611203eSToby Isaac #if 0 && defined(PETSC_USE_DEBUG)
27910611203eSToby Isaac       {
27920611203eSToby Isaac         PetscReal maxAbs = 0.;
27930611203eSToby Isaac 
27940611203eSToby Isaac         for (l = 0; l < dimC; l++) {
27950611203eSToby Isaac           maxAbs = PetscMax(maxAbs,PetscAbsReal(resNeg[l]));
27960611203eSToby Isaac         }
2797087ef6b2SMatthew G. Knepley         ierr = PetscInfo4(dm,"cell %D, point %D, iter %D: res %g\n",cell,j,i,(double) maxAbs);CHKERRQ(ierr);
27980611203eSToby Isaac       }
27990611203eSToby Isaac #endif
28009d150b73SToby Isaac 
28019d150b73SToby Isaac       ierr = DMPlexCoordinatesToReference_NewtonUpdate(dimC,dimR,J,invJ,work,resNeg,guess);CHKERRQ(ierr);
28029d150b73SToby Isaac     }
28039d150b73SToby Isaac   }
280469291d52SBarry Smith   ierr = DMRestoreWorkArray(dm, 3 * dimR * dimC, MPIU_SCALAR, &J);CHKERRQ(ierr);
280569291d52SBarry Smith   ierr = DMRestoreWorkArray(dm, 2 * coordSize + dimR + dimC, MPIU_REAL, &cellData);CHKERRQ(ierr);
28069d150b73SToby Isaac   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
28079d150b73SToby Isaac   PetscFunctionReturn(0);
28089d150b73SToby Isaac }
28099d150b73SToby Isaac 
28109d150b73SToby Isaac static PetscErrorCode DMPlexReferenceToCoordinates_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[], Vec coords, PetscInt dimC, PetscInt dimR)
28119d150b73SToby Isaac {
28129d150b73SToby Isaac   PetscInt       coordSize, i, j, k, l, numV = (1 << dimR);
28139d150b73SToby Isaac   PetscScalar    *coordsScalar = NULL;
28149d150b73SToby Isaac   PetscReal      *cellData, *cellCoords, *cellCoeffs;
28159d150b73SToby Isaac   PetscErrorCode ierr;
28169d150b73SToby Isaac 
28179d150b73SToby Isaac   PetscFunctionBegin;
28189d150b73SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
28199d150b73SToby Isaac   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
282013903a91SSatish Balay   if (coordSize < dimC * numV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expecting at least %D coordinates, got %D",dimC * (1 << dimR), coordSize);
282169291d52SBarry Smith   ierr = DMGetWorkArray(dm, 2 * coordSize, MPIU_REAL, &cellData);CHKERRQ(ierr);
28229d150b73SToby Isaac   cellCoords = &cellData[0];
28239d150b73SToby Isaac   cellCoeffs = &cellData[coordSize];
28249d150b73SToby Isaac   if (dimR == 2) {
28259d150b73SToby Isaac     const PetscInt zToPlex[4] = {0, 1, 3, 2};
28269d150b73SToby Isaac 
28279d150b73SToby Isaac     for (i = 0; i < 4; i++) {
28289d150b73SToby Isaac       PetscInt plexI = zToPlex[i];
28299d150b73SToby Isaac 
28309d150b73SToby Isaac       for (j = 0; j < dimC; j++) {
28319d150b73SToby Isaac         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
28329d150b73SToby Isaac       }
28339d150b73SToby Isaac     }
28349d150b73SToby Isaac   } else if (dimR == 3) {
28359d150b73SToby Isaac     const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6};
28369d150b73SToby Isaac 
28379d150b73SToby Isaac     for (i = 0; i < 8; i++) {
28389d150b73SToby Isaac       PetscInt plexI = zToPlex[i];
28399d150b73SToby Isaac 
28409d150b73SToby Isaac       for (j = 0; j < dimC; j++) {
28419d150b73SToby Isaac         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
28429d150b73SToby Isaac       }
28439d150b73SToby Isaac     }
28449d150b73SToby Isaac   } else {
28459d150b73SToby Isaac     for (i = 0; i < coordSize; i++) {cellCoords[i] = PetscRealPart(coordsScalar[i]);}
28469d150b73SToby Isaac   }
28479d150b73SToby Isaac   /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */
28489d150b73SToby Isaac   for (i = 0; i < dimR; i++) {
28499d150b73SToby Isaac     PetscReal *swap;
28509d150b73SToby Isaac 
28519d150b73SToby Isaac     for (j = 0; j < (numV / 2); j++) {
28529d150b73SToby Isaac       for (k = 0; k < dimC; k++) {
28539d150b73SToby Isaac         cellCoeffs[dimC * j + k]                = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]);
28549d150b73SToby Isaac         cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]);
28559d150b73SToby Isaac       }
28569d150b73SToby Isaac     }
28579d150b73SToby Isaac 
28589d150b73SToby Isaac     if (i < dimR - 1) {
28599d150b73SToby Isaac       swap = cellCoeffs;
28609d150b73SToby Isaac       cellCoeffs = cellCoords;
28619d150b73SToby Isaac       cellCoords = swap;
28629d150b73SToby Isaac     }
28639d150b73SToby Isaac   }
2864580bdb30SBarry Smith   ierr = PetscArrayzero(realCoords,numPoints * dimC);CHKERRQ(ierr);
28659d150b73SToby Isaac   for (j = 0; j < numPoints; j++) {
28669d150b73SToby Isaac     const PetscReal *guess  = &refCoords[dimR * j];
28679d150b73SToby Isaac     PetscReal       *mapped = &realCoords[dimC * j];
28689d150b73SToby Isaac 
28699d150b73SToby Isaac     for (k = 0; k < numV; k++) {
28709d150b73SToby Isaac       PetscReal extCoord = 1.;
28719d150b73SToby Isaac       for (l = 0; l < dimR; l++) {
28729d150b73SToby Isaac         PetscReal coord = guess[l];
28739d150b73SToby Isaac         PetscInt  dep   = (k & (1 << l)) >> l;
28749d150b73SToby Isaac 
28759d150b73SToby Isaac         extCoord *= dep * coord + !dep;
28769d150b73SToby Isaac       }
28779d150b73SToby Isaac       for (l = 0; l < dimC; l++) {
28789d150b73SToby Isaac         PetscReal coeff = cellCoeffs[dimC * k + l];
28799d150b73SToby Isaac 
28809d150b73SToby Isaac         mapped[l] += coeff * extCoord;
28819d150b73SToby Isaac       }
28829d150b73SToby Isaac     }
28839d150b73SToby Isaac   }
288469291d52SBarry Smith   ierr = DMRestoreWorkArray(dm, 2 * coordSize, MPIU_REAL, &cellData);CHKERRQ(ierr);
28859d150b73SToby Isaac   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
28869d150b73SToby Isaac   PetscFunctionReturn(0);
28879d150b73SToby Isaac }
28889d150b73SToby Isaac 
28899c3cf19fSMatthew G. Knepley /* TODO: TOBY please fix this for Nc > 1 */
28909c3cf19fSMatthew G. Knepley static PetscErrorCode DMPlexCoordinatesToReference_FE(DM dm, PetscFE fe, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[], Vec coords, PetscInt Nc, PetscInt dimR)
28919d150b73SToby Isaac {
28929c3cf19fSMatthew G. Knepley   PetscInt       numComp, pdim, i, j, k, l, m, maxIter = 7, coordSize;
2893c6e120d1SToby Isaac   PetscScalar    *nodes = NULL;
2894c6e120d1SToby Isaac   PetscReal      *invV, *modes;
2895c6e120d1SToby Isaac   PetscReal      *B, *D, *resNeg;
2896c6e120d1SToby Isaac   PetscScalar    *J, *invJ, *work;
28979d150b73SToby Isaac   PetscErrorCode ierr;
28989d150b73SToby Isaac 
28999d150b73SToby Isaac   PetscFunctionBegin;
29009c3cf19fSMatthew G. Knepley   ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr);
29019d150b73SToby Isaac   ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
29029c3cf19fSMatthew G. Knepley   if (numComp != Nc) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"coordinate discretization must have as many components (%D) as embedding dimension (!= %D)",numComp,Nc);
29039d150b73SToby Isaac   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
29049d150b73SToby Isaac   /* convert nodes to values in the stable evaluation basis */
290569291d52SBarry Smith   ierr = DMGetWorkArray(dm,pdim,MPIU_REAL,&modes);CHKERRQ(ierr);
29069d150b73SToby Isaac   invV = fe->invV;
2907012b7cc6SMatthew G. Knepley   for (i = 0; i < pdim; ++i) {
2908012b7cc6SMatthew G. Knepley     modes[i] = 0.;
2909012b7cc6SMatthew G. Knepley     for (j = 0; j < pdim; ++j) {
2910012b7cc6SMatthew G. Knepley       modes[i] += invV[i * pdim + j] * PetscRealPart(nodes[j]);
29119d150b73SToby Isaac     }
29129d150b73SToby Isaac   }
291369291d52SBarry Smith   ierr   = DMGetWorkArray(dm,pdim * Nc + pdim * Nc * dimR + Nc,MPIU_REAL,&B);CHKERRQ(ierr);
29149c3cf19fSMatthew G. Knepley   D      = &B[pdim*Nc];
29159c3cf19fSMatthew G. Knepley   resNeg = &D[pdim*Nc * dimR];
291669291d52SBarry Smith   ierr = DMGetWorkArray(dm,3 * Nc * dimR,MPIU_SCALAR,&J);CHKERRQ(ierr);
29179c3cf19fSMatthew G. Knepley   invJ = &J[Nc * dimR];
29189c3cf19fSMatthew G. Knepley   work = &invJ[Nc * dimR];
29199d150b73SToby Isaac   for (i = 0; i < numPoints * dimR; i++) {refCoords[i] = 0.;}
29209d150b73SToby Isaac   for (j = 0; j < numPoints; j++) {
29219b1f03cbSToby 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 */
29229d150b73SToby Isaac       PetscReal *guess = &refCoords[j * dimR];
29239d150b73SToby Isaac       ierr = PetscSpaceEvaluate(fe->basisSpace, 1, guess, B, D, NULL);CHKERRQ(ierr);
29249c3cf19fSMatthew G. Knepley       for (k = 0; k < Nc; k++) {resNeg[k] = realCoords[j * Nc + k];}
29259c3cf19fSMatthew G. Knepley       for (k = 0; k < Nc * dimR; k++) {J[k] = 0.;}
29269c3cf19fSMatthew G. Knepley       for (k = 0; k < pdim; k++) {
29279c3cf19fSMatthew G. Knepley         for (l = 0; l < Nc; l++) {
2928012b7cc6SMatthew G. Knepley           resNeg[l] -= modes[k] * B[k * Nc + l];
29299d150b73SToby Isaac           for (m = 0; m < dimR; m++) {
2930012b7cc6SMatthew G. Knepley             J[l * dimR + m] += modes[k] * D[(k * Nc + l) * dimR + m];
29319d150b73SToby Isaac           }
29329d150b73SToby Isaac         }
29339d150b73SToby Isaac       }
29340611203eSToby Isaac #if 0 && defined(PETSC_USE_DEBUG)
29350611203eSToby Isaac       {
29360611203eSToby Isaac         PetscReal maxAbs = 0.;
29370611203eSToby Isaac 
29389c3cf19fSMatthew G. Knepley         for (l = 0; l < Nc; l++) {
29390611203eSToby Isaac           maxAbs = PetscMax(maxAbs,PetscAbsReal(resNeg[l]));
29400611203eSToby Isaac         }
2941087ef6b2SMatthew G. Knepley         ierr = PetscInfo4(dm,"cell %D, point %D, iter %D: res %g\n",cell,j,i,(double) maxAbs);CHKERRQ(ierr);
29420611203eSToby Isaac       }
29430611203eSToby Isaac #endif
29449c3cf19fSMatthew G. Knepley       ierr = DMPlexCoordinatesToReference_NewtonUpdate(Nc,dimR,J,invJ,work,resNeg,guess);CHKERRQ(ierr);
29459d150b73SToby Isaac     }
29469d150b73SToby Isaac   }
294769291d52SBarry Smith   ierr = DMRestoreWorkArray(dm,3 * Nc * dimR,MPIU_SCALAR,&J);CHKERRQ(ierr);
294869291d52SBarry Smith   ierr = DMRestoreWorkArray(dm,pdim * Nc + pdim * Nc * dimR + Nc,MPIU_REAL,&B);CHKERRQ(ierr);
294969291d52SBarry Smith   ierr = DMRestoreWorkArray(dm,pdim,MPIU_REAL,&modes);CHKERRQ(ierr);
29509d150b73SToby Isaac   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
29519d150b73SToby Isaac   PetscFunctionReturn(0);
29529d150b73SToby Isaac }
29539d150b73SToby Isaac 
29549c3cf19fSMatthew G. Knepley /* TODO: TOBY please fix this for Nc > 1 */
29559c3cf19fSMatthew G. Knepley static PetscErrorCode DMPlexReferenceToCoordinates_FE(DM dm, PetscFE fe, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[], Vec coords, PetscInt Nc, PetscInt dimR)
29569d150b73SToby Isaac {
29579c3cf19fSMatthew G. Knepley   PetscInt       numComp, pdim, i, j, k, l, coordSize;
2958c6e120d1SToby Isaac   PetscScalar    *nodes = NULL;
2959c6e120d1SToby Isaac   PetscReal      *invV, *modes;
29609d150b73SToby Isaac   PetscReal      *B;
29619d150b73SToby Isaac   PetscErrorCode ierr;
29629d150b73SToby Isaac 
29639d150b73SToby Isaac   PetscFunctionBegin;
29649c3cf19fSMatthew G. Knepley   ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr);
29659d150b73SToby Isaac   ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
29669c3cf19fSMatthew G. Knepley   if (numComp != Nc) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"coordinate discretization must have as many components (%D) as embedding dimension (!= %D)",numComp,Nc);
29679d150b73SToby Isaac   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
29689d150b73SToby Isaac   /* convert nodes to values in the stable evaluation basis */
296969291d52SBarry Smith   ierr = DMGetWorkArray(dm,pdim,MPIU_REAL,&modes);CHKERRQ(ierr);
29709d150b73SToby Isaac   invV = fe->invV;
2971012b7cc6SMatthew G. Knepley   for (i = 0; i < pdim; ++i) {
2972012b7cc6SMatthew G. Knepley     modes[i] = 0.;
2973012b7cc6SMatthew G. Knepley     for (j = 0; j < pdim; ++j) {
2974012b7cc6SMatthew G. Knepley       modes[i] += invV[i * pdim + j] * PetscRealPart(nodes[j]);
29759d150b73SToby Isaac     }
29769d150b73SToby Isaac   }
297769291d52SBarry Smith   ierr = DMGetWorkArray(dm,numPoints * pdim * Nc,MPIU_REAL,&B);CHKERRQ(ierr);
2978012b7cc6SMatthew G. Knepley   ierr = PetscSpaceEvaluate(fe->basisSpace, numPoints, refCoords, B, NULL, NULL);CHKERRQ(ierr);
29799c3cf19fSMatthew G. Knepley   for (i = 0; i < numPoints * Nc; i++) {realCoords[i] = 0.;}
29809d150b73SToby Isaac   for (j = 0; j < numPoints; j++) {
29819c3cf19fSMatthew G. Knepley     PetscReal *mapped = &realCoords[j * Nc];
29829d150b73SToby Isaac 
29839c3cf19fSMatthew G. Knepley     for (k = 0; k < pdim; k++) {
29849c3cf19fSMatthew G. Knepley       for (l = 0; l < Nc; l++) {
298540cf36b3SToby Isaac         mapped[l] += modes[k] * B[(j * pdim + k) * Nc + l];
29869d150b73SToby Isaac       }
29879d150b73SToby Isaac     }
29889d150b73SToby Isaac   }
298969291d52SBarry Smith   ierr = DMRestoreWorkArray(dm,numPoints * pdim * Nc,MPIU_REAL,&B);CHKERRQ(ierr);
299069291d52SBarry Smith   ierr = DMRestoreWorkArray(dm,pdim,MPIU_REAL,&modes);CHKERRQ(ierr);
29919d150b73SToby Isaac   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
29929d150b73SToby Isaac   PetscFunctionReturn(0);
29939d150b73SToby Isaac }
29949d150b73SToby Isaac 
2995d6143a4eSToby Isaac /*@
2996d6143a4eSToby Isaac   DMPlexCoordinatesToReference - Pull coordinates back from the mesh to the reference element using a single element
2997d6143a4eSToby Isaac   map.  This inversion will be accurate inside the reference element, but may be inaccurate for mappings that do not
2998d6143a4eSToby Isaac   extend uniquely outside the reference cell (e.g, most non-affine maps)
2999d6143a4eSToby Isaac 
3000d6143a4eSToby Isaac   Not collective
3001d6143a4eSToby Isaac 
3002d6143a4eSToby Isaac   Input Parameters:
3003d6143a4eSToby Isaac + dm         - The mesh, with coordinate maps defined either by a PetscDS for the coordinate DM (see DMGetCoordinateDM()) or
3004d6143a4eSToby Isaac                implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or
3005d6143a4eSToby Isaac                as a multilinear map for tensor-product elements
3006d6143a4eSToby Isaac . cell       - the cell whose map is used.
3007d6143a4eSToby Isaac . numPoints  - the number of points to locate
30081b266c99SBarry Smith - realCoords - (numPoints x coordinate dimension) array of coordinates (see DMGetCoordinateDim())
3009d6143a4eSToby Isaac 
3010d6143a4eSToby Isaac   Output Parameters:
3011d6143a4eSToby Isaac . refCoords  - (numPoints x dimension) array of reference coordinates (see DMGetDimension())
30121b266c99SBarry Smith 
30131b266c99SBarry Smith   Level: intermediate
301473c9229bSMatthew Knepley 
301573c9229bSMatthew Knepley .seealso: DMPlexReferenceToCoordinates()
3016d6143a4eSToby Isaac @*/
3017d6143a4eSToby Isaac PetscErrorCode DMPlexCoordinatesToReference(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[])
3018d6143a4eSToby Isaac {
3019485ad865SMatthew G. Knepley   PetscInt       dimC, dimR, depth, cStart, cEnd, i;
30209d150b73SToby Isaac   DM             coordDM = NULL;
30219d150b73SToby Isaac   Vec            coords;
30229d150b73SToby Isaac   PetscFE        fe = NULL;
30239d150b73SToby Isaac   PetscErrorCode ierr;
30249d150b73SToby Isaac 
3025d6143a4eSToby Isaac   PetscFunctionBegin;
30269d150b73SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
30279d150b73SToby Isaac   ierr = DMGetDimension(dm,&dimR);CHKERRQ(ierr);
30289d150b73SToby Isaac   ierr = DMGetCoordinateDim(dm,&dimC);CHKERRQ(ierr);
30299d150b73SToby Isaac   if (dimR <= 0 || dimC <= 0 || numPoints <= 0) PetscFunctionReturn(0);
30309d150b73SToby Isaac   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
30319d150b73SToby Isaac   ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
30329d150b73SToby Isaac   ierr = DMGetCoordinateDM(dm,&coordDM);CHKERRQ(ierr);
30339d150b73SToby Isaac   if (coordDM) {
30349d150b73SToby Isaac     PetscInt coordFields;
30359d150b73SToby Isaac 
30369d150b73SToby Isaac     ierr = DMGetNumFields(coordDM,&coordFields);CHKERRQ(ierr);
30379d150b73SToby Isaac     if (coordFields) {
30389d150b73SToby Isaac       PetscClassId id;
30399d150b73SToby Isaac       PetscObject  disc;
30409d150b73SToby Isaac 
304144a7f3ddSMatthew G. Knepley       ierr = DMGetField(coordDM,0,NULL,&disc);CHKERRQ(ierr);
30429d150b73SToby Isaac       ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
30439d150b73SToby Isaac       if (id == PETSCFE_CLASSID) {
30449d150b73SToby Isaac         fe = (PetscFE) disc;
30459d150b73SToby Isaac       }
30469d150b73SToby Isaac     }
30479d150b73SToby Isaac   }
3048485ad865SMatthew G. Knepley   ierr = DMPlexGetInteriorCellStratum(dm,&cStart,&cEnd);CHKERRQ(ierr);
304913903a91SSatish Balay   if (cell < cStart || cell >= cEnd) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D not in cell range [%D,%D)",cell,cStart,cEnd);
30509d150b73SToby Isaac   if (!fe) { /* implicit discretization: affine or multilinear */
30519d150b73SToby Isaac     PetscInt  coneSize;
30529d150b73SToby Isaac     PetscBool isSimplex, isTensor;
30539d150b73SToby Isaac 
30549d150b73SToby Isaac     ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr);
30559d150b73SToby Isaac     isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE;
30569d150b73SToby Isaac     isTensor  = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE;
30579d150b73SToby Isaac     if (isSimplex) {
30589d150b73SToby Isaac       PetscReal detJ, *v0, *J, *invJ;
30599d150b73SToby Isaac 
306069291d52SBarry Smith       ierr = DMGetWorkArray(dm,dimC + 2 * dimC * dimC, MPIU_REAL, &v0);CHKERRQ(ierr);
30619d150b73SToby Isaac       J    = &v0[dimC];
30629d150b73SToby Isaac       invJ = &J[dimC * dimC];
30639d150b73SToby Isaac       ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, &detJ);CHKERRQ(ierr);
30649d150b73SToby Isaac       for (i = 0; i < numPoints; i++) { /* Apply the inverse affine transformation for each point */
3065c330f8ffSToby Isaac         const PetscReal x0[3] = {-1.,-1.,-1.};
3066c330f8ffSToby Isaac 
3067c330f8ffSToby Isaac         CoordinatesRealToRef(dimC, dimR, x0, v0, invJ, &realCoords[dimC * i], &refCoords[dimR * i]);
30689d150b73SToby Isaac       }
306969291d52SBarry Smith       ierr = DMRestoreWorkArray(dm,dimC + 2 * dimC * dimC, MPIU_REAL, &v0);CHKERRQ(ierr);
30709d150b73SToby Isaac     } else if (isTensor) {
30719d150b73SToby Isaac       ierr = DMPlexCoordinatesToReference_Tensor(coordDM, cell, numPoints, realCoords, refCoords, coords, dimC, dimR);CHKERRQ(ierr);
30729d150b73SToby Isaac     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unrecognized cone size %D",coneSize);
30739d150b73SToby Isaac   } else {
30749d150b73SToby Isaac     ierr = DMPlexCoordinatesToReference_FE(coordDM, fe, cell, numPoints, realCoords, refCoords, coords, dimC, dimR);CHKERRQ(ierr);
30759d150b73SToby Isaac   }
30769d150b73SToby Isaac   PetscFunctionReturn(0);
30779d150b73SToby Isaac }
30789d150b73SToby Isaac 
30799d150b73SToby Isaac /*@
30809d150b73SToby Isaac   DMPlexReferenceToCoordinates - Map references coordinates to coordinates in the the mesh for a single element map.
30819d150b73SToby Isaac 
30829d150b73SToby Isaac   Not collective
30839d150b73SToby Isaac 
30849d150b73SToby Isaac   Input Parameters:
30859d150b73SToby Isaac + dm         - The mesh, with coordinate maps defined either by a PetscDS for the coordinate DM (see DMGetCoordinateDM()) or
30869d150b73SToby Isaac                implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or
30879d150b73SToby Isaac                as a multilinear map for tensor-product elements
30889d150b73SToby Isaac . cell       - the cell whose map is used.
30899d150b73SToby Isaac . numPoints  - the number of points to locate
3090a2b725a8SWilliam Gropp - refCoords  - (numPoints x dimension) array of reference coordinates (see DMGetDimension())
30919d150b73SToby Isaac 
30929d150b73SToby Isaac   Output Parameters:
30939d150b73SToby Isaac . realCoords - (numPoints x coordinate dimension) array of coordinates (see DMGetCoordinateDim())
30941b266c99SBarry Smith 
30951b266c99SBarry Smith    Level: intermediate
309673c9229bSMatthew Knepley 
309773c9229bSMatthew Knepley .seealso: DMPlexCoordinatesToReference()
30989d150b73SToby Isaac @*/
30999d150b73SToby Isaac PetscErrorCode DMPlexReferenceToCoordinates(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[])
31009d150b73SToby Isaac {
3101485ad865SMatthew G. Knepley   PetscInt       dimC, dimR, depth, cStart, cEnd, i;
31029d150b73SToby Isaac   DM             coordDM = NULL;
31039d150b73SToby Isaac   Vec            coords;
31049d150b73SToby Isaac   PetscFE        fe = NULL;
31059d150b73SToby Isaac   PetscErrorCode ierr;
31069d150b73SToby Isaac 
31079d150b73SToby Isaac   PetscFunctionBegin;
31089d150b73SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
31099d150b73SToby Isaac   ierr = DMGetDimension(dm,&dimR);CHKERRQ(ierr);
31109d150b73SToby Isaac   ierr = DMGetCoordinateDim(dm,&dimC);CHKERRQ(ierr);
31119d150b73SToby Isaac   if (dimR <= 0 || dimC <= 0 || numPoints <= 0) PetscFunctionReturn(0);
31129d150b73SToby Isaac   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
31139d150b73SToby Isaac   ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
31149d150b73SToby Isaac   ierr = DMGetCoordinateDM(dm,&coordDM);CHKERRQ(ierr);
31159d150b73SToby Isaac   if (coordDM) {
31169d150b73SToby Isaac     PetscInt coordFields;
31179d150b73SToby Isaac 
31189d150b73SToby Isaac     ierr = DMGetNumFields(coordDM,&coordFields);CHKERRQ(ierr);
31199d150b73SToby Isaac     if (coordFields) {
31209d150b73SToby Isaac       PetscClassId id;
31219d150b73SToby Isaac       PetscObject  disc;
31229d150b73SToby Isaac 
312344a7f3ddSMatthew G. Knepley       ierr = DMGetField(coordDM,0,NULL,&disc);CHKERRQ(ierr);
31249d150b73SToby Isaac       ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
31259d150b73SToby Isaac       if (id == PETSCFE_CLASSID) {
31269d150b73SToby Isaac         fe = (PetscFE) disc;
31279d150b73SToby Isaac       }
31289d150b73SToby Isaac     }
31299d150b73SToby Isaac   }
3130485ad865SMatthew G. Knepley   ierr = DMPlexGetInteriorCellStratum(dm,&cStart,&cEnd);CHKERRQ(ierr);
313113903a91SSatish Balay   if (cell < cStart || cell >= cEnd) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D not in cell range [%D,%D)",cell,cStart,cEnd);
31329d150b73SToby Isaac   if (!fe) { /* implicit discretization: affine or multilinear */
31339d150b73SToby Isaac     PetscInt  coneSize;
31349d150b73SToby Isaac     PetscBool isSimplex, isTensor;
31359d150b73SToby Isaac 
31369d150b73SToby Isaac     ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr);
31379d150b73SToby Isaac     isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE;
31389d150b73SToby Isaac     isTensor  = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE;
31399d150b73SToby Isaac     if (isSimplex) {
31409d150b73SToby Isaac       PetscReal detJ, *v0, *J;
31419d150b73SToby Isaac 
314269291d52SBarry Smith       ierr = DMGetWorkArray(dm,dimC + 2 * dimC * dimC, MPIU_REAL, &v0);CHKERRQ(ierr);
31439d150b73SToby Isaac       J    = &v0[dimC];
31449d150b73SToby Isaac       ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, NULL, &detJ);CHKERRQ(ierr);
3145c330f8ffSToby Isaac       for (i = 0; i < numPoints; i++) { /* Apply the affine transformation for each point */
3146c330f8ffSToby Isaac         const PetscReal xi0[3] = {-1.,-1.,-1.};
3147c330f8ffSToby Isaac 
3148c330f8ffSToby Isaac         CoordinatesRefToReal(dimC, dimR, xi0, v0, J, &refCoords[dimR * i], &realCoords[dimC * i]);
31499d150b73SToby Isaac       }
315069291d52SBarry Smith       ierr = DMRestoreWorkArray(dm,dimC + 2 * dimC * dimC, MPIU_REAL, &v0);CHKERRQ(ierr);
31519d150b73SToby Isaac     } else if (isTensor) {
31529d150b73SToby Isaac       ierr = DMPlexReferenceToCoordinates_Tensor(coordDM, cell, numPoints, refCoords, realCoords, coords, dimC, dimR);CHKERRQ(ierr);
31539d150b73SToby Isaac     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unrecognized cone size %D",coneSize);
31549d150b73SToby Isaac   } else {
31559d150b73SToby Isaac     ierr = DMPlexReferenceToCoordinates_FE(coordDM, fe, cell, numPoints, refCoords, realCoords, coords, dimC, dimR);CHKERRQ(ierr);
31569d150b73SToby Isaac   }
3157d6143a4eSToby Isaac   PetscFunctionReturn(0);
3158d6143a4eSToby Isaac }
31590139fca9SMatthew G. Knepley 
31600139fca9SMatthew G. Knepley /*@C
31610139fca9SMatthew G. Knepley   DMPlexRemapGeometry - This function maps the original DM coordinates to new coordinates.
31620139fca9SMatthew G. Knepley 
31630139fca9SMatthew G. Knepley   Not collective
31640139fca9SMatthew G. Knepley 
31650139fca9SMatthew G. Knepley   Input Parameters:
31660139fca9SMatthew G. Knepley + dm      - The DM
31670139fca9SMatthew G. Knepley . time    - The time
31680139fca9SMatthew G. Knepley - func    - The function transforming current coordinates to new coordaintes
31690139fca9SMatthew G. Knepley 
31700139fca9SMatthew G. Knepley    Calling sequence of func:
31710139fca9SMatthew G. Knepley $    func(PetscInt dim, PetscInt Nf, PetscInt NfAux,
31720139fca9SMatthew G. Knepley $         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
31730139fca9SMatthew G. Knepley $         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
31740139fca9SMatthew G. Knepley $         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]);
31750139fca9SMatthew G. Knepley 
31760139fca9SMatthew G. Knepley +  dim          - The spatial dimension
31770139fca9SMatthew G. Knepley .  Nf           - The number of input fields (here 1)
31780139fca9SMatthew G. Knepley .  NfAux        - The number of input auxiliary fields
31790139fca9SMatthew G. Knepley .  uOff         - The offset of the coordinates in u[] (here 0)
31800139fca9SMatthew G. Knepley .  uOff_x       - The offset of the coordinates in u_x[] (here 0)
31810139fca9SMatthew G. Knepley .  u            - The coordinate values at this point in space
31820139fca9SMatthew G. Knepley .  u_t          - The coordinate time derivative at this point in space (here NULL)
31830139fca9SMatthew G. Knepley .  u_x          - The coordinate derivatives at this point in space
31840139fca9SMatthew G. Knepley .  aOff         - The offset of each auxiliary field in u[]
31850139fca9SMatthew G. Knepley .  aOff_x       - The offset of each auxiliary field in u_x[]
31860139fca9SMatthew G. Knepley .  a            - The auxiliary field values at this point in space
31870139fca9SMatthew G. Knepley .  a_t          - The auxiliary field time derivative at this point in space (or NULL)
31880139fca9SMatthew G. Knepley .  a_x          - The auxiliary field derivatives at this point in space
31890139fca9SMatthew G. Knepley .  t            - The current time
31900139fca9SMatthew G. Knepley .  x            - The coordinates of this point (here not used)
31910139fca9SMatthew G. Knepley .  numConstants - The number of constants
31920139fca9SMatthew G. Knepley .  constants    - The value of each constant
31930139fca9SMatthew G. Knepley -  f            - The new coordinates at this point in space
31940139fca9SMatthew G. Knepley 
31950139fca9SMatthew G. Knepley   Level: intermediate
31960139fca9SMatthew G. Knepley 
31970139fca9SMatthew G. Knepley .seealso: DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMProjectFieldLocal(), DMProjectFieldLabelLocal()
31980139fca9SMatthew G. Knepley @*/
31990139fca9SMatthew G. Knepley PetscErrorCode DMPlexRemapGeometry(DM dm, PetscReal time,
32000139fca9SMatthew G. Knepley                                    void (*func)(PetscInt, PetscInt, PetscInt,
32010139fca9SMatthew G. Knepley                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
32020139fca9SMatthew G. Knepley                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
32030139fca9SMatthew G. Knepley                                                 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
32040139fca9SMatthew G. Knepley {
32050139fca9SMatthew G. Knepley   DM             cdm;
32060139fca9SMatthew G. Knepley   Vec            lCoords, tmpCoords;
32070139fca9SMatthew G. Knepley   PetscErrorCode ierr;
32080139fca9SMatthew G. Knepley 
32090139fca9SMatthew G. Knepley   PetscFunctionBegin;
32100139fca9SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
32110139fca9SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &lCoords);CHKERRQ(ierr);
32120139fca9SMatthew G. Knepley   ierr = DMGetLocalVector(cdm, &tmpCoords);CHKERRQ(ierr);
32130139fca9SMatthew G. Knepley   ierr = VecCopy(lCoords, tmpCoords);CHKERRQ(ierr);
32140139fca9SMatthew G. Knepley   ierr = DMProjectFieldLocal(cdm, time, tmpCoords, &func, INSERT_VALUES, lCoords);CHKERRQ(ierr);
32150139fca9SMatthew G. Knepley   ierr = DMRestoreLocalVector(cdm, &tmpCoords);CHKERRQ(ierr);
32160139fca9SMatthew G. Knepley   PetscFunctionReturn(0);
32170139fca9SMatthew G. Knepley }
32180139fca9SMatthew G. Knepley 
32190139fca9SMatthew G. Knepley /* Shear applies the transformation, assuming we fix z,
32200139fca9SMatthew G. Knepley   / 1  0  m_0 \
32210139fca9SMatthew G. Knepley   | 0  1  m_1 |
32220139fca9SMatthew G. Knepley   \ 0  0   1  /
32230139fca9SMatthew G. Knepley */
32240139fca9SMatthew G. Knepley static void f0_shear(PetscInt dim, PetscInt Nf, PetscInt NfAux,
32250139fca9SMatthew G. Knepley                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
32260139fca9SMatthew G. Knepley                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
32270139fca9SMatthew G. Knepley                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar coords[])
32280139fca9SMatthew G. Knepley {
32290139fca9SMatthew G. Knepley   const PetscInt Nc = uOff[1]-uOff[0];
3230c1f1bd54SMatthew G. Knepley   const PetscInt ax = (PetscInt) PetscRealPart(constants[0]);
32310139fca9SMatthew G. Knepley   PetscInt       c;
32320139fca9SMatthew G. Knepley 
32330139fca9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
32340139fca9SMatthew G. Knepley     coords[c] = u[c] + constants[c+1]*u[ax];
32350139fca9SMatthew G. Knepley   }
32360139fca9SMatthew G. Knepley }
32370139fca9SMatthew G. Knepley 
32380139fca9SMatthew G. Knepley /*@C
32390139fca9SMatthew G. Knepley   DMPlexShearGeometry - This shears the domain, meaning adds a multiple of the shear coordinate to all other coordinates.
32400139fca9SMatthew G. Knepley 
32410139fca9SMatthew G. Knepley   Not collective
32420139fca9SMatthew G. Knepley 
32430139fca9SMatthew G. Knepley   Input Parameters:
32440139fca9SMatthew G. Knepley + dm          - The DM
32453ee9839eSMatthew G. Knepley . direction   - The shear coordinate direction, e.g. 0 is the x-axis
32460139fca9SMatthew G. Knepley - multipliers - The multiplier m for each direction which is not the shear direction
32470139fca9SMatthew G. Knepley 
32480139fca9SMatthew G. Knepley   Level: intermediate
32490139fca9SMatthew G. Knepley 
32500139fca9SMatthew G. Knepley .seealso: DMPlexRemapGeometry()
32510139fca9SMatthew G. Knepley @*/
32523ee9839eSMatthew G. Knepley PetscErrorCode DMPlexShearGeometry(DM dm, DMDirection direction, PetscReal multipliers[])
32530139fca9SMatthew G. Knepley {
32540139fca9SMatthew G. Knepley   DM             cdm;
32550139fca9SMatthew G. Knepley   PetscDS        cds;
32560139fca9SMatthew G. Knepley   PetscObject    obj;
32570139fca9SMatthew G. Knepley   PetscClassId   id;
32580139fca9SMatthew G. Knepley   PetscScalar   *moduli;
32593ee9839eSMatthew G. Knepley   const PetscInt dir = (PetscInt) direction;
32600139fca9SMatthew G. Knepley   PetscInt       dE, d, e;
32610139fca9SMatthew G. Knepley   PetscErrorCode ierr;
32620139fca9SMatthew G. Knepley 
32630139fca9SMatthew G. Knepley   PetscFunctionBegin;
32640139fca9SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
32650139fca9SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dE);CHKERRQ(ierr);
32660139fca9SMatthew G. Knepley   ierr = PetscMalloc1(dE+1, &moduli);CHKERRQ(ierr);
32670139fca9SMatthew G. Knepley   moduli[0] = dir;
32680139fca9SMatthew G. Knepley   for (d = 0, e = 0; d < dE; ++d) moduli[d] = d == dir ? 0.0 : (multipliers ? multipliers[e++] : 1.0);
32690139fca9SMatthew G. Knepley   ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr);
32700139fca9SMatthew G. Knepley   ierr = PetscDSGetDiscretization(cds, 0, &obj);CHKERRQ(ierr);
32710139fca9SMatthew G. Knepley   ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
32720139fca9SMatthew G. Knepley   if (id != PETSCFE_CLASSID) {
32730139fca9SMatthew G. Knepley     Vec           lCoords;
32740139fca9SMatthew G. Knepley     PetscSection  cSection;
32750139fca9SMatthew G. Knepley     PetscScalar  *coords;
32760139fca9SMatthew G. Knepley     PetscInt      vStart, vEnd, v;
32770139fca9SMatthew G. Knepley 
32780139fca9SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
32790139fca9SMatthew G. Knepley     ierr = DMGetCoordinateSection(dm, &cSection);CHKERRQ(ierr);
32800139fca9SMatthew G. Knepley     ierr = DMGetCoordinatesLocal(dm, &lCoords);CHKERRQ(ierr);
32810139fca9SMatthew G. Knepley     ierr = VecGetArray(lCoords, &coords);CHKERRQ(ierr);
32820139fca9SMatthew G. Knepley     for (v = vStart; v < vEnd; ++v) {
32830139fca9SMatthew G. Knepley       PetscReal ds;
32840139fca9SMatthew G. Knepley       PetscInt  off, c;
32850139fca9SMatthew G. Knepley 
32860139fca9SMatthew G. Knepley       ierr = PetscSectionGetOffset(cSection, v, &off);CHKERRQ(ierr);
32870139fca9SMatthew G. Knepley       ds   = PetscRealPart(coords[off+dir]);
32880139fca9SMatthew G. Knepley       for (c = 0; c < dE; ++c) coords[off+c] += moduli[c]*ds;
32890139fca9SMatthew G. Knepley     }
32900139fca9SMatthew G. Knepley     ierr = VecRestoreArray(lCoords, &coords);CHKERRQ(ierr);
32910139fca9SMatthew G. Knepley   } else {
32920139fca9SMatthew G. Knepley     ierr = PetscDSSetConstants(cds, dE+1, moduli);CHKERRQ(ierr);
32930139fca9SMatthew G. Knepley     ierr = DMPlexRemapGeometry(dm, 0.0, f0_shear);CHKERRQ(ierr);
32940139fca9SMatthew G. Knepley   }
32950139fca9SMatthew G. Knepley   ierr = PetscFree(moduli);CHKERRQ(ierr);
32960139fca9SMatthew G. Knepley   PetscFunctionReturn(0);
32970139fca9SMatthew G. Knepley }
3298