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