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 /*@C 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 /*@C 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 /*@C 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 PetscInt i; 947 948 PetscUseTypeMethod(label, getstratumis, v, &is); 949 PetscCall(ISLocate(is, point, &i)); 950 PetscCall(ISDestroy(&is)); 951 *contains = (PetscBool)(i >= 0); 952 } else { 953 PetscCall(PetscHSetIHas(label->ht[v], point, contains)); 954 } 955 } else { // value is not present 956 *contains = PETSC_FALSE; 957 } 958 } 959 PetscFunctionReturn(PETSC_SUCCESS); 960 } 961 962 /*@ 963 DMLabelGetDefaultValue - Get the default value returned by `DMLabelGetValue()` if a point has not been explicitly given a value. 964 When a label is created, it is initialized to -1. 965 966 Not Collective 967 968 Input Parameter: 969 . label - a `DMLabel` object 970 971 Output Parameter: 972 . defaultValue - the default value 973 974 Level: beginner 975 976 .seealso: `DMLabel`, `DM`, `DMLabelSetDefaultValue()`, `DMLabelGetValue()`, `DMLabelSetValue()` 977 @*/ 978 PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue) 979 { 980 PetscFunctionBegin; 981 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 982 *defaultValue = label->defaultValue; 983 PetscFunctionReturn(PETSC_SUCCESS); 984 } 985 986 /*@ 987 DMLabelSetDefaultValue - Set the default value returned by `DMLabelGetValue()` if a point has not been explicitly given a value. 988 When a label is created, it is initialized to -1. 989 990 Not Collective 991 992 Input Parameter: 993 . label - a `DMLabel` object 994 995 Output Parameter: 996 . defaultValue - the default value 997 998 Level: beginner 999 1000 .seealso: `DMLabel`, `DM`, `DMLabelGetDefaultValue()`, `DMLabelGetValue()`, `DMLabelSetValue()` 1001 @*/ 1002 PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue) 1003 { 1004 PetscFunctionBegin; 1005 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1006 label->defaultValue = defaultValue; 1007 PetscFunctionReturn(PETSC_SUCCESS); 1008 } 1009 1010 /*@ 1011 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 1012 `DMLabelSetDefaultValue()`) 1013 1014 Not Collective 1015 1016 Input Parameters: 1017 + label - the `DMLabel` 1018 - point - the point 1019 1020 Output Parameter: 1021 . value - The point value, or the default value (-1 by default) 1022 1023 Level: intermediate 1024 1025 Note: 1026 A label may assign multiple values to a point. No guarantees are made about which value is returned in that case. 1027 Use `DMLabelStratumHasPoint()` to check for inclusion in a specific value stratum. 1028 1029 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelSetValue()`, `DMLabelClearValue()`, `DMLabelGetDefaultValue()`, `DMLabelSetDefaultValue()` 1030 @*/ 1031 PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value) 1032 { 1033 PetscInt v; 1034 1035 PetscFunctionBeginHot; 1036 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1037 PetscAssertPointer(value, 3); 1038 *value = label->defaultValue; 1039 for (v = 0; v < label->numStrata; ++v) { 1040 if (label->validIS[v] || label->readonly) { 1041 IS is; 1042 PetscInt i; 1043 1044 PetscUseTypeMethod(label, getstratumis, v, &is); 1045 PetscCall(ISLocate(label->points[v], point, &i)); 1046 PetscCall(ISDestroy(&is)); 1047 if (i >= 0) { 1048 *value = label->stratumValues[v]; 1049 break; 1050 } 1051 } else { 1052 PetscBool has; 1053 1054 PetscCall(PetscHSetIHas(label->ht[v], point, &has)); 1055 if (has) { 1056 *value = label->stratumValues[v]; 1057 break; 1058 } 1059 } 1060 } 1061 PetscFunctionReturn(PETSC_SUCCESS); 1062 } 1063 1064 /*@ 1065 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 1066 be changed with `DMLabelSetDefaultValue()` to something different), then this function will do nothing. 1067 1068 Not Collective 1069 1070 Input Parameters: 1071 + label - the `DMLabel` 1072 . point - the point 1073 - value - The point value 1074 1075 Level: intermediate 1076 1077 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelClearValue()`, `DMLabelGetDefaultValue()`, `DMLabelSetDefaultValue()` 1078 @*/ 1079 PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value) 1080 { 1081 PetscInt v; 1082 1083 PetscFunctionBegin; 1084 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1085 /* Find label value, add new entry if needed */ 1086 if (value == label->defaultValue) PetscFunctionReturn(PETSC_SUCCESS); 1087 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1088 PetscCall(DMLabelLookupAddStratum(label, value, &v)); 1089 /* Set key */ 1090 PetscCall(DMLabelMakeInvalid_Private(label, v)); 1091 PetscCall(PetscHSetIAdd(label->ht[v], point)); 1092 PetscFunctionReturn(PETSC_SUCCESS); 1093 } 1094 1095 /*@ 1096 DMLabelClearValue - Clear the value a label assigns to a point 1097 1098 Not Collective 1099 1100 Input Parameters: 1101 + label - the `DMLabel` 1102 . point - the point 1103 - value - The point value 1104 1105 Level: intermediate 1106 1107 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()` 1108 @*/ 1109 PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value) 1110 { 1111 PetscInt v; 1112 1113 PetscFunctionBegin; 1114 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1115 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1116 /* Find label value */ 1117 PetscCall(DMLabelLookupStratum(label, value, &v)); 1118 if (v < 0) PetscFunctionReturn(PETSC_SUCCESS); 1119 1120 if (label->bt) { 1121 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); 1122 PetscCall(PetscBTClear(label->bt, point - label->pStart)); 1123 } 1124 1125 /* Delete key */ 1126 PetscCall(DMLabelMakeInvalid_Private(label, v)); 1127 PetscCall(PetscHSetIDel(label->ht[v], point)); 1128 PetscFunctionReturn(PETSC_SUCCESS); 1129 } 1130 1131 /*@ 1132 DMLabelInsertIS - Set all points in the `IS` to a value 1133 1134 Not Collective 1135 1136 Input Parameters: 1137 + label - the `DMLabel` 1138 . is - the point `IS` 1139 - value - The point value 1140 1141 Level: intermediate 1142 1143 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1144 @*/ 1145 PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value) 1146 { 1147 PetscInt v, n, p; 1148 const PetscInt *points; 1149 1150 PetscFunctionBegin; 1151 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1152 PetscValidHeaderSpecific(is, IS_CLASSID, 2); 1153 /* Find label value, add new entry if needed */ 1154 if (value == label->defaultValue) PetscFunctionReturn(PETSC_SUCCESS); 1155 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1156 PetscCall(DMLabelLookupAddStratum(label, value, &v)); 1157 /* Set keys */ 1158 PetscCall(DMLabelMakeInvalid_Private(label, v)); 1159 PetscCall(ISGetLocalSize(is, &n)); 1160 PetscCall(ISGetIndices(is, &points)); 1161 for (p = 0; p < n; ++p) PetscCall(PetscHSetIAdd(label->ht[v], points[p])); 1162 PetscCall(ISRestoreIndices(is, &points)); 1163 PetscFunctionReturn(PETSC_SUCCESS); 1164 } 1165 1166 /*@ 1167 DMLabelGetNumValues - Get the number of values that the `DMLabel` takes 1168 1169 Not Collective 1170 1171 Input Parameter: 1172 . label - the `DMLabel` 1173 1174 Output Parameter: 1175 . numValues - the number of values 1176 1177 Level: intermediate 1178 1179 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1180 @*/ 1181 PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues) 1182 { 1183 PetscFunctionBegin; 1184 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1185 PetscAssertPointer(numValues, 2); 1186 *numValues = label->numStrata; 1187 PetscFunctionReturn(PETSC_SUCCESS); 1188 } 1189 1190 /*@ 1191 DMLabelGetValueIS - Get an `IS` of all values that the `DMlabel` takes 1192 1193 Not Collective 1194 1195 Input Parameter: 1196 . label - the `DMLabel` 1197 1198 Output Parameter: 1199 . values - the value `IS` 1200 1201 Level: intermediate 1202 1203 Notes: 1204 The output `IS` should be destroyed when no longer needed. 1205 Strata which are allocated but empty [`DMLabelGetStratumSize()` yields 0] are counted. 1206 If you need to count only nonempty strata, use `DMLabelGetNonEmptyStratumValuesIS()`. 1207 1208 .seealso: `DMLabel`, `DM`, `DMLabelGetNonEmptyStratumValuesIS()`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1209 @*/ 1210 PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values) 1211 { 1212 PetscFunctionBegin; 1213 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1214 PetscAssertPointer(values, 2); 1215 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values)); 1216 PetscFunctionReturn(PETSC_SUCCESS); 1217 } 1218 1219 /*@ 1220 DMLabelGetNonEmptyStratumValuesIS - Get an `IS` of all values that the `DMlabel` takes 1221 1222 Not Collective 1223 1224 Input Parameter: 1225 . label - the `DMLabel` 1226 1227 Output Parameter: 1228 . values - the value `IS` 1229 1230 Level: intermediate 1231 1232 Notes: 1233 The output `IS` should be destroyed when no longer needed. 1234 This is similar to `DMLabelGetValueIS()` but counts only nonempty strata. 1235 1236 .seealso: `DMLabel`, `DM`, `DMLabelGetValueIS()`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1237 @*/ 1238 PetscErrorCode DMLabelGetNonEmptyStratumValuesIS(DMLabel label, IS *values) 1239 { 1240 PetscInt i, j; 1241 PetscInt *valuesArr; 1242 1243 PetscFunctionBegin; 1244 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1245 PetscAssertPointer(values, 2); 1246 PetscCall(PetscMalloc1(label->numStrata, &valuesArr)); 1247 for (i = 0, j = 0; i < label->numStrata; i++) { 1248 PetscInt n; 1249 1250 PetscCall(DMLabelGetStratumSize_Private(label, i, &n)); 1251 if (n) valuesArr[j++] = label->stratumValues[i]; 1252 } 1253 if (j == label->numStrata) { 1254 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values)); 1255 } else { 1256 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, j, valuesArr, PETSC_COPY_VALUES, values)); 1257 } 1258 PetscCall(PetscFree(valuesArr)); 1259 PetscFunctionReturn(PETSC_SUCCESS); 1260 } 1261 1262 /*@ 1263 DMLabelGetValueIndex - Get the index of a given value in the list of values for the `DMlabel`, or -1 if it is not present 1264 1265 Not Collective 1266 1267 Input Parameters: 1268 + label - the `DMLabel` 1269 - value - the value 1270 1271 Output Parameter: 1272 . index - the index of value in the list of values 1273 1274 Level: intermediate 1275 1276 .seealso: `DMLabel`, `DM`, `DMLabelGetValueIS()`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1277 @*/ 1278 PetscErrorCode DMLabelGetValueIndex(DMLabel label, PetscInt value, PetscInt *index) 1279 { 1280 PetscInt v; 1281 1282 PetscFunctionBegin; 1283 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1284 PetscAssertPointer(index, 3); 1285 /* Do not assume they are sorted */ 1286 for (v = 0; v < label->numStrata; ++v) 1287 if (label->stratumValues[v] == value) break; 1288 if (v >= label->numStrata) *index = -1; 1289 else *index = v; 1290 PetscFunctionReturn(PETSC_SUCCESS); 1291 } 1292 1293 /*@ 1294 DMLabelHasStratum - Determine whether points exist with the given value 1295 1296 Not Collective 1297 1298 Input Parameters: 1299 + label - the `DMLabel` 1300 - value - the stratum value 1301 1302 Output Parameter: 1303 . exists - Flag saying whether points exist 1304 1305 Level: intermediate 1306 1307 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1308 @*/ 1309 PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists) 1310 { 1311 PetscInt v; 1312 1313 PetscFunctionBegin; 1314 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1315 PetscAssertPointer(exists, 3); 1316 PetscCall(DMLabelLookupStratum(label, value, &v)); 1317 *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE; 1318 PetscFunctionReturn(PETSC_SUCCESS); 1319 } 1320 1321 /*@ 1322 DMLabelGetStratumSize - Get the size of a stratum 1323 1324 Not Collective 1325 1326 Input Parameters: 1327 + label - the `DMLabel` 1328 - value - the stratum value 1329 1330 Output Parameter: 1331 . size - The number of points in the stratum 1332 1333 Level: intermediate 1334 1335 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1336 @*/ 1337 PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size) 1338 { 1339 PetscInt v; 1340 1341 PetscFunctionBegin; 1342 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1343 PetscAssertPointer(size, 3); 1344 PetscCall(DMLabelLookupStratum(label, value, &v)); 1345 PetscCall(DMLabelGetStratumSize_Private(label, v, size)); 1346 PetscFunctionReturn(PETSC_SUCCESS); 1347 } 1348 1349 /*@ 1350 DMLabelGetStratumBounds - Get the largest and smallest point of a stratum 1351 1352 Not Collective 1353 1354 Input Parameters: 1355 + label - the `DMLabel` 1356 - value - the stratum value 1357 1358 Output Parameters: 1359 + start - the smallest point in the stratum 1360 - end - the largest point in the stratum 1361 1362 Level: intermediate 1363 1364 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1365 @*/ 1366 PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end) 1367 { 1368 IS is; 1369 PetscInt v, min, max; 1370 1371 PetscFunctionBegin; 1372 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1373 if (start) { 1374 PetscAssertPointer(start, 3); 1375 *start = -1; 1376 } 1377 if (end) { 1378 PetscAssertPointer(end, 4); 1379 *end = -1; 1380 } 1381 PetscCall(DMLabelLookupStratum(label, value, &v)); 1382 if (v < 0) PetscFunctionReturn(PETSC_SUCCESS); 1383 PetscCall(DMLabelMakeValid_Private(label, v)); 1384 if (label->stratumSizes[v] <= 0) PetscFunctionReturn(PETSC_SUCCESS); 1385 PetscUseTypeMethod(label, getstratumis, v, &is); 1386 PetscCall(ISGetMinMax(is, &min, &max)); 1387 PetscCall(ISDestroy(&is)); 1388 if (start) *start = min; 1389 if (end) *end = max + 1; 1390 PetscFunctionReturn(PETSC_SUCCESS); 1391 } 1392 1393 static PetscErrorCode DMLabelGetStratumIS_Concrete(DMLabel label, PetscInt v, IS *pointIS) 1394 { 1395 PetscFunctionBegin; 1396 PetscCall(PetscObjectReference((PetscObject)label->points[v])); 1397 *pointIS = label->points[v]; 1398 PetscFunctionReturn(PETSC_SUCCESS); 1399 } 1400 1401 /*@ 1402 DMLabelGetStratumIS - Get an `IS` with the stratum points 1403 1404 Not Collective 1405 1406 Input Parameters: 1407 + label - the `DMLabel` 1408 - value - the stratum value 1409 1410 Output Parameter: 1411 . points - The stratum points 1412 1413 Level: intermediate 1414 1415 Notes: 1416 The output `IS` should be destroyed when no longer needed. 1417 Returns `NULL` if the stratum is empty. 1418 1419 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1420 @*/ 1421 PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points) 1422 { 1423 PetscInt v; 1424 1425 PetscFunctionBegin; 1426 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1427 PetscAssertPointer(points, 3); 1428 *points = NULL; 1429 PetscCall(DMLabelLookupStratum(label, value, &v)); 1430 if (v < 0) PetscFunctionReturn(PETSC_SUCCESS); 1431 PetscCall(DMLabelMakeValid_Private(label, v)); 1432 PetscUseTypeMethod(label, getstratumis, v, points); 1433 PetscFunctionReturn(PETSC_SUCCESS); 1434 } 1435 1436 /*@ 1437 DMLabelSetStratumIS - Set the stratum points using an `IS` 1438 1439 Not Collective 1440 1441 Input Parameters: 1442 + label - the `DMLabel` 1443 . value - the stratum value 1444 - is - The stratum points 1445 1446 Level: intermediate 1447 1448 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1449 @*/ 1450 PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is) 1451 { 1452 PetscInt v; 1453 1454 PetscFunctionBegin; 1455 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1456 PetscValidHeaderSpecific(is, IS_CLASSID, 3); 1457 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1458 PetscCall(DMLabelLookupAddStratum(label, value, &v)); 1459 if (is == label->points[v]) PetscFunctionReturn(PETSC_SUCCESS); 1460 PetscCall(DMLabelClearStratum(label, value)); 1461 PetscCall(ISGetLocalSize(is, &label->stratumSizes[v])); 1462 PetscCall(PetscObjectReference((PetscObject)is)); 1463 PetscCall(ISDestroy(&label->points[v])); 1464 label->points[v] = is; 1465 label->validIS[v] = PETSC_TRUE; 1466 PetscCall(PetscObjectStateIncrease((PetscObject)label)); 1467 if (label->bt) { 1468 const PetscInt *points; 1469 PetscInt p; 1470 1471 PetscCall(ISGetIndices(is, &points)); 1472 for (p = 0; p < label->stratumSizes[v]; ++p) { 1473 const PetscInt point = points[p]; 1474 1475 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); 1476 PetscCall(PetscBTSet(label->bt, point - label->pStart)); 1477 } 1478 } 1479 PetscFunctionReturn(PETSC_SUCCESS); 1480 } 1481 1482 /*@ 1483 DMLabelClearStratum - Remove a stratum 1484 1485 Not Collective 1486 1487 Input Parameters: 1488 + label - the `DMLabel` 1489 - value - the stratum value 1490 1491 Level: intermediate 1492 1493 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1494 @*/ 1495 PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value) 1496 { 1497 PetscInt v; 1498 1499 PetscFunctionBegin; 1500 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1501 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1502 PetscCall(DMLabelLookupStratum(label, value, &v)); 1503 if (v < 0) PetscFunctionReturn(PETSC_SUCCESS); 1504 if (label->validIS[v]) { 1505 if (label->bt) { 1506 PetscInt i; 1507 const PetscInt *points; 1508 1509 PetscCall(ISGetIndices(label->points[v], &points)); 1510 for (i = 0; i < label->stratumSizes[v]; ++i) { 1511 const PetscInt point = points[i]; 1512 1513 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); 1514 PetscCall(PetscBTClear(label->bt, point - label->pStart)); 1515 } 1516 PetscCall(ISRestoreIndices(label->points[v], &points)); 1517 } 1518 label->stratumSizes[v] = 0; 1519 PetscCall(ISDestroy(&label->points[v])); 1520 PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &label->points[v])); 1521 PetscCall(PetscObjectSetName((PetscObject)label->points[v], "indices")); 1522 PetscCall(PetscObjectStateIncrease((PetscObject)label)); 1523 } else { 1524 PetscCall(PetscHSetIClear(label->ht[v])); 1525 } 1526 PetscFunctionReturn(PETSC_SUCCESS); 1527 } 1528 1529 /*@ 1530 DMLabelSetStratumBounds - Efficiently give a contiguous set of points a given label value 1531 1532 Not Collective 1533 1534 Input Parameters: 1535 + label - The `DMLabel` 1536 . value - The label value for all points 1537 . pStart - The first point 1538 - pEnd - A point beyond all marked points 1539 1540 Level: intermediate 1541 1542 Note: 1543 The marks points are [`pStart`, `pEnd`), and only the bounds are stored. 1544 1545 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelSetStratumIS()`, `DMLabelGetStratumIS()` 1546 @*/ 1547 PetscErrorCode DMLabelSetStratumBounds(DMLabel label, PetscInt value, PetscInt pStart, PetscInt pEnd) 1548 { 1549 IS pIS; 1550 1551 PetscFunctionBegin; 1552 PetscCall(ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &pIS)); 1553 PetscCall(DMLabelSetStratumIS(label, value, pIS)); 1554 PetscCall(ISDestroy(&pIS)); 1555 PetscFunctionReturn(PETSC_SUCCESS); 1556 } 1557 1558 /*@ 1559 DMLabelGetStratumPointIndex - Get the index of a point in a given stratum 1560 1561 Not Collective 1562 1563 Input Parameters: 1564 + label - The `DMLabel` 1565 . value - The label value 1566 - p - A point with this value 1567 1568 Output Parameter: 1569 . 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 1570 1571 Level: intermediate 1572 1573 .seealso: `DMLabel`, `DM`, `DMLabelGetValueIndex()`, `DMLabelGetStratumIS()`, `DMLabelCreate()` 1574 @*/ 1575 PetscErrorCode DMLabelGetStratumPointIndex(DMLabel label, PetscInt value, PetscInt p, PetscInt *index) 1576 { 1577 IS pointIS; 1578 const PetscInt *indices; 1579 PetscInt v; 1580 1581 PetscFunctionBegin; 1582 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1583 PetscAssertPointer(index, 4); 1584 *index = -1; 1585 PetscCall(DMLabelLookupStratum(label, value, &v)); 1586 if (v < 0) PetscFunctionReturn(PETSC_SUCCESS); 1587 PetscCall(DMLabelMakeValid_Private(label, v)); 1588 PetscUseTypeMethod(label, getstratumis, v, &pointIS); 1589 PetscCall(ISGetIndices(pointIS, &indices)); 1590 PetscCall(PetscFindInt(p, label->stratumSizes[v], indices, index)); 1591 PetscCall(ISRestoreIndices(pointIS, &indices)); 1592 PetscCall(ISDestroy(&pointIS)); 1593 PetscFunctionReturn(PETSC_SUCCESS); 1594 } 1595 1596 /*@ 1597 DMLabelFilter - Remove all points outside of [`start`, `end`) 1598 1599 Not Collective 1600 1601 Input Parameters: 1602 + label - the `DMLabel` 1603 . start - the first point kept 1604 - end - one more than the last point kept 1605 1606 Level: intermediate 1607 1608 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1609 @*/ 1610 PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end) 1611 { 1612 PetscInt v; 1613 1614 PetscFunctionBegin; 1615 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1616 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1617 PetscCall(DMLabelDestroyIndex(label)); 1618 PetscCall(DMLabelMakeAllValid_Private(label)); 1619 for (v = 0; v < label->numStrata; ++v) { 1620 PetscCall(ISGeneralFilter(label->points[v], start, end)); 1621 PetscCall(ISGetLocalSize(label->points[v], &label->stratumSizes[v])); 1622 } 1623 PetscCall(DMLabelCreateIndex(label, start, end)); 1624 PetscFunctionReturn(PETSC_SUCCESS); 1625 } 1626 1627 /*@ 1628 DMLabelPermute - Create a new label with permuted points 1629 1630 Not Collective 1631 1632 Input Parameters: 1633 + label - the `DMLabel` 1634 - permutation - the point permutation 1635 1636 Output Parameter: 1637 . labelNew - the new label containing the permuted points 1638 1639 Level: intermediate 1640 1641 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1642 @*/ 1643 PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew) 1644 { 1645 const PetscInt *perm; 1646 PetscInt numValues, numPoints, v, q; 1647 1648 PetscFunctionBegin; 1649 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1650 PetscValidHeaderSpecific(permutation, IS_CLASSID, 2); 1651 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1652 PetscCall(DMLabelMakeAllValid_Private(label)); 1653 PetscCall(DMLabelDuplicate(label, labelNew)); 1654 PetscCall(DMLabelGetNumValues(*labelNew, &numValues)); 1655 PetscCall(ISGetLocalSize(permutation, &numPoints)); 1656 PetscCall(ISGetIndices(permutation, &perm)); 1657 for (v = 0; v < numValues; ++v) { 1658 const PetscInt size = (*labelNew)->stratumSizes[v]; 1659 const PetscInt *points; 1660 PetscInt *pointsNew; 1661 1662 PetscCall(ISGetIndices((*labelNew)->points[v], &points)); 1663 PetscCall(PetscCalloc1(size, &pointsNew)); 1664 for (q = 0; q < size; ++q) { 1665 const PetscInt point = points[q]; 1666 1667 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); 1668 pointsNew[q] = perm[point]; 1669 } 1670 PetscCall(ISRestoreIndices((*labelNew)->points[v], &points)); 1671 PetscCall(PetscSortInt(size, pointsNew)); 1672 PetscCall(ISDestroy(&((*labelNew)->points[v]))); 1673 if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) { 1674 PetscCall(ISCreateStride(PETSC_COMM_SELF, size, pointsNew[0], 1, &((*labelNew)->points[v]))); 1675 PetscCall(PetscFree(pointsNew)); 1676 } else { 1677 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size, pointsNew, PETSC_OWN_POINTER, &((*labelNew)->points[v]))); 1678 } 1679 PetscCall(PetscObjectSetName((PetscObject)((*labelNew)->points[v]), "indices")); 1680 } 1681 PetscCall(ISRestoreIndices(permutation, &perm)); 1682 if (label->bt) { 1683 PetscCall(PetscBTDestroy(&label->bt)); 1684 PetscCall(DMLabelCreateIndex(label, label->pStart, label->pEnd)); 1685 } 1686 PetscFunctionReturn(PETSC_SUCCESS); 1687 } 1688 1689 static PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata) 1690 { 1691 MPI_Comm comm; 1692 PetscInt s, l, nroots, nleaves, offset, size; 1693 PetscInt *remoteOffsets, *rootStrata, *rootIdx; 1694 PetscSection rootSection; 1695 PetscSF labelSF; 1696 1697 PetscFunctionBegin; 1698 if (label) PetscCall(DMLabelMakeAllValid_Private(label)); 1699 PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 1700 /* Build a section of stratum values per point, generate the according SF 1701 and distribute point-wise stratum values to leaves. */ 1702 PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL)); 1703 PetscCall(PetscSectionCreate(comm, &rootSection)); 1704 PetscCall(PetscSectionSetChart(rootSection, 0, nroots)); 1705 if (label) { 1706 for (s = 0; s < label->numStrata; ++s) { 1707 const PetscInt *points; 1708 1709 PetscCall(ISGetIndices(label->points[s], &points)); 1710 for (l = 0; l < label->stratumSizes[s]; l++) PetscCall(PetscSectionAddDof(rootSection, points[l], 1)); 1711 PetscCall(ISRestoreIndices(label->points[s], &points)); 1712 } 1713 } 1714 PetscCall(PetscSectionSetUp(rootSection)); 1715 /* Create a point-wise array of stratum values */ 1716 PetscCall(PetscSectionGetStorageSize(rootSection, &size)); 1717 PetscCall(PetscMalloc1(size, &rootStrata)); 1718 PetscCall(PetscCalloc1(nroots, &rootIdx)); 1719 if (label) { 1720 for (s = 0; s < label->numStrata; ++s) { 1721 const PetscInt *points; 1722 1723 PetscCall(ISGetIndices(label->points[s], &points)); 1724 for (l = 0; l < label->stratumSizes[s]; l++) { 1725 const PetscInt p = points[l]; 1726 PetscCall(PetscSectionGetOffset(rootSection, p, &offset)); 1727 rootStrata[offset + rootIdx[p]++] = label->stratumValues[s]; 1728 } 1729 PetscCall(ISRestoreIndices(label->points[s], &points)); 1730 } 1731 } 1732 /* Build SF that maps label points to remote processes */ 1733 PetscCall(PetscSectionCreate(comm, leafSection)); 1734 PetscCall(PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection)); 1735 PetscCall(PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF)); 1736 PetscCall(PetscFree(remoteOffsets)); 1737 /* Send the strata for each point over the derived SF */ 1738 PetscCall(PetscSectionGetStorageSize(*leafSection, &size)); 1739 PetscCall(PetscMalloc1(size, leafStrata)); 1740 PetscCall(PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata, MPI_REPLACE)); 1741 PetscCall(PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata, MPI_REPLACE)); 1742 /* Clean up */ 1743 PetscCall(PetscFree(rootStrata)); 1744 PetscCall(PetscFree(rootIdx)); 1745 PetscCall(PetscSectionDestroy(&rootSection)); 1746 PetscCall(PetscSFDestroy(&labelSF)); 1747 PetscFunctionReturn(PETSC_SUCCESS); 1748 } 1749 1750 /*@ 1751 DMLabelDistribute - Create a new label pushed forward over the `PetscSF` 1752 1753 Collective 1754 1755 Input Parameters: 1756 + label - the `DMLabel` 1757 - sf - the map from old to new distribution 1758 1759 Output Parameter: 1760 . labelNew - the new redistributed label 1761 1762 Level: intermediate 1763 1764 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()` 1765 @*/ 1766 PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew) 1767 { 1768 MPI_Comm comm; 1769 PetscSection leafSection; 1770 PetscInt p, pStart, pEnd, s, size, dof, offset, stratum; 1771 PetscInt *leafStrata, *strataIdx; 1772 PetscInt **points; 1773 const char *lname = NULL; 1774 char *name; 1775 PetscInt nameSize; 1776 PetscHSetI stratumHash; 1777 size_t len = 0; 1778 PetscMPIInt rank; 1779 1780 PetscFunctionBegin; 1781 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 1782 if (label) { 1783 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1784 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1785 PetscCall(DMLabelMakeAllValid_Private(label)); 1786 } 1787 PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 1788 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1789 /* Bcast name */ 1790 if (rank == 0) { 1791 PetscCall(PetscObjectGetName((PetscObject)label, &lname)); 1792 PetscCall(PetscStrlen(lname, &len)); 1793 } 1794 nameSize = len; 1795 PetscCallMPI(MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm)); 1796 PetscCall(PetscMalloc1(nameSize + 1, &name)); 1797 if (rank == 0) PetscCall(PetscArraycpy(name, lname, nameSize + 1)); 1798 PetscCallMPI(MPI_Bcast(name, nameSize + 1, MPI_CHAR, 0, comm)); 1799 PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, labelNew)); 1800 PetscCall(PetscFree(name)); 1801 /* Bcast defaultValue */ 1802 if (rank == 0) (*labelNew)->defaultValue = label->defaultValue; 1803 PetscCallMPI(MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm)); 1804 /* Distribute stratum values over the SF and get the point mapping on the receiver */ 1805 PetscCall(DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata)); 1806 /* Determine received stratum values and initialise new label*/ 1807 PetscCall(PetscHSetICreate(&stratumHash)); 1808 PetscCall(PetscSectionGetStorageSize(leafSection, &size)); 1809 for (p = 0; p < size; ++p) PetscCall(PetscHSetIAdd(stratumHash, leafStrata[p])); 1810 PetscCall(PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata)); 1811 PetscCall(PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS)); 1812 for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE; 1813 PetscCall(PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues)); 1814 /* Turn leafStrata into indices rather than stratum values */ 1815 offset = 0; 1816 PetscCall(PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues)); 1817 PetscCall(PetscSortInt((*labelNew)->numStrata, (*labelNew)->stratumValues)); 1818 for (s = 0; s < (*labelNew)->numStrata; ++s) PetscCall(PetscHMapISet((*labelNew)->hmap, (*labelNew)->stratumValues[s], s)); 1819 for (p = 0; p < size; ++p) { 1820 for (s = 0; s < (*labelNew)->numStrata; ++s) { 1821 if (leafStrata[p] == (*labelNew)->stratumValues[s]) { 1822 leafStrata[p] = s; 1823 break; 1824 } 1825 } 1826 } 1827 /* Rebuild the point strata on the receiver */ 1828 PetscCall(PetscCalloc1((*labelNew)->numStrata, &(*labelNew)->stratumSizes)); 1829 PetscCall(PetscSectionGetChart(leafSection, &pStart, &pEnd)); 1830 for (p = pStart; p < pEnd; p++) { 1831 PetscCall(PetscSectionGetDof(leafSection, p, &dof)); 1832 PetscCall(PetscSectionGetOffset(leafSection, p, &offset)); 1833 for (s = 0; s < dof; s++) (*labelNew)->stratumSizes[leafStrata[offset + s]]++; 1834 } 1835 PetscCall(PetscCalloc1((*labelNew)->numStrata, &(*labelNew)->ht)); 1836 PetscCall(PetscCalloc1((*labelNew)->numStrata, &(*labelNew)->points)); 1837 PetscCall(PetscCalloc1((*labelNew)->numStrata, &points)); 1838 for (s = 0; s < (*labelNew)->numStrata; ++s) { 1839 PetscCall(PetscHSetICreate(&(*labelNew)->ht[s])); 1840 PetscCall(PetscMalloc1((*labelNew)->stratumSizes[s], &points[s])); 1841 } 1842 /* Insert points into new strata */ 1843 PetscCall(PetscCalloc1((*labelNew)->numStrata, &strataIdx)); 1844 PetscCall(PetscSectionGetChart(leafSection, &pStart, &pEnd)); 1845 for (p = pStart; p < pEnd; p++) { 1846 PetscCall(PetscSectionGetDof(leafSection, p, &dof)); 1847 PetscCall(PetscSectionGetOffset(leafSection, p, &offset)); 1848 for (s = 0; s < dof; s++) { 1849 stratum = leafStrata[offset + s]; 1850 points[stratum][strataIdx[stratum]++] = p; 1851 } 1852 } 1853 for (s = 0; s < (*labelNew)->numStrata; s++) { 1854 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, (*labelNew)->stratumSizes[s], &points[s][0], PETSC_OWN_POINTER, &((*labelNew)->points[s]))); 1855 PetscCall(PetscObjectSetName((PetscObject)((*labelNew)->points[s]), "indices")); 1856 } 1857 PetscCall(PetscFree(points)); 1858 PetscCall(PetscHSetIDestroy(&stratumHash)); 1859 PetscCall(PetscFree(leafStrata)); 1860 PetscCall(PetscFree(strataIdx)); 1861 PetscCall(PetscSectionDestroy(&leafSection)); 1862 PetscFunctionReturn(PETSC_SUCCESS); 1863 } 1864 1865 /*@ 1866 DMLabelGather - Gather all label values from leafs into roots 1867 1868 Collective 1869 1870 Input Parameters: 1871 + label - the `DMLabel` 1872 - sf - the `PetscSF` communication map 1873 1874 Output Parameter: 1875 . labelNew - the new `DMLabel` with localised leaf values 1876 1877 Level: developer 1878 1879 Note: 1880 This is the inverse operation to `DMLabelDistribute()`. 1881 1882 .seealso: `DMLabel`, `DM`, `DMLabelDistribute()` 1883 @*/ 1884 PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew) 1885 { 1886 MPI_Comm comm; 1887 PetscSection rootSection; 1888 PetscSF sfLabel; 1889 PetscSFNode *rootPoints, *leafPoints; 1890 PetscInt p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset; 1891 const PetscInt *rootDegree, *ilocal; 1892 PetscInt *rootStrata; 1893 const char *lname; 1894 char *name; 1895 PetscInt nameSize; 1896 size_t len = 0; 1897 PetscMPIInt rank, size; 1898 1899 PetscFunctionBegin; 1900 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1901 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 1902 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 1903 PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 1904 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1905 PetscCallMPI(MPI_Comm_size(comm, &size)); 1906 /* Bcast name */ 1907 if (rank == 0) { 1908 PetscCall(PetscObjectGetName((PetscObject)label, &lname)); 1909 PetscCall(PetscStrlen(lname, &len)); 1910 } 1911 nameSize = len; 1912 PetscCallMPI(MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm)); 1913 PetscCall(PetscMalloc1(nameSize + 1, &name)); 1914 if (rank == 0) PetscCall(PetscArraycpy(name, lname, nameSize + 1)); 1915 PetscCallMPI(MPI_Bcast(name, nameSize + 1, MPI_CHAR, 0, comm)); 1916 PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, labelNew)); 1917 PetscCall(PetscFree(name)); 1918 /* Gather rank/index pairs of leaves into local roots to build 1919 an inverse, multi-rooted SF. Note that this ignores local leaf 1920 indexing due to the use of the multiSF in PetscSFGather. */ 1921 PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL)); 1922 PetscCall(PetscMalloc1(nroots, &leafPoints)); 1923 for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1; 1924 for (p = 0; p < nleaves; p++) { 1925 PetscInt ilp = ilocal ? ilocal[p] : p; 1926 1927 leafPoints[ilp].index = ilp; 1928 leafPoints[ilp].rank = rank; 1929 } 1930 PetscCall(PetscSFComputeDegreeBegin(sf, &rootDegree)); 1931 PetscCall(PetscSFComputeDegreeEnd(sf, &rootDegree)); 1932 for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p]; 1933 PetscCall(PetscMalloc1(nmultiroots, &rootPoints)); 1934 PetscCall(PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints)); 1935 PetscCall(PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints)); 1936 PetscCall(PetscSFCreate(comm, &sfLabel)); 1937 PetscCall(PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER)); 1938 /* Migrate label over inverted SF to pull stratum values at leaves into roots. */ 1939 PetscCall(DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata)); 1940 /* Rebuild the point strata on the receiver */ 1941 for (p = 0, idx = 0; p < nroots; p++) { 1942 for (d = 0; d < rootDegree[p]; d++) { 1943 PetscCall(PetscSectionGetDof(rootSection, idx + d, &dof)); 1944 PetscCall(PetscSectionGetOffset(rootSection, idx + d, &offset)); 1945 for (s = 0; s < dof; s++) PetscCall(DMLabelSetValue(*labelNew, p, rootStrata[offset + s])); 1946 } 1947 idx += rootDegree[p]; 1948 } 1949 PetscCall(PetscFree(leafPoints)); 1950 PetscCall(PetscFree(rootStrata)); 1951 PetscCall(PetscSectionDestroy(&rootSection)); 1952 PetscCall(PetscSFDestroy(&sfLabel)); 1953 PetscFunctionReturn(PETSC_SUCCESS); 1954 } 1955 1956 static PetscErrorCode DMLabelPropagateInit_Internal(DMLabel label, PetscSF pointSF, PetscInt valArray[]) 1957 { 1958 const PetscInt *degree; 1959 const PetscInt *points; 1960 PetscInt Nr, r, Nl, l, val, defVal; 1961 1962 PetscFunctionBegin; 1963 PetscCall(DMLabelGetDefaultValue(label, &defVal)); 1964 /* Add in leaves */ 1965 PetscCall(PetscSFGetGraph(pointSF, &Nr, &Nl, &points, NULL)); 1966 for (l = 0; l < Nl; ++l) { 1967 PetscCall(DMLabelGetValue(label, points[l], &val)); 1968 if (val != defVal) valArray[points[l]] = val; 1969 } 1970 /* Add in shared roots */ 1971 PetscCall(PetscSFComputeDegreeBegin(pointSF, °ree)); 1972 PetscCall(PetscSFComputeDegreeEnd(pointSF, °ree)); 1973 for (r = 0; r < Nr; ++r) { 1974 if (degree[r]) { 1975 PetscCall(DMLabelGetValue(label, r, &val)); 1976 if (val != defVal) valArray[r] = val; 1977 } 1978 } 1979 PetscFunctionReturn(PETSC_SUCCESS); 1980 } 1981 1982 static PetscErrorCode DMLabelPropagateFini_Internal(DMLabel label, PetscSF pointSF, PetscInt valArray[], PetscErrorCode (*markPoint)(DMLabel, PetscInt, PetscInt, void *), void *ctx) 1983 { 1984 const PetscInt *degree; 1985 const PetscInt *points; 1986 PetscInt Nr, r, Nl, l, val, defVal; 1987 1988 PetscFunctionBegin; 1989 PetscCall(DMLabelGetDefaultValue(label, &defVal)); 1990 /* Read out leaves */ 1991 PetscCall(PetscSFGetGraph(pointSF, &Nr, &Nl, &points, NULL)); 1992 for (l = 0; l < Nl; ++l) { 1993 const PetscInt p = points[l]; 1994 const PetscInt cval = valArray[p]; 1995 1996 if (cval != defVal) { 1997 PetscCall(DMLabelGetValue(label, p, &val)); 1998 if (val == defVal) { 1999 PetscCall(DMLabelSetValue(label, p, cval)); 2000 if (markPoint) PetscCall((*markPoint)(label, p, cval, ctx)); 2001 } 2002 } 2003 } 2004 /* Read out shared roots */ 2005 PetscCall(PetscSFComputeDegreeBegin(pointSF, °ree)); 2006 PetscCall(PetscSFComputeDegreeEnd(pointSF, °ree)); 2007 for (r = 0; r < Nr; ++r) { 2008 if (degree[r]) { 2009 const PetscInt cval = valArray[r]; 2010 2011 if (cval != defVal) { 2012 PetscCall(DMLabelGetValue(label, r, &val)); 2013 if (val == defVal) { 2014 PetscCall(DMLabelSetValue(label, r, cval)); 2015 if (markPoint) PetscCall((*markPoint)(label, r, cval, ctx)); 2016 } 2017 } 2018 } 2019 } 2020 PetscFunctionReturn(PETSC_SUCCESS); 2021 } 2022 2023 /*@ 2024 DMLabelPropagateBegin - Setup a cycle of label propagation 2025 2026 Collective 2027 2028 Input Parameters: 2029 + label - The `DMLabel` to propagate across processes 2030 - sf - The `PetscSF` describing parallel layout of the label points 2031 2032 Level: intermediate 2033 2034 .seealso: `DMLabel`, `DM`, `DMLabelPropagateEnd()`, `DMLabelPropagatePush()` 2035 @*/ 2036 PetscErrorCode DMLabelPropagateBegin(DMLabel label, PetscSF sf) 2037 { 2038 PetscInt Nr, r, defVal; 2039 PetscMPIInt size; 2040 2041 PetscFunctionBegin; 2042 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 2043 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sf), &size)); 2044 if (size > 1) { 2045 PetscCall(DMLabelGetDefaultValue(label, &defVal)); 2046 PetscCall(PetscSFGetGraph(sf, &Nr, NULL, NULL, NULL)); 2047 if (Nr >= 0) PetscCall(PetscMalloc1(Nr, &label->propArray)); 2048 for (r = 0; r < Nr; ++r) label->propArray[r] = defVal; 2049 } 2050 PetscFunctionReturn(PETSC_SUCCESS); 2051 } 2052 2053 /*@ 2054 DMLabelPropagateEnd - Tear down a cycle of label propagation 2055 2056 Collective 2057 2058 Input Parameters: 2059 + label - The `DMLabel` to propagate across processes 2060 - pointSF - The `PetscSF` describing parallel layout of the label points 2061 2062 Level: intermediate 2063 2064 .seealso: `DMLabel`, `DM`, `DMLabelPropagateBegin()`, `DMLabelPropagatePush()` 2065 @*/ 2066 PetscErrorCode DMLabelPropagateEnd(DMLabel label, PetscSF pointSF) 2067 { 2068 PetscFunctionBegin; 2069 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 2070 PetscCall(PetscFree(label->propArray)); 2071 label->propArray = NULL; 2072 PetscFunctionReturn(PETSC_SUCCESS); 2073 } 2074 2075 /*@C 2076 DMLabelPropagatePush - Tear down a cycle of label propagation 2077 2078 Collective 2079 2080 Input Parameters: 2081 + label - The `DMLabel` to propagate across processes 2082 . pointSF - The `PetscSF` describing parallel layout of the label points 2083 . markPoint - An optional callback that is called when a point is marked, or `NULL` 2084 - ctx - An optional user context for the callback, or `NULL` 2085 2086 Calling sequence of `markPoint`: 2087 + label - The `DMLabel` 2088 . p - The point being marked 2089 . val - The label value for `p` 2090 - ctx - An optional user context 2091 2092 Level: intermediate 2093 2094 .seealso: `DMLabel`, `DM`, `DMLabelPropagateBegin()`, `DMLabelPropagateEnd()` 2095 @*/ 2096 PetscErrorCode DMLabelPropagatePush(DMLabel label, PetscSF pointSF, PetscErrorCode (*markPoint)(DMLabel label, PetscInt p, PetscInt val, void *ctx), void *ctx) 2097 { 2098 PetscInt *valArray = label->propArray, Nr; 2099 PetscMPIInt size; 2100 2101 PetscFunctionBegin; 2102 PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered"); 2103 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pointSF), &size)); 2104 PetscCall(PetscSFGetGraph(pointSF, &Nr, NULL, NULL, NULL)); 2105 if (size > 1 && Nr >= 0) { 2106 /* Communicate marked edges 2107 The current implementation allocates an array the size of the number of root. We put the label values into the 2108 array, and then call PetscSFReduce()+PetscSFBcast() to make the marks consistent. 2109 2110 TODO: We could use in-place communication with a different SF 2111 We use MPI_SUM for the Reduce, and check the result against the rootdegree. If sum >= rootdegree+1, then the edge has 2112 already been marked. If not, it might have been handled on the process in this round, but we add it anyway. 2113 2114 In order to update the queue with the new edges from the label communication, we use BcastAnOp(MPI_SUM), so that new 2115 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 2116 edge to the queue. 2117 */ 2118 PetscCall(DMLabelPropagateInit_Internal(label, pointSF, valArray)); 2119 PetscCall(PetscSFReduceBegin(pointSF, MPIU_INT, valArray, valArray, MPI_MAX)); 2120 PetscCall(PetscSFReduceEnd(pointSF, MPIU_INT, valArray, valArray, MPI_MAX)); 2121 PetscCall(PetscSFBcastBegin(pointSF, MPIU_INT, valArray, valArray, MPI_REPLACE)); 2122 PetscCall(PetscSFBcastEnd(pointSF, MPIU_INT, valArray, valArray, MPI_REPLACE)); 2123 PetscCall(DMLabelPropagateFini_Internal(label, pointSF, valArray, markPoint, ctx)); 2124 } 2125 PetscFunctionReturn(PETSC_SUCCESS); 2126 } 2127 2128 /*@ 2129 DMLabelConvertToSection - Make a `PetscSection`/`IS` pair that encodes the label 2130 2131 Not Collective 2132 2133 Input Parameter: 2134 . label - the `DMLabel` 2135 2136 Output Parameters: 2137 + section - the section giving offsets for each stratum 2138 - is - An `IS` containing all the label points 2139 2140 Level: developer 2141 2142 .seealso: `DMLabel`, `DM`, `DMLabelDistribute()` 2143 @*/ 2144 PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is) 2145 { 2146 IS vIS; 2147 const PetscInt *values; 2148 PetscInt *points; 2149 PetscInt nV, vS = 0, vE = 0, v, N; 2150 2151 PetscFunctionBegin; 2152 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 2153 PetscCall(DMLabelGetNumValues(label, &nV)); 2154 PetscCall(DMLabelGetValueIS(label, &vIS)); 2155 PetscCall(ISGetIndices(vIS, &values)); 2156 if (nV) { 2157 vS = values[0]; 2158 vE = values[0] + 1; 2159 } 2160 for (v = 1; v < nV; ++v) { 2161 vS = PetscMin(vS, values[v]); 2162 vE = PetscMax(vE, values[v] + 1); 2163 } 2164 PetscCall(PetscSectionCreate(PETSC_COMM_SELF, section)); 2165 PetscCall(PetscSectionSetChart(*section, vS, vE)); 2166 for (v = 0; v < nV; ++v) { 2167 PetscInt n; 2168 2169 PetscCall(DMLabelGetStratumSize(label, values[v], &n)); 2170 PetscCall(PetscSectionSetDof(*section, values[v], n)); 2171 } 2172 PetscCall(PetscSectionSetUp(*section)); 2173 PetscCall(PetscSectionGetStorageSize(*section, &N)); 2174 PetscCall(PetscMalloc1(N, &points)); 2175 for (v = 0; v < nV; ++v) { 2176 IS is; 2177 const PetscInt *spoints; 2178 PetscInt dof, off, p; 2179 2180 PetscCall(PetscSectionGetDof(*section, values[v], &dof)); 2181 PetscCall(PetscSectionGetOffset(*section, values[v], &off)); 2182 PetscCall(DMLabelGetStratumIS(label, values[v], &is)); 2183 PetscCall(ISGetIndices(is, &spoints)); 2184 for (p = 0; p < dof; ++p) points[off + p] = spoints[p]; 2185 PetscCall(ISRestoreIndices(is, &spoints)); 2186 PetscCall(ISDestroy(&is)); 2187 } 2188 PetscCall(ISRestoreIndices(vIS, &values)); 2189 PetscCall(ISDestroy(&vIS)); 2190 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is)); 2191 PetscFunctionReturn(PETSC_SUCCESS); 2192 } 2193 2194 /*@C 2195 DMLabelRegister - Adds a new label component implementation 2196 2197 Not Collective 2198 2199 Input Parameters: 2200 + name - The name of a new user-defined creation routine 2201 - create_func - The creation routine itself 2202 2203 Notes: 2204 `DMLabelRegister()` may be called multiple times to add several user-defined labels 2205 2206 Example Usage: 2207 .vb 2208 DMLabelRegister("my_label", MyLabelCreate); 2209 .ve 2210 2211 Then, your label type can be chosen with the procedural interface via 2212 .vb 2213 DMLabelCreate(MPI_Comm, DMLabel *); 2214 DMLabelSetType(DMLabel, "my_label"); 2215 .ve 2216 or at runtime via the option 2217 .vb 2218 -dm_label_type my_label 2219 .ve 2220 2221 Level: advanced 2222 2223 .seealso: `DMLabel`, `DM`, `DMLabelType`, `DMLabelRegisterAll()`, `DMLabelRegisterDestroy()` 2224 @*/ 2225 PetscErrorCode DMLabelRegister(const char name[], PetscErrorCode (*create_func)(DMLabel)) 2226 { 2227 PetscFunctionBegin; 2228 PetscCall(DMInitializePackage()); 2229 PetscCall(PetscFunctionListAdd(&DMLabelList, name, create_func)); 2230 PetscFunctionReturn(PETSC_SUCCESS); 2231 } 2232 2233 PETSC_EXTERN PetscErrorCode DMLabelCreate_Concrete(DMLabel); 2234 PETSC_EXTERN PetscErrorCode DMLabelCreate_Ephemeral(DMLabel); 2235 2236 /*@C 2237 DMLabelRegisterAll - Registers all of the `DMLabel` implementations in the `DM` package. 2238 2239 Not Collective 2240 2241 Level: advanced 2242 2243 .seealso: `DMLabel`, `DM`, `DMRegisterAll()`, `DMLabelRegisterDestroy()` 2244 @*/ 2245 PetscErrorCode DMLabelRegisterAll(void) 2246 { 2247 PetscFunctionBegin; 2248 if (DMLabelRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 2249 DMLabelRegisterAllCalled = PETSC_TRUE; 2250 2251 PetscCall(DMLabelRegister(DMLABELCONCRETE, DMLabelCreate_Concrete)); 2252 PetscCall(DMLabelRegister(DMLABELEPHEMERAL, DMLabelCreate_Ephemeral)); 2253 PetscFunctionReturn(PETSC_SUCCESS); 2254 } 2255 2256 /*@C 2257 DMLabelRegisterDestroy - This function destroys the `DMLabel` registry. It is called from `PetscFinalize()`. 2258 2259 Level: developer 2260 2261 .seealso: `DMLabel`, `DM`, `PetscInitialize()` 2262 @*/ 2263 PetscErrorCode DMLabelRegisterDestroy(void) 2264 { 2265 PetscFunctionBegin; 2266 PetscCall(PetscFunctionListDestroy(&DMLabelList)); 2267 DMLabelRegisterAllCalled = PETSC_FALSE; 2268 PetscFunctionReturn(PETSC_SUCCESS); 2269 } 2270 2271 /*@C 2272 DMLabelSetType - Sets the particular implementation for a label. 2273 2274 Collective 2275 2276 Input Parameters: 2277 + label - The label 2278 - method - The name of the label type 2279 2280 Options Database Key: 2281 . -dm_label_type <type> - Sets the label type; use -help for a list of available types or see `DMLabelType` 2282 2283 Level: intermediate 2284 2285 .seealso: `DMLabel`, `DM`, `DMLabelGetType()`, `DMLabelCreate()` 2286 @*/ 2287 PetscErrorCode DMLabelSetType(DMLabel label, DMLabelType method) 2288 { 2289 PetscErrorCode (*r)(DMLabel); 2290 PetscBool match; 2291 2292 PetscFunctionBegin; 2293 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 2294 PetscCall(PetscObjectTypeCompare((PetscObject)label, method, &match)); 2295 if (match) PetscFunctionReturn(PETSC_SUCCESS); 2296 2297 PetscCall(DMLabelRegisterAll()); 2298 PetscCall(PetscFunctionListFind(DMLabelList, method, &r)); 2299 PetscCheck(r, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DMLabel type: %s", method); 2300 2301 PetscTryTypeMethod(label, destroy); 2302 PetscCall(PetscMemzero(label->ops, sizeof(*label->ops))); 2303 PetscCall(PetscObjectChangeTypeName((PetscObject)label, method)); 2304 PetscCall((*r)(label)); 2305 PetscFunctionReturn(PETSC_SUCCESS); 2306 } 2307 2308 /*@C 2309 DMLabelGetType - Gets the type name (as a string) from the label. 2310 2311 Not Collective 2312 2313 Input Parameter: 2314 . label - The `DMLabel` 2315 2316 Output Parameter: 2317 . type - The `DMLabel` type name 2318 2319 Level: intermediate 2320 2321 .seealso: `DMLabel`, `DM`, `DMLabelSetType()`, `DMLabelCreate()` 2322 @*/ 2323 PetscErrorCode DMLabelGetType(DMLabel label, DMLabelType *type) 2324 { 2325 PetscFunctionBegin; 2326 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 2327 PetscAssertPointer(type, 2); 2328 PetscCall(DMLabelRegisterAll()); 2329 *type = ((PetscObject)label)->type_name; 2330 PetscFunctionReturn(PETSC_SUCCESS); 2331 } 2332 2333 static PetscErrorCode DMLabelInitialize_Concrete(DMLabel label) 2334 { 2335 PetscFunctionBegin; 2336 label->ops->view = DMLabelView_Concrete; 2337 label->ops->setup = NULL; 2338 label->ops->duplicate = DMLabelDuplicate_Concrete; 2339 label->ops->getstratumis = DMLabelGetStratumIS_Concrete; 2340 PetscFunctionReturn(PETSC_SUCCESS); 2341 } 2342 2343 PETSC_EXTERN PetscErrorCode DMLabelCreate_Concrete(DMLabel label) 2344 { 2345 PetscFunctionBegin; 2346 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 2347 PetscCall(DMLabelInitialize_Concrete(label)); 2348 PetscFunctionReturn(PETSC_SUCCESS); 2349 } 2350 2351 /*@ 2352 PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using 2353 the local section and an `PetscSF` describing the section point overlap. 2354 2355 Collective 2356 2357 Input Parameters: 2358 + s - The `PetscSection` for the local field layout 2359 . sf - The `PetscSF` describing parallel layout of the section points 2360 . includeConstraints - By default this is `PETSC_FALSE`, meaning that the global field vector will not possess constrained dofs 2361 . label - The label specifying the points 2362 - labelValue - The label stratum specifying the points 2363 2364 Output Parameter: 2365 . gsection - The `PetscSection` for the global field layout 2366 2367 Level: developer 2368 2369 Note: 2370 This gives negative sizes and offsets to points not owned by this process 2371 2372 .seealso: `DMLabel`, `DM`, `PetscSectionCreate()` 2373 @*/ 2374 PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection) 2375 { 2376 PetscInt *neg = NULL, *tmpOff = NULL; 2377 PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots; 2378 2379 PetscFunctionBegin; 2380 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2381 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 2382 PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4); 2383 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)s), gsection)); 2384 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 2385 PetscCall(PetscSectionSetChart(*gsection, pStart, pEnd)); 2386 PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 2387 if (nroots >= 0) { 2388 PetscCheck(nroots >= pEnd - pStart, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %" PetscInt_FMT " < %" PetscInt_FMT " section size", nroots, pEnd - pStart); 2389 PetscCall(PetscCalloc1(nroots, &neg)); 2390 if (nroots > pEnd - pStart) { 2391 PetscCall(PetscCalloc1(nroots, &tmpOff)); 2392 } else { 2393 tmpOff = &(*gsection)->atlasDof[-pStart]; 2394 } 2395 } 2396 /* Mark ghost points with negative dof */ 2397 for (p = pStart; p < pEnd; ++p) { 2398 PetscInt value; 2399 2400 PetscCall(DMLabelGetValue(label, p, &value)); 2401 if (value != labelValue) continue; 2402 PetscCall(PetscSectionGetDof(s, p, &dof)); 2403 PetscCall(PetscSectionSetDof(*gsection, p, dof)); 2404 PetscCall(PetscSectionGetConstraintDof(s, p, &cdof)); 2405 if (!includeConstraints && cdof > 0) PetscCall(PetscSectionSetConstraintDof(*gsection, p, cdof)); 2406 if (neg) neg[p] = -(dof + 1); 2407 } 2408 PetscCall(PetscSectionSetUpBC(*gsection)); 2409 if (nroots >= 0) { 2410 PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE)); 2411 PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE)); 2412 if (nroots > pEnd - pStart) { 2413 for (p = pStart; p < pEnd; ++p) { 2414 if (tmpOff[p] < 0) (*gsection)->atlasDof[p - pStart] = tmpOff[p]; 2415 } 2416 } 2417 } 2418 /* Calculate new sizes, get process offset, and calculate point offsets */ 2419 for (p = 0, off = 0; p < pEnd - pStart; ++p) { 2420 cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0; 2421 (*gsection)->atlasOff[p] = off; 2422 off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p] - cdof : 0; 2423 } 2424 PetscCallMPI(MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)s))); 2425 globalOff -= off; 2426 for (p = 0, off = 0; p < pEnd - pStart; ++p) { 2427 (*gsection)->atlasOff[p] += globalOff; 2428 if (neg) neg[p] = -((*gsection)->atlasOff[p] + 1); 2429 } 2430 /* Put in negative offsets for ghost points */ 2431 if (nroots >= 0) { 2432 PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE)); 2433 PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE)); 2434 if (nroots > pEnd - pStart) { 2435 for (p = pStart; p < pEnd; ++p) { 2436 if (tmpOff[p] < 0) (*gsection)->atlasOff[p - pStart] = tmpOff[p]; 2437 } 2438 } 2439 } 2440 if (nroots >= 0 && nroots > pEnd - pStart) PetscCall(PetscFree(tmpOff)); 2441 PetscCall(PetscFree(neg)); 2442 PetscFunctionReturn(PETSC_SUCCESS); 2443 } 2444 2445 typedef struct _n_PetscSectionSym_Label { 2446 DMLabel label; 2447 PetscCopyMode *modes; 2448 PetscInt *sizes; 2449 const PetscInt ***perms; 2450 const PetscScalar ***rots; 2451 PetscInt (*minMaxOrients)[2]; 2452 PetscInt numStrata; /* numStrata is only increasing, functions as a state */ 2453 } PetscSectionSym_Label; 2454 2455 static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym) 2456 { 2457 PetscInt i, j; 2458 PetscSectionSym_Label *sl = (PetscSectionSym_Label *)sym->data; 2459 2460 PetscFunctionBegin; 2461 for (i = 0; i <= sl->numStrata; i++) { 2462 if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) { 2463 for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 2464 if (sl->perms[i]) PetscCall(PetscFree(sl->perms[i][j])); 2465 if (sl->rots[i]) PetscCall(PetscFree(sl->rots[i][j])); 2466 } 2467 if (sl->perms[i]) { 2468 const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]]; 2469 2470 PetscCall(PetscFree(perms)); 2471 } 2472 if (sl->rots[i]) { 2473 const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]]; 2474 2475 PetscCall(PetscFree(rots)); 2476 } 2477 } 2478 } 2479 PetscCall(PetscFree5(sl->modes, sl->sizes, sl->perms, sl->rots, sl->minMaxOrients)); 2480 PetscCall(DMLabelDestroy(&sl->label)); 2481 sl->numStrata = 0; 2482 PetscFunctionReturn(PETSC_SUCCESS); 2483 } 2484 2485 static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym) 2486 { 2487 PetscFunctionBegin; 2488 PetscCall(PetscSectionSymLabelReset(sym)); 2489 PetscCall(PetscFree(sym->data)); 2490 PetscFunctionReturn(PETSC_SUCCESS); 2491 } 2492 2493 static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer) 2494 { 2495 PetscSectionSym_Label *sl = (PetscSectionSym_Label *)sym->data; 2496 PetscBool isAscii; 2497 DMLabel label = sl->label; 2498 const char *name; 2499 2500 PetscFunctionBegin; 2501 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isAscii)); 2502 if (isAscii) { 2503 PetscInt i, j, k; 2504 PetscViewerFormat format; 2505 2506 PetscCall(PetscViewerGetFormat(viewer, &format)); 2507 if (label) { 2508 PetscCall(PetscViewerGetFormat(viewer, &format)); 2509 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 2510 PetscCall(PetscViewerASCIIPushTab(viewer)); 2511 PetscCall(DMLabelView(label, viewer)); 2512 PetscCall(PetscViewerASCIIPopTab(viewer)); 2513 } else { 2514 PetscCall(PetscObjectGetName((PetscObject)sl->label, &name)); 2515 PetscCall(PetscViewerASCIIPrintf(viewer, " Label '%s'\n", name)); 2516 } 2517 } else { 2518 PetscCall(PetscViewerASCIIPrintf(viewer, "No label given\n")); 2519 } 2520 PetscCall(PetscViewerASCIIPushTab(viewer)); 2521 for (i = 0; i <= sl->numStrata; i++) { 2522 PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue; 2523 2524 if (!(sl->perms[i] || sl->rots[i])) { 2525 PetscCall(PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %" PetscInt_FMT " (%" PetscInt_FMT " dofs per point): no symmetries\n", value, sl->sizes[i])); 2526 } else { 2527 PetscCall(PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %" PetscInt_FMT " (%" PetscInt_FMT " dofs per point):\n", value, sl->sizes[i])); 2528 PetscCall(PetscViewerASCIIPushTab(viewer)); 2529 PetscCall(PetscViewerASCIIPrintf(viewer, "Orientation range: [%" PetscInt_FMT ", %" PetscInt_FMT ")\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1])); 2530 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 2531 PetscCall(PetscViewerASCIIPushTab(viewer)); 2532 for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 2533 if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) { 2534 PetscCall(PetscViewerASCIIPrintf(viewer, "Orientation %" PetscInt_FMT ": identity\n", j)); 2535 } else { 2536 PetscInt tab; 2537 2538 PetscCall(PetscViewerASCIIPrintf(viewer, "Orientation %" PetscInt_FMT ":\n", j)); 2539 PetscCall(PetscViewerASCIIPushTab(viewer)); 2540 PetscCall(PetscViewerASCIIGetTab(viewer, &tab)); 2541 if (sl->perms[i] && sl->perms[i][j]) { 2542 PetscCall(PetscViewerASCIIPrintf(viewer, "Permutation:")); 2543 PetscCall(PetscViewerASCIISetTab(viewer, 0)); 2544 for (k = 0; k < sl->sizes[i]; k++) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, sl->perms[i][j][k])); 2545 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 2546 PetscCall(PetscViewerASCIISetTab(viewer, tab)); 2547 } 2548 if (sl->rots[i] && sl->rots[i][j]) { 2549 PetscCall(PetscViewerASCIIPrintf(viewer, "Rotations: ")); 2550 PetscCall(PetscViewerASCIISetTab(viewer, 0)); 2551 #if defined(PETSC_USE_COMPLEX) 2552 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]))); 2553 #else 2554 for (k = 0; k < sl->sizes[i]; k++) PetscCall(PetscViewerASCIIPrintf(viewer, " %+g", (double)sl->rots[i][j][k])); 2555 #endif 2556 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 2557 PetscCall(PetscViewerASCIISetTab(viewer, tab)); 2558 } 2559 PetscCall(PetscViewerASCIIPopTab(viewer)); 2560 } 2561 } 2562 PetscCall(PetscViewerASCIIPopTab(viewer)); 2563 } 2564 PetscCall(PetscViewerASCIIPopTab(viewer)); 2565 } 2566 } 2567 PetscCall(PetscViewerASCIIPopTab(viewer)); 2568 } 2569 PetscFunctionReturn(PETSC_SUCCESS); 2570 } 2571 2572 /*@ 2573 PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries 2574 2575 Logically 2576 2577 Input Parameters: 2578 + sym - the section symmetries 2579 - label - the `DMLabel` describing the types of points 2580 2581 Level: developer: 2582 2583 .seealso: `DMLabel`, `DM`, `PetscSectionSymLabelSetStratum()`, `PetscSectionSymCreateLabel()`, `PetscSectionGetPointSyms()` 2584 @*/ 2585 PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label) 2586 { 2587 PetscSectionSym_Label *sl; 2588 2589 PetscFunctionBegin; 2590 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1); 2591 sl = (PetscSectionSym_Label *)sym->data; 2592 if (sl->label && sl->label != label) PetscCall(PetscSectionSymLabelReset(sym)); 2593 if (label) { 2594 sl->label = label; 2595 PetscCall(PetscObjectReference((PetscObject)label)); 2596 PetscCall(DMLabelGetNumValues(label, &sl->numStrata)); 2597 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)); 2598 PetscCall(PetscMemzero((void *)sl->modes, (sl->numStrata + 1) * sizeof(PetscCopyMode))); 2599 PetscCall(PetscMemzero((void *)sl->sizes, (sl->numStrata + 1) * sizeof(PetscInt))); 2600 PetscCall(PetscMemzero((void *)sl->perms, (sl->numStrata + 1) * sizeof(const PetscInt **))); 2601 PetscCall(PetscMemzero((void *)sl->rots, (sl->numStrata + 1) * sizeof(const PetscScalar **))); 2602 PetscCall(PetscMemzero((void *)sl->minMaxOrients, (sl->numStrata + 1) * sizeof(PetscInt[2]))); 2603 } 2604 PetscFunctionReturn(PETSC_SUCCESS); 2605 } 2606 2607 /*@C 2608 PetscSectionSymLabelGetStratum - get the symmetries for the orientations of a stratum 2609 2610 Logically Collective 2611 2612 Input Parameters: 2613 + sym - the section symmetries 2614 - stratum - the stratum value in the label that we are assigning symmetries for 2615 2616 Output Parameters: 2617 + size - the number of dofs for points in the `stratum` of the label 2618 . minOrient - the smallest orientation for a point in this `stratum` 2619 . maxOrient - one greater than the largest orientation for a ppoint in this `stratum` (i.e., orientations are in the range [`minOrient`, `maxOrient`)) 2620 . perms - `NULL` if there are no permutations, or (`maxOrient` - `minOrient`) permutations, one for each orientation. A `NULL` permutation is the identity 2621 - 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 2622 2623 Level: developer 2624 2625 .seealso: `DMLabel`, `DM`, `PetscSectionSymLabelSetStratum()`, `PetscSectionSymCreate()`, `PetscSectionSetSym()`, `PetscSectionGetPointSyms()`, `PetscSectionSymCreateLabel()` 2626 @*/ 2627 PetscErrorCode PetscSectionSymLabelGetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt *size, PetscInt *minOrient, PetscInt *maxOrient, const PetscInt ***perms, const PetscScalar ***rots) 2628 { 2629 PetscSectionSym_Label *sl; 2630 const char *name; 2631 PetscInt i; 2632 2633 PetscFunctionBegin; 2634 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1); 2635 sl = (PetscSectionSym_Label *)sym->data; 2636 PetscCheck(sl->label, PetscObjectComm((PetscObject)sym), PETSC_ERR_ARG_WRONGSTATE, "No label set yet"); 2637 for (i = 0; i <= sl->numStrata; i++) { 2638 PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue; 2639 2640 if (stratum == value) break; 2641 } 2642 PetscCall(PetscObjectGetName((PetscObject)sl->label, &name)); 2643 PetscCheck(i <= sl->numStrata, PetscObjectComm((PetscObject)sym), PETSC_ERR_ARG_OUTOFRANGE, "Stratum %" PetscInt_FMT " not found in label %s", stratum, name); 2644 if (size) { 2645 PetscAssertPointer(size, 3); 2646 *size = sl->sizes[i]; 2647 } 2648 if (minOrient) { 2649 PetscAssertPointer(minOrient, 4); 2650 *minOrient = sl->minMaxOrients[i][0]; 2651 } 2652 if (maxOrient) { 2653 PetscAssertPointer(maxOrient, 5); 2654 *maxOrient = sl->minMaxOrients[i][1]; 2655 } 2656 if (perms) { 2657 PetscAssertPointer(perms, 6); 2658 *perms = PetscSafePointerPlusOffset(sl->perms[i], sl->minMaxOrients[i][0]); 2659 } 2660 if (rots) { 2661 PetscAssertPointer(rots, 7); 2662 *rots = PetscSafePointerPlusOffset(sl->rots[i], sl->minMaxOrients[i][0]); 2663 } 2664 PetscFunctionReturn(PETSC_SUCCESS); 2665 } 2666 2667 /*@C 2668 PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum 2669 2670 Logically 2671 2672 Input Parameters: 2673 + sym - the section symmetries 2674 . stratum - the stratum value in the label that we are assigning symmetries for 2675 . size - the number of dofs for points in the `stratum` of the label 2676 . minOrient - the smallest orientation for a point in this `stratum` 2677 . maxOrient - one greater than the largest orientation for a point in this `stratum` (i.e., orientations are in the range [`minOrient`, `maxOrient`)) 2678 . mode - how `sym` should copy the `perms` and `rots` arrays 2679 . perms - `NULL` if there are no permutations, or (`maxOrient` - `minOrient`) permutations, one for each orientation. A `NULL` permutation is the identity 2680 - 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 2681 2682 Level: developer 2683 2684 .seealso: `DMLabel`, `DM`, `PetscSectionSymLabelGetStratum()`, `PetscSectionSymCreate()`, `PetscSectionSetSym()`, `PetscSectionGetPointSyms()`, `PetscSectionSymCreateLabel()` 2685 @*/ 2686 PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots) 2687 { 2688 PetscSectionSym_Label *sl; 2689 const char *name; 2690 PetscInt i, j, k; 2691 2692 PetscFunctionBegin; 2693 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1); 2694 sl = (PetscSectionSym_Label *)sym->data; 2695 PetscCheck(sl->label, PetscObjectComm((PetscObject)sym), PETSC_ERR_ARG_WRONGSTATE, "No label set yet"); 2696 for (i = 0; i <= sl->numStrata; i++) { 2697 PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue; 2698 2699 if (stratum == value) break; 2700 } 2701 PetscCall(PetscObjectGetName((PetscObject)sl->label, &name)); 2702 PetscCheck(i <= sl->numStrata, PetscObjectComm((PetscObject)sym), PETSC_ERR_ARG_OUTOFRANGE, "Stratum %" PetscInt_FMT " not found in label %s", stratum, name); 2703 sl->sizes[i] = size; 2704 sl->modes[i] = mode; 2705 sl->minMaxOrients[i][0] = minOrient; 2706 sl->minMaxOrients[i][1] = maxOrient; 2707 if (mode == PETSC_COPY_VALUES) { 2708 if (perms) { 2709 PetscInt **ownPerms; 2710 2711 PetscCall(PetscCalloc1(maxOrient - minOrient, &ownPerms)); 2712 for (j = 0; j < maxOrient - minOrient; j++) { 2713 if (perms[j]) { 2714 PetscCall(PetscMalloc1(size, &ownPerms[j])); 2715 for (k = 0; k < size; k++) ownPerms[j][k] = perms[j][k]; 2716 } 2717 } 2718 sl->perms[i] = (const PetscInt **)&ownPerms[-minOrient]; 2719 } 2720 if (rots) { 2721 PetscScalar **ownRots; 2722 2723 PetscCall(PetscCalloc1(maxOrient - minOrient, &ownRots)); 2724 for (j = 0; j < maxOrient - minOrient; j++) { 2725 if (rots[j]) { 2726 PetscCall(PetscMalloc1(size, &ownRots[j])); 2727 for (k = 0; k < size; k++) ownRots[j][k] = rots[j][k]; 2728 } 2729 } 2730 sl->rots[i] = (const PetscScalar **)&ownRots[-minOrient]; 2731 } 2732 } else { 2733 sl->perms[i] = PetscSafePointerPlusOffset(perms, -minOrient); 2734 sl->rots[i] = PetscSafePointerPlusOffset(rots, -minOrient); 2735 } 2736 PetscFunctionReturn(PETSC_SUCCESS); 2737 } 2738 2739 static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots) 2740 { 2741 PetscInt i, j, numStrata; 2742 PetscSectionSym_Label *sl; 2743 DMLabel label; 2744 2745 PetscFunctionBegin; 2746 sl = (PetscSectionSym_Label *)sym->data; 2747 numStrata = sl->numStrata; 2748 label = sl->label; 2749 for (i = 0; i < numPoints; i++) { 2750 PetscInt point = points[2 * i]; 2751 PetscInt ornt = points[2 * i + 1]; 2752 2753 for (j = 0; j < numStrata; j++) { 2754 if (label->validIS[j]) { 2755 PetscInt k; 2756 2757 PetscCall(ISLocate(label->points[j], point, &k)); 2758 if (k >= 0) break; 2759 } else { 2760 PetscBool has; 2761 2762 PetscCall(PetscHSetIHas(label->ht[j], point, &has)); 2763 if (has) break; 2764 } 2765 } 2766 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], 2767 j < numStrata ? label->stratumValues[j] : label->defaultValue); 2768 if (perms) perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL; 2769 if (rots) rots[i] = sl->rots[j] ? sl->rots[j][ornt] : NULL; 2770 } 2771 PetscFunctionReturn(PETSC_SUCCESS); 2772 } 2773 2774 static PetscErrorCode PetscSectionSymCopy_Label(PetscSectionSym sym, PetscSectionSym nsym) 2775 { 2776 PetscSectionSym_Label *sl = (PetscSectionSym_Label *)nsym->data; 2777 IS valIS; 2778 const PetscInt *values; 2779 PetscInt Nv, v; 2780 2781 PetscFunctionBegin; 2782 PetscCall(DMLabelGetNumValues(sl->label, &Nv)); 2783 PetscCall(DMLabelGetValueIS(sl->label, &valIS)); 2784 PetscCall(ISGetIndices(valIS, &values)); 2785 for (v = 0; v < Nv; ++v) { 2786 const PetscInt val = values[v]; 2787 PetscInt size, minOrient, maxOrient; 2788 const PetscInt **perms; 2789 const PetscScalar **rots; 2790 2791 PetscCall(PetscSectionSymLabelGetStratum(sym, val, &size, &minOrient, &maxOrient, &perms, &rots)); 2792 PetscCall(PetscSectionSymLabelSetStratum(nsym, val, size, minOrient, maxOrient, PETSC_COPY_VALUES, perms, rots)); 2793 } 2794 PetscCall(ISDestroy(&valIS)); 2795 PetscFunctionReturn(PETSC_SUCCESS); 2796 } 2797 2798 static PetscErrorCode PetscSectionSymDistribute_Label(PetscSectionSym sym, PetscSF migrationSF, PetscSectionSym *dsym) 2799 { 2800 PetscSectionSym_Label *sl = (PetscSectionSym_Label *)sym->data; 2801 DMLabel dlabel; 2802 2803 PetscFunctionBegin; 2804 PetscCall(DMLabelDistribute(sl->label, migrationSF, &dlabel)); 2805 PetscCall(PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)sym), dlabel, dsym)); 2806 PetscCall(DMLabelDestroy(&dlabel)); 2807 PetscCall(PetscSectionSymCopy(sym, *dsym)); 2808 PetscFunctionReturn(PETSC_SUCCESS); 2809 } 2810 2811 PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym) 2812 { 2813 PetscSectionSym_Label *sl; 2814 2815 PetscFunctionBegin; 2816 PetscCall(PetscNew(&sl)); 2817 sym->ops->getpoints = PetscSectionSymGetPoints_Label; 2818 sym->ops->distribute = PetscSectionSymDistribute_Label; 2819 sym->ops->copy = PetscSectionSymCopy_Label; 2820 sym->ops->view = PetscSectionSymView_Label; 2821 sym->ops->destroy = PetscSectionSymDestroy_Label; 2822 sym->data = (void *)sl; 2823 PetscFunctionReturn(PETSC_SUCCESS); 2824 } 2825 2826 /*@ 2827 PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label 2828 2829 Collective 2830 2831 Input Parameters: 2832 + comm - the MPI communicator for the new symmetry 2833 - label - the label defining the strata 2834 2835 Output Parameter: 2836 . sym - the section symmetries 2837 2838 Level: developer 2839 2840 .seealso: `DMLabel`, `DM`, `PetscSectionSymCreate()`, `PetscSectionSetSym()`, `PetscSectionGetSym()`, `PetscSectionSymLabelSetStratum()`, `PetscSectionGetPointSyms()` 2841 @*/ 2842 PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym) 2843 { 2844 PetscFunctionBegin; 2845 PetscCall(DMInitializePackage()); 2846 PetscCall(PetscSectionSymCreate(comm, sym)); 2847 PetscCall(PetscSectionSymSetType(*sym, PETSCSECTIONSYMLABEL)); 2848 PetscCall(PetscSectionSymLabelSetLabel(*sym, label)); 2849 PetscFunctionReturn(PETSC_SUCCESS); 2850 } 2851