xref: /petsc/src/dm/impls/plex/plexsubmesh.c (revision efa14ee051fe8a2bfbd76140cdcdff4f22172457)
1e6ccafaeSMatthew G Knepley #include <petsc-private/pleximpl.h>   /*I      "petscdmplex.h"   I*/
2e6ccafaeSMatthew G Knepley 
3e6ccafaeSMatthew G Knepley extern PetscErrorCode DMPlexGetNumFaceVertices_Internal(DM, PetscInt, PetscInt *);
4e6ccafaeSMatthew G Knepley 
5e6ccafaeSMatthew G Knepley #undef __FUNCT__
6*efa14ee0SMatthew G Knepley #define __FUNCT__ "DMPlexMarkSubmesh_Uninterpolated"
7*efa14ee0SMatthew G Knepley /* Here we need the explicit assumption that:
8*efa14ee0SMatthew G Knepley 
9*efa14ee0SMatthew G Knepley      For any marked cell, the marked vertices constitute a single face
10*efa14ee0SMatthew G Knepley */
11*efa14ee0SMatthew G Knepley static PetscErrorCode DMPlexMarkSubmesh_Uninterpolated(DM dm, DMLabel vertexLabel, DMLabel subpointMap, PetscInt *numFaces, PetscInt *nFV, DM subdm)
12*efa14ee0SMatthew G Knepley {
13*efa14ee0SMatthew G Knepley   IS               subvertexIS;
14*efa14ee0SMatthew G Knepley   const PetscInt  *subvertices;
15*efa14ee0SMatthew G Knepley   PetscInt        *pStart, *pEnd, *pMax;
16*efa14ee0SMatthew G Knepley   PetscInt         depth, dim, d, numSubVerticesInitial = 0, v;
17*efa14ee0SMatthew G Knepley   PetscErrorCode   ierr;
18*efa14ee0SMatthew G Knepley 
19*efa14ee0SMatthew G Knepley   PetscFunctionBegin;
20*efa14ee0SMatthew G Knepley   *numFaces = 0;
21*efa14ee0SMatthew G Knepley   *nFV      = 0;
22*efa14ee0SMatthew G Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
23*efa14ee0SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
24*efa14ee0SMatthew G Knepley   ierr = PetscMalloc3(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd,dim+1,PetscInt,&pMax);CHKERRQ(ierr);
25*efa14ee0SMatthew G Knepley   ierr = DMPlexGetHybridBounds(dm, &pMax[depth], depth>1 ? &pMax[depth-1] : PETSC_NULL, depth > 2 ? &pMax[1] : PETSC_NULL, &pMax[0]);CHKERRQ(ierr);
26*efa14ee0SMatthew G Knepley   for (d = 0; d <= depth; ++d) {
27*efa14ee0SMatthew G Knepley     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
28*efa14ee0SMatthew G Knepley     if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]);
29*efa14ee0SMatthew G Knepley   }
30*efa14ee0SMatthew G Knepley   /* Loop over initial vertices and mark all faces in the collective star() */
31*efa14ee0SMatthew G Knepley   ierr = DMLabelGetStratumIS(vertexLabel, 1, &subvertexIS);CHKERRQ(ierr);
32*efa14ee0SMatthew G Knepley   if (subvertexIS) {
33*efa14ee0SMatthew G Knepley     ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr);
34*efa14ee0SMatthew G Knepley     ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
35*efa14ee0SMatthew G Knepley   }
36*efa14ee0SMatthew G Knepley   for (v = 0; v < numSubVerticesInitial; ++v) {
37*efa14ee0SMatthew G Knepley     const PetscInt vertex = subvertices[v];
38*efa14ee0SMatthew G Knepley     PetscInt      *star   = PETSC_NULL;
39*efa14ee0SMatthew G Knepley     PetscInt       starSize, s, numCells = 0, c;
40*efa14ee0SMatthew G Knepley 
41*efa14ee0SMatthew G Knepley     ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
42*efa14ee0SMatthew G Knepley     for (s = 0; s < starSize*2; s += 2) {
43*efa14ee0SMatthew G Knepley       const PetscInt point = star[s];
44*efa14ee0SMatthew G Knepley       if ((point >= pStart[depth]) && (point < pEnd[depth])) star[numCells++] = point;
45*efa14ee0SMatthew G Knepley     }
46*efa14ee0SMatthew G Knepley     for (c = 0; c < numCells; ++c) {
47*efa14ee0SMatthew G Knepley       const PetscInt cell    = star[c];
48*efa14ee0SMatthew G Knepley       PetscInt      *closure = PETSC_NULL;
49*efa14ee0SMatthew G Knepley       PetscInt       closureSize, cl;
50*efa14ee0SMatthew G Knepley       PetscInt       cellLoc, numCorners = 0, faceSize = 0;
51*efa14ee0SMatthew G Knepley 
52*efa14ee0SMatthew G Knepley       ierr = DMLabelGetValue(subpointMap, cell, &cellLoc);CHKERRQ(ierr);
53*efa14ee0SMatthew G Knepley       if (cellLoc == dim) continue;
54*efa14ee0SMatthew G Knepley       if (cellLoc >= 0) SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "Cell %d has dimension %d in the surface label", cell, cellLoc);
55*efa14ee0SMatthew G Knepley       ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
56*efa14ee0SMatthew G Knepley       for (cl = 0; cl < closureSize*2; cl += 2) {
57*efa14ee0SMatthew G Knepley         const PetscInt point = closure[cl];
58*efa14ee0SMatthew G Knepley         PetscInt       vertexLoc;
59*efa14ee0SMatthew G Knepley 
60*efa14ee0SMatthew G Knepley         if ((point >= pStart[0]) && (point < pEnd[0])) {
61*efa14ee0SMatthew G Knepley           ++numCorners;
62*efa14ee0SMatthew G Knepley           ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr);
63*efa14ee0SMatthew G Knepley           if (vertexLoc >= 0) closure[faceSize++] = point;
64*efa14ee0SMatthew G Knepley         }
65*efa14ee0SMatthew G Knepley       }
66*efa14ee0SMatthew G Knepley       if (!(*nFV)) {ierr = DMPlexGetNumFaceVertices_Internal(dm, numCorners, nFV);CHKERRQ(ierr);}
67*efa14ee0SMatthew G Knepley       if (faceSize > *nFV) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize);
68*efa14ee0SMatthew G Knepley       if (faceSize == *nFV) {
69*efa14ee0SMatthew G Knepley         ++(*numFaces);
70*efa14ee0SMatthew G Knepley         for (cl = 0; cl < faceSize; ++cl) {
71*efa14ee0SMatthew G Knepley           ierr = DMLabelSetValue(subpointMap, closure[cl], 0);CHKERRQ(ierr);
72*efa14ee0SMatthew G Knepley         }
73*efa14ee0SMatthew G Knepley         ierr = DMLabelSetValue(subpointMap, cell, 2);CHKERRQ(ierr);
74*efa14ee0SMatthew G Knepley       }
75*efa14ee0SMatthew G Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
76*efa14ee0SMatthew G Knepley     }
77*efa14ee0SMatthew G Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
78*efa14ee0SMatthew G Knepley   }
79*efa14ee0SMatthew G Knepley   if (subvertexIS) {
80*efa14ee0SMatthew G Knepley     ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
81*efa14ee0SMatthew G Knepley   }
82*efa14ee0SMatthew G Knepley   ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr);
83*efa14ee0SMatthew G Knepley   ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr);
84*efa14ee0SMatthew G Knepley   PetscFunctionReturn(0);
85*efa14ee0SMatthew G Knepley }
86*efa14ee0SMatthew G Knepley 
87*efa14ee0SMatthew G Knepley #undef __FUNCT__
88*efa14ee0SMatthew G Knepley #define __FUNCT__ "DMPlexMarkSubmesh_Interpolated"
89*efa14ee0SMatthew G Knepley static PetscErrorCode DMPlexMarkSubmesh_Interpolated(DM dm, DMLabel vertexLabel, DMLabel subpointMap, DM subdm)
90*efa14ee0SMatthew G Knepley {
91*efa14ee0SMatthew G Knepley   IS               subvertexIS;
92*efa14ee0SMatthew G Knepley   const PetscInt  *subvertices;
93*efa14ee0SMatthew G Knepley   PetscInt        *pStart, *pEnd, *pMax;
94*efa14ee0SMatthew G Knepley   PetscInt         dim, d, numSubVerticesInitial = 0, v;
95*efa14ee0SMatthew G Knepley   PetscErrorCode   ierr;
96*efa14ee0SMatthew G Knepley 
97*efa14ee0SMatthew G Knepley   PetscFunctionBegin;
98*efa14ee0SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
99*efa14ee0SMatthew G Knepley   ierr = PetscMalloc3(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd,dim+1,PetscInt,&pMax);CHKERRQ(ierr);
100*efa14ee0SMatthew G Knepley   ierr = DMPlexGetHybridBounds(dm, &pMax[dim], dim>1 ? &pMax[dim-1] : PETSC_NULL, dim > 2 ? &pMax[1] : PETSC_NULL, &pMax[0]);CHKERRQ(ierr);
101*efa14ee0SMatthew G Knepley   for (d = 0; d <= dim; ++d) {
102*efa14ee0SMatthew G Knepley     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
103*efa14ee0SMatthew G Knepley     if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]);
104*efa14ee0SMatthew G Knepley   }
105*efa14ee0SMatthew G Knepley   /* Loop over initial vertices and mark all faces in the collective star() */
106*efa14ee0SMatthew G Knepley   ierr = DMLabelGetStratumIS(vertexLabel, 1, &subvertexIS);CHKERRQ(ierr);
107*efa14ee0SMatthew G Knepley   if (subvertexIS) {
108*efa14ee0SMatthew G Knepley     ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr);
109*efa14ee0SMatthew G Knepley     ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
110*efa14ee0SMatthew G Knepley   }
111*efa14ee0SMatthew G Knepley   for (v = 0; v < numSubVerticesInitial; ++v) {
112*efa14ee0SMatthew G Knepley     const PetscInt vertex = subvertices[v];
113*efa14ee0SMatthew G Knepley     PetscInt      *star   = PETSC_NULL;
114*efa14ee0SMatthew G Knepley     PetscInt       starSize, s, numFaces = 0, f;
115*efa14ee0SMatthew G Knepley 
116*efa14ee0SMatthew G Knepley     ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
117*efa14ee0SMatthew G Knepley     for (s = 0; s < starSize*2; s += 2) {
118*efa14ee0SMatthew G Knepley       const PetscInt point = star[s];
119*efa14ee0SMatthew G Knepley       if ((point >= pStart[dim-1]) && (point < pEnd[dim-1])) star[numFaces++] = point;
120*efa14ee0SMatthew G Knepley     }
121*efa14ee0SMatthew G Knepley     for (f = 0; f < numFaces; ++f) {
122*efa14ee0SMatthew G Knepley       const PetscInt face    = star[f];
123*efa14ee0SMatthew G Knepley       PetscInt      *closure = PETSC_NULL;
124*efa14ee0SMatthew G Knepley       PetscInt       closureSize, c;
125*efa14ee0SMatthew G Knepley       PetscInt       faceLoc;
126*efa14ee0SMatthew G Knepley 
127*efa14ee0SMatthew G Knepley       ierr = DMLabelGetValue(subpointMap, face, &faceLoc);CHKERRQ(ierr);
128*efa14ee0SMatthew G Knepley       if (faceLoc == dim-1) continue;
129*efa14ee0SMatthew G Knepley       if (faceLoc >= 0) SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "Face %d has dimension %d in the surface label", face, faceLoc);
130*efa14ee0SMatthew G Knepley       ierr = DMPlexGetTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
131*efa14ee0SMatthew G Knepley       for (c = 0; c < closureSize*2; c += 2) {
132*efa14ee0SMatthew G Knepley         const PetscInt point = closure[c];
133*efa14ee0SMatthew G Knepley         PetscInt       vertexLoc;
134*efa14ee0SMatthew G Knepley 
135*efa14ee0SMatthew G Knepley         if ((point >= pStart[0]) && (point < pEnd[0])) {
136*efa14ee0SMatthew G Knepley           ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr);
137*efa14ee0SMatthew G Knepley           if (vertexLoc < 0) break;
138*efa14ee0SMatthew G Knepley         }
139*efa14ee0SMatthew G Knepley       }
140*efa14ee0SMatthew G Knepley       if (c == closureSize*2) {
141*efa14ee0SMatthew G Knepley         const PetscInt *support;
142*efa14ee0SMatthew G Knepley         PetscInt        supportSize, s;
143*efa14ee0SMatthew G Knepley 
144*efa14ee0SMatthew G Knepley         for (c = 0; c < closureSize*2; c += 2) {
145*efa14ee0SMatthew G Knepley           const PetscInt point = closure[c];
146*efa14ee0SMatthew G Knepley 
147*efa14ee0SMatthew G Knepley           for (d = 0; d < dim; ++d) {
148*efa14ee0SMatthew G Knepley             if ((point >= pStart[d]) && (point < pEnd[d])) {
149*efa14ee0SMatthew G Knepley               ierr = DMLabelSetValue(subpointMap, point, d);CHKERRQ(ierr);
150*efa14ee0SMatthew G Knepley               break;
151*efa14ee0SMatthew G Knepley             }
152*efa14ee0SMatthew G Knepley           }
153*efa14ee0SMatthew G Knepley         }
154*efa14ee0SMatthew G Knepley         ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr);
155*efa14ee0SMatthew G Knepley         ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr);
156*efa14ee0SMatthew G Knepley         for (s = 0; s < supportSize; ++s) {
157*efa14ee0SMatthew G Knepley           ierr = DMLabelSetValue(subpointMap, support[s], dim);CHKERRQ(ierr);
158*efa14ee0SMatthew G Knepley         }
159*efa14ee0SMatthew G Knepley       }
160*efa14ee0SMatthew G Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
161*efa14ee0SMatthew G Knepley     }
162*efa14ee0SMatthew G Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
163*efa14ee0SMatthew G Knepley   }
164*efa14ee0SMatthew G Knepley   if (subvertexIS) {
165*efa14ee0SMatthew G Knepley     ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr);
166*efa14ee0SMatthew G Knepley   }
167*efa14ee0SMatthew G Knepley   ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr);
168*efa14ee0SMatthew G Knepley   ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr);
169*efa14ee0SMatthew G Knepley   PetscFunctionReturn(0);
170*efa14ee0SMatthew G Knepley }
171*efa14ee0SMatthew G Knepley 
172*efa14ee0SMatthew G Knepley #undef __FUNCT__
173e6ccafaeSMatthew G Knepley #define __FUNCT__ "DMPlexGetFaceOrientation"
174e6ccafaeSMatthew G Knepley PetscErrorCode DMPlexGetFaceOrientation(DM dm, PetscInt cell, PetscInt numCorners, PetscInt indices[], PetscInt oppositeVertex, PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented)
175e6ccafaeSMatthew G Knepley {
176e6ccafaeSMatthew G Knepley   MPI_Comm       comm      = ((PetscObject) dm)->comm;
177e6ccafaeSMatthew G Knepley   PetscBool      posOrient = PETSC_FALSE;
178e6ccafaeSMatthew G Knepley   const PetscInt debug     = 0;
179e6ccafaeSMatthew G Knepley   PetscInt       cellDim, faceSize, f;
180e6ccafaeSMatthew G Knepley   PetscErrorCode ierr;
181e6ccafaeSMatthew G Knepley 
182e6ccafaeSMatthew G Knepley   ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr);
183e6ccafaeSMatthew G Knepley   if (debug) {PetscPrintf(comm, "cellDim: %d numCorners: %d\n", cellDim, numCorners);CHKERRQ(ierr);}
184e6ccafaeSMatthew G Knepley 
185e6ccafaeSMatthew G Knepley   if (cellDim == numCorners-1) {
186e6ccafaeSMatthew G Knepley     /* Simplices */
187e6ccafaeSMatthew G Knepley     faceSize  = numCorners-1;
188e6ccafaeSMatthew G Knepley     posOrient = !(oppositeVertex%2) ? PETSC_TRUE : PETSC_FALSE;
189e6ccafaeSMatthew G Knepley   } else if (cellDim == 1 && numCorners == 3) {
190e6ccafaeSMatthew G Knepley     /* Quadratic line */
191e6ccafaeSMatthew G Knepley     faceSize  = 1;
192e6ccafaeSMatthew G Knepley     posOrient = PETSC_TRUE;
193e6ccafaeSMatthew G Knepley   } else if (cellDim == 2 && numCorners == 4) {
194e6ccafaeSMatthew G Knepley     /* Quads */
195e6ccafaeSMatthew G Knepley     faceSize = 2;
196e6ccafaeSMatthew G Knepley     if ((indices[1] > indices[0]) && (indices[1] - indices[0] == 1)) {
197e6ccafaeSMatthew G Knepley       posOrient = PETSC_TRUE;
198e6ccafaeSMatthew G Knepley     } else if ((indices[0] == 3) && (indices[1] == 0)) {
199e6ccafaeSMatthew G Knepley       posOrient = PETSC_TRUE;
200e6ccafaeSMatthew G Knepley     } else {
201e6ccafaeSMatthew G Knepley       if (((indices[0] > indices[1]) && (indices[0] - indices[1] == 1)) || ((indices[0] == 0) && (indices[1] == 3))) {
202e6ccafaeSMatthew G Knepley         posOrient = PETSC_FALSE;
203e6ccafaeSMatthew G Knepley       } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossedge");
204e6ccafaeSMatthew G Knepley     }
205e6ccafaeSMatthew G Knepley   } else if (cellDim == 2 && numCorners == 6) {
206e6ccafaeSMatthew G Knepley     /* Quadratic triangle (I hate this) */
207e6ccafaeSMatthew G Knepley     /* Edges are determined by the first 2 vertices (corners of edges) */
208e6ccafaeSMatthew G Knepley     const PetscInt faceSizeTri = 3;
209e6ccafaeSMatthew G Knepley     PetscInt       sortedIndices[3], i, iFace;
210e6ccafaeSMatthew G Knepley     PetscBool      found                    = PETSC_FALSE;
211e6ccafaeSMatthew G Knepley     PetscInt       faceVerticesTriSorted[9] = {
212e6ccafaeSMatthew G Knepley       0, 3,  4, /* bottom */
213e6ccafaeSMatthew G Knepley       1, 4,  5, /* right */
214e6ccafaeSMatthew G Knepley       2, 3,  5, /* left */
215e6ccafaeSMatthew G Knepley     };
216e6ccafaeSMatthew G Knepley     PetscInt       faceVerticesTri[9] = {
217e6ccafaeSMatthew G Knepley       0, 3,  4, /* bottom */
218e6ccafaeSMatthew G Knepley       1, 4,  5, /* right */
219e6ccafaeSMatthew G Knepley       2, 5,  3, /* left */
220e6ccafaeSMatthew G Knepley     };
221e6ccafaeSMatthew G Knepley 
222e6ccafaeSMatthew G Knepley     faceSize = faceSizeTri;
223e6ccafaeSMatthew G Knepley     for (i = 0; i < faceSizeTri; ++i) sortedIndices[i] = indices[i];
224e6ccafaeSMatthew G Knepley     ierr = PetscSortInt(faceSizeTri, sortedIndices);CHKERRQ(ierr);
225e6ccafaeSMatthew G Knepley     for (iFace = 0; iFace < 3; ++iFace) {
226e6ccafaeSMatthew G Knepley       const PetscInt ii = iFace*faceSizeTri;
227e6ccafaeSMatthew G Knepley       PetscInt       fVertex, cVertex;
228e6ccafaeSMatthew G Knepley 
229e6ccafaeSMatthew G Knepley       if ((sortedIndices[0] == faceVerticesTriSorted[ii+0]) &&
230e6ccafaeSMatthew G Knepley           (sortedIndices[1] == faceVerticesTriSorted[ii+1])) {
231e6ccafaeSMatthew G Knepley         for (fVertex = 0; fVertex < faceSizeTri; ++fVertex) {
232e6ccafaeSMatthew G Knepley           for (cVertex = 0; cVertex < faceSizeTri; ++cVertex) {
233e6ccafaeSMatthew G Knepley             if (indices[cVertex] == faceVerticesTri[ii+fVertex]) {
234e6ccafaeSMatthew G Knepley               faceVertices[fVertex] = origVertices[cVertex];
235e6ccafaeSMatthew G Knepley               break;
236e6ccafaeSMatthew G Knepley             }
237e6ccafaeSMatthew G Knepley           }
238e6ccafaeSMatthew G Knepley         }
239e6ccafaeSMatthew G Knepley         found = PETSC_TRUE;
240e6ccafaeSMatthew G Knepley         break;
241e6ccafaeSMatthew G Knepley       }
242e6ccafaeSMatthew G Knepley     }
243e6ccafaeSMatthew G Knepley     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tri crossface");
244e6ccafaeSMatthew G Knepley     if (posOriented) *posOriented = PETSC_TRUE;
245e6ccafaeSMatthew G Knepley     PetscFunctionReturn(0);
246e6ccafaeSMatthew G Knepley   } else if (cellDim == 2 && numCorners == 9) {
247e6ccafaeSMatthew G Knepley     /* Quadratic quad (I hate this) */
248e6ccafaeSMatthew G Knepley     /* Edges are determined by the first 2 vertices (corners of edges) */
249e6ccafaeSMatthew G Knepley     const PetscInt faceSizeQuad = 3;
250e6ccafaeSMatthew G Knepley     PetscInt       sortedIndices[3], i, iFace;
251e6ccafaeSMatthew G Knepley     PetscBool      found                      = PETSC_FALSE;
252e6ccafaeSMatthew G Knepley     PetscInt       faceVerticesQuadSorted[12] = {
253e6ccafaeSMatthew G Knepley       0, 1,  4, /* bottom */
254e6ccafaeSMatthew G Knepley       1, 2,  5, /* right */
255e6ccafaeSMatthew G Knepley       2, 3,  6, /* top */
256e6ccafaeSMatthew G Knepley       0, 3,  7, /* left */
257e6ccafaeSMatthew G Knepley     };
258e6ccafaeSMatthew G Knepley     PetscInt       faceVerticesQuad[12] = {
259e6ccafaeSMatthew G Knepley       0, 1,  4, /* bottom */
260e6ccafaeSMatthew G Knepley       1, 2,  5, /* right */
261e6ccafaeSMatthew G Knepley       2, 3,  6, /* top */
262e6ccafaeSMatthew G Knepley       3, 0,  7, /* left */
263e6ccafaeSMatthew G Knepley     };
264e6ccafaeSMatthew G Knepley 
265e6ccafaeSMatthew G Knepley     faceSize = faceSizeQuad;
266e6ccafaeSMatthew G Knepley     for (i = 0; i < faceSizeQuad; ++i) sortedIndices[i] = indices[i];
267e6ccafaeSMatthew G Knepley     ierr = PetscSortInt(faceSizeQuad, sortedIndices);CHKERRQ(ierr);
268e6ccafaeSMatthew G Knepley     for (iFace = 0; iFace < 4; ++iFace) {
269e6ccafaeSMatthew G Knepley       const PetscInt ii = iFace*faceSizeQuad;
270e6ccafaeSMatthew G Knepley       PetscInt       fVertex, cVertex;
271e6ccafaeSMatthew G Knepley 
272e6ccafaeSMatthew G Knepley       if ((sortedIndices[0] == faceVerticesQuadSorted[ii+0]) &&
273e6ccafaeSMatthew G Knepley           (sortedIndices[1] == faceVerticesQuadSorted[ii+1])) {
274e6ccafaeSMatthew G Knepley         for (fVertex = 0; fVertex < faceSizeQuad; ++fVertex) {
275e6ccafaeSMatthew G Knepley           for (cVertex = 0; cVertex < faceSizeQuad; ++cVertex) {
276e6ccafaeSMatthew G Knepley             if (indices[cVertex] == faceVerticesQuad[ii+fVertex]) {
277e6ccafaeSMatthew G Knepley               faceVertices[fVertex] = origVertices[cVertex];
278e6ccafaeSMatthew G Knepley               break;
279e6ccafaeSMatthew G Knepley             }
280e6ccafaeSMatthew G Knepley           }
281e6ccafaeSMatthew G Knepley         }
282e6ccafaeSMatthew G Knepley         found = PETSC_TRUE;
283e6ccafaeSMatthew G Knepley         break;
284e6ccafaeSMatthew G Knepley       }
285e6ccafaeSMatthew G Knepley     }
286e6ccafaeSMatthew G Knepley     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossface");
287e6ccafaeSMatthew G Knepley     if (posOriented) *posOriented = PETSC_TRUE;
288e6ccafaeSMatthew G Knepley     PetscFunctionReturn(0);
289e6ccafaeSMatthew G Knepley   } else if (cellDim == 3 && numCorners == 8) {
290e6ccafaeSMatthew G Knepley     /* Hexes
291e6ccafaeSMatthew G Knepley        A hex is two oriented quads with the normal of the first
292e6ccafaeSMatthew G Knepley        pointing up at the second.
293e6ccafaeSMatthew G Knepley 
294e6ccafaeSMatthew G Knepley           7---6
295e6ccafaeSMatthew G Knepley          /|  /|
296e6ccafaeSMatthew G Knepley         4---5 |
297e6ccafaeSMatthew G Knepley         | 3-|-2
298e6ccafaeSMatthew G Knepley         |/  |/
299e6ccafaeSMatthew G Knepley         0---1
300e6ccafaeSMatthew G Knepley 
301e6ccafaeSMatthew G Knepley         Faces are determined by the first 4 vertices (corners of faces) */
302e6ccafaeSMatthew G Knepley     const PetscInt faceSizeHex = 4;
303e6ccafaeSMatthew G Knepley     PetscInt       sortedIndices[4], i, iFace;
304e6ccafaeSMatthew G Knepley     PetscBool      found                     = PETSC_FALSE;
305e6ccafaeSMatthew G Knepley     PetscInt       faceVerticesHexSorted[24] = {
306e6ccafaeSMatthew G Knepley       0, 1, 2, 3,  /* bottom */
307e6ccafaeSMatthew G Knepley       4, 5, 6, 7,  /* top */
308e6ccafaeSMatthew G Knepley       0, 1, 4, 5,  /* front */
309e6ccafaeSMatthew G Knepley       1, 2, 5, 6,  /* right */
310e6ccafaeSMatthew G Knepley       2, 3, 6, 7,  /* back */
311e6ccafaeSMatthew G Knepley       0, 3, 4, 7,  /* left */
312e6ccafaeSMatthew G Knepley     };
313e6ccafaeSMatthew G Knepley     PetscInt       faceVerticesHex[24] = {
314e6ccafaeSMatthew G Knepley       3, 2, 1, 0,  /* bottom */
315e6ccafaeSMatthew G Knepley       4, 5, 6, 7,  /* top */
316e6ccafaeSMatthew G Knepley       0, 1, 5, 4,  /* front */
317e6ccafaeSMatthew G Knepley       1, 2, 6, 5,  /* right */
318e6ccafaeSMatthew G Knepley       2, 3, 7, 6,  /* back */
319e6ccafaeSMatthew G Knepley       3, 0, 4, 7,  /* left */
320e6ccafaeSMatthew G Knepley     };
321e6ccafaeSMatthew G Knepley 
322e6ccafaeSMatthew G Knepley     faceSize = faceSizeHex;
323e6ccafaeSMatthew G Knepley     for (i = 0; i < faceSizeHex; ++i) sortedIndices[i] = indices[i];
324e6ccafaeSMatthew G Knepley     ierr = PetscSortInt(faceSizeHex, sortedIndices);CHKERRQ(ierr);
325e6ccafaeSMatthew G Knepley     for (iFace = 0; iFace < 6; ++iFace) {
326e6ccafaeSMatthew G Knepley       const PetscInt ii = iFace*faceSizeHex;
327e6ccafaeSMatthew G Knepley       PetscInt       fVertex, cVertex;
328e6ccafaeSMatthew G Knepley 
329e6ccafaeSMatthew G Knepley       if ((sortedIndices[0] == faceVerticesHexSorted[ii+0]) &&
330e6ccafaeSMatthew G Knepley           (sortedIndices[1] == faceVerticesHexSorted[ii+1]) &&
331e6ccafaeSMatthew G Knepley           (sortedIndices[2] == faceVerticesHexSorted[ii+2]) &&
332e6ccafaeSMatthew G Knepley           (sortedIndices[3] == faceVerticesHexSorted[ii+3])) {
333e6ccafaeSMatthew G Knepley         for (fVertex = 0; fVertex < faceSizeHex; ++fVertex) {
334e6ccafaeSMatthew G Knepley           for (cVertex = 0; cVertex < faceSizeHex; ++cVertex) {
335e6ccafaeSMatthew G Knepley             if (indices[cVertex] == faceVerticesHex[ii+fVertex]) {
336e6ccafaeSMatthew G Knepley               faceVertices[fVertex] = origVertices[cVertex];
337e6ccafaeSMatthew G Knepley               break;
338e6ccafaeSMatthew G Knepley             }
339e6ccafaeSMatthew G Knepley           }
340e6ccafaeSMatthew G Knepley         }
341e6ccafaeSMatthew G Knepley         found = PETSC_TRUE;
342e6ccafaeSMatthew G Knepley         break;
343e6ccafaeSMatthew G Knepley       }
344e6ccafaeSMatthew G Knepley     }
345e6ccafaeSMatthew G Knepley     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface");
346e6ccafaeSMatthew G Knepley     if (posOriented) *posOriented = PETSC_TRUE;
347e6ccafaeSMatthew G Knepley     PetscFunctionReturn(0);
348e6ccafaeSMatthew G Knepley   } else if (cellDim == 3 && numCorners == 10) {
349e6ccafaeSMatthew G Knepley     /* Quadratic tet */
350e6ccafaeSMatthew G Knepley     /* Faces are determined by the first 3 vertices (corners of faces) */
351e6ccafaeSMatthew G Knepley     const PetscInt faceSizeTet = 6;
352e6ccafaeSMatthew G Knepley     PetscInt       sortedIndices[6], i, iFace;
353e6ccafaeSMatthew G Knepley     PetscBool      found                     = PETSC_FALSE;
354e6ccafaeSMatthew G Knepley     PetscInt       faceVerticesTetSorted[24] = {
355e6ccafaeSMatthew G Knepley       0, 1, 2,  6, 7, 8, /* bottom */
356e6ccafaeSMatthew G Knepley       0, 3, 4,  6, 7, 9,  /* front */
357e6ccafaeSMatthew G Knepley       1, 4, 5,  7, 8, 9,  /* right */
358e6ccafaeSMatthew G Knepley       2, 3, 5,  6, 8, 9,  /* left */
359e6ccafaeSMatthew G Knepley     };
360e6ccafaeSMatthew G Knepley     PetscInt       faceVerticesTet[24] = {
361e6ccafaeSMatthew G Knepley       0, 1, 2,  6, 7, 8, /* bottom */
362e6ccafaeSMatthew G Knepley       0, 4, 3,  6, 7, 9,  /* front */
363e6ccafaeSMatthew G Knepley       1, 5, 4,  7, 8, 9,  /* right */
364e6ccafaeSMatthew G Knepley       2, 3, 5,  8, 6, 9,  /* left */
365e6ccafaeSMatthew G Knepley     };
366e6ccafaeSMatthew G Knepley 
367e6ccafaeSMatthew G Knepley     faceSize = faceSizeTet;
368e6ccafaeSMatthew G Knepley     for (i = 0; i < faceSizeTet; ++i) sortedIndices[i] = indices[i];
369e6ccafaeSMatthew G Knepley     ierr = PetscSortInt(faceSizeTet, sortedIndices);CHKERRQ(ierr);
370e6ccafaeSMatthew G Knepley     for (iFace=0; iFace < 4; ++iFace) {
371e6ccafaeSMatthew G Knepley       const PetscInt ii = iFace*faceSizeTet;
372e6ccafaeSMatthew G Knepley       PetscInt       fVertex, cVertex;
373e6ccafaeSMatthew G Knepley 
374e6ccafaeSMatthew G Knepley       if ((sortedIndices[0] == faceVerticesTetSorted[ii+0]) &&
375e6ccafaeSMatthew G Knepley           (sortedIndices[1] == faceVerticesTetSorted[ii+1]) &&
376e6ccafaeSMatthew G Knepley           (sortedIndices[2] == faceVerticesTetSorted[ii+2]) &&
377e6ccafaeSMatthew G Knepley           (sortedIndices[3] == faceVerticesTetSorted[ii+3])) {
378e6ccafaeSMatthew G Knepley         for (fVertex = 0; fVertex < faceSizeTet; ++fVertex) {
379e6ccafaeSMatthew G Knepley           for (cVertex = 0; cVertex < faceSizeTet; ++cVertex) {
380e6ccafaeSMatthew G Knepley             if (indices[cVertex] == faceVerticesTet[ii+fVertex]) {
381e6ccafaeSMatthew G Knepley               faceVertices[fVertex] = origVertices[cVertex];
382e6ccafaeSMatthew G Knepley               break;
383e6ccafaeSMatthew G Knepley             }
384e6ccafaeSMatthew G Knepley           }
385e6ccafaeSMatthew G Knepley         }
386e6ccafaeSMatthew G Knepley         found = PETSC_TRUE;
387e6ccafaeSMatthew G Knepley         break;
388e6ccafaeSMatthew G Knepley       }
389e6ccafaeSMatthew G Knepley     }
390e6ccafaeSMatthew G Knepley     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tet crossface");
391e6ccafaeSMatthew G Knepley     if (posOriented) *posOriented = PETSC_TRUE;
392e6ccafaeSMatthew G Knepley     PetscFunctionReturn(0);
393e6ccafaeSMatthew G Knepley   } else if (cellDim == 3 && numCorners == 27) {
394e6ccafaeSMatthew G Knepley     /* Quadratic hexes (I hate this)
395e6ccafaeSMatthew G Knepley        A hex is two oriented quads with the normal of the first
396e6ccafaeSMatthew G Knepley        pointing up at the second.
397e6ccafaeSMatthew G Knepley 
398e6ccafaeSMatthew G Knepley          7---6
399e6ccafaeSMatthew G Knepley         /|  /|
400e6ccafaeSMatthew G Knepley        4---5 |
401e6ccafaeSMatthew G Knepley        | 3-|-2
402e6ccafaeSMatthew G Knepley        |/  |/
403e6ccafaeSMatthew G Knepley        0---1
404e6ccafaeSMatthew G Knepley 
405e6ccafaeSMatthew G Knepley        Faces are determined by the first 4 vertices (corners of faces) */
406e6ccafaeSMatthew G Knepley     const PetscInt faceSizeQuadHex = 9;
407e6ccafaeSMatthew G Knepley     PetscInt       sortedIndices[9], i, iFace;
408e6ccafaeSMatthew G Knepley     PetscBool      found                         = PETSC_FALSE;
409e6ccafaeSMatthew G Knepley     PetscInt       faceVerticesQuadHexSorted[54] = {
410e6ccafaeSMatthew G Knepley       0, 1, 2, 3,  8, 9, 10, 11,  24, /* bottom */
411e6ccafaeSMatthew G Knepley       4, 5, 6, 7,  12, 13, 14, 15,  25, /* top */
412e6ccafaeSMatthew G Knepley       0, 1, 4, 5,  8, 12, 16, 17,  22, /* front */
413e6ccafaeSMatthew G Knepley       1, 2, 5, 6,  9, 13, 17, 18,  21, /* right */
414e6ccafaeSMatthew G Knepley       2, 3, 6, 7,  10, 14, 18, 19,  23, /* back */
415e6ccafaeSMatthew G Knepley       0, 3, 4, 7,  11, 15, 16, 19,  20, /* left */
416e6ccafaeSMatthew G Knepley     };
417e6ccafaeSMatthew G Knepley     PetscInt       faceVerticesQuadHex[54] = {
418e6ccafaeSMatthew G Knepley       3, 2, 1, 0,  10, 9, 8, 11,  24, /* bottom */
419e6ccafaeSMatthew G Knepley       4, 5, 6, 7,  12, 13, 14, 15,  25, /* top */
420e6ccafaeSMatthew G Knepley       0, 1, 5, 4,  8, 17, 12, 16,  22, /* front */
421e6ccafaeSMatthew G Knepley       1, 2, 6, 5,  9, 18, 13, 17,  21, /* right */
422e6ccafaeSMatthew G Knepley       2, 3, 7, 6,  10, 19, 14, 18,  23, /* back */
423e6ccafaeSMatthew G Knepley       3, 0, 4, 7,  11, 16, 15, 19,  20 /* left */
424e6ccafaeSMatthew G Knepley     };
425e6ccafaeSMatthew G Knepley 
426e6ccafaeSMatthew G Knepley     faceSize = faceSizeQuadHex;
427e6ccafaeSMatthew G Knepley     for (i = 0; i < faceSizeQuadHex; ++i) sortedIndices[i] = indices[i];
428e6ccafaeSMatthew G Knepley     ierr = PetscSortInt(faceSizeQuadHex, sortedIndices);CHKERRQ(ierr);
429e6ccafaeSMatthew G Knepley     for (iFace = 0; iFace < 6; ++iFace) {
430e6ccafaeSMatthew G Knepley       const PetscInt ii = iFace*faceSizeQuadHex;
431e6ccafaeSMatthew G Knepley       PetscInt       fVertex, cVertex;
432e6ccafaeSMatthew G Knepley 
433e6ccafaeSMatthew G Knepley       if ((sortedIndices[0] == faceVerticesQuadHexSorted[ii+0]) &&
434e6ccafaeSMatthew G Knepley           (sortedIndices[1] == faceVerticesQuadHexSorted[ii+1]) &&
435e6ccafaeSMatthew G Knepley           (sortedIndices[2] == faceVerticesQuadHexSorted[ii+2]) &&
436e6ccafaeSMatthew G Knepley           (sortedIndices[3] == faceVerticesQuadHexSorted[ii+3])) {
437e6ccafaeSMatthew G Knepley         for (fVertex = 0; fVertex < faceSizeQuadHex; ++fVertex) {
438e6ccafaeSMatthew G Knepley           for (cVertex = 0; cVertex < faceSizeQuadHex; ++cVertex) {
439e6ccafaeSMatthew G Knepley             if (indices[cVertex] == faceVerticesQuadHex[ii+fVertex]) {
440e6ccafaeSMatthew G Knepley               faceVertices[fVertex] = origVertices[cVertex];
441e6ccafaeSMatthew G Knepley               break;
442e6ccafaeSMatthew G Knepley             }
443e6ccafaeSMatthew G Knepley           }
444e6ccafaeSMatthew G Knepley         }
445e6ccafaeSMatthew G Knepley         found = PETSC_TRUE;
446e6ccafaeSMatthew G Knepley         break;
447e6ccafaeSMatthew G Knepley       }
448e6ccafaeSMatthew G Knepley     }
449e6ccafaeSMatthew G Knepley     if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface");
450e6ccafaeSMatthew G Knepley     if (posOriented) *posOriented = PETSC_TRUE;
451e6ccafaeSMatthew G Knepley     PetscFunctionReturn(0);
452e6ccafaeSMatthew G Knepley   } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Unknown cell type for faceOrientation().");
453e6ccafaeSMatthew G Knepley   if (!posOrient) {
454e6ccafaeSMatthew G Knepley     if (debug) {ierr = PetscPrintf(comm, "  Reversing initial face orientation\n");CHKERRQ(ierr);}
455e6ccafaeSMatthew G Knepley     for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[faceSize-1 - f];
456e6ccafaeSMatthew G Knepley   } else {
457e6ccafaeSMatthew G Knepley     if (debug) {ierr = PetscPrintf(comm, "  Keeping initial face orientation\n");CHKERRQ(ierr);}
458e6ccafaeSMatthew G Knepley     for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[f];
459e6ccafaeSMatthew G Knepley   }
460e6ccafaeSMatthew G Knepley   if (posOriented) *posOriented = posOrient;
461e6ccafaeSMatthew G Knepley   PetscFunctionReturn(0);
462e6ccafaeSMatthew G Knepley }
463e6ccafaeSMatthew G Knepley 
464e6ccafaeSMatthew G Knepley #undef __FUNCT__
465e6ccafaeSMatthew G Knepley #define __FUNCT__ "DMPlexGetOrientedFace"
466e6ccafaeSMatthew G Knepley /*
467e6ccafaeSMatthew G Knepley     Given a cell and a face, as a set of vertices,
468e6ccafaeSMatthew G Knepley       return the oriented face, as a set of vertices, in faceVertices
469e6ccafaeSMatthew G Knepley     The orientation is such that the face normal points out of the cell
470e6ccafaeSMatthew G Knepley */
471e6ccafaeSMatthew G Knepley PetscErrorCode DMPlexGetOrientedFace(DM dm, PetscInt cell, PetscInt faceSize, const PetscInt face[], PetscInt numCorners, PetscInt indices[], PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented)
472e6ccafaeSMatthew G Knepley {
473e6ccafaeSMatthew G Knepley   const PetscInt *cone = PETSC_NULL;
474e6ccafaeSMatthew G Knepley   PetscInt        coneSize, v, f, v2;
475e6ccafaeSMatthew G Knepley   PetscInt        oppositeVertex = -1;
476e6ccafaeSMatthew G Knepley   PetscErrorCode  ierr;
477e6ccafaeSMatthew G Knepley 
478e6ccafaeSMatthew G Knepley   PetscFunctionBegin;
479e6ccafaeSMatthew G Knepley   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
480e6ccafaeSMatthew G Knepley   ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
481e6ccafaeSMatthew G Knepley   for (v = 0, v2 = 0; v < coneSize; ++v) {
482e6ccafaeSMatthew G Knepley     PetscBool found = PETSC_FALSE;
483e6ccafaeSMatthew G Knepley 
484e6ccafaeSMatthew G Knepley     for (f = 0; f < faceSize; ++f) {
485e6ccafaeSMatthew G Knepley       if (face[f] == cone[v]) {
486e6ccafaeSMatthew G Knepley         found = PETSC_TRUE; break;
487e6ccafaeSMatthew G Knepley       }
488e6ccafaeSMatthew G Knepley     }
489e6ccafaeSMatthew G Knepley     if (found) {
490e6ccafaeSMatthew G Knepley       indices[v2]      = v;
491e6ccafaeSMatthew G Knepley       origVertices[v2] = cone[v];
492e6ccafaeSMatthew G Knepley       ++v2;
493e6ccafaeSMatthew G Knepley     } else {
494e6ccafaeSMatthew G Knepley       oppositeVertex = v;
495e6ccafaeSMatthew G Knepley     }
496e6ccafaeSMatthew G Knepley   }
497e6ccafaeSMatthew G Knepley   ierr = DMPlexGetFaceOrientation(dm, cell, numCorners, indices, oppositeVertex, origVertices, faceVertices, posOriented);CHKERRQ(ierr);
498e6ccafaeSMatthew G Knepley   PetscFunctionReturn(0);
499e6ccafaeSMatthew G Knepley }
500e6ccafaeSMatthew G Knepley 
501e6ccafaeSMatthew G Knepley #undef __FUNCT__
502e6ccafaeSMatthew G Knepley #define __FUNCT__ "DMPlexInsertFace_Private"
503e6ccafaeSMatthew G Knepley /*
504e6ccafaeSMatthew G Knepley   DMPlexInsertFace_Private - Puts a face into the mesh
505e6ccafaeSMatthew G Knepley 
506e6ccafaeSMatthew G Knepley   Not collective
507e6ccafaeSMatthew G Knepley 
508e6ccafaeSMatthew G Knepley   Input Parameters:
509e6ccafaeSMatthew G Knepley   + dm              - The DMPlex
510e6ccafaeSMatthew G Knepley   . numFaceVertex   - The number of vertices in the face
511e6ccafaeSMatthew G Knepley   . faceVertices    - The vertices in the face for dm
512e6ccafaeSMatthew G Knepley   . subfaceVertices - The vertices in the face for subdm
513e6ccafaeSMatthew G Knepley   . numCorners      - The number of vertices in the cell
514e6ccafaeSMatthew G Knepley   . cell            - A cell in dm containing the face
515e6ccafaeSMatthew G Knepley   . subcell         - A cell in subdm containing the face
516e6ccafaeSMatthew G Knepley   . firstFace       - First face in the mesh
517e6ccafaeSMatthew G Knepley   - newFacePoint    - Next face in the mesh
518e6ccafaeSMatthew G Knepley 
519e6ccafaeSMatthew G Knepley   Output Parameters:
520e6ccafaeSMatthew G Knepley   . newFacePoint - Contains next face point number on input, updated on output
521e6ccafaeSMatthew G Knepley 
522e6ccafaeSMatthew G Knepley   Level: developer
523e6ccafaeSMatthew G Knepley */
524e6ccafaeSMatthew G Knepley PetscErrorCode DMPlexInsertFace_Private(DM dm, DM subdm, PetscInt numFaceVertices, const PetscInt faceVertices[], const PetscInt subfaceVertices[], PetscInt numCorners, PetscInt cell, PetscInt subcell, PetscInt firstFace, PetscInt *newFacePoint)
525e6ccafaeSMatthew G Knepley {
526e6ccafaeSMatthew G Knepley   MPI_Comm        comm    = ((PetscObject) dm)->comm;
527e6ccafaeSMatthew G Knepley   DM_Plex        *submesh = (DM_Plex*) subdm->data;
528e6ccafaeSMatthew G Knepley   const PetscInt *faces;
529e6ccafaeSMatthew G Knepley   PetscInt        numFaces, coneSize;
530e6ccafaeSMatthew G Knepley   PetscErrorCode  ierr;
531e6ccafaeSMatthew G Knepley 
532e6ccafaeSMatthew G Knepley   PetscFunctionBegin;
533e6ccafaeSMatthew G Knepley   ierr = DMPlexGetConeSize(subdm, subcell, &coneSize);CHKERRQ(ierr);
534e6ccafaeSMatthew G Knepley   if (coneSize != 1) SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cone size of cell %d is %d != 1", cell, coneSize);
535e6ccafaeSMatthew G Knepley #if 0
536e6ccafaeSMatthew G Knepley   /* Cannot use this because support() has not been constructed yet */
537e6ccafaeSMatthew G Knepley   ierr = DMPlexGetJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr);
538e6ccafaeSMatthew G Knepley #else
539e6ccafaeSMatthew G Knepley   {
540e6ccafaeSMatthew G Knepley     PetscInt f;
541e6ccafaeSMatthew G Knepley 
542e6ccafaeSMatthew G Knepley     numFaces = 0;
543e6ccafaeSMatthew G Knepley     ierr     = DMGetWorkArray(subdm, 1, PETSC_INT, (void**) &faces);CHKERRQ(ierr);
544e6ccafaeSMatthew G Knepley     for (f = firstFace; f < *newFacePoint; ++f) {
545e6ccafaeSMatthew G Knepley       PetscInt dof, off, d;
546e6ccafaeSMatthew G Knepley 
547e6ccafaeSMatthew G Knepley       ierr = PetscSectionGetDof(submesh->coneSection, f, &dof);CHKERRQ(ierr);
548e6ccafaeSMatthew G Knepley       ierr = PetscSectionGetOffset(submesh->coneSection, f, &off);CHKERRQ(ierr);
549e6ccafaeSMatthew G Knepley       /* Yes, I know this is quadratic, but I expect the sizes to be <5 */
550e6ccafaeSMatthew G Knepley       for (d = 0; d < dof; ++d) {
551e6ccafaeSMatthew G Knepley         const PetscInt p = submesh->cones[off+d];
552e6ccafaeSMatthew G Knepley         PetscInt       v;
553e6ccafaeSMatthew G Knepley 
554e6ccafaeSMatthew G Knepley         for (v = 0; v < numFaceVertices; ++v) {
555e6ccafaeSMatthew G Knepley           if (subfaceVertices[v] == p) break;
556e6ccafaeSMatthew G Knepley         }
557e6ccafaeSMatthew G Knepley         if (v == numFaceVertices) break;
558e6ccafaeSMatthew G Knepley       }
559e6ccafaeSMatthew G Knepley       if (d == dof) {
560e6ccafaeSMatthew G Knepley         numFaces               = 1;
561e6ccafaeSMatthew G Knepley         ((PetscInt*) faces)[0] = f;
562e6ccafaeSMatthew G Knepley       }
563e6ccafaeSMatthew G Knepley     }
564e6ccafaeSMatthew G Knepley   }
565e6ccafaeSMatthew G Knepley #endif
566e6ccafaeSMatthew G Knepley   if (numFaces > 1) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Vertex set had %d faces, not one", numFaces);
567e6ccafaeSMatthew G Knepley   else if (numFaces == 1) {
568e6ccafaeSMatthew G Knepley     /* Add the other cell neighbor for this face */
569e6ccafaeSMatthew G Knepley     ierr = DMPlexSetCone(subdm, cell, faces);CHKERRQ(ierr);
570e6ccafaeSMatthew G Knepley   } else {
571e6ccafaeSMatthew G Knepley     PetscInt *indices, *origVertices, *orientedVertices, *orientedSubVertices, v, ov;
572e6ccafaeSMatthew G Knepley     PetscBool posOriented;
573e6ccafaeSMatthew G Knepley 
574e6ccafaeSMatthew G Knepley     ierr                = DMGetWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr);
575e6ccafaeSMatthew G Knepley     origVertices        = &orientedVertices[numFaceVertices];
576e6ccafaeSMatthew G Knepley     indices             = &orientedVertices[numFaceVertices*2];
577e6ccafaeSMatthew G Knepley     orientedSubVertices = &orientedVertices[numFaceVertices*3];
578e6ccafaeSMatthew G Knepley     ierr                = DMPlexGetOrientedFace(dm, cell, numFaceVertices, faceVertices, numCorners, indices, origVertices, orientedVertices, &posOriented);CHKERRQ(ierr);
579e6ccafaeSMatthew G Knepley     /* TODO: I know that routine should return a permutation, not the indices */
580e6ccafaeSMatthew G Knepley     for (v = 0; v < numFaceVertices; ++v) {
581e6ccafaeSMatthew G Knepley       const PetscInt vertex = faceVertices[v], subvertex = subfaceVertices[v];
582e6ccafaeSMatthew G Knepley       for (ov = 0; ov < numFaceVertices; ++ov) {
583e6ccafaeSMatthew G Knepley         if (orientedVertices[ov] == vertex) {
584e6ccafaeSMatthew G Knepley           orientedSubVertices[ov] = subvertex;
585e6ccafaeSMatthew G Knepley           break;
586e6ccafaeSMatthew G Knepley         }
587e6ccafaeSMatthew G Knepley       }
588e6ccafaeSMatthew G Knepley       if (ov == numFaceVertices) SETERRQ1(comm, PETSC_ERR_PLIB, "Could not find face vertex %d in orientated set", vertex);
589e6ccafaeSMatthew G Knepley     }
590e6ccafaeSMatthew G Knepley     ierr = DMPlexSetCone(subdm, *newFacePoint, orientedSubVertices);CHKERRQ(ierr);
591e6ccafaeSMatthew G Knepley     ierr = DMPlexSetCone(subdm, subcell, newFacePoint);CHKERRQ(ierr);
592e6ccafaeSMatthew G Knepley     ierr = DMRestoreWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr);
593e6ccafaeSMatthew G Knepley     ++(*newFacePoint);
594e6ccafaeSMatthew G Knepley   }
595e6ccafaeSMatthew G Knepley   ierr = DMPlexRestoreJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr);
596e6ccafaeSMatthew G Knepley   PetscFunctionReturn(0);
597e6ccafaeSMatthew G Knepley }
598e6ccafaeSMatthew G Knepley 
599e6ccafaeSMatthew G Knepley #undef __FUNCT__
600e6ccafaeSMatthew G Knepley #define __FUNCT__ "DMPlexCreateSubmesh_Uninterpolated"
601*efa14ee0SMatthew G Knepley static PetscErrorCode DMPlexCreateSubmesh_Uninterpolated(DM dm, const char vertexLabelName[], DM subdm)
602e6ccafaeSMatthew G Knepley {
603e6ccafaeSMatthew G Knepley   MPI_Comm        comm = ((PetscObject) dm)->comm;
604*efa14ee0SMatthew G Knepley   DMLabel         vertexLabel, subpointMap;
605*efa14ee0SMatthew G Knepley   IS              subvertexIS,  subcellIS;
606*efa14ee0SMatthew G Knepley   const PetscInt *subVertices, *subCells;
607*efa14ee0SMatthew G Knepley   PetscInt        numSubVertices, firstSubVertex, numSubCells;
608*efa14ee0SMatthew G Knepley   PetscInt       *subface, maxConeSize, numSubFaces, firstSubFace, newFacePoint, nFV;
609*efa14ee0SMatthew G Knepley   PetscInt        vStart, vEnd, c, f;
610e6ccafaeSMatthew G Knepley   PetscErrorCode  ierr;
611e6ccafaeSMatthew G Knepley 
612e6ccafaeSMatthew G Knepley   PetscFunctionBegin;
613*efa14ee0SMatthew G Knepley   /* Create subpointMap which marks the submesh */
614*efa14ee0SMatthew G Knepley   ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr);
615*efa14ee0SMatthew G Knepley   ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr);
616*efa14ee0SMatthew G Knepley   ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr);
617*efa14ee0SMatthew G Knepley   ierr = DMPlexGetLabel(dm, vertexLabelName, &vertexLabel);CHKERRQ(ierr);
618*efa14ee0SMatthew G Knepley   ierr = DMPlexMarkSubmesh_Uninterpolated(dm, vertexLabel, subpointMap, &numSubFaces, &nFV, subdm);CHKERRQ(ierr);
619*efa14ee0SMatthew G Knepley   /* Setup chart */
620*efa14ee0SMatthew G Knepley   ierr = DMLabelGetStratumSize(subpointMap, 0, &numSubVertices);CHKERRQ(ierr);
621*efa14ee0SMatthew G Knepley   ierr = DMLabelGetStratumSize(subpointMap, 2, &numSubCells);CHKERRQ(ierr);
622*efa14ee0SMatthew G Knepley   ierr = DMPlexSetChart(subdm, 0, numSubCells+numSubFaces+numSubVertices);CHKERRQ(ierr);
623*efa14ee0SMatthew G Knepley   ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr);
624e6ccafaeSMatthew G Knepley   /* Set cone sizes */
625e6ccafaeSMatthew G Knepley   firstSubVertex = numSubCells;
626*efa14ee0SMatthew G Knepley   firstSubFace   = numSubCells+numSubVertices;
627e6ccafaeSMatthew G Knepley   newFacePoint   = firstSubFace;
628*efa14ee0SMatthew G Knepley   ierr = DMLabelGetStratumIS(subpointMap, 0, &subvertexIS);CHKERRQ(ierr);
629*efa14ee0SMatthew G Knepley   if (subvertexIS) {ierr = ISGetIndices(subvertexIS, &subVertices);CHKERRQ(ierr);}
630*efa14ee0SMatthew G Knepley   ierr = DMLabelGetStratumIS(subpointMap, 2, &subcellIS);CHKERRQ(ierr);
631*efa14ee0SMatthew G Knepley   if (subcellIS) {ierr = ISGetIndices(subcellIS, &subCells);CHKERRQ(ierr);}
632e6ccafaeSMatthew G Knepley   for (c = 0; c < numSubCells; ++c) {
633e6ccafaeSMatthew G Knepley     ierr = DMPlexSetConeSize(subdm, c, 1);CHKERRQ(ierr);
634e6ccafaeSMatthew G Knepley   }
635e6ccafaeSMatthew G Knepley   for (f = firstSubFace; f < firstSubFace+numSubFaces; ++f) {
636e6ccafaeSMatthew G Knepley     ierr = DMPlexSetConeSize(subdm, f, nFV);CHKERRQ(ierr);
637e6ccafaeSMatthew G Knepley   }
638e6ccafaeSMatthew G Knepley   ierr = DMSetUp(subdm);CHKERRQ(ierr);
639e6ccafaeSMatthew G Knepley   /* Create face cones */
640*efa14ee0SMatthew G Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
641*efa14ee0SMatthew G Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, PETSC_NULL);CHKERRQ(ierr);
642*efa14ee0SMatthew G Knepley   ierr = DMGetWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr);
643e6ccafaeSMatthew G Knepley   for (c = 0; c < numSubCells; ++c) {
644e6ccafaeSMatthew G Knepley     const PetscInt cell    = subCells[c];
645*efa14ee0SMatthew G Knepley     const PetscInt subcell = c;
646e6ccafaeSMatthew G Knepley     PetscInt      *closure = PETSC_NULL;
647*efa14ee0SMatthew G Knepley     PetscInt       closureSize, cl, numCorners = 0, faceSize = 0;
648e6ccafaeSMatthew G Knepley 
649e6ccafaeSMatthew G Knepley     ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
650*efa14ee0SMatthew G Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
651*efa14ee0SMatthew G Knepley       const PetscInt point = closure[cl];
652e6ccafaeSMatthew G Knepley       PetscInt       subVertex;
653e6ccafaeSMatthew G Knepley 
654*efa14ee0SMatthew G Knepley       if ((point >= vStart) && (point < vEnd)) {
655*efa14ee0SMatthew G Knepley         ++numCorners;
656*efa14ee0SMatthew G Knepley         ierr = PetscFindInt(point, numSubVertices, subVertices, &subVertex);CHKERRQ(ierr);
657*efa14ee0SMatthew G Knepley         if (subVertex >= 0) {
658*efa14ee0SMatthew G Knepley           closure[faceSize] = point;
659*efa14ee0SMatthew G Knepley           subface[faceSize] = subVertex;
660e6ccafaeSMatthew G Knepley           ++faceSize;
661e6ccafaeSMatthew G Knepley         }
662e6ccafaeSMatthew G Knepley       }
663e6ccafaeSMatthew G Knepley     }
664e6ccafaeSMatthew G Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
665*efa14ee0SMatthew G Knepley     if (faceSize > nFV) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize);
666*efa14ee0SMatthew G Knepley     if (faceSize == nFV) {
667*efa14ee0SMatthew G Knepley       ierr = DMPlexInsertFace_Private(dm, subdm, faceSize, closure, subface, numCorners, cell, subcell, firstSubFace, &newFacePoint);CHKERRQ(ierr);
668e6ccafaeSMatthew G Knepley     }
669e6ccafaeSMatthew G Knepley   }
670*efa14ee0SMatthew G Knepley   ierr = DMRestoreWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr);
671e6ccafaeSMatthew G Knepley   ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr);
672e6ccafaeSMatthew G Knepley   ierr = DMPlexStratify(subdm);CHKERRQ(ierr);
673e6ccafaeSMatthew G Knepley   /* Build coordinates */
674*efa14ee0SMatthew G Knepley   {
675*efa14ee0SMatthew G Knepley     PetscSection coordSection, subCoordSection;
676*efa14ee0SMatthew G Knepley     Vec          coordinates, subCoordinates;
677*efa14ee0SMatthew G Knepley     PetscScalar *coords, *subCoords;
678*efa14ee0SMatthew G Knepley     PetscInt     coordSize, v;
679*efa14ee0SMatthew G Knepley 
680e6ccafaeSMatthew G Knepley     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
681e6ccafaeSMatthew G Knepley     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
682e6ccafaeSMatthew G Knepley     ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr);
683*efa14ee0SMatthew G Knepley     ierr = PetscSectionSetChart(subCoordSection, firstSubVertex, firstSubVertex+numSubVertices);CHKERRQ(ierr);
684*efa14ee0SMatthew G Knepley     for (v = 0; v < numSubVertices; ++v) {
685*efa14ee0SMatthew G Knepley       const PetscInt vertex    = subVertices[v];
686*efa14ee0SMatthew G Knepley       const PetscInt subvertex = firstSubVertex+v;
687*efa14ee0SMatthew G Knepley       PetscInt       dof;
688*efa14ee0SMatthew G Knepley 
689*efa14ee0SMatthew G Knepley       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
690*efa14ee0SMatthew G Knepley       ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr);
691e6ccafaeSMatthew G Knepley     }
692e6ccafaeSMatthew G Knepley     ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr);
693e6ccafaeSMatthew G Knepley     ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr);
694e6ccafaeSMatthew G Knepley     ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr);
695e6ccafaeSMatthew G Knepley     ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
696e6ccafaeSMatthew G Knepley     ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr);
697e6ccafaeSMatthew G Knepley     ierr = VecGetArray(coordinates,    &coords);CHKERRQ(ierr);
698e6ccafaeSMatthew G Knepley     ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr);
699*efa14ee0SMatthew G Knepley     for (v = 0; v < numSubVertices; ++v) {
700*efa14ee0SMatthew G Knepley       const PetscInt vertex    = subVertices[v];
701*efa14ee0SMatthew G Knepley       const PetscInt subvertex = firstSubVertex+v;
702*efa14ee0SMatthew G Knepley       PetscInt       dof, off, sdof, soff, d;
703e6ccafaeSMatthew G Knepley 
704e6ccafaeSMatthew G Knepley       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
705e6ccafaeSMatthew G Knepley       ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr);
706*efa14ee0SMatthew G Knepley       ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr);
707*efa14ee0SMatthew G Knepley       ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr);
708*efa14ee0SMatthew G Knepley       if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof);
709e6ccafaeSMatthew G Knepley       for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d];
710e6ccafaeSMatthew G Knepley     }
711e6ccafaeSMatthew G Knepley     ierr = VecRestoreArray(coordinates,    &coords);CHKERRQ(ierr);
712e6ccafaeSMatthew G Knepley     ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr);
713e6ccafaeSMatthew G Knepley     ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr);
714e6ccafaeSMatthew G Knepley     ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr);
715e6ccafaeSMatthew G Knepley   }
716*efa14ee0SMatthew G Knepley   /* Cleanup */
717*efa14ee0SMatthew G Knepley   if (subvertexIS) {ierr = ISRestoreIndices(subvertexIS, &subVertices);CHKERRQ(ierr);}
718*efa14ee0SMatthew G Knepley   ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr);
719*efa14ee0SMatthew G Knepley   if (subcellIS) {ierr = ISRestoreIndices(subcellIS, &subCells);CHKERRQ(ierr);}
720*efa14ee0SMatthew G Knepley   ierr = ISDestroy(&subcellIS);CHKERRQ(ierr);
721e6ccafaeSMatthew G Knepley   PetscFunctionReturn(0);
722e6ccafaeSMatthew G Knepley }
723e6ccafaeSMatthew G Knepley 
724e6ccafaeSMatthew G Knepley #undef __FUNCT__
725e6ccafaeSMatthew G Knepley #define __FUNCT__ "DMPlexCreateSubmesh_Interpolated"
726e6ccafaeSMatthew G Knepley static PetscErrorCode DMPlexCreateSubmesh_Interpolated(DM dm, const char vertexLabelName[], DM subdm)
727e6ccafaeSMatthew G Knepley {
728e6ccafaeSMatthew G Knepley   MPI_Comm         comm = ((PetscObject) dm)->comm;
729e6ccafaeSMatthew G Knepley   DMLabel          subpointMap, vertexLabel;
730*efa14ee0SMatthew G Knepley   IS              *subpointIS;
731*efa14ee0SMatthew G Knepley   const PetscInt **subpoints;
732*efa14ee0SMatthew G Knepley   PetscInt        *numSubPoints, *firstSubPoint, *coneNew;
733*efa14ee0SMatthew G Knepley   PetscInt         totSubPoints = 0, maxConeSize, dim, p, d, v;
734e6ccafaeSMatthew G Knepley   PetscErrorCode   ierr;
735e6ccafaeSMatthew G Knepley 
736e6ccafaeSMatthew G Knepley   PetscFunctionBegin;
737*efa14ee0SMatthew G Knepley   /* Create subpointMap which marks the submesh */
738e6ccafaeSMatthew G Knepley   ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr);
739*efa14ee0SMatthew G Knepley   ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr);
740*efa14ee0SMatthew G Knepley   ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr);
741e6ccafaeSMatthew G Knepley   ierr = DMPlexGetLabel(dm, vertexLabelName, &vertexLabel);CHKERRQ(ierr);
742*efa14ee0SMatthew G Knepley   ierr = DMPlexMarkSubmesh_Interpolated(dm, vertexLabel, subpointMap, subdm);CHKERRQ(ierr);
743*efa14ee0SMatthew G Knepley   /* Setup chart */
744*efa14ee0SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
745*efa14ee0SMatthew G Knepley   ierr = PetscMalloc4(dim+1,PetscInt,&numSubPoints,dim+1,PetscInt,&firstSubPoint,dim+1,IS,&subpointIS,dim+1,const PetscInt *,&subpoints);CHKERRQ(ierr);
746e6ccafaeSMatthew G Knepley   for (d = 0; d <= dim; ++d) {
747e6ccafaeSMatthew G Knepley     ierr = DMLabelGetStratumSize(subpointMap, d, &numSubPoints[d]);CHKERRQ(ierr);
748e6ccafaeSMatthew G Knepley     totSubPoints += numSubPoints[d];
749e6ccafaeSMatthew G Knepley   }
750e6ccafaeSMatthew G Knepley   ierr = DMPlexSetChart(subdm, 0, totSubPoints);CHKERRQ(ierr);
751e6ccafaeSMatthew G Knepley   ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr);
752e6ccafaeSMatthew G Knepley   /* Set cone sizes */
753e6ccafaeSMatthew G Knepley   firstSubPoint[dim] = 0;
754e6ccafaeSMatthew G Knepley   firstSubPoint[0]   = firstSubPoint[dim] + numSubPoints[dim];
755e6ccafaeSMatthew G Knepley   if (dim > 1) {firstSubPoint[dim-1] = firstSubPoint[0]     + numSubPoints[0];}
756e6ccafaeSMatthew G Knepley   if (dim > 2) {firstSubPoint[dim-2] = firstSubPoint[dim-1] + numSubPoints[dim-1];}
757e6ccafaeSMatthew G Knepley   for (d = 0; d <= dim; ++d) {
758e6ccafaeSMatthew G Knepley     ierr = DMLabelGetStratumIS(subpointMap, d, &subpointIS[d]);CHKERRQ(ierr);
759e6ccafaeSMatthew G Knepley     ierr = ISGetIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr);
760e6ccafaeSMatthew G Knepley   }
761e6ccafaeSMatthew G Knepley   for (d = 0; d <= dim; ++d) {
762e6ccafaeSMatthew G Knepley     for (p = 0; p < numSubPoints[d]; ++p) {
763e6ccafaeSMatthew G Knepley       const PetscInt  point    = subpoints[d][p];
764e6ccafaeSMatthew G Knepley       const PetscInt  subpoint = firstSubPoint[d] + p;
765e6ccafaeSMatthew G Knepley       const PetscInt *cone;
766e6ccafaeSMatthew G Knepley       PetscInt        coneSize, coneSizeNew, c, val;
767e6ccafaeSMatthew G Knepley 
768e6ccafaeSMatthew G Knepley       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
769e6ccafaeSMatthew G Knepley       ierr = DMPlexSetConeSize(subdm, subpoint, coneSize);CHKERRQ(ierr);
770e6ccafaeSMatthew G Knepley       if (d == dim) {
771e6ccafaeSMatthew G Knepley         ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
772e6ccafaeSMatthew G Knepley         for (c = 0, coneSizeNew = 0; c < coneSize; ++c) {
773e6ccafaeSMatthew G Knepley           ierr = DMLabelGetValue(subpointMap, cone[c], &val);CHKERRQ(ierr);
774e6ccafaeSMatthew G Knepley           if (val >= 0) coneSizeNew++;
775e6ccafaeSMatthew G Knepley         }
776e6ccafaeSMatthew G Knepley         ierr = DMPlexSetConeSize(subdm, subpoint, coneSizeNew);CHKERRQ(ierr);
777e6ccafaeSMatthew G Knepley       }
778e6ccafaeSMatthew G Knepley     }
779e6ccafaeSMatthew G Knepley   }
780e6ccafaeSMatthew G Knepley   ierr = DMSetUp(subdm);CHKERRQ(ierr);
781e6ccafaeSMatthew G Knepley   /* Set cones */
782e6ccafaeSMatthew G Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, PETSC_NULL);CHKERRQ(ierr);
783e6ccafaeSMatthew G Knepley   ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &coneNew);CHKERRQ(ierr);
784e6ccafaeSMatthew G Knepley   for (d = 0; d <= dim; ++d) {
785e6ccafaeSMatthew G Knepley     for (p = 0; p < numSubPoints[d]; ++p) {
786e6ccafaeSMatthew G Knepley       const PetscInt  point    = subpoints[d][p];
787e6ccafaeSMatthew G Knepley       const PetscInt  subpoint = firstSubPoint[d] + p;
788e6ccafaeSMatthew G Knepley       const PetscInt *cone;
789e6ccafaeSMatthew G Knepley       PetscInt        coneSize, subconeSize, coneSizeNew, c, subc;
790e6ccafaeSMatthew G Knepley 
791e6ccafaeSMatthew G Knepley       ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
792e6ccafaeSMatthew G Knepley       ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr);
793e6ccafaeSMatthew G Knepley       ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
794e6ccafaeSMatthew G Knepley       for (c = 0, coneSizeNew = 0; c < coneSize; ++c) {
795e6ccafaeSMatthew G Knepley         ierr = PetscFindInt(cone[c], numSubPoints[d-1], subpoints[d-1], &subc);CHKERRQ(ierr);
796e6ccafaeSMatthew G Knepley         if (subc >= 0) coneNew[coneSizeNew++] = firstSubPoint[d-1] + subc;
797e6ccafaeSMatthew G Knepley       }
798e6ccafaeSMatthew G Knepley       if (coneSizeNew != subconeSize) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of cone points located %d does not match subcone size %d", coneSizeNew, subconeSize);
799e6ccafaeSMatthew G Knepley       ierr = DMPlexSetCone(subdm, subpoint, coneNew);CHKERRQ(ierr);
800e6ccafaeSMatthew G Knepley     }
801e6ccafaeSMatthew G Knepley   }
802e6ccafaeSMatthew G Knepley   ierr = PetscFree(coneNew);CHKERRQ(ierr);
803e6ccafaeSMatthew G Knepley   ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr);
804e6ccafaeSMatthew G Knepley   ierr = DMPlexStratify(subdm);CHKERRQ(ierr);
805e6ccafaeSMatthew G Knepley   /* Build coordinates */
806e6ccafaeSMatthew G Knepley   {
807e6ccafaeSMatthew G Knepley     PetscSection coordSection, subCoordSection;
808e6ccafaeSMatthew G Knepley     Vec          coordinates, subCoordinates;
809e6ccafaeSMatthew G Knepley     PetscScalar *coords, *subCoords;
810e6ccafaeSMatthew G Knepley     PetscInt     coordSize;
811e6ccafaeSMatthew G Knepley 
812e6ccafaeSMatthew G Knepley     ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
813e6ccafaeSMatthew G Knepley     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
814e6ccafaeSMatthew G Knepley     ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr);
815e6ccafaeSMatthew G Knepley     ierr = PetscSectionSetChart(subCoordSection, firstSubPoint[0], firstSubPoint[0]+numSubPoints[0]);CHKERRQ(ierr);
816e6ccafaeSMatthew G Knepley     for (v = 0; v < numSubPoints[0]; ++v) {
817e6ccafaeSMatthew G Knepley       const PetscInt vertex    = subpoints[0][v];
818e6ccafaeSMatthew G Knepley       const PetscInt subvertex = firstSubPoint[0]+v;
819e6ccafaeSMatthew G Knepley       PetscInt       dof;
820e6ccafaeSMatthew G Knepley 
821e6ccafaeSMatthew G Knepley       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
822e6ccafaeSMatthew G Knepley       ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr);
823e6ccafaeSMatthew G Knepley     }
824e6ccafaeSMatthew G Knepley     ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr);
825e6ccafaeSMatthew G Knepley     ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr);
826e6ccafaeSMatthew G Knepley     ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr);
827e6ccafaeSMatthew G Knepley     ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
828e6ccafaeSMatthew G Knepley     ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr);
829e6ccafaeSMatthew G Knepley     ierr = VecGetArray(coordinates,    &coords);CHKERRQ(ierr);
830e6ccafaeSMatthew G Knepley     ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr);
831e6ccafaeSMatthew G Knepley     for (v = 0; v < numSubPoints[0]; ++v) {
832e6ccafaeSMatthew G Knepley       const PetscInt vertex    = subpoints[0][v];
833e6ccafaeSMatthew G Knepley       const PetscInt subvertex = firstSubPoint[0]+v;
834e6ccafaeSMatthew G Knepley       PetscInt dof, off, sdof, soff, d;
835e6ccafaeSMatthew G Knepley 
836e6ccafaeSMatthew G Knepley       ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr);
837e6ccafaeSMatthew G Knepley       ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr);
838e6ccafaeSMatthew G Knepley       ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr);
839e6ccafaeSMatthew G Knepley       ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr);
840e6ccafaeSMatthew G Knepley       if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof);
841*efa14ee0SMatthew G Knepley       for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d];
842e6ccafaeSMatthew G Knepley     }
843e6ccafaeSMatthew G Knepley     ierr = VecRestoreArray(coordinates,    &coords);CHKERRQ(ierr);
844e6ccafaeSMatthew G Knepley     ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr);
845e6ccafaeSMatthew G Knepley     ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr);
846e6ccafaeSMatthew G Knepley     ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr);
847e6ccafaeSMatthew G Knepley   }
848*efa14ee0SMatthew G Knepley   /* Cleanup */
849e6ccafaeSMatthew G Knepley   for (d = 0; d <= dim; ++d) {
850e6ccafaeSMatthew G Knepley     ierr = ISRestoreIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr);
851e6ccafaeSMatthew G Knepley     ierr = ISDestroy(&subpointIS[d]);CHKERRQ(ierr);
852e6ccafaeSMatthew G Knepley   }
853*efa14ee0SMatthew G Knepley   ierr = PetscFree4(numSubPoints,firstSubPoint,subpointIS,subpoints);CHKERRQ(ierr);
854e6ccafaeSMatthew G Knepley   PetscFunctionReturn(0);
855e6ccafaeSMatthew G Knepley }
856e6ccafaeSMatthew G Knepley 
857e6ccafaeSMatthew G Knepley #undef __FUNCT__
858e6ccafaeSMatthew G Knepley #define __FUNCT__ "DMPlexCreateSubmesh"
859e6ccafaeSMatthew G Knepley /*
860e6ccafaeSMatthew G Knepley   DMPlexCreateSubmesh - Extract a hypersurface from the mesh using vertices defined by a label
861e6ccafaeSMatthew G Knepley 
862e6ccafaeSMatthew G Knepley   Input Parameters:
863e6ccafaeSMatthew G Knepley + dm           - The original mesh
864e6ccafaeSMatthew G Knepley - vertexLabel  - The DMLabel marking vertices contained in the surface
865e6ccafaeSMatthew G Knepley 
866e6ccafaeSMatthew G Knepley   Output Parameter:
867e6ccafaeSMatthew G Knepley . subdm - The surface mesh
868e6ccafaeSMatthew G Knepley 
869e6ccafaeSMatthew G Knepley   Note: This function produces a DMLabel mapping original points in the submesh to their depth. This can be obtained using DMPlexGetSubpointMap().
870e6ccafaeSMatthew G Knepley 
871e6ccafaeSMatthew G Knepley   Level: developer
872e6ccafaeSMatthew G Knepley 
873e6ccafaeSMatthew G Knepley .seealso: DMPlexGetSubpointMap(), DMPlexGetLabel(), DMLabelSetValue()
874e6ccafaeSMatthew G Knepley */
875e6ccafaeSMatthew G Knepley PetscErrorCode DMPlexCreateSubmesh(DM dm, const char vertexLabel[], DM *subdm)
876e6ccafaeSMatthew G Knepley {
877e6ccafaeSMatthew G Knepley   PetscInt       dim, depth;
878e6ccafaeSMatthew G Knepley   PetscErrorCode ierr;
879e6ccafaeSMatthew G Knepley 
880e6ccafaeSMatthew G Knepley   PetscFunctionBegin;
881e6ccafaeSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
882e6ccafaeSMatthew G Knepley   PetscValidCharPointer(vertexLabel, 2);
883e6ccafaeSMatthew G Knepley   PetscValidPointer(subdm, 4);
884e6ccafaeSMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
885e6ccafaeSMatthew G Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
886e6ccafaeSMatthew G Knepley   ierr = DMCreate(((PetscObject) dm)->comm, subdm);CHKERRQ(ierr);
887e6ccafaeSMatthew G Knepley   ierr = DMSetType(*subdm, DMPLEX);CHKERRQ(ierr);
888e6ccafaeSMatthew G Knepley   ierr = DMPlexSetDimension(*subdm, dim-1);CHKERRQ(ierr);
889e6ccafaeSMatthew G Knepley   if (depth == dim) {
890e6ccafaeSMatthew G Knepley     ierr = DMPlexCreateSubmesh_Interpolated(dm, vertexLabel, *subdm);CHKERRQ(ierr);
891e6ccafaeSMatthew G Knepley   } else {
892e6ccafaeSMatthew G Knepley     ierr = DMPlexCreateSubmesh_Uninterpolated(dm, vertexLabel, *subdm);CHKERRQ(ierr);
893e6ccafaeSMatthew G Knepley   }
894e6ccafaeSMatthew G Knepley   PetscFunctionReturn(0);
895e6ccafaeSMatthew G Knepley }
896e6ccafaeSMatthew G Knepley 
897e6ccafaeSMatthew G Knepley #undef __FUNCT__
898e6ccafaeSMatthew G Knepley #define __FUNCT__ "DMPlexGetSubpointMap"
899e6ccafaeSMatthew G Knepley PetscErrorCode DMPlexGetSubpointMap(DM dm, DMLabel *subpointMap)
900e6ccafaeSMatthew G Knepley {
901e6ccafaeSMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
902e6ccafaeSMatthew G Knepley 
903e6ccafaeSMatthew G Knepley   PetscFunctionBegin;
904e6ccafaeSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
905e6ccafaeSMatthew G Knepley   PetscValidPointer(subpointMap, 2);
906e6ccafaeSMatthew G Knepley   *subpointMap = mesh->subpointMap;
907e6ccafaeSMatthew G Knepley   PetscFunctionReturn(0);
908e6ccafaeSMatthew G Knepley }
909e6ccafaeSMatthew G Knepley 
910e6ccafaeSMatthew G Knepley #undef __FUNCT__
911e6ccafaeSMatthew G Knepley #define __FUNCT__ "DMPlexSetSubpointMap"
912e6ccafaeSMatthew G Knepley /* Note: Should normally not be called by the user, since it is set in DMPlexCreateSubmesh() */
913e6ccafaeSMatthew G Knepley PetscErrorCode DMPlexSetSubpointMap(DM dm, DMLabel subpointMap)
914e6ccafaeSMatthew G Knepley {
915e6ccafaeSMatthew G Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
916e6ccafaeSMatthew G Knepley   PetscErrorCode ierr;
917e6ccafaeSMatthew G Knepley 
918e6ccafaeSMatthew G Knepley   PetscFunctionBegin;
919e6ccafaeSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
920e6ccafaeSMatthew G Knepley   ierr = DMLabelDestroy(&mesh->subpointMap);CHKERRQ(ierr);
921e6ccafaeSMatthew G Knepley   mesh->subpointMap = subpointMap;
922e6ccafaeSMatthew G Knepley   ++mesh->subpointMap->refct;CHKERRQ(ierr);
923e6ccafaeSMatthew G Knepley   PetscFunctionReturn(0);
924e6ccafaeSMatthew G Knepley }
925e6ccafaeSMatthew G Knepley 
926e6ccafaeSMatthew G Knepley #undef __FUNCT__
927e6ccafaeSMatthew G Knepley #define __FUNCT__ "DMPlexCreateSubpointIS"
928e6ccafaeSMatthew G Knepley /*
929e6ccafaeSMatthew G Knepley   DMPlexCreateSubpointIS - Creates an IS covering the entire subdm chart with the original points as data
930e6ccafaeSMatthew G Knepley 
931e6ccafaeSMatthew G Knepley   Input Parameter:
932e6ccafaeSMatthew G Knepley . dm - The submesh DM
933e6ccafaeSMatthew G Knepley 
934e6ccafaeSMatthew G Knepley   Output Parameter:
935e6ccafaeSMatthew G Knepley . subpointIS - The IS of all the points from the original mesh in this submesh, or PETSC_NULL if this is not a submesh
936e6ccafaeSMatthew G Knepley 
937e6ccafaeSMatthew G Knepley   Note: This is IS is guaranteed to be sorted by the construction of the submesh
938e6ccafaeSMatthew G Knepley */
939e6ccafaeSMatthew G Knepley PetscErrorCode DMPlexCreateSubpointIS(DM dm, IS *subpointIS)
940e6ccafaeSMatthew G Knepley {
941e6ccafaeSMatthew G Knepley   MPI_Comm        comm = ((PetscObject) dm)->comm;
942e6ccafaeSMatthew G Knepley   DMLabel         subpointMap;
943e6ccafaeSMatthew G Knepley   IS              is;
944e6ccafaeSMatthew G Knepley   const PetscInt *opoints;
945e6ccafaeSMatthew G Knepley   PetscInt       *points, *depths;
946e6ccafaeSMatthew G Knepley   PetscInt        depth, depStart, depEnd, d, pStart, pEnd, p, n, off;
947e6ccafaeSMatthew G Knepley   PetscErrorCode  ierr;
948e6ccafaeSMatthew G Knepley 
949e6ccafaeSMatthew G Knepley   PetscFunctionBegin;
950e6ccafaeSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
951e6ccafaeSMatthew G Knepley   PetscValidPointer(subpointIS, 2);
952e6ccafaeSMatthew G Knepley   *subpointIS = PETSC_NULL;
953e6ccafaeSMatthew G Knepley   ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr);
954e6ccafaeSMatthew G Knepley   if (subpointMap) {
955e6ccafaeSMatthew G Knepley     ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
956e6ccafaeSMatthew G Knepley     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
957e6ccafaeSMatthew G Knepley     if (pStart) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Submeshes must start the point numbering at 0, not %d", pStart);
958e6ccafaeSMatthew G Knepley     ierr = DMGetWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr);
959e6ccafaeSMatthew G Knepley     depths[0] = depth;
960e6ccafaeSMatthew G Knepley     depths[1] = 0;
961e6ccafaeSMatthew G Knepley     for(d = 2; d <= depth; ++d) {depths[d] = depth+1 - d;}
962e6ccafaeSMatthew G Knepley     ierr = PetscMalloc(pEnd * sizeof(PetscInt), &points);CHKERRQ(ierr);
963e6ccafaeSMatthew G Knepley     for(d = 0, off = 0; d <= depth; ++d) {
964e6ccafaeSMatthew G Knepley       const PetscInt dep = depths[d];
965e6ccafaeSMatthew G Knepley 
966e6ccafaeSMatthew G Knepley       ierr = DMPlexGetDepthStratum(dm, dep, &depStart, &depEnd);CHKERRQ(ierr);
967e6ccafaeSMatthew G Knepley       ierr = DMLabelGetStratumSize(subpointMap, dep, &n);CHKERRQ(ierr);
968e6ccafaeSMatthew G Knepley       if (d < 2) { /* Only check vertices and cells for now since the map is broken for others */
969e6ccafaeSMatthew G Knepley         if (n != depEnd-depStart) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "The number of mapped submesh points %d at depth %d should be %d", n, dep, depEnd-depStart);
970e6ccafaeSMatthew G Knepley       } else {
971e6ccafaeSMatthew G Knepley         if (!n) for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = PETSC_MAX_INT;
972e6ccafaeSMatthew G Knepley       }
973e6ccafaeSMatthew G Knepley       if (n) {
974e6ccafaeSMatthew G Knepley         ierr = DMLabelGetStratumIS(subpointMap, dep, &is);CHKERRQ(ierr);
975e6ccafaeSMatthew G Knepley         ierr = ISGetIndices(is, &opoints);CHKERRQ(ierr);
976e6ccafaeSMatthew G Knepley         for(p = 0; p < n; ++p, ++off) points[off] = opoints[p];
977e6ccafaeSMatthew G Knepley         ierr = ISRestoreIndices(is, &opoints);CHKERRQ(ierr);
978e6ccafaeSMatthew G Knepley         ierr = ISDestroy(&is);CHKERRQ(ierr);
979e6ccafaeSMatthew G Knepley       }
980e6ccafaeSMatthew G Knepley     }
981e6ccafaeSMatthew G Knepley     ierr = DMRestoreWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr);
982e6ccafaeSMatthew G Knepley     if (off != pEnd) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "The number of mapped submesh points %d should be %d", off, pEnd);
983e6ccafaeSMatthew G Knepley     ierr = ISCreateGeneral(comm, pEnd, points, PETSC_OWN_POINTER, subpointIS);CHKERRQ(ierr);
984e6ccafaeSMatthew G Knepley   }
985e6ccafaeSMatthew G Knepley   PetscFunctionReturn(0);
986e6ccafaeSMatthew G Knepley }
987