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