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