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