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, *supportNew; 544 PetscInt shift = 100, depth = 0, dep, dim, d, numSP = 0, sp, maxConeSize, maxSupportSize, numLabels, 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 = PetscMalloc5(depth+1,PetscInt,&depthShift,depth+1,PetscInt,&depthOffset,depth+1,PetscInt,&pMaxNew,maxConeSize*3,PetscInt,&coneNew,maxSupportSize,PetscInt,&supportNew);CHKERRQ(ierr); 559 ierr = PetscMalloc3(depth+1,IS,&pointIS,depth+1,PetscInt,&numSplitPoints,depth+1,const PetscInt*,&splitPoints);CHKERRQ(ierr); 560 ierr = PetscMemzero(depthShift, (depth+1) * sizeof(PetscInt));CHKERRQ(ierr); 561 for (d = 0; d <= depth; ++d) { 562 ierr = DMPlexGetDepthStratum(dm, d, NULL, &pMaxNew[d]);CHKERRQ(ierr); 563 numSplitPoints[d] = 0; 564 splitPoints[d] = NULL; 565 pointIS[d] = NULL; 566 } 567 for (sp = 0; sp < numSP; ++sp) { 568 const PetscInt dep = values[sp]; 569 570 if ((dep < 0) || (dep > depth)) continue; 571 ierr = DMLabelGetStratumSize(label, dep, &depthShift[dep]);CHKERRQ(ierr); 572 ierr = DMLabelGetStratumIS(label, dep, &pointIS[dep]);CHKERRQ(ierr); 573 if (pointIS[dep]) { 574 ierr = ISGetLocalSize(pointIS[dep], &numSplitPoints[dep]);CHKERRQ(ierr); 575 ierr = ISGetIndices(pointIS[dep], &splitPoints[dep]);CHKERRQ(ierr); 576 } 577 } 578 if (depth >= 0) { 579 /* Calculate number of additional points */ 580 depthShift[depth] = depthShift[depth-1]; /* There is a cohesive cell for every split face */ 581 depthShift[1] += depthShift[0]; /* There is a cohesive edge for every split vertex */ 582 /* Calculate hybrid bound for each dimension */ 583 pMaxNew[0] += depthShift[depth]; 584 if (depth > 1) pMaxNew[dim-1] += depthShift[depth] + depthShift[0]; 585 if (depth > 2) pMaxNew[1] += depthShift[depth] + depthShift[0] + depthShift[dim-1]; 586 587 /* Calculate point offset for each dimension */ 588 depthOffset[depth] = 0; 589 depthOffset[0] = depthOffset[depth] + depthShift[depth]; 590 if (depth > 1) depthOffset[dim-1] = depthOffset[0] + depthShift[0]; 591 if (depth > 2) depthOffset[1] = depthOffset[dim-1] + depthShift[dim-1]; 592 } 593 ierr = DMPlexShiftSizes_Internal(dm, depthShift, sdm);CHKERRQ(ierr); 594 /* Step 3: Set cone/support sizes for new points */ 595 for (dep = 0; dep <= depth; ++dep) { 596 for (p = 0; p < numSplitPoints[dep]; ++p) { 597 const PetscInt oldp = splitPoints[dep][p]; 598 const PetscInt newp = depthOffset[dep] + oldp; 599 const PetscInt splitp = pMaxNew[dep] + p; 600 const PetscInt *support; 601 PetscInt coneSize, supportSize, q, e; 602 603 ierr = DMPlexGetConeSize(dm, oldp, &coneSize);CHKERRQ(ierr); 604 ierr = DMPlexSetConeSize(sdm, splitp, coneSize);CHKERRQ(ierr); 605 ierr = DMPlexGetSupportSize(dm, oldp, &supportSize);CHKERRQ(ierr); 606 ierr = DMPlexSetSupportSize(sdm, splitp, supportSize);CHKERRQ(ierr); 607 if (dep == depth-1) { 608 const PetscInt ccell = pMaxNew[depth] + p; 609 /* Add cohesive cells, they are prisms */ 610 ierr = DMPlexSetConeSize(sdm, ccell, 2 + coneSize);CHKERRQ(ierr); 611 } else if (dep == 0) { 612 const PetscInt cedge = pMaxNew[1] + (depthShift[1] - depthShift[0]) + p; 613 614 ierr = DMPlexGetSupport(dm, oldp, &support);CHKERRQ(ierr); 615 /* Split old vertex: Edges in old split faces and new cohesive edge */ 616 for (e = 0, q = 0; e < supportSize; ++e) { 617 PetscInt val; 618 619 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 620 if ((val == 1) || (val == (shift + 1))) ++q; 621 } 622 ierr = DMPlexSetSupportSize(sdm, newp, q+1);CHKERRQ(ierr); 623 /* Split new vertex: Edges in new split faces and new cohesive edge */ 624 for (e = 0, q = 0; e < supportSize; ++e) { 625 PetscInt val; 626 627 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 628 if ((val == 1) || (val == -(shift + 1))) ++q; 629 } 630 ierr = DMPlexSetSupportSize(sdm, splitp, q+1);CHKERRQ(ierr); 631 /* Add cohesive edges */ 632 ierr = DMPlexSetConeSize(sdm, cedge, 2);CHKERRQ(ierr); 633 /* Punt for now on support, you loop over closure, extract faces, check which ones are in the label */ 634 } else if (dep == dim-2) { 635 ierr = DMPlexGetSupport(dm, oldp, &support);CHKERRQ(ierr); 636 /* Split old edge: Faces in positive side cells and old split faces */ 637 for (e = 0, q = 0; e < supportSize; ++e) { 638 PetscInt val; 639 640 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 641 if ((val == dim-1) || (val == (shift + dim-1))) ++q; 642 } 643 ierr = DMPlexSetSupportSize(sdm, newp, q);CHKERRQ(ierr); 644 /* Split new edge: Faces in negative side cells and new split faces */ 645 for (e = 0, q = 0; e < supportSize; ++e) { 646 PetscInt val; 647 648 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 649 if ((val == dim-1) || (val == -(shift + dim-1))) ++q; 650 } 651 ierr = DMPlexSetSupportSize(sdm, splitp, q);CHKERRQ(ierr); 652 } 653 } 654 } 655 /* Step 4: Setup split DM */ 656 ierr = DMSetUp(sdm);CHKERRQ(ierr); 657 ierr = DMPlexShiftPoints_Internal(dm, depthShift, sdm);CHKERRQ(ierr); 658 /* Step 6: Set cones and supports for new points */ 659 for (dep = 0; dep <= depth; ++dep) { 660 for (p = 0; p < numSplitPoints[dep]; ++p) { 661 const PetscInt oldp = splitPoints[dep][p]; 662 const PetscInt newp = depthOffset[dep] + oldp; 663 const PetscInt splitp = pMaxNew[dep] + p; 664 const PetscInt *cone, *support, *ornt; 665 PetscInt coneSize, supportSize, q, v, e, s; 666 667 ierr = DMPlexGetConeSize(dm, oldp, &coneSize);CHKERRQ(ierr); 668 ierr = DMPlexGetCone(dm, oldp, &cone);CHKERRQ(ierr); 669 ierr = DMPlexGetConeOrientation(dm, oldp, &ornt);CHKERRQ(ierr); 670 ierr = DMPlexGetSupportSize(dm, oldp, &supportSize);CHKERRQ(ierr); 671 ierr = DMPlexGetSupport(dm, oldp, &support);CHKERRQ(ierr); 672 if (dep == depth-1) { 673 const PetscInt ccell = pMaxNew[depth] + p; 674 const PetscInt *supportF; 675 676 /* Split face: copy in old face to new face to start */ 677 ierr = DMPlexGetSupport(sdm, newp, &supportF);CHKERRQ(ierr); 678 ierr = DMPlexSetSupport(sdm, splitp, supportF);CHKERRQ(ierr); 679 /* Split old face: old vertices/edges in cone so no change */ 680 /* Split new face: new vertices/edges in cone */ 681 for (q = 0; q < coneSize; ++q) { 682 ierr = PetscFindInt(cone[q], numSplitPoints[dim-2], splitPoints[dim-2], &v);CHKERRQ(ierr); 683 684 coneNew[2+q] = pMaxNew[dim-2] + v; 685 } 686 ierr = DMPlexSetCone(sdm, splitp, &coneNew[2]);CHKERRQ(ierr); 687 ierr = DMPlexSetConeOrientation(sdm, splitp, ornt);CHKERRQ(ierr); 688 /* Cohesive cell: Old and new split face, then new cohesive edges */ 689 coneNew[0] = newp; 690 coneNew[1] = splitp; 691 for (q = 0; q < coneSize; ++q) { 692 coneNew[2+q] = (pMaxNew[1] - pMaxNew[dim-2]) + (depthShift[1] - depthShift[0]) + coneNew[2+q]; 693 } 694 ierr = DMPlexSetCone(sdm, ccell, coneNew);CHKERRQ(ierr); 695 696 697 for (s = 0; s < supportSize; ++s) { 698 PetscInt val; 699 700 ierr = DMLabelGetValue(label, support[s], &val);CHKERRQ(ierr); 701 if (val < 0) { 702 /* Split old face: Replace negative side cell with cohesive cell */ 703 ierr = DMPlexInsertSupport(sdm, newp, s, ccell);CHKERRQ(ierr); 704 } else { 705 /* Split new face: Replace positive side cell with cohesive cell */ 706 ierr = DMPlexInsertSupport(sdm, splitp, s, ccell);CHKERRQ(ierr); 707 } 708 } 709 } else if (dep == 0) { 710 const PetscInt cedge = pMaxNew[1] + (depthShift[1] - depthShift[0]) + p; 711 712 /* Split old vertex: Edges in old split faces and new cohesive edge */ 713 for (e = 0, q = 0; e < supportSize; ++e) { 714 PetscInt val; 715 716 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 717 if ((val == 1) || (val == (shift + 1))) { 718 supportNew[q++] = depthOffset[1] + support[e]; 719 } 720 } 721 supportNew[q] = cedge; 722 723 ierr = DMPlexSetSupport(sdm, newp, supportNew);CHKERRQ(ierr); 724 /* Split new vertex: Edges in new split faces and new cohesive edge */ 725 for (e = 0, q = 0; e < supportSize; ++e) { 726 PetscInt val, edge; 727 728 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 729 if (val == 1) { 730 ierr = PetscFindInt(support[e], numSplitPoints[1], splitPoints[1], &edge);CHKERRQ(ierr); 731 if (edge < 0) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Edge %d is not a split edge", support[e]); 732 supportNew[q++] = pMaxNew[1] + edge; 733 } else if (val == -(shift + 1)) { 734 supportNew[q++] = depthOffset[1] + support[e]; 735 } 736 } 737 supportNew[q] = cedge; 738 ierr = DMPlexSetSupport(sdm, splitp, supportNew);CHKERRQ(ierr); 739 /* Cohesive edge: Old and new split vertex, punting on support */ 740 coneNew[0] = newp; 741 coneNew[1] = splitp; 742 ierr = DMPlexSetCone(sdm, cedge, coneNew);CHKERRQ(ierr); 743 } else if (dep == dim-2) { 744 /* Split old edge: old vertices in cone so no change */ 745 /* Split new edge: new vertices in cone */ 746 for (q = 0; q < coneSize; ++q) { 747 ierr = PetscFindInt(cone[q], numSplitPoints[dim-3], splitPoints[dim-3], &v);CHKERRQ(ierr); 748 749 coneNew[q] = pMaxNew[dim-3] + v; 750 } 751 ierr = DMPlexSetCone(sdm, splitp, coneNew);CHKERRQ(ierr); 752 /* Split old edge: Faces in positive side cells and old split faces */ 753 for (e = 0, q = 0; e < supportSize; ++e) { 754 PetscInt val; 755 756 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 757 if ((val == dim-1) || (val == (shift + dim-1))) { 758 supportNew[q++] = depthOffset[dim-1] + support[e]; 759 } 760 } 761 ierr = DMPlexSetSupport(sdm, newp, supportNew);CHKERRQ(ierr); 762 /* Split new edge: Faces in negative side cells and new split faces */ 763 for (e = 0, q = 0; e < supportSize; ++e) { 764 PetscInt val, face; 765 766 ierr = DMLabelGetValue(label, support[e], &val);CHKERRQ(ierr); 767 if (val == dim-1) { 768 ierr = PetscFindInt(support[e], numSplitPoints[dim-1], splitPoints[dim-1], &face);CHKERRQ(ierr); 769 if (face < 0) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Face %d is not a split face", support[e]); 770 supportNew[q++] = pMaxNew[dim-1] + face; 771 } else if (val == -(shift + dim-1)) { 772 supportNew[q++] = depthOffset[dim-1] + support[e]; 773 } 774 } 775 ierr = DMPlexSetSupport(sdm, splitp, supportNew);CHKERRQ(ierr); 776 } 777 } 778 } 779 /* Step 6b: Replace split points in negative side cones */ 780 for (sp = 0; sp < numSP; ++sp) { 781 PetscInt dep = values[sp]; 782 IS pIS; 783 PetscInt numPoints; 784 const PetscInt *points; 785 786 if (dep >= 0) continue; 787 ierr = DMLabelGetStratumIS(label, dep, &pIS);CHKERRQ(ierr); 788 if (!pIS) continue; 789 dep = -dep - shift; 790 ierr = ISGetLocalSize(pIS, &numPoints);CHKERRQ(ierr); 791 ierr = ISGetIndices(pIS, &points);CHKERRQ(ierr); 792 for (p = 0; p < numPoints; ++p) { 793 const PetscInt oldp = points[p]; 794 const PetscInt newp = depthOffset[dep] + oldp; 795 const PetscInt *cone; 796 PetscInt coneSize, c; 797 PetscBool replaced = PETSC_FALSE; 798 799 /* Negative edge: replace split vertex */ 800 /* Negative cell: replace split face */ 801 ierr = DMPlexGetConeSize(sdm, newp, &coneSize);CHKERRQ(ierr); 802 ierr = DMPlexGetCone(sdm, newp, &cone);CHKERRQ(ierr); 803 for (c = 0; c < coneSize; ++c) { 804 const PetscInt coldp = cone[c] - depthOffset[dep-1]; 805 PetscInt csplitp, cp, val; 806 807 ierr = DMLabelGetValue(label, coldp, &val);CHKERRQ(ierr); 808 if (val == dep-1) { 809 ierr = PetscFindInt(coldp, numSplitPoints[dep-1], splitPoints[dep-1], &cp);CHKERRQ(ierr); 810 if (cp < 0) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "Point %d is not a split point of dimension %d", oldp, dep-1); 811 csplitp = pMaxNew[dep-1] + cp; 812 ierr = DMPlexInsertCone(sdm, newp, c, csplitp);CHKERRQ(ierr); 813 replaced = PETSC_TRUE; 814 } 815 } 816 if (!replaced) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "The cone of point %d does not contain split points", oldp); 817 } 818 ierr = ISRestoreIndices(pIS, &points);CHKERRQ(ierr); 819 ierr = ISDestroy(&pIS);CHKERRQ(ierr); 820 } 821 /* Step 7: Stratify */ 822 ierr = DMPlexStratify(sdm);CHKERRQ(ierr); 823 /* Step 8: Coordinates */ 824 ierr = DMPlexShiftCoordinates_Internal(dm, depthShift, sdm);CHKERRQ(ierr); 825 ierr = DMPlexGetCoordinateSection(sdm, &coordSection);CHKERRQ(ierr); 826 ierr = DMGetCoordinatesLocal(sdm, &coordinates);CHKERRQ(ierr); 827 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 828 for (v = 0; v < (numSplitPoints ? numSplitPoints[0] : 0); ++v) { 829 const PetscInt newp = depthOffset[0] + splitPoints[0][v]; 830 const PetscInt splitp = pMaxNew[0] + v; 831 PetscInt dof, off, soff, d; 832 833 ierr = PetscSectionGetDof(coordSection, newp, &dof);CHKERRQ(ierr); 834 ierr = PetscSectionGetOffset(coordSection, newp, &off);CHKERRQ(ierr); 835 ierr = PetscSectionGetOffset(coordSection, splitp, &soff);CHKERRQ(ierr); 836 for (d = 0; d < dof; ++d) coords[soff+d] = coords[off+d]; 837 } 838 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 839 /* Step 9: SF, if I can figure this out we can split the mesh in parallel */ 840 ierr = DMPlexShiftSF_Internal(dm, depthShift, sdm);CHKERRQ(ierr); 841 /* Step 10: Labels */ 842 ierr = DMPlexShiftLabels_Internal(dm, depthShift, sdm);CHKERRQ(ierr); 843 ierr = DMPlexGetNumLabels(sdm, &numLabels);CHKERRQ(ierr); 844 for (dep = 0; dep <= depth; ++dep) { 845 for (p = 0; p < numSplitPoints[dep]; ++p) { 846 const PetscInt newp = depthOffset[dep] + splitPoints[dep][p]; 847 const PetscInt splitp = pMaxNew[dep] + p; 848 PetscInt l; 849 850 for (l = 0; l < numLabels; ++l) { 851 DMLabel mlabel; 852 const char *lname; 853 PetscInt val; 854 855 ierr = DMPlexGetLabelName(sdm, l, &lname);CHKERRQ(ierr); 856 ierr = DMPlexGetLabel(sdm, lname, &mlabel);CHKERRQ(ierr); 857 ierr = DMLabelGetValue(mlabel, newp, &val);CHKERRQ(ierr); 858 if (val >= 0) { 859 ierr = DMLabelSetValue(mlabel, splitp, val);CHKERRQ(ierr); 860 if (dep == 0) { 861 const PetscInt cedge = pMaxNew[1] + (depthShift[1] - depthShift[0]) + p; 862 ierr = DMLabelSetValue(mlabel, cedge, val);CHKERRQ(ierr); 863 } 864 } 865 } 866 } 867 } 868 for (sp = 0; sp < numSP; ++sp) { 869 const PetscInt dep = values[sp]; 870 871 if ((dep < 0) || (dep > depth)) continue; 872 if (pointIS[dep]) {ierr = ISRestoreIndices(pointIS[dep], &splitPoints[dep]);CHKERRQ(ierr);} 873 ierr = ISDestroy(&pointIS[dep]);CHKERRQ(ierr); 874 } 875 if (label) { 876 ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr); 877 ierr = ISDestroy(&valueIS);CHKERRQ(ierr); 878 } 879 ierr = PetscFree5(depthShift, depthOffset, pMaxNew, coneNew, supportNew);CHKERRQ(ierr); 880 ierr = PetscFree3(pointIS, numSplitPoints, splitPoints);CHKERRQ(ierr); 881 PetscFunctionReturn(0); 882 } 883 884 #undef __FUNCT__ 885 #define __FUNCT__ "DMPlexConstructCohesiveCells" 886 /*@C 887 DMPlexConstructCohesiveCells - Construct cohesive cells which split the face along an internal interface 888 889 Collective on dm 890 891 Input Parameters: 892 + dm - The original DM 893 - labelName - The label specifying the boundary faces (this could be auto-generated) 894 895 Output Parameters: 896 - dmSplit - The new DM 897 898 Level: developer 899 900 .seealso: DMCreate(), DMPlexLabelCohesiveComplete() 901 @*/ 902 PetscErrorCode DMPlexConstructCohesiveCells(DM dm, DMLabel label, DM *dmSplit) 903 { 904 DM sdm; 905 PetscInt dim; 906 PetscErrorCode ierr; 907 908 PetscFunctionBegin; 909 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 910 PetscValidPointer(dmSplit, 4); 911 ierr = DMCreate(PetscObjectComm((PetscObject)dm), &sdm);CHKERRQ(ierr); 912 ierr = DMSetType(sdm, DMPLEX);CHKERRQ(ierr); 913 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 914 ierr = DMPlexSetDimension(sdm, dim);CHKERRQ(ierr); 915 switch (dim) { 916 case 2: 917 case 3: 918 ierr = DMPlexConstructCohesiveCells_Internal(dm, label, sdm);CHKERRQ(ierr); 919 break; 920 default: 921 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct cohesive cells for dimension %d", dim); 922 } 923 *dmSplit = sdm; 924 PetscFunctionReturn(0); 925 } 926 927 #undef __FUNCT__ 928 #define __FUNCT__ "DMPlexLabelCohesiveComplete" 929 /*@ 930 DMPlexLabelCohesiveComplete - Starting with a label marking vertices on an internal surface, we add all other mesh pieces 931 to complete the surface 932 933 Input Parameters: 934 + dm - The DM 935 - label - A DMLabel marking the surface vertices 936 937 Output Parameter: 938 . label - A DMLabel marking all surface points 939 940 Level: developer 941 942 .seealso: DMPlexConstructCohesiveCells(), DMPlexLabelComplete() 943 @*/ 944 PetscErrorCode DMPlexLabelCohesiveComplete(DM dm, DMLabel label) 945 { 946 IS dimIS; 947 const PetscInt *points; 948 PetscInt shift = 100, dim, dep, cStart, cEnd, numPoints, p, val; 949 PetscErrorCode ierr; 950 951 PetscFunctionBegin; 952 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 953 /* Cell orientation for face gives the side of the fault */ 954 ierr = DMLabelGetStratumIS(label, dim-1, &dimIS);CHKERRQ(ierr); 955 if (!dimIS) PetscFunctionReturn(0); 956 ierr = ISGetLocalSize(dimIS, &numPoints);CHKERRQ(ierr); 957 ierr = ISGetIndices(dimIS, &points);CHKERRQ(ierr); 958 for (p = 0; p < numPoints; ++p) { 959 const PetscInt *support; 960 PetscInt supportSize, s; 961 962 ierr = DMPlexGetSupportSize(dm, points[p], &supportSize);CHKERRQ(ierr); 963 if (supportSize != 2) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Split face %d has %d != 2 supports", points[p], supportSize); 964 ierr = DMPlexGetSupport(dm, points[p], &support);CHKERRQ(ierr); 965 for (s = 0; s < supportSize; ++s) { 966 const PetscInt *cone, *ornt; 967 PetscInt coneSize, c; 968 PetscBool pos = PETSC_TRUE; 969 970 ierr = DMPlexGetConeSize(dm, support[s], &coneSize);CHKERRQ(ierr); 971 ierr = DMPlexGetCone(dm, support[s], &cone);CHKERRQ(ierr); 972 ierr = DMPlexGetConeOrientation(dm, support[s], &ornt);CHKERRQ(ierr); 973 for (c = 0; c < coneSize; ++c) { 974 if (cone[c] == points[p]) { 975 if (ornt[c] >= 0) { 976 ierr = DMLabelSetValue(label, support[s], shift+dim);CHKERRQ(ierr); 977 } else { 978 ierr = DMLabelSetValue(label, support[s], -(shift+dim));CHKERRQ(ierr); 979 pos = PETSC_FALSE; 980 } 981 break; 982 } 983 } 984 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]); 985 /* Put faces touching the fault in the label */ 986 for (c = 0; c < coneSize; ++c) { 987 const PetscInt point = cone[c]; 988 989 ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr); 990 if (val == -1) { 991 PetscInt *closure = NULL; 992 PetscInt closureSize, cl; 993 994 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 995 for (cl = 0; cl < closureSize*2; cl += 2) { 996 const PetscInt clp = closure[cl]; 997 998 ierr = DMLabelGetValue(label, clp, &val);CHKERRQ(ierr); 999 if ((val >= 0) && (val < dim-1)) { 1000 ierr = DMLabelSetValue(label, point, pos == PETSC_TRUE ? shift+dim-1 : -(shift+dim-1));CHKERRQ(ierr); 1001 break; 1002 } 1003 } 1004 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1005 } 1006 } 1007 } 1008 } 1009 ierr = ISRestoreIndices(dimIS, &points);CHKERRQ(ierr); 1010 ierr = ISDestroy(&dimIS);CHKERRQ(ierr); 1011 /* Search for other cells/faces/edges connected to the fault by a vertex */ 1012 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1013 ierr = DMLabelGetStratumIS(label, 0, &dimIS);CHKERRQ(ierr); 1014 if (!dimIS) PetscFunctionReturn(0); 1015 ierr = ISGetLocalSize(dimIS, &numPoints);CHKERRQ(ierr); 1016 ierr = ISGetIndices(dimIS, &points);CHKERRQ(ierr); 1017 for (p = 0; p < numPoints; ++p) { 1018 PetscInt *star = NULL; 1019 PetscInt starSize, s; 1020 PetscInt again = 1; /* 0: Finished 1: Keep iterating after a change 2: No change */ 1021 1022 /* First mark cells connected to the fault */ 1023 ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1024 while (again) { 1025 if (again > 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Could not classify all cells connected to the fault"); 1026 again = 0; 1027 for (s = 0; s < starSize*2; s += 2) { 1028 const PetscInt point = star[s]; 1029 const PetscInt *cone; 1030 PetscInt coneSize, c; 1031 1032 if ((point < cStart) || (point >= cEnd)) continue; 1033 ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr); 1034 if (val != -1) continue; 1035 again = 2; 1036 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 1037 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 1038 for (c = 0; c < coneSize; ++c) { 1039 ierr = DMLabelGetValue(label, cone[c], &val);CHKERRQ(ierr); 1040 if (val != -1) { 1041 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); 1042 if (val > 0) { 1043 ierr = DMLabelSetValue(label, point, shift+dim);CHKERRQ(ierr); 1044 } else { 1045 ierr = DMLabelSetValue(label, point, -(shift+dim));CHKERRQ(ierr); 1046 } 1047 again = 1; 1048 break; 1049 } 1050 } 1051 } 1052 } 1053 /* Classify the rest by cell membership */ 1054 for (s = 0; s < starSize*2; s += 2) { 1055 const PetscInt point = star[s]; 1056 1057 ierr = DMLabelGetValue(label, point, &val);CHKERRQ(ierr); 1058 if (val == -1) { 1059 PetscInt *sstar = NULL; 1060 PetscInt sstarSize, ss; 1061 PetscBool marked = PETSC_FALSE; 1062 1063 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &sstarSize, &sstar);CHKERRQ(ierr); 1064 for (ss = 0; ss < sstarSize*2; ss += 2) { 1065 const PetscInt spoint = sstar[ss]; 1066 1067 if ((spoint < cStart) || (spoint >= cEnd)) continue; 1068 ierr = DMLabelGetValue(label, spoint, &val);CHKERRQ(ierr); 1069 if (val == -1) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Cell %d in star of %d does not have a valid label", spoint, point); 1070 ierr = DMPlexGetLabelValue(dm, "depth", point, &dep);CHKERRQ(ierr); 1071 if (val > 0) { 1072 ierr = DMLabelSetValue(label, point, shift+dep);CHKERRQ(ierr); 1073 } else { 1074 ierr = DMLabelSetValue(label, point, -(shift+dep));CHKERRQ(ierr); 1075 } 1076 marked = PETSC_TRUE; 1077 break; 1078 } 1079 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &sstarSize, &sstar);CHKERRQ(ierr); 1080 if (!marked) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Point %d could not be classified", point); 1081 } 1082 } 1083 ierr = DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1084 } 1085 ierr = ISRestoreIndices(dimIS, &points);CHKERRQ(ierr); 1086 ierr = ISDestroy(&dimIS);CHKERRQ(ierr); 1087 PetscFunctionReturn(0); 1088 } 1089 1090 #undef __FUNCT__ 1091 #define __FUNCT__ "DMPlexMarkSubmesh_Uninterpolated" 1092 /* Here we need the explicit assumption that: 1093 1094 For any marked cell, the marked vertices constitute a single face 1095 */ 1096 static PetscErrorCode DMPlexMarkSubmesh_Uninterpolated(DM dm, DMLabel vertexLabel, DMLabel subpointMap, PetscInt *numFaces, PetscInt *nFV, DM subdm) 1097 { 1098 IS subvertexIS; 1099 const PetscInt *subvertices; 1100 PetscInt *pStart, *pEnd, *pMax; 1101 PetscInt depth, dim, d, numSubVerticesInitial = 0, v; 1102 PetscErrorCode ierr; 1103 1104 PetscFunctionBegin; 1105 *numFaces = 0; 1106 *nFV = 0; 1107 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1108 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1109 ierr = PetscMalloc3(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd,dim+1,PetscInt,&pMax);CHKERRQ(ierr); 1110 ierr = DMPlexGetHybridBounds(dm, &pMax[depth], depth>1 ? &pMax[depth-1] : NULL, depth > 2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr); 1111 for (d = 0; d <= depth; ++d) { 1112 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1113 if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]); 1114 } 1115 /* Loop over initial vertices and mark all faces in the collective star() */ 1116 ierr = DMLabelGetStratumIS(vertexLabel, 1, &subvertexIS);CHKERRQ(ierr); 1117 if (subvertexIS) { 1118 ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr); 1119 ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1120 } 1121 for (v = 0; v < numSubVerticesInitial; ++v) { 1122 const PetscInt vertex = subvertices[v]; 1123 PetscInt *star = NULL; 1124 PetscInt starSize, s, numCells = 0, c; 1125 1126 ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1127 for (s = 0; s < starSize*2; s += 2) { 1128 const PetscInt point = star[s]; 1129 if ((point >= pStart[depth]) && (point < pEnd[depth])) star[numCells++] = point; 1130 } 1131 for (c = 0; c < numCells; ++c) { 1132 const PetscInt cell = star[c]; 1133 PetscInt *closure = NULL; 1134 PetscInt closureSize, cl; 1135 PetscInt cellLoc, numCorners = 0, faceSize = 0; 1136 1137 ierr = DMLabelGetValue(subpointMap, cell, &cellLoc);CHKERRQ(ierr); 1138 if (cellLoc == 2) continue; 1139 if (cellLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Cell %d has dimension %d in the surface label", cell, cellLoc); 1140 ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1141 for (cl = 0; cl < closureSize*2; cl += 2) { 1142 const PetscInt point = closure[cl]; 1143 PetscInt vertexLoc; 1144 1145 if ((point >= pStart[0]) && (point < pEnd[0])) { 1146 ++numCorners; 1147 ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr); 1148 if (vertexLoc >= 0) closure[faceSize++] = point; 1149 } 1150 } 1151 if (!(*nFV)) {ierr = DMPlexGetNumFaceVertices(dm, dim, numCorners, nFV);CHKERRQ(ierr);} 1152 if (faceSize > *nFV) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize); 1153 if (faceSize == *nFV) { 1154 ++(*numFaces); 1155 for (cl = 0; cl < faceSize; ++cl) { 1156 ierr = DMLabelSetValue(subpointMap, closure[cl], 0);CHKERRQ(ierr); 1157 } 1158 ierr = DMLabelSetValue(subpointMap, cell, 2);CHKERRQ(ierr); 1159 } 1160 ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1161 } 1162 ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1163 } 1164 if (subvertexIS) { 1165 ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1166 } 1167 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 1168 ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr); 1169 PetscFunctionReturn(0); 1170 } 1171 1172 #undef __FUNCT__ 1173 #define __FUNCT__ "DMPlexMarkSubmesh_Interpolated" 1174 static PetscErrorCode DMPlexMarkSubmesh_Interpolated(DM dm, DMLabel vertexLabel, DMLabel subpointMap, DM subdm) 1175 { 1176 IS subvertexIS; 1177 const PetscInt *subvertices; 1178 PetscInt *pStart, *pEnd, *pMax; 1179 PetscInt dim, d, numSubVerticesInitial = 0, v; 1180 PetscErrorCode ierr; 1181 1182 PetscFunctionBegin; 1183 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1184 ierr = PetscMalloc3(dim+1,PetscInt,&pStart,dim+1,PetscInt,&pEnd,dim+1,PetscInt,&pMax);CHKERRQ(ierr); 1185 ierr = DMPlexGetHybridBounds(dm, &pMax[dim], dim>1 ? &pMax[dim-1] : NULL, dim > 2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr); 1186 for (d = 0; d <= dim; ++d) { 1187 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 1188 if (pMax[d] >= 0) pEnd[d] = PetscMin(pEnd[d], pMax[d]); 1189 } 1190 /* Loop over initial vertices and mark all faces in the collective star() */ 1191 ierr = DMLabelGetStratumIS(vertexLabel, 1, &subvertexIS);CHKERRQ(ierr); 1192 if (subvertexIS) { 1193 ierr = ISGetSize(subvertexIS, &numSubVerticesInitial);CHKERRQ(ierr); 1194 ierr = ISGetIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1195 } 1196 for (v = 0; v < numSubVerticesInitial; ++v) { 1197 const PetscInt vertex = subvertices[v]; 1198 PetscInt *star = NULL; 1199 PetscInt starSize, s, numFaces = 0, f; 1200 1201 ierr = DMPlexGetTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1202 for (s = 0; s < starSize*2; s += 2) { 1203 const PetscInt point = star[s]; 1204 if ((point >= pStart[dim-1]) && (point < pEnd[dim-1])) star[numFaces++] = point; 1205 } 1206 for (f = 0; f < numFaces; ++f) { 1207 const PetscInt face = star[f]; 1208 PetscInt *closure = NULL; 1209 PetscInt closureSize, c; 1210 PetscInt faceLoc; 1211 1212 ierr = DMLabelGetValue(subpointMap, face, &faceLoc);CHKERRQ(ierr); 1213 if (faceLoc == dim-1) continue; 1214 if (faceLoc >= 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Face %d has dimension %d in the surface label", face, faceLoc); 1215 ierr = DMPlexGetTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1216 for (c = 0; c < closureSize*2; c += 2) { 1217 const PetscInt point = closure[c]; 1218 PetscInt vertexLoc; 1219 1220 if ((point >= pStart[0]) && (point < pEnd[0])) { 1221 ierr = DMLabelGetValue(vertexLabel, point, &vertexLoc);CHKERRQ(ierr); 1222 if (vertexLoc < 0) break; 1223 } 1224 } 1225 if (c == closureSize*2) { 1226 const PetscInt *support; 1227 PetscInt supportSize, s; 1228 1229 for (c = 0; c < closureSize*2; c += 2) { 1230 const PetscInt point = closure[c]; 1231 1232 for (d = 0; d < dim; ++d) { 1233 if ((point >= pStart[d]) && (point < pEnd[d])) { 1234 ierr = DMLabelSetValue(subpointMap, point, d);CHKERRQ(ierr); 1235 break; 1236 } 1237 } 1238 } 1239 ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); 1240 ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); 1241 for (s = 0; s < supportSize; ++s) { 1242 ierr = DMLabelSetValue(subpointMap, support[s], dim);CHKERRQ(ierr); 1243 } 1244 } 1245 ierr = DMPlexRestoreTransitiveClosure(dm, face, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1246 } 1247 ierr = DMPlexRestoreTransitiveClosure(dm, vertex, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1248 } 1249 if (subvertexIS) { 1250 ierr = ISRestoreIndices(subvertexIS, &subvertices);CHKERRQ(ierr); 1251 } 1252 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 1253 ierr = PetscFree3(pStart,pEnd,pMax);CHKERRQ(ierr); 1254 PetscFunctionReturn(0); 1255 } 1256 1257 #undef __FUNCT__ 1258 #define __FUNCT__ "DMPlexGetFaceOrientation" 1259 PetscErrorCode DMPlexGetFaceOrientation(DM dm, PetscInt cell, PetscInt numCorners, PetscInt indices[], PetscInt oppositeVertex, PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented) 1260 { 1261 MPI_Comm comm; 1262 PetscBool posOrient = PETSC_FALSE; 1263 const PetscInt debug = 0; 1264 PetscInt cellDim, faceSize, f; 1265 PetscErrorCode ierr; 1266 1267 PetscFunctionBegin; 1268 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1269 ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr); 1270 if (debug) {PetscPrintf(comm, "cellDim: %d numCorners: %d\n", cellDim, numCorners);CHKERRQ(ierr);} 1271 1272 if (cellDim == numCorners-1) { 1273 /* Simplices */ 1274 faceSize = numCorners-1; 1275 posOrient = !(oppositeVertex%2) ? PETSC_TRUE : PETSC_FALSE; 1276 } else if (cellDim == 1 && numCorners == 3) { 1277 /* Quadratic line */ 1278 faceSize = 1; 1279 posOrient = PETSC_TRUE; 1280 } else if (cellDim == 2 && numCorners == 4) { 1281 /* Quads */ 1282 faceSize = 2; 1283 if ((indices[1] > indices[0]) && (indices[1] - indices[0] == 1)) { 1284 posOrient = PETSC_TRUE; 1285 } else if ((indices[0] == 3) && (indices[1] == 0)) { 1286 posOrient = PETSC_TRUE; 1287 } else { 1288 if (((indices[0] > indices[1]) && (indices[0] - indices[1] == 1)) || ((indices[0] == 0) && (indices[1] == 3))) { 1289 posOrient = PETSC_FALSE; 1290 } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossedge"); 1291 } 1292 } else if (cellDim == 2 && numCorners == 6) { 1293 /* Quadratic triangle (I hate this) */ 1294 /* Edges are determined by the first 2 vertices (corners of edges) */ 1295 const PetscInt faceSizeTri = 3; 1296 PetscInt sortedIndices[3], i, iFace; 1297 PetscBool found = PETSC_FALSE; 1298 PetscInt faceVerticesTriSorted[9] = { 1299 0, 3, 4, /* bottom */ 1300 1, 4, 5, /* right */ 1301 2, 3, 5, /* left */ 1302 }; 1303 PetscInt faceVerticesTri[9] = { 1304 0, 3, 4, /* bottom */ 1305 1, 4, 5, /* right */ 1306 2, 5, 3, /* left */ 1307 }; 1308 1309 faceSize = faceSizeTri; 1310 for (i = 0; i < faceSizeTri; ++i) sortedIndices[i] = indices[i]; 1311 ierr = PetscSortInt(faceSizeTri, sortedIndices);CHKERRQ(ierr); 1312 for (iFace = 0; iFace < 3; ++iFace) { 1313 const PetscInt ii = iFace*faceSizeTri; 1314 PetscInt fVertex, cVertex; 1315 1316 if ((sortedIndices[0] == faceVerticesTriSorted[ii+0]) && 1317 (sortedIndices[1] == faceVerticesTriSorted[ii+1])) { 1318 for (fVertex = 0; fVertex < faceSizeTri; ++fVertex) { 1319 for (cVertex = 0; cVertex < faceSizeTri; ++cVertex) { 1320 if (indices[cVertex] == faceVerticesTri[ii+fVertex]) { 1321 faceVertices[fVertex] = origVertices[cVertex]; 1322 break; 1323 } 1324 } 1325 } 1326 found = PETSC_TRUE; 1327 break; 1328 } 1329 } 1330 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tri crossface"); 1331 if (posOriented) *posOriented = PETSC_TRUE; 1332 PetscFunctionReturn(0); 1333 } else if (cellDim == 2 && numCorners == 9) { 1334 /* Quadratic quad (I hate this) */ 1335 /* Edges are determined by the first 2 vertices (corners of edges) */ 1336 const PetscInt faceSizeQuad = 3; 1337 PetscInt sortedIndices[3], i, iFace; 1338 PetscBool found = PETSC_FALSE; 1339 PetscInt faceVerticesQuadSorted[12] = { 1340 0, 1, 4, /* bottom */ 1341 1, 2, 5, /* right */ 1342 2, 3, 6, /* top */ 1343 0, 3, 7, /* left */ 1344 }; 1345 PetscInt faceVerticesQuad[12] = { 1346 0, 1, 4, /* bottom */ 1347 1, 2, 5, /* right */ 1348 2, 3, 6, /* top */ 1349 3, 0, 7, /* left */ 1350 }; 1351 1352 faceSize = faceSizeQuad; 1353 for (i = 0; i < faceSizeQuad; ++i) sortedIndices[i] = indices[i]; 1354 ierr = PetscSortInt(faceSizeQuad, sortedIndices);CHKERRQ(ierr); 1355 for (iFace = 0; iFace < 4; ++iFace) { 1356 const PetscInt ii = iFace*faceSizeQuad; 1357 PetscInt fVertex, cVertex; 1358 1359 if ((sortedIndices[0] == faceVerticesQuadSorted[ii+0]) && 1360 (sortedIndices[1] == faceVerticesQuadSorted[ii+1])) { 1361 for (fVertex = 0; fVertex < faceSizeQuad; ++fVertex) { 1362 for (cVertex = 0; cVertex < faceSizeQuad; ++cVertex) { 1363 if (indices[cVertex] == faceVerticesQuad[ii+fVertex]) { 1364 faceVertices[fVertex] = origVertices[cVertex]; 1365 break; 1366 } 1367 } 1368 } 1369 found = PETSC_TRUE; 1370 break; 1371 } 1372 } 1373 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid quad crossface"); 1374 if (posOriented) *posOriented = PETSC_TRUE; 1375 PetscFunctionReturn(0); 1376 } else if (cellDim == 3 && numCorners == 8) { 1377 /* Hexes 1378 A hex is two oriented quads with the normal of the first 1379 pointing up at the second. 1380 1381 7---6 1382 /| /| 1383 4---5 | 1384 | 3-|-2 1385 |/ |/ 1386 0---1 1387 1388 Faces are determined by the first 4 vertices (corners of faces) */ 1389 const PetscInt faceSizeHex = 4; 1390 PetscInt sortedIndices[4], i, iFace; 1391 PetscBool found = PETSC_FALSE; 1392 PetscInt faceVerticesHexSorted[24] = { 1393 0, 1, 2, 3, /* bottom */ 1394 4, 5, 6, 7, /* top */ 1395 0, 1, 4, 5, /* front */ 1396 1, 2, 5, 6, /* right */ 1397 2, 3, 6, 7, /* back */ 1398 0, 3, 4, 7, /* left */ 1399 }; 1400 PetscInt faceVerticesHex[24] = { 1401 3, 2, 1, 0, /* bottom */ 1402 4, 5, 6, 7, /* top */ 1403 0, 1, 5, 4, /* front */ 1404 1, 2, 6, 5, /* right */ 1405 2, 3, 7, 6, /* back */ 1406 3, 0, 4, 7, /* left */ 1407 }; 1408 1409 faceSize = faceSizeHex; 1410 for (i = 0; i < faceSizeHex; ++i) sortedIndices[i] = indices[i]; 1411 ierr = PetscSortInt(faceSizeHex, sortedIndices);CHKERRQ(ierr); 1412 for (iFace = 0; iFace < 6; ++iFace) { 1413 const PetscInt ii = iFace*faceSizeHex; 1414 PetscInt fVertex, cVertex; 1415 1416 if ((sortedIndices[0] == faceVerticesHexSorted[ii+0]) && 1417 (sortedIndices[1] == faceVerticesHexSorted[ii+1]) && 1418 (sortedIndices[2] == faceVerticesHexSorted[ii+2]) && 1419 (sortedIndices[3] == faceVerticesHexSorted[ii+3])) { 1420 for (fVertex = 0; fVertex < faceSizeHex; ++fVertex) { 1421 for (cVertex = 0; cVertex < faceSizeHex; ++cVertex) { 1422 if (indices[cVertex] == faceVerticesHex[ii+fVertex]) { 1423 faceVertices[fVertex] = origVertices[cVertex]; 1424 break; 1425 } 1426 } 1427 } 1428 found = PETSC_TRUE; 1429 break; 1430 } 1431 } 1432 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface"); 1433 if (posOriented) *posOriented = PETSC_TRUE; 1434 PetscFunctionReturn(0); 1435 } else if (cellDim == 3 && numCorners == 10) { 1436 /* Quadratic tet */ 1437 /* Faces are determined by the first 3 vertices (corners of faces) */ 1438 const PetscInt faceSizeTet = 6; 1439 PetscInt sortedIndices[6], i, iFace; 1440 PetscBool found = PETSC_FALSE; 1441 PetscInt faceVerticesTetSorted[24] = { 1442 0, 1, 2, 6, 7, 8, /* bottom */ 1443 0, 3, 4, 6, 7, 9, /* front */ 1444 1, 4, 5, 7, 8, 9, /* right */ 1445 2, 3, 5, 6, 8, 9, /* left */ 1446 }; 1447 PetscInt faceVerticesTet[24] = { 1448 0, 1, 2, 6, 7, 8, /* bottom */ 1449 0, 4, 3, 6, 7, 9, /* front */ 1450 1, 5, 4, 7, 8, 9, /* right */ 1451 2, 3, 5, 8, 6, 9, /* left */ 1452 }; 1453 1454 faceSize = faceSizeTet; 1455 for (i = 0; i < faceSizeTet; ++i) sortedIndices[i] = indices[i]; 1456 ierr = PetscSortInt(faceSizeTet, sortedIndices);CHKERRQ(ierr); 1457 for (iFace=0; iFace < 4; ++iFace) { 1458 const PetscInt ii = iFace*faceSizeTet; 1459 PetscInt fVertex, cVertex; 1460 1461 if ((sortedIndices[0] == faceVerticesTetSorted[ii+0]) && 1462 (sortedIndices[1] == faceVerticesTetSorted[ii+1]) && 1463 (sortedIndices[2] == faceVerticesTetSorted[ii+2]) && 1464 (sortedIndices[3] == faceVerticesTetSorted[ii+3])) { 1465 for (fVertex = 0; fVertex < faceSizeTet; ++fVertex) { 1466 for (cVertex = 0; cVertex < faceSizeTet; ++cVertex) { 1467 if (indices[cVertex] == faceVerticesTet[ii+fVertex]) { 1468 faceVertices[fVertex] = origVertices[cVertex]; 1469 break; 1470 } 1471 } 1472 } 1473 found = PETSC_TRUE; 1474 break; 1475 } 1476 } 1477 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid tet crossface"); 1478 if (posOriented) *posOriented = PETSC_TRUE; 1479 PetscFunctionReturn(0); 1480 } else if (cellDim == 3 && numCorners == 27) { 1481 /* Quadratic hexes (I hate this) 1482 A hex is two oriented quads with the normal of the first 1483 pointing up at the second. 1484 1485 7---6 1486 /| /| 1487 4---5 | 1488 | 3-|-2 1489 |/ |/ 1490 0---1 1491 1492 Faces are determined by the first 4 vertices (corners of faces) */ 1493 const PetscInt faceSizeQuadHex = 9; 1494 PetscInt sortedIndices[9], i, iFace; 1495 PetscBool found = PETSC_FALSE; 1496 PetscInt faceVerticesQuadHexSorted[54] = { 1497 0, 1, 2, 3, 8, 9, 10, 11, 24, /* bottom */ 1498 4, 5, 6, 7, 12, 13, 14, 15, 25, /* top */ 1499 0, 1, 4, 5, 8, 12, 16, 17, 22, /* front */ 1500 1, 2, 5, 6, 9, 13, 17, 18, 21, /* right */ 1501 2, 3, 6, 7, 10, 14, 18, 19, 23, /* back */ 1502 0, 3, 4, 7, 11, 15, 16, 19, 20, /* left */ 1503 }; 1504 PetscInt faceVerticesQuadHex[54] = { 1505 3, 2, 1, 0, 10, 9, 8, 11, 24, /* bottom */ 1506 4, 5, 6, 7, 12, 13, 14, 15, 25, /* top */ 1507 0, 1, 5, 4, 8, 17, 12, 16, 22, /* front */ 1508 1, 2, 6, 5, 9, 18, 13, 17, 21, /* right */ 1509 2, 3, 7, 6, 10, 19, 14, 18, 23, /* back */ 1510 3, 0, 4, 7, 11, 16, 15, 19, 20 /* left */ 1511 }; 1512 1513 faceSize = faceSizeQuadHex; 1514 for (i = 0; i < faceSizeQuadHex; ++i) sortedIndices[i] = indices[i]; 1515 ierr = PetscSortInt(faceSizeQuadHex, sortedIndices);CHKERRQ(ierr); 1516 for (iFace = 0; iFace < 6; ++iFace) { 1517 const PetscInt ii = iFace*faceSizeQuadHex; 1518 PetscInt fVertex, cVertex; 1519 1520 if ((sortedIndices[0] == faceVerticesQuadHexSorted[ii+0]) && 1521 (sortedIndices[1] == faceVerticesQuadHexSorted[ii+1]) && 1522 (sortedIndices[2] == faceVerticesQuadHexSorted[ii+2]) && 1523 (sortedIndices[3] == faceVerticesQuadHexSorted[ii+3])) { 1524 for (fVertex = 0; fVertex < faceSizeQuadHex; ++fVertex) { 1525 for (cVertex = 0; cVertex < faceSizeQuadHex; ++cVertex) { 1526 if (indices[cVertex] == faceVerticesQuadHex[ii+fVertex]) { 1527 faceVertices[fVertex] = origVertices[cVertex]; 1528 break; 1529 } 1530 } 1531 } 1532 found = PETSC_TRUE; 1533 break; 1534 } 1535 } 1536 if (!found) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid hex crossface"); 1537 if (posOriented) *posOriented = PETSC_TRUE; 1538 PetscFunctionReturn(0); 1539 } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Unknown cell type for faceOrientation()."); 1540 if (!posOrient) { 1541 if (debug) {ierr = PetscPrintf(comm, " Reversing initial face orientation\n");CHKERRQ(ierr);} 1542 for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[faceSize-1 - f]; 1543 } else { 1544 if (debug) {ierr = PetscPrintf(comm, " Keeping initial face orientation\n");CHKERRQ(ierr);} 1545 for (f = 0; f < faceSize; ++f) faceVertices[f] = origVertices[f]; 1546 } 1547 if (posOriented) *posOriented = posOrient; 1548 PetscFunctionReturn(0); 1549 } 1550 1551 #undef __FUNCT__ 1552 #define __FUNCT__ "DMPlexGetOrientedFace" 1553 /* 1554 Given a cell and a face, as a set of vertices, 1555 return the oriented face, as a set of vertices, in faceVertices 1556 The orientation is such that the face normal points out of the cell 1557 */ 1558 PetscErrorCode DMPlexGetOrientedFace(DM dm, PetscInt cell, PetscInt faceSize, const PetscInt face[], PetscInt numCorners, PetscInt indices[], PetscInt origVertices[], PetscInt faceVertices[], PetscBool *posOriented) 1559 { 1560 const PetscInt *cone = NULL; 1561 PetscInt coneSize, v, f, v2; 1562 PetscInt oppositeVertex = -1; 1563 PetscErrorCode ierr; 1564 1565 PetscFunctionBegin; 1566 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 1567 ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 1568 for (v = 0, v2 = 0; v < coneSize; ++v) { 1569 PetscBool found = PETSC_FALSE; 1570 1571 for (f = 0; f < faceSize; ++f) { 1572 if (face[f] == cone[v]) { 1573 found = PETSC_TRUE; break; 1574 } 1575 } 1576 if (found) { 1577 indices[v2] = v; 1578 origVertices[v2] = cone[v]; 1579 ++v2; 1580 } else { 1581 oppositeVertex = v; 1582 } 1583 } 1584 ierr = DMPlexGetFaceOrientation(dm, cell, numCorners, indices, oppositeVertex, origVertices, faceVertices, posOriented);CHKERRQ(ierr); 1585 PetscFunctionReturn(0); 1586 } 1587 1588 #undef __FUNCT__ 1589 #define __FUNCT__ "DMPlexInsertFace_Internal" 1590 /* 1591 DMPlexInsertFace_Internal - Puts a face into the mesh 1592 1593 Not collective 1594 1595 Input Parameters: 1596 + dm - The DMPlex 1597 . numFaceVertex - The number of vertices in the face 1598 . faceVertices - The vertices in the face for dm 1599 . subfaceVertices - The vertices in the face for subdm 1600 . numCorners - The number of vertices in the cell 1601 . cell - A cell in dm containing the face 1602 . subcell - A cell in subdm containing the face 1603 . firstFace - First face in the mesh 1604 - newFacePoint - Next face in the mesh 1605 1606 Output Parameters: 1607 . newFacePoint - Contains next face point number on input, updated on output 1608 1609 Level: developer 1610 */ 1611 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) 1612 { 1613 MPI_Comm comm; 1614 DM_Plex *submesh = (DM_Plex*) subdm->data; 1615 const PetscInt *faces; 1616 PetscInt numFaces, coneSize; 1617 PetscErrorCode ierr; 1618 1619 PetscFunctionBegin; 1620 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1621 ierr = DMPlexGetConeSize(subdm, subcell, &coneSize);CHKERRQ(ierr); 1622 if (coneSize != 1) SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cone size of cell %d is %d != 1", cell, coneSize); 1623 #if 0 1624 /* Cannot use this because support() has not been constructed yet */ 1625 ierr = DMPlexGetJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr); 1626 #else 1627 { 1628 PetscInt f; 1629 1630 numFaces = 0; 1631 ierr = DMGetWorkArray(subdm, 1, PETSC_INT, (void**) &faces);CHKERRQ(ierr); 1632 for (f = firstFace; f < *newFacePoint; ++f) { 1633 PetscInt dof, off, d; 1634 1635 ierr = PetscSectionGetDof(submesh->coneSection, f, &dof);CHKERRQ(ierr); 1636 ierr = PetscSectionGetOffset(submesh->coneSection, f, &off);CHKERRQ(ierr); 1637 /* Yes, I know this is quadratic, but I expect the sizes to be <5 */ 1638 for (d = 0; d < dof; ++d) { 1639 const PetscInt p = submesh->cones[off+d]; 1640 PetscInt v; 1641 1642 for (v = 0; v < numFaceVertices; ++v) { 1643 if (subfaceVertices[v] == p) break; 1644 } 1645 if (v == numFaceVertices) break; 1646 } 1647 if (d == dof) { 1648 numFaces = 1; 1649 ((PetscInt*) faces)[0] = f; 1650 } 1651 } 1652 } 1653 #endif 1654 if (numFaces > 1) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Vertex set had %d faces, not one", numFaces); 1655 else if (numFaces == 1) { 1656 /* Add the other cell neighbor for this face */ 1657 ierr = DMPlexSetCone(subdm, cell, faces);CHKERRQ(ierr); 1658 } else { 1659 PetscInt *indices, *origVertices, *orientedVertices, *orientedSubVertices, v, ov; 1660 PetscBool posOriented; 1661 1662 ierr = DMGetWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr); 1663 origVertices = &orientedVertices[numFaceVertices]; 1664 indices = &orientedVertices[numFaceVertices*2]; 1665 orientedSubVertices = &orientedVertices[numFaceVertices*3]; 1666 ierr = DMPlexGetOrientedFace(dm, cell, numFaceVertices, faceVertices, numCorners, indices, origVertices, orientedVertices, &posOriented);CHKERRQ(ierr); 1667 /* TODO: I know that routine should return a permutation, not the indices */ 1668 for (v = 0; v < numFaceVertices; ++v) { 1669 const PetscInt vertex = faceVertices[v], subvertex = subfaceVertices[v]; 1670 for (ov = 0; ov < numFaceVertices; ++ov) { 1671 if (orientedVertices[ov] == vertex) { 1672 orientedSubVertices[ov] = subvertex; 1673 break; 1674 } 1675 } 1676 if (ov == numFaceVertices) SETERRQ1(comm, PETSC_ERR_PLIB, "Could not find face vertex %d in orientated set", vertex); 1677 } 1678 ierr = DMPlexSetCone(subdm, *newFacePoint, orientedSubVertices);CHKERRQ(ierr); 1679 ierr = DMPlexSetCone(subdm, subcell, newFacePoint);CHKERRQ(ierr); 1680 ierr = DMRestoreWorkArray(subdm, 4*numFaceVertices * sizeof(PetscInt), PETSC_INT, &orientedVertices);CHKERRQ(ierr); 1681 ++(*newFacePoint); 1682 } 1683 ierr = DMPlexRestoreJoin(subdm, numFaceVertices, subfaceVertices, &numFaces, &faces);CHKERRQ(ierr); 1684 PetscFunctionReturn(0); 1685 } 1686 1687 #undef __FUNCT__ 1688 #define __FUNCT__ "DMPlexCreateSubmesh_Uninterpolated" 1689 static PetscErrorCode DMPlexCreateSubmesh_Uninterpolated(DM dm, const char vertexLabelName[], DM subdm) 1690 { 1691 MPI_Comm comm; 1692 DMLabel vertexLabel, subpointMap; 1693 IS subvertexIS, subcellIS; 1694 const PetscInt *subVertices, *subCells; 1695 PetscInt numSubVertices, firstSubVertex, numSubCells; 1696 PetscInt *subface, maxConeSize, numSubFaces, firstSubFace, newFacePoint, nFV; 1697 PetscInt vStart, vEnd, c, f; 1698 PetscErrorCode ierr; 1699 1700 PetscFunctionBegin; 1701 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1702 /* Create subpointMap which marks the submesh */ 1703 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 1704 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 1705 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 1706 ierr = DMPlexGetLabel(dm, vertexLabelName, &vertexLabel);CHKERRQ(ierr); 1707 ierr = DMPlexMarkSubmesh_Uninterpolated(dm, vertexLabel, subpointMap, &numSubFaces, &nFV, subdm);CHKERRQ(ierr); 1708 /* Setup chart */ 1709 ierr = DMLabelGetStratumSize(subpointMap, 0, &numSubVertices);CHKERRQ(ierr); 1710 ierr = DMLabelGetStratumSize(subpointMap, 2, &numSubCells);CHKERRQ(ierr); 1711 ierr = DMPlexSetChart(subdm, 0, numSubCells+numSubFaces+numSubVertices);CHKERRQ(ierr); 1712 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 1713 /* Set cone sizes */ 1714 firstSubVertex = numSubCells; 1715 firstSubFace = numSubCells+numSubVertices; 1716 newFacePoint = firstSubFace; 1717 ierr = DMLabelGetStratumIS(subpointMap, 0, &subvertexIS);CHKERRQ(ierr); 1718 if (subvertexIS) {ierr = ISGetIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 1719 ierr = DMLabelGetStratumIS(subpointMap, 2, &subcellIS);CHKERRQ(ierr); 1720 if (subcellIS) {ierr = ISGetIndices(subcellIS, &subCells);CHKERRQ(ierr);} 1721 for (c = 0; c < numSubCells; ++c) { 1722 ierr = DMPlexSetConeSize(subdm, c, 1);CHKERRQ(ierr); 1723 } 1724 for (f = firstSubFace; f < firstSubFace+numSubFaces; ++f) { 1725 ierr = DMPlexSetConeSize(subdm, f, nFV);CHKERRQ(ierr); 1726 } 1727 ierr = DMSetUp(subdm);CHKERRQ(ierr); 1728 /* Create face cones */ 1729 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1730 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 1731 ierr = DMGetWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 1732 for (c = 0; c < numSubCells; ++c) { 1733 const PetscInt cell = subCells[c]; 1734 const PetscInt subcell = c; 1735 PetscInt *closure = NULL; 1736 PetscInt closureSize, cl, numCorners = 0, faceSize = 0; 1737 1738 ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1739 for (cl = 0; cl < closureSize*2; cl += 2) { 1740 const PetscInt point = closure[cl]; 1741 PetscInt subVertex; 1742 1743 if ((point >= vStart) && (point < vEnd)) { 1744 ++numCorners; 1745 ierr = PetscFindInt(point, numSubVertices, subVertices, &subVertex);CHKERRQ(ierr); 1746 if (subVertex >= 0) { 1747 closure[faceSize] = point; 1748 subface[faceSize] = firstSubVertex+subVertex; 1749 ++faceSize; 1750 } 1751 } 1752 } 1753 if (faceSize > nFV) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Invalid submesh: Too many vertices %d of an element on the surface", faceSize); 1754 if (faceSize == nFV) { 1755 ierr = DMPlexInsertFace_Internal(dm, subdm, faceSize, closure, subface, numCorners, cell, subcell, firstSubFace, &newFacePoint);CHKERRQ(ierr); 1756 } 1757 ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1758 } 1759 ierr = DMRestoreWorkArray(subdm, maxConeSize, PETSC_INT, (void**) &subface);CHKERRQ(ierr); 1760 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 1761 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 1762 /* Build coordinates */ 1763 { 1764 PetscSection coordSection, subCoordSection; 1765 Vec coordinates, subCoordinates; 1766 PetscScalar *coords, *subCoords; 1767 PetscInt coordSize, v; 1768 1769 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1770 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1771 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 1772 ierr = PetscSectionSetChart(subCoordSection, firstSubVertex, firstSubVertex+numSubVertices);CHKERRQ(ierr); 1773 for (v = 0; v < numSubVertices; ++v) { 1774 const PetscInt vertex = subVertices[v]; 1775 const PetscInt subvertex = firstSubVertex+v; 1776 PetscInt dof; 1777 1778 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 1779 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 1780 } 1781 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 1782 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 1783 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 1784 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1785 ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr); 1786 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1787 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 1788 for (v = 0; v < numSubVertices; ++v) { 1789 const PetscInt vertex = subVertices[v]; 1790 const PetscInt subvertex = firstSubVertex+v; 1791 PetscInt dof, off, sdof, soff, d; 1792 1793 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 1794 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 1795 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 1796 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 1797 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 1798 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 1799 } 1800 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1801 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 1802 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 1803 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 1804 } 1805 /* Cleanup */ 1806 if (subvertexIS) {ierr = ISRestoreIndices(subvertexIS, &subVertices);CHKERRQ(ierr);} 1807 ierr = ISDestroy(&subvertexIS);CHKERRQ(ierr); 1808 if (subcellIS) {ierr = ISRestoreIndices(subcellIS, &subCells);CHKERRQ(ierr);} 1809 ierr = ISDestroy(&subcellIS);CHKERRQ(ierr); 1810 PetscFunctionReturn(0); 1811 } 1812 1813 #undef __FUNCT__ 1814 #define __FUNCT__ "DMPlexCreateSubmesh_Interpolated" 1815 static PetscErrorCode DMPlexCreateSubmesh_Interpolated(DM dm, const char vertexLabelName[], DM subdm) 1816 { 1817 MPI_Comm comm; 1818 DMLabel subpointMap, vertexLabel; 1819 IS *subpointIS; 1820 const PetscInt **subpoints; 1821 PetscInt *numSubPoints, *firstSubPoint, *coneNew; 1822 PetscInt totSubPoints = 0, maxConeSize, dim, p, d, v; 1823 PetscErrorCode ierr; 1824 1825 PetscFunctionBegin; 1826 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1827 /* Create subpointMap which marks the submesh */ 1828 ierr = DMLabelCreate("subpoint_map", &subpointMap);CHKERRQ(ierr); 1829 ierr = DMPlexSetSubpointMap(subdm, subpointMap);CHKERRQ(ierr); 1830 ierr = DMLabelDestroy(&subpointMap);CHKERRQ(ierr); 1831 ierr = DMPlexGetLabel(dm, vertexLabelName, &vertexLabel);CHKERRQ(ierr); 1832 ierr = DMPlexMarkSubmesh_Interpolated(dm, vertexLabel, subpointMap, subdm);CHKERRQ(ierr); 1833 /* Setup chart */ 1834 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1835 ierr = PetscMalloc4(dim+1,PetscInt,&numSubPoints,dim+1,PetscInt,&firstSubPoint,dim+1,IS,&subpointIS,dim+1,const PetscInt *,&subpoints);CHKERRQ(ierr); 1836 for (d = 0; d <= dim; ++d) { 1837 ierr = DMLabelGetStratumSize(subpointMap, d, &numSubPoints[d]);CHKERRQ(ierr); 1838 totSubPoints += numSubPoints[d]; 1839 } 1840 ierr = DMPlexSetChart(subdm, 0, totSubPoints);CHKERRQ(ierr); 1841 ierr = DMPlexSetVTKCellHeight(subdm, 1);CHKERRQ(ierr); 1842 /* Set cone sizes */ 1843 firstSubPoint[dim] = 0; 1844 firstSubPoint[0] = firstSubPoint[dim] + numSubPoints[dim]; 1845 if (dim > 1) {firstSubPoint[dim-1] = firstSubPoint[0] + numSubPoints[0];} 1846 if (dim > 2) {firstSubPoint[dim-2] = firstSubPoint[dim-1] + numSubPoints[dim-1];} 1847 for (d = 0; d <= dim; ++d) { 1848 ierr = DMLabelGetStratumIS(subpointMap, d, &subpointIS[d]);CHKERRQ(ierr); 1849 ierr = ISGetIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 1850 } 1851 for (d = 0; d <= dim; ++d) { 1852 for (p = 0; p < numSubPoints[d]; ++p) { 1853 const PetscInt point = subpoints[d][p]; 1854 const PetscInt subpoint = firstSubPoint[d] + p; 1855 const PetscInt *cone; 1856 PetscInt coneSize, coneSizeNew, c, val; 1857 1858 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 1859 ierr = DMPlexSetConeSize(subdm, subpoint, coneSize);CHKERRQ(ierr); 1860 if (d == dim) { 1861 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 1862 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 1863 ierr = DMLabelGetValue(subpointMap, cone[c], &val);CHKERRQ(ierr); 1864 if (val >= 0) coneSizeNew++; 1865 } 1866 ierr = DMPlexSetConeSize(subdm, subpoint, coneSizeNew);CHKERRQ(ierr); 1867 } 1868 } 1869 } 1870 ierr = DMSetUp(subdm);CHKERRQ(ierr); 1871 /* Set cones */ 1872 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 1873 ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &coneNew);CHKERRQ(ierr); 1874 for (d = 0; d <= dim; ++d) { 1875 for (p = 0; p < numSubPoints[d]; ++p) { 1876 const PetscInt point = subpoints[d][p]; 1877 const PetscInt subpoint = firstSubPoint[d] + p; 1878 const PetscInt *cone; 1879 PetscInt coneSize, subconeSize, coneSizeNew, c, subc; 1880 1881 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 1882 ierr = DMPlexGetConeSize(subdm, subpoint, &subconeSize);CHKERRQ(ierr); 1883 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 1884 for (c = 0, coneSizeNew = 0; c < coneSize; ++c) { 1885 ierr = PetscFindInt(cone[c], numSubPoints[d-1], subpoints[d-1], &subc);CHKERRQ(ierr); 1886 if (subc >= 0) coneNew[coneSizeNew++] = firstSubPoint[d-1] + subc; 1887 } 1888 if (coneSizeNew != subconeSize) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of cone points located %d does not match subcone size %d", coneSizeNew, subconeSize); 1889 ierr = DMPlexSetCone(subdm, subpoint, coneNew);CHKERRQ(ierr); 1890 } 1891 } 1892 ierr = PetscFree(coneNew);CHKERRQ(ierr); 1893 ierr = DMPlexSymmetrize(subdm);CHKERRQ(ierr); 1894 ierr = DMPlexStratify(subdm);CHKERRQ(ierr); 1895 /* Build coordinates */ 1896 { 1897 PetscSection coordSection, subCoordSection; 1898 Vec coordinates, subCoordinates; 1899 PetscScalar *coords, *subCoords; 1900 PetscInt coordSize; 1901 1902 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1903 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1904 ierr = DMPlexGetCoordinateSection(subdm, &subCoordSection);CHKERRQ(ierr); 1905 ierr = PetscSectionSetChart(subCoordSection, firstSubPoint[0], firstSubPoint[0]+numSubPoints[0]);CHKERRQ(ierr); 1906 for (v = 0; v < numSubPoints[0]; ++v) { 1907 const PetscInt vertex = subpoints[0][v]; 1908 const PetscInt subvertex = firstSubPoint[0]+v; 1909 PetscInt dof; 1910 1911 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 1912 ierr = PetscSectionSetDof(subCoordSection, subvertex, dof);CHKERRQ(ierr); 1913 } 1914 ierr = PetscSectionSetUp(subCoordSection);CHKERRQ(ierr); 1915 ierr = PetscSectionGetStorageSize(subCoordSection, &coordSize);CHKERRQ(ierr); 1916 ierr = VecCreate(comm, &subCoordinates);CHKERRQ(ierr); 1917 ierr = VecSetSizes(subCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1918 ierr = VecSetFromOptions(subCoordinates);CHKERRQ(ierr); 1919 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1920 ierr = VecGetArray(subCoordinates, &subCoords);CHKERRQ(ierr); 1921 for (v = 0; v < numSubPoints[0]; ++v) { 1922 const PetscInt vertex = subpoints[0][v]; 1923 const PetscInt subvertex = firstSubPoint[0]+v; 1924 PetscInt dof, off, sdof, soff, d; 1925 1926 ierr = PetscSectionGetDof(coordSection, vertex, &dof);CHKERRQ(ierr); 1927 ierr = PetscSectionGetOffset(coordSection, vertex, &off);CHKERRQ(ierr); 1928 ierr = PetscSectionGetDof(subCoordSection, subvertex, &sdof);CHKERRQ(ierr); 1929 ierr = PetscSectionGetOffset(subCoordSection, subvertex, &soff);CHKERRQ(ierr); 1930 if (dof != sdof) SETERRQ4(comm, PETSC_ERR_PLIB, "Coordinate dimension %d on subvertex %d, vertex %d should be %d", sdof, subvertex, vertex, dof); 1931 for (d = 0; d < dof; ++d) subCoords[soff+d] = coords[off+d]; 1932 } 1933 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1934 ierr = VecRestoreArray(subCoordinates, &subCoords);CHKERRQ(ierr); 1935 ierr = DMSetCoordinatesLocal(subdm, subCoordinates);CHKERRQ(ierr); 1936 ierr = VecDestroy(&subCoordinates);CHKERRQ(ierr); 1937 } 1938 /* Cleanup */ 1939 for (d = 0; d <= dim; ++d) { 1940 ierr = ISRestoreIndices(subpointIS[d], &subpoints[d]);CHKERRQ(ierr); 1941 ierr = ISDestroy(&subpointIS[d]);CHKERRQ(ierr); 1942 } 1943 ierr = PetscFree4(numSubPoints,firstSubPoint,subpointIS,subpoints);CHKERRQ(ierr); 1944 PetscFunctionReturn(0); 1945 } 1946 1947 #undef __FUNCT__ 1948 #define __FUNCT__ "DMPlexCreateSubmesh" 1949 /* 1950 DMPlexCreateSubmesh - Extract a hypersurface from the mesh using vertices defined by a label 1951 1952 Input Parameters: 1953 + dm - The original mesh 1954 - vertexLabel - The DMLabel marking vertices contained in the surface 1955 1956 Output Parameter: 1957 . subdm - The surface mesh 1958 1959 Note: This function produces a DMLabel mapping original points in the submesh to their depth. This can be obtained using DMPlexGetSubpointMap(). 1960 1961 Level: developer 1962 1963 .seealso: DMPlexGetSubpointMap(), DMPlexGetLabel(), DMLabelSetValue() 1964 */ 1965 PetscErrorCode DMPlexCreateSubmesh(DM dm, const char vertexLabel[], DM *subdm) 1966 { 1967 PetscInt dim, depth; 1968 PetscErrorCode ierr; 1969 1970 PetscFunctionBegin; 1971 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1972 PetscValidCharPointer(vertexLabel, 2); 1973 PetscValidPointer(subdm, 4); 1974 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1975 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1976 ierr = DMCreate(PetscObjectComm((PetscObject)dm), subdm);CHKERRQ(ierr); 1977 ierr = DMSetType(*subdm, DMPLEX);CHKERRQ(ierr); 1978 ierr = DMPlexSetDimension(*subdm, dim-1);CHKERRQ(ierr); 1979 if (depth == dim) { 1980 ierr = DMPlexCreateSubmesh_Interpolated(dm, vertexLabel, *subdm);CHKERRQ(ierr); 1981 } else { 1982 ierr = DMPlexCreateSubmesh_Uninterpolated(dm, vertexLabel, *subdm);CHKERRQ(ierr); 1983 } 1984 PetscFunctionReturn(0); 1985 } 1986 1987 #undef __FUNCT__ 1988 #define __FUNCT__ "DMPlexGetSubpointMap" 1989 PetscErrorCode DMPlexGetSubpointMap(DM dm, DMLabel *subpointMap) 1990 { 1991 DM_Plex *mesh = (DM_Plex*) dm->data; 1992 1993 PetscFunctionBegin; 1994 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1995 PetscValidPointer(subpointMap, 2); 1996 *subpointMap = mesh->subpointMap; 1997 PetscFunctionReturn(0); 1998 } 1999 2000 #undef __FUNCT__ 2001 #define __FUNCT__ "DMPlexSetSubpointMap" 2002 /* Note: Should normally not be called by the user, since it is set in DMPlexCreateSubmesh() */ 2003 PetscErrorCode DMPlexSetSubpointMap(DM dm, DMLabel subpointMap) 2004 { 2005 DM_Plex *mesh = (DM_Plex *) dm->data; 2006 PetscErrorCode ierr; 2007 2008 PetscFunctionBegin; 2009 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2010 ierr = DMLabelDestroy(&mesh->subpointMap);CHKERRQ(ierr); 2011 mesh->subpointMap = subpointMap; 2012 ++mesh->subpointMap->refct;CHKERRQ(ierr); 2013 PetscFunctionReturn(0); 2014 } 2015 2016 #undef __FUNCT__ 2017 #define __FUNCT__ "DMPlexCreateSubpointIS" 2018 /* 2019 DMPlexCreateSubpointIS - Creates an IS covering the entire subdm chart with the original points as data 2020 2021 Input Parameter: 2022 . dm - The submesh DM 2023 2024 Output Parameter: 2025 . subpointIS - The IS of all the points from the original mesh in this submesh, or NULL if this is not a submesh 2026 2027 Note: This is IS is guaranteed to be sorted by the construction of the submesh 2028 */ 2029 PetscErrorCode DMPlexCreateSubpointIS(DM dm, IS *subpointIS) 2030 { 2031 MPI_Comm comm; 2032 DMLabel subpointMap; 2033 IS is; 2034 const PetscInt *opoints; 2035 PetscInt *points, *depths; 2036 PetscInt depth, depStart, depEnd, d, pStart, pEnd, p, n, off; 2037 PetscErrorCode ierr; 2038 2039 PetscFunctionBegin; 2040 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2041 PetscValidPointer(subpointIS, 2); 2042 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2043 *subpointIS = NULL; 2044 ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 2045 if (subpointMap) { 2046 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 2047 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 2048 if (pStart) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Submeshes must start the point numbering at 0, not %d", pStart); 2049 ierr = DMGetWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr); 2050 depths[0] = depth; 2051 depths[1] = 0; 2052 for(d = 2; d <= depth; ++d) {depths[d] = depth+1 - d;} 2053 ierr = PetscMalloc(pEnd * sizeof(PetscInt), &points);CHKERRQ(ierr); 2054 for(d = 0, off = 0; d <= depth; ++d) { 2055 const PetscInt dep = depths[d]; 2056 2057 ierr = DMPlexGetDepthStratum(dm, dep, &depStart, &depEnd);CHKERRQ(ierr); 2058 ierr = DMLabelGetStratumSize(subpointMap, dep, &n);CHKERRQ(ierr); 2059 if (((d < 2) && (depth > 1)) || (d == 1)) { /* Only check vertices and cells for now since the map is broken for others */ 2060 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); 2061 } else { 2062 if (!n) { 2063 if (d == 0) { 2064 /* Missing cells */ 2065 for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = -1; 2066 } else { 2067 /* Missing faces */ 2068 for(p = 0; p < depEnd-depStart; ++p, ++off) points[off] = PETSC_MAX_INT; 2069 } 2070 } 2071 } 2072 if (n) { 2073 ierr = DMLabelGetStratumIS(subpointMap, dep, &is);CHKERRQ(ierr); 2074 ierr = ISGetIndices(is, &opoints);CHKERRQ(ierr); 2075 for(p = 0; p < n; ++p, ++off) points[off] = opoints[p]; 2076 ierr = ISRestoreIndices(is, &opoints);CHKERRQ(ierr); 2077 ierr = ISDestroy(&is);CHKERRQ(ierr); 2078 } 2079 } 2080 ierr = DMRestoreWorkArray(dm, depth+1, PETSC_INT, &depths);CHKERRQ(ierr); 2081 if (off != pEnd) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "The number of mapped submesh points %d should be %d", off, pEnd); 2082 ierr = ISCreateGeneral(comm, pEnd, points, PETSC_OWN_POINTER, subpointIS);CHKERRQ(ierr); 2083 } 2084 PetscFunctionReturn(0); 2085 } 2086