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