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