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