1 #include <petsc/private/dmlabelimpl.h> /*I "petscdmlabel.h" I*/ 2 #include <petsc/private/isimpl.h> /*I "petscis.h" I*/ 3 #include <petscsf.h> 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "DMLabelCreate" 7 /*@C 8 DMLabelCreate - Create a DMLabel object, which is a multimap 9 10 Input parameter: 11 . name - The label name 12 13 Output parameter: 14 . label - The DMLabel 15 16 Level: beginner 17 18 .seealso: DMLabelDestroy() 19 @*/ 20 PetscErrorCode DMLabelCreate(const char name[], DMLabel *label) 21 { 22 PetscErrorCode ierr; 23 24 PetscFunctionBegin; 25 ierr = PetscNew(label);CHKERRQ(ierr); 26 ierr = PetscStrallocpy(name, &(*label)->name);CHKERRQ(ierr); 27 28 (*label)->refct = 1; 29 (*label)->state = -1; 30 (*label)->numStrata = 0; 31 (*label)->defaultValue = -1; 32 (*label)->stratumValues = NULL; 33 (*label)->arrayValid = NULL; 34 (*label)->stratumSizes = NULL; 35 (*label)->points = NULL; 36 (*label)->ht = NULL; 37 (*label)->pStart = -1; 38 (*label)->pEnd = -1; 39 (*label)->bt = NULL; 40 PetscFunctionReturn(0); 41 } 42 43 #undef __FUNCT__ 44 #define __FUNCT__ "DMLabelMakeValid_Private" 45 /* 46 DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format 47 48 Input parameter: 49 + label - The DMLabel 50 - v - The stratum value 51 52 Output parameter: 53 . label - The DMLabel with stratum in sorted list format 54 55 Level: developer 56 57 .seealso: DMLabelCreate() 58 */ 59 static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v) 60 { 61 PetscInt off; 62 PetscErrorCode ierr; 63 64 if (label->arrayValid[v]) return 0; 65 if (v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeValid_Private\n", v); 66 PetscFunctionBegin; 67 PetscHashISize(label->ht[v], label->stratumSizes[v]); 68 69 ierr = PetscMalloc1(label->stratumSizes[v], &label->points[v]);CHKERRQ(ierr); 70 off = 0; 71 ierr = PetscHashIGetKeys(label->ht[v], &off, &(label->points[v][0]));CHKERRQ(ierr); 72 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]); 73 PetscHashIClear(label->ht[v]); 74 ierr = PetscSortInt(label->stratumSizes[v], label->points[v]);CHKERRQ(ierr); 75 if (label->bt) { 76 PetscInt p; 77 78 for (p = 0; p < label->stratumSizes[v]; ++p) { 79 const PetscInt point = label->points[v][p]; 80 81 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); 82 ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 83 } 84 } 85 label->arrayValid[v] = PETSC_TRUE; 86 ++label->state; 87 PetscFunctionReturn(0); 88 } 89 90 #undef __FUNCT__ 91 #define __FUNCT__ "DMLabelMakeAllValid_Private" 92 /* 93 DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format 94 95 Input parameter: 96 . label - The DMLabel 97 98 Output parameter: 99 . label - The DMLabel with all strata in sorted list format 100 101 Level: developer 102 103 .seealso: DMLabelCreate() 104 */ 105 static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label) 106 { 107 PetscInt v; 108 PetscErrorCode ierr; 109 110 PetscFunctionBegin; 111 for (v = 0; v < label->numStrata; v++){ 112 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 113 } 114 PetscFunctionReturn(0); 115 } 116 117 #undef __FUNCT__ 118 #define __FUNCT__ "DMLabelMakeInvalid_Private" 119 /* 120 DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format 121 122 Input parameter: 123 + label - The DMLabel 124 - v - The stratum value 125 126 Output parameter: 127 . label - The DMLabel with stratum in hash format 128 129 Level: developer 130 131 .seealso: DMLabelCreate() 132 */ 133 static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v) 134 { 135 PETSC_UNUSED PetscHashIIter ret, iter; 136 PetscInt p; 137 PetscErrorCode ierr; 138 139 PetscFunctionBegin; 140 if (!label->arrayValid[v]) PetscFunctionReturn(0); 141 for (p = 0; p < label->stratumSizes[v]; ++p) PetscHashIPut(label->ht[v], label->points[v][p], ret, iter); 142 ierr = PetscFree(label->points[v]);CHKERRQ(ierr); 143 label->arrayValid[v] = PETSC_FALSE; 144 PetscFunctionReturn(0); 145 } 146 147 #undef __FUNCT__ 148 #define __FUNCT__ "DMLabelGetState" 149 PetscErrorCode DMLabelGetState(DMLabel label, PetscObjectState *state) 150 { 151 PetscFunctionBegin; 152 PetscValidPointer(state, 2); 153 *state = label->state; 154 PetscFunctionReturn(0); 155 } 156 157 #undef __FUNCT__ 158 #define __FUNCT__ "DMLabelAddStratum" 159 PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value) 160 { 161 PetscInt v, *tmpV, *tmpS, **tmpP; 162 PetscHashI *tmpH; 163 PetscBool *tmpB; 164 PetscErrorCode ierr; 165 166 PetscFunctionBegin; 167 168 ierr = PetscMalloc1((label->numStrata+1), &tmpV);CHKERRQ(ierr); 169 ierr = PetscMalloc1((label->numStrata+1), &tmpS);CHKERRQ(ierr); 170 ierr = PetscMalloc1((label->numStrata+1), &tmpH);CHKERRQ(ierr); 171 ierr = PetscMalloc1((label->numStrata+1), &tmpP);CHKERRQ(ierr); 172 ierr = PetscMalloc1((label->numStrata+1), &tmpB);CHKERRQ(ierr); 173 for (v = 0; v < label->numStrata; ++v) { 174 tmpV[v] = label->stratumValues[v]; 175 tmpS[v] = label->stratumSizes[v]; 176 tmpH[v] = label->ht[v]; 177 tmpP[v] = label->points[v]; 178 tmpB[v] = label->arrayValid[v]; 179 } 180 tmpV[v] = value; 181 tmpS[v] = 0; 182 PetscHashICreate(tmpH[v]); 183 tmpP[v] = NULL; 184 tmpB[v] = PETSC_TRUE; 185 ++label->numStrata; 186 ierr = PetscFree(label->stratumValues);CHKERRQ(ierr); 187 ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr); 188 ierr = PetscFree(label->ht);CHKERRQ(ierr); 189 ierr = PetscFree(label->points);CHKERRQ(ierr); 190 ierr = PetscFree(label->arrayValid);CHKERRQ(ierr); 191 label->stratumValues = tmpV; 192 label->stratumSizes = tmpS; 193 label->ht = tmpH; 194 label->points = tmpP; 195 label->arrayValid = tmpB; 196 197 PetscFunctionReturn(0); 198 } 199 200 #undef __FUNCT__ 201 #define __FUNCT__ "DMLabelGetName" 202 /*@C 203 DMLabelGetName - Return the name of a DMLabel object 204 205 Input parameter: 206 . label - The DMLabel 207 208 Output parameter: 209 . name - The label name 210 211 Level: beginner 212 213 .seealso: DMLabelCreate() 214 @*/ 215 PetscErrorCode DMLabelGetName(DMLabel label, const char **name) 216 { 217 PetscFunctionBegin; 218 PetscValidPointer(name, 2); 219 *name = label->name; 220 PetscFunctionReturn(0); 221 } 222 223 #undef __FUNCT__ 224 #define __FUNCT__ "DMLabelView_Ascii" 225 static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer) 226 { 227 PetscInt v; 228 PetscMPIInt rank; 229 PetscErrorCode ierr; 230 231 PetscFunctionBegin; 232 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr); 233 ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 234 if (label) { 235 ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", label->name);CHKERRQ(ierr); 236 if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, " Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);} 237 for (v = 0; v < label->numStrata; ++v) { 238 const PetscInt value = label->stratumValues[v]; 239 PetscInt p; 240 241 for (p = 0; p < label->stratumSizes[v]; ++p) { 242 ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, label->points[v][p], value);CHKERRQ(ierr); 243 } 244 } 245 } 246 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 247 ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 248 PetscFunctionReturn(0); 249 } 250 251 #undef __FUNCT__ 252 #define __FUNCT__ "DMLabelView" 253 /*@C 254 DMLabelView - View the label 255 256 Input Parameters: 257 + label - The DMLabel 258 - viewer - The PetscViewer 259 260 Level: intermediate 261 262 .seealso: DMLabelCreate(), DMLabelDestroy() 263 @*/ 264 PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer) 265 { 266 PetscBool iascii; 267 PetscErrorCode ierr; 268 269 PetscFunctionBegin; 270 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 271 if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 272 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 273 if (iascii) { 274 ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr); 275 } 276 PetscFunctionReturn(0); 277 } 278 279 #undef __FUNCT__ 280 #define __FUNCT__ "DMLabelDestroy" 281 PetscErrorCode DMLabelDestroy(DMLabel *label) 282 { 283 PetscInt v; 284 PetscErrorCode ierr; 285 286 PetscFunctionBegin; 287 if (!(*label)) PetscFunctionReturn(0); 288 if (--(*label)->refct > 0) PetscFunctionReturn(0); 289 ierr = PetscFree((*label)->name);CHKERRQ(ierr); 290 ierr = PetscFree((*label)->stratumValues);CHKERRQ(ierr); 291 ierr = PetscFree((*label)->stratumSizes);CHKERRQ(ierr); 292 for (v = 0; v < (*label)->numStrata; ++v) {ierr = PetscFree((*label)->points[v]);CHKERRQ(ierr);} 293 ierr = PetscFree((*label)->points);CHKERRQ(ierr); 294 ierr = PetscFree((*label)->arrayValid);CHKERRQ(ierr); 295 if ((*label)->ht) { 296 for (v = 0; v < (*label)->numStrata; ++v) {PetscHashIDestroy((*label)->ht[v]);} 297 ierr = PetscFree((*label)->ht);CHKERRQ(ierr); 298 } 299 ierr = PetscBTDestroy(&(*label)->bt);CHKERRQ(ierr); 300 ierr = PetscFree(*label);CHKERRQ(ierr); 301 PetscFunctionReturn(0); 302 } 303 304 #undef __FUNCT__ 305 #define __FUNCT__ "DMLabelDuplicate" 306 PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew) 307 { 308 PetscInt v, q; 309 PetscErrorCode ierr; 310 311 PetscFunctionBegin; 312 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 313 ierr = PetscNew(labelnew);CHKERRQ(ierr); 314 ierr = PetscStrallocpy(label->name, &(*labelnew)->name);CHKERRQ(ierr); 315 316 (*labelnew)->refct = 1; 317 (*labelnew)->numStrata = label->numStrata; 318 (*labelnew)->defaultValue = label->defaultValue; 319 if (label->numStrata) { 320 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr); 321 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr); 322 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr); 323 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr); 324 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->arrayValid);CHKERRQ(ierr); 325 /* Could eliminate unused space here */ 326 for (v = 0; v < label->numStrata; ++v) { 327 ierr = PetscMalloc1(label->stratumSizes[v], &(*labelnew)->points[v]);CHKERRQ(ierr); 328 PetscHashICreate((*labelnew)->ht[v]); 329 (*labelnew)->arrayValid[v] = PETSC_TRUE; 330 (*labelnew)->stratumValues[v] = label->stratumValues[v]; 331 (*labelnew)->stratumSizes[v] = label->stratumSizes[v]; 332 for (q = 0; q < label->stratumSizes[v]; ++q) { 333 (*labelnew)->points[v][q] = label->points[v][q]; 334 } 335 } 336 } 337 (*labelnew)->pStart = -1; 338 (*labelnew)->pEnd = -1; 339 (*labelnew)->bt = NULL; 340 PetscFunctionReturn(0); 341 } 342 343 #undef __FUNCT__ 344 #define __FUNCT__ "DMLabelCreateIndex" 345 /* This can be hooked into SetValue(), ClearValue(), etc. for updating */ 346 PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd) 347 { 348 PetscInt v; 349 PetscErrorCode ierr; 350 351 PetscFunctionBegin; 352 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 353 if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 354 label->pStart = pStart; 355 label->pEnd = pEnd; 356 ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr); 357 ierr = PetscBTMemzero(pEnd - pStart, label->bt);CHKERRQ(ierr); 358 for (v = 0; v < label->numStrata; ++v) { 359 PetscInt i; 360 361 for (i = 0; i < label->stratumSizes[v]; ++i) { 362 const PetscInt point = label->points[v][i]; 363 364 if ((point < pStart) || (point >= pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, pStart, pEnd); 365 ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr); 366 } 367 } 368 PetscFunctionReturn(0); 369 } 370 371 #undef __FUNCT__ 372 #define __FUNCT__ "DMLabelDestroyIndex" 373 PetscErrorCode DMLabelDestroyIndex(DMLabel label) 374 { 375 PetscErrorCode ierr; 376 377 PetscFunctionBegin; 378 label->pStart = -1; 379 label->pEnd = -1; 380 if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 381 PetscFunctionReturn(0); 382 } 383 384 #undef __FUNCT__ 385 #define __FUNCT__ "DMLabelHasValue" 386 /*@ 387 DMLabelHasValue - Determine whether a label assigns the value to any point 388 389 Input Parameters: 390 + label - the DMLabel 391 - value - the value 392 393 Output Parameter: 394 . contains - Flag indicating whether the label maps this value to any point 395 396 Level: developer 397 398 .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue() 399 @*/ 400 PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains) 401 { 402 PetscInt v; 403 404 PetscFunctionBegin; 405 PetscValidPointer(contains, 3); 406 for (v = 0; v < label->numStrata; ++v) { 407 if (value == label->stratumValues[v]) break; 408 } 409 *contains = (v < label->numStrata ? PETSC_TRUE : PETSC_FALSE); 410 PetscFunctionReturn(0); 411 } 412 413 #undef __FUNCT__ 414 #define __FUNCT__ "DMLabelHasPoint" 415 /*@ 416 DMLabelHasPoint - Determine whether a label assigns a value to a point 417 418 Input Parameters: 419 + label - the DMLabel 420 - point - the point 421 422 Output Parameter: 423 . contains - Flag indicating whether the label maps this point to a value 424 425 Note: The user must call DMLabelCreateIndex() before this function. 426 427 Level: developer 428 429 .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue() 430 @*/ 431 PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains) 432 { 433 PetscErrorCode ierr; 434 435 PetscFunctionBeginHot; 436 PetscValidPointer(contains, 3); 437 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 438 #if defined(PETSC_USE_DEBUG) 439 if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()"); 440 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); 441 #endif 442 *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE; 443 PetscFunctionReturn(0); 444 } 445 446 #undef __FUNCT__ 447 #define __FUNCT__ "DMLabelStratumHasPoint" 448 /*@ 449 DMLabelStratumHasPoint - Return true if the stratum contains a point 450 451 Input Parameters: 452 + label - the DMLabel 453 . value - the stratum value 454 - point - the point 455 456 Output Parameter: 457 . contains - true if the stratum contains the point 458 459 Level: intermediate 460 461 .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue() 462 @*/ 463 PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains) 464 { 465 PetscInt v; 466 PetscErrorCode ierr; 467 468 PetscFunctionBegin; 469 PetscValidPointer(contains, 4); 470 *contains = PETSC_FALSE; 471 for (v = 0; v < label->numStrata; ++v) { 472 if (label->stratumValues[v] == value) { 473 if (label->arrayValid[v]) { 474 PetscInt i; 475 476 ierr = PetscFindInt(point, label->stratumSizes[v], &label->points[v][0], &i);CHKERRQ(ierr); 477 if (i >= 0) { 478 *contains = PETSC_TRUE; 479 break; 480 } 481 } else { 482 PetscBool has; 483 484 PetscHashIHasKey(label->ht[v], point, has); 485 if (has) { 486 *contains = PETSC_TRUE; 487 break; 488 } 489 } 490 } 491 } 492 PetscFunctionReturn(0); 493 } 494 495 #undef __FUNCT__ 496 #define __FUNCT__ "DMLabelGetDefaultValue" 497 /* 498 DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 499 When a label is created, it is initialized to -1. 500 501 Input parameter: 502 . label - a DMLabel object 503 504 Output parameter: 505 . defaultValue - the default value 506 507 Level: beginner 508 509 .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 510 */ 511 PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue) 512 { 513 PetscFunctionBegin; 514 *defaultValue = label->defaultValue; 515 PetscFunctionReturn(0); 516 } 517 518 #undef __FUNCT__ 519 #define __FUNCT__ "DMLabelSetDefaultValue" 520 /* 521 DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 522 When a label is created, it is initialized to -1. 523 524 Input parameter: 525 . label - a DMLabel object 526 527 Output parameter: 528 . defaultValue - the default value 529 530 Level: beginner 531 532 .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 533 */ 534 PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue) 535 { 536 PetscFunctionBegin; 537 label->defaultValue = defaultValue; 538 PetscFunctionReturn(0); 539 } 540 541 #undef __FUNCT__ 542 #define __FUNCT__ "DMLabelGetValue" 543 /*@ 544 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()) 545 546 Input Parameters: 547 + label - the DMLabel 548 - point - the point 549 550 Output Parameter: 551 . value - The point value, or -1 552 553 Level: intermediate 554 555 .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 556 @*/ 557 PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value) 558 { 559 PetscInt v; 560 PetscErrorCode ierr; 561 562 PetscFunctionBegin; 563 PetscValidPointer(value, 3); 564 *value = label->defaultValue; 565 for (v = 0; v < label->numStrata; ++v) { 566 if (label->arrayValid[v]) { 567 PetscInt i; 568 569 ierr = PetscFindInt(point, label->stratumSizes[v], &label->points[v][0], &i);CHKERRQ(ierr); 570 if (i >= 0) { 571 *value = label->stratumValues[v]; 572 break; 573 } 574 } else { 575 PetscBool has; 576 577 PetscHashIHasKey(label->ht[v], point, has); 578 if (has) { 579 *value = label->stratumValues[v]; 580 break; 581 } 582 } 583 } 584 PetscFunctionReturn(0); 585 } 586 587 #undef __FUNCT__ 588 #define __FUNCT__ "DMLabelSetValue" 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 #undef __FUNCT__ 622 #define __FUNCT__ "DMLabelClearValue" 623 /*@ 624 DMLabelClearValue - Clear the value a label assigns to a point 625 626 Input Parameters: 627 + label - the DMLabel 628 . point - the point 629 - value - The point value 630 631 Level: intermediate 632 633 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue() 634 @*/ 635 PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value) 636 { 637 PetscInt v, p; 638 PetscErrorCode ierr; 639 640 PetscFunctionBegin; 641 /* Find label value */ 642 for (v = 0; v < label->numStrata; ++v) { 643 if (label->stratumValues[v] == value) break; 644 } 645 if (v >= label->numStrata) PetscFunctionReturn(0); 646 if (label->arrayValid[v]) { 647 /* Check whether point exists */ 648 ierr = PetscFindInt(point, label->stratumSizes[v], &label->points[v][0], &p);CHKERRQ(ierr); 649 if (p >= 0) { 650 ierr = PetscMemmove(&label->points[v][p], &label->points[v][p+1], (label->stratumSizes[v]-p-1) * sizeof(PetscInt));CHKERRQ(ierr); 651 --label->stratumSizes[v]; 652 if (label->bt) { 653 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); 654 ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 655 } 656 } 657 } else { 658 ierr = PetscHashIDelKey(label->ht[v], point);CHKERRQ(ierr); 659 } 660 PetscFunctionReturn(0); 661 } 662 663 #undef __FUNCT__ 664 #define __FUNCT__ "DMLabelInsertIS" 665 /*@ 666 DMLabelInsertIS - Set all points in the IS to a value 667 668 Input Parameters: 669 + label - the DMLabel 670 . is - the point IS 671 - value - The point value 672 673 Level: intermediate 674 675 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 676 @*/ 677 PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value) 678 { 679 const PetscInt *points; 680 PetscInt n, p; 681 PetscErrorCode ierr; 682 683 PetscFunctionBegin; 684 PetscValidHeaderSpecific(is, IS_CLASSID, 2); 685 ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 686 ierr = ISGetIndices(is, &points);CHKERRQ(ierr); 687 for (p = 0; p < n; ++p) {ierr = DMLabelSetValue(label, points[p], value);CHKERRQ(ierr);} 688 ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr); 689 PetscFunctionReturn(0); 690 } 691 692 #undef __FUNCT__ 693 #define __FUNCT__ "DMLabelGetNumValues" 694 PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues) 695 { 696 PetscFunctionBegin; 697 PetscValidPointer(numValues, 2); 698 *numValues = label->numStrata; 699 PetscFunctionReturn(0); 700 } 701 702 #undef __FUNCT__ 703 #define __FUNCT__ "DMLabelGetValueIS" 704 PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values) 705 { 706 PetscErrorCode ierr; 707 708 PetscFunctionBegin; 709 PetscValidPointer(values, 2); 710 ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr); 711 PetscFunctionReturn(0); 712 } 713 714 #undef __FUNCT__ 715 #define __FUNCT__ "DMLabelGetStratumSize" 716 PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size) 717 { 718 PetscInt v; 719 PetscErrorCode ierr; 720 721 PetscFunctionBegin; 722 PetscValidPointer(size, 3); 723 *size = 0; 724 for (v = 0; v < label->numStrata; ++v) { 725 if (label->stratumValues[v] == value) { 726 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 727 *size = label->stratumSizes[v]; 728 break; 729 } 730 } 731 PetscFunctionReturn(0); 732 } 733 734 #undef __FUNCT__ 735 #define __FUNCT__ "DMLabelGetStratumBounds" 736 PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end) 737 { 738 PetscInt v; 739 PetscErrorCode ierr; 740 741 PetscFunctionBegin; 742 if (start) {PetscValidPointer(start, 3); *start = 0;} 743 if (end) {PetscValidPointer(end, 4); *end = 0;} 744 for (v = 0; v < label->numStrata; ++v) { 745 if (label->stratumValues[v] != value) continue; 746 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 747 if (label->stratumSizes[v] <= 0) break; 748 if (start) *start = label->points[v][0]; 749 if (end) *end = label->points[v][label->stratumSizes[v]-1]+1; 750 break; 751 } 752 PetscFunctionReturn(0); 753 } 754 755 #undef __FUNCT__ 756 #define __FUNCT__ "DMLabelGetStratumIS" 757 PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points) 758 { 759 PetscInt v; 760 PetscErrorCode ierr; 761 762 PetscFunctionBegin; 763 PetscValidPointer(points, 3); 764 *points = NULL; 765 for (v = 0; v < label->numStrata; ++v) { 766 if (label->stratumValues[v] == value) { 767 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 768 if (label->arrayValid[v]) { 769 ierr = ISCreateGeneral(PETSC_COMM_SELF, label->stratumSizes[v], &label->points[v][0], PETSC_COPY_VALUES, points);CHKERRQ(ierr); 770 ierr = PetscObjectSetName((PetscObject) *points, "indices");CHKERRQ(ierr); 771 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to implement this to speedup Stratify"); 772 break; 773 } 774 } 775 PetscFunctionReturn(0); 776 } 777 778 #undef __FUNCT__ 779 #define __FUNCT__ "DMLabelClearStratum" 780 PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value) 781 { 782 PetscInt v; 783 PetscErrorCode ierr; 784 785 PetscFunctionBegin; 786 for (v = 0; v < label->numStrata; ++v) { 787 if (label->stratumValues[v] == value) break; 788 } 789 if (v >= label->numStrata) PetscFunctionReturn(0); 790 if (label->bt) { 791 PetscInt i; 792 793 for (i = 0; i < label->stratumSizes[v]; ++i) { 794 const PetscInt point = label->points[v][i]; 795 796 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); 797 ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 798 } 799 } 800 if (label->arrayValid[v]) { 801 label->stratumSizes[v] = 0; 802 } else { 803 PetscHashIClear(label->ht[v]); 804 } 805 PetscFunctionReturn(0); 806 } 807 808 #undef __FUNCT__ 809 #define __FUNCT__ "DMLabelFilter" 810 PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end) 811 { 812 PetscInt v; 813 PetscErrorCode ierr; 814 815 PetscFunctionBegin; 816 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 817 label->pStart = start; 818 label->pEnd = end; 819 if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 820 /* Could squish offsets, but would only make sense if I reallocate the storage */ 821 for (v = 0; v < label->numStrata; ++v) { 822 PetscInt off, q; 823 824 for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) { 825 const PetscInt point = label->points[v][q]; 826 827 if ((point < start) || (point >= end)) continue; 828 label->points[v][off++] = point; 829 } 830 label->stratumSizes[v] = off; 831 } 832 ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr); 833 PetscFunctionReturn(0); 834 } 835 836 #undef __FUNCT__ 837 #define __FUNCT__ "DMLabelPermute" 838 PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew) 839 { 840 const PetscInt *perm; 841 PetscInt numValues, numPoints, v, q; 842 PetscErrorCode ierr; 843 844 PetscFunctionBegin; 845 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 846 ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr); 847 ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr); 848 ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr); 849 ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr); 850 for (v = 0; v < numValues; ++v) { 851 const PetscInt size = (*labelNew)->stratumSizes[v]; 852 853 for (q = 0; q < size; ++q) { 854 const PetscInt point = (*labelNew)->points[v][q]; 855 856 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); 857 (*labelNew)->points[v][q] = perm[point]; 858 } 859 ierr = PetscSortInt(size, &(*labelNew)->points[v][0]);CHKERRQ(ierr); 860 } 861 ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr); 862 if (label->bt) { 863 ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 864 ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr); 865 } 866 PetscFunctionReturn(0); 867 } 868 869 #undef __FUNCT__ 870 #define __FUNCT__ "DMLabelDistribute_Internal" 871 PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata) 872 { 873 MPI_Comm comm; 874 PetscInt s, l, nroots, nleaves, dof, offset, size; 875 PetscInt *remoteOffsets, *rootStrata, *rootIdx; 876 PetscSection rootSection; 877 PetscSF labelSF; 878 PetscErrorCode ierr; 879 880 PetscFunctionBegin; 881 if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 882 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 883 /* Build a section of stratum values per point, generate the according SF 884 and distribute point-wise stratum values to leaves. */ 885 ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 886 ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 887 ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr); 888 if (label) { 889 for (s = 0; s < label->numStrata; ++s) { 890 for (l = 0; l < label->stratumSizes[s]; l++) { 891 ierr = PetscSectionGetDof(rootSection, label->points[s][l], &dof);CHKERRQ(ierr); 892 ierr = PetscSectionSetDof(rootSection, label->points[s][l], dof+1);CHKERRQ(ierr); 893 } 894 } 895 } 896 ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 897 /* Create a point-wise array of stratum values */ 898 ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr); 899 ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr); 900 ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr); 901 if (label) { 902 for (s = 0; s < label->numStrata; ++s) { 903 for (l = 0; l < label->stratumSizes[s]; l++) { 904 const PetscInt p = label->points[s][l]; 905 ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr); 906 rootStrata[offset+rootIdx[p]++] = label->stratumValues[s]; 907 } 908 } 909 } 910 /* Build SF that maps label points to remote processes */ 911 ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr); 912 ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr); 913 ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr); 914 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 915 /* Send the strata for each point over the derived SF */ 916 ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr); 917 ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr); 918 ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 919 ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 920 /* Clean up */ 921 ierr = PetscFree(rootStrata);CHKERRQ(ierr); 922 ierr = PetscFree(rootIdx);CHKERRQ(ierr); 923 ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 924 ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr); 925 PetscFunctionReturn(0); 926 } 927 928 #undef __FUNCT__ 929 #define __FUNCT__ "DMLabelDistribute" 930 PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew) 931 { 932 MPI_Comm comm; 933 PetscSection leafSection; 934 PetscInt p, pStart, pEnd, s, size, dof, offset, stratum; 935 PetscInt *leafStrata, *strataIdx; 936 char *name; 937 PetscInt nameSize; 938 PetscHashI stratumHash; 939 PETSC_UNUSED PetscHashIIter ret, iter; 940 size_t len = 0; 941 PetscMPIInt rank; 942 PetscErrorCode ierr; 943 944 PetscFunctionBegin; 945 if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 946 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 947 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 948 /* Bcast name */ 949 if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 950 nameSize = len; 951 ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 952 ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 953 if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 954 ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 955 ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 956 ierr = PetscFree(name);CHKERRQ(ierr); 957 /* Bcast defaultValue */ 958 if (!rank) (*labelNew)->defaultValue = label->defaultValue; 959 ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 960 /* Distribute stratum values over the SF and get the point mapping on the receiver */ 961 ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr); 962 /* Determine received stratum values and initialise new label*/ 963 PetscHashICreate(stratumHash); 964 ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr); 965 for (p = 0; p < size; ++p) PetscHashIPut(stratumHash, leafStrata[p], ret, iter); 966 PetscHashISize(stratumHash, (*labelNew)->numStrata); 967 ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->arrayValid);CHKERRQ(ierr); 968 for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->arrayValid[s] = PETSC_TRUE; 969 ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr); 970 /* Turn leafStrata into indices rather than stratum values */ 971 offset = 0; 972 ierr = PetscHashIGetKeys(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr); 973 for (s = 0; s < (*labelNew)->numStrata; ++s) { 974 for (p = 0; p < size; ++p) { 975 if (leafStrata[p] == (*labelNew)->stratumValues[s]) leafStrata[p] = s; 976 } 977 } 978 /* Rebuild the point strata on the receiver */ 979 ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr); 980 ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 981 for (p=pStart; p<pEnd; p++) { 982 ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 983 ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 984 for (s=0; s<dof; s++) { 985 (*labelNew)->stratumSizes[leafStrata[offset+s]]++; 986 } 987 } 988 ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr); 989 ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr); 990 for (s = 0; s < (*labelNew)->numStrata; ++s) { 991 PetscHashICreate((*labelNew)->ht[s]); 992 ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(*labelNew)->points[s]);CHKERRQ(ierr); 993 } 994 /* Insert points into new strata */ 995 ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr); 996 ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 997 for (p=pStart; p<pEnd; p++) { 998 ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 999 ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1000 for (s=0; s<dof; s++) { 1001 stratum = leafStrata[offset+s]; 1002 (*labelNew)->points[stratum][strataIdx[stratum]++] = p; 1003 } 1004 } 1005 PetscHashIDestroy(stratumHash); 1006 ierr = PetscFree(leafStrata);CHKERRQ(ierr); 1007 ierr = PetscFree(strataIdx);CHKERRQ(ierr); 1008 ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 1009 PetscFunctionReturn(0); 1010 } 1011 1012 #undef __FUNCT__ 1013 #define __FUNCT__ "DMLabelGather" 1014 /*@ 1015 DMLabelGather - Gather all label values from leafs into roots 1016 1017 Input Parameters: 1018 + label - the DMLabel 1019 . point - the Star Forest point communication map 1020 1021 Input Parameters: 1022 + label - the new DMLabel with localised leaf values 1023 1024 Level: developer 1025 1026 Note: This is the inverse operation to DMLabelDistribute. 1027 1028 .seealso: DMLabelDistribute() 1029 @*/ 1030 PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew) 1031 { 1032 MPI_Comm comm; 1033 PetscSection rootSection; 1034 PetscSF sfLabel; 1035 PetscSFNode *rootPoints, *leafPoints; 1036 PetscInt p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset; 1037 const PetscInt *rootDegree, *ilocal; 1038 PetscInt *rootStrata; 1039 char *name; 1040 PetscInt nameSize; 1041 size_t len = 0; 1042 PetscMPIInt rank, numProcs; 1043 PetscErrorCode ierr; 1044 1045 PetscFunctionBegin; 1046 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 1047 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1048 ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 1049 /* Bcast name */ 1050 if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 1051 nameSize = len; 1052 ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1053 ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1054 if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 1055 ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1056 ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 1057 ierr = PetscFree(name);CHKERRQ(ierr); 1058 /* Gather rank/index pairs of leaves into local roots to build 1059 an inverse, multi-rooted SF. Note that this ignores local leaf 1060 indexing due to the use of the multiSF in PetscSFGather. */ 1061 ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr); 1062 ierr = PetscMalloc1(nleaves, &leafPoints);CHKERRQ(ierr); 1063 for (p = 0; p < nleaves; p++) { 1064 leafPoints[p].index = ilocal[p]; 1065 leafPoints[p].rank = rank; 1066 } 1067 ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr); 1068 ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr); 1069 for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p]; 1070 ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr); 1071 ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 1072 ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 1073 ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr); 1074 ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 1075 /* Migrate label over inverted SF to pull stratum values at leaves into roots. */ 1076 ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr); 1077 /* Rebuild the point strata on the receiver */ 1078 for (p = 0, idx = 0; p < nroots; p++) { 1079 for (d = 0; d < rootDegree[p]; d++) { 1080 ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr); 1081 ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr); 1082 for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);} 1083 } 1084 idx += rootDegree[p]; 1085 } 1086 ierr = PetscFree(leafPoints);CHKERRQ(ierr); 1087 ierr = PetscFree(rootStrata);CHKERRQ(ierr); 1088 ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 1089 ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr); 1090 PetscFunctionReturn(0); 1091 } 1092 1093 #undef __FUNCT__ 1094 #define __FUNCT__ "DMLabelConvertToSection" 1095 PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is) 1096 { 1097 IS vIS; 1098 const PetscInt *values; 1099 PetscInt *points; 1100 PetscInt nV, vS = 0, vE = 0, v, N; 1101 PetscErrorCode ierr; 1102 1103 PetscFunctionBegin; 1104 ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr); 1105 ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr); 1106 ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr); 1107 if (nV) {vS = values[0]; vE = values[0]+1;} 1108 for (v = 1; v < nV; ++v) { 1109 vS = PetscMin(vS, values[v]); 1110 vE = PetscMax(vE, values[v]+1); 1111 } 1112 ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr); 1113 ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr); 1114 for (v = 0; v < nV; ++v) { 1115 PetscInt n; 1116 1117 ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr); 1118 ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr); 1119 } 1120 ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 1121 ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr); 1122 ierr = PetscMalloc1(N, &points);CHKERRQ(ierr); 1123 for (v = 0; v < nV; ++v) { 1124 IS is; 1125 const PetscInt *spoints; 1126 PetscInt dof, off, p; 1127 1128 ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr); 1129 ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr); 1130 ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr); 1131 ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr); 1132 for (p = 0; p < dof; ++p) points[off+p] = spoints[p]; 1133 ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr); 1134 ierr = ISDestroy(&is);CHKERRQ(ierr); 1135 } 1136 ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr); 1137 ierr = ISDestroy(&vIS);CHKERRQ(ierr); 1138 ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr); 1139 PetscFunctionReturn(0); 1140 } 1141 1142 #undef __FUNCT__ 1143 #define __FUNCT__ "PetscSectionCreateGlobalSectionLabel" 1144 /*@C 1145 PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using 1146 the local section and an SF describing the section point overlap. 1147 1148 Input Parameters: 1149 + s - The PetscSection for the local field layout 1150 . sf - The SF describing parallel layout of the section points 1151 . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs 1152 . label - The label specifying the points 1153 - labelValue - The label stratum specifying the points 1154 1155 Output Parameter: 1156 . gsection - The PetscSection for the global field layout 1157 1158 Note: This gives negative sizes and offsets to points not owned by this process 1159 1160 Level: developer 1161 1162 .seealso: PetscSectionCreate() 1163 @*/ 1164 PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection) 1165 { 1166 PetscInt *neg = NULL, *tmpOff = NULL; 1167 PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots; 1168 PetscErrorCode ierr; 1169 1170 PetscFunctionBegin; 1171 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr); 1172 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 1173 ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr); 1174 ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 1175 if (nroots >= 0) { 1176 if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart); 1177 ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr); 1178 if (nroots > pEnd-pStart) { 1179 ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr); 1180 } else { 1181 tmpOff = &(*gsection)->atlasDof[-pStart]; 1182 } 1183 } 1184 /* Mark ghost points with negative dof */ 1185 for (p = pStart; p < pEnd; ++p) { 1186 PetscInt value; 1187 1188 ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr); 1189 if (value != labelValue) continue; 1190 ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr); 1191 ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr); 1192 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 1193 if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);} 1194 if (neg) neg[p] = -(dof+1); 1195 } 1196 ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr); 1197 if (nroots >= 0) { 1198 ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1199 ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1200 if (nroots > pEnd-pStart) { 1201 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];} 1202 } 1203 } 1204 /* Calculate new sizes, get proccess offset, and calculate point offsets */ 1205 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1206 cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0; 1207 (*gsection)->atlasOff[p] = off; 1208 off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0; 1209 } 1210 ierr = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr); 1211 globalOff -= off; 1212 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1213 (*gsection)->atlasOff[p] += globalOff; 1214 if (neg) neg[p] = -((*gsection)->atlasOff[p]+1); 1215 } 1216 /* Put in negative offsets for ghost points */ 1217 if (nroots >= 0) { 1218 ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1219 ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1220 if (nroots > pEnd-pStart) { 1221 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];} 1222 } 1223 } 1224 if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);} 1225 ierr = PetscFree(neg);CHKERRQ(ierr); 1226 PetscFunctionReturn(0); 1227 } 1228 1229