xref: /petsc/src/dm/impls/plex/plexgeometry.c (revision 4f99dae567c24cd8db2a416c01ed058f644d10b4)
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 
3337900f7dSMatthew G. Knepley .seealso: DMPlexCreate(), DMGetCoordinatesLocal()
343985bb02SVaclav Hapla @*/
353985bb02SVaclav Hapla PetscErrorCode DMPlexFindVertices(DM dm, PetscInt npoints, const PetscReal coord[], PetscReal eps, PetscInt dagPoints[])
363985bb02SVaclav Hapla {
3737900f7dSMatthew G. Knepley   PetscInt          c, cdim, 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;
4537900f7dSMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &cdim);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);
4976bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
50335ef845SVaclav Hapla     /* check coordinate section is consistent with DM dimension */
51335ef845SVaclav Hapla     PetscSection      cs;
52335ef845SVaclav Hapla     PetscInt          ndof;
53335ef845SVaclav Hapla 
54335ef845SVaclav Hapla     ierr = DMGetCoordinateSection(dm, &cs);CHKERRQ(ierr);
553985bb02SVaclav Hapla     for (p = vStart; p < vEnd; p++) {
563985bb02SVaclav Hapla       ierr = PetscSectionGetDof(cs, p, &ndof);CHKERRQ(ierr);
5737900f7dSMatthew G. Knepley       if (PetscUnlikely(ndof != cdim)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "point %D: ndof = %D != %D = cdim", p, ndof, cdim);
58335ef845SVaclav Hapla     }
59335ef845SVaclav Hapla   }
60eca9f518SVaclav Hapla   if (eps == 0.0) {
6137900f7dSMatthew G. Knepley     for (i=0,j=0; i < npoints; i++,j+=cdim) {
62eca9f518SVaclav Hapla       dagPoints[i] = -1;
6337900f7dSMatthew G. Knepley       for (p = vStart,o=0; p < vEnd; p++,o+=cdim) {
6437900f7dSMatthew G. Knepley         for (c = 0; c < cdim; c++) {
65eca9f518SVaclav Hapla           if (coord[j+c] != PetscRealPart(allCoords[o+c])) break;
66eca9f518SVaclav Hapla         }
6737900f7dSMatthew G. Knepley         if (c == cdim) {
68eca9f518SVaclav Hapla           dagPoints[i] = p;
69eca9f518SVaclav Hapla           break;
70eca9f518SVaclav Hapla         }
71eca9f518SVaclav Hapla       }
72eca9f518SVaclav Hapla     }
73eca9f518SVaclav Hapla     ierr = VecRestoreArrayRead(allCoordsVec, &allCoords);CHKERRQ(ierr);
74eca9f518SVaclav Hapla     PetscFunctionReturn(0);
75eca9f518SVaclav Hapla   }
7637900f7dSMatthew G. Knepley   for (i=0,j=0; i < npoints; i++,j+=cdim) {
77335ef845SVaclav Hapla     dagPoints[i] = -1;
7837900f7dSMatthew G. Knepley     for (p = vStart,o=0; p < vEnd; p++,o+=cdim) {
793985bb02SVaclav Hapla       norm = 0.0;
8037900f7dSMatthew G. Knepley       for (c = 0; c < cdim; c++) {
813985bb02SVaclav Hapla         norm += PetscSqr(coord[j+c] - PetscRealPart(allCoords[o+c]));
823985bb02SVaclav Hapla       }
833985bb02SVaclav Hapla       norm = PetscSqrtReal(norm);
843985bb02SVaclav Hapla       if (norm <= eps) {
853985bb02SVaclav Hapla         dagPoints[i] = p;
863985bb02SVaclav Hapla         break;
873985bb02SVaclav Hapla       }
883985bb02SVaclav Hapla     }
893985bb02SVaclav Hapla   }
903985bb02SVaclav Hapla   ierr = VecRestoreArrayRead(allCoordsVec, &allCoords);CHKERRQ(ierr);
913985bb02SVaclav Hapla   PetscFunctionReturn(0);
923985bb02SVaclav Hapla }
933985bb02SVaclav Hapla 
94fea14342SMatthew G. Knepley static PetscErrorCode DMPlexGetLineIntersection_2D_Internal(const PetscReal segmentA[], const PetscReal segmentB[], PetscReal intersection[], PetscBool *hasIntersection)
95fea14342SMatthew G. Knepley {
96fea14342SMatthew G. Knepley   const PetscReal p0_x  = segmentA[0*2+0];
97fea14342SMatthew G. Knepley   const PetscReal p0_y  = segmentA[0*2+1];
98fea14342SMatthew G. Knepley   const PetscReal p1_x  = segmentA[1*2+0];
99fea14342SMatthew G. Knepley   const PetscReal p1_y  = segmentA[1*2+1];
100fea14342SMatthew G. Knepley   const PetscReal p2_x  = segmentB[0*2+0];
101fea14342SMatthew G. Knepley   const PetscReal p2_y  = segmentB[0*2+1];
102fea14342SMatthew G. Knepley   const PetscReal p3_x  = segmentB[1*2+0];
103fea14342SMatthew G. Knepley   const PetscReal p3_y  = segmentB[1*2+1];
104fea14342SMatthew G. Knepley   const PetscReal s1_x  = p1_x - p0_x;
105fea14342SMatthew G. Knepley   const PetscReal s1_y  = p1_y - p0_y;
106fea14342SMatthew G. Knepley   const PetscReal s2_x  = p3_x - p2_x;
107fea14342SMatthew G. Knepley   const PetscReal s2_y  = p3_y - p2_y;
108fea14342SMatthew G. Knepley   const PetscReal denom = (-s2_x * s1_y + s1_x * s2_y);
109fea14342SMatthew G. Knepley 
110fea14342SMatthew G. Knepley   PetscFunctionBegin;
111fea14342SMatthew G. Knepley   *hasIntersection = PETSC_FALSE;
112fea14342SMatthew G. Knepley   /* Non-parallel lines */
113fea14342SMatthew G. Knepley   if (denom != 0.0) {
114fea14342SMatthew G. Knepley     const PetscReal s = (-s1_y * (p0_x - p2_x) + s1_x * (p0_y - p2_y)) / denom;
115fea14342SMatthew G. Knepley     const PetscReal t = ( s2_x * (p0_y - p2_y) - s2_y * (p0_x - p2_x)) / denom;
116fea14342SMatthew G. Knepley 
117fea14342SMatthew G. Knepley     if (s >= 0 && s <= 1 && t >= 0 && t <= 1) {
118fea14342SMatthew G. Knepley       *hasIntersection = PETSC_TRUE;
119fea14342SMatthew G. Knepley       if (intersection) {
120fea14342SMatthew G. Knepley         intersection[0] = p0_x + (t * s1_x);
121fea14342SMatthew G. Knepley         intersection[1] = p0_y + (t * s1_y);
122fea14342SMatthew G. Knepley       }
123fea14342SMatthew G. Knepley     }
124fea14342SMatthew G. Knepley   }
125fea14342SMatthew G. Knepley   PetscFunctionReturn(0);
126fea14342SMatthew G. Knepley }
127fea14342SMatthew G. Knepley 
128ddce0771SMatthew G. Knepley /* The plane is segmentB x segmentC: https://en.wikipedia.org/wiki/Line%E2%80%93plane_intersection */
129ddce0771SMatthew G. Knepley static PetscErrorCode DMPlexGetLinePlaneIntersection_3D_Internal(const PetscReal segmentA[], const PetscReal segmentB[], const PetscReal segmentC[], PetscReal intersection[], PetscBool *hasIntersection)
130ddce0771SMatthew G. Knepley {
131ddce0771SMatthew G. Knepley   const PetscReal p0_x  = segmentA[0*3+0];
132ddce0771SMatthew G. Knepley   const PetscReal p0_y  = segmentA[0*3+1];
133ddce0771SMatthew G. Knepley   const PetscReal p0_z  = segmentA[0*3+2];
134ddce0771SMatthew G. Knepley   const PetscReal p1_x  = segmentA[1*3+0];
135ddce0771SMatthew G. Knepley   const PetscReal p1_y  = segmentA[1*3+1];
136ddce0771SMatthew G. Knepley   const PetscReal p1_z  = segmentA[1*3+2];
137ddce0771SMatthew G. Knepley   const PetscReal q0_x  = segmentB[0*3+0];
138ddce0771SMatthew G. Knepley   const PetscReal q0_y  = segmentB[0*3+1];
139ddce0771SMatthew G. Knepley   const PetscReal q0_z  = segmentB[0*3+2];
140ddce0771SMatthew G. Knepley   const PetscReal q1_x  = segmentB[1*3+0];
141ddce0771SMatthew G. Knepley   const PetscReal q1_y  = segmentB[1*3+1];
142ddce0771SMatthew G. Knepley   const PetscReal q1_z  = segmentB[1*3+2];
143ddce0771SMatthew G. Knepley   const PetscReal r0_x  = segmentC[0*3+0];
144ddce0771SMatthew G. Knepley   const PetscReal r0_y  = segmentC[0*3+1];
145ddce0771SMatthew G. Knepley   const PetscReal r0_z  = segmentC[0*3+2];
146ddce0771SMatthew G. Knepley   const PetscReal r1_x  = segmentC[1*3+0];
147ddce0771SMatthew G. Knepley   const PetscReal r1_y  = segmentC[1*3+1];
148ddce0771SMatthew G. Knepley   const PetscReal r1_z  = segmentC[1*3+2];
149ddce0771SMatthew G. Knepley   const PetscReal s0_x  = p1_x - p0_x;
150ddce0771SMatthew G. Knepley   const PetscReal s0_y  = p1_y - p0_y;
151ddce0771SMatthew G. Knepley   const PetscReal s0_z  = p1_z - p0_z;
152ddce0771SMatthew G. Knepley   const PetscReal s1_x  = q1_x - q0_x;
153ddce0771SMatthew G. Knepley   const PetscReal s1_y  = q1_y - q0_y;
154ddce0771SMatthew G. Knepley   const PetscReal s1_z  = q1_z - q0_z;
155ddce0771SMatthew G. Knepley   const PetscReal s2_x  = r1_x - r0_x;
156ddce0771SMatthew G. Knepley   const PetscReal s2_y  = r1_y - r0_y;
157ddce0771SMatthew G. Knepley   const PetscReal s2_z  = r1_z - r0_z;
158ddce0771SMatthew G. Knepley   const PetscReal s3_x  = s1_y*s2_z - s1_z*s2_y; /* s1 x s2 */
159ddce0771SMatthew G. Knepley   const PetscReal s3_y  = s1_z*s2_x - s1_x*s2_z;
160ddce0771SMatthew G. Knepley   const PetscReal s3_z  = s1_x*s2_y - s1_y*s2_x;
161ddce0771SMatthew G. Knepley   const PetscReal s4_x  = s0_y*s2_z - s0_z*s2_y; /* s0 x s2 */
162ddce0771SMatthew G. Knepley   const PetscReal s4_y  = s0_z*s2_x - s0_x*s2_z;
163ddce0771SMatthew G. Knepley   const PetscReal s4_z  = s0_x*s2_y - s0_y*s2_x;
164ddce0771SMatthew G. Knepley   const PetscReal s5_x  = s1_y*s0_z - s1_z*s0_y; /* s1 x s0 */
165ddce0771SMatthew G. Knepley   const PetscReal s5_y  = s1_z*s0_x - s1_x*s0_z;
166ddce0771SMatthew G. Knepley   const PetscReal s5_z  = s1_x*s0_y - s1_y*s0_x;
167ddce0771SMatthew G. Knepley   const PetscReal denom = -(s0_x*s3_x + s0_y*s3_y + s0_z*s3_z); /* -s0 . (s1 x s2) */
168ddce0771SMatthew G. Knepley 
169ddce0771SMatthew G. Knepley   PetscFunctionBegin;
170ddce0771SMatthew G. Knepley   *hasIntersection = PETSC_FALSE;
171ddce0771SMatthew G. Knepley   /* Line not parallel to plane */
172ddce0771SMatthew G. Knepley   if (denom != 0.0) {
173ddce0771SMatthew G. Knepley     const PetscReal t = (s3_x * (p0_x - q0_x) + s3_y * (p0_y - q0_y) + s3_z * (p0_z - q0_z)) / denom;
174ddce0771SMatthew G. Knepley     const PetscReal u = (s4_x * (p0_x - q0_x) + s4_y * (p0_y - q0_y) + s4_z * (p0_z - q0_z)) / denom;
175ddce0771SMatthew G. Knepley     const PetscReal v = (s5_x * (p0_x - q0_x) + s5_y * (p0_y - q0_y) + s5_z * (p0_z - q0_z)) / denom;
176ddce0771SMatthew G. Knepley 
177ddce0771SMatthew G. Knepley     if (t >= 0 && t <= 1 && u >= 0 && u <= 1 && v >= 0 && v <= 1) {
178ddce0771SMatthew G. Knepley       *hasIntersection = PETSC_TRUE;
179ddce0771SMatthew G. Knepley       if (intersection) {
180ddce0771SMatthew G. Knepley         intersection[0] = p0_x + (t * s0_x);
181ddce0771SMatthew G. Knepley         intersection[1] = p0_y + (t * s0_y);
182ddce0771SMatthew G. Knepley         intersection[2] = p0_z + (t * s0_z);
183ddce0771SMatthew G. Knepley       }
184ddce0771SMatthew G. Knepley     }
185ddce0771SMatthew G. Knepley   }
186ddce0771SMatthew G. Knepley   PetscFunctionReturn(0);
187ddce0771SMatthew G. Knepley }
188ddce0771SMatthew G. Knepley 
18914bbb9f0SLawrence Mitchell static PetscErrorCode DMPlexLocatePoint_Simplex_1D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
19014bbb9f0SLawrence Mitchell {
19114bbb9f0SLawrence Mitchell   const PetscReal eps = PETSC_SQRT_MACHINE_EPSILON;
19214bbb9f0SLawrence Mitchell   const PetscReal x   = PetscRealPart(point[0]);
19314bbb9f0SLawrence Mitchell   PetscReal       v0, J, invJ, detJ;
19414bbb9f0SLawrence Mitchell   PetscReal       xi;
19514bbb9f0SLawrence Mitchell   PetscErrorCode  ierr;
19614bbb9f0SLawrence Mitchell 
19714bbb9f0SLawrence Mitchell   PetscFunctionBegin;
19814bbb9f0SLawrence Mitchell   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, &v0, &J, &invJ, &detJ);CHKERRQ(ierr);
19914bbb9f0SLawrence Mitchell   xi   = invJ*(x - v0);
20014bbb9f0SLawrence Mitchell 
20114bbb9f0SLawrence Mitchell   if ((xi >= -eps) && (xi <= 2.+eps)) *cell = c;
20214bbb9f0SLawrence Mitchell   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
20314bbb9f0SLawrence Mitchell   PetscFunctionReturn(0);
20414bbb9f0SLawrence Mitchell }
20514bbb9f0SLawrence Mitchell 
206ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
207ccd2543fSMatthew G Knepley {
208ccd2543fSMatthew G Knepley   const PetscInt  embedDim = 2;
209f5ebc837SMatthew G. Knepley   const PetscReal eps      = PETSC_SQRT_MACHINE_EPSILON;
210ccd2543fSMatthew G Knepley   PetscReal       x        = PetscRealPart(point[0]);
211ccd2543fSMatthew G Knepley   PetscReal       y        = PetscRealPart(point[1]);
212ccd2543fSMatthew G Knepley   PetscReal       v0[2], J[4], invJ[4], detJ;
213ccd2543fSMatthew G Knepley   PetscReal       xi, eta;
214ccd2543fSMatthew G Knepley   PetscErrorCode  ierr;
215ccd2543fSMatthew G Knepley 
216ccd2543fSMatthew G Knepley   PetscFunctionBegin;
2178e0841e0SMatthew G. Knepley   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
218ccd2543fSMatthew G Knepley   xi  = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]);
219ccd2543fSMatthew G Knepley   eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]);
220ccd2543fSMatthew G Knepley 
221f5ebc837SMatthew G. Knepley   if ((xi >= -eps) && (eta >= -eps) && (xi + eta <= 2.0+eps)) *cell = c;
222c1496c66SMatthew G. Knepley   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
223ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
224ccd2543fSMatthew G Knepley }
225ccd2543fSMatthew G Knepley 
22662a38674SMatthew G. Knepley static PetscErrorCode DMPlexClosestPoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscReal cpoint[])
22762a38674SMatthew G. Knepley {
22862a38674SMatthew G. Knepley   const PetscInt  embedDim = 2;
22962a38674SMatthew G. Knepley   PetscReal       x        = PetscRealPart(point[0]);
23062a38674SMatthew G. Knepley   PetscReal       y        = PetscRealPart(point[1]);
23162a38674SMatthew G. Knepley   PetscReal       v0[2], J[4], invJ[4], detJ;
23262a38674SMatthew G. Knepley   PetscReal       xi, eta, r;
23362a38674SMatthew G. Knepley   PetscErrorCode  ierr;
23462a38674SMatthew G. Knepley 
23562a38674SMatthew G. Knepley   PetscFunctionBegin;
23662a38674SMatthew G. Knepley   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
23762a38674SMatthew G. Knepley   xi  = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]);
23862a38674SMatthew G. Knepley   eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]);
23962a38674SMatthew G. Knepley 
24062a38674SMatthew G. Knepley   xi  = PetscMax(xi,  0.0);
24162a38674SMatthew G. Knepley   eta = PetscMax(eta, 0.0);
24262a38674SMatthew G. Knepley   if (xi + eta > 2.0) {
24362a38674SMatthew G. Knepley     r    = (xi + eta)/2.0;
24462a38674SMatthew G. Knepley     xi  /= r;
24562a38674SMatthew G. Knepley     eta /= r;
24662a38674SMatthew G. Knepley   }
24762a38674SMatthew G. Knepley   cpoint[0] = J[0*embedDim+0]*xi + J[0*embedDim+1]*eta + v0[0];
24862a38674SMatthew G. Knepley   cpoint[1] = J[1*embedDim+0]*xi + J[1*embedDim+1]*eta + v0[1];
24962a38674SMatthew G. Knepley   PetscFunctionReturn(0);
25062a38674SMatthew G. Knepley }
25162a38674SMatthew G. Knepley 
252ba2698f1SMatthew G. Knepley static PetscErrorCode DMPlexLocatePoint_Quad_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
253ccd2543fSMatthew G Knepley {
254ccd2543fSMatthew G Knepley   PetscSection       coordSection;
255ccd2543fSMatthew G Knepley   Vec             coordsLocal;
256a1e44745SMatthew G. Knepley   PetscScalar    *coords = NULL;
257ccd2543fSMatthew G Knepley   const PetscInt  faces[8]  = {0, 1, 1, 2, 2, 3, 3, 0};
258ccd2543fSMatthew G Knepley   PetscReal       x         = PetscRealPart(point[0]);
259ccd2543fSMatthew G Knepley   PetscReal       y         = PetscRealPart(point[1]);
260ccd2543fSMatthew G Knepley   PetscInt        crossings = 0, f;
261ccd2543fSMatthew G Knepley   PetscErrorCode  ierr;
262ccd2543fSMatthew G Knepley 
263ccd2543fSMatthew G Knepley   PetscFunctionBegin;
264ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
26569d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
266ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
267ccd2543fSMatthew G Knepley   for (f = 0; f < 4; ++f) {
268ccd2543fSMatthew G Knepley     PetscReal x_i   = PetscRealPart(coords[faces[2*f+0]*2+0]);
269ccd2543fSMatthew G Knepley     PetscReal y_i   = PetscRealPart(coords[faces[2*f+0]*2+1]);
270ccd2543fSMatthew G Knepley     PetscReal x_j   = PetscRealPart(coords[faces[2*f+1]*2+0]);
271ccd2543fSMatthew G Knepley     PetscReal y_j   = PetscRealPart(coords[faces[2*f+1]*2+1]);
272ccd2543fSMatthew G Knepley     PetscReal slope = (y_j - y_i) / (x_j - x_i);
273ccd2543fSMatthew G Knepley     PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE;
274ccd2543fSMatthew G Knepley     PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE;
275ccd2543fSMatthew G Knepley     PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE;
276ccd2543fSMatthew G Knepley     if ((cond1 || cond2)  && above) ++crossings;
277ccd2543fSMatthew G Knepley   }
278ccd2543fSMatthew G Knepley   if (crossings % 2) *cell = c;
279c1496c66SMatthew G. Knepley   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
280ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
281ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
282ccd2543fSMatthew G Knepley }
283ccd2543fSMatthew G Knepley 
284ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
285ccd2543fSMatthew G Knepley {
286ccd2543fSMatthew G Knepley   const PetscInt  embedDim = 3;
28737900f7dSMatthew G. Knepley   const PetscReal eps      = PETSC_SQRT_MACHINE_EPSILON;
288ccd2543fSMatthew G Knepley   PetscReal       v0[3], J[9], invJ[9], detJ;
289ccd2543fSMatthew G Knepley   PetscReal       x = PetscRealPart(point[0]);
290ccd2543fSMatthew G Knepley   PetscReal       y = PetscRealPart(point[1]);
291ccd2543fSMatthew G Knepley   PetscReal       z = PetscRealPart(point[2]);
292ccd2543fSMatthew G Knepley   PetscReal       xi, eta, zeta;
293ccd2543fSMatthew G Knepley   PetscErrorCode  ierr;
294ccd2543fSMatthew G Knepley 
295ccd2543fSMatthew G Knepley   PetscFunctionBegin;
2968e0841e0SMatthew G. Knepley   ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
297ccd2543fSMatthew G Knepley   xi   = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]);
298ccd2543fSMatthew G Knepley   eta  = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]);
299ccd2543fSMatthew G Knepley   zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]);
300ccd2543fSMatthew G Knepley 
30137900f7dSMatthew G. Knepley   if ((xi >= -eps) && (eta >= -eps) && (zeta >= -eps) && (xi + eta + zeta <= 2.0+eps)) *cell = c;
302c1496c66SMatthew G. Knepley   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
303ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
304ccd2543fSMatthew G Knepley }
305ccd2543fSMatthew G Knepley 
306ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
307ccd2543fSMatthew G Knepley {
308ccd2543fSMatthew G Knepley   PetscSection   coordSection;
309ccd2543fSMatthew G Knepley   Vec            coordsLocal;
310872a9804SMatthew G. Knepley   PetscScalar   *coords = NULL;
311fb150da6SMatthew G. Knepley   const PetscInt faces[24] = {0, 3, 2, 1,  5, 4, 7, 6,  3, 0, 4, 5,
312fb150da6SMatthew G. Knepley                               1, 2, 6, 7,  3, 5, 6, 2,  0, 1, 7, 4};
313ccd2543fSMatthew G Knepley   PetscBool      found = PETSC_TRUE;
314ccd2543fSMatthew G Knepley   PetscInt       f;
315ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
316ccd2543fSMatthew G Knepley 
317ccd2543fSMatthew G Knepley   PetscFunctionBegin;
318ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
31969d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
320ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
321ccd2543fSMatthew G Knepley   for (f = 0; f < 6; ++f) {
322ccd2543fSMatthew G Knepley     /* Check the point is under plane */
323ccd2543fSMatthew G Knepley     /*   Get face normal */
324ccd2543fSMatthew G Knepley     PetscReal v_i[3];
325ccd2543fSMatthew G Knepley     PetscReal v_j[3];
326ccd2543fSMatthew G Knepley     PetscReal normal[3];
327ccd2543fSMatthew G Knepley     PetscReal pp[3];
328ccd2543fSMatthew G Knepley     PetscReal dot;
329ccd2543fSMatthew G Knepley 
330ccd2543fSMatthew G Knepley     v_i[0]    = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]);
331ccd2543fSMatthew G Knepley     v_i[1]    = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]);
332ccd2543fSMatthew G Knepley     v_i[2]    = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]);
333ccd2543fSMatthew G Knepley     v_j[0]    = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]);
334ccd2543fSMatthew G Knepley     v_j[1]    = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]);
335ccd2543fSMatthew G Knepley     v_j[2]    = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]);
336ccd2543fSMatthew G Knepley     normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1];
337ccd2543fSMatthew G Knepley     normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2];
338ccd2543fSMatthew G Knepley     normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0];
339ccd2543fSMatthew G Knepley     pp[0]     = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]);
340ccd2543fSMatthew G Knepley     pp[1]     = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]);
341ccd2543fSMatthew G Knepley     pp[2]     = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]);
342ccd2543fSMatthew G Knepley     dot       = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2];
343ccd2543fSMatthew G Knepley 
344ccd2543fSMatthew G Knepley     /* Check that projected point is in face (2D location problem) */
345ccd2543fSMatthew G Knepley     if (dot < 0.0) {
346ccd2543fSMatthew G Knepley       found = PETSC_FALSE;
347ccd2543fSMatthew G Knepley       break;
348ccd2543fSMatthew G Knepley     }
349ccd2543fSMatthew G Knepley   }
350ccd2543fSMatthew G Knepley   if (found) *cell = c;
351c1496c66SMatthew G. Knepley   else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
352ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr);
353ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
354ccd2543fSMatthew G Knepley }
355ccd2543fSMatthew G Knepley 
356c4eade1cSMatthew G. Knepley static PetscErrorCode PetscGridHashInitialize_Internal(PetscGridHash box, PetscInt dim, const PetscScalar point[])
357c4eade1cSMatthew G. Knepley {
358c4eade1cSMatthew G. Knepley   PetscInt d;
359c4eade1cSMatthew G. Knepley 
360c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
361c4eade1cSMatthew G. Knepley   box->dim = dim;
362c4eade1cSMatthew G. Knepley   for (d = 0; d < dim; ++d) box->lower[d] = box->upper[d] = PetscRealPart(point[d]);
363c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
364c4eade1cSMatthew G. Knepley }
365c4eade1cSMatthew G. Knepley 
366c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashCreate(MPI_Comm comm, PetscInt dim, const PetscScalar point[], PetscGridHash *box)
367c4eade1cSMatthew G. Knepley {
368c4eade1cSMatthew G. Knepley   PetscErrorCode ierr;
369c4eade1cSMatthew G. Knepley 
370c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
371c4eade1cSMatthew G. Knepley   ierr = PetscMalloc1(1, box);CHKERRQ(ierr);
372c4eade1cSMatthew G. Knepley   ierr = PetscGridHashInitialize_Internal(*box, dim, point);CHKERRQ(ierr);
373c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
374c4eade1cSMatthew G. Knepley }
375c4eade1cSMatthew G. Knepley 
376c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashEnlarge(PetscGridHash box, const PetscScalar point[])
377c4eade1cSMatthew G. Knepley {
378c4eade1cSMatthew G. Knepley   PetscInt d;
379c4eade1cSMatthew G. Knepley 
380c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
381c4eade1cSMatthew G. Knepley   for (d = 0; d < box->dim; ++d) {
382c4eade1cSMatthew G. Knepley     box->lower[d] = PetscMin(box->lower[d], PetscRealPart(point[d]));
383c4eade1cSMatthew G. Knepley     box->upper[d] = PetscMax(box->upper[d], PetscRealPart(point[d]));
384c4eade1cSMatthew G. Knepley   }
385c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
386c4eade1cSMatthew G. Knepley }
387c4eade1cSMatthew G. Knepley 
38862a38674SMatthew G. Knepley /*
38962a38674SMatthew G. Knepley   PetscGridHashSetGrid - Divide the grid into boxes
39062a38674SMatthew G. Knepley 
39162a38674SMatthew G. Knepley   Not collective
39262a38674SMatthew G. Knepley 
39362a38674SMatthew G. Knepley   Input Parameters:
39462a38674SMatthew G. Knepley + box - The grid hash object
39562a38674SMatthew G. Knepley . n   - The number of boxes in each dimension, or PETSC_DETERMINE
39662a38674SMatthew G. Knepley - h   - The box size in each dimension, only used if n[d] == PETSC_DETERMINE
39762a38674SMatthew G. Knepley 
39862a38674SMatthew G. Knepley   Level: developer
39962a38674SMatthew G. Knepley 
40062a38674SMatthew G. Knepley .seealso: PetscGridHashCreate()
40162a38674SMatthew G. Knepley */
402c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashSetGrid(PetscGridHash box, const PetscInt n[], const PetscReal h[])
403c4eade1cSMatthew G. Knepley {
404c4eade1cSMatthew G. Knepley   PetscInt d;
405c4eade1cSMatthew G. Knepley 
406c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
407c4eade1cSMatthew G. Knepley   for (d = 0; d < box->dim; ++d) {
408c4eade1cSMatthew G. Knepley     box->extent[d] = box->upper[d] - box->lower[d];
409c4eade1cSMatthew G. Knepley     if (n[d] == PETSC_DETERMINE) {
410c4eade1cSMatthew G. Knepley       box->h[d] = h[d];
411c4eade1cSMatthew G. Knepley       box->n[d] = PetscCeilReal(box->extent[d]/h[d]);
412c4eade1cSMatthew G. Knepley     } else {
413c4eade1cSMatthew G. Knepley       box->n[d] = n[d];
414c4eade1cSMatthew G. Knepley       box->h[d] = box->extent[d]/n[d];
415c4eade1cSMatthew G. Knepley     }
416c4eade1cSMatthew G. Knepley   }
417c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
418c4eade1cSMatthew G. Knepley }
419c4eade1cSMatthew G. Knepley 
42062a38674SMatthew G. Knepley /*
42162a38674SMatthew G. Knepley   PetscGridHashGetEnclosingBox - Find the grid boxes containing each input point
42262a38674SMatthew G. Knepley 
42362a38674SMatthew G. Knepley   Not collective
42462a38674SMatthew G. Knepley 
42562a38674SMatthew G. Knepley   Input Parameters:
42662a38674SMatthew G. Knepley + box       - The grid hash object
42762a38674SMatthew G. Knepley . numPoints - The number of input points
42862a38674SMatthew G. Knepley - points    - The input point coordinates
42962a38674SMatthew G. Knepley 
43062a38674SMatthew G. Knepley   Output Parameters:
43162a38674SMatthew G. Knepley + dboxes    - An array of numPoints*dim integers expressing the enclosing box as (i_0, i_1, ..., i_dim)
43262a38674SMatthew G. Knepley - boxes     - An array of numPoints integers expressing the enclosing box as single number, or NULL
43362a38674SMatthew G. Knepley 
43462a38674SMatthew G. Knepley   Level: developer
43562a38674SMatthew G. Knepley 
43662a38674SMatthew G. Knepley .seealso: PetscGridHashCreate()
43762a38674SMatthew G. Knepley */
4381c6dfc3eSMatthew G. Knepley PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[])
439c4eade1cSMatthew G. Knepley {
440c4eade1cSMatthew G. Knepley   const PetscReal *lower = box->lower;
441c4eade1cSMatthew G. Knepley   const PetscReal *upper = box->upper;
442c4eade1cSMatthew G. Knepley   const PetscReal *h     = box->h;
443c4eade1cSMatthew G. Knepley   const PetscInt  *n     = box->n;
444c4eade1cSMatthew G. Knepley   const PetscInt   dim   = box->dim;
445c4eade1cSMatthew G. Knepley   PetscInt         d, p;
446c4eade1cSMatthew G. Knepley 
447c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
448c4eade1cSMatthew G. Knepley   for (p = 0; p < numPoints; ++p) {
449c4eade1cSMatthew G. Knepley     for (d = 0; d < dim; ++d) {
4501c6dfc3eSMatthew G. Knepley       PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]);
451c4eade1cSMatthew G. Knepley 
4521c6dfc3eSMatthew G. Knepley       if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1;
4532a705cacSMatthew G. Knepley       if (dbox == -1   && PetscAbsReal(PetscRealPart(points[p*dim+d]) - lower[d]) < 1.0e-9) dbox = 0;
454c4eade1cSMatthew 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",
455087ef6b2SMatthew 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);
456c4eade1cSMatthew G. Knepley       dboxes[p*dim+d] = dbox;
457c4eade1cSMatthew G. Knepley     }
458ddce0771SMatthew G. Knepley     if (boxes) for (d = dim-2, boxes[p] = dboxes[p*dim+dim-1]; d >= 0; --d) boxes[p] = boxes[p]*n[d] + dboxes[p*dim+d];
459c4eade1cSMatthew G. Knepley   }
460c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
461c4eade1cSMatthew G. Knepley }
462c4eade1cSMatthew G. Knepley 
463af74b616SDave May /*
464af74b616SDave May  PetscGridHashGetEnclosingBoxQuery - Find the grid boxes containing each input point
465af74b616SDave May 
466af74b616SDave May  Not collective
467af74b616SDave May 
468af74b616SDave May   Input Parameters:
469af74b616SDave May + box       - The grid hash object
470af74b616SDave May . numPoints - The number of input points
471af74b616SDave May - points    - The input point coordinates
472af74b616SDave May 
473af74b616SDave May   Output Parameters:
474af74b616SDave May + dboxes    - An array of numPoints*dim integers expressing the enclosing box as (i_0, i_1, ..., i_dim)
475af74b616SDave May . boxes     - An array of numPoints integers expressing the enclosing box as single number, or NULL
476af74b616SDave May - found     - Flag indicating if point was located within a box
477af74b616SDave May 
478af74b616SDave May   Level: developer
479af74b616SDave May 
480af74b616SDave May .seealso: PetscGridHashGetEnclosingBox()
481af74b616SDave May */
482af74b616SDave May PetscErrorCode PetscGridHashGetEnclosingBoxQuery(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[],PetscBool *found)
483af74b616SDave May {
484af74b616SDave May   const PetscReal *lower = box->lower;
485af74b616SDave May   const PetscReal *upper = box->upper;
486af74b616SDave May   const PetscReal *h     = box->h;
487af74b616SDave May   const PetscInt  *n     = box->n;
488af74b616SDave May   const PetscInt   dim   = box->dim;
489af74b616SDave May   PetscInt         d, p;
490af74b616SDave May 
491af74b616SDave May   PetscFunctionBegin;
492af74b616SDave May   *found = PETSC_FALSE;
493af74b616SDave May   for (p = 0; p < numPoints; ++p) {
494af74b616SDave May     for (d = 0; d < dim; ++d) {
495af74b616SDave May       PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]);
496af74b616SDave May 
497af74b616SDave May       if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1;
498af74b616SDave May       if (dbox < 0 || dbox >= n[d]) {
499af74b616SDave May         PetscFunctionReturn(0);
500af74b616SDave May       }
501af74b616SDave May       dboxes[p*dim+d] = dbox;
502af74b616SDave May     }
503ddce0771SMatthew G. Knepley     if (boxes) for (d = dim-2, boxes[p] = dboxes[p*dim+dim-1]; d >= 0; --d) boxes[p] = boxes[p]*n[d] + dboxes[p*dim+d];
504af74b616SDave May   }
505af74b616SDave May   *found = PETSC_TRUE;
506af74b616SDave May   PetscFunctionReturn(0);
507af74b616SDave May }
508af74b616SDave May 
509c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashDestroy(PetscGridHash *box)
510c4eade1cSMatthew G. Knepley {
511c4eade1cSMatthew G. Knepley   PetscErrorCode ierr;
512c4eade1cSMatthew G. Knepley 
513c4eade1cSMatthew G. Knepley   PetscFunctionBegin;
514c4eade1cSMatthew G. Knepley   if (*box) {
515c4eade1cSMatthew G. Knepley     ierr = PetscSectionDestroy(&(*box)->cellSection);CHKERRQ(ierr);
516c4eade1cSMatthew G. Knepley     ierr = ISDestroy(&(*box)->cells);CHKERRQ(ierr);
517c4eade1cSMatthew G. Knepley     ierr = DMLabelDestroy(&(*box)->cellsSparse);CHKERRQ(ierr);
518c4eade1cSMatthew G. Knepley   }
519c4eade1cSMatthew G. Knepley   ierr = PetscFree(*box);CHKERRQ(ierr);
520c4eade1cSMatthew G. Knepley   PetscFunctionReturn(0);
521c4eade1cSMatthew G. Knepley }
522c4eade1cSMatthew G. Knepley 
523cafe43deSMatthew G. Knepley PetscErrorCode DMPlexLocatePoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cellStart, PetscInt *cell)
524cafe43deSMatthew G. Knepley {
525ba2698f1SMatthew G. Knepley   DMPolytopeType ct;
526cafe43deSMatthew G. Knepley   PetscErrorCode ierr;
527cafe43deSMatthew G. Knepley 
528cafe43deSMatthew G. Knepley   PetscFunctionBegin;
529ba2698f1SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, cellStart, &ct);CHKERRQ(ierr);
530ba2698f1SMatthew G. Knepley   switch (ct) {
53114bbb9f0SLawrence Mitchell     case DM_POLYTOPE_SEGMENT:
53214bbb9f0SLawrence Mitchell     ierr = DMPlexLocatePoint_Simplex_1D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);break;
533ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:
534ba2698f1SMatthew G. Knepley     ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);break;
535ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_QUADRILATERAL:
536ba2698f1SMatthew G. Knepley     ierr = DMPlexLocatePoint_Quad_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);break;
537ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_TETRAHEDRON:
538ba2698f1SMatthew G. Knepley     ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);break;
539ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_HEXAHEDRON:
540ba2698f1SMatthew G. Knepley     ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr);break;
541ba2698f1SMatthew G. Knepley     default: SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell %D with type %s", cellStart, DMPolytopeTypes[ct]);
542cafe43deSMatthew G. Knepley   }
543cafe43deSMatthew G. Knepley   PetscFunctionReturn(0);
544cafe43deSMatthew G. Knepley }
545cafe43deSMatthew G. Knepley 
54662a38674SMatthew G. Knepley /*
54762a38674SMatthew G. Knepley   DMPlexClosestPoint_Internal - Returns the closest point in the cell to the given point
54862a38674SMatthew G. Knepley */
54962a38674SMatthew G. Knepley PetscErrorCode DMPlexClosestPoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cell, PetscReal cpoint[])
55062a38674SMatthew G. Knepley {
551ba2698f1SMatthew G. Knepley   DMPolytopeType ct;
55262a38674SMatthew G. Knepley   PetscErrorCode ierr;
55362a38674SMatthew G. Knepley 
55462a38674SMatthew G. Knepley   PetscFunctionBegin;
555ba2698f1SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, cell, &ct);CHKERRQ(ierr);
556ba2698f1SMatthew G. Knepley   switch (ct) {
557ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:
558ba2698f1SMatthew G. Knepley     ierr = DMPlexClosestPoint_Simplex_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);break;
55962a38674SMatthew G. Knepley #if 0
560ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_QUADRILATERAL:
561ba2698f1SMatthew G. Knepley     ierr = DMPlexClosestPoint_General_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);break;
562ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_TETRAHEDRON:
563ba2698f1SMatthew G. Knepley     ierr = DMPlexClosestPoint_Simplex_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);break;
564ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_HEXAHEDRON:
565ba2698f1SMatthew G. Knepley     ierr = DMPlexClosestPoint_General_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr);break;
56662a38674SMatthew G. Knepley #endif
567ba2698f1SMatthew G. Knepley     default: SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell %D with type %s", cell, DMPolytopeTypes[ct]);
56862a38674SMatthew G. Knepley   }
56962a38674SMatthew G. Knepley   PetscFunctionReturn(0);
57062a38674SMatthew G. Knepley }
57162a38674SMatthew G. Knepley 
57262a38674SMatthew G. Knepley /*
57362a38674SMatthew G. Knepley   DMPlexComputeGridHash_Internal - Create a grid hash structure covering the Plex
57462a38674SMatthew G. Knepley 
575d083f849SBarry Smith   Collective on dm
57662a38674SMatthew G. Knepley 
57762a38674SMatthew G. Knepley   Input Parameter:
57862a38674SMatthew G. Knepley . dm - The Plex
57962a38674SMatthew G. Knepley 
58062a38674SMatthew G. Knepley   Output Parameter:
58162a38674SMatthew G. Knepley . localBox - The grid hash object
58262a38674SMatthew G. Knepley 
58362a38674SMatthew G. Knepley   Level: developer
58462a38674SMatthew G. Knepley 
58562a38674SMatthew G. Knepley .seealso: PetscGridHashCreate(), PetscGridHashGetEnclosingBox()
58662a38674SMatthew G. Knepley */
587cafe43deSMatthew G. Knepley PetscErrorCode DMPlexComputeGridHash_Internal(DM dm, PetscGridHash *localBox)
588cafe43deSMatthew G. Knepley {
589ddce0771SMatthew G. Knepley   const PetscInt     debug = 0;
590cafe43deSMatthew G. Knepley   MPI_Comm           comm;
591cafe43deSMatthew G. Knepley   PetscGridHash      lbox;
592cafe43deSMatthew G. Knepley   Vec                coordinates;
593cafe43deSMatthew G. Knepley   PetscSection       coordSection;
594cafe43deSMatthew G. Knepley   Vec                coordsLocal;
595cafe43deSMatthew G. Knepley   const PetscScalar *coords;
596ddce0771SMatthew G. Knepley   PetscScalar       *edgeCoords;
597722d0f5cSMatthew G. Knepley   PetscInt          *dboxes, *boxes;
598ddce0771SMatthew G. Knepley   PetscInt           n[3] = {2, 2, 2};
599ddce0771SMatthew G. Knepley   PetscInt           dim, N, maxConeSize, cStart, cEnd, c, eStart, eEnd, i;
600ddce0771SMatthew G. Knepley   PetscBool          flg;
601cafe43deSMatthew G. Knepley   PetscErrorCode     ierr;
602cafe43deSMatthew G. Knepley 
603cafe43deSMatthew G. Knepley   PetscFunctionBegin;
604cafe43deSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
605cafe43deSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
606cafe43deSMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
607ddce0771SMatthew G. Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr);
608ddce0771SMatthew G. Knepley   ierr = DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
609cafe43deSMatthew G. Knepley   ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr);
610cafe43deSMatthew G. Knepley   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
611cafe43deSMatthew G. Knepley   ierr = PetscGridHashCreate(comm, dim, coords, &lbox);CHKERRQ(ierr);
612cafe43deSMatthew G. Knepley   for (i = 0; i < N; i += dim) {ierr = PetscGridHashEnlarge(lbox, &coords[i]);CHKERRQ(ierr);}
613cafe43deSMatthew G. Knepley   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
614ddce0771SMatthew G. Knepley   c    = dim;
615ddce0771SMatthew G. Knepley   ierr = PetscOptionsGetIntArray(NULL, ((PetscObject) dm)->prefix, "-dm_plex_hash_box_faces", n, &c, &flg);CHKERRQ(ierr);
616ddce0771SMatthew G. Knepley   if (flg) {for (i = c; i < dim; ++i) n[i] = n[c-1];}
617ddce0771SMatthew G. Knepley   else     {for (i = 0; i < dim; ++i) n[i] = PetscMax(2, PetscFloorReal(PetscPowReal((PetscReal) (cEnd - cStart), 1.0/dim) * 0.8));}
618cafe43deSMatthew G. Knepley   ierr = PetscGridHashSetGrid(lbox, n, NULL);CHKERRQ(ierr);
619cafe43deSMatthew G. Knepley #if 0
620cafe43deSMatthew G. Knepley   /* Could define a custom reduction to merge these */
621820f2d46SBarry Smith   ierr = MPIU_Allreduce(lbox->lower, gbox->lower, 3, MPIU_REAL, MPI_MIN, comm);CHKERRMPI(ierr);
622820f2d46SBarry Smith   ierr = MPIU_Allreduce(lbox->upper, gbox->upper, 3, MPIU_REAL, MPI_MAX, comm);CHKERRMPI(ierr);
623cafe43deSMatthew G. Knepley #endif
624cafe43deSMatthew G. Knepley   /* Is there a reason to snap the local bounding box to a division of the global box? */
625cafe43deSMatthew G. Knepley   /* Should we compute all overlaps of local boxes? We could do this with a rendevouz scheme partitioning the global box */
626cafe43deSMatthew G. Knepley   /* Create label */
627ddce0771SMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr);
628d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, "cells", &lbox->cellsSparse);CHKERRQ(ierr);
629cafe43deSMatthew G. Knepley   ierr = DMLabelCreateIndex(lbox->cellsSparse, cStart, cEnd);CHKERRQ(ierr);
630a8d69d7bSBarry Smith   /* Compute boxes which overlap each cell: https://stackoverflow.com/questions/13790208/triangle-square-intersection-test-in-2d */
631cafe43deSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
632cafe43deSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
633ddce0771SMatthew G. Knepley   ierr = PetscCalloc3(16 * dim, &dboxes, 16, &boxes, PetscPowInt(maxConeSize, dim) * dim, &edgeCoords);CHKERRQ(ierr);
634cafe43deSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
635cafe43deSMatthew G. Knepley     const PetscReal *h       = lbox->h;
636cafe43deSMatthew G. Knepley     PetscScalar     *ccoords = NULL;
63738353de4SMatthew G. Knepley     PetscInt         csize   = 0;
638ddce0771SMatthew G. Knepley     PetscInt        *closure = NULL;
639ddce0771SMatthew G. Knepley     PetscInt         Ncl, cl, Ne = 0;
640cafe43deSMatthew G. Knepley     PetscScalar      point[3];
641cafe43deSMatthew G. Knepley     PetscInt         dlim[6], d, e, i, j, k;
642cafe43deSMatthew G. Knepley 
643ddce0771SMatthew G. Knepley     /* Get all edges in cell */
644ddce0771SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr);
645ddce0771SMatthew G. Knepley     for (cl = 0; cl < Ncl*2; ++cl) {
646ddce0771SMatthew G. Knepley       if ((closure[cl] >= eStart) && (closure[cl] < eEnd)) {
647ddce0771SMatthew G. Knepley         PetscScalar *ecoords = &edgeCoords[Ne*dim*2];
648ddce0771SMatthew G. Knepley         PetscInt     ecsize  = dim*2;
649ddce0771SMatthew G. Knepley 
650e032ffc8SMatthew G. Knepley         ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, closure[cl], &ecsize, &ecoords);CHKERRQ(ierr);
651ddce0771SMatthew G. Knepley         if (ecsize != dim*2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Got %D coords for edge, instead of %D", ecsize, dim*2);
652ddce0771SMatthew G. Knepley         ++Ne;
653ddce0771SMatthew G. Knepley       }
654ddce0771SMatthew G. Knepley     }
655ddce0771SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &Ncl, &closure);CHKERRQ(ierr);
656cafe43deSMatthew G. Knepley     /* Find boxes enclosing each vertex */
65738353de4SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, &csize, &ccoords);CHKERRQ(ierr);
65838353de4SMatthew G. Knepley     ierr = PetscGridHashGetEnclosingBox(lbox, csize/dim, ccoords, dboxes, boxes);CHKERRQ(ierr);
659722d0f5cSMatthew G. Knepley     /* Mark cells containing the vertices */
660ddce0771SMatthew G. Knepley     for (e = 0; e < csize/dim; ++e) {
661ddce0771SMatthew G. Knepley       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D has vertex in box %D (%D, %D, %D)\n", c, boxes[e], dboxes[e*dim+0], dim > 1 ? dboxes[e*dim+1] : -1, dim > 2 ? dboxes[e*dim+2] : -1);CHKERRQ(ierr);}
662ddce0771SMatthew G. Knepley       ierr = DMLabelSetValue(lbox->cellsSparse, c, boxes[e]);CHKERRQ(ierr);
663ddce0771SMatthew G. Knepley     }
664cafe43deSMatthew G. Knepley     /* Get grid of boxes containing these */
665cafe43deSMatthew G. Knepley     for (d = 0;   d < dim; ++d) {dlim[d*2+0] = dlim[d*2+1] = dboxes[d];}
6662291669eSMatthew G. Knepley     for (d = dim; d < 3;   ++d) {dlim[d*2+0] = dlim[d*2+1] = 0;}
667cafe43deSMatthew G. Knepley     for (e = 1; e < dim+1; ++e) {
668cafe43deSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
669cafe43deSMatthew G. Knepley         dlim[d*2+0] = PetscMin(dlim[d*2+0], dboxes[e*dim+d]);
670cafe43deSMatthew G. Knepley         dlim[d*2+1] = PetscMax(dlim[d*2+1], dboxes[e*dim+d]);
671cafe43deSMatthew G. Knepley       }
672cafe43deSMatthew G. Knepley     }
673fea14342SMatthew G. Knepley     /* Check for intersection of box with cell */
674cafe43deSMatthew 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]) {
675cafe43deSMatthew 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]) {
676cafe43deSMatthew 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]) {
677cafe43deSMatthew G. Knepley           const PetscInt box = (k*lbox->n[1] + j)*lbox->n[0] + i;
678cafe43deSMatthew G. Knepley           PetscScalar    cpoint[3];
679fea14342SMatthew G. Knepley           PetscInt       cell, edge, ii, jj, kk;
680cafe43deSMatthew G. Knepley 
681ddce0771SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "Box %D: (%.2g, %.2g, %.2g) -- (%.2g, %.2g, %.2g)\n", box, point[0], point[1], point[2], point[0] + h[0], point[1] + h[1], point[2] + h[2]);CHKERRQ(ierr);}
682ddce0771SMatthew G. Knepley           /* Check whether cell contains any vertex of this subbox TODO vectorize this */
683cafe43deSMatthew G. Knepley           for (kk = 0, cpoint[2] = point[2]; kk < (dim > 2 ? 2 : 1); ++kk, cpoint[2] += h[2]) {
684cafe43deSMatthew G. Knepley             for (jj = 0, cpoint[1] = point[1]; jj < (dim > 1 ? 2 : 1); ++jj, cpoint[1] += h[1]) {
685cafe43deSMatthew G. Knepley               for (ii = 0, cpoint[0] = point[0]; ii < 2; ++ii, cpoint[0] += h[0]) {
686cafe43deSMatthew G. Knepley 
687cafe43deSMatthew G. Knepley                 ierr = DMPlexLocatePoint_Internal(dm, dim, cpoint, c, &cell);CHKERRQ(ierr);
6880b6bfacdSStefano Zampini                 if (cell >= 0) {
689ddce0771SMatthew G. Knepley                   if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  Cell %D contains vertex (%.2g, %.2g, %.2g) of box %D\n", c, cpoint[0], cpoint[1], cpoint[2], box);CHKERRQ(ierr);}
6900b6bfacdSStefano Zampini                   ierr = DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr);
6910b6bfacdSStefano Zampini                   jj = kk = 2;
6920b6bfacdSStefano Zampini                   break;
6930b6bfacdSStefano Zampini                 }
694cafe43deSMatthew G. Knepley               }
695cafe43deSMatthew G. Knepley             }
696cafe43deSMatthew G. Knepley           }
697ddce0771SMatthew G. Knepley           /* Check whether cell edge intersects any face of these subboxes TODO vectorize this */
698ddce0771SMatthew G. Knepley           for (edge = 0; edge < Ne; ++edge) {
699a5cae605SSatish Balay             PetscReal segA[6] = {0.,0.,0.,0.,0.,0.};
700a5cae605SSatish Balay             PetscReal segB[6] = {0.,0.,0.,0.,0.,0.};
701a5cae605SSatish Balay             PetscReal segC[6] = {0.,0.,0.,0.,0.,0.};
702fea14342SMatthew G. Knepley 
703aab5bcd8SJed Brown             if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected dim %d > 3",dim);
704ddce0771SMatthew G. Knepley             for (d = 0; d < dim*2; ++d) segA[d] = PetscRealPart(edgeCoords[edge*dim*2+d]);
705ddce0771SMatthew G. Knepley             /* 1D: (x) -- (x+h)               0 -- 1
706ddce0771SMatthew G. Knepley                2D: (x,   y)   -- (x,   y+h)   (0, 0) -- (0, 1)
707ddce0771SMatthew G. Knepley                    (x+h, y)   -- (x+h, y+h)   (1, 0) -- (1, 1)
708ddce0771SMatthew G. Knepley                    (x,   y)   -- (x+h, y)     (0, 0) -- (1, 0)
709ddce0771SMatthew G. Knepley                    (x,   y+h) -- (x+h, y+h)   (0, 1) -- (1, 1)
710ddce0771SMatthew G. Knepley                3D: (x,   y,   z)   -- (x,   y+h, z),   (x,   y,   z)   -- (x,   y,   z+h) (0, 0, 0) -- (0, 1, 0), (0, 0, 0) -- (0, 0, 1)
711ddce0771SMatthew G. Knepley                    (x+h, y,   z)   -- (x+h, y+h, z),   (x+h, y,   z)   -- (x+h, y,   z+h) (1, 0, 0) -- (1, 1, 0), (1, 0, 0) -- (1, 0, 1)
712ddce0771SMatthew G. Knepley                    (x,   y,   z)   -- (x+h, y,   z),   (x,   y,   z)   -- (x,   y,   z+h) (0, 0, 0) -- (1, 0, 0), (0, 0, 0) -- (0, 0, 1)
713ddce0771SMatthew G. Knepley                    (x,   y+h, z)   -- (x+h, y+h, z),   (x,   y+h, z)   -- (x,   y+h, z+h) (0, 1, 0) -- (1, 1, 0), (0, 1, 0) -- (0, 1, 1)
714ddce0771SMatthew G. Knepley                    (x,   y,   z)   -- (x+h, y,   z),   (x,   y,   z)   -- (x,   y+h, z)   (0, 0, 0) -- (1, 0, 0), (0, 0, 0) -- (0, 1, 0)
715ddce0771SMatthew G. Knepley                    (x,   y,   z+h) -- (x+h, y,   z+h), (x,   y,   z+h) -- (x,   y+h, z+h) (0, 0, 1) -- (1, 0, 1), (0, 0, 1) -- (0, 1, 1)
716ddce0771SMatthew G. Knepley              */
717ddce0771SMatthew G. Knepley             /* Loop over faces with normal in direction d */
718ddce0771SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
719ddce0771SMatthew G. Knepley               PetscBool intersects = PETSC_FALSE;
720ddce0771SMatthew G. Knepley               PetscInt  e = (d+1)%dim;
721ddce0771SMatthew G. Knepley               PetscInt  f = (d+2)%dim;
722ddce0771SMatthew G. Knepley 
723ddce0771SMatthew G. Knepley               /* There are two faces in each dimension */
724ddce0771SMatthew G. Knepley               for (ii = 0; ii < 2; ++ii) {
725ddce0771SMatthew G. Knepley                 segB[d]     = PetscRealPart(point[d] + ii*h[d]);
726ddce0771SMatthew G. Knepley                 segB[dim+d] = PetscRealPart(point[d] + ii*h[d]);
727ddce0771SMatthew G. Knepley                 segC[d]     = PetscRealPart(point[d] + ii*h[d]);
728ddce0771SMatthew G. Knepley                 segC[dim+d] = PetscRealPart(point[d] + ii*h[d]);
729ddce0771SMatthew G. Knepley                 if (dim > 1) {
730ddce0771SMatthew G. Knepley                   segB[e]     = PetscRealPart(point[e] + 0*h[e]);
731ddce0771SMatthew G. Knepley                   segB[dim+e] = PetscRealPart(point[e] + 1*h[e]);
732ddce0771SMatthew G. Knepley                   segC[e]     = PetscRealPart(point[e] + 0*h[e]);
733ddce0771SMatthew G. Knepley                   segC[dim+e] = PetscRealPart(point[e] + 0*h[e]);
734ddce0771SMatthew G. Knepley                 }
735ddce0771SMatthew G. Knepley                 if (dim > 2) {
736ddce0771SMatthew G. Knepley                   segB[f]     = PetscRealPart(point[f] + 0*h[f]);
737ddce0771SMatthew G. Knepley                   segB[dim+f] = PetscRealPart(point[f] + 0*h[f]);
738ddce0771SMatthew G. Knepley                   segC[f]     = PetscRealPart(point[f] + 0*h[f]);
739ddce0771SMatthew G. Knepley                   segC[dim+f] = PetscRealPart(point[f] + 1*h[f]);
740ddce0771SMatthew G. Knepley                 }
741ddce0771SMatthew G. Knepley                 if (dim == 2) {
742ddce0771SMatthew G. Knepley                   ierr = DMPlexGetLineIntersection_2D_Internal(segA, segB, NULL, &intersects);CHKERRQ(ierr);
743ddce0771SMatthew G. Knepley                 } else if (dim == 3) {
744ddce0771SMatthew G. Knepley                   ierr = DMPlexGetLinePlaneIntersection_3D_Internal(segA, segB, segC, NULL, &intersects);CHKERRQ(ierr);
745ddce0771SMatthew G. Knepley                 }
746ddce0771SMatthew G. Knepley                 if (intersects) {
747ddce0771SMatthew G. Knepley                   if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  Cell %D edge %D (%.2g, %.2g, %.2g)--(%.2g, %.2g, %.2g) intersects box %D, face (%.2g, %.2g, %.2g)--(%.2g, %.2g, %.2g) (%.2g, %.2g, %.2g)--(%.2g, %.2g, %.2g)\n", c, edge, segA[0], segA[1], segA[2], segA[3], segA[4], segA[5], box, segB[0], segB[1], segB[2], segB[3], segB[4], segB[5], segC[0], segC[1], segC[2], segC[3], segC[4], segC[5]);CHKERRQ(ierr);}
748ddce0771SMatthew G. Knepley                   ierr = DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); edge = Ne; break;
749ddce0771SMatthew G. Knepley                 }
750ddce0771SMatthew G. Knepley               }
751ddce0771SMatthew G. Knepley             }
752cafe43deSMatthew G. Knepley           }
753fea14342SMatthew G. Knepley         }
754fea14342SMatthew G. Knepley       }
755fea14342SMatthew G. Knepley     }
756fea14342SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &ccoords);CHKERRQ(ierr);
757fea14342SMatthew G. Knepley   }
758ddce0771SMatthew G. Knepley   ierr = PetscFree3(dboxes, boxes, edgeCoords);CHKERRQ(ierr);
759ddce0771SMatthew G. Knepley   if (debug) {ierr = DMLabelView(lbox->cellsSparse, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
760cafe43deSMatthew G. Knepley   ierr = DMLabelConvertToSection(lbox->cellsSparse, &lbox->cellSection, &lbox->cells);CHKERRQ(ierr);
761cafe43deSMatthew G. Knepley   ierr = DMLabelDestroy(&lbox->cellsSparse);CHKERRQ(ierr);
762cafe43deSMatthew G. Knepley   *localBox = lbox;
763cafe43deSMatthew G. Knepley   PetscFunctionReturn(0);
764cafe43deSMatthew G. Knepley }
765cafe43deSMatthew G. Knepley 
76662a38674SMatthew G. Knepley PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, DMPointLocationType ltype, PetscSF cellSF)
767ccd2543fSMatthew G Knepley {
768ddce0771SMatthew G. Knepley   const PetscInt  debug = 0;
769cafe43deSMatthew G. Knepley   DM_Plex        *mesh = (DM_Plex *) dm->data;
770af74b616SDave May   PetscBool       hash = mesh->useHashLocation, reuse = PETSC_FALSE;
7713a93e3b7SToby Isaac   PetscInt        bs, numPoints, p, numFound, *found = NULL;
772412e9a14SMatthew G. Knepley   PetscInt        dim, cStart, cEnd, numCells, c, d;
773cafe43deSMatthew G. Knepley   const PetscInt *boxCells;
7743a93e3b7SToby Isaac   PetscSFNode    *cells;
775ccd2543fSMatthew G Knepley   PetscScalar    *a;
7763a93e3b7SToby Isaac   PetscMPIInt     result;
777af74b616SDave May   PetscLogDouble  t0,t1;
7789cb35068SDave May   PetscReal       gmin[3],gmax[3];
7799cb35068SDave May   PetscInt        terminating_query_type[] = { 0, 0, 0 };
780ccd2543fSMatthew G Knepley   PetscErrorCode  ierr;
781ccd2543fSMatthew G Knepley 
782ccd2543fSMatthew G Knepley   PetscFunctionBegin;
783cadf77a0SMark Adams   ierr = PetscLogEventBegin(DMPLEX_LocatePoints,0,0,0,0);CHKERRQ(ierr);
784af74b616SDave May   ierr = PetscTime(&t0);CHKERRQ(ierr);
785080342d1SMatthew 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.");
786cafe43deSMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
787cafe43deSMatthew G. Knepley   ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr);
788ffc4695bSBarry Smith   ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)cellSF),PETSC_COMM_SELF,&result);CHKERRMPI(ierr);
7893a93e3b7SToby Isaac   if (result != MPI_IDENT && result != MPI_CONGRUENT) SETERRQ(PetscObjectComm((PetscObject)cellSF),PETSC_ERR_SUP, "Trying parallel point location: only local point location supported");
790cafe43deSMatthew 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);
791412e9a14SMatthew G. Knepley   ierr = DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
792ccd2543fSMatthew G Knepley   ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr);
793ccd2543fSMatthew G Knepley   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
794ccd2543fSMatthew G Knepley   numPoints /= bs;
795af74b616SDave May   {
796af74b616SDave May     const PetscSFNode *sf_cells;
797af74b616SDave May 
798af74b616SDave May     ierr = PetscSFGetGraph(cellSF,NULL,NULL,NULL,&sf_cells);CHKERRQ(ierr);
799af74b616SDave May     if (sf_cells) {
800af74b616SDave May       ierr = PetscInfo(dm,"[DMLocatePoints_Plex] Re-using existing StarForest node list\n");CHKERRQ(ierr);
801af74b616SDave May       cells = (PetscSFNode*)sf_cells;
802af74b616SDave May       reuse = PETSC_TRUE;
803af74b616SDave May     } else {
804af74b616SDave May       ierr = PetscInfo(dm,"[DMLocatePoints_Plex] Creating and initializing new StarForest node list\n");CHKERRQ(ierr);
805785e854fSJed Brown       ierr = PetscMalloc1(numPoints, &cells);CHKERRQ(ierr);
806af74b616SDave May       /* initialize cells if created */
807af74b616SDave May       for (p=0; p<numPoints; p++) {
808af74b616SDave May         cells[p].rank  = 0;
809af74b616SDave May         cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
810af74b616SDave May       }
811af74b616SDave May     }
812af74b616SDave May   }
8139cb35068SDave May   /* define domain bounding box */
8149cb35068SDave May   {
8159cb35068SDave May     Vec coorglobal;
8169cb35068SDave May 
8179cb35068SDave May     ierr = DMGetCoordinates(dm,&coorglobal);CHKERRQ(ierr);
8189cb35068SDave May     ierr = VecStrideMaxAll(coorglobal,NULL,gmax);CHKERRQ(ierr);
8199cb35068SDave May     ierr = VecStrideMinAll(coorglobal,NULL,gmin);CHKERRQ(ierr);
8209cb35068SDave May   }
821953fc75cSMatthew G. Knepley   if (hash) {
822ac6ec2abSMatthew G. Knepley     if (!mesh->lbox) {ierr = PetscInfo(dm, "Initializing grid hashing");CHKERRQ(ierr);ierr = DMPlexComputeGridHash_Internal(dm, &mesh->lbox);CHKERRQ(ierr);}
823cafe43deSMatthew G. Knepley     /* Designate the local box for each point */
824cafe43deSMatthew G. Knepley     /* Send points to correct process */
825cafe43deSMatthew G. Knepley     /* Search cells that lie in each subbox */
826cafe43deSMatthew G. Knepley     /*   Should we bin points before doing search? */
827cafe43deSMatthew G. Knepley     ierr = ISGetIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr);
828953fc75cSMatthew G. Knepley   }
8293a93e3b7SToby Isaac   for (p = 0, numFound = 0; p < numPoints; ++p) {
830ccd2543fSMatthew G Knepley     const PetscScalar *point = &a[p*bs];
831e56f9228SJed Brown     PetscInt           dbin[3] = {-1,-1,-1}, bin, cell = -1, cellOffset;
8329cb35068SDave May     PetscBool          point_outside_domain = PETSC_FALSE;
833ccd2543fSMatthew G Knepley 
8349cb35068SDave May     /* check bounding box of domain */
8359cb35068SDave May     for (d=0; d<dim; d++) {
836a5f152d1SDave May       if (PetscRealPart(point[d]) < gmin[d]) { point_outside_domain = PETSC_TRUE; break; }
837a5f152d1SDave May       if (PetscRealPart(point[d]) > gmax[d]) { point_outside_domain = PETSC_TRUE; break; }
8389cb35068SDave May     }
8399cb35068SDave May     if (point_outside_domain) {
840e9b685f5SMatthew G. Knepley       cells[p].rank = 0;
841e9b685f5SMatthew G. Knepley       cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
8429cb35068SDave May       terminating_query_type[0]++;
8439cb35068SDave May       continue;
8449cb35068SDave May     }
845ccd2543fSMatthew G Knepley 
846af74b616SDave May     /* check initial values in cells[].index - abort early if found */
847af74b616SDave May     if (cells[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
848af74b616SDave May       c = cells[p].index;
8493a93e3b7SToby Isaac       cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
850af74b616SDave May       ierr = DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);CHKERRQ(ierr);
851af74b616SDave May       if (cell >= 0) {
852af74b616SDave May         cells[p].rank = 0;
853af74b616SDave May         cells[p].index = cell;
854af74b616SDave May         numFound++;
855af74b616SDave May       }
856af74b616SDave May     }
8579cb35068SDave May     if (cells[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
8589cb35068SDave May       terminating_query_type[1]++;
8599cb35068SDave May       continue;
8609cb35068SDave May     }
861af74b616SDave May 
862953fc75cSMatthew G. Knepley     if (hash) {
863af74b616SDave May       PetscBool found_box;
864af74b616SDave May 
865ddce0771SMatthew G. Knepley       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "Checking point %D (%.2g, %.2g, %.2g)\n", p, point[0], point[1], point[2]);CHKERRQ(ierr);}
866af74b616SDave May       /* allow for case that point is outside box - abort early */
867af74b616SDave May       ierr = PetscGridHashGetEnclosingBoxQuery(mesh->lbox, 1, point, dbin, &bin,&found_box);CHKERRQ(ierr);
868af74b616SDave May       if (found_box) {
869ddce0771SMatthew G. Knepley         if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  Found point in box %D (%D, %D, %D)\n", bin, dbin[0], dbin[1], dbin[2]);CHKERRQ(ierr);}
870cafe43deSMatthew G. Knepley         /* TODO Lay an interface over this so we can switch between Section (dense) and Label (sparse) */
871cafe43deSMatthew G. Knepley         ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr);
872cafe43deSMatthew G. Knepley         ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr);
873cafe43deSMatthew G. Knepley         for (c = cellOffset; c < cellOffset + numCells; ++c) {
874ddce0771SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    Checking for point in cell %D\n", boxCells[c]);CHKERRQ(ierr);}
875cafe43deSMatthew G. Knepley           ierr = DMPlexLocatePoint_Internal(dm, dim, point, boxCells[c], &cell);CHKERRQ(ierr);
8763a93e3b7SToby Isaac           if (cell >= 0) {
877ddce0771SMatthew G. Knepley             if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "      FOUND in cell %D\n", cell);CHKERRQ(ierr);}
8783a93e3b7SToby Isaac             cells[p].rank = 0;
8793a93e3b7SToby Isaac             cells[p].index = cell;
8803a93e3b7SToby Isaac             numFound++;
8819cb35068SDave May             terminating_query_type[2]++;
8823a93e3b7SToby Isaac             break;
883ccd2543fSMatthew G Knepley           }
8843a93e3b7SToby Isaac         }
885af74b616SDave May       }
886953fc75cSMatthew G. Knepley     } else {
887953fc75cSMatthew G. Knepley       for (c = cStart; c < cEnd; ++c) {
888953fc75cSMatthew G. Knepley         ierr = DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);CHKERRQ(ierr);
8893a93e3b7SToby Isaac         if (cell >= 0) {
8903a93e3b7SToby Isaac           cells[p].rank = 0;
8913a93e3b7SToby Isaac           cells[p].index = cell;
8923a93e3b7SToby Isaac           numFound++;
8939cb35068SDave May           terminating_query_type[2]++;
8943a93e3b7SToby Isaac           break;
895953fc75cSMatthew G. Knepley         }
896953fc75cSMatthew G. Knepley       }
8973a93e3b7SToby Isaac     }
898ccd2543fSMatthew G Knepley   }
899953fc75cSMatthew G. Knepley   if (hash) {ierr = ISRestoreIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr);}
90062a38674SMatthew G. Knepley   if (ltype == DM_POINTLOCATION_NEAREST && hash && numFound < numPoints) {
90162a38674SMatthew G. Knepley     for (p = 0; p < numPoints; p++) {
90262a38674SMatthew G. Knepley       const PetscScalar *point = &a[p*bs];
90362a38674SMatthew G. Knepley       PetscReal          cpoint[3], diff[3], dist, distMax = PETSC_MAX_REAL;
904e56f9228SJed Brown       PetscInt           dbin[3] = {-1,-1,-1}, bin, cellOffset, d;
90562a38674SMatthew G. Knepley 
906e9b685f5SMatthew G. Knepley       if (cells[p].index < 0) {
90762a38674SMatthew G. Knepley         ++numFound;
90862a38674SMatthew G. Knepley         ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr);
90962a38674SMatthew G. Knepley         ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr);
91062a38674SMatthew G. Knepley         ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr);
91162a38674SMatthew G. Knepley         for (c = cellOffset; c < cellOffset + numCells; ++c) {
91262a38674SMatthew G. Knepley           ierr = DMPlexClosestPoint_Internal(dm, dim, point, boxCells[c], cpoint);CHKERRQ(ierr);
913b716b415SMatthew G. Knepley           for (d = 0; d < dim; ++d) diff[d] = cpoint[d] - PetscRealPart(point[d]);
91462a38674SMatthew G. Knepley           dist = DMPlex_NormD_Internal(dim, diff);
91562a38674SMatthew G. Knepley           if (dist < distMax) {
91662a38674SMatthew G. Knepley             for (d = 0; d < dim; ++d) a[p*bs+d] = cpoint[d];
91762a38674SMatthew G. Knepley             cells[p].rank  = 0;
91862a38674SMatthew G. Knepley             cells[p].index = boxCells[c];
91962a38674SMatthew G. Knepley             distMax = dist;
92062a38674SMatthew G. Knepley             break;
92162a38674SMatthew G. Knepley           }
92262a38674SMatthew G. Knepley         }
92362a38674SMatthew G. Knepley       }
92462a38674SMatthew G. Knepley     }
92562a38674SMatthew G. Knepley   }
92662a38674SMatthew G. Knepley   /* This code is only be relevant when interfaced to parallel point location */
927cafe43deSMatthew G. Knepley   /* Check for highest numbered proc that claims a point (do we care?) */
9282d1fa6caSMatthew G. Knepley   if (ltype == DM_POINTLOCATION_REMOVE && numFound < numPoints) {
9293a93e3b7SToby Isaac     ierr = PetscMalloc1(numFound,&found);CHKERRQ(ierr);
9303a93e3b7SToby Isaac     for (p = 0, numFound = 0; p < numPoints; p++) {
9313a93e3b7SToby Isaac       if (cells[p].rank >= 0 && cells[p].index >= 0) {
9323a93e3b7SToby Isaac         if (numFound < p) {
9333a93e3b7SToby Isaac           cells[numFound] = cells[p];
9343a93e3b7SToby Isaac         }
9353a93e3b7SToby Isaac         found[numFound++] = p;
9363a93e3b7SToby Isaac       }
9373a93e3b7SToby Isaac     }
9383a93e3b7SToby Isaac   }
93962a38674SMatthew G. Knepley   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
940af74b616SDave May   if (!reuse) {
9413a93e3b7SToby Isaac     ierr = PetscSFSetGraph(cellSF, cEnd - cStart, numFound, found, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER);CHKERRQ(ierr);
942af74b616SDave May   }
943af74b616SDave May   ierr = PetscTime(&t1);CHKERRQ(ierr);
9449cb35068SDave May   if (hash) {
9452d4ee042Sprj-     ierr = PetscInfo3(dm,"[DMLocatePoints_Plex] terminating_query_type : %D [outside domain] : %D [inside initial cell] : %D [hash]\n",terminating_query_type[0],terminating_query_type[1],terminating_query_type[2]);CHKERRQ(ierr);
9469cb35068SDave May   } else {
9472d4ee042Sprj-     ierr = PetscInfo3(dm,"[DMLocatePoints_Plex] terminating_query_type : %D [outside domain] : %D [inside initial cell] : %D [brute-force]\n",terminating_query_type[0],terminating_query_type[1],terminating_query_type[2]);CHKERRQ(ierr);
9489cb35068SDave May   }
949af74b616SDave 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);
950cadf77a0SMark Adams   ierr = PetscLogEventEnd(DMPLEX_LocatePoints,0,0,0,0);CHKERRQ(ierr);
951ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
952ccd2543fSMatthew G Knepley }
953ccd2543fSMatthew G Knepley 
954741bfc07SMatthew G. Knepley /*@C
955741bfc07SMatthew G. Knepley   DMPlexComputeProjection2Dto1D - Rewrite coordinates to be the 1D projection of the 2D coordinates
956741bfc07SMatthew G. Knepley 
957741bfc07SMatthew G. Knepley   Not collective
958741bfc07SMatthew G. Knepley 
9596b867d5aSJose E. Roman   Input/Output Parameter:
9606b867d5aSJose E. Roman . coords - The coordinates of a segment, on output the new y-coordinate, and 0 for x
961741bfc07SMatthew G. Knepley 
9626b867d5aSJose E. Roman   Output Parameter:
9636b867d5aSJose E. Roman . R - The rotation which accomplishes the projection
964741bfc07SMatthew G. Knepley 
965741bfc07SMatthew G. Knepley   Level: developer
966741bfc07SMatthew G. Knepley 
967741bfc07SMatthew G. Knepley .seealso: DMPlexComputeProjection3Dto1D(), DMPlexComputeProjection3Dto2D()
968741bfc07SMatthew G. Knepley @*/
969741bfc07SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar coords[], PetscReal R[])
97017fe8556SMatthew G. Knepley {
97117fe8556SMatthew G. Knepley   const PetscReal x = PetscRealPart(coords[2] - coords[0]);
97217fe8556SMatthew G. Knepley   const PetscReal y = PetscRealPart(coords[3] - coords[1]);
9738b49ba18SBarry Smith   const PetscReal r = PetscSqrtReal(x*x + y*y), c = x/r, s = y/r;
97417fe8556SMatthew G. Knepley 
97517fe8556SMatthew G. Knepley   PetscFunctionBegin;
9761c99cf0cSGeoffrey Irving   R[0] = c; R[1] = -s;
9771c99cf0cSGeoffrey Irving   R[2] = s; R[3] =  c;
97817fe8556SMatthew G. Knepley   coords[0] = 0.0;
9797f07f362SMatthew G. Knepley   coords[1] = r;
98017fe8556SMatthew G. Knepley   PetscFunctionReturn(0);
98117fe8556SMatthew G. Knepley }
98217fe8556SMatthew G. Knepley 
983741bfc07SMatthew G. Knepley /*@C
984741bfc07SMatthew G. Knepley   DMPlexComputeProjection3Dto1D - Rewrite coordinates to be the 1D projection of the 3D coordinates
98528dbe442SToby Isaac 
986741bfc07SMatthew G. Knepley   Not collective
98728dbe442SToby Isaac 
9886b867d5aSJose E. Roman   Input/Output Parameter:
9896b867d5aSJose E. Roman . coords - The coordinates of a segment; on output, the new y-coordinate, and 0 for x and z
990741bfc07SMatthew G. Knepley 
9916b867d5aSJose E. Roman   Output Parameter:
9926b867d5aSJose E. Roman . R - The rotation which accomplishes the projection
993741bfc07SMatthew G. Knepley 
994741bfc07SMatthew 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
995741bfc07SMatthew G. Knepley 
996741bfc07SMatthew G. Knepley   Level: developer
997741bfc07SMatthew G. Knepley 
998741bfc07SMatthew G. Knepley .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto2D()
999741bfc07SMatthew G. Knepley @*/
1000741bfc07SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar coords[], PetscReal R[])
100128dbe442SToby Isaac {
100228dbe442SToby Isaac   PetscReal      x    = PetscRealPart(coords[3] - coords[0]);
100328dbe442SToby Isaac   PetscReal      y    = PetscRealPart(coords[4] - coords[1]);
100428dbe442SToby Isaac   PetscReal      z    = PetscRealPart(coords[5] - coords[2]);
100528dbe442SToby Isaac   PetscReal      r    = PetscSqrtReal(x*x + y*y + z*z);
100628dbe442SToby Isaac   PetscReal      rinv = 1. / r;
100728dbe442SToby Isaac   PetscFunctionBegin;
100828dbe442SToby Isaac 
100928dbe442SToby Isaac   x *= rinv; y *= rinv; z *= rinv;
101028dbe442SToby Isaac   if (x > 0.) {
101128dbe442SToby Isaac     PetscReal inv1pX   = 1./ (1. + x);
101228dbe442SToby Isaac 
101328dbe442SToby Isaac     R[0] = x; R[1] = -y;              R[2] = -z;
101428dbe442SToby Isaac     R[3] = y; R[4] = 1. - y*y*inv1pX; R[5] =     -y*z*inv1pX;
101528dbe442SToby Isaac     R[6] = z; R[7] =     -y*z*inv1pX; R[8] = 1. - z*z*inv1pX;
101628dbe442SToby Isaac   }
101728dbe442SToby Isaac   else {
101828dbe442SToby Isaac     PetscReal inv1mX   = 1./ (1. - x);
101928dbe442SToby Isaac 
102028dbe442SToby Isaac     R[0] = x; R[1] = z;               R[2] = y;
102128dbe442SToby Isaac     R[3] = y; R[4] =     -y*z*inv1mX; R[5] = 1. - y*y*inv1mX;
102228dbe442SToby Isaac     R[6] = z; R[7] = 1. - z*z*inv1mX; R[8] =     -y*z*inv1mX;
102328dbe442SToby Isaac   }
102428dbe442SToby Isaac   coords[0] = 0.0;
102528dbe442SToby Isaac   coords[1] = r;
102628dbe442SToby Isaac   PetscFunctionReturn(0);
102728dbe442SToby Isaac }
102828dbe442SToby Isaac 
1029741bfc07SMatthew G. Knepley /*@
1030c871b86eSJed Brown   DMPlexComputeProjection3Dto2D - Rewrite coordinates of 3 or more coplanar 3D points to a common 2D basis for the
1031c871b86eSJed Brown     plane.  The normal is defined by positive orientation of the first 3 points.
1032741bfc07SMatthew G. Knepley 
1033741bfc07SMatthew G. Knepley   Not collective
1034741bfc07SMatthew G. Knepley 
1035741bfc07SMatthew G. Knepley   Input Parameter:
10366b867d5aSJose E. Roman . coordSize - Length of coordinate array (3x number of points); must be at least 9 (3 points)
1037741bfc07SMatthew G. Knepley 
10386b867d5aSJose E. Roman   Input/Output Parameter:
10396b867d5aSJose E. Roman . coords - The interlaced coordinates of each coplanar 3D point; on output the first
10406b867d5aSJose E. Roman            2*coordSize/3 entries contain interlaced 2D points, with the rest undefined
10416b867d5aSJose E. Roman 
10426b867d5aSJose E. Roman   Output Parameter:
10436b867d5aSJose E. Roman . R - 3x3 row-major rotation matrix whose columns are the tangent basis [t1, t2, n].  Multiplying by R^T transforms from original frame to tangent frame.
1044741bfc07SMatthew G. Knepley 
1045741bfc07SMatthew G. Knepley   Level: developer
1046741bfc07SMatthew G. Knepley 
1047741bfc07SMatthew G. Knepley .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto1D()
1048741bfc07SMatthew G. Knepley @*/
1049741bfc07SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt coordSize, PetscScalar coords[], PetscReal R[])
1050ccd2543fSMatthew G Knepley {
1051c871b86eSJed Brown   PetscReal x1[3], x2[3], n[3], c[3], norm;
1052ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
1053c871b86eSJed Brown   PetscInt       d, p;
1054ccd2543fSMatthew G Knepley 
1055ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1056ccd2543fSMatthew G Knepley   /* 0) Calculate normal vector */
1057ccd2543fSMatthew G Knepley   for (d = 0; d < dim; ++d) {
10581ee9d5ecSMatthew G. Knepley     x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]);
10591ee9d5ecSMatthew G. Knepley     x2[d] = PetscRealPart(coords[2*dim+d] - coords[0*dim+d]);
1060ccd2543fSMatthew G Knepley   }
1061c871b86eSJed Brown   // n = x1 \otimes x2
1062ccd2543fSMatthew G Knepley   n[0] = x1[1]*x2[2] - x1[2]*x2[1];
1063ccd2543fSMatthew G Knepley   n[1] = x1[2]*x2[0] - x1[0]*x2[2];
1064ccd2543fSMatthew G Knepley   n[2] = x1[0]*x2[1] - x1[1]*x2[0];
10658b49ba18SBarry Smith   norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]);
1066c871b86eSJed Brown   for (d = 0; d < dim; d++) n[d] /= norm;
1067c871b86eSJed Brown   norm = PetscSqrtReal(x1[0] * x1[0] + x1[1] * x1[1] + x1[2] * x1[2]);
1068c871b86eSJed Brown   for (d = 0; d < dim; d++) x1[d] /= norm;
1069c871b86eSJed Brown   // x2 = n \otimes x1
1070c871b86eSJed Brown   x2[0] = n[1] * x1[2] - n[2] * x1[1];
1071c871b86eSJed Brown   x2[1] = n[2] * x1[0] - n[0] * x1[2];
1072c871b86eSJed Brown   x2[2] = n[0] * x1[1] - n[1] * x1[0];
1073c871b86eSJed Brown   for (d=0; d<dim; d++) {
1074c871b86eSJed Brown     R[d * dim + 0] = x1[d];
1075c871b86eSJed Brown     R[d * dim + 1] = x2[d];
1076c871b86eSJed Brown     R[d * dim + 2] = n[d];
1077c871b86eSJed Brown     c[d] = PetscRealPart(coords[0*dim + d]);
107873868372SMatthew G. Knepley   }
1079c871b86eSJed Brown   for (p=0; p<coordSize/dim; p++) {
1080c871b86eSJed Brown     PetscReal y[3];
1081c871b86eSJed Brown     for (d=0; d<dim; d++) y[d] = PetscRealPart(coords[p*dim + d]) - c[d];
1082c871b86eSJed Brown     for (d=0; d<2; d++) coords[p*2+d] = R[0*dim + d] * y[0] + R[1*dim + d] * y[1] + R[2*dim + d] * y[2];
10837f07f362SMatthew G. Knepley   }
1084ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1085ccd2543fSMatthew G Knepley }
1086ccd2543fSMatthew G Knepley 
10876322fe33SJed Brown PETSC_UNUSED
1088834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[])
1089834e62ceSMatthew G. Knepley {
1090834e62ceSMatthew G. Knepley   /* Signed volume is 1/2 the determinant
1091834e62ceSMatthew G. Knepley 
1092834e62ceSMatthew G. Knepley    |  1  1  1 |
1093834e62ceSMatthew G. Knepley    | x0 x1 x2 |
1094834e62ceSMatthew G. Knepley    | y0 y1 y2 |
1095834e62ceSMatthew G. Knepley 
1096834e62ceSMatthew G. Knepley      but if x0,y0 is the origin, we have
1097834e62ceSMatthew G. Knepley 
1098834e62ceSMatthew G. Knepley    | x1 x2 |
1099834e62ceSMatthew G. Knepley    | y1 y2 |
1100834e62ceSMatthew G. Knepley   */
1101834e62ceSMatthew G. Knepley   const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1];
1102834e62ceSMatthew G. Knepley   const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1];
1103834e62ceSMatthew G. Knepley   PetscReal       M[4], detM;
1104834e62ceSMatthew G. Knepley   M[0] = x1; M[1] = x2;
110586623015SMatthew G. Knepley   M[2] = y1; M[3] = y2;
1106923591dfSMatthew G. Knepley   DMPlex_Det2D_Internal(&detM, M);
1107834e62ceSMatthew G. Knepley   *vol = 0.5*detM;
11083bc0b13bSBarry Smith   (void)PetscLogFlops(5.0);
1109834e62ceSMatthew G. Knepley }
1110834e62ceSMatthew G. Knepley 
1111834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[])
1112834e62ceSMatthew G. Knepley {
1113923591dfSMatthew G. Knepley   DMPlex_Det2D_Internal(vol, coords);
1114834e62ceSMatthew G. Knepley   *vol *= 0.5;
1115834e62ceSMatthew G. Knepley }
1116834e62ceSMatthew G. Knepley 
11176322fe33SJed Brown PETSC_UNUSED
1118834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[])
1119834e62ceSMatthew G. Knepley {
1120834e62ceSMatthew G. Knepley   /* Signed volume is 1/6th of the determinant
1121834e62ceSMatthew G. Knepley 
1122834e62ceSMatthew G. Knepley    |  1  1  1  1 |
1123834e62ceSMatthew G. Knepley    | x0 x1 x2 x3 |
1124834e62ceSMatthew G. Knepley    | y0 y1 y2 y3 |
1125834e62ceSMatthew G. Knepley    | z0 z1 z2 z3 |
1126834e62ceSMatthew G. Knepley 
1127834e62ceSMatthew G. Knepley      but if x0,y0,z0 is the origin, we have
1128834e62ceSMatthew G. Knepley 
1129834e62ceSMatthew G. Knepley    | x1 x2 x3 |
1130834e62ceSMatthew G. Knepley    | y1 y2 y3 |
1131834e62ceSMatthew G. Knepley    | z1 z2 z3 |
1132834e62ceSMatthew G. Knepley   */
1133834e62ceSMatthew G. Knepley   const PetscReal x1 = coords[3] - coords[0], y1 = coords[4]  - coords[1], z1 = coords[5]  - coords[2];
1134834e62ceSMatthew G. Knepley   const PetscReal x2 = coords[6] - coords[0], y2 = coords[7]  - coords[1], z2 = coords[8]  - coords[2];
1135834e62ceSMatthew G. Knepley   const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2];
11360a3da2c2SToby Isaac   const PetscReal onesixth = ((PetscReal)1./(PetscReal)6.);
1137834e62ceSMatthew G. Knepley   PetscReal       M[9], detM;
1138834e62ceSMatthew G. Knepley   M[0] = x1; M[1] = x2; M[2] = x3;
1139834e62ceSMatthew G. Knepley   M[3] = y1; M[4] = y2; M[5] = y3;
1140834e62ceSMatthew G. Knepley   M[6] = z1; M[7] = z2; M[8] = z3;
1141923591dfSMatthew G. Knepley   DMPlex_Det3D_Internal(&detM, M);
11420a3da2c2SToby Isaac   *vol = -onesixth*detM;
11433bc0b13bSBarry Smith   (void)PetscLogFlops(10.0);
1144834e62ceSMatthew G. Knepley }
1145834e62ceSMatthew G. Knepley 
11460ec8681fSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[])
11470ec8681fSMatthew G. Knepley {
11480a3da2c2SToby Isaac   const PetscReal onesixth = ((PetscReal)1./(PetscReal)6.);
1149923591dfSMatthew G. Knepley   DMPlex_Det3D_Internal(vol, coords);
11500a3da2c2SToby Isaac   *vol *= -onesixth;
11510ec8681fSMatthew G. Knepley }
11520ec8681fSMatthew G. Knepley 
1153cb92db44SToby Isaac static PetscErrorCode DMPlexComputePointGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1154cb92db44SToby Isaac {
1155cb92db44SToby Isaac   PetscSection   coordSection;
1156cb92db44SToby Isaac   Vec            coordinates;
1157cb92db44SToby Isaac   const PetscScalar *coords;
1158cb92db44SToby Isaac   PetscInt       dim, d, off;
1159cb92db44SToby Isaac   PetscErrorCode ierr;
1160cb92db44SToby Isaac 
1161cb92db44SToby Isaac   PetscFunctionBegin;
1162cb92db44SToby Isaac   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1163cb92db44SToby Isaac   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1164cb92db44SToby Isaac   ierr = PetscSectionGetDof(coordSection,e,&dim);CHKERRQ(ierr);
1165cb92db44SToby Isaac   if (!dim) PetscFunctionReturn(0);
1166cb92db44SToby Isaac   ierr = PetscSectionGetOffset(coordSection,e,&off);CHKERRQ(ierr);
1167cb92db44SToby Isaac   ierr = VecGetArrayRead(coordinates,&coords);CHKERRQ(ierr);
1168cb92db44SToby Isaac   if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[off + d]);}
1169cb92db44SToby Isaac   ierr = VecRestoreArrayRead(coordinates,&coords);CHKERRQ(ierr);
1170cb92db44SToby Isaac   *detJ = 1.;
1171cb92db44SToby Isaac   if (J) {
1172cb92db44SToby Isaac     for (d = 0; d < dim * dim; d++) J[d] = 0.;
1173cb92db44SToby Isaac     for (d = 0; d < dim; d++) J[d * dim + d] = 1.;
1174cb92db44SToby Isaac     if (invJ) {
1175cb92db44SToby Isaac       for (d = 0; d < dim * dim; d++) invJ[d] = 0.;
1176cb92db44SToby Isaac       for (d = 0; d < dim; d++) invJ[d * dim + d] = 1.;
1177cb92db44SToby Isaac     }
1178cb92db44SToby Isaac   }
1179cb92db44SToby Isaac   PetscFunctionReturn(0);
1180cb92db44SToby Isaac }
1181cb92db44SToby Isaac 
118217fe8556SMatthew G. Knepley static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
118317fe8556SMatthew G. Knepley {
118417fe8556SMatthew G. Knepley   PetscSection   coordSection;
118517fe8556SMatthew G. Knepley   Vec            coordinates;
1186a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
11878bf5c034SToby Isaac   PetscInt       numCoords, d, pStart, pEnd, numSelfCoords = 0;
118817fe8556SMatthew G. Knepley   PetscErrorCode ierr;
118917fe8556SMatthew G. Knepley 
119017fe8556SMatthew G. Knepley   PetscFunctionBegin;
119117fe8556SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
119269d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
11938bf5c034SToby Isaac   ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr);
11948bf5c034SToby Isaac   if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);}
119517fe8556SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
11968bf5c034SToby Isaac   numCoords = numSelfCoords ? numSelfCoords : numCoords;
1197adac9986SMatthew G. Knepley   if (invJ && !J) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL");
11987f07f362SMatthew G. Knepley   *detJ = 0.0;
119928dbe442SToby Isaac   if (numCoords == 6) {
120028dbe442SToby Isaac     const PetscInt dim = 3;
120128dbe442SToby Isaac     PetscReal      R[9], J0;
120228dbe442SToby Isaac 
120328dbe442SToby Isaac     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1204741bfc07SMatthew G. Knepley     ierr = DMPlexComputeProjection3Dto1D(coords, R);CHKERRQ(ierr);
120528dbe442SToby Isaac     if (J)    {
120628dbe442SToby Isaac       J0   = 0.5*PetscRealPart(coords[1]);
120728dbe442SToby Isaac       J[0] = R[0]*J0; J[1] = R[1]; J[2] = R[2];
120828dbe442SToby Isaac       J[3] = R[3]*J0; J[4] = R[4]; J[5] = R[5];
120928dbe442SToby Isaac       J[6] = R[6]*J0; J[7] = R[7]; J[8] = R[8];
121028dbe442SToby Isaac       DMPlex_Det3D_Internal(detJ, J);
121128dbe442SToby Isaac       if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1212adac9986SMatthew G. Knepley     }
121328dbe442SToby Isaac   } else if (numCoords == 4) {
12147f07f362SMatthew G. Knepley     const PetscInt dim = 2;
12157f07f362SMatthew G. Knepley     PetscReal      R[4], J0;
12167f07f362SMatthew G. Knepley 
12177f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1218741bfc07SMatthew G. Knepley     ierr = DMPlexComputeProjection2Dto1D(coords, R);CHKERRQ(ierr);
121917fe8556SMatthew G. Knepley     if (J)    {
12207f07f362SMatthew G. Knepley       J0   = 0.5*PetscRealPart(coords[1]);
12217f07f362SMatthew G. Knepley       J[0] = R[0]*J0; J[1] = R[1];
12227f07f362SMatthew G. Knepley       J[2] = R[2]*J0; J[3] = R[3];
1223923591dfSMatthew G. Knepley       DMPlex_Det2D_Internal(detJ, J);
1224923591dfSMatthew G. Knepley       if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1225adac9986SMatthew G. Knepley     }
12267f07f362SMatthew G. Knepley   } else if (numCoords == 2) {
12277f07f362SMatthew G. Knepley     const PetscInt dim = 1;
12287f07f362SMatthew G. Knepley 
12297f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
12307f07f362SMatthew G. Knepley     if (J)    {
12317f07f362SMatthew G. Knepley       J[0]  = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0]));
123217fe8556SMatthew G. Knepley       *detJ = J[0];
12333bc0b13bSBarry Smith       ierr = PetscLogFlops(2.0);CHKERRQ(ierr);
12343bc0b13bSBarry Smith       if (invJ) {invJ[0] = 1.0/J[0]; ierr = PetscLogFlops(1.0);CHKERRQ(ierr);}
1235adac9986SMatthew G. Knepley     }
1236796f034aSJed Brown   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %D != 2", numCoords);
123717fe8556SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
123817fe8556SMatthew G. Knepley   PetscFunctionReturn(0);
123917fe8556SMatthew G. Knepley }
124017fe8556SMatthew G. Knepley 
1241ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1242ccd2543fSMatthew G Knepley {
1243ccd2543fSMatthew G Knepley   PetscSection   coordSection;
1244ccd2543fSMatthew G Knepley   Vec            coordinates;
1245a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
124603d7ed2eSStefano Zampini   PetscInt       numCoords, numSelfCoords = 0, d, f, g, pStart, pEnd;
1247ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
1248ccd2543fSMatthew G Knepley 
1249ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1250ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
125169d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
125203d7ed2eSStefano Zampini   ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr);
125303d7ed2eSStefano Zampini   if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);}
1254ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
125503d7ed2eSStefano Zampini   numCoords = numSelfCoords ? numSelfCoords : numCoords;
12567f07f362SMatthew G. Knepley   *detJ = 0.0;
1257ccd2543fSMatthew G Knepley   if (numCoords == 9) {
12587f07f362SMatthew G. Knepley     const PetscInt dim = 3;
12597f07f362SMatthew 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};
12607f07f362SMatthew G. Knepley 
12617f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1262741bfc07SMatthew G. Knepley     ierr = DMPlexComputeProjection3Dto2D(numCoords, coords, R);CHKERRQ(ierr);
12637f07f362SMatthew G. Knepley     if (J)    {
1264b7ad821dSMatthew G. Knepley       const PetscInt pdim = 2;
1265b7ad821dSMatthew G. Knepley 
1266b7ad821dSMatthew G. Knepley       for (d = 0; d < pdim; d++) {
1267b7ad821dSMatthew G. Knepley         for (f = 0; f < pdim; f++) {
1268b7ad821dSMatthew G. Knepley           J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
1269ccd2543fSMatthew G Knepley         }
12707f07f362SMatthew G. Knepley       }
12713bc0b13bSBarry Smith       ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1272923591dfSMatthew G. Knepley       DMPlex_Det3D_Internal(detJ, J0);
12737f07f362SMatthew G. Knepley       for (d = 0; d < dim; d++) {
12747f07f362SMatthew G. Knepley         for (f = 0; f < dim; f++) {
12757f07f362SMatthew G. Knepley           J[d*dim+f] = 0.0;
12767f07f362SMatthew G. Knepley           for (g = 0; g < dim; g++) {
12777f07f362SMatthew G. Knepley             J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
12787f07f362SMatthew G. Knepley           }
12797f07f362SMatthew G. Knepley         }
12807f07f362SMatthew G. Knepley       }
12813bc0b13bSBarry Smith       ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
12827f07f362SMatthew G. Knepley     }
1283923591dfSMatthew G. Knepley     if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
12847f07f362SMatthew G. Knepley   } else if (numCoords == 6) {
12857f07f362SMatthew G. Knepley     const PetscInt dim = 2;
12867f07f362SMatthew G. Knepley 
12877f07f362SMatthew G. Knepley     if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1288ccd2543fSMatthew G Knepley     if (J)    {
1289ccd2543fSMatthew G Knepley       for (d = 0; d < dim; d++) {
1290ccd2543fSMatthew G Knepley         for (f = 0; f < dim; f++) {
1291ccd2543fSMatthew G Knepley           J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
1292ccd2543fSMatthew G Knepley         }
1293ccd2543fSMatthew G Knepley       }
12943bc0b13bSBarry Smith       ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1295923591dfSMatthew G. Knepley       DMPlex_Det2D_Internal(detJ, J);
1296ccd2543fSMatthew G Knepley     }
1297923591dfSMatthew G. Knepley     if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1298796f034aSJed Brown   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %D != 6 or 9", numCoords);
1299ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
1300ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1301ccd2543fSMatthew G Knepley }
1302ccd2543fSMatthew G Knepley 
1303412e9a14SMatthew G. Knepley static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscBool isTensor, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1304ccd2543fSMatthew G Knepley {
1305ccd2543fSMatthew G Knepley   PetscSection   coordSection;
1306ccd2543fSMatthew G Knepley   Vec            coordinates;
1307a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
13080d29256aSToby Isaac   PetscInt       numCoords, numSelfCoords = 0, d, f, g, pStart, pEnd;
1309ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
1310ccd2543fSMatthew G Knepley 
1311ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1312ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
131369d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
13140d29256aSToby Isaac   ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr);
13150d29256aSToby Isaac   if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);}
131699dec3a6SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
131771f58de1SToby Isaac   numCoords = numSelfCoords ? numSelfCoords : numCoords;
1318dfccc68fSToby Isaac   if (!Nq) {
1319412e9a14SMatthew G. Knepley     PetscInt vorder[4] = {0, 1, 2, 3};
1320412e9a14SMatthew G. Knepley 
1321412e9a14SMatthew G. Knepley     if (isTensor) {vorder[2] = 3; vorder[3] = 2;}
13227f07f362SMatthew G. Knepley     *detJ = 0.0;
132399dec3a6SMatthew G. Knepley     if (numCoords == 12) {
132499dec3a6SMatthew G. Knepley       const PetscInt dim = 3;
132599dec3a6SMatthew 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};
132699dec3a6SMatthew G. Knepley 
1327dfccc68fSToby Isaac       if (v)   {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);}
1328741bfc07SMatthew G. Knepley       ierr = DMPlexComputeProjection3Dto2D(numCoords, coords, R);CHKERRQ(ierr);
132999dec3a6SMatthew G. Knepley       if (J)    {
133099dec3a6SMatthew G. Knepley         const PetscInt pdim = 2;
133199dec3a6SMatthew G. Knepley 
133299dec3a6SMatthew G. Knepley         for (d = 0; d < pdim; d++) {
1333412e9a14SMatthew G. Knepley           J0[d*dim+0] = 0.5*(PetscRealPart(coords[vorder[1]*pdim+d]) - PetscRealPart(coords[vorder[0]*pdim+d]));
1334412e9a14SMatthew G. Knepley           J0[d*dim+1] = 0.5*(PetscRealPart(coords[vorder[2]*pdim+d]) - PetscRealPart(coords[vorder[1]*pdim+d]));
133599dec3a6SMatthew G. Knepley         }
13363bc0b13bSBarry Smith         ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1337923591dfSMatthew G. Knepley         DMPlex_Det3D_Internal(detJ, J0);
133899dec3a6SMatthew G. Knepley         for (d = 0; d < dim; d++) {
133999dec3a6SMatthew G. Knepley           for (f = 0; f < dim; f++) {
134099dec3a6SMatthew G. Knepley             J[d*dim+f] = 0.0;
134199dec3a6SMatthew G. Knepley             for (g = 0; g < dim; g++) {
134299dec3a6SMatthew G. Knepley               J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
134399dec3a6SMatthew G. Knepley             }
134499dec3a6SMatthew G. Knepley           }
134599dec3a6SMatthew G. Knepley         }
13463bc0b13bSBarry Smith         ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
134799dec3a6SMatthew G. Knepley       }
1348923591dfSMatthew G. Knepley       if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
134971f58de1SToby Isaac     } else if (numCoords == 8) {
135099dec3a6SMatthew G. Knepley       const PetscInt dim = 2;
135199dec3a6SMatthew G. Knepley 
1352dfccc68fSToby Isaac       if (v)   {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);}
1353ccd2543fSMatthew G Knepley       if (J)    {
1354ccd2543fSMatthew G Knepley         for (d = 0; d < dim; d++) {
1355412e9a14SMatthew G. Knepley           J[d*dim+0] = 0.5*(PetscRealPart(coords[vorder[1]*dim+d]) - PetscRealPart(coords[vorder[0]*dim+d]));
1356412e9a14SMatthew G. Knepley           J[d*dim+1] = 0.5*(PetscRealPart(coords[vorder[3]*dim+d]) - PetscRealPart(coords[vorder[0]*dim+d]));
1357ccd2543fSMatthew G Knepley         }
13583bc0b13bSBarry Smith         ierr = PetscLogFlops(8.0);CHKERRQ(ierr);
1359923591dfSMatthew G. Knepley         DMPlex_Det2D_Internal(detJ, J);
1360ccd2543fSMatthew G Knepley       }
1361923591dfSMatthew G. Knepley       if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1362796f034aSJed Brown     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords);
1363dfccc68fSToby Isaac   } else {
1364dfccc68fSToby Isaac     const PetscInt Nv = 4;
1365dfccc68fSToby Isaac     const PetscInt dimR = 2;
1366412e9a14SMatthew G. Knepley     PetscInt  zToPlex[4] = {0, 1, 3, 2};
1367dfccc68fSToby Isaac     PetscReal zOrder[12];
1368dfccc68fSToby Isaac     PetscReal zCoeff[12];
1369dfccc68fSToby Isaac     PetscInt  i, j, k, l, dim;
1370dfccc68fSToby Isaac 
1371412e9a14SMatthew G. Knepley     if (isTensor) {zToPlex[2] = 2; zToPlex[3] = 3;}
1372dfccc68fSToby Isaac     if (numCoords == 12) {
1373dfccc68fSToby Isaac       dim = 3;
1374dfccc68fSToby Isaac     } else if (numCoords == 8) {
1375dfccc68fSToby Isaac       dim = 2;
1376dfccc68fSToby Isaac     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords);
1377dfccc68fSToby Isaac     for (i = 0; i < Nv; i++) {
1378dfccc68fSToby Isaac       PetscInt zi = zToPlex[i];
1379dfccc68fSToby Isaac 
1380dfccc68fSToby Isaac       for (j = 0; j < dim; j++) {
1381dfccc68fSToby Isaac         zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]);
1382dfccc68fSToby Isaac       }
1383dfccc68fSToby Isaac     }
1384dfccc68fSToby Isaac     for (j = 0; j < dim; j++) {
1385dfccc68fSToby Isaac       zCoeff[dim * 0 + j] = 0.25 * (  zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1386dfccc68fSToby Isaac       zCoeff[dim * 1 + j] = 0.25 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1387dfccc68fSToby Isaac       zCoeff[dim * 2 + j] = 0.25 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1388dfccc68fSToby Isaac       zCoeff[dim * 3 + j] = 0.25 * (  zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1389dfccc68fSToby Isaac     }
1390dfccc68fSToby Isaac     for (i = 0; i < Nq; i++) {
1391dfccc68fSToby Isaac       PetscReal xi = points[dimR * i], eta = points[dimR * i + 1];
1392dfccc68fSToby Isaac 
1393dfccc68fSToby Isaac       if (v) {
1394dfccc68fSToby Isaac         PetscReal extPoint[4];
1395dfccc68fSToby Isaac 
1396dfccc68fSToby Isaac         extPoint[0] = 1.;
1397dfccc68fSToby Isaac         extPoint[1] = xi;
1398dfccc68fSToby Isaac         extPoint[2] = eta;
1399dfccc68fSToby Isaac         extPoint[3] = xi * eta;
1400dfccc68fSToby Isaac         for (j = 0; j < dim; j++) {
1401dfccc68fSToby Isaac           PetscReal val = 0.;
1402dfccc68fSToby Isaac 
1403dfccc68fSToby Isaac           for (k = 0; k < Nv; k++) {
1404dfccc68fSToby Isaac             val += extPoint[k] * zCoeff[dim * k + j];
1405dfccc68fSToby Isaac           }
1406dfccc68fSToby Isaac           v[i * dim + j] = val;
1407dfccc68fSToby Isaac         }
1408dfccc68fSToby Isaac       }
1409dfccc68fSToby Isaac       if (J) {
1410dfccc68fSToby Isaac         PetscReal extJ[8];
1411dfccc68fSToby Isaac 
1412dfccc68fSToby Isaac         extJ[0] = 0.;
1413dfccc68fSToby Isaac         extJ[1] = 0.;
1414dfccc68fSToby Isaac         extJ[2] = 1.;
1415dfccc68fSToby Isaac         extJ[3] = 0.;
1416dfccc68fSToby Isaac         extJ[4] = 0.;
1417dfccc68fSToby Isaac         extJ[5] = 1.;
1418dfccc68fSToby Isaac         extJ[6] = eta;
1419dfccc68fSToby Isaac         extJ[7] = xi;
1420dfccc68fSToby Isaac         for (j = 0; j < dim; j++) {
1421dfccc68fSToby Isaac           for (k = 0; k < dimR; k++) {
1422dfccc68fSToby Isaac             PetscReal val = 0.;
1423dfccc68fSToby Isaac 
1424dfccc68fSToby Isaac             for (l = 0; l < Nv; l++) {
1425dfccc68fSToby Isaac               val += zCoeff[dim * l + j] * extJ[dimR * l + k];
1426dfccc68fSToby Isaac             }
1427dfccc68fSToby Isaac             J[i * dim * dim + dim * j + k] = val;
1428dfccc68fSToby Isaac           }
1429dfccc68fSToby Isaac         }
1430dfccc68fSToby Isaac         if (dim == 3) { /* put the cross product in the third component of the Jacobian */
1431dfccc68fSToby Isaac           PetscReal x, y, z;
1432dfccc68fSToby Isaac           PetscReal *iJ = &J[i * dim * dim];
1433dfccc68fSToby Isaac           PetscReal norm;
1434dfccc68fSToby Isaac 
1435dfccc68fSToby Isaac           x = iJ[1 * dim + 0] * iJ[2 * dim + 1] - iJ[1 * dim + 1] * iJ[2 * dim + 0];
1436dfccc68fSToby Isaac           y = iJ[0 * dim + 1] * iJ[2 * dim + 0] - iJ[0 * dim + 0] * iJ[2 * dim + 1];
1437dfccc68fSToby Isaac           z = iJ[0 * dim + 0] * iJ[1 * dim + 1] - iJ[0 * dim + 1] * iJ[1 * dim + 0];
1438dfccc68fSToby Isaac           norm = PetscSqrtReal(x * x + y * y + z * z);
1439dfccc68fSToby Isaac           iJ[2] = x / norm;
1440dfccc68fSToby Isaac           iJ[5] = y / norm;
1441dfccc68fSToby Isaac           iJ[8] = z / norm;
1442dfccc68fSToby Isaac           DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]);
1443dfccc68fSToby Isaac           if (invJ) {DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);}
1444dfccc68fSToby Isaac         } else {
1445dfccc68fSToby Isaac           DMPlex_Det2D_Internal(&detJ[i], &J[i * dim * dim]);
1446dfccc68fSToby Isaac           if (invJ) {DMPlex_Invert2D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);}
1447dfccc68fSToby Isaac         }
1448dfccc68fSToby Isaac       }
1449dfccc68fSToby Isaac     }
1450dfccc68fSToby Isaac   }
145199dec3a6SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr);
1452ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1453ccd2543fSMatthew G Knepley }
1454ccd2543fSMatthew G Knepley 
1455ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1456ccd2543fSMatthew G Knepley {
1457ccd2543fSMatthew G Knepley   PetscSection   coordSection;
1458ccd2543fSMatthew G Knepley   Vec            coordinates;
1459a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
1460ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
146199dec3a6SMatthew G. Knepley   PetscInt       d;
1462ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
1463ccd2543fSMatthew G Knepley 
1464ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1465ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
146669d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1467ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
14687f07f362SMatthew G. Knepley   *detJ = 0.0;
14697f07f362SMatthew G. Knepley   if (v0)   {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1470ccd2543fSMatthew G Knepley   if (J)    {
1471ccd2543fSMatthew G Knepley     for (d = 0; d < dim; d++) {
1472f0df753eSMatthew G. Knepley       /* I orient with outward face normals */
1473f0df753eSMatthew G. Knepley       J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d]));
1474f0df753eSMatthew G. Knepley       J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
1475f0df753eSMatthew G. Knepley       J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1476ccd2543fSMatthew G Knepley     }
14773bc0b13bSBarry Smith     ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
1478923591dfSMatthew G. Knepley     DMPlex_Det3D_Internal(detJ, J);
1479ccd2543fSMatthew G Knepley   }
1480923591dfSMatthew G. Knepley   if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1481ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1482ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1483ccd2543fSMatthew G Knepley }
1484ccd2543fSMatthew G Knepley 
1485dfccc68fSToby Isaac static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1486ccd2543fSMatthew G Knepley {
1487ccd2543fSMatthew G Knepley   PetscSection   coordSection;
1488ccd2543fSMatthew G Knepley   Vec            coordinates;
1489a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
1490ccd2543fSMatthew G Knepley   const PetscInt dim = 3;
1491ccd2543fSMatthew G Knepley   PetscInt       d;
1492ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
1493ccd2543fSMatthew G Knepley 
1494ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1495ccd2543fSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
149669d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1497ccd2543fSMatthew G Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1498dfccc68fSToby Isaac   if (!Nq) {
14997f07f362SMatthew G. Knepley     *detJ = 0.0;
1500dfccc68fSToby Isaac     if (v)   {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);}
1501ccd2543fSMatthew G Knepley     if (J)    {
1502ccd2543fSMatthew G Knepley       for (d = 0; d < dim; d++) {
1503f0df753eSMatthew G. Knepley         J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1504f0df753eSMatthew G. Knepley         J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
1505f0df753eSMatthew G. Knepley         J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d]));
1506ccd2543fSMatthew G Knepley       }
15073bc0b13bSBarry Smith       ierr = PetscLogFlops(18.0);CHKERRQ(ierr);
1508923591dfSMatthew G. Knepley       DMPlex_Det3D_Internal(detJ, J);
1509ccd2543fSMatthew G Knepley     }
1510923591dfSMatthew G. Knepley     if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1511dfccc68fSToby Isaac   } else {
1512dfccc68fSToby Isaac     const PetscInt Nv = 8;
1513dfccc68fSToby Isaac     const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6};
1514dfccc68fSToby Isaac     const PetscInt dim = 3;
1515dfccc68fSToby Isaac     const PetscInt dimR = 3;
1516dfccc68fSToby Isaac     PetscReal zOrder[24];
1517dfccc68fSToby Isaac     PetscReal zCoeff[24];
1518dfccc68fSToby Isaac     PetscInt  i, j, k, l;
1519dfccc68fSToby Isaac 
1520dfccc68fSToby Isaac     for (i = 0; i < Nv; i++) {
1521dfccc68fSToby Isaac       PetscInt zi = zToPlex[i];
1522dfccc68fSToby Isaac 
1523dfccc68fSToby Isaac       for (j = 0; j < dim; j++) {
1524dfccc68fSToby Isaac         zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]);
1525dfccc68fSToby Isaac       }
1526dfccc68fSToby Isaac     }
1527dfccc68fSToby Isaac     for (j = 0; j < dim; j++) {
1528dfccc68fSToby 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]);
1529dfccc68fSToby 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]);
1530dfccc68fSToby 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]);
1531dfccc68fSToby 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]);
1532dfccc68fSToby 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]);
1533dfccc68fSToby 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]);
1534dfccc68fSToby 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]);
1535dfccc68fSToby 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]);
1536dfccc68fSToby Isaac     }
1537dfccc68fSToby Isaac     for (i = 0; i < Nq; i++) {
1538dfccc68fSToby Isaac       PetscReal xi = points[dimR * i], eta = points[dimR * i + 1], theta = points[dimR * i + 2];
1539dfccc68fSToby Isaac 
1540dfccc68fSToby Isaac       if (v) {
154191d2b7ceSToby Isaac         PetscReal extPoint[8];
1542dfccc68fSToby Isaac 
1543dfccc68fSToby Isaac         extPoint[0] = 1.;
1544dfccc68fSToby Isaac         extPoint[1] = xi;
1545dfccc68fSToby Isaac         extPoint[2] = eta;
1546dfccc68fSToby Isaac         extPoint[3] = xi * eta;
1547dfccc68fSToby Isaac         extPoint[4] = theta;
1548dfccc68fSToby Isaac         extPoint[5] = theta * xi;
1549dfccc68fSToby Isaac         extPoint[6] = theta * eta;
1550dfccc68fSToby Isaac         extPoint[7] = theta * eta * xi;
1551dfccc68fSToby Isaac         for (j = 0; j < dim; j++) {
1552dfccc68fSToby Isaac           PetscReal val = 0.;
1553dfccc68fSToby Isaac 
1554dfccc68fSToby Isaac           for (k = 0; k < Nv; k++) {
1555dfccc68fSToby Isaac             val += extPoint[k] * zCoeff[dim * k + j];
1556dfccc68fSToby Isaac           }
1557dfccc68fSToby Isaac           v[i * dim + j] = val;
1558dfccc68fSToby Isaac         }
1559dfccc68fSToby Isaac       }
1560dfccc68fSToby Isaac       if (J) {
1561dfccc68fSToby Isaac         PetscReal extJ[24];
1562dfccc68fSToby Isaac 
1563dfccc68fSToby Isaac         extJ[0]  = 0.         ; extJ[1]  = 0.        ; extJ[2]  = 0.      ;
1564dfccc68fSToby Isaac         extJ[3]  = 1.         ; extJ[4]  = 0.        ; extJ[5]  = 0.      ;
1565dfccc68fSToby Isaac         extJ[6]  = 0.         ; extJ[7]  = 1.        ; extJ[8]  = 0.      ;
1566dfccc68fSToby Isaac         extJ[9]  = eta        ; extJ[10] = xi        ; extJ[11] = 0.      ;
1567dfccc68fSToby Isaac         extJ[12] = 0.         ; extJ[13] = 0.        ; extJ[14] = 1.      ;
1568dfccc68fSToby Isaac         extJ[15] = theta      ; extJ[16] = 0.        ; extJ[17] = xi      ;
1569dfccc68fSToby Isaac         extJ[18] = 0.         ; extJ[19] = theta     ; extJ[20] = eta     ;
1570dfccc68fSToby Isaac         extJ[21] = theta * eta; extJ[22] = theta * xi; extJ[23] = eta * xi;
1571dfccc68fSToby Isaac 
1572dfccc68fSToby Isaac         for (j = 0; j < dim; j++) {
1573dfccc68fSToby Isaac           for (k = 0; k < dimR; k++) {
1574dfccc68fSToby Isaac             PetscReal val = 0.;
1575dfccc68fSToby Isaac 
1576dfccc68fSToby Isaac             for (l = 0; l < Nv; l++) {
1577dfccc68fSToby Isaac               val += zCoeff[dim * l + j] * extJ[dimR * l + k];
1578dfccc68fSToby Isaac             }
1579dfccc68fSToby Isaac             J[i * dim * dim + dim * j + k] = val;
1580dfccc68fSToby Isaac           }
1581dfccc68fSToby Isaac         }
1582dfccc68fSToby Isaac         DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]);
1583dfccc68fSToby Isaac         if (invJ) {DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);}
1584dfccc68fSToby Isaac       }
1585dfccc68fSToby Isaac     }
1586dfccc68fSToby Isaac   }
1587ccd2543fSMatthew G Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr);
1588ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1589ccd2543fSMatthew G Knepley }
1590ccd2543fSMatthew G Knepley 
1591dfccc68fSToby Isaac static PetscErrorCode DMPlexComputeCellGeometryFEM_Implicit(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
1592dfccc68fSToby Isaac {
1593ba2698f1SMatthew G. Knepley   DMPolytopeType  ct;
1594dfccc68fSToby Isaac   PetscInt        depth, dim, coordDim, coneSize, i;
1595dfccc68fSToby Isaac   PetscInt        Nq = 0;
1596dfccc68fSToby Isaac   const PetscReal *points = NULL;
1597dfccc68fSToby Isaac   DMLabel         depthLabel;
1598c330f8ffSToby Isaac   PetscReal       xi0[3] = {-1.,-1.,-1.}, v0[3], J0[9], detJ0;
1599dfccc68fSToby Isaac   PetscBool       isAffine = PETSC_TRUE;
1600dfccc68fSToby Isaac   PetscErrorCode  ierr;
1601dfccc68fSToby Isaac 
1602dfccc68fSToby Isaac   PetscFunctionBegin;
1603dfccc68fSToby Isaac   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1604dfccc68fSToby Isaac   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
1605dfccc68fSToby Isaac   ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr);
1606dfccc68fSToby Isaac   ierr = DMLabelGetValue(depthLabel, cell, &dim);CHKERRQ(ierr);
1607dfccc68fSToby Isaac   if (depth == 1 && dim == 1) {
1608dfccc68fSToby Isaac     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1609dfccc68fSToby Isaac   }
1610dfccc68fSToby Isaac   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
1611dfccc68fSToby Isaac   if (coordDim > 3) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported coordinate dimension %D > 3", coordDim);
16129c3cf19fSMatthew G. Knepley   if (quad) {ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, &points, NULL);CHKERRQ(ierr);}
1613ba2698f1SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, cell, &ct);CHKERRQ(ierr);
1614ba2698f1SMatthew G. Knepley   switch (ct) {
1615ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_POINT:
1616dfccc68fSToby Isaac     ierr = DMPlexComputePointGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr);
1617dfccc68fSToby Isaac     isAffine = PETSC_FALSE;
1618dfccc68fSToby Isaac     break;
1619ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_SEGMENT:
1620412e9a14SMatthew G. Knepley     case DM_POLYTOPE_POINT_PRISM_TENSOR:
1621ba2698f1SMatthew G. Knepley     if (Nq) {ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr);}
1622ba2698f1SMatthew G. Knepley     else    {ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v,  J,  invJ,  detJ);CHKERRQ(ierr);}
1623dfccc68fSToby Isaac     break;
1624ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:
1625ba2698f1SMatthew G. Knepley     if (Nq) {ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr);}
1626ba2698f1SMatthew G. Knepley     else    {ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v,  J,  invJ,  detJ);CHKERRQ(ierr);}
1627dfccc68fSToby Isaac     break;
1628ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_QUADRILATERAL:
1629412e9a14SMatthew G. Knepley     ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, PETSC_FALSE, Nq, points, v, J, invJ, detJ);CHKERRQ(ierr);
1630412e9a14SMatthew G. Knepley     isAffine = PETSC_FALSE;
1631412e9a14SMatthew G. Knepley     break;
1632412e9a14SMatthew G. Knepley     case DM_POLYTOPE_SEG_PRISM_TENSOR:
1633412e9a14SMatthew G. Knepley     ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, PETSC_TRUE, Nq, points, v, J, invJ, detJ);CHKERRQ(ierr);
1634dfccc68fSToby Isaac     isAffine = PETSC_FALSE;
1635dfccc68fSToby Isaac     break;
1636ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_TETRAHEDRON:
1637ba2698f1SMatthew G. Knepley     if (Nq) {ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr);}
1638ba2698f1SMatthew G. Knepley     else    {ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v,  J,  invJ,  detJ);CHKERRQ(ierr);}
1639dfccc68fSToby Isaac     break;
1640ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_HEXAHEDRON:
1641dfccc68fSToby Isaac     ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, Nq, points, v, J, invJ, detJ);CHKERRQ(ierr);
1642dfccc68fSToby Isaac     isAffine = PETSC_FALSE;
1643dfccc68fSToby Isaac     break;
1644ba2698f1SMatthew G. Knepley     default: SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No element geometry for cell %D with type %s", cell, DMPolytopeTypes[ct]);
1645dfccc68fSToby Isaac   }
16467318780aSToby Isaac   if (isAffine && Nq) {
1647dfccc68fSToby Isaac     if (v) {
1648dfccc68fSToby Isaac       for (i = 0; i < Nq; i++) {
1649c330f8ffSToby Isaac         CoordinatesRefToReal(coordDim, dim, xi0, v0, J0, &points[dim * i], &v[coordDim * i]);
1650dfccc68fSToby Isaac       }
1651dfccc68fSToby Isaac     }
16527318780aSToby Isaac     if (detJ) {
16537318780aSToby Isaac       for (i = 0; i < Nq; i++) {
16547318780aSToby Isaac         detJ[i] = detJ0;
1655dfccc68fSToby Isaac       }
16567318780aSToby Isaac     }
16577318780aSToby Isaac     if (J) {
16587318780aSToby Isaac       PetscInt k;
16597318780aSToby Isaac 
16607318780aSToby Isaac       for (i = 0, k = 0; i < Nq; i++) {
1661dfccc68fSToby Isaac         PetscInt j;
1662dfccc68fSToby Isaac 
16637318780aSToby Isaac         for (j = 0; j < coordDim * coordDim; j++, k++) {
16647318780aSToby Isaac           J[k] = J0[j];
16657318780aSToby Isaac         }
16667318780aSToby Isaac       }
16677318780aSToby Isaac     }
16687318780aSToby Isaac     if (invJ) {
16697318780aSToby Isaac       PetscInt k;
16707318780aSToby Isaac       switch (coordDim) {
16717318780aSToby Isaac       case 0:
16727318780aSToby Isaac         break;
16737318780aSToby Isaac       case 1:
16747318780aSToby Isaac         invJ[0] = 1./J0[0];
16757318780aSToby Isaac         break;
16767318780aSToby Isaac       case 2:
16777318780aSToby Isaac         DMPlex_Invert2D_Internal(invJ, J0, detJ0);
16787318780aSToby Isaac         break;
16797318780aSToby Isaac       case 3:
16807318780aSToby Isaac         DMPlex_Invert3D_Internal(invJ, J0, detJ0);
16817318780aSToby Isaac         break;
16827318780aSToby Isaac       }
16837318780aSToby Isaac       for (i = 1, k = coordDim * coordDim; i < Nq; i++) {
16847318780aSToby Isaac         PetscInt j;
16857318780aSToby Isaac 
16867318780aSToby Isaac         for (j = 0; j < coordDim * coordDim; j++, k++) {
16877318780aSToby Isaac           invJ[k] = invJ[j];
16887318780aSToby Isaac         }
1689dfccc68fSToby Isaac       }
1690dfccc68fSToby Isaac     }
1691dfccc68fSToby Isaac   }
1692dfccc68fSToby Isaac   PetscFunctionReturn(0);
1693dfccc68fSToby Isaac }
1694dfccc68fSToby Isaac 
1695ccd2543fSMatthew G Knepley /*@C
16968e0841e0SMatthew G. Knepley   DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell
1697ccd2543fSMatthew G Knepley 
1698d083f849SBarry Smith   Collective on dm
1699ccd2543fSMatthew G Knepley 
17004165533cSJose E. Roman   Input Parameters:
1701ccd2543fSMatthew G Knepley + dm   - the DM
1702ccd2543fSMatthew G Knepley - cell - the cell
1703ccd2543fSMatthew G Knepley 
17044165533cSJose E. Roman   Output Parameters:
1705ccd2543fSMatthew G Knepley + v0   - the translation part of this affine transform
1706ccd2543fSMatthew G Knepley . J    - the Jacobian of the transform from the reference element
1707ccd2543fSMatthew G Knepley . invJ - the inverse of the Jacobian
1708ccd2543fSMatthew G Knepley - detJ - the Jacobian determinant
1709ccd2543fSMatthew G Knepley 
1710ccd2543fSMatthew G Knepley   Level: advanced
1711ccd2543fSMatthew G Knepley 
1712ccd2543fSMatthew G Knepley   Fortran Notes:
1713ccd2543fSMatthew G Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
1714ccd2543fSMatthew G Knepley   include petsc.h90 in your code.
1715ccd2543fSMatthew G Knepley 
1716e8964c0aSStefano Zampini .seealso: DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinates()
1717ccd2543fSMatthew G Knepley @*/
17188e0841e0SMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
1719ccd2543fSMatthew G Knepley {
1720ccd2543fSMatthew G Knepley   PetscErrorCode ierr;
1721ccd2543fSMatthew G Knepley 
1722ccd2543fSMatthew G Knepley   PetscFunctionBegin;
1723dfccc68fSToby Isaac   ierr = DMPlexComputeCellGeometryFEM_Implicit(dm,cell,NULL,v0,J,invJ,detJ);CHKERRQ(ierr);
17248e0841e0SMatthew G. Knepley   PetscFunctionReturn(0);
17258e0841e0SMatthew G. Knepley }
17268e0841e0SMatthew G. Knepley 
1727dfccc68fSToby Isaac static PetscErrorCode DMPlexComputeCellGeometryFEM_FE(DM dm, PetscFE fe, PetscInt point, PetscQuadrature quad, PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
17288e0841e0SMatthew G. Knepley {
1729dfccc68fSToby Isaac   PetscQuadrature   feQuad;
17308e0841e0SMatthew G. Knepley   PetscSection      coordSection;
17318e0841e0SMatthew G. Knepley   Vec               coordinates;
17328e0841e0SMatthew G. Knepley   PetscScalar      *coords = NULL;
17338e0841e0SMatthew G. Knepley   const PetscReal  *quadPoints;
1734ef0bb6c7SMatthew G. Knepley   PetscTabulation T;
1735f960e424SToby Isaac   PetscInt          dim, cdim, pdim, qdim, Nq, numCoords, q;
17368e0841e0SMatthew G. Knepley   PetscErrorCode    ierr;
17378e0841e0SMatthew G. Knepley 
17388e0841e0SMatthew G. Knepley   PetscFunctionBegin;
17398e0841e0SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
17408e0841e0SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
17418e0841e0SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr);
17428e0841e0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
17438e0841e0SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
1744dfccc68fSToby Isaac   if (!quad) { /* use the first point of the first functional of the dual space */
1745dfccc68fSToby Isaac     PetscDualSpace dsp;
1746dfccc68fSToby Isaac 
1747dfccc68fSToby Isaac     ierr = PetscFEGetDualSpace(fe, &dsp);CHKERRQ(ierr);
1748dfccc68fSToby Isaac     ierr = PetscDualSpaceGetFunctional(dsp, 0, &quad);CHKERRQ(ierr);
17499c3cf19fSMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, &qdim, NULL, &Nq, &quadPoints, NULL);CHKERRQ(ierr);
1750dfccc68fSToby Isaac     Nq = 1;
1751dfccc68fSToby Isaac   } else {
17529c3cf19fSMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, &qdim, NULL, &Nq, &quadPoints, NULL);CHKERRQ(ierr);
1753dfccc68fSToby Isaac   }
175491d2b7ceSToby Isaac   ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr);
1755dfccc68fSToby Isaac   ierr = PetscFEGetQuadrature(fe, &feQuad);CHKERRQ(ierr);
1756dfccc68fSToby Isaac   if (feQuad == quad) {
1757f9244615SMatthew G. Knepley     ierr = PetscFEGetCellTabulation(fe, J ? 1 : 0, &T);CHKERRQ(ierr);
17588e0841e0SMatthew 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);
1759dfccc68fSToby Isaac   } else {
1760ef0bb6c7SMatthew G. Knepley     ierr = PetscFECreateTabulation(fe, 1, Nq, quadPoints, J ? 1 : 0, &T);CHKERRQ(ierr);
1761dfccc68fSToby Isaac   }
1762dfccc68fSToby Isaac   if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %d != quadrature dimension %d", dim, qdim);
1763ef0bb6c7SMatthew G. Knepley   {
1764ef0bb6c7SMatthew G. Knepley     const PetscReal *basis    = T->T[0];
1765ef0bb6c7SMatthew G. Knepley     const PetscReal *basisDer = T->T[1];
1766ef0bb6c7SMatthew G. Knepley     PetscReal        detJt;
1767ef0bb6c7SMatthew G. Knepley 
1768a2a9e04cSMatthew G. Knepley #if defined(PETSC_USE_DEBUG)
1769a2a9e04cSMatthew G. Knepley     if (Nq != T->Np)     SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Np %D != %D", Nq, T->Np);
1770a2a9e04cSMatthew G. Knepley     if (pdim != T->Nb)   SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nb %D != %D", pdim, T->Nb);
1771a2a9e04cSMatthew G. Knepley     if (dim != T->Nc)    SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nc %D != %D", dim, T->Nc);
1772a2a9e04cSMatthew G. Knepley     if (cdim != T->cdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "cdim %D != %D", cdim, T->cdim);
1773a2a9e04cSMatthew G. Knepley #endif
1774dfccc68fSToby Isaac     if (v) {
1775580bdb30SBarry Smith       ierr = PetscArrayzero(v, Nq*cdim);CHKERRQ(ierr);
1776f960e424SToby Isaac       for (q = 0; q < Nq; ++q) {
1777f960e424SToby Isaac         PetscInt i, k;
1778f960e424SToby Isaac 
1779301b184aSMatthew G. Knepley         for (k = 0; k < pdim; ++k) {
1780301b184aSMatthew G. Knepley           const PetscInt vertex = k/cdim;
1781301b184aSMatthew G. Knepley           for (i = 0; i < cdim; ++i) {
1782301b184aSMatthew G. Knepley             v[q*cdim + i] += basis[(q*pdim + k)*cdim + i] * PetscRealPart(coords[vertex*cdim + i]);
1783301b184aSMatthew G. Knepley           }
1784301b184aSMatthew G. Knepley         }
1785f960e424SToby Isaac         ierr = PetscLogFlops(2.0*pdim*cdim);CHKERRQ(ierr);
1786f960e424SToby Isaac       }
1787f960e424SToby Isaac     }
17888e0841e0SMatthew G. Knepley     if (J) {
1789580bdb30SBarry Smith       ierr = PetscArrayzero(J, Nq*cdim*cdim);CHKERRQ(ierr);
17908e0841e0SMatthew G. Knepley       for (q = 0; q < Nq; ++q) {
17918e0841e0SMatthew G. Knepley         PetscInt i, j, k, c, r;
17928e0841e0SMatthew G. Knepley 
17938e0841e0SMatthew G. Knepley         /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */
1794301b184aSMatthew G. Knepley         for (k = 0; k < pdim; ++k) {
1795301b184aSMatthew G. Knepley           const PetscInt vertex = k/cdim;
1796301b184aSMatthew G. Knepley           for (j = 0; j < dim; ++j) {
1797301b184aSMatthew G. Knepley             for (i = 0; i < cdim; ++i) {
1798301b184aSMatthew G. Knepley               J[(q*cdim + i)*cdim + j] += basisDer[((q*pdim + k)*cdim + i)*dim + j] * PetscRealPart(coords[vertex*cdim + i]);
1799301b184aSMatthew G. Knepley             }
1800301b184aSMatthew G. Knepley           }
1801301b184aSMatthew G. Knepley         }
18023bc0b13bSBarry Smith         ierr = PetscLogFlops(2.0*pdim*dim*cdim);CHKERRQ(ierr);
18038e0841e0SMatthew G. Knepley         if (cdim > dim) {
18048e0841e0SMatthew G. Knepley           for (c = dim; c < cdim; ++c)
18058e0841e0SMatthew G. Knepley             for (r = 0; r < cdim; ++r)
18068e0841e0SMatthew G. Knepley               J[r*cdim+c] = r == c ? 1.0 : 0.0;
18078e0841e0SMatthew G. Knepley         }
1808f960e424SToby Isaac         if (!detJ && !invJ) continue;
1809a63b72c6SToby Isaac         detJt = 0.;
18108e0841e0SMatthew G. Knepley         switch (cdim) {
18118e0841e0SMatthew G. Knepley         case 3:
1812037dc194SToby Isaac           DMPlex_Det3D_Internal(&detJt, &J[q*cdim*dim]);
1813037dc194SToby Isaac           if (invJ) {DMPlex_Invert3D_Internal(&invJ[q*cdim*dim], &J[q*cdim*dim], detJt);}
181417fe8556SMatthew G. Knepley           break;
181549dc4407SMatthew G. Knepley         case 2:
18169f328543SToby Isaac           DMPlex_Det2D_Internal(&detJt, &J[q*cdim*dim]);
1817037dc194SToby Isaac           if (invJ) {DMPlex_Invert2D_Internal(&invJ[q*cdim*dim], &J[q*cdim*dim], detJt);}
181849dc4407SMatthew G. Knepley           break;
18198e0841e0SMatthew G. Knepley         case 1:
1820037dc194SToby Isaac           detJt = J[q*cdim*dim];
1821037dc194SToby Isaac           if (invJ) invJ[q*cdim*dim] = 1.0/detJt;
182249dc4407SMatthew G. Knepley         }
1823f960e424SToby Isaac         if (detJ) detJ[q] = detJt;
182449dc4407SMatthew G. Knepley       }
1825ef0bb6c7SMatthew G. Knepley     } else if (detJ || invJ) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Need J to compute invJ or detJ");
182649dc4407SMatthew G. Knepley   }
1827ef0bb6c7SMatthew G. Knepley   if (feQuad != quad) {ierr = PetscTabulationDestroy(&T);CHKERRQ(ierr);}
18288e0841e0SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr);
18298e0841e0SMatthew G. Knepley   PetscFunctionReturn(0);
18308e0841e0SMatthew G. Knepley }
18318e0841e0SMatthew G. Knepley 
18328e0841e0SMatthew G. Knepley /*@C
18338e0841e0SMatthew G. Knepley   DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell
18348e0841e0SMatthew G. Knepley 
1835d083f849SBarry Smith   Collective on dm
18368e0841e0SMatthew G. Knepley 
18374165533cSJose E. Roman   Input Parameters:
18388e0841e0SMatthew G. Knepley + dm   - the DM
18398e0841e0SMatthew G. Knepley . cell - the cell
1840dfccc68fSToby Isaac - quad - the quadrature containing the points in the reference element where the geometry will be evaluated.  If quad == NULL, geometry will be
1841dfccc68fSToby Isaac          evaluated at the first vertex of the reference element
18428e0841e0SMatthew G. Knepley 
18434165533cSJose E. Roman   Output Parameters:
1844dfccc68fSToby Isaac + v    - the image of the transformed quadrature points, otherwise the image of the first vertex in the closure of the reference element
18458e0841e0SMatthew G. Knepley . J    - the Jacobian of the transform from the reference element at each quadrature point
18468e0841e0SMatthew G. Knepley . invJ - the inverse of the Jacobian at each quadrature point
18478e0841e0SMatthew G. Knepley - detJ - the Jacobian determinant at each quadrature point
18488e0841e0SMatthew G. Knepley 
18498e0841e0SMatthew G. Knepley   Level: advanced
18508e0841e0SMatthew G. Knepley 
18518e0841e0SMatthew G. Knepley   Fortran Notes:
18528e0841e0SMatthew G. Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
18538e0841e0SMatthew G. Knepley   include petsc.h90 in your code.
18548e0841e0SMatthew G. Knepley 
1855e8964c0aSStefano Zampini .seealso: DMGetCoordinateSection(), DMGetCoordinates()
18568e0841e0SMatthew G. Knepley @*/
1857dfccc68fSToby Isaac PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
18588e0841e0SMatthew G. Knepley {
1859bb4a5db5SMatthew G. Knepley   DM             cdm;
1860dfccc68fSToby Isaac   PetscFE        fe = NULL;
18618e0841e0SMatthew G. Knepley   PetscErrorCode ierr;
18628e0841e0SMatthew G. Knepley 
18638e0841e0SMatthew G. Knepley   PetscFunctionBegin;
18642d89661fSMatthew G. Knepley   PetscValidPointer(detJ, 7);
1865bb4a5db5SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
1866bb4a5db5SMatthew G. Knepley   if (cdm) {
1867dfccc68fSToby Isaac     PetscClassId id;
1868dfccc68fSToby Isaac     PetscInt     numFields;
1869e5e52638SMatthew G. Knepley     PetscDS      prob;
1870dfccc68fSToby Isaac     PetscObject  disc;
1871dfccc68fSToby Isaac 
1872bb4a5db5SMatthew G. Knepley     ierr = DMGetNumFields(cdm, &numFields);CHKERRQ(ierr);
1873dfccc68fSToby Isaac     if (numFields) {
1874bb4a5db5SMatthew G. Knepley       ierr = DMGetDS(cdm, &prob);CHKERRQ(ierr);
1875dfccc68fSToby Isaac       ierr = PetscDSGetDiscretization(prob,0,&disc);CHKERRQ(ierr);
1876dfccc68fSToby Isaac       ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
1877dfccc68fSToby Isaac       if (id == PETSCFE_CLASSID) {
1878dfccc68fSToby Isaac         fe = (PetscFE) disc;
1879dfccc68fSToby Isaac       }
1880dfccc68fSToby Isaac     }
1881dfccc68fSToby Isaac   }
1882dfccc68fSToby Isaac   if (!fe) {ierr = DMPlexComputeCellGeometryFEM_Implicit(dm, cell, quad, v, J, invJ, detJ);CHKERRQ(ierr);}
1883dfccc68fSToby Isaac   else     {ierr = DMPlexComputeCellGeometryFEM_FE(dm, fe, cell, quad, v, J, invJ, detJ);CHKERRQ(ierr);}
1884ccd2543fSMatthew G Knepley   PetscFunctionReturn(0);
1885ccd2543fSMatthew G Knepley }
1886834e62ceSMatthew G. Knepley 
1887011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1888cc08537eSMatthew G. Knepley {
1889cc08537eSMatthew G. Knepley   PetscSection   coordSection;
1890cc08537eSMatthew G. Knepley   Vec            coordinates;
1891a1e44745SMatthew G. Knepley   PetscScalar   *coords = NULL;
189206e2781eSMatthew G. Knepley   PetscScalar    tmp[2];
1893714b99b6SMatthew G. Knepley   PetscInt       coordSize, d;
1894cc08537eSMatthew G. Knepley   PetscErrorCode ierr;
1895cc08537eSMatthew G. Knepley 
1896cc08537eSMatthew G. Knepley   PetscFunctionBegin;
1897cc08537eSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
189869d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1899cc08537eSMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
19002e17dfb7SMatthew G. Knepley   ierr = DMLocalizeCoordinate_Internal(dm, dim, coords, &coords[dim], tmp);CHKERRQ(ierr);
1901cc08537eSMatthew G. Knepley   if (centroid) {
1902714b99b6SMatthew G. Knepley     for (d = 0; d < dim; ++d) centroid[d] = 0.5*PetscRealPart(coords[d] + tmp[d]);
1903cc08537eSMatthew G. Knepley   }
1904cc08537eSMatthew G. Knepley   if (normal) {
1905a60a936bSMatthew G. Knepley     PetscReal norm;
1906a60a936bSMatthew G. Knepley 
1907714b99b6SMatthew G. Knepley     if (dim != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "We only support 2D edges right now");
190806e2781eSMatthew G. Knepley     normal[0]  = -PetscRealPart(coords[1] - tmp[1]);
190906e2781eSMatthew G. Knepley     normal[1]  =  PetscRealPart(coords[0] - tmp[0]);
1910714b99b6SMatthew G. Knepley     norm       = DMPlex_NormD_Internal(dim, normal);
1911714b99b6SMatthew G. Knepley     for (d = 0; d < dim; ++d) normal[d] /= norm;
1912cc08537eSMatthew G. Knepley   }
1913cc08537eSMatthew G. Knepley   if (vol) {
1914714b99b6SMatthew G. Knepley     *vol = 0.0;
1915714b99b6SMatthew G. Knepley     for (d = 0; d < dim; ++d) *vol += PetscSqr(PetscRealPart(coords[d] - tmp[d]));
1916714b99b6SMatthew G. Knepley     *vol = PetscSqrtReal(*vol);
1917cc08537eSMatthew G. Knepley   }
1918cc08537eSMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1919cc08537eSMatthew G. Knepley   PetscFunctionReturn(0);
1920cc08537eSMatthew G. Knepley }
1921cc08537eSMatthew G. Knepley 
1922cc08537eSMatthew G. Knepley /* Centroid_i = (\sum_n A_n Cn_i) / A */
1923011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1924cc08537eSMatthew G. Knepley {
1925412e9a14SMatthew G. Knepley   DMPolytopeType ct;
1926cc08537eSMatthew G. Knepley   PetscSection   coordSection;
1927cc08537eSMatthew G. Knepley   Vec            coordinates;
1928cc08537eSMatthew G. Knepley   PetscScalar   *coords = NULL;
1929793a2a13SMatthew G. Knepley   PetscInt       fv[4] = {0, 1, 2, 3};
1930*4f99dae5SMatthew G. Knepley   PetscInt       cdim, coordSize, numCorners, p, d;
1931cc08537eSMatthew G. Knepley   PetscErrorCode ierr;
1932cc08537eSMatthew G. Knepley 
1933cc08537eSMatthew G. Knepley   PetscFunctionBegin;
1934793a2a13SMatthew G. Knepley   /* Must check for hybrid cells because prisms have a different orientation scheme */
1935412e9a14SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, cell, &ct);CHKERRQ(ierr);
1936412e9a14SMatthew G. Knepley   switch (ct) {
1937*4f99dae5SMatthew G. Knepley     case DM_POLYTOPE_SEG_PRISM_TENSOR: fv[2] = 3; fv[3] = 2;break;
1938412e9a14SMatthew G. Knepley     default: break;
1939412e9a14SMatthew G. Knepley   }
1940cc08537eSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
19410a1d6728SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr);
194269d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1943cc08537eSMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1944*4f99dae5SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
1945793a2a13SMatthew G. Knepley 
1946*4f99dae5SMatthew G. Knepley   if (cdim > 2) {
1947*4f99dae5SMatthew G. Knepley     PetscReal c[3] = {0., 0., 0.}, n[3] = {0., 0., 0.}, norm;
1948*4f99dae5SMatthew G. Knepley 
1949*4f99dae5SMatthew G. Knepley     for (p = 0; p < numCorners-2; ++p) {
1950*4f99dae5SMatthew G. Knepley       const PetscReal x0 = PetscRealPart(coords[cdim*fv[p+1]+0] - coords[0]), x1 = PetscRealPart(coords[cdim*fv[p+2]+0] - coords[0]);
1951*4f99dae5SMatthew G. Knepley       const PetscReal y0 = PetscRealPart(coords[cdim*fv[p+1]+1] - coords[1]), y1 = PetscRealPart(coords[cdim*fv[p+2]+1] - coords[1]);
1952*4f99dae5SMatthew G. Knepley       const PetscReal z0 = PetscRealPart(coords[cdim*fv[p+1]+2] - coords[2]), z1 = PetscRealPart(coords[cdim*fv[p+2]+2] - coords[2]);
1953*4f99dae5SMatthew G. Knepley       const PetscReal dx = y0*z1 - z0*y1;
1954*4f99dae5SMatthew G. Knepley       const PetscReal dy = z0*x1 - x0*z1;
1955*4f99dae5SMatthew G. Knepley       const PetscReal dz = x0*y1 - y0*x1;
1956*4f99dae5SMatthew G. Knepley       PetscReal       a  = PetscSqrtReal(dx*dx + dy*dy + dz*dz);
1957*4f99dae5SMatthew G. Knepley 
1958*4f99dae5SMatthew G. Knepley       n[0] += dx;
1959*4f99dae5SMatthew G. Knepley       n[1] += dy;
1960*4f99dae5SMatthew G. Knepley       n[2] += dz;
1961*4f99dae5SMatthew G. Knepley       c[0] += a * PetscRealPart(coords[0] + coords[cdim*fv[p+1]+0] + coords[cdim*fv[p+2]+0])/3.;
1962*4f99dae5SMatthew G. Knepley       c[1] += a * PetscRealPart(coords[1] + coords[cdim*fv[p+1]+1] + coords[cdim*fv[p+2]+1])/3.;
1963*4f99dae5SMatthew G. Knepley       c[2] += a * PetscRealPart(coords[2] + coords[cdim*fv[p+1]+2] + coords[cdim*fv[p+2]+2])/3.;
1964ceee4971SMatthew G. Knepley     }
1965*4f99dae5SMatthew G. Knepley     norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]);
1966*4f99dae5SMatthew G. Knepley     n[0] /= norm;
1967*4f99dae5SMatthew G. Knepley     n[1] /= norm;
1968*4f99dae5SMatthew G. Knepley     n[2] /= norm;
1969*4f99dae5SMatthew G. Knepley     c[0] /= norm;
1970*4f99dae5SMatthew G. Knepley     c[1] /= norm;
1971*4f99dae5SMatthew G. Knepley     c[2] /= norm;
1972*4f99dae5SMatthew G. Knepley     if (vol) *vol = 0.5*norm;
1973*4f99dae5SMatthew G. Knepley     if (centroid) for (d = 0; d < cdim; ++d) centroid[d] = c[d];
1974*4f99dae5SMatthew G. Knepley     if (normal) for (d = 0; d < cdim; ++d) normal[d] = n[d];
1975011ea5d8SMatthew G. Knepley   } else {
1976*4f99dae5SMatthew G. Knepley     PetscReal vsum = 0.0, csum[2] = {0.0, 0.0}, vtmp, ctmp[4] = {0., 0., 0., 0.};
1977*4f99dae5SMatthew G. Knepley 
19780a1d6728SMatthew G. Knepley     for (p = 0; p < numCorners; ++p) {
1979793a2a13SMatthew G. Knepley       const PetscInt pi  = p < 4 ? fv[p] : p;
1980793a2a13SMatthew G. Knepley       const PetscInt pin = p < 3 ? fv[(p+1)%numCorners] : (p+1)%numCorners;
19810a1d6728SMatthew G. Knepley       /* Need to do this copy to get types right */
1982*4f99dae5SMatthew G. Knepley       for (d = 0; d < cdim; ++d) {
1983*4f99dae5SMatthew G. Knepley         ctmp[d]      = PetscRealPart(coords[pi*cdim+d]);
1984*4f99dae5SMatthew G. Knepley         ctmp[cdim+d] = PetscRealPart(coords[pin*cdim+d]);
19850a1d6728SMatthew G. Knepley       }
19860a1d6728SMatthew G. Knepley       Volume_Triangle_Origin_Internal(&vtmp, ctmp);
19870a1d6728SMatthew G. Knepley       vsum += vtmp;
1988*4f99dae5SMatthew G. Knepley       for (d = 0; d < cdim; ++d) csum[d] += (ctmp[d] + ctmp[cdim+d])*vtmp;
19890a1d6728SMatthew G. Knepley     }
1990*4f99dae5SMatthew G. Knepley     if (vol) *vol = PetscAbsReal(vsum);
1991*4f99dae5SMatthew G. Knepley     if (centroid) for (d = 0; d < cdim; ++d) centroid[d] = csum[d] / ((cdim+1)*vsum);
1992*4f99dae5SMatthew G. Knepley     if (normal) for (d = 0; d < cdim; ++d) normal[d] = 0.0;
19930a1d6728SMatthew G. Knepley   }
19940a1d6728SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
1995cc08537eSMatthew G. Knepley   PetscFunctionReturn(0);
1996cc08537eSMatthew G. Knepley }
1997cc08537eSMatthew G. Knepley 
19980ec8681fSMatthew G. Knepley /* Centroid_i = (\sum_n V_n Cn_i) / V */
1999011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
20000ec8681fSMatthew G. Knepley {
2001412e9a14SMatthew G. Knepley   DMPolytopeType  ct;
20020ec8681fSMatthew G. Knepley   PetscSection    coordSection;
20030ec8681fSMatthew G. Knepley   Vec             coordinates;
20040ec8681fSMatthew G. Knepley   PetscScalar    *coords = NULL;
200586623015SMatthew G. Knepley   PetscReal       vsum = 0.0, vtmp, coordsTmp[3*3];
2006a7df9edeSMatthew G. Knepley   const PetscInt *faces, *facesO;
2007793a2a13SMatthew G. Knepley   PetscBool       isHybrid = PETSC_FALSE;
2008412e9a14SMatthew G. Knepley   PetscInt        numFaces, f, coordSize, p, d;
20090ec8681fSMatthew G. Knepley   PetscErrorCode  ierr;
20100ec8681fSMatthew G. Knepley 
20110ec8681fSMatthew G. Knepley   PetscFunctionBegin;
2012f6dae198SJed Brown   if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No support for dim %D > 3",dim);
2013793a2a13SMatthew G. Knepley   /* Must check for hybrid cells because prisms have a different orientation scheme */
2014412e9a14SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, cell, &ct);CHKERRQ(ierr);
2015412e9a14SMatthew G. Knepley   switch (ct) {
2016412e9a14SMatthew G. Knepley     case DM_POLYTOPE_POINT_PRISM_TENSOR:
2017412e9a14SMatthew G. Knepley     case DM_POLYTOPE_SEG_PRISM_TENSOR:
2018412e9a14SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM_TENSOR:
2019412e9a14SMatthew G. Knepley     case DM_POLYTOPE_QUAD_PRISM_TENSOR:
2020412e9a14SMatthew G. Knepley       isHybrid = PETSC_TRUE;
2021412e9a14SMatthew G. Knepley     default: break;
2022412e9a14SMatthew G. Knepley   }
2023793a2a13SMatthew G. Knepley 
20240ec8681fSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
202569d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
20260ec8681fSMatthew G. Knepley 
2027d9a81ebdSMatthew G. Knepley   if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0;
20280ec8681fSMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr);
20290ec8681fSMatthew G. Knepley   ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr);
2030a7df9edeSMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr);
20310ec8681fSMatthew G. Knepley   for (f = 0; f < numFaces; ++f) {
2032793a2a13SMatthew G. Knepley     PetscBool      flip = isHybrid && f == 0 ? PETSC_TRUE : PETSC_FALSE; /* The first hybrid face is reversed */
2033ba2698f1SMatthew G. Knepley     DMPolytopeType ct;
2034793a2a13SMatthew G. Knepley 
2035011ea5d8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr);
2036ba2698f1SMatthew G. Knepley     ierr = DMPlexGetCellType(dm, faces[f], &ct);CHKERRQ(ierr);
2037ba2698f1SMatthew G. Knepley     switch (ct) {
2038ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:
20390ec8681fSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
20401ee9d5ecSMatthew G. Knepley         coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]);
20411ee9d5ecSMatthew G. Knepley         coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]);
20421ee9d5ecSMatthew G. Knepley         coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]);
20430ec8681fSMatthew G. Knepley       }
20440ec8681fSMatthew G. Knepley       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
2045793a2a13SMatthew G. Knepley       if (facesO[f] < 0 || flip) vtmp = -vtmp;
20460ec8681fSMatthew G. Knepley       vsum += vtmp;
20474f25033aSJed Brown       if (centroid) {           /* Centroid of OABC = (a+b+c)/4 */
20480ec8681fSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
20491ee9d5ecSMatthew G. Knepley           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
20500ec8681fSMatthew G. Knepley         }
20510ec8681fSMatthew G. Knepley       }
20520ec8681fSMatthew G. Knepley       break;
2053ba2698f1SMatthew G. Knepley     case DM_POLYTOPE_QUADRILATERAL:
2054412e9a14SMatthew G. Knepley     case DM_POLYTOPE_SEG_PRISM_TENSOR:
2055793a2a13SMatthew G. Knepley     {
2056793a2a13SMatthew G. Knepley       PetscInt fv[4] = {0, 1, 2, 3};
2057793a2a13SMatthew G. Knepley 
2058793a2a13SMatthew G. Knepley       /* Side faces for hybrid cells are are stored as tensor products */
2059793a2a13SMatthew G. Knepley       if (isHybrid && f > 1) {fv[2] = 3; fv[3] = 2;}
20600ec8681fSMatthew G. Knepley       /* DO FOR PYRAMID */
20610ec8681fSMatthew G. Knepley       /* First tet */
20620ec8681fSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
2063793a2a13SMatthew G. Knepley         coordsTmp[0*dim+d] = PetscRealPart(coords[fv[0]*dim+d]);
2064793a2a13SMatthew G. Knepley         coordsTmp[1*dim+d] = PetscRealPart(coords[fv[1]*dim+d]);
2065793a2a13SMatthew G. Knepley         coordsTmp[2*dim+d] = PetscRealPart(coords[fv[3]*dim+d]);
20660ec8681fSMatthew G. Knepley       }
20670ec8681fSMatthew G. Knepley       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
2068793a2a13SMatthew G. Knepley       if (facesO[f] < 0 || flip) vtmp = -vtmp;
20690ec8681fSMatthew G. Knepley       vsum += vtmp;
20700ec8681fSMatthew G. Knepley       if (centroid) {
20710ec8681fSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
20720ec8681fSMatthew G. Knepley           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
20730ec8681fSMatthew G. Knepley         }
20740ec8681fSMatthew G. Knepley       }
20750ec8681fSMatthew G. Knepley       /* Second tet */
20760ec8681fSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
2077793a2a13SMatthew G. Knepley         coordsTmp[0*dim+d] = PetscRealPart(coords[fv[1]*dim+d]);
2078793a2a13SMatthew G. Knepley         coordsTmp[1*dim+d] = PetscRealPart(coords[fv[2]*dim+d]);
2079793a2a13SMatthew G. Knepley         coordsTmp[2*dim+d] = PetscRealPart(coords[fv[3]*dim+d]);
20800ec8681fSMatthew G. Knepley       }
20810ec8681fSMatthew G. Knepley       Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
2082793a2a13SMatthew G. Knepley       if (facesO[f] < 0 || flip) vtmp = -vtmp;
20830ec8681fSMatthew G. Knepley       vsum += vtmp;
20840ec8681fSMatthew G. Knepley       if (centroid) {
20850ec8681fSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
20860ec8681fSMatthew G. Knepley           for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
20870ec8681fSMatthew G. Knepley         }
20880ec8681fSMatthew G. Knepley       }
20890ec8681fSMatthew G. Knepley       break;
2090793a2a13SMatthew G. Knepley     }
20910ec8681fSMatthew G. Knepley     default:
2092ba2698f1SMatthew G. Knepley       SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle face %D of type %s", faces[f], DMPolytopeTypes[ct]);
20930ec8681fSMatthew G. Knepley     }
20944f25033aSJed Brown     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr);
20950ec8681fSMatthew G. Knepley   }
20968763be8eSMatthew G. Knepley   if (vol)     *vol = PetscAbsReal(vsum);
20970ec8681fSMatthew G. Knepley   if (normal)   for (d = 0; d < dim; ++d) normal[d]    = 0.0;
2098d9a81ebdSMatthew G. Knepley   if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4);
20990ec8681fSMatthew G. Knepley   PetscFunctionReturn(0);
21000ec8681fSMatthew G. Knepley }
21010ec8681fSMatthew G. Knepley 
2102834e62ceSMatthew G. Knepley /*@C
2103834e62ceSMatthew G. Knepley   DMPlexComputeCellGeometryFVM - Compute the volume for a given cell
2104834e62ceSMatthew G. Knepley 
2105d083f849SBarry Smith   Collective on dm
2106834e62ceSMatthew G. Knepley 
21074165533cSJose E. Roman   Input Parameters:
2108834e62ceSMatthew G. Knepley + dm   - the DM
2109834e62ceSMatthew G. Knepley - cell - the cell
2110834e62ceSMatthew G. Knepley 
21114165533cSJose E. Roman   Output Parameters:
2112834e62ceSMatthew G. Knepley + volume   - the cell volume
2113cc08537eSMatthew G. Knepley . centroid - the cell centroid
2114cc08537eSMatthew G. Knepley - normal - the cell normal, if appropriate
2115834e62ceSMatthew G. Knepley 
2116834e62ceSMatthew G. Knepley   Level: advanced
2117834e62ceSMatthew G. Knepley 
2118834e62ceSMatthew G. Knepley   Fortran Notes:
2119834e62ceSMatthew G. Knepley   Since it returns arrays, this routine is only available in Fortran 90, and you must
2120834e62ceSMatthew G. Knepley   include petsc.h90 in your code.
2121834e62ceSMatthew G. Knepley 
2122e8964c0aSStefano Zampini .seealso: DMGetCoordinateSection(), DMGetCoordinates()
2123834e62ceSMatthew G. Knepley @*/
2124cc08537eSMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
2125834e62ceSMatthew G. Knepley {
21260ec8681fSMatthew G. Knepley   PetscInt       depth, dim;
2127834e62ceSMatthew G. Knepley   PetscErrorCode ierr;
2128834e62ceSMatthew G. Knepley 
2129834e62ceSMatthew G. Knepley   PetscFunctionBegin;
2130834e62ceSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2131c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2132834e62ceSMatthew G. Knepley   if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated");
2133ba2698f1SMatthew G. Knepley   ierr = DMPlexGetPointDepth(dm, cell, &depth);CHKERRQ(ierr);
2134011ea5d8SMatthew G. Knepley   switch (depth) {
2135cc08537eSMatthew G. Knepley   case 1:
2136011ea5d8SMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
2137cc08537eSMatthew G. Knepley     break;
2138834e62ceSMatthew G. Knepley   case 2:
2139011ea5d8SMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
2140834e62ceSMatthew G. Knepley     break;
2141834e62ceSMatthew G. Knepley   case 3:
2142011ea5d8SMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr);
2143834e62ceSMatthew G. Knepley     break;
2144834e62ceSMatthew G. Knepley   default:
2145b81cf158SMatthew G. Knepley     SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D (depth %D) for element geometry computation", dim, depth);
2146834e62ceSMatthew G. Knepley   }
2147834e62ceSMatthew G. Knepley   PetscFunctionReturn(0);
2148834e62ceSMatthew G. Knepley }
2149113c68e6SMatthew G. Knepley 
2150c501906fSMatthew G. Knepley /*@
2151c501906fSMatthew G. Knepley   DMPlexComputeGeometryFEM - Precompute cell geometry for the entire mesh
2152c501906fSMatthew G. Knepley 
2153c501906fSMatthew G. Knepley   Collective on dm
2154c501906fSMatthew G. Knepley 
2155c501906fSMatthew G. Knepley   Input Parameter:
2156c501906fSMatthew G. Knepley . dm - The DMPlex
2157c501906fSMatthew G. Knepley 
2158c501906fSMatthew G. Knepley   Output Parameter:
2159c501906fSMatthew G. Knepley . cellgeom - A vector with the cell geometry data for each cell
2160c501906fSMatthew G. Knepley 
2161c501906fSMatthew G. Knepley   Level: beginner
2162c501906fSMatthew G. Knepley 
2163c501906fSMatthew G. Knepley @*/
2164c0d900a5SMatthew G. Knepley PetscErrorCode DMPlexComputeGeometryFEM(DM dm, Vec *cellgeom)
2165c0d900a5SMatthew G. Knepley {
2166c0d900a5SMatthew G. Knepley   DM             dmCell;
2167c0d900a5SMatthew G. Knepley   Vec            coordinates;
2168c0d900a5SMatthew G. Knepley   PetscSection   coordSection, sectionCell;
2169c0d900a5SMatthew G. Knepley   PetscScalar   *cgeom;
2170412e9a14SMatthew G. Knepley   PetscInt       cStart, cEnd, c;
2171c0d900a5SMatthew G. Knepley   PetscErrorCode ierr;
2172c0d900a5SMatthew G. Knepley 
2173c0d900a5SMatthew G. Knepley   PetscFunctionBegin;
2174c0d900a5SMatthew G. Knepley   ierr = DMClone(dm, &dmCell);CHKERRQ(ierr);
2175c0d900a5SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2176c0d900a5SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2177c0d900a5SMatthew G. Knepley   ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr);
2178c0d900a5SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr);
2179c0d900a5SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionCell);CHKERRQ(ierr);
2180412e9a14SMatthew G. Knepley   ierr = DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2181c0d900a5SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr);
2182c0d900a5SMatthew G. Knepley   /* TODO This needs to be multiplied by Nq for non-affine */
2183cf0b7c11SKarl Rupp   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFEGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);}
2184c0d900a5SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr);
218592fd8e1eSJed Brown   ierr = DMSetLocalSection(dmCell, sectionCell);CHKERRQ(ierr);
2186c0d900a5SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionCell);CHKERRQ(ierr);
2187c0d900a5SMatthew G. Knepley   ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr);
2188c0d900a5SMatthew G. Knepley   ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr);
2189c0d900a5SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
2190cf0b7c11SKarl Rupp     PetscFEGeom *cg;
2191c0d900a5SMatthew G. Knepley 
2192c0d900a5SMatthew G. Knepley     ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
2193580bdb30SBarry Smith     ierr = PetscArrayzero(cg, 1);CHKERRQ(ierr);
2194cf0b7c11SKarl Rupp     ierr = DMPlexComputeCellGeometryFEM(dmCell, c, NULL, cg->v, cg->J, cg->invJ, cg->detJ);CHKERRQ(ierr);
2195087ef6b2SMatthew 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);
2196c0d900a5SMatthew G. Knepley   }
2197c0d900a5SMatthew G. Knepley   PetscFunctionReturn(0);
2198c0d900a5SMatthew G. Knepley }
2199c0d900a5SMatthew G. Knepley 
2200891a9168SMatthew G. Knepley /*@
2201891a9168SMatthew G. Knepley   DMPlexComputeGeometryFVM - Computes the cell and face geometry for a finite volume method
2202891a9168SMatthew G. Knepley 
2203891a9168SMatthew G. Knepley   Input Parameter:
2204891a9168SMatthew G. Knepley . dm - The DM
2205891a9168SMatthew G. Knepley 
2206891a9168SMatthew G. Knepley   Output Parameters:
2207891a9168SMatthew G. Knepley + cellgeom - A Vec of PetscFVCellGeom data
2208a2b725a8SWilliam Gropp - facegeom - A Vec of PetscFVFaceGeom data
2209891a9168SMatthew G. Knepley 
2210891a9168SMatthew G. Knepley   Level: developer
2211891a9168SMatthew G. Knepley 
2212891a9168SMatthew G. Knepley .seealso: PetscFVFaceGeom, PetscFVCellGeom, DMPlexComputeGeometryFEM()
2213891a9168SMatthew G. Knepley @*/
2214113c68e6SMatthew G. Knepley PetscErrorCode DMPlexComputeGeometryFVM(DM dm, Vec *cellgeom, Vec *facegeom)
2215113c68e6SMatthew G. Knepley {
2216113c68e6SMatthew G. Knepley   DM             dmFace, dmCell;
2217113c68e6SMatthew G. Knepley   DMLabel        ghostLabel;
2218113c68e6SMatthew G. Knepley   PetscSection   sectionFace, sectionCell;
2219113c68e6SMatthew G. Knepley   PetscSection   coordSection;
2220113c68e6SMatthew G. Knepley   Vec            coordinates;
2221113c68e6SMatthew G. Knepley   PetscScalar   *fgeom, *cgeom;
2222113c68e6SMatthew G. Knepley   PetscReal      minradius, gminradius;
2223113c68e6SMatthew G. Knepley   PetscInt       dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f;
2224113c68e6SMatthew G. Knepley   PetscErrorCode ierr;
2225113c68e6SMatthew G. Knepley 
2226113c68e6SMatthew G. Knepley   PetscFunctionBegin;
2227113c68e6SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2228113c68e6SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
2229113c68e6SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
2230113c68e6SMatthew G. Knepley   /* Make cell centroids and volumes */
2231113c68e6SMatthew G. Knepley   ierr = DMClone(dm, &dmCell);CHKERRQ(ierr);
2232113c68e6SMatthew G. Knepley   ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr);
2233113c68e6SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr);
2234113c68e6SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionCell);CHKERRQ(ierr);
2235113c68e6SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2236485ad865SMatthew G. Knepley   ierr = DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);CHKERRQ(ierr);
2237113c68e6SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr);
22389e5edeeeSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVCellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);}
2239113c68e6SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr);
224092fd8e1eSJed Brown   ierr = DMSetLocalSection(dmCell, sectionCell);CHKERRQ(ierr);
2241113c68e6SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionCell);CHKERRQ(ierr);
2242113c68e6SMatthew G. Knepley   ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr);
2243485ad865SMatthew G. Knepley   if (cEndInterior < 0) cEndInterior = cEnd;
2244113c68e6SMatthew G. Knepley   ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr);
2245113c68e6SMatthew G. Knepley   for (c = cStart; c < cEndInterior; ++c) {
2246113c68e6SMatthew G. Knepley     PetscFVCellGeom *cg;
2247113c68e6SMatthew G. Knepley 
2248113c68e6SMatthew G. Knepley     ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
2249580bdb30SBarry Smith     ierr = PetscArrayzero(cg, 1);CHKERRQ(ierr);
2250113c68e6SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr);
2251113c68e6SMatthew G. Knepley   }
2252113c68e6SMatthew G. Knepley   /* Compute face normals and minimum cell radius */
2253113c68e6SMatthew G. Knepley   ierr = DMClone(dm, &dmFace);CHKERRQ(ierr);
2254113c68e6SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionFace);CHKERRQ(ierr);
2255113c68e6SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
2256113c68e6SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr);
22579e5edeeeSMatthew G. Knepley   for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVFaceGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);}
2258113c68e6SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr);
225992fd8e1eSJed Brown   ierr = DMSetLocalSection(dmFace, sectionFace);CHKERRQ(ierr);
2260113c68e6SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionFace);CHKERRQ(ierr);
2261113c68e6SMatthew G. Knepley   ierr = DMCreateLocalVector(dmFace, facegeom);CHKERRQ(ierr);
2262113c68e6SMatthew G. Knepley   ierr = VecGetArray(*facegeom, &fgeom);CHKERRQ(ierr);
2263c58f1c22SToby Isaac   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
2264113c68e6SMatthew G. Knepley   minradius = PETSC_MAX_REAL;
2265113c68e6SMatthew G. Knepley   for (f = fStart; f < fEnd; ++f) {
2266113c68e6SMatthew G. Knepley     PetscFVFaceGeom *fg;
2267113c68e6SMatthew G. Knepley     PetscReal        area;
2268412e9a14SMatthew G. Knepley     const PetscInt  *cells;
2269412e9a14SMatthew G. Knepley     PetscInt         ncells, ghost = -1, d, numChildren;
2270113c68e6SMatthew G. Knepley 
22719ac3fadcSMatthew G. Knepley     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);}
227250d63984SToby Isaac     ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr);
2273412e9a14SMatthew G. Knepley     ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr);
2274412e9a14SMatthew G. Knepley     ierr = DMPlexGetSupportSize(dm, f, &ncells);CHKERRQ(ierr);
2275412e9a14SMatthew G. Knepley     /* It is possible to get a face with no support when using partition overlap */
2276412e9a14SMatthew G. Knepley     if (!ncells || ghost >= 0 || numChildren) continue;
2277113c68e6SMatthew G. Knepley     ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr);
2278113c68e6SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr);
2279113c68e6SMatthew G. Knepley     for (d = 0; d < dim; ++d) fg->normal[d] *= area;
2280113c68e6SMatthew G. Knepley     /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */
2281113c68e6SMatthew G. Knepley     {
2282113c68e6SMatthew G. Knepley       PetscFVCellGeom *cL, *cR;
2283113c68e6SMatthew G. Knepley       PetscReal       *lcentroid, *rcentroid;
22840453c0cdSMatthew G. Knepley       PetscReal        l[3], r[3], v[3];
2285113c68e6SMatthew G. Knepley 
2286113c68e6SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr);
2287113c68e6SMatthew G. Knepley       lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid;
228806348e87SToby Isaac       if (ncells > 1) {
228906348e87SToby Isaac         ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr);
2290113c68e6SMatthew G. Knepley         rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid;
229106348e87SToby Isaac       }
229206348e87SToby Isaac       else {
229306348e87SToby Isaac         rcentroid = fg->centroid;
229406348e87SToby Isaac       }
22952e17dfb7SMatthew G. Knepley       ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, lcentroid, l);CHKERRQ(ierr);
22962e17dfb7SMatthew G. Knepley       ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, rcentroid, r);CHKERRQ(ierr);
22970453c0cdSMatthew G. Knepley       DMPlex_WaxpyD_Internal(dim, -1, l, r, v);
2298113c68e6SMatthew G. Knepley       if (DMPlex_DotRealD_Internal(dim, fg->normal, v) < 0) {
2299113c68e6SMatthew G. Knepley         for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d];
2300113c68e6SMatthew G. Knepley       }
2301113c68e6SMatthew G. Knepley       if (DMPlex_DotRealD_Internal(dim, fg->normal, v) <= 0) {
2302113c68e6SMatthew 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]);
2303113c68e6SMatthew 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]);
2304113c68e6SMatthew G. Knepley         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f);
2305113c68e6SMatthew G. Knepley       }
2306113c68e6SMatthew G. Knepley       if (cells[0] < cEndInterior) {
2307113c68e6SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cL->centroid, v);
2308113c68e6SMatthew G. Knepley         minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v));
2309113c68e6SMatthew G. Knepley       }
231006348e87SToby Isaac       if (ncells > 1 && cells[1] < cEndInterior) {
2311113c68e6SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cR->centroid, v);
2312113c68e6SMatthew G. Knepley         minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v));
2313113c68e6SMatthew G. Knepley       }
2314113c68e6SMatthew G. Knepley     }
2315113c68e6SMatthew G. Knepley   }
2316820f2d46SBarry Smith   ierr = MPIU_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm));CHKERRMPI(ierr);
2317113c68e6SMatthew G. Knepley   ierr = DMPlexSetMinRadius(dm, gminradius);CHKERRQ(ierr);
2318113c68e6SMatthew G. Knepley   /* Compute centroids of ghost cells */
2319113c68e6SMatthew G. Knepley   for (c = cEndInterior; c < cEnd; ++c) {
2320113c68e6SMatthew G. Knepley     PetscFVFaceGeom *fg;
2321113c68e6SMatthew G. Knepley     const PetscInt  *cone,    *support;
2322113c68e6SMatthew G. Knepley     PetscInt         coneSize, supportSize, s;
2323113c68e6SMatthew G. Knepley 
2324113c68e6SMatthew G. Knepley     ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr);
2325113c68e6SMatthew G. Knepley     if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize);
2326113c68e6SMatthew G. Knepley     ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr);
2327113c68e6SMatthew G. Knepley     ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr);
232850d63984SToby Isaac     if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 2", cone[0], supportSize);
2329113c68e6SMatthew G. Knepley     ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr);
2330113c68e6SMatthew G. Knepley     ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr);
2331113c68e6SMatthew G. Knepley     for (s = 0; s < 2; ++s) {
2332113c68e6SMatthew G. Knepley       /* Reflect ghost centroid across plane of face */
2333113c68e6SMatthew G. Knepley       if (support[s] == c) {
2334640bce14SSatish Balay         PetscFVCellGeom       *ci;
2335113c68e6SMatthew G. Knepley         PetscFVCellGeom       *cg;
2336113c68e6SMatthew G. Knepley         PetscReal              c2f[3], a;
2337113c68e6SMatthew G. Knepley 
2338113c68e6SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr);
2339113c68e6SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */
2340113c68e6SMatthew G. Knepley         a    = DMPlex_DotRealD_Internal(dim, c2f, fg->normal)/DMPlex_DotRealD_Internal(dim, fg->normal, fg->normal);
2341113c68e6SMatthew G. Knepley         ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr);
2342113c68e6SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, 2*a, fg->normal, ci->centroid, cg->centroid);
2343113c68e6SMatthew G. Knepley         cg->volume = ci->volume;
2344113c68e6SMatthew G. Knepley       }
2345113c68e6SMatthew G. Knepley     }
2346113c68e6SMatthew G. Knepley   }
2347113c68e6SMatthew G. Knepley   ierr = VecRestoreArray(*facegeom, &fgeom);CHKERRQ(ierr);
2348113c68e6SMatthew G. Knepley   ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr);
2349113c68e6SMatthew G. Knepley   ierr = DMDestroy(&dmCell);CHKERRQ(ierr);
2350113c68e6SMatthew G. Knepley   ierr = DMDestroy(&dmFace);CHKERRQ(ierr);
2351113c68e6SMatthew G. Knepley   PetscFunctionReturn(0);
2352113c68e6SMatthew G. Knepley }
2353113c68e6SMatthew G. Knepley 
2354113c68e6SMatthew G. Knepley /*@C
2355113c68e6SMatthew G. Knepley   DMPlexGetMinRadius - Returns the minimum distance from any cell centroid to a face
2356113c68e6SMatthew G. Knepley 
2357113c68e6SMatthew G. Knepley   Not collective
2358113c68e6SMatthew G. Knepley 
23594165533cSJose E. Roman   Input Parameter:
2360113c68e6SMatthew G. Knepley . dm - the DM
2361113c68e6SMatthew G. Knepley 
23624165533cSJose E. Roman   Output Parameter:
2363a5b23f4aSJose E. Roman . minradius - the minimum cell radius
2364113c68e6SMatthew G. Knepley 
2365113c68e6SMatthew G. Knepley   Level: developer
2366113c68e6SMatthew G. Knepley 
2367113c68e6SMatthew G. Knepley .seealso: DMGetCoordinates()
2368113c68e6SMatthew G. Knepley @*/
2369113c68e6SMatthew G. Knepley PetscErrorCode DMPlexGetMinRadius(DM dm, PetscReal *minradius)
2370113c68e6SMatthew G. Knepley {
2371113c68e6SMatthew G. Knepley   PetscFunctionBegin;
2372113c68e6SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2373113c68e6SMatthew G. Knepley   PetscValidPointer(minradius,2);
2374113c68e6SMatthew G. Knepley   *minradius = ((DM_Plex*) dm->data)->minradius;
2375113c68e6SMatthew G. Knepley   PetscFunctionReturn(0);
2376113c68e6SMatthew G. Knepley }
2377113c68e6SMatthew G. Knepley 
2378113c68e6SMatthew G. Knepley /*@C
2379113c68e6SMatthew G. Knepley   DMPlexSetMinRadius - Sets the minimum distance from the cell centroid to a face
2380113c68e6SMatthew G. Knepley 
2381113c68e6SMatthew G. Knepley   Logically collective
2382113c68e6SMatthew G. Knepley 
23834165533cSJose E. Roman   Input Parameters:
2384113c68e6SMatthew G. Knepley + dm - the DM
2385a5b23f4aSJose E. Roman - minradius - the minimum cell radius
2386113c68e6SMatthew G. Knepley 
2387113c68e6SMatthew G. Knepley   Level: developer
2388113c68e6SMatthew G. Knepley 
2389113c68e6SMatthew G. Knepley .seealso: DMSetCoordinates()
2390113c68e6SMatthew G. Knepley @*/
2391113c68e6SMatthew G. Knepley PetscErrorCode DMPlexSetMinRadius(DM dm, PetscReal minradius)
2392113c68e6SMatthew G. Knepley {
2393113c68e6SMatthew G. Knepley   PetscFunctionBegin;
2394113c68e6SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2395113c68e6SMatthew G. Knepley   ((DM_Plex*) dm->data)->minradius = minradius;
2396113c68e6SMatthew G. Knepley   PetscFunctionReturn(0);
2397113c68e6SMatthew G. Knepley }
2398856ac710SMatthew G. Knepley 
2399856ac710SMatthew G. Knepley static PetscErrorCode BuildGradientReconstruction_Internal(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
2400856ac710SMatthew G. Knepley {
2401856ac710SMatthew G. Knepley   DMLabel        ghostLabel;
2402856ac710SMatthew G. Knepley   PetscScalar   *dx, *grad, **gref;
2403856ac710SMatthew G. Knepley   PetscInt       dim, cStart, cEnd, c, cEndInterior, maxNumFaces;
2404856ac710SMatthew G. Knepley   PetscErrorCode ierr;
2405856ac710SMatthew G. Knepley 
2406856ac710SMatthew G. Knepley   PetscFunctionBegin;
2407856ac710SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2408856ac710SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2409485ad865SMatthew G. Knepley   ierr = DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);CHKERRQ(ierr);
2410089217ebSMatthew G. Knepley   cEndInterior = cEndInterior < 0 ? cEnd : cEndInterior;
2411856ac710SMatthew G. Knepley   ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr);
2412856ac710SMatthew G. Knepley   ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr);
2413c58f1c22SToby Isaac   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
2414856ac710SMatthew G. Knepley   ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr);
2415856ac710SMatthew G. Knepley   for (c = cStart; c < cEndInterior; c++) {
2416856ac710SMatthew G. Knepley     const PetscInt        *faces;
2417856ac710SMatthew G. Knepley     PetscInt               numFaces, usedFaces, f, d;
2418640bce14SSatish Balay     PetscFVCellGeom        *cg;
2419856ac710SMatthew G. Knepley     PetscBool              boundary;
2420856ac710SMatthew G. Knepley     PetscInt               ghost;
2421856ac710SMatthew G. Knepley 
2422856ac710SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
2423856ac710SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr);
2424856ac710SMatthew G. Knepley     ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr);
2425856ac710SMatthew 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);
2426856ac710SMatthew G. Knepley     for (f = 0, usedFaces = 0; f < numFaces; ++f) {
2427640bce14SSatish Balay       PetscFVCellGeom       *cg1;
2428856ac710SMatthew G. Knepley       PetscFVFaceGeom       *fg;
2429856ac710SMatthew G. Knepley       const PetscInt        *fcells;
2430856ac710SMatthew G. Knepley       PetscInt               ncell, side;
2431856ac710SMatthew G. Knepley 
2432856ac710SMatthew G. Knepley       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
2433a6ba4734SToby Isaac       ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
2434856ac710SMatthew G. Knepley       if ((ghost >= 0) || boundary) continue;
2435856ac710SMatthew G. Knepley       ierr  = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr);
2436856ac710SMatthew G. Knepley       side  = (c != fcells[0]); /* c is on left=0 or right=1 of face */
2437856ac710SMatthew G. Knepley       ncell = fcells[!side];    /* the neighbor */
2438856ac710SMatthew G. Knepley       ierr  = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr);
2439856ac710SMatthew G. Knepley       ierr  = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr);
2440856ac710SMatthew G. Knepley       for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d];
2441856ac710SMatthew G. Knepley       gref[usedFaces++] = fg->grad[side];  /* Gradient reconstruction term will go here */
2442856ac710SMatthew G. Knepley     }
2443856ac710SMatthew G. Knepley     if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?");
2444856ac710SMatthew G. Knepley     ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr);
2445856ac710SMatthew G. Knepley     for (f = 0, usedFaces = 0; f < numFaces; ++f) {
2446856ac710SMatthew G. Knepley       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
2447a6ba4734SToby Isaac       ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
2448856ac710SMatthew G. Knepley       if ((ghost >= 0) || boundary) continue;
2449856ac710SMatthew G. Knepley       for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d];
2450856ac710SMatthew G. Knepley       ++usedFaces;
2451856ac710SMatthew G. Knepley     }
2452856ac710SMatthew G. Knepley   }
2453856ac710SMatthew G. Knepley   ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr);
2454856ac710SMatthew G. Knepley   PetscFunctionReturn(0);
2455856ac710SMatthew G. Knepley }
2456856ac710SMatthew G. Knepley 
2457b81db932SToby Isaac static PetscErrorCode BuildGradientReconstruction_Internal_Tree(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
2458b81db932SToby Isaac {
2459b81db932SToby Isaac   DMLabel        ghostLabel;
2460b81db932SToby Isaac   PetscScalar   *dx, *grad, **gref;
2461b81db932SToby Isaac   PetscInt       dim, cStart, cEnd, c, cEndInterior, fStart, fEnd, f, nStart, nEnd, maxNumFaces = 0;
2462b81db932SToby Isaac   PetscSection   neighSec;
2463b81db932SToby Isaac   PetscInt     (*neighbors)[2];
2464b81db932SToby Isaac   PetscInt      *counter;
2465b81db932SToby Isaac   PetscErrorCode ierr;
2466b81db932SToby Isaac 
2467b81db932SToby Isaac   PetscFunctionBegin;
2468b81db932SToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2469b81db932SToby Isaac   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2470485ad865SMatthew G. Knepley   ierr = DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);CHKERRQ(ierr);
2471485ad865SMatthew G. Knepley   if (cEndInterior < 0) cEndInterior = cEnd;
2472b81db932SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&neighSec);CHKERRQ(ierr);
2473b81db932SToby Isaac   ierr = PetscSectionSetChart(neighSec,cStart,cEndInterior);CHKERRQ(ierr);
2474b81db932SToby Isaac   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
2475c58f1c22SToby Isaac   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
2476b81db932SToby Isaac   for (f = fStart; f < fEnd; f++) {
2477b81db932SToby Isaac     const PetscInt        *fcells;
2478b81db932SToby Isaac     PetscBool              boundary;
24795bc680faSToby Isaac     PetscInt               ghost = -1;
2480b81db932SToby Isaac     PetscInt               numChildren, numCells, c;
2481b81db932SToby Isaac 
248206348e87SToby Isaac     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);}
2483a6ba4734SToby Isaac     ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr);
2484b81db932SToby Isaac     ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr);
2485b81db932SToby Isaac     if ((ghost >= 0) || boundary || numChildren) continue;
2486b81db932SToby Isaac     ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr);
248706348e87SToby Isaac     if (numCells == 2) {
2488b81db932SToby Isaac       ierr = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr);
2489b81db932SToby Isaac       for (c = 0; c < 2; c++) {
2490b81db932SToby Isaac         PetscInt cell = fcells[c];
2491b81db932SToby Isaac 
2492e6885bbbSToby Isaac         if (cell >= cStart && cell < cEndInterior) {
2493b81db932SToby Isaac           ierr = PetscSectionAddDof(neighSec,cell,1);CHKERRQ(ierr);
2494b81db932SToby Isaac         }
2495b81db932SToby Isaac       }
2496b81db932SToby Isaac     }
249706348e87SToby Isaac   }
2498b81db932SToby Isaac   ierr = PetscSectionSetUp(neighSec);CHKERRQ(ierr);
2499b81db932SToby Isaac   ierr = PetscSectionGetMaxDof(neighSec,&maxNumFaces);CHKERRQ(ierr);
2500b81db932SToby Isaac   ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr);
2501b81db932SToby Isaac   nStart = 0;
2502b81db932SToby Isaac   ierr = PetscSectionGetStorageSize(neighSec,&nEnd);CHKERRQ(ierr);
2503b81db932SToby Isaac   ierr = PetscMalloc1((nEnd-nStart),&neighbors);CHKERRQ(ierr);
2504b81db932SToby Isaac   ierr = PetscCalloc1((cEndInterior-cStart),&counter);CHKERRQ(ierr);
2505b81db932SToby Isaac   for (f = fStart; f < fEnd; f++) {
2506b81db932SToby Isaac     const PetscInt        *fcells;
2507b81db932SToby Isaac     PetscBool              boundary;
25085bc680faSToby Isaac     PetscInt               ghost = -1;
2509b81db932SToby Isaac     PetscInt               numChildren, numCells, c;
2510b81db932SToby Isaac 
251106348e87SToby Isaac     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);}
2512a6ba4734SToby Isaac     ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr);
2513b81db932SToby Isaac     ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr);
2514b81db932SToby Isaac     if ((ghost >= 0) || boundary || numChildren) continue;
2515b81db932SToby Isaac     ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr);
251606348e87SToby Isaac     if (numCells == 2) {
2517b81db932SToby Isaac       ierr  = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr);
2518b81db932SToby Isaac       for (c = 0; c < 2; c++) {
2519b81db932SToby Isaac         PetscInt cell = fcells[c], off;
2520b81db932SToby Isaac 
2521e6885bbbSToby Isaac         if (cell >= cStart && cell < cEndInterior) {
2522b81db932SToby Isaac           ierr = PetscSectionGetOffset(neighSec,cell,&off);CHKERRQ(ierr);
2523b81db932SToby Isaac           off += counter[cell - cStart]++;
2524b81db932SToby Isaac           neighbors[off][0] = f;
2525b81db932SToby Isaac           neighbors[off][1] = fcells[1 - c];
2526b81db932SToby Isaac         }
2527b81db932SToby Isaac       }
2528b81db932SToby Isaac     }
252906348e87SToby Isaac   }
2530b81db932SToby Isaac   ierr = PetscFree(counter);CHKERRQ(ierr);
2531b81db932SToby Isaac   ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr);
2532b81db932SToby Isaac   for (c = cStart; c < cEndInterior; c++) {
2533317218b9SToby Isaac     PetscInt               numFaces, f, d, off, ghost = -1;
2534640bce14SSatish Balay     PetscFVCellGeom        *cg;
2535b81db932SToby Isaac 
2536b81db932SToby Isaac     ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
2537b81db932SToby Isaac     ierr = PetscSectionGetDof(neighSec, c, &numFaces);CHKERRQ(ierr);
2538b81db932SToby Isaac     ierr = PetscSectionGetOffset(neighSec, c, &off);CHKERRQ(ierr);
2539317218b9SToby Isaac     if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, c, &ghost);CHKERRQ(ierr);}
2540317218b9SToby 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);
2541b81db932SToby Isaac     for (f = 0; f < numFaces; ++f) {
2542640bce14SSatish Balay       PetscFVCellGeom       *cg1;
2543b81db932SToby Isaac       PetscFVFaceGeom       *fg;
2544b81db932SToby Isaac       const PetscInt        *fcells;
2545b81db932SToby Isaac       PetscInt               ncell, side, nface;
2546b81db932SToby Isaac 
2547b81db932SToby Isaac       nface = neighbors[off + f][0];
2548b81db932SToby Isaac       ncell = neighbors[off + f][1];
2549b81db932SToby Isaac       ierr  = DMPlexGetSupport(dm,nface,&fcells);CHKERRQ(ierr);
2550b81db932SToby Isaac       side  = (c != fcells[0]);
2551b81db932SToby Isaac       ierr  = DMPlexPointLocalRef(dmFace, nface, fgeom, &fg);CHKERRQ(ierr);
2552b81db932SToby Isaac       ierr  = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr);
2553b81db932SToby Isaac       for (d = 0; d < dim; ++d) dx[f*dim+d] = cg1->centroid[d] - cg->centroid[d];
2554b81db932SToby Isaac       gref[f] = fg->grad[side];  /* Gradient reconstruction term will go here */
2555b81db932SToby Isaac     }
2556b81db932SToby Isaac     ierr = PetscFVComputeGradient(fvm, numFaces, dx, grad);CHKERRQ(ierr);
2557b81db932SToby Isaac     for (f = 0; f < numFaces; ++f) {
2558b81db932SToby Isaac       for (d = 0; d < dim; ++d) gref[f][d] = grad[f*dim+d];
2559b81db932SToby Isaac     }
2560b81db932SToby Isaac   }
2561b81db932SToby Isaac   ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr);
25625fe94518SToby Isaac   ierr = PetscSectionDestroy(&neighSec);CHKERRQ(ierr);
2563b81db932SToby Isaac   ierr = PetscFree(neighbors);CHKERRQ(ierr);
2564b81db932SToby Isaac   PetscFunctionReturn(0);
2565b81db932SToby Isaac }
2566b81db932SToby Isaac 
2567856ac710SMatthew G. Knepley /*@
2568856ac710SMatthew G. Knepley   DMPlexComputeGradientFVM - Compute geometric factors for gradient reconstruction, which are stored in the geometry data, and compute layout for gradient data
2569856ac710SMatthew G. Knepley 
2570d083f849SBarry Smith   Collective on dm
2571856ac710SMatthew G. Knepley 
25724165533cSJose E. Roman   Input Parameters:
2573856ac710SMatthew G. Knepley + dm  - The DM
2574856ac710SMatthew G. Knepley . fvm - The PetscFV
25758f9f38e3SMatthew G. Knepley - cellGeometry - The face geometry from DMPlexComputeCellGeometryFVM()
2576856ac710SMatthew G. Knepley 
25776b867d5aSJose E. Roman   Input/Output Parameter:
25786b867d5aSJose E. Roman . faceGeometry - The face geometry from DMPlexComputeFaceGeometryFVM(); on output
25796b867d5aSJose E. Roman                  the geometric factors for gradient calculation are inserted
25806b867d5aSJose E. Roman 
25816b867d5aSJose E. Roman   Output Parameter:
25826b867d5aSJose E. Roman . dmGrad - The DM describing the layout of gradient data
2583856ac710SMatthew G. Knepley 
2584856ac710SMatthew G. Knepley   Level: developer
2585856ac710SMatthew G. Knepley 
2586856ac710SMatthew G. Knepley .seealso: DMPlexGetFaceGeometryFVM(), DMPlexGetCellGeometryFVM()
2587856ac710SMatthew G. Knepley @*/
2588856ac710SMatthew G. Knepley PetscErrorCode DMPlexComputeGradientFVM(DM dm, PetscFV fvm, Vec faceGeometry, Vec cellGeometry, DM *dmGrad)
2589856ac710SMatthew G. Knepley {
2590856ac710SMatthew G. Knepley   DM             dmFace, dmCell;
2591856ac710SMatthew G. Knepley   PetscScalar   *fgeom, *cgeom;
2592b81db932SToby Isaac   PetscSection   sectionGrad, parentSection;
2593856ac710SMatthew G. Knepley   PetscInt       dim, pdim, cStart, cEnd, cEndInterior, c;
2594856ac710SMatthew G. Knepley   PetscErrorCode ierr;
2595856ac710SMatthew G. Knepley 
2596856ac710SMatthew G. Knepley   PetscFunctionBegin;
2597856ac710SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2598856ac710SMatthew G. Knepley   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
2599856ac710SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2600485ad865SMatthew G. Knepley   ierr = DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);CHKERRQ(ierr);
2601856ac710SMatthew G. Knepley   /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */
2602856ac710SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
2603856ac710SMatthew G. Knepley   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
2604856ac710SMatthew G. Knepley   ierr = VecGetArray(faceGeometry, &fgeom);CHKERRQ(ierr);
2605856ac710SMatthew G. Knepley   ierr = VecGetArray(cellGeometry, &cgeom);CHKERRQ(ierr);
2606b81db932SToby Isaac   ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
2607b81db932SToby Isaac   if (!parentSection) {
2608856ac710SMatthew G. Knepley     ierr = BuildGradientReconstruction_Internal(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr);
2609b5a3613cSMatthew G. Knepley   } else {
2610b81db932SToby Isaac     ierr = BuildGradientReconstruction_Internal_Tree(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr);
2611b81db932SToby Isaac   }
2612856ac710SMatthew G. Knepley   ierr = VecRestoreArray(faceGeometry, &fgeom);CHKERRQ(ierr);
2613856ac710SMatthew G. Knepley   ierr = VecRestoreArray(cellGeometry, &cgeom);CHKERRQ(ierr);
2614856ac710SMatthew G. Knepley   /* Create storage for gradients */
2615856ac710SMatthew G. Knepley   ierr = DMClone(dm, dmGrad);CHKERRQ(ierr);
2616856ac710SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionGrad);CHKERRQ(ierr);
2617856ac710SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr);
2618856ac710SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);}
2619856ac710SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr);
262092fd8e1eSJed Brown   ierr = DMSetLocalSection(*dmGrad, sectionGrad);CHKERRQ(ierr);
2621856ac710SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionGrad);CHKERRQ(ierr);
2622856ac710SMatthew G. Knepley   PetscFunctionReturn(0);
2623856ac710SMatthew G. Knepley }
2624b27d5b9eSToby Isaac 
2625c501906fSMatthew G. Knepley /*@
2626c501906fSMatthew G. Knepley   DMPlexGetDataFVM - Retrieve precomputed cell geometry
2627c501906fSMatthew G. Knepley 
2628d083f849SBarry Smith   Collective on dm
2629c501906fSMatthew G. Knepley 
26304165533cSJose E. Roman   Input Parameters:
2631c501906fSMatthew G. Knepley + dm  - The DM
26326b867d5aSJose E. Roman - fv  - The PetscFV
2633c501906fSMatthew G. Knepley 
2634c501906fSMatthew G. Knepley   Output Parameters:
2635c501906fSMatthew G. Knepley + cellGeometry - The cell geometry
2636c501906fSMatthew G. Knepley . faceGeometry - The face geometry
26376b867d5aSJose E. Roman - gradDM       - The gradient matrices
2638c501906fSMatthew G. Knepley 
2639c501906fSMatthew G. Knepley   Level: developer
2640c501906fSMatthew G. Knepley 
2641c501906fSMatthew G. Knepley .seealso: DMPlexComputeGeometryFVM()
2642c501906fSMatthew G. Knepley @*/
2643b27d5b9eSToby Isaac PetscErrorCode DMPlexGetDataFVM(DM dm, PetscFV fv, Vec *cellgeom, Vec *facegeom, DM *gradDM)
2644b27d5b9eSToby Isaac {
2645b27d5b9eSToby Isaac   PetscObject    cellgeomobj, facegeomobj;
2646b27d5b9eSToby Isaac   PetscErrorCode ierr;
2647b27d5b9eSToby Isaac 
2648b27d5b9eSToby Isaac   PetscFunctionBegin;
2649b27d5b9eSToby Isaac   ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);CHKERRQ(ierr);
2650b27d5b9eSToby Isaac   if (!cellgeomobj) {
2651b27d5b9eSToby Isaac     Vec cellgeomInt, facegeomInt;
2652b27d5b9eSToby Isaac 
2653b27d5b9eSToby Isaac     ierr = DMPlexComputeGeometryFVM(dm, &cellgeomInt, &facegeomInt);CHKERRQ(ierr);
2654b27d5b9eSToby Isaac     ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_cellgeom_fvm",(PetscObject)cellgeomInt);CHKERRQ(ierr);
2655b27d5b9eSToby Isaac     ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_facegeom_fvm",(PetscObject)facegeomInt);CHKERRQ(ierr);
2656b27d5b9eSToby Isaac     ierr = VecDestroy(&cellgeomInt);CHKERRQ(ierr);
2657b27d5b9eSToby Isaac     ierr = VecDestroy(&facegeomInt);CHKERRQ(ierr);
2658b27d5b9eSToby Isaac     ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);CHKERRQ(ierr);
2659b27d5b9eSToby Isaac   }
2660b27d5b9eSToby Isaac   ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_facegeom_fvm", &facegeomobj);CHKERRQ(ierr);
2661b27d5b9eSToby Isaac   if (cellgeom) *cellgeom = (Vec) cellgeomobj;
2662b27d5b9eSToby Isaac   if (facegeom) *facegeom = (Vec) facegeomobj;
2663b27d5b9eSToby Isaac   if (gradDM) {
2664b27d5b9eSToby Isaac     PetscObject gradobj;
2665b27d5b9eSToby Isaac     PetscBool   computeGradients;
2666b27d5b9eSToby Isaac 
2667b27d5b9eSToby Isaac     ierr = PetscFVGetComputeGradients(fv,&computeGradients);CHKERRQ(ierr);
2668b27d5b9eSToby Isaac     if (!computeGradients) {
2669b27d5b9eSToby Isaac       *gradDM = NULL;
2670b27d5b9eSToby Isaac       PetscFunctionReturn(0);
2671b27d5b9eSToby Isaac     }
2672b27d5b9eSToby Isaac     ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);CHKERRQ(ierr);
2673b27d5b9eSToby Isaac     if (!gradobj) {
2674b27d5b9eSToby Isaac       DM dmGradInt;
2675b27d5b9eSToby Isaac 
2676b27d5b9eSToby Isaac       ierr = DMPlexComputeGradientFVM(dm,fv,(Vec) facegeomobj,(Vec) cellgeomobj,&dmGradInt);CHKERRQ(ierr);
2677b27d5b9eSToby Isaac       ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_dmgrad_fvm", (PetscObject)dmGradInt);CHKERRQ(ierr);
2678b27d5b9eSToby Isaac       ierr = DMDestroy(&dmGradInt);CHKERRQ(ierr);
2679b27d5b9eSToby Isaac       ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);CHKERRQ(ierr);
2680b27d5b9eSToby Isaac     }
2681b27d5b9eSToby Isaac     *gradDM = (DM) gradobj;
2682b27d5b9eSToby Isaac   }
2683b27d5b9eSToby Isaac   PetscFunctionReturn(0);
2684b27d5b9eSToby Isaac }
2685d6143a4eSToby Isaac 
26869d150b73SToby Isaac static PetscErrorCode DMPlexCoordinatesToReference_NewtonUpdate(PetscInt dimC, PetscInt dimR, PetscScalar *J, PetscScalar *invJ, PetscScalar *work,  PetscReal *resNeg, PetscReal *guess)
26879d150b73SToby Isaac {
26889d150b73SToby Isaac   PetscInt l, m;
26899d150b73SToby Isaac 
2690cd345991SToby Isaac   PetscFunctionBeginHot;
26919d150b73SToby Isaac   if (dimC == dimR && dimR <= 3) {
26929d150b73SToby Isaac     /* invert Jacobian, multiply */
26939d150b73SToby Isaac     PetscScalar det, idet;
26949d150b73SToby Isaac 
26959d150b73SToby Isaac     switch (dimR) {
26969d150b73SToby Isaac     case 1:
26979d150b73SToby Isaac       invJ[0] = 1./ J[0];
26989d150b73SToby Isaac       break;
26999d150b73SToby Isaac     case 2:
27009d150b73SToby Isaac       det = J[0] * J[3] - J[1] * J[2];
27019d150b73SToby Isaac       idet = 1./det;
27029d150b73SToby Isaac       invJ[0] =  J[3] * idet;
27039d150b73SToby Isaac       invJ[1] = -J[1] * idet;
27049d150b73SToby Isaac       invJ[2] = -J[2] * idet;
27059d150b73SToby Isaac       invJ[3] =  J[0] * idet;
27069d150b73SToby Isaac       break;
27079d150b73SToby Isaac     case 3:
27089d150b73SToby Isaac       {
27099d150b73SToby Isaac         invJ[0] = J[4] * J[8] - J[5] * J[7];
27109d150b73SToby Isaac         invJ[1] = J[2] * J[7] - J[1] * J[8];
27119d150b73SToby Isaac         invJ[2] = J[1] * J[5] - J[2] * J[4];
27129d150b73SToby Isaac         det = invJ[0] * J[0] + invJ[1] * J[3] + invJ[2] * J[6];
27139d150b73SToby Isaac         idet = 1./det;
27149d150b73SToby Isaac         invJ[0] *= idet;
27159d150b73SToby Isaac         invJ[1] *= idet;
27169d150b73SToby Isaac         invJ[2] *= idet;
27179d150b73SToby Isaac         invJ[3]  = idet * (J[5] * J[6] - J[3] * J[8]);
27189d150b73SToby Isaac         invJ[4]  = idet * (J[0] * J[8] - J[2] * J[6]);
27199d150b73SToby Isaac         invJ[5]  = idet * (J[2] * J[3] - J[0] * J[5]);
27209d150b73SToby Isaac         invJ[6]  = idet * (J[3] * J[7] - J[4] * J[6]);
27219d150b73SToby Isaac         invJ[7]  = idet * (J[1] * J[6] - J[0] * J[7]);
27229d150b73SToby Isaac         invJ[8]  = idet * (J[0] * J[4] - J[1] * J[3]);
27239d150b73SToby Isaac       }
27249d150b73SToby Isaac       break;
27259d150b73SToby Isaac     }
27269d150b73SToby Isaac     for (l = 0; l < dimR; l++) {
27279d150b73SToby Isaac       for (m = 0; m < dimC; m++) {
2728c6e120d1SToby Isaac         guess[l] += PetscRealPart(invJ[l * dimC + m]) * resNeg[m];
27299d150b73SToby Isaac       }
27309d150b73SToby Isaac     }
27319d150b73SToby Isaac   } else {
27329d150b73SToby Isaac #if defined(PETSC_USE_COMPLEX)
27339d150b73SToby Isaac     char transpose = 'C';
27349d150b73SToby Isaac #else
27359d150b73SToby Isaac     char transpose = 'T';
27369d150b73SToby Isaac #endif
27379d150b73SToby Isaac     PetscBLASInt m = dimR;
27389d150b73SToby Isaac     PetscBLASInt n = dimC;
27399d150b73SToby Isaac     PetscBLASInt one = 1;
27409d150b73SToby Isaac     PetscBLASInt worksize = dimR * dimC, info;
27419d150b73SToby Isaac 
27429d150b73SToby Isaac     for (l = 0; l < dimC; l++) {invJ[l] = resNeg[l];}
27439d150b73SToby Isaac 
27449d150b73SToby Isaac     PetscStackCallBLAS("LAPACKgels",LAPACKgels_(&transpose,&m,&n,&one,J,&m,invJ,&n,work,&worksize, &info));
27459d150b73SToby Isaac     if (info != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELS");
27469d150b73SToby Isaac 
2747c6e120d1SToby Isaac     for (l = 0; l < dimR; l++) {guess[l] += PetscRealPart(invJ[l]);}
27489d150b73SToby Isaac   }
27499d150b73SToby Isaac   PetscFunctionReturn(0);
27509d150b73SToby Isaac }
27519d150b73SToby Isaac 
27529d150b73SToby Isaac static PetscErrorCode DMPlexCoordinatesToReference_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[], Vec coords, PetscInt dimC, PetscInt dimR)
27539d150b73SToby Isaac {
2754c0cbe899SToby Isaac   PetscInt       coordSize, i, j, k, l, m, maxIts = 7, numV = (1 << dimR);
27559d150b73SToby Isaac   PetscScalar    *coordsScalar = NULL;
27569d150b73SToby Isaac   PetscReal      *cellData, *cellCoords, *cellCoeffs, *extJ, *resNeg;
27579d150b73SToby Isaac   PetscScalar    *J, *invJ, *work;
27589d150b73SToby Isaac   PetscErrorCode ierr;
27599d150b73SToby Isaac 
27609d150b73SToby Isaac   PetscFunctionBegin;
27619d150b73SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
27629d150b73SToby Isaac   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
276313903a91SSatish Balay   if (coordSize < dimC * numV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expecting at least %D coordinates, got %D",dimC * (1 << dimR), coordSize);
276469291d52SBarry Smith   ierr = DMGetWorkArray(dm, 2 * coordSize + dimR + dimC, MPIU_REAL, &cellData);CHKERRQ(ierr);
276569291d52SBarry Smith   ierr = DMGetWorkArray(dm, 3 * dimR * dimC, MPIU_SCALAR, &J);CHKERRQ(ierr);
27669d150b73SToby Isaac   cellCoords = &cellData[0];
27679d150b73SToby Isaac   cellCoeffs = &cellData[coordSize];
27689d150b73SToby Isaac   extJ       = &cellData[2 * coordSize];
27699d150b73SToby Isaac   resNeg     = &cellData[2 * coordSize + dimR];
27709d150b73SToby Isaac   invJ       = &J[dimR * dimC];
27719d150b73SToby Isaac   work       = &J[2 * dimR * dimC];
27729d150b73SToby Isaac   if (dimR == 2) {
27739d150b73SToby Isaac     const PetscInt zToPlex[4] = {0, 1, 3, 2};
27749d150b73SToby Isaac 
27759d150b73SToby Isaac     for (i = 0; i < 4; i++) {
27769d150b73SToby Isaac       PetscInt plexI = zToPlex[i];
27779d150b73SToby Isaac 
27789d150b73SToby Isaac       for (j = 0; j < dimC; j++) {
27799d150b73SToby Isaac         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
27809d150b73SToby Isaac       }
27819d150b73SToby Isaac     }
27829d150b73SToby Isaac   } else if (dimR == 3) {
27839d150b73SToby Isaac     const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6};
27849d150b73SToby Isaac 
27859d150b73SToby Isaac     for (i = 0; i < 8; i++) {
27869d150b73SToby Isaac       PetscInt plexI = zToPlex[i];
27879d150b73SToby Isaac 
27889d150b73SToby Isaac       for (j = 0; j < dimC; j++) {
27899d150b73SToby Isaac         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
27909d150b73SToby Isaac       }
27919d150b73SToby Isaac     }
27929d150b73SToby Isaac   } else {
27939d150b73SToby Isaac     for (i = 0; i < coordSize; i++) {cellCoords[i] = PetscRealPart(coordsScalar[i]);}
27949d150b73SToby Isaac   }
27959d150b73SToby Isaac   /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */
27969d150b73SToby Isaac   for (i = 0; i < dimR; i++) {
27979d150b73SToby Isaac     PetscReal *swap;
27989d150b73SToby Isaac 
27999d150b73SToby Isaac     for (j = 0; j < (numV / 2); j++) {
28009d150b73SToby Isaac       for (k = 0; k < dimC; k++) {
28019d150b73SToby Isaac         cellCoeffs[dimC * j + k]                = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]);
28029d150b73SToby Isaac         cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]);
28039d150b73SToby Isaac       }
28049d150b73SToby Isaac     }
28059d150b73SToby Isaac 
28069d150b73SToby Isaac     if (i < dimR - 1) {
28079d150b73SToby Isaac       swap = cellCoeffs;
28089d150b73SToby Isaac       cellCoeffs = cellCoords;
28099d150b73SToby Isaac       cellCoords = swap;
28109d150b73SToby Isaac     }
28119d150b73SToby Isaac   }
2812580bdb30SBarry Smith   ierr = PetscArrayzero(refCoords,numPoints * dimR);CHKERRQ(ierr);
28139d150b73SToby Isaac   for (j = 0; j < numPoints; j++) {
28149d150b73SToby Isaac     for (i = 0; i < maxIts; i++) {
28159d150b73SToby Isaac       PetscReal *guess = &refCoords[dimR * j];
28169d150b73SToby Isaac 
28179d150b73SToby Isaac       /* compute -residual and Jacobian */
28189d150b73SToby Isaac       for (k = 0; k < dimC; k++) {resNeg[k] = realCoords[dimC * j + k];}
28199d150b73SToby Isaac       for (k = 0; k < dimC * dimR; k++) {J[k] = 0.;}
28209d150b73SToby Isaac       for (k = 0; k < numV; k++) {
28219d150b73SToby Isaac         PetscReal extCoord = 1.;
28229d150b73SToby Isaac         for (l = 0; l < dimR; l++) {
28239d150b73SToby Isaac           PetscReal coord = guess[l];
28249d150b73SToby Isaac           PetscInt  dep   = (k & (1 << l)) >> l;
28259d150b73SToby Isaac 
28269d150b73SToby Isaac           extCoord *= dep * coord + !dep;
28279d150b73SToby Isaac           extJ[l] = dep;
28289d150b73SToby Isaac 
28299d150b73SToby Isaac           for (m = 0; m < dimR; m++) {
28309d150b73SToby Isaac             PetscReal coord = guess[m];
28319d150b73SToby Isaac             PetscInt  dep   = ((k & (1 << m)) >> m) && (m != l);
28329d150b73SToby Isaac             PetscReal mult  = dep * coord + !dep;
28339d150b73SToby Isaac 
28349d150b73SToby Isaac             extJ[l] *= mult;
28359d150b73SToby Isaac           }
28369d150b73SToby Isaac         }
28379d150b73SToby Isaac         for (l = 0; l < dimC; l++) {
28389d150b73SToby Isaac           PetscReal coeff = cellCoeffs[dimC * k + l];
28399d150b73SToby Isaac 
28409d150b73SToby Isaac           resNeg[l] -= coeff * extCoord;
28419d150b73SToby Isaac           for (m = 0; m < dimR; m++) {
28429d150b73SToby Isaac             J[dimR * l + m] += coeff * extJ[m];
28439d150b73SToby Isaac           }
28449d150b73SToby Isaac         }
28459d150b73SToby Isaac       }
284676bd3646SJed Brown       if (0 && PetscDefined(USE_DEBUG)) {
28470611203eSToby Isaac         PetscReal maxAbs = 0.;
28480611203eSToby Isaac 
28490611203eSToby Isaac         for (l = 0; l < dimC; l++) {
28500611203eSToby Isaac           maxAbs = PetscMax(maxAbs,PetscAbsReal(resNeg[l]));
28510611203eSToby Isaac         }
2852087ef6b2SMatthew G. Knepley         ierr = PetscInfo4(dm,"cell %D, point %D, iter %D: res %g\n",cell,j,i,(double) maxAbs);CHKERRQ(ierr);
28530611203eSToby Isaac       }
28549d150b73SToby Isaac 
28559d150b73SToby Isaac       ierr = DMPlexCoordinatesToReference_NewtonUpdate(dimC,dimR,J,invJ,work,resNeg,guess);CHKERRQ(ierr);
28569d150b73SToby Isaac     }
28579d150b73SToby Isaac   }
285869291d52SBarry Smith   ierr = DMRestoreWorkArray(dm, 3 * dimR * dimC, MPIU_SCALAR, &J);CHKERRQ(ierr);
285969291d52SBarry Smith   ierr = DMRestoreWorkArray(dm, 2 * coordSize + dimR + dimC, MPIU_REAL, &cellData);CHKERRQ(ierr);
28609d150b73SToby Isaac   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
28619d150b73SToby Isaac   PetscFunctionReturn(0);
28629d150b73SToby Isaac }
28639d150b73SToby Isaac 
28649d150b73SToby Isaac static PetscErrorCode DMPlexReferenceToCoordinates_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[], Vec coords, PetscInt dimC, PetscInt dimR)
28659d150b73SToby Isaac {
28669d150b73SToby Isaac   PetscInt       coordSize, i, j, k, l, numV = (1 << dimR);
28679d150b73SToby Isaac   PetscScalar    *coordsScalar = NULL;
28689d150b73SToby Isaac   PetscReal      *cellData, *cellCoords, *cellCoeffs;
28699d150b73SToby Isaac   PetscErrorCode ierr;
28709d150b73SToby Isaac 
28719d150b73SToby Isaac   PetscFunctionBegin;
28729d150b73SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
28739d150b73SToby Isaac   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
287413903a91SSatish Balay   if (coordSize < dimC * numV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expecting at least %D coordinates, got %D",dimC * (1 << dimR), coordSize);
287569291d52SBarry Smith   ierr = DMGetWorkArray(dm, 2 * coordSize, MPIU_REAL, &cellData);CHKERRQ(ierr);
28769d150b73SToby Isaac   cellCoords = &cellData[0];
28779d150b73SToby Isaac   cellCoeffs = &cellData[coordSize];
28789d150b73SToby Isaac   if (dimR == 2) {
28799d150b73SToby Isaac     const PetscInt zToPlex[4] = {0, 1, 3, 2};
28809d150b73SToby Isaac 
28819d150b73SToby Isaac     for (i = 0; i < 4; i++) {
28829d150b73SToby Isaac       PetscInt plexI = zToPlex[i];
28839d150b73SToby Isaac 
28849d150b73SToby Isaac       for (j = 0; j < dimC; j++) {
28859d150b73SToby Isaac         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
28869d150b73SToby Isaac       }
28879d150b73SToby Isaac     }
28889d150b73SToby Isaac   } else if (dimR == 3) {
28899d150b73SToby Isaac     const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6};
28909d150b73SToby Isaac 
28919d150b73SToby Isaac     for (i = 0; i < 8; i++) {
28929d150b73SToby Isaac       PetscInt plexI = zToPlex[i];
28939d150b73SToby Isaac 
28949d150b73SToby Isaac       for (j = 0; j < dimC; j++) {
28959d150b73SToby Isaac         cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
28969d150b73SToby Isaac       }
28979d150b73SToby Isaac     }
28989d150b73SToby Isaac   } else {
28999d150b73SToby Isaac     for (i = 0; i < coordSize; i++) {cellCoords[i] = PetscRealPart(coordsScalar[i]);}
29009d150b73SToby Isaac   }
29019d150b73SToby Isaac   /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */
29029d150b73SToby Isaac   for (i = 0; i < dimR; i++) {
29039d150b73SToby Isaac     PetscReal *swap;
29049d150b73SToby Isaac 
29059d150b73SToby Isaac     for (j = 0; j < (numV / 2); j++) {
29069d150b73SToby Isaac       for (k = 0; k < dimC; k++) {
29079d150b73SToby Isaac         cellCoeffs[dimC * j + k]                = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]);
29089d150b73SToby Isaac         cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]);
29099d150b73SToby Isaac       }
29109d150b73SToby Isaac     }
29119d150b73SToby Isaac 
29129d150b73SToby Isaac     if (i < dimR - 1) {
29139d150b73SToby Isaac       swap = cellCoeffs;
29149d150b73SToby Isaac       cellCoeffs = cellCoords;
29159d150b73SToby Isaac       cellCoords = swap;
29169d150b73SToby Isaac     }
29179d150b73SToby Isaac   }
2918580bdb30SBarry Smith   ierr = PetscArrayzero(realCoords,numPoints * dimC);CHKERRQ(ierr);
29199d150b73SToby Isaac   for (j = 0; j < numPoints; j++) {
29209d150b73SToby Isaac     const PetscReal *guess  = &refCoords[dimR * j];
29219d150b73SToby Isaac     PetscReal       *mapped = &realCoords[dimC * j];
29229d150b73SToby Isaac 
29239d150b73SToby Isaac     for (k = 0; k < numV; k++) {
29249d150b73SToby Isaac       PetscReal extCoord = 1.;
29259d150b73SToby Isaac       for (l = 0; l < dimR; l++) {
29269d150b73SToby Isaac         PetscReal coord = guess[l];
29279d150b73SToby Isaac         PetscInt  dep   = (k & (1 << l)) >> l;
29289d150b73SToby Isaac 
29299d150b73SToby Isaac         extCoord *= dep * coord + !dep;
29309d150b73SToby Isaac       }
29319d150b73SToby Isaac       for (l = 0; l < dimC; l++) {
29329d150b73SToby Isaac         PetscReal coeff = cellCoeffs[dimC * k + l];
29339d150b73SToby Isaac 
29349d150b73SToby Isaac         mapped[l] += coeff * extCoord;
29359d150b73SToby Isaac       }
29369d150b73SToby Isaac     }
29379d150b73SToby Isaac   }
293869291d52SBarry Smith   ierr = DMRestoreWorkArray(dm, 2 * coordSize, MPIU_REAL, &cellData);CHKERRQ(ierr);
29399d150b73SToby Isaac   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr);
29409d150b73SToby Isaac   PetscFunctionReturn(0);
29419d150b73SToby Isaac }
29429d150b73SToby Isaac 
29439c3cf19fSMatthew G. Knepley /* TODO: TOBY please fix this for Nc > 1 */
29449c3cf19fSMatthew 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)
29459d150b73SToby Isaac {
29469c3cf19fSMatthew G. Knepley   PetscInt       numComp, pdim, i, j, k, l, m, maxIter = 7, coordSize;
2947c6e120d1SToby Isaac   PetscScalar    *nodes = NULL;
2948c6e120d1SToby Isaac   PetscReal      *invV, *modes;
2949c6e120d1SToby Isaac   PetscReal      *B, *D, *resNeg;
2950c6e120d1SToby Isaac   PetscScalar    *J, *invJ, *work;
29519d150b73SToby Isaac   PetscErrorCode ierr;
29529d150b73SToby Isaac 
29539d150b73SToby Isaac   PetscFunctionBegin;
29549c3cf19fSMatthew G. Knepley   ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr);
29559d150b73SToby Isaac   ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
29569c3cf19fSMatthew 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);
29579d150b73SToby Isaac   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
29589d150b73SToby Isaac   /* convert nodes to values in the stable evaluation basis */
295969291d52SBarry Smith   ierr = DMGetWorkArray(dm,pdim,MPIU_REAL,&modes);CHKERRQ(ierr);
29609d150b73SToby Isaac   invV = fe->invV;
2961012b7cc6SMatthew G. Knepley   for (i = 0; i < pdim; ++i) {
2962012b7cc6SMatthew G. Knepley     modes[i] = 0.;
2963012b7cc6SMatthew G. Knepley     for (j = 0; j < pdim; ++j) {
2964012b7cc6SMatthew G. Knepley       modes[i] += invV[i * pdim + j] * PetscRealPart(nodes[j]);
29659d150b73SToby Isaac     }
29669d150b73SToby Isaac   }
296769291d52SBarry Smith   ierr   = DMGetWorkArray(dm,pdim * Nc + pdim * Nc * dimR + Nc,MPIU_REAL,&B);CHKERRQ(ierr);
29689c3cf19fSMatthew G. Knepley   D      = &B[pdim*Nc];
29699c3cf19fSMatthew G. Knepley   resNeg = &D[pdim*Nc * dimR];
297069291d52SBarry Smith   ierr = DMGetWorkArray(dm,3 * Nc * dimR,MPIU_SCALAR,&J);CHKERRQ(ierr);
29719c3cf19fSMatthew G. Knepley   invJ = &J[Nc * dimR];
29729c3cf19fSMatthew G. Knepley   work = &invJ[Nc * dimR];
29739d150b73SToby Isaac   for (i = 0; i < numPoints * dimR; i++) {refCoords[i] = 0.;}
29749d150b73SToby Isaac   for (j = 0; j < numPoints; j++) {
29759b1f03cbSToby 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 */
29769d150b73SToby Isaac       PetscReal *guess = &refCoords[j * dimR];
29779d150b73SToby Isaac       ierr = PetscSpaceEvaluate(fe->basisSpace, 1, guess, B, D, NULL);CHKERRQ(ierr);
29789c3cf19fSMatthew G. Knepley       for (k = 0; k < Nc; k++) {resNeg[k] = realCoords[j * Nc + k];}
29799c3cf19fSMatthew G. Knepley       for (k = 0; k < Nc * dimR; k++) {J[k] = 0.;}
29809c3cf19fSMatthew G. Knepley       for (k = 0; k < pdim; k++) {
29819c3cf19fSMatthew G. Knepley         for (l = 0; l < Nc; l++) {
2982012b7cc6SMatthew G. Knepley           resNeg[l] -= modes[k] * B[k * Nc + l];
29839d150b73SToby Isaac           for (m = 0; m < dimR; m++) {
2984012b7cc6SMatthew G. Knepley             J[l * dimR + m] += modes[k] * D[(k * Nc + l) * dimR + m];
29859d150b73SToby Isaac           }
29869d150b73SToby Isaac         }
29879d150b73SToby Isaac       }
298876bd3646SJed Brown       if (0 && PetscDefined(USE_DEBUG)) {
29890611203eSToby Isaac         PetscReal maxAbs = 0.;
29900611203eSToby Isaac 
29919c3cf19fSMatthew G. Knepley         for (l = 0; l < Nc; l++) {
29920611203eSToby Isaac           maxAbs = PetscMax(maxAbs,PetscAbsReal(resNeg[l]));
29930611203eSToby Isaac         }
2994087ef6b2SMatthew G. Knepley         ierr = PetscInfo4(dm,"cell %D, point %D, iter %D: res %g\n",cell,j,i,(double) maxAbs);CHKERRQ(ierr);
29950611203eSToby Isaac       }
29969c3cf19fSMatthew G. Knepley       ierr = DMPlexCoordinatesToReference_NewtonUpdate(Nc,dimR,J,invJ,work,resNeg,guess);CHKERRQ(ierr);
29979d150b73SToby Isaac     }
29989d150b73SToby Isaac   }
299969291d52SBarry Smith   ierr = DMRestoreWorkArray(dm,3 * Nc * dimR,MPIU_SCALAR,&J);CHKERRQ(ierr);
300069291d52SBarry Smith   ierr = DMRestoreWorkArray(dm,pdim * Nc + pdim * Nc * dimR + Nc,MPIU_REAL,&B);CHKERRQ(ierr);
300169291d52SBarry Smith   ierr = DMRestoreWorkArray(dm,pdim,MPIU_REAL,&modes);CHKERRQ(ierr);
30029d150b73SToby Isaac   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
30039d150b73SToby Isaac   PetscFunctionReturn(0);
30049d150b73SToby Isaac }
30059d150b73SToby Isaac 
30069c3cf19fSMatthew G. Knepley /* TODO: TOBY please fix this for Nc > 1 */
30079c3cf19fSMatthew 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)
30089d150b73SToby Isaac {
30099c3cf19fSMatthew G. Knepley   PetscInt       numComp, pdim, i, j, k, l, coordSize;
3010c6e120d1SToby Isaac   PetscScalar    *nodes = NULL;
3011c6e120d1SToby Isaac   PetscReal      *invV, *modes;
30129d150b73SToby Isaac   PetscReal      *B;
30139d150b73SToby Isaac   PetscErrorCode ierr;
30149d150b73SToby Isaac 
30159d150b73SToby Isaac   PetscFunctionBegin;
30169c3cf19fSMatthew G. Knepley   ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr);
30179d150b73SToby Isaac   ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
30189c3cf19fSMatthew 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);
30199d150b73SToby Isaac   ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
30209d150b73SToby Isaac   /* convert nodes to values in the stable evaluation basis */
302169291d52SBarry Smith   ierr = DMGetWorkArray(dm,pdim,MPIU_REAL,&modes);CHKERRQ(ierr);
30229d150b73SToby Isaac   invV = fe->invV;
3023012b7cc6SMatthew G. Knepley   for (i = 0; i < pdim; ++i) {
3024012b7cc6SMatthew G. Knepley     modes[i] = 0.;
3025012b7cc6SMatthew G. Knepley     for (j = 0; j < pdim; ++j) {
3026012b7cc6SMatthew G. Knepley       modes[i] += invV[i * pdim + j] * PetscRealPart(nodes[j]);
30279d150b73SToby Isaac     }
30289d150b73SToby Isaac   }
302969291d52SBarry Smith   ierr = DMGetWorkArray(dm,numPoints * pdim * Nc,MPIU_REAL,&B);CHKERRQ(ierr);
3030012b7cc6SMatthew G. Knepley   ierr = PetscSpaceEvaluate(fe->basisSpace, numPoints, refCoords, B, NULL, NULL);CHKERRQ(ierr);
30319c3cf19fSMatthew G. Knepley   for (i = 0; i < numPoints * Nc; i++) {realCoords[i] = 0.;}
30329d150b73SToby Isaac   for (j = 0; j < numPoints; j++) {
30339c3cf19fSMatthew G. Knepley     PetscReal *mapped = &realCoords[j * Nc];
30349d150b73SToby Isaac 
30359c3cf19fSMatthew G. Knepley     for (k = 0; k < pdim; k++) {
30369c3cf19fSMatthew G. Knepley       for (l = 0; l < Nc; l++) {
303740cf36b3SToby Isaac         mapped[l] += modes[k] * B[(j * pdim + k) * Nc + l];
30389d150b73SToby Isaac       }
30399d150b73SToby Isaac     }
30409d150b73SToby Isaac   }
304169291d52SBarry Smith   ierr = DMRestoreWorkArray(dm,numPoints * pdim * Nc,MPIU_REAL,&B);CHKERRQ(ierr);
304269291d52SBarry Smith   ierr = DMRestoreWorkArray(dm,pdim,MPIU_REAL,&modes);CHKERRQ(ierr);
30439d150b73SToby Isaac   ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr);
30449d150b73SToby Isaac   PetscFunctionReturn(0);
30459d150b73SToby Isaac }
30469d150b73SToby Isaac 
3047d6143a4eSToby Isaac /*@
3048d6143a4eSToby Isaac   DMPlexCoordinatesToReference - Pull coordinates back from the mesh to the reference element using a single element
3049d6143a4eSToby Isaac   map.  This inversion will be accurate inside the reference element, but may be inaccurate for mappings that do not
3050d6143a4eSToby Isaac   extend uniquely outside the reference cell (e.g, most non-affine maps)
3051d6143a4eSToby Isaac 
3052d6143a4eSToby Isaac   Not collective
3053d6143a4eSToby Isaac 
3054d6143a4eSToby Isaac   Input Parameters:
3055d6143a4eSToby Isaac + dm         - The mesh, with coordinate maps defined either by a PetscDS for the coordinate DM (see DMGetCoordinateDM()) or
3056d6143a4eSToby Isaac                implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or
3057d6143a4eSToby Isaac                as a multilinear map for tensor-product elements
3058d6143a4eSToby Isaac . cell       - the cell whose map is used.
3059d6143a4eSToby Isaac . numPoints  - the number of points to locate
30601b266c99SBarry Smith - realCoords - (numPoints x coordinate dimension) array of coordinates (see DMGetCoordinateDim())
3061d6143a4eSToby Isaac 
3062d6143a4eSToby Isaac   Output Parameters:
3063d6143a4eSToby Isaac . refCoords  - (numPoints x dimension) array of reference coordinates (see DMGetDimension())
30641b266c99SBarry Smith 
30651b266c99SBarry Smith   Level: intermediate
306673c9229bSMatthew Knepley 
306773c9229bSMatthew Knepley .seealso: DMPlexReferenceToCoordinates()
3068d6143a4eSToby Isaac @*/
3069d6143a4eSToby Isaac PetscErrorCode DMPlexCoordinatesToReference(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[])
3070d6143a4eSToby Isaac {
3071485ad865SMatthew G. Knepley   PetscInt       dimC, dimR, depth, cStart, cEnd, i;
30729d150b73SToby Isaac   DM             coordDM = NULL;
30739d150b73SToby Isaac   Vec            coords;
30749d150b73SToby Isaac   PetscFE        fe = NULL;
30759d150b73SToby Isaac   PetscErrorCode ierr;
30769d150b73SToby Isaac 
3077d6143a4eSToby Isaac   PetscFunctionBegin;
30789d150b73SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
30799d150b73SToby Isaac   ierr = DMGetDimension(dm,&dimR);CHKERRQ(ierr);
30809d150b73SToby Isaac   ierr = DMGetCoordinateDim(dm,&dimC);CHKERRQ(ierr);
30819d150b73SToby Isaac   if (dimR <= 0 || dimC <= 0 || numPoints <= 0) PetscFunctionReturn(0);
30829d150b73SToby Isaac   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
30839d150b73SToby Isaac   ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
30849d150b73SToby Isaac   ierr = DMGetCoordinateDM(dm,&coordDM);CHKERRQ(ierr);
30859d150b73SToby Isaac   if (coordDM) {
30869d150b73SToby Isaac     PetscInt coordFields;
30879d150b73SToby Isaac 
30889d150b73SToby Isaac     ierr = DMGetNumFields(coordDM,&coordFields);CHKERRQ(ierr);
30899d150b73SToby Isaac     if (coordFields) {
30909d150b73SToby Isaac       PetscClassId id;
30919d150b73SToby Isaac       PetscObject  disc;
30929d150b73SToby Isaac 
309344a7f3ddSMatthew G. Knepley       ierr = DMGetField(coordDM,0,NULL,&disc);CHKERRQ(ierr);
30949d150b73SToby Isaac       ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
30959d150b73SToby Isaac       if (id == PETSCFE_CLASSID) {
30969d150b73SToby Isaac         fe = (PetscFE) disc;
30979d150b73SToby Isaac       }
30989d150b73SToby Isaac     }
30999d150b73SToby Isaac   }
3100412e9a14SMatthew G. Knepley   ierr = DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
310113903a91SSatish 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);
31029d150b73SToby Isaac   if (!fe) { /* implicit discretization: affine or multilinear */
31039d150b73SToby Isaac     PetscInt  coneSize;
31049d150b73SToby Isaac     PetscBool isSimplex, isTensor;
31059d150b73SToby Isaac 
31069d150b73SToby Isaac     ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr);
31079d150b73SToby Isaac     isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE;
31089d150b73SToby Isaac     isTensor  = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE;
31099d150b73SToby Isaac     if (isSimplex) {
31109d150b73SToby Isaac       PetscReal detJ, *v0, *J, *invJ;
31119d150b73SToby Isaac 
311269291d52SBarry Smith       ierr = DMGetWorkArray(dm,dimC + 2 * dimC * dimC, MPIU_REAL, &v0);CHKERRQ(ierr);
31139d150b73SToby Isaac       J    = &v0[dimC];
31149d150b73SToby Isaac       invJ = &J[dimC * dimC];
31159d150b73SToby Isaac       ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, &detJ);CHKERRQ(ierr);
31169d150b73SToby Isaac       for (i = 0; i < numPoints; i++) { /* Apply the inverse affine transformation for each point */
3117c330f8ffSToby Isaac         const PetscReal x0[3] = {-1.,-1.,-1.};
3118c330f8ffSToby Isaac 
3119c330f8ffSToby Isaac         CoordinatesRealToRef(dimC, dimR, x0, v0, invJ, &realCoords[dimC * i], &refCoords[dimR * i]);
31209d150b73SToby Isaac       }
312169291d52SBarry Smith       ierr = DMRestoreWorkArray(dm,dimC + 2 * dimC * dimC, MPIU_REAL, &v0);CHKERRQ(ierr);
31229d150b73SToby Isaac     } else if (isTensor) {
31239d150b73SToby Isaac       ierr = DMPlexCoordinatesToReference_Tensor(coordDM, cell, numPoints, realCoords, refCoords, coords, dimC, dimR);CHKERRQ(ierr);
31249d150b73SToby Isaac     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unrecognized cone size %D",coneSize);
31259d150b73SToby Isaac   } else {
31269d150b73SToby Isaac     ierr = DMPlexCoordinatesToReference_FE(coordDM, fe, cell, numPoints, realCoords, refCoords, coords, dimC, dimR);CHKERRQ(ierr);
31279d150b73SToby Isaac   }
31289d150b73SToby Isaac   PetscFunctionReturn(0);
31299d150b73SToby Isaac }
31309d150b73SToby Isaac 
31319d150b73SToby Isaac /*@
31329d150b73SToby Isaac   DMPlexReferenceToCoordinates - Map references coordinates to coordinates in the the mesh for a single element map.
31339d150b73SToby Isaac 
31349d150b73SToby Isaac   Not collective
31359d150b73SToby Isaac 
31369d150b73SToby Isaac   Input Parameters:
31379d150b73SToby Isaac + dm         - The mesh, with coordinate maps defined either by a PetscDS for the coordinate DM (see DMGetCoordinateDM()) or
31389d150b73SToby Isaac                implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or
31399d150b73SToby Isaac                as a multilinear map for tensor-product elements
31409d150b73SToby Isaac . cell       - the cell whose map is used.
31419d150b73SToby Isaac . numPoints  - the number of points to locate
3142a2b725a8SWilliam Gropp - refCoords  - (numPoints x dimension) array of reference coordinates (see DMGetDimension())
31439d150b73SToby Isaac 
31449d150b73SToby Isaac   Output Parameters:
31459d150b73SToby Isaac . realCoords - (numPoints x coordinate dimension) array of coordinates (see DMGetCoordinateDim())
31461b266c99SBarry Smith 
31471b266c99SBarry Smith    Level: intermediate
314873c9229bSMatthew Knepley 
314973c9229bSMatthew Knepley .seealso: DMPlexCoordinatesToReference()
31509d150b73SToby Isaac @*/
31519d150b73SToby Isaac PetscErrorCode DMPlexReferenceToCoordinates(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[])
31529d150b73SToby Isaac {
3153485ad865SMatthew G. Knepley   PetscInt       dimC, dimR, depth, cStart, cEnd, i;
31549d150b73SToby Isaac   DM             coordDM = NULL;
31559d150b73SToby Isaac   Vec            coords;
31569d150b73SToby Isaac   PetscFE        fe = NULL;
31579d150b73SToby Isaac   PetscErrorCode ierr;
31589d150b73SToby Isaac 
31599d150b73SToby Isaac   PetscFunctionBegin;
31609d150b73SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
31619d150b73SToby Isaac   ierr = DMGetDimension(dm,&dimR);CHKERRQ(ierr);
31629d150b73SToby Isaac   ierr = DMGetCoordinateDim(dm,&dimC);CHKERRQ(ierr);
31639d150b73SToby Isaac   if (dimR <= 0 || dimC <= 0 || numPoints <= 0) PetscFunctionReturn(0);
31649d150b73SToby Isaac   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
31659d150b73SToby Isaac   ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
31669d150b73SToby Isaac   ierr = DMGetCoordinateDM(dm,&coordDM);CHKERRQ(ierr);
31679d150b73SToby Isaac   if (coordDM) {
31689d150b73SToby Isaac     PetscInt coordFields;
31699d150b73SToby Isaac 
31709d150b73SToby Isaac     ierr = DMGetNumFields(coordDM,&coordFields);CHKERRQ(ierr);
31719d150b73SToby Isaac     if (coordFields) {
31729d150b73SToby Isaac       PetscClassId id;
31739d150b73SToby Isaac       PetscObject  disc;
31749d150b73SToby Isaac 
317544a7f3ddSMatthew G. Knepley       ierr = DMGetField(coordDM,0,NULL,&disc);CHKERRQ(ierr);
31769d150b73SToby Isaac       ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
31779d150b73SToby Isaac       if (id == PETSCFE_CLASSID) {
31789d150b73SToby Isaac         fe = (PetscFE) disc;
31799d150b73SToby Isaac       }
31809d150b73SToby Isaac     }
31819d150b73SToby Isaac   }
3182412e9a14SMatthew G. Knepley   ierr = DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
318313903a91SSatish 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);
31849d150b73SToby Isaac   if (!fe) { /* implicit discretization: affine or multilinear */
31859d150b73SToby Isaac     PetscInt  coneSize;
31869d150b73SToby Isaac     PetscBool isSimplex, isTensor;
31879d150b73SToby Isaac 
31889d150b73SToby Isaac     ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr);
31899d150b73SToby Isaac     isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE;
31909d150b73SToby Isaac     isTensor  = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE;
31919d150b73SToby Isaac     if (isSimplex) {
31929d150b73SToby Isaac       PetscReal detJ, *v0, *J;
31939d150b73SToby Isaac 
319469291d52SBarry Smith       ierr = DMGetWorkArray(dm,dimC + 2 * dimC * dimC, MPIU_REAL, &v0);CHKERRQ(ierr);
31959d150b73SToby Isaac       J    = &v0[dimC];
31969d150b73SToby Isaac       ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, NULL, &detJ);CHKERRQ(ierr);
3197c330f8ffSToby Isaac       for (i = 0; i < numPoints; i++) { /* Apply the affine transformation for each point */
3198c330f8ffSToby Isaac         const PetscReal xi0[3] = {-1.,-1.,-1.};
3199c330f8ffSToby Isaac 
3200c330f8ffSToby Isaac         CoordinatesRefToReal(dimC, dimR, xi0, v0, J, &refCoords[dimR * i], &realCoords[dimC * i]);
32019d150b73SToby Isaac       }
320269291d52SBarry Smith       ierr = DMRestoreWorkArray(dm,dimC + 2 * dimC * dimC, MPIU_REAL, &v0);CHKERRQ(ierr);
32039d150b73SToby Isaac     } else if (isTensor) {
32049d150b73SToby Isaac       ierr = DMPlexReferenceToCoordinates_Tensor(coordDM, cell, numPoints, refCoords, realCoords, coords, dimC, dimR);CHKERRQ(ierr);
32059d150b73SToby Isaac     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unrecognized cone size %D",coneSize);
32069d150b73SToby Isaac   } else {
32079d150b73SToby Isaac     ierr = DMPlexReferenceToCoordinates_FE(coordDM, fe, cell, numPoints, refCoords, realCoords, coords, dimC, dimR);CHKERRQ(ierr);
32089d150b73SToby Isaac   }
3209d6143a4eSToby Isaac   PetscFunctionReturn(0);
3210d6143a4eSToby Isaac }
32110139fca9SMatthew G. Knepley 
32120139fca9SMatthew G. Knepley /*@C
32130139fca9SMatthew G. Knepley   DMPlexRemapGeometry - This function maps the original DM coordinates to new coordinates.
32140139fca9SMatthew G. Knepley 
32150139fca9SMatthew G. Knepley   Not collective
32160139fca9SMatthew G. Knepley 
32170139fca9SMatthew G. Knepley   Input Parameters:
32180139fca9SMatthew G. Knepley + dm      - The DM
32190139fca9SMatthew G. Knepley . time    - The time
32200139fca9SMatthew G. Knepley - func    - The function transforming current coordinates to new coordaintes
32210139fca9SMatthew G. Knepley 
32220139fca9SMatthew G. Knepley    Calling sequence of func:
32230139fca9SMatthew G. Knepley $    func(PetscInt dim, PetscInt Nf, PetscInt NfAux,
32240139fca9SMatthew G. Knepley $         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
32250139fca9SMatthew G. Knepley $         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
32260139fca9SMatthew G. Knepley $         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]);
32270139fca9SMatthew G. Knepley 
32280139fca9SMatthew G. Knepley +  dim          - The spatial dimension
32290139fca9SMatthew G. Knepley .  Nf           - The number of input fields (here 1)
32300139fca9SMatthew G. Knepley .  NfAux        - The number of input auxiliary fields
32310139fca9SMatthew G. Knepley .  uOff         - The offset of the coordinates in u[] (here 0)
32320139fca9SMatthew G. Knepley .  uOff_x       - The offset of the coordinates in u_x[] (here 0)
32330139fca9SMatthew G. Knepley .  u            - The coordinate values at this point in space
32340139fca9SMatthew G. Knepley .  u_t          - The coordinate time derivative at this point in space (here NULL)
32350139fca9SMatthew G. Knepley .  u_x          - The coordinate derivatives at this point in space
32360139fca9SMatthew G. Knepley .  aOff         - The offset of each auxiliary field in u[]
32370139fca9SMatthew G. Knepley .  aOff_x       - The offset of each auxiliary field in u_x[]
32380139fca9SMatthew G. Knepley .  a            - The auxiliary field values at this point in space
32390139fca9SMatthew G. Knepley .  a_t          - The auxiliary field time derivative at this point in space (or NULL)
32400139fca9SMatthew G. Knepley .  a_x          - The auxiliary field derivatives at this point in space
32410139fca9SMatthew G. Knepley .  t            - The current time
32420139fca9SMatthew G. Knepley .  x            - The coordinates of this point (here not used)
32430139fca9SMatthew G. Knepley .  numConstants - The number of constants
32440139fca9SMatthew G. Knepley .  constants    - The value of each constant
32450139fca9SMatthew G. Knepley -  f            - The new coordinates at this point in space
32460139fca9SMatthew G. Knepley 
32470139fca9SMatthew G. Knepley   Level: intermediate
32480139fca9SMatthew G. Knepley 
32490139fca9SMatthew G. Knepley .seealso: DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMProjectFieldLocal(), DMProjectFieldLabelLocal()
32500139fca9SMatthew G. Knepley @*/
32510139fca9SMatthew G. Knepley PetscErrorCode DMPlexRemapGeometry(DM dm, PetscReal time,
32520139fca9SMatthew G. Knepley                                    void (*func)(PetscInt, PetscInt, PetscInt,
32530139fca9SMatthew G. Knepley                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
32540139fca9SMatthew G. Knepley                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
32550139fca9SMatthew G. Knepley                                                 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
32560139fca9SMatthew G. Knepley {
32570139fca9SMatthew G. Knepley   DM             cdm;
32588bf1a49fSMatthew G. Knepley   DMField        cf;
32590139fca9SMatthew G. Knepley   Vec            lCoords, tmpCoords;
32600139fca9SMatthew G. Knepley   PetscErrorCode ierr;
32610139fca9SMatthew G. Knepley 
32620139fca9SMatthew G. Knepley   PetscFunctionBegin;
32630139fca9SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
32640139fca9SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &lCoords);CHKERRQ(ierr);
32650139fca9SMatthew G. Knepley   ierr = DMGetLocalVector(cdm, &tmpCoords);CHKERRQ(ierr);
32660139fca9SMatthew G. Knepley   ierr = VecCopy(lCoords, tmpCoords);CHKERRQ(ierr);
32678bf1a49fSMatthew G. Knepley   /* We have to do the coordinate field manually right now since the coordinate DM will not have its own */
32688bf1a49fSMatthew G. Knepley   ierr = DMGetCoordinateField(dm, &cf);CHKERRQ(ierr);
32698bf1a49fSMatthew G. Knepley   cdm->coordinateField = cf;
32700139fca9SMatthew G. Knepley   ierr = DMProjectFieldLocal(cdm, time, tmpCoords, &func, INSERT_VALUES, lCoords);CHKERRQ(ierr);
32718bf1a49fSMatthew G. Knepley   cdm->coordinateField = NULL;
32720139fca9SMatthew G. Knepley   ierr = DMRestoreLocalVector(cdm, &tmpCoords);CHKERRQ(ierr);
3273cdaaecf7SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dm, lCoords);CHKERRQ(ierr);
32740139fca9SMatthew G. Knepley   PetscFunctionReturn(0);
32750139fca9SMatthew G. Knepley }
32760139fca9SMatthew G. Knepley 
32770139fca9SMatthew G. Knepley /* Shear applies the transformation, assuming we fix z,
32780139fca9SMatthew G. Knepley   / 1  0  m_0 \
32790139fca9SMatthew G. Knepley   | 0  1  m_1 |
32800139fca9SMatthew G. Knepley   \ 0  0   1  /
32810139fca9SMatthew G. Knepley */
32820139fca9SMatthew G. Knepley static void f0_shear(PetscInt dim, PetscInt Nf, PetscInt NfAux,
32830139fca9SMatthew G. Knepley                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
32840139fca9SMatthew G. Knepley                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
32850139fca9SMatthew G. Knepley                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar coords[])
32860139fca9SMatthew G. Knepley {
32870139fca9SMatthew G. Knepley   const PetscInt Nc = uOff[1]-uOff[0];
3288c1f1bd54SMatthew G. Knepley   const PetscInt ax = (PetscInt) PetscRealPart(constants[0]);
32890139fca9SMatthew G. Knepley   PetscInt       c;
32900139fca9SMatthew G. Knepley 
32910139fca9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
32920139fca9SMatthew G. Knepley     coords[c] = u[c] + constants[c+1]*u[ax];
32930139fca9SMatthew G. Knepley   }
32940139fca9SMatthew G. Knepley }
32950139fca9SMatthew G. Knepley 
32960139fca9SMatthew G. Knepley /*@C
32970139fca9SMatthew G. Knepley   DMPlexShearGeometry - This shears the domain, meaning adds a multiple of the shear coordinate to all other coordinates.
32980139fca9SMatthew G. Knepley 
32990139fca9SMatthew G. Knepley   Not collective
33000139fca9SMatthew G. Knepley 
33010139fca9SMatthew G. Knepley   Input Parameters:
33020139fca9SMatthew G. Knepley + dm          - The DM
33033ee9839eSMatthew G. Knepley . direction   - The shear coordinate direction, e.g. 0 is the x-axis
33040139fca9SMatthew G. Knepley - multipliers - The multiplier m for each direction which is not the shear direction
33050139fca9SMatthew G. Knepley 
33060139fca9SMatthew G. Knepley   Level: intermediate
33070139fca9SMatthew G. Knepley 
33080139fca9SMatthew G. Knepley .seealso: DMPlexRemapGeometry()
33090139fca9SMatthew G. Knepley @*/
33103ee9839eSMatthew G. Knepley PetscErrorCode DMPlexShearGeometry(DM dm, DMDirection direction, PetscReal multipliers[])
33110139fca9SMatthew G. Knepley {
33120139fca9SMatthew G. Knepley   DM             cdm;
33130139fca9SMatthew G. Knepley   PetscDS        cds;
33140139fca9SMatthew G. Knepley   PetscObject    obj;
33150139fca9SMatthew G. Knepley   PetscClassId   id;
33160139fca9SMatthew G. Knepley   PetscScalar   *moduli;
33173ee9839eSMatthew G. Knepley   const PetscInt dir = (PetscInt) direction;
33180139fca9SMatthew G. Knepley   PetscInt       dE, d, e;
33190139fca9SMatthew G. Knepley   PetscErrorCode ierr;
33200139fca9SMatthew G. Knepley 
33210139fca9SMatthew G. Knepley   PetscFunctionBegin;
33220139fca9SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
33230139fca9SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dE);CHKERRQ(ierr);
33240139fca9SMatthew G. Knepley   ierr = PetscMalloc1(dE+1, &moduli);CHKERRQ(ierr);
33250139fca9SMatthew G. Knepley   moduli[0] = dir;
3326cdaaecf7SMatthew G. Knepley   for (d = 0, e = 0; d < dE; ++d) moduli[d+1] = d == dir ? 0.0 : (multipliers ? multipliers[e++] : 1.0);
33270139fca9SMatthew G. Knepley   ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr);
33280139fca9SMatthew G. Knepley   ierr = PetscDSGetDiscretization(cds, 0, &obj);CHKERRQ(ierr);
33290139fca9SMatthew G. Knepley   ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
33300139fca9SMatthew G. Knepley   if (id != PETSCFE_CLASSID) {
33310139fca9SMatthew G. Knepley     Vec           lCoords;
33320139fca9SMatthew G. Knepley     PetscSection  cSection;
33330139fca9SMatthew G. Knepley     PetscScalar  *coords;
33340139fca9SMatthew G. Knepley     PetscInt      vStart, vEnd, v;
33350139fca9SMatthew G. Knepley 
33360139fca9SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
33370139fca9SMatthew G. Knepley     ierr = DMGetCoordinateSection(dm, &cSection);CHKERRQ(ierr);
33380139fca9SMatthew G. Knepley     ierr = DMGetCoordinatesLocal(dm, &lCoords);CHKERRQ(ierr);
33390139fca9SMatthew G. Knepley     ierr = VecGetArray(lCoords, &coords);CHKERRQ(ierr);
33400139fca9SMatthew G. Knepley     for (v = vStart; v < vEnd; ++v) {
33410139fca9SMatthew G. Knepley       PetscReal ds;
33420139fca9SMatthew G. Knepley       PetscInt  off, c;
33430139fca9SMatthew G. Knepley 
33440139fca9SMatthew G. Knepley       ierr = PetscSectionGetOffset(cSection, v, &off);CHKERRQ(ierr);
33450139fca9SMatthew G. Knepley       ds   = PetscRealPart(coords[off+dir]);
33460139fca9SMatthew G. Knepley       for (c = 0; c < dE; ++c) coords[off+c] += moduli[c]*ds;
33470139fca9SMatthew G. Knepley     }
33480139fca9SMatthew G. Knepley     ierr = VecRestoreArray(lCoords, &coords);CHKERRQ(ierr);
33490139fca9SMatthew G. Knepley   } else {
33500139fca9SMatthew G. Knepley     ierr = PetscDSSetConstants(cds, dE+1, moduli);CHKERRQ(ierr);
33510139fca9SMatthew G. Knepley     ierr = DMPlexRemapGeometry(dm, 0.0, f0_shear);CHKERRQ(ierr);
33520139fca9SMatthew G. Knepley   }
33530139fca9SMatthew G. Knepley   ierr = PetscFree(moduli);CHKERRQ(ierr);
33540139fca9SMatthew G. Knepley   PetscFunctionReturn(0);
33550139fca9SMatthew G. Knepley }
3356