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 - subdm - The subDM associated with the label, or NULL 943 944 Output Parameter: 945 . label - A DMLabel marking all surface points 946 947 Level: developer 948 949 .seealso: DMPlexConstructCohesiveCells(), DMPlexLabelComplete() 950 @*/ 951 PetscErrorCode DMPlexLabelCohesiveComplete(DM dm, DMLabel label, DM subdm) 952 { 953 IS dimIS, subpointIS; 954 const PetscInt *points, *subpoints; 955 PetscInt shift = 100, dim, dep, cStart, cEnd, numPoints, numSubpoints, p, val; 956 PetscErrorCode ierr; 957 958 PetscFunctionBegin; 959 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 960 if (subdm) { 961 ierr = DMPlexCreateSubpointIS(subdm, &subpointIS);CHKERRQ(ierr); 962 if (subpointIS) { 963 ierr = ISGetLocalSize(subpointIS, &numSubpoints);CHKERRQ(ierr); 964 ierr = ISGetIndices(subpointIS, &subpoints);CHKERRQ(ierr); 965 } 966 } 967 /* Cell orientation for face gives the side of the fault */ 968 ierr = DMLabelGetStratumIS(label, dim-1, &dimIS);CHKERRQ(ierr); 969 if (!dimIS) PetscFunctionReturn(0); 970 ierr = ISGetLocalSize(dimIS, &numPoints);CHKERRQ(ierr); 971 ierr = ISGetIndices(dimIS, &points);CHKERRQ(ierr); 972 for (p = 0; p < numPoints; ++p) { 973 const PetscInt *support; 974 PetscInt supportSize, s; 975 976 ierr = DMPlexGetSupportSize(dm, points[p], &supportSize);CHKERRQ(ierr); 977 if (supportSize != 2) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Split face %d has %d != 2 supports", points[p], supportSize); 978 ierr = DMPlexGetSupport(dm, points[p], &support);CHKERRQ(ierr); 979 for (s = 0; s < supportSize; ++s) { 980 const PetscInt *cone, *ornt; 981 PetscInt coneSize, c; 982 PetscBool pos = PETSC_TRUE; 983 984 ierr = DMPlexGetConeSize(dm, support[s], &coneSize);CHKERRQ(ierr); 985 ierr = DMPlexGetCone(dm, support[s], &cone);CHKERRQ(ierr); 986 ierr = DMPlexGetConeOrientation(dm, support[s], &ornt);CHKERRQ(ierr); 987 for (c = 0; c < coneSize; ++c) { 988 if (cone[c] == points[p]) { 989 PetscInt o = ornt[c]; 990 991 if (subdm) { 992 const PetscInt *subcone, *subornt; 993 PetscInt subpoint, subface, subconeSize, sc; 994 995 ierr = PetscFindInt(support[s], numSubpoints, subpoints, &subpoint);CHKERRQ(ierr); 996 ierr = PetscFindInt(points[p], numSubpoints, subpoints, &subface);CHKERRQ(ierr); 997 ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr); 998 ierr = DMPlexGetCone(subdm, subpoint, &subcone);CHKERRQ(ierr); 999 ierr = DMPlexGetConeOrientation(subdm, subpoint, &subornt);CHKERRQ(ierr); 1000 for (sc = 0; sc < subconeSize; ++sc) { 1001 if (subcone[sc] == subface) { 1002 o = subornt[0]; 1003 break; 1004 } 1005 } 1006 if (sc >= subconeSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find point %d in cone for subpoint %d", points[p], subpoint); 1007 } 1008 #if 1 1009 if (o >= 0) { 1010 #else 1011 if (o < 0) { 1012 #endif 1013 ierr = DMLabelSetValue(label, support[s], shift+dim);CHKERRQ(ierr); 1014 } else { 1015 ierr = DMLabelSetValue(label, support[s], -(shift+dim));CHKERRQ(ierr); 1016 pos = PETSC_FALSE; 1017 } 1018 break; 1019 } 1020 } 1021 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]); 1022 /* Put faces touching the fault in the label */ 1023 for (c = 0; c < coneSize; ++c) { 1024 const PetscInt point = cone[c]; 1025 1026 ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr); 1027 if (val == -1) { 1028 PetscInt *closure = NULL; 1029 PetscInt closureSize, cl; 1030 1031 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1032 for (cl = 0; cl < closureSize*2; cl += 2) { 1033 const PetscInt clp = closure[cl]; 1034 1035 ierr = DMLabelGetValue(label, clp, &val);CHKERRQ(ierr); 1036 if ((val >= 0) && (val < dim-1)) { 1037 ierr = DMLabelSetValue(label, point, pos == PETSC_TRUE ? shift+dim-1 : -(shift+dim-1));CHKERRQ(ierr); 1038 break; 1039 } 1040 } 1041 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1042 } 1043 } 1044 } 1045 } 1046 if (subdm) { 1047 if (subpointIS) {ierr = ISRestoreIndices(subpointIS, &subpoints);CHKERRQ(ierr);} 1048 ierr = ISDestroy(&subpointIS);CHKERRQ(ierr); 1049 } 1050 ierr = ISRestoreIndices(dimIS, &points);CHKERRQ(ierr); 1051 ierr = ISDestroy(&dimIS);CHKERRQ(ierr); 1052 /* Search for other cells/faces/edges connected to the fault by a vertex */ 1053 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1054 ierr = DMLabelGetStratumIS(label, 0, &dimIS);CHKERRQ(ierr); 1055 if (!dimIS) PetscFunctionReturn(0); 1056 ierr = ISGetLocalSize(dimIS, &numPoints);CHKERRQ(ierr); 1057 ierr = ISGetIndices(dimIS, &points);CHKERRQ(ierr); 1058 for (p = 0; p < numPoints; ++p) { 1059 PetscInt *star = NULL; 1060 PetscInt starSize, s; 1061 PetscInt again = 1; /* 0: Finished 1: Keep iterating after a change 2: No change */ 1062 1063 /* First mark cells connected to the fault */ 1064 ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1065 while (again) { 1066 if (again > 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Could not classify all cells connected to the fault"); 1067 again = 0; 1068 for (s = 0; s < starSize*2; s += 2) { 1069 const PetscInt point = star[s]; 1070 const PetscInt *cone; 1071 PetscInt coneSize, c; 1072 1073 if ((point < cStart) || (point >= cEnd)) continue; 1074 ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr); 1075 if (val != -1) continue; 1076 again = 2; 1077 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 1078 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 1079 for (c = 0; c < coneSize; ++c) { 1080 ierr = DMLabelGetValue(label, cone[c], &val);CHKERRQ(ierr); 1081 if (val != -1) { 1082 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); 1083 if (val > 0) { 1084 ierr = DMLabelSetValue(label, point, shift+dim);CHKERRQ(ierr); 1085 } else { 1086 ierr = DMLabelSetValue(label, point, -(shift+dim));CHKERRQ(ierr); 1087 } 1088 again = 1; 1089 break; 1090 } 1091 } 1092 } 1093 } 1094 /* Classify the rest by cell membership */ 1095 for (s = 0; s < starSize*2; s += 2) { 1096 const PetscInt point = star[s]; 1097 1098 ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr); 1099 if (val == -1) { 1100 PetscInt *sstar = NULL; 1101 PetscInt sstarSize, ss; 1102 PetscBool marked = PETSC_FALSE; 1103 1104 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &sstarSize, &sstar);CHKERRQ(ierr); 1105 for (ss = 0; ss < sstarSize*2; ss += 2) { 1106 const PetscInt spoint = sstar[ss]; 1107 1108 if ((spoint < cStart) || (spoint >= cEnd)) continue; 1109 ierr = DMLabelGetValue(label, spoint, &val);CHKERRQ(ierr); 1110 if (val == -1) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Cell %d in star of %d does not have a valid label", spoint, point); 1111 ierr = DMPlexGetLabelValue(dm, "depth", point, &dep);CHKERRQ(ierr); 1112 if (val > 0) { 1113 ierr = DMLabelSetValue(label, point, shift+dep);CHKERRQ(ierr); 1114 } else { 1115 ierr = DMLabelSetValue(label, point, -(shift+dep));CHKERRQ(ierr); 1116 } 1117 marked = PETSC_TRUE; 1118 break; 1119 } 1120 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &sstarSize, &sstar);CHKERRQ(ierr); 1121 if (!marked) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Point %d could not be classified", point); 1122 } 1123 } 1124 ierr = DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1125 } 1126 ierr = ISRestoreIndices(dimIS, &points);CHKERRQ(ierr); 1127 ierr = ISDestroy(&dimIS);CHKERRQ(ierr); 1128 PetscFunctionReturn(0); 1129 } 1130 1131 #undef __FUNCT__ 1132 #define __FUNCT__ "DMPlexMarkSubmesh_Uninterpolated" 1133 /* Here we need the explicit assumption that: 1134 1135 For any marked cell, the marked vertices constitute a single face 1136 */ 1137 static PetscErrorCode DMPlexMarkSubmesh_Uninterpolated(DM dm, DMLabel vertexLabel, PetscInt value, DMLabel subpointMap, PetscInt *numFaces, PetscInt *nFV, DM subdm) 1138 { 1139 IS subvertexIS; 1140 const PetscInt *subvertices; 1141 PetscInt *pStart, *pEnd, *pMax, pSize; 1142 PetscInt depth, dim, d, numSubVerticesInitial = 0, v; 1143 PetscErrorCode ierr; 1144 1145 PetscFunctionBegin; 1146 *numFaces = 0; 1147 *nFV = 0; 1148 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1149 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1150 pSize = PetscMax(depth, dim) + 1; 1151 ierr = PetscMalloc3(pSize,PetscInt,&pStart,pSize,PetscInt,&pEnd,pSize,PetscInt,&pMax);CHKERRQ(ierr); 1152 ierr = DMPlexGetHybridBounds(dm, &pMax[depth], depth>1 ? &pMax[depth-1] : NULL, depth > 2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr); 1153 for (d = 0; d <= depth; ++d) { 1154 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1155 if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]); 1156 } 1157 /* Loop over initial vertices and mark all faces in the collective star() */ 1158 ierr = DMLabelGetStratumIS(vertexLabel, value, &subvertexIS);CHKERRQ(ierr); 1159 if (subvertexIS) { 1160 ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr); 1161 ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1162 } 1163 for (v = 0; v < numSubVerticesInitial; ++v) { 1164 const PetscInt vertex = subvertices[v]; 1165 PetscInt *star = NULL; 1166 PetscInt starSize, s, numCells = 0, c; 1167 1168 ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1169 for (s = 0; s < starSize*2; s += 2) { 1170 const PetscInt point = star[s]; 1171 if ((point >= pStart[depth]) && (point < pEnd[depth])) star[numCells++] = point; 1172 } 1173 for (c = 0; c < numCells; ++c) { 1174 const PetscInt cell = star[c]; 1175 PetscInt *closure = NULL; 1176 PetscInt closureSize, cl; 1177 PetscInt cellLoc, numCorners = 0, faceSize = 0; 1178 1179 ierr = DMLabelGetValue(subpointMap, cell, &cellLoc);CHKERRQ(ierr); 1180 if (cellLoc == 2) continue; 1181 if (cellLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Cell %d has dimension %d in the surface label", cell, cellLoc); 1182 ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1183 for (cl = 0; cl < closureSize*2; cl += 2) { 1184 const PetscInt point = closure[cl]; 1185 PetscInt vertexLoc; 1186 1187 if ((point >= pStart[0]) && (point < pEnd[0])) { 1188 ++numCorners; 1189 ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr); 1190 if (vertexLoc == value) closure[faceSize++] = point; 1191 } 1192 } 1193 if (!(*nFV)) {ierr = DMPlexGetNumFaceVertices(dm, dim, numCorners, nFV);CHKERRQ(ierr);} 1194 if (faceSize > *nFV) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize); 1195 if (faceSize == *nFV) { 1196 const PetscInt *cells = NULL; 1197 PetscInt numCells, nc; 1198 1199 ++(*numFaces); 1200 for (cl = 0; cl < faceSize; ++cl) { 1201 ierr = DMLabelSetValue(subpointMap, closure[cl], 0);CHKERRQ(ierr); 1202 } 1203 ierr = DMPlexGetJoin(dm, faceSize, closure, &numCells, &cells);CHKERRQ(ierr); 1204 for (nc = 0; nc < numCells; ++nc) { 1205 ierr = DMLabelSetValue(subpointMap, cells[nc], 2);CHKERRQ(ierr); 1206 } 1207 ierr = DMPlexRestoreJoin(dm, faceSize, closure, &numCells, &cells);CHKERRQ(ierr); 1208 } 1209 ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1210 } 1211 ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1212 } 1213 if (subvertexIS) { 1214 ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1215 } 1216 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 1217 ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr); 1218 PetscFunctionReturn(0); 1219 } 1220 1221 #undef __FUNCT__ 1222 #define __FUNCT__ "DMPlexMarkSubmesh_Interpolated" 1223 static PetscErrorCode DMPlexMarkSubmesh_Interpolated(DM dm, DMLabel vertexLabel, PetscInt value, DMLabel subpointMap, DM subdm) 1224 { 1225 IS subvertexIS; 1226 const PetscInt *subvertices; 1227 PetscInt *pStart, *pEnd, *pMax; 1228 PetscInt dim, d, numSubVerticesInitial = 0, v; 1229 PetscErrorCode ierr; 1230 1231 PetscFunctionBegin; 1232 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1233 ierr = PetscMalloc3(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd,dim+1,PetscInt,&pMax);CHKERRQ(ierr); 1234 ierr = DMPlexGetHybridBounds(dm, &pMax[dim], dim>1 ? &pMax[dim-1] : NULL, dim > 2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr); 1235 for (d = 0; d <= dim; ++d) { 1236 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1237 if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]); 1238 } 1239 /* Loop over initial vertices and mark all faces in the collective star() */ 1240 ierr = DMLabelGetStratumIS(vertexLabel, value, &subvertexIS);CHKERRQ(ierr); 1241 if (subvertexIS) { 1242 ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr); 1243 ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1244 } 1245 for (v = 0; v < numSubVerticesInitial; ++v) { 1246 const PetscInt vertex = subvertices[v]; 1247 PetscInt *star = NULL; 1248 PetscInt starSize, s, numFaces = 0, f; 1249 1250 ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1251 for (s = 0; s < starSize*2; s += 2) { 1252 const PetscInt point = star[s]; 1253 if ((point >= pStart[dim-1]) && (point < pEnd[dim-1])) star[numFaces++] = point; 1254 } 1255 for (f = 0; f < numFaces; ++f) { 1256 const PetscInt face = star[f]; 1257 PetscInt *closure = NULL; 1258 PetscInt closureSize, c; 1259 PetscInt faceLoc; 1260 1261 ierr = DMLabelGetValue(subpointMap, face, &faceLoc);CHKERRQ(ierr); 1262 if (faceLoc == dim-1) continue; 1263 if (faceLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Face %d has dimension %d in the surface label", face, faceLoc); 1264 ierr = DMPlexGetTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1265 for (c = 0; c < closureSize*2; c += 2) { 1266 const PetscInt point = closure[c]; 1267 PetscInt vertexLoc; 1268 1269 if ((point >= pStart[0]) && (point < pEnd[0])) { 1270 ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr); 1271 if (vertexLoc != value) break; 1272 } 1273 } 1274 if (c == closureSize*2) { 1275 const PetscInt *support; 1276 PetscInt supportSize, s; 1277 1278 for (c = 0; c < closureSize*2; c += 2) { 1279 const PetscInt point = closure[c]; 1280 1281 for (d = 0; d < dim; ++d) { 1282 if ((point >= pStart[d]) && (point < pEnd[d])) { 1283 ierr = DMLabelSetValue(subpointMap, point, d);CHKERRQ(ierr); 1284 break; 1285 } 1286 } 1287 } 1288 ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 1289 ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); 1290 for (s = 0; s < supportSize; ++s) { 1291 ierr = DMLabelSetValue(subpointMap, support[s], dim);CHKERRQ(ierr); 1292 } 1293 } 1294 ierr = DMPlexRestoreTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1295 } 1296 ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1297 } 1298 if (subvertexIS) { 1299 ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1300 } 1301 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 1302 ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr); 1303 PetscFunctionReturn(0); 1304 } 1305 1306 #undef __FUNCT__ 1307 #define __FUNCT__ "DMPlexMarkCohesiveSubmesh_Uninterpolated" 1308 static PetscErrorCode DMPlexMarkCohesiveSubmesh_Uninterpolated(DM dm, PetscBool hasLagrange, DMLabel subpointMap, PetscInt *numFaces, PetscInt *nFV, PetscInt *subCells[], DM subdm) 1309 { 1310 const PetscInt *cone; 1311 PetscInt dim, cMax, cEnd, c, p, coneSize; 1312 PetscErrorCode ierr; 1313 1314 PetscFunctionBegin; 1315 *numFaces = 0; 1316 *nFV = 0; 1317 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1318 ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr); 1319 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1320 if (cMax < 0) PetscFunctionReturn(0); 1321 ierr = DMPlexGetConeSize(dm, cMax, &coneSize);CHKERRQ(ierr); 1322 *numFaces = cEnd - cMax; 1323 *nFV = hasLagrange ? coneSize/3 : coneSize/2; 1324 ierr = PetscMalloc(*numFaces *2 * sizeof(PetscInt), subCells);CHKERRQ(ierr); 1325 for (c = cMax; c < cEnd; ++c) { 1326 const PetscInt *cells; 1327 PetscInt numCells; 1328 1329 ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 1330 for (p = 0; p < *nFV; ++p) { 1331 ierr = DMLabelSetValue(subpointMap, cone[p], 0);CHKERRQ(ierr); 1332 } 1333 /* Negative face */ 1334 ierr = DMPlexGetJoin(dm, *nFV, cone, &numCells, &cells);CHKERRQ(ierr); 1335 if (numCells != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells"); 1336 for (p = 0; p < numCells; ++p) { 1337 ierr = DMLabelSetValue(subpointMap, cells[p], 2);CHKERRQ(ierr); 1338 (*subCells)[(c-cMax)*2+p] = cells[p]; 1339 } 1340 ierr = DMPlexRestoreJoin(dm, *nFV, cone, &numCells, &cells);CHKERRQ(ierr); 1341 /* Positive face is not included */ 1342 } 1343 PetscFunctionReturn(0); 1344 } 1345 1346 #undef __FUNCT__ 1347 #define __FUNCT__ "DMPlexMarkCohesiveSubmesh_Interpolated" 1348 static PetscErrorCode DMPlexMarkCohesiveSubmesh_Interpolated(DM dm, PetscBool hasLagrange, DMLabel subpointMap, DM subdm) 1349 { 1350 PetscInt *pStart, *pEnd; 1351 PetscInt dim, cMax, cEnd, c, d; 1352 PetscErrorCode ierr; 1353 1354 PetscFunctionBegin; 1355 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1356 ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr); 1357 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1358 if (cMax < 0) PetscFunctionReturn(0); 1359 ierr = PetscMalloc2(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd);CHKERRQ(ierr); 1360 for (d = 0; d <= dim; ++d) { 1361 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1362 } 1363 for (c = cMax; c < cEnd; ++c) { 1364 const PetscInt *cone; 1365 PetscInt *closure = NULL; 1366 PetscInt coneSize, closureSize, cl; 1367 1368 ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 1369 if (hasLagrange) { 1370 if (coneSize != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells"); 1371 } else { 1372 if (coneSize != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells"); 1373 } 1374 /* Negative face */ 1375 ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 1376 ierr = DMPlexGetTransitiveClosure(dm, cone[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1377 for (cl = 0; cl < closureSize*2; cl += 2) { 1378 const PetscInt point = closure[cl]; 1379 1380 for (d = 0; d <= dim; ++d) { 1381 if ((point >= pStart[d]) && (point < pEnd[d])) { 1382 ierr = DMLabelSetValue(subpointMap, point, d);CHKERRQ(ierr); 1383 break; 1384 } 1385 } 1386 } 1387 ierr = DMPlexRestoreTransitiveClosure(dm, cone[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1388 /* Cells -- positive face is not included */ 1389 for (cl = 0; cl < 1; ++cl) { 1390 const PetscInt *support; 1391 PetscInt supportSize, s; 1392 1393 ierr = DMPlexGetSupportSize(dm, cone[cl], &supportSize);CHKERRQ(ierr); 1394 if (supportSize != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive faces should separate two cells"); 1395 ierr = DMPlexGetSupport(dm, cone[cl], &support);CHKERRQ(ierr); 1396 for (s = 0; s < supportSize; ++s) { 1397 ierr = DMLabelSetValue(subpointMap, support[s], dim);CHKERRQ(ierr); 1398 } 1399 } 1400 } 1401 ierr = PetscFree2(pStart, pEnd);CHKERRQ(ierr); 1402 PetscFunctionReturn(0); 1403 } 1404 1405 #undef __FUNCT__ 1406 #define __FUNCT__ "DMPlexGetFaceOrientation" 1407 PetscErrorCode DMPlexGetFaceOrientation(DM dm, PetscInt cell, PetscInt numCorners, PetscInt indices[], PetscInt oppositeVertex, PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented) 1408 { 1409 MPI_Comm comm; 1410 PetscBool posOrient = PETSC_FALSE; 1411 const PetscInt debug = 0; 1412 PetscInt cellDim, faceSize, f; 1413 PetscErrorCode ierr; 1414 1415 PetscFunctionBegin; 1416 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1417 ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr); 1418 if (debug) {PetscPrintf(comm, "cellDim: %d numCorners: %d\n", cellDim, numCorners);CHKERRQ(ierr);} 1419 1420 if (cellDim == numCorners-1) { 1421 /* Simplices */ 1422 faceSize = numCorners-1; 1423 posOrient = !(oppositeVertex%2) ? PETSC_TRUE : PETSC_FALSE; 1424 } else if (cellDim == 1 && numCorners == 3) { 1425 /* Quadratic line */ 1426 faceSize = 1; 1427 posOrient = PETSC_TRUE; 1428 } else if (cellDim == 2 && numCorners == 4) { 1429 /* Quads */ 1430 faceSize = 2; 1431 if ((indices[1] > indices[0]) && (indices[1] - indices[0] == 1)) { 1432 posOrient = PETSC_TRUE; 1433 } else if ((indices[0] == 3) && (indices[1] == 0)) { 1434 posOrient = PETSC_TRUE; 1435 } else { 1436 if (((indices[0] > indices[1]) && (indices[0] - indices[1] == 1)) || ((indices[0] == 0) && (indices[1] == 3))) { 1437 posOrient = PETSC_FALSE; 1438 } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossedge"); 1439 } 1440 } else if (cellDim == 2 && numCorners == 6) { 1441 /* Quadratic triangle (I hate this) */ 1442 /* Edges are determined by the first 2 vertices (corners of edges) */ 1443 const PetscInt faceSizeTri = 3; 1444 PetscInt sortedIndices[3], i, iFace; 1445 PetscBool found = PETSC_FALSE; 1446 PetscInt faceVerticesTriSorted[9] = { 1447 0, 3, 4, /* bottom */ 1448 1, 4, 5, /* right */ 1449 2, 3, 5, /* left */ 1450 }; 1451 PetscInt faceVerticesTri[9] = { 1452 0, 3, 4, /* bottom */ 1453 1, 4, 5, /* right */ 1454 2, 5, 3, /* left */ 1455 }; 1456 1457 faceSize = faceSizeTri; 1458 for (i = 0; i < faceSizeTri; ++i) sortedIndices[i] = indices[i]; 1459 ierr = PetscSortInt(faceSizeTri, sortedIndices);CHKERRQ(ierr); 1460 for (iFace = 0; iFace < 3; ++iFace) { 1461 const PetscInt ii = iFace*faceSizeTri; 1462 PetscInt fVertex, cVertex; 1463 1464 if ((sortedIndices[0] == faceVerticesTriSorted[ii+0]) && 1465 (sortedIndices[1] == faceVerticesTriSorted[ii+1])) { 1466 for (fVertex = 0; fVertex < faceSizeTri; ++fVertex) { 1467 for (cVertex = 0; cVertex < faceSizeTri; ++cVertex) { 1468 if (indices[cVertex] == faceVerticesTri[ii+fVertex]) { 1469 faceVertices[fVertex] = origVertices[cVertex]; 1470 break; 1471 } 1472 } 1473 } 1474 found = PETSC_TRUE; 1475 break; 1476 } 1477 } 1478 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tri crossface"); 1479 if (posOriented) *posOriented = PETSC_TRUE; 1480 PetscFunctionReturn(0); 1481 } else if (cellDim == 2 && numCorners == 9) { 1482 /* Quadratic quad (I hate this) */ 1483 /* Edges are determined by the first 2 vertices (corners of edges) */ 1484 const PetscInt faceSizeQuad = 3; 1485 PetscInt sortedIndices[3], i, iFace; 1486 PetscBool found = PETSC_FALSE; 1487 PetscInt faceVerticesQuadSorted[12] = { 1488 0, 1, 4, /* bottom */ 1489 1, 2, 5, /* right */ 1490 2, 3, 6, /* top */ 1491 0, 3, 7, /* left */ 1492 }; 1493 PetscInt faceVerticesQuad[12] = { 1494 0, 1, 4, /* bottom */ 1495 1, 2, 5, /* right */ 1496 2, 3, 6, /* top */ 1497 3, 0, 7, /* left */ 1498 }; 1499 1500 faceSize = faceSizeQuad; 1501 for (i = 0; i < faceSizeQuad; ++i) sortedIndices[i] = indices[i]; 1502 ierr = PetscSortInt(faceSizeQuad, sortedIndices);CHKERRQ(ierr); 1503 for (iFace = 0; iFace < 4; ++iFace) { 1504 const PetscInt ii = iFace*faceSizeQuad; 1505 PetscInt fVertex, cVertex; 1506 1507 if ((sortedIndices[0] == faceVerticesQuadSorted[ii+0]) && 1508 (sortedIndices[1] == faceVerticesQuadSorted[ii+1])) { 1509 for (fVertex = 0; fVertex < faceSizeQuad; ++fVertex) { 1510 for (cVertex = 0; cVertex < faceSizeQuad; ++cVertex) { 1511 if (indices[cVertex] == faceVerticesQuad[ii+fVertex]) { 1512 faceVertices[fVertex] = origVertices[cVertex]; 1513 break; 1514 } 1515 } 1516 } 1517 found = PETSC_TRUE; 1518 break; 1519 } 1520 } 1521 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossface"); 1522 if (posOriented) *posOriented = PETSC_TRUE; 1523 PetscFunctionReturn(0); 1524 } else if (cellDim == 3 && numCorners == 8) { 1525 /* Hexes 1526 A hex is two oriented quads with the normal of the first 1527 pointing up at the second. 1528 1529 7---6 1530 /| /| 1531 4---5 | 1532 | 3-|-2 1533 |/ |/ 1534 0---1 1535 1536 Faces are determined by the first 4 vertices (corners of faces) */ 1537 const PetscInt faceSizeHex = 4; 1538 PetscInt sortedIndices[4], i, iFace; 1539 PetscBool found = PETSC_FALSE; 1540 PetscInt faceVerticesHexSorted[24] = { 1541 0, 1, 2, 3, /* bottom */ 1542 4, 5, 6, 7, /* top */ 1543 0, 1, 4, 5, /* front */ 1544 1, 2, 5, 6, /* right */ 1545 2, 3, 6, 7, /* back */ 1546 0, 3, 4, 7, /* left */ 1547 }; 1548 PetscInt faceVerticesHex[24] = { 1549 3, 2, 1, 0, /* bottom */ 1550 4, 5, 6, 7, /* top */ 1551 0, 1, 5, 4, /* front */ 1552 1, 2, 6, 5, /* right */ 1553 2, 3, 7, 6, /* back */ 1554 3, 0, 4, 7, /* left */ 1555 }; 1556 1557 faceSize = faceSizeHex; 1558 for (i = 0; i < faceSizeHex; ++i) sortedIndices[i] = indices[i]; 1559 ierr = PetscSortInt(faceSizeHex, sortedIndices);CHKERRQ(ierr); 1560 for (iFace = 0; iFace < 6; ++iFace) { 1561 const PetscInt ii = iFace*faceSizeHex; 1562 PetscInt fVertex, cVertex; 1563 1564 if ((sortedIndices[0] == faceVerticesHexSorted[ii+0]) && 1565 (sortedIndices[1] == faceVerticesHexSorted[ii+1]) && 1566 (sortedIndices[2] == faceVerticesHexSorted[ii+2]) && 1567 (sortedIndices[3] == faceVerticesHexSorted[ii+3])) { 1568 for (fVertex = 0; fVertex < faceSizeHex; ++fVertex) { 1569 for (cVertex = 0; cVertex < faceSizeHex; ++cVertex) { 1570 if (indices[cVertex] == faceVerticesHex[ii+fVertex]) { 1571 faceVertices[fVertex] = origVertices[cVertex]; 1572 break; 1573 } 1574 } 1575 } 1576 found = PETSC_TRUE; 1577 break; 1578 } 1579 } 1580 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface"); 1581 if (posOriented) *posOriented = PETSC_TRUE; 1582 PetscFunctionReturn(0); 1583 } else if (cellDim == 3 && numCorners == 10) { 1584 /* Quadratic tet */ 1585 /* Faces are determined by the first 3 vertices (corners of faces) */ 1586 const PetscInt faceSizeTet = 6; 1587 PetscInt sortedIndices[6], i, iFace; 1588 PetscBool found = PETSC_FALSE; 1589 PetscInt faceVerticesTetSorted[24] = { 1590 0, 1, 2, 6, 7, 8, /* bottom */ 1591 0, 3, 4, 6, 7, 9, /* front */ 1592 1, 4, 5, 7, 8, 9, /* right */ 1593 2, 3, 5, 6, 8, 9, /* left */ 1594 }; 1595 PetscInt faceVerticesTet[24] = { 1596 0, 1, 2, 6, 7, 8, /* bottom */ 1597 0, 4, 3, 6, 7, 9, /* front */ 1598 1, 5, 4, 7, 8, 9, /* right */ 1599 2, 3, 5, 8, 6, 9, /* left */ 1600 }; 1601 1602 faceSize = faceSizeTet; 1603 for (i = 0; i < faceSizeTet; ++i) sortedIndices[i] = indices[i]; 1604 ierr = PetscSortInt(faceSizeTet, sortedIndices);CHKERRQ(ierr); 1605 for (iFace=0; iFace < 4; ++iFace) { 1606 const PetscInt ii = iFace*faceSizeTet; 1607 PetscInt fVertex, cVertex; 1608 1609 if ((sortedIndices[0] == faceVerticesTetSorted[ii+0]) && 1610 (sortedIndices[1] == faceVerticesTetSorted[ii+1]) && 1611 (sortedIndices[2] == faceVerticesTetSorted[ii+2]) && 1612 (sortedIndices[3] == faceVerticesTetSorted[ii+3])) { 1613 for (fVertex = 0; fVertex < faceSizeTet; ++fVertex) { 1614 for (cVertex = 0; cVertex < faceSizeTet; ++cVertex) { 1615 if (indices[cVertex] == faceVerticesTet[ii+fVertex]) { 1616 faceVertices[fVertex] = origVertices[cVertex]; 1617 break; 1618 } 1619 } 1620 } 1621 found = PETSC_TRUE; 1622 break; 1623 } 1624 } 1625 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tet crossface"); 1626 if (posOriented) *posOriented = PETSC_TRUE; 1627 PetscFunctionReturn(0); 1628 } else if (cellDim == 3 && numCorners == 27) { 1629 /* Quadratic hexes (I hate this) 1630 A hex is two oriented quads with the normal of the first 1631 pointing up at the second. 1632 1633 7---6 1634 /| /| 1635 4---5 | 1636 | 3-|-2 1637 |/ |/ 1638 0---1 1639 1640 Faces are determined by the first 4 vertices (corners of faces) */ 1641 const PetscInt faceSizeQuadHex = 9; 1642 PetscInt sortedIndices[9], i, iFace; 1643 PetscBool found = PETSC_FALSE; 1644 PetscInt faceVerticesQuadHexSorted[54] = { 1645 0, 1, 2, 3, 8, 9, 10, 11, 24, /* bottom */ 1646 4, 5, 6, 7, 12, 13, 14, 15, 25, /* top */ 1647 0, 1, 4, 5, 8, 12, 16, 17, 22, /* front */ 1648 1, 2, 5, 6, 9, 13, 17, 18, 21, /* right */ 1649 2, 3, 6, 7, 10, 14, 18, 19, 23, /* back */ 1650 0, 3, 4, 7, 11, 15, 16, 19, 20, /* left */ 1651 }; 1652 PetscInt faceVerticesQuadHex[54] = { 1653 3, 2, 1, 0, 10, 9, 8, 11, 24, /* bottom */ 1654 4, 5, 6, 7, 12, 13, 14, 15, 25, /* top */ 1655 0, 1, 5, 4, 8, 17, 12, 16, 22, /* front */ 1656 1, 2, 6, 5, 9, 18, 13, 17, 21, /* right */ 1657 2, 3, 7, 6, 10, 19, 14, 18, 23, /* back */ 1658 3, 0, 4, 7, 11, 16, 15, 19, 20 /* left */ 1659 }; 1660 1661 faceSize = faceSizeQuadHex; 1662 for (i = 0; i < faceSizeQuadHex; ++i) sortedIndices[i] = indices[i]; 1663 ierr = PetscSortInt(faceSizeQuadHex, sortedIndices);CHKERRQ(ierr); 1664 for (iFace = 0; iFace < 6; ++iFace) { 1665 const PetscInt ii = iFace*faceSizeQuadHex; 1666 PetscInt fVertex, cVertex; 1667 1668 if ((sortedIndices[0] == faceVerticesQuadHexSorted[ii+0]) && 1669 (sortedIndices[1] == faceVerticesQuadHexSorted[ii+1]) && 1670 (sortedIndices[2] == faceVerticesQuadHexSorted[ii+2]) && 1671 (sortedIndices[3] == faceVerticesQuadHexSorted[ii+3])) { 1672 for (fVertex = 0; fVertex < faceSizeQuadHex; ++fVertex) { 1673 for (cVertex = 0; cVertex < faceSizeQuadHex; ++cVertex) { 1674 if (indices[cVertex] == faceVerticesQuadHex[ii+fVertex]) { 1675 faceVertices[fVertex] = origVertices[cVertex]; 1676 break; 1677 } 1678 } 1679 } 1680 found = PETSC_TRUE; 1681 break; 1682 } 1683 } 1684 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface"); 1685 if (posOriented) *posOriented = PETSC_TRUE; 1686 PetscFunctionReturn(0); 1687 } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Unknown cell type for faceOrientation()."); 1688 if (!posOrient) { 1689 if (debug) {ierr = PetscPrintf(comm, " Reversing initial face orientation\n");CHKERRQ(ierr);} 1690 for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[faceSize-1 - f]; 1691 } else { 1692 if (debug) {ierr = PetscPrintf(comm, " Keeping initial face orientation\n");CHKERRQ(ierr);} 1693 for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[f]; 1694 } 1695 if (posOriented) *posOriented = posOrient; 1696 PetscFunctionReturn(0); 1697 } 1698 1699 #undef __FUNCT__ 1700 #define __FUNCT__ "DMPlexGetOrientedFace" 1701 /* 1702 Given a cell and a face, as a set of vertices, 1703 return the oriented face, as a set of vertices, in faceVertices 1704 The orientation is such that the face normal points out of the cell 1705 */ 1706 PetscErrorCode DMPlexGetOrientedFace(DM dm, PetscInt cell, PetscInt faceSize, const PetscInt face[], PetscInt numCorners, PetscInt indices[], PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented) 1707 { 1708 const PetscInt *cone = NULL; 1709 PetscInt coneSize, v, f, v2; 1710 PetscInt oppositeVertex = -1; 1711 PetscErrorCode ierr; 1712 1713 PetscFunctionBegin; 1714 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 1715 ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 1716 for (v = 0, v2 = 0; v < coneSize; ++v) { 1717 PetscBool found = PETSC_FALSE; 1718 1719 for (f = 0; f < faceSize; ++f) { 1720 if (face[f] == cone[v]) { 1721 found = PETSC_TRUE; break; 1722 } 1723 } 1724 if (found) { 1725 indices[v2] = v; 1726 origVertices[v2] = cone[v]; 1727 ++v2; 1728 } else { 1729 oppositeVertex = v; 1730 } 1731 } 1732 ierr = DMPlexGetFaceOrientation(dm, cell, numCorners, indices, oppositeVertex, origVertices, faceVertices, posOriented);CHKERRQ(ierr); 1733 PetscFunctionReturn(0); 1734 } 1735 1736 #undef __FUNCT__ 1737 #define __FUNCT__ "DMPlexInsertFace_Internal" 1738 /* 1739 DMPlexInsertFace_Internal - Puts a face into the mesh 1740 1741 Not collective 1742 1743 Input Parameters: 1744 + dm - The DMPlex 1745 . numFaceVertex - The number of vertices in the face 1746 . faceVertices - The vertices in the face for dm 1747 . subfaceVertices - The vertices in the face for subdm 1748 . numCorners - The number of vertices in the cell 1749 . cell - A cell in dm containing the face 1750 . subcell - A cell in subdm containing the face 1751 . firstFace - First face in the mesh 1752 - newFacePoint - Next face in the mesh 1753 1754 Output Parameters: 1755 . newFacePoint - Contains next face point number on input, updated on output 1756 1757 Level: developer 1758 */ 1759 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) 1760 { 1761 MPI_Comm comm; 1762 DM_Plex *submesh = (DM_Plex*) subdm->data; 1763 const PetscInt *faces; 1764 PetscInt numFaces, coneSize; 1765 PetscErrorCode ierr; 1766 1767 PetscFunctionBegin; 1768 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1769 ierr = DMPlexGetConeSize(subdm, subcell, &coneSize);CHKERRQ(ierr); 1770 if (coneSize != 1) SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cone size of cell %d is %d != 1", cell, coneSize); 1771 #if 0 1772 /* Cannot use this because support() has not been constructed yet */ 1773 ierr = DMPlexGetJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr); 1774 #else 1775 { 1776 PetscInt f; 1777 1778 numFaces = 0; 1779 ierr = DMGetWorkArray(subdm, 1, PETSC_INT, (void **) &faces);CHKERRQ(ierr); 1780 for (f = firstFace; f < *newFacePoint; ++f) { 1781 PetscInt dof, off, d; 1782 1783 ierr = PetscSectionGetDof(submesh->coneSection, f, &dof);CHKERRQ(ierr); 1784 ierr = PetscSectionGetOffset(submesh->coneSection, f, &off);CHKERRQ(ierr); 1785 /* Yes, I know this is quadratic, but I expect the sizes to be <5 */ 1786 for (d = 0; d < dof; ++d) { 1787 const PetscInt p = submesh->cones[off+d]; 1788 PetscInt v; 1789 1790 for (v = 0; v < numFaceVertices; ++v) { 1791 if (subfaceVertices[v] == p) break; 1792 } 1793 if (v == numFaceVertices) break; 1794 } 1795 if (d == dof) { 1796 numFaces = 1; 1797 ((PetscInt*) faces)[0] = f; 1798 } 1799 } 1800 } 1801 #endif 1802 if (numFaces > 1) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Vertex set had %d faces, not one", numFaces); 1803 else if (numFaces == 1) { 1804 /* Add the other cell neighbor for this face */ 1805 ierr = DMPlexSetCone(subdm, subcell, faces);CHKERRQ(ierr); 1806 } else { 1807 PetscInt *indices, *origVertices, *orientedVertices, *orientedSubVertices, v, ov; 1808 PetscBool posOriented; 1809 1810 ierr = DMGetWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr); 1811 origVertices = &orientedVertices[numFaceVertices]; 1812 indices = &orientedVertices[numFaceVertices*2]; 1813 orientedSubVertices = &orientedVertices[numFaceVertices*3]; 1814 ierr = DMPlexGetOrientedFace(dm, cell, numFaceVertices, faceVertices, numCorners, indices, origVertices, orientedVertices, &posOriented);CHKERRQ(ierr); 1815 /* TODO: I know that routine should return a permutation, not the indices */ 1816 for (v = 0; v < numFaceVertices; ++v) { 1817 const PetscInt vertex = faceVertices[v], subvertex = subfaceVertices[v]; 1818 for (ov = 0; ov < numFaceVertices; ++ov) { 1819 if (orientedVertices[ov] == vertex) { 1820 orientedSubVertices[ov] = subvertex; 1821 break; 1822 } 1823 } 1824 if (ov == numFaceVertices) SETERRQ1(comm, PETSC_ERR_PLIB, "Could not find face vertex %d in orientated set", vertex); 1825 } 1826 ierr = DMPlexSetCone(subdm, *newFacePoint, orientedSubVertices);CHKERRQ(ierr); 1827 ierr = DMPlexSetCone(subdm, subcell, newFacePoint);CHKERRQ(ierr); 1828 ierr = DMRestoreWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr); 1829 ++(*newFacePoint); 1830 } 1831 #if 0 1832 ierr = DMPlexRestoreJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr); 1833 #else 1834 ierr = DMRestoreWorkArray(subdm, 1, PETSC_INT, (void **) &faces);CHKERRQ(ierr); 1835 #endif 1836 PetscFunctionReturn(0); 1837 } 1838 1839 #undef __FUNCT__ 1840 #define __FUNCT__ "DMPlexCreateSubmesh_Uninterpolated" 1841 static PetscErrorCode DMPlexCreateSubmesh_Uninterpolated(DM dm, const char vertexLabelName[], PetscInt value, DM subdm) 1842 { 1843 MPI_Comm comm; 1844 DMLabel vertexLabel, subpointMap; 1845 IS subvertexIS, subcellIS; 1846 const PetscInt *subVertices, *subCells; 1847 PetscInt numSubVertices, firstSubVertex, numSubCells; 1848 PetscInt *subface, maxConeSize, numSubFaces, firstSubFace, newFacePoint, nFV; 1849 PetscInt vStart, vEnd, c, f; 1850 PetscErrorCode ierr; 1851 1852 PetscFunctionBegin; 1853 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1854 /* Create subpointMap which marks the submesh */ 1855 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 1856 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 1857 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 1858 ierr = DMPlexGetLabel(dm, vertexLabelName, &vertexLabel);CHKERRQ(ierr); 1859 ierr = DMPlexMarkSubmesh_Uninterpolated(dm, vertexLabel, value, subpointMap, &numSubFaces, &nFV, subdm);CHKERRQ(ierr); 1860 /* Setup chart */ 1861 ierr = DMLabelGetStratumSize(subpointMap, 0, &numSubVertices);CHKERRQ(ierr); 1862 ierr = DMLabelGetStratumSize(subpointMap, 2, &numSubCells);CHKERRQ(ierr); 1863 ierr = DMPlexSetChart(subdm, 0, numSubCells+numSubFaces+numSubVertices);CHKERRQ(ierr); 1864 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 1865 /* Set cone sizes */ 1866 firstSubVertex = numSubCells; 1867 firstSubFace = numSubCells+numSubVertices; 1868 newFacePoint = firstSubFace; 1869 ierr = DMLabelGetStratumIS(subpointMap, 0, &subvertexIS);CHKERRQ(ierr); 1870 if (subvertexIS) {ierr = ISGetIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 1871 ierr = DMLabelGetStratumIS(subpointMap, 2, &subcellIS);CHKERRQ(ierr); 1872 if (subcellIS) {ierr = ISGetIndices(subcellIS, &subCells);CHKERRQ(ierr);} 1873 for (c = 0; c < numSubCells; ++c) { 1874 ierr = DMPlexSetConeSize(subdm, c, 1);CHKERRQ(ierr); 1875 } 1876 for (f = firstSubFace; f < firstSubFace+numSubFaces; ++f) { 1877 ierr = DMPlexSetConeSize(subdm, f, nFV);CHKERRQ(ierr); 1878 } 1879 ierr = DMSetUp(subdm);CHKERRQ(ierr); 1880 /* Create face cones */ 1881 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1882 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 1883 ierr = DMGetWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 1884 for (c = 0; c < numSubCells; ++c) { 1885 const PetscInt cell = subCells[c]; 1886 const PetscInt subcell = c; 1887 PetscInt *closure = NULL; 1888 PetscInt closureSize, cl, numCorners = 0, faceSize = 0; 1889 1890 ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1891 for (cl = 0; cl < closureSize*2; cl += 2) { 1892 const PetscInt point = closure[cl]; 1893 PetscInt subVertex; 1894 1895 if ((point >= vStart) && (point < vEnd)) { 1896 ++numCorners; 1897 ierr = PetscFindInt(point, numSubVertices, subVertices, &subVertex);CHKERRQ(ierr); 1898 if (subVertex >= 0) { 1899 closure[faceSize] = point; 1900 subface[faceSize] = firstSubVertex+subVertex; 1901 ++faceSize; 1902 } 1903 } 1904 } 1905 if (faceSize > nFV) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize); 1906 if (faceSize == nFV) { 1907 ierr = DMPlexInsertFace_Internal(dm, subdm, faceSize, closure, subface, numCorners, cell, subcell, firstSubFace, &newFacePoint);CHKERRQ(ierr); 1908 } 1909 ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1910 } 1911 ierr = DMRestoreWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 1912 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 1913 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 1914 /* Build coordinates */ 1915 { 1916 PetscSection coordSection, subCoordSection; 1917 Vec coordinates, subCoordinates; 1918 PetscScalar *coords, *subCoords; 1919 PetscInt numComp, coordSize, v; 1920 1921 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1922 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1923 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 1924 ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr); 1925 ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr); 1926 ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr); 1927 ierr = PetscSectionSetChart(subCoordSection, firstSubVertex, firstSubVertex+numSubVertices);CHKERRQ(ierr); 1928 for (v = 0; v < numSubVertices; ++v) { 1929 const PetscInt vertex = subVertices[v]; 1930 const PetscInt subvertex = firstSubVertex+v; 1931 PetscInt dof; 1932 1933 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 1934 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 1935 ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr); 1936 } 1937 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 1938 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 1939 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 1940 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1941 ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr); 1942 if (coordSize) { 1943 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1944 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 1945 for (v = 0; v < numSubVertices; ++v) { 1946 const PetscInt vertex = subVertices[v]; 1947 const PetscInt subvertex = firstSubVertex+v; 1948 PetscInt dof, off, sdof, soff, d; 1949 1950 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 1951 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 1952 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 1953 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 1954 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 1955 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 1956 } 1957 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1958 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 1959 } 1960 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 1961 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 1962 } 1963 /* Cleanup */ 1964 if (subvertexIS) {ierr = ISRestoreIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 1965 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 1966 if (subcellIS) {ierr = ISRestoreIndices(subcellIS, &subCells);CHKERRQ(ierr);} 1967 ierr = ISDestroy(&subcellIS);CHKERRQ(ierr); 1968 PetscFunctionReturn(0); 1969 } 1970 1971 #undef __FUNCT__ 1972 #define __FUNCT__ "DMPlexCreateSubmesh_Interpolated" 1973 static PetscErrorCode DMPlexCreateSubmesh_Interpolated(DM dm, const char vertexLabelName[], PetscInt value, DM subdm) 1974 { 1975 MPI_Comm comm; 1976 DMLabel subpointMap, vertexLabel; 1977 IS *subpointIS; 1978 const PetscInt **subpoints; 1979 PetscInt *numSubPoints, *firstSubPoint, *coneNew; 1980 PetscInt totSubPoints = 0, maxConeSize, dim, p, d, v; 1981 PetscErrorCode ierr; 1982 1983 PetscFunctionBegin; 1984 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1985 /* Create subpointMap which marks the submesh */ 1986 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 1987 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 1988 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 1989 ierr = DMPlexGetLabel(dm, vertexLabelName, &vertexLabel);CHKERRQ(ierr); 1990 ierr = DMPlexMarkSubmesh_Interpolated(dm, vertexLabel, value, subpointMap, subdm);CHKERRQ(ierr); 1991 /* Setup chart */ 1992 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1993 ierr = PetscMalloc4(dim+1,PetscInt,&numSubPoints,dim+1,PetscInt,&firstSubPoint,dim+1,IS,&subpointIS,dim+1,const PetscInt *,&subpoints);CHKERRQ(ierr); 1994 for (d = 0; d <= dim; ++d) { 1995 ierr = DMLabelGetStratumSize(subpointMap, d, &numSubPoints[d]);CHKERRQ(ierr); 1996 totSubPoints += numSubPoints[d]; 1997 } 1998 ierr = DMPlexSetChart(subdm, 0, totSubPoints);CHKERRQ(ierr); 1999 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 2000 /* Set cone sizes */ 2001 firstSubPoint[dim] = 0; 2002 firstSubPoint[0] = firstSubPoint[dim] + numSubPoints[dim]; 2003 if (dim > 1) {firstSubPoint[dim-1] = firstSubPoint[0] + numSubPoints[0];} 2004 if (dim > 2) {firstSubPoint[dim-2] = firstSubPoint[dim-1] + numSubPoints[dim-1];} 2005 for (d = 0; d <= dim; ++d) { 2006 ierr = DMLabelGetStratumIS(subpointMap, d, &subpointIS[d]);CHKERRQ(ierr); 2007 ierr = ISGetIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 2008 } 2009 for (d = 0; d <= dim; ++d) { 2010 for (p = 0; p < numSubPoints[d]; ++p) { 2011 const PetscInt point = subpoints[d][p]; 2012 const PetscInt subpoint = firstSubPoint[d] + p; 2013 const PetscInt *cone; 2014 PetscInt coneSize, coneSizeNew, c, val; 2015 2016 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 2017 ierr = DMPlexSetConeSize(subdm, subpoint, coneSize);CHKERRQ(ierr); 2018 if (d == dim) { 2019 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 2020 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 2021 ierr = DMLabelGetValue(subpointMap, cone[c], &val);CHKERRQ(ierr); 2022 if (val >= 0) coneSizeNew++; 2023 } 2024 ierr = DMPlexSetConeSize(subdm, subpoint, coneSizeNew);CHKERRQ(ierr); 2025 } 2026 } 2027 } 2028 ierr = DMSetUp(subdm);CHKERRQ(ierr); 2029 /* Set cones */ 2030 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 2031 ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &coneNew);CHKERRQ(ierr); 2032 for (d = 0; d <= dim; ++d) { 2033 for (p = 0; p < numSubPoints[d]; ++p) { 2034 const PetscInt point = subpoints[d][p]; 2035 const PetscInt subpoint = firstSubPoint[d] + p; 2036 const PetscInt *cone, *ornt; 2037 PetscInt coneSize, subconeSize, coneSizeNew, c, subc; 2038 2039 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 2040 ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr); 2041 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 2042 ierr = DMPlexGetConeOrientation(dm, point, &ornt);CHKERRQ(ierr); 2043 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 2044 ierr = PetscFindInt(cone[c], numSubPoints[d-1], subpoints[d-1], &subc);CHKERRQ(ierr); 2045 if (subc >= 0) coneNew[coneSizeNew++] = firstSubPoint[d-1] + subc; 2046 } 2047 if (coneSizeNew != subconeSize) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of cone points located %d does not match subcone size %d", coneSizeNew, subconeSize); 2048 ierr = DMPlexSetCone(subdm, subpoint, coneNew);CHKERRQ(ierr); 2049 ierr = DMPlexSetConeOrientation(subdm, subpoint, ornt);CHKERRQ(ierr); 2050 } 2051 } 2052 ierr = PetscFree(coneNew);CHKERRQ(ierr); 2053 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 2054 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 2055 /* Build coordinates */ 2056 { 2057 PetscSection coordSection, subCoordSection; 2058 Vec coordinates, subCoordinates; 2059 PetscScalar *coords, *subCoords; 2060 PetscInt numComp, coordSize; 2061 2062 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2063 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 2064 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 2065 ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr); 2066 ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr); 2067 ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr); 2068 ierr = PetscSectionSetChart(subCoordSection, firstSubPoint[0], firstSubPoint[0]+numSubPoints[0]);CHKERRQ(ierr); 2069 for (v = 0; v < numSubPoints[0]; ++v) { 2070 const PetscInt vertex = subpoints[0][v]; 2071 const PetscInt subvertex = firstSubPoint[0]+v; 2072 PetscInt dof; 2073 2074 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2075 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 2076 ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr); 2077 } 2078 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 2079 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 2080 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 2081 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2082 ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr); 2083 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2084 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2085 for (v = 0; v < numSubPoints[0]; ++v) { 2086 const PetscInt vertex = subpoints[0][v]; 2087 const PetscInt subvertex = firstSubPoint[0]+v; 2088 PetscInt dof, off, sdof, soff, d; 2089 2090 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2091 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 2092 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 2093 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 2094 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 2095 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 2096 } 2097 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2098 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2099 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 2100 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 2101 } 2102 /* Cleanup */ 2103 for (d = 0; d <= dim; ++d) { 2104 ierr = ISRestoreIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 2105 ierr = ISDestroy(&subpointIS[d]);CHKERRQ(ierr); 2106 } 2107 ierr = PetscFree4(numSubPoints,firstSubPoint,subpointIS,subpoints);CHKERRQ(ierr); 2108 PetscFunctionReturn(0); 2109 } 2110 2111 #undef __FUNCT__ 2112 #define __FUNCT__ "DMPlexCreateSubmesh" 2113 /*@C 2114 DMPlexCreateSubmesh - Extract a hypersurface from the mesh using vertices defined by a label 2115 2116 Input Parameters: 2117 + dm - The original mesh 2118 . vertexLabel - The DMLabel marking vertices contained in the surface 2119 - value - The label value to use 2120 2121 Output Parameter: 2122 . subdm - The surface mesh 2123 2124 Note: This function produces a DMLabel mapping original points in the submesh to their depth. This can be obtained using DMPlexGetSubpointMap(). 2125 2126 Level: developer 2127 2128 .seealso: DMPlexGetSubpointMap(), DMPlexGetLabel(), DMLabelSetValue() 2129 @*/ 2130 PetscErrorCode DMPlexCreateSubmesh(DM dm, const char vertexLabel[], PetscInt value, DM *subdm) 2131 { 2132 PetscInt dim, depth; 2133 PetscErrorCode ierr; 2134 2135 PetscFunctionBegin; 2136 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2137 PetscValidCharPointer(vertexLabel, 2); 2138 PetscValidPointer(subdm, 3); 2139 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 2140 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 2141 ierr = DMCreate(PetscObjectComm((PetscObject)dm), subdm);CHKERRQ(ierr); 2142 ierr = DMSetType(*subdm, DMPLEX);CHKERRQ(ierr); 2143 ierr = DMPlexSetDimension(*subdm, dim-1);CHKERRQ(ierr); 2144 if (depth == dim) { 2145 ierr = DMPlexCreateSubmesh_Interpolated(dm, vertexLabel, value, *subdm);CHKERRQ(ierr); 2146 } else { 2147 ierr = DMPlexCreateSubmesh_Uninterpolated(dm, vertexLabel, value, *subdm);CHKERRQ(ierr); 2148 } 2149 PetscFunctionReturn(0); 2150 } 2151 2152 #undef __FUNCT__ 2153 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh_Uninterpolated" 2154 static PetscErrorCode DMPlexCreateCohesiveSubmesh_Uninterpolated(DM dm, PetscBool hasLagrange, DM subdm) 2155 { 2156 MPI_Comm comm; 2157 DMLabel subpointMap; 2158 IS subvertexIS; 2159 const PetscInt *subVertices; 2160 PetscInt numSubVertices, firstSubVertex, numSubCells, *subCells; 2161 PetscInt *subface, maxConeSize, numSubFaces, firstSubFace, newFacePoint, nFV; 2162 PetscInt cMax, c, f; 2163 PetscErrorCode ierr; 2164 2165 PetscFunctionBegin; 2166 ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 2167 /* Create subpointMap which marks the submesh */ 2168 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 2169 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 2170 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 2171 ierr = DMPlexMarkCohesiveSubmesh_Uninterpolated(dm, hasLagrange, subpointMap, &numSubFaces, &nFV, &subCells, subdm);CHKERRQ(ierr); 2172 /* Setup chart */ 2173 ierr = DMLabelGetStratumSize(subpointMap, 0, &numSubVertices);CHKERRQ(ierr); 2174 ierr = DMLabelGetStratumSize(subpointMap, 2, &numSubCells);CHKERRQ(ierr); 2175 ierr = DMPlexSetChart(subdm, 0, numSubCells+numSubFaces+numSubVertices);CHKERRQ(ierr); 2176 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 2177 /* Set cone sizes */ 2178 firstSubVertex = numSubCells; 2179 firstSubFace = numSubCells+numSubVertices; 2180 newFacePoint = firstSubFace; 2181 ierr = DMLabelGetStratumIS(subpointMap, 0, &subvertexIS);CHKERRQ(ierr); 2182 if (subvertexIS) {ierr = ISGetIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 2183 for (c = 0; c < numSubCells; ++c) { 2184 ierr = DMPlexSetConeSize(subdm, c, 1);CHKERRQ(ierr); 2185 } 2186 for (f = firstSubFace; f < firstSubFace+numSubFaces; ++f) { 2187 ierr = DMPlexSetConeSize(subdm, f, nFV);CHKERRQ(ierr); 2188 } 2189 ierr = DMSetUp(subdm);CHKERRQ(ierr); 2190 /* Create face cones */ 2191 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 2192 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 2193 ierr = DMGetWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 2194 for (c = 0; c < numSubCells; ++c) { 2195 const PetscInt cell = subCells[c]; 2196 const PetscInt subcell = c; 2197 const PetscInt *cone, *cells; 2198 PetscInt numCells, subVertex, p, v; 2199 2200 if (cell < cMax) continue; 2201 ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 2202 for (v = 0; v < nFV; ++v) { 2203 ierr = PetscFindInt(cone[v], numSubVertices, subVertices, &subVertex);CHKERRQ(ierr); 2204 subface[v] = firstSubVertex+subVertex; 2205 } 2206 ierr = DMPlexSetCone(subdm, newFacePoint, subface);CHKERRQ(ierr); 2207 ierr = DMPlexSetCone(subdm, subcell, &newFacePoint);CHKERRQ(ierr); 2208 ierr = DMPlexGetJoin(dm, nFV, cone, &numCells, &cells);CHKERRQ(ierr); 2209 if (numCells != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cells should separate two cells"); 2210 for (p = 0; p < numCells; ++p) { 2211 PetscInt negsubcell; 2212 2213 if (cells[p] >= cMax) continue; 2214 /* I know this is a crap search */ 2215 for (negsubcell = 0; negsubcell < numSubCells; ++negsubcell) { 2216 if (subCells[negsubcell] == cells[p]) break; 2217 } 2218 if (negsubcell == numSubCells) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find negative face neighbor for cohesive cell %d", cell); 2219 ierr = DMPlexSetCone(subdm, negsubcell, &newFacePoint);CHKERRQ(ierr); 2220 } 2221 ierr = DMPlexRestoreJoin(dm, nFV, cone, &numCells, &cells);CHKERRQ(ierr); 2222 ++newFacePoint; 2223 } 2224 ierr = DMRestoreWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 2225 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 2226 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 2227 /* Build coordinates */ 2228 { 2229 PetscSection coordSection, subCoordSection; 2230 Vec coordinates, subCoordinates; 2231 PetscScalar *coords, *subCoords; 2232 PetscInt numComp, coordSize, v; 2233 2234 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2235 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 2236 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 2237 ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr); 2238 ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr); 2239 ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr); 2240 ierr = PetscSectionSetChart(subCoordSection, firstSubVertex, firstSubVertex+numSubVertices);CHKERRQ(ierr); 2241 for (v = 0; v < numSubVertices; ++v) { 2242 const PetscInt vertex = subVertices[v]; 2243 const PetscInt subvertex = firstSubVertex+v; 2244 PetscInt dof; 2245 2246 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2247 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 2248 ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr); 2249 } 2250 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 2251 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 2252 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 2253 if (coordSize) { 2254 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2255 ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr); 2256 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2257 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2258 for (v = 0; v < numSubVertices; ++v) { 2259 const PetscInt vertex = subVertices[v]; 2260 const PetscInt subvertex = firstSubVertex+v; 2261 PetscInt dof, off, sdof, soff, d; 2262 2263 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2264 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 2265 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 2266 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 2267 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 2268 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 2269 } 2270 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2271 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2272 } 2273 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 2274 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 2275 } 2276 /* Cleanup */ 2277 if (subvertexIS) {ierr = ISRestoreIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 2278 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 2279 ierr = PetscFree(subCells);CHKERRQ(ierr); 2280 PetscFunctionReturn(0); 2281 } 2282 2283 #undef __FUNCT__ 2284 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh_Interpolated" 2285 static PetscErrorCode DMPlexCreateCohesiveSubmesh_Interpolated(DM dm, PetscBool hasLagrange, DM subdm) 2286 { 2287 MPI_Comm comm; 2288 DMLabel subpointMap; 2289 IS *subpointIS; 2290 const PetscInt **subpoints; 2291 PetscInt *numSubPoints, *firstSubPoint, *coneNew; 2292 PetscInt totSubPoints = 0, maxConeSize, dim, p, d, v; 2293 PetscErrorCode ierr; 2294 2295 PetscFunctionBegin; 2296 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2297 /* Create subpointMap which marks the submesh */ 2298 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 2299 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 2300 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 2301 ierr = DMPlexMarkCohesiveSubmesh_Interpolated(dm, hasLagrange, subpointMap, subdm);CHKERRQ(ierr); 2302 /* Setup chart */ 2303 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 2304 ierr = PetscMalloc4(dim+1,PetscInt,&numSubPoints,dim+1,PetscInt,&firstSubPoint,dim+1,IS,&subpointIS,dim+1,const PetscInt *,&subpoints);CHKERRQ(ierr); 2305 for (d = 0; d <= dim; ++d) { 2306 ierr = DMLabelGetStratumSize(subpointMap, d, &numSubPoints[d]);CHKERRQ(ierr); 2307 totSubPoints += numSubPoints[d]; 2308 } 2309 ierr = DMPlexSetChart(subdm, 0, totSubPoints);CHKERRQ(ierr); 2310 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 2311 /* Set cone sizes */ 2312 firstSubPoint[dim] = 0; 2313 firstSubPoint[0] = firstSubPoint[dim] + numSubPoints[dim]; 2314 if (dim > 1) {firstSubPoint[dim-1] = firstSubPoint[0] + numSubPoints[0];} 2315 if (dim > 2) {firstSubPoint[dim-2] = firstSubPoint[dim-1] + numSubPoints[dim-1];} 2316 for (d = 0; d <= dim; ++d) { 2317 ierr = DMLabelGetStratumIS(subpointMap, d, &subpointIS[d]);CHKERRQ(ierr); 2318 ierr = ISGetIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 2319 } 2320 for (d = 0; d <= dim; ++d) { 2321 for (p = 0; p < numSubPoints[d]; ++p) { 2322 const PetscInt point = subpoints[d][p]; 2323 const PetscInt subpoint = firstSubPoint[d] + p; 2324 const PetscInt *cone; 2325 PetscInt coneSize, coneSizeNew, c, val; 2326 2327 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 2328 ierr = DMPlexSetConeSize(subdm, subpoint, coneSize);CHKERRQ(ierr); 2329 if (d == dim) { 2330 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 2331 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 2332 ierr = DMLabelGetValue(subpointMap, cone[c], &val);CHKERRQ(ierr); 2333 if (val >= 0) coneSizeNew++; 2334 } 2335 ierr = DMPlexSetConeSize(subdm, subpoint, coneSizeNew);CHKERRQ(ierr); 2336 } 2337 } 2338 } 2339 ierr = DMSetUp(subdm);CHKERRQ(ierr); 2340 /* Set cones */ 2341 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 2342 ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &coneNew);CHKERRQ(ierr); 2343 for (d = 0; d <= dim; ++d) { 2344 for (p = 0; p < numSubPoints[d]; ++p) { 2345 const PetscInt point = subpoints[d][p]; 2346 const PetscInt subpoint = firstSubPoint[d] + p; 2347 const PetscInt *cone; 2348 PetscInt coneSize, subconeSize, coneSizeNew, c, subc; 2349 2350 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 2351 ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr); 2352 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 2353 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 2354 ierr = PetscFindInt(cone[c], numSubPoints[d-1], subpoints[d-1], &subc);CHKERRQ(ierr); 2355 if (subc >= 0) coneNew[coneSizeNew++] = firstSubPoint[d-1] + subc; 2356 } 2357 if (coneSizeNew != subconeSize) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of cone points located %d does not match subcone size %d", coneSizeNew, subconeSize); 2358 ierr = DMPlexSetCone(subdm, subpoint, coneNew);CHKERRQ(ierr); 2359 } 2360 } 2361 ierr = PetscFree(coneNew);CHKERRQ(ierr); 2362 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 2363 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 2364 /* Build coordinates */ 2365 { 2366 PetscSection coordSection, subCoordSection; 2367 Vec coordinates, subCoordinates; 2368 PetscScalar *coords, *subCoords; 2369 PetscInt numComp, coordSize; 2370 2371 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2372 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 2373 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 2374 ierr = PetscSectionSetNumFields(subCoordSection, 1);CHKERRQ(ierr); 2375 ierr = PetscSectionGetFieldComponents(coordSection, 0, &numComp);CHKERRQ(ierr); 2376 ierr = PetscSectionSetFieldComponents(subCoordSection, 0, numComp);CHKERRQ(ierr); 2377 ierr = PetscSectionSetChart(subCoordSection, firstSubPoint[0], firstSubPoint[0]+numSubPoints[0]);CHKERRQ(ierr); 2378 for (v = 0; v < numSubPoints[0]; ++v) { 2379 const PetscInt vertex = subpoints[0][v]; 2380 const PetscInt subvertex = firstSubPoint[0]+v; 2381 PetscInt dof; 2382 2383 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2384 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 2385 ierr = PetscSectionSetFieldDof(subCoordSection, subvertex, 0, dof);CHKERRQ(ierr); 2386 } 2387 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 2388 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 2389 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 2390 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2391 ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr); 2392 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2393 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2394 for (v = 0; v < numSubPoints[0]; ++v) { 2395 const PetscInt vertex = subpoints[0][v]; 2396 const PetscInt subvertex = firstSubPoint[0]+v; 2397 PetscInt dof, off, sdof, soff, d; 2398 2399 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 2400 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 2401 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 2402 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 2403 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 2404 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 2405 } 2406 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2407 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 2408 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 2409 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 2410 } 2411 /* Cleanup */ 2412 for (d = 0; d <= dim; ++d) { 2413 ierr = ISRestoreIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 2414 ierr = ISDestroy(&subpointIS[d]);CHKERRQ(ierr); 2415 } 2416 ierr = PetscFree4(numSubPoints,firstSubPoint,subpointIS,subpoints);CHKERRQ(ierr); 2417 PetscFunctionReturn(0); 2418 } 2419 2420 #undef __FUNCT__ 2421 #define __FUNCT__ "DMPlexCreateCohesiveSubmesh" 2422 /* 2423 DMPlexCreateCohesiveSubmesh - Extract from a mesh with cohesive cells the hypersurface defined by one face of the cells. 2424 2425 Input Parameters: 2426 + dm - The original mesh 2427 - hasLagrange - The mesh has Lagrange unknowns in the cohesive cells 2428 2429 Output Parameter: 2430 . subdm - The surface mesh 2431 2432 Note: This function produces a DMLabel mapping original points in the submesh to their depth. This can be obtained using DMPlexGetSubpointMap(). 2433 2434 Level: developer 2435 2436 .seealso: DMPlexGetSubpointMap(), DMPlexCreateSubmesh() 2437 */ 2438 PetscErrorCode DMPlexCreateCohesiveSubmesh(DM dm, PetscBool hasLagrange, DM *subdm) 2439 { 2440 PetscInt dim, depth; 2441 PetscErrorCode ierr; 2442 2443 PetscFunctionBegin; 2444 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2445 PetscValidPointer(subdm, 3); 2446 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 2447 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 2448 ierr = DMCreate(PetscObjectComm((PetscObject)dm), subdm);CHKERRQ(ierr); 2449 ierr = DMSetType(*subdm, DMPLEX);CHKERRQ(ierr); 2450 ierr = DMPlexSetDimension(*subdm, dim-1);CHKERRQ(ierr); 2451 if (depth == dim) { 2452 ierr = DMPlexCreateCohesiveSubmesh_Interpolated(dm, hasLagrange, *subdm);CHKERRQ(ierr); 2453 } else { 2454 ierr = DMPlexCreateCohesiveSubmesh_Uninterpolated(dm, hasLagrange, *subdm);CHKERRQ(ierr); 2455 } 2456 PetscFunctionReturn(0); 2457 } 2458 2459 #undef __FUNCT__ 2460 #define __FUNCT__ "DMPlexGetSubpointMap" 2461 /*@ 2462 DMPlexGetSubpointMap - Returns a DMLabel with point dimension as values 2463 2464 Input Parameter: 2465 . dm - The submesh DM 2466 2467 Output Parameter: 2468 . subpointMap - The DMLabel of all the points from the original mesh in this submesh, or NULL if this is not a submesh 2469 2470 Level: developer 2471 2472 .seealso: DMPlexCreateSubmesh(), DMPlexCreateSubpointIS() 2473 @*/ 2474 PetscErrorCode DMPlexGetSubpointMap(DM dm, DMLabel *subpointMap) 2475 { 2476 DM_Plex *mesh = (DM_Plex*) dm->data; 2477 2478 PetscFunctionBegin; 2479 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2480 PetscValidPointer(subpointMap, 2); 2481 *subpointMap = mesh->subpointMap; 2482 PetscFunctionReturn(0); 2483 } 2484 2485 #undef __FUNCT__ 2486 #define __FUNCT__ "DMPlexSetSubpointMap" 2487 /* Note: Should normally not be called by the user, since it is set in DMPlexCreateSubmesh() */ 2488 PetscErrorCode DMPlexSetSubpointMap(DM dm, DMLabel subpointMap) 2489 { 2490 DM_Plex *mesh = (DM_Plex *) dm->data; 2491 DMLabel tmp; 2492 PetscErrorCode ierr; 2493 2494 PetscFunctionBegin; 2495 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2496 tmp = mesh->subpointMap; 2497 mesh->subpointMap = subpointMap; 2498 ++mesh->subpointMap->refct; 2499 ierr = DMLabelDestroy(&tmp);CHKERRQ(ierr); 2500 PetscFunctionReturn(0); 2501 } 2502 2503 #undef __FUNCT__ 2504 #define __FUNCT__ "DMPlexCreateSubpointIS" 2505 /*@ 2506 DMPlexCreateSubpointIS - Creates an IS covering the entire subdm chart with the original points as data 2507 2508 Input Parameter: 2509 . dm - The submesh DM 2510 2511 Output Parameter: 2512 . subpointIS - The IS of all the points from the original mesh in this submesh, or NULL if this is not a submesh 2513 2514 Note: This is IS is guaranteed to be sorted by the construction of the submesh 2515 2516 Level: developer 2517 2518 .seealso: DMPlexCreateSubmesh(), DMPlexGetSubpointMap() 2519 @*/ 2520 PetscErrorCode DMPlexCreateSubpointIS(DM dm, IS *subpointIS) 2521 { 2522 MPI_Comm comm; 2523 DMLabel subpointMap; 2524 IS is; 2525 const PetscInt *opoints; 2526 PetscInt *points, *depths; 2527 PetscInt depth, depStart, depEnd, d, pStart, pEnd, p, n, off; 2528 PetscErrorCode ierr; 2529 2530 PetscFunctionBegin; 2531 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2532 PetscValidPointer(subpointIS, 2); 2533 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2534 *subpointIS = NULL; 2535 ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 2536 if (subpointMap) { 2537 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 2538 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 2539 if (pStart) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Submeshes must start the point numbering at 0, not %d", pStart); 2540 ierr = DMGetWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr); 2541 depths[0] = depth; 2542 depths[1] = 0; 2543 for(d = 2; d <= depth; ++d) {depths[d] = depth+1 - d;} 2544 ierr = PetscMalloc(pEnd * sizeof(PetscInt), &points);CHKERRQ(ierr); 2545 for(d = 0, off = 0; d <= depth; ++d) { 2546 const PetscInt dep = depths[d]; 2547 2548 ierr = DMPlexGetDepthStratum(dm, dep, &depStart, &depEnd);CHKERRQ(ierr); 2549 ierr = DMLabelGetStratumSize(subpointMap, dep, &n);CHKERRQ(ierr); 2550 if (((d < 2) && (depth > 1)) || (d == 1)) { /* Only check vertices and cells for now since the map is broken for others */ 2551 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); 2552 } else { 2553 if (!n) { 2554 if (d == 0) { 2555 /* Missing cells */ 2556 for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = -1; 2557 } else { 2558 /* Missing faces */ 2559 for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = PETSC_MAX_INT; 2560 } 2561 } 2562 } 2563 if (n) { 2564 ierr = DMLabelGetStratumIS(subpointMap, dep, &is);CHKERRQ(ierr); 2565 ierr = ISGetIndices(is, &opoints);CHKERRQ(ierr); 2566 for(p = 0; p < n; ++p, ++off) points[off] = opoints[p]; 2567 ierr = ISRestoreIndices(is, &opoints);CHKERRQ(ierr); 2568 ierr = ISDestroy(&is);CHKERRQ(ierr); 2569 } 2570 } 2571 ierr = DMRestoreWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr); 2572 if (off != pEnd) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "The number of mapped submesh points %d should be %d", off, pEnd); 2573 ierr = ISCreateGeneral(comm, pEnd, points, PETSC_OWN_POINTER, subpointIS);CHKERRQ(ierr); 2574 } 2575 PetscFunctionReturn(0); 2576 } 2577