1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <petsc/private/matorderimpl.h> /*I "petscmat.h" I*/ 3 #include <petsc/private/dmlabelimpl.h> 4 5 static PetscErrorCode DMPlexCreateOrderingClosure_Static(DM dm, PetscInt numPoints, const PetscInt pperm[], PetscInt **clperm, PetscInt **invclperm) 6 { 7 PetscInt *perm, *iperm; 8 PetscInt depth, d, pStart, pEnd, fStart, fMax, fEnd, p; 9 10 PetscFunctionBegin; 11 PetscCall(DMPlexGetDepth(dm, &depth)); 12 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 13 PetscCall(PetscMalloc1(pEnd - pStart, &perm)); 14 PetscCall(PetscMalloc1(pEnd - pStart, &iperm)); 15 for (p = pStart; p < pEnd; ++p) iperm[p] = -1; 16 for (d = depth; d > 0; --d) { 17 PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd)); 18 PetscCall(DMPlexGetDepthStratum(dm, d - 1, &fStart, &fEnd)); 19 fMax = fStart; 20 for (p = pStart; p < pEnd; ++p) { 21 const PetscInt *cone; 22 PetscInt point, coneSize, c; 23 24 if (d == depth) { 25 perm[p] = pperm[p]; 26 iperm[pperm[p]] = p; 27 } 28 point = perm[p]; 29 PetscCall(DMPlexGetConeSize(dm, point, &coneSize)); 30 PetscCall(DMPlexGetCone(dm, point, &cone)); 31 for (c = 0; c < coneSize; ++c) { 32 const PetscInt oldc = cone[c]; 33 const PetscInt newc = iperm[oldc]; 34 35 if (newc < 0) { 36 perm[fMax] = oldc; 37 iperm[oldc] = fMax++; 38 } 39 } 40 } 41 PetscCheck(fMax == fEnd, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of depth %" PetscInt_FMT " faces %" PetscInt_FMT " does not match permuted number %" PetscInt_FMT, d, fEnd - fStart, fMax - fStart); 42 } 43 *clperm = perm; 44 *invclperm = iperm; 45 PetscFunctionReturn(PETSC_SUCCESS); 46 } 47 48 /*@ 49 DMPlexGetOrdering - Calculate a reordering of the mesh 50 51 Collective on dm 52 53 Input Parameters: 54 + dm - The DMPlex object 55 . otype - type of reordering, one of the following: 56 $ MATORDERINGNATURAL - Natural 57 $ MATORDERINGND - Nested Dissection 58 $ MATORDERING1WD - One-way Dissection 59 $ MATORDERINGRCM - Reverse Cuthill-McKee 60 $ MATORDERINGQMD - Quotient Minimum Degree 61 - label - [Optional] Label used to segregate ordering into sets, or NULL 62 63 Output Parameter: 64 . perm - The point permutation as an IS, perm[old point number] = new point number 65 66 Note: The label is used to group sets of points together by label value. This makes it easy to reorder a mesh which 67 has different types of cells, and then loop over each set of reordered cells for assembly. 68 69 Level: intermediate 70 71 .seealso: `DMPlexPermute()`, `MatGetOrdering()` 72 @*/ 73 PetscErrorCode DMPlexGetOrdering(DM dm, MatOrderingType otype, DMLabel label, IS *perm) 74 { 75 PetscInt numCells = 0; 76 PetscInt *start = NULL, *adjacency = NULL, *cperm, *clperm = NULL, *invclperm = NULL, *mask, *xls, pStart, pEnd, c, i; 77 78 PetscFunctionBegin; 79 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 80 PetscValidPointer(perm, 4); 81 PetscCall(DMPlexCreateNeighborCSR(dm, 0, &numCells, &start, &adjacency)); 82 PetscCall(PetscMalloc3(numCells, &cperm, numCells, &mask, numCells * 2, &xls)); 83 if (numCells) { 84 /* Shift for Fortran numbering */ 85 for (i = 0; i < start[numCells]; ++i) ++adjacency[i]; 86 for (i = 0; i <= numCells; ++i) ++start[i]; 87 PetscCall(SPARSEPACKgenrcm(&numCells, start, adjacency, cperm, mask, xls)); 88 } 89 PetscCall(PetscFree(start)); 90 PetscCall(PetscFree(adjacency)); 91 /* Shift for Fortran numbering */ 92 for (c = 0; c < numCells; ++c) --cperm[c]; 93 /* Segregate */ 94 if (label) { 95 IS valueIS; 96 const PetscInt *valuesTmp; 97 PetscInt *values; 98 PetscInt numValues, numPoints = 0; 99 PetscInt *sperm, *vsize, *voff, v; 100 101 // Can't directly sort the valueIS, since it is a view into the DMLabel 102 PetscCall(DMLabelGetValueIS(label, &valueIS)); 103 PetscCall(ISGetLocalSize(valueIS, &numValues)); 104 PetscCall(ISGetIndices(valueIS, &valuesTmp)); 105 PetscCall(PetscCalloc4(numCells, &sperm, numValues, &values, numValues, &vsize, numValues + 1, &voff)); 106 PetscCall(PetscArraycpy(values, valuesTmp, numValues)); 107 PetscCall(PetscSortInt(numValues, values)); 108 PetscCall(ISRestoreIndices(valueIS, &valuesTmp)); 109 PetscCall(ISDestroy(&valueIS)); 110 for (v = 0; v < numValues; ++v) { 111 PetscCall(DMLabelGetStratumSize(label, values[v], &vsize[v])); 112 if (v < numValues - 1) voff[v + 2] += vsize[v] + voff[v + 1]; 113 numPoints += vsize[v]; 114 } 115 PetscCheck(numPoints == numCells, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label only covers %" PetscInt_FMT " cells != %" PetscInt_FMT " total cells", numPoints, numCells); 116 for (c = 0; c < numCells; ++c) { 117 const PetscInt oldc = cperm[c]; 118 PetscInt val, vloc; 119 120 PetscCall(DMLabelGetValue(label, oldc, &val)); 121 PetscCheck(val != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cell %" PetscInt_FMT " not present in label", oldc); 122 PetscCall(PetscFindInt(val, numValues, values, &vloc)); 123 PetscCheck(vloc >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %" PetscInt_FMT " not present label", val); 124 sperm[voff[vloc + 1]++] = oldc; 125 } 126 for (v = 0; v < numValues; ++v) PetscCheck(voff[v + 1] - voff[v] == vsize[v], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of %" PetscInt_FMT " values found is %" PetscInt_FMT " != %" PetscInt_FMT, values[v], voff[v + 1] - voff[v], vsize[v]); 127 PetscCall(PetscArraycpy(cperm, sperm, numCells)); 128 PetscCall(PetscFree4(sperm, values, vsize, voff)); 129 } 130 /* Construct closure */ 131 PetscCall(DMPlexCreateOrderingClosure_Static(dm, numCells, cperm, &clperm, &invclperm)); 132 PetscCall(PetscFree3(cperm, mask, xls)); 133 PetscCall(PetscFree(clperm)); 134 /* Invert permutation */ 135 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 136 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), pEnd - pStart, invclperm, PETSC_OWN_POINTER, perm)); 137 PetscFunctionReturn(PETSC_SUCCESS); 138 } 139 140 /*@ 141 DMPlexGetOrdering1D - Reorder the vertices so that the mesh is in a line 142 143 Collective on dm 144 145 Input Parameter: 146 . dm - The DMPlex object 147 148 Output Parameter: 149 . perm - The point permutation as an IS, perm[old point number] = new point number 150 151 Level: intermediate 152 153 .seealso: `DMPlexGetOrdering()`, `DMPlexPermute()`, `MatGetOrdering()` 154 @*/ 155 PetscErrorCode DMPlexGetOrdering1D(DM dm, IS *perm) 156 { 157 PetscInt *points; 158 const PetscInt *support, *cone; 159 PetscInt dim, pStart, pEnd, cStart, cEnd, c, vStart, vEnd, v, suppSize, lastCell = 0; 160 161 PetscFunctionBegin; 162 PetscCall(DMGetDimension(dm, &dim)); 163 PetscCheck(dim == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Input mesh must be one dimensional, not %" PetscInt_FMT, dim); 164 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 165 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 166 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 167 PetscCall(PetscMalloc1(pEnd - pStart, &points)); 168 for (c = cStart; c < cEnd; ++c) points[c] = c; 169 for (v = vStart; v < vEnd; ++v) points[v] = v; 170 for (v = vStart; v < vEnd; ++v) { 171 PetscCall(DMPlexGetSupportSize(dm, v, &suppSize)); 172 PetscCall(DMPlexGetSupport(dm, v, &support)); 173 if (suppSize == 1) { 174 lastCell = support[0]; 175 break; 176 } 177 } 178 if (v < vEnd) { 179 PetscInt pos = cEnd; 180 181 points[v] = pos++; 182 while (lastCell >= cStart) { 183 PetscCall(DMPlexGetCone(dm, lastCell, &cone)); 184 if (cone[0] == v) v = cone[1]; 185 else v = cone[0]; 186 PetscCall(DMPlexGetSupport(dm, v, &support)); 187 PetscCall(DMPlexGetSupportSize(dm, v, &suppSize)); 188 if (suppSize == 1) { 189 lastCell = -1; 190 } else { 191 if (support[0] == lastCell) lastCell = support[1]; 192 else lastCell = support[0]; 193 } 194 points[v] = pos++; 195 } 196 PetscCheck(pos == pEnd, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Last vertex was %" PetscInt_FMT ", not %" PetscInt_FMT, pos, pEnd); 197 } 198 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), pEnd - pStart, points, PETSC_OWN_POINTER, perm)); 199 PetscFunctionReturn(PETSC_SUCCESS); 200 } 201 202 static PetscErrorCode DMPlexRemapCoordinates_Private(IS perm, PetscSection cs, Vec coordinates, PetscSection *csNew, Vec *coordinatesNew) 203 { 204 PetscScalar *coords, *coordsNew; 205 const PetscInt *pperm; 206 PetscInt pStart, pEnd, p; 207 const char *name; 208 209 PetscFunctionBegin; 210 PetscCall(PetscSectionPermute(cs, perm, csNew)); 211 PetscCall(VecDuplicate(coordinates, coordinatesNew)); 212 PetscCall(PetscObjectGetName((PetscObject)coordinates, &name)); 213 PetscCall(PetscObjectSetName((PetscObject)*coordinatesNew, name)); 214 PetscCall(VecGetArray(coordinates, &coords)); 215 PetscCall(VecGetArray(*coordinatesNew, &coordsNew)); 216 PetscCall(PetscSectionGetChart(*csNew, &pStart, &pEnd)); 217 PetscCall(ISGetIndices(perm, &pperm)); 218 for (p = pStart; p < pEnd; ++p) { 219 PetscInt dof, off, offNew, d; 220 221 PetscCall(PetscSectionGetDof(*csNew, p, &dof)); 222 PetscCall(PetscSectionGetOffset(cs, p, &off)); 223 PetscCall(PetscSectionGetOffset(*csNew, pperm[p], &offNew)); 224 for (d = 0; d < dof; ++d) coordsNew[offNew + d] = coords[off + d]; 225 } 226 PetscCall(ISRestoreIndices(perm, &pperm)); 227 PetscCall(VecRestoreArray(coordinates, &coords)); 228 PetscCall(VecRestoreArray(*coordinatesNew, &coordsNew)); 229 PetscFunctionReturn(PETSC_SUCCESS); 230 } 231 232 /*@ 233 DMPlexPermute - Reorder the mesh according to the input permutation 234 235 Collective on dm 236 237 Input Parameters: 238 + dm - The DMPlex object 239 - perm - The point permutation, perm[old point number] = new point number 240 241 Output Parameter: 242 . pdm - The permuted DM 243 244 Level: intermediate 245 246 .seealso: `MatPermute()` 247 @*/ 248 PetscErrorCode DMPlexPermute(DM dm, IS perm, DM *pdm) 249 { 250 DM_Plex *plex = (DM_Plex *)dm->data, *plexNew; 251 PetscInt dim, cdim; 252 const char *name; 253 254 PetscFunctionBegin; 255 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 256 PetscValidHeaderSpecific(perm, IS_CLASSID, 2); 257 PetscValidPointer(pdm, 3); 258 PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), pdm)); 259 PetscCall(DMSetType(*pdm, DMPLEX)); 260 PetscCall(PetscObjectGetName((PetscObject)dm, &name)); 261 PetscCall(PetscObjectSetName((PetscObject)*pdm, name)); 262 PetscCall(DMGetDimension(dm, &dim)); 263 PetscCall(DMSetDimension(*pdm, dim)); 264 PetscCall(DMGetCoordinateDim(dm, &cdim)); 265 PetscCall(DMSetCoordinateDim(*pdm, cdim)); 266 PetscCall(DMCopyDisc(dm, *pdm)); 267 if (dm->localSection) { 268 PetscSection section, sectionNew; 269 270 PetscCall(DMGetLocalSection(dm, §ion)); 271 PetscCall(PetscSectionPermute(section, perm, §ionNew)); 272 PetscCall(DMSetLocalSection(*pdm, sectionNew)); 273 PetscCall(PetscSectionDestroy(§ionNew)); 274 } 275 plexNew = (DM_Plex *)(*pdm)->data; 276 /* Ignore ltogmap, ltogmapb */ 277 /* Ignore sf, sectionSF */ 278 /* Ignore globalVertexNumbers, globalCellNumbers */ 279 /* Reorder labels */ 280 { 281 PetscInt numLabels, l; 282 DMLabel label, labelNew; 283 284 PetscCall(DMGetNumLabels(dm, &numLabels)); 285 for (l = 0; l < numLabels; ++l) { 286 PetscCall(DMGetLabelByNum(dm, l, &label)); 287 PetscCall(DMLabelPermute(label, perm, &labelNew)); 288 PetscCall(DMAddLabel(*pdm, labelNew)); 289 PetscCall(DMLabelDestroy(&labelNew)); 290 } 291 PetscCall(DMGetLabel(*pdm, "depth", &(*pdm)->depthLabel)); 292 if (plex->subpointMap) PetscCall(DMLabelPermute(plex->subpointMap, perm, &plexNew->subpointMap)); 293 } 294 if ((*pdm)->celltypeLabel) { 295 DMLabel ctLabel; 296 297 // Reset label for fast lookup 298 PetscCall(DMPlexGetCellTypeLabel(*pdm, &ctLabel)); 299 PetscCall(DMLabelMakeAllInvalid_Internal(ctLabel)); 300 } 301 /* Reorder topology */ 302 { 303 const PetscInt *pperm; 304 PetscInt n, pStart, pEnd, p; 305 306 PetscCall(PetscSectionDestroy(&plexNew->coneSection)); 307 PetscCall(PetscSectionPermute(plex->coneSection, perm, &plexNew->coneSection)); 308 PetscCall(PetscSectionGetStorageSize(plexNew->coneSection, &n)); 309 PetscCall(PetscMalloc1(n, &plexNew->cones)); 310 PetscCall(PetscMalloc1(n, &plexNew->coneOrientations)); 311 PetscCall(ISGetIndices(perm, &pperm)); 312 PetscCall(PetscSectionGetChart(plex->coneSection, &pStart, &pEnd)); 313 for (p = pStart; p < pEnd; ++p) { 314 PetscInt dof, off, offNew, d; 315 316 PetscCall(PetscSectionGetDof(plexNew->coneSection, pperm[p], &dof)); 317 PetscCall(PetscSectionGetOffset(plex->coneSection, p, &off)); 318 PetscCall(PetscSectionGetOffset(plexNew->coneSection, pperm[p], &offNew)); 319 for (d = 0; d < dof; ++d) { 320 plexNew->cones[offNew + d] = pperm[plex->cones[off + d]]; 321 plexNew->coneOrientations[offNew + d] = plex->coneOrientations[off + d]; 322 } 323 } 324 PetscCall(PetscSectionDestroy(&plexNew->supportSection)); 325 PetscCall(PetscSectionPermute(plex->supportSection, perm, &plexNew->supportSection)); 326 PetscCall(PetscSectionGetStorageSize(plexNew->supportSection, &n)); 327 PetscCall(PetscMalloc1(n, &plexNew->supports)); 328 PetscCall(PetscSectionGetChart(plex->supportSection, &pStart, &pEnd)); 329 for (p = pStart; p < pEnd; ++p) { 330 PetscInt dof, off, offNew, d; 331 332 PetscCall(PetscSectionGetDof(plexNew->supportSection, pperm[p], &dof)); 333 PetscCall(PetscSectionGetOffset(plex->supportSection, p, &off)); 334 PetscCall(PetscSectionGetOffset(plexNew->supportSection, pperm[p], &offNew)); 335 for (d = 0; d < dof; ++d) plexNew->supports[offNew + d] = pperm[plex->supports[off + d]]; 336 } 337 PetscCall(ISRestoreIndices(perm, &pperm)); 338 } 339 /* Remap coordinates */ 340 { 341 DM cdm, cdmNew; 342 PetscSection cs, csNew; 343 Vec coordinates, coordinatesNew; 344 345 PetscCall(DMGetCoordinateSection(dm, &cs)); 346 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 347 PetscCall(DMPlexRemapCoordinates_Private(perm, cs, coordinates, &csNew, &coordinatesNew)); 348 PetscCall(DMSetCoordinateSection(*pdm, PETSC_DETERMINE, csNew)); 349 PetscCall(DMSetCoordinatesLocal(*pdm, coordinatesNew)); 350 PetscCall(PetscSectionDestroy(&csNew)); 351 PetscCall(VecDestroy(&coordinatesNew)); 352 353 PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 354 if (cdm) { 355 PetscCall(DMGetCoordinateDM(*pdm, &cdm)); 356 PetscCall(DMClone(cdm, &cdmNew)); 357 PetscCall(DMSetCellCoordinateDM(*pdm, cdmNew)); 358 PetscCall(DMDestroy(&cdmNew)); 359 PetscCall(DMGetCellCoordinateSection(dm, &cs)); 360 PetscCall(DMGetCellCoordinatesLocal(dm, &coordinates)); 361 PetscCall(DMPlexRemapCoordinates_Private(perm, cs, coordinates, &csNew, &coordinatesNew)); 362 PetscCall(DMSetCellCoordinateSection(*pdm, PETSC_DETERMINE, csNew)); 363 PetscCall(DMSetCellCoordinatesLocal(*pdm, coordinatesNew)); 364 PetscCall(PetscSectionDestroy(&csNew)); 365 PetscCall(VecDestroy(&coordinatesNew)); 366 } 367 } 368 PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_TRUE, *pdm)); 369 (*pdm)->setupcalled = PETSC_TRUE; 370 PetscFunctionReturn(PETSC_SUCCESS); 371 } 372 373 PetscErrorCode DMPlexReorderSetDefault_Plex(DM dm, DMPlexReorderDefaultFlag reorder) 374 { 375 DM_Plex *mesh = (DM_Plex *)dm->data; 376 377 PetscFunctionBegin; 378 mesh->reorderDefault = reorder; 379 PetscFunctionReturn(PETSC_SUCCESS); 380 } 381 382 /*@ 383 DMPlexReorderSetDefault - Set flag indicating whether the DM should be reordered by default 384 385 Logically collective 386 387 Input Parameters: 388 + dm - The DM 389 - reorder - Flag for reordering 390 391 Level: intermediate 392 393 .seealso: `DMPlexReorderGetDefault()` 394 @*/ 395 PetscErrorCode DMPlexReorderSetDefault(DM dm, DMPlexReorderDefaultFlag reorder) 396 { 397 PetscFunctionBegin; 398 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 399 PetscTryMethod(dm, "DMPlexReorderSetDefault_C", (DM, DMPlexReorderDefaultFlag), (dm, reorder)); 400 PetscFunctionReturn(PETSC_SUCCESS); 401 } 402 403 PetscErrorCode DMPlexReorderGetDefault_Plex(DM dm, DMPlexReorderDefaultFlag *reorder) 404 { 405 DM_Plex *mesh = (DM_Plex *)dm->data; 406 407 PetscFunctionBegin; 408 *reorder = mesh->reorderDefault; 409 PetscFunctionReturn(PETSC_SUCCESS); 410 } 411 412 /*@ 413 DMPlexReorderGetDefault - Get flag indicating whether the DM should be reordered by default 414 415 Not collective 416 417 Input Parameter: 418 . dm - The DM 419 420 Output Parameter: 421 . reorder - Flag for reordering 422 423 Level: intermediate 424 425 .seealso: `DMPlexReorderSetDefault()` 426 @*/ 427 PetscErrorCode DMPlexReorderGetDefault(DM dm, DMPlexReorderDefaultFlag *reorder) 428 { 429 PetscFunctionBegin; 430 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 431 PetscValidPointer(reorder, 2); 432 PetscUseMethod(dm, "DMPlexReorderGetDefault_C", (DM, DMPlexReorderDefaultFlag *), (dm, reorder)); 433 PetscFunctionReturn(PETSC_SUCCESS); 434 } 435