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 = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&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 PetscErrorCode DMLabelDestroy(DMLabel *label) 278 { 279 PetscInt v; 280 PetscErrorCode ierr; 281 282 PetscFunctionBegin; 283 if (!(*label)) PetscFunctionReturn(0); 284 if (--(*label)->refct > 0) PetscFunctionReturn(0); 285 ierr = PetscFree((*label)->name);CHKERRQ(ierr); 286 ierr = PetscFree((*label)->stratumValues);CHKERRQ(ierr); 287 ierr = PetscFree((*label)->stratumSizes);CHKERRQ(ierr); 288 for (v = 0; v < (*label)->numStrata; ++v) {ierr = ISDestroy(&((*label)->points[v]));CHKERRQ(ierr);} 289 ierr = PetscFree((*label)->points);CHKERRQ(ierr); 290 ierr = PetscFree((*label)->validIS);CHKERRQ(ierr); 291 if ((*label)->ht) { 292 for (v = 0; v < (*label)->numStrata; ++v) {PetscHashIDestroy((*label)->ht[v]);} 293 ierr = PetscFree((*label)->ht);CHKERRQ(ierr); 294 } 295 ierr = PetscBTDestroy(&(*label)->bt);CHKERRQ(ierr); 296 ierr = PetscFree(*label);CHKERRQ(ierr); 297 PetscFunctionReturn(0); 298 } 299 300 PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew) 301 { 302 PetscInt v; 303 PetscErrorCode ierr; 304 305 PetscFunctionBegin; 306 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 307 ierr = PetscNew(labelnew);CHKERRQ(ierr); 308 ierr = PetscStrallocpy(label->name, &(*labelnew)->name);CHKERRQ(ierr); 309 310 (*labelnew)->refct = 1; 311 (*labelnew)->numStrata = label->numStrata; 312 (*labelnew)->defaultValue = label->defaultValue; 313 if (label->numStrata) { 314 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr); 315 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr); 316 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr); 317 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr); 318 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr); 319 /* Could eliminate unused space here */ 320 for (v = 0; v < label->numStrata; ++v) { 321 PetscHashICreate((*labelnew)->ht[v]); 322 (*labelnew)->validIS[v] = PETSC_TRUE; 323 (*labelnew)->stratumValues[v] = label->stratumValues[v]; 324 (*labelnew)->stratumSizes[v] = label->stratumSizes[v]; 325 ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr); 326 (*labelnew)->points[v] = label->points[v]; 327 } 328 } 329 (*labelnew)->pStart = -1; 330 (*labelnew)->pEnd = -1; 331 (*labelnew)->bt = NULL; 332 PetscFunctionReturn(0); 333 } 334 335 /* This can be hooked into SetValue(), ClearValue(), etc. for updating */ 336 PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd) 337 { 338 PetscInt v; 339 PetscErrorCode ierr; 340 341 PetscFunctionBegin; 342 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 343 if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 344 label->pStart = pStart; 345 label->pEnd = pEnd; 346 ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr); 347 ierr = PetscBTMemzero(pEnd - pStart, label->bt);CHKERRQ(ierr); 348 for (v = 0; v < label->numStrata; ++v) { 349 const PetscInt *points; 350 PetscInt i; 351 352 ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 353 for (i = 0; i < label->stratumSizes[v]; ++i) { 354 const PetscInt point = points[i]; 355 356 if ((point < pStart) || (point >= pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, pStart, pEnd); 357 ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr); 358 } 359 ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 360 } 361 PetscFunctionReturn(0); 362 } 363 364 PetscErrorCode DMLabelDestroyIndex(DMLabel label) 365 { 366 PetscErrorCode ierr; 367 368 PetscFunctionBegin; 369 label->pStart = -1; 370 label->pEnd = -1; 371 if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 372 PetscFunctionReturn(0); 373 } 374 375 /*@ 376 DMLabelHasValue - Determine whether a label assigns the value to any point 377 378 Input Parameters: 379 + label - the DMLabel 380 - value - the value 381 382 Output Parameter: 383 . contains - Flag indicating whether the label maps this value to any point 384 385 Level: developer 386 387 .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue() 388 @*/ 389 PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains) 390 { 391 PetscInt v; 392 393 PetscFunctionBegin; 394 PetscValidPointer(contains, 3); 395 for (v = 0; v < label->numStrata; ++v) { 396 if (value == label->stratumValues[v]) break; 397 } 398 *contains = (v < label->numStrata ? PETSC_TRUE : PETSC_FALSE); 399 PetscFunctionReturn(0); 400 } 401 402 /*@ 403 DMLabelHasPoint - Determine whether a label assigns a value to a point 404 405 Input Parameters: 406 + label - the DMLabel 407 - point - the point 408 409 Output Parameter: 410 . contains - Flag indicating whether the label maps this point to a value 411 412 Note: The user must call DMLabelCreateIndex() before this function. 413 414 Level: developer 415 416 .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue() 417 @*/ 418 PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains) 419 { 420 PetscErrorCode ierr; 421 422 PetscFunctionBeginHot; 423 PetscValidPointer(contains, 3); 424 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 425 #if defined(PETSC_USE_DEBUG) 426 if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()"); 427 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); 428 #endif 429 *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE; 430 PetscFunctionReturn(0); 431 } 432 433 /*@ 434 DMLabelStratumHasPoint - Return true if the stratum contains a point 435 436 Input Parameters: 437 + label - the DMLabel 438 . value - the stratum value 439 - point - the point 440 441 Output Parameter: 442 . contains - true if the stratum contains the point 443 444 Level: intermediate 445 446 .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue() 447 @*/ 448 PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains) 449 { 450 PetscInt v; 451 PetscErrorCode ierr; 452 453 PetscFunctionBegin; 454 PetscValidPointer(contains, 4); 455 *contains = PETSC_FALSE; 456 for (v = 0; v < label->numStrata; ++v) { 457 if (label->stratumValues[v] == value) { 458 if (label->validIS[v]) { 459 PetscInt i; 460 461 ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr); 462 if (i >= 0) { 463 *contains = PETSC_TRUE; 464 break; 465 } 466 } else { 467 PetscBool has; 468 469 PetscHashIHasKey(label->ht[v], point, has); 470 if (has) { 471 *contains = PETSC_TRUE; 472 break; 473 } 474 } 475 } 476 } 477 PetscFunctionReturn(0); 478 } 479 480 /* 481 DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 482 When a label is created, it is initialized to -1. 483 484 Input parameter: 485 . label - a DMLabel object 486 487 Output parameter: 488 . defaultValue - the default value 489 490 Level: beginner 491 492 .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 493 */ 494 PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue) 495 { 496 PetscFunctionBegin; 497 *defaultValue = label->defaultValue; 498 PetscFunctionReturn(0); 499 } 500 501 /* 502 DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 503 When a label is created, it is initialized to -1. 504 505 Input parameter: 506 . label - a DMLabel object 507 508 Output parameter: 509 . defaultValue - the default value 510 511 Level: beginner 512 513 .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 514 */ 515 PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue) 516 { 517 PetscFunctionBegin; 518 label->defaultValue = defaultValue; 519 PetscFunctionReturn(0); 520 } 521 522 /*@ 523 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()) 524 525 Input Parameters: 526 + label - the DMLabel 527 - point - the point 528 529 Output Parameter: 530 . value - The point value, or -1 531 532 Level: intermediate 533 534 .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 535 @*/ 536 PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value) 537 { 538 PetscInt v; 539 PetscErrorCode ierr; 540 541 PetscFunctionBegin; 542 PetscValidPointer(value, 3); 543 *value = label->defaultValue; 544 for (v = 0; v < label->numStrata; ++v) { 545 if (label->validIS[v]) { 546 PetscInt i; 547 548 ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr); 549 if (i >= 0) { 550 *value = label->stratumValues[v]; 551 break; 552 } 553 } else { 554 PetscBool has; 555 556 PetscHashIHasKey(label->ht[v], point, has); 557 if (has) { 558 *value = label->stratumValues[v]; 559 break; 560 } 561 } 562 } 563 PetscFunctionReturn(0); 564 } 565 566 /*@ 567 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. 568 569 Input Parameters: 570 + label - the DMLabel 571 . point - the point 572 - value - The point value 573 574 Level: intermediate 575 576 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 577 @*/ 578 PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value) 579 { 580 PETSC_UNUSED PetscHashIIter iter, ret; 581 PetscInt v; 582 PetscErrorCode ierr; 583 584 PetscFunctionBegin; 585 /* Find, or add, label value */ 586 if (value == label->defaultValue) PetscFunctionReturn(0); 587 for (v = 0; v < label->numStrata; ++v) { 588 if (label->stratumValues[v] == value) break; 589 } 590 /* Create new table */ 591 if (v >= label->numStrata) {ierr = DMLabelAddStratum(label, value);CHKERRQ(ierr);} 592 ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 593 /* Set key */ 594 PetscHashIPut(label->ht[v], point, ret, iter); 595 PetscFunctionReturn(0); 596 } 597 598 /*@ 599 DMLabelClearValue - Clear the value a label assigns to a point 600 601 Input Parameters: 602 + label - the DMLabel 603 . point - the point 604 - value - The point value 605 606 Level: intermediate 607 608 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue() 609 @*/ 610 PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value) 611 { 612 PetscInt v; 613 PetscErrorCode ierr; 614 615 PetscFunctionBegin; 616 /* Find label value */ 617 for (v = 0; v < label->numStrata; ++v) { 618 if (label->stratumValues[v] == value) break; 619 } 620 if (v >= label->numStrata) PetscFunctionReturn(0); 621 if (label->validIS[v]) { 622 ierr = DMLabelMakeInvalid_Private(label,v);CHKERRQ(ierr); 623 } 624 if (label->bt) { 625 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); 626 ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 627 } 628 ierr = PetscHashIDelKey(label->ht[v], point);CHKERRQ(ierr); 629 PetscFunctionReturn(0); 630 } 631 632 /*@ 633 DMLabelInsertIS - Set all points in the IS to a value 634 635 Input Parameters: 636 + label - the DMLabel 637 . is - the point IS 638 - value - The point value 639 640 Level: intermediate 641 642 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 643 @*/ 644 PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value) 645 { 646 const PetscInt *points; 647 PetscInt n, p; 648 PetscErrorCode ierr; 649 650 PetscFunctionBegin; 651 PetscValidHeaderSpecific(is, IS_CLASSID, 2); 652 ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 653 ierr = ISGetIndices(is, &points);CHKERRQ(ierr); 654 for (p = 0; p < n; ++p) {ierr = DMLabelSetValue(label, points[p], value);CHKERRQ(ierr);} 655 ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr); 656 PetscFunctionReturn(0); 657 } 658 659 PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues) 660 { 661 PetscFunctionBegin; 662 PetscValidPointer(numValues, 2); 663 *numValues = label->numStrata; 664 PetscFunctionReturn(0); 665 } 666 667 PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values) 668 { 669 PetscErrorCode ierr; 670 671 PetscFunctionBegin; 672 PetscValidPointer(values, 2); 673 ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr); 674 PetscFunctionReturn(0); 675 } 676 677 PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists) 678 { 679 PetscInt v; 680 681 PetscFunctionBegin; 682 PetscValidPointer(exists, 3); 683 *exists = PETSC_FALSE; 684 for (v = 0; v < label->numStrata; ++v) { 685 if (label->stratumValues[v] == value) { 686 *exists = PETSC_TRUE; 687 break; 688 } 689 } 690 PetscFunctionReturn(0); 691 } 692 693 PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size) 694 { 695 PetscInt v; 696 PetscErrorCode ierr; 697 698 PetscFunctionBegin; 699 PetscValidPointer(size, 3); 700 *size = 0; 701 for (v = 0; v < label->numStrata; ++v) { 702 if (label->stratumValues[v] == value) { 703 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 704 *size = label->stratumSizes[v]; 705 break; 706 } 707 } 708 PetscFunctionReturn(0); 709 } 710 711 PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end) 712 { 713 PetscInt v; 714 PetscErrorCode ierr; 715 716 PetscFunctionBegin; 717 if (start) {PetscValidPointer(start, 3); *start = 0;} 718 if (end) {PetscValidPointer(end, 4); *end = 0;} 719 for (v = 0; v < label->numStrata; ++v) { 720 PetscInt min, max; 721 722 if (label->stratumValues[v] != value) continue; 723 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 724 if (label->stratumSizes[v] <= 0) break; 725 ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr); 726 if (start) *start = min; 727 if (end) *end = max+1; 728 break; 729 } 730 PetscFunctionReturn(0); 731 } 732 733 PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points) 734 { 735 PetscInt v; 736 PetscErrorCode ierr; 737 738 PetscFunctionBegin; 739 PetscValidPointer(points, 3); 740 *points = NULL; 741 for (v = 0; v < label->numStrata; ++v) { 742 if (label->stratumValues[v] == value) { 743 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 744 if (label->validIS[v]) { 745 ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr); 746 *points = label->points[v]; 747 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to implement this to speedup Stratify"); 748 break; 749 } 750 } 751 PetscFunctionReturn(0); 752 } 753 754 PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is) 755 { 756 PetscInt v, numStrata; 757 PetscErrorCode ierr; 758 759 PetscFunctionBegin; 760 numStrata = label->numStrata; 761 for (v = 0; v < numStrata; v++) { 762 if (label->stratumValues[v] == value) break; 763 } 764 if (v >= numStrata) {ierr = DMLabelAddStratum(label,value);CHKERRQ(ierr);} 765 if (is == label->points[v]) PetscFunctionReturn(0); 766 ierr = DMLabelClearStratum(label,value);CHKERRQ(ierr); 767 ierr = ISGetLocalSize(is,&(label->stratumSizes[v]));CHKERRQ(ierr); 768 label->stratumValues[v] = value; 769 label->validIS[v] = PETSC_TRUE; 770 ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 771 ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 772 if (label->bt) { 773 const PetscInt *points; 774 PetscInt p; 775 776 ierr = ISGetIndices(is,&points);CHKERRQ(ierr); 777 for (p = 0; p < label->stratumSizes[v]; ++p) { 778 const PetscInt point = points[p]; 779 780 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); 781 ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 782 } 783 } 784 label->points[v] = is; 785 PetscFunctionReturn(0); 786 } 787 788 789 PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value) 790 { 791 PetscInt v; 792 PetscErrorCode ierr; 793 794 PetscFunctionBegin; 795 for (v = 0; v < label->numStrata; ++v) { 796 if (label->stratumValues[v] == value) break; 797 } 798 if (v >= label->numStrata) PetscFunctionReturn(0); 799 if (label->validIS[v]) { 800 if (label->bt) { 801 PetscInt i; 802 const PetscInt *points; 803 804 ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 805 for (i = 0; i < label->stratumSizes[v]; ++i) { 806 const PetscInt point = points[i]; 807 808 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); 809 ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 810 } 811 ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr); 812 } 813 ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 814 label->stratumSizes[v] = 0; 815 ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr); 816 ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr); 817 } else { 818 PetscHashIClear(label->ht[v]); 819 } 820 PetscFunctionReturn(0); 821 } 822 823 PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end) 824 { 825 PetscInt v; 826 PetscErrorCode ierr; 827 828 PetscFunctionBegin; 829 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 830 label->pStart = start; 831 label->pEnd = end; 832 if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 833 /* Could squish offsets, but would only make sense if I reallocate the storage */ 834 for (v = 0; v < label->numStrata; ++v) { 835 PetscInt off, q; 836 const PetscInt *points; 837 PetscInt *pointsNew = NULL; 838 839 ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 840 for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) { 841 const PetscInt point = points[q]; 842 843 if ((point < start) || (point >= end)) { 844 if (!pointsNew) { 845 ierr = PetscMalloc1(label->stratumSizes[v],&pointsNew);CHKERRQ(ierr); 846 ierr = PetscMemcpy(pointsNew,points,(size_t) off * sizeof(PetscInt));CHKERRQ(ierr); 847 } 848 continue; 849 } 850 if (pointsNew) { 851 pointsNew[off++] = point; 852 } 853 } 854 ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 855 if (pointsNew) { 856 ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 857 ierr = ISCreateGeneral(PETSC_COMM_SELF,off,pointsNew,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr); 858 ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr); 859 } 860 label->stratumSizes[v] = off; 861 } 862 ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr); 863 PetscFunctionReturn(0); 864 } 865 866 PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew) 867 { 868 const PetscInt *perm; 869 PetscInt numValues, numPoints, v, q; 870 PetscErrorCode ierr; 871 872 PetscFunctionBegin; 873 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 874 ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr); 875 ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr); 876 ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr); 877 ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr); 878 for (v = 0; v < numValues; ++v) { 879 const PetscInt size = (*labelNew)->stratumSizes[v]; 880 const PetscInt *points; 881 PetscInt *pointsNew; 882 883 ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 884 ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr); 885 for (q = 0; q < size; ++q) { 886 const PetscInt point = points[q]; 887 888 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); 889 pointsNew[q] = perm[point]; 890 } 891 ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 892 ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr); 893 ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr); 894 ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr); 895 ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr); 896 } 897 ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr); 898 if (label->bt) { 899 ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 900 ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr); 901 } 902 PetscFunctionReturn(0); 903 } 904 905 PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata) 906 { 907 MPI_Comm comm; 908 PetscInt s, l, nroots, nleaves, dof, offset, size; 909 PetscInt *remoteOffsets, *rootStrata, *rootIdx; 910 PetscSection rootSection; 911 PetscSF labelSF; 912 PetscErrorCode ierr; 913 914 PetscFunctionBegin; 915 if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 916 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 917 /* Build a section of stratum values per point, generate the according SF 918 and distribute point-wise stratum values to leaves. */ 919 ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 920 ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 921 ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr); 922 if (label) { 923 for (s = 0; s < label->numStrata; ++s) { 924 const PetscInt *points; 925 926 ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 927 for (l = 0; l < label->stratumSizes[s]; l++) { 928 ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr); 929 ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr); 930 } 931 ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 932 } 933 } 934 ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 935 /* Create a point-wise array of stratum values */ 936 ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr); 937 ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr); 938 ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr); 939 if (label) { 940 for (s = 0; s < label->numStrata; ++s) { 941 const PetscInt *points; 942 943 ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 944 for (l = 0; l < label->stratumSizes[s]; l++) { 945 const PetscInt p = points[l]; 946 ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr); 947 rootStrata[offset+rootIdx[p]++] = label->stratumValues[s]; 948 } 949 ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 950 } 951 } 952 /* Build SF that maps label points to remote processes */ 953 ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr); 954 ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr); 955 ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr); 956 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 957 /* Send the strata for each point over the derived SF */ 958 ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr); 959 ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr); 960 ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 961 ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 962 /* Clean up */ 963 ierr = PetscFree(rootStrata);CHKERRQ(ierr); 964 ierr = PetscFree(rootIdx);CHKERRQ(ierr); 965 ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 966 ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr); 967 PetscFunctionReturn(0); 968 } 969 970 PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew) 971 { 972 MPI_Comm comm; 973 PetscSection leafSection; 974 PetscInt p, pStart, pEnd, s, size, dof, offset, stratum; 975 PetscInt *leafStrata, *strataIdx; 976 PetscInt **points; 977 char *name; 978 PetscInt nameSize; 979 PetscHashI stratumHash; 980 PETSC_UNUSED PetscHashIIter ret, iter; 981 size_t len = 0; 982 PetscMPIInt rank; 983 PetscErrorCode ierr; 984 985 PetscFunctionBegin; 986 if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 987 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 988 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 989 /* Bcast name */ 990 if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 991 nameSize = len; 992 ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 993 ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 994 if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 995 ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 996 ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 997 ierr = PetscFree(name);CHKERRQ(ierr); 998 /* Bcast defaultValue */ 999 if (!rank) (*labelNew)->defaultValue = label->defaultValue; 1000 ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1001 /* Distribute stratum values over the SF and get the point mapping on the receiver */ 1002 ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr); 1003 /* Determine received stratum values and initialise new label*/ 1004 PetscHashICreate(stratumHash); 1005 ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr); 1006 for (p = 0; p < size; ++p) PetscHashIPut(stratumHash, leafStrata[p], ret, iter); 1007 PetscHashISize(stratumHash, (*labelNew)->numStrata); 1008 ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr); 1009 for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE; 1010 ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr); 1011 /* Turn leafStrata into indices rather than stratum values */ 1012 offset = 0; 1013 ierr = PetscHashIGetKeys(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr); 1014 for (p = 0; p < size; ++p) { 1015 for (s = 0; s < (*labelNew)->numStrata; ++s) { 1016 if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;} 1017 } 1018 } 1019 /* Rebuild the point strata on the receiver */ 1020 ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr); 1021 ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1022 for (p=pStart; p<pEnd; p++) { 1023 ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1024 ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1025 for (s=0; s<dof; s++) { 1026 (*labelNew)->stratumSizes[leafStrata[offset+s]]++; 1027 } 1028 } 1029 ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr); 1030 ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr); 1031 ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr); 1032 for (s = 0; s < (*labelNew)->numStrata; ++s) { 1033 PetscHashICreate((*labelNew)->ht[s]); 1034 ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr); 1035 } 1036 /* Insert points into new strata */ 1037 ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr); 1038 ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1039 for (p=pStart; p<pEnd; p++) { 1040 ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1041 ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1042 for (s=0; s<dof; s++) { 1043 stratum = leafStrata[offset+s]; 1044 points[stratum][strataIdx[stratum]++] = p; 1045 } 1046 } 1047 for (s = 0; s < (*labelNew)->numStrata; s++) { 1048 ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr); 1049 ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr); 1050 } 1051 ierr = PetscFree(points);CHKERRQ(ierr); 1052 PetscHashIDestroy(stratumHash); 1053 ierr = PetscFree(leafStrata);CHKERRQ(ierr); 1054 ierr = PetscFree(strataIdx);CHKERRQ(ierr); 1055 ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 1056 PetscFunctionReturn(0); 1057 } 1058 1059 /*@ 1060 DMLabelGather - Gather all label values from leafs into roots 1061 1062 Input Parameters: 1063 + label - the DMLabel 1064 . point - the Star Forest point communication map 1065 1066 Input Parameters: 1067 + label - the new DMLabel with localised leaf values 1068 1069 Level: developer 1070 1071 Note: This is the inverse operation to DMLabelDistribute. 1072 1073 .seealso: DMLabelDistribute() 1074 @*/ 1075 PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew) 1076 { 1077 MPI_Comm comm; 1078 PetscSection rootSection; 1079 PetscSF sfLabel; 1080 PetscSFNode *rootPoints, *leafPoints; 1081 PetscInt p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset; 1082 const PetscInt *rootDegree, *ilocal; 1083 PetscInt *rootStrata; 1084 char *name; 1085 PetscInt nameSize; 1086 size_t len = 0; 1087 PetscMPIInt rank, size; 1088 PetscErrorCode ierr; 1089 1090 PetscFunctionBegin; 1091 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 1092 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1093 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 1094 /* Bcast name */ 1095 if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 1096 nameSize = len; 1097 ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1098 ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1099 if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 1100 ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1101 ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 1102 ierr = PetscFree(name);CHKERRQ(ierr); 1103 /* Gather rank/index pairs of leaves into local roots to build 1104 an inverse, multi-rooted SF. Note that this ignores local leaf 1105 indexing due to the use of the multiSF in PetscSFGather. */ 1106 ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr); 1107 ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr); 1108 for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1; 1109 for (p = 0; p < nleaves; p++) { 1110 leafPoints[ilocal[p]].index = ilocal[p]; 1111 leafPoints[ilocal[p]].rank = rank; 1112 } 1113 ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr); 1114 ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr); 1115 for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p]; 1116 ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr); 1117 ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 1118 ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 1119 ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr); 1120 ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 1121 /* Migrate label over inverted SF to pull stratum values at leaves into roots. */ 1122 ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr); 1123 /* Rebuild the point strata on the receiver */ 1124 for (p = 0, idx = 0; p < nroots; p++) { 1125 for (d = 0; d < rootDegree[p]; d++) { 1126 ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr); 1127 ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr); 1128 for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);} 1129 } 1130 idx += rootDegree[p]; 1131 } 1132 ierr = PetscFree(leafPoints);CHKERRQ(ierr); 1133 ierr = PetscFree(rootStrata);CHKERRQ(ierr); 1134 ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 1135 ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr); 1136 PetscFunctionReturn(0); 1137 } 1138 1139 PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is) 1140 { 1141 IS vIS; 1142 const PetscInt *values; 1143 PetscInt *points; 1144 PetscInt nV, vS = 0, vE = 0, v, N; 1145 PetscErrorCode ierr; 1146 1147 PetscFunctionBegin; 1148 ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr); 1149 ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr); 1150 ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr); 1151 if (nV) {vS = values[0]; vE = values[0]+1;} 1152 for (v = 1; v < nV; ++v) { 1153 vS = PetscMin(vS, values[v]); 1154 vE = PetscMax(vE, values[v]+1); 1155 } 1156 ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr); 1157 ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr); 1158 for (v = 0; v < nV; ++v) { 1159 PetscInt n; 1160 1161 ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr); 1162 ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr); 1163 } 1164 ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 1165 ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr); 1166 ierr = PetscMalloc1(N, &points);CHKERRQ(ierr); 1167 for (v = 0; v < nV; ++v) { 1168 IS is; 1169 const PetscInt *spoints; 1170 PetscInt dof, off, p; 1171 1172 ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr); 1173 ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr); 1174 ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr); 1175 ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr); 1176 for (p = 0; p < dof; ++p) points[off+p] = spoints[p]; 1177 ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr); 1178 ierr = ISDestroy(&is);CHKERRQ(ierr); 1179 } 1180 ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr); 1181 ierr = ISDestroy(&vIS);CHKERRQ(ierr); 1182 ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr); 1183 PetscFunctionReturn(0); 1184 } 1185 1186 /*@C 1187 PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using 1188 the local section and an SF describing the section point overlap. 1189 1190 Input Parameters: 1191 + s - The PetscSection for the local field layout 1192 . sf - The SF describing parallel layout of the section points 1193 . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs 1194 . label - The label specifying the points 1195 - labelValue - The label stratum specifying the points 1196 1197 Output Parameter: 1198 . gsection - The PetscSection for the global field layout 1199 1200 Note: This gives negative sizes and offsets to points not owned by this process 1201 1202 Level: developer 1203 1204 .seealso: PetscSectionCreate() 1205 @*/ 1206 PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection) 1207 { 1208 PetscInt *neg = NULL, *tmpOff = NULL; 1209 PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots; 1210 PetscErrorCode ierr; 1211 1212 PetscFunctionBegin; 1213 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr); 1214 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 1215 ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr); 1216 ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 1217 if (nroots >= 0) { 1218 if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart); 1219 ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr); 1220 if (nroots > pEnd-pStart) { 1221 ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr); 1222 } else { 1223 tmpOff = &(*gsection)->atlasDof[-pStart]; 1224 } 1225 } 1226 /* Mark ghost points with negative dof */ 1227 for (p = pStart; p < pEnd; ++p) { 1228 PetscInt value; 1229 1230 ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr); 1231 if (value != labelValue) continue; 1232 ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr); 1233 ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr); 1234 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 1235 if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);} 1236 if (neg) neg[p] = -(dof+1); 1237 } 1238 ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr); 1239 if (nroots >= 0) { 1240 ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1241 ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1242 if (nroots > pEnd-pStart) { 1243 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];} 1244 } 1245 } 1246 /* Calculate new sizes, get proccess offset, and calculate point offsets */ 1247 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1248 cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0; 1249 (*gsection)->atlasOff[p] = off; 1250 off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0; 1251 } 1252 ierr = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr); 1253 globalOff -= off; 1254 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1255 (*gsection)->atlasOff[p] += globalOff; 1256 if (neg) neg[p] = -((*gsection)->atlasOff[p]+1); 1257 } 1258 /* Put in negative offsets for ghost points */ 1259 if (nroots >= 0) { 1260 ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1261 ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1262 if (nroots > pEnd-pStart) { 1263 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];} 1264 } 1265 } 1266 if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);} 1267 ierr = PetscFree(neg);CHKERRQ(ierr); 1268 PetscFunctionReturn(0); 1269 } 1270 1271 typedef struct _n_PetscSectionSym_Label 1272 { 1273 DMLabel label; 1274 PetscCopyMode *modes; 1275 PetscInt *sizes; 1276 const PetscInt ***perms; 1277 const PetscScalar ***rots; 1278 PetscInt (*minMaxOrients)[2]; 1279 PetscInt numStrata; /* numStrata is only increasing, functions as a state */ 1280 } PetscSectionSym_Label; 1281 1282 static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym) 1283 { 1284 PetscInt i, j; 1285 PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data; 1286 PetscErrorCode ierr; 1287 1288 PetscFunctionBegin; 1289 for (i = 0; i <= sl->numStrata; i++) { 1290 if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) { 1291 for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 1292 if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);} 1293 if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);} 1294 } 1295 if (sl->perms[i]) { 1296 const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]]; 1297 1298 ierr = PetscFree(perms);CHKERRQ(ierr); 1299 } 1300 if (sl->rots[i]) { 1301 const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]]; 1302 1303 ierr = PetscFree(rots);CHKERRQ(ierr); 1304 } 1305 } 1306 } 1307 ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr); 1308 ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr); 1309 sl->numStrata = 0; 1310 PetscFunctionReturn(0); 1311 } 1312 1313 static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym) 1314 { 1315 PetscErrorCode ierr; 1316 1317 PetscFunctionBegin; 1318 ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr); 1319 ierr = PetscFree(sym->data);CHKERRQ(ierr); 1320 PetscFunctionReturn(0); 1321 } 1322 1323 static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer) 1324 { 1325 PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data; 1326 PetscBool isAscii; 1327 DMLabel label = sl->label; 1328 PetscErrorCode ierr; 1329 1330 PetscFunctionBegin; 1331 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr); 1332 if (isAscii) { 1333 PetscInt i, j, k; 1334 PetscViewerFormat format; 1335 1336 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1337 if (label) { 1338 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1339 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 1340 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1341 ierr = DMLabelView(label, viewer);CHKERRQ(ierr); 1342 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1343 } else { 1344 ierr = PetscViewerASCIIPrintf(viewer," Label '%s'\n",sl->label->name);CHKERRQ(ierr); 1345 } 1346 } else { 1347 ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr); 1348 } 1349 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1350 for (i = 0; i <= sl->numStrata; i++) { 1351 PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue; 1352 1353 if (!(sl->perms[i] || sl->rots[i])) { 1354 ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr); 1355 } else { 1356 ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr); 1357 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1358 ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr); 1359 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 1360 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1361 for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 1362 if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) { 1363 ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr); 1364 } else { 1365 PetscInt tab; 1366 1367 ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr); 1368 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1369 ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr); 1370 if (sl->perms[i] && sl->perms[i][j]) { 1371 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr); 1372 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr); 1373 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);} 1374 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1375 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr); 1376 } 1377 if (sl->rots[i] && sl->rots[i][j]) { 1378 ierr = PetscViewerASCIIPrintf(viewer,"Rotations: ");CHKERRQ(ierr); 1379 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr); 1380 #if defined(PETSC_USE_COMPLEX) 1381 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);} 1382 #else 1383 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);} 1384 #endif 1385 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1386 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr); 1387 } 1388 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1389 } 1390 } 1391 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1392 } 1393 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1394 } 1395 } 1396 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1397 } 1398 PetscFunctionReturn(0); 1399 } 1400 1401 /*@ 1402 PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries 1403 1404 Logically collective on sym 1405 1406 Input parameters: 1407 + sym - the section symmetries 1408 - label - the DMLabel describing the types of points 1409 1410 Level: developer: 1411 1412 .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms() 1413 @*/ 1414 PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label) 1415 { 1416 PetscSectionSym_Label *sl; 1417 PetscErrorCode ierr; 1418 1419 PetscFunctionBegin; 1420 PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1); 1421 sl = (PetscSectionSym_Label *) sym->data; 1422 if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);} 1423 if (label) { 1424 label->refct++; 1425 sl->label = label; 1426 ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr); 1427 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); 1428 ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr); 1429 ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr); 1430 ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr); 1431 ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr); 1432 ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr); 1433 } 1434 PetscFunctionReturn(0); 1435 } 1436 1437 /*@C 1438 PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum 1439 1440 Logically collective on PetscSectionSym 1441 1442 InputParameters: 1443 + sys - the section symmetries 1444 . stratum - the stratum value in the label that we are assigning symmetries for 1445 . size - the number of dofs for points in the stratum of the label 1446 . minOrient - the smallest orientation for a point in this stratum 1447 . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient)) 1448 . mode - how sym should copy the perms and rots arrays 1449 . perms - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation. A NULL permutation is the identity 1450 + 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 1451 1452 Level: developer 1453 1454 .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel() 1455 @*/ 1456 PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots) 1457 { 1458 PetscInt i, j, k; 1459 PetscSectionSym_Label *sl; 1460 PetscErrorCode ierr; 1461 1462 PetscFunctionBegin; 1463 PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1); 1464 sl = (PetscSectionSym_Label *) sym->data; 1465 if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet"); 1466 for (i = 0; i <= sl->numStrata; i++) { 1467 PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue; 1468 1469 if (stratum == value) break; 1470 } 1471 if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,sl->label->name); 1472 sl->sizes[i] = size; 1473 sl->modes[i] = mode; 1474 sl->minMaxOrients[i][0] = minOrient; 1475 sl->minMaxOrients[i][1] = maxOrient; 1476 if (mode == PETSC_COPY_VALUES) { 1477 if (perms) { 1478 PetscInt **ownPerms; 1479 1480 ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr); 1481 for (j = 0; j < maxOrient-minOrient; j++) { 1482 if (perms[j]) { 1483 ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr); 1484 for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];} 1485 } 1486 } 1487 sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient]; 1488 } 1489 if (rots) { 1490 PetscScalar **ownRots; 1491 1492 ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr); 1493 for (j = 0; j < maxOrient-minOrient; j++) { 1494 if (rots[j]) { 1495 ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr); 1496 for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];} 1497 } 1498 } 1499 sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient]; 1500 } 1501 } else { 1502 sl->perms[i] = perms ? &perms[-minOrient] : NULL; 1503 sl->rots[i] = rots ? &rots[-minOrient] : NULL; 1504 } 1505 PetscFunctionReturn(0); 1506 } 1507 1508 static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots) 1509 { 1510 PetscInt i, j, numStrata; 1511 PetscSectionSym_Label *sl; 1512 DMLabel label; 1513 PetscErrorCode ierr; 1514 1515 PetscFunctionBegin; 1516 sl = (PetscSectionSym_Label *) sym->data; 1517 numStrata = sl->numStrata; 1518 label = sl->label; 1519 for (i = 0; i < numPoints; i++) { 1520 PetscInt point = points[2*i]; 1521 PetscInt ornt = points[2*i+1]; 1522 1523 for (j = 0; j < numStrata; j++) { 1524 if (label->validIS[j]) { 1525 PetscInt k; 1526 1527 ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr); 1528 if (k >= 0) break; 1529 } else { 1530 PetscBool has; 1531 1532 PetscHashIHasKey(label->ht[j], point, has); 1533 if (has) break; 1534 } 1535 } 1536 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); 1537 if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;} 1538 if (rots) {rots[i] = sl->rots[j] ? sl->rots[j][ornt] : NULL;} 1539 } 1540 PetscFunctionReturn(0); 1541 } 1542 1543 PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym) 1544 { 1545 PetscSectionSym_Label *sl; 1546 PetscErrorCode ierr; 1547 1548 PetscFunctionBegin; 1549 ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr); 1550 sym->ops->getpoints = PetscSectionSymGetPoints_Label; 1551 sym->ops->view = PetscSectionSymView_Label; 1552 sym->ops->destroy = PetscSectionSymDestroy_Label; 1553 sym->data = (void *) sl; 1554 PetscFunctionReturn(0); 1555 } 1556 1557 /*@ 1558 PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label 1559 1560 Collective 1561 1562 Input Parameters: 1563 + comm - the MPI communicator for the new symmetry 1564 - label - the label defining the strata 1565 1566 Output Parameters: 1567 . sym - the section symmetries 1568 1569 Level: developer 1570 1571 .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms() 1572 @*/ 1573 PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym) 1574 { 1575 PetscErrorCode ierr; 1576 1577 PetscFunctionBegin; 1578 ierr = DMInitializePackage();CHKERRQ(ierr); 1579 ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr); 1580 ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr); 1581 ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr); 1582 PetscFunctionReturn(0); 1583 } 1584