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