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,VECSTANDARD);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 /* All points connected to the fault are inside a cell, so at the top level we will only check cells */ 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 = again == 1 ? 1 : 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 PetscInt side, closureSize, cl, st, *closure = NULL; 1155 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); 1156 if (val > 0) side = 1; 1157 else side = -1; 1158 ierr = DMLabelSetValue(label, point, side*(shift+dim));CHKERRQ(ierr); 1159 /* Mark all other cell parts */ 1160 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1161 for (cl = 0; cl < closureSize*2; cl += 2) { 1162 const PetscInt clpoint = closure[cl]; 1163 for (st = 0; st < starSize*2; st += 2) { 1164 if (star[st] == clpoint) { 1165 ierr = DMLabelGetValue(label, clpoint, &val);CHKERRQ(ierr); 1166 if (val != -1) continue; 1167 ierr = DMLabelGetValue(depthLabel, clpoint, &dep);CHKERRQ(ierr); 1168 ierr = DMLabelSetValue(label, clpoint, side*(shift+dep));CHKERRQ(ierr); 1169 break; 1170 } 1171 } 1172 } 1173 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1174 again = 1; 1175 break; 1176 } 1177 } 1178 } 1179 } 1180 /* Classify the rest by cell membership */ 1181 for (s = 0; s < starSize*2; s += 2) { 1182 const PetscInt point = star[s]; 1183 1184 ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr); 1185 if (val == -1) { 1186 PetscInt *sstar = NULL; 1187 PetscInt sstarSize, ss; 1188 PetscBool marked = PETSC_FALSE; 1189 1190 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &sstarSize, &sstar);CHKERRQ(ierr); 1191 for (ss = 0; ss < sstarSize*2; ss += 2) { 1192 const PetscInt spoint = sstar[ss]; 1193 1194 if ((spoint < cStart) || (spoint >= cEnd)) continue; 1195 ierr = DMLabelGetValue(label, spoint, &val);CHKERRQ(ierr); 1196 if (val == -1) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Cell %d in star of %d does not have a valid label", spoint, point); 1197 ierr = DMLabelGetValue(depthLabel, point, &dep);CHKERRQ(ierr); 1198 if (val > 0) { 1199 ierr = DMLabelSetValue(label, point, shift+dep);CHKERRQ(ierr); 1200 } else { 1201 ierr = DMLabelSetValue(label, point, -(shift+dep));CHKERRQ(ierr); 1202 } 1203 marked = PETSC_TRUE; 1204 break; 1205 } 1206 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &sstarSize, &sstar);CHKERRQ(ierr); 1207 if (!marked) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Point %d could not be classified", point); 1208 } 1209 } 1210 ierr = DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1211 } 1212 ierr = ISRestoreIndices(dimIS, &points);CHKERRQ(ierr); 1213 ierr = ISDestroy(&dimIS);CHKERRQ(ierr); 1214 PetscFunctionReturn(0); 1215 } 1216 1217 #undef __FUNCT__ 1218 #define __FUNCT__ "DMPlexCreateHybridMesh" 1219 /*@C 1220 DMPlexCreateHybridMesh - Create a mesh with hybrid cells along an internal interface 1221 1222 Collective on dm 1223 1224 Input Parameters: 1225 + dm - The original DM 1226 - labelName - The label specifying the interface vertices 1227 1228 Output Parameters: 1229 + hybridLabel - The label fully marking the interface 1230 - dmHybrid - The new DM 1231 1232 Level: developer 1233 1234 .seealso: DMPlexConstructCohesiveCells(), DMPlexLabelCohesiveComplete(), DMCreate() 1235 @*/ 1236 PetscErrorCode DMPlexCreateHybridMesh(DM dm, DMLabel label, DMLabel *hybridLabel, DM *dmHybrid) 1237 { 1238 DM idm; 1239 DMLabel subpointMap, hlabel; 1240 PetscInt dim; 1241 PetscErrorCode ierr; 1242 1243 PetscFunctionBegin; 1244 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1245 if (hybridLabel) PetscValidPointer(hybridLabel, 3); 1246 PetscValidPointer(dmHybrid, 4); 1247 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1248 ierr = DMPlexCreateSubmesh(dm, label, 1, &idm);CHKERRQ(ierr); 1249 ierr = DMPlexOrient(idm);CHKERRQ(ierr); 1250 ierr = DMPlexGetSubpointMap(idm, &subpointMap);CHKERRQ(ierr); 1251 ierr = DMLabelDuplicate(subpointMap, &hlabel);CHKERRQ(ierr); 1252 ierr = DMLabelClearStratum(hlabel, dim);CHKERRQ(ierr); 1253 ierr = DMPlexLabelCohesiveComplete(dm, hlabel, PETSC_FALSE, idm);CHKERRQ(ierr); 1254 ierr = DMDestroy(&idm);CHKERRQ(ierr); 1255 ierr = DMPlexConstructCohesiveCells(dm, hlabel, dmHybrid);CHKERRQ(ierr); 1256 if (hybridLabel) *hybridLabel = hlabel; 1257 else {ierr = DMLabelDestroy(&hlabel);CHKERRQ(ierr);} 1258 PetscFunctionReturn(0); 1259 } 1260 1261 #undef __FUNCT__ 1262 #define __FUNCT__ "DMPlexMarkSubmesh_Uninterpolated" 1263 /* Here we need the explicit assumption that: 1264 1265 For any marked cell, the marked vertices constitute a single face 1266 */ 1267 static PetscErrorCode DMPlexMarkSubmesh_Uninterpolated(DM dm, DMLabel vertexLabel, PetscInt value, DMLabel subpointMap, PetscInt *numFaces, PetscInt *nFV, DM subdm) 1268 { 1269 IS subvertexIS = NULL; 1270 const PetscInt *subvertices; 1271 PetscInt *pStart, *pEnd, *pMax, pSize; 1272 PetscInt depth, dim, d, numSubVerticesInitial = 0, v; 1273 PetscErrorCode ierr; 1274 1275 PetscFunctionBegin; 1276 *numFaces = 0; 1277 *nFV = 0; 1278 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1279 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1280 pSize = PetscMax(depth, dim) + 1; 1281 ierr = PetscMalloc3(pSize,PetscInt,&pStart,pSize,PetscInt,&pEnd,pSize,PetscInt,&pMax);CHKERRQ(ierr); 1282 ierr = DMPlexGetHybridBounds(dm, depth >= 0 ? &pMax[depth] : NULL, depth>1 ? &pMax[depth-1] : NULL, depth>2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr); 1283 for (d = 0; d <= depth; ++d) { 1284 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1285 if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]); 1286 } 1287 /* Loop over initial vertices and mark all faces in the collective star() */ 1288 if (vertexLabel) {ierr = DMLabelGetStratumIS(vertexLabel, value, &subvertexIS);CHKERRQ(ierr);} 1289 if (subvertexIS) { 1290 ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr); 1291 ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1292 } 1293 for (v = 0; v < numSubVerticesInitial; ++v) { 1294 const PetscInt vertex = subvertices[v]; 1295 PetscInt *star = NULL; 1296 PetscInt starSize, s, numCells = 0, c; 1297 1298 ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1299 for (s = 0; s < starSize*2; s += 2) { 1300 const PetscInt point = star[s]; 1301 if ((point >= pStart[depth]) && (point < pEnd[depth])) star[numCells++] = point; 1302 } 1303 for (c = 0; c < numCells; ++c) { 1304 const PetscInt cell = star[c]; 1305 PetscInt *closure = NULL; 1306 PetscInt closureSize, cl; 1307 PetscInt cellLoc, numCorners = 0, faceSize = 0; 1308 1309 ierr = DMLabelGetValue(subpointMap, cell, &cellLoc);CHKERRQ(ierr); 1310 if (cellLoc == 2) continue; 1311 if (cellLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Cell %d has dimension %d in the surface label", cell, cellLoc); 1312 ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1313 for (cl = 0; cl < closureSize*2; cl += 2) { 1314 const PetscInt point = closure[cl]; 1315 PetscInt vertexLoc; 1316 1317 if ((point >= pStart[0]) && (point < pEnd[0])) { 1318 ++numCorners; 1319 ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr); 1320 if (vertexLoc == value) closure[faceSize++] = point; 1321 } 1322 } 1323 if (!(*nFV)) {ierr = DMPlexGetNumFaceVertices(dm, dim, numCorners, nFV);CHKERRQ(ierr);} 1324 if (faceSize > *nFV) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize); 1325 if (faceSize == *nFV) { 1326 const PetscInt *cells = NULL; 1327 PetscInt numCells, nc; 1328 1329 ++(*numFaces); 1330 for (cl = 0; cl < faceSize; ++cl) { 1331 ierr = DMLabelSetValue(subpointMap, closure[cl], 0);CHKERRQ(ierr); 1332 } 1333 ierr = DMPlexGetJoin(dm, faceSize, closure, &numCells, &cells);CHKERRQ(ierr); 1334 for (nc = 0; nc < numCells; ++nc) { 1335 ierr = DMLabelSetValue(subpointMap, cells[nc], 2);CHKERRQ(ierr); 1336 } 1337 ierr = DMPlexRestoreJoin(dm, faceSize, closure, &numCells, &cells);CHKERRQ(ierr); 1338 } 1339 ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1340 } 1341 ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1342 } 1343 if (subvertexIS) { 1344 ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1345 } 1346 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 1347 ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr); 1348 PetscFunctionReturn(0); 1349 } 1350 1351 #undef __FUNCT__ 1352 #define __FUNCT__ "DMPlexMarkSubmesh_Interpolated" 1353 static PetscErrorCode DMPlexMarkSubmesh_Interpolated(DM dm, DMLabel vertexLabel, PetscInt value, DMLabel subpointMap, DM subdm) 1354 { 1355 IS subvertexIS; 1356 const PetscInt *subvertices; 1357 PetscInt *pStart, *pEnd, *pMax; 1358 PetscInt dim, d, numSubVerticesInitial = 0, v; 1359 PetscErrorCode ierr; 1360 1361 PetscFunctionBegin; 1362 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1363 ierr = PetscMalloc3(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd,dim+1,PetscInt,&pMax);CHKERRQ(ierr); 1364 ierr = DMPlexGetHybridBounds(dm, &pMax[dim], dim>1 ? &pMax[dim-1] : NULL, dim > 2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr); 1365 for (d = 0; d <= dim; ++d) { 1366 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1367 if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]); 1368 } 1369 /* Loop over initial vertices and mark all faces in the collective star() */ 1370 ierr = DMLabelGetStratumIS(vertexLabel, value, &subvertexIS);CHKERRQ(ierr); 1371 if (subvertexIS) { 1372 ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr); 1373 ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1374 } 1375 for (v = 0; v < numSubVerticesInitial; ++v) { 1376 const PetscInt vertex = subvertices[v]; 1377 PetscInt *star = NULL; 1378 PetscInt starSize, s, numFaces = 0, f; 1379 1380 ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1381 for (s = 0; s < starSize*2; s += 2) { 1382 const PetscInt point = star[s]; 1383 if ((point >= pStart[dim-1]) && (point < pEnd[dim-1])) star[numFaces++] = point; 1384 } 1385 for (f = 0; f < numFaces; ++f) { 1386 const PetscInt face = star[f]; 1387 PetscInt *closure = NULL; 1388 PetscInt closureSize, c; 1389 PetscInt faceLoc; 1390 1391 ierr = DMLabelGetValue(subpointMap, face, &faceLoc);CHKERRQ(ierr); 1392 if (faceLoc == dim-1) continue; 1393 if (faceLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Face %d has dimension %d in the surface label", face, faceLoc); 1394 ierr = DMPlexGetTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1395 for (c = 0; c < closureSize*2; c += 2) { 1396 const PetscInt point = closure[c]; 1397 PetscInt vertexLoc; 1398 1399 if ((point >= pStart[0]) && (point < pEnd[0])) { 1400 ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr); 1401 if (vertexLoc != value) break; 1402 } 1403 } 1404 if (c == closureSize*2) { 1405 const PetscInt *support; 1406 PetscInt supportSize, s; 1407 1408 for (c = 0; c < closureSize*2; c += 2) { 1409 const PetscInt point = closure[c]; 1410 1411 for (d = 0; d < dim; ++d) { 1412 if ((point >= pStart[d]) && (point < pEnd[d])) { 1413 ierr = DMLabelSetValue(subpointMap, point, d);CHKERRQ(ierr); 1414 break; 1415 } 1416 } 1417 } 1418 ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 1419 ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); 1420 for (s = 0; s < supportSize; ++s) { 1421 ierr = DMLabelSetValue(subpointMap, support[s], dim);CHKERRQ(ierr); 1422 } 1423 } 1424 ierr = DMPlexRestoreTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1425 } 1426 ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1427 } 1428 if (subvertexIS) { 1429 ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1430 } 1431 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 1432 ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr); 1433 PetscFunctionReturn(0); 1434 } 1435 1436 #undef __FUNCT__ 1437 #define __FUNCT__ "DMPlexMarkCohesiveSubmesh_Uninterpolated" 1438 static PetscErrorCode DMPlexMarkCohesiveSubmesh_Uninterpolated(DM dm, PetscBool hasLagrange, const char labelname[], PetscInt value, DMLabel subpointMap, PetscInt *numFaces, PetscInt *nFV, PetscInt *subCells[], DM subdm) 1439 { 1440 DMLabel label = NULL; 1441 const PetscInt *cone; 1442 PetscInt dim, cMax, cEnd, c, subc = 0, p, coneSize; 1443 PetscErrorCode ierr; 1444 1445 PetscFunctionBegin; 1446 *numFaces = 0; 1447 *nFV = 0; 1448 if (labelname) {ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);} 1449 *subCells = NULL; 1450 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1451 ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr); 1452 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1453 if (cMax < 0) PetscFunctionReturn(0); 1454 if (label) { 1455 for (c = cMax; c < cEnd; ++c) { 1456 PetscInt val; 1457 1458 ierr = DMLabelGetValue(label, c, &val);CHKERRQ(ierr); 1459 if (val == value) { 1460 ++(*numFaces); 1461 ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 1462 } 1463 } 1464 } else { 1465 *numFaces = cEnd - cMax; 1466 ierr = DMPlexGetConeSize(dm, cMax, &coneSize);CHKERRQ(ierr); 1467 } 1468 *nFV = hasLagrange ? coneSize/3 : coneSize/2; 1469 ierr = PetscMalloc(*numFaces *2 * sizeof(PetscInt), subCells);CHKERRQ(ierr); 1470 for (c = cMax; c < cEnd; ++c) { 1471 const PetscInt *cells; 1472 PetscInt numCells; 1473 1474 if (label) { 1475 PetscInt val; 1476 1477 ierr = DMLabelGetValue(label, c, &val);CHKERRQ(ierr); 1478 if (val != value) continue; 1479 } 1480 ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 1481 for (p = 0; p < *nFV; ++p) { 1482 ierr = DMLabelSetValue(subpointMap, cone[p], 0);CHKERRQ(ierr); 1483 } 1484 /* Negative face */ 1485 ierr = DMPlexGetJoin(dm, *nFV, cone, &numCells, &cells);CHKERRQ(ierr); 1486 /* Not true in parallel 1487 if (numCells != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells"); */ 1488 for (p = 0; p < numCells; ++p) { 1489 ierr = DMLabelSetValue(subpointMap, cells[p], 2);CHKERRQ(ierr); 1490 (*subCells)[subc++] = cells[p]; 1491 } 1492 ierr = DMPlexRestoreJoin(dm, *nFV, cone, &numCells, &cells);CHKERRQ(ierr); 1493 /* Positive face is not included */ 1494 } 1495 PetscFunctionReturn(0); 1496 } 1497 1498 #undef __FUNCT__ 1499 #define __FUNCT__ "DMPlexMarkCohesiveSubmesh_Interpolated" 1500 static PetscErrorCode DMPlexMarkCohesiveSubmesh_Interpolated(DM dm, PetscBool hasLagrange, const char labelname[], PetscInt value, DMLabel subpointMap, DM subdm) 1501 { 1502 DMLabel label = NULL; 1503 PetscInt *pStart, *pEnd; 1504 PetscInt dim, cMax, cEnd, c, d; 1505 PetscErrorCode ierr; 1506 1507 PetscFunctionBegin; 1508 if (labelname) {ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);} 1509 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1510 ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr); 1511 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1512 if (cMax < 0) PetscFunctionReturn(0); 1513 ierr = PetscMalloc2(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd);CHKERRQ(ierr); 1514 for (d = 0; d <= dim; ++d) { 1515 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1516 } 1517 for (c = cMax; c < cEnd; ++c) { 1518 const PetscInt *cone; 1519 PetscInt *closure = NULL; 1520 PetscInt coneSize, closureSize, cl; 1521 1522 if (label) { 1523 PetscInt val; 1524 1525 ierr = DMLabelGetValue(label, c, &val);CHKERRQ(ierr); 1526 if (val != value) continue; 1527 } 1528 ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 1529 if (hasLagrange) { 1530 if (coneSize != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells"); 1531 } else { 1532 if (coneSize != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells"); 1533 } 1534 /* Negative face */ 1535 ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 1536 ierr = DMPlexGetTransitiveClosure(dm, cone[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1537 for (cl = 0; cl < closureSize*2; cl += 2) { 1538 const PetscInt point = closure[cl]; 1539 1540 for (d = 0; d <= dim; ++d) { 1541 if ((point >= pStart[d]) && (point < pEnd[d])) { 1542 ierr = DMLabelSetValue(subpointMap, point, d);CHKERRQ(ierr); 1543 break; 1544 } 1545 } 1546 } 1547 ierr = DMPlexRestoreTransitiveClosure(dm, cone[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1548 /* Cells -- positive face is not included */ 1549 for (cl = 0; cl < 1; ++cl) { 1550 const PetscInt *support; 1551 PetscInt supportSize, s; 1552 1553 ierr = DMPlexGetSupportSize(dm, cone[cl], &supportSize);CHKERRQ(ierr); 1554 if (supportSize != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive faces should separate two cells"); 1555 ierr = DMPlexGetSupport(dm, cone[cl], &support);CHKERRQ(ierr); 1556 for (s = 0; s < supportSize; ++s) { 1557 ierr = DMLabelSetValue(subpointMap, support[s], dim);CHKERRQ(ierr); 1558 } 1559 } 1560 } 1561 ierr = PetscFree2(pStart, pEnd);CHKERRQ(ierr); 1562 PetscFunctionReturn(0); 1563 } 1564 1565 #undef __FUNCT__ 1566 #define __FUNCT__ "DMPlexGetFaceOrientation" 1567 PetscErrorCode DMPlexGetFaceOrientation(DM dm, PetscInt cell, PetscInt numCorners, PetscInt indices[], PetscInt oppositeVertex, PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented) 1568 { 1569 MPI_Comm comm; 1570 PetscBool posOrient = PETSC_FALSE; 1571 const PetscInt debug = 0; 1572 PetscInt cellDim, faceSize, f; 1573 PetscErrorCode ierr; 1574 1575 PetscFunctionBegin; 1576 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1577 ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr); 1578 if (debug) {PetscPrintf(comm, "cellDim: %d numCorners: %d\n", cellDim, numCorners);CHKERRQ(ierr);} 1579 1580 if (cellDim == 1 && numCorners == 2) { 1581 /* Triangle */ 1582 faceSize = numCorners-1; 1583 posOrient = !(oppositeVertex%2) ? PETSC_TRUE : PETSC_FALSE; 1584 } else if (cellDim == 2 && numCorners == 3) { 1585 /* Triangle */ 1586 faceSize = numCorners-1; 1587 posOrient = !(oppositeVertex%2) ? PETSC_TRUE : PETSC_FALSE; 1588 } else if (cellDim == 3 && numCorners == 4) { 1589 /* Tetrahedron */ 1590 faceSize = numCorners-1; 1591 posOrient = (oppositeVertex%2) ? PETSC_TRUE : PETSC_FALSE; 1592 } else if (cellDim == 1 && numCorners == 3) { 1593 /* Quadratic line */ 1594 faceSize = 1; 1595 posOrient = PETSC_TRUE; 1596 } else if (cellDim == 2 && numCorners == 4) { 1597 /* Quads */ 1598 faceSize = 2; 1599 if ((indices[1] > indices[0]) && (indices[1] - indices[0] == 1)) { 1600 posOrient = PETSC_TRUE; 1601 } else if ((indices[0] == 3) && (indices[1] == 0)) { 1602 posOrient = PETSC_TRUE; 1603 } else { 1604 if (((indices[0] > indices[1]) && (indices[0] - indices[1] == 1)) || ((indices[0] == 0) && (indices[1] == 3))) { 1605 posOrient = PETSC_FALSE; 1606 } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossedge"); 1607 } 1608 } else if (cellDim == 2 && numCorners == 6) { 1609 /* Quadratic triangle (I hate this) */ 1610 /* Edges are determined by the first 2 vertices (corners of edges) */ 1611 const PetscInt faceSizeTri = 3; 1612 PetscInt sortedIndices[3], i, iFace; 1613 PetscBool found = PETSC_FALSE; 1614 PetscInt faceVerticesTriSorted[9] = { 1615 0, 3, 4, /* bottom */ 1616 1, 4, 5, /* right */ 1617 2, 3, 5, /* left */ 1618 }; 1619 PetscInt faceVerticesTri[9] = { 1620 0, 3, 4, /* bottom */ 1621 1, 4, 5, /* right */ 1622 2, 5, 3, /* left */ 1623 }; 1624 1625 faceSize = faceSizeTri; 1626 for (i = 0; i < faceSizeTri; ++i) sortedIndices[i] = indices[i]; 1627 ierr = PetscSortInt(faceSizeTri, sortedIndices);CHKERRQ(ierr); 1628 for (iFace = 0; iFace < 3; ++iFace) { 1629 const PetscInt ii = iFace*faceSizeTri; 1630 PetscInt fVertex, cVertex; 1631 1632 if ((sortedIndices[0] == faceVerticesTriSorted[ii+0]) && 1633 (sortedIndices[1] == faceVerticesTriSorted[ii+1])) { 1634 for (fVertex = 0; fVertex < faceSizeTri; ++fVertex) { 1635 for (cVertex = 0; cVertex < faceSizeTri; ++cVertex) { 1636 if (indices[cVertex] == faceVerticesTri[ii+fVertex]) { 1637 faceVertices[fVertex] = origVertices[cVertex]; 1638 break; 1639 } 1640 } 1641 } 1642 found = PETSC_TRUE; 1643 break; 1644 } 1645 } 1646 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tri crossface"); 1647 if (posOriented) *posOriented = PETSC_TRUE; 1648 PetscFunctionReturn(0); 1649 } else if (cellDim == 2 && numCorners == 9) { 1650 /* Quadratic quad (I hate this) */ 1651 /* Edges are determined by the first 2 vertices (corners of edges) */ 1652 const PetscInt faceSizeQuad = 3; 1653 PetscInt sortedIndices[3], i, iFace; 1654 PetscBool found = PETSC_FALSE; 1655 PetscInt faceVerticesQuadSorted[12] = { 1656 0, 1, 4, /* bottom */ 1657 1, 2, 5, /* right */ 1658 2, 3, 6, /* top */ 1659 0, 3, 7, /* left */ 1660 }; 1661 PetscInt faceVerticesQuad[12] = { 1662 0, 1, 4, /* bottom */ 1663 1, 2, 5, /* right */ 1664 2, 3, 6, /* top */ 1665 3, 0, 7, /* left */ 1666 }; 1667 1668 faceSize = faceSizeQuad; 1669 for (i = 0; i < faceSizeQuad; ++i) sortedIndices[i] = indices[i]; 1670 ierr = PetscSortInt(faceSizeQuad, sortedIndices);CHKERRQ(ierr); 1671 for (iFace = 0; iFace < 4; ++iFace) { 1672 const PetscInt ii = iFace*faceSizeQuad; 1673 PetscInt fVertex, cVertex; 1674 1675 if ((sortedIndices[0] == faceVerticesQuadSorted[ii+0]) && 1676 (sortedIndices[1] == faceVerticesQuadSorted[ii+1])) { 1677 for (fVertex = 0; fVertex < faceSizeQuad; ++fVertex) { 1678 for (cVertex = 0; cVertex < faceSizeQuad; ++cVertex) { 1679 if (indices[cVertex] == faceVerticesQuad[ii+fVertex]) { 1680 faceVertices[fVertex] = origVertices[cVertex]; 1681 break; 1682 } 1683 } 1684 } 1685 found = PETSC_TRUE; 1686 break; 1687 } 1688 } 1689 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossface"); 1690 if (posOriented) *posOriented = PETSC_TRUE; 1691 PetscFunctionReturn(0); 1692 } else if (cellDim == 3 && numCorners == 8) { 1693 /* Hexes 1694 A hex is two oriented quads with the normal of the first 1695 pointing up at the second. 1696 1697 7---6 1698 /| /| 1699 4---5 | 1700 | 1-|-2 1701 |/ |/ 1702 0---3 1703 1704 Faces are determined by the first 4 vertices (corners of faces) */ 1705 const PetscInt faceSizeHex = 4; 1706 PetscInt sortedIndices[4], i, iFace; 1707 PetscBool found = PETSC_FALSE; 1708 PetscInt faceVerticesHexSorted[24] = { 1709 0, 1, 2, 3, /* bottom */ 1710 4, 5, 6, 7, /* top */ 1711 0, 3, 4, 5, /* front */ 1712 2, 3, 5, 6, /* right */ 1713 1, 2, 6, 7, /* back */ 1714 0, 1, 4, 7, /* left */ 1715 }; 1716 PetscInt faceVerticesHex[24] = { 1717 1, 2, 3, 0, /* bottom */ 1718 4, 5, 6, 7, /* top */ 1719 0, 3, 5, 4, /* front */ 1720 3, 2, 6, 5, /* right */ 1721 2, 1, 7, 6, /* back */ 1722 1, 0, 4, 7, /* left */ 1723 }; 1724 1725 faceSize = faceSizeHex; 1726 for (i = 0; i < faceSizeHex; ++i) sortedIndices[i] = indices[i]; 1727 ierr = PetscSortInt(faceSizeHex, sortedIndices);CHKERRQ(ierr); 1728 for (iFace = 0; iFace < 6; ++iFace) { 1729 const PetscInt ii = iFace*faceSizeHex; 1730 PetscInt fVertex, cVertex; 1731 1732 if ((sortedIndices[0] == faceVerticesHexSorted[ii+0]) && 1733 (sortedIndices[1] == faceVerticesHexSorted[ii+1]) && 1734 (sortedIndices[2] == faceVerticesHexSorted[ii+2]) && 1735 (sortedIndices[3] == faceVerticesHexSorted[ii+3])) { 1736 for (fVertex = 0; fVertex < faceSizeHex; ++fVertex) { 1737 for (cVertex = 0; cVertex < faceSizeHex; ++cVertex) { 1738 if (indices[cVertex] == faceVerticesHex[ii+fVertex]) { 1739 faceVertices[fVertex] = origVertices[cVertex]; 1740 break; 1741 } 1742 } 1743 } 1744 found = PETSC_TRUE; 1745 break; 1746 } 1747 } 1748 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface"); 1749 if (posOriented) *posOriented = PETSC_TRUE; 1750 PetscFunctionReturn(0); 1751 } else if (cellDim == 3 && numCorners == 10) { 1752 /* Quadratic tet */ 1753 /* Faces are determined by the first 3 vertices (corners of faces) */ 1754 const PetscInt faceSizeTet = 6; 1755 PetscInt sortedIndices[6], i, iFace; 1756 PetscBool found = PETSC_FALSE; 1757 PetscInt faceVerticesTetSorted[24] = { 1758 0, 1, 2, 6, 7, 8, /* bottom */ 1759 0, 3, 4, 6, 7, 9, /* front */ 1760 1, 4, 5, 7, 8, 9, /* right */ 1761 2, 3, 5, 6, 8, 9, /* left */ 1762 }; 1763 PetscInt faceVerticesTet[24] = { 1764 0, 1, 2, 6, 7, 8, /* bottom */ 1765 0, 4, 3, 6, 7, 9, /* front */ 1766 1, 5, 4, 7, 8, 9, /* right */ 1767 2, 3, 5, 8, 6, 9, /* left */ 1768 }; 1769 1770 faceSize = faceSizeTet; 1771 for (i = 0; i < faceSizeTet; ++i) sortedIndices[i] = indices[i]; 1772 ierr = PetscSortInt(faceSizeTet, sortedIndices);CHKERRQ(ierr); 1773 for (iFace=0; iFace < 4; ++iFace) { 1774 const PetscInt ii = iFace*faceSizeTet; 1775 PetscInt fVertex, cVertex; 1776 1777 if ((sortedIndices[0] == faceVerticesTetSorted[ii+0]) && 1778 (sortedIndices[1] == faceVerticesTetSorted[ii+1]) && 1779 (sortedIndices[2] == faceVerticesTetSorted[ii+2]) && 1780 (sortedIndices[3] == faceVerticesTetSorted[ii+3])) { 1781 for (fVertex = 0; fVertex < faceSizeTet; ++fVertex) { 1782 for (cVertex = 0; cVertex < faceSizeTet; ++cVertex) { 1783 if (indices[cVertex] == faceVerticesTet[ii+fVertex]) { 1784 faceVertices[fVertex] = origVertices[cVertex]; 1785 break; 1786 } 1787 } 1788 } 1789 found = PETSC_TRUE; 1790 break; 1791 } 1792 } 1793 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tet crossface"); 1794 if (posOriented) *posOriented = PETSC_TRUE; 1795 PetscFunctionReturn(0); 1796 } else if (cellDim == 3 && numCorners == 27) { 1797 /* Quadratic hexes (I hate this) 1798 A hex is two oriented quads with the normal of the first 1799 pointing up at the second. 1800 1801 7---6 1802 /| /| 1803 4---5 | 1804 | 3-|-2 1805 |/ |/ 1806 0---1 1807 1808 Faces are determined by the first 4 vertices (corners of faces) */ 1809 const PetscInt faceSizeQuadHex = 9; 1810 PetscInt sortedIndices[9], i, iFace; 1811 PetscBool found = PETSC_FALSE; 1812 PetscInt faceVerticesQuadHexSorted[54] = { 1813 0, 1, 2, 3, 8, 9, 10, 11, 24, /* bottom */ 1814 4, 5, 6, 7, 12, 13, 14, 15, 25, /* top */ 1815 0, 1, 4, 5, 8, 12, 16, 17, 22, /* front */ 1816 1, 2, 5, 6, 9, 13, 17, 18, 21, /* right */ 1817 2, 3, 6, 7, 10, 14, 18, 19, 23, /* back */ 1818 0, 3, 4, 7, 11, 15, 16, 19, 20, /* left */ 1819 }; 1820 PetscInt faceVerticesQuadHex[54] = { 1821 3, 2, 1, 0, 10, 9, 8, 11, 24, /* bottom */ 1822 4, 5, 6, 7, 12, 13, 14, 15, 25, /* top */ 1823 0, 1, 5, 4, 8, 17, 12, 16, 22, /* front */ 1824 1, 2, 6, 5, 9, 18, 13, 17, 21, /* right */ 1825 2, 3, 7, 6, 10, 19, 14, 18, 23, /* back */ 1826 3, 0, 4, 7, 11, 16, 15, 19, 20 /* left */ 1827 }; 1828 1829 faceSize = faceSizeQuadHex; 1830 for (i = 0; i < faceSizeQuadHex; ++i) sortedIndices[i] = indices[i]; 1831 ierr = PetscSortInt(faceSizeQuadHex, sortedIndices);CHKERRQ(ierr); 1832 for (iFace = 0; iFace < 6; ++iFace) { 1833 const PetscInt ii = iFace*faceSizeQuadHex; 1834 PetscInt fVertex, cVertex; 1835 1836 if ((sortedIndices[0] == faceVerticesQuadHexSorted[ii+0]) && 1837 (sortedIndices[1] == faceVerticesQuadHexSorted[ii+1]) && 1838 (sortedIndices[2] == faceVerticesQuadHexSorted[ii+2]) && 1839 (sortedIndices[3] == faceVerticesQuadHexSorted[ii+3])) { 1840 for (fVertex = 0; fVertex < faceSizeQuadHex; ++fVertex) { 1841 for (cVertex = 0; cVertex < faceSizeQuadHex; ++cVertex) { 1842 if (indices[cVertex] == faceVerticesQuadHex[ii+fVertex]) { 1843 faceVertices[fVertex] = origVertices[cVertex]; 1844 break; 1845 } 1846 } 1847 } 1848 found = PETSC_TRUE; 1849 break; 1850 } 1851 } 1852 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface"); 1853 if (posOriented) *posOriented = PETSC_TRUE; 1854 PetscFunctionReturn(0); 1855 } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Unknown cell type for faceOrientation()."); 1856 if (!posOrient) { 1857 if (debug) {ierr = PetscPrintf(comm, " Reversing initial face orientation\n");CHKERRQ(ierr);} 1858 for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[faceSize-1 - f]; 1859 } else { 1860 if (debug) {ierr = PetscPrintf(comm, " Keeping initial face orientation\n");CHKERRQ(ierr);} 1861 for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[f]; 1862 } 1863 if (posOriented) *posOriented = posOrient; 1864 PetscFunctionReturn(0); 1865 } 1866 1867 #undef __FUNCT__ 1868 #define __FUNCT__ "DMPlexGetOrientedFace" 1869 /* 1870 Given a cell and a face, as a set of vertices, 1871 return the oriented face, as a set of vertices, in faceVertices 1872 The orientation is such that the face normal points out of the cell 1873 */ 1874 PetscErrorCode DMPlexGetOrientedFace(DM dm, PetscInt cell, PetscInt faceSize, const PetscInt face[], PetscInt numCorners, PetscInt indices[], PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented) 1875 { 1876 const PetscInt *cone = NULL; 1877 PetscInt coneSize, v, f, v2; 1878 PetscInt oppositeVertex = -1; 1879 PetscErrorCode ierr; 1880 1881 PetscFunctionBegin; 1882 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 1883 ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 1884 for (v = 0, v2 = 0; v < coneSize; ++v) { 1885 PetscBool found = PETSC_FALSE; 1886 1887 for (f = 0; f < faceSize; ++f) { 1888 if (face[f] == cone[v]) { 1889 found = PETSC_TRUE; break; 1890 } 1891 } 1892 if (found) { 1893 indices[v2] = v; 1894 origVertices[v2] = cone[v]; 1895 ++v2; 1896 } else { 1897 oppositeVertex = v; 1898 } 1899 } 1900 ierr = DMPlexGetFaceOrientation(dm, cell, numCorners, indices, oppositeVertex, origVertices, faceVertices, posOriented);CHKERRQ(ierr); 1901 PetscFunctionReturn(0); 1902 } 1903 1904 #undef __FUNCT__ 1905 #define __FUNCT__ "DMPlexInsertFace_Internal" 1906 /* 1907 DMPlexInsertFace_Internal - Puts a face into the mesh 1908 1909 Not collective 1910 1911 Input Parameters: 1912 + dm - The DMPlex 1913 . numFaceVertex - The number of vertices in the face 1914 . faceVertices - The vertices in the face for dm 1915 . subfaceVertices - The vertices in the face for subdm 1916 . numCorners - The number of vertices in the cell 1917 . cell - A cell in dm containing the face 1918 . subcell - A cell in subdm containing the face 1919 . firstFace - First face in the mesh 1920 - newFacePoint - Next face in the mesh 1921 1922 Output Parameters: 1923 . newFacePoint - Contains next face point number on input, updated on output 1924 1925 Level: developer 1926 */ 1927 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) 1928 { 1929 MPI_Comm comm; 1930 DM_Plex *submesh = (DM_Plex*) subdm->data; 1931 const PetscInt *faces; 1932 PetscInt numFaces, coneSize; 1933 PetscErrorCode ierr; 1934 1935 PetscFunctionBegin; 1936 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1937 ierr = DMPlexGetConeSize(subdm, subcell, &coneSize);CHKERRQ(ierr); 1938 if (coneSize != 1) SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cone size of cell %d is %d != 1", cell, coneSize); 1939 #if 0 1940 /* Cannot use this because support() has not been constructed yet */ 1941 ierr = DMPlexGetJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr); 1942 #else 1943 { 1944 PetscInt f; 1945 1946 numFaces = 0; 1947 ierr = DMGetWorkArray(subdm, 1, PETSC_INT, (void **) &faces);CHKERRQ(ierr); 1948 for (f = firstFace; f < *newFacePoint; ++f) { 1949 PetscInt dof, off, d; 1950 1951 ierr = PetscSectionGetDof(submesh->coneSection, f, &dof);CHKERRQ(ierr); 1952 ierr = PetscSectionGetOffset(submesh->coneSection, f, &off);CHKERRQ(ierr); 1953 /* Yes, I know this is quadratic, but I expect the sizes to be <5 */ 1954 for (d = 0; d < dof; ++d) { 1955 const PetscInt p = submesh->cones[off+d]; 1956 PetscInt v; 1957 1958 for (v = 0; v < numFaceVertices; ++v) { 1959 if (subfaceVertices[v] == p) break; 1960 } 1961 if (v == numFaceVertices) break; 1962 } 1963 if (d == dof) { 1964 numFaces = 1; 1965 ((PetscInt*) faces)[0] = f; 1966 } 1967 } 1968 } 1969 #endif 1970 if (numFaces > 1) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Vertex set had %d faces, not one", numFaces); 1971 else if (numFaces == 1) { 1972 /* Add the other cell neighbor for this face */ 1973 ierr = DMPlexSetCone(subdm, subcell, faces);CHKERRQ(ierr); 1974 } else { 1975 PetscInt *indices, *origVertices, *orientedVertices, *orientedSubVertices, v, ov; 1976 PetscBool posOriented; 1977 1978 ierr = DMGetWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr); 1979 origVertices = &orientedVertices[numFaceVertices]; 1980 indices = &orientedVertices[numFaceVertices*2]; 1981 orientedSubVertices = &orientedVertices[numFaceVertices*3]; 1982 ierr = DMPlexGetOrientedFace(dm, cell, numFaceVertices, faceVertices, numCorners, indices, origVertices, orientedVertices, &posOriented);CHKERRQ(ierr); 1983 /* TODO: I know that routine should return a permutation, not the indices */ 1984 for (v = 0; v < numFaceVertices; ++v) { 1985 const PetscInt vertex = faceVertices[v], subvertex = subfaceVertices[v]; 1986 for (ov = 0; ov < numFaceVertices; ++ov) { 1987 if (orientedVertices[ov] == vertex) { 1988 orientedSubVertices[ov] = subvertex; 1989 break; 1990 } 1991 } 1992 if (ov == numFaceVertices) SETERRQ1(comm, PETSC_ERR_PLIB, "Could not find face vertex %d in orientated set", vertex); 1993 } 1994 ierr = DMPlexSetCone(subdm, *newFacePoint, orientedSubVertices);CHKERRQ(ierr); 1995 ierr = DMPlexSetCone(subdm, subcell, newFacePoint);CHKERRQ(ierr); 1996 ierr = DMRestoreWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr); 1997 ++(*newFacePoint); 1998 } 1999 #if 0 2000 ierr = DMPlexRestoreJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr); 2001 #else 2002 ierr = DMRestoreWorkArray(subdm, 1, PETSC_INT, (void **) &faces);CHKERRQ(ierr); 2003 #endif 2004 PetscFunctionReturn(0); 2005 } 2006 2007 #undef __FUNCT__ 2008 #define __FUNCT__ "DMPlexCreateSubmesh_Uninterpolated" 2009 static PetscErrorCode DMPlexCreateSubmesh_Uninterpolated(DM dm, DMLabel vertexLabel, PetscInt value, DM subdm) 2010 { 2011 MPI_Comm comm; 2012 DMLabel subpointMap; 2013 IS subvertexIS, subcellIS; 2014 const PetscInt *subVertices, *subCells; 2015 PetscInt numSubVertices, firstSubVertex, numSubCells; 2016 PetscInt *subface, maxConeSize, numSubFaces = 0, firstSubFace, newFacePoint, nFV = 0; 2017 PetscInt vStart, vEnd, c, f; 2018 PetscErrorCode ierr; 2019 2020 PetscFunctionBegin; 2021 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2022 /* Create subpointMap which marks the submesh */ 2023 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 2024 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 2025 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 2026 if (vertexLabel) {ierr = DMPlexMarkSubmesh_Uninterpolated(dm, vertexLabel, value, subpointMap, &numSubFaces, &nFV, subdm);CHKERRQ(ierr);} 2027 /* Setup chart */ 2028 ierr = DMLabelGetStratumSize(subpointMap, 0, &numSubVertices);CHKERRQ(ierr); 2029 ierr = DMLabelGetStratumSize(subpointMap, 2, &numSubCells);CHKERRQ(ierr); 2030 ierr = DMPlexSetChart(subdm, 0, numSubCells+numSubFaces+numSubVertices);CHKERRQ(ierr); 2031 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 2032 /* Set cone sizes */ 2033 firstSubVertex = numSubCells; 2034 firstSubFace = numSubCells+numSubVertices; 2035 newFacePoint = firstSubFace; 2036 ierr = DMLabelGetStratumIS(subpointMap, 0, &subvertexIS);CHKERRQ(ierr); 2037 if (subvertexIS) {ierr = ISGetIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 2038 ierr = DMLabelGetStratumIS(subpointMap, 2, &subcellIS);CHKERRQ(ierr); 2039 if (subcellIS) {ierr = ISGetIndices(subcellIS, &subCells);CHKERRQ(ierr);} 2040 for (c = 0; c < numSubCells; ++c) { 2041 ierr = DMPlexSetConeSize(subdm, c, 1);CHKERRQ(ierr); 2042 } 2043 for (f = firstSubFace; f < firstSubFace+numSubFaces; ++f) { 2044 ierr = DMPlexSetConeSize(subdm, f, nFV);CHKERRQ(ierr); 2045 } 2046 ierr = DMSetUp(subdm);CHKERRQ(ierr); 2047 /* Create face cones */ 2048 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 2049 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 2050 ierr = DMGetWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 2051 for (c = 0; c < numSubCells; ++c) { 2052 const PetscInt cell = subCells[c]; 2053 const PetscInt subcell = c; 2054 PetscInt *closure = NULL; 2055 PetscInt closureSize, cl, numCorners = 0, faceSize = 0; 2056 2057 ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 2058 for (cl = 0; cl < closureSize*2; cl += 2) { 2059 const PetscInt point = closure[cl]; 2060 PetscInt subVertex; 2061 2062 if ((point >= vStart) && (point < vEnd)) { 2063 ++numCorners; 2064 ierr = PetscFindInt(point, numSubVertices, subVertices, &subVertex);CHKERRQ(ierr); 2065 if (subVertex >= 0) { 2066 closure[faceSize] = point; 2067 subface[faceSize] = firstSubVertex+subVertex; 2068 ++faceSize; 2069 } 2070 } 2071 } 2072 if (faceSize > nFV) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize); 2073 if (faceSize == nFV) { 2074 ierr = DMPlexInsertFace_Internal(dm, subdm, faceSize, closure, subface, numCorners, cell, subcell, firstSubFace, &newFacePoint);CHKERRQ(ierr); 2075 } 2076 ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 2077 } 2078 ierr = DMRestoreWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 2079 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 2080 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 2081 /* Build coordinates */ 2082 { 2083 PetscSection coordSection, subCoordSection; 2084 Vec coordinates, subCoordinates; 2085 PetscScalar *coords, *subCoords; 2086 PetscInt numComp, coordSize, v; 2087 2088 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2089 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 2090 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 2091 ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr); 2092 ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr); 2093 ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr); 2094 ierr = PetscSectionSetChart(subCoordSection, firstSubVertex, firstSubVertex+numSubVertices);CHKERRQ(ierr); 2095 for (v = 0; v < numSubVertices; ++v) { 2096 const PetscInt vertex = subVertices[v]; 2097 const PetscInt subvertex = firstSubVertex+v; 2098 PetscInt dof; 2099 2100 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2101 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 2102 ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr); 2103 } 2104 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 2105 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 2106 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 2107 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2108 ierr = VecSetType(subCoordinates,VECSTANDARD);CHKERRQ(ierr); 2109 if (coordSize) { 2110 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2111 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2112 for (v = 0; v < numSubVertices; ++v) { 2113 const PetscInt vertex = subVertices[v]; 2114 const PetscInt subvertex = firstSubVertex+v; 2115 PetscInt dof, off, sdof, soff, d; 2116 2117 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2118 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 2119 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 2120 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 2121 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 2122 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 2123 } 2124 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2125 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2126 } 2127 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 2128 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 2129 } 2130 /* Cleanup */ 2131 if (subvertexIS) {ierr = ISRestoreIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 2132 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 2133 if (subcellIS) {ierr = ISRestoreIndices(subcellIS, &subCells);CHKERRQ(ierr);} 2134 ierr = ISDestroy(&subcellIS);CHKERRQ(ierr); 2135 PetscFunctionReturn(0); 2136 } 2137 2138 #undef __FUNCT__ 2139 #define __FUNCT__ "DMPlexCreateSubmesh_Interpolated" 2140 static PetscErrorCode DMPlexCreateSubmesh_Interpolated(DM dm, DMLabel vertexLabel, PetscInt value, DM subdm) 2141 { 2142 MPI_Comm comm; 2143 DMLabel subpointMap; 2144 IS *subpointIS; 2145 const PetscInt **subpoints; 2146 PetscInt *numSubPoints, *firstSubPoint, *coneNew, *coneONew; 2147 PetscInt totSubPoints = 0, maxConeSize, dim, p, d, v; 2148 PetscErrorCode ierr; 2149 2150 PetscFunctionBegin; 2151 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2152 /* Create subpointMap which marks the submesh */ 2153 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 2154 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 2155 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 2156 if (vertexLabel) {ierr = DMPlexMarkSubmesh_Interpolated(dm, vertexLabel, value, subpointMap, subdm);CHKERRQ(ierr);} 2157 /* Setup chart */ 2158 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 2159 ierr = PetscMalloc4(dim+1,PetscInt,&numSubPoints,dim+1,PetscInt,&firstSubPoint,dim+1,IS,&subpointIS,dim+1,const PetscInt *,&subpoints);CHKERRQ(ierr); 2160 for (d = 0; d <= dim; ++d) { 2161 ierr = DMLabelGetStratumSize(subpointMap, d, &numSubPoints[d]);CHKERRQ(ierr); 2162 totSubPoints += numSubPoints[d]; 2163 } 2164 ierr = DMPlexSetChart(subdm, 0, totSubPoints);CHKERRQ(ierr); 2165 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 2166 /* Set cone sizes */ 2167 firstSubPoint[dim] = 0; 2168 firstSubPoint[0] = firstSubPoint[dim] + numSubPoints[dim]; 2169 if (dim > 1) {firstSubPoint[dim-1] = firstSubPoint[0] + numSubPoints[0];} 2170 if (dim > 2) {firstSubPoint[dim-2] = firstSubPoint[dim-1] + numSubPoints[dim-1];} 2171 for (d = 0; d <= dim; ++d) { 2172 ierr = DMLabelGetStratumIS(subpointMap, d, &subpointIS[d]);CHKERRQ(ierr); 2173 ierr = ISGetIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 2174 } 2175 for (d = 0; d <= dim; ++d) { 2176 for (p = 0; p < numSubPoints[d]; ++p) { 2177 const PetscInt point = subpoints[d][p]; 2178 const PetscInt subpoint = firstSubPoint[d] + p; 2179 const PetscInt *cone; 2180 PetscInt coneSize, coneSizeNew, c, val; 2181 2182 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 2183 ierr = DMPlexSetConeSize(subdm, subpoint, coneSize);CHKERRQ(ierr); 2184 if (d == dim) { 2185 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 2186 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 2187 ierr = DMLabelGetValue(subpointMap, cone[c], &val);CHKERRQ(ierr); 2188 if (val >= 0) coneSizeNew++; 2189 } 2190 ierr = DMPlexSetConeSize(subdm, subpoint, coneSizeNew);CHKERRQ(ierr); 2191 } 2192 } 2193 } 2194 ierr = DMSetUp(subdm);CHKERRQ(ierr); 2195 /* Set cones */ 2196 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 2197 ierr = PetscMalloc2(maxConeSize,PetscInt,&coneNew,maxConeSize,PetscInt,&coneONew);CHKERRQ(ierr); 2198 for (d = 0; d <= dim; ++d) { 2199 for (p = 0; p < numSubPoints[d]; ++p) { 2200 const PetscInt point = subpoints[d][p]; 2201 const PetscInt subpoint = firstSubPoint[d] + p; 2202 const PetscInt *cone, *ornt; 2203 PetscInt coneSize, subconeSize, coneSizeNew, c, subc; 2204 2205 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 2206 ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr); 2207 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 2208 ierr = DMPlexGetConeOrientation(dm, point, &ornt);CHKERRQ(ierr); 2209 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 2210 ierr = PetscFindInt(cone[c], numSubPoints[d-1], subpoints[d-1], &subc);CHKERRQ(ierr); 2211 if (subc >= 0) { 2212 coneNew[coneSizeNew] = firstSubPoint[d-1] + subc; 2213 coneONew[coneSizeNew] = ornt[c]; 2214 ++coneSizeNew; 2215 } 2216 } 2217 if (coneSizeNew != subconeSize) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of cone points located %d does not match subcone size %d", coneSizeNew, subconeSize); 2218 ierr = DMPlexSetCone(subdm, subpoint, coneNew);CHKERRQ(ierr); 2219 ierr = DMPlexSetConeOrientation(subdm, subpoint, coneONew);CHKERRQ(ierr); 2220 } 2221 } 2222 ierr = PetscFree2(coneNew,coneONew);CHKERRQ(ierr); 2223 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 2224 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 2225 /* Build coordinates */ 2226 { 2227 PetscSection coordSection, subCoordSection; 2228 Vec coordinates, subCoordinates; 2229 PetscScalar *coords, *subCoords; 2230 PetscInt numComp, coordSize; 2231 2232 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2233 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 2234 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 2235 ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr); 2236 ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr); 2237 ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr); 2238 ierr = PetscSectionSetChart(subCoordSection, firstSubPoint[0], firstSubPoint[0]+numSubPoints[0]);CHKERRQ(ierr); 2239 for (v = 0; v < numSubPoints[0]; ++v) { 2240 const PetscInt vertex = subpoints[0][v]; 2241 const PetscInt subvertex = firstSubPoint[0]+v; 2242 PetscInt dof; 2243 2244 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2245 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 2246 ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr); 2247 } 2248 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 2249 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 2250 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 2251 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2252 ierr = VecSetType(subCoordinates,VECSTANDARD);CHKERRQ(ierr); 2253 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2254 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2255 for (v = 0; v < numSubPoints[0]; ++v) { 2256 const PetscInt vertex = subpoints[0][v]; 2257 const PetscInt subvertex = firstSubPoint[0]+v; 2258 PetscInt dof, off, sdof, soff, d; 2259 2260 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2261 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 2262 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 2263 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 2264 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 2265 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 2266 } 2267 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2268 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2269 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 2270 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 2271 } 2272 /* Cleanup */ 2273 for (d = 0; d <= dim; ++d) { 2274 ierr = ISRestoreIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 2275 ierr = ISDestroy(&subpointIS[d]);CHKERRQ(ierr); 2276 } 2277 ierr = PetscFree4(numSubPoints,firstSubPoint,subpointIS,subpoints);CHKERRQ(ierr); 2278 PetscFunctionReturn(0); 2279 } 2280 2281 #undef __FUNCT__ 2282 #define __FUNCT__ "DMPlexCreateSubmesh" 2283 /*@ 2284 DMPlexCreateSubmesh - Extract a hypersurface from the mesh using vertices defined by a label 2285 2286 Input Parameters: 2287 + dm - The original mesh 2288 . vertexLabel - The DMLabel marking vertices contained in the surface 2289 - value - The label value to use 2290 2291 Output Parameter: 2292 . subdm - The surface mesh 2293 2294 Note: This function produces a DMLabel mapping original points in the submesh to their depth. This can be obtained using DMPlexGetSubpointMap(). 2295 2296 Level: developer 2297 2298 .seealso: DMPlexGetSubpointMap(), DMPlexGetLabel(), DMLabelSetValue() 2299 @*/ 2300 PetscErrorCode DMPlexCreateSubmesh(DM dm, DMLabel vertexLabel, PetscInt value, DM *subdm) 2301 { 2302 PetscInt dim, depth; 2303 PetscErrorCode ierr; 2304 2305 PetscFunctionBegin; 2306 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2307 PetscValidPointer(subdm, 3); 2308 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 2309 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 2310 ierr = DMCreate(PetscObjectComm((PetscObject)dm), subdm);CHKERRQ(ierr); 2311 ierr = DMSetType(*subdm, DMPLEX);CHKERRQ(ierr); 2312 ierr = DMPlexSetDimension(*subdm, dim-1);CHKERRQ(ierr); 2313 if (depth == dim) { 2314 ierr = DMPlexCreateSubmesh_Interpolated(dm, vertexLabel, value, *subdm);CHKERRQ(ierr); 2315 } else { 2316 ierr = DMPlexCreateSubmesh_Uninterpolated(dm, vertexLabel, value, *subdm);CHKERRQ(ierr); 2317 } 2318 PetscFunctionReturn(0); 2319 } 2320 2321 #undef __FUNCT__ 2322 #define __FUNCT__ "DMPlexFilterPoint_Internal" 2323 PETSC_STATIC_INLINE PetscInt DMPlexFilterPoint_Internal(PetscInt point, PetscInt firstSubPoint, PetscInt numSubPoints, const PetscInt subPoints[]) 2324 { 2325 PetscInt subPoint; 2326 PetscErrorCode ierr; 2327 2328 ierr = PetscFindInt(point, numSubPoints, subPoints, &subPoint); if (ierr < 0) return ierr; 2329 return subPoint < 0 ? subPoint : firstSubPoint+subPoint; 2330 } 2331 2332 #undef __FUNCT__ 2333 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh_Uninterpolated" 2334 static PetscErrorCode DMPlexCreateCohesiveSubmesh_Uninterpolated(DM dm, PetscBool hasLagrange, const char label[], PetscInt value, DM subdm) 2335 { 2336 MPI_Comm comm; 2337 DMLabel subpointMap; 2338 IS subvertexIS; 2339 const PetscInt *subVertices; 2340 PetscInt numSubVertices, firstSubVertex, numSubCells, *subCells = NULL; 2341 PetscInt *subface, maxConeSize, numSubFaces, firstSubFace, newFacePoint, nFV; 2342 PetscInt cMax, c, f; 2343 PetscErrorCode ierr; 2344 2345 PetscFunctionBegin; 2346 ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 2347 /* Create subpointMap which marks the submesh */ 2348 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 2349 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 2350 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 2351 ierr = DMPlexMarkCohesiveSubmesh_Uninterpolated(dm, hasLagrange, label, value, subpointMap, &numSubFaces, &nFV, &subCells, subdm);CHKERRQ(ierr); 2352 /* Setup chart */ 2353 ierr = DMLabelGetStratumSize(subpointMap, 0, &numSubVertices);CHKERRQ(ierr); 2354 ierr = DMLabelGetStratumSize(subpointMap, 2, &numSubCells);CHKERRQ(ierr); 2355 ierr = DMPlexSetChart(subdm, 0, numSubCells+numSubFaces+numSubVertices);CHKERRQ(ierr); 2356 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 2357 /* Set cone sizes */ 2358 firstSubVertex = numSubCells; 2359 firstSubFace = numSubCells+numSubVertices; 2360 newFacePoint = firstSubFace; 2361 ierr = DMLabelGetStratumIS(subpointMap, 0, &subvertexIS);CHKERRQ(ierr); 2362 if (subvertexIS) {ierr = ISGetIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 2363 for (c = 0; c < numSubCells; ++c) { 2364 ierr = DMPlexSetConeSize(subdm, c, 1);CHKERRQ(ierr); 2365 } 2366 for (f = firstSubFace; f < firstSubFace+numSubFaces; ++f) { 2367 ierr = DMPlexSetConeSize(subdm, f, nFV);CHKERRQ(ierr); 2368 } 2369 ierr = DMSetUp(subdm);CHKERRQ(ierr); 2370 /* Create face cones */ 2371 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 2372 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 2373 ierr = DMGetWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 2374 for (c = 0; c < numSubCells; ++c) { 2375 const PetscInt cell = subCells[c]; 2376 const PetscInt subcell = c; 2377 const PetscInt *cone, *cells; 2378 PetscInt numCells, subVertex, p, v; 2379 2380 if (cell < cMax) continue; 2381 ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 2382 for (v = 0; v < nFV; ++v) { 2383 ierr = PetscFindInt(cone[v], numSubVertices, subVertices, &subVertex);CHKERRQ(ierr); 2384 subface[v] = firstSubVertex+subVertex; 2385 } 2386 ierr = DMPlexSetCone(subdm, newFacePoint, subface);CHKERRQ(ierr); 2387 ierr = DMPlexSetCone(subdm, subcell, &newFacePoint);CHKERRQ(ierr); 2388 ierr = DMPlexGetJoin(dm, nFV, cone, &numCells, &cells);CHKERRQ(ierr); 2389 /* Not true in parallel 2390 if (numCells != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells"); */ 2391 for (p = 0; p < numCells; ++p) { 2392 PetscInt negsubcell; 2393 2394 if (cells[p] >= cMax) continue; 2395 /* I know this is a crap search */ 2396 for (negsubcell = 0; negsubcell < numSubCells; ++negsubcell) { 2397 if (subCells[negsubcell] == cells[p]) break; 2398 } 2399 if (negsubcell == numSubCells) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find negative face neighbor for cohesive cell %d", cell); 2400 ierr = DMPlexSetCone(subdm, negsubcell, &newFacePoint);CHKERRQ(ierr); 2401 } 2402 ierr = DMPlexRestoreJoin(dm, nFV, cone, &numCells, &cells);CHKERRQ(ierr); 2403 ++newFacePoint; 2404 } 2405 ierr = DMRestoreWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 2406 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 2407 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 2408 /* Build coordinates */ 2409 { 2410 PetscSection coordSection, subCoordSection; 2411 Vec coordinates, subCoordinates; 2412 PetscScalar *coords, *subCoords; 2413 PetscInt numComp, coordSize, v; 2414 2415 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2416 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 2417 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 2418 ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr); 2419 ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr); 2420 ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr); 2421 ierr = PetscSectionSetChart(subCoordSection, firstSubVertex, firstSubVertex+numSubVertices);CHKERRQ(ierr); 2422 for (v = 0; v < numSubVertices; ++v) { 2423 const PetscInt vertex = subVertices[v]; 2424 const PetscInt subvertex = firstSubVertex+v; 2425 PetscInt dof; 2426 2427 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2428 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 2429 ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr); 2430 } 2431 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 2432 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 2433 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 2434 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2435 ierr = VecSetType(subCoordinates,VECSTANDARD);CHKERRQ(ierr); 2436 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2437 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2438 for (v = 0; v < numSubVertices; ++v) { 2439 const PetscInt vertex = subVertices[v]; 2440 const PetscInt subvertex = firstSubVertex+v; 2441 PetscInt dof, off, sdof, soff, d; 2442 2443 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2444 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 2445 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 2446 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 2447 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 2448 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 2449 } 2450 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2451 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2452 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 2453 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 2454 } 2455 /* Build SF */ 2456 CHKMEMQ; 2457 { 2458 PetscSF sfPoint, sfPointSub; 2459 const PetscSFNode *remotePoints; 2460 PetscSFNode *sremotePoints, *newLocalPoints, *newOwners; 2461 const PetscInt *localPoints; 2462 PetscInt *slocalPoints; 2463 PetscInt numRoots, numLeaves, numSubRoots = numSubCells+numSubFaces+numSubVertices, numSubLeaves = 0, l, sl, ll, pStart, pEnd, p, vStart, vEnd; 2464 PetscMPIInt rank; 2465 2466 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 2467 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 2468 ierr = DMGetPointSF(subdm, &sfPointSub);CHKERRQ(ierr); 2469 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 2470 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 2471 ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 2472 if (numRoots >= 0) { 2473 /* Only vertices should be shared */ 2474 ierr = PetscMalloc2(pEnd-pStart,PetscSFNode,&newLocalPoints,numRoots,PetscSFNode,&newOwners);CHKERRQ(ierr); 2475 for (p = 0; p < pEnd-pStart; ++p) { 2476 newLocalPoints[p].rank = -2; 2477 newLocalPoints[p].index = -2; 2478 } 2479 /* Set subleaves */ 2480 for (l = 0; l < numLeaves; ++l) { 2481 const PetscInt point = localPoints[l]; 2482 const PetscInt subPoint = DMPlexFilterPoint_Internal(point, firstSubVertex, numSubVertices, subVertices); 2483 2484 if ((point < vStart) && (point >= vEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Should not be mapping anything but vertices, %d", point); 2485 if (subPoint < 0) continue; 2486 newLocalPoints[point-pStart].rank = rank; 2487 newLocalPoints[point-pStart].index = subPoint; 2488 ++numSubLeaves; 2489 } 2490 /* Must put in owned subpoints */ 2491 for (p = pStart; p < pEnd; ++p) { 2492 const PetscInt subPoint = DMPlexFilterPoint_Internal(p, firstSubVertex, numSubVertices, subVertices); 2493 2494 if (subPoint < 0) { 2495 newOwners[p-pStart].rank = -3; 2496 newOwners[p-pStart].index = -3; 2497 } else { 2498 newOwners[p-pStart].rank = rank; 2499 newOwners[p-pStart].index = subPoint; 2500 } 2501 } 2502 ierr = PetscSFReduceBegin(sfPoint, MPIU_2INT, newLocalPoints, newOwners, MPI_MAXLOC);CHKERRQ(ierr); 2503 ierr = PetscSFReduceEnd(sfPoint, MPIU_2INT, newLocalPoints, newOwners, MPI_MAXLOC);CHKERRQ(ierr); 2504 ierr = PetscSFBcastBegin(sfPoint, MPIU_2INT, newOwners, newLocalPoints);CHKERRQ(ierr); 2505 ierr = PetscSFBcastEnd(sfPoint, MPIU_2INT, newOwners, newLocalPoints);CHKERRQ(ierr); 2506 ierr = PetscMalloc(numSubLeaves * sizeof(PetscInt), &slocalPoints);CHKERRQ(ierr); 2507 ierr = PetscMalloc(numSubLeaves * sizeof(PetscSFNode), &sremotePoints);CHKERRQ(ierr); 2508 for (l = 0, sl = 0, ll = 0; l < numLeaves; ++l) { 2509 const PetscInt point = localPoints[l]; 2510 const PetscInt subPoint = DMPlexFilterPoint_Internal(point, firstSubVertex, numSubVertices, subVertices); 2511 2512 if (subPoint < 0) continue; 2513 if (newLocalPoints[point].rank == rank) {++ll; continue;} 2514 slocalPoints[sl] = subPoint; 2515 sremotePoints[sl].rank = newLocalPoints[point].rank; 2516 sremotePoints[sl].index = newLocalPoints[point].index; 2517 if (sremotePoints[sl].rank < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid remote rank for local point %d", point); 2518 if (sremotePoints[sl].index < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid remote subpoint for local point %d", point); 2519 ++sl; 2520 } 2521 ierr = PetscFree2(newLocalPoints,newOwners);CHKERRQ(ierr); 2522 if (sl + ll != numSubLeaves) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Mismatch in number of subleaves %d + %d != %d", sl, ll, numSubLeaves); 2523 ierr = PetscSFSetGraph(sfPointSub, numSubRoots, sl, slocalPoints, PETSC_OWN_POINTER, sremotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 2524 } 2525 } 2526 CHKMEMQ; 2527 /* Cleanup */ 2528 if (subvertexIS) {ierr = ISRestoreIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 2529 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 2530 ierr = PetscFree(subCells);CHKERRQ(ierr); 2531 PetscFunctionReturn(0); 2532 } 2533 2534 #undef __FUNCT__ 2535 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh_Interpolated" 2536 static PetscErrorCode DMPlexCreateCohesiveSubmesh_Interpolated(DM dm, PetscBool hasLagrange, const char label[], PetscInt value, DM subdm) 2537 { 2538 MPI_Comm comm; 2539 DMLabel subpointMap; 2540 IS *subpointIS; 2541 const PetscInt **subpoints; 2542 PetscInt *numSubPoints, *firstSubPoint, *coneNew; 2543 PetscInt totSubPoints = 0, maxConeSize, dim, p, d, v; 2544 PetscErrorCode ierr; 2545 2546 PetscFunctionBegin; 2547 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2548 /* Create subpointMap which marks the submesh */ 2549 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 2550 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 2551 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 2552 ierr = DMPlexMarkCohesiveSubmesh_Interpolated(dm, hasLagrange, label, value, subpointMap, subdm);CHKERRQ(ierr); 2553 /* Setup chart */ 2554 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 2555 ierr = PetscMalloc4(dim+1,PetscInt,&numSubPoints,dim+1,PetscInt,&firstSubPoint,dim+1,IS,&subpointIS,dim+1,const PetscInt *,&subpoints);CHKERRQ(ierr); 2556 for (d = 0; d <= dim; ++d) { 2557 ierr = DMLabelGetStratumSize(subpointMap, d, &numSubPoints[d]);CHKERRQ(ierr); 2558 totSubPoints += numSubPoints[d]; 2559 } 2560 ierr = DMPlexSetChart(subdm, 0, totSubPoints);CHKERRQ(ierr); 2561 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 2562 /* Set cone sizes */ 2563 firstSubPoint[dim] = 0; 2564 firstSubPoint[0] = firstSubPoint[dim] + numSubPoints[dim]; 2565 if (dim > 1) {firstSubPoint[dim-1] = firstSubPoint[0] + numSubPoints[0];} 2566 if (dim > 2) {firstSubPoint[dim-2] = firstSubPoint[dim-1] + numSubPoints[dim-1];} 2567 for (d = 0; d <= dim; ++d) { 2568 ierr = DMLabelGetStratumIS(subpointMap, d, &subpointIS[d]);CHKERRQ(ierr); 2569 ierr = ISGetIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 2570 } 2571 for (d = 0; d <= dim; ++d) { 2572 for (p = 0; p < numSubPoints[d]; ++p) { 2573 const PetscInt point = subpoints[d][p]; 2574 const PetscInt subpoint = firstSubPoint[d] + p; 2575 const PetscInt *cone; 2576 PetscInt coneSize, coneSizeNew, c, val; 2577 2578 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 2579 ierr = DMPlexSetConeSize(subdm, subpoint, coneSize);CHKERRQ(ierr); 2580 if (d == dim) { 2581 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 2582 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 2583 ierr = DMLabelGetValue(subpointMap, cone[c], &val);CHKERRQ(ierr); 2584 if (val >= 0) coneSizeNew++; 2585 } 2586 ierr = DMPlexSetConeSize(subdm, subpoint, coneSizeNew);CHKERRQ(ierr); 2587 } 2588 } 2589 } 2590 ierr = DMSetUp(subdm);CHKERRQ(ierr); 2591 /* Set cones */ 2592 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 2593 ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &coneNew);CHKERRQ(ierr); 2594 for (d = 0; d <= dim; ++d) { 2595 for (p = 0; p < numSubPoints[d]; ++p) { 2596 const PetscInt point = subpoints[d][p]; 2597 const PetscInt subpoint = firstSubPoint[d] + p; 2598 const PetscInt *cone; 2599 PetscInt coneSize, subconeSize, coneSizeNew, c, subc; 2600 2601 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 2602 ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr); 2603 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 2604 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 2605 ierr = PetscFindInt(cone[c], numSubPoints[d-1], subpoints[d-1], &subc);CHKERRQ(ierr); 2606 if (subc >= 0) coneNew[coneSizeNew++] = firstSubPoint[d-1] + subc; 2607 } 2608 if (coneSizeNew != subconeSize) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of cone points located %d does not match subcone size %d", coneSizeNew, subconeSize); 2609 ierr = DMPlexSetCone(subdm, subpoint, coneNew);CHKERRQ(ierr); 2610 } 2611 } 2612 ierr = PetscFree(coneNew);CHKERRQ(ierr); 2613 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 2614 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 2615 /* Build coordinates */ 2616 { 2617 PetscSection coordSection, subCoordSection; 2618 Vec coordinates, subCoordinates; 2619 PetscScalar *coords, *subCoords; 2620 PetscInt numComp, coordSize; 2621 2622 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2623 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 2624 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 2625 ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr); 2626 ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr); 2627 ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr); 2628 ierr = PetscSectionSetChart(subCoordSection, firstSubPoint[0], firstSubPoint[0]+numSubPoints[0]);CHKERRQ(ierr); 2629 for (v = 0; v < numSubPoints[0]; ++v) { 2630 const PetscInt vertex = subpoints[0][v]; 2631 const PetscInt subvertex = firstSubPoint[0]+v; 2632 PetscInt dof; 2633 2634 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2635 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 2636 ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr); 2637 } 2638 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 2639 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 2640 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 2641 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2642 ierr = VecSetType(subCoordinates,VECSTANDARD);CHKERRQ(ierr); 2643 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2644 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2645 for (v = 0; v < numSubPoints[0]; ++v) { 2646 const PetscInt vertex = subpoints[0][v]; 2647 const PetscInt subvertex = firstSubPoint[0]+v; 2648 PetscInt dof, off, sdof, soff, d; 2649 2650 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2651 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 2652 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 2653 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 2654 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 2655 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 2656 } 2657 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2658 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2659 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 2660 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 2661 } 2662 /* Cleanup */ 2663 for (d = 0; d <= dim; ++d) { 2664 ierr = ISRestoreIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 2665 ierr = ISDestroy(&subpointIS[d]);CHKERRQ(ierr); 2666 } 2667 ierr = PetscFree4(numSubPoints,firstSubPoint,subpointIS,subpoints);CHKERRQ(ierr); 2668 PetscFunctionReturn(0); 2669 } 2670 2671 #undef __FUNCT__ 2672 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh" 2673 /* 2674 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. 2675 2676 Input Parameters: 2677 + dm - The original mesh 2678 . hasLagrange - The mesh has Lagrange unknowns in the cohesive cells 2679 . label - A label name, or PETSC_NULL 2680 - value - A label value 2681 2682 Output Parameter: 2683 . subdm - The surface mesh 2684 2685 Note: This function produces a DMLabel mapping original points in the submesh to their depth. This can be obtained using DMPlexGetSubpointMap(). 2686 2687 Level: developer 2688 2689 .seealso: DMPlexGetSubpointMap(), DMPlexCreateSubmesh() 2690 */ 2691 PetscErrorCode DMPlexCreateCohesiveSubmesh(DM dm, PetscBool hasLagrange, const char label[], PetscInt value, DM *subdm) 2692 { 2693 PetscInt dim, depth; 2694 PetscErrorCode ierr; 2695 2696 PetscFunctionBegin; 2697 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2698 PetscValidPointer(subdm, 5); 2699 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 2700 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 2701 ierr = DMCreate(PetscObjectComm((PetscObject)dm), subdm);CHKERRQ(ierr); 2702 ierr = DMSetType(*subdm, DMPLEX);CHKERRQ(ierr); 2703 ierr = DMPlexSetDimension(*subdm, dim-1);CHKERRQ(ierr); 2704 if (depth == dim) { 2705 ierr = DMPlexCreateCohesiveSubmesh_Interpolated(dm, hasLagrange, label, value, *subdm);CHKERRQ(ierr); 2706 } else { 2707 ierr = DMPlexCreateCohesiveSubmesh_Uninterpolated(dm, hasLagrange, label, value, *subdm);CHKERRQ(ierr); 2708 } 2709 PetscFunctionReturn(0); 2710 } 2711 2712 #undef __FUNCT__ 2713 #define __FUNCT__ "DMPlexGetSubpointMap" 2714 /*@ 2715 DMPlexGetSubpointMap - Returns a DMLabel with point dimension as values 2716 2717 Input Parameter: 2718 . dm - The submesh DM 2719 2720 Output Parameter: 2721 . subpointMap - The DMLabel of all the points from the original mesh in this submesh, or NULL if this is not a submesh 2722 2723 Level: developer 2724 2725 .seealso: DMPlexCreateSubmesh(), DMPlexCreateSubpointIS() 2726 @*/ 2727 PetscErrorCode DMPlexGetSubpointMap(DM dm, DMLabel *subpointMap) 2728 { 2729 DM_Plex *mesh = (DM_Plex*) dm->data; 2730 2731 PetscFunctionBegin; 2732 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2733 PetscValidPointer(subpointMap, 2); 2734 *subpointMap = mesh->subpointMap; 2735 PetscFunctionReturn(0); 2736 } 2737 2738 #undef __FUNCT__ 2739 #define __FUNCT__ "DMPlexSetSubpointMap" 2740 /* Note: Should normally not be called by the user, since it is set in DMPlexCreateSubmesh() */ 2741 PetscErrorCode DMPlexSetSubpointMap(DM dm, DMLabel subpointMap) 2742 { 2743 DM_Plex *mesh = (DM_Plex *) dm->data; 2744 DMLabel tmp; 2745 PetscErrorCode ierr; 2746 2747 PetscFunctionBegin; 2748 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2749 tmp = mesh->subpointMap; 2750 mesh->subpointMap = subpointMap; 2751 ++mesh->subpointMap->refct; 2752 ierr = DMLabelDestroy(&tmp);CHKERRQ(ierr); 2753 PetscFunctionReturn(0); 2754 } 2755 2756 #undef __FUNCT__ 2757 #define __FUNCT__ "DMPlexCreateSubpointIS" 2758 /*@ 2759 DMPlexCreateSubpointIS - Creates an IS covering the entire subdm chart with the original points as data 2760 2761 Input Parameter: 2762 . dm - The submesh DM 2763 2764 Output Parameter: 2765 . subpointIS - The IS of all the points from the original mesh in this submesh, or NULL if this is not a submesh 2766 2767 Note: This is IS is guaranteed to be sorted by the construction of the submesh 2768 2769 Level: developer 2770 2771 .seealso: DMPlexCreateSubmesh(), DMPlexGetSubpointMap() 2772 @*/ 2773 PetscErrorCode DMPlexCreateSubpointIS(DM dm, IS *subpointIS) 2774 { 2775 MPI_Comm comm; 2776 DMLabel subpointMap; 2777 IS is; 2778 const PetscInt *opoints; 2779 PetscInt *points, *depths; 2780 PetscInt depth, depStart, depEnd, d, pStart, pEnd, p, n, off; 2781 PetscErrorCode ierr; 2782 2783 PetscFunctionBegin; 2784 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2785 PetscValidPointer(subpointIS, 2); 2786 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2787 *subpointIS = NULL; 2788 ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 2789 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 2790 if (subpointMap && depth >= 0) { 2791 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 2792 if (pStart) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Submeshes must start the point numbering at 0, not %d", pStart); 2793 ierr = DMGetWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr); 2794 depths[0] = depth; 2795 depths[1] = 0; 2796 for(d = 2; d <= depth; ++d) {depths[d] = depth+1 - d;} 2797 ierr = PetscMalloc(pEnd * sizeof(PetscInt), &points);CHKERRQ(ierr); 2798 for(d = 0, off = 0; d <= depth; ++d) { 2799 const PetscInt dep = depths[d]; 2800 2801 ierr = DMPlexGetDepthStratum(dm, dep, &depStart, &depEnd);CHKERRQ(ierr); 2802 ierr = DMLabelGetStratumSize(subpointMap, dep, &n);CHKERRQ(ierr); 2803 if (((d < 2) && (depth > 1)) || (d == 1)) { /* Only check vertices and cells for now since the map is broken for others */ 2804 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); 2805 } else { 2806 if (!n) { 2807 if (d == 0) { 2808 /* Missing cells */ 2809 for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = -1; 2810 } else { 2811 /* Missing faces */ 2812 for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = PETSC_MAX_INT; 2813 } 2814 } 2815 } 2816 if (n) { 2817 ierr = DMLabelGetStratumIS(subpointMap, dep, &is);CHKERRQ(ierr); 2818 ierr = ISGetIndices(is, &opoints);CHKERRQ(ierr); 2819 for(p = 0; p < n; ++p, ++off) points[off] = opoints[p]; 2820 ierr = ISRestoreIndices(is, &opoints);CHKERRQ(ierr); 2821 ierr = ISDestroy(&is);CHKERRQ(ierr); 2822 } 2823 } 2824 ierr = DMRestoreWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr); 2825 if (off != pEnd) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "The number of mapped submesh points %d should be %d", off, pEnd); 2826 ierr = ISCreateGeneral(PETSC_COMM_SELF, pEnd, points, PETSC_OWN_POINTER, subpointIS);CHKERRQ(ierr); 2827 } 2828 PetscFunctionReturn(0); 2829 } 2830