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(0); 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 *values; 97 PetscInt numValues, numPoints = 0; 98 PetscInt *sperm, *vsize, *voff, v; 99 100 PetscCall(DMLabelGetValueIS(label, &valueIS)); 101 PetscCall(ISSort(valueIS)); 102 PetscCall(ISGetLocalSize(valueIS, &numValues)); 103 PetscCall(ISGetIndices(valueIS, &values)); 104 PetscCall(PetscCalloc3(numCells,&sperm,numValues,&vsize,numValues+1,&voff)); 105 for (v = 0; v < numValues; ++v) { 106 PetscCall(DMLabelGetStratumSize(label, values[v], &vsize[v])); 107 if (v < numValues-1) voff[v+2] += vsize[v] + voff[v+1]; 108 numPoints += vsize[v]; 109 } 110 PetscCheck(numPoints == numCells, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label only covers %" PetscInt_FMT " cells != %" PetscInt_FMT " total cells", numPoints, numCells); 111 for (c = 0; c < numCells; ++c) { 112 const PetscInt oldc = cperm[c]; 113 PetscInt val, vloc; 114 115 PetscCall(DMLabelGetValue(label, oldc, &val)); 116 PetscCheck(val != -1,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cell %" PetscInt_FMT " not present in label", oldc); 117 PetscCall(PetscFindInt(val, numValues, values, &vloc)); 118 PetscCheck(vloc >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %" PetscInt_FMT " not present label", val); 119 sperm[voff[vloc+1]++] = oldc; 120 } 121 for (v = 0; v < numValues; ++v) { 122 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]); 123 } 124 PetscCall(ISRestoreIndices(valueIS, &values)); 125 PetscCall(ISDestroy(&valueIS)); 126 PetscCall(PetscArraycpy(cperm, sperm, numCells)); 127 PetscCall(PetscFree3(sperm, vsize, voff)); 128 } 129 /* Construct closure */ 130 PetscCall(DMPlexCreateOrderingClosure_Static(dm, numCells, cperm, &clperm, &invclperm)); 131 PetscCall(PetscFree3(cperm,mask,xls)); 132 PetscCall(PetscFree(clperm)); 133 /* Invert permutation */ 134 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 135 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd-pStart, invclperm, PETSC_OWN_POINTER, perm)); 136 PetscFunctionReturn(0); 137 } 138 139 /*@ 140 DMPlexGetOrdering1D - Reorder the vertices so that the mesh is in a line 141 142 Collective on dm 143 144 Input Parameter: 145 . dm - The DMPlex object 146 147 Output Parameter: 148 . perm - The point permutation as an IS, perm[old point number] = new point number 149 150 Level: intermediate 151 152 .seealso: `DMPlexGetOrdering()`, `DMPlexPermute()`, `MatGetOrdering()` 153 @*/ 154 PetscErrorCode DMPlexGetOrdering1D(DM dm, IS *perm) 155 { 156 PetscInt *points; 157 const PetscInt *support, *cone; 158 PetscInt dim, pStart, pEnd, cStart, cEnd, c, vStart, vEnd, v, suppSize, lastCell = 0; 159 160 PetscFunctionBegin; 161 PetscCall(DMGetDimension(dm, &dim)); 162 PetscCheck(dim == 1, PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Input mesh must be one dimensional, not %" PetscInt_FMT, dim); 163 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 164 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 165 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 166 PetscCall(PetscMalloc1(pEnd-pStart, &points)); 167 for (c = cStart; c < cEnd; ++c) points[c] = c; 168 for (v = vStart; v < vEnd; ++v) points[v] = v; 169 for (v = vStart; v < vEnd; ++v) { 170 PetscCall(DMPlexGetSupportSize(dm, v, &suppSize)); 171 PetscCall(DMPlexGetSupport(dm, v, &support)); 172 if (suppSize == 1) {lastCell = support[0]; break;} 173 } 174 if (v < vEnd) { 175 PetscInt pos = cEnd; 176 177 points[v] = pos++; 178 while (lastCell >= cStart) { 179 PetscCall(DMPlexGetCone(dm, lastCell, &cone)); 180 if (cone[0] == v) v = cone[1]; 181 else v = cone[0]; 182 PetscCall(DMPlexGetSupport(dm, v, &support)); 183 PetscCall(DMPlexGetSupportSize(dm, v, &suppSize)); 184 if (suppSize == 1) {lastCell = -1;} 185 else { 186 if (support[0] == lastCell) lastCell = support[1]; 187 else lastCell = support[0]; 188 } 189 points[v] = pos++; 190 } 191 PetscCheck(pos == pEnd, PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Last vertex was %" PetscInt_FMT ", not %" PetscInt_FMT, pos, pEnd); 192 } 193 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd-pStart, points, PETSC_OWN_POINTER, perm)); 194 PetscFunctionReturn(0); 195 } 196 197 static PetscErrorCode DMPlexRemapCoordinates_Private(IS perm, PetscSection cs, Vec coordinates, PetscSection *csNew, Vec *coordinatesNew) 198 { 199 PetscScalar *coords, *coordsNew; 200 const PetscInt *pperm; 201 PetscInt pStart, pEnd, p; 202 const char *name; 203 204 PetscFunctionBegin; 205 PetscCall(PetscSectionPermute(cs, perm, csNew)); 206 PetscCall(VecDuplicate(coordinates, coordinatesNew)); 207 PetscCall(PetscObjectGetName((PetscObject) coordinates, &name)); 208 PetscCall(PetscObjectSetName((PetscObject) *coordinatesNew, name)); 209 PetscCall(VecGetArray(coordinates, &coords)); 210 PetscCall(VecGetArray(*coordinatesNew, &coordsNew)); 211 PetscCall(PetscSectionGetChart(*csNew, &pStart, &pEnd)); 212 PetscCall(ISGetIndices(perm, &pperm)); 213 for (p = pStart; p < pEnd; ++p) { 214 PetscInt dof, off, offNew, d; 215 216 PetscCall(PetscSectionGetDof(*csNew, p, &dof)); 217 PetscCall(PetscSectionGetOffset(cs, p, &off)); 218 PetscCall(PetscSectionGetOffset(*csNew, pperm[p], &offNew)); 219 for (d = 0; d < dof; ++d) coordsNew[offNew+d] = coords[off+d]; 220 } 221 PetscCall(ISRestoreIndices(perm, &pperm)); 222 PetscCall(VecRestoreArray(coordinates, &coords)); 223 PetscCall(VecRestoreArray(*coordinatesNew, &coordsNew)); 224 PetscFunctionReturn(0); 225 } 226 227 /*@ 228 DMPlexPermute - Reorder the mesh according to the input permutation 229 230 Collective on dm 231 232 Input Parameters: 233 + dm - The DMPlex object 234 - perm - The point permutation, perm[old point number] = new point number 235 236 Output Parameter: 237 . pdm - The permuted DM 238 239 Level: intermediate 240 241 .seealso: `MatPermute()` 242 @*/ 243 PetscErrorCode DMPlexPermute(DM dm, IS perm, DM *pdm) 244 { 245 DM_Plex *plex = (DM_Plex *) dm->data, *plexNew; 246 PetscInt dim, cdim; 247 const char *name; 248 249 PetscFunctionBegin; 250 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 251 PetscValidHeaderSpecific(perm, IS_CLASSID, 2); 252 PetscValidPointer(pdm, 3); 253 PetscCall(DMCreate(PetscObjectComm((PetscObject) dm), pdm)); 254 PetscCall(DMSetType(*pdm, DMPLEX)); 255 PetscCall(PetscObjectGetName((PetscObject) dm, &name)); 256 PetscCall(PetscObjectSetName((PetscObject) *pdm, name)); 257 PetscCall(DMGetDimension(dm, &dim)); 258 PetscCall(DMSetDimension(*pdm, dim)); 259 PetscCall(DMGetCoordinateDim(dm, &cdim)); 260 PetscCall(DMSetCoordinateDim(*pdm, cdim)); 261 PetscCall(DMCopyDisc(dm, *pdm)); 262 if (dm->localSection) { 263 PetscSection section, sectionNew; 264 265 PetscCall(DMGetLocalSection(dm, §ion)); 266 PetscCall(PetscSectionPermute(section, perm, §ionNew)); 267 PetscCall(DMSetLocalSection(*pdm, sectionNew)); 268 PetscCall(PetscSectionDestroy(§ionNew)); 269 } 270 plexNew = (DM_Plex *) (*pdm)->data; 271 /* Ignore ltogmap, ltogmapb */ 272 /* Ignore sf, sectionSF */ 273 /* Ignore globalVertexNumbers, globalCellNumbers */ 274 /* Reorder labels */ 275 { 276 PetscInt numLabels, l; 277 DMLabel label, labelNew; 278 279 PetscCall(DMGetNumLabels(dm, &numLabels)); 280 for (l = 0; l < numLabels; ++l) { 281 PetscCall(DMGetLabelByNum(dm, l, &label)); 282 PetscCall(DMLabelPermute(label, perm, &labelNew)); 283 PetscCall(DMAddLabel(*pdm, labelNew)); 284 PetscCall(DMLabelDestroy(&labelNew)); 285 } 286 PetscCall(DMGetLabel(*pdm, "depth", &(*pdm)->depthLabel)); 287 if (plex->subpointMap) PetscCall(DMLabelPermute(plex->subpointMap, perm, &plexNew->subpointMap)); 288 } 289 if ((*pdm)->celltypeLabel) { 290 DMLabel ctLabel; 291 292 // Reset label for fast lookup 293 PetscCall(DMPlexGetCellTypeLabel(*pdm, &ctLabel)); 294 PetscCall(DMLabelMakeAllInvalid_Internal(ctLabel)); 295 } 296 /* Reorder topology */ 297 { 298 const PetscInt *pperm; 299 PetscInt n, pStart, pEnd, p; 300 301 PetscCall(PetscSectionDestroy(&plexNew->coneSection)); 302 PetscCall(PetscSectionPermute(plex->coneSection, perm, &plexNew->coneSection)); 303 PetscCall(PetscSectionGetStorageSize(plexNew->coneSection, &n)); 304 PetscCall(PetscMalloc1(n, &plexNew->cones)); 305 PetscCall(PetscMalloc1(n, &plexNew->coneOrientations)); 306 PetscCall(ISGetIndices(perm, &pperm)); 307 PetscCall(PetscSectionGetChart(plex->coneSection, &pStart, &pEnd)); 308 for (p = pStart; p < pEnd; ++p) { 309 PetscInt dof, off, offNew, d; 310 311 PetscCall(PetscSectionGetDof(plexNew->coneSection, pperm[p], &dof)); 312 PetscCall(PetscSectionGetOffset(plex->coneSection, p, &off)); 313 PetscCall(PetscSectionGetOffset(plexNew->coneSection, pperm[p], &offNew)); 314 for (d = 0; d < dof; ++d) { 315 plexNew->cones[offNew+d] = pperm[plex->cones[off+d]]; 316 plexNew->coneOrientations[offNew+d] = plex->coneOrientations[off+d]; 317 } 318 } 319 PetscCall(PetscSectionDestroy(&plexNew->supportSection)); 320 PetscCall(PetscSectionPermute(plex->supportSection, perm, &plexNew->supportSection)); 321 PetscCall(PetscSectionGetStorageSize(plexNew->supportSection, &n)); 322 PetscCall(PetscMalloc1(n, &plexNew->supports)); 323 PetscCall(PetscSectionGetChart(plex->supportSection, &pStart, &pEnd)); 324 for (p = pStart; p < pEnd; ++p) { 325 PetscInt dof, off, offNew, d; 326 327 PetscCall(PetscSectionGetDof(plexNew->supportSection, pperm[p], &dof)); 328 PetscCall(PetscSectionGetOffset(plex->supportSection, p, &off)); 329 PetscCall(PetscSectionGetOffset(plexNew->supportSection, pperm[p], &offNew)); 330 for (d = 0; d < dof; ++d) { 331 plexNew->supports[offNew+d] = pperm[plex->supports[off+d]]; 332 } 333 } 334 PetscCall(ISRestoreIndices(perm, &pperm)); 335 } 336 /* Remap coordinates */ 337 { 338 DM cdm, cdmNew; 339 PetscSection cs, csNew; 340 Vec coordinates, coordinatesNew; 341 342 PetscCall(DMGetCoordinateSection(dm, &cs)); 343 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 344 PetscCall(DMPlexRemapCoordinates_Private(perm, cs, coordinates, &csNew, &coordinatesNew)); 345 PetscCall(DMSetCoordinateSection(*pdm, PETSC_DETERMINE, csNew)); 346 PetscCall(DMSetCoordinatesLocal(*pdm, coordinatesNew)); 347 PetscCall(PetscSectionDestroy(&csNew)); 348 PetscCall(VecDestroy(&coordinatesNew)); 349 350 PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 351 if (cdm) { 352 PetscCall(DMGetCoordinateDM(*pdm, &cdm)); 353 PetscCall(DMClone(cdm, &cdmNew)); 354 PetscCall(DMSetCellCoordinateDM(*pdm, cdmNew)); 355 PetscCall(DMDestroy(&cdmNew)); 356 PetscCall(DMGetCellCoordinateSection(dm, &cs)); 357 PetscCall(DMGetCellCoordinatesLocal(dm, &coordinates)); 358 PetscCall(DMPlexRemapCoordinates_Private(perm, cs, coordinates, &csNew, &coordinatesNew)); 359 PetscCall(DMSetCellCoordinateSection(*pdm, PETSC_DETERMINE, csNew)); 360 PetscCall(DMSetCellCoordinatesLocal(*pdm, coordinatesNew)); 361 PetscCall(PetscSectionDestroy(&csNew)); 362 PetscCall(VecDestroy(&coordinatesNew)); 363 } 364 } 365 PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_TRUE, *pdm)); 366 (*pdm)->setupcalled = PETSC_TRUE; 367 PetscFunctionReturn(0); 368 } 369 370 PetscErrorCode DMPlexReorderSetDefault_Plex(DM dm, DMPlexReorderDefaultFlag reorder) 371 { 372 DM_Plex *mesh = (DM_Plex *) dm->data; 373 374 PetscFunctionBegin; 375 mesh->reorderDefault = reorder; 376 PetscFunctionReturn(0); 377 } 378 379 /*@ 380 DMPlexReorderSetDefault - Set flag indicating whether the DM should be reordered by default 381 382 Logically collective 383 384 Input Parameters: 385 + dm - The DM 386 - reorder - Flag for reordering 387 388 Level: intermediate 389 390 .seealso: `DMPlexReorderGetDefault()` 391 @*/ 392 PetscErrorCode DMPlexReorderSetDefault(DM dm, DMPlexReorderDefaultFlag reorder) 393 { 394 PetscFunctionBegin; 395 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 396 PetscTryMethod(dm,"DMPlexReorderSetDefault_C",(DM,DMPlexReorderDefaultFlag),(dm,reorder)); 397 PetscFunctionReturn(0); 398 } 399 400 PetscErrorCode DMPlexReorderGetDefault_Plex(DM dm, DMPlexReorderDefaultFlag *reorder) 401 { 402 DM_Plex *mesh = (DM_Plex *) dm->data; 403 404 PetscFunctionBegin; 405 *reorder = mesh->reorderDefault; 406 PetscFunctionReturn(0); 407 } 408 409 /*@ 410 DMPlexReorderGetDefault - Get flag indicating whether the DM should be reordered by default 411 412 Not collective 413 414 Input Parameter: 415 . dm - The DM 416 417 Output Parameter: 418 . reorder - Flag for reordering 419 420 Level: intermediate 421 422 .seealso: `DMPlexReorderSetDefault()` 423 @*/ 424 PetscErrorCode DMPlexReorderGetDefault(DM dm, DMPlexReorderDefaultFlag *reorder) 425 { 426 PetscFunctionBegin; 427 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 428 PetscValidPointer(reorder, 2); 429 PetscUseMethod(dm,"DMPlexReorderGetDefault_C",(DM,DMPlexReorderDefaultFlag*),(dm,reorder)); 430 PetscFunctionReturn(0); 431 } 432