1 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <petscsf.h> 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "DMPlexMarkBoundaryFaces" 6 /*@ 7 DMPlexMarkBoundaryFaces - Mark all faces on the boundary 8 9 Not Collective 10 11 Input Parameter: 12 . dm - The original DM 13 14 Output Parameter: 15 . label - The DMLabel marking boundary faces with value 1 16 17 Level: developer 18 19 .seealso: DMLabelCreate(), DMPlexCreateLabel() 20 @*/ 21 PetscErrorCode DMPlexMarkBoundaryFaces(DM dm, DMLabel label) 22 { 23 PetscInt fStart, fEnd, f; 24 PetscErrorCode ierr; 25 26 PetscFunctionBegin; 27 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 28 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 29 for (f = fStart; f < fEnd; ++f) { 30 PetscInt supportSize; 31 32 ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr); 33 if (supportSize == 1) { 34 ierr = DMLabelSetValue(label, f, 1);CHKERRQ(ierr); 35 } 36 } 37 PetscFunctionReturn(0); 38 } 39 40 #undef __FUNCT__ 41 #define __FUNCT__ "DMPlexLabelComplete" 42 /*@ 43 DMPlexLabelComplete - Starting with a label marking points on a surface, we add the transitive closure to the surface 44 45 Input Parameters: 46 + dm - The DM 47 - label - A DMLabel marking the surface points 48 49 Output Parameter: 50 . label - A DMLabel marking all surface points in the transitive closure 51 52 Level: developer 53 54 .seealso: DMPlexLabelCohesiveComplete() 55 @*/ 56 PetscErrorCode DMPlexLabelComplete(DM dm, DMLabel label) 57 { 58 IS valueIS; 59 const PetscInt *values; 60 PetscInt numValues, v; 61 PetscErrorCode ierr; 62 63 PetscFunctionBegin; 64 ierr = DMLabelGetNumValues(label, &numValues);CHKERRQ(ierr); 65 ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr); 66 ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr); 67 for (v = 0; v < numValues; ++v) { 68 IS pointIS; 69 const PetscInt *points; 70 PetscInt numPoints, p; 71 72 ierr = DMLabelGetStratumSize(label, values[v], &numPoints);CHKERRQ(ierr); 73 ierr = DMLabelGetStratumIS(label, values[v], &pointIS);CHKERRQ(ierr); 74 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 75 for (p = 0; p < numPoints; ++p) { 76 PetscInt *closure = NULL; 77 PetscInt closureSize, c; 78 79 ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 80 for (c = 0; c < closureSize*2; c += 2) { 81 ierr = DMLabelSetValue(label, closure[c], values[v]);CHKERRQ(ierr); 82 } 83 ierr = DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 84 } 85 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 86 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 87 } 88 ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr); 89 ierr = ISDestroy(&valueIS);CHKERRQ(ierr); 90 PetscFunctionReturn(0); 91 } 92 93 #undef __FUNCT__ 94 #define __FUNCT__ "DMPlexShiftPoint_Internal" 95 PETSC_STATIC_INLINE PetscInt DMPlexShiftPoint_Internal(PetscInt p, PetscInt depth, PetscInt depthEnd[], PetscInt depthShift[]) 96 { 97 if (depth < 0) return p; 98 /* Cells */ if (p < depthEnd[depth]) return p; 99 /* Vertices */ if (p < depthEnd[0]) return p + depthShift[depth]; 100 /* Faces */ if (p < depthEnd[depth-1]) return p + depthShift[depth] + depthShift[0]; 101 /* Edges */ return p + depthShift[depth] + depthShift[0] + depthShift[depth-1]; 102 } 103 104 #undef __FUNCT__ 105 #define __FUNCT__ "DMPlexShiftSizes_Internal" 106 static PetscErrorCode DMPlexShiftSizes_Internal(DM dm, PetscInt depthShift[], DM dmNew) 107 { 108 PetscInt *depthEnd; 109 PetscInt depth = 0, d, pStart, pEnd, p; 110 PetscErrorCode ierr; 111 112 PetscFunctionBegin; 113 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 114 if (depth < 0) PetscFunctionReturn(0); 115 ierr = PetscMalloc((depth+1) * sizeof(PetscInt), &depthEnd);CHKERRQ(ierr); 116 /* Step 1: Expand chart */ 117 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 118 for (d = 0; d <= depth; ++d) { 119 pEnd += depthShift[d]; 120 ierr = DMPlexGetDepthStratum(dm, d, NULL, &depthEnd[d]);CHKERRQ(ierr); 121 } 122 ierr = DMPlexSetChart(dmNew, pStart, pEnd);CHKERRQ(ierr); 123 /* Step 2: Set cone and support sizes */ 124 for (d = 0; d <= depth; ++d) { 125 ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 126 for (p = pStart; p < pEnd; ++p) { 127 PetscInt newp = DMPlexShiftPoint_Internal(p, depth, depthEnd, depthShift); 128 PetscInt size; 129 130 ierr = DMPlexGetConeSize(dm, p, &size);CHKERRQ(ierr); 131 ierr = DMPlexSetConeSize(dmNew, newp, size);CHKERRQ(ierr); 132 ierr = DMPlexGetSupportSize(dm, p, &size);CHKERRQ(ierr); 133 ierr = DMPlexSetSupportSize(dmNew, newp, size);CHKERRQ(ierr); 134 } 135 } 136 ierr = PetscFree(depthEnd);CHKERRQ(ierr); 137 PetscFunctionReturn(0); 138 } 139 140 #undef __FUNCT__ 141 #define __FUNCT__ "DMPlexShiftPoints_Internal" 142 static PetscErrorCode DMPlexShiftPoints_Internal(DM dm, PetscInt depthShift[], DM dmNew) 143 { 144 PetscInt *depthEnd, *newpoints; 145 PetscInt depth = 0, d, maxConeSize, maxSupportSize, pStart, pEnd, p; 146 PetscErrorCode ierr; 147 148 PetscFunctionBegin; 149 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 150 if (depth < 0) PetscFunctionReturn(0); 151 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 152 ierr = PetscMalloc2(depth+1,PetscInt,&depthEnd,PetscMax(maxConeSize, maxSupportSize),PetscInt,&newpoints);CHKERRQ(ierr); 153 for (d = 0; d <= depth; ++d) { 154 ierr = DMPlexGetDepthStratum(dm, d, NULL, &depthEnd[d]);CHKERRQ(ierr); 155 } 156 /* Step 5: Set cones and supports */ 157 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 158 for (p = pStart; p < pEnd; ++p) { 159 const PetscInt *points = NULL, *orientations = NULL; 160 PetscInt size, i, newp = DMPlexShiftPoint_Internal(p, depth, depthEnd, depthShift); 161 162 ierr = DMPlexGetConeSize(dm, p, &size);CHKERRQ(ierr); 163 ierr = DMPlexGetCone(dm, p, &points);CHKERRQ(ierr); 164 ierr = DMPlexGetConeOrientation(dm, p, &orientations);CHKERRQ(ierr); 165 for (i = 0; i < size; ++i) { 166 newpoints[i] = DMPlexShiftPoint_Internal(points[i], depth, depthEnd, depthShift); 167 } 168 ierr = DMPlexSetCone(dmNew, newp, newpoints);CHKERRQ(ierr); 169 ierr = DMPlexSetConeOrientation(dmNew, newp, orientations);CHKERRQ(ierr); 170 ierr = DMPlexGetSupportSize(dm, p, &size);CHKERRQ(ierr); 171 ierr = DMPlexGetSupport(dm, p, &points);CHKERRQ(ierr); 172 for (i = 0; i < size; ++i) { 173 newpoints[i] = DMPlexShiftPoint_Internal(points[i], depth, depthEnd, depthShift); 174 } 175 ierr = DMPlexSetSupport(dmNew, newp, newpoints);CHKERRQ(ierr); 176 } 177 ierr = PetscFree2(depthEnd,newpoints);CHKERRQ(ierr); 178 PetscFunctionReturn(0); 179 } 180 181 #undef __FUNCT__ 182 #define __FUNCT__ "DMPlexShiftCoordinates_Internal" 183 static PetscErrorCode DMPlexShiftCoordinates_Internal(DM dm, PetscInt depthShift[], DM dmNew) 184 { 185 PetscSection coordSection, newCoordSection; 186 Vec coordinates, newCoordinates; 187 PetscScalar *coords, *newCoords; 188 PetscInt *depthEnd, coordSize; 189 PetscInt dim, depth = 0, d, vStart, vEnd, vStartNew, vEndNew, v; 190 PetscErrorCode ierr; 191 192 PetscFunctionBegin; 193 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 194 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 195 ierr = PetscMalloc((depth+1) * sizeof(PetscInt), &depthEnd);CHKERRQ(ierr); 196 for (d = 0; d <= depth; ++d) { 197 ierr = DMPlexGetDepthStratum(dm, d, NULL, &depthEnd[d]);CHKERRQ(ierr); 198 } 199 /* Step 8: Convert coordinates */ 200 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 201 ierr = DMPlexGetDepthStratum(dmNew, 0, &vStartNew, &vEndNew);CHKERRQ(ierr); 202 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 203 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &newCoordSection);CHKERRQ(ierr); 204 ierr = PetscSectionSetNumFields(newCoordSection, 1);CHKERRQ(ierr); 205 ierr = PetscSectionSetFieldComponents(newCoordSection, 0, dim);CHKERRQ(ierr); 206 ierr = PetscSectionSetChart(newCoordSection, vStartNew, vEndNew);CHKERRQ(ierr); 207 for (v = vStartNew; v < vEndNew; ++v) { 208 ierr = PetscSectionSetDof(newCoordSection, v, dim);CHKERRQ(ierr); 209 ierr = PetscSectionSetFieldDof(newCoordSection, v, 0, dim);CHKERRQ(ierr); 210 } 211 ierr = PetscSectionSetUp(newCoordSection);CHKERRQ(ierr); 212 ierr = DMPlexSetCoordinateSection(dmNew, newCoordSection);CHKERRQ(ierr); 213 ierr = PetscSectionGetStorageSize(newCoordSection, &coordSize);CHKERRQ(ierr); 214 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &newCoordinates);CHKERRQ(ierr); 215 ierr = PetscObjectSetName((PetscObject) newCoordinates, "coordinates");CHKERRQ(ierr); 216 ierr = VecSetSizes(newCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 217 ierr = VecSetType(newCoordinates,dm->vectype);CHKERRQ(ierr); 218 ierr = DMSetCoordinatesLocal(dmNew, newCoordinates);CHKERRQ(ierr); 219 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 220 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 221 ierr = VecGetArray(newCoordinates, &newCoords);CHKERRQ(ierr); 222 for (v = vStart; v < vEnd; ++v) { 223 PetscInt dof, off, noff, d; 224 225 ierr = PetscSectionGetDof(coordSection, v, &dof);CHKERRQ(ierr); 226 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 227 ierr = PetscSectionGetOffset(newCoordSection, DMPlexShiftPoint_Internal(v, depth, depthEnd, depthShift), &noff);CHKERRQ(ierr); 228 for (d = 0; d < dof; ++d) { 229 newCoords[noff+d] = coords[off+d]; 230 } 231 } 232 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 233 ierr = VecRestoreArray(newCoordinates, &newCoords);CHKERRQ(ierr); 234 ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr); 235 ierr = PetscSectionDestroy(&newCoordSection);CHKERRQ(ierr); 236 ierr = PetscFree(depthEnd);CHKERRQ(ierr); 237 PetscFunctionReturn(0); 238 } 239 240 #undef __FUNCT__ 241 #define __FUNCT__ "DMPlexShiftSF_Internal" 242 static PetscErrorCode DMPlexShiftSF_Internal(DM dm, PetscInt depthShift[], DM dmNew) 243 { 244 PetscInt *depthEnd; 245 PetscInt depth = 0, d; 246 PetscSF sfPoint, sfPointNew; 247 const PetscSFNode *remotePoints; 248 PetscSFNode *gremotePoints; 249 const PetscInt *localPoints; 250 PetscInt *glocalPoints, *newLocation, *newRemoteLocation; 251 PetscInt numRoots, numLeaves, l, pStart, pEnd, totShift = 0; 252 PetscErrorCode ierr; 253 254 PetscFunctionBegin; 255 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 256 ierr = PetscMalloc((depth+1) * sizeof(PetscInt), &depthEnd);CHKERRQ(ierr); 257 for (d = 0; d <= depth; ++d) { 258 totShift += depthShift[d]; 259 ierr = DMPlexGetDepthStratum(dm, d, NULL, &depthEnd[d]);CHKERRQ(ierr); 260 } 261 /* Step 9: Convert pointSF */ 262 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 263 ierr = DMGetPointSF(dmNew, &sfPointNew);CHKERRQ(ierr); 264 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 265 ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 266 if (numRoots >= 0) { 267 ierr = PetscMalloc2(numRoots,PetscInt,&newLocation,pEnd-pStart,PetscInt,&newRemoteLocation);CHKERRQ(ierr); 268 for (l=0; l<numRoots; l++) newLocation[l] = DMPlexShiftPoint_Internal(l, depth, depthEnd, depthShift); 269 ierr = PetscSFBcastBegin(sfPoint, MPIU_INT, newLocation, newRemoteLocation);CHKERRQ(ierr); 270 ierr = PetscSFBcastEnd(sfPoint, MPIU_INT, newLocation, newRemoteLocation);CHKERRQ(ierr); 271 ierr = PetscMalloc(numLeaves * sizeof(PetscInt), &glocalPoints);CHKERRQ(ierr); 272 ierr = PetscMalloc(numLeaves * sizeof(PetscSFNode), &gremotePoints);CHKERRQ(ierr); 273 for (l = 0; l < numLeaves; ++l) { 274 glocalPoints[l] = DMPlexShiftPoint_Internal(localPoints[l], depth, depthEnd, depthShift); 275 gremotePoints[l].rank = remotePoints[l].rank; 276 gremotePoints[l].index = newRemoteLocation[localPoints[l]]; 277 } 278 ierr = PetscFree2(newLocation,newRemoteLocation);CHKERRQ(ierr); 279 ierr = PetscSFSetGraph(sfPointNew, numRoots + totShift, numLeaves, glocalPoints, PETSC_OWN_POINTER, gremotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 280 } 281 ierr = PetscFree(depthEnd);CHKERRQ(ierr); 282 PetscFunctionReturn(0); 283 } 284 285 #undef __FUNCT__ 286 #define __FUNCT__ "DMPlexShiftLabels_Internal" 287 static PetscErrorCode DMPlexShiftLabels_Internal(DM dm, PetscInt depthShift[], DM dmNew) 288 { 289 PetscSF sfPoint; 290 DMLabel vtkLabel, ghostLabel; 291 PetscInt *depthEnd; 292 const PetscSFNode *leafRemote; 293 const PetscInt *leafLocal; 294 PetscInt depth = 0, d, numLeaves, numLabels, l, cStart, cEnd, c, fStart, fEnd, f; 295 PetscMPIInt rank; 296 PetscErrorCode ierr; 297 298 PetscFunctionBegin; 299 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 300 ierr = PetscMalloc((depth+1) * sizeof(PetscInt), &depthEnd);CHKERRQ(ierr); 301 for (d = 0; d <= depth; ++d) { 302 ierr = DMPlexGetDepthStratum(dm, d, NULL, &depthEnd[d]);CHKERRQ(ierr); 303 } 304 /* Step 10: Convert labels */ 305 ierr = DMPlexGetNumLabels(dm, &numLabels);CHKERRQ(ierr); 306 for (l = 0; l < numLabels; ++l) { 307 DMLabel label, newlabel; 308 const char *lname; 309 PetscBool isDepth; 310 IS valueIS; 311 const PetscInt *values; 312 PetscInt numValues, val; 313 314 ierr = DMPlexGetLabelName(dm, l, &lname);CHKERRQ(ierr); 315 ierr = PetscStrcmp(lname, "depth", &isDepth);CHKERRQ(ierr); 316 if (isDepth) continue; 317 ierr = DMPlexCreateLabel(dmNew, lname);CHKERRQ(ierr); 318 ierr = DMPlexGetLabel(dm, lname, &label);CHKERRQ(ierr); 319 ierr = DMPlexGetLabel(dmNew, lname, &newlabel);CHKERRQ(ierr); 320 ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr); 321 ierr = ISGetLocalSize(valueIS, &numValues);CHKERRQ(ierr); 322 ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr); 323 for (val = 0; val < numValues; ++val) { 324 IS pointIS; 325 const PetscInt *points; 326 PetscInt numPoints, p; 327 328 ierr = DMLabelGetStratumIS(label, values[val], &pointIS);CHKERRQ(ierr); 329 ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr); 330 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 331 for (p = 0; p < numPoints; ++p) { 332 const PetscInt newpoint = DMPlexShiftPoint_Internal(points[p], depth, depthEnd, depthShift); 333 334 ierr = DMLabelSetValue(newlabel, newpoint, values[val]);CHKERRQ(ierr); 335 } 336 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 337 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 338 } 339 ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr); 340 ierr = ISDestroy(&valueIS);CHKERRQ(ierr); 341 } 342 ierr = PetscFree(depthEnd);CHKERRQ(ierr); 343 /* Step 11: Make label for output (vtk) and to mark ghost points (ghost) */ 344 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 345 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 346 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 347 ierr = PetscSFGetGraph(sfPoint, NULL, &numLeaves, &leafLocal, &leafRemote);CHKERRQ(ierr); 348 ierr = DMPlexCreateLabel(dmNew, "vtk");CHKERRQ(ierr); 349 ierr = DMPlexCreateLabel(dmNew, "ghost");CHKERRQ(ierr); 350 ierr = DMPlexGetLabel(dmNew, "vtk", &vtkLabel);CHKERRQ(ierr); 351 ierr = DMPlexGetLabel(dmNew, "ghost", &ghostLabel);CHKERRQ(ierr); 352 for (l = 0, c = cStart; l < numLeaves && c < cEnd; ++l, ++c) { 353 for (; c < leafLocal[l] && c < cEnd; ++c) { 354 ierr = DMLabelSetValue(vtkLabel, c, 1);CHKERRQ(ierr); 355 } 356 if (leafLocal[l] >= cEnd) break; 357 if (leafRemote[l].rank == rank) { 358 ierr = DMLabelSetValue(vtkLabel, c, 1);CHKERRQ(ierr); 359 } else { 360 ierr = DMLabelSetValue(ghostLabel, c, 2);CHKERRQ(ierr); 361 } 362 } 363 for (; c < cEnd; ++c) { 364 ierr = DMLabelSetValue(vtkLabel, c, 1);CHKERRQ(ierr); 365 } 366 if (0) { 367 ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr); 368 ierr = DMLabelView(vtkLabel, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 369 ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 370 } 371 ierr = DMPlexGetHeightStratum(dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr); 372 for (f = fStart; f < fEnd; ++f) { 373 PetscInt numCells; 374 375 ierr = DMPlexGetSupportSize(dmNew, f, &numCells);CHKERRQ(ierr); 376 if (numCells < 2) { 377 ierr = DMLabelSetValue(ghostLabel, f, 1);CHKERRQ(ierr); 378 } else { 379 const PetscInt *cells = NULL; 380 PetscInt vA, vB; 381 382 ierr = DMPlexGetSupport(dmNew, f, &cells);CHKERRQ(ierr); 383 ierr = DMLabelGetValue(vtkLabel, cells[0], &vA);CHKERRQ(ierr); 384 ierr = DMLabelGetValue(vtkLabel, cells[1], &vB);CHKERRQ(ierr); 385 if (!vA && !vB) {ierr = DMLabelSetValue(ghostLabel, f, 1);CHKERRQ(ierr);} 386 } 387 } 388 if (0) { 389 ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr); 390 ierr = DMLabelView(ghostLabel, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 391 ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 392 } 393 PetscFunctionReturn(0); 394 } 395 396 #undef __FUNCT__ 397 #define __FUNCT__ "DMPlexConstructGhostCells_Internal" 398 static PetscErrorCode DMPlexConstructGhostCells_Internal(DM dm, DMLabel label, PetscInt *numGhostCells, DM gdm) 399 { 400 IS valueIS; 401 const PetscInt *values; 402 PetscInt *depthShift; 403 PetscInt depth = 0, numFS, fs, ghostCell, cEnd, c; 404 PetscErrorCode ierr; 405 406 PetscFunctionBegin; 407 /* Count ghost cells */ 408 ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr); 409 ierr = ISGetLocalSize(valueIS, &numFS);CHKERRQ(ierr); 410 ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr); 411 412 *numGhostCells = 0; 413 for (fs = 0; fs < numFS; ++fs) { 414 PetscInt numBdFaces; 415 416 ierr = DMLabelGetStratumSize(label, values[fs], &numBdFaces);CHKERRQ(ierr); 417 418 *numGhostCells += numBdFaces; 419 } 420 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 421 ierr = PetscMalloc((depth+1) * sizeof(PetscInt), &depthShift);CHKERRQ(ierr); 422 ierr = PetscMemzero(depthShift, (depth+1) * sizeof(PetscInt));CHKERRQ(ierr); 423 if (depth >= 0) depthShift[depth] = *numGhostCells; 424 ierr = DMPlexShiftSizes_Internal(dm, depthShift, gdm);CHKERRQ(ierr); 425 /* Step 3: Set cone/support sizes for new points */ 426 ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr); 427 for (c = cEnd; c < cEnd + *numGhostCells; ++c) { 428 ierr = DMPlexSetConeSize(gdm, c, 1);CHKERRQ(ierr); 429 } 430 for (fs = 0; fs < numFS; ++fs) { 431 IS faceIS; 432 const PetscInt *faces; 433 PetscInt numFaces, f; 434 435 ierr = DMLabelGetStratumIS(label, values[fs], &faceIS);CHKERRQ(ierr); 436 ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr); 437 ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr); 438 for (f = 0; f < numFaces; ++f) { 439 PetscInt size; 440 441 ierr = DMPlexGetSupportSize(dm, faces[f], &size);CHKERRQ(ierr); 442 if (size != 1) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM has boundary face %d with %d support cells", faces[f], size); 443 ierr = DMPlexSetSupportSize(gdm, faces[f] + *numGhostCells, 2);CHKERRQ(ierr); 444 } 445 ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 446 ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 447 } 448 /* Step 4: Setup ghosted DM */ 449 ierr = DMSetUp(gdm);CHKERRQ(ierr); 450 ierr = DMPlexShiftPoints_Internal(dm, depthShift, gdm);CHKERRQ(ierr); 451 /* Step 6: Set cones and supports for new points */ 452 ghostCell = cEnd; 453 for (fs = 0; fs < numFS; ++fs) { 454 IS faceIS; 455 const PetscInt *faces; 456 PetscInt numFaces, f; 457 458 ierr = DMLabelGetStratumIS(label, values[fs], &faceIS);CHKERRQ(ierr); 459 ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr); 460 ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr); 461 for (f = 0; f < numFaces; ++f, ++ghostCell) { 462 PetscInt newFace = faces[f] + *numGhostCells; 463 464 ierr = DMPlexSetCone(gdm, ghostCell, &newFace);CHKERRQ(ierr); 465 ierr = DMPlexInsertSupport(gdm, newFace, 1, ghostCell);CHKERRQ(ierr); 466 } 467 ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 468 ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 469 } 470 ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr); 471 ierr = ISDestroy(&valueIS);CHKERRQ(ierr); 472 /* Step 7: Stratify */ 473 ierr = DMPlexStratify(gdm);CHKERRQ(ierr); 474 ierr = DMPlexShiftCoordinates_Internal(dm, depthShift, gdm);CHKERRQ(ierr); 475 ierr = DMPlexShiftSF_Internal(dm, depthShift, gdm);CHKERRQ(ierr); 476 ierr = DMPlexShiftLabels_Internal(dm, depthShift, gdm);CHKERRQ(ierr); 477 ierr = PetscFree(depthShift);CHKERRQ(ierr); 478 PetscFunctionReturn(0); 479 } 480 481 #undef __FUNCT__ 482 #define __FUNCT__ "DMPlexConstructGhostCells" 483 /*@C 484 DMPlexConstructGhostCells - Construct ghost cells which connect to every boundary face 485 486 Collective on dm 487 488 Input Parameters: 489 + dm - The original DM 490 - labelName - The label specifying the boundary faces, or "Face Sets" if this is NULL 491 492 Output Parameters: 493 + numGhostCells - The number of ghost cells added to the DM 494 - dmGhosted - The new DM 495 496 Note: If no label exists of that name, one will be created marking all boundary faces 497 498 Level: developer 499 500 .seealso: DMCreate() 501 */ 502 PetscErrorCode DMPlexConstructGhostCells(DM dm, const char labelName[], PetscInt *numGhostCells, DM *dmGhosted) 503 { 504 DM gdm; 505 DMLabel label; 506 const char *name = labelName ? labelName : "Face Sets"; 507 PetscInt dim; 508 PetscErrorCode ierr; 509 510 PetscFunctionBegin; 511 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 512 PetscValidPointer(numGhostCells, 3); 513 PetscValidPointer(dmGhosted, 4); 514 ierr = DMCreate(PetscObjectComm((PetscObject)dm), &gdm);CHKERRQ(ierr); 515 ierr = DMSetType(gdm, DMPLEX);CHKERRQ(ierr); 516 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 517 ierr = DMPlexSetDimension(gdm, dim);CHKERRQ(ierr); 518 ierr = DMPlexGetLabel(dm, name, &label);CHKERRQ(ierr); 519 if (!label) { 520 /* Get label for boundary faces */ 521 ierr = DMPlexCreateLabel(dm, name);CHKERRQ(ierr); 522 ierr = DMPlexGetLabel(dm, name, &label);CHKERRQ(ierr); 523 ierr = DMPlexMarkBoundaryFaces(dm, label);CHKERRQ(ierr); 524 } 525 ierr = DMPlexConstructGhostCells_Internal(dm, label, numGhostCells, gdm);CHKERRQ(ierr); 526 ierr = DMSetFromOptions(gdm);CHKERRQ(ierr); 527 *dmGhosted = gdm; 528 PetscFunctionReturn(0); 529 } 530 531 #undef __FUNCT__ 532 #define __FUNCT__ "DMPlexConstructCohesiveCells_Internal" 533 /* 534 We are adding two kinds of points here: 535 Replicated: Copies of points which exist in the mesh, such as vertices identified across a fault 536 Hybrid: Entirely new points, such as cohesive cells 537 */ 538 static PetscErrorCode DMPlexConstructCohesiveCells_Internal(DM dm, DMLabel label, DM sdm) 539 { 540 MPI_Comm comm; 541 IS valueIS; 542 PetscInt numSP = 0; /* The number of depths for which we have replicated points */ 543 const PetscInt *values; /* List of depths for which we have replicated points */ 544 IS *pointIS; 545 PetscInt *numSplitPoints; /* The number of replicated points at each depth */ 546 PetscInt *numHybridPoints; /* The number of hybrid points at each depth */ 547 const PetscInt **splitPoints; /* Replicated points for each depth */ 548 PetscSection coordSection; 549 Vec coordinates; 550 PetscScalar *coords; 551 PetscInt depths[4]; /* Depths in the order that plex points are numbered */ 552 PetscInt *depthShift; /* Number of replicated+hybrid points at each depth */ 553 PetscInt *depthOffset; /* Prefix sums of depthShift */ 554 PetscInt *pMaxNew; /* The first replicated point at each depth in the new mesh, hybrids come after this */ 555 PetscInt *coneNew, *coneONew, *supportNew; 556 PetscInt shift = 100, depth = 0, dep, dim, d, sp, maxConeSize, maxSupportSize, numLabels, vStart, vEnd, pEnd, p, v; 557 PetscErrorCode ierr; 558 559 PetscFunctionBegin; 560 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 561 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 562 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 563 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 564 depths[0] = depth; 565 depths[1] = 0; 566 depths[2] = depth-1; 567 depths[3] = 1; 568 /* Count split points and add cohesive cells */ 569 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 570 ierr = PetscMalloc6(depth+1,PetscInt,&depthShift,depth+1,PetscInt,&depthOffset,depth+1,PetscInt,&pMaxNew,maxConeSize*3,PetscInt,&coneNew,maxConeSize*3,PetscInt,&coneONew,maxSupportSize,PetscInt,&supportNew);CHKERRQ(ierr); 571 ierr = PetscMalloc4(depth+1,IS,&pointIS,depth+1,PetscInt,&numSplitPoints,depth+1,PetscInt,&numHybridPoints,depth+1,const PetscInt*,&splitPoints);CHKERRQ(ierr); 572 ierr = PetscMemzero(depthShift, (depth+1) * sizeof(PetscInt));CHKERRQ(ierr); 573 ierr = PetscMemzero(depthOffset, (depth+1) * sizeof(PetscInt));CHKERRQ(ierr); 574 for (d = 0; d <= depth; ++d) { 575 ierr = DMPlexGetDepthStratum(dm, d, NULL, &pMaxNew[d]);CHKERRQ(ierr); 576 numSplitPoints[d] = 0; 577 numHybridPoints[d] = 0; 578 splitPoints[d] = NULL; 579 pointIS[d] = NULL; 580 } 581 if (label) { 582 ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr); 583 ierr = ISGetLocalSize(valueIS, &numSP);CHKERRQ(ierr); 584 ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr); 585 } 586 for (sp = 0; sp < numSP; ++sp) { 587 const PetscInt dep = values[sp]; 588 589 if ((dep < 0) || (dep > depth)) continue; 590 ierr = DMLabelGetStratumIS(label, dep, &pointIS[dep]);CHKERRQ(ierr); 591 if (pointIS[dep]) { 592 ierr = ISGetLocalSize(pointIS[dep], &numSplitPoints[dep]);CHKERRQ(ierr); 593 ierr = ISGetIndices(pointIS[dep], &splitPoints[dep]);CHKERRQ(ierr); 594 } 595 } 596 /* Calculate number of hybrid points */ 597 for (d = 1; d <= depth; ++d) numHybridPoints[d] = numSplitPoints[d-1]; /* There is a hybrid cell/face/edge for every split face/edge/vertex */ 598 for (d = 0; d <= depth; ++d) depthShift[d] = numSplitPoints[d] + numHybridPoints[d]; 599 for (d = 1; d <= depth; ++d) depthOffset[depths[d]] = depthOffset[depths[d-1]] + depthShift[depths[d-1]]; 600 for (d = 0; d <= depth; ++d) pMaxNew[d] += depthOffset[d]; 601 ierr = DMPlexShiftSizes_Internal(dm, depthShift, sdm);CHKERRQ(ierr); 602 /* Step 3: Set cone/support sizes for new points */ 603 for (dep = 0; dep <= depth; ++dep) { 604 for (p = 0; p < numSplitPoints[dep]; ++p) { 605 const PetscInt oldp = splitPoints[dep][p]; 606 const PetscInt newp = oldp + depthOffset[dep]; 607 const PetscInt splitp = p + pMaxNew[dep]; 608 const PetscInt *support; 609 PetscInt coneSize, supportSize, qf, qn, qp, e; 610 611 ierr = DMPlexGetConeSize(dm, oldp, &coneSize);CHKERRQ(ierr); 612 ierr = DMPlexSetConeSize(sdm, splitp, coneSize);CHKERRQ(ierr); 613 ierr = DMPlexGetSupportSize(dm, oldp, &supportSize);CHKERRQ(ierr); 614 ierr = DMPlexSetSupportSize(sdm, splitp, supportSize);CHKERRQ(ierr); 615 if (dep == depth-1) { 616 const PetscInt hybcell = p + pMaxNew[dep+1] + numSplitPoints[dep+1]; 617 618 /* Add cohesive cells, they are prisms */ 619 ierr = DMPlexSetConeSize(sdm, hybcell, 2 + coneSize);CHKERRQ(ierr); 620 } else if (dep == 0) { 621 const PetscInt hybedge = p + pMaxNew[dep+1] + numSplitPoints[dep+1]; 622 623 ierr = DMPlexGetSupport(dm, oldp, &support);CHKERRQ(ierr); 624 for (e = 0, qn = 0, qp = 0, qf = 0; e < supportSize; ++e) { 625 PetscInt val; 626 627 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 628 if (val == 1) ++qf; 629 if ((val == 1) || (val == (shift + 1))) ++qn; 630 if ((val == 1) || (val == -(shift + 1))) ++qp; 631 } 632 /* Split old vertex: Edges into original vertex and new cohesive edge */ 633 ierr = DMPlexSetSupportSize(sdm, newp, qn+1);CHKERRQ(ierr); 634 /* Split new vertex: Edges into split vertex and new cohesive edge */ 635 ierr = DMPlexSetSupportSize(sdm, splitp, qp+1);CHKERRQ(ierr); 636 /* Add hybrid edge */ 637 ierr = DMPlexSetConeSize(sdm, hybedge, 2);CHKERRQ(ierr); 638 ierr = DMPlexSetSupportSize(sdm, hybedge, qf);CHKERRQ(ierr); 639 } else if (dep == dim-2) { 640 const PetscInt hybface = p + pMaxNew[dep+1] + numSplitPoints[dep+1]; 641 642 ierr = DMPlexGetSupport(dm, oldp, &support);CHKERRQ(ierr); 643 for (e = 0, qn = 0, qp = 0, qf = 0; e < supportSize; ++e) { 644 PetscInt val; 645 646 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 647 if (val == dim-1) ++qf; 648 if ((val == dim-1) || (val == (shift + dim-1))) ++qn; 649 if ((val == dim-1) || (val == -(shift + dim-1))) ++qp; 650 } 651 /* Split old edge: Faces into original edge and cohesive face (positive side?) */ 652 ierr = DMPlexSetSupportSize(sdm, newp, qn+1);CHKERRQ(ierr); 653 /* Split new edge: Faces into split edge and cohesive face (negative side?) */ 654 ierr = DMPlexSetSupportSize(sdm, splitp, qp+1);CHKERRQ(ierr); 655 /* Add hybrid face */ 656 ierr = DMPlexSetConeSize(sdm, hybface, 4);CHKERRQ(ierr); 657 ierr = DMPlexSetSupportSize(sdm, hybface, qf);CHKERRQ(ierr); 658 } 659 } 660 } 661 /* Step 4: Setup split DM */ 662 ierr = DMSetUp(sdm);CHKERRQ(ierr); 663 ierr = DMPlexShiftPoints_Internal(dm, depthShift, sdm);CHKERRQ(ierr); 664 /* Step 6: Set cones and supports for new points */ 665 for (dep = 0; dep <= depth; ++dep) { 666 for (p = 0; p < numSplitPoints[dep]; ++p) { 667 const PetscInt oldp = splitPoints[dep][p]; 668 const PetscInt newp = oldp + depthOffset[dep]; 669 const PetscInt splitp = p + pMaxNew[dep]; 670 const PetscInt *cone, *support, *ornt; 671 PetscInt coneSize, supportSize, q, qf, qn, qp, v, e, s; 672 673 ierr = DMPlexGetConeSize(dm, oldp, &coneSize);CHKERRQ(ierr); 674 ierr = DMPlexGetCone(dm, oldp, &cone);CHKERRQ(ierr); 675 ierr = DMPlexGetConeOrientation(dm, oldp, &ornt);CHKERRQ(ierr); 676 ierr = DMPlexGetSupportSize(dm, oldp, &supportSize);CHKERRQ(ierr); 677 ierr = DMPlexGetSupport(dm, oldp, &support);CHKERRQ(ierr); 678 if (dep == depth-1) { 679 const PetscInt hybcell = p + pMaxNew[dep+1] + numSplitPoints[dep+1]; 680 const PetscInt *supportF; 681 682 /* Split face: copy in old face to new face to start */ 683 ierr = DMPlexGetSupport(sdm, newp, &supportF);CHKERRQ(ierr); 684 ierr = DMPlexSetSupport(sdm, splitp, supportF);CHKERRQ(ierr); 685 /* Split old face: old vertices/edges in cone so no change */ 686 /* Split new face: new vertices/edges in cone */ 687 for (q = 0; q < coneSize; ++q) { 688 ierr = PetscFindInt(cone[q], numSplitPoints[dep-1], splitPoints[dep-1], &v);CHKERRQ(ierr); 689 690 coneNew[2+q] = v + pMaxNew[dep-1]; 691 } 692 ierr = DMPlexSetCone(sdm, splitp, &coneNew[2]);CHKERRQ(ierr); 693 ierr = DMPlexSetConeOrientation(sdm, splitp, ornt);CHKERRQ(ierr); 694 /* Face support */ 695 for (s = 0; s < supportSize; ++s) { 696 PetscInt val; 697 698 ierr = DMLabelGetValue(label, support[s], &val);CHKERRQ(ierr); 699 if (val < 0) { 700 /* Split old face: Replace negative side cell with cohesive cell */ 701 ierr = DMPlexInsertSupport(sdm, newp, s, hybcell);CHKERRQ(ierr); 702 } else { 703 /* Split new face: Replace positive side cell with cohesive cell */ 704 ierr = DMPlexInsertSupport(sdm, splitp, s, hybcell);CHKERRQ(ierr); 705 /* Get orientation for cohesive face */ 706 { 707 const PetscInt *ncone, *nconeO; 708 PetscInt nconeSize, nc; 709 710 ierr = DMPlexGetConeSize(dm, support[s], &nconeSize);CHKERRQ(ierr); 711 ierr = DMPlexGetCone(dm, support[s], &ncone);CHKERRQ(ierr); 712 ierr = DMPlexGetConeOrientation(dm, support[s], &nconeO);CHKERRQ(ierr); 713 for (nc = 0; nc < nconeSize; ++nc) { 714 if (ncone[nc] == oldp) { 715 coneONew[0] = nconeO[nc]; 716 break; 717 } 718 } 719 if (nc >= nconeSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not locate face %d in neighboring cell %d", oldp, support[s]); 720 } 721 } 722 } 723 /* Cohesive cell: Old and new split face, then new cohesive faces */ 724 coneNew[0] = newp; /* Extracted negative side orientation above */ 725 coneNew[1] = splitp; 726 coneONew[1] = coneONew[0]; 727 for (q = 0; q < coneSize; ++q) { 728 coneNew[2+q] += (pMaxNew[dep] - pMaxNew[dep-1]) + numSplitPoints[dep]; 729 coneONew[2+q] = 0; 730 } 731 ierr = DMPlexSetCone(sdm, hybcell, coneNew);CHKERRQ(ierr); 732 ierr = DMPlexSetConeOrientation(sdm, hybcell, coneONew);CHKERRQ(ierr); 733 } else if (dep == 0) { 734 const PetscInt hybedge = p + pMaxNew[dep+1] + numSplitPoints[dep+1]; 735 736 /* Split old vertex: Edges in old split faces and new cohesive edge */ 737 for (e = 0, qn = 0; e < supportSize; ++e) { 738 PetscInt val; 739 740 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 741 if ((val == 1) || (val == (shift + 1))) { 742 supportNew[qn++] = support[e] + depthOffset[dep+1]; 743 } 744 } 745 supportNew[qn] = hybedge; 746 ierr = DMPlexSetSupport(sdm, newp, supportNew);CHKERRQ(ierr); 747 /* Split new vertex: Edges in new split faces and new cohesive edge */ 748 for (e = 0, qp = 0; e < supportSize; ++e) { 749 PetscInt val, edge; 750 751 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 752 if (val == 1) { 753 ierr = PetscFindInt(support[e], numSplitPoints[dep+1], splitPoints[dep+1], &edge);CHKERRQ(ierr); 754 if (edge < 0) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Edge %d is not a split edge", support[e]); 755 supportNew[qp++] = edge + pMaxNew[dep+1]; 756 } else if (val == -(shift + 1)) { 757 supportNew[qp++] = support[e] + depthOffset[dep+1]; 758 } 759 } 760 supportNew[qp] = hybedge; 761 ierr = DMPlexSetSupport(sdm, splitp, supportNew);CHKERRQ(ierr); 762 /* Hybrid edge: Old and new split vertex */ 763 coneNew[0] = newp; 764 coneNew[1] = splitp; 765 ierr = DMPlexSetCone(sdm, hybedge, coneNew);CHKERRQ(ierr); 766 for (e = 0, qf = 0; e < supportSize; ++e) { 767 PetscInt val, edge; 768 769 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 770 if (val == 1) { 771 ierr = PetscFindInt(support[e], numSplitPoints[dep+1], splitPoints[dep+1], &edge);CHKERRQ(ierr); 772 if (edge < 0) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Edge %d is not a split edge", support[e]); 773 supportNew[qf++] = edge + pMaxNew[dep+2] + numSplitPoints[dep+2]; 774 } 775 } 776 ierr = DMPlexSetSupport(sdm, hybedge, supportNew);CHKERRQ(ierr); 777 } else if (dep == dim-2) { 778 const PetscInt hybface = p + pMaxNew[dep+1] + numSplitPoints[dep+1]; 779 780 /* Split old edge: old vertices in cone so no change */ 781 /* Split new edge: new vertices in cone */ 782 for (q = 0; q < coneSize; ++q) { 783 ierr = PetscFindInt(cone[q], numSplitPoints[dep-1], splitPoints[dep-1], &v);CHKERRQ(ierr); 784 785 coneNew[q] = v + pMaxNew[dep-1]; 786 } 787 ierr = DMPlexSetCone(sdm, splitp, coneNew);CHKERRQ(ierr); 788 /* Split old edge: Faces in positive side cells and old split faces */ 789 for (e = 0, q = 0; e < supportSize; ++e) { 790 PetscInt val; 791 792 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 793 if (val == dim-1) { 794 supportNew[q++] = support[e] + depthOffset[dep+1]; 795 } else if (val == (shift + dim-1)) { 796 supportNew[q++] = support[e] + depthOffset[dep+1]; 797 } 798 } 799 supportNew[q++] = p + pMaxNew[dep+1] + numSplitPoints[dep+1]; 800 ierr = DMPlexSetSupport(sdm, newp, supportNew);CHKERRQ(ierr); 801 /* Split new edge: Faces in negative side cells and new split faces */ 802 for (e = 0, q = 0; e < supportSize; ++e) { 803 PetscInt val, face; 804 805 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 806 if (val == dim-1) { 807 ierr = PetscFindInt(support[e], numSplitPoints[dep+1], splitPoints[dep+1], &face);CHKERRQ(ierr); 808 if (face < 0) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Face %d is not a split face", support[e]); 809 supportNew[q++] = face + pMaxNew[dep+1]; 810 } else if (val == -(shift + dim-1)) { 811 supportNew[q++] = support[e] + depthOffset[dep+1]; 812 } 813 } 814 supportNew[q++] = p + pMaxNew[dep+1] + numSplitPoints[dep+1]; 815 ierr = DMPlexSetSupport(sdm, splitp, supportNew);CHKERRQ(ierr); 816 /* Hybrid face */ 817 coneNew[0] = newp; 818 coneNew[1] = splitp; 819 for (v = 0; v < coneSize; ++v) { 820 PetscInt vertex; 821 ierr = PetscFindInt(cone[v], numSplitPoints[dep-1], splitPoints[dep-1], &vertex);CHKERRQ(ierr); 822 if (vertex < 0) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Vertex %d is not a split vertex", support[e]); 823 coneNew[2+v] = vertex + pMaxNew[dep] + numSplitPoints[dep]; 824 } 825 ierr = DMPlexSetCone(sdm, hybface, coneNew);CHKERRQ(ierr); 826 for (e = 0, qf = 0; e < supportSize; ++e) { 827 PetscInt val, face; 828 829 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 830 if (val == dim-1) { 831 ierr = PetscFindInt(support[e], numSplitPoints[dep+1], splitPoints[dep+1], &face);CHKERRQ(ierr); 832 if (face < 0) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Face %d is not a split face", support[e]); 833 supportNew[qf++] = face + pMaxNew[dep+2] + numSplitPoints[dep+2]; 834 } 835 } 836 ierr = DMPlexSetSupport(sdm, hybface, supportNew);CHKERRQ(ierr); 837 } 838 } 839 } 840 /* Step 6b: Replace split points in negative side cones */ 841 for (sp = 0; sp < numSP; ++sp) { 842 PetscInt dep = values[sp]; 843 IS pIS; 844 PetscInt numPoints; 845 const PetscInt *points; 846 847 if (dep >= 0) continue; 848 ierr = DMLabelGetStratumIS(label, dep, &pIS);CHKERRQ(ierr); 849 if (!pIS) continue; 850 dep = -dep - shift; 851 ierr = ISGetLocalSize(pIS, &numPoints);CHKERRQ(ierr); 852 ierr = ISGetIndices(pIS, &points);CHKERRQ(ierr); 853 for (p = 0; p < numPoints; ++p) { 854 const PetscInt oldp = points[p]; 855 const PetscInt newp = depthOffset[dep] + oldp; 856 const PetscInt *cone; 857 PetscInt coneSize, c; 858 /* PetscBool replaced = PETSC_FALSE; */ 859 860 /* Negative edge: replace split vertex */ 861 /* Negative cell: replace split face */ 862 ierr = DMPlexGetConeSize(sdm, newp, &coneSize);CHKERRQ(ierr); 863 ierr = DMPlexGetCone(sdm, newp, &cone);CHKERRQ(ierr); 864 for (c = 0; c < coneSize; ++c) { 865 const PetscInt coldp = cone[c] - depthOffset[dep-1]; 866 PetscInt csplitp, cp, val; 867 868 ierr = DMLabelGetValue(label, coldp, &val);CHKERRQ(ierr); 869 if (val == dep-1) { 870 ierr = PetscFindInt(coldp, numSplitPoints[dep-1], splitPoints[dep-1], &cp);CHKERRQ(ierr); 871 if (cp < 0) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "Point %d is not a split point of dimension %d", oldp, dep-1); 872 csplitp = pMaxNew[dep-1] + cp; 873 ierr = DMPlexInsertCone(sdm, newp, c, csplitp);CHKERRQ(ierr); 874 /* replaced = PETSC_TRUE; */ 875 } 876 } 877 /* Cells with only a vertex or edge on the submesh have no replacement */ 878 /* if (!replaced) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "The cone of point %d does not contain split points", oldp); */ 879 } 880 ierr = ISRestoreIndices(pIS, &points);CHKERRQ(ierr); 881 ierr = ISDestroy(&pIS);CHKERRQ(ierr); 882 } 883 /* Step 7: Stratify */ 884 ierr = DMPlexStratify(sdm);CHKERRQ(ierr); 885 /* Step 8: Coordinates */ 886 ierr = DMPlexShiftCoordinates_Internal(dm, depthShift, sdm);CHKERRQ(ierr); 887 ierr = DMPlexGetCoordinateSection(sdm, &coordSection);CHKERRQ(ierr); 888 ierr = DMGetCoordinatesLocal(sdm, &coordinates);CHKERRQ(ierr); 889 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 890 for (v = 0; v < (numSplitPoints ? numSplitPoints[0] : 0); ++v) { 891 const PetscInt newp = depthOffset[0] + splitPoints[0][v]; 892 const PetscInt splitp = pMaxNew[0] + v; 893 PetscInt dof, off, soff, d; 894 895 ierr = PetscSectionGetDof(coordSection, newp, &dof);CHKERRQ(ierr); 896 ierr = PetscSectionGetOffset(coordSection, newp, &off);CHKERRQ(ierr); 897 ierr = PetscSectionGetOffset(coordSection, splitp, &soff);CHKERRQ(ierr); 898 for (d = 0; d < dof; ++d) coords[soff+d] = coords[off+d]; 899 } 900 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 901 /* Step 9: SF, if I can figure this out we can split the mesh in parallel */ 902 ierr = DMPlexShiftSF_Internal(dm, depthShift, sdm);CHKERRQ(ierr); 903 /* Step 10: Labels */ 904 ierr = DMPlexShiftLabels_Internal(dm, depthShift, sdm);CHKERRQ(ierr); 905 ierr = DMPlexGetNumLabels(sdm, &numLabels);CHKERRQ(ierr); 906 for (dep = 0; dep <= depth; ++dep) { 907 for (p = 0; p < numSplitPoints[dep]; ++p) { 908 const PetscInt newp = depthOffset[dep] + splitPoints[dep][p]; 909 const PetscInt splitp = pMaxNew[dep] + p; 910 PetscInt l; 911 912 for (l = 0; l < numLabels; ++l) { 913 DMLabel mlabel; 914 const char *lname; 915 PetscInt val; 916 PetscBool isDepth; 917 918 ierr = DMPlexGetLabelName(sdm, l, &lname);CHKERRQ(ierr); 919 ierr = PetscStrcmp(lname, "depth", &isDepth);CHKERRQ(ierr); 920 if (isDepth) continue; 921 ierr = DMPlexGetLabel(sdm, lname, &mlabel);CHKERRQ(ierr); 922 ierr = DMLabelGetValue(mlabel, newp, &val);CHKERRQ(ierr); 923 if (val >= 0) { 924 ierr = DMLabelSetValue(mlabel, splitp, val);CHKERRQ(ierr); 925 #if 0 926 /* Do not put cohesive edges into the label */ 927 if (dep == 0) { 928 const PetscInt cedge = p + pMaxNew[dep+1] + numSplitPoints[dep+1]; 929 ierr = DMLabelSetValue(mlabel, cedge, val);CHKERRQ(ierr); 930 } else if (dep == dim-2) { 931 const PetscInt cface = p + pMaxNew[dep+1] + numSplitPoints[dep+1]; 932 ierr = DMLabelSetValue(mlabel, cface, val);CHKERRQ(ierr); 933 } 934 /* Do not put cohesive faces into the label */ 935 #endif 936 } 937 } 938 } 939 } 940 for (sp = 0; sp < numSP; ++sp) { 941 const PetscInt dep = values[sp]; 942 943 if ((dep < 0) || (dep > depth)) continue; 944 if (pointIS[dep]) {ierr = ISRestoreIndices(pointIS[dep], &splitPoints[dep]);CHKERRQ(ierr);} 945 ierr = ISDestroy(&pointIS[dep]);CHKERRQ(ierr); 946 } 947 if (label) { 948 ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr); 949 ierr = ISDestroy(&valueIS);CHKERRQ(ierr); 950 } 951 for (d = 0; d <= depth; ++d) { 952 ierr = DMPlexGetDepthStratum(sdm, d, NULL, &pEnd);CHKERRQ(ierr); 953 pMaxNew[d] = pEnd - numHybridPoints[d]; 954 } 955 ierr = DMPlexSetHybridBounds(sdm, depth >= 0 ? pMaxNew[depth] : PETSC_DETERMINE, depth>1 ? pMaxNew[depth-1] : PETSC_DETERMINE, depth>2 ? pMaxNew[1] : PETSC_DETERMINE, depth >= 0 ? pMaxNew[0] : PETSC_DETERMINE);CHKERRQ(ierr); 956 ierr = PetscFree6(depthShift, depthOffset, pMaxNew, coneNew, coneONew, supportNew);CHKERRQ(ierr); 957 ierr = PetscFree4(pointIS, numSplitPoints, numHybridPoints, splitPoints);CHKERRQ(ierr); 958 PetscFunctionReturn(0); 959 } 960 961 #undef __FUNCT__ 962 #define __FUNCT__ "DMPlexConstructCohesiveCells" 963 /*@C 964 DMPlexConstructCohesiveCells - Construct cohesive cells which split the face along an internal interface 965 966 Collective on dm 967 968 Input Parameters: 969 + dm - The original DM 970 - label - The label specifying the boundary faces (this could be auto-generated) 971 972 Output Parameters: 973 - dmSplit - The new DM 974 975 Level: developer 976 977 .seealso: DMCreate(), DMPlexLabelCohesiveComplete() 978 @*/ 979 PetscErrorCode DMPlexConstructCohesiveCells(DM dm, DMLabel label, DM *dmSplit) 980 { 981 DM sdm; 982 PetscInt dim; 983 PetscErrorCode ierr; 984 985 PetscFunctionBegin; 986 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 987 PetscValidPointer(dmSplit, 3); 988 ierr = DMCreate(PetscObjectComm((PetscObject)dm), &sdm);CHKERRQ(ierr); 989 ierr = DMSetType(sdm, DMPLEX);CHKERRQ(ierr); 990 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 991 ierr = DMPlexSetDimension(sdm, dim);CHKERRQ(ierr); 992 switch (dim) { 993 case 2: 994 case 3: 995 ierr = DMPlexConstructCohesiveCells_Internal(dm, label, sdm);CHKERRQ(ierr); 996 break; 997 default: 998 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct cohesive cells for dimension %d", dim); 999 } 1000 *dmSplit = sdm; 1001 PetscFunctionReturn(0); 1002 } 1003 1004 #undef __FUNCT__ 1005 #define __FUNCT__ "DMPlexLabelCohesiveComplete" 1006 /*@ 1007 DMPlexLabelCohesiveComplete - Starting with a label marking vertices on an internal surface, we add all other mesh pieces 1008 to complete the surface 1009 1010 Input Parameters: 1011 + dm - The DM 1012 . label - A DMLabel marking the surface vertices 1013 . flip - Flag to flip the submesh normal and replace points on the other side 1014 - subdm - The subDM associated with the label, or NULL 1015 1016 Output Parameter: 1017 . label - A DMLabel marking all surface points 1018 1019 Level: developer 1020 1021 .seealso: DMPlexConstructCohesiveCells(), DMPlexLabelComplete() 1022 @*/ 1023 PetscErrorCode DMPlexLabelCohesiveComplete(DM dm, DMLabel label, PetscBool flip, DM subdm) 1024 { 1025 DMLabel depthLabel; 1026 IS dimIS, subpointIS; 1027 const PetscInt *points, *subpoints; 1028 const PetscInt rev = flip ? -1 : 1; 1029 PetscInt shift = 100, dim, dep, cStart, cEnd, numPoints, numSubpoints, p, val; 1030 PetscErrorCode ierr; 1031 1032 PetscFunctionBegin; 1033 ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 1034 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1035 if (subdm) { 1036 ierr = DMPlexCreateSubpointIS(subdm, &subpointIS);CHKERRQ(ierr); 1037 if (subpointIS) { 1038 ierr = ISGetLocalSize(subpointIS, &numSubpoints);CHKERRQ(ierr); 1039 ierr = ISGetIndices(subpointIS, &subpoints);CHKERRQ(ierr); 1040 } 1041 } 1042 /* Cell orientation for face gives the side of the fault */ 1043 ierr = DMLabelGetStratumIS(label, dim-1, &dimIS);CHKERRQ(ierr); 1044 if (!dimIS) PetscFunctionReturn(0); 1045 ierr = ISGetLocalSize(dimIS, &numPoints);CHKERRQ(ierr); 1046 ierr = ISGetIndices(dimIS, &points);CHKERRQ(ierr); 1047 for (p = 0; p < numPoints; ++p) { 1048 const PetscInt *support; 1049 PetscInt supportSize, s; 1050 1051 ierr = DMPlexGetSupportSize(dm, points[p], &supportSize);CHKERRQ(ierr); 1052 if (supportSize != 2) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Split face %d has %d != 2 supports", points[p], supportSize); 1053 ierr = DMPlexGetSupport(dm, points[p], &support);CHKERRQ(ierr); 1054 for (s = 0; s < supportSize; ++s) { 1055 const PetscInt *cone, *ornt; 1056 PetscInt coneSize, c; 1057 PetscBool pos = PETSC_TRUE; 1058 1059 ierr = DMPlexGetConeSize(dm, support[s], &coneSize);CHKERRQ(ierr); 1060 ierr = DMPlexGetCone(dm, support[s], &cone);CHKERRQ(ierr); 1061 ierr = DMPlexGetConeOrientation(dm, support[s], &ornt);CHKERRQ(ierr); 1062 for (c = 0; c < coneSize; ++c) { 1063 if (cone[c] == points[p]) { 1064 PetscInt o = ornt[c]; 1065 1066 if (subdm) { 1067 const PetscInt *subcone, *subornt; 1068 PetscInt subpoint, subface, subconeSize, sc; 1069 1070 ierr = PetscFindInt(support[s], numSubpoints, subpoints, &subpoint);CHKERRQ(ierr); 1071 ierr = PetscFindInt(points[p], numSubpoints, subpoints, &subface);CHKERRQ(ierr); 1072 ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr); 1073 ierr = DMPlexGetCone(subdm, subpoint, &subcone);CHKERRQ(ierr); 1074 ierr = DMPlexGetConeOrientation(subdm, subpoint, &subornt);CHKERRQ(ierr); 1075 for (sc = 0; sc < subconeSize; ++sc) { 1076 if (subcone[sc] == subface) { 1077 o = subornt[0]; 1078 break; 1079 } 1080 } 1081 if (sc >= subconeSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find point %d in cone for subpoint %d", points[p], subpoint); 1082 } 1083 if (o >= 0) { 1084 ierr = DMLabelSetValue(label, support[s], rev*(shift+dim));CHKERRQ(ierr); 1085 pos = rev > 0 ? PETSC_TRUE : PETSC_FALSE; 1086 } else { 1087 ierr = DMLabelSetValue(label, support[s], -rev*(shift+dim));CHKERRQ(ierr); 1088 pos = rev > 0 ? PETSC_FALSE : PETSC_TRUE; 1089 } 1090 break; 1091 } 1092 } 1093 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]); 1094 /* Put faces touching the fault in the label */ 1095 for (c = 0; c < coneSize; ++c) { 1096 const PetscInt point = cone[c]; 1097 1098 ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr); 1099 if (val == -1) { 1100 PetscInt *closure = NULL; 1101 PetscInt closureSize, cl; 1102 1103 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1104 for (cl = 0; cl < closureSize*2; cl += 2) { 1105 const PetscInt clp = closure[cl]; 1106 1107 ierr = DMLabelGetValue(label, clp, &val);CHKERRQ(ierr); 1108 if ((val >= 0) && (val < dim-1)) { 1109 ierr = DMLabelSetValue(label, point, pos == PETSC_TRUE ? shift+dim-1 : -(shift+dim-1));CHKERRQ(ierr); 1110 break; 1111 } 1112 } 1113 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1114 } 1115 } 1116 } 1117 } 1118 if (subdm) { 1119 if (subpointIS) {ierr = ISRestoreIndices(subpointIS, &subpoints);CHKERRQ(ierr);} 1120 ierr = ISDestroy(&subpointIS);CHKERRQ(ierr); 1121 } 1122 ierr = ISRestoreIndices(dimIS, &points);CHKERRQ(ierr); 1123 ierr = ISDestroy(&dimIS);CHKERRQ(ierr); 1124 /* Search for other cells/faces/edges connected to the fault by a vertex */ 1125 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1126 ierr = DMLabelGetStratumIS(label, 0, &dimIS);CHKERRQ(ierr); 1127 if (!dimIS) PetscFunctionReturn(0); 1128 ierr = ISGetLocalSize(dimIS, &numPoints);CHKERRQ(ierr); 1129 ierr = ISGetIndices(dimIS, &points);CHKERRQ(ierr); 1130 for (p = 0; p < numPoints; ++p) { 1131 PetscInt *star = NULL; 1132 PetscInt starSize, s; 1133 PetscInt again = 1; /* 0: Finished 1: Keep iterating after a change 2: No change */ 1134 1135 /* First mark cells connected to the fault */ 1136 ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1137 while (again) { 1138 if (again > 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Could not classify all cells connected to the fault"); 1139 again = 0; 1140 for (s = 0; s < starSize*2; s += 2) { 1141 const PetscInt point = star[s]; 1142 const PetscInt *cone; 1143 PetscInt coneSize, c; 1144 1145 if ((point < cStart) || (point >= cEnd)) continue; 1146 ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr); 1147 if (val != -1) continue; 1148 again = 2; 1149 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 1150 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 1151 for (c = 0; c < coneSize; ++c) { 1152 ierr = DMLabelGetValue(label, cone[c], &val);CHKERRQ(ierr); 1153 if (val != -1) { 1154 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); 1155 if (val > 0) { 1156 ierr = DMLabelSetValue(label, point, shift+dim);CHKERRQ(ierr); 1157 } else { 1158 ierr = DMLabelSetValue(label, point, -(shift+dim));CHKERRQ(ierr); 1159 } 1160 again = 1; 1161 break; 1162 } 1163 } 1164 } 1165 } 1166 /* Classify the rest by cell membership */ 1167 for (s = 0; s < starSize*2; s += 2) { 1168 const PetscInt point = star[s]; 1169 1170 ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr); 1171 if (val == -1) { 1172 PetscInt *sstar = NULL; 1173 PetscInt sstarSize, ss; 1174 PetscBool marked = PETSC_FALSE; 1175 1176 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &sstarSize, &sstar);CHKERRQ(ierr); 1177 for (ss = 0; ss < sstarSize*2; ss += 2) { 1178 const PetscInt spoint = sstar[ss]; 1179 1180 if ((spoint < cStart) || (spoint >= cEnd)) continue; 1181 ierr = DMLabelGetValue(label, spoint, &val);CHKERRQ(ierr); 1182 if (val == -1) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Cell %d in star of %d does not have a valid label", spoint, point); 1183 ierr = DMLabelGetValue(depthLabel, point, &dep);CHKERRQ(ierr); 1184 if (val > 0) { 1185 ierr = DMLabelSetValue(label, point, shift+dep);CHKERRQ(ierr); 1186 } else { 1187 ierr = DMLabelSetValue(label, point, -(shift+dep));CHKERRQ(ierr); 1188 } 1189 marked = PETSC_TRUE; 1190 break; 1191 } 1192 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &sstarSize, &sstar);CHKERRQ(ierr); 1193 if (!marked) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Point %d could not be classified", point); 1194 } 1195 } 1196 ierr = DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1197 } 1198 ierr = ISRestoreIndices(dimIS, &points);CHKERRQ(ierr); 1199 ierr = ISDestroy(&dimIS);CHKERRQ(ierr); 1200 PetscFunctionReturn(0); 1201 } 1202 1203 #undef __FUNCT__ 1204 #define __FUNCT__ "DMPlexCreateHybridMesh" 1205 /*@C 1206 DMPlexCreateHybridMesh - Create a mesh with hybrid cells along an internal interface 1207 1208 Collective on dm 1209 1210 Input Parameters: 1211 + dm - The original DM 1212 - labelName - The label specifying the interface vertices 1213 1214 Output Parameters: 1215 + hybridLabel - The label fully marking the interface 1216 - dmHybrid - The new DM 1217 1218 Level: developer 1219 1220 .seealso: DMPlexConstructCohesiveCells(), DMPlexLabelCohesiveComplete(), DMCreate() 1221 @*/ 1222 PetscErrorCode DMPlexCreateHybridMesh(DM dm, DMLabel label, DMLabel *hybridLabel, DM *dmHybrid) 1223 { 1224 DM idm; 1225 DMLabel subpointMap, hlabel; 1226 PetscInt dim; 1227 PetscErrorCode ierr; 1228 1229 PetscFunctionBegin; 1230 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1231 if (hybridLabel) PetscValidPointer(hybridLabel, 3); 1232 PetscValidPointer(dmHybrid, 4); 1233 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1234 ierr = DMPlexCreateSubmesh(dm, label, 1, &idm);CHKERRQ(ierr); 1235 ierr = DMPlexOrient(idm);CHKERRQ(ierr); 1236 ierr = DMPlexGetSubpointMap(idm, &subpointMap);CHKERRQ(ierr); 1237 ierr = DMLabelDuplicate(subpointMap, &hlabel);CHKERRQ(ierr); 1238 ierr = DMLabelClearStratum(hlabel, dim);CHKERRQ(ierr); 1239 ierr = DMPlexLabelCohesiveComplete(dm, hlabel, PETSC_FALSE, idm);CHKERRQ(ierr); 1240 ierr = DMDestroy(&idm);CHKERRQ(ierr); 1241 ierr = DMPlexConstructCohesiveCells(dm, hlabel, dmHybrid);CHKERRQ(ierr); 1242 if (hybridLabel) *hybridLabel = hlabel; 1243 else {ierr = DMLabelDestroy(&hlabel);CHKERRQ(ierr);} 1244 PetscFunctionReturn(0); 1245 } 1246 1247 #undef __FUNCT__ 1248 #define __FUNCT__ "DMPlexMarkSubmesh_Uninterpolated" 1249 /* Here we need the explicit assumption that: 1250 1251 For any marked cell, the marked vertices constitute a single face 1252 */ 1253 static PetscErrorCode DMPlexMarkSubmesh_Uninterpolated(DM dm, DMLabel vertexLabel, PetscInt value, DMLabel subpointMap, PetscInt *numFaces, PetscInt *nFV, DM subdm) 1254 { 1255 IS subvertexIS = NULL; 1256 const PetscInt *subvertices; 1257 PetscInt *pStart, *pEnd, *pMax, pSize; 1258 PetscInt depth, dim, d, numSubVerticesInitial = 0, v; 1259 PetscErrorCode ierr; 1260 1261 PetscFunctionBegin; 1262 *numFaces = 0; 1263 *nFV = 0; 1264 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1265 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1266 pSize = PetscMax(depth, dim) + 1; 1267 ierr = PetscMalloc3(pSize,PetscInt,&pStart,pSize,PetscInt,&pEnd,pSize,PetscInt,&pMax);CHKERRQ(ierr); 1268 ierr = DMPlexGetHybridBounds(dm, depth >= 0 ? &pMax[depth] : NULL, depth>1 ? &pMax[depth-1] : NULL, depth>2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr); 1269 for (d = 0; d <= depth; ++d) { 1270 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1271 if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]); 1272 } 1273 /* Loop over initial vertices and mark all faces in the collective star() */ 1274 if (vertexLabel) {ierr = DMLabelGetStratumIS(vertexLabel, value, &subvertexIS);CHKERRQ(ierr);} 1275 if (subvertexIS) { 1276 ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr); 1277 ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1278 } 1279 for (v = 0; v < numSubVerticesInitial; ++v) { 1280 const PetscInt vertex = subvertices[v]; 1281 PetscInt *star = NULL; 1282 PetscInt starSize, s, numCells = 0, c; 1283 1284 ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1285 for (s = 0; s < starSize*2; s += 2) { 1286 const PetscInt point = star[s]; 1287 if ((point >= pStart[depth]) && (point < pEnd[depth])) star[numCells++] = point; 1288 } 1289 for (c = 0; c < numCells; ++c) { 1290 const PetscInt cell = star[c]; 1291 PetscInt *closure = NULL; 1292 PetscInt closureSize, cl; 1293 PetscInt cellLoc, numCorners = 0, faceSize = 0; 1294 1295 ierr = DMLabelGetValue(subpointMap, cell, &cellLoc);CHKERRQ(ierr); 1296 if (cellLoc == 2) continue; 1297 if (cellLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Cell %d has dimension %d in the surface label", cell, cellLoc); 1298 ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1299 for (cl = 0; cl < closureSize*2; cl += 2) { 1300 const PetscInt point = closure[cl]; 1301 PetscInt vertexLoc; 1302 1303 if ((point >= pStart[0]) && (point < pEnd[0])) { 1304 ++numCorners; 1305 ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr); 1306 if (vertexLoc == value) closure[faceSize++] = point; 1307 } 1308 } 1309 if (!(*nFV)) {ierr = DMPlexGetNumFaceVertices(dm, dim, numCorners, nFV);CHKERRQ(ierr);} 1310 if (faceSize > *nFV) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize); 1311 if (faceSize == *nFV) { 1312 const PetscInt *cells = NULL; 1313 PetscInt numCells, nc; 1314 1315 ++(*numFaces); 1316 for (cl = 0; cl < faceSize; ++cl) { 1317 ierr = DMLabelSetValue(subpointMap, closure[cl], 0);CHKERRQ(ierr); 1318 } 1319 ierr = DMPlexGetJoin(dm, faceSize, closure, &numCells, &cells);CHKERRQ(ierr); 1320 for (nc = 0; nc < numCells; ++nc) { 1321 ierr = DMLabelSetValue(subpointMap, cells[nc], 2);CHKERRQ(ierr); 1322 } 1323 ierr = DMPlexRestoreJoin(dm, faceSize, closure, &numCells, &cells);CHKERRQ(ierr); 1324 } 1325 ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1326 } 1327 ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1328 } 1329 if (subvertexIS) { 1330 ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1331 } 1332 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 1333 ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr); 1334 PetscFunctionReturn(0); 1335 } 1336 1337 #undef __FUNCT__ 1338 #define __FUNCT__ "DMPlexMarkSubmesh_Interpolated" 1339 static PetscErrorCode DMPlexMarkSubmesh_Interpolated(DM dm, DMLabel vertexLabel, PetscInt value, DMLabel subpointMap, DM subdm) 1340 { 1341 IS subvertexIS; 1342 const PetscInt *subvertices; 1343 PetscInt *pStart, *pEnd, *pMax; 1344 PetscInt dim, d, numSubVerticesInitial = 0, v; 1345 PetscErrorCode ierr; 1346 1347 PetscFunctionBegin; 1348 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1349 ierr = PetscMalloc3(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd,dim+1,PetscInt,&pMax);CHKERRQ(ierr); 1350 ierr = DMPlexGetHybridBounds(dm, &pMax[dim], dim>1 ? &pMax[dim-1] : NULL, dim > 2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr); 1351 for (d = 0; d <= dim; ++d) { 1352 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1353 if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]); 1354 } 1355 /* Loop over initial vertices and mark all faces in the collective star() */ 1356 ierr = DMLabelGetStratumIS(vertexLabel, value, &subvertexIS);CHKERRQ(ierr); 1357 if (subvertexIS) { 1358 ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr); 1359 ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1360 } 1361 for (v = 0; v < numSubVerticesInitial; ++v) { 1362 const PetscInt vertex = subvertices[v]; 1363 PetscInt *star = NULL; 1364 PetscInt starSize, s, numFaces = 0, f; 1365 1366 ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1367 for (s = 0; s < starSize*2; s += 2) { 1368 const PetscInt point = star[s]; 1369 if ((point >= pStart[dim-1]) && (point < pEnd[dim-1])) star[numFaces++] = point; 1370 } 1371 for (f = 0; f < numFaces; ++f) { 1372 const PetscInt face = star[f]; 1373 PetscInt *closure = NULL; 1374 PetscInt closureSize, c; 1375 PetscInt faceLoc; 1376 1377 ierr = DMLabelGetValue(subpointMap, face, &faceLoc);CHKERRQ(ierr); 1378 if (faceLoc == dim-1) continue; 1379 if (faceLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Face %d has dimension %d in the surface label", face, faceLoc); 1380 ierr = DMPlexGetTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1381 for (c = 0; c < closureSize*2; c += 2) { 1382 const PetscInt point = closure[c]; 1383 PetscInt vertexLoc; 1384 1385 if ((point >= pStart[0]) && (point < pEnd[0])) { 1386 ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr); 1387 if (vertexLoc != value) break; 1388 } 1389 } 1390 if (c == closureSize*2) { 1391 const PetscInt *support; 1392 PetscInt supportSize, s; 1393 1394 for (c = 0; c < closureSize*2; c += 2) { 1395 const PetscInt point = closure[c]; 1396 1397 for (d = 0; d < dim; ++d) { 1398 if ((point >= pStart[d]) && (point < pEnd[d])) { 1399 ierr = DMLabelSetValue(subpointMap, point, d);CHKERRQ(ierr); 1400 break; 1401 } 1402 } 1403 } 1404 ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 1405 ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); 1406 for (s = 0; s < supportSize; ++s) { 1407 ierr = DMLabelSetValue(subpointMap, support[s], dim);CHKERRQ(ierr); 1408 } 1409 } 1410 ierr = DMPlexRestoreTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1411 } 1412 ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1413 } 1414 if (subvertexIS) { 1415 ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1416 } 1417 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 1418 ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr); 1419 PetscFunctionReturn(0); 1420 } 1421 1422 #undef __FUNCT__ 1423 #define __FUNCT__ "DMPlexMarkCohesiveSubmesh_Uninterpolated" 1424 static PetscErrorCode DMPlexMarkCohesiveSubmesh_Uninterpolated(DM dm, PetscBool hasLagrange, const char labelname[], PetscInt value, DMLabel subpointMap, PetscInt *numFaces, PetscInt *nFV, PetscInt *subCells[], DM subdm) 1425 { 1426 DMLabel label = NULL; 1427 const PetscInt *cone; 1428 PetscInt dim, cMax, cEnd, c, subc = 0, p, coneSize; 1429 PetscErrorCode ierr; 1430 1431 PetscFunctionBegin; 1432 *numFaces = 0; 1433 *nFV = 0; 1434 if (labelname) {ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);} 1435 *subCells = NULL; 1436 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1437 ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr); 1438 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1439 if (cMax < 0) PetscFunctionReturn(0); 1440 if (label) { 1441 for (c = cMax; c < cEnd; ++c) { 1442 PetscInt val; 1443 1444 ierr = DMLabelGetValue(label, c, &val);CHKERRQ(ierr); 1445 if (val == value) { 1446 ++(*numFaces); 1447 ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 1448 } 1449 } 1450 } else { 1451 *numFaces = cEnd - cMax; 1452 ierr = DMPlexGetConeSize(dm, cMax, &coneSize);CHKERRQ(ierr); 1453 } 1454 *nFV = hasLagrange ? coneSize/3 : coneSize/2; 1455 ierr = PetscMalloc(*numFaces *2 * sizeof(PetscInt), subCells);CHKERRQ(ierr); 1456 for (c = cMax; c < cEnd; ++c) { 1457 const PetscInt *cells; 1458 PetscInt numCells; 1459 1460 if (label) { 1461 PetscInt val; 1462 1463 ierr = DMLabelGetValue(label, c, &val);CHKERRQ(ierr); 1464 if (val != value) continue; 1465 } 1466 ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 1467 for (p = 0; p < *nFV; ++p) { 1468 ierr = DMLabelSetValue(subpointMap, cone[p], 0);CHKERRQ(ierr); 1469 } 1470 /* Negative face */ 1471 ierr = DMPlexGetJoin(dm, *nFV, cone, &numCells, &cells);CHKERRQ(ierr); 1472 /* Not true in parallel 1473 if (numCells != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells"); */ 1474 for (p = 0; p < numCells; ++p) { 1475 ierr = DMLabelSetValue(subpointMap, cells[p], 2);CHKERRQ(ierr); 1476 (*subCells)[subc++] = cells[p]; 1477 } 1478 ierr = DMPlexRestoreJoin(dm, *nFV, cone, &numCells, &cells);CHKERRQ(ierr); 1479 /* Positive face is not included */ 1480 } 1481 PetscFunctionReturn(0); 1482 } 1483 1484 #undef __FUNCT__ 1485 #define __FUNCT__ "DMPlexMarkCohesiveSubmesh_Interpolated" 1486 static PetscErrorCode DMPlexMarkCohesiveSubmesh_Interpolated(DM dm, PetscBool hasLagrange, const char labelname[], PetscInt value, DMLabel subpointMap, DM subdm) 1487 { 1488 DMLabel label = NULL; 1489 PetscInt *pStart, *pEnd; 1490 PetscInt dim, cMax, cEnd, c, d; 1491 PetscErrorCode ierr; 1492 1493 PetscFunctionBegin; 1494 if (labelname) {ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);} 1495 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1496 ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr); 1497 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1498 if (cMax < 0) PetscFunctionReturn(0); 1499 ierr = PetscMalloc2(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd);CHKERRQ(ierr); 1500 for (d = 0; d <= dim; ++d) { 1501 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1502 } 1503 for (c = cMax; c < cEnd; ++c) { 1504 const PetscInt *cone; 1505 PetscInt *closure = NULL; 1506 PetscInt coneSize, closureSize, cl; 1507 1508 if (label) { 1509 PetscInt val; 1510 1511 ierr = DMLabelGetValue(label, c, &val);CHKERRQ(ierr); 1512 if (val != value) continue; 1513 } 1514 ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 1515 if (hasLagrange) { 1516 if (coneSize != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells"); 1517 } else { 1518 if (coneSize != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells"); 1519 } 1520 /* Negative face */ 1521 ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 1522 ierr = DMPlexGetTransitiveClosure(dm, cone[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1523 for (cl = 0; cl < closureSize*2; cl += 2) { 1524 const PetscInt point = closure[cl]; 1525 1526 for (d = 0; d <= dim; ++d) { 1527 if ((point >= pStart[d]) && (point < pEnd[d])) { 1528 ierr = DMLabelSetValue(subpointMap, point, d);CHKERRQ(ierr); 1529 break; 1530 } 1531 } 1532 } 1533 ierr = DMPlexRestoreTransitiveClosure(dm, cone[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1534 /* Cells -- positive face is not included */ 1535 for (cl = 0; cl < 1; ++cl) { 1536 const PetscInt *support; 1537 PetscInt supportSize, s; 1538 1539 ierr = DMPlexGetSupportSize(dm, cone[cl], &supportSize);CHKERRQ(ierr); 1540 if (supportSize != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive faces should separate two cells"); 1541 ierr = DMPlexGetSupport(dm, cone[cl], &support);CHKERRQ(ierr); 1542 for (s = 0; s < supportSize; ++s) { 1543 ierr = DMLabelSetValue(subpointMap, support[s], dim);CHKERRQ(ierr); 1544 } 1545 } 1546 } 1547 ierr = PetscFree2(pStart, pEnd);CHKERRQ(ierr); 1548 PetscFunctionReturn(0); 1549 } 1550 1551 #undef __FUNCT__ 1552 #define __FUNCT__ "DMPlexGetFaceOrientation" 1553 PetscErrorCode DMPlexGetFaceOrientation(DM dm, PetscInt cell, PetscInt numCorners, PetscInt indices[], PetscInt oppositeVertex, PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented) 1554 { 1555 MPI_Comm comm; 1556 PetscBool posOrient = PETSC_FALSE; 1557 const PetscInt debug = 0; 1558 PetscInt cellDim, faceSize, f; 1559 PetscErrorCode ierr; 1560 1561 PetscFunctionBegin; 1562 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1563 ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr); 1564 if (debug) {PetscPrintf(comm, "cellDim: %d numCorners: %d\n", cellDim, numCorners);CHKERRQ(ierr);} 1565 1566 if (cellDim == 1 && numCorners == 2) { 1567 /* Triangle */ 1568 faceSize = numCorners-1; 1569 posOrient = !(oppositeVertex%2) ? PETSC_TRUE : PETSC_FALSE; 1570 } else if (cellDim == 2 && numCorners == 3) { 1571 /* Triangle */ 1572 faceSize = numCorners-1; 1573 posOrient = !(oppositeVertex%2) ? PETSC_TRUE : PETSC_FALSE; 1574 } else if (cellDim == 3 && numCorners == 4) { 1575 /* Tetrahedron */ 1576 faceSize = numCorners-1; 1577 posOrient = (oppositeVertex%2) ? PETSC_TRUE : PETSC_FALSE; 1578 } else if (cellDim == 1 && numCorners == 3) { 1579 /* Quadratic line */ 1580 faceSize = 1; 1581 posOrient = PETSC_TRUE; 1582 } else if (cellDim == 2 && numCorners == 4) { 1583 /* Quads */ 1584 faceSize = 2; 1585 if ((indices[1] > indices[0]) && (indices[1] - indices[0] == 1)) { 1586 posOrient = PETSC_TRUE; 1587 } else if ((indices[0] == 3) && (indices[1] == 0)) { 1588 posOrient = PETSC_TRUE; 1589 } else { 1590 if (((indices[0] > indices[1]) && (indices[0] - indices[1] == 1)) || ((indices[0] == 0) && (indices[1] == 3))) { 1591 posOrient = PETSC_FALSE; 1592 } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossedge"); 1593 } 1594 } else if (cellDim == 2 && numCorners == 6) { 1595 /* Quadratic triangle (I hate this) */ 1596 /* Edges are determined by the first 2 vertices (corners of edges) */ 1597 const PetscInt faceSizeTri = 3; 1598 PetscInt sortedIndices[3], i, iFace; 1599 PetscBool found = PETSC_FALSE; 1600 PetscInt faceVerticesTriSorted[9] = { 1601 0, 3, 4, /* bottom */ 1602 1, 4, 5, /* right */ 1603 2, 3, 5, /* left */ 1604 }; 1605 PetscInt faceVerticesTri[9] = { 1606 0, 3, 4, /* bottom */ 1607 1, 4, 5, /* right */ 1608 2, 5, 3, /* left */ 1609 }; 1610 1611 faceSize = faceSizeTri; 1612 for (i = 0; i < faceSizeTri; ++i) sortedIndices[i] = indices[i]; 1613 ierr = PetscSortInt(faceSizeTri, sortedIndices);CHKERRQ(ierr); 1614 for (iFace = 0; iFace < 3; ++iFace) { 1615 const PetscInt ii = iFace*faceSizeTri; 1616 PetscInt fVertex, cVertex; 1617 1618 if ((sortedIndices[0] == faceVerticesTriSorted[ii+0]) && 1619 (sortedIndices[1] == faceVerticesTriSorted[ii+1])) { 1620 for (fVertex = 0; fVertex < faceSizeTri; ++fVertex) { 1621 for (cVertex = 0; cVertex < faceSizeTri; ++cVertex) { 1622 if (indices[cVertex] == faceVerticesTri[ii+fVertex]) { 1623 faceVertices[fVertex] = origVertices[cVertex]; 1624 break; 1625 } 1626 } 1627 } 1628 found = PETSC_TRUE; 1629 break; 1630 } 1631 } 1632 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tri crossface"); 1633 if (posOriented) *posOriented = PETSC_TRUE; 1634 PetscFunctionReturn(0); 1635 } else if (cellDim == 2 && numCorners == 9) { 1636 /* Quadratic quad (I hate this) */ 1637 /* Edges are determined by the first 2 vertices (corners of edges) */ 1638 const PetscInt faceSizeQuad = 3; 1639 PetscInt sortedIndices[3], i, iFace; 1640 PetscBool found = PETSC_FALSE; 1641 PetscInt faceVerticesQuadSorted[12] = { 1642 0, 1, 4, /* bottom */ 1643 1, 2, 5, /* right */ 1644 2, 3, 6, /* top */ 1645 0, 3, 7, /* left */ 1646 }; 1647 PetscInt faceVerticesQuad[12] = { 1648 0, 1, 4, /* bottom */ 1649 1, 2, 5, /* right */ 1650 2, 3, 6, /* top */ 1651 3, 0, 7, /* left */ 1652 }; 1653 1654 faceSize = faceSizeQuad; 1655 for (i = 0; i < faceSizeQuad; ++i) sortedIndices[i] = indices[i]; 1656 ierr = PetscSortInt(faceSizeQuad, sortedIndices);CHKERRQ(ierr); 1657 for (iFace = 0; iFace < 4; ++iFace) { 1658 const PetscInt ii = iFace*faceSizeQuad; 1659 PetscInt fVertex, cVertex; 1660 1661 if ((sortedIndices[0] == faceVerticesQuadSorted[ii+0]) && 1662 (sortedIndices[1] == faceVerticesQuadSorted[ii+1])) { 1663 for (fVertex = 0; fVertex < faceSizeQuad; ++fVertex) { 1664 for (cVertex = 0; cVertex < faceSizeQuad; ++cVertex) { 1665 if (indices[cVertex] == faceVerticesQuad[ii+fVertex]) { 1666 faceVertices[fVertex] = origVertices[cVertex]; 1667 break; 1668 } 1669 } 1670 } 1671 found = PETSC_TRUE; 1672 break; 1673 } 1674 } 1675 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossface"); 1676 if (posOriented) *posOriented = PETSC_TRUE; 1677 PetscFunctionReturn(0); 1678 } else if (cellDim == 3 && numCorners == 8) { 1679 /* Hexes 1680 A hex is two oriented quads with the normal of the first 1681 pointing up at the second. 1682 1683 7---6 1684 /| /| 1685 4---5 | 1686 | 1-|-2 1687 |/ |/ 1688 0---3 1689 1690 Faces are determined by the first 4 vertices (corners of faces) */ 1691 const PetscInt faceSizeHex = 4; 1692 PetscInt sortedIndices[4], i, iFace; 1693 PetscBool found = PETSC_FALSE; 1694 PetscInt faceVerticesHexSorted[24] = { 1695 0, 1, 2, 3, /* bottom */ 1696 4, 5, 6, 7, /* top */ 1697 0, 3, 4, 5, /* front */ 1698 2, 3, 5, 6, /* right */ 1699 1, 2, 6, 7, /* back */ 1700 0, 1, 4, 7, /* left */ 1701 }; 1702 PetscInt faceVerticesHex[24] = { 1703 1, 2, 3, 0, /* bottom */ 1704 4, 5, 6, 7, /* top */ 1705 0, 3, 5, 4, /* front */ 1706 3, 2, 6, 5, /* right */ 1707 2, 1, 7, 6, /* back */ 1708 1, 0, 4, 7, /* left */ 1709 }; 1710 1711 faceSize = faceSizeHex; 1712 for (i = 0; i < faceSizeHex; ++i) sortedIndices[i] = indices[i]; 1713 ierr = PetscSortInt(faceSizeHex, sortedIndices);CHKERRQ(ierr); 1714 for (iFace = 0; iFace < 6; ++iFace) { 1715 const PetscInt ii = iFace*faceSizeHex; 1716 PetscInt fVertex, cVertex; 1717 1718 if ((sortedIndices[0] == faceVerticesHexSorted[ii+0]) && 1719 (sortedIndices[1] == faceVerticesHexSorted[ii+1]) && 1720 (sortedIndices[2] == faceVerticesHexSorted[ii+2]) && 1721 (sortedIndices[3] == faceVerticesHexSorted[ii+3])) { 1722 for (fVertex = 0; fVertex < faceSizeHex; ++fVertex) { 1723 for (cVertex = 0; cVertex < faceSizeHex; ++cVertex) { 1724 if (indices[cVertex] == faceVerticesHex[ii+fVertex]) { 1725 faceVertices[fVertex] = origVertices[cVertex]; 1726 break; 1727 } 1728 } 1729 } 1730 found = PETSC_TRUE; 1731 break; 1732 } 1733 } 1734 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface"); 1735 if (posOriented) *posOriented = PETSC_TRUE; 1736 PetscFunctionReturn(0); 1737 } else if (cellDim == 3 && numCorners == 10) { 1738 /* Quadratic tet */ 1739 /* Faces are determined by the first 3 vertices (corners of faces) */ 1740 const PetscInt faceSizeTet = 6; 1741 PetscInt sortedIndices[6], i, iFace; 1742 PetscBool found = PETSC_FALSE; 1743 PetscInt faceVerticesTetSorted[24] = { 1744 0, 1, 2, 6, 7, 8, /* bottom */ 1745 0, 3, 4, 6, 7, 9, /* front */ 1746 1, 4, 5, 7, 8, 9, /* right */ 1747 2, 3, 5, 6, 8, 9, /* left */ 1748 }; 1749 PetscInt faceVerticesTet[24] = { 1750 0, 1, 2, 6, 7, 8, /* bottom */ 1751 0, 4, 3, 6, 7, 9, /* front */ 1752 1, 5, 4, 7, 8, 9, /* right */ 1753 2, 3, 5, 8, 6, 9, /* left */ 1754 }; 1755 1756 faceSize = faceSizeTet; 1757 for (i = 0; i < faceSizeTet; ++i) sortedIndices[i] = indices[i]; 1758 ierr = PetscSortInt(faceSizeTet, sortedIndices);CHKERRQ(ierr); 1759 for (iFace=0; iFace < 4; ++iFace) { 1760 const PetscInt ii = iFace*faceSizeTet; 1761 PetscInt fVertex, cVertex; 1762 1763 if ((sortedIndices[0] == faceVerticesTetSorted[ii+0]) && 1764 (sortedIndices[1] == faceVerticesTetSorted[ii+1]) && 1765 (sortedIndices[2] == faceVerticesTetSorted[ii+2]) && 1766 (sortedIndices[3] == faceVerticesTetSorted[ii+3])) { 1767 for (fVertex = 0; fVertex < faceSizeTet; ++fVertex) { 1768 for (cVertex = 0; cVertex < faceSizeTet; ++cVertex) { 1769 if (indices[cVertex] == faceVerticesTet[ii+fVertex]) { 1770 faceVertices[fVertex] = origVertices[cVertex]; 1771 break; 1772 } 1773 } 1774 } 1775 found = PETSC_TRUE; 1776 break; 1777 } 1778 } 1779 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tet crossface"); 1780 if (posOriented) *posOriented = PETSC_TRUE; 1781 PetscFunctionReturn(0); 1782 } else if (cellDim == 3 && numCorners == 27) { 1783 /* Quadratic hexes (I hate this) 1784 A hex is two oriented quads with the normal of the first 1785 pointing up at the second. 1786 1787 7---6 1788 /| /| 1789 4---5 | 1790 | 3-|-2 1791 |/ |/ 1792 0---1 1793 1794 Faces are determined by the first 4 vertices (corners of faces) */ 1795 const PetscInt faceSizeQuadHex = 9; 1796 PetscInt sortedIndices[9], i, iFace; 1797 PetscBool found = PETSC_FALSE; 1798 PetscInt faceVerticesQuadHexSorted[54] = { 1799 0, 1, 2, 3, 8, 9, 10, 11, 24, /* bottom */ 1800 4, 5, 6, 7, 12, 13, 14, 15, 25, /* top */ 1801 0, 1, 4, 5, 8, 12, 16, 17, 22, /* front */ 1802 1, 2, 5, 6, 9, 13, 17, 18, 21, /* right */ 1803 2, 3, 6, 7, 10, 14, 18, 19, 23, /* back */ 1804 0, 3, 4, 7, 11, 15, 16, 19, 20, /* left */ 1805 }; 1806 PetscInt faceVerticesQuadHex[54] = { 1807 3, 2, 1, 0, 10, 9, 8, 11, 24, /* bottom */ 1808 4, 5, 6, 7, 12, 13, 14, 15, 25, /* top */ 1809 0, 1, 5, 4, 8, 17, 12, 16, 22, /* front */ 1810 1, 2, 6, 5, 9, 18, 13, 17, 21, /* right */ 1811 2, 3, 7, 6, 10, 19, 14, 18, 23, /* back */ 1812 3, 0, 4, 7, 11, 16, 15, 19, 20 /* left */ 1813 }; 1814 1815 faceSize = faceSizeQuadHex; 1816 for (i = 0; i < faceSizeQuadHex; ++i) sortedIndices[i] = indices[i]; 1817 ierr = PetscSortInt(faceSizeQuadHex, sortedIndices);CHKERRQ(ierr); 1818 for (iFace = 0; iFace < 6; ++iFace) { 1819 const PetscInt ii = iFace*faceSizeQuadHex; 1820 PetscInt fVertex, cVertex; 1821 1822 if ((sortedIndices[0] == faceVerticesQuadHexSorted[ii+0]) && 1823 (sortedIndices[1] == faceVerticesQuadHexSorted[ii+1]) && 1824 (sortedIndices[2] == faceVerticesQuadHexSorted[ii+2]) && 1825 (sortedIndices[3] == faceVerticesQuadHexSorted[ii+3])) { 1826 for (fVertex = 0; fVertex < faceSizeQuadHex; ++fVertex) { 1827 for (cVertex = 0; cVertex < faceSizeQuadHex; ++cVertex) { 1828 if (indices[cVertex] == faceVerticesQuadHex[ii+fVertex]) { 1829 faceVertices[fVertex] = origVertices[cVertex]; 1830 break; 1831 } 1832 } 1833 } 1834 found = PETSC_TRUE; 1835 break; 1836 } 1837 } 1838 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface"); 1839 if (posOriented) *posOriented = PETSC_TRUE; 1840 PetscFunctionReturn(0); 1841 } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Unknown cell type for faceOrientation()."); 1842 if (!posOrient) { 1843 if (debug) {ierr = PetscPrintf(comm, " Reversing initial face orientation\n");CHKERRQ(ierr);} 1844 for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[faceSize-1 - f]; 1845 } else { 1846 if (debug) {ierr = PetscPrintf(comm, " Keeping initial face orientation\n");CHKERRQ(ierr);} 1847 for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[f]; 1848 } 1849 if (posOriented) *posOriented = posOrient; 1850 PetscFunctionReturn(0); 1851 } 1852 1853 #undef __FUNCT__ 1854 #define __FUNCT__ "DMPlexGetOrientedFace" 1855 /* 1856 Given a cell and a face, as a set of vertices, 1857 return the oriented face, as a set of vertices, in faceVertices 1858 The orientation is such that the face normal points out of the cell 1859 */ 1860 PetscErrorCode DMPlexGetOrientedFace(DM dm, PetscInt cell, PetscInt faceSize, const PetscInt face[], PetscInt numCorners, PetscInt indices[], PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented) 1861 { 1862 const PetscInt *cone = NULL; 1863 PetscInt coneSize, v, f, v2; 1864 PetscInt oppositeVertex = -1; 1865 PetscErrorCode ierr; 1866 1867 PetscFunctionBegin; 1868 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 1869 ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 1870 for (v = 0, v2 = 0; v < coneSize; ++v) { 1871 PetscBool found = PETSC_FALSE; 1872 1873 for (f = 0; f < faceSize; ++f) { 1874 if (face[f] == cone[v]) { 1875 found = PETSC_TRUE; break; 1876 } 1877 } 1878 if (found) { 1879 indices[v2] = v; 1880 origVertices[v2] = cone[v]; 1881 ++v2; 1882 } else { 1883 oppositeVertex = v; 1884 } 1885 } 1886 ierr = DMPlexGetFaceOrientation(dm, cell, numCorners, indices, oppositeVertex, origVertices, faceVertices, posOriented);CHKERRQ(ierr); 1887 PetscFunctionReturn(0); 1888 } 1889 1890 #undef __FUNCT__ 1891 #define __FUNCT__ "DMPlexInsertFace_Internal" 1892 /* 1893 DMPlexInsertFace_Internal - Puts a face into the mesh 1894 1895 Not collective 1896 1897 Input Parameters: 1898 + dm - The DMPlex 1899 . numFaceVertex - The number of vertices in the face 1900 . faceVertices - The vertices in the face for dm 1901 . subfaceVertices - The vertices in the face for subdm 1902 . numCorners - The number of vertices in the cell 1903 . cell - A cell in dm containing the face 1904 . subcell - A cell in subdm containing the face 1905 . firstFace - First face in the mesh 1906 - newFacePoint - Next face in the mesh 1907 1908 Output Parameters: 1909 . newFacePoint - Contains next face point number on input, updated on output 1910 1911 Level: developer 1912 */ 1913 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) 1914 { 1915 MPI_Comm comm; 1916 DM_Plex *submesh = (DM_Plex*) subdm->data; 1917 const PetscInt *faces; 1918 PetscInt numFaces, coneSize; 1919 PetscErrorCode ierr; 1920 1921 PetscFunctionBegin; 1922 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1923 ierr = DMPlexGetConeSize(subdm, subcell, &coneSize);CHKERRQ(ierr); 1924 if (coneSize != 1) SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cone size of cell %d is %d != 1", cell, coneSize); 1925 #if 0 1926 /* Cannot use this because support() has not been constructed yet */ 1927 ierr = DMPlexGetJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr); 1928 #else 1929 { 1930 PetscInt f; 1931 1932 numFaces = 0; 1933 ierr = DMGetWorkArray(subdm, 1, PETSC_INT, (void **) &faces);CHKERRQ(ierr); 1934 for (f = firstFace; f < *newFacePoint; ++f) { 1935 PetscInt dof, off, d; 1936 1937 ierr = PetscSectionGetDof(submesh->coneSection, f, &dof);CHKERRQ(ierr); 1938 ierr = PetscSectionGetOffset(submesh->coneSection, f, &off);CHKERRQ(ierr); 1939 /* Yes, I know this is quadratic, but I expect the sizes to be <5 */ 1940 for (d = 0; d < dof; ++d) { 1941 const PetscInt p = submesh->cones[off+d]; 1942 PetscInt v; 1943 1944 for (v = 0; v < numFaceVertices; ++v) { 1945 if (subfaceVertices[v] == p) break; 1946 } 1947 if (v == numFaceVertices) break; 1948 } 1949 if (d == dof) { 1950 numFaces = 1; 1951 ((PetscInt*) faces)[0] = f; 1952 } 1953 } 1954 } 1955 #endif 1956 if (numFaces > 1) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Vertex set had %d faces, not one", numFaces); 1957 else if (numFaces == 1) { 1958 /* Add the other cell neighbor for this face */ 1959 ierr = DMPlexSetCone(subdm, subcell, faces);CHKERRQ(ierr); 1960 } else { 1961 PetscInt *indices, *origVertices, *orientedVertices, *orientedSubVertices, v, ov; 1962 PetscBool posOriented; 1963 1964 ierr = DMGetWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr); 1965 origVertices = &orientedVertices[numFaceVertices]; 1966 indices = &orientedVertices[numFaceVertices*2]; 1967 orientedSubVertices = &orientedVertices[numFaceVertices*3]; 1968 ierr = DMPlexGetOrientedFace(dm, cell, numFaceVertices, faceVertices, numCorners, indices, origVertices, orientedVertices, &posOriented);CHKERRQ(ierr); 1969 /* TODO: I know that routine should return a permutation, not the indices */ 1970 for (v = 0; v < numFaceVertices; ++v) { 1971 const PetscInt vertex = faceVertices[v], subvertex = subfaceVertices[v]; 1972 for (ov = 0; ov < numFaceVertices; ++ov) { 1973 if (orientedVertices[ov] == vertex) { 1974 orientedSubVertices[ov] = subvertex; 1975 break; 1976 } 1977 } 1978 if (ov == numFaceVertices) SETERRQ1(comm, PETSC_ERR_PLIB, "Could not find face vertex %d in orientated set", vertex); 1979 } 1980 ierr = DMPlexSetCone(subdm, *newFacePoint, orientedSubVertices);CHKERRQ(ierr); 1981 ierr = DMPlexSetCone(subdm, subcell, newFacePoint);CHKERRQ(ierr); 1982 ierr = DMRestoreWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr); 1983 ++(*newFacePoint); 1984 } 1985 #if 0 1986 ierr = DMPlexRestoreJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr); 1987 #else 1988 ierr = DMRestoreWorkArray(subdm, 1, PETSC_INT, (void **) &faces);CHKERRQ(ierr); 1989 #endif 1990 PetscFunctionReturn(0); 1991 } 1992 1993 #undef __FUNCT__ 1994 #define __FUNCT__ "DMPlexCreateSubmesh_Uninterpolated" 1995 static PetscErrorCode DMPlexCreateSubmesh_Uninterpolated(DM dm, DMLabel vertexLabel, PetscInt value, DM subdm) 1996 { 1997 MPI_Comm comm; 1998 DMLabel subpointMap; 1999 IS subvertexIS, subcellIS; 2000 const PetscInt *subVertices, *subCells; 2001 PetscInt numSubVertices, firstSubVertex, numSubCells; 2002 PetscInt *subface, maxConeSize, numSubFaces = 0, firstSubFace, newFacePoint, nFV = 0; 2003 PetscInt vStart, vEnd, c, f; 2004 PetscErrorCode ierr; 2005 2006 PetscFunctionBegin; 2007 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2008 /* Create subpointMap which marks the submesh */ 2009 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 2010 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 2011 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 2012 if (vertexLabel) {ierr = DMPlexMarkSubmesh_Uninterpolated(dm, vertexLabel, value, subpointMap, &numSubFaces, &nFV, subdm);CHKERRQ(ierr);} 2013 /* Setup chart */ 2014 ierr = DMLabelGetStratumSize(subpointMap, 0, &numSubVertices);CHKERRQ(ierr); 2015 ierr = DMLabelGetStratumSize(subpointMap, 2, &numSubCells);CHKERRQ(ierr); 2016 ierr = DMPlexSetChart(subdm, 0, numSubCells+numSubFaces+numSubVertices);CHKERRQ(ierr); 2017 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 2018 /* Set cone sizes */ 2019 firstSubVertex = numSubCells; 2020 firstSubFace = numSubCells+numSubVertices; 2021 newFacePoint = firstSubFace; 2022 ierr = DMLabelGetStratumIS(subpointMap, 0, &subvertexIS);CHKERRQ(ierr); 2023 if (subvertexIS) {ierr = ISGetIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 2024 ierr = DMLabelGetStratumIS(subpointMap, 2, &subcellIS);CHKERRQ(ierr); 2025 if (subcellIS) {ierr = ISGetIndices(subcellIS, &subCells);CHKERRQ(ierr);} 2026 for (c = 0; c < numSubCells; ++c) { 2027 ierr = DMPlexSetConeSize(subdm, c, 1);CHKERRQ(ierr); 2028 } 2029 for (f = firstSubFace; f < firstSubFace+numSubFaces; ++f) { 2030 ierr = DMPlexSetConeSize(subdm, f, nFV);CHKERRQ(ierr); 2031 } 2032 ierr = DMSetUp(subdm);CHKERRQ(ierr); 2033 /* Create face cones */ 2034 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 2035 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 2036 ierr = DMGetWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 2037 for (c = 0; c < numSubCells; ++c) { 2038 const PetscInt cell = subCells[c]; 2039 const PetscInt subcell = c; 2040 PetscInt *closure = NULL; 2041 PetscInt closureSize, cl, numCorners = 0, faceSize = 0; 2042 2043 ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 2044 for (cl = 0; cl < closureSize*2; cl += 2) { 2045 const PetscInt point = closure[cl]; 2046 PetscInt subVertex; 2047 2048 if ((point >= vStart) && (point < vEnd)) { 2049 ++numCorners; 2050 ierr = PetscFindInt(point, numSubVertices, subVertices, &subVertex);CHKERRQ(ierr); 2051 if (subVertex >= 0) { 2052 closure[faceSize] = point; 2053 subface[faceSize] = firstSubVertex+subVertex; 2054 ++faceSize; 2055 } 2056 } 2057 } 2058 if (faceSize > nFV) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize); 2059 if (faceSize == nFV) { 2060 ierr = DMPlexInsertFace_Internal(dm, subdm, faceSize, closure, subface, numCorners, cell, subcell, firstSubFace, &newFacePoint);CHKERRQ(ierr); 2061 } 2062 ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 2063 } 2064 ierr = DMRestoreWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 2065 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 2066 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 2067 /* Build coordinates */ 2068 { 2069 PetscSection coordSection, subCoordSection; 2070 Vec coordinates, subCoordinates; 2071 PetscScalar *coords, *subCoords; 2072 PetscInt numComp, coordSize, v; 2073 2074 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2075 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 2076 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 2077 ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr); 2078 ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr); 2079 ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr); 2080 ierr = PetscSectionSetChart(subCoordSection, firstSubVertex, firstSubVertex+numSubVertices);CHKERRQ(ierr); 2081 for (v = 0; v < numSubVertices; ++v) { 2082 const PetscInt vertex = subVertices[v]; 2083 const PetscInt subvertex = firstSubVertex+v; 2084 PetscInt dof; 2085 2086 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2087 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 2088 ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr); 2089 } 2090 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 2091 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 2092 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 2093 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2094 ierr = VecSetType(subCoordinates,dm->vectype);CHKERRQ(ierr); 2095 if (coordSize) { 2096 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2097 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2098 for (v = 0; v < numSubVertices; ++v) { 2099 const PetscInt vertex = subVertices[v]; 2100 const PetscInt subvertex = firstSubVertex+v; 2101 PetscInt dof, off, sdof, soff, d; 2102 2103 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2104 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 2105 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 2106 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 2107 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 2108 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 2109 } 2110 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2111 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2112 } 2113 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 2114 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 2115 } 2116 /* Cleanup */ 2117 if (subvertexIS) {ierr = ISRestoreIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 2118 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 2119 if (subcellIS) {ierr = ISRestoreIndices(subcellIS, &subCells);CHKERRQ(ierr);} 2120 ierr = ISDestroy(&subcellIS);CHKERRQ(ierr); 2121 PetscFunctionReturn(0); 2122 } 2123 2124 #undef __FUNCT__ 2125 #define __FUNCT__ "DMPlexCreateSubmesh_Interpolated" 2126 static PetscErrorCode DMPlexCreateSubmesh_Interpolated(DM dm, DMLabel vertexLabel, PetscInt value, DM subdm) 2127 { 2128 MPI_Comm comm; 2129 DMLabel subpointMap; 2130 IS *subpointIS; 2131 const PetscInt **subpoints; 2132 PetscInt *numSubPoints, *firstSubPoint, *coneNew, *coneONew; 2133 PetscInt totSubPoints = 0, maxConeSize, dim, p, d, v; 2134 PetscErrorCode ierr; 2135 2136 PetscFunctionBegin; 2137 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2138 /* Create subpointMap which marks the submesh */ 2139 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 2140 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 2141 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 2142 if (vertexLabel) {ierr = DMPlexMarkSubmesh_Interpolated(dm, vertexLabel, value, subpointMap, subdm);CHKERRQ(ierr);} 2143 /* Setup chart */ 2144 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 2145 ierr = PetscMalloc4(dim+1,PetscInt,&numSubPoints,dim+1,PetscInt,&firstSubPoint,dim+1,IS,&subpointIS,dim+1,const PetscInt *,&subpoints);CHKERRQ(ierr); 2146 for (d = 0; d <= dim; ++d) { 2147 ierr = DMLabelGetStratumSize(subpointMap, d, &numSubPoints[d]);CHKERRQ(ierr); 2148 totSubPoints += numSubPoints[d]; 2149 } 2150 ierr = DMPlexSetChart(subdm, 0, totSubPoints);CHKERRQ(ierr); 2151 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 2152 /* Set cone sizes */ 2153 firstSubPoint[dim] = 0; 2154 firstSubPoint[0] = firstSubPoint[dim] + numSubPoints[dim]; 2155 if (dim > 1) {firstSubPoint[dim-1] = firstSubPoint[0] + numSubPoints[0];} 2156 if (dim > 2) {firstSubPoint[dim-2] = firstSubPoint[dim-1] + numSubPoints[dim-1];} 2157 for (d = 0; d <= dim; ++d) { 2158 ierr = DMLabelGetStratumIS(subpointMap, d, &subpointIS[d]);CHKERRQ(ierr); 2159 ierr = ISGetIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 2160 } 2161 for (d = 0; d <= dim; ++d) { 2162 for (p = 0; p < numSubPoints[d]; ++p) { 2163 const PetscInt point = subpoints[d][p]; 2164 const PetscInt subpoint = firstSubPoint[d] + p; 2165 const PetscInt *cone; 2166 PetscInt coneSize, coneSizeNew, c, val; 2167 2168 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 2169 ierr = DMPlexSetConeSize(subdm, subpoint, coneSize);CHKERRQ(ierr); 2170 if (d == dim) { 2171 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 2172 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 2173 ierr = DMLabelGetValue(subpointMap, cone[c], &val);CHKERRQ(ierr); 2174 if (val >= 0) coneSizeNew++; 2175 } 2176 ierr = DMPlexSetConeSize(subdm, subpoint, coneSizeNew);CHKERRQ(ierr); 2177 } 2178 } 2179 } 2180 ierr = DMSetUp(subdm);CHKERRQ(ierr); 2181 /* Set cones */ 2182 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 2183 ierr = PetscMalloc2(maxConeSize,PetscInt,&coneNew,maxConeSize,PetscInt,&coneONew);CHKERRQ(ierr); 2184 for (d = 0; d <= dim; ++d) { 2185 for (p = 0; p < numSubPoints[d]; ++p) { 2186 const PetscInt point = subpoints[d][p]; 2187 const PetscInt subpoint = firstSubPoint[d] + p; 2188 const PetscInt *cone, *ornt; 2189 PetscInt coneSize, subconeSize, coneSizeNew, c, subc; 2190 2191 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 2192 ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr); 2193 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 2194 ierr = DMPlexGetConeOrientation(dm, point, &ornt);CHKERRQ(ierr); 2195 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 2196 ierr = PetscFindInt(cone[c], numSubPoints[d-1], subpoints[d-1], &subc);CHKERRQ(ierr); 2197 if (subc >= 0) { 2198 coneNew[coneSizeNew] = firstSubPoint[d-1] + subc; 2199 coneONew[coneSizeNew] = ornt[c]; 2200 ++coneSizeNew; 2201 } 2202 } 2203 if (coneSizeNew != subconeSize) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of cone points located %d does not match subcone size %d", coneSizeNew, subconeSize); 2204 ierr = DMPlexSetCone(subdm, subpoint, coneNew);CHKERRQ(ierr); 2205 ierr = DMPlexSetConeOrientation(subdm, subpoint, coneONew);CHKERRQ(ierr); 2206 } 2207 } 2208 ierr = PetscFree2(coneNew,coneONew);CHKERRQ(ierr); 2209 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 2210 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 2211 /* Build coordinates */ 2212 { 2213 PetscSection coordSection, subCoordSection; 2214 Vec coordinates, subCoordinates; 2215 PetscScalar *coords, *subCoords; 2216 PetscInt numComp, coordSize; 2217 2218 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2219 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 2220 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 2221 ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr); 2222 ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr); 2223 ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr); 2224 ierr = PetscSectionSetChart(subCoordSection, firstSubPoint[0], firstSubPoint[0]+numSubPoints[0]);CHKERRQ(ierr); 2225 for (v = 0; v < numSubPoints[0]; ++v) { 2226 const PetscInt vertex = subpoints[0][v]; 2227 const PetscInt subvertex = firstSubPoint[0]+v; 2228 PetscInt dof; 2229 2230 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2231 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 2232 ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr); 2233 } 2234 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 2235 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 2236 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 2237 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2238 ierr = VecSetType(subCoordinates,dm->vectype);CHKERRQ(ierr); 2239 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2240 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2241 for (v = 0; v < numSubPoints[0]; ++v) { 2242 const PetscInt vertex = subpoints[0][v]; 2243 const PetscInt subvertex = firstSubPoint[0]+v; 2244 PetscInt dof, off, sdof, soff, d; 2245 2246 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2247 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 2248 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 2249 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 2250 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 2251 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 2252 } 2253 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2254 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2255 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 2256 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 2257 } 2258 /* Cleanup */ 2259 for (d = 0; d <= dim; ++d) { 2260 ierr = ISRestoreIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 2261 ierr = ISDestroy(&subpointIS[d]);CHKERRQ(ierr); 2262 } 2263 ierr = PetscFree4(numSubPoints,firstSubPoint,subpointIS,subpoints);CHKERRQ(ierr); 2264 PetscFunctionReturn(0); 2265 } 2266 2267 #undef __FUNCT__ 2268 #define __FUNCT__ "DMPlexCreateSubmesh" 2269 /*@ 2270 DMPlexCreateSubmesh - Extract a hypersurface from the mesh using vertices defined by a label 2271 2272 Input Parameters: 2273 + dm - The original mesh 2274 . vertexLabel - The DMLabel marking vertices contained in the surface 2275 - value - The label value to use 2276 2277 Output Parameter: 2278 . subdm - The surface mesh 2279 2280 Note: This function produces a DMLabel mapping original points in the submesh to their depth. This can be obtained using DMPlexGetSubpointMap(). 2281 2282 Level: developer 2283 2284 .seealso: DMPlexGetSubpointMap(), DMPlexGetLabel(), DMLabelSetValue() 2285 @*/ 2286 PetscErrorCode DMPlexCreateSubmesh(DM dm, DMLabel vertexLabel, PetscInt value, DM *subdm) 2287 { 2288 PetscInt dim, depth; 2289 PetscErrorCode ierr; 2290 2291 PetscFunctionBegin; 2292 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2293 PetscValidPointer(subdm, 3); 2294 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 2295 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 2296 ierr = DMCreate(PetscObjectComm((PetscObject)dm), subdm);CHKERRQ(ierr); 2297 ierr = DMSetType(*subdm, DMPLEX);CHKERRQ(ierr); 2298 ierr = DMPlexSetDimension(*subdm, dim-1);CHKERRQ(ierr); 2299 if (depth == dim) { 2300 ierr = DMPlexCreateSubmesh_Interpolated(dm, vertexLabel, value, *subdm);CHKERRQ(ierr); 2301 } else { 2302 ierr = DMPlexCreateSubmesh_Uninterpolated(dm, vertexLabel, value, *subdm);CHKERRQ(ierr); 2303 } 2304 PetscFunctionReturn(0); 2305 } 2306 2307 #undef __FUNCT__ 2308 #define __FUNCT__ "DMPlexFilterPoint_Internal" 2309 PETSC_STATIC_INLINE PetscInt DMPlexFilterPoint_Internal(PetscInt point, PetscInt firstSubPoint, PetscInt numSubPoints, const PetscInt subPoints[]) 2310 { 2311 PetscInt subPoint; 2312 PetscErrorCode ierr; 2313 2314 ierr = PetscFindInt(point, numSubPoints, subPoints, &subPoint); if (ierr < 0) return ierr; 2315 return subPoint < 0 ? subPoint : firstSubPoint+subPoint; 2316 } 2317 2318 #undef __FUNCT__ 2319 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh_Uninterpolated" 2320 static PetscErrorCode DMPlexCreateCohesiveSubmesh_Uninterpolated(DM dm, PetscBool hasLagrange, const char label[], PetscInt value, DM subdm) 2321 { 2322 MPI_Comm comm; 2323 DMLabel subpointMap; 2324 IS subvertexIS; 2325 const PetscInt *subVertices; 2326 PetscInt numSubVertices, firstSubVertex, numSubCells, *subCells = NULL; 2327 PetscInt *subface, maxConeSize, numSubFaces, firstSubFace, newFacePoint, nFV; 2328 PetscInt cMax, c, f; 2329 PetscErrorCode ierr; 2330 2331 PetscFunctionBegin; 2332 ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 2333 /* Create subpointMap which marks the submesh */ 2334 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 2335 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 2336 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 2337 ierr = DMPlexMarkCohesiveSubmesh_Uninterpolated(dm, hasLagrange, label, value, subpointMap, &numSubFaces, &nFV, &subCells, subdm);CHKERRQ(ierr); 2338 /* Setup chart */ 2339 ierr = DMLabelGetStratumSize(subpointMap, 0, &numSubVertices);CHKERRQ(ierr); 2340 ierr = DMLabelGetStratumSize(subpointMap, 2, &numSubCells);CHKERRQ(ierr); 2341 ierr = DMPlexSetChart(subdm, 0, numSubCells+numSubFaces+numSubVertices);CHKERRQ(ierr); 2342 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 2343 /* Set cone sizes */ 2344 firstSubVertex = numSubCells; 2345 firstSubFace = numSubCells+numSubVertices; 2346 newFacePoint = firstSubFace; 2347 ierr = DMLabelGetStratumIS(subpointMap, 0, &subvertexIS);CHKERRQ(ierr); 2348 if (subvertexIS) {ierr = ISGetIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 2349 for (c = 0; c < numSubCells; ++c) { 2350 ierr = DMPlexSetConeSize(subdm, c, 1);CHKERRQ(ierr); 2351 } 2352 for (f = firstSubFace; f < firstSubFace+numSubFaces; ++f) { 2353 ierr = DMPlexSetConeSize(subdm, f, nFV);CHKERRQ(ierr); 2354 } 2355 ierr = DMSetUp(subdm);CHKERRQ(ierr); 2356 /* Create face cones */ 2357 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 2358 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 2359 ierr = DMGetWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 2360 for (c = 0; c < numSubCells; ++c) { 2361 const PetscInt cell = subCells[c]; 2362 const PetscInt subcell = c; 2363 const PetscInt *cone, *cells; 2364 PetscInt numCells, subVertex, p, v; 2365 2366 if (cell < cMax) continue; 2367 ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 2368 for (v = 0; v < nFV; ++v) { 2369 ierr = PetscFindInt(cone[v], numSubVertices, subVertices, &subVertex);CHKERRQ(ierr); 2370 subface[v] = firstSubVertex+subVertex; 2371 } 2372 ierr = DMPlexSetCone(subdm, newFacePoint, subface);CHKERRQ(ierr); 2373 ierr = DMPlexSetCone(subdm, subcell, &newFacePoint);CHKERRQ(ierr); 2374 ierr = DMPlexGetJoin(dm, nFV, cone, &numCells, &cells);CHKERRQ(ierr); 2375 /* Not true in parallel 2376 if (numCells != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells"); */ 2377 for (p = 0; p < numCells; ++p) { 2378 PetscInt negsubcell; 2379 2380 if (cells[p] >= cMax) continue; 2381 /* I know this is a crap search */ 2382 for (negsubcell = 0; negsubcell < numSubCells; ++negsubcell) { 2383 if (subCells[negsubcell] == cells[p]) break; 2384 } 2385 if (negsubcell == numSubCells) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find negative face neighbor for cohesive cell %d", cell); 2386 ierr = DMPlexSetCone(subdm, negsubcell, &newFacePoint);CHKERRQ(ierr); 2387 } 2388 ierr = DMPlexRestoreJoin(dm, nFV, cone, &numCells, &cells);CHKERRQ(ierr); 2389 ++newFacePoint; 2390 } 2391 ierr = DMRestoreWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 2392 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 2393 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 2394 /* Build coordinates */ 2395 { 2396 PetscSection coordSection, subCoordSection; 2397 Vec coordinates, subCoordinates; 2398 PetscScalar *coords, *subCoords; 2399 PetscInt numComp, coordSize, v; 2400 2401 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2402 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 2403 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 2404 ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr); 2405 ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr); 2406 ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr); 2407 ierr = PetscSectionSetChart(subCoordSection, firstSubVertex, firstSubVertex+numSubVertices);CHKERRQ(ierr); 2408 for (v = 0; v < numSubVertices; ++v) { 2409 const PetscInt vertex = subVertices[v]; 2410 const PetscInt subvertex = firstSubVertex+v; 2411 PetscInt dof; 2412 2413 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2414 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 2415 ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr); 2416 } 2417 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 2418 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 2419 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 2420 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2421 ierr = VecSetType(subCoordinates,dm->vectype);CHKERRQ(ierr); 2422 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2423 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2424 for (v = 0; v < numSubVertices; ++v) { 2425 const PetscInt vertex = subVertices[v]; 2426 const PetscInt subvertex = firstSubVertex+v; 2427 PetscInt dof, off, sdof, soff, d; 2428 2429 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2430 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 2431 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 2432 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 2433 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 2434 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 2435 } 2436 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2437 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2438 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 2439 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 2440 } 2441 /* Build SF */ 2442 CHKMEMQ; 2443 { 2444 PetscSF sfPoint, sfPointSub; 2445 const PetscSFNode *remotePoints; 2446 PetscSFNode *sremotePoints, *newLocalPoints, *newOwners; 2447 const PetscInt *localPoints; 2448 PetscInt *slocalPoints; 2449 PetscInt numRoots, numLeaves, numSubRoots = numSubCells+numSubFaces+numSubVertices, numSubLeaves = 0, l, sl, ll, pStart, pEnd, p, vStart, vEnd; 2450 PetscMPIInt rank; 2451 2452 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 2453 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 2454 ierr = DMGetPointSF(subdm, &sfPointSub);CHKERRQ(ierr); 2455 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 2456 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 2457 ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 2458 if (numRoots >= 0) { 2459 /* Only vertices should be shared */ 2460 ierr = PetscMalloc2(pEnd-pStart,PetscSFNode,&newLocalPoints,numRoots,PetscSFNode,&newOwners);CHKERRQ(ierr); 2461 for (p = 0; p < pEnd-pStart; ++p) { 2462 newLocalPoints[p].rank = -2; 2463 newLocalPoints[p].index = -2; 2464 } 2465 /* Set subleaves */ 2466 for (l = 0; l < numLeaves; ++l) { 2467 const PetscInt point = localPoints[l]; 2468 const PetscInt subPoint = DMPlexFilterPoint_Internal(point, firstSubVertex, numSubVertices, subVertices); 2469 2470 if ((point < vStart) && (point >= vEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Should not be mapping anything but vertices, %d", point); 2471 if (subPoint < 0) continue; 2472 newLocalPoints[point-pStart].rank = rank; 2473 newLocalPoints[point-pStart].index = subPoint; 2474 ++numSubLeaves; 2475 } 2476 /* Must put in owned subpoints */ 2477 for (p = pStart; p < pEnd; ++p) { 2478 const PetscInt subPoint = DMPlexFilterPoint_Internal(p, firstSubVertex, numSubVertices, subVertices); 2479 2480 if (subPoint < 0) { 2481 newOwners[p-pStart].rank = -3; 2482 newOwners[p-pStart].index = -3; 2483 } else { 2484 newOwners[p-pStart].rank = rank; 2485 newOwners[p-pStart].index = subPoint; 2486 } 2487 } 2488 ierr = PetscSFReduceBegin(sfPoint, MPIU_2INT, newLocalPoints, newOwners, MPI_MAXLOC);CHKERRQ(ierr); 2489 ierr = PetscSFReduceEnd(sfPoint, MPIU_2INT, newLocalPoints, newOwners, MPI_MAXLOC);CHKERRQ(ierr); 2490 ierr = PetscSFBcastBegin(sfPoint, MPIU_2INT, newOwners, newLocalPoints);CHKERRQ(ierr); 2491 ierr = PetscSFBcastEnd(sfPoint, MPIU_2INT, newOwners, newLocalPoints);CHKERRQ(ierr); 2492 ierr = PetscMalloc(numSubLeaves * sizeof(PetscInt), &slocalPoints);CHKERRQ(ierr); 2493 ierr = PetscMalloc(numSubLeaves * sizeof(PetscSFNode), &sremotePoints);CHKERRQ(ierr); 2494 for (l = 0, sl = 0, ll = 0; l < numLeaves; ++l) { 2495 const PetscInt point = localPoints[l]; 2496 const PetscInt subPoint = DMPlexFilterPoint_Internal(point, firstSubVertex, numSubVertices, subVertices); 2497 2498 if (subPoint < 0) continue; 2499 if (newLocalPoints[point].rank == rank) {++ll; continue;} 2500 slocalPoints[sl] = subPoint; 2501 sremotePoints[sl].rank = newLocalPoints[point].rank; 2502 sremotePoints[sl].index = newLocalPoints[point].index; 2503 if (sremotePoints[sl].rank < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid remote rank for local point %d", point); 2504 if (sremotePoints[sl].index < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid remote subpoint for local point %d", point); 2505 ++sl; 2506 } 2507 ierr = PetscFree2(newLocalPoints,newOwners);CHKERRQ(ierr); 2508 if (sl + ll != numSubLeaves) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Mismatch in number of subleaves %d + %d != %d", sl, ll, numSubLeaves); 2509 ierr = PetscSFSetGraph(sfPointSub, numSubRoots, sl, slocalPoints, PETSC_OWN_POINTER, sremotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 2510 } 2511 } 2512 CHKMEMQ; 2513 /* Cleanup */ 2514 if (subvertexIS) {ierr = ISRestoreIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 2515 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 2516 ierr = PetscFree(subCells);CHKERRQ(ierr); 2517 PetscFunctionReturn(0); 2518 } 2519 2520 #undef __FUNCT__ 2521 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh_Interpolated" 2522 static PetscErrorCode DMPlexCreateCohesiveSubmesh_Interpolated(DM dm, PetscBool hasLagrange, const char label[], PetscInt value, DM subdm) 2523 { 2524 MPI_Comm comm; 2525 DMLabel subpointMap; 2526 IS *subpointIS; 2527 const PetscInt **subpoints; 2528 PetscInt *numSubPoints, *firstSubPoint, *coneNew; 2529 PetscInt totSubPoints = 0, maxConeSize, dim, p, d, v; 2530 PetscErrorCode ierr; 2531 2532 PetscFunctionBegin; 2533 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2534 /* Create subpointMap which marks the submesh */ 2535 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 2536 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 2537 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 2538 ierr = DMPlexMarkCohesiveSubmesh_Interpolated(dm, hasLagrange, label, value, subpointMap, subdm);CHKERRQ(ierr); 2539 /* Setup chart */ 2540 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 2541 ierr = PetscMalloc4(dim+1,PetscInt,&numSubPoints,dim+1,PetscInt,&firstSubPoint,dim+1,IS,&subpointIS,dim+1,const PetscInt *,&subpoints);CHKERRQ(ierr); 2542 for (d = 0; d <= dim; ++d) { 2543 ierr = DMLabelGetStratumSize(subpointMap, d, &numSubPoints[d]);CHKERRQ(ierr); 2544 totSubPoints += numSubPoints[d]; 2545 } 2546 ierr = DMPlexSetChart(subdm, 0, totSubPoints);CHKERRQ(ierr); 2547 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 2548 /* Set cone sizes */ 2549 firstSubPoint[dim] = 0; 2550 firstSubPoint[0] = firstSubPoint[dim] + numSubPoints[dim]; 2551 if (dim > 1) {firstSubPoint[dim-1] = firstSubPoint[0] + numSubPoints[0];} 2552 if (dim > 2) {firstSubPoint[dim-2] = firstSubPoint[dim-1] + numSubPoints[dim-1];} 2553 for (d = 0; d <= dim; ++d) { 2554 ierr = DMLabelGetStratumIS(subpointMap, d, &subpointIS[d]);CHKERRQ(ierr); 2555 ierr = ISGetIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 2556 } 2557 for (d = 0; d <= dim; ++d) { 2558 for (p = 0; p < numSubPoints[d]; ++p) { 2559 const PetscInt point = subpoints[d][p]; 2560 const PetscInt subpoint = firstSubPoint[d] + p; 2561 const PetscInt *cone; 2562 PetscInt coneSize, coneSizeNew, c, val; 2563 2564 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 2565 ierr = DMPlexSetConeSize(subdm, subpoint, coneSize);CHKERRQ(ierr); 2566 if (d == dim) { 2567 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 2568 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 2569 ierr = DMLabelGetValue(subpointMap, cone[c], &val);CHKERRQ(ierr); 2570 if (val >= 0) coneSizeNew++; 2571 } 2572 ierr = DMPlexSetConeSize(subdm, subpoint, coneSizeNew);CHKERRQ(ierr); 2573 } 2574 } 2575 } 2576 ierr = DMSetUp(subdm);CHKERRQ(ierr); 2577 /* Set cones */ 2578 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 2579 ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &coneNew);CHKERRQ(ierr); 2580 for (d = 0; d <= dim; ++d) { 2581 for (p = 0; p < numSubPoints[d]; ++p) { 2582 const PetscInt point = subpoints[d][p]; 2583 const PetscInt subpoint = firstSubPoint[d] + p; 2584 const PetscInt *cone; 2585 PetscInt coneSize, subconeSize, coneSizeNew, c, subc; 2586 2587 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 2588 ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr); 2589 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 2590 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 2591 ierr = PetscFindInt(cone[c], numSubPoints[d-1], subpoints[d-1], &subc);CHKERRQ(ierr); 2592 if (subc >= 0) coneNew[coneSizeNew++] = firstSubPoint[d-1] + subc; 2593 } 2594 if (coneSizeNew != subconeSize) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of cone points located %d does not match subcone size %d", coneSizeNew, subconeSize); 2595 ierr = DMPlexSetCone(subdm, subpoint, coneNew);CHKERRQ(ierr); 2596 } 2597 } 2598 ierr = PetscFree(coneNew);CHKERRQ(ierr); 2599 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 2600 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 2601 /* Build coordinates */ 2602 { 2603 PetscSection coordSection, subCoordSection; 2604 Vec coordinates, subCoordinates; 2605 PetscScalar *coords, *subCoords; 2606 PetscInt numComp, coordSize; 2607 2608 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2609 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 2610 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 2611 ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr); 2612 ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr); 2613 ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr); 2614 ierr = PetscSectionSetChart(subCoordSection, firstSubPoint[0], firstSubPoint[0]+numSubPoints[0]);CHKERRQ(ierr); 2615 for (v = 0; v < numSubPoints[0]; ++v) { 2616 const PetscInt vertex = subpoints[0][v]; 2617 const PetscInt subvertex = firstSubPoint[0]+v; 2618 PetscInt dof; 2619 2620 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2621 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 2622 ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr); 2623 } 2624 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 2625 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 2626 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 2627 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2628 ierr = VecSetType(subCoordinates,dm->vectype);CHKERRQ(ierr); 2629 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2630 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2631 for (v = 0; v < numSubPoints[0]; ++v) { 2632 const PetscInt vertex = subpoints[0][v]; 2633 const PetscInt subvertex = firstSubPoint[0]+v; 2634 PetscInt dof, off, sdof, soff, d; 2635 2636 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2637 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 2638 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 2639 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 2640 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 2641 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 2642 } 2643 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2644 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2645 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 2646 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 2647 } 2648 /* Cleanup */ 2649 for (d = 0; d <= dim; ++d) { 2650 ierr = ISRestoreIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 2651 ierr = ISDestroy(&subpointIS[d]);CHKERRQ(ierr); 2652 } 2653 ierr = PetscFree4(numSubPoints,firstSubPoint,subpointIS,subpoints);CHKERRQ(ierr); 2654 PetscFunctionReturn(0); 2655 } 2656 2657 #undef __FUNCT__ 2658 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh" 2659 /* 2660 DMPlexCreateCohesiveSubmesh - Extract from a mesh with cohesive cells the hypersurface defined by one face of the cells. Optionally, a Label an be given to restrict the cells. 2661 2662 Input Parameters: 2663 + dm - The original mesh 2664 . hasLagrange - The mesh has Lagrange unknowns in the cohesive cells 2665 . label - A label name, or PETSC_NULL 2666 - value - A label value 2667 2668 Output Parameter: 2669 . subdm - The surface mesh 2670 2671 Note: This function produces a DMLabel mapping original points in the submesh to their depth. This can be obtained using DMPlexGetSubpointMap(). 2672 2673 Level: developer 2674 2675 .seealso: DMPlexGetSubpointMap(), DMPlexCreateSubmesh() 2676 */ 2677 PetscErrorCode DMPlexCreateCohesiveSubmesh(DM dm, PetscBool hasLagrange, const char label[], PetscInt value, DM *subdm) 2678 { 2679 PetscInt dim, depth; 2680 PetscErrorCode ierr; 2681 2682 PetscFunctionBegin; 2683 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2684 PetscValidPointer(subdm, 5); 2685 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 2686 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 2687 ierr = DMCreate(PetscObjectComm((PetscObject)dm), subdm);CHKERRQ(ierr); 2688 ierr = DMSetType(*subdm, DMPLEX);CHKERRQ(ierr); 2689 ierr = DMPlexSetDimension(*subdm, dim-1);CHKERRQ(ierr); 2690 if (depth == dim) { 2691 ierr = DMPlexCreateCohesiveSubmesh_Interpolated(dm, hasLagrange, label, value, *subdm);CHKERRQ(ierr); 2692 } else { 2693 ierr = DMPlexCreateCohesiveSubmesh_Uninterpolated(dm, hasLagrange, label, value, *subdm);CHKERRQ(ierr); 2694 } 2695 PetscFunctionReturn(0); 2696 } 2697 2698 #undef __FUNCT__ 2699 #define __FUNCT__ "DMPlexGetSubpointMap" 2700 /*@ 2701 DMPlexGetSubpointMap - Returns a DMLabel with point dimension as values 2702 2703 Input Parameter: 2704 . dm - The submesh DM 2705 2706 Output Parameter: 2707 . subpointMap - The DMLabel of all the points from the original mesh in this submesh, or NULL if this is not a submesh 2708 2709 Level: developer 2710 2711 .seealso: DMPlexCreateSubmesh(), DMPlexCreateSubpointIS() 2712 @*/ 2713 PetscErrorCode DMPlexGetSubpointMap(DM dm, DMLabel *subpointMap) 2714 { 2715 DM_Plex *mesh = (DM_Plex*) dm->data; 2716 2717 PetscFunctionBegin; 2718 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2719 PetscValidPointer(subpointMap, 2); 2720 *subpointMap = mesh->subpointMap; 2721 PetscFunctionReturn(0); 2722 } 2723 2724 #undef __FUNCT__ 2725 #define __FUNCT__ "DMPlexSetSubpointMap" 2726 /* Note: Should normally not be called by the user, since it is set in DMPlexCreateSubmesh() */ 2727 PetscErrorCode DMPlexSetSubpointMap(DM dm, DMLabel subpointMap) 2728 { 2729 DM_Plex *mesh = (DM_Plex *) dm->data; 2730 DMLabel tmp; 2731 PetscErrorCode ierr; 2732 2733 PetscFunctionBegin; 2734 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2735 tmp = mesh->subpointMap; 2736 mesh->subpointMap = subpointMap; 2737 ++mesh->subpointMap->refct; 2738 ierr = DMLabelDestroy(&tmp);CHKERRQ(ierr); 2739 PetscFunctionReturn(0); 2740 } 2741 2742 #undef __FUNCT__ 2743 #define __FUNCT__ "DMPlexCreateSubpointIS" 2744 /*@ 2745 DMPlexCreateSubpointIS - Creates an IS covering the entire subdm chart with the original points as data 2746 2747 Input Parameter: 2748 . dm - The submesh DM 2749 2750 Output Parameter: 2751 . subpointIS - The IS of all the points from the original mesh in this submesh, or NULL if this is not a submesh 2752 2753 Note: This is IS is guaranteed to be sorted by the construction of the submesh 2754 2755 Level: developer 2756 2757 .seealso: DMPlexCreateSubmesh(), DMPlexGetSubpointMap() 2758 @*/ 2759 PetscErrorCode DMPlexCreateSubpointIS(DM dm, IS *subpointIS) 2760 { 2761 MPI_Comm comm; 2762 DMLabel subpointMap; 2763 IS is; 2764 const PetscInt *opoints; 2765 PetscInt *points, *depths; 2766 PetscInt depth, depStart, depEnd, d, pStart, pEnd, p, n, off; 2767 PetscErrorCode ierr; 2768 2769 PetscFunctionBegin; 2770 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2771 PetscValidPointer(subpointIS, 2); 2772 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2773 *subpointIS = NULL; 2774 ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 2775 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 2776 if (subpointMap && depth >= 0) { 2777 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 2778 if (pStart) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Submeshes must start the point numbering at 0, not %d", pStart); 2779 ierr = DMGetWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr); 2780 depths[0] = depth; 2781 depths[1] = 0; 2782 for(d = 2; d <= depth; ++d) {depths[d] = depth+1 - d;} 2783 ierr = PetscMalloc(pEnd * sizeof(PetscInt), &points);CHKERRQ(ierr); 2784 for(d = 0, off = 0; d <= depth; ++d) { 2785 const PetscInt dep = depths[d]; 2786 2787 ierr = DMPlexGetDepthStratum(dm, dep, &depStart, &depEnd);CHKERRQ(ierr); 2788 ierr = DMLabelGetStratumSize(subpointMap, dep, &n);CHKERRQ(ierr); 2789 if (((d < 2) && (depth > 1)) || (d == 1)) { /* Only check vertices and cells for now since the map is broken for others */ 2790 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); 2791 } else { 2792 if (!n) { 2793 if (d == 0) { 2794 /* Missing cells */ 2795 for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = -1; 2796 } else { 2797 /* Missing faces */ 2798 for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = PETSC_MAX_INT; 2799 } 2800 } 2801 } 2802 if (n) { 2803 ierr = DMLabelGetStratumIS(subpointMap, dep, &is);CHKERRQ(ierr); 2804 ierr = ISGetIndices(is, &opoints);CHKERRQ(ierr); 2805 for(p = 0; p < n; ++p, ++off) points[off] = opoints[p]; 2806 ierr = ISRestoreIndices(is, &opoints);CHKERRQ(ierr); 2807 ierr = ISDestroy(&is);CHKERRQ(ierr); 2808 } 2809 } 2810 ierr = DMRestoreWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr); 2811 if (off != pEnd) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "The number of mapped submesh points %d should be %d", off, pEnd); 2812 ierr = ISCreateGeneral(PETSC_COMM_SELF, pEnd, points, PETSC_OWN_POINTER, subpointIS);CHKERRQ(ierr); 2813 } 2814 PetscFunctionReturn(0); 2815 } 2816