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