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