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