1 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 extern PetscErrorCode DMPlexGetNumFaceVertices_Internal(DM, PetscInt, PetscInt, PetscInt *); 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "DMPlexMarkBoundaryFaces" 7 /*@ 8 DMPlexMarkBoundaryFaces - Mark all faces on the boundary 9 10 Not Collective 11 12 Input Parameter: 13 . dm - The original DM 14 15 Output Parameter: 16 . label - The DMLabel marking boundary faces with value 1 17 18 Level: developer 19 20 .seealso: DMLabelCreate(), DMPlexCreateLabel() 21 @*/ 22 PetscErrorCode DMPlexMarkBoundaryFaces(DM dm, DMLabel label) 23 { 24 PetscInt fStart, fEnd, f; 25 PetscErrorCode ierr; 26 27 PetscFunctionBegin; 28 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 29 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 30 for (f = fStart; f < fEnd; ++f) { 31 PetscInt supportSize; 32 33 ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr); 34 if (supportSize == 1) { 35 ierr = DMLabelSetValue(label, f, 1);CHKERRQ(ierr); 36 } 37 } 38 PetscFunctionReturn(0); 39 } 40 41 #undef __FUNCT__ 42 #define __FUNCT__ "DMPlexShiftPoint_Internal" 43 PETSC_STATIC_INLINE PetscInt DMPlexShiftPoint_Internal(PetscInt p, PetscInt depth, PetscInt depthEnd[], PetscInt depthShift[]) 44 { 45 if (depth < 0) return p; 46 /* Cells */ if (p < depthEnd[depth]) return p; 47 /* Vertices */ if (p < depthEnd[0]) return p + depthShift[depth]; 48 /* Faces */ if (p < depthEnd[depth-1]) return p + depthShift[depth] + depthShift[0]; 49 /* Edges */ return p + depthShift[depth] + depthShift[0] + depthShift[depth-1]; 50 } 51 52 #undef __FUNCT__ 53 #define __FUNCT__ "DMPlexShiftSizes_Internal" 54 static PetscErrorCode DMPlexShiftSizes_Internal(DM dm, PetscInt depthShift[], DM dmNew) 55 { 56 PetscInt *depthEnd; 57 PetscInt depth = 0, d, pStart, pEnd, p; 58 PetscErrorCode ierr; 59 60 PetscFunctionBegin; 61 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 62 if (depth < 0) PetscFunctionReturn(0); 63 ierr = PetscMalloc((depth+1) * sizeof(PetscInt), &depthEnd);CHKERRQ(ierr); 64 /* Step 1: Expand chart */ 65 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 66 for (d = 0; d <= depth; ++d) { 67 pEnd += depthShift[d]; 68 ierr = DMPlexGetDepthStratum(dm, d, NULL, &depthEnd[d]);CHKERRQ(ierr); 69 } 70 ierr = DMPlexSetChart(dmNew, pStart, pEnd);CHKERRQ(ierr); 71 /* Step 2: Set cone and support sizes */ 72 for (d = 0; d <= depth; ++d) { 73 ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 74 for (p = pStart; p < pEnd; ++p) { 75 PetscInt newp = DMPlexShiftPoint_Internal(p, depth, depthEnd, depthShift); 76 PetscInt size; 77 78 ierr = DMPlexGetConeSize(dm, p, &size);CHKERRQ(ierr); 79 ierr = DMPlexSetConeSize(dmNew, newp, size);CHKERRQ(ierr); 80 ierr = DMPlexGetSupportSize(dm, p, &size);CHKERRQ(ierr); 81 ierr = DMPlexSetSupportSize(dmNew, newp, size);CHKERRQ(ierr); 82 } 83 } 84 ierr = PetscFree(depthEnd);CHKERRQ(ierr); 85 PetscFunctionReturn(0); 86 } 87 88 #undef __FUNCT__ 89 #define __FUNCT__ "DMPlexShiftPoints_Internal" 90 static PetscErrorCode DMPlexShiftPoints_Internal(DM dm, PetscInt depthShift[], DM dmNew) 91 { 92 PetscInt *depthEnd, *newpoints; 93 PetscInt depth = 0, d, maxConeSize, maxSupportSize, pStart, pEnd, p; 94 PetscErrorCode ierr; 95 96 PetscFunctionBegin; 97 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 98 if (depth < 0) PetscFunctionReturn(0); 99 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 100 ierr = PetscMalloc2(depth+1,PetscInt,&depthEnd,PetscMax(maxConeSize, maxSupportSize),PetscInt,&newpoints);CHKERRQ(ierr); 101 for (d = 0; d <= depth; ++d) { 102 ierr = DMPlexGetDepthStratum(dm, d, NULL, &depthEnd[d]);CHKERRQ(ierr); 103 } 104 /* Step 5: Set cones and supports */ 105 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 106 for (p = pStart; p < pEnd; ++p) { 107 const PetscInt *points = NULL, *orientations = NULL; 108 PetscInt size, i, newp = DMPlexShiftPoint_Internal(p, depth, depthEnd, depthShift); 109 110 ierr = DMPlexGetConeSize(dm, p, &size);CHKERRQ(ierr); 111 ierr = DMPlexGetCone(dm, p, &points);CHKERRQ(ierr); 112 ierr = DMPlexGetConeOrientation(dm, p, &orientations);CHKERRQ(ierr); 113 for (i = 0; i < size; ++i) { 114 newpoints[i] = DMPlexShiftPoint_Internal(points[i], depth, depthEnd, depthShift); 115 } 116 ierr = DMPlexSetCone(dmNew, newp, newpoints);CHKERRQ(ierr); 117 ierr = DMPlexSetConeOrientation(dmNew, newp, orientations);CHKERRQ(ierr); 118 ierr = DMPlexGetSupportSize(dm, p, &size);CHKERRQ(ierr); 119 ierr = DMPlexGetSupport(dm, p, &points);CHKERRQ(ierr); 120 for (i = 0; i < size; ++i) { 121 newpoints[i] = DMPlexShiftPoint_Internal(points[i], depth, depthEnd, depthShift); 122 } 123 ierr = DMPlexSetSupport(dmNew, newp, newpoints);CHKERRQ(ierr); 124 } 125 ierr = PetscFree2(depthEnd,newpoints);CHKERRQ(ierr); 126 PetscFunctionReturn(0); 127 } 128 129 #undef __FUNCT__ 130 #define __FUNCT__ "DMPlexShiftCoordinates_Internal" 131 static PetscErrorCode DMPlexShiftCoordinates_Internal(DM dm, PetscInt depthShift[], DM dmNew) 132 { 133 PetscSection coordSection, newCoordSection; 134 Vec coordinates, newCoordinates; 135 PetscScalar *coords, *newCoords; 136 PetscInt *depthEnd, coordSize; 137 PetscInt dim, depth = 0, d, vStart, vEnd, vStartNew, vEndNew, v; 138 PetscErrorCode ierr; 139 140 PetscFunctionBegin; 141 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 142 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 143 ierr = PetscMalloc((depth+1) * sizeof(PetscInt), &depthEnd);CHKERRQ(ierr); 144 for (d = 0; d <= depth; ++d) { 145 ierr = DMPlexGetDepthStratum(dm, d, NULL, &depthEnd[d]);CHKERRQ(ierr); 146 } 147 /* Step 8: Convert coordinates */ 148 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 149 ierr = DMPlexGetDepthStratum(dmNew, 0, &vStartNew, &vEndNew);CHKERRQ(ierr); 150 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 151 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &newCoordSection);CHKERRQ(ierr); 152 ierr = PetscSectionSetNumFields(newCoordSection, 1);CHKERRQ(ierr); 153 ierr = PetscSectionSetFieldComponents(newCoordSection, 0, dim);CHKERRQ(ierr); 154 ierr = PetscSectionSetChart(newCoordSection, vStartNew, vEndNew);CHKERRQ(ierr); 155 for (v = vStartNew; v < vEndNew; ++v) { 156 ierr = PetscSectionSetDof(newCoordSection, v, dim);CHKERRQ(ierr); 157 ierr = PetscSectionSetFieldDof(newCoordSection, v, 0, dim);CHKERRQ(ierr); 158 } 159 ierr = PetscSectionSetUp(newCoordSection);CHKERRQ(ierr); 160 ierr = DMPlexSetCoordinateSection(dmNew, newCoordSection);CHKERRQ(ierr); 161 ierr = PetscSectionGetStorageSize(newCoordSection, &coordSize);CHKERRQ(ierr); 162 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &newCoordinates);CHKERRQ(ierr); 163 ierr = PetscObjectSetName((PetscObject) newCoordinates, "coordinates");CHKERRQ(ierr); 164 ierr = VecSetSizes(newCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 165 ierr = VecSetFromOptions(newCoordinates);CHKERRQ(ierr); 166 ierr = DMSetCoordinatesLocal(dmNew, newCoordinates);CHKERRQ(ierr); 167 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 168 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 169 ierr = VecGetArray(newCoordinates, &newCoords);CHKERRQ(ierr); 170 for (v = vStart; v < vEnd; ++v) { 171 PetscInt dof, off, noff, d; 172 173 ierr = PetscSectionGetDof(coordSection, v, &dof);CHKERRQ(ierr); 174 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 175 ierr = PetscSectionGetOffset(newCoordSection, DMPlexShiftPoint_Internal(v, depth, depthEnd, depthShift), &noff);CHKERRQ(ierr); 176 for (d = 0; d < dof; ++d) { 177 newCoords[noff+d] = coords[off+d]; 178 } 179 } 180 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 181 ierr = VecRestoreArray(newCoordinates, &newCoords);CHKERRQ(ierr); 182 ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr); 183 ierr = PetscSectionDestroy(&newCoordSection);CHKERRQ(ierr); 184 ierr = PetscFree(depthEnd);CHKERRQ(ierr); 185 PetscFunctionReturn(0); 186 } 187 188 #undef __FUNCT__ 189 #define __FUNCT__ "DMPlexShiftSF_Internal" 190 static PetscErrorCode DMPlexShiftSF_Internal(DM dm, PetscInt depthShift[], DM dmNew) 191 { 192 PetscInt *depthEnd; 193 PetscInt depth = 0, d; 194 PetscSF sfPoint, sfPointNew; 195 const PetscSFNode *remotePoints; 196 PetscSFNode *gremotePoints; 197 const PetscInt *localPoints; 198 PetscInt *glocalPoints, *newLocation, *newRemoteLocation; 199 PetscInt numRoots, numLeaves, l, pStart, pEnd, totShift = 0; 200 PetscMPIInt numProcs; 201 PetscErrorCode ierr; 202 203 PetscFunctionBegin; 204 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 205 ierr = PetscMalloc((depth+1) * sizeof(PetscInt), &depthEnd);CHKERRQ(ierr); 206 for (d = 0; d <= depth; ++d) { 207 totShift += depthShift[d]; 208 ierr = DMPlexGetDepthStratum(dm, d, NULL, &depthEnd[d]);CHKERRQ(ierr); 209 } 210 /* Step 9: Convert pointSF */ 211 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr); 212 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 213 ierr = DMGetPointSF(dmNew, &sfPointNew);CHKERRQ(ierr); 214 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 215 ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 216 if (numRoots >= 0) { 217 ierr = PetscMalloc2(numRoots,PetscInt,&newLocation,pEnd-pStart,PetscInt,&newRemoteLocation);CHKERRQ(ierr); 218 for (l=0; l<numRoots; l++) newLocation[l] = DMPlexShiftPoint_Internal(l, depth, depthEnd, depthShift); 219 ierr = PetscSFBcastBegin(sfPoint, MPIU_INT, newLocation, newRemoteLocation);CHKERRQ(ierr); 220 ierr = PetscSFBcastEnd(sfPoint, MPIU_INT, newLocation, newRemoteLocation);CHKERRQ(ierr); 221 ierr = PetscMalloc(numLeaves * sizeof(PetscInt), &glocalPoints);CHKERRQ(ierr); 222 ierr = PetscMalloc(numLeaves * sizeof(PetscSFNode), &gremotePoints);CHKERRQ(ierr); 223 for (l = 0; l < numLeaves; ++l) { 224 glocalPoints[l] = DMPlexShiftPoint_Internal(localPoints[l], depth, depthEnd, depthShift); 225 gremotePoints[l].rank = remotePoints[l].rank; 226 gremotePoints[l].index = newRemoteLocation[localPoints[l]]; 227 } 228 ierr = PetscFree2(newLocation,newRemoteLocation);CHKERRQ(ierr); 229 ierr = PetscSFSetGraph(sfPointNew, numRoots + totShift, numLeaves, glocalPoints, PETSC_OWN_POINTER, gremotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 230 } 231 ierr = PetscFree(depthEnd);CHKERRQ(ierr); 232 PetscFunctionReturn(0); 233 } 234 235 #undef __FUNCT__ 236 #define __FUNCT__ "DMPlexShiftLabels_Internal" 237 static PetscErrorCode DMPlexShiftLabels_Internal(DM dm, PetscInt depthShift[], DM dmNew) 238 { 239 PetscSF sfPoint; 240 DMLabel vtkLabel, ghostLabel; 241 PetscInt *depthEnd; 242 const PetscSFNode *leafRemote; 243 const PetscInt *leafLocal; 244 PetscInt depth = 0, d, numLeaves, numLabels, l, cStart, cEnd, c, fStart, fEnd, f; 245 PetscMPIInt rank; 246 PetscErrorCode ierr; 247 248 PetscFunctionBegin; 249 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 250 ierr = PetscMalloc((depth+1) * sizeof(PetscInt), &depthEnd);CHKERRQ(ierr); 251 for (d = 0; d <= depth; ++d) { 252 ierr = DMPlexGetDepthStratum(dm, d, NULL, &depthEnd[d]);CHKERRQ(ierr); 253 } 254 /* Step 10: Convert labels */ 255 ierr = DMPlexGetNumLabels(dm, &numLabels);CHKERRQ(ierr); 256 for (l = 0; l < numLabels; ++l) { 257 DMLabel label, newlabel; 258 const char *lname; 259 PetscBool isDepth; 260 IS valueIS; 261 const PetscInt *values; 262 PetscInt numValues, val; 263 264 ierr = DMPlexGetLabelName(dm, l, &lname);CHKERRQ(ierr); 265 ierr = PetscStrcmp(lname, "depth", &isDepth);CHKERRQ(ierr); 266 if (isDepth) continue; 267 ierr = DMPlexCreateLabel(dmNew, lname);CHKERRQ(ierr); 268 ierr = DMPlexGetLabel(dm, lname, &label);CHKERRQ(ierr); 269 ierr = DMPlexGetLabel(dmNew, lname, &newlabel);CHKERRQ(ierr); 270 ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr); 271 ierr = ISGetLocalSize(valueIS, &numValues);CHKERRQ(ierr); 272 ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr); 273 for (val = 0; val < numValues; ++val) { 274 IS pointIS; 275 const PetscInt *points; 276 PetscInt numPoints, p; 277 278 ierr = DMLabelGetStratumIS(label, values[val], &pointIS);CHKERRQ(ierr); 279 ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr); 280 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 281 for (p = 0; p < numPoints; ++p) { 282 const PetscInt newpoint = DMPlexShiftPoint_Internal(points[p], depth, depthEnd, depthShift); 283 284 ierr = DMLabelSetValue(newlabel, newpoint, values[val]);CHKERRQ(ierr); 285 } 286 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 287 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 288 } 289 ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr); 290 ierr = ISDestroy(&valueIS);CHKERRQ(ierr); 291 } 292 ierr = PetscFree(depthEnd);CHKERRQ(ierr); 293 /* Step 11: Make label for output (vtk) and to mark ghost points (ghost) */ 294 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 295 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 296 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 297 ierr = PetscSFGetGraph(sfPoint, NULL, &numLeaves, &leafLocal, &leafRemote);CHKERRQ(ierr); 298 ierr = DMPlexCreateLabel(dmNew, "vtk");CHKERRQ(ierr); 299 ierr = DMPlexCreateLabel(dmNew, "ghost");CHKERRQ(ierr); 300 ierr = DMPlexGetLabel(dmNew, "vtk", &vtkLabel);CHKERRQ(ierr); 301 ierr = DMPlexGetLabel(dmNew, "ghost", &ghostLabel);CHKERRQ(ierr); 302 for (l = 0, c = cStart; l < numLeaves && c < cEnd; ++l, ++c) { 303 for (; c < leafLocal[l] && c < cEnd; ++c) { 304 ierr = DMLabelSetValue(vtkLabel, c, 1);CHKERRQ(ierr); 305 } 306 if (leafLocal[l] >= cEnd) break; 307 if (leafRemote[l].rank == rank) { 308 ierr = DMLabelSetValue(vtkLabel, c, 1);CHKERRQ(ierr); 309 } else { 310 ierr = DMLabelSetValue(ghostLabel, c, 2);CHKERRQ(ierr); 311 } 312 } 313 for (; c < cEnd; ++c) { 314 ierr = DMLabelSetValue(vtkLabel, c, 1);CHKERRQ(ierr); 315 } 316 if (0) { 317 ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr); 318 ierr = DMLabelView(vtkLabel, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 319 ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 320 } 321 ierr = DMPlexGetHeightStratum(dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr); 322 for (f = fStart; f < fEnd; ++f) { 323 PetscInt numCells; 324 325 ierr = DMPlexGetSupportSize(dmNew, f, &numCells);CHKERRQ(ierr); 326 if (numCells < 2) { 327 ierr = DMLabelSetValue(ghostLabel, f, 1);CHKERRQ(ierr); 328 } else { 329 const PetscInt *cells = NULL; 330 PetscInt vA, vB; 331 332 ierr = DMPlexGetSupport(dmNew, f, &cells);CHKERRQ(ierr); 333 ierr = DMLabelGetValue(vtkLabel, cells[0], &vA);CHKERRQ(ierr); 334 ierr = DMLabelGetValue(vtkLabel, cells[1], &vB);CHKERRQ(ierr); 335 if (!vA && !vB) {ierr = DMLabelSetValue(ghostLabel, f, 1);CHKERRQ(ierr);} 336 } 337 } 338 if (0) { 339 ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr); 340 ierr = DMLabelView(ghostLabel, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 341 ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 342 } 343 PetscFunctionReturn(0); 344 } 345 346 #undef __FUNCT__ 347 #define __FUNCT__ "DMPlexConstructGhostCells_Internal" 348 static PetscErrorCode DMPlexConstructGhostCells_Internal(DM dm, DMLabel label, PetscInt *numGhostCells, DM gdm) 349 { 350 IS valueIS; 351 const PetscInt *values; 352 PetscInt *depthShift; 353 PetscInt depth = 0, numFS, fs, ghostCell, cEnd, c; 354 PetscErrorCode ierr; 355 356 PetscFunctionBegin; 357 /* Count ghost cells */ 358 ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr); 359 ierr = ISGetLocalSize(valueIS, &numFS);CHKERRQ(ierr); 360 ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr); 361 362 *numGhostCells = 0; 363 for (fs = 0; fs < numFS; ++fs) { 364 PetscInt numBdFaces; 365 366 ierr = DMLabelGetStratumSize(label, values[fs], &numBdFaces);CHKERRQ(ierr); 367 368 *numGhostCells += numBdFaces; 369 } 370 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 371 ierr = PetscMalloc((depth+1) * sizeof(PetscInt), &depthShift);CHKERRQ(ierr); 372 ierr = PetscMemzero(depthShift, (depth+1) * sizeof(PetscInt));CHKERRQ(ierr); 373 if (depth >= 0) depthShift[depth] = *numGhostCells; 374 ierr = DMPlexShiftSizes_Internal(dm, depthShift, gdm);CHKERRQ(ierr); 375 /* Step 3: Set cone/support sizes for new points */ 376 ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr); 377 for (c = cEnd; c < cEnd + *numGhostCells; ++c) { 378 ierr = DMPlexSetConeSize(gdm, c, 1);CHKERRQ(ierr); 379 } 380 for (fs = 0; fs < numFS; ++fs) { 381 IS faceIS; 382 const PetscInt *faces; 383 PetscInt numFaces, f; 384 385 ierr = DMLabelGetStratumIS(label, values[fs], &faceIS);CHKERRQ(ierr); 386 ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr); 387 ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr); 388 for (f = 0; f < numFaces; ++f) { 389 PetscInt size; 390 391 ierr = DMPlexGetSupportSize(dm, faces[f], &size);CHKERRQ(ierr); 392 if (size != 1) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM has boundary face %d with %d support cells", faces[f], size); 393 ierr = DMPlexSetSupportSize(gdm, faces[f] + *numGhostCells, 2);CHKERRQ(ierr); 394 } 395 ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 396 ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 397 } 398 /* Step 4: Setup ghosted DM */ 399 ierr = DMSetUp(gdm);CHKERRQ(ierr); 400 ierr = DMPlexShiftPoints_Internal(dm, depthShift, gdm);CHKERRQ(ierr); 401 /* Step 6: Set cones and supports for new points */ 402 ghostCell = cEnd; 403 for (fs = 0; fs < numFS; ++fs) { 404 IS faceIS; 405 const PetscInt *faces; 406 PetscInt numFaces, f; 407 408 ierr = DMLabelGetStratumIS(label, values[fs], &faceIS);CHKERRQ(ierr); 409 ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr); 410 ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr); 411 for (f = 0; f < numFaces; ++f, ++ghostCell) { 412 PetscInt newFace = faces[f] + *numGhostCells; 413 414 ierr = DMPlexSetCone(gdm, ghostCell, &newFace);CHKERRQ(ierr); 415 ierr = DMPlexInsertSupport(gdm, newFace, 1, ghostCell);CHKERRQ(ierr); 416 } 417 ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 418 ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 419 } 420 ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr); 421 ierr = ISDestroy(&valueIS);CHKERRQ(ierr); 422 /* Step 7: Stratify */ 423 ierr = DMPlexStratify(gdm);CHKERRQ(ierr); 424 ierr = DMPlexShiftCoordinates_Internal(dm, depthShift, gdm);CHKERRQ(ierr); 425 ierr = DMPlexShiftSF_Internal(dm, depthShift, gdm);CHKERRQ(ierr); 426 ierr = DMPlexShiftLabels_Internal(dm, depthShift, gdm);CHKERRQ(ierr); 427 ierr = PetscFree(depthShift);CHKERRQ(ierr); 428 PetscFunctionReturn(0); 429 } 430 431 #undef __FUNCT__ 432 #define __FUNCT__ "DMPlexConstructGhostCells" 433 /*@C 434 DMPlexConstructGhostCells - Construct ghost cells which connect to every boundary face 435 436 Collective on dm 437 438 Input Parameters: 439 + dm - The original DM 440 - labelName - The label specifying the boundary faces, or "Face Sets" if this is NULL 441 442 Output Parameters: 443 + numGhostCells - The number of ghost cells added to the DM 444 - dmGhosted - The new DM 445 446 Note: If no label exists of that name, one will be created marking all boundary faces 447 448 Level: developer 449 450 .seealso: DMCreate() 451 */ 452 PetscErrorCode DMPlexConstructGhostCells(DM dm, const char labelName[], PetscInt *numGhostCells, DM *dmGhosted) 453 { 454 DM gdm; 455 DMLabel label; 456 const char *name = labelName ? labelName : "Face Sets"; 457 PetscInt dim; 458 PetscErrorCode ierr; 459 460 PetscFunctionBegin; 461 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 462 PetscValidPointer(numGhostCells, 3); 463 PetscValidPointer(dmGhosted, 4); 464 ierr = DMCreate(PetscObjectComm((PetscObject)dm), &gdm);CHKERRQ(ierr); 465 ierr = DMSetType(gdm, DMPLEX);CHKERRQ(ierr); 466 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 467 ierr = DMPlexSetDimension(gdm, dim);CHKERRQ(ierr); 468 ierr = DMPlexGetLabel(dm, name, &label);CHKERRQ(ierr); 469 if (!label) { 470 /* Get label for boundary faces */ 471 ierr = DMPlexCreateLabel(dm, name);CHKERRQ(ierr); 472 ierr = DMPlexGetLabel(dm, name, &label);CHKERRQ(ierr); 473 ierr = DMPlexMarkBoundaryFaces(dm, label);CHKERRQ(ierr); 474 } 475 ierr = DMPlexConstructGhostCells_Internal(dm, label, numGhostCells, gdm);CHKERRQ(ierr); 476 ierr = DMSetFromOptions(gdm);CHKERRQ(ierr); 477 *dmGhosted = gdm; 478 PetscFunctionReturn(0); 479 } 480 481 #undef __FUNCT__ 482 #define __FUNCT__ "DMPlexConstructCohesiveCells_Internal" 483 static PetscErrorCode DMPlexConstructCohesiveCells_Internal(DM dm, DMLabel label, DM sdm) 484 { 485 MPI_Comm comm; 486 IS valueIS, *pointIS; 487 const PetscInt *values, **splitPoints; 488 PetscSection coordSection; 489 Vec coordinates; 490 PetscScalar *coords; 491 PetscInt *depthShift, *depthOffset, *pMaxNew, *numSplitPoints, *coneNew, *supportNew; 492 PetscInt shift = 100, depth = 0, dep, dim, d, numSP = 0, sp, maxConeSize, maxSupportSize, numLabels, p, v; 493 PetscErrorCode ierr; 494 495 PetscFunctionBegin; 496 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 497 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 498 /* Count split points and add cohesive cells */ 499 if (label) { 500 ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr); 501 ierr = ISGetLocalSize(valueIS, &numSP);CHKERRQ(ierr); 502 ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr); 503 } 504 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 505 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 506 ierr = PetscMalloc5(depth+1,PetscInt,&depthShift,depth+1,PetscInt,&depthOffset,depth+1,PetscInt,&pMaxNew,maxConeSize*3,PetscInt,&coneNew,maxSupportSize,PetscInt,&supportNew);CHKERRQ(ierr); 507 ierr = PetscMalloc3(depth+1,IS,&pointIS,depth+1,PetscInt,&numSplitPoints,depth+1,const PetscInt*,&splitPoints);CHKERRQ(ierr); 508 ierr = PetscMemzero(depthShift, (depth+1) * sizeof(PetscInt));CHKERRQ(ierr); 509 for (d = 0; d <= depth; ++d) { 510 ierr = DMPlexGetDepthStratum(dm, d, NULL, &pMaxNew[d]);CHKERRQ(ierr); 511 numSplitPoints[d] = 0; 512 splitPoints[d] = NULL; 513 pointIS[d] = NULL; 514 } 515 for (sp = 0; sp < numSP; ++sp) { 516 const PetscInt dep = values[sp]; 517 518 if ((dep < 0) || (dep > depth)) continue; 519 ierr = DMLabelGetStratumSize(label, dep, &depthShift[dep]);CHKERRQ(ierr); 520 ierr = DMLabelGetStratumIS(label, dep, &pointIS[dep]);CHKERRQ(ierr); 521 if (pointIS[dep]) { 522 ierr = ISGetLocalSize(pointIS[dep], &numSplitPoints[dep]);CHKERRQ(ierr); 523 ierr = ISGetIndices(pointIS[dep], &splitPoints[dep]);CHKERRQ(ierr); 524 } 525 } 526 if (depth >= 0) { 527 /* Calculate number of additional points */ 528 depthShift[depth] = depthShift[depth-1]; /* There is a cohesive cell for every split face */ 529 depthShift[1] += depthShift[0]; /* There is a cohesive edge for every split vertex */ 530 /* Calculate hybrid bound for each dimension */ 531 pMaxNew[0] += depthShift[depth]; 532 if (depth > 1) pMaxNew[dim-1] += depthShift[depth] + depthShift[0]; 533 if (depth > 2) pMaxNew[1] += depthShift[depth] + depthShift[0] + depthShift[dim-1]; 534 535 /* Calculate point offset for each dimension */ 536 depthOffset[depth] = 0; 537 depthOffset[0] = depthOffset[depth] + depthShift[depth]; 538 if (depth > 1) depthOffset[dim-1] = depthOffset[0] + depthShift[0]; 539 if (depth > 2) depthOffset[1] = depthOffset[dim-1] + depthShift[dim-1]; 540 } 541 ierr = DMPlexShiftSizes_Internal(dm, depthShift, sdm);CHKERRQ(ierr); 542 /* Step 3: Set cone/support sizes for new points */ 543 for (dep = 0; dep <= depth; ++dep) { 544 for (p = 0; p < numSplitPoints[dep]; ++p) { 545 const PetscInt oldp = splitPoints[dep][p]; 546 const PetscInt newp = depthOffset[dep] + oldp; 547 const PetscInt splitp = pMaxNew[dep] + p; 548 const PetscInt *support; 549 PetscInt coneSize, supportSize, q, e; 550 551 ierr = DMPlexGetConeSize(dm, oldp, &coneSize);CHKERRQ(ierr); 552 ierr = DMPlexSetConeSize(sdm, splitp, coneSize);CHKERRQ(ierr); 553 ierr = DMPlexGetSupportSize(dm, oldp, &supportSize);CHKERRQ(ierr); 554 ierr = DMPlexSetSupportSize(sdm, splitp, supportSize);CHKERRQ(ierr); 555 if (dep == depth-1) { 556 const PetscInt ccell = pMaxNew[depth] + p; 557 /* Add cohesive cells, they are prisms */ 558 ierr = DMPlexSetConeSize(sdm, ccell, 2 + coneSize);CHKERRQ(ierr); 559 } else if (dep == 0) { 560 const PetscInt cedge = pMaxNew[1] + (depthShift[1] - depthShift[0]) + p; 561 562 ierr = DMPlexGetSupport(dm, oldp, &support);CHKERRQ(ierr); 563 /* Split old vertex: Edges in old split faces and new cohesive edge */ 564 for (e = 0, q = 0; e < supportSize; ++e) { 565 PetscInt val; 566 567 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 568 if ((val == 1) || (val == (shift + 1))) ++q; 569 } 570 ierr = DMPlexSetSupportSize(sdm, newp, q+1);CHKERRQ(ierr); 571 /* Split new vertex: Edges in new split faces and new cohesive edge */ 572 for (e = 0, q = 0; e < supportSize; ++e) { 573 PetscInt val; 574 575 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 576 if ((val == 1) || (val == -(shift + 1))) ++q; 577 } 578 ierr = DMPlexSetSupportSize(sdm, splitp, q+1);CHKERRQ(ierr); 579 /* Add cohesive edges */ 580 ierr = DMPlexSetConeSize(sdm, cedge, 2);CHKERRQ(ierr); 581 /* Punt for now on support, you loop over closure, extract faces, check which ones are in the label */ 582 } else if (dep == dim-2) { 583 ierr = DMPlexGetSupport(dm, oldp, &support);CHKERRQ(ierr); 584 /* Split old edge: Faces in positive side cells and old split faces */ 585 for (e = 0, q = 0; e < supportSize; ++e) { 586 PetscInt val; 587 588 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 589 if ((val == dim-1) || (val == (shift + dim-1))) ++q; 590 } 591 ierr = DMPlexSetSupportSize(sdm, newp, q);CHKERRQ(ierr); 592 /* Split new edge: Faces in negative side cells and new split faces */ 593 for (e = 0, q = 0; e < supportSize; ++e) { 594 PetscInt val; 595 596 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 597 if ((val == dim-1) || (val == -(shift + dim-1))) ++q; 598 } 599 ierr = DMPlexSetSupportSize(sdm, splitp, q);CHKERRQ(ierr); 600 } 601 } 602 } 603 /* Step 4: Setup split DM */ 604 ierr = DMSetUp(sdm);CHKERRQ(ierr); 605 ierr = DMPlexShiftPoints_Internal(dm, depthShift, sdm);CHKERRQ(ierr); 606 /* Step 6: Set cones and supports for new points */ 607 for (dep = 0; dep <= depth; ++dep) { 608 for (p = 0; p < numSplitPoints[dep]; ++p) { 609 const PetscInt oldp = splitPoints[dep][p]; 610 const PetscInt newp = depthOffset[dep] + oldp; 611 const PetscInt splitp = pMaxNew[dep] + p; 612 const PetscInt *cone, *support, *ornt; 613 PetscInt coneSize, supportSize, q, v, e, s; 614 615 ierr = DMPlexGetConeSize(dm, oldp, &coneSize);CHKERRQ(ierr); 616 ierr = DMPlexGetCone(dm, oldp, &cone);CHKERRQ(ierr); 617 ierr = DMPlexGetConeOrientation(dm, oldp, &ornt);CHKERRQ(ierr); 618 ierr = DMPlexGetSupportSize(dm, oldp, &supportSize);CHKERRQ(ierr); 619 ierr = DMPlexGetSupport(dm, oldp, &support);CHKERRQ(ierr); 620 if (dep == depth-1) { 621 const PetscInt ccell = pMaxNew[depth] + p; 622 const PetscInt *supportF; 623 624 /* Split face: copy in old face to new face to start */ 625 ierr = DMPlexGetSupport(sdm, newp, &supportF);CHKERRQ(ierr); 626 ierr = DMPlexSetSupport(sdm, splitp, supportF);CHKERRQ(ierr); 627 /* Split old face: old vertices/edges in cone so no change */ 628 /* Split new face: new vertices/edges in cone */ 629 for (q = 0; q < coneSize; ++q) { 630 ierr = PetscFindInt(cone[q], numSplitPoints[dim-2], splitPoints[dim-2], &v);CHKERRQ(ierr); 631 632 coneNew[2+q] = pMaxNew[dim-2] + v; 633 } 634 ierr = DMPlexSetCone(sdm, splitp, &coneNew[2]);CHKERRQ(ierr); 635 ierr = DMPlexSetConeOrientation(sdm, splitp, ornt);CHKERRQ(ierr); 636 /* Cohesive cell: Old and new split face, then new cohesive edges */ 637 coneNew[0] = newp; 638 coneNew[1] = splitp; 639 for (q = 0; q < coneSize; ++q) { 640 coneNew[2+q] = (pMaxNew[1] - pMaxNew[dim-2]) + (depthShift[1] - depthShift[0]) + coneNew[2+q]; 641 } 642 ierr = DMPlexSetCone(sdm, ccell, coneNew);CHKERRQ(ierr); 643 644 645 for (s = 0; s < supportSize; ++s) { 646 PetscInt val; 647 648 ierr = DMLabelGetValue(label, support[s], &val);CHKERRQ(ierr); 649 if (val < 0) { 650 /* Split old face: Replace negative side cell with cohesive cell */ 651 ierr = DMPlexInsertSupport(sdm, newp, s, ccell);CHKERRQ(ierr); 652 } else { 653 /* Split new face: Replace positive side cell with cohesive cell */ 654 ierr = DMPlexInsertSupport(sdm, splitp, s, ccell);CHKERRQ(ierr); 655 } 656 } 657 } else if (dep == 0) { 658 const PetscInt cedge = pMaxNew[1] + (depthShift[1] - depthShift[0]) + p; 659 660 /* Split old vertex: Edges in old split faces and new cohesive edge */ 661 for (e = 0, q = 0; e < supportSize; ++e) { 662 PetscInt val; 663 664 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 665 if ((val == 1) || (val == (shift + 1))) { 666 supportNew[q++] = depthOffset[1] + support[e]; 667 } 668 } 669 supportNew[q] = cedge; 670 671 ierr = DMPlexSetSupport(sdm, newp, supportNew);CHKERRQ(ierr); 672 /* Split new vertex: Edges in new split faces and new cohesive edge */ 673 for (e = 0, q = 0; e < supportSize; ++e) { 674 PetscInt val, edge; 675 676 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 677 if (val == 1) { 678 ierr = PetscFindInt(support[e], numSplitPoints[1], splitPoints[1], &edge);CHKERRQ(ierr); 679 if (edge < 0) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Edge %d is not a split edge", support[e]); 680 supportNew[q++] = pMaxNew[1] + edge; 681 } else if (val == -(shift + 1)) { 682 supportNew[q++] = depthOffset[1] + support[e]; 683 } 684 } 685 supportNew[q] = cedge; 686 ierr = DMPlexSetSupport(sdm, splitp, supportNew);CHKERRQ(ierr); 687 /* Cohesive edge: Old and new split vertex, punting on support */ 688 coneNew[0] = newp; 689 coneNew[1] = splitp; 690 ierr = DMPlexSetCone(sdm, cedge, coneNew);CHKERRQ(ierr); 691 } else if (dep == dim-2) { 692 /* Split old edge: old vertices in cone so no change */ 693 /* Split new edge: new vertices in cone */ 694 for (q = 0; q < coneSize; ++q) { 695 ierr = PetscFindInt(cone[q], numSplitPoints[dim-3], splitPoints[dim-3], &v);CHKERRQ(ierr); 696 697 coneNew[q] = pMaxNew[dim-3] + v; 698 } 699 ierr = DMPlexSetCone(sdm, splitp, coneNew);CHKERRQ(ierr); 700 /* Split old edge: Faces in positive side cells and old split faces */ 701 for (e = 0, q = 0; e < supportSize; ++e) { 702 PetscInt val; 703 704 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 705 if ((val == dim-1) || (val == (shift + dim-1))) { 706 supportNew[q++] = depthOffset[dim-1] + support[e]; 707 } 708 } 709 ierr = DMPlexSetSupport(sdm, newp, supportNew);CHKERRQ(ierr); 710 /* Split new edge: Faces in negative side cells and new split faces */ 711 for (e = 0, q = 0; e < supportSize; ++e) { 712 PetscInt val, face; 713 714 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 715 if (val == dim-1) { 716 ierr = PetscFindInt(support[e], numSplitPoints[dim-1], splitPoints[dim-1], &face);CHKERRQ(ierr); 717 if (face < 0) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Face %d is not a split face", support[e]); 718 supportNew[q++] = pMaxNew[dim-1] + face; 719 } else if (val == -(shift + dim-1)) { 720 supportNew[q++] = depthOffset[dim-1] + support[e]; 721 } 722 } 723 ierr = DMPlexSetSupport(sdm, splitp, supportNew);CHKERRQ(ierr); 724 } 725 } 726 } 727 /* Step 6b: Replace split points in negative side cones */ 728 for (sp = 0; sp < numSP; ++sp) { 729 PetscInt dep = values[sp]; 730 IS pIS; 731 PetscInt numPoints; 732 const PetscInt *points; 733 734 if (dep >= 0) continue; 735 ierr = DMLabelGetStratumIS(label, dep, &pIS);CHKERRQ(ierr); 736 if (!pIS) continue; 737 dep = -dep - shift; 738 ierr = ISGetLocalSize(pIS, &numPoints);CHKERRQ(ierr); 739 ierr = ISGetIndices(pIS, &points);CHKERRQ(ierr); 740 for (p = 0; p < numPoints; ++p) { 741 const PetscInt oldp = points[p]; 742 const PetscInt newp = depthOffset[dep] + oldp; 743 const PetscInt *cone; 744 PetscInt coneSize, c; 745 PetscBool replaced = PETSC_FALSE; 746 747 /* Negative edge: replace split vertex */ 748 /* Negative cell: replace split face */ 749 ierr = DMPlexGetConeSize(sdm, newp, &coneSize);CHKERRQ(ierr); 750 ierr = DMPlexGetCone(sdm, newp, &cone);CHKERRQ(ierr); 751 for (c = 0; c < coneSize; ++c) { 752 const PetscInt coldp = cone[c] - depthOffset[dep-1]; 753 PetscInt csplitp, cp, val; 754 755 ierr = DMLabelGetValue(label, coldp, &val);CHKERRQ(ierr); 756 if (val == dep-1) { 757 ierr = PetscFindInt(coldp, numSplitPoints[dep-1], splitPoints[dep-1], &cp);CHKERRQ(ierr); 758 if (cp < 0) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "Point %d is not a split point of dimension %d", oldp, dep-1); 759 csplitp = pMaxNew[dep-1] + cp; 760 ierr = DMPlexInsertCone(sdm, newp, c, csplitp);CHKERRQ(ierr); 761 replaced = PETSC_TRUE; 762 } 763 } 764 if (!replaced) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "The cone of point %d does not contain split points", oldp); 765 } 766 ierr = ISRestoreIndices(pIS, &points);CHKERRQ(ierr); 767 ierr = ISDestroy(&pIS);CHKERRQ(ierr); 768 } 769 /* Step 7: Stratify */ 770 ierr = DMPlexStratify(sdm);CHKERRQ(ierr); 771 /* Step 8: Coordinates */ 772 ierr = DMPlexShiftCoordinates_Internal(dm, depthShift, sdm);CHKERRQ(ierr); 773 ierr = DMPlexGetCoordinateSection(sdm, &coordSection);CHKERRQ(ierr); 774 ierr = DMGetCoordinatesLocal(sdm, &coordinates);CHKERRQ(ierr); 775 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 776 for (v = 0; v < (numSplitPoints ? numSplitPoints[0] : 0); ++v) { 777 const PetscInt newp = depthOffset[0] + splitPoints[0][v]; 778 const PetscInt splitp = pMaxNew[0] + v; 779 PetscInt dof, off, soff, d; 780 781 ierr = PetscSectionGetDof(coordSection, newp, &dof);CHKERRQ(ierr); 782 ierr = PetscSectionGetOffset(coordSection, newp, &off);CHKERRQ(ierr); 783 ierr = PetscSectionGetOffset(coordSection, splitp, &soff);CHKERRQ(ierr); 784 for (d = 0; d < dof; ++d) coords[soff+d] = coords[off+d]; 785 } 786 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 787 /* Step 9: SF, if I can figure this out we can split the mesh in parallel */ 788 ierr = DMPlexShiftSF_Internal(dm, depthShift, sdm);CHKERRQ(ierr); 789 /* Step 10: Labels */ 790 ierr = DMPlexShiftLabels_Internal(dm, depthShift, sdm);CHKERRQ(ierr); 791 ierr = DMPlexGetNumLabels(sdm, &numLabels);CHKERRQ(ierr); 792 for (dep = 0; dep <= depth; ++dep) { 793 for (p = 0; p < numSplitPoints[dep]; ++p) { 794 const PetscInt newp = depthOffset[dep] + splitPoints[dep][p]; 795 const PetscInt splitp = pMaxNew[dep] + p; 796 PetscInt l; 797 798 for (l = 0; l < numLabels; ++l) { 799 DMLabel mlabel; 800 const char *lname; 801 PetscInt val; 802 803 ierr = DMPlexGetLabelName(sdm, l, &lname);CHKERRQ(ierr); 804 ierr = DMPlexGetLabel(sdm, lname, &mlabel);CHKERRQ(ierr); 805 ierr = DMLabelGetValue(mlabel, newp, &val);CHKERRQ(ierr); 806 if (val >= 0) { 807 ierr = DMLabelSetValue(mlabel, splitp, val);CHKERRQ(ierr); 808 if (dep == 0) { 809 const PetscInt cedge = pMaxNew[1] + (depthShift[1] - depthShift[0]) + p; 810 ierr = DMLabelSetValue(mlabel, cedge, val);CHKERRQ(ierr); 811 } 812 } 813 } 814 } 815 } 816 for (sp = 0; sp < numSP; ++sp) { 817 const PetscInt dep = values[sp]; 818 819 if ((dep < 0) || (dep > depth)) continue; 820 if (pointIS[dep]) {ierr = ISRestoreIndices(pointIS[dep], &splitPoints[dep]);CHKERRQ(ierr);} 821 ierr = ISDestroy(&pointIS[dep]);CHKERRQ(ierr); 822 } 823 if (label) { 824 ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr); 825 ierr = ISDestroy(&valueIS);CHKERRQ(ierr); 826 } 827 ierr = PetscFree5(depthShift, depthOffset, pMaxNew, coneNew, supportNew);CHKERRQ(ierr); 828 ierr = PetscFree3(pointIS, numSplitPoints, splitPoints);CHKERRQ(ierr); 829 PetscFunctionReturn(0); 830 } 831 832 #undef __FUNCT__ 833 #define __FUNCT__ "DMPlexConstructCohesiveCells" 834 /*@C 835 DMPlexConstructCohesiveCells - Construct cohesive cells which split the face along an internal interface 836 837 Collective on dm 838 839 Input Parameters: 840 + dm - The original DM 841 - labelName - The label specifying the boundary faces (this could be auto-generated) 842 843 Output Parameters: 844 - dmSplit - The new DM 845 846 Level: developer 847 848 .seealso: DMCreate(), DMLabelCohesiveComplete() 849 @*/ 850 PetscErrorCode DMPlexConstructCohesiveCells(DM dm, DMLabel label, DM *dmSplit) 851 { 852 DM sdm; 853 PetscInt dim; 854 PetscErrorCode ierr; 855 856 PetscFunctionBegin; 857 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 858 PetscValidPointer(dmSplit, 4); 859 ierr = DMCreate(PetscObjectComm((PetscObject)dm), &sdm);CHKERRQ(ierr); 860 ierr = DMSetType(sdm, DMPLEX);CHKERRQ(ierr); 861 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 862 ierr = DMPlexSetDimension(sdm, dim);CHKERRQ(ierr); 863 switch (dim) { 864 case 2: 865 case 3: 866 ierr = DMPlexConstructCohesiveCells_Internal(dm, label, sdm);CHKERRQ(ierr); 867 break; 868 default: 869 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct cohesive cells for dimension %d", dim); 870 } 871 *dmSplit = sdm; 872 PetscFunctionReturn(0); 873 } 874 875 #undef __FUNCT__ 876 #define __FUNCT__ "DMLabelCohesiveComplete" 877 /*@ 878 DMLabelCohesiveComplete - Starting with a label marking vertices on an internal surface, we add all other mesh pieces 879 to complete the surface 880 881 Input Parameters: 882 + dm - The DM 883 - label - A DMLabel marking the surface vertices 884 885 Output Parameter: 886 . label - A DMLabel marking all surface points 887 888 Level: developer 889 890 .seealso: DMPlexConstructCohesiveCells() 891 @*/ 892 PetscErrorCode DMLabelCohesiveComplete(DM dm, DMLabel label) 893 { 894 IS dimIS; 895 const PetscInt *points; 896 PetscInt shift = 100, dim, dep, cStart, cEnd, numPoints, p, val; 897 PetscErrorCode ierr; 898 899 PetscFunctionBegin; 900 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 901 /* Cell orientation for face gives the side of the fault */ 902 ierr = DMLabelGetStratumIS(label, dim-1, &dimIS);CHKERRQ(ierr); 903 if (!dimIS) PetscFunctionReturn(0); 904 ierr = ISGetLocalSize(dimIS, &numPoints);CHKERRQ(ierr); 905 ierr = ISGetIndices(dimIS, &points);CHKERRQ(ierr); 906 for (p = 0; p < numPoints; ++p) { 907 const PetscInt *support; 908 PetscInt supportSize, s; 909 910 ierr = DMPlexGetSupportSize(dm, points[p], &supportSize);CHKERRQ(ierr); 911 if (supportSize != 2) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Split face %d has %d != 2 supports", points[p], supportSize); 912 ierr = DMPlexGetSupport(dm, points[p], &support);CHKERRQ(ierr); 913 for (s = 0; s < supportSize; ++s) { 914 const PetscInt *cone, *ornt; 915 PetscInt coneSize, c; 916 PetscBool pos = PETSC_TRUE; 917 918 ierr = DMPlexGetConeSize(dm, support[s], &coneSize);CHKERRQ(ierr); 919 ierr = DMPlexGetCone(dm, support[s], &cone);CHKERRQ(ierr); 920 ierr = DMPlexGetConeOrientation(dm, support[s], &ornt);CHKERRQ(ierr); 921 for (c = 0; c < coneSize; ++c) { 922 if (cone[c] == points[p]) { 923 if (ornt[c] >= 0) { 924 ierr = DMLabelSetValue(label, support[s], shift+dim);CHKERRQ(ierr); 925 } else { 926 ierr = DMLabelSetValue(label, support[s], -(shift+dim));CHKERRQ(ierr); 927 pos = PETSC_FALSE; 928 } 929 break; 930 } 931 } 932 if (c == coneSize) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Cell split face %d support does not have it in the cone", points[p]); 933 /* Put faces touching the fault in the label */ 934 for (c = 0; c < coneSize; ++c) { 935 const PetscInt point = cone[c]; 936 937 ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr); 938 if (val == -1) { 939 PetscInt *closure = NULL; 940 PetscInt closureSize, cl; 941 942 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 943 for (cl = 0; cl < closureSize*2; cl += 2) { 944 const PetscInt clp = closure[cl]; 945 946 ierr = DMLabelGetValue(label, clp, &val);CHKERRQ(ierr); 947 if ((val >= 0) && (val < dim-1)) { 948 ierr = DMLabelSetValue(label, point, pos == PETSC_TRUE ? shift+dim-1 : -(shift+dim-1));CHKERRQ(ierr); 949 break; 950 } 951 } 952 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 953 } 954 } 955 } 956 } 957 ierr = ISRestoreIndices(dimIS, &points);CHKERRQ(ierr); 958 ierr = ISDestroy(&dimIS);CHKERRQ(ierr); 959 /* Search for other cells/faces/edges connected to the fault by a vertex */ 960 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 961 ierr = DMLabelGetStratumIS(label, 0, &dimIS);CHKERRQ(ierr); 962 if (!dimIS) PetscFunctionReturn(0); 963 ierr = ISGetLocalSize(dimIS, &numPoints);CHKERRQ(ierr); 964 ierr = ISGetIndices(dimIS, &points);CHKERRQ(ierr); 965 for (p = 0; p < numPoints; ++p) { 966 PetscInt *star = NULL; 967 PetscInt starSize, s; 968 PetscInt again = 1; /* 0: Finished 1: Keep iterating after a change 2: No change */ 969 970 /* First mark cells connected to the fault */ 971 ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 972 while (again) { 973 if (again > 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Could not classify all cells connected to the fault"); 974 again = 0; 975 for (s = 0; s < starSize*2; s += 2) { 976 const PetscInt point = star[s]; 977 const PetscInt *cone; 978 PetscInt coneSize, c; 979 980 if ((point < cStart) || (point >= cEnd)) continue; 981 ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr); 982 if (val != -1) continue; 983 again = 2; 984 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 985 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 986 for (c = 0; c < coneSize; ++c) { 987 ierr = DMLabelGetValue(label, cone[c], &val);CHKERRQ(ierr); 988 if (val != -1) { 989 if (abs(val) < shift) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Face %d on cell %d has an invalid label %d", cone[c], point, val); 990 if (val > 0) { 991 ierr = DMLabelSetValue(label, point, shift+dim);CHKERRQ(ierr); 992 } else { 993 ierr = DMLabelSetValue(label, point, -(shift+dim));CHKERRQ(ierr); 994 } 995 again = 1; 996 break; 997 } 998 } 999 } 1000 } 1001 /* Classify the rest by cell membership */ 1002 for (s = 0; s < starSize*2; s += 2) { 1003 const PetscInt point = star[s]; 1004 1005 ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr); 1006 if (val == -1) { 1007 PetscInt *sstar = NULL; 1008 PetscInt sstarSize, ss; 1009 PetscBool marked = PETSC_FALSE; 1010 1011 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &sstarSize, &sstar);CHKERRQ(ierr); 1012 for (ss = 0; ss < sstarSize*2; ss += 2) { 1013 const PetscInt spoint = sstar[ss]; 1014 1015 if ((spoint < cStart) || (spoint >= cEnd)) continue; 1016 ierr = DMLabelGetValue(label, spoint, &val);CHKERRQ(ierr); 1017 if (val == -1) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Cell %d in star of %d does not have a valid label", spoint, point); 1018 ierr = DMPlexGetLabelValue(dm, "depth", point, &dep);CHKERRQ(ierr); 1019 if (val > 0) { 1020 ierr = DMLabelSetValue(label, point, shift+dep);CHKERRQ(ierr); 1021 } else { 1022 ierr = DMLabelSetValue(label, point, -(shift+dep));CHKERRQ(ierr); 1023 } 1024 marked = PETSC_TRUE; 1025 break; 1026 } 1027 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &sstarSize, &sstar);CHKERRQ(ierr); 1028 if (!marked) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Point %d could not be classified", point); 1029 } 1030 } 1031 ierr = DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1032 } 1033 ierr = ISRestoreIndices(dimIS, &points);CHKERRQ(ierr); 1034 ierr = ISDestroy(&dimIS);CHKERRQ(ierr); 1035 PetscFunctionReturn(0); 1036 } 1037 1038 #undef __FUNCT__ 1039 #define __FUNCT__ "DMPlexMarkSubmesh_Uninterpolated" 1040 /* Here we need the explicit assumption that: 1041 1042 For any marked cell, the marked vertices constitute a single face 1043 */ 1044 static PetscErrorCode DMPlexMarkSubmesh_Uninterpolated(DM dm, DMLabel vertexLabel, DMLabel subpointMap, PetscInt *numFaces, PetscInt *nFV, DM subdm) 1045 { 1046 IS subvertexIS; 1047 const PetscInt *subvertices; 1048 PetscInt *pStart, *pEnd, *pMax; 1049 PetscInt depth, dim, d, numSubVerticesInitial = 0, v; 1050 PetscErrorCode ierr; 1051 1052 PetscFunctionBegin; 1053 *numFaces = 0; 1054 *nFV = 0; 1055 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1056 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1057 ierr = PetscMalloc3(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd,dim+1,PetscInt,&pMax);CHKERRQ(ierr); 1058 ierr = DMPlexGetHybridBounds(dm, &pMax[depth], depth>1 ? &pMax[depth-1] : NULL, depth > 2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr); 1059 for (d = 0; d <= depth; ++d) { 1060 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1061 if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]); 1062 } 1063 /* Loop over initial vertices and mark all faces in the collective star() */ 1064 ierr = DMLabelGetStratumIS(vertexLabel, 1, &subvertexIS);CHKERRQ(ierr); 1065 if (subvertexIS) { 1066 ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr); 1067 ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1068 } 1069 for (v = 0; v < numSubVerticesInitial; ++v) { 1070 const PetscInt vertex = subvertices[v]; 1071 PetscInt *star = NULL; 1072 PetscInt starSize, s, numCells = 0, c; 1073 1074 ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1075 for (s = 0; s < starSize*2; s += 2) { 1076 const PetscInt point = star[s]; 1077 if ((point >= pStart[depth]) && (point < pEnd[depth])) star[numCells++] = point; 1078 } 1079 for (c = 0; c < numCells; ++c) { 1080 const PetscInt cell = star[c]; 1081 PetscInt *closure = NULL; 1082 PetscInt closureSize, cl; 1083 PetscInt cellLoc, numCorners = 0, faceSize = 0; 1084 1085 ierr = DMLabelGetValue(subpointMap, cell, &cellLoc);CHKERRQ(ierr); 1086 if (cellLoc == 2) continue; 1087 if (cellLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Cell %d has dimension %d in the surface label", cell, cellLoc); 1088 ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1089 for (cl = 0; cl < closureSize*2; cl += 2) { 1090 const PetscInt point = closure[cl]; 1091 PetscInt vertexLoc; 1092 1093 if ((point >= pStart[0]) && (point < pEnd[0])) { 1094 ++numCorners; 1095 ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr); 1096 if (vertexLoc >= 0) closure[faceSize++] = point; 1097 } 1098 } 1099 if (!(*nFV)) {ierr = DMPlexGetNumFaceVertices_Internal(dm, dim, numCorners, nFV);CHKERRQ(ierr);} 1100 if (faceSize > *nFV) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize); 1101 if (faceSize == *nFV) { 1102 ++(*numFaces); 1103 for (cl = 0; cl < faceSize; ++cl) { 1104 ierr = DMLabelSetValue(subpointMap, closure[cl], 0);CHKERRQ(ierr); 1105 } 1106 ierr = DMLabelSetValue(subpointMap, cell, 2);CHKERRQ(ierr); 1107 } 1108 ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1109 } 1110 ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1111 } 1112 if (subvertexIS) { 1113 ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1114 } 1115 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 1116 ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr); 1117 PetscFunctionReturn(0); 1118 } 1119 1120 #undef __FUNCT__ 1121 #define __FUNCT__ "DMPlexMarkSubmesh_Interpolated" 1122 static PetscErrorCode DMPlexMarkSubmesh_Interpolated(DM dm, DMLabel vertexLabel, DMLabel subpointMap, DM subdm) 1123 { 1124 IS subvertexIS; 1125 const PetscInt *subvertices; 1126 PetscInt *pStart, *pEnd, *pMax; 1127 PetscInt dim, d, numSubVerticesInitial = 0, v; 1128 PetscErrorCode ierr; 1129 1130 PetscFunctionBegin; 1131 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1132 ierr = PetscMalloc3(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd,dim+1,PetscInt,&pMax);CHKERRQ(ierr); 1133 ierr = DMPlexGetHybridBounds(dm, &pMax[dim], dim>1 ? &pMax[dim-1] : NULL, dim > 2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr); 1134 for (d = 0; d <= dim; ++d) { 1135 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1136 if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]); 1137 } 1138 /* Loop over initial vertices and mark all faces in the collective star() */ 1139 ierr = DMLabelGetStratumIS(vertexLabel, 1, &subvertexIS);CHKERRQ(ierr); 1140 if (subvertexIS) { 1141 ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr); 1142 ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1143 } 1144 for (v = 0; v < numSubVerticesInitial; ++v) { 1145 const PetscInt vertex = subvertices[v]; 1146 PetscInt *star = NULL; 1147 PetscInt starSize, s, numFaces = 0, f; 1148 1149 ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1150 for (s = 0; s < starSize*2; s += 2) { 1151 const PetscInt point = star[s]; 1152 if ((point >= pStart[dim-1]) && (point < pEnd[dim-1])) star[numFaces++] = point; 1153 } 1154 for (f = 0; f < numFaces; ++f) { 1155 const PetscInt face = star[f]; 1156 PetscInt *closure = NULL; 1157 PetscInt closureSize, c; 1158 PetscInt faceLoc; 1159 1160 ierr = DMLabelGetValue(subpointMap, face, &faceLoc);CHKERRQ(ierr); 1161 if (faceLoc == dim-1) continue; 1162 if (faceLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Face %d has dimension %d in the surface label", face, faceLoc); 1163 ierr = DMPlexGetTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1164 for (c = 0; c < closureSize*2; c += 2) { 1165 const PetscInt point = closure[c]; 1166 PetscInt vertexLoc; 1167 1168 if ((point >= pStart[0]) && (point < pEnd[0])) { 1169 ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr); 1170 if (vertexLoc < 0) break; 1171 } 1172 } 1173 if (c == closureSize*2) { 1174 const PetscInt *support; 1175 PetscInt supportSize, s; 1176 1177 for (c = 0; c < closureSize*2; c += 2) { 1178 const PetscInt point = closure[c]; 1179 1180 for (d = 0; d < dim; ++d) { 1181 if ((point >= pStart[d]) && (point < pEnd[d])) { 1182 ierr = DMLabelSetValue(subpointMap, point, d);CHKERRQ(ierr); 1183 break; 1184 } 1185 } 1186 } 1187 ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 1188 ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); 1189 for (s = 0; s < supportSize; ++s) { 1190 ierr = DMLabelSetValue(subpointMap, support[s], dim);CHKERRQ(ierr); 1191 } 1192 } 1193 ierr = DMPlexRestoreTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1194 } 1195 ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1196 } 1197 if (subvertexIS) { 1198 ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1199 } 1200 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 1201 ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr); 1202 PetscFunctionReturn(0); 1203 } 1204 1205 #undef __FUNCT__ 1206 #define __FUNCT__ "DMPlexGetFaceOrientation" 1207 PetscErrorCode DMPlexGetFaceOrientation(DM dm, PetscInt cell, PetscInt numCorners, PetscInt indices[], PetscInt oppositeVertex, PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented) 1208 { 1209 MPI_Comm comm; 1210 PetscBool posOrient = PETSC_FALSE; 1211 const PetscInt debug = 0; 1212 PetscInt cellDim, faceSize, f; 1213 PetscErrorCode ierr; 1214 1215 PetscFunctionBegin; 1216 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1217 ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr); 1218 if (debug) {PetscPrintf(comm, "cellDim: %d numCorners: %d\n", cellDim, numCorners);CHKERRQ(ierr);} 1219 1220 if (cellDim == numCorners-1) { 1221 /* Simplices */ 1222 faceSize = numCorners-1; 1223 posOrient = !(oppositeVertex%2) ? PETSC_TRUE : PETSC_FALSE; 1224 } else if (cellDim == 1 && numCorners == 3) { 1225 /* Quadratic line */ 1226 faceSize = 1; 1227 posOrient = PETSC_TRUE; 1228 } else if (cellDim == 2 && numCorners == 4) { 1229 /* Quads */ 1230 faceSize = 2; 1231 if ((indices[1] > indices[0]) && (indices[1] - indices[0] == 1)) { 1232 posOrient = PETSC_TRUE; 1233 } else if ((indices[0] == 3) && (indices[1] == 0)) { 1234 posOrient = PETSC_TRUE; 1235 } else { 1236 if (((indices[0] > indices[1]) && (indices[0] - indices[1] == 1)) || ((indices[0] == 0) && (indices[1] == 3))) { 1237 posOrient = PETSC_FALSE; 1238 } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossedge"); 1239 } 1240 } else if (cellDim == 2 && numCorners == 6) { 1241 /* Quadratic triangle (I hate this) */ 1242 /* Edges are determined by the first 2 vertices (corners of edges) */ 1243 const PetscInt faceSizeTri = 3; 1244 PetscInt sortedIndices[3], i, iFace; 1245 PetscBool found = PETSC_FALSE; 1246 PetscInt faceVerticesTriSorted[9] = { 1247 0, 3, 4, /* bottom */ 1248 1, 4, 5, /* right */ 1249 2, 3, 5, /* left */ 1250 }; 1251 PetscInt faceVerticesTri[9] = { 1252 0, 3, 4, /* bottom */ 1253 1, 4, 5, /* right */ 1254 2, 5, 3, /* left */ 1255 }; 1256 1257 faceSize = faceSizeTri; 1258 for (i = 0; i < faceSizeTri; ++i) sortedIndices[i] = indices[i]; 1259 ierr = PetscSortInt(faceSizeTri, sortedIndices);CHKERRQ(ierr); 1260 for (iFace = 0; iFace < 3; ++iFace) { 1261 const PetscInt ii = iFace*faceSizeTri; 1262 PetscInt fVertex, cVertex; 1263 1264 if ((sortedIndices[0] == faceVerticesTriSorted[ii+0]) && 1265 (sortedIndices[1] == faceVerticesTriSorted[ii+1])) { 1266 for (fVertex = 0; fVertex < faceSizeTri; ++fVertex) { 1267 for (cVertex = 0; cVertex < faceSizeTri; ++cVertex) { 1268 if (indices[cVertex] == faceVerticesTri[ii+fVertex]) { 1269 faceVertices[fVertex] = origVertices[cVertex]; 1270 break; 1271 } 1272 } 1273 } 1274 found = PETSC_TRUE; 1275 break; 1276 } 1277 } 1278 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tri crossface"); 1279 if (posOriented) *posOriented = PETSC_TRUE; 1280 PetscFunctionReturn(0); 1281 } else if (cellDim == 2 && numCorners == 9) { 1282 /* Quadratic quad (I hate this) */ 1283 /* Edges are determined by the first 2 vertices (corners of edges) */ 1284 const PetscInt faceSizeQuad = 3; 1285 PetscInt sortedIndices[3], i, iFace; 1286 PetscBool found = PETSC_FALSE; 1287 PetscInt faceVerticesQuadSorted[12] = { 1288 0, 1, 4, /* bottom */ 1289 1, 2, 5, /* right */ 1290 2, 3, 6, /* top */ 1291 0, 3, 7, /* left */ 1292 }; 1293 PetscInt faceVerticesQuad[12] = { 1294 0, 1, 4, /* bottom */ 1295 1, 2, 5, /* right */ 1296 2, 3, 6, /* top */ 1297 3, 0, 7, /* left */ 1298 }; 1299 1300 faceSize = faceSizeQuad; 1301 for (i = 0; i < faceSizeQuad; ++i) sortedIndices[i] = indices[i]; 1302 ierr = PetscSortInt(faceSizeQuad, sortedIndices);CHKERRQ(ierr); 1303 for (iFace = 0; iFace < 4; ++iFace) { 1304 const PetscInt ii = iFace*faceSizeQuad; 1305 PetscInt fVertex, cVertex; 1306 1307 if ((sortedIndices[0] == faceVerticesQuadSorted[ii+0]) && 1308 (sortedIndices[1] == faceVerticesQuadSorted[ii+1])) { 1309 for (fVertex = 0; fVertex < faceSizeQuad; ++fVertex) { 1310 for (cVertex = 0; cVertex < faceSizeQuad; ++cVertex) { 1311 if (indices[cVertex] == faceVerticesQuad[ii+fVertex]) { 1312 faceVertices[fVertex] = origVertices[cVertex]; 1313 break; 1314 } 1315 } 1316 } 1317 found = PETSC_TRUE; 1318 break; 1319 } 1320 } 1321 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossface"); 1322 if (posOriented) *posOriented = PETSC_TRUE; 1323 PetscFunctionReturn(0); 1324 } else if (cellDim == 3 && numCorners == 8) { 1325 /* Hexes 1326 A hex is two oriented quads with the normal of the first 1327 pointing up at the second. 1328 1329 7---6 1330 /| /| 1331 4---5 | 1332 | 3-|-2 1333 |/ |/ 1334 0---1 1335 1336 Faces are determined by the first 4 vertices (corners of faces) */ 1337 const PetscInt faceSizeHex = 4; 1338 PetscInt sortedIndices[4], i, iFace; 1339 PetscBool found = PETSC_FALSE; 1340 PetscInt faceVerticesHexSorted[24] = { 1341 0, 1, 2, 3, /* bottom */ 1342 4, 5, 6, 7, /* top */ 1343 0, 1, 4, 5, /* front */ 1344 1, 2, 5, 6, /* right */ 1345 2, 3, 6, 7, /* back */ 1346 0, 3, 4, 7, /* left */ 1347 }; 1348 PetscInt faceVerticesHex[24] = { 1349 3, 2, 1, 0, /* bottom */ 1350 4, 5, 6, 7, /* top */ 1351 0, 1, 5, 4, /* front */ 1352 1, 2, 6, 5, /* right */ 1353 2, 3, 7, 6, /* back */ 1354 3, 0, 4, 7, /* left */ 1355 }; 1356 1357 faceSize = faceSizeHex; 1358 for (i = 0; i < faceSizeHex; ++i) sortedIndices[i] = indices[i]; 1359 ierr = PetscSortInt(faceSizeHex, sortedIndices);CHKERRQ(ierr); 1360 for (iFace = 0; iFace < 6; ++iFace) { 1361 const PetscInt ii = iFace*faceSizeHex; 1362 PetscInt fVertex, cVertex; 1363 1364 if ((sortedIndices[0] == faceVerticesHexSorted[ii+0]) && 1365 (sortedIndices[1] == faceVerticesHexSorted[ii+1]) && 1366 (sortedIndices[2] == faceVerticesHexSorted[ii+2]) && 1367 (sortedIndices[3] == faceVerticesHexSorted[ii+3])) { 1368 for (fVertex = 0; fVertex < faceSizeHex; ++fVertex) { 1369 for (cVertex = 0; cVertex < faceSizeHex; ++cVertex) { 1370 if (indices[cVertex] == faceVerticesHex[ii+fVertex]) { 1371 faceVertices[fVertex] = origVertices[cVertex]; 1372 break; 1373 } 1374 } 1375 } 1376 found = PETSC_TRUE; 1377 break; 1378 } 1379 } 1380 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface"); 1381 if (posOriented) *posOriented = PETSC_TRUE; 1382 PetscFunctionReturn(0); 1383 } else if (cellDim == 3 && numCorners == 10) { 1384 /* Quadratic tet */ 1385 /* Faces are determined by the first 3 vertices (corners of faces) */ 1386 const PetscInt faceSizeTet = 6; 1387 PetscInt sortedIndices[6], i, iFace; 1388 PetscBool found = PETSC_FALSE; 1389 PetscInt faceVerticesTetSorted[24] = { 1390 0, 1, 2, 6, 7, 8, /* bottom */ 1391 0, 3, 4, 6, 7, 9, /* front */ 1392 1, 4, 5, 7, 8, 9, /* right */ 1393 2, 3, 5, 6, 8, 9, /* left */ 1394 }; 1395 PetscInt faceVerticesTet[24] = { 1396 0, 1, 2, 6, 7, 8, /* bottom */ 1397 0, 4, 3, 6, 7, 9, /* front */ 1398 1, 5, 4, 7, 8, 9, /* right */ 1399 2, 3, 5, 8, 6, 9, /* left */ 1400 }; 1401 1402 faceSize = faceSizeTet; 1403 for (i = 0; i < faceSizeTet; ++i) sortedIndices[i] = indices[i]; 1404 ierr = PetscSortInt(faceSizeTet, sortedIndices);CHKERRQ(ierr); 1405 for (iFace=0; iFace < 4; ++iFace) { 1406 const PetscInt ii = iFace*faceSizeTet; 1407 PetscInt fVertex, cVertex; 1408 1409 if ((sortedIndices[0] == faceVerticesTetSorted[ii+0]) && 1410 (sortedIndices[1] == faceVerticesTetSorted[ii+1]) && 1411 (sortedIndices[2] == faceVerticesTetSorted[ii+2]) && 1412 (sortedIndices[3] == faceVerticesTetSorted[ii+3])) { 1413 for (fVertex = 0; fVertex < faceSizeTet; ++fVertex) { 1414 for (cVertex = 0; cVertex < faceSizeTet; ++cVertex) { 1415 if (indices[cVertex] == faceVerticesTet[ii+fVertex]) { 1416 faceVertices[fVertex] = origVertices[cVertex]; 1417 break; 1418 } 1419 } 1420 } 1421 found = PETSC_TRUE; 1422 break; 1423 } 1424 } 1425 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tet crossface"); 1426 if (posOriented) *posOriented = PETSC_TRUE; 1427 PetscFunctionReturn(0); 1428 } else if (cellDim == 3 && numCorners == 27) { 1429 /* Quadratic hexes (I hate this) 1430 A hex is two oriented quads with the normal of the first 1431 pointing up at the second. 1432 1433 7---6 1434 /| /| 1435 4---5 | 1436 | 3-|-2 1437 |/ |/ 1438 0---1 1439 1440 Faces are determined by the first 4 vertices (corners of faces) */ 1441 const PetscInt faceSizeQuadHex = 9; 1442 PetscInt sortedIndices[9], i, iFace; 1443 PetscBool found = PETSC_FALSE; 1444 PetscInt faceVerticesQuadHexSorted[54] = { 1445 0, 1, 2, 3, 8, 9, 10, 11, 24, /* bottom */ 1446 4, 5, 6, 7, 12, 13, 14, 15, 25, /* top */ 1447 0, 1, 4, 5, 8, 12, 16, 17, 22, /* front */ 1448 1, 2, 5, 6, 9, 13, 17, 18, 21, /* right */ 1449 2, 3, 6, 7, 10, 14, 18, 19, 23, /* back */ 1450 0, 3, 4, 7, 11, 15, 16, 19, 20, /* left */ 1451 }; 1452 PetscInt faceVerticesQuadHex[54] = { 1453 3, 2, 1, 0, 10, 9, 8, 11, 24, /* bottom */ 1454 4, 5, 6, 7, 12, 13, 14, 15, 25, /* top */ 1455 0, 1, 5, 4, 8, 17, 12, 16, 22, /* front */ 1456 1, 2, 6, 5, 9, 18, 13, 17, 21, /* right */ 1457 2, 3, 7, 6, 10, 19, 14, 18, 23, /* back */ 1458 3, 0, 4, 7, 11, 16, 15, 19, 20 /* left */ 1459 }; 1460 1461 faceSize = faceSizeQuadHex; 1462 for (i = 0; i < faceSizeQuadHex; ++i) sortedIndices[i] = indices[i]; 1463 ierr = PetscSortInt(faceSizeQuadHex, sortedIndices);CHKERRQ(ierr); 1464 for (iFace = 0; iFace < 6; ++iFace) { 1465 const PetscInt ii = iFace*faceSizeQuadHex; 1466 PetscInt fVertex, cVertex; 1467 1468 if ((sortedIndices[0] == faceVerticesQuadHexSorted[ii+0]) && 1469 (sortedIndices[1] == faceVerticesQuadHexSorted[ii+1]) && 1470 (sortedIndices[2] == faceVerticesQuadHexSorted[ii+2]) && 1471 (sortedIndices[3] == faceVerticesQuadHexSorted[ii+3])) { 1472 for (fVertex = 0; fVertex < faceSizeQuadHex; ++fVertex) { 1473 for (cVertex = 0; cVertex < faceSizeQuadHex; ++cVertex) { 1474 if (indices[cVertex] == faceVerticesQuadHex[ii+fVertex]) { 1475 faceVertices[fVertex] = origVertices[cVertex]; 1476 break; 1477 } 1478 } 1479 } 1480 found = PETSC_TRUE; 1481 break; 1482 } 1483 } 1484 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface"); 1485 if (posOriented) *posOriented = PETSC_TRUE; 1486 PetscFunctionReturn(0); 1487 } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Unknown cell type for faceOrientation()."); 1488 if (!posOrient) { 1489 if (debug) {ierr = PetscPrintf(comm, " Reversing initial face orientation\n");CHKERRQ(ierr);} 1490 for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[faceSize-1 - f]; 1491 } else { 1492 if (debug) {ierr = PetscPrintf(comm, " Keeping initial face orientation\n");CHKERRQ(ierr);} 1493 for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[f]; 1494 } 1495 if (posOriented) *posOriented = posOrient; 1496 PetscFunctionReturn(0); 1497 } 1498 1499 #undef __FUNCT__ 1500 #define __FUNCT__ "DMPlexGetOrientedFace" 1501 /* 1502 Given a cell and a face, as a set of vertices, 1503 return the oriented face, as a set of vertices, in faceVertices 1504 The orientation is such that the face normal points out of the cell 1505 */ 1506 PetscErrorCode DMPlexGetOrientedFace(DM dm, PetscInt cell, PetscInt faceSize, const PetscInt face[], PetscInt numCorners, PetscInt indices[], PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented) 1507 { 1508 const PetscInt *cone = NULL; 1509 PetscInt coneSize, v, f, v2; 1510 PetscInt oppositeVertex = -1; 1511 PetscErrorCode ierr; 1512 1513 PetscFunctionBegin; 1514 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 1515 ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 1516 for (v = 0, v2 = 0; v < coneSize; ++v) { 1517 PetscBool found = PETSC_FALSE; 1518 1519 for (f = 0; f < faceSize; ++f) { 1520 if (face[f] == cone[v]) { 1521 found = PETSC_TRUE; break; 1522 } 1523 } 1524 if (found) { 1525 indices[v2] = v; 1526 origVertices[v2] = cone[v]; 1527 ++v2; 1528 } else { 1529 oppositeVertex = v; 1530 } 1531 } 1532 ierr = DMPlexGetFaceOrientation(dm, cell, numCorners, indices, oppositeVertex, origVertices, faceVertices, posOriented);CHKERRQ(ierr); 1533 PetscFunctionReturn(0); 1534 } 1535 1536 #undef __FUNCT__ 1537 #define __FUNCT__ "DMPlexInsertFace_Internal" 1538 /* 1539 DMPlexInsertFace_Internal - Puts a face into the mesh 1540 1541 Not collective 1542 1543 Input Parameters: 1544 + dm - The DMPlex 1545 . numFaceVertex - The number of vertices in the face 1546 . faceVertices - The vertices in the face for dm 1547 . subfaceVertices - The vertices in the face for subdm 1548 . numCorners - The number of vertices in the cell 1549 . cell - A cell in dm containing the face 1550 . subcell - A cell in subdm containing the face 1551 . firstFace - First face in the mesh 1552 - newFacePoint - Next face in the mesh 1553 1554 Output Parameters: 1555 . newFacePoint - Contains next face point number on input, updated on output 1556 1557 Level: developer 1558 */ 1559 static PetscErrorCode DMPlexInsertFace_Internal(DM dm, DM subdm, PetscInt numFaceVertices, const PetscInt faceVertices[], const PetscInt subfaceVertices[], PetscInt numCorners, PetscInt cell, PetscInt subcell, PetscInt firstFace, PetscInt *newFacePoint) 1560 { 1561 MPI_Comm comm; 1562 DM_Plex *submesh = (DM_Plex*) subdm->data; 1563 const PetscInt *faces; 1564 PetscInt numFaces, coneSize; 1565 PetscErrorCode ierr; 1566 1567 PetscFunctionBegin; 1568 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1569 ierr = DMPlexGetConeSize(subdm, subcell, &coneSize);CHKERRQ(ierr); 1570 if (coneSize != 1) SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cone size of cell %d is %d != 1", cell, coneSize); 1571 #if 0 1572 /* Cannot use this because support() has not been constructed yet */ 1573 ierr = DMPlexGetJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr); 1574 #else 1575 { 1576 PetscInt f; 1577 1578 numFaces = 0; 1579 ierr = DMGetWorkArray(subdm, 1, PETSC_INT, (void**) &faces);CHKERRQ(ierr); 1580 for (f = firstFace; f < *newFacePoint; ++f) { 1581 PetscInt dof, off, d; 1582 1583 ierr = PetscSectionGetDof(submesh->coneSection, f, &dof);CHKERRQ(ierr); 1584 ierr = PetscSectionGetOffset(submesh->coneSection, f, &off);CHKERRQ(ierr); 1585 /* Yes, I know this is quadratic, but I expect the sizes to be <5 */ 1586 for (d = 0; d < dof; ++d) { 1587 const PetscInt p = submesh->cones[off+d]; 1588 PetscInt v; 1589 1590 for (v = 0; v < numFaceVertices; ++v) { 1591 if (subfaceVertices[v] == p) break; 1592 } 1593 if (v == numFaceVertices) break; 1594 } 1595 if (d == dof) { 1596 numFaces = 1; 1597 ((PetscInt*) faces)[0] = f; 1598 } 1599 } 1600 } 1601 #endif 1602 if (numFaces > 1) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Vertex set had %d faces, not one", numFaces); 1603 else if (numFaces == 1) { 1604 /* Add the other cell neighbor for this face */ 1605 ierr = DMPlexSetCone(subdm, cell, faces);CHKERRQ(ierr); 1606 } else { 1607 PetscInt *indices, *origVertices, *orientedVertices, *orientedSubVertices, v, ov; 1608 PetscBool posOriented; 1609 1610 ierr = DMGetWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr); 1611 origVertices = &orientedVertices[numFaceVertices]; 1612 indices = &orientedVertices[numFaceVertices*2]; 1613 orientedSubVertices = &orientedVertices[numFaceVertices*3]; 1614 ierr = DMPlexGetOrientedFace(dm, cell, numFaceVertices, faceVertices, numCorners, indices, origVertices, orientedVertices, &posOriented);CHKERRQ(ierr); 1615 /* TODO: I know that routine should return a permutation, not the indices */ 1616 for (v = 0; v < numFaceVertices; ++v) { 1617 const PetscInt vertex = faceVertices[v], subvertex = subfaceVertices[v]; 1618 for (ov = 0; ov < numFaceVertices; ++ov) { 1619 if (orientedVertices[ov] == vertex) { 1620 orientedSubVertices[ov] = subvertex; 1621 break; 1622 } 1623 } 1624 if (ov == numFaceVertices) SETERRQ1(comm, PETSC_ERR_PLIB, "Could not find face vertex %d in orientated set", vertex); 1625 } 1626 ierr = DMPlexSetCone(subdm, *newFacePoint, orientedSubVertices);CHKERRQ(ierr); 1627 ierr = DMPlexSetCone(subdm, subcell, newFacePoint);CHKERRQ(ierr); 1628 ierr = DMRestoreWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr); 1629 ++(*newFacePoint); 1630 } 1631 ierr = DMPlexRestoreJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr); 1632 PetscFunctionReturn(0); 1633 } 1634 1635 #undef __FUNCT__ 1636 #define __FUNCT__ "DMPlexCreateSubmesh_Uninterpolated" 1637 static PetscErrorCode DMPlexCreateSubmesh_Uninterpolated(DM dm, const char vertexLabelName[], DM subdm) 1638 { 1639 MPI_Comm comm; 1640 DMLabel vertexLabel, subpointMap; 1641 IS subvertexIS, subcellIS; 1642 const PetscInt *subVertices, *subCells; 1643 PetscInt numSubVertices, firstSubVertex, numSubCells; 1644 PetscInt *subface, maxConeSize, numSubFaces, firstSubFace, newFacePoint, nFV; 1645 PetscInt vStart, vEnd, c, f; 1646 PetscErrorCode ierr; 1647 1648 PetscFunctionBegin; 1649 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1650 /* Create subpointMap which marks the submesh */ 1651 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 1652 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 1653 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 1654 ierr = DMPlexGetLabel(dm, vertexLabelName, &vertexLabel);CHKERRQ(ierr); 1655 ierr = DMPlexMarkSubmesh_Uninterpolated(dm, vertexLabel, subpointMap, &numSubFaces, &nFV, subdm);CHKERRQ(ierr); 1656 /* Setup chart */ 1657 ierr = DMLabelGetStratumSize(subpointMap, 0, &numSubVertices);CHKERRQ(ierr); 1658 ierr = DMLabelGetStratumSize(subpointMap, 2, &numSubCells);CHKERRQ(ierr); 1659 ierr = DMPlexSetChart(subdm, 0, numSubCells+numSubFaces+numSubVertices);CHKERRQ(ierr); 1660 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 1661 /* Set cone sizes */ 1662 firstSubVertex = numSubCells; 1663 firstSubFace = numSubCells+numSubVertices; 1664 newFacePoint = firstSubFace; 1665 ierr = DMLabelGetStratumIS(subpointMap, 0, &subvertexIS);CHKERRQ(ierr); 1666 if (subvertexIS) {ierr = ISGetIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 1667 ierr = DMLabelGetStratumIS(subpointMap, 2, &subcellIS);CHKERRQ(ierr); 1668 if (subcellIS) {ierr = ISGetIndices(subcellIS, &subCells);CHKERRQ(ierr);} 1669 for (c = 0; c < numSubCells; ++c) { 1670 ierr = DMPlexSetConeSize(subdm, c, 1);CHKERRQ(ierr); 1671 } 1672 for (f = firstSubFace; f < firstSubFace+numSubFaces; ++f) { 1673 ierr = DMPlexSetConeSize(subdm, f, nFV);CHKERRQ(ierr); 1674 } 1675 ierr = DMSetUp(subdm);CHKERRQ(ierr); 1676 /* Create face cones */ 1677 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1678 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 1679 ierr = DMGetWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 1680 for (c = 0; c < numSubCells; ++c) { 1681 const PetscInt cell = subCells[c]; 1682 const PetscInt subcell = c; 1683 PetscInt *closure = NULL; 1684 PetscInt closureSize, cl, numCorners = 0, faceSize = 0; 1685 1686 ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1687 for (cl = 0; cl < closureSize*2; cl += 2) { 1688 const PetscInt point = closure[cl]; 1689 PetscInt subVertex; 1690 1691 if ((point >= vStart) && (point < vEnd)) { 1692 ++numCorners; 1693 ierr = PetscFindInt(point, numSubVertices, subVertices, &subVertex);CHKERRQ(ierr); 1694 if (subVertex >= 0) { 1695 closure[faceSize] = point; 1696 subface[faceSize] = firstSubVertex+subVertex; 1697 ++faceSize; 1698 } 1699 } 1700 } 1701 if (faceSize > nFV) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize); 1702 if (faceSize == nFV) { 1703 ierr = DMPlexInsertFace_Internal(dm, subdm, faceSize, closure, subface, numCorners, cell, subcell, firstSubFace, &newFacePoint);CHKERRQ(ierr); 1704 } 1705 ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1706 } 1707 ierr = DMRestoreWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 1708 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 1709 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 1710 /* Build coordinates */ 1711 { 1712 PetscSection coordSection, subCoordSection; 1713 Vec coordinates, subCoordinates; 1714 PetscScalar *coords, *subCoords; 1715 PetscInt coordSize, v; 1716 1717 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1718 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1719 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 1720 ierr = PetscSectionSetChart(subCoordSection, firstSubVertex, firstSubVertex+numSubVertices);CHKERRQ(ierr); 1721 for (v = 0; v < numSubVertices; ++v) { 1722 const PetscInt vertex = subVertices[v]; 1723 const PetscInt subvertex = firstSubVertex+v; 1724 PetscInt dof; 1725 1726 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 1727 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 1728 } 1729 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 1730 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 1731 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 1732 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1733 ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr); 1734 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1735 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 1736 for (v = 0; v < numSubVertices; ++v) { 1737 const PetscInt vertex = subVertices[v]; 1738 const PetscInt subvertex = firstSubVertex+v; 1739 PetscInt dof, off, sdof, soff, d; 1740 1741 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 1742 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 1743 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 1744 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 1745 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 1746 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 1747 } 1748 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1749 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 1750 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 1751 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 1752 } 1753 /* Cleanup */ 1754 if (subvertexIS) {ierr = ISRestoreIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 1755 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 1756 if (subcellIS) {ierr = ISRestoreIndices(subcellIS, &subCells);CHKERRQ(ierr);} 1757 ierr = ISDestroy(&subcellIS);CHKERRQ(ierr); 1758 PetscFunctionReturn(0); 1759 } 1760 1761 #undef __FUNCT__ 1762 #define __FUNCT__ "DMPlexCreateSubmesh_Interpolated" 1763 static PetscErrorCode DMPlexCreateSubmesh_Interpolated(DM dm, const char vertexLabelName[], DM subdm) 1764 { 1765 MPI_Comm comm; 1766 DMLabel subpointMap, vertexLabel; 1767 IS *subpointIS; 1768 const PetscInt **subpoints; 1769 PetscInt *numSubPoints, *firstSubPoint, *coneNew; 1770 PetscInt totSubPoints = 0, maxConeSize, dim, p, d, v; 1771 PetscErrorCode ierr; 1772 1773 PetscFunctionBegin; 1774 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1775 /* Create subpointMap which marks the submesh */ 1776 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 1777 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 1778 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 1779 ierr = DMPlexGetLabel(dm, vertexLabelName, &vertexLabel);CHKERRQ(ierr); 1780 ierr = DMPlexMarkSubmesh_Interpolated(dm, vertexLabel, subpointMap, subdm);CHKERRQ(ierr); 1781 /* Setup chart */ 1782 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1783 ierr = PetscMalloc4(dim+1,PetscInt,&numSubPoints,dim+1,PetscInt,&firstSubPoint,dim+1,IS,&subpointIS,dim+1,const PetscInt *,&subpoints);CHKERRQ(ierr); 1784 for (d = 0; d <= dim; ++d) { 1785 ierr = DMLabelGetStratumSize(subpointMap, d, &numSubPoints[d]);CHKERRQ(ierr); 1786 totSubPoints += numSubPoints[d]; 1787 } 1788 ierr = DMPlexSetChart(subdm, 0, totSubPoints);CHKERRQ(ierr); 1789 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 1790 /* Set cone sizes */ 1791 firstSubPoint[dim] = 0; 1792 firstSubPoint[0] = firstSubPoint[dim] + numSubPoints[dim]; 1793 if (dim > 1) {firstSubPoint[dim-1] = firstSubPoint[0] + numSubPoints[0];} 1794 if (dim > 2) {firstSubPoint[dim-2] = firstSubPoint[dim-1] + numSubPoints[dim-1];} 1795 for (d = 0; d <= dim; ++d) { 1796 ierr = DMLabelGetStratumIS(subpointMap, d, &subpointIS[d]);CHKERRQ(ierr); 1797 ierr = ISGetIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 1798 } 1799 for (d = 0; d <= dim; ++d) { 1800 for (p = 0; p < numSubPoints[d]; ++p) { 1801 const PetscInt point = subpoints[d][p]; 1802 const PetscInt subpoint = firstSubPoint[d] + p; 1803 const PetscInt *cone; 1804 PetscInt coneSize, coneSizeNew, c, val; 1805 1806 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 1807 ierr = DMPlexSetConeSize(subdm, subpoint, coneSize);CHKERRQ(ierr); 1808 if (d == dim) { 1809 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 1810 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 1811 ierr = DMLabelGetValue(subpointMap, cone[c], &val);CHKERRQ(ierr); 1812 if (val >= 0) coneSizeNew++; 1813 } 1814 ierr = DMPlexSetConeSize(subdm, subpoint, coneSizeNew);CHKERRQ(ierr); 1815 } 1816 } 1817 } 1818 ierr = DMSetUp(subdm);CHKERRQ(ierr); 1819 /* Set cones */ 1820 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 1821 ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &coneNew);CHKERRQ(ierr); 1822 for (d = 0; d <= dim; ++d) { 1823 for (p = 0; p < numSubPoints[d]; ++p) { 1824 const PetscInt point = subpoints[d][p]; 1825 const PetscInt subpoint = firstSubPoint[d] + p; 1826 const PetscInt *cone; 1827 PetscInt coneSize, subconeSize, coneSizeNew, c, subc; 1828 1829 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 1830 ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr); 1831 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 1832 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 1833 ierr = PetscFindInt(cone[c], numSubPoints[d-1], subpoints[d-1], &subc);CHKERRQ(ierr); 1834 if (subc >= 0) coneNew[coneSizeNew++] = firstSubPoint[d-1] + subc; 1835 } 1836 if (coneSizeNew != subconeSize) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of cone points located %d does not match subcone size %d", coneSizeNew, subconeSize); 1837 ierr = DMPlexSetCone(subdm, subpoint, coneNew);CHKERRQ(ierr); 1838 } 1839 } 1840 ierr = PetscFree(coneNew);CHKERRQ(ierr); 1841 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 1842 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 1843 /* Build coordinates */ 1844 { 1845 PetscSection coordSection, subCoordSection; 1846 Vec coordinates, subCoordinates; 1847 PetscScalar *coords, *subCoords; 1848 PetscInt coordSize; 1849 1850 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1851 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1852 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 1853 ierr = PetscSectionSetChart(subCoordSection, firstSubPoint[0], firstSubPoint[0]+numSubPoints[0]);CHKERRQ(ierr); 1854 for (v = 0; v < numSubPoints[0]; ++v) { 1855 const PetscInt vertex = subpoints[0][v]; 1856 const PetscInt subvertex = firstSubPoint[0]+v; 1857 PetscInt dof; 1858 1859 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 1860 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 1861 } 1862 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 1863 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 1864 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 1865 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1866 ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr); 1867 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1868 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 1869 for (v = 0; v < numSubPoints[0]; ++v) { 1870 const PetscInt vertex = subpoints[0][v]; 1871 const PetscInt subvertex = firstSubPoint[0]+v; 1872 PetscInt dof, off, sdof, soff, d; 1873 1874 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 1875 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 1876 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 1877 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 1878 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 1879 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 1880 } 1881 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1882 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 1883 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 1884 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 1885 } 1886 /* Cleanup */ 1887 for (d = 0; d <= dim; ++d) { 1888 ierr = ISRestoreIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 1889 ierr = ISDestroy(&subpointIS[d]);CHKERRQ(ierr); 1890 } 1891 ierr = PetscFree4(numSubPoints,firstSubPoint,subpointIS,subpoints);CHKERRQ(ierr); 1892 PetscFunctionReturn(0); 1893 } 1894 1895 #undef __FUNCT__ 1896 #define __FUNCT__ "DMPlexCreateSubmesh" 1897 /* 1898 DMPlexCreateSubmesh - Extract a hypersurface from the mesh using vertices defined by a label 1899 1900 Input Parameters: 1901 + dm - The original mesh 1902 - vertexLabel - The DMLabel marking vertices contained in the surface 1903 1904 Output Parameter: 1905 . subdm - The surface mesh 1906 1907 Note: This function produces a DMLabel mapping original points in the submesh to their depth. This can be obtained using DMPlexGetSubpointMap(). 1908 1909 Level: developer 1910 1911 .seealso: DMPlexGetSubpointMap(), DMPlexGetLabel(), DMLabelSetValue() 1912 */ 1913 PetscErrorCode DMPlexCreateSubmesh(DM dm, const char vertexLabel[], DM *subdm) 1914 { 1915 PetscInt dim, depth; 1916 PetscErrorCode ierr; 1917 1918 PetscFunctionBegin; 1919 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1920 PetscValidCharPointer(vertexLabel, 2); 1921 PetscValidPointer(subdm, 4); 1922 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1923 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1924 ierr = DMCreate(PetscObjectComm((PetscObject)dm), subdm);CHKERRQ(ierr); 1925 ierr = DMSetType(*subdm, DMPLEX);CHKERRQ(ierr); 1926 ierr = DMPlexSetDimension(*subdm, dim-1);CHKERRQ(ierr); 1927 if (depth == dim) { 1928 ierr = DMPlexCreateSubmesh_Interpolated(dm, vertexLabel, *subdm);CHKERRQ(ierr); 1929 } else { 1930 ierr = DMPlexCreateSubmesh_Uninterpolated(dm, vertexLabel, *subdm);CHKERRQ(ierr); 1931 } 1932 PetscFunctionReturn(0); 1933 } 1934 1935 #undef __FUNCT__ 1936 #define __FUNCT__ "DMPlexGetSubpointMap" 1937 PetscErrorCode DMPlexGetSubpointMap(DM dm, DMLabel *subpointMap) 1938 { 1939 DM_Plex *mesh = (DM_Plex*) dm->data; 1940 1941 PetscFunctionBegin; 1942 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1943 PetscValidPointer(subpointMap, 2); 1944 *subpointMap = mesh->subpointMap; 1945 PetscFunctionReturn(0); 1946 } 1947 1948 #undef __FUNCT__ 1949 #define __FUNCT__ "DMPlexSetSubpointMap" 1950 /* Note: Should normally not be called by the user, since it is set in DMPlexCreateSubmesh() */ 1951 PetscErrorCode DMPlexSetSubpointMap(DM dm, DMLabel subpointMap) 1952 { 1953 DM_Plex *mesh = (DM_Plex *) dm->data; 1954 PetscErrorCode ierr; 1955 1956 PetscFunctionBegin; 1957 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1958 ierr = DMLabelDestroy(&mesh->subpointMap);CHKERRQ(ierr); 1959 mesh->subpointMap = subpointMap; 1960 ++mesh->subpointMap->refct;CHKERRQ(ierr); 1961 PetscFunctionReturn(0); 1962 } 1963 1964 #undef __FUNCT__ 1965 #define __FUNCT__ "DMPlexCreateSubpointIS" 1966 /* 1967 DMPlexCreateSubpointIS - Creates an IS covering the entire subdm chart with the original points as data 1968 1969 Input Parameter: 1970 . dm - The submesh DM 1971 1972 Output Parameter: 1973 . subpointIS - The IS of all the points from the original mesh in this submesh, or NULL if this is not a submesh 1974 1975 Note: This is IS is guaranteed to be sorted by the construction of the submesh 1976 */ 1977 PetscErrorCode DMPlexCreateSubpointIS(DM dm, IS *subpointIS) 1978 { 1979 MPI_Comm comm; 1980 DMLabel subpointMap; 1981 IS is; 1982 const PetscInt *opoints; 1983 PetscInt *points, *depths; 1984 PetscInt depth, depStart, depEnd, d, pStart, pEnd, p, n, off; 1985 PetscErrorCode ierr; 1986 1987 PetscFunctionBegin; 1988 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1989 PetscValidPointer(subpointIS, 2); 1990 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1991 *subpointIS = NULL; 1992 ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 1993 if (subpointMap) { 1994 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1995 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 1996 if (pStart) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Submeshes must start the point numbering at 0, not %d", pStart); 1997 ierr = DMGetWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr); 1998 depths[0] = depth; 1999 depths[1] = 0; 2000 for(d = 2; d <= depth; ++d) {depths[d] = depth+1 - d;} 2001 ierr = PetscMalloc(pEnd * sizeof(PetscInt), &points);CHKERRQ(ierr); 2002 for(d = 0, off = 0; d <= depth; ++d) { 2003 const PetscInt dep = depths[d]; 2004 2005 ierr = DMPlexGetDepthStratum(dm, dep, &depStart, &depEnd);CHKERRQ(ierr); 2006 ierr = DMLabelGetStratumSize(subpointMap, dep, &n);CHKERRQ(ierr); 2007 if (((d < 2) && (depth > 1)) || (d == 1)) { /* Only check vertices and cells for now since the map is broken for others */ 2008 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); 2009 } else { 2010 if (!n) { 2011 if (d == 0) { 2012 /* Missing cells */ 2013 for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = -1; 2014 } else { 2015 /* Missing faces */ 2016 for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = PETSC_MAX_INT; 2017 } 2018 } 2019 } 2020 if (n) { 2021 ierr = DMLabelGetStratumIS(subpointMap, dep, &is);CHKERRQ(ierr); 2022 ierr = ISGetIndices(is, &opoints);CHKERRQ(ierr); 2023 for(p = 0; p < n; ++p, ++off) points[off] = opoints[p]; 2024 ierr = ISRestoreIndices(is, &opoints);CHKERRQ(ierr); 2025 ierr = ISDestroy(&is);CHKERRQ(ierr); 2026 } 2027 } 2028 ierr = DMRestoreWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr); 2029 if (off != pEnd) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "The number of mapped submesh points %d should be %d", off, pEnd); 2030 ierr = ISCreateGeneral(comm, pEnd, points, PETSC_OWN_POINTER, subpointIS);CHKERRQ(ierr); 2031 } 2032 PetscFunctionReturn(0); 2033 } 2034