1 #include <petscdm.h> 2 #include <petsc/private/dmlabelimpl.h> /*I "petscdmlabel.h" I*/ 3 #include <petsc/private/isimpl.h> /*I "petscis.h" I*/ 4 #include <petscsf.h> 5 6 /*@C 7 DMLabelCreate - Create a DMLabel object, which is a multimap 8 9 Input parameter: 10 . name - The label name 11 12 Output parameter: 13 . label - The DMLabel 14 15 Level: beginner 16 17 .seealso: DMLabelDestroy() 18 @*/ 19 PetscErrorCode DMLabelCreate(const char name[], DMLabel *label) 20 { 21 PetscErrorCode ierr; 22 23 PetscFunctionBegin; 24 ierr = PetscNew(label);CHKERRQ(ierr); 25 ierr = PetscStrallocpy(name, &(*label)->name);CHKERRQ(ierr); 26 27 (*label)->refct = 1; 28 (*label)->state = -1; 29 (*label)->numStrata = 0; 30 (*label)->defaultValue = -1; 31 (*label)->stratumValues = NULL; 32 (*label)->validIS = NULL; 33 (*label)->stratumSizes = NULL; 34 (*label)->points = NULL; 35 (*label)->ht = NULL; 36 (*label)->pStart = -1; 37 (*label)->pEnd = -1; 38 (*label)->bt = NULL; 39 PetscFunctionReturn(0); 40 } 41 42 /* 43 DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format 44 45 Input parameter: 46 + label - The DMLabel 47 - v - The stratum value 48 49 Output parameter: 50 . label - The DMLabel with stratum in sorted list format 51 52 Level: developer 53 54 .seealso: DMLabelCreate() 55 */ 56 static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v) 57 { 58 PetscInt off = 0; 59 PetscInt *pointArray; 60 PetscErrorCode ierr; 61 62 if (label->validIS[v]) return 0; 63 PetscFunctionBegin; 64 if (v < 0 || v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeValid_Private\n", v); 65 PetscHashISize(label->ht[v], label->stratumSizes[v]); 66 67 ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr); 68 ierr = PetscHashIGetKeys(label->ht[v], &off, pointArray);CHKERRQ(ierr); 69 if (off != label->stratumSizes[v]) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of contributed points %D from value %D should be %D", off, label->stratumValues[v], label->stratumSizes[v]); 70 PetscHashIClear(label->ht[v]); 71 ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr); 72 if (label->bt) { 73 PetscInt p; 74 75 for (p = 0; p < label->stratumSizes[v]; ++p) { 76 const PetscInt point = pointArray[p]; 77 78 if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, label->pStart, label->pEnd); 79 ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 80 } 81 } 82 ierr = ISCreateGeneral(PETSC_COMM_SELF,label->stratumSizes[v],pointArray,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr); 83 ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr); 84 label->validIS[v] = PETSC_TRUE; 85 ++label->state; 86 PetscFunctionReturn(0); 87 } 88 89 /* 90 DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format 91 92 Input parameter: 93 . label - The DMLabel 94 95 Output parameter: 96 . label - The DMLabel with all strata in sorted list format 97 98 Level: developer 99 100 .seealso: DMLabelCreate() 101 */ 102 static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label) 103 { 104 PetscInt v; 105 PetscErrorCode ierr; 106 107 PetscFunctionBegin; 108 for (v = 0; v < label->numStrata; v++){ 109 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 110 } 111 PetscFunctionReturn(0); 112 } 113 114 /* 115 DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format 116 117 Input parameter: 118 + label - The DMLabel 119 - v - The stratum value 120 121 Output parameter: 122 . label - The DMLabel with stratum in hash format 123 124 Level: developer 125 126 .seealso: DMLabelCreate() 127 */ 128 static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v) 129 { 130 PETSC_UNUSED PetscHashIIter ret, iter; 131 PetscInt p; 132 const PetscInt *points; 133 PetscErrorCode ierr; 134 135 if (!label->validIS[v]) return 0; 136 PetscFunctionBegin; 137 if (v < 0 || v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeInvalid_Private\n", v); 138 if (label->points[v]) { 139 ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 140 for (p = 0; p < label->stratumSizes[v]; ++p) PetscHashIPut(label->ht[v], points[p], ret, iter); 141 ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 142 ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 143 } 144 label->validIS[v] = PETSC_FALSE; 145 PetscFunctionReturn(0); 146 } 147 148 PetscErrorCode DMLabelGetState(DMLabel label, PetscObjectState *state) 149 { 150 PetscFunctionBegin; 151 PetscValidPointer(state, 2); 152 *state = label->state; 153 PetscFunctionReturn(0); 154 } 155 156 PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index) 157 { 158 PetscInt v; 159 PetscFunctionBegin; 160 for (*index = -1, v = 0; v < label->numStrata; ++v) 161 if (label->stratumValues[v] == value) { *index = v; break; } 162 PetscFunctionReturn(0); 163 } 164 165 static PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index) 166 { 167 PetscInt v, *tmpV, *tmpS; 168 IS *tmpP; 169 PetscHashI *tmpH; 170 PetscBool *tmpB; 171 PetscErrorCode ierr; 172 173 PetscFunctionBegin; 174 ierr = PetscMalloc1((label->numStrata+1), &tmpV);CHKERRQ(ierr); 175 ierr = PetscMalloc1((label->numStrata+1), &tmpS);CHKERRQ(ierr); 176 ierr = PetscMalloc1((label->numStrata+1), &tmpH);CHKERRQ(ierr); 177 ierr = PetscMalloc1((label->numStrata+1), &tmpP);CHKERRQ(ierr); 178 ierr = PetscMalloc1((label->numStrata+1), &tmpB);CHKERRQ(ierr); 179 for (v = 0; v < label->numStrata; ++v) { 180 tmpV[v] = label->stratumValues[v]; 181 tmpS[v] = label->stratumSizes[v]; 182 tmpH[v] = label->ht[v]; 183 tmpP[v] = label->points[v]; 184 tmpB[v] = label->validIS[v]; 185 } 186 tmpV[v] = value; 187 tmpS[v] = 0; 188 PetscHashICreate(tmpH[v]); 189 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&tmpP[v]);CHKERRQ(ierr); 190 tmpB[v] = PETSC_TRUE; 191 192 ++label->numStrata; 193 ierr = PetscFree(label->stratumValues);CHKERRQ(ierr); 194 ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr); 195 ierr = PetscFree(label->ht);CHKERRQ(ierr); 196 ierr = PetscFree(label->points);CHKERRQ(ierr); 197 ierr = PetscFree(label->validIS);CHKERRQ(ierr); 198 label->stratumValues = tmpV; 199 label->stratumSizes = tmpS; 200 label->ht = tmpH; 201 label->points = tmpP; 202 label->validIS = tmpB; 203 204 *index = v; 205 PetscFunctionReturn(0); 206 } 207 208 PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value) 209 { 210 PetscInt v; 211 PetscErrorCode ierr; 212 213 PetscFunctionBegin; 214 ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 215 if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);} 216 PetscFunctionReturn(0); 217 } 218 219 /*@C 220 DMLabelGetName - Return the name of a DMLabel object 221 222 Input parameter: 223 . label - The DMLabel 224 225 Output parameter: 226 . name - The label name 227 228 Level: beginner 229 230 .seealso: DMLabelCreate() 231 @*/ 232 PetscErrorCode DMLabelGetName(DMLabel label, const char **name) 233 { 234 PetscFunctionBegin; 235 PetscValidPointer(name, 2); 236 *name = label->name; 237 PetscFunctionReturn(0); 238 } 239 240 static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer) 241 { 242 PetscInt v; 243 PetscMPIInt rank; 244 PetscErrorCode ierr; 245 246 PetscFunctionBegin; 247 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr); 248 ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 249 if (label) { 250 ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", label->name);CHKERRQ(ierr); 251 if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, " Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);} 252 for (v = 0; v < label->numStrata; ++v) { 253 const PetscInt value = label->stratumValues[v]; 254 const PetscInt *points; 255 PetscInt p; 256 257 ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 258 for (p = 0; p < label->stratumSizes[v]; ++p) { 259 ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr); 260 } 261 ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 262 } 263 } 264 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 265 ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 266 PetscFunctionReturn(0); 267 } 268 269 /*@C 270 DMLabelView - View the label 271 272 Input Parameters: 273 + label - The DMLabel 274 - viewer - The PetscViewer 275 276 Level: intermediate 277 278 .seealso: DMLabelCreate(), DMLabelDestroy() 279 @*/ 280 PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer) 281 { 282 PetscBool iascii; 283 PetscErrorCode ierr; 284 285 PetscFunctionBegin; 286 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 287 if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 288 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 289 if (iascii) { 290 ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr); 291 } 292 PetscFunctionReturn(0); 293 } 294 295 /*@ 296 DMLabelDestroy - Destroys a DMLabel 297 298 Input Parameter: 299 . label - The DMLabel 300 301 Level: beginner 302 303 .seealso: DMLabelCreate() 304 @*/ 305 PetscErrorCode DMLabelDestroy(DMLabel *label) 306 { 307 PetscInt v; 308 PetscErrorCode ierr; 309 310 PetscFunctionBegin; 311 if (!(*label)) PetscFunctionReturn(0); 312 if (--(*label)->refct > 0) PetscFunctionReturn(0); 313 ierr = PetscFree((*label)->name);CHKERRQ(ierr); 314 ierr = PetscFree((*label)->stratumValues);CHKERRQ(ierr); 315 ierr = PetscFree((*label)->stratumSizes);CHKERRQ(ierr); 316 for (v = 0; v < (*label)->numStrata; ++v) { 317 PetscHashIDestroy((*label)->ht[v]); 318 ierr = ISDestroy(&((*label)->points[v]));CHKERRQ(ierr); 319 } 320 ierr = PetscFree((*label)->ht);CHKERRQ(ierr); 321 ierr = PetscFree((*label)->points);CHKERRQ(ierr); 322 ierr = PetscFree((*label)->validIS);CHKERRQ(ierr); 323 ierr = PetscBTDestroy(&(*label)->bt);CHKERRQ(ierr); 324 ierr = PetscFree(*label);CHKERRQ(ierr); 325 PetscFunctionReturn(0); 326 } 327 328 /*@ 329 DMLabelDuplicate - Duplicates a DMLabel 330 331 Input Parameter: 332 . label - The DMLabel 333 334 Output Parameter: 335 . labelnew - location to put new vector 336 337 Level: intermediate 338 339 .seealso: DMLabelCreate(), DMLabelDestroy() 340 @*/ 341 PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew) 342 { 343 PetscInt v; 344 PetscErrorCode ierr; 345 346 PetscFunctionBegin; 347 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 348 ierr = PetscNew(labelnew);CHKERRQ(ierr); 349 ierr = PetscStrallocpy(label->name, &(*labelnew)->name);CHKERRQ(ierr); 350 351 (*labelnew)->refct = 1; 352 (*labelnew)->numStrata = label->numStrata; 353 (*labelnew)->defaultValue = label->defaultValue; 354 if (label->numStrata) { 355 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr); 356 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr); 357 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr); 358 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr); 359 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr); 360 /* Could eliminate unused space here */ 361 for (v = 0; v < label->numStrata; ++v) { 362 PetscHashICreate((*labelnew)->ht[v]); 363 (*labelnew)->validIS[v] = PETSC_TRUE; 364 (*labelnew)->stratumValues[v] = label->stratumValues[v]; 365 (*labelnew)->stratumSizes[v] = label->stratumSizes[v]; 366 ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr); 367 (*labelnew)->points[v] = label->points[v]; 368 } 369 } 370 (*labelnew)->pStart = -1; 371 (*labelnew)->pEnd = -1; 372 (*labelnew)->bt = NULL; 373 PetscFunctionReturn(0); 374 } 375 376 /* This can be hooked into SetValue(), ClearValue(), etc. for updating */ 377 PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd) 378 { 379 PetscInt v; 380 PetscErrorCode ierr; 381 382 PetscFunctionBegin; 383 ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr); 384 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 385 label->pStart = pStart; 386 label->pEnd = pEnd; 387 ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr); 388 for (v = 0; v < label->numStrata; ++v) { 389 const PetscInt *points; 390 PetscInt i; 391 392 ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 393 for (i = 0; i < label->stratumSizes[v]; ++i) { 394 const PetscInt point = points[i]; 395 396 if ((point < pStart) || (point >= pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, pStart, pEnd); 397 ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr); 398 } 399 ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr); 400 } 401 PetscFunctionReturn(0); 402 } 403 404 PetscErrorCode DMLabelDestroyIndex(DMLabel label) 405 { 406 PetscErrorCode ierr; 407 408 PetscFunctionBegin; 409 label->pStart = -1; 410 label->pEnd = -1; 411 ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 412 PetscFunctionReturn(0); 413 } 414 415 /*@ 416 DMLabelHasValue - Determine whether a label assigns the value to any point 417 418 Input Parameters: 419 + label - the DMLabel 420 - value - the value 421 422 Output Parameter: 423 . contains - Flag indicating whether the label maps this value to any point 424 425 Level: developer 426 427 .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue() 428 @*/ 429 PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains) 430 { 431 PetscInt v; 432 PetscErrorCode ierr; 433 434 PetscFunctionBegin; 435 PetscValidPointer(contains, 3); 436 ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 437 *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE; 438 PetscFunctionReturn(0); 439 } 440 441 /*@ 442 DMLabelHasPoint - Determine whether a label assigns a value to a point 443 444 Input Parameters: 445 + label - the DMLabel 446 - point - the point 447 448 Output Parameter: 449 . contains - Flag indicating whether the label maps this point to a value 450 451 Note: The user must call DMLabelCreateIndex() before this function. 452 453 Level: developer 454 455 .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue() 456 @*/ 457 PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains) 458 { 459 PetscErrorCode ierr; 460 461 PetscFunctionBeginHot; 462 PetscValidPointer(contains, 3); 463 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 464 #if defined(PETSC_USE_DEBUG) 465 if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()"); 466 if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, label->pStart, label->pEnd); 467 #endif 468 *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE; 469 PetscFunctionReturn(0); 470 } 471 472 /*@ 473 DMLabelStratumHasPoint - Return true if the stratum contains a point 474 475 Input Parameters: 476 + label - the DMLabel 477 . value - the stratum value 478 - point - the point 479 480 Output Parameter: 481 . contains - true if the stratum contains the point 482 483 Level: intermediate 484 485 .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue() 486 @*/ 487 PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains) 488 { 489 PetscInt v; 490 PetscErrorCode ierr; 491 492 PetscFunctionBeginHot; 493 PetscValidPointer(contains, 4); 494 *contains = PETSC_FALSE; 495 ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 496 if (v < 0) PetscFunctionReturn(0); 497 498 if (label->validIS[v]) { 499 PetscInt i; 500 501 ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr); 502 if (i >= 0) *contains = PETSC_TRUE; 503 } else { 504 PetscBool has; 505 506 PetscHashIHasKey(label->ht[v], point, has); 507 if (has) *contains = PETSC_TRUE; 508 } 509 PetscFunctionReturn(0); 510 } 511 512 /*@ 513 DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 514 When a label is created, it is initialized to -1. 515 516 Input parameter: 517 . label - a DMLabel object 518 519 Output parameter: 520 . defaultValue - the default value 521 522 Level: beginner 523 524 .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 525 @*/ 526 PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue) 527 { 528 PetscFunctionBegin; 529 *defaultValue = label->defaultValue; 530 PetscFunctionReturn(0); 531 } 532 533 /*@ 534 DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 535 When a label is created, it is initialized to -1. 536 537 Input parameter: 538 . label - a DMLabel object 539 540 Output parameter: 541 . defaultValue - the default value 542 543 Level: beginner 544 545 .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 546 @*/ 547 PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue) 548 { 549 PetscFunctionBegin; 550 label->defaultValue = defaultValue; 551 PetscFunctionReturn(0); 552 } 553 554 /*@ 555 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()) 556 557 Input Parameters: 558 + label - the DMLabel 559 - point - the point 560 561 Output Parameter: 562 . value - The point value, or -1 563 564 Level: intermediate 565 566 .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 567 @*/ 568 PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value) 569 { 570 PetscInt v; 571 PetscErrorCode ierr; 572 573 PetscFunctionBeginHot; 574 PetscValidPointer(value, 3); 575 *value = label->defaultValue; 576 for (v = 0; v < label->numStrata; ++v) { 577 if (label->validIS[v]) { 578 PetscInt i; 579 580 ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr); 581 if (i >= 0) { 582 *value = label->stratumValues[v]; 583 break; 584 } 585 } else { 586 PetscBool has; 587 588 PetscHashIHasKey(label->ht[v], point, has); 589 if (has) { 590 *value = label->stratumValues[v]; 591 break; 592 } 593 } 594 } 595 PetscFunctionReturn(0); 596 } 597 598 /*@ 599 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 somethingg different), then this function will do nothing. 600 601 Input Parameters: 602 + label - the DMLabel 603 . point - the point 604 - value - The point value 605 606 Level: intermediate 607 608 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 609 @*/ 610 PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value) 611 { 612 PETSC_UNUSED PetscHashIIter iter, ret; 613 PetscInt v; 614 PetscErrorCode ierr; 615 616 PetscFunctionBegin; 617 /* Find label value, add new entry if needed */ 618 if (value == label->defaultValue) PetscFunctionReturn(0); 619 ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 620 if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);} 621 /* Set key */ 622 ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 623 PetscHashIPut(label->ht[v], point, ret, iter); 624 PetscFunctionReturn(0); 625 } 626 627 /*@ 628 DMLabelClearValue - Clear the value a label assigns to a point 629 630 Input Parameters: 631 + label - the DMLabel 632 . point - the point 633 - value - The point value 634 635 Level: intermediate 636 637 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue() 638 @*/ 639 PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value) 640 { 641 PetscInt v; 642 PetscErrorCode ierr; 643 644 PetscFunctionBegin; 645 /* Find label value */ 646 ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 647 if (v < 0) PetscFunctionReturn(0); 648 649 if (label->bt) { 650 if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, label->pStart, label->pEnd); 651 ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 652 } 653 654 /* Delete key */ 655 ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 656 ierr = PetscHashIDelKey(label->ht[v], point);CHKERRQ(ierr); 657 PetscFunctionReturn(0); 658 } 659 660 /*@ 661 DMLabelInsertIS - Set all points in the IS to a value 662 663 Input Parameters: 664 + label - the DMLabel 665 . is - the point IS 666 - value - The point value 667 668 Level: intermediate 669 670 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 671 @*/ 672 PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value) 673 { 674 PETSC_UNUSED PetscHashIIter iter, ret; 675 PetscInt v, n, p; 676 const PetscInt *points; 677 PetscErrorCode ierr; 678 679 PetscFunctionBegin; 680 PetscValidHeaderSpecific(is, IS_CLASSID, 2); 681 /* Find label value, add new entry if needed */ 682 if (value == label->defaultValue) PetscFunctionReturn(0); 683 ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 684 if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);} 685 /* Set keys */ 686 ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 687 ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 688 ierr = ISGetIndices(is, &points);CHKERRQ(ierr); 689 for (p = 0; p < n; ++p) PetscHashIPut(label->ht[v], points[p], ret, iter); 690 ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr); 691 PetscFunctionReturn(0); 692 } 693 694 /*@ 695 DMLabelGetNumValues - Get the number of values that the DMLabel takes 696 697 Input Parameter: 698 . label - the DMLabel 699 700 Output Paramater: 701 . numValues - the number of values 702 703 Level: intermediate 704 705 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 706 @*/ 707 PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues) 708 { 709 PetscFunctionBegin; 710 PetscValidPointer(numValues, 2); 711 *numValues = label->numStrata; 712 PetscFunctionReturn(0); 713 } 714 715 /*@ 716 DMLabelGetValueIS - Get an IS of all values that the DMlabel takes 717 718 Input Parameter: 719 . label - the DMLabel 720 721 Output Paramater: 722 . is - the value IS 723 724 Level: intermediate 725 726 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 727 @*/ 728 PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values) 729 { 730 PetscErrorCode ierr; 731 732 PetscFunctionBegin; 733 PetscValidPointer(values, 2); 734 ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr); 735 PetscFunctionReturn(0); 736 } 737 738 /*@ 739 DMLabelHasStratum - Determine whether points exist with the given value 740 741 Input Parameters: 742 + label - the DMLabel 743 - value - the stratum value 744 745 Output Paramater: 746 . exists - Flag saying whether points exist 747 748 Level: intermediate 749 750 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 751 @*/ 752 PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists) 753 { 754 PetscInt v; 755 PetscErrorCode ierr; 756 757 PetscFunctionBegin; 758 PetscValidPointer(exists, 3); 759 ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 760 *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE; 761 PetscFunctionReturn(0); 762 } 763 764 /*@ 765 DMLabelGetStratumSize - Get the size of a stratum 766 767 Input Parameters: 768 + label - the DMLabel 769 - value - the stratum value 770 771 Output Paramater: 772 . size - The number of points in the stratum 773 774 Level: intermediate 775 776 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 777 @*/ 778 PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size) 779 { 780 PetscInt v; 781 PetscErrorCode ierr; 782 783 PetscFunctionBegin; 784 PetscValidPointer(size, 3); 785 *size = 0; 786 ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 787 if (v < 0) PetscFunctionReturn(0); 788 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 789 *size = label->stratumSizes[v]; 790 PetscFunctionReturn(0); 791 } 792 793 /*@ 794 DMLabelGetStratumBounds - Get the largest and smallest point of a stratum 795 796 Input Parameters: 797 + label - the DMLabel 798 - value - the stratum value 799 800 Output Paramaters: 801 + start - the smallest point in the stratum 802 - end - the largest point in the stratum 803 804 Level: intermediate 805 806 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 807 @*/ 808 PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end) 809 { 810 PetscInt v, min, max; 811 PetscErrorCode ierr; 812 813 PetscFunctionBegin; 814 if (start) {PetscValidPointer(start, 3); *start = 0;} 815 if (end) {PetscValidPointer(end, 4); *end = 0;} 816 ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 817 if (v < 0) PetscFunctionReturn(0); 818 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 819 if (label->stratumSizes[v] <= 0) PetscFunctionReturn(0); 820 ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr); 821 if (start) *start = min; 822 if (end) *end = max+1; 823 PetscFunctionReturn(0); 824 } 825 826 /*@ 827 DMLabelGetStratumIS - Get an IS with the stratum points 828 829 Input Parameters: 830 + label - the DMLabel 831 - value - the stratum value 832 833 Output Paramater: 834 . points - The stratum points 835 836 Level: intermediate 837 838 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 839 @*/ 840 PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points) 841 { 842 PetscInt v; 843 PetscErrorCode ierr; 844 845 PetscFunctionBegin; 846 PetscValidPointer(points, 3); 847 *points = NULL; 848 ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 849 if (v < 0) PetscFunctionReturn(0); 850 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 851 ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr); 852 *points = label->points[v]; 853 PetscFunctionReturn(0); 854 } 855 856 /*@ 857 DMLabelSetStratumIS - Set the stratum points using an IS 858 859 Input Parameters: 860 + label - the DMLabel 861 . value - the stratum value 862 - points - The stratum points 863 864 Level: intermediate 865 866 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 867 @*/ 868 PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is) 869 { 870 PetscInt v; 871 PetscErrorCode ierr; 872 873 PetscFunctionBegin; 874 ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 875 if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);} 876 if (is == label->points[v]) PetscFunctionReturn(0); 877 ierr = DMLabelClearStratum(label, value);CHKERRQ(ierr); 878 ierr = ISGetLocalSize(is, &(label->stratumSizes[v]));CHKERRQ(ierr); 879 ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 880 ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 881 label->points[v] = is; 882 label->validIS[v] = PETSC_TRUE; 883 if (label->bt) { 884 const PetscInt *points; 885 PetscInt p; 886 887 ierr = ISGetIndices(is,&points);CHKERRQ(ierr); 888 for (p = 0; p < label->stratumSizes[v]; ++p) { 889 const PetscInt point = points[p]; 890 891 if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, label->pStart, label->pEnd); 892 ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 893 } 894 } 895 PetscFunctionReturn(0); 896 } 897 898 /*@ 899 DMLabelClearStratum - Remove a stratum 900 901 Input Parameters: 902 + label - the DMLabel 903 - value - the stratum value 904 905 Level: intermediate 906 907 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 908 @*/ 909 PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value) 910 { 911 PetscInt v; 912 PetscErrorCode ierr; 913 914 PetscFunctionBegin; 915 ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 916 if (v < 0) PetscFunctionReturn(0); 917 if (label->validIS[v]) { 918 if (label->bt) { 919 PetscInt i; 920 const PetscInt *points; 921 922 ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 923 for (i = 0; i < label->stratumSizes[v]; ++i) { 924 const PetscInt point = points[i]; 925 926 if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, label->pStart, label->pEnd); 927 ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 928 } 929 ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr); 930 } 931 label->stratumSizes[v] = 0; 932 ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr); 933 ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr); 934 ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr); 935 } else { 936 PetscHashIClear(label->ht[v]); 937 } 938 PetscFunctionReturn(0); 939 } 940 941 /*@ 942 DMLabelFilter - Remove all points outside of [start, end) 943 944 Input Parameters: 945 + label - the DMLabel 946 . start - the first point 947 - end - the last point 948 949 Level: intermediate 950 951 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 952 @*/ 953 PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end) 954 { 955 PetscInt v; 956 PetscErrorCode ierr; 957 958 PetscFunctionBegin; 959 ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr); 960 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 961 for (v = 0; v < label->numStrata; ++v) { 962 PetscInt off, q; 963 const PetscInt *points; 964 PetscInt numPointsNew = 0, *pointsNew = NULL; 965 966 ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 967 for (q = 0; q < label->stratumSizes[v]; ++q) 968 if (points[q] >= start && points[q] < end) 969 numPointsNew++; 970 ierr = PetscMalloc1(numPointsNew, &pointsNew);CHKERRQ(ierr); 971 for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) { 972 if (points[q] >= start && points[q] < end) 973 pointsNew[off++] = points[q]; 974 } 975 ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 976 977 label->stratumSizes[v] = numPointsNew; 978 ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr); 979 ierr = ISCreateGeneral(PETSC_COMM_SELF,numPointsNew, pointsNew, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr); 980 ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr); 981 } 982 ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr); 983 PetscFunctionReturn(0); 984 } 985 986 /*@ 987 DMLabelPermute - Create a new label with permuted points 988 989 Input Parameters: 990 + label - the DMLabel 991 - permutation - the point permutation 992 993 Output Parameter: 994 . labelnew - the new label containing the permuted points 995 996 Level: intermediate 997 998 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 999 @*/ 1000 PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew) 1001 { 1002 const PetscInt *perm; 1003 PetscInt numValues, numPoints, v, q; 1004 PetscErrorCode ierr; 1005 1006 PetscFunctionBegin; 1007 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 1008 ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr); 1009 ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr); 1010 ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr); 1011 ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr); 1012 for (v = 0; v < numValues; ++v) { 1013 const PetscInt size = (*labelNew)->stratumSizes[v]; 1014 const PetscInt *points; 1015 PetscInt *pointsNew; 1016 1017 ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 1018 ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr); 1019 for (q = 0; q < size; ++q) { 1020 const PetscInt point = points[q]; 1021 1022 if ((point < 0) || (point >= numPoints)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [0, %D) for the remapping", point, numPoints); 1023 pointsNew[q] = perm[point]; 1024 } 1025 ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 1026 ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr); 1027 ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr); 1028 if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) { 1029 ierr = ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));CHKERRQ(ierr); 1030 ierr = PetscFree(pointsNew);CHKERRQ(ierr); 1031 } else { 1032 ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr); 1033 } 1034 ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr); 1035 } 1036 ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr); 1037 if (label->bt) { 1038 ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 1039 ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr); 1040 } 1041 PetscFunctionReturn(0); 1042 } 1043 1044 PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata) 1045 { 1046 MPI_Comm comm; 1047 PetscInt s, l, nroots, nleaves, dof, offset, size; 1048 PetscInt *remoteOffsets, *rootStrata, *rootIdx; 1049 PetscSection rootSection; 1050 PetscSF labelSF; 1051 PetscErrorCode ierr; 1052 1053 PetscFunctionBegin; 1054 if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 1055 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 1056 /* Build a section of stratum values per point, generate the according SF 1057 and distribute point-wise stratum values to leaves. */ 1058 ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 1059 ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 1060 ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr); 1061 if (label) { 1062 for (s = 0; s < label->numStrata; ++s) { 1063 const PetscInt *points; 1064 1065 ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 1066 for (l = 0; l < label->stratumSizes[s]; l++) { 1067 ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr); 1068 ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr); 1069 } 1070 ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 1071 } 1072 } 1073 ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 1074 /* Create a point-wise array of stratum values */ 1075 ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr); 1076 ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr); 1077 ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr); 1078 if (label) { 1079 for (s = 0; s < label->numStrata; ++s) { 1080 const PetscInt *points; 1081 1082 ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 1083 for (l = 0; l < label->stratumSizes[s]; l++) { 1084 const PetscInt p = points[l]; 1085 ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr); 1086 rootStrata[offset+rootIdx[p]++] = label->stratumValues[s]; 1087 } 1088 ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 1089 } 1090 } 1091 /* Build SF that maps label points to remote processes */ 1092 ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr); 1093 ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr); 1094 ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr); 1095 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 1096 /* Send the strata for each point over the derived SF */ 1097 ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr); 1098 ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr); 1099 ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 1100 ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 1101 /* Clean up */ 1102 ierr = PetscFree(rootStrata);CHKERRQ(ierr); 1103 ierr = PetscFree(rootIdx);CHKERRQ(ierr); 1104 ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 1105 ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr); 1106 PetscFunctionReturn(0); 1107 } 1108 1109 /*@ 1110 DMLabelDistribute - Create a new label pushed forward over the PetscSF 1111 1112 Input Parameters: 1113 + label - the DMLabel 1114 - sf - the map from old to new distribution 1115 1116 Output Parameter: 1117 . labelnew - the new redistributed label 1118 1119 Level: intermediate 1120 1121 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 1122 @*/ 1123 PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew) 1124 { 1125 MPI_Comm comm; 1126 PetscSection leafSection; 1127 PetscInt p, pStart, pEnd, s, size, dof, offset, stratum; 1128 PetscInt *leafStrata, *strataIdx; 1129 PetscInt **points; 1130 char *name; 1131 PetscInt nameSize; 1132 PetscHashI stratumHash; 1133 PETSC_UNUSED PetscHashIIter ret, iter; 1134 size_t len = 0; 1135 PetscMPIInt rank; 1136 PetscErrorCode ierr; 1137 1138 PetscFunctionBegin; 1139 if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 1140 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 1141 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1142 /* Bcast name */ 1143 if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 1144 nameSize = len; 1145 ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1146 ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1147 if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 1148 ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1149 ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 1150 ierr = PetscFree(name);CHKERRQ(ierr); 1151 /* Bcast defaultValue */ 1152 if (!rank) (*labelNew)->defaultValue = label->defaultValue; 1153 ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1154 /* Distribute stratum values over the SF and get the point mapping on the receiver */ 1155 ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr); 1156 /* Determine received stratum values and initialise new label*/ 1157 PetscHashICreate(stratumHash); 1158 ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr); 1159 for (p = 0; p < size; ++p) PetscHashIPut(stratumHash, leafStrata[p], ret, iter); 1160 PetscHashISize(stratumHash, (*labelNew)->numStrata); 1161 ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr); 1162 for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE; 1163 ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr); 1164 /* Turn leafStrata into indices rather than stratum values */ 1165 offset = 0; 1166 ierr = PetscHashIGetKeys(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr); 1167 ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr); 1168 for (p = 0; p < size; ++p) { 1169 for (s = 0; s < (*labelNew)->numStrata; ++s) { 1170 if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;} 1171 } 1172 } 1173 /* Rebuild the point strata on the receiver */ 1174 ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr); 1175 ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1176 for (p=pStart; p<pEnd; p++) { 1177 ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1178 ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1179 for (s=0; s<dof; s++) { 1180 (*labelNew)->stratumSizes[leafStrata[offset+s]]++; 1181 } 1182 } 1183 ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr); 1184 ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr); 1185 ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr); 1186 for (s = 0; s < (*labelNew)->numStrata; ++s) { 1187 PetscHashICreate((*labelNew)->ht[s]); 1188 ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr); 1189 } 1190 /* Insert points into new strata */ 1191 ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr); 1192 ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1193 for (p=pStart; p<pEnd; p++) { 1194 ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1195 ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1196 for (s=0; s<dof; s++) { 1197 stratum = leafStrata[offset+s]; 1198 points[stratum][strataIdx[stratum]++] = p; 1199 } 1200 } 1201 for (s = 0; s < (*labelNew)->numStrata; s++) { 1202 ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr); 1203 ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr); 1204 } 1205 ierr = PetscFree(points);CHKERRQ(ierr); 1206 PetscHashIDestroy(stratumHash); 1207 ierr = PetscFree(leafStrata);CHKERRQ(ierr); 1208 ierr = PetscFree(strataIdx);CHKERRQ(ierr); 1209 ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 1210 PetscFunctionReturn(0); 1211 } 1212 1213 /*@ 1214 DMLabelGather - Gather all label values from leafs into roots 1215 1216 Input Parameters: 1217 + label - the DMLabel 1218 - sf - the Star Forest point communication map 1219 1220 Output Parameters: 1221 . labelNew - the new DMLabel with localised leaf values 1222 1223 Level: developer 1224 1225 Note: This is the inverse operation to DMLabelDistribute. 1226 1227 .seealso: DMLabelDistribute() 1228 @*/ 1229 PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew) 1230 { 1231 MPI_Comm comm; 1232 PetscSection rootSection; 1233 PetscSF sfLabel; 1234 PetscSFNode *rootPoints, *leafPoints; 1235 PetscInt p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset; 1236 const PetscInt *rootDegree, *ilocal; 1237 PetscInt *rootStrata; 1238 char *name; 1239 PetscInt nameSize; 1240 size_t len = 0; 1241 PetscMPIInt rank, size; 1242 PetscErrorCode ierr; 1243 1244 PetscFunctionBegin; 1245 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 1246 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1247 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 1248 /* Bcast name */ 1249 if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 1250 nameSize = len; 1251 ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1252 ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1253 if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 1254 ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1255 ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 1256 ierr = PetscFree(name);CHKERRQ(ierr); 1257 /* Gather rank/index pairs of leaves into local roots to build 1258 an inverse, multi-rooted SF. Note that this ignores local leaf 1259 indexing due to the use of the multiSF in PetscSFGather. */ 1260 ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr); 1261 ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr); 1262 for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1; 1263 for (p = 0; p < nleaves; p++) { 1264 leafPoints[ilocal[p]].index = ilocal[p]; 1265 leafPoints[ilocal[p]].rank = rank; 1266 } 1267 ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr); 1268 ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr); 1269 for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p]; 1270 ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr); 1271 ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 1272 ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 1273 ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr); 1274 ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 1275 /* Migrate label over inverted SF to pull stratum values at leaves into roots. */ 1276 ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr); 1277 /* Rebuild the point strata on the receiver */ 1278 for (p = 0, idx = 0; p < nroots; p++) { 1279 for (d = 0; d < rootDegree[p]; d++) { 1280 ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr); 1281 ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr); 1282 for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);} 1283 } 1284 idx += rootDegree[p]; 1285 } 1286 ierr = PetscFree(leafPoints);CHKERRQ(ierr); 1287 ierr = PetscFree(rootStrata);CHKERRQ(ierr); 1288 ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 1289 ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr); 1290 PetscFunctionReturn(0); 1291 } 1292 1293 /*@ 1294 DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label 1295 1296 Input Parameter: 1297 . label - the DMLabel 1298 1299 Output Parameters: 1300 + section - the section giving offsets for each stratum 1301 - is - An IS containing all the label points 1302 1303 Level: developer 1304 1305 .seealso: DMLabelDistribute() 1306 @*/ 1307 PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is) 1308 { 1309 IS vIS; 1310 const PetscInt *values; 1311 PetscInt *points; 1312 PetscInt nV, vS = 0, vE = 0, v, N; 1313 PetscErrorCode ierr; 1314 1315 PetscFunctionBegin; 1316 ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr); 1317 ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr); 1318 ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr); 1319 if (nV) {vS = values[0]; vE = values[0]+1;} 1320 for (v = 1; v < nV; ++v) { 1321 vS = PetscMin(vS, values[v]); 1322 vE = PetscMax(vE, values[v]+1); 1323 } 1324 ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr); 1325 ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr); 1326 for (v = 0; v < nV; ++v) { 1327 PetscInt n; 1328 1329 ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr); 1330 ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr); 1331 } 1332 ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 1333 ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr); 1334 ierr = PetscMalloc1(N, &points);CHKERRQ(ierr); 1335 for (v = 0; v < nV; ++v) { 1336 IS is; 1337 const PetscInt *spoints; 1338 PetscInt dof, off, p; 1339 1340 ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr); 1341 ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr); 1342 ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr); 1343 ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr); 1344 for (p = 0; p < dof; ++p) points[off+p] = spoints[p]; 1345 ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr); 1346 ierr = ISDestroy(&is);CHKERRQ(ierr); 1347 } 1348 ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr); 1349 ierr = ISDestroy(&vIS);CHKERRQ(ierr); 1350 ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr); 1351 PetscFunctionReturn(0); 1352 } 1353 1354 /*@ 1355 PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using 1356 the local section and an SF describing the section point overlap. 1357 1358 Input Parameters: 1359 + s - The PetscSection for the local field layout 1360 . sf - The SF describing parallel layout of the section points 1361 . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs 1362 . label - The label specifying the points 1363 - labelValue - The label stratum specifying the points 1364 1365 Output Parameter: 1366 . gsection - The PetscSection for the global field layout 1367 1368 Note: This gives negative sizes and offsets to points not owned by this process 1369 1370 Level: developer 1371 1372 .seealso: PetscSectionCreate() 1373 @*/ 1374 PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection) 1375 { 1376 PetscInt *neg = NULL, *tmpOff = NULL; 1377 PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots; 1378 PetscErrorCode ierr; 1379 1380 PetscFunctionBegin; 1381 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr); 1382 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 1383 ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr); 1384 ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 1385 if (nroots >= 0) { 1386 if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart); 1387 ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr); 1388 if (nroots > pEnd-pStart) { 1389 ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr); 1390 } else { 1391 tmpOff = &(*gsection)->atlasDof[-pStart]; 1392 } 1393 } 1394 /* Mark ghost points with negative dof */ 1395 for (p = pStart; p < pEnd; ++p) { 1396 PetscInt value; 1397 1398 ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr); 1399 if (value != labelValue) continue; 1400 ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr); 1401 ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr); 1402 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 1403 if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);} 1404 if (neg) neg[p] = -(dof+1); 1405 } 1406 ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr); 1407 if (nroots >= 0) { 1408 ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1409 ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1410 if (nroots > pEnd-pStart) { 1411 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];} 1412 } 1413 } 1414 /* Calculate new sizes, get proccess offset, and calculate point offsets */ 1415 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1416 cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0; 1417 (*gsection)->atlasOff[p] = off; 1418 off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0; 1419 } 1420 ierr = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr); 1421 globalOff -= off; 1422 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1423 (*gsection)->atlasOff[p] += globalOff; 1424 if (neg) neg[p] = -((*gsection)->atlasOff[p]+1); 1425 } 1426 /* Put in negative offsets for ghost points */ 1427 if (nroots >= 0) { 1428 ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1429 ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1430 if (nroots > pEnd-pStart) { 1431 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];} 1432 } 1433 } 1434 if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);} 1435 ierr = PetscFree(neg);CHKERRQ(ierr); 1436 PetscFunctionReturn(0); 1437 } 1438 1439 typedef struct _n_PetscSectionSym_Label 1440 { 1441 DMLabel label; 1442 PetscCopyMode *modes; 1443 PetscInt *sizes; 1444 const PetscInt ***perms; 1445 const PetscScalar ***rots; 1446 PetscInt (*minMaxOrients)[2]; 1447 PetscInt numStrata; /* numStrata is only increasing, functions as a state */ 1448 } PetscSectionSym_Label; 1449 1450 static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym) 1451 { 1452 PetscInt i, j; 1453 PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data; 1454 PetscErrorCode ierr; 1455 1456 PetscFunctionBegin; 1457 for (i = 0; i <= sl->numStrata; i++) { 1458 if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) { 1459 for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 1460 if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);} 1461 if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);} 1462 } 1463 if (sl->perms[i]) { 1464 const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]]; 1465 1466 ierr = PetscFree(perms);CHKERRQ(ierr); 1467 } 1468 if (sl->rots[i]) { 1469 const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]]; 1470 1471 ierr = PetscFree(rots);CHKERRQ(ierr); 1472 } 1473 } 1474 } 1475 ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr); 1476 ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr); 1477 sl->numStrata = 0; 1478 PetscFunctionReturn(0); 1479 } 1480 1481 static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym) 1482 { 1483 PetscErrorCode ierr; 1484 1485 PetscFunctionBegin; 1486 ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr); 1487 ierr = PetscFree(sym->data);CHKERRQ(ierr); 1488 PetscFunctionReturn(0); 1489 } 1490 1491 static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer) 1492 { 1493 PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data; 1494 PetscBool isAscii; 1495 DMLabel label = sl->label; 1496 PetscErrorCode ierr; 1497 1498 PetscFunctionBegin; 1499 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr); 1500 if (isAscii) { 1501 PetscInt i, j, k; 1502 PetscViewerFormat format; 1503 1504 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1505 if (label) { 1506 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1507 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 1508 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1509 ierr = DMLabelView(label, viewer);CHKERRQ(ierr); 1510 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1511 } else { 1512 ierr = PetscViewerASCIIPrintf(viewer," Label '%s'\n",sl->label->name);CHKERRQ(ierr); 1513 } 1514 } else { 1515 ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr); 1516 } 1517 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1518 for (i = 0; i <= sl->numStrata; i++) { 1519 PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue; 1520 1521 if (!(sl->perms[i] || sl->rots[i])) { 1522 ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr); 1523 } else { 1524 ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr); 1525 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1526 ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr); 1527 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 1528 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1529 for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 1530 if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) { 1531 ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr); 1532 } else { 1533 PetscInt tab; 1534 1535 ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr); 1536 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1537 ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr); 1538 if (sl->perms[i] && sl->perms[i][j]) { 1539 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr); 1540 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr); 1541 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);} 1542 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1543 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr); 1544 } 1545 if (sl->rots[i] && sl->rots[i][j]) { 1546 ierr = PetscViewerASCIIPrintf(viewer,"Rotations: ");CHKERRQ(ierr); 1547 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr); 1548 #if defined(PETSC_USE_COMPLEX) 1549 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f+i*%+f",PetscRealPart(sl->rots[i][j][k]),PetscImaginaryPart(sl->rots[i][j][k]));CHKERRQ(ierr);} 1550 #else 1551 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);} 1552 #endif 1553 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1554 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr); 1555 } 1556 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1557 } 1558 } 1559 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1560 } 1561 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1562 } 1563 } 1564 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1565 } 1566 PetscFunctionReturn(0); 1567 } 1568 1569 /*@ 1570 PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries 1571 1572 Logically collective on sym 1573 1574 Input parameters: 1575 + sym - the section symmetries 1576 - label - the DMLabel describing the types of points 1577 1578 Level: developer: 1579 1580 .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms() 1581 @*/ 1582 PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label) 1583 { 1584 PetscSectionSym_Label *sl; 1585 PetscErrorCode ierr; 1586 1587 PetscFunctionBegin; 1588 PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1); 1589 sl = (PetscSectionSym_Label *) sym->data; 1590 if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);} 1591 if (label) { 1592 label->refct++; 1593 sl->label = label; 1594 ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr); 1595 ierr = 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);CHKERRQ(ierr); 1596 ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr); 1597 ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr); 1598 ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr); 1599 ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr); 1600 ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr); 1601 } 1602 PetscFunctionReturn(0); 1603 } 1604 1605 /*@C 1606 PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum 1607 1608 Logically collective on PetscSectionSym 1609 1610 InputParameters: 1611 + sys - the section symmetries 1612 . stratum - the stratum value in the label that we are assigning symmetries for 1613 . size - the number of dofs for points in the stratum of the label 1614 . minOrient - the smallest orientation for a point in this stratum 1615 . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient)) 1616 . mode - how sym should copy the perms and rots arrays 1617 . perms - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation. A NULL permutation is the identity 1618 + 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 1619 1620 Level: developer 1621 1622 .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel() 1623 @*/ 1624 PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots) 1625 { 1626 PetscInt i, j, k; 1627 PetscSectionSym_Label *sl; 1628 PetscErrorCode ierr; 1629 1630 PetscFunctionBegin; 1631 PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1); 1632 sl = (PetscSectionSym_Label *) sym->data; 1633 if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet"); 1634 for (i = 0; i <= sl->numStrata; i++) { 1635 PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue; 1636 1637 if (stratum == value) break; 1638 } 1639 if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,sl->label->name); 1640 sl->sizes[i] = size; 1641 sl->modes[i] = mode; 1642 sl->minMaxOrients[i][0] = minOrient; 1643 sl->minMaxOrients[i][1] = maxOrient; 1644 if (mode == PETSC_COPY_VALUES) { 1645 if (perms) { 1646 PetscInt **ownPerms; 1647 1648 ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr); 1649 for (j = 0; j < maxOrient-minOrient; j++) { 1650 if (perms[j]) { 1651 ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr); 1652 for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];} 1653 } 1654 } 1655 sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient]; 1656 } 1657 if (rots) { 1658 PetscScalar **ownRots; 1659 1660 ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr); 1661 for (j = 0; j < maxOrient-minOrient; j++) { 1662 if (rots[j]) { 1663 ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr); 1664 for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];} 1665 } 1666 } 1667 sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient]; 1668 } 1669 } else { 1670 sl->perms[i] = perms ? &perms[-minOrient] : NULL; 1671 sl->rots[i] = rots ? &rots[-minOrient] : NULL; 1672 } 1673 PetscFunctionReturn(0); 1674 } 1675 1676 static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots) 1677 { 1678 PetscInt i, j, numStrata; 1679 PetscSectionSym_Label *sl; 1680 DMLabel label; 1681 PetscErrorCode ierr; 1682 1683 PetscFunctionBegin; 1684 sl = (PetscSectionSym_Label *) sym->data; 1685 numStrata = sl->numStrata; 1686 label = sl->label; 1687 for (i = 0; i < numPoints; i++) { 1688 PetscInt point = points[2*i]; 1689 PetscInt ornt = points[2*i+1]; 1690 1691 for (j = 0; j < numStrata; j++) { 1692 if (label->validIS[j]) { 1693 PetscInt k; 1694 1695 ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr); 1696 if (k >= 0) break; 1697 } else { 1698 PetscBool has; 1699 1700 PetscHashIHasKey(label->ht[j], point, has); 1701 if (has) break; 1702 } 1703 } 1704 if ((sl->minMaxOrients[j][1] > sl->minMaxOrients[j][0]) && (ornt < sl->minMaxOrients[j][0] || ornt >= sl->minMaxOrients[j][1])) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D orientation %D not in range [%D, %D) for stratum %D",point,ornt,sl->minMaxOrients[j][0],sl->minMaxOrients[j][1],j < numStrata ? label->stratumValues[j] : label->defaultValue); 1705 if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;} 1706 if (rots) {rots[i] = sl->rots[j] ? sl->rots[j][ornt] : NULL;} 1707 } 1708 PetscFunctionReturn(0); 1709 } 1710 1711 PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym) 1712 { 1713 PetscSectionSym_Label *sl; 1714 PetscErrorCode ierr; 1715 1716 PetscFunctionBegin; 1717 ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr); 1718 sym->ops->getpoints = PetscSectionSymGetPoints_Label; 1719 sym->ops->view = PetscSectionSymView_Label; 1720 sym->ops->destroy = PetscSectionSymDestroy_Label; 1721 sym->data = (void *) sl; 1722 PetscFunctionReturn(0); 1723 } 1724 1725 /*@ 1726 PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label 1727 1728 Collective 1729 1730 Input Parameters: 1731 + comm - the MPI communicator for the new symmetry 1732 - label - the label defining the strata 1733 1734 Output Parameters: 1735 . sym - the section symmetries 1736 1737 Level: developer 1738 1739 .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms() 1740 @*/ 1741 PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym) 1742 { 1743 PetscErrorCode ierr; 1744 1745 PetscFunctionBegin; 1746 ierr = DMInitializePackage();CHKERRQ(ierr); 1747 ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr); 1748 ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr); 1749 ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr); 1750 PetscFunctionReturn(0); 1751 } 1752