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 PetscInt i; 491 492 ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr); 493 if (i >= 0) { 494 *contains = PETSC_TRUE; 495 break; 496 } 497 } else { 498 PetscBool has; 499 500 PetscHashIHasKey(label->ht[v], point, has); 501 if (has) { 502 *contains = PETSC_TRUE; 503 break; 504 } 505 } 506 } 507 } 508 PetscFunctionReturn(0); 509 } 510 511 #undef __FUNCT__ 512 #define __FUNCT__ "DMLabelGetDefaultValue" 513 /* 514 DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 515 When a label is created, it is initialized to -1. 516 517 Input parameter: 518 . label - a DMLabel object 519 520 Output parameter: 521 . defaultValue - the default value 522 523 Level: beginner 524 525 .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 526 */ 527 PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue) 528 { 529 PetscFunctionBegin; 530 *defaultValue = label->defaultValue; 531 PetscFunctionReturn(0); 532 } 533 534 #undef __FUNCT__ 535 #define __FUNCT__ "DMLabelSetDefaultValue" 536 /* 537 DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 538 When a label is created, it is initialized to -1. 539 540 Input parameter: 541 . label - a DMLabel object 542 543 Output parameter: 544 . defaultValue - the default value 545 546 Level: beginner 547 548 .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 549 */ 550 PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue) 551 { 552 PetscFunctionBegin; 553 label->defaultValue = defaultValue; 554 PetscFunctionReturn(0); 555 } 556 557 #undef __FUNCT__ 558 #define __FUNCT__ "DMLabelGetValue" 559 /*@ 560 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()) 561 562 Input Parameters: 563 + label - the DMLabel 564 - point - the point 565 566 Output Parameter: 567 . value - The point value, or -1 568 569 Level: intermediate 570 571 .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 572 @*/ 573 PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value) 574 { 575 PetscInt v; 576 PetscErrorCode ierr; 577 578 PetscFunctionBegin; 579 PetscValidPointer(value, 3); 580 *value = label->defaultValue; 581 for (v = 0; v < label->numStrata; ++v) { 582 if (label->validIS[v]) { 583 PetscInt i; 584 585 ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr); 586 if (i >= 0) { 587 *value = label->stratumValues[v]; 588 break; 589 } 590 } else { 591 PetscBool has; 592 593 PetscHashIHasKey(label->ht[v], point, has); 594 if (has) { 595 *value = label->stratumValues[v]; 596 break; 597 } 598 } 599 } 600 PetscFunctionReturn(0); 601 } 602 603 #undef __FUNCT__ 604 #define __FUNCT__ "DMLabelSetValue" 605 /*@ 606 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. 607 608 Input Parameters: 609 + label - the DMLabel 610 . point - the point 611 - value - The point value 612 613 Level: intermediate 614 615 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 616 @*/ 617 PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value) 618 { 619 PETSC_UNUSED PetscHashIIter iter, ret; 620 PetscInt v; 621 PetscErrorCode ierr; 622 623 PetscFunctionBegin; 624 /* Find, or add, label value */ 625 if (value == label->defaultValue) PetscFunctionReturn(0); 626 for (v = 0; v < label->numStrata; ++v) { 627 if (label->stratumValues[v] == value) break; 628 } 629 /* Create new table */ 630 if (v >= label->numStrata) {ierr = DMLabelAddStratum(label, value);CHKERRQ(ierr);} 631 ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 632 /* Set key */ 633 PetscHashIPut(label->ht[v], point, ret, iter); 634 PetscFunctionReturn(0); 635 } 636 637 #undef __FUNCT__ 638 #define __FUNCT__ "DMLabelClearValue" 639 /*@ 640 DMLabelClearValue - Clear the value a label assigns to a point 641 642 Input Parameters: 643 + label - the DMLabel 644 . point - the point 645 - value - The point value 646 647 Level: intermediate 648 649 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue() 650 @*/ 651 PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value) 652 { 653 PetscInt v; 654 PetscErrorCode ierr; 655 656 PetscFunctionBegin; 657 /* Find label value */ 658 for (v = 0; v < label->numStrata; ++v) { 659 if (label->stratumValues[v] == value) break; 660 } 661 if (v >= label->numStrata) PetscFunctionReturn(0); 662 if (label->validIS[v]) { 663 ierr = DMLabelMakeInvalid_Private(label,v);CHKERRQ(ierr); 664 } 665 if (label->bt) { 666 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); 667 ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 668 } 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__ "DMLabelSetStratumIS" 812 PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is) 813 { 814 PetscInt v, numStrata; 815 PetscErrorCode ierr; 816 817 PetscFunctionBegin; 818 numStrata = label->numStrata; 819 for (v = 0; v < numStrata; v++) { 820 if (label->stratumValues[v] == value) break; 821 } 822 if (v >= numStrata) {ierr = DMLabelAddStratum(label,value);CHKERRQ(ierr);} 823 if (is == label->points[v]) PetscFunctionReturn(0); 824 ierr = DMLabelClearStratum(label,value);CHKERRQ(ierr); 825 ierr = ISGetLocalSize(is,&(label->stratumSizes[v]));CHKERRQ(ierr); 826 label->stratumValues[v] = value; 827 label->validIS[v] = PETSC_TRUE; 828 ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 829 ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 830 if (label->bt) { 831 const PetscInt *points; 832 PetscInt p; 833 834 ierr = ISGetIndices(is,&points);CHKERRQ(ierr); 835 for (p = 0; p < label->stratumSizes[v]; ++p) { 836 const PetscInt point = points[p]; 837 838 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); 839 ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 840 } 841 } 842 label->points[v] = is; 843 PetscFunctionReturn(0); 844 } 845 846 847 #undef __FUNCT__ 848 #define __FUNCT__ "DMLabelClearStratum" 849 PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value) 850 { 851 PetscInt v; 852 PetscErrorCode ierr; 853 854 PetscFunctionBegin; 855 for (v = 0; v < label->numStrata; ++v) { 856 if (label->stratumValues[v] == value) break; 857 } 858 if (v >= label->numStrata) PetscFunctionReturn(0); 859 if (label->validIS[v]) { 860 if (label->bt) { 861 PetscInt i; 862 const PetscInt *points; 863 864 ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 865 for (i = 0; i < label->stratumSizes[v]; ++i) { 866 const PetscInt point = points[i]; 867 868 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); 869 ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 870 } 871 ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr); 872 } 873 ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 874 label->stratumSizes[v] = 0; 875 ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr); 876 ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr); 877 } else { 878 PetscHashIClear(label->ht[v]); 879 } 880 PetscFunctionReturn(0); 881 } 882 883 #undef __FUNCT__ 884 #define __FUNCT__ "DMLabelFilter" 885 PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end) 886 { 887 PetscInt v; 888 PetscErrorCode ierr; 889 890 PetscFunctionBegin; 891 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 892 label->pStart = start; 893 label->pEnd = end; 894 if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 895 /* Could squish offsets, but would only make sense if I reallocate the storage */ 896 for (v = 0; v < label->numStrata; ++v) { 897 PetscInt off, q; 898 const PetscInt *points; 899 PetscInt *pointsNew = NULL; 900 901 ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 902 for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) { 903 const PetscInt point = points[q]; 904 905 if ((point < start) || (point >= end)) { 906 if (!pointsNew) { 907 ierr = PetscMalloc1(label->stratumSizes[v],&pointsNew);CHKERRQ(ierr); 908 ierr = PetscMemcpy(pointsNew,points,(size_t) off * sizeof(PetscInt));CHKERRQ(ierr); 909 } 910 continue; 911 } 912 if (pointsNew) { 913 pointsNew[off++] = point; 914 } 915 } 916 ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 917 if (pointsNew) { 918 ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 919 ierr = ISCreateGeneral(PETSC_COMM_SELF,off,pointsNew,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr); 920 ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr); 921 } 922 label->stratumSizes[v] = off; 923 } 924 ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr); 925 PetscFunctionReturn(0); 926 } 927 928 #undef __FUNCT__ 929 #define __FUNCT__ "DMLabelPermute" 930 PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew) 931 { 932 const PetscInt *perm; 933 PetscInt numValues, numPoints, v, q; 934 PetscErrorCode ierr; 935 936 PetscFunctionBegin; 937 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 938 ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr); 939 ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr); 940 ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr); 941 ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr); 942 for (v = 0; v < numValues; ++v) { 943 const PetscInt size = (*labelNew)->stratumSizes[v]; 944 const PetscInt *points; 945 PetscInt *pointsNew; 946 947 ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 948 ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr); 949 for (q = 0; q < size; ++q) { 950 const PetscInt point = points[q]; 951 952 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); 953 pointsNew[q] = perm[point]; 954 } 955 ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 956 ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr); 957 ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr); 958 ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr); 959 ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr); 960 } 961 ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr); 962 if (label->bt) { 963 ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 964 ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr); 965 } 966 PetscFunctionReturn(0); 967 } 968 969 #undef __FUNCT__ 970 #define __FUNCT__ "DMLabelDistribute_Internal" 971 PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata) 972 { 973 MPI_Comm comm; 974 PetscInt s, l, nroots, nleaves, dof, offset, size; 975 PetscInt *remoteOffsets, *rootStrata, *rootIdx; 976 PetscSection rootSection; 977 PetscSF labelSF; 978 PetscErrorCode ierr; 979 980 PetscFunctionBegin; 981 if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 982 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 983 /* Build a section of stratum values per point, generate the according SF 984 and distribute point-wise stratum values to leaves. */ 985 ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 986 ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 987 ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr); 988 if (label) { 989 for (s = 0; s < label->numStrata; ++s) { 990 const PetscInt *points; 991 992 ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 993 for (l = 0; l < label->stratumSizes[s]; l++) { 994 ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr); 995 ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr); 996 } 997 ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 998 } 999 } 1000 ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 1001 /* Create a point-wise array of stratum values */ 1002 ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr); 1003 ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr); 1004 ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr); 1005 if (label) { 1006 for (s = 0; s < label->numStrata; ++s) { 1007 const PetscInt *points; 1008 1009 ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 1010 for (l = 0; l < label->stratumSizes[s]; l++) { 1011 const PetscInt p = points[l]; 1012 ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr); 1013 rootStrata[offset+rootIdx[p]++] = label->stratumValues[s]; 1014 } 1015 ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 1016 } 1017 } 1018 /* Build SF that maps label points to remote processes */ 1019 ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr); 1020 ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr); 1021 ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr); 1022 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 1023 /* Send the strata for each point over the derived SF */ 1024 ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr); 1025 ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr); 1026 ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 1027 ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 1028 /* Clean up */ 1029 ierr = PetscFree(rootStrata);CHKERRQ(ierr); 1030 ierr = PetscFree(rootIdx);CHKERRQ(ierr); 1031 ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 1032 ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr); 1033 PetscFunctionReturn(0); 1034 } 1035 1036 #undef __FUNCT__ 1037 #define __FUNCT__ "DMLabelDistribute" 1038 PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew) 1039 { 1040 MPI_Comm comm; 1041 PetscSection leafSection; 1042 PetscInt p, pStart, pEnd, s, size, dof, offset, stratum; 1043 PetscInt *leafStrata, *strataIdx; 1044 PetscInt **points; 1045 char *name; 1046 PetscInt nameSize; 1047 PetscHashI stratumHash; 1048 PETSC_UNUSED PetscHashIIter ret, iter; 1049 size_t len = 0; 1050 PetscMPIInt rank; 1051 PetscErrorCode ierr; 1052 1053 PetscFunctionBegin; 1054 if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 1055 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 1056 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1057 /* Bcast name */ 1058 if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 1059 nameSize = len; 1060 ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1061 ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1062 if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 1063 ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1064 ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 1065 ierr = PetscFree(name);CHKERRQ(ierr); 1066 /* Bcast defaultValue */ 1067 if (!rank) (*labelNew)->defaultValue = label->defaultValue; 1068 ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1069 /* Distribute stratum values over the SF and get the point mapping on the receiver */ 1070 ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr); 1071 /* Determine received stratum values and initialise new label*/ 1072 PetscHashICreate(stratumHash); 1073 ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr); 1074 for (p = 0; p < size; ++p) PetscHashIPut(stratumHash, leafStrata[p], ret, iter); 1075 PetscHashISize(stratumHash, (*labelNew)->numStrata); 1076 ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr); 1077 for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE; 1078 ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr); 1079 /* Turn leafStrata into indices rather than stratum values */ 1080 offset = 0; 1081 ierr = PetscHashIGetKeys(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr); 1082 for (p = 0; p < size; ++p) { 1083 for (s = 0; s < (*labelNew)->numStrata; ++s) { 1084 if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;} 1085 } 1086 } 1087 /* Rebuild the point strata on the receiver */ 1088 ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr); 1089 ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1090 for (p=pStart; p<pEnd; p++) { 1091 ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1092 ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1093 for (s=0; s<dof; s++) { 1094 (*labelNew)->stratumSizes[leafStrata[offset+s]]++; 1095 } 1096 } 1097 ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr); 1098 ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr); 1099 ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr); 1100 for (s = 0; s < (*labelNew)->numStrata; ++s) { 1101 PetscHashICreate((*labelNew)->ht[s]); 1102 ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr); 1103 } 1104 /* Insert points into new strata */ 1105 ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr); 1106 ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1107 for (p=pStart; p<pEnd; p++) { 1108 ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1109 ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1110 for (s=0; s<dof; s++) { 1111 stratum = leafStrata[offset+s]; 1112 points[stratum][strataIdx[stratum]++] = p; 1113 } 1114 } 1115 for (s = 0; s < (*labelNew)->numStrata; s++) { 1116 ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr); 1117 ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr); 1118 } 1119 ierr = PetscFree(points);CHKERRQ(ierr); 1120 PetscHashIDestroy(stratumHash); 1121 ierr = PetscFree(leafStrata);CHKERRQ(ierr); 1122 ierr = PetscFree(strataIdx);CHKERRQ(ierr); 1123 ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 1124 PetscFunctionReturn(0); 1125 } 1126 1127 #undef __FUNCT__ 1128 #define __FUNCT__ "DMLabelGather" 1129 /*@ 1130 DMLabelGather - Gather all label values from leafs into roots 1131 1132 Input Parameters: 1133 + label - the DMLabel 1134 . point - the Star Forest point communication map 1135 1136 Input Parameters: 1137 + label - the new DMLabel with localised leaf values 1138 1139 Level: developer 1140 1141 Note: This is the inverse operation to DMLabelDistribute. 1142 1143 .seealso: DMLabelDistribute() 1144 @*/ 1145 PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew) 1146 { 1147 MPI_Comm comm; 1148 PetscSection rootSection; 1149 PetscSF sfLabel; 1150 PetscSFNode *rootPoints, *leafPoints; 1151 PetscInt p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset; 1152 const PetscInt *rootDegree, *ilocal; 1153 PetscInt *rootStrata; 1154 char *name; 1155 PetscInt nameSize; 1156 size_t len = 0; 1157 PetscMPIInt rank, numProcs; 1158 PetscErrorCode ierr; 1159 1160 PetscFunctionBegin; 1161 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 1162 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1163 ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 1164 /* Bcast name */ 1165 if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 1166 nameSize = len; 1167 ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1168 ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1169 if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 1170 ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1171 ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 1172 ierr = PetscFree(name);CHKERRQ(ierr); 1173 /* Gather rank/index pairs of leaves into local roots to build 1174 an inverse, multi-rooted SF. Note that this ignores local leaf 1175 indexing due to the use of the multiSF in PetscSFGather. */ 1176 ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr); 1177 ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr); 1178 for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1; 1179 for (p = 0; p < nleaves; p++) { 1180 leafPoints[ilocal[p]].index = ilocal[p]; 1181 leafPoints[ilocal[p]].rank = rank; 1182 } 1183 ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr); 1184 ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr); 1185 for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p]; 1186 ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr); 1187 ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 1188 ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 1189 ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr); 1190 ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 1191 /* Migrate label over inverted SF to pull stratum values at leaves into roots. */ 1192 ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr); 1193 /* Rebuild the point strata on the receiver */ 1194 for (p = 0, idx = 0; p < nroots; p++) { 1195 for (d = 0; d < rootDegree[p]; d++) { 1196 ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr); 1197 ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr); 1198 for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);} 1199 } 1200 idx += rootDegree[p]; 1201 } 1202 ierr = PetscFree(leafPoints);CHKERRQ(ierr); 1203 ierr = PetscFree(rootStrata);CHKERRQ(ierr); 1204 ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 1205 ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr); 1206 PetscFunctionReturn(0); 1207 } 1208 1209 #undef __FUNCT__ 1210 #define __FUNCT__ "DMLabelConvertToSection" 1211 PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is) 1212 { 1213 IS vIS; 1214 const PetscInt *values; 1215 PetscInt *points; 1216 PetscInt nV, vS = 0, vE = 0, v, N; 1217 PetscErrorCode ierr; 1218 1219 PetscFunctionBegin; 1220 ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr); 1221 ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr); 1222 ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr); 1223 if (nV) {vS = values[0]; vE = values[0]+1;} 1224 for (v = 1; v < nV; ++v) { 1225 vS = PetscMin(vS, values[v]); 1226 vE = PetscMax(vE, values[v]+1); 1227 } 1228 ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr); 1229 ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr); 1230 for (v = 0; v < nV; ++v) { 1231 PetscInt n; 1232 1233 ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr); 1234 ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr); 1235 } 1236 ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 1237 ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr); 1238 ierr = PetscMalloc1(N, &points);CHKERRQ(ierr); 1239 for (v = 0; v < nV; ++v) { 1240 IS is; 1241 const PetscInt *spoints; 1242 PetscInt dof, off, p; 1243 1244 ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr); 1245 ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr); 1246 ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr); 1247 ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr); 1248 for (p = 0; p < dof; ++p) points[off+p] = spoints[p]; 1249 ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr); 1250 ierr = ISDestroy(&is);CHKERRQ(ierr); 1251 } 1252 ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr); 1253 ierr = ISDestroy(&vIS);CHKERRQ(ierr); 1254 ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr); 1255 PetscFunctionReturn(0); 1256 } 1257 1258 #undef __FUNCT__ 1259 #define __FUNCT__ "PetscSectionCreateGlobalSectionLabel" 1260 /*@C 1261 PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using 1262 the local section and an SF describing the section point overlap. 1263 1264 Input Parameters: 1265 + s - The PetscSection for the local field layout 1266 . sf - The SF describing parallel layout of the section points 1267 . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs 1268 . label - The label specifying the points 1269 - labelValue - The label stratum specifying the points 1270 1271 Output Parameter: 1272 . gsection - The PetscSection for the global field layout 1273 1274 Note: This gives negative sizes and offsets to points not owned by this process 1275 1276 Level: developer 1277 1278 .seealso: PetscSectionCreate() 1279 @*/ 1280 PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection) 1281 { 1282 PetscInt *neg = NULL, *tmpOff = NULL; 1283 PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots; 1284 PetscErrorCode ierr; 1285 1286 PetscFunctionBegin; 1287 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr); 1288 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 1289 ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr); 1290 ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 1291 if (nroots >= 0) { 1292 if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart); 1293 ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr); 1294 if (nroots > pEnd-pStart) { 1295 ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr); 1296 } else { 1297 tmpOff = &(*gsection)->atlasDof[-pStart]; 1298 } 1299 } 1300 /* Mark ghost points with negative dof */ 1301 for (p = pStart; p < pEnd; ++p) { 1302 PetscInt value; 1303 1304 ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr); 1305 if (value != labelValue) continue; 1306 ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr); 1307 ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr); 1308 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 1309 if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);} 1310 if (neg) neg[p] = -(dof+1); 1311 } 1312 ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr); 1313 if (nroots >= 0) { 1314 ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1315 ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1316 if (nroots > pEnd-pStart) { 1317 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];} 1318 } 1319 } 1320 /* Calculate new sizes, get proccess offset, and calculate point offsets */ 1321 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1322 cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0; 1323 (*gsection)->atlasOff[p] = off; 1324 off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0; 1325 } 1326 ierr = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr); 1327 globalOff -= off; 1328 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1329 (*gsection)->atlasOff[p] += globalOff; 1330 if (neg) neg[p] = -((*gsection)->atlasOff[p]+1); 1331 } 1332 /* Put in negative offsets for ghost points */ 1333 if (nroots >= 0) { 1334 ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1335 ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1336 if (nroots > pEnd-pStart) { 1337 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];} 1338 } 1339 } 1340 if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);} 1341 ierr = PetscFree(neg);CHKERRQ(ierr); 1342 PetscFunctionReturn(0); 1343 } 1344 1345