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