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