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