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