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]) {ierr = DMLabelMakeInvalid_Private(label,v);CHKERRQ(ierr);} 663 ierr = PetscHashIDelKey(label->ht[v], point);CHKERRQ(ierr); 664 PetscFunctionReturn(0); 665 } 666 667 #undef __FUNCT__ 668 #define __FUNCT__ "DMLabelInsertIS" 669 /*@ 670 DMLabelInsertIS - Set all points in the IS to a value 671 672 Input Parameters: 673 + label - the DMLabel 674 . is - the point IS 675 - value - The point value 676 677 Level: intermediate 678 679 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 680 @*/ 681 PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value) 682 { 683 const PetscInt *points; 684 PetscInt n, p; 685 PetscErrorCode ierr; 686 687 PetscFunctionBegin; 688 PetscValidHeaderSpecific(is, IS_CLASSID, 2); 689 ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 690 ierr = ISGetIndices(is, &points);CHKERRQ(ierr); 691 for (p = 0; p < n; ++p) {ierr = DMLabelSetValue(label, points[p], value);CHKERRQ(ierr);} 692 ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr); 693 PetscFunctionReturn(0); 694 } 695 696 #undef __FUNCT__ 697 #define __FUNCT__ "DMLabelGetNumValues" 698 PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues) 699 { 700 PetscFunctionBegin; 701 PetscValidPointer(numValues, 2); 702 *numValues = label->numStrata; 703 PetscFunctionReturn(0); 704 } 705 706 #undef __FUNCT__ 707 #define __FUNCT__ "DMLabelGetValueIS" 708 PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values) 709 { 710 PetscErrorCode ierr; 711 712 PetscFunctionBegin; 713 PetscValidPointer(values, 2); 714 ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr); 715 PetscFunctionReturn(0); 716 } 717 718 #undef __FUNCT__ 719 #define __FUNCT__ "DMLabelHasStratum" 720 PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists) 721 { 722 PetscInt v; 723 724 PetscFunctionBegin; 725 PetscValidPointer(exists, 3); 726 *exists = PETSC_FALSE; 727 for (v = 0; v < label->numStrata; ++v) { 728 if (label->stratumValues[v] == value) { 729 *exists = PETSC_TRUE; 730 break; 731 } 732 } 733 PetscFunctionReturn(0); 734 } 735 736 #undef __FUNCT__ 737 #define __FUNCT__ "DMLabelGetStratumSize" 738 PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size) 739 { 740 PetscInt v; 741 PetscErrorCode ierr; 742 743 PetscFunctionBegin; 744 PetscValidPointer(size, 3); 745 *size = 0; 746 for (v = 0; v < label->numStrata; ++v) { 747 if (label->stratumValues[v] == value) { 748 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 749 *size = label->stratumSizes[v]; 750 break; 751 } 752 } 753 PetscFunctionReturn(0); 754 } 755 756 #undef __FUNCT__ 757 #define __FUNCT__ "DMLabelGetStratumBounds" 758 PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end) 759 { 760 PetscInt v; 761 PetscErrorCode ierr; 762 763 PetscFunctionBegin; 764 if (start) {PetscValidPointer(start, 3); *start = 0;} 765 if (end) {PetscValidPointer(end, 4); *end = 0;} 766 for (v = 0; v < label->numStrata; ++v) { 767 const PetscInt *points; 768 769 if (label->stratumValues[v] != value) continue; 770 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 771 if (label->stratumSizes[v] <= 0) break; 772 ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 773 if (start) *start = points[0]; 774 if (end) *end = points[label->stratumSizes[v]-1]+1; 775 ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 776 break; 777 } 778 PetscFunctionReturn(0); 779 } 780 781 #undef __FUNCT__ 782 #define __FUNCT__ "DMLabelGetStratumIS" 783 PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points) 784 { 785 PetscInt v; 786 PetscErrorCode ierr; 787 788 PetscFunctionBegin; 789 PetscValidPointer(points, 3); 790 *points = NULL; 791 for (v = 0; v < label->numStrata; ++v) { 792 if (label->stratumValues[v] == value) { 793 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 794 if (label->validIS[v]) { 795 ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr); 796 *points = label->points[v]; 797 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to implement this to speedup Stratify"); 798 break; 799 } 800 } 801 PetscFunctionReturn(0); 802 } 803 804 #undef __FUNCT__ 805 #define __FUNCT__ "DMLabelSetStratumIS" 806 PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is) 807 { 808 PetscInt v, numStrata; 809 PetscErrorCode ierr; 810 811 PetscFunctionBegin; 812 numStrata = label->numStrata; 813 for (v = 0; v < numStrata; v++) { 814 if (label->stratumValues[v] == value) break; 815 } 816 if (v >= numStrata) {ierr = DMLabelAddStratum(label,value);CHKERRQ(ierr);} 817 if (is == label->points[v]) PetscFunctionReturn(0); 818 ierr = DMLabelClearStratum(label,value);CHKERRQ(ierr); 819 ierr = ISGetLocalSize(is,&(label->stratumSizes[v]));CHKERRQ(ierr); 820 label->stratumValues[v] = value; 821 label->validIS[v] = PETSC_TRUE; 822 ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 823 ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 824 if (label->bt) { 825 const PetscInt *points; 826 PetscInt p; 827 828 ierr = ISGetIndices(is,&points);CHKERRQ(ierr); 829 for (p = 0; p < label->stratumSizes[v]; ++p) { 830 const PetscInt point = points[p]; 831 832 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); 833 ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 834 } 835 } 836 label->points[v] = is; 837 PetscFunctionReturn(0); 838 } 839 840 841 #undef __FUNCT__ 842 #define __FUNCT__ "DMLabelClearStratum" 843 PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value) 844 { 845 PetscInt v; 846 PetscErrorCode ierr; 847 848 PetscFunctionBegin; 849 for (v = 0; v < label->numStrata; ++v) { 850 if (label->stratumValues[v] == value) break; 851 } 852 if (v >= label->numStrata) PetscFunctionReturn(0); 853 if (label->validIS[v]) { 854 if (label->bt) { 855 PetscInt i; 856 const PetscInt *points; 857 858 ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 859 for (i = 0; i < label->stratumSizes[v]; ++i) { 860 const PetscInt point = points[i]; 861 862 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); 863 ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 864 } 865 ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr); 866 } 867 ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 868 label->stratumSizes[v] = 0; 869 ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr); 870 ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr); 871 } else { 872 PetscHashIClear(label->ht[v]); 873 } 874 PetscFunctionReturn(0); 875 } 876 877 #undef __FUNCT__ 878 #define __FUNCT__ "DMLabelFilter" 879 PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end) 880 { 881 PetscInt v; 882 PetscErrorCode ierr; 883 884 PetscFunctionBegin; 885 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 886 label->pStart = start; 887 label->pEnd = end; 888 if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 889 /* Could squish offsets, but would only make sense if I reallocate the storage */ 890 for (v = 0; v < label->numStrata; ++v) { 891 PetscInt off, q; 892 const PetscInt *points; 893 PetscInt *pointsNew = NULL; 894 895 ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 896 for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) { 897 const PetscInt point = points[q]; 898 899 if ((point < start) || (point >= end)) { 900 if (!pointsNew) { 901 ierr = PetscMalloc1(label->stratumSizes[v],&pointsNew);CHKERRQ(ierr); 902 ierr = PetscMemcpy(pointsNew,points,(size_t) off * sizeof(PetscInt));CHKERRQ(ierr); 903 } 904 continue; 905 } 906 if (pointsNew) { 907 pointsNew[off++] = point; 908 } 909 } 910 ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 911 if (pointsNew) { 912 ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 913 ierr = ISCreateGeneral(PETSC_COMM_SELF,off,pointsNew,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr); 914 ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr); 915 } 916 label->stratumSizes[v] = off; 917 } 918 ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr); 919 PetscFunctionReturn(0); 920 } 921 922 #undef __FUNCT__ 923 #define __FUNCT__ "DMLabelPermute" 924 PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew) 925 { 926 const PetscInt *perm; 927 PetscInt numValues, numPoints, v, q; 928 PetscErrorCode ierr; 929 930 PetscFunctionBegin; 931 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 932 ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr); 933 ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr); 934 ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr); 935 ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr); 936 for (v = 0; v < numValues; ++v) { 937 const PetscInt size = (*labelNew)->stratumSizes[v]; 938 const PetscInt *points; 939 PetscInt *pointsNew; 940 941 ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 942 ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr); 943 for (q = 0; q < size; ++q) { 944 const PetscInt point = points[q]; 945 946 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); 947 pointsNew[q] = perm[point]; 948 } 949 ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 950 ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr); 951 ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr); 952 ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr); 953 ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr); 954 } 955 ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr); 956 if (label->bt) { 957 ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 958 ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr); 959 } 960 PetscFunctionReturn(0); 961 } 962 963 #undef __FUNCT__ 964 #define __FUNCT__ "DMLabelDistribute_Internal" 965 PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata) 966 { 967 MPI_Comm comm; 968 PetscInt s, l, nroots, nleaves, dof, offset, size; 969 PetscInt *remoteOffsets, *rootStrata, *rootIdx; 970 PetscSection rootSection; 971 PetscSF labelSF; 972 PetscErrorCode ierr; 973 974 PetscFunctionBegin; 975 if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 976 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 977 /* Build a section of stratum values per point, generate the according SF 978 and distribute point-wise stratum values to leaves. */ 979 ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 980 ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 981 ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr); 982 if (label) { 983 for (s = 0; s < label->numStrata; ++s) { 984 const PetscInt *points; 985 986 ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 987 for (l = 0; l < label->stratumSizes[s]; l++) { 988 ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr); 989 ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr); 990 } 991 ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 992 } 993 } 994 ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 995 /* Create a point-wise array of stratum values */ 996 ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr); 997 ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr); 998 ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr); 999 if (label) { 1000 for (s = 0; s < label->numStrata; ++s) { 1001 const PetscInt *points; 1002 1003 ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 1004 for (l = 0; l < label->stratumSizes[s]; l++) { 1005 const PetscInt p = points[l]; 1006 ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr); 1007 rootStrata[offset+rootIdx[p]++] = label->stratumValues[s]; 1008 } 1009 ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 1010 } 1011 } 1012 /* Build SF that maps label points to remote processes */ 1013 ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr); 1014 ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr); 1015 ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr); 1016 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 1017 /* Send the strata for each point over the derived SF */ 1018 ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr); 1019 ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr); 1020 ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 1021 ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 1022 /* Clean up */ 1023 ierr = PetscFree(rootStrata);CHKERRQ(ierr); 1024 ierr = PetscFree(rootIdx);CHKERRQ(ierr); 1025 ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 1026 ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr); 1027 PetscFunctionReturn(0); 1028 } 1029 1030 #undef __FUNCT__ 1031 #define __FUNCT__ "DMLabelDistribute" 1032 PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew) 1033 { 1034 MPI_Comm comm; 1035 PetscSection leafSection; 1036 PetscInt p, pStart, pEnd, s, size, dof, offset, stratum; 1037 PetscInt *leafStrata, *strataIdx; 1038 PetscInt **points; 1039 char *name; 1040 PetscInt nameSize; 1041 PetscHashI stratumHash; 1042 PETSC_UNUSED PetscHashIIter ret, iter; 1043 size_t len = 0; 1044 PetscMPIInt rank; 1045 PetscErrorCode ierr; 1046 1047 PetscFunctionBegin; 1048 if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 1049 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 1050 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1051 /* Bcast name */ 1052 if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 1053 nameSize = len; 1054 ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1055 ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1056 if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 1057 ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1058 ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 1059 ierr = PetscFree(name);CHKERRQ(ierr); 1060 /* Bcast defaultValue */ 1061 if (!rank) (*labelNew)->defaultValue = label->defaultValue; 1062 ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1063 /* Distribute stratum values over the SF and get the point mapping on the receiver */ 1064 ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr); 1065 /* Determine received stratum values and initialise new label*/ 1066 PetscHashICreate(stratumHash); 1067 ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr); 1068 for (p = 0; p < size; ++p) PetscHashIPut(stratumHash, leafStrata[p], ret, iter); 1069 PetscHashISize(stratumHash, (*labelNew)->numStrata); 1070 ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr); 1071 for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE; 1072 ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr); 1073 /* Turn leafStrata into indices rather than stratum values */ 1074 offset = 0; 1075 ierr = PetscHashIGetKeys(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr); 1076 for (p = 0; p < size; ++p) { 1077 for (s = 0; s < (*labelNew)->numStrata; ++s) { 1078 if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;} 1079 } 1080 } 1081 /* Rebuild the point strata on the receiver */ 1082 ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr); 1083 ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1084 for (p=pStart; p<pEnd; p++) { 1085 ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1086 ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1087 for (s=0; s<dof; s++) { 1088 (*labelNew)->stratumSizes[leafStrata[offset+s]]++; 1089 } 1090 } 1091 ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr); 1092 ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr); 1093 ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr); 1094 for (s = 0; s < (*labelNew)->numStrata; ++s) { 1095 PetscHashICreate((*labelNew)->ht[s]); 1096 ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr); 1097 } 1098 /* Insert points into new strata */ 1099 ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr); 1100 ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1101 for (p=pStart; p<pEnd; p++) { 1102 ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1103 ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1104 for (s=0; s<dof; s++) { 1105 stratum = leafStrata[offset+s]; 1106 points[stratum][strataIdx[stratum]++] = p; 1107 } 1108 } 1109 for (s = 0; s < (*labelNew)->numStrata; s++) { 1110 ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr); 1111 ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr); 1112 } 1113 ierr = PetscFree(points);CHKERRQ(ierr); 1114 PetscHashIDestroy(stratumHash); 1115 ierr = PetscFree(leafStrata);CHKERRQ(ierr); 1116 ierr = PetscFree(strataIdx);CHKERRQ(ierr); 1117 ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 1118 PetscFunctionReturn(0); 1119 } 1120 1121 #undef __FUNCT__ 1122 #define __FUNCT__ "DMLabelGather" 1123 /*@ 1124 DMLabelGather - Gather all label values from leafs into roots 1125 1126 Input Parameters: 1127 + label - the DMLabel 1128 . point - the Star Forest point communication map 1129 1130 Input Parameters: 1131 + label - the new DMLabel with localised leaf values 1132 1133 Level: developer 1134 1135 Note: This is the inverse operation to DMLabelDistribute. 1136 1137 .seealso: DMLabelDistribute() 1138 @*/ 1139 PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew) 1140 { 1141 MPI_Comm comm; 1142 PetscSection rootSection; 1143 PetscSF sfLabel; 1144 PetscSFNode *rootPoints, *leafPoints; 1145 PetscInt p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset; 1146 const PetscInt *rootDegree, *ilocal; 1147 PetscInt *rootStrata; 1148 char *name; 1149 PetscInt nameSize; 1150 size_t len = 0; 1151 PetscMPIInt rank, numProcs; 1152 PetscErrorCode ierr; 1153 1154 PetscFunctionBegin; 1155 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 1156 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1157 ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 1158 /* Bcast name */ 1159 if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 1160 nameSize = len; 1161 ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1162 ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1163 if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 1164 ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1165 ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 1166 ierr = PetscFree(name);CHKERRQ(ierr); 1167 /* Gather rank/index pairs of leaves into local roots to build 1168 an inverse, multi-rooted SF. Note that this ignores local leaf 1169 indexing due to the use of the multiSF in PetscSFGather. */ 1170 ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr); 1171 ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr); 1172 for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1; 1173 for (p = 0; p < nleaves; p++) { 1174 leafPoints[ilocal[p]].index = ilocal[p]; 1175 leafPoints[ilocal[p]].rank = rank; 1176 } 1177 ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr); 1178 ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr); 1179 for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p]; 1180 ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr); 1181 ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 1182 ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 1183 ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr); 1184 ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 1185 /* Migrate label over inverted SF to pull stratum values at leaves into roots. */ 1186 ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr); 1187 /* Rebuild the point strata on the receiver */ 1188 for (p = 0, idx = 0; p < nroots; p++) { 1189 for (d = 0; d < rootDegree[p]; d++) { 1190 ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr); 1191 ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr); 1192 for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);} 1193 } 1194 idx += rootDegree[p]; 1195 } 1196 ierr = PetscFree(leafPoints);CHKERRQ(ierr); 1197 ierr = PetscFree(rootStrata);CHKERRQ(ierr); 1198 ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 1199 ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr); 1200 PetscFunctionReturn(0); 1201 } 1202 1203 #undef __FUNCT__ 1204 #define __FUNCT__ "DMLabelConvertToSection" 1205 PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is) 1206 { 1207 IS vIS; 1208 const PetscInt *values; 1209 PetscInt *points; 1210 PetscInt nV, vS = 0, vE = 0, v, N; 1211 PetscErrorCode ierr; 1212 1213 PetscFunctionBegin; 1214 ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr); 1215 ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr); 1216 ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr); 1217 if (nV) {vS = values[0]; vE = values[0]+1;} 1218 for (v = 1; v < nV; ++v) { 1219 vS = PetscMin(vS, values[v]); 1220 vE = PetscMax(vE, values[v]+1); 1221 } 1222 ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr); 1223 ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr); 1224 for (v = 0; v < nV; ++v) { 1225 PetscInt n; 1226 1227 ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr); 1228 ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr); 1229 } 1230 ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 1231 ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr); 1232 ierr = PetscMalloc1(N, &points);CHKERRQ(ierr); 1233 for (v = 0; v < nV; ++v) { 1234 IS is; 1235 const PetscInt *spoints; 1236 PetscInt dof, off, p; 1237 1238 ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr); 1239 ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr); 1240 ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr); 1241 ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr); 1242 for (p = 0; p < dof; ++p) points[off+p] = spoints[p]; 1243 ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr); 1244 ierr = ISDestroy(&is);CHKERRQ(ierr); 1245 } 1246 ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr); 1247 ierr = ISDestroy(&vIS);CHKERRQ(ierr); 1248 ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr); 1249 PetscFunctionReturn(0); 1250 } 1251 1252 #undef __FUNCT__ 1253 #define __FUNCT__ "PetscSectionCreateGlobalSectionLabel" 1254 /*@C 1255 PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using 1256 the local section and an SF describing the section point overlap. 1257 1258 Input Parameters: 1259 + s - The PetscSection for the local field layout 1260 . sf - The SF describing parallel layout of the section points 1261 . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs 1262 . label - The label specifying the points 1263 - labelValue - The label stratum specifying the points 1264 1265 Output Parameter: 1266 . gsection - The PetscSection for the global field layout 1267 1268 Note: This gives negative sizes and offsets to points not owned by this process 1269 1270 Level: developer 1271 1272 .seealso: PetscSectionCreate() 1273 @*/ 1274 PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection) 1275 { 1276 PetscInt *neg = NULL, *tmpOff = NULL; 1277 PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots; 1278 PetscErrorCode ierr; 1279 1280 PetscFunctionBegin; 1281 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr); 1282 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 1283 ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr); 1284 ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 1285 if (nroots >= 0) { 1286 if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart); 1287 ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr); 1288 if (nroots > pEnd-pStart) { 1289 ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr); 1290 } else { 1291 tmpOff = &(*gsection)->atlasDof[-pStart]; 1292 } 1293 } 1294 /* Mark ghost points with negative dof */ 1295 for (p = pStart; p < pEnd; ++p) { 1296 PetscInt value; 1297 1298 ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr); 1299 if (value != labelValue) continue; 1300 ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr); 1301 ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr); 1302 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 1303 if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);} 1304 if (neg) neg[p] = -(dof+1); 1305 } 1306 ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr); 1307 if (nroots >= 0) { 1308 ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1309 ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1310 if (nroots > pEnd-pStart) { 1311 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];} 1312 } 1313 } 1314 /* Calculate new sizes, get proccess offset, and calculate point offsets */ 1315 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1316 cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0; 1317 (*gsection)->atlasOff[p] = off; 1318 off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0; 1319 } 1320 ierr = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr); 1321 globalOff -= off; 1322 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1323 (*gsection)->atlasOff[p] += globalOff; 1324 if (neg) neg[p] = -((*gsection)->atlasOff[p]+1); 1325 } 1326 /* Put in negative offsets for ghost points */ 1327 if (nroots >= 0) { 1328 ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1329 ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1330 if (nroots > pEnd-pStart) { 1331 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];} 1332 } 1333 } 1334 if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);} 1335 ierr = PetscFree(neg);CHKERRQ(ierr); 1336 PetscFunctionReturn(0); 1337 } 1338 1339