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