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