1 #include <petscdm.h> 2 #include <petsc/private/dmlabelimpl.h> /*I "petscdmlabel.h" I*/ 3 #include <petsc/private/sectionimpl.h> /*I "petscsection.h" I*/ 4 #include <petscsf.h> 5 #include <petscsection.h> 6 7 PetscFunctionList DMLabelList = NULL; 8 PetscBool DMLabelRegisterAllCalled = PETSC_FALSE; 9 10 /*@ 11 DMLabelCreate - Create a `DMLabel` object, which is a multimap 12 13 Collective 14 15 Input Parameters: 16 + comm - The communicator, usually `PETSC_COMM_SELF` 17 - name - The label name 18 19 Output Parameter: 20 . label - The `DMLabel` 21 22 Level: beginner 23 24 Notes: 25 The label name is actually usually the `PetscObject` name. 26 One can get/set it with `PetscObjectGetName()`/`PetscObjectSetName()`. 27 28 .seealso: `DMLabel`, `DM`, `DMLabelDestroy()` 29 @*/ 30 PetscErrorCode DMLabelCreate(MPI_Comm comm, const char name[], DMLabel *label) 31 { 32 PetscFunctionBegin; 33 PetscAssertPointer(label, 3); 34 PetscCall(DMInitializePackage()); 35 36 PetscCall(PetscHeaderCreate(*label, DMLABEL_CLASSID, "DMLabel", "DMLabel", "DM", comm, DMLabelDestroy, DMLabelView)); 37 38 (*label)->numStrata = 0; 39 (*label)->defaultValue = -1; 40 (*label)->stratumValues = NULL; 41 (*label)->validIS = NULL; 42 (*label)->stratumSizes = NULL; 43 (*label)->points = NULL; 44 (*label)->ht = NULL; 45 (*label)->pStart = -1; 46 (*label)->pEnd = -1; 47 (*label)->bt = NULL; 48 PetscCall(PetscHMapICreate(&(*label)->hmap)); 49 PetscCall(PetscObjectSetName((PetscObject)*label, name)); 50 PetscCall(DMLabelSetType(*label, DMLABELCONCRETE)); 51 PetscFunctionReturn(PETSC_SUCCESS); 52 } 53 54 /*@ 55 DMLabelSetUp - SetUp a `DMLabel` object 56 57 Collective 58 59 Input Parameters: 60 . label - The `DMLabel` 61 62 Level: intermediate 63 64 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelDestroy()` 65 @*/ 66 PetscErrorCode DMLabelSetUp(DMLabel label) 67 { 68 PetscFunctionBegin; 69 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 70 PetscTryTypeMethod(label, setup); 71 PetscFunctionReturn(PETSC_SUCCESS); 72 } 73 74 /* 75 DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format 76 77 Not collective 78 79 Input parameter: 80 + label - The `DMLabel` 81 - v - The stratum value 82 83 Output parameter: 84 . label - The `DMLabel` with stratum in sorted list format 85 86 Level: developer 87 88 .seealso: `DMLabel`, `DM`, `DMLabelCreate()` 89 */ 90 static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v) 91 { 92 IS is; 93 PetscInt off = 0, *pointArray, p; 94 95 PetscFunctionBegin; 96 if ((PetscLikely(v >= 0 && v < label->numStrata) && label->validIS[v]) || label->readonly) PetscFunctionReturn(PETSC_SUCCESS); 97 PetscCheck(v >= 0 && v < label->numStrata, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %" PetscInt_FMT " in DMLabelMakeValid_Private", v); 98 PetscCall(PetscHSetIGetSize(label->ht[v], &label->stratumSizes[v])); 99 PetscCall(PetscMalloc1(label->stratumSizes[v], &pointArray)); 100 PetscCall(PetscHSetIGetElems(label->ht[v], &off, pointArray)); 101 PetscCall(PetscHSetIClear(label->ht[v])); 102 PetscCall(PetscSortInt(label->stratumSizes[v], pointArray)); 103 if (label->bt) { 104 for (p = 0; p < label->stratumSizes[v]; ++p) { 105 const PetscInt point = pointArray[p]; 106 PetscCheck(!(point < label->pStart) && !(point >= label->pEnd), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %" PetscInt_FMT " is not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, label->pStart, label->pEnd); 107 PetscCall(PetscBTSet(label->bt, point - label->pStart)); 108 } 109 } 110 if (label->stratumSizes[v] > 0 && pointArray[label->stratumSizes[v] - 1] == pointArray[0] + label->stratumSizes[v] - 1) { 111 PetscCall(ISCreateStride(PETSC_COMM_SELF, label->stratumSizes[v], pointArray[0], 1, &is)); 112 PetscCall(PetscFree(pointArray)); 113 } else { 114 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, label->stratumSizes[v], pointArray, PETSC_OWN_POINTER, &is)); 115 } 116 PetscCall(PetscObjectSetName((PetscObject)is, "indices")); 117 label->points[v] = is; 118 label->validIS[v] = PETSC_TRUE; 119 PetscCall(PetscObjectStateIncrease((PetscObject)label)); 120 PetscFunctionReturn(PETSC_SUCCESS); 121 } 122 123 /* 124 DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format 125 126 Not Collective 127 128 Input parameter: 129 . label - The `DMLabel` 130 131 Output parameter: 132 . label - The `DMLabel` with all strata in sorted list format 133 134 Level: developer 135 136 .seealso: `DMLabel`, `DM`, `DMLabelCreate()` 137 */ 138 static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label) 139 { 140 PetscInt v; 141 142 PetscFunctionBegin; 143 for (v = 0; v < label->numStrata; v++) PetscCall(DMLabelMakeValid_Private(label, v)); 144 PetscFunctionReturn(PETSC_SUCCESS); 145 } 146 147 /* 148 DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format 149 150 Not Collective 151 152 Input parameter: 153 + label - The `DMLabel` 154 - v - The stratum value 155 156 Output parameter: 157 . label - The `DMLabel` with stratum in hash format 158 159 Level: developer 160 161 .seealso: `DMLabel`, `DM`, `DMLabelCreate()` 162 */ 163 static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v) 164 { 165 PetscInt p; 166 const PetscInt *points; 167 168 PetscFunctionBegin; 169 if ((PetscLikely(v >= 0 && v < label->numStrata) && !label->validIS[v]) || label->readonly) PetscFunctionReturn(PETSC_SUCCESS); 170 PetscCheck(v >= 0 && v < label->numStrata, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %" PetscInt_FMT " in DMLabelMakeInvalid_Private", v); 171 if (label->points[v]) { 172 PetscCall(ISGetIndices(label->points[v], &points)); 173 for (p = 0; p < label->stratumSizes[v]; ++p) PetscCall(PetscHSetIAdd(label->ht[v], points[p])); 174 PetscCall(ISRestoreIndices(label->points[v], &points)); 175 PetscCall(ISDestroy(&label->points[v])); 176 } 177 label->validIS[v] = PETSC_FALSE; 178 PetscFunctionReturn(PETSC_SUCCESS); 179 } 180 181 PetscErrorCode DMLabelMakeAllInvalid_Internal(DMLabel label) 182 { 183 PetscInt v; 184 185 PetscFunctionBegin; 186 for (v = 0; v < label->numStrata; v++) PetscCall(DMLabelMakeInvalid_Private(label, v)); 187 PetscFunctionReturn(PETSC_SUCCESS); 188 } 189 190 #if !defined(DMLABEL_LOOKUP_THRESHOLD) 191 #define DMLABEL_LOOKUP_THRESHOLD 16 192 #endif 193 194 PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index) 195 { 196 PetscInt v; 197 198 PetscFunctionBegin; 199 *index = -1; 200 if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD || label->readonly) { 201 for (v = 0; v < label->numStrata; ++v) 202 if (label->stratumValues[v] == value) { 203 *index = v; 204 break; 205 } 206 } else { 207 PetscCall(PetscHMapIGet(label->hmap, value, index)); 208 } 209 if (PetscDefined(USE_DEBUG) && !label->readonly) { /* Check strata hash map consistency */ 210 PetscInt len, loc = -1; 211 PetscCall(PetscHMapIGetSize(label->hmap, &len)); 212 PetscCheck(len == label->numStrata, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map size"); 213 if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) { 214 PetscCall(PetscHMapIGet(label->hmap, value, &loc)); 215 } else { 216 for (v = 0; v < label->numStrata; ++v) 217 if (label->stratumValues[v] == value) { 218 loc = v; 219 break; 220 } 221 } 222 PetscCheck(loc == *index, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map lookup"); 223 } 224 PetscFunctionReturn(PETSC_SUCCESS); 225 } 226 227 static inline PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index) 228 { 229 PetscInt v; 230 PetscInt *tmpV; 231 PetscInt *tmpS; 232 PetscHSetI *tmpH, ht; 233 IS *tmpP, is; 234 PetscBool *tmpB; 235 PetscHMapI hmap = label->hmap; 236 237 PetscFunctionBegin; 238 v = label->numStrata; 239 tmpV = label->stratumValues; 240 tmpS = label->stratumSizes; 241 tmpH = label->ht; 242 tmpP = label->points; 243 tmpB = label->validIS; 244 { /* TODO: PetscRealloc() is broken, use malloc+memcpy+free */ 245 PetscInt *oldV = tmpV; 246 PetscInt *oldS = tmpS; 247 PetscHSetI *oldH = tmpH; 248 IS *oldP = tmpP; 249 PetscBool *oldB = tmpB; 250 PetscCall(PetscMalloc((v + 1) * sizeof(*tmpV), &tmpV)); 251 PetscCall(PetscMalloc((v + 1) * sizeof(*tmpS), &tmpS)); 252 PetscCall(PetscCalloc((v + 1) * sizeof(*tmpH), &tmpH)); 253 PetscCall(PetscCalloc((v + 1) * sizeof(*tmpP), &tmpP)); 254 PetscCall(PetscMalloc((v + 1) * sizeof(*tmpB), &tmpB)); 255 PetscCall(PetscArraycpy(tmpV, oldV, v)); 256 PetscCall(PetscArraycpy(tmpS, oldS, v)); 257 PetscCall(PetscArraycpy(tmpH, oldH, v)); 258 PetscCall(PetscArraycpy(tmpP, oldP, v)); 259 PetscCall(PetscArraycpy(tmpB, oldB, v)); 260 PetscCall(PetscFree(oldV)); 261 PetscCall(PetscFree(oldS)); 262 PetscCall(PetscFree(oldH)); 263 PetscCall(PetscFree(oldP)); 264 PetscCall(PetscFree(oldB)); 265 } 266 label->numStrata = v + 1; 267 label->stratumValues = tmpV; 268 label->stratumSizes = tmpS; 269 label->ht = tmpH; 270 label->points = tmpP; 271 label->validIS = tmpB; 272 PetscCall(PetscHSetICreate(&ht)); 273 PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &is)); 274 PetscCall(PetscHMapISet(hmap, value, v)); 275 tmpV[v] = value; 276 tmpS[v] = 0; 277 tmpH[v] = ht; 278 tmpP[v] = is; 279 tmpB[v] = PETSC_TRUE; 280 PetscCall(PetscObjectStateIncrease((PetscObject)label)); 281 *index = v; 282 PetscFunctionReturn(PETSC_SUCCESS); 283 } 284 285 static inline PetscErrorCode DMLabelLookupAddStratum(DMLabel label, PetscInt value, PetscInt *index) 286 { 287 PetscFunctionBegin; 288 PetscCall(DMLabelLookupStratum(label, value, index)); 289 if (*index < 0) PetscCall(DMLabelNewStratum(label, value, index)); 290 PetscFunctionReturn(PETSC_SUCCESS); 291 } 292 293 PetscErrorCode DMLabelGetStratumSize_Private(DMLabel label, PetscInt v, PetscInt *size) 294 { 295 PetscFunctionBegin; 296 *size = 0; 297 if (v < 0) PetscFunctionReturn(PETSC_SUCCESS); 298 if (label->readonly || label->validIS[v]) { 299 *size = label->stratumSizes[v]; 300 } else { 301 PetscCall(PetscHSetIGetSize(label->ht[v], size)); 302 } 303 PetscFunctionReturn(PETSC_SUCCESS); 304 } 305 306 /*@ 307 DMLabelAddStratum - Adds a new stratum value in a `DMLabel` 308 309 Input Parameters: 310 + label - The `DMLabel` 311 - value - The stratum value 312 313 Level: beginner 314 315 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelDestroy()` 316 @*/ 317 PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value) 318 { 319 PetscInt v; 320 321 PetscFunctionBegin; 322 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 323 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 324 PetscCall(DMLabelLookupAddStratum(label, value, &v)); 325 PetscFunctionReturn(PETSC_SUCCESS); 326 } 327 328 /*@ 329 DMLabelAddStrata - Adds new stratum values in a `DMLabel` 330 331 Not Collective 332 333 Input Parameters: 334 + label - The `DMLabel` 335 . numStrata - The number of stratum values 336 - stratumValues - The stratum values 337 338 Level: beginner 339 340 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelDestroy()` 341 @*/ 342 PetscErrorCode DMLabelAddStrata(DMLabel label, PetscInt numStrata, const PetscInt stratumValues[]) 343 { 344 PetscInt *values, v; 345 346 PetscFunctionBegin; 347 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 348 if (numStrata) PetscAssertPointer(stratumValues, 3); 349 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 350 PetscCall(PetscMalloc1(numStrata, &values)); 351 PetscCall(PetscArraycpy(values, stratumValues, numStrata)); 352 PetscCall(PetscSortRemoveDupsInt(&numStrata, values)); 353 if (!label->numStrata) { /* Fast preallocation */ 354 PetscInt *tmpV; 355 PetscInt *tmpS; 356 PetscHSetI *tmpH, ht; 357 IS *tmpP, is; 358 PetscBool *tmpB; 359 PetscHMapI hmap = label->hmap; 360 361 PetscCall(PetscMalloc1(numStrata, &tmpV)); 362 PetscCall(PetscMalloc1(numStrata, &tmpS)); 363 PetscCall(PetscCalloc1(numStrata, &tmpH)); 364 PetscCall(PetscCalloc1(numStrata, &tmpP)); 365 PetscCall(PetscMalloc1(numStrata, &tmpB)); 366 label->numStrata = numStrata; 367 label->stratumValues = tmpV; 368 label->stratumSizes = tmpS; 369 label->ht = tmpH; 370 label->points = tmpP; 371 label->validIS = tmpB; 372 for (v = 0; v < numStrata; ++v) { 373 PetscCall(PetscHSetICreate(&ht)); 374 PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &is)); 375 PetscCall(PetscHMapISet(hmap, values[v], v)); 376 tmpV[v] = values[v]; 377 tmpS[v] = 0; 378 tmpH[v] = ht; 379 tmpP[v] = is; 380 tmpB[v] = PETSC_TRUE; 381 } 382 PetscCall(PetscObjectStateIncrease((PetscObject)label)); 383 } else { 384 for (v = 0; v < numStrata; ++v) PetscCall(DMLabelAddStratum(label, values[v])); 385 } 386 PetscCall(PetscFree(values)); 387 PetscFunctionReturn(PETSC_SUCCESS); 388 } 389 390 /*@ 391 DMLabelAddStrataIS - Adds new stratum values in a `DMLabel` 392 393 Not Collective 394 395 Input Parameters: 396 + label - The `DMLabel` 397 - valueIS - Index set with stratum values 398 399 Level: beginner 400 401 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelDestroy()` 402 @*/ 403 PetscErrorCode DMLabelAddStrataIS(DMLabel label, IS valueIS) 404 { 405 PetscInt numStrata; 406 const PetscInt *stratumValues; 407 408 PetscFunctionBegin; 409 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 410 PetscValidHeaderSpecific(valueIS, IS_CLASSID, 2); 411 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 412 PetscCall(ISGetLocalSize(valueIS, &numStrata)); 413 PetscCall(ISGetIndices(valueIS, &stratumValues)); 414 PetscCall(DMLabelAddStrata(label, numStrata, stratumValues)); 415 PetscFunctionReturn(PETSC_SUCCESS); 416 } 417 418 static PetscErrorCode DMLabelView_Concrete_Ascii(DMLabel label, PetscViewer viewer) 419 { 420 PetscInt v; 421 PetscMPIInt rank; 422 423 PetscFunctionBegin; 424 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank)); 425 PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 426 if (label) { 427 const char *name; 428 429 PetscCall(PetscObjectGetName((PetscObject)label, &name)); 430 PetscCall(PetscViewerASCIIPrintf(viewer, "Label '%s':\n", name)); 431 if (label->bt) PetscCall(PetscViewerASCIIPrintf(viewer, " Index has been calculated in [%" PetscInt_FMT ", %" PetscInt_FMT ")\n", label->pStart, label->pEnd)); 432 for (v = 0; v < label->numStrata; ++v) { 433 const PetscInt value = label->stratumValues[v]; 434 const PetscInt *points; 435 PetscInt p; 436 437 PetscCall(ISGetIndices(label->points[v], &points)); 438 for (p = 0; p < label->stratumSizes[v]; ++p) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %" PetscInt_FMT " (%" PetscInt_FMT ")\n", rank, points[p], value)); 439 PetscCall(ISRestoreIndices(label->points[v], &points)); 440 } 441 } 442 PetscCall(PetscViewerFlush(viewer)); 443 PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 444 PetscFunctionReturn(PETSC_SUCCESS); 445 } 446 447 static PetscErrorCode DMLabelView_Concrete(DMLabel label, PetscViewer viewer) 448 { 449 PetscBool iascii; 450 451 PetscFunctionBegin; 452 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 453 if (iascii) PetscCall(DMLabelView_Concrete_Ascii(label, viewer)); 454 PetscFunctionReturn(PETSC_SUCCESS); 455 } 456 457 /*@ 458 DMLabelView - View the label 459 460 Collective 461 462 Input Parameters: 463 + label - The `DMLabel` 464 - viewer - The `PetscViewer` 465 466 Level: intermediate 467 468 .seealso: `DMLabel`, `PetscViewer`, `DM`, `DMLabelCreate()`, `DMLabelDestroy()` 469 @*/ 470 PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer) 471 { 472 PetscFunctionBegin; 473 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 474 if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)label), &viewer)); 475 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 476 PetscCall(DMLabelMakeAllValid_Private(label)); 477 PetscUseTypeMethod(label, view, viewer); 478 PetscFunctionReturn(PETSC_SUCCESS); 479 } 480 481 /*@ 482 DMLabelReset - Destroys internal data structures in a `DMLabel` 483 484 Not Collective 485 486 Input Parameter: 487 . label - The `DMLabel` 488 489 Level: beginner 490 491 .seealso: `DMLabel`, `DM`, `DMLabelDestroy()`, `DMLabelCreate()` 492 @*/ 493 PetscErrorCode DMLabelReset(DMLabel label) 494 { 495 PetscInt v; 496 497 PetscFunctionBegin; 498 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 499 for (v = 0; v < label->numStrata; ++v) { 500 if (label->ht[v]) PetscCall(PetscHSetIDestroy(&label->ht[v])); 501 PetscCall(ISDestroy(&label->points[v])); 502 } 503 label->numStrata = 0; 504 PetscCall(PetscFree(label->stratumValues)); 505 PetscCall(PetscFree(label->stratumSizes)); 506 PetscCall(PetscFree(label->ht)); 507 PetscCall(PetscFree(label->points)); 508 PetscCall(PetscFree(label->validIS)); 509 PetscCall(PetscHMapIReset(label->hmap)); 510 label->pStart = -1; 511 label->pEnd = -1; 512 PetscCall(PetscBTDestroy(&label->bt)); 513 PetscFunctionReturn(PETSC_SUCCESS); 514 } 515 516 /*@ 517 DMLabelDestroy - Destroys a `DMLabel` 518 519 Collective 520 521 Input Parameter: 522 . label - The `DMLabel` 523 524 Level: beginner 525 526 .seealso: `DMLabel`, `DM`, `DMLabelReset()`, `DMLabelCreate()` 527 @*/ 528 PetscErrorCode DMLabelDestroy(DMLabel *label) 529 { 530 PetscFunctionBegin; 531 if (!*label) PetscFunctionReturn(PETSC_SUCCESS); 532 PetscValidHeaderSpecific(*label, DMLABEL_CLASSID, 1); 533 if (--((PetscObject)*label)->refct > 0) { 534 *label = NULL; 535 PetscFunctionReturn(PETSC_SUCCESS); 536 } 537 PetscCall(DMLabelReset(*label)); 538 PetscCall(PetscHMapIDestroy(&(*label)->hmap)); 539 PetscCall(PetscHeaderDestroy(label)); 540 PetscFunctionReturn(PETSC_SUCCESS); 541 } 542 543 static PetscErrorCode DMLabelDuplicate_Concrete(DMLabel label, DMLabel *labelnew) 544 { 545 PetscFunctionBegin; 546 for (PetscInt v = 0; v < label->numStrata; ++v) { 547 PetscCall(PetscHSetICreate(&(*labelnew)->ht[v])); 548 PetscCall(PetscObjectReference((PetscObject)label->points[v])); 549 (*labelnew)->points[v] = label->points[v]; 550 } 551 PetscCall(PetscHMapIDestroy(&(*labelnew)->hmap)); 552 PetscCall(PetscHMapIDuplicate(label->hmap, &(*labelnew)->hmap)); 553 PetscFunctionReturn(PETSC_SUCCESS); 554 } 555 556 /*@ 557 DMLabelDuplicate - Duplicates a `DMLabel` 558 559 Collective 560 561 Input Parameter: 562 . label - The `DMLabel` 563 564 Output Parameter: 565 . labelnew - new label 566 567 Level: intermediate 568 569 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelDestroy()` 570 @*/ 571 PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew) 572 { 573 const char *name; 574 575 PetscFunctionBegin; 576 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 577 PetscCall(DMLabelMakeAllValid_Private(label)); 578 PetscCall(PetscObjectGetName((PetscObject)label, &name)); 579 PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)label), name, labelnew)); 580 581 (*labelnew)->numStrata = label->numStrata; 582 (*labelnew)->defaultValue = label->defaultValue; 583 (*labelnew)->readonly = label->readonly; 584 PetscCall(PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues)); 585 PetscCall(PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes)); 586 PetscCall(PetscCalloc1(label->numStrata, &(*labelnew)->ht)); 587 PetscCall(PetscCalloc1(label->numStrata, &(*labelnew)->points)); 588 PetscCall(PetscMalloc1(label->numStrata, &(*labelnew)->validIS)); 589 for (PetscInt v = 0; v < label->numStrata; ++v) { 590 (*labelnew)->stratumValues[v] = label->stratumValues[v]; 591 (*labelnew)->stratumSizes[v] = label->stratumSizes[v]; 592 (*labelnew)->validIS[v] = PETSC_TRUE; 593 } 594 (*labelnew)->pStart = -1; 595 (*labelnew)->pEnd = -1; 596 (*labelnew)->bt = NULL; 597 PetscUseTypeMethod(label, duplicate, labelnew); 598 PetscFunctionReturn(PETSC_SUCCESS); 599 } 600 601 /*@C 602 DMLabelCompare - Compare two `DMLabel` objects 603 604 Collective; No Fortran Support 605 606 Input Parameters: 607 + comm - Comm over which to compare labels 608 . l0 - First `DMLabel` 609 - l1 - Second `DMLabel` 610 611 Output Parameters: 612 + equal - (Optional) Flag whether the two labels are equal 613 - message - (Optional) Message describing the difference 614 615 Level: intermediate 616 617 Notes: 618 The output flag equal is the same on all processes. 619 If it is passed as `NULL` and difference is found, an error is thrown on all processes. 620 Make sure to pass `NULL` on all processes. 621 622 The output message is set independently on each rank. 623 It is set to `NULL` if no difference was found on the current rank. It must be freed by user. 624 If message is passed as `NULL` and difference is found, the difference description is printed to stderr in synchronized manner. 625 Make sure to pass `NULL` on all processes. 626 627 For the comparison, we ignore the order of stratum values, and strata with no points. 628 629 The communicator needs to be specified because currently `DMLabel` can live on `PETSC_COMM_SELF` even if the underlying `DM` is parallel. 630 631 Developer Note: 632 Fortran stub cannot be generated automatically because `message` must be freed with `PetscFree()` 633 634 .seealso: `DMLabel`, `DM`, `DMCompareLabels()`, `DMLabelGetNumValues()`, `DMLabelGetDefaultValue()`, `DMLabelGetNonEmptyStratumValuesIS()`, `DMLabelGetStratumIS()` 635 @*/ 636 PetscErrorCode DMLabelCompare(MPI_Comm comm, DMLabel l0, DMLabel l1, PetscBool *equal, char **message) 637 { 638 const char *name0, *name1; 639 char msg[PETSC_MAX_PATH_LEN] = ""; 640 PetscBool eq; 641 PetscMPIInt rank; 642 643 PetscFunctionBegin; 644 PetscValidHeaderSpecific(l0, DMLABEL_CLASSID, 2); 645 PetscValidHeaderSpecific(l1, DMLABEL_CLASSID, 3); 646 if (equal) PetscAssertPointer(equal, 4); 647 if (message) PetscAssertPointer(message, 5); 648 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 649 PetscCall(PetscObjectGetName((PetscObject)l0, &name0)); 650 PetscCall(PetscObjectGetName((PetscObject)l1, &name1)); 651 { 652 PetscInt v0, v1; 653 654 PetscCall(DMLabelGetDefaultValue(l0, &v0)); 655 PetscCall(DMLabelGetDefaultValue(l1, &v1)); 656 eq = (PetscBool)(v0 == v1); 657 if (!eq) PetscCall(PetscSNPrintf(msg, sizeof(msg), "Default value of DMLabel l0 \"%s\" = %" PetscInt_FMT " != %" PetscInt_FMT " = Default value of DMLabel l1 \"%s\"", name0, v0, v1, name1)); 658 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPIU_BOOL, MPI_LAND, comm)); 659 if (!eq) goto finish; 660 } 661 { 662 IS is0, is1; 663 664 PetscCall(DMLabelGetNonEmptyStratumValuesIS(l0, &is0)); 665 PetscCall(DMLabelGetNonEmptyStratumValuesIS(l1, &is1)); 666 PetscCall(ISEqual(is0, is1, &eq)); 667 PetscCall(ISDestroy(&is0)); 668 PetscCall(ISDestroy(&is1)); 669 if (!eq) PetscCall(PetscSNPrintf(msg, sizeof(msg), "Stratum values in DMLabel l0 \"%s\" are different than in DMLabel l1 \"%s\"", name0, name1)); 670 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPIU_BOOL, MPI_LAND, comm)); 671 if (!eq) goto finish; 672 } 673 { 674 PetscInt i, nValues; 675 676 PetscCall(DMLabelGetNumValues(l0, &nValues)); 677 for (i = 0; i < nValues; i++) { 678 const PetscInt v = l0->stratumValues[i]; 679 PetscInt n; 680 IS is0, is1; 681 682 PetscCall(DMLabelGetStratumSize_Private(l0, i, &n)); 683 if (!n) continue; 684 PetscCall(DMLabelGetStratumIS(l0, v, &is0)); 685 PetscCall(DMLabelGetStratumIS(l1, v, &is1)); 686 PetscCall(ISEqualUnsorted(is0, is1, &eq)); 687 PetscCall(ISDestroy(&is0)); 688 PetscCall(ISDestroy(&is1)); 689 if (!eq) { 690 PetscCall(PetscSNPrintf(msg, sizeof(msg), "Stratum #%" PetscInt_FMT " with value %" PetscInt_FMT " contains different points in DMLabel l0 \"%s\" and DMLabel l1 \"%s\"", i, v, name0, name1)); 691 break; 692 } 693 } 694 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPIU_BOOL, MPI_LAND, comm)); 695 } 696 finish: 697 /* If message output arg not set, print to stderr */ 698 if (message) { 699 *message = NULL; 700 if (msg[0]) PetscCall(PetscStrallocpy(msg, message)); 701 } else { 702 if (msg[0]) PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] %s\n", rank, msg)); 703 PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR)); 704 } 705 /* If same output arg not ser and labels are not equal, throw error */ 706 if (equal) *equal = eq; 707 else PetscCheck(eq, comm, PETSC_ERR_ARG_INCOMP, "DMLabels l0 \"%s\" and l1 \"%s\" are not equal", name0, name1); 708 PetscFunctionReturn(PETSC_SUCCESS); 709 } 710 711 /*@ 712 DMLabelComputeIndex - Create an index structure for membership determination, automatically determining the bounds 713 714 Not Collective 715 716 Input Parameter: 717 . label - The `DMLabel` 718 719 Level: intermediate 720 721 .seealso: `DMLabel`, `DM`, `DMLabelHasPoint()`, `DMLabelCreateIndex()`, `DMLabelDestroyIndex()`, `DMLabelGetValue()`, `DMLabelSetValue()` 722 @*/ 723 PetscErrorCode DMLabelComputeIndex(DMLabel label) 724 { 725 PetscInt pStart = PETSC_MAX_INT, pEnd = -1, v; 726 727 PetscFunctionBegin; 728 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 729 PetscCall(DMLabelMakeAllValid_Private(label)); 730 for (v = 0; v < label->numStrata; ++v) { 731 const PetscInt *points; 732 PetscInt i; 733 734 PetscCall(ISGetIndices(label->points[v], &points)); 735 for (i = 0; i < label->stratumSizes[v]; ++i) { 736 const PetscInt point = points[i]; 737 738 pStart = PetscMin(point, pStart); 739 pEnd = PetscMax(point + 1, pEnd); 740 } 741 PetscCall(ISRestoreIndices(label->points[v], &points)); 742 } 743 label->pStart = pStart == PETSC_MAX_INT ? -1 : pStart; 744 label->pEnd = pEnd; 745 PetscCall(DMLabelCreateIndex(label, label->pStart, label->pEnd)); 746 PetscFunctionReturn(PETSC_SUCCESS); 747 } 748 749 /*@ 750 DMLabelCreateIndex - Create an index structure for membership determination 751 752 Not Collective 753 754 Input Parameters: 755 + label - The `DMLabel` 756 . pStart - The smallest point 757 - pEnd - The largest point + 1 758 759 Level: intermediate 760 761 .seealso: `DMLabel`, `DM`, `DMLabelHasPoint()`, `DMLabelComputeIndex()`, `DMLabelDestroyIndex()`, `DMLabelGetValue()`, `DMLabelSetValue()` 762 @*/ 763 PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd) 764 { 765 PetscInt v; 766 767 PetscFunctionBegin; 768 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 769 PetscCall(DMLabelDestroyIndex(label)); 770 PetscCall(DMLabelMakeAllValid_Private(label)); 771 label->pStart = pStart; 772 label->pEnd = pEnd; 773 /* This can be hooked into SetValue(), ClearValue(), etc. for updating */ 774 PetscCall(PetscBTCreate(pEnd - pStart, &label->bt)); 775 for (v = 0; v < label->numStrata; ++v) { 776 IS pointIS; 777 const PetscInt *points; 778 PetscInt i; 779 780 PetscUseTypeMethod(label, getstratumis, v, &pointIS); 781 PetscCall(ISGetIndices(pointIS, &points)); 782 for (i = 0; i < label->stratumSizes[v]; ++i) { 783 const PetscInt point = points[i]; 784 785 PetscCheck(!(point < pStart) && !(point >= pEnd), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %" PetscInt_FMT " in stratum %" PetscInt_FMT " is not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, label->stratumValues[v], pStart, pEnd); 786 PetscCall(PetscBTSet(label->bt, point - pStart)); 787 } 788 PetscCall(ISRestoreIndices(label->points[v], &points)); 789 PetscCall(ISDestroy(&pointIS)); 790 } 791 PetscFunctionReturn(PETSC_SUCCESS); 792 } 793 794 /*@ 795 DMLabelDestroyIndex - Destroy the index structure 796 797 Not Collective 798 799 Input Parameter: 800 . label - the `DMLabel` 801 802 Level: intermediate 803 804 .seealso: `DMLabel`, `DM`, `DMLabelHasPoint()`, `DMLabelCreateIndex()`, `DMLabelGetValue()`, `DMLabelSetValue()` 805 @*/ 806 PetscErrorCode DMLabelDestroyIndex(DMLabel label) 807 { 808 PetscFunctionBegin; 809 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 810 label->pStart = -1; 811 label->pEnd = -1; 812 PetscCall(PetscBTDestroy(&label->bt)); 813 PetscFunctionReturn(PETSC_SUCCESS); 814 } 815 816 /*@ 817 DMLabelGetBounds - Return the smallest and largest point in the label 818 819 Not Collective 820 821 Input Parameter: 822 . label - the `DMLabel` 823 824 Output Parameters: 825 + pStart - The smallest point 826 - pEnd - The largest point + 1 827 828 Level: intermediate 829 830 Note: 831 This will compute an index for the label if one does not exist. 832 833 .seealso: `DMLabel`, `DM`, `DMLabelHasPoint()`, `DMLabelCreateIndex()`, `DMLabelGetValue()`, `DMLabelSetValue()` 834 @*/ 835 PetscErrorCode DMLabelGetBounds(DMLabel label, PetscInt *pStart, PetscInt *pEnd) 836 { 837 PetscFunctionBegin; 838 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 839 if ((label->pStart == -1) && (label->pEnd == -1)) PetscCall(DMLabelComputeIndex(label)); 840 if (pStart) { 841 PetscAssertPointer(pStart, 2); 842 *pStart = label->pStart; 843 } 844 if (pEnd) { 845 PetscAssertPointer(pEnd, 3); 846 *pEnd = label->pEnd; 847 } 848 PetscFunctionReturn(PETSC_SUCCESS); 849 } 850 851 /*@ 852 DMLabelHasValue - Determine whether a label assigns the value to any point 853 854 Not Collective 855 856 Input Parameters: 857 + label - the `DMLabel` 858 - value - the value 859 860 Output Parameter: 861 . contains - Flag indicating whether the label maps this value to any point 862 863 Level: developer 864 865 .seealso: `DMLabel`, `DM`, `DMLabelHasPoint()`, `DMLabelGetValue()`, `DMLabelSetValue()` 866 @*/ 867 PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains) 868 { 869 PetscInt v; 870 871 PetscFunctionBegin; 872 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 873 PetscAssertPointer(contains, 3); 874 PetscCall(DMLabelLookupStratum(label, value, &v)); 875 *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE; 876 PetscFunctionReturn(PETSC_SUCCESS); 877 } 878 879 /*@ 880 DMLabelHasPoint - Determine whether a label assigns a value to a point 881 882 Not Collective 883 884 Input Parameters: 885 + label - the `DMLabel` 886 - point - the point 887 888 Output Parameter: 889 . contains - Flag indicating whether the label maps this point to a value 890 891 Level: developer 892 893 Note: 894 The user must call `DMLabelCreateIndex()` before this function. 895 896 .seealso: `DMLabel`, `DM`, `DMLabelCreateIndex()`, `DMLabelGetValue()`, `DMLabelSetValue()` 897 @*/ 898 PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains) 899 { 900 PetscFunctionBeginHot; 901 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 902 PetscAssertPointer(contains, 3); 903 PetscCall(DMLabelMakeAllValid_Private(label)); 904 if (PetscDefined(USE_DEBUG)) { 905 PetscCheck(label->bt, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()"); 906 PetscCheck(!(point < label->pStart) && !(point >= label->pEnd), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %" PetscInt_FMT " is not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, label->pStart, label->pEnd); 907 } 908 *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE; 909 PetscFunctionReturn(PETSC_SUCCESS); 910 } 911 912 /*@ 913 DMLabelStratumHasPoint - Return true if the stratum contains a point 914 915 Not Collective 916 917 Input Parameters: 918 + label - the `DMLabel` 919 . value - the stratum value 920 - point - the point 921 922 Output Parameter: 923 . contains - true if the stratum contains the point 924 925 Level: intermediate 926 927 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelSetValue()`, `DMLabelClearValue()` 928 @*/ 929 PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains) 930 { 931 PetscFunctionBeginHot; 932 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 933 PetscAssertPointer(contains, 4); 934 if (value == label->defaultValue) { 935 PetscInt pointVal; 936 937 PetscCall(DMLabelGetValue(label, point, &pointVal)); 938 *contains = (PetscBool)(pointVal == value); 939 } else { 940 PetscInt v; 941 942 PetscCall(DMLabelLookupStratum(label, value, &v)); 943 if (v >= 0) { 944 if (label->validIS[v] || label->readonly) { 945 IS is; 946 const PetscInt *indices; 947 PetscInt i; 948 949 PetscUseTypeMethod(label, getstratumis, v, &is); 950 // For some reason, ISLocate() is insanely slow here 951 PetscCall(ISGetIndices(is, &indices)); 952 PetscCall(PetscFindInt(point, label->stratumSizes[v], indices, &i)); 953 PetscCall(ISRestoreIndices(is, &indices)); 954 PetscCall(ISDestroy(&is)); 955 *contains = (PetscBool)(i >= 0); 956 } else { 957 PetscCall(PetscHSetIHas(label->ht[v], point, contains)); 958 } 959 } else { // value is not present 960 *contains = PETSC_FALSE; 961 } 962 } 963 PetscFunctionReturn(PETSC_SUCCESS); 964 } 965 966 /*@ 967 DMLabelGetDefaultValue - Get the default value returned by `DMLabelGetValue()` if a point has not been explicitly given a value. 968 When a label is created, it is initialized to -1. 969 970 Not Collective 971 972 Input Parameter: 973 . label - a `DMLabel` object 974 975 Output Parameter: 976 . defaultValue - the default value 977 978 Level: beginner 979 980 .seealso: `DMLabel`, `DM`, `DMLabelSetDefaultValue()`, `DMLabelGetValue()`, `DMLabelSetValue()` 981 @*/ 982 PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue) 983 { 984 PetscFunctionBegin; 985 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 986 *defaultValue = label->defaultValue; 987 PetscFunctionReturn(PETSC_SUCCESS); 988 } 989 990 /*@ 991 DMLabelSetDefaultValue - Set the default value returned by `DMLabelGetValue()` if a point has not been explicitly given a value. 992 When a label is created, it is initialized to -1. 993 994 Not Collective 995 996 Input Parameter: 997 . label - a `DMLabel` object 998 999 Output Parameter: 1000 . defaultValue - the default value 1001 1002 Level: beginner 1003 1004 .seealso: `DMLabel`, `DM`, `DMLabelGetDefaultValue()`, `DMLabelGetValue()`, `DMLabelSetValue()` 1005 @*/ 1006 PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue) 1007 { 1008 PetscFunctionBegin; 1009 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1010 label->defaultValue = defaultValue; 1011 PetscFunctionReturn(PETSC_SUCCESS); 1012 } 1013 1014 /*@ 1015 DMLabelGetValue - Return the value a label assigns to a point, or the label's default value (which is initially -1, and can be changed with 1016 `DMLabelSetDefaultValue()`) 1017 1018 Not Collective 1019 1020 Input Parameters: 1021 + label - the `DMLabel` 1022 - point - the point 1023 1024 Output Parameter: 1025 . value - The point value, or the default value (-1 by default) 1026 1027 Level: intermediate 1028 1029 Note: 1030 A label may assign multiple values to a point. No guarantees are made about which value is returned in that case. 1031 Use `DMLabelStratumHasPoint()` to check for inclusion in a specific value stratum. 1032 1033 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelSetValue()`, `DMLabelClearValue()`, `DMLabelGetDefaultValue()`, `DMLabelSetDefaultValue()` 1034 @*/ 1035 PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value) 1036 { 1037 PetscInt v; 1038 1039 PetscFunctionBeginHot; 1040 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1041 PetscAssertPointer(value, 3); 1042 *value = label->defaultValue; 1043 for (v = 0; v < label->numStrata; ++v) { 1044 if (label->validIS[v] || label->readonly) { 1045 IS is; 1046 const PetscInt *indices; 1047 PetscInt i; 1048 1049 PetscUseTypeMethod(label, getstratumis, v, &is); 1050 // For some reason, ISLocate() is insanely slow here 1051 PetscCall(ISGetIndices(is, &indices)); 1052 PetscCall(PetscFindInt(point, label->stratumSizes[v], indices, &i)); 1053 PetscCall(ISRestoreIndices(is, &indices)); 1054 PetscCall(ISDestroy(&is)); 1055 if (i >= 0) { 1056 *value = label->stratumValues[v]; 1057 break; 1058 } 1059 } else { 1060 PetscBool has; 1061 1062 PetscCall(PetscHSetIHas(label->ht[v], point, &has)); 1063 if (has) { 1064 *value = label->stratumValues[v]; 1065 break; 1066 } 1067 } 1068 } 1069 PetscFunctionReturn(PETSC_SUCCESS); 1070 } 1071 1072 /*@ 1073 DMLabelSetValue - Set the value a label assigns to a point. If the value is the same as the label's default value (which is initially -1, and can 1074 be changed with `DMLabelSetDefaultValue()` to something different), then this function will do nothing. 1075 1076 Not Collective 1077 1078 Input Parameters: 1079 + label - the `DMLabel` 1080 . point - the point 1081 - value - The point value 1082 1083 Level: intermediate 1084 1085 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelClearValue()`, `DMLabelGetDefaultValue()`, `DMLabelSetDefaultValue()` 1086 @*/ 1087 PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value) 1088 { 1089 PetscInt v; 1090 1091 PetscFunctionBegin; 1092 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1093 /* Find label value, add new entry if needed */ 1094 if (value == label->defaultValue) PetscFunctionReturn(PETSC_SUCCESS); 1095 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1096 PetscCall(DMLabelLookupAddStratum(label, value, &v)); 1097 /* Set key */ 1098 PetscCall(DMLabelMakeInvalid_Private(label, v)); 1099 PetscCall(PetscHSetIAdd(label->ht[v], point)); 1100 PetscFunctionReturn(PETSC_SUCCESS); 1101 } 1102 1103 /*@ 1104 DMLabelClearValue - Clear the value a label assigns to a point 1105 1106 Not Collective 1107 1108 Input Parameters: 1109 + label - the `DMLabel` 1110 . point - the point 1111 - value - The point value 1112 1113 Level: intermediate 1114 1115 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()` 1116 @*/ 1117 PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value) 1118 { 1119 PetscInt v; 1120 1121 PetscFunctionBegin; 1122 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1123 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1124 /* Find label value */ 1125 PetscCall(DMLabelLookupStratum(label, value, &v)); 1126 if (v < 0) PetscFunctionReturn(PETSC_SUCCESS); 1127 1128 if (label->bt) { 1129 PetscCheck(!(point < label->pStart) && !(point >= label->pEnd), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %" PetscInt_FMT " is not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, label->pStart, label->pEnd); 1130 PetscCall(PetscBTClear(label->bt, point - label->pStart)); 1131 } 1132 1133 /* Delete key */ 1134 PetscCall(DMLabelMakeInvalid_Private(label, v)); 1135 PetscCall(PetscHSetIDel(label->ht[v], point)); 1136 PetscFunctionReturn(PETSC_SUCCESS); 1137 } 1138 1139 /*@ 1140 DMLabelInsertIS - Set all points in the `IS` to a value 1141 1142 Not Collective 1143 1144 Input Parameters: 1145 + label - the `DMLabel` 1146 . is - the point `IS` 1147 - value - The point value 1148 1149 Level: intermediate 1150 1151 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1152 @*/ 1153 PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value) 1154 { 1155 PetscInt v, n, p; 1156 const PetscInt *points; 1157 1158 PetscFunctionBegin; 1159 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1160 PetscValidHeaderSpecific(is, IS_CLASSID, 2); 1161 /* Find label value, add new entry if needed */ 1162 if (value == label->defaultValue) PetscFunctionReturn(PETSC_SUCCESS); 1163 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1164 PetscCall(DMLabelLookupAddStratum(label, value, &v)); 1165 /* Set keys */ 1166 PetscCall(DMLabelMakeInvalid_Private(label, v)); 1167 PetscCall(ISGetLocalSize(is, &n)); 1168 PetscCall(ISGetIndices(is, &points)); 1169 for (p = 0; p < n; ++p) PetscCall(PetscHSetIAdd(label->ht[v], points[p])); 1170 PetscCall(ISRestoreIndices(is, &points)); 1171 PetscFunctionReturn(PETSC_SUCCESS); 1172 } 1173 1174 /*@ 1175 DMLabelGetNumValues - Get the number of values that the `DMLabel` takes 1176 1177 Not Collective 1178 1179 Input Parameter: 1180 . label - the `DMLabel` 1181 1182 Output Parameter: 1183 . numValues - the number of values 1184 1185 Level: intermediate 1186 1187 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1188 @*/ 1189 PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues) 1190 { 1191 PetscFunctionBegin; 1192 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1193 PetscAssertPointer(numValues, 2); 1194 *numValues = label->numStrata; 1195 PetscFunctionReturn(PETSC_SUCCESS); 1196 } 1197 1198 /*@ 1199 DMLabelGetValueIS - Get an `IS` of all values that the `DMlabel` takes 1200 1201 Not Collective 1202 1203 Input Parameter: 1204 . label - the `DMLabel` 1205 1206 Output Parameter: 1207 . values - the value `IS` 1208 1209 Level: intermediate 1210 1211 Notes: 1212 The output `IS` should be destroyed when no longer needed. 1213 Strata which are allocated but empty [`DMLabelGetStratumSize()` yields 0] are counted. 1214 If you need to count only nonempty strata, use `DMLabelGetNonEmptyStratumValuesIS()`. 1215 1216 .seealso: `DMLabel`, `DM`, `DMLabelGetNonEmptyStratumValuesIS()`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1217 @*/ 1218 PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values) 1219 { 1220 PetscFunctionBegin; 1221 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1222 PetscAssertPointer(values, 2); 1223 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values)); 1224 PetscFunctionReturn(PETSC_SUCCESS); 1225 } 1226 1227 /*@ 1228 DMLabelGetValueBounds - Return the smallest and largest value in the label 1229 1230 Not Collective 1231 1232 Input Parameter: 1233 . label - the `DMLabel` 1234 1235 Output Parameters: 1236 + minValue - The smallest value 1237 - maxValue - The largest value 1238 1239 Level: intermediate 1240 1241 .seealso: `DMLabel`, `DM`, `DMLabelGetBounds()`, `DMLabelGetValue()`, `DMLabelSetValue()` 1242 @*/ 1243 PetscErrorCode DMLabelGetValueBounds(DMLabel label, PetscInt *minValue, PetscInt *maxValue) 1244 { 1245 PetscInt min = PETSC_MAX_INT, max = PETSC_MIN_INT; 1246 1247 PetscFunctionBegin; 1248 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1249 for (PetscInt v = 0; v < label->numStrata; ++v) { 1250 min = PetscMin(min, label->stratumValues[v]); 1251 max = PetscMax(max, label->stratumValues[v]); 1252 } 1253 if (minValue) { 1254 PetscAssertPointer(minValue, 2); 1255 *minValue = min; 1256 } 1257 if (maxValue) { 1258 PetscAssertPointer(maxValue, 3); 1259 *maxValue = max; 1260 } 1261 PetscFunctionReturn(PETSC_SUCCESS); 1262 } 1263 1264 /*@ 1265 DMLabelGetNonEmptyStratumValuesIS - Get an `IS` of all values that the `DMlabel` takes 1266 1267 Not Collective 1268 1269 Input Parameter: 1270 . label - the `DMLabel` 1271 1272 Output Parameter: 1273 . values - the value `IS` 1274 1275 Level: intermediate 1276 1277 Notes: 1278 The output `IS` should be destroyed when no longer needed. 1279 This is similar to `DMLabelGetValueIS()` but counts only nonempty strata. 1280 1281 .seealso: `DMLabel`, `DM`, `DMLabelGetValueIS()`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1282 @*/ 1283 PetscErrorCode DMLabelGetNonEmptyStratumValuesIS(DMLabel label, IS *values) 1284 { 1285 PetscInt i, j; 1286 PetscInt *valuesArr; 1287 1288 PetscFunctionBegin; 1289 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1290 PetscAssertPointer(values, 2); 1291 PetscCall(PetscMalloc1(label->numStrata, &valuesArr)); 1292 for (i = 0, j = 0; i < label->numStrata; i++) { 1293 PetscInt n; 1294 1295 PetscCall(DMLabelGetStratumSize_Private(label, i, &n)); 1296 if (n) valuesArr[j++] = label->stratumValues[i]; 1297 } 1298 if (j == label->numStrata) { 1299 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values)); 1300 } else { 1301 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, j, valuesArr, PETSC_COPY_VALUES, values)); 1302 } 1303 PetscCall(PetscFree(valuesArr)); 1304 PetscFunctionReturn(PETSC_SUCCESS); 1305 } 1306 1307 /*@ 1308 DMLabelGetValueIndex - Get the index of a given value in the list of values for the `DMlabel`, or -1 if it is not present 1309 1310 Not Collective 1311 1312 Input Parameters: 1313 + label - the `DMLabel` 1314 - value - the value 1315 1316 Output Parameter: 1317 . index - the index of value in the list of values 1318 1319 Level: intermediate 1320 1321 .seealso: `DMLabel`, `DM`, `DMLabelGetValueIS()`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1322 @*/ 1323 PetscErrorCode DMLabelGetValueIndex(DMLabel label, PetscInt value, PetscInt *index) 1324 { 1325 PetscInt v; 1326 1327 PetscFunctionBegin; 1328 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1329 PetscAssertPointer(index, 3); 1330 /* Do not assume they are sorted */ 1331 for (v = 0; v < label->numStrata; ++v) 1332 if (label->stratumValues[v] == value) break; 1333 if (v >= label->numStrata) *index = -1; 1334 else *index = v; 1335 PetscFunctionReturn(PETSC_SUCCESS); 1336 } 1337 1338 /*@ 1339 DMLabelHasStratum - Determine whether points exist with the given value 1340 1341 Not Collective 1342 1343 Input Parameters: 1344 + label - the `DMLabel` 1345 - value - the stratum value 1346 1347 Output Parameter: 1348 . exists - Flag saying whether points exist 1349 1350 Level: intermediate 1351 1352 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1353 @*/ 1354 PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists) 1355 { 1356 PetscInt v; 1357 1358 PetscFunctionBegin; 1359 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1360 PetscAssertPointer(exists, 3); 1361 PetscCall(DMLabelLookupStratum(label, value, &v)); 1362 *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE; 1363 PetscFunctionReturn(PETSC_SUCCESS); 1364 } 1365 1366 /*@ 1367 DMLabelGetStratumSize - Get the size of a stratum 1368 1369 Not Collective 1370 1371 Input Parameters: 1372 + label - the `DMLabel` 1373 - value - the stratum value 1374 1375 Output Parameter: 1376 . size - The number of points in the stratum 1377 1378 Level: intermediate 1379 1380 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1381 @*/ 1382 PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size) 1383 { 1384 PetscInt v; 1385 1386 PetscFunctionBegin; 1387 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1388 PetscAssertPointer(size, 3); 1389 PetscCall(DMLabelLookupStratum(label, value, &v)); 1390 PetscCall(DMLabelGetStratumSize_Private(label, v, size)); 1391 PetscFunctionReturn(PETSC_SUCCESS); 1392 } 1393 1394 /*@ 1395 DMLabelGetStratumBounds - Get the largest and smallest point of a stratum 1396 1397 Not Collective 1398 1399 Input Parameters: 1400 + label - the `DMLabel` 1401 - value - the stratum value 1402 1403 Output Parameters: 1404 + start - the smallest point in the stratum 1405 - end - the largest point in the stratum 1406 1407 Level: intermediate 1408 1409 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1410 @*/ 1411 PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end) 1412 { 1413 IS is; 1414 PetscInt v, min, max; 1415 1416 PetscFunctionBegin; 1417 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1418 if (start) { 1419 PetscAssertPointer(start, 3); 1420 *start = -1; 1421 } 1422 if (end) { 1423 PetscAssertPointer(end, 4); 1424 *end = -1; 1425 } 1426 PetscCall(DMLabelLookupStratum(label, value, &v)); 1427 if (v < 0) PetscFunctionReturn(PETSC_SUCCESS); 1428 PetscCall(DMLabelMakeValid_Private(label, v)); 1429 if (label->stratumSizes[v] <= 0) PetscFunctionReturn(PETSC_SUCCESS); 1430 PetscUseTypeMethod(label, getstratumis, v, &is); 1431 PetscCall(ISGetMinMax(is, &min, &max)); 1432 PetscCall(ISDestroy(&is)); 1433 if (start) *start = min; 1434 if (end) *end = max + 1; 1435 PetscFunctionReturn(PETSC_SUCCESS); 1436 } 1437 1438 static PetscErrorCode DMLabelGetStratumIS_Concrete(DMLabel label, PetscInt v, IS *pointIS) 1439 { 1440 PetscFunctionBegin; 1441 PetscCall(PetscObjectReference((PetscObject)label->points[v])); 1442 *pointIS = label->points[v]; 1443 PetscFunctionReturn(PETSC_SUCCESS); 1444 } 1445 1446 /*@ 1447 DMLabelGetStratumIS - Get an `IS` with the stratum points 1448 1449 Not Collective 1450 1451 Input Parameters: 1452 + label - the `DMLabel` 1453 - value - the stratum value 1454 1455 Output Parameter: 1456 . points - The stratum points 1457 1458 Level: intermediate 1459 1460 Notes: 1461 The output `IS` should be destroyed when no longer needed. 1462 Returns `NULL` if the stratum is empty. 1463 1464 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1465 @*/ 1466 PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points) 1467 { 1468 PetscInt v; 1469 1470 PetscFunctionBegin; 1471 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1472 PetscAssertPointer(points, 3); 1473 *points = NULL; 1474 PetscCall(DMLabelLookupStratum(label, value, &v)); 1475 if (v < 0) PetscFunctionReturn(PETSC_SUCCESS); 1476 PetscCall(DMLabelMakeValid_Private(label, v)); 1477 PetscUseTypeMethod(label, getstratumis, v, points); 1478 PetscFunctionReturn(PETSC_SUCCESS); 1479 } 1480 1481 /*@ 1482 DMLabelSetStratumIS - Set the stratum points using an `IS` 1483 1484 Not Collective 1485 1486 Input Parameters: 1487 + label - the `DMLabel` 1488 . value - the stratum value 1489 - is - The stratum points 1490 1491 Level: intermediate 1492 1493 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1494 @*/ 1495 PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is) 1496 { 1497 PetscInt v; 1498 1499 PetscFunctionBegin; 1500 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1501 PetscValidHeaderSpecific(is, IS_CLASSID, 3); 1502 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1503 PetscCall(DMLabelLookupAddStratum(label, value, &v)); 1504 if (is == label->points[v]) PetscFunctionReturn(PETSC_SUCCESS); 1505 PetscCall(DMLabelClearStratum(label, value)); 1506 PetscCall(ISGetLocalSize(is, &label->stratumSizes[v])); 1507 PetscCall(PetscObjectReference((PetscObject)is)); 1508 PetscCall(ISDestroy(&label->points[v])); 1509 label->points[v] = is; 1510 label->validIS[v] = PETSC_TRUE; 1511 PetscCall(PetscObjectStateIncrease((PetscObject)label)); 1512 if (label->bt) { 1513 const PetscInt *points; 1514 PetscInt p; 1515 1516 PetscCall(ISGetIndices(is, &points)); 1517 for (p = 0; p < label->stratumSizes[v]; ++p) { 1518 const PetscInt point = points[p]; 1519 1520 PetscCheck(!(point < label->pStart) && !(point >= label->pEnd), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %" PetscInt_FMT " is not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, label->pStart, label->pEnd); 1521 PetscCall(PetscBTSet(label->bt, point - label->pStart)); 1522 } 1523 } 1524 PetscFunctionReturn(PETSC_SUCCESS); 1525 } 1526 1527 /*@ 1528 DMLabelClearStratum - Remove a stratum 1529 1530 Not Collective 1531 1532 Input Parameters: 1533 + label - the `DMLabel` 1534 - value - the stratum value 1535 1536 Level: intermediate 1537 1538 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1539 @*/ 1540 PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value) 1541 { 1542 PetscInt v; 1543 1544 PetscFunctionBegin; 1545 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1546 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1547 PetscCall(DMLabelLookupStratum(label, value, &v)); 1548 if (v < 0) PetscFunctionReturn(PETSC_SUCCESS); 1549 if (label->validIS[v]) { 1550 if (label->bt) { 1551 PetscInt i; 1552 const PetscInt *points; 1553 1554 PetscCall(ISGetIndices(label->points[v], &points)); 1555 for (i = 0; i < label->stratumSizes[v]; ++i) { 1556 const PetscInt point = points[i]; 1557 1558 PetscCheck(!(point < label->pStart) && !(point >= label->pEnd), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %" PetscInt_FMT " is not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, label->pStart, label->pEnd); 1559 PetscCall(PetscBTClear(label->bt, point - label->pStart)); 1560 } 1561 PetscCall(ISRestoreIndices(label->points[v], &points)); 1562 } 1563 label->stratumSizes[v] = 0; 1564 PetscCall(ISDestroy(&label->points[v])); 1565 PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &label->points[v])); 1566 PetscCall(PetscObjectSetName((PetscObject)label->points[v], "indices")); 1567 PetscCall(PetscObjectStateIncrease((PetscObject)label)); 1568 } else { 1569 PetscCall(PetscHSetIClear(label->ht[v])); 1570 } 1571 PetscFunctionReturn(PETSC_SUCCESS); 1572 } 1573 1574 /*@ 1575 DMLabelSetStratumBounds - Efficiently give a contiguous set of points a given label value 1576 1577 Not Collective 1578 1579 Input Parameters: 1580 + label - The `DMLabel` 1581 . value - The label value for all points 1582 . pStart - The first point 1583 - pEnd - A point beyond all marked points 1584 1585 Level: intermediate 1586 1587 Note: 1588 The marks points are [`pStart`, `pEnd`), and only the bounds are stored. 1589 1590 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelSetStratumIS()`, `DMLabelGetStratumIS()` 1591 @*/ 1592 PetscErrorCode DMLabelSetStratumBounds(DMLabel label, PetscInt value, PetscInt pStart, PetscInt pEnd) 1593 { 1594 IS pIS; 1595 1596 PetscFunctionBegin; 1597 PetscCall(ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &pIS)); 1598 PetscCall(DMLabelSetStratumIS(label, value, pIS)); 1599 PetscCall(ISDestroy(&pIS)); 1600 PetscFunctionReturn(PETSC_SUCCESS); 1601 } 1602 1603 /*@ 1604 DMLabelGetStratumPointIndex - Get the index of a point in a given stratum 1605 1606 Not Collective 1607 1608 Input Parameters: 1609 + label - The `DMLabel` 1610 . value - The label value 1611 - p - A point with this value 1612 1613 Output Parameter: 1614 . index - The index of this point in the stratum, or -1 if the point is not in the stratum or the stratum does not exist 1615 1616 Level: intermediate 1617 1618 .seealso: `DMLabel`, `DM`, `DMLabelGetValueIndex()`, `DMLabelGetStratumIS()`, `DMLabelCreate()` 1619 @*/ 1620 PetscErrorCode DMLabelGetStratumPointIndex(DMLabel label, PetscInt value, PetscInt p, PetscInt *index) 1621 { 1622 IS pointIS; 1623 const PetscInt *indices; 1624 PetscInt v; 1625 1626 PetscFunctionBegin; 1627 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1628 PetscAssertPointer(index, 4); 1629 *index = -1; 1630 PetscCall(DMLabelLookupStratum(label, value, &v)); 1631 if (v < 0) PetscFunctionReturn(PETSC_SUCCESS); 1632 PetscCall(DMLabelMakeValid_Private(label, v)); 1633 PetscUseTypeMethod(label, getstratumis, v, &pointIS); 1634 PetscCall(ISGetIndices(pointIS, &indices)); 1635 PetscCall(PetscFindInt(p, label->stratumSizes[v], indices, index)); 1636 PetscCall(ISRestoreIndices(pointIS, &indices)); 1637 PetscCall(ISDestroy(&pointIS)); 1638 PetscFunctionReturn(PETSC_SUCCESS); 1639 } 1640 1641 /*@ 1642 DMLabelFilter - Remove all points outside of [`start`, `end`) 1643 1644 Not Collective 1645 1646 Input Parameters: 1647 + label - the `DMLabel` 1648 . start - the first point kept 1649 - end - one more than the last point kept 1650 1651 Level: intermediate 1652 1653 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1654 @*/ 1655 PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end) 1656 { 1657 PetscInt v; 1658 1659 PetscFunctionBegin; 1660 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1661 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1662 PetscCall(DMLabelDestroyIndex(label)); 1663 PetscCall(DMLabelMakeAllValid_Private(label)); 1664 for (v = 0; v < label->numStrata; ++v) { 1665 PetscCall(ISGeneralFilter(label->points[v], start, end)); 1666 PetscCall(ISGetLocalSize(label->points[v], &label->stratumSizes[v])); 1667 } 1668 PetscCall(DMLabelCreateIndex(label, start, end)); 1669 PetscFunctionReturn(PETSC_SUCCESS); 1670 } 1671 1672 /*@ 1673 DMLabelPermute - Create a new label with permuted points 1674 1675 Not Collective 1676 1677 Input Parameters: 1678 + label - the `DMLabel` 1679 - permutation - the point permutation 1680 1681 Output Parameter: 1682 . labelNew - the new label containing the permuted points 1683 1684 Level: intermediate 1685 1686 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1687 @*/ 1688 PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew) 1689 { 1690 const PetscInt *perm; 1691 PetscInt numValues, numPoints, v, q; 1692 1693 PetscFunctionBegin; 1694 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1695 PetscValidHeaderSpecific(permutation, IS_CLASSID, 2); 1696 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1697 PetscCall(DMLabelMakeAllValid_Private(label)); 1698 PetscCall(DMLabelDuplicate(label, labelNew)); 1699 PetscCall(DMLabelGetNumValues(*labelNew, &numValues)); 1700 PetscCall(ISGetLocalSize(permutation, &numPoints)); 1701 PetscCall(ISGetIndices(permutation, &perm)); 1702 for (v = 0; v < numValues; ++v) { 1703 const PetscInt size = (*labelNew)->stratumSizes[v]; 1704 const PetscInt *points; 1705 PetscInt *pointsNew; 1706 1707 PetscCall(ISGetIndices((*labelNew)->points[v], &points)); 1708 PetscCall(PetscCalloc1(size, &pointsNew)); 1709 for (q = 0; q < size; ++q) { 1710 const PetscInt point = points[q]; 1711 1712 PetscCheck(!(point < 0) && !(point >= numPoints), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ") for the remapping", point, numPoints); 1713 pointsNew[q] = perm[point]; 1714 } 1715 PetscCall(ISRestoreIndices((*labelNew)->points[v], &points)); 1716 PetscCall(PetscSortInt(size, pointsNew)); 1717 PetscCall(ISDestroy(&((*labelNew)->points[v]))); 1718 if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) { 1719 PetscCall(ISCreateStride(PETSC_COMM_SELF, size, pointsNew[0], 1, &((*labelNew)->points[v]))); 1720 PetscCall(PetscFree(pointsNew)); 1721 } else { 1722 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size, pointsNew, PETSC_OWN_POINTER, &((*labelNew)->points[v]))); 1723 } 1724 PetscCall(PetscObjectSetName((PetscObject)((*labelNew)->points[v]), "indices")); 1725 } 1726 PetscCall(ISRestoreIndices(permutation, &perm)); 1727 if (label->bt) { 1728 PetscCall(PetscBTDestroy(&label->bt)); 1729 PetscCall(DMLabelCreateIndex(label, label->pStart, label->pEnd)); 1730 } 1731 PetscFunctionReturn(PETSC_SUCCESS); 1732 } 1733 1734 static PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata) 1735 { 1736 MPI_Comm comm; 1737 PetscInt s, l, nroots, nleaves, offset, size; 1738 PetscInt *remoteOffsets, *rootStrata, *rootIdx; 1739 PetscSection rootSection; 1740 PetscSF labelSF; 1741 1742 PetscFunctionBegin; 1743 if (label) PetscCall(DMLabelMakeAllValid_Private(label)); 1744 PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 1745 /* Build a section of stratum values per point, generate the according SF 1746 and distribute point-wise stratum values to leaves. */ 1747 PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL)); 1748 PetscCall(PetscSectionCreate(comm, &rootSection)); 1749 PetscCall(PetscSectionSetChart(rootSection, 0, nroots)); 1750 if (label) { 1751 for (s = 0; s < label->numStrata; ++s) { 1752 const PetscInt *points; 1753 1754 PetscCall(ISGetIndices(label->points[s], &points)); 1755 for (l = 0; l < label->stratumSizes[s]; l++) PetscCall(PetscSectionAddDof(rootSection, points[l], 1)); 1756 PetscCall(ISRestoreIndices(label->points[s], &points)); 1757 } 1758 } 1759 PetscCall(PetscSectionSetUp(rootSection)); 1760 /* Create a point-wise array of stratum values */ 1761 PetscCall(PetscSectionGetStorageSize(rootSection, &size)); 1762 PetscCall(PetscMalloc1(size, &rootStrata)); 1763 PetscCall(PetscCalloc1(nroots, &rootIdx)); 1764 if (label) { 1765 for (s = 0; s < label->numStrata; ++s) { 1766 const PetscInt *points; 1767 1768 PetscCall(ISGetIndices(label->points[s], &points)); 1769 for (l = 0; l < label->stratumSizes[s]; l++) { 1770 const PetscInt p = points[l]; 1771 PetscCall(PetscSectionGetOffset(rootSection, p, &offset)); 1772 rootStrata[offset + rootIdx[p]++] = label->stratumValues[s]; 1773 } 1774 PetscCall(ISRestoreIndices(label->points[s], &points)); 1775 } 1776 } 1777 /* Build SF that maps label points to remote processes */ 1778 PetscCall(PetscSectionCreate(comm, leafSection)); 1779 PetscCall(PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection)); 1780 PetscCall(PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF)); 1781 PetscCall(PetscFree(remoteOffsets)); 1782 /* Send the strata for each point over the derived SF */ 1783 PetscCall(PetscSectionGetStorageSize(*leafSection, &size)); 1784 PetscCall(PetscMalloc1(size, leafStrata)); 1785 PetscCall(PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata, MPI_REPLACE)); 1786 PetscCall(PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata, MPI_REPLACE)); 1787 /* Clean up */ 1788 PetscCall(PetscFree(rootStrata)); 1789 PetscCall(PetscFree(rootIdx)); 1790 PetscCall(PetscSectionDestroy(&rootSection)); 1791 PetscCall(PetscSFDestroy(&labelSF)); 1792 PetscFunctionReturn(PETSC_SUCCESS); 1793 } 1794 1795 /*@ 1796 DMLabelDistribute - Create a new label pushed forward over the `PetscSF` 1797 1798 Collective 1799 1800 Input Parameters: 1801 + label - the `DMLabel` 1802 - sf - the map from old to new distribution 1803 1804 Output Parameter: 1805 . labelNew - the new redistributed label 1806 1807 Level: intermediate 1808 1809 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1810 @*/ 1811 PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew) 1812 { 1813 MPI_Comm comm; 1814 PetscSection leafSection; 1815 PetscInt p, pStart, pEnd, s, size, dof, offset, stratum; 1816 PetscInt *leafStrata, *strataIdx; 1817 PetscInt **points; 1818 const char *lname = NULL; 1819 char *name; 1820 PetscInt nameSize; 1821 PetscHSetI stratumHash; 1822 size_t len = 0; 1823 PetscMPIInt rank; 1824 1825 PetscFunctionBegin; 1826 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 1827 if (label) { 1828 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1829 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1830 PetscCall(DMLabelMakeAllValid_Private(label)); 1831 } 1832 PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 1833 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1834 /* Bcast name */ 1835 if (rank == 0) { 1836 PetscCall(PetscObjectGetName((PetscObject)label, &lname)); 1837 PetscCall(PetscStrlen(lname, &len)); 1838 } 1839 nameSize = len; 1840 PetscCallMPI(MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm)); 1841 PetscCall(PetscMalloc1(nameSize + 1, &name)); 1842 if (rank == 0) PetscCall(PetscArraycpy(name, lname, nameSize + 1)); 1843 PetscCallMPI(MPI_Bcast(name, nameSize + 1, MPI_CHAR, 0, comm)); 1844 PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, labelNew)); 1845 PetscCall(PetscFree(name)); 1846 /* Bcast defaultValue */ 1847 if (rank == 0) (*labelNew)->defaultValue = label->defaultValue; 1848 PetscCallMPI(MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm)); 1849 /* Distribute stratum values over the SF and get the point mapping on the receiver */ 1850 PetscCall(DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata)); 1851 /* Determine received stratum values and initialise new label*/ 1852 PetscCall(PetscHSetICreate(&stratumHash)); 1853 PetscCall(PetscSectionGetStorageSize(leafSection, &size)); 1854 for (p = 0; p < size; ++p) PetscCall(PetscHSetIAdd(stratumHash, leafStrata[p])); 1855 PetscCall(PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata)); 1856 PetscCall(PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS)); 1857 for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE; 1858 PetscCall(PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues)); 1859 /* Turn leafStrata into indices rather than stratum values */ 1860 offset = 0; 1861 PetscCall(PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues)); 1862 PetscCall(PetscSortInt((*labelNew)->numStrata, (*labelNew)->stratumValues)); 1863 for (s = 0; s < (*labelNew)->numStrata; ++s) PetscCall(PetscHMapISet((*labelNew)->hmap, (*labelNew)->stratumValues[s], s)); 1864 for (p = 0; p < size; ++p) { 1865 for (s = 0; s < (*labelNew)->numStrata; ++s) { 1866 if (leafStrata[p] == (*labelNew)->stratumValues[s]) { 1867 leafStrata[p] = s; 1868 break; 1869 } 1870 } 1871 } 1872 /* Rebuild the point strata on the receiver */ 1873 PetscCall(PetscCalloc1((*labelNew)->numStrata, &(*labelNew)->stratumSizes)); 1874 PetscCall(PetscSectionGetChart(leafSection, &pStart, &pEnd)); 1875 for (p = pStart; p < pEnd; p++) { 1876 PetscCall(PetscSectionGetDof(leafSection, p, &dof)); 1877 PetscCall(PetscSectionGetOffset(leafSection, p, &offset)); 1878 for (s = 0; s < dof; s++) (*labelNew)->stratumSizes[leafStrata[offset + s]]++; 1879 } 1880 PetscCall(PetscCalloc1((*labelNew)->numStrata, &(*labelNew)->ht)); 1881 PetscCall(PetscCalloc1((*labelNew)->numStrata, &(*labelNew)->points)); 1882 PetscCall(PetscCalloc1((*labelNew)->numStrata, &points)); 1883 for (s = 0; s < (*labelNew)->numStrata; ++s) { 1884 PetscCall(PetscHSetICreate(&(*labelNew)->ht[s])); 1885 PetscCall(PetscMalloc1((*labelNew)->stratumSizes[s], &points[s])); 1886 } 1887 /* Insert points into new strata */ 1888 PetscCall(PetscCalloc1((*labelNew)->numStrata, &strataIdx)); 1889 PetscCall(PetscSectionGetChart(leafSection, &pStart, &pEnd)); 1890 for (p = pStart; p < pEnd; p++) { 1891 PetscCall(PetscSectionGetDof(leafSection, p, &dof)); 1892 PetscCall(PetscSectionGetOffset(leafSection, p, &offset)); 1893 for (s = 0; s < dof; s++) { 1894 stratum = leafStrata[offset + s]; 1895 points[stratum][strataIdx[stratum]++] = p; 1896 } 1897 } 1898 for (s = 0; s < (*labelNew)->numStrata; s++) { 1899 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, (*labelNew)->stratumSizes[s], &points[s][0], PETSC_OWN_POINTER, &((*labelNew)->points[s]))); 1900 PetscCall(PetscObjectSetName((PetscObject)((*labelNew)->points[s]), "indices")); 1901 } 1902 PetscCall(PetscFree(points)); 1903 PetscCall(PetscHSetIDestroy(&stratumHash)); 1904 PetscCall(PetscFree(leafStrata)); 1905 PetscCall(PetscFree(strataIdx)); 1906 PetscCall(PetscSectionDestroy(&leafSection)); 1907 PetscFunctionReturn(PETSC_SUCCESS); 1908 } 1909 1910 /*@ 1911 DMLabelGather - Gather all label values from leafs into roots 1912 1913 Collective 1914 1915 Input Parameters: 1916 + label - the `DMLabel` 1917 - sf - the `PetscSF` communication map 1918 1919 Output Parameter: 1920 . labelNew - the new `DMLabel` with localised leaf values 1921 1922 Level: developer 1923 1924 Note: 1925 This is the inverse operation to `DMLabelDistribute()`. 1926 1927 .seealso: `DMLabel`, `DM`, `DMLabelDistribute()` 1928 @*/ 1929 PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew) 1930 { 1931 MPI_Comm comm; 1932 PetscSection rootSection; 1933 PetscSF sfLabel; 1934 PetscSFNode *rootPoints, *leafPoints; 1935 PetscInt p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset; 1936 const PetscInt *rootDegree, *ilocal; 1937 PetscInt *rootStrata; 1938 const char *lname; 1939 char *name; 1940 PetscInt nameSize; 1941 size_t len = 0; 1942 PetscMPIInt rank, size; 1943 1944 PetscFunctionBegin; 1945 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1946 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 1947 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1948 PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 1949 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1950 PetscCallMPI(MPI_Comm_size(comm, &size)); 1951 /* Bcast name */ 1952 if (rank == 0) { 1953 PetscCall(PetscObjectGetName((PetscObject)label, &lname)); 1954 PetscCall(PetscStrlen(lname, &len)); 1955 } 1956 nameSize = len; 1957 PetscCallMPI(MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm)); 1958 PetscCall(PetscMalloc1(nameSize + 1, &name)); 1959 if (rank == 0) PetscCall(PetscArraycpy(name, lname, nameSize + 1)); 1960 PetscCallMPI(MPI_Bcast(name, nameSize + 1, MPI_CHAR, 0, comm)); 1961 PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, labelNew)); 1962 PetscCall(PetscFree(name)); 1963 /* Gather rank/index pairs of leaves into local roots to build 1964 an inverse, multi-rooted SF. Note that this ignores local leaf 1965 indexing due to the use of the multiSF in PetscSFGather. */ 1966 PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL)); 1967 PetscCall(PetscMalloc1(nroots, &leafPoints)); 1968 for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1; 1969 for (p = 0; p < nleaves; p++) { 1970 PetscInt ilp = ilocal ? ilocal[p] : p; 1971 1972 leafPoints[ilp].index = ilp; 1973 leafPoints[ilp].rank = rank; 1974 } 1975 PetscCall(PetscSFComputeDegreeBegin(sf, &rootDegree)); 1976 PetscCall(PetscSFComputeDegreeEnd(sf, &rootDegree)); 1977 for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p]; 1978 PetscCall(PetscMalloc1(nmultiroots, &rootPoints)); 1979 PetscCall(PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints)); 1980 PetscCall(PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints)); 1981 PetscCall(PetscSFCreate(comm, &sfLabel)); 1982 PetscCall(PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER)); 1983 /* Migrate label over inverted SF to pull stratum values at leaves into roots. */ 1984 PetscCall(DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata)); 1985 /* Rebuild the point strata on the receiver */ 1986 for (p = 0, idx = 0; p < nroots; p++) { 1987 for (d = 0; d < rootDegree[p]; d++) { 1988 PetscCall(PetscSectionGetDof(rootSection, idx + d, &dof)); 1989 PetscCall(PetscSectionGetOffset(rootSection, idx + d, &offset)); 1990 for (s = 0; s < dof; s++) PetscCall(DMLabelSetValue(*labelNew, p, rootStrata[offset + s])); 1991 } 1992 idx += rootDegree[p]; 1993 } 1994 PetscCall(PetscFree(leafPoints)); 1995 PetscCall(PetscFree(rootStrata)); 1996 PetscCall(PetscSectionDestroy(&rootSection)); 1997 PetscCall(PetscSFDestroy(&sfLabel)); 1998 PetscFunctionReturn(PETSC_SUCCESS); 1999 } 2000 2001 static PetscErrorCode DMLabelPropagateInit_Internal(DMLabel label, PetscSF pointSF, PetscInt valArray[]) 2002 { 2003 const PetscInt *degree; 2004 const PetscInt *points; 2005 PetscInt Nr, r, Nl, l, val, defVal; 2006 2007 PetscFunctionBegin; 2008 PetscCall(DMLabelGetDefaultValue(label, &defVal)); 2009 /* Add in leaves */ 2010 PetscCall(PetscSFGetGraph(pointSF, &Nr, &Nl, &points, NULL)); 2011 for (l = 0; l < Nl; ++l) { 2012 PetscCall(DMLabelGetValue(label, points[l], &val)); 2013 if (val != defVal) valArray[points[l]] = val; 2014 } 2015 /* Add in shared roots */ 2016 PetscCall(PetscSFComputeDegreeBegin(pointSF, °ree)); 2017 PetscCall(PetscSFComputeDegreeEnd(pointSF, °ree)); 2018 for (r = 0; r < Nr; ++r) { 2019 if (degree[r]) { 2020 PetscCall(DMLabelGetValue(label, r, &val)); 2021 if (val != defVal) valArray[r] = val; 2022 } 2023 } 2024 PetscFunctionReturn(PETSC_SUCCESS); 2025 } 2026 2027 static PetscErrorCode DMLabelPropagateFini_Internal(DMLabel label, PetscSF pointSF, PetscInt valArray[], PetscErrorCode (*markPoint)(DMLabel, PetscInt, PetscInt, void *), void *ctx) 2028 { 2029 const PetscInt *degree; 2030 const PetscInt *points; 2031 PetscInt Nr, r, Nl, l, val, defVal; 2032 2033 PetscFunctionBegin; 2034 PetscCall(DMLabelGetDefaultValue(label, &defVal)); 2035 /* Read out leaves */ 2036 PetscCall(PetscSFGetGraph(pointSF, &Nr, &Nl, &points, NULL)); 2037 for (l = 0; l < Nl; ++l) { 2038 const PetscInt p = points[l]; 2039 const PetscInt cval = valArray[p]; 2040 2041 if (cval != defVal) { 2042 PetscCall(DMLabelGetValue(label, p, &val)); 2043 if (val == defVal) { 2044 PetscCall(DMLabelSetValue(label, p, cval)); 2045 if (markPoint) PetscCall((*markPoint)(label, p, cval, ctx)); 2046 } 2047 } 2048 } 2049 /* Read out shared roots */ 2050 PetscCall(PetscSFComputeDegreeBegin(pointSF, °ree)); 2051 PetscCall(PetscSFComputeDegreeEnd(pointSF, °ree)); 2052 for (r = 0; r < Nr; ++r) { 2053 if (degree[r]) { 2054 const PetscInt cval = valArray[r]; 2055 2056 if (cval != defVal) { 2057 PetscCall(DMLabelGetValue(label, r, &val)); 2058 if (val == defVal) { 2059 PetscCall(DMLabelSetValue(label, r, cval)); 2060 if (markPoint) PetscCall((*markPoint)(label, r, cval, ctx)); 2061 } 2062 } 2063 } 2064 } 2065 PetscFunctionReturn(PETSC_SUCCESS); 2066 } 2067 2068 /*@ 2069 DMLabelPropagateBegin - Setup a cycle of label propagation 2070 2071 Collective 2072 2073 Input Parameters: 2074 + label - The `DMLabel` to propagate across processes 2075 - sf - The `PetscSF` describing parallel layout of the label points 2076 2077 Level: intermediate 2078 2079 .seealso: `DMLabel`, `DM`, `DMLabelPropagateEnd()`, `DMLabelPropagatePush()` 2080 @*/ 2081 PetscErrorCode DMLabelPropagateBegin(DMLabel label, PetscSF sf) 2082 { 2083 PetscInt Nr, r, defVal; 2084 PetscMPIInt size; 2085 2086 PetscFunctionBegin; 2087 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 2088 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sf), &size)); 2089 if (size > 1) { 2090 PetscCall(DMLabelGetDefaultValue(label, &defVal)); 2091 PetscCall(PetscSFGetGraph(sf, &Nr, NULL, NULL, NULL)); 2092 if (Nr >= 0) PetscCall(PetscMalloc1(Nr, &label->propArray)); 2093 for (r = 0; r < Nr; ++r) label->propArray[r] = defVal; 2094 } 2095 PetscFunctionReturn(PETSC_SUCCESS); 2096 } 2097 2098 /*@ 2099 DMLabelPropagateEnd - Tear down a cycle of label propagation 2100 2101 Collective 2102 2103 Input Parameters: 2104 + label - The `DMLabel` to propagate across processes 2105 - pointSF - The `PetscSF` describing parallel layout of the label points 2106 2107 Level: intermediate 2108 2109 .seealso: `DMLabel`, `DM`, `DMLabelPropagateBegin()`, `DMLabelPropagatePush()` 2110 @*/ 2111 PetscErrorCode DMLabelPropagateEnd(DMLabel label, PetscSF pointSF) 2112 { 2113 PetscFunctionBegin; 2114 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 2115 PetscCall(PetscFree(label->propArray)); 2116 label->propArray = NULL; 2117 PetscFunctionReturn(PETSC_SUCCESS); 2118 } 2119 2120 /*@C 2121 DMLabelPropagatePush - Tear down a cycle of label propagation 2122 2123 Collective 2124 2125 Input Parameters: 2126 + label - The `DMLabel` to propagate across processes 2127 . pointSF - The `PetscSF` describing parallel layout of the label points 2128 . markPoint - An optional callback that is called when a point is marked, or `NULL` 2129 - ctx - An optional user context for the callback, or `NULL` 2130 2131 Calling sequence of `markPoint`: 2132 + label - The `DMLabel` 2133 . p - The point being marked 2134 . val - The label value for `p` 2135 - ctx - An optional user context 2136 2137 Level: intermediate 2138 2139 .seealso: `DMLabel`, `DM`, `DMLabelPropagateBegin()`, `DMLabelPropagateEnd()` 2140 @*/ 2141 PetscErrorCode DMLabelPropagatePush(DMLabel label, PetscSF pointSF, PetscErrorCode (*markPoint)(DMLabel label, PetscInt p, PetscInt val, void *ctx), void *ctx) 2142 { 2143 PetscInt *valArray = label->propArray, Nr; 2144 PetscMPIInt size; 2145 2146 PetscFunctionBegin; 2147 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 2148 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pointSF), &size)); 2149 PetscCall(PetscSFGetGraph(pointSF, &Nr, NULL, NULL, NULL)); 2150 if (size > 1 && Nr >= 0) { 2151 /* Communicate marked edges 2152 The current implementation allocates an array the size of the number of root. We put the label values into the 2153 array, and then call PetscSFReduce()+PetscSFBcast() to make the marks consistent. 2154 2155 TODO: We could use in-place communication with a different SF 2156 We use MPI_SUM for the Reduce, and check the result against the rootdegree. If sum >= rootdegree+1, then the edge has 2157 already been marked. If not, it might have been handled on the process in this round, but we add it anyway. 2158 2159 In order to update the queue with the new edges from the label communication, we use BcastAnOp(MPI_SUM), so that new 2160 values will have 1+0=1 and old values will have 1+1=2. Loop over these, resetting the values to 1, and adding any new 2161 edge to the queue. 2162 */ 2163 PetscCall(DMLabelPropagateInit_Internal(label, pointSF, valArray)); 2164 PetscCall(PetscSFReduceBegin(pointSF, MPIU_INT, valArray, valArray, MPI_MAX)); 2165 PetscCall(PetscSFReduceEnd(pointSF, MPIU_INT, valArray, valArray, MPI_MAX)); 2166 PetscCall(PetscSFBcastBegin(pointSF, MPIU_INT, valArray, valArray, MPI_REPLACE)); 2167 PetscCall(PetscSFBcastEnd(pointSF, MPIU_INT, valArray, valArray, MPI_REPLACE)); 2168 PetscCall(DMLabelPropagateFini_Internal(label, pointSF, valArray, markPoint, ctx)); 2169 } 2170 PetscFunctionReturn(PETSC_SUCCESS); 2171 } 2172 2173 /*@ 2174 DMLabelConvertToSection - Make a `PetscSection`/`IS` pair that encodes the label 2175 2176 Not Collective 2177 2178 Input Parameter: 2179 . label - the `DMLabel` 2180 2181 Output Parameters: 2182 + section - the section giving offsets for each stratum 2183 - is - An `IS` containing all the label points 2184 2185 Level: developer 2186 2187 .seealso: `DMLabel`, `DM`, `DMLabelDistribute()` 2188 @*/ 2189 PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is) 2190 { 2191 IS vIS; 2192 const PetscInt *values; 2193 PetscInt *points; 2194 PetscInt nV, vS = 0, vE = 0, v, N; 2195 2196 PetscFunctionBegin; 2197 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 2198 PetscCall(DMLabelGetNumValues(label, &nV)); 2199 PetscCall(DMLabelGetValueIS(label, &vIS)); 2200 PetscCall(ISGetIndices(vIS, &values)); 2201 if (nV) { 2202 vS = values[0]; 2203 vE = values[0] + 1; 2204 } 2205 for (v = 1; v < nV; ++v) { 2206 vS = PetscMin(vS, values[v]); 2207 vE = PetscMax(vE, values[v] + 1); 2208 } 2209 PetscCall(PetscSectionCreate(PETSC_COMM_SELF, section)); 2210 PetscCall(PetscSectionSetChart(*section, vS, vE)); 2211 for (v = 0; v < nV; ++v) { 2212 PetscInt n; 2213 2214 PetscCall(DMLabelGetStratumSize(label, values[v], &n)); 2215 PetscCall(PetscSectionSetDof(*section, values[v], n)); 2216 } 2217 PetscCall(PetscSectionSetUp(*section)); 2218 PetscCall(PetscSectionGetStorageSize(*section, &N)); 2219 PetscCall(PetscMalloc1(N, &points)); 2220 for (v = 0; v < nV; ++v) { 2221 IS is; 2222 const PetscInt *spoints; 2223 PetscInt dof, off, p; 2224 2225 PetscCall(PetscSectionGetDof(*section, values[v], &dof)); 2226 PetscCall(PetscSectionGetOffset(*section, values[v], &off)); 2227 PetscCall(DMLabelGetStratumIS(label, values[v], &is)); 2228 PetscCall(ISGetIndices(is, &spoints)); 2229 for (p = 0; p < dof; ++p) points[off + p] = spoints[p]; 2230 PetscCall(ISRestoreIndices(is, &spoints)); 2231 PetscCall(ISDestroy(&is)); 2232 } 2233 PetscCall(ISRestoreIndices(vIS, &values)); 2234 PetscCall(ISDestroy(&vIS)); 2235 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is)); 2236 PetscFunctionReturn(PETSC_SUCCESS); 2237 } 2238 2239 /*@C 2240 DMLabelRegister - Adds a new label component implementation 2241 2242 Not Collective 2243 2244 Input Parameters: 2245 + name - The name of a new user-defined creation routine 2246 - create_func - The creation routine itself 2247 2248 Notes: 2249 `DMLabelRegister()` may be called multiple times to add several user-defined labels 2250 2251 Example Usage: 2252 .vb 2253 DMLabelRegister("my_label", MyLabelCreate); 2254 .ve 2255 2256 Then, your label type can be chosen with the procedural interface via 2257 .vb 2258 DMLabelCreate(MPI_Comm, DMLabel *); 2259 DMLabelSetType(DMLabel, "my_label"); 2260 .ve 2261 or at runtime via the option 2262 .vb 2263 -dm_label_type my_label 2264 .ve 2265 2266 Level: advanced 2267 2268 .seealso: `DMLabel`, `DM`, `DMLabelType`, `DMLabelRegisterAll()`, `DMLabelRegisterDestroy()` 2269 @*/ 2270 PetscErrorCode DMLabelRegister(const char name[], PetscErrorCode (*create_func)(DMLabel)) 2271 { 2272 PetscFunctionBegin; 2273 PetscCall(DMInitializePackage()); 2274 PetscCall(PetscFunctionListAdd(&DMLabelList, name, create_func)); 2275 PetscFunctionReturn(PETSC_SUCCESS); 2276 } 2277 2278 PETSC_EXTERN PetscErrorCode DMLabelCreate_Concrete(DMLabel); 2279 PETSC_EXTERN PetscErrorCode DMLabelCreate_Ephemeral(DMLabel); 2280 2281 /*@C 2282 DMLabelRegisterAll - Registers all of the `DMLabel` implementations in the `DM` package. 2283 2284 Not Collective 2285 2286 Level: advanced 2287 2288 .seealso: `DMLabel`, `DM`, `DMRegisterAll()`, `DMLabelRegisterDestroy()` 2289 @*/ 2290 PetscErrorCode DMLabelRegisterAll(void) 2291 { 2292 PetscFunctionBegin; 2293 if (DMLabelRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 2294 DMLabelRegisterAllCalled = PETSC_TRUE; 2295 2296 PetscCall(DMLabelRegister(DMLABELCONCRETE, DMLabelCreate_Concrete)); 2297 PetscCall(DMLabelRegister(DMLABELEPHEMERAL, DMLabelCreate_Ephemeral)); 2298 PetscFunctionReturn(PETSC_SUCCESS); 2299 } 2300 2301 /*@C 2302 DMLabelRegisterDestroy - This function destroys the `DMLabel` registry. It is called from `PetscFinalize()`. 2303 2304 Level: developer 2305 2306 .seealso: `DMLabel`, `DM`, `PetscInitialize()` 2307 @*/ 2308 PetscErrorCode DMLabelRegisterDestroy(void) 2309 { 2310 PetscFunctionBegin; 2311 PetscCall(PetscFunctionListDestroy(&DMLabelList)); 2312 DMLabelRegisterAllCalled = PETSC_FALSE; 2313 PetscFunctionReturn(PETSC_SUCCESS); 2314 } 2315 2316 /*@ 2317 DMLabelSetType - Sets the particular implementation for a label. 2318 2319 Collective 2320 2321 Input Parameters: 2322 + label - The label 2323 - method - The name of the label type 2324 2325 Options Database Key: 2326 . -dm_label_type <type> - Sets the label type; use -help for a list of available types or see `DMLabelType` 2327 2328 Level: intermediate 2329 2330 .seealso: `DMLabel`, `DM`, `DMLabelGetType()`, `DMLabelCreate()` 2331 @*/ 2332 PetscErrorCode DMLabelSetType(DMLabel label, DMLabelType method) 2333 { 2334 PetscErrorCode (*r)(DMLabel); 2335 PetscBool match; 2336 2337 PetscFunctionBegin; 2338 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 2339 PetscCall(PetscObjectTypeCompare((PetscObject)label, method, &match)); 2340 if (match) PetscFunctionReturn(PETSC_SUCCESS); 2341 2342 PetscCall(DMLabelRegisterAll()); 2343 PetscCall(PetscFunctionListFind(DMLabelList, method, &r)); 2344 PetscCheck(r, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DMLabel type: %s", method); 2345 2346 PetscTryTypeMethod(label, destroy); 2347 PetscCall(PetscMemzero(label->ops, sizeof(*label->ops))); 2348 PetscCall(PetscObjectChangeTypeName((PetscObject)label, method)); 2349 PetscCall((*r)(label)); 2350 PetscFunctionReturn(PETSC_SUCCESS); 2351 } 2352 2353 /*@ 2354 DMLabelGetType - Gets the type name (as a string) from the label. 2355 2356 Not Collective 2357 2358 Input Parameter: 2359 . label - The `DMLabel` 2360 2361 Output Parameter: 2362 . type - The `DMLabel` type name 2363 2364 Level: intermediate 2365 2366 .seealso: `DMLabel`, `DM`, `DMLabelSetType()`, `DMLabelCreate()` 2367 @*/ 2368 PetscErrorCode DMLabelGetType(DMLabel label, DMLabelType *type) 2369 { 2370 PetscFunctionBegin; 2371 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 2372 PetscAssertPointer(type, 2); 2373 PetscCall(DMLabelRegisterAll()); 2374 *type = ((PetscObject)label)->type_name; 2375 PetscFunctionReturn(PETSC_SUCCESS); 2376 } 2377 2378 static PetscErrorCode DMLabelInitialize_Concrete(DMLabel label) 2379 { 2380 PetscFunctionBegin; 2381 label->ops->view = DMLabelView_Concrete; 2382 label->ops->setup = NULL; 2383 label->ops->duplicate = DMLabelDuplicate_Concrete; 2384 label->ops->getstratumis = DMLabelGetStratumIS_Concrete; 2385 PetscFunctionReturn(PETSC_SUCCESS); 2386 } 2387 2388 PETSC_EXTERN PetscErrorCode DMLabelCreate_Concrete(DMLabel label) 2389 { 2390 PetscFunctionBegin; 2391 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 2392 PetscCall(DMLabelInitialize_Concrete(label)); 2393 PetscFunctionReturn(PETSC_SUCCESS); 2394 } 2395 2396 /*@ 2397 PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using 2398 the local section and an `PetscSF` describing the section point overlap. 2399 2400 Collective 2401 2402 Input Parameters: 2403 + s - The `PetscSection` for the local field layout 2404 . sf - The `PetscSF` describing parallel layout of the section points 2405 . includeConstraints - By default this is `PETSC_FALSE`, meaning that the global field vector will not possess constrained dofs 2406 . label - The label specifying the points 2407 - labelValue - The label stratum specifying the points 2408 2409 Output Parameter: 2410 . gsection - The `PetscSection` for the global field layout 2411 2412 Level: developer 2413 2414 Note: 2415 This gives negative sizes and offsets to points not owned by this process 2416 2417 .seealso: `DMLabel`, `DM`, `PetscSectionCreate()` 2418 @*/ 2419 PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection) 2420 { 2421 PetscInt *neg = NULL, *tmpOff = NULL; 2422 PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots; 2423 2424 PetscFunctionBegin; 2425 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2426 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 2427 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4); 2428 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)s), gsection)); 2429 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 2430 PetscCall(PetscSectionSetChart(*gsection, pStart, pEnd)); 2431 PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 2432 if (nroots >= 0) { 2433 PetscCheck(nroots >= pEnd - pStart, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %" PetscInt_FMT " < %" PetscInt_FMT " section size", nroots, pEnd - pStart); 2434 PetscCall(PetscCalloc1(nroots, &neg)); 2435 if (nroots > pEnd - pStart) { 2436 PetscCall(PetscCalloc1(nroots, &tmpOff)); 2437 } else { 2438 tmpOff = &(*gsection)->atlasDof[-pStart]; 2439 } 2440 } 2441 /* Mark ghost points with negative dof */ 2442 for (p = pStart; p < pEnd; ++p) { 2443 PetscInt value; 2444 2445 PetscCall(DMLabelGetValue(label, p, &value)); 2446 if (value != labelValue) continue; 2447 PetscCall(PetscSectionGetDof(s, p, &dof)); 2448 PetscCall(PetscSectionSetDof(*gsection, p, dof)); 2449 PetscCall(PetscSectionGetConstraintDof(s, p, &cdof)); 2450 if (!includeConstraints && cdof > 0) PetscCall(PetscSectionSetConstraintDof(*gsection, p, cdof)); 2451 if (neg) neg[p] = -(dof + 1); 2452 } 2453 PetscCall(PetscSectionSetUpBC(*gsection)); 2454 if (nroots >= 0) { 2455 PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE)); 2456 PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE)); 2457 if (nroots > pEnd - pStart) { 2458 for (p = pStart; p < pEnd; ++p) { 2459 if (tmpOff[p] < 0) (*gsection)->atlasDof[p - pStart] = tmpOff[p]; 2460 } 2461 } 2462 } 2463 /* Calculate new sizes, get process offset, and calculate point offsets */ 2464 for (p = 0, off = 0; p < pEnd - pStart; ++p) { 2465 cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0; 2466 (*gsection)->atlasOff[p] = off; 2467 off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p] - cdof : 0; 2468 } 2469 PetscCallMPI(MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)s))); 2470 globalOff -= off; 2471 for (p = 0, off = 0; p < pEnd - pStart; ++p) { 2472 (*gsection)->atlasOff[p] += globalOff; 2473 if (neg) neg[p] = -((*gsection)->atlasOff[p] + 1); 2474 } 2475 /* Put in negative offsets for ghost points */ 2476 if (nroots >= 0) { 2477 PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE)); 2478 PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE)); 2479 if (nroots > pEnd - pStart) { 2480 for (p = pStart; p < pEnd; ++p) { 2481 if (tmpOff[p] < 0) (*gsection)->atlasOff[p - pStart] = tmpOff[p]; 2482 } 2483 } 2484 } 2485 if (nroots >= 0 && nroots > pEnd - pStart) PetscCall(PetscFree(tmpOff)); 2486 PetscCall(PetscFree(neg)); 2487 PetscFunctionReturn(PETSC_SUCCESS); 2488 } 2489 2490 typedef struct _n_PetscSectionSym_Label { 2491 DMLabel label; 2492 PetscCopyMode *modes; 2493 PetscInt *sizes; 2494 const PetscInt ***perms; 2495 const PetscScalar ***rots; 2496 PetscInt (*minMaxOrients)[2]; 2497 PetscInt numStrata; /* numStrata is only increasing, functions as a state */ 2498 } PetscSectionSym_Label; 2499 2500 static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym) 2501 { 2502 PetscInt i, j; 2503 PetscSectionSym_Label *sl = (PetscSectionSym_Label *)sym->data; 2504 2505 PetscFunctionBegin; 2506 for (i = 0; i <= sl->numStrata; i++) { 2507 if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) { 2508 for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 2509 if (sl->perms[i]) PetscCall(PetscFree(sl->perms[i][j])); 2510 if (sl->rots[i]) PetscCall(PetscFree(sl->rots[i][j])); 2511 } 2512 if (sl->perms[i]) { 2513 const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]]; 2514 2515 PetscCall(PetscFree(perms)); 2516 } 2517 if (sl->rots[i]) { 2518 const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]]; 2519 2520 PetscCall(PetscFree(rots)); 2521 } 2522 } 2523 } 2524 PetscCall(PetscFree5(sl->modes, sl->sizes, sl->perms, sl->rots, sl->minMaxOrients)); 2525 PetscCall(DMLabelDestroy(&sl->label)); 2526 sl->numStrata = 0; 2527 PetscFunctionReturn(PETSC_SUCCESS); 2528 } 2529 2530 static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym) 2531 { 2532 PetscFunctionBegin; 2533 PetscCall(PetscSectionSymLabelReset(sym)); 2534 PetscCall(PetscFree(sym->data)); 2535 PetscFunctionReturn(PETSC_SUCCESS); 2536 } 2537 2538 static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer) 2539 { 2540 PetscSectionSym_Label *sl = (PetscSectionSym_Label *)sym->data; 2541 PetscBool isAscii; 2542 DMLabel label = sl->label; 2543 const char *name; 2544 2545 PetscFunctionBegin; 2546 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isAscii)); 2547 if (isAscii) { 2548 PetscInt i, j, k; 2549 PetscViewerFormat format; 2550 2551 PetscCall(PetscViewerGetFormat(viewer, &format)); 2552 if (label) { 2553 PetscCall(PetscViewerGetFormat(viewer, &format)); 2554 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 2555 PetscCall(PetscViewerASCIIPushTab(viewer)); 2556 PetscCall(DMLabelView(label, viewer)); 2557 PetscCall(PetscViewerASCIIPopTab(viewer)); 2558 } else { 2559 PetscCall(PetscObjectGetName((PetscObject)sl->label, &name)); 2560 PetscCall(PetscViewerASCIIPrintf(viewer, " Label '%s'\n", name)); 2561 } 2562 } else { 2563 PetscCall(PetscViewerASCIIPrintf(viewer, "No label given\n")); 2564 } 2565 PetscCall(PetscViewerASCIIPushTab(viewer)); 2566 for (i = 0; i <= sl->numStrata; i++) { 2567 PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue; 2568 2569 if (!(sl->perms[i] || sl->rots[i])) { 2570 PetscCall(PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %" PetscInt_FMT " (%" PetscInt_FMT " dofs per point): no symmetries\n", value, sl->sizes[i])); 2571 } else { 2572 PetscCall(PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %" PetscInt_FMT " (%" PetscInt_FMT " dofs per point):\n", value, sl->sizes[i])); 2573 PetscCall(PetscViewerASCIIPushTab(viewer)); 2574 PetscCall(PetscViewerASCIIPrintf(viewer, "Orientation range: [%" PetscInt_FMT ", %" PetscInt_FMT ")\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1])); 2575 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 2576 PetscCall(PetscViewerASCIIPushTab(viewer)); 2577 for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 2578 if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) { 2579 PetscCall(PetscViewerASCIIPrintf(viewer, "Orientation %" PetscInt_FMT ": identity\n", j)); 2580 } else { 2581 PetscInt tab; 2582 2583 PetscCall(PetscViewerASCIIPrintf(viewer, "Orientation %" PetscInt_FMT ":\n", j)); 2584 PetscCall(PetscViewerASCIIPushTab(viewer)); 2585 PetscCall(PetscViewerASCIIGetTab(viewer, &tab)); 2586 if (sl->perms[i] && sl->perms[i][j]) { 2587 PetscCall(PetscViewerASCIIPrintf(viewer, "Permutation:")); 2588 PetscCall(PetscViewerASCIISetTab(viewer, 0)); 2589 for (k = 0; k < sl->sizes[i]; k++) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, sl->perms[i][j][k])); 2590 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 2591 PetscCall(PetscViewerASCIISetTab(viewer, tab)); 2592 } 2593 if (sl->rots[i] && sl->rots[i][j]) { 2594 PetscCall(PetscViewerASCIIPrintf(viewer, "Rotations: ")); 2595 PetscCall(PetscViewerASCIISetTab(viewer, 0)); 2596 #if defined(PETSC_USE_COMPLEX) 2597 for (k = 0; k < sl->sizes[i]; k++) PetscCall(PetscViewerASCIIPrintf(viewer, " %+g+i*%+g", (double)PetscRealPart(sl->rots[i][j][k]), (double)PetscImaginaryPart(sl->rots[i][j][k]))); 2598 #else 2599 for (k = 0; k < sl->sizes[i]; k++) PetscCall(PetscViewerASCIIPrintf(viewer, " %+g", (double)sl->rots[i][j][k])); 2600 #endif 2601 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 2602 PetscCall(PetscViewerASCIISetTab(viewer, tab)); 2603 } 2604 PetscCall(PetscViewerASCIIPopTab(viewer)); 2605 } 2606 } 2607 PetscCall(PetscViewerASCIIPopTab(viewer)); 2608 } 2609 PetscCall(PetscViewerASCIIPopTab(viewer)); 2610 } 2611 } 2612 PetscCall(PetscViewerASCIIPopTab(viewer)); 2613 } 2614 PetscFunctionReturn(PETSC_SUCCESS); 2615 } 2616 2617 /*@ 2618 PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries 2619 2620 Logically 2621 2622 Input Parameters: 2623 + sym - the section symmetries 2624 - label - the `DMLabel` describing the types of points 2625 2626 Level: developer: 2627 2628 .seealso: `DMLabel`, `DM`, `PetscSectionSymLabelSetStratum()`, `PetscSectionSymCreateLabel()`, `PetscSectionGetPointSyms()` 2629 @*/ 2630 PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label) 2631 { 2632 PetscSectionSym_Label *sl; 2633 2634 PetscFunctionBegin; 2635 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1); 2636 sl = (PetscSectionSym_Label *)sym->data; 2637 if (sl->label && sl->label != label) PetscCall(PetscSectionSymLabelReset(sym)); 2638 if (label) { 2639 sl->label = label; 2640 PetscCall(PetscObjectReference((PetscObject)label)); 2641 PetscCall(DMLabelGetNumValues(label, &sl->numStrata)); 2642 PetscCall(PetscMalloc5(sl->numStrata + 1, &sl->modes, sl->numStrata + 1, &sl->sizes, sl->numStrata + 1, &sl->perms, sl->numStrata + 1, &sl->rots, sl->numStrata + 1, &sl->minMaxOrients)); 2643 PetscCall(PetscMemzero((void *)sl->modes, (sl->numStrata + 1) * sizeof(PetscCopyMode))); 2644 PetscCall(PetscMemzero((void *)sl->sizes, (sl->numStrata + 1) * sizeof(PetscInt))); 2645 PetscCall(PetscMemzero((void *)sl->perms, (sl->numStrata + 1) * sizeof(const PetscInt **))); 2646 PetscCall(PetscMemzero((void *)sl->rots, (sl->numStrata + 1) * sizeof(const PetscScalar **))); 2647 PetscCall(PetscMemzero((void *)sl->minMaxOrients, (sl->numStrata + 1) * sizeof(PetscInt[2]))); 2648 } 2649 PetscFunctionReturn(PETSC_SUCCESS); 2650 } 2651 2652 /*@C 2653 PetscSectionSymLabelGetStratum - get the symmetries for the orientations of a stratum 2654 2655 Logically Collective 2656 2657 Input Parameters: 2658 + sym - the section symmetries 2659 - stratum - the stratum value in the label that we are assigning symmetries for 2660 2661 Output Parameters: 2662 + size - the number of dofs for points in the `stratum` of the label 2663 . minOrient - the smallest orientation for a point in this `stratum` 2664 . maxOrient - one greater than the largest orientation for a ppoint in this `stratum` (i.e., orientations are in the range [`minOrient`, `maxOrient`)) 2665 . perms - `NULL` if there are no permutations, or (`maxOrient` - `minOrient`) permutations, one for each orientation. A `NULL` permutation is the identity 2666 - rots - `NULL` if there are no rotations, or (`maxOrient` - `minOrient`) sets of rotations, one for each orientation. A `NULL` set of orientations is the identity 2667 2668 Level: developer 2669 2670 .seealso: `DMLabel`, `DM`, `PetscSectionSymLabelSetStratum()`, `PetscSectionSymCreate()`, `PetscSectionSetSym()`, `PetscSectionGetPointSyms()`, `PetscSectionSymCreateLabel()` 2671 @*/ 2672 PetscErrorCode PetscSectionSymLabelGetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt *size, PetscInt *minOrient, PetscInt *maxOrient, const PetscInt ***perms, const PetscScalar ***rots) 2673 { 2674 PetscSectionSym_Label *sl; 2675 const char *name; 2676 PetscInt i; 2677 2678 PetscFunctionBegin; 2679 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1); 2680 sl = (PetscSectionSym_Label *)sym->data; 2681 PetscCheck(sl->label, PetscObjectComm((PetscObject)sym), PETSC_ERR_ARG_WRONGSTATE, "No label set yet"); 2682 for (i = 0; i <= sl->numStrata; i++) { 2683 PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue; 2684 2685 if (stratum == value) break; 2686 } 2687 PetscCall(PetscObjectGetName((PetscObject)sl->label, &name)); 2688 PetscCheck(i <= sl->numStrata, PetscObjectComm((PetscObject)sym), PETSC_ERR_ARG_OUTOFRANGE, "Stratum %" PetscInt_FMT " not found in label %s", stratum, name); 2689 if (size) { 2690 PetscAssertPointer(size, 3); 2691 *size = sl->sizes[i]; 2692 } 2693 if (minOrient) { 2694 PetscAssertPointer(minOrient, 4); 2695 *minOrient = sl->minMaxOrients[i][0]; 2696 } 2697 if (maxOrient) { 2698 PetscAssertPointer(maxOrient, 5); 2699 *maxOrient = sl->minMaxOrients[i][1]; 2700 } 2701 if (perms) { 2702 PetscAssertPointer(perms, 6); 2703 *perms = PetscSafePointerPlusOffset(sl->perms[i], sl->minMaxOrients[i][0]); 2704 } 2705 if (rots) { 2706 PetscAssertPointer(rots, 7); 2707 *rots = PetscSafePointerPlusOffset(sl->rots[i], sl->minMaxOrients[i][0]); 2708 } 2709 PetscFunctionReturn(PETSC_SUCCESS); 2710 } 2711 2712 /*@C 2713 PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum 2714 2715 Logically 2716 2717 Input Parameters: 2718 + sym - the section symmetries 2719 . stratum - the stratum value in the label that we are assigning symmetries for 2720 . size - the number of dofs for points in the `stratum` of the label 2721 . minOrient - the smallest orientation for a point in this `stratum` 2722 . maxOrient - one greater than the largest orientation for a point in this `stratum` (i.e., orientations are in the range [`minOrient`, `maxOrient`)) 2723 . mode - how `sym` should copy the `perms` and `rots` arrays 2724 . perms - `NULL` if there are no permutations, or (`maxOrient` - `minOrient`) permutations, one for each orientation. A `NULL` permutation is the identity 2725 - rots - `NULL` if there are no rotations, or (`maxOrient` - `minOrient`) sets of rotations, one for each orientation. A `NULL` set of orientations is the identity 2726 2727 Level: developer 2728 2729 .seealso: `DMLabel`, `DM`, `PetscSectionSymLabelGetStratum()`, `PetscSectionSymCreate()`, `PetscSectionSetSym()`, `PetscSectionGetPointSyms()`, `PetscSectionSymCreateLabel()` 2730 @*/ 2731 PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots) 2732 { 2733 PetscSectionSym_Label *sl; 2734 const char *name; 2735 PetscInt i, j, k; 2736 2737 PetscFunctionBegin; 2738 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1); 2739 sl = (PetscSectionSym_Label *)sym->data; 2740 PetscCheck(sl->label, PetscObjectComm((PetscObject)sym), PETSC_ERR_ARG_WRONGSTATE, "No label set yet"); 2741 for (i = 0; i <= sl->numStrata; i++) { 2742 PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue; 2743 2744 if (stratum == value) break; 2745 } 2746 PetscCall(PetscObjectGetName((PetscObject)sl->label, &name)); 2747 PetscCheck(i <= sl->numStrata, PetscObjectComm((PetscObject)sym), PETSC_ERR_ARG_OUTOFRANGE, "Stratum %" PetscInt_FMT " not found in label %s", stratum, name); 2748 sl->sizes[i] = size; 2749 sl->modes[i] = mode; 2750 sl->minMaxOrients[i][0] = minOrient; 2751 sl->minMaxOrients[i][1] = maxOrient; 2752 if (mode == PETSC_COPY_VALUES) { 2753 if (perms) { 2754 PetscInt **ownPerms; 2755 2756 PetscCall(PetscCalloc1(maxOrient - minOrient, &ownPerms)); 2757 for (j = 0; j < maxOrient - minOrient; j++) { 2758 if (perms[j]) { 2759 PetscCall(PetscMalloc1(size, &ownPerms[j])); 2760 for (k = 0; k < size; k++) ownPerms[j][k] = perms[j][k]; 2761 } 2762 } 2763 sl->perms[i] = (const PetscInt **)&ownPerms[-minOrient]; 2764 } 2765 if (rots) { 2766 PetscScalar **ownRots; 2767 2768 PetscCall(PetscCalloc1(maxOrient - minOrient, &ownRots)); 2769 for (j = 0; j < maxOrient - minOrient; j++) { 2770 if (rots[j]) { 2771 PetscCall(PetscMalloc1(size, &ownRots[j])); 2772 for (k = 0; k < size; k++) ownRots[j][k] = rots[j][k]; 2773 } 2774 } 2775 sl->rots[i] = (const PetscScalar **)&ownRots[-minOrient]; 2776 } 2777 } else { 2778 sl->perms[i] = PetscSafePointerPlusOffset(perms, -minOrient); 2779 sl->rots[i] = PetscSafePointerPlusOffset(rots, -minOrient); 2780 } 2781 PetscFunctionReturn(PETSC_SUCCESS); 2782 } 2783 2784 static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots) 2785 { 2786 PetscInt i, j, numStrata; 2787 PetscSectionSym_Label *sl; 2788 DMLabel label; 2789 2790 PetscFunctionBegin; 2791 sl = (PetscSectionSym_Label *)sym->data; 2792 numStrata = sl->numStrata; 2793 label = sl->label; 2794 for (i = 0; i < numPoints; i++) { 2795 PetscInt point = points[2 * i]; 2796 PetscInt ornt = points[2 * i + 1]; 2797 2798 for (j = 0; j < numStrata; j++) { 2799 if (label->validIS[j]) { 2800 const PetscInt *indices; 2801 PetscInt k; 2802 2803 // For some reason, ISLocate() is insanely slow here 2804 PetscCall(ISGetIndices(label->points[j], &indices)); 2805 PetscCall(PetscFindInt(point, label->stratumSizes[j], indices, &k)); 2806 PetscCall(ISRestoreIndices(label->points[j], &indices)); 2807 if (k >= 0) break; 2808 } else { 2809 PetscBool has; 2810 2811 PetscCall(PetscHSetIHas(label->ht[j], point, &has)); 2812 if (has) break; 2813 } 2814 } 2815 PetscCheck(!(sl->minMaxOrients[j][1] > sl->minMaxOrients[j][0]) || !(ornt < sl->minMaxOrients[j][0] || ornt >= sl->minMaxOrients[j][1]), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "point %" PetscInt_FMT " orientation %" PetscInt_FMT " not in range [%" PetscInt_FMT ", %" PetscInt_FMT ") for stratum %" PetscInt_FMT, point, ornt, sl->minMaxOrients[j][0], sl->minMaxOrients[j][1], 2816 j < numStrata ? label->stratumValues[j] : label->defaultValue); 2817 if (perms) perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL; 2818 if (rots) rots[i] = sl->rots[j] ? sl->rots[j][ornt] : NULL; 2819 } 2820 PetscFunctionReturn(PETSC_SUCCESS); 2821 } 2822 2823 static PetscErrorCode PetscSectionSymCopy_Label(PetscSectionSym sym, PetscSectionSym nsym) 2824 { 2825 PetscSectionSym_Label *sl = (PetscSectionSym_Label *)nsym->data; 2826 IS valIS; 2827 const PetscInt *values; 2828 PetscInt Nv, v; 2829 2830 PetscFunctionBegin; 2831 PetscCall(DMLabelGetNumValues(sl->label, &Nv)); 2832 PetscCall(DMLabelGetValueIS(sl->label, &valIS)); 2833 PetscCall(ISGetIndices(valIS, &values)); 2834 for (v = 0; v < Nv; ++v) { 2835 const PetscInt val = values[v]; 2836 PetscInt size, minOrient, maxOrient; 2837 const PetscInt **perms; 2838 const PetscScalar **rots; 2839 2840 PetscCall(PetscSectionSymLabelGetStratum(sym, val, &size, &minOrient, &maxOrient, &perms, &rots)); 2841 PetscCall(PetscSectionSymLabelSetStratum(nsym, val, size, minOrient, maxOrient, PETSC_COPY_VALUES, perms, rots)); 2842 } 2843 PetscCall(ISDestroy(&valIS)); 2844 PetscFunctionReturn(PETSC_SUCCESS); 2845 } 2846 2847 static PetscErrorCode PetscSectionSymDistribute_Label(PetscSectionSym sym, PetscSF migrationSF, PetscSectionSym *dsym) 2848 { 2849 PetscSectionSym_Label *sl = (PetscSectionSym_Label *)sym->data; 2850 DMLabel dlabel; 2851 2852 PetscFunctionBegin; 2853 PetscCall(DMLabelDistribute(sl->label, migrationSF, &dlabel)); 2854 PetscCall(PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)sym), dlabel, dsym)); 2855 PetscCall(DMLabelDestroy(&dlabel)); 2856 PetscCall(PetscSectionSymCopy(sym, *dsym)); 2857 PetscFunctionReturn(PETSC_SUCCESS); 2858 } 2859 2860 PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym) 2861 { 2862 PetscSectionSym_Label *sl; 2863 2864 PetscFunctionBegin; 2865 PetscCall(PetscNew(&sl)); 2866 sym->ops->getpoints = PetscSectionSymGetPoints_Label; 2867 sym->ops->distribute = PetscSectionSymDistribute_Label; 2868 sym->ops->copy = PetscSectionSymCopy_Label; 2869 sym->ops->view = PetscSectionSymView_Label; 2870 sym->ops->destroy = PetscSectionSymDestroy_Label; 2871 sym->data = (void *)sl; 2872 PetscFunctionReturn(PETSC_SUCCESS); 2873 } 2874 2875 /*@ 2876 PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label 2877 2878 Collective 2879 2880 Input Parameters: 2881 + comm - the MPI communicator for the new symmetry 2882 - label - the label defining the strata 2883 2884 Output Parameter: 2885 . sym - the section symmetries 2886 2887 Level: developer 2888 2889 .seealso: `DMLabel`, `DM`, `PetscSectionSymCreate()`, `PetscSectionSetSym()`, `PetscSectionGetSym()`, `PetscSectionSymLabelSetStratum()`, `PetscSectionGetPointSyms()` 2890 @*/ 2891 PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym) 2892 { 2893 PetscFunctionBegin; 2894 PetscCall(DMInitializePackage()); 2895 PetscCall(PetscSectionSymCreate(comm, sym)); 2896 PetscCall(PetscSectionSymSetType(*sym, PETSCSECTIONSYMLABEL)); 2897 PetscCall(PetscSectionSymLabelSetLabel(*sym, label)); 2898 PetscFunctionReturn(PETSC_SUCCESS); 2899 } 2900