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