1 #include <petscdm.h> 2 #include <petsc/private/dmlabelimpl.h> /*I "petscdmlabel.h" I*/ 3 #include <petsc/private/isimpl.h> /*I "petscis.h" I*/ 4 #include <petscsf.h> 5 6 #undef __FUNCT__ 7 #define __FUNCT__ "DMLabelCreate" 8 /*@C 9 DMLabelCreate - Create a DMLabel object, which is a multimap 10 11 Input parameter: 12 . name - The label name 13 14 Output parameter: 15 . label - The DMLabel 16 17 Level: beginner 18 19 .seealso: DMLabelDestroy() 20 @*/ 21 PetscErrorCode DMLabelCreate(const char name[], DMLabel *label) 22 { 23 PetscErrorCode ierr; 24 25 PetscFunctionBegin; 26 ierr = PetscNew(label);CHKERRQ(ierr); 27 ierr = PetscStrallocpy(name, &(*label)->name);CHKERRQ(ierr); 28 29 (*label)->refct = 1; 30 (*label)->state = -1; 31 (*label)->numStrata = 0; 32 (*label)->defaultValue = -1; 33 (*label)->stratumValues = NULL; 34 (*label)->validIS = NULL; 35 (*label)->stratumSizes = NULL; 36 (*label)->points = NULL; 37 (*label)->ht = NULL; 38 (*label)->pStart = -1; 39 (*label)->pEnd = -1; 40 (*label)->bt = NULL; 41 PetscFunctionReturn(0); 42 } 43 44 #undef __FUNCT__ 45 #define __FUNCT__ "DMLabelMakeValid_Private" 46 /* 47 DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format 48 49 Input parameter: 50 + label - The DMLabel 51 - v - The stratum value 52 53 Output parameter: 54 . label - The DMLabel with stratum in sorted list format 55 56 Level: developer 57 58 .seealso: DMLabelCreate() 59 */ 60 static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v) 61 { 62 PetscInt off; 63 PetscInt *pointArray; 64 PetscErrorCode ierr; 65 66 if (label->validIS[v]) return 0; 67 if (v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeValid_Private\n", v); 68 PetscFunctionBegin; 69 PetscHashISize(label->ht[v], label->stratumSizes[v]); 70 71 ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr); 72 off = 0; 73 ierr = PetscHashIGetKeys(label->ht[v], &off, pointArray);CHKERRQ(ierr); 74 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]); 75 PetscHashIClear(label->ht[v]); 76 ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr); 77 if (label->bt) { 78 PetscInt p; 79 80 for (p = 0; p < label->stratumSizes[v]; ++p) { 81 const PetscInt point = pointArray[p]; 82 83 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); 84 ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 85 } 86 } 87 ierr = ISCreateGeneral(PETSC_COMM_SELF,label->stratumSizes[v],pointArray,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr); 88 ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr); 89 label->validIS[v] = PETSC_TRUE; 90 ++label->state; 91 PetscFunctionReturn(0); 92 } 93 94 #undef __FUNCT__ 95 #define __FUNCT__ "DMLabelMakeAllValid_Private" 96 /* 97 DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format 98 99 Input parameter: 100 . label - The DMLabel 101 102 Output parameter: 103 . label - The DMLabel with all strata in sorted list format 104 105 Level: developer 106 107 .seealso: DMLabelCreate() 108 */ 109 static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label) 110 { 111 PetscInt v; 112 PetscErrorCode ierr; 113 114 PetscFunctionBegin; 115 for (v = 0; v < label->numStrata; v++){ 116 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 117 } 118 PetscFunctionReturn(0); 119 } 120 121 #undef __FUNCT__ 122 #define __FUNCT__ "DMLabelMakeInvalid_Private" 123 /* 124 DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format 125 126 Input parameter: 127 + label - The DMLabel 128 - v - The stratum value 129 130 Output parameter: 131 . label - The DMLabel with stratum in hash format 132 133 Level: developer 134 135 .seealso: DMLabelCreate() 136 */ 137 static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v) 138 { 139 PETSC_UNUSED PetscHashIIter ret, iter; 140 PetscInt p; 141 const PetscInt *points; 142 PetscErrorCode ierr; 143 144 PetscFunctionBegin; 145 if (!label->validIS[v]) PetscFunctionReturn(0); 146 if (label->points[v]) { 147 ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 148 for (p = 0; p < label->stratumSizes[v]; ++p) PetscHashIPut(label->ht[v], points[p], ret, iter); 149 ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 150 ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 151 } 152 label->validIS[v] = PETSC_FALSE; 153 PetscFunctionReturn(0); 154 } 155 156 #undef __FUNCT__ 157 #define __FUNCT__ "DMLabelGetState" 158 PetscErrorCode DMLabelGetState(DMLabel label, PetscObjectState *state) 159 { 160 PetscFunctionBegin; 161 PetscValidPointer(state, 2); 162 *state = label->state; 163 PetscFunctionReturn(0); 164 } 165 166 #undef __FUNCT__ 167 #define __FUNCT__ "DMLabelAddStratum" 168 PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value) 169 { 170 PetscInt v, *tmpV, *tmpS; 171 IS *tmpP; 172 PetscHashI *tmpH; 173 PetscBool *tmpB; 174 PetscErrorCode ierr; 175 176 PetscFunctionBegin; 177 178 for (v = 0; v < label->numStrata; v++) { 179 if (label->stratumValues[v] == value) PetscFunctionReturn(0); 180 } 181 ierr = PetscMalloc1((label->numStrata+1), &tmpV);CHKERRQ(ierr); 182 ierr = PetscMalloc1((label->numStrata+1), &tmpS);CHKERRQ(ierr); 183 ierr = PetscMalloc1((label->numStrata+1), &tmpH);CHKERRQ(ierr); 184 ierr = PetscMalloc1((label->numStrata+1), &tmpP);CHKERRQ(ierr); 185 ierr = PetscMalloc1((label->numStrata+1), &tmpB);CHKERRQ(ierr); 186 for (v = 0; v < label->numStrata; ++v) { 187 tmpV[v] = label->stratumValues[v]; 188 tmpS[v] = label->stratumSizes[v]; 189 tmpH[v] = label->ht[v]; 190 tmpP[v] = label->points[v]; 191 tmpB[v] = label->validIS[v]; 192 } 193 tmpV[v] = value; 194 tmpS[v] = 0; 195 PetscHashICreate(tmpH[v]); 196 ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&tmpP[v]);CHKERRQ(ierr); 197 tmpB[v] = PETSC_TRUE; 198 ++label->numStrata; 199 ierr = PetscFree(label->stratumValues);CHKERRQ(ierr); 200 ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr); 201 ierr = PetscFree(label->ht);CHKERRQ(ierr); 202 ierr = PetscFree(label->points);CHKERRQ(ierr); 203 ierr = PetscFree(label->validIS);CHKERRQ(ierr); 204 label->stratumValues = tmpV; 205 label->stratumSizes = tmpS; 206 label->ht = tmpH; 207 label->points = tmpP; 208 label->validIS = tmpB; 209 210 PetscFunctionReturn(0); 211 } 212 213 #undef __FUNCT__ 214 #define __FUNCT__ "DMLabelGetName" 215 /*@C 216 DMLabelGetName - Return the name of a DMLabel object 217 218 Input parameter: 219 . label - The DMLabel 220 221 Output parameter: 222 . name - The label name 223 224 Level: beginner 225 226 .seealso: DMLabelCreate() 227 @*/ 228 PetscErrorCode DMLabelGetName(DMLabel label, const char **name) 229 { 230 PetscFunctionBegin; 231 PetscValidPointer(name, 2); 232 *name = label->name; 233 PetscFunctionReturn(0); 234 } 235 236 #undef __FUNCT__ 237 #define __FUNCT__ "DMLabelView_Ascii" 238 static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer) 239 { 240 PetscInt v; 241 PetscMPIInt rank; 242 PetscErrorCode ierr; 243 244 PetscFunctionBegin; 245 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr); 246 ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 247 if (label) { 248 ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", label->name);CHKERRQ(ierr); 249 if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, " Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);} 250 for (v = 0; v < label->numStrata; ++v) { 251 const PetscInt value = label->stratumValues[v]; 252 const PetscInt *points; 253 PetscInt p; 254 255 ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 256 for (p = 0; p < label->stratumSizes[v]; ++p) { 257 ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr); 258 } 259 ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 260 } 261 } 262 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 263 ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 264 PetscFunctionReturn(0); 265 } 266 267 #undef __FUNCT__ 268 #define __FUNCT__ "DMLabelView" 269 /*@C 270 DMLabelView - View the label 271 272 Input Parameters: 273 + label - The DMLabel 274 - viewer - The PetscViewer 275 276 Level: intermediate 277 278 .seealso: DMLabelCreate(), DMLabelDestroy() 279 @*/ 280 PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer) 281 { 282 PetscBool iascii; 283 PetscErrorCode ierr; 284 285 PetscFunctionBegin; 286 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 287 if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 288 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 289 if (iascii) { 290 ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr); 291 } 292 PetscFunctionReturn(0); 293 } 294 295 #undef __FUNCT__ 296 #define __FUNCT__ "DMLabelDestroy" 297 PetscErrorCode DMLabelDestroy(DMLabel *label) 298 { 299 PetscInt v; 300 PetscErrorCode ierr; 301 302 PetscFunctionBegin; 303 if (!(*label)) PetscFunctionReturn(0); 304 if (--(*label)->refct > 0) PetscFunctionReturn(0); 305 ierr = PetscFree((*label)->name);CHKERRQ(ierr); 306 ierr = PetscFree((*label)->stratumValues);CHKERRQ(ierr); 307 ierr = PetscFree((*label)->stratumSizes);CHKERRQ(ierr); 308 for (v = 0; v < (*label)->numStrata; ++v) {ierr = ISDestroy(&((*label)->points[v]));CHKERRQ(ierr);} 309 ierr = PetscFree((*label)->points);CHKERRQ(ierr); 310 ierr = PetscFree((*label)->validIS);CHKERRQ(ierr); 311 if ((*label)->ht) { 312 for (v = 0; v < (*label)->numStrata; ++v) {PetscHashIDestroy((*label)->ht[v]);} 313 ierr = PetscFree((*label)->ht);CHKERRQ(ierr); 314 } 315 ierr = PetscBTDestroy(&(*label)->bt);CHKERRQ(ierr); 316 ierr = PetscFree(*label);CHKERRQ(ierr); 317 PetscFunctionReturn(0); 318 } 319 320 #undef __FUNCT__ 321 #define __FUNCT__ "DMLabelDuplicate" 322 PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew) 323 { 324 PetscInt v; 325 PetscErrorCode ierr; 326 327 PetscFunctionBegin; 328 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 329 ierr = PetscNew(labelnew);CHKERRQ(ierr); 330 ierr = PetscStrallocpy(label->name, &(*labelnew)->name);CHKERRQ(ierr); 331 332 (*labelnew)->refct = 1; 333 (*labelnew)->numStrata = label->numStrata; 334 (*labelnew)->defaultValue = label->defaultValue; 335 if (label->numStrata) { 336 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr); 337 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr); 338 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr); 339 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr); 340 ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr); 341 /* Could eliminate unused space here */ 342 for (v = 0; v < label->numStrata; ++v) { 343 PetscHashICreate((*labelnew)->ht[v]); 344 (*labelnew)->validIS[v] = PETSC_TRUE; 345 (*labelnew)->stratumValues[v] = label->stratumValues[v]; 346 (*labelnew)->stratumSizes[v] = label->stratumSizes[v]; 347 ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr); 348 (*labelnew)->points[v] = label->points[v]; 349 } 350 } 351 (*labelnew)->pStart = -1; 352 (*labelnew)->pEnd = -1; 353 (*labelnew)->bt = NULL; 354 PetscFunctionReturn(0); 355 } 356 357 #undef __FUNCT__ 358 #define __FUNCT__ "DMLabelCreateIndex" 359 /* This can be hooked into SetValue(), ClearValue(), etc. for updating */ 360 PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd) 361 { 362 PetscInt v; 363 PetscErrorCode ierr; 364 365 PetscFunctionBegin; 366 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 367 if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 368 label->pStart = pStart; 369 label->pEnd = pEnd; 370 ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr); 371 ierr = PetscBTMemzero(pEnd - pStart, label->bt);CHKERRQ(ierr); 372 for (v = 0; v < label->numStrata; ++v) { 373 const PetscInt *points; 374 PetscInt i; 375 376 ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 377 for (i = 0; i < label->stratumSizes[v]; ++i) { 378 const PetscInt point = points[i]; 379 380 if ((point < pStart) || (point >= pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, pStart, pEnd); 381 ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr); 382 } 383 ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 384 } 385 PetscFunctionReturn(0); 386 } 387 388 #undef __FUNCT__ 389 #define __FUNCT__ "DMLabelDestroyIndex" 390 PetscErrorCode DMLabelDestroyIndex(DMLabel label) 391 { 392 PetscErrorCode ierr; 393 394 PetscFunctionBegin; 395 label->pStart = -1; 396 label->pEnd = -1; 397 if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 398 PetscFunctionReturn(0); 399 } 400 401 #undef __FUNCT__ 402 #define __FUNCT__ "DMLabelHasValue" 403 /*@ 404 DMLabelHasValue - Determine whether a label assigns the value to any point 405 406 Input Parameters: 407 + label - the DMLabel 408 - value - the value 409 410 Output Parameter: 411 . contains - Flag indicating whether the label maps this value to any point 412 413 Level: developer 414 415 .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue() 416 @*/ 417 PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains) 418 { 419 PetscInt v; 420 421 PetscFunctionBegin; 422 PetscValidPointer(contains, 3); 423 for (v = 0; v < label->numStrata; ++v) { 424 if (value == label->stratumValues[v]) break; 425 } 426 *contains = (v < label->numStrata ? PETSC_TRUE : PETSC_FALSE); 427 PetscFunctionReturn(0); 428 } 429 430 #undef __FUNCT__ 431 #define __FUNCT__ "DMLabelHasPoint" 432 /*@ 433 DMLabelHasPoint - Determine whether a label assigns a value to a point 434 435 Input Parameters: 436 + label - the DMLabel 437 - point - the point 438 439 Output Parameter: 440 . contains - Flag indicating whether the label maps this point to a value 441 442 Note: The user must call DMLabelCreateIndex() before this function. 443 444 Level: developer 445 446 .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue() 447 @*/ 448 PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains) 449 { 450 PetscErrorCode ierr; 451 452 PetscFunctionBeginHot; 453 PetscValidPointer(contains, 3); 454 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 455 #if defined(PETSC_USE_DEBUG) 456 if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()"); 457 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); 458 #endif 459 *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE; 460 PetscFunctionReturn(0); 461 } 462 463 #undef __FUNCT__ 464 #define __FUNCT__ "DMLabelStratumHasPoint" 465 /*@ 466 DMLabelStratumHasPoint - Return true if the stratum contains a point 467 468 Input Parameters: 469 + label - the DMLabel 470 . value - the stratum value 471 - point - the point 472 473 Output Parameter: 474 . contains - true if the stratum contains the point 475 476 Level: intermediate 477 478 .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue() 479 @*/ 480 PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains) 481 { 482 PetscInt v; 483 PetscErrorCode ierr; 484 485 PetscFunctionBegin; 486 PetscValidPointer(contains, 4); 487 *contains = PETSC_FALSE; 488 for (v = 0; v < label->numStrata; ++v) { 489 if (label->stratumValues[v] == value) { 490 if (label->validIS[v]) { 491 PetscInt i; 492 493 ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr); 494 if (i >= 0) { 495 *contains = PETSC_TRUE; 496 break; 497 } 498 } else { 499 PetscBool has; 500 501 PetscHashIHasKey(label->ht[v], point, has); 502 if (has) { 503 *contains = PETSC_TRUE; 504 break; 505 } 506 } 507 } 508 } 509 PetscFunctionReturn(0); 510 } 511 512 #undef __FUNCT__ 513 #define __FUNCT__ "DMLabelGetDefaultValue" 514 /* 515 DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 516 When a label is created, it is initialized to -1. 517 518 Input parameter: 519 . label - a DMLabel object 520 521 Output parameter: 522 . defaultValue - the default value 523 524 Level: beginner 525 526 .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 527 */ 528 PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue) 529 { 530 PetscFunctionBegin; 531 *defaultValue = label->defaultValue; 532 PetscFunctionReturn(0); 533 } 534 535 #undef __FUNCT__ 536 #define __FUNCT__ "DMLabelSetDefaultValue" 537 /* 538 DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 539 When a label is created, it is initialized to -1. 540 541 Input parameter: 542 . label - a DMLabel object 543 544 Output parameter: 545 . defaultValue - the default value 546 547 Level: beginner 548 549 .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 550 */ 551 PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue) 552 { 553 PetscFunctionBegin; 554 label->defaultValue = defaultValue; 555 PetscFunctionReturn(0); 556 } 557 558 #undef __FUNCT__ 559 #define __FUNCT__ "DMLabelGetValue" 560 /*@ 561 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()) 562 563 Input Parameters: 564 + label - the DMLabel 565 - point - the point 566 567 Output Parameter: 568 . value - The point value, or -1 569 570 Level: intermediate 571 572 .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 573 @*/ 574 PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value) 575 { 576 PetscInt v; 577 PetscErrorCode ierr; 578 579 PetscFunctionBegin; 580 PetscValidPointer(value, 3); 581 *value = label->defaultValue; 582 for (v = 0; v < label->numStrata; ++v) { 583 if (label->validIS[v]) { 584 PetscInt i; 585 586 ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr); 587 if (i >= 0) { 588 *value = label->stratumValues[v]; 589 break; 590 } 591 } else { 592 PetscBool has; 593 594 PetscHashIHasKey(label->ht[v], point, has); 595 if (has) { 596 *value = label->stratumValues[v]; 597 break; 598 } 599 } 600 } 601 PetscFunctionReturn(0); 602 } 603 604 #undef __FUNCT__ 605 #define __FUNCT__ "DMLabelSetValue" 606 /*@ 607 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. 608 609 Input Parameters: 610 + label - the DMLabel 611 . point - the point 612 - value - The point value 613 614 Level: intermediate 615 616 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 617 @*/ 618 PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value) 619 { 620 PETSC_UNUSED PetscHashIIter iter, ret; 621 PetscInt v; 622 PetscErrorCode ierr; 623 624 PetscFunctionBegin; 625 /* Find, or add, label value */ 626 if (value == label->defaultValue) PetscFunctionReturn(0); 627 for (v = 0; v < label->numStrata; ++v) { 628 if (label->stratumValues[v] == value) break; 629 } 630 /* Create new table */ 631 if (v >= label->numStrata) {ierr = DMLabelAddStratum(label, value);CHKERRQ(ierr);} 632 ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 633 /* Set key */ 634 PetscHashIPut(label->ht[v], point, ret, iter); 635 PetscFunctionReturn(0); 636 } 637 638 #undef __FUNCT__ 639 #define __FUNCT__ "DMLabelClearValue" 640 /*@ 641 DMLabelClearValue - Clear the value a label assigns to a point 642 643 Input Parameters: 644 + label - the DMLabel 645 . point - the point 646 - value - The point value 647 648 Level: intermediate 649 650 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue() 651 @*/ 652 PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value) 653 { 654 PetscInt v; 655 PetscErrorCode ierr; 656 657 PetscFunctionBegin; 658 /* Find label value */ 659 for (v = 0; v < label->numStrata; ++v) { 660 if (label->stratumValues[v] == value) break; 661 } 662 if (v >= label->numStrata) PetscFunctionReturn(0); 663 if (label->validIS[v]) { 664 ierr = DMLabelMakeInvalid_Private(label,v);CHKERRQ(ierr); 665 } 666 if (label->bt) { 667 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); 668 ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 669 } 670 ierr = PetscHashIDelKey(label->ht[v], point);CHKERRQ(ierr); 671 PetscFunctionReturn(0); 672 } 673 674 #undef __FUNCT__ 675 #define __FUNCT__ "DMLabelInsertIS" 676 /*@ 677 DMLabelInsertIS - Set all points in the IS to a value 678 679 Input Parameters: 680 + label - the DMLabel 681 . is - the point IS 682 - value - The point value 683 684 Level: intermediate 685 686 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 687 @*/ 688 PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value) 689 { 690 const PetscInt *points; 691 PetscInt n, p; 692 PetscErrorCode ierr; 693 694 PetscFunctionBegin; 695 PetscValidHeaderSpecific(is, IS_CLASSID, 2); 696 ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 697 ierr = ISGetIndices(is, &points);CHKERRQ(ierr); 698 for (p = 0; p < n; ++p) {ierr = DMLabelSetValue(label, points[p], value);CHKERRQ(ierr);} 699 ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr); 700 PetscFunctionReturn(0); 701 } 702 703 #undef __FUNCT__ 704 #define __FUNCT__ "DMLabelGetNumValues" 705 PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues) 706 { 707 PetscFunctionBegin; 708 PetscValidPointer(numValues, 2); 709 *numValues = label->numStrata; 710 PetscFunctionReturn(0); 711 } 712 713 #undef __FUNCT__ 714 #define __FUNCT__ "DMLabelGetValueIS" 715 PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values) 716 { 717 PetscErrorCode ierr; 718 719 PetscFunctionBegin; 720 PetscValidPointer(values, 2); 721 ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr); 722 PetscFunctionReturn(0); 723 } 724 725 #undef __FUNCT__ 726 #define __FUNCT__ "DMLabelHasStratum" 727 PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists) 728 { 729 PetscInt v; 730 731 PetscFunctionBegin; 732 PetscValidPointer(exists, 3); 733 *exists = PETSC_FALSE; 734 for (v = 0; v < label->numStrata; ++v) { 735 if (label->stratumValues[v] == value) { 736 *exists = PETSC_TRUE; 737 break; 738 } 739 } 740 PetscFunctionReturn(0); 741 } 742 743 #undef __FUNCT__ 744 #define __FUNCT__ "DMLabelGetStratumSize" 745 PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size) 746 { 747 PetscInt v; 748 PetscErrorCode ierr; 749 750 PetscFunctionBegin; 751 PetscValidPointer(size, 3); 752 *size = 0; 753 for (v = 0; v < label->numStrata; ++v) { 754 if (label->stratumValues[v] == value) { 755 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 756 *size = label->stratumSizes[v]; 757 break; 758 } 759 } 760 PetscFunctionReturn(0); 761 } 762 763 #undef __FUNCT__ 764 #define __FUNCT__ "DMLabelGetStratumBounds" 765 PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end) 766 { 767 PetscInt v; 768 PetscErrorCode ierr; 769 770 PetscFunctionBegin; 771 if (start) {PetscValidPointer(start, 3); *start = 0;} 772 if (end) {PetscValidPointer(end, 4); *end = 0;} 773 for (v = 0; v < label->numStrata; ++v) { 774 const PetscInt *points; 775 776 if (label->stratumValues[v] != value) continue; 777 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 778 if (label->stratumSizes[v] <= 0) break; 779 ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 780 if (start) *start = points[0]; 781 if (end) *end = points[label->stratumSizes[v]-1]+1; 782 ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 783 break; 784 } 785 PetscFunctionReturn(0); 786 } 787 788 #undef __FUNCT__ 789 #define __FUNCT__ "DMLabelGetStratumIS" 790 PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points) 791 { 792 PetscInt v; 793 PetscErrorCode ierr; 794 795 PetscFunctionBegin; 796 PetscValidPointer(points, 3); 797 *points = NULL; 798 for (v = 0; v < label->numStrata; ++v) { 799 if (label->stratumValues[v] == value) { 800 ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 801 if (label->validIS[v]) { 802 ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr); 803 *points = label->points[v]; 804 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to implement this to speedup Stratify"); 805 break; 806 } 807 } 808 PetscFunctionReturn(0); 809 } 810 811 #undef __FUNCT__ 812 #define __FUNCT__ "DMLabelSetStratumIS" 813 PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is) 814 { 815 PetscInt v, numStrata; 816 PetscErrorCode ierr; 817 818 PetscFunctionBegin; 819 numStrata = label->numStrata; 820 for (v = 0; v < numStrata; v++) { 821 if (label->stratumValues[v] == value) break; 822 } 823 if (v >= numStrata) {ierr = DMLabelAddStratum(label,value);CHKERRQ(ierr);} 824 if (is == label->points[v]) PetscFunctionReturn(0); 825 ierr = DMLabelClearStratum(label,value);CHKERRQ(ierr); 826 ierr = ISGetLocalSize(is,&(label->stratumSizes[v]));CHKERRQ(ierr); 827 label->stratumValues[v] = value; 828 label->validIS[v] = PETSC_TRUE; 829 ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 830 ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 831 if (label->bt) { 832 const PetscInt *points; 833 PetscInt p; 834 835 ierr = ISGetIndices(is,&points);CHKERRQ(ierr); 836 for (p = 0; p < label->stratumSizes[v]; ++p) { 837 const PetscInt point = points[p]; 838 839 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); 840 ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 841 } 842 } 843 label->points[v] = is; 844 PetscFunctionReturn(0); 845 } 846 847 848 #undef __FUNCT__ 849 #define __FUNCT__ "DMLabelClearStratum" 850 PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value) 851 { 852 PetscInt v; 853 PetscErrorCode ierr; 854 855 PetscFunctionBegin; 856 for (v = 0; v < label->numStrata; ++v) { 857 if (label->stratumValues[v] == value) break; 858 } 859 if (v >= label->numStrata) PetscFunctionReturn(0); 860 if (label->validIS[v]) { 861 if (label->bt) { 862 PetscInt i; 863 const PetscInt *points; 864 865 ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 866 for (i = 0; i < label->stratumSizes[v]; ++i) { 867 const PetscInt point = points[i]; 868 869 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); 870 ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 871 } 872 ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr); 873 } 874 ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 875 label->stratumSizes[v] = 0; 876 ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr); 877 ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr); 878 } else { 879 PetscHashIClear(label->ht[v]); 880 } 881 PetscFunctionReturn(0); 882 } 883 884 #undef __FUNCT__ 885 #define __FUNCT__ "DMLabelFilter" 886 PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end) 887 { 888 PetscInt v; 889 PetscErrorCode ierr; 890 891 PetscFunctionBegin; 892 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 893 label->pStart = start; 894 label->pEnd = end; 895 if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 896 /* Could squish offsets, but would only make sense if I reallocate the storage */ 897 for (v = 0; v < label->numStrata; ++v) { 898 PetscInt off, q; 899 const PetscInt *points; 900 PetscInt *pointsNew = NULL; 901 902 ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 903 for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) { 904 const PetscInt point = points[q]; 905 906 if ((point < start) || (point >= end)) { 907 if (!pointsNew) { 908 ierr = PetscMalloc1(label->stratumSizes[v],&pointsNew);CHKERRQ(ierr); 909 ierr = PetscMemcpy(pointsNew,points,(size_t) off * sizeof(PetscInt));CHKERRQ(ierr); 910 } 911 continue; 912 } 913 if (pointsNew) { 914 pointsNew[off++] = point; 915 } 916 } 917 ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 918 if (pointsNew) { 919 ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 920 ierr = ISCreateGeneral(PETSC_COMM_SELF,off,pointsNew,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr); 921 ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr); 922 } 923 label->stratumSizes[v] = off; 924 } 925 ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr); 926 PetscFunctionReturn(0); 927 } 928 929 #undef __FUNCT__ 930 #define __FUNCT__ "DMLabelPermute" 931 PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew) 932 { 933 const PetscInt *perm; 934 PetscInt numValues, numPoints, v, q; 935 PetscErrorCode ierr; 936 937 PetscFunctionBegin; 938 ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 939 ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr); 940 ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr); 941 ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr); 942 ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr); 943 for (v = 0; v < numValues; ++v) { 944 const PetscInt size = (*labelNew)->stratumSizes[v]; 945 const PetscInt *points; 946 PetscInt *pointsNew; 947 948 ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 949 ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr); 950 for (q = 0; q < size; ++q) { 951 const PetscInt point = points[q]; 952 953 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); 954 pointsNew[q] = perm[point]; 955 } 956 ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 957 ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr); 958 ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr); 959 ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr); 960 ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr); 961 } 962 ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr); 963 if (label->bt) { 964 ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 965 ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr); 966 } 967 PetscFunctionReturn(0); 968 } 969 970 #undef __FUNCT__ 971 #define __FUNCT__ "DMLabelDistribute_Internal" 972 PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata) 973 { 974 MPI_Comm comm; 975 PetscInt s, l, nroots, nleaves, dof, offset, size; 976 PetscInt *remoteOffsets, *rootStrata, *rootIdx; 977 PetscSection rootSection; 978 PetscSF labelSF; 979 PetscErrorCode ierr; 980 981 PetscFunctionBegin; 982 if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 983 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 984 /* Build a section of stratum values per point, generate the according SF 985 and distribute point-wise stratum values to leaves. */ 986 ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 987 ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 988 ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr); 989 if (label) { 990 for (s = 0; s < label->numStrata; ++s) { 991 const PetscInt *points; 992 993 ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 994 for (l = 0; l < label->stratumSizes[s]; l++) { 995 ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr); 996 ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr); 997 } 998 ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 999 } 1000 } 1001 ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 1002 /* Create a point-wise array of stratum values */ 1003 ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr); 1004 ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr); 1005 ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr); 1006 if (label) { 1007 for (s = 0; s < label->numStrata; ++s) { 1008 const PetscInt *points; 1009 1010 ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 1011 for (l = 0; l < label->stratumSizes[s]; l++) { 1012 const PetscInt p = points[l]; 1013 ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr); 1014 rootStrata[offset+rootIdx[p]++] = label->stratumValues[s]; 1015 } 1016 ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 1017 } 1018 } 1019 /* Build SF that maps label points to remote processes */ 1020 ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr); 1021 ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr); 1022 ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr); 1023 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 1024 /* Send the strata for each point over the derived SF */ 1025 ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr); 1026 ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr); 1027 ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 1028 ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 1029 /* Clean up */ 1030 ierr = PetscFree(rootStrata);CHKERRQ(ierr); 1031 ierr = PetscFree(rootIdx);CHKERRQ(ierr); 1032 ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 1033 ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr); 1034 PetscFunctionReturn(0); 1035 } 1036 1037 #undef __FUNCT__ 1038 #define __FUNCT__ "DMLabelDistribute" 1039 PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew) 1040 { 1041 MPI_Comm comm; 1042 PetscSection leafSection; 1043 PetscInt p, pStart, pEnd, s, size, dof, offset, stratum; 1044 PetscInt *leafStrata, *strataIdx; 1045 PetscInt **points; 1046 char *name; 1047 PetscInt nameSize; 1048 PetscHashI stratumHash; 1049 PETSC_UNUSED PetscHashIIter ret, iter; 1050 size_t len = 0; 1051 PetscMPIInt rank; 1052 PetscErrorCode ierr; 1053 1054 PetscFunctionBegin; 1055 if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 1056 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 1057 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1058 /* Bcast name */ 1059 if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 1060 nameSize = len; 1061 ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1062 ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1063 if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 1064 ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1065 ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 1066 ierr = PetscFree(name);CHKERRQ(ierr); 1067 /* Bcast defaultValue */ 1068 if (!rank) (*labelNew)->defaultValue = label->defaultValue; 1069 ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1070 /* Distribute stratum values over the SF and get the point mapping on the receiver */ 1071 ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr); 1072 /* Determine received stratum values and initialise new label*/ 1073 PetscHashICreate(stratumHash); 1074 ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr); 1075 for (p = 0; p < size; ++p) PetscHashIPut(stratumHash, leafStrata[p], ret, iter); 1076 PetscHashISize(stratumHash, (*labelNew)->numStrata); 1077 ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr); 1078 for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE; 1079 ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr); 1080 /* Turn leafStrata into indices rather than stratum values */ 1081 offset = 0; 1082 ierr = PetscHashIGetKeys(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr); 1083 for (p = 0; p < size; ++p) { 1084 for (s = 0; s < (*labelNew)->numStrata; ++s) { 1085 if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;} 1086 } 1087 } 1088 /* Rebuild the point strata on the receiver */ 1089 ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr); 1090 ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1091 for (p=pStart; p<pEnd; p++) { 1092 ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1093 ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1094 for (s=0; s<dof; s++) { 1095 (*labelNew)->stratumSizes[leafStrata[offset+s]]++; 1096 } 1097 } 1098 ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr); 1099 ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr); 1100 ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr); 1101 for (s = 0; s < (*labelNew)->numStrata; ++s) { 1102 PetscHashICreate((*labelNew)->ht[s]); 1103 ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr); 1104 } 1105 /* Insert points into new strata */ 1106 ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr); 1107 ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1108 for (p=pStart; p<pEnd; p++) { 1109 ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1110 ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1111 for (s=0; s<dof; s++) { 1112 stratum = leafStrata[offset+s]; 1113 points[stratum][strataIdx[stratum]++] = p; 1114 } 1115 } 1116 for (s = 0; s < (*labelNew)->numStrata; s++) { 1117 ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr); 1118 ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr); 1119 } 1120 ierr = PetscFree(points);CHKERRQ(ierr); 1121 PetscHashIDestroy(stratumHash); 1122 ierr = PetscFree(leafStrata);CHKERRQ(ierr); 1123 ierr = PetscFree(strataIdx);CHKERRQ(ierr); 1124 ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 1125 PetscFunctionReturn(0); 1126 } 1127 1128 #undef __FUNCT__ 1129 #define __FUNCT__ "DMLabelGather" 1130 /*@ 1131 DMLabelGather - Gather all label values from leafs into roots 1132 1133 Input Parameters: 1134 + label - the DMLabel 1135 . point - the Star Forest point communication map 1136 1137 Input Parameters: 1138 + label - the new DMLabel with localised leaf values 1139 1140 Level: developer 1141 1142 Note: This is the inverse operation to DMLabelDistribute. 1143 1144 .seealso: DMLabelDistribute() 1145 @*/ 1146 PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew) 1147 { 1148 MPI_Comm comm; 1149 PetscSection rootSection; 1150 PetscSF sfLabel; 1151 PetscSFNode *rootPoints, *leafPoints; 1152 PetscInt p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset; 1153 const PetscInt *rootDegree, *ilocal; 1154 PetscInt *rootStrata; 1155 char *name; 1156 PetscInt nameSize; 1157 size_t len = 0; 1158 PetscMPIInt rank, numProcs; 1159 PetscErrorCode ierr; 1160 1161 PetscFunctionBegin; 1162 ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 1163 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1164 ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 1165 /* Bcast name */ 1166 if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 1167 nameSize = len; 1168 ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1169 ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1170 if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 1171 ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1172 ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 1173 ierr = PetscFree(name);CHKERRQ(ierr); 1174 /* Gather rank/index pairs of leaves into local roots to build 1175 an inverse, multi-rooted SF. Note that this ignores local leaf 1176 indexing due to the use of the multiSF in PetscSFGather. */ 1177 ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr); 1178 ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr); 1179 for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1; 1180 for (p = 0; p < nleaves; p++) { 1181 leafPoints[ilocal[p]].index = ilocal[p]; 1182 leafPoints[ilocal[p]].rank = rank; 1183 } 1184 ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr); 1185 ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr); 1186 for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p]; 1187 ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr); 1188 ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 1189 ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 1190 ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr); 1191 ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 1192 /* Migrate label over inverted SF to pull stratum values at leaves into roots. */ 1193 ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr); 1194 /* Rebuild the point strata on the receiver */ 1195 for (p = 0, idx = 0; p < nroots; p++) { 1196 for (d = 0; d < rootDegree[p]; d++) { 1197 ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr); 1198 ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr); 1199 for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);} 1200 } 1201 idx += rootDegree[p]; 1202 } 1203 ierr = PetscFree(leafPoints);CHKERRQ(ierr); 1204 ierr = PetscFree(rootStrata);CHKERRQ(ierr); 1205 ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 1206 ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr); 1207 PetscFunctionReturn(0); 1208 } 1209 1210 #undef __FUNCT__ 1211 #define __FUNCT__ "DMLabelConvertToSection" 1212 PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is) 1213 { 1214 IS vIS; 1215 const PetscInt *values; 1216 PetscInt *points; 1217 PetscInt nV, vS = 0, vE = 0, v, N; 1218 PetscErrorCode ierr; 1219 1220 PetscFunctionBegin; 1221 ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr); 1222 ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr); 1223 ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr); 1224 if (nV) {vS = values[0]; vE = values[0]+1;} 1225 for (v = 1; v < nV; ++v) { 1226 vS = PetscMin(vS, values[v]); 1227 vE = PetscMax(vE, values[v]+1); 1228 } 1229 ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr); 1230 ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr); 1231 for (v = 0; v < nV; ++v) { 1232 PetscInt n; 1233 1234 ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr); 1235 ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr); 1236 } 1237 ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 1238 ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr); 1239 ierr = PetscMalloc1(N, &points);CHKERRQ(ierr); 1240 for (v = 0; v < nV; ++v) { 1241 IS is; 1242 const PetscInt *spoints; 1243 PetscInt dof, off, p; 1244 1245 ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr); 1246 ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr); 1247 ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr); 1248 ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr); 1249 for (p = 0; p < dof; ++p) points[off+p] = spoints[p]; 1250 ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr); 1251 ierr = ISDestroy(&is);CHKERRQ(ierr); 1252 } 1253 ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr); 1254 ierr = ISDestroy(&vIS);CHKERRQ(ierr); 1255 ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr); 1256 PetscFunctionReturn(0); 1257 } 1258 1259 #undef __FUNCT__ 1260 #define __FUNCT__ "PetscSectionCreateGlobalSectionLabel" 1261 /*@C 1262 PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using 1263 the local section and an SF describing the section point overlap. 1264 1265 Input Parameters: 1266 + s - The PetscSection for the local field layout 1267 . sf - The SF describing parallel layout of the section points 1268 . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs 1269 . label - The label specifying the points 1270 - labelValue - The label stratum specifying the points 1271 1272 Output Parameter: 1273 . gsection - The PetscSection for the global field layout 1274 1275 Note: This gives negative sizes and offsets to points not owned by this process 1276 1277 Level: developer 1278 1279 .seealso: PetscSectionCreate() 1280 @*/ 1281 PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection) 1282 { 1283 PetscInt *neg = NULL, *tmpOff = NULL; 1284 PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots; 1285 PetscErrorCode ierr; 1286 1287 PetscFunctionBegin; 1288 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr); 1289 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 1290 ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr); 1291 ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 1292 if (nroots >= 0) { 1293 if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart); 1294 ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr); 1295 if (nroots > pEnd-pStart) { 1296 ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr); 1297 } else { 1298 tmpOff = &(*gsection)->atlasDof[-pStart]; 1299 } 1300 } 1301 /* Mark ghost points with negative dof */ 1302 for (p = pStart; p < pEnd; ++p) { 1303 PetscInt value; 1304 1305 ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr); 1306 if (value != labelValue) continue; 1307 ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr); 1308 ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr); 1309 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 1310 if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);} 1311 if (neg) neg[p] = -(dof+1); 1312 } 1313 ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr); 1314 if (nroots >= 0) { 1315 ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1316 ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1317 if (nroots > pEnd-pStart) { 1318 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];} 1319 } 1320 } 1321 /* Calculate new sizes, get proccess offset, and calculate point offsets */ 1322 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1323 cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0; 1324 (*gsection)->atlasOff[p] = off; 1325 off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0; 1326 } 1327 ierr = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr); 1328 globalOff -= off; 1329 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1330 (*gsection)->atlasOff[p] += globalOff; 1331 if (neg) neg[p] = -((*gsection)->atlasOff[p]+1); 1332 } 1333 /* Put in negative offsets for ghost points */ 1334 if (nroots >= 0) { 1335 ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1336 ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1337 if (nroots > pEnd-pStart) { 1338 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];} 1339 } 1340 } 1341 if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);} 1342 ierr = PetscFree(neg);CHKERRQ(ierr); 1343 PetscFunctionReturn(0); 1344 } 1345 1346 typedef struct _n_PetscSectionSym_Label 1347 { 1348 DMLabel label; 1349 PetscCopyMode *modes; 1350 PetscInt *sizes; 1351 const PetscInt ***perms; 1352 const PetscScalar ***rots; 1353 PetscInt (*minMaxOrients)[2]; 1354 PetscInt numStrata; /* numStrata is only increasing, functions as a state */ 1355 } PetscSectionSym_Label; 1356 1357 #undef __FUNCT__ 1358 #define __FUNCT__ "PetscSectionSymLabelReset" 1359 static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym) 1360 { 1361 PetscInt i, j; 1362 PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data; 1363 PetscErrorCode ierr; 1364 1365 PetscFunctionBegin; 1366 for (i = 0; i <= sl->numStrata; i++) { 1367 if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) { 1368 for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 1369 if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);} 1370 if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);} 1371 } 1372 if (sl->perms[i]) { 1373 const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]]; 1374 1375 ierr = PetscFree(perms);CHKERRQ(ierr); 1376 } 1377 if (sl->rots[i]) { 1378 const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]]; 1379 1380 ierr = PetscFree(rots);CHKERRQ(ierr); 1381 } 1382 } 1383 } 1384 ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr); 1385 ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr); 1386 sl->numStrata = 0; 1387 PetscFunctionReturn(0); 1388 } 1389 1390 #undef __FUNCT__ 1391 #define __FUNCT__ "PetscSectionSymDestroy_Label" 1392 static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym) 1393 { 1394 PetscErrorCode ierr; 1395 1396 PetscFunctionBegin; 1397 ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr); 1398 ierr = PetscFree(sym->data);CHKERRQ(ierr); 1399 PetscFunctionReturn(0); 1400 } 1401 1402 #undef __FUNCT__ 1403 #define __FUNCT__ "PetscSectionSymView_Label" 1404 static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer) 1405 { 1406 PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data; 1407 PetscBool isAscii; 1408 DMLabel label = sl->label; 1409 PetscErrorCode ierr; 1410 1411 PetscFunctionBegin; 1412 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr); 1413 if (isAscii) { 1414 PetscInt i, j, k; 1415 PetscViewerFormat format; 1416 1417 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1418 if (label) { 1419 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1420 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 1421 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1422 ierr = DMLabelView(label, viewer);CHKERRQ(ierr); 1423 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1424 } else { 1425 ierr = PetscViewerASCIIPrintf(viewer," Label '%s'\n",sl->label->name);CHKERRQ(ierr); 1426 } 1427 } else { 1428 ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr); 1429 } 1430 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1431 for (i = 0; i <= sl->numStrata; i++) { 1432 PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue; 1433 1434 if (!(sl->perms[i] || sl->rots[i])) { 1435 ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr); 1436 } else { 1437 ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr); 1438 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1439 ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr); 1440 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 1441 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1442 for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 1443 if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) { 1444 ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr); 1445 } else { 1446 PetscInt tab; 1447 1448 ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr); 1449 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1450 ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr); 1451 if (sl->perms[i] && sl->perms[i][j]) { 1452 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr); 1453 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr); 1454 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);} 1455 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1456 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr); 1457 } 1458 if (sl->rots[i] && sl->rots[i][j]) { 1459 ierr = PetscViewerASCIIPrintf(viewer,"Rotations: ");CHKERRQ(ierr); 1460 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr); 1461 #if defined(PETSC_USE_COMPLEX) 1462 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f+i*%+f",PetscRealPart(sl->rots[i][j][k]),PetscImaginaryPart(sl->rots[i][j][k]));CHKERRQ(ierr);} 1463 #else 1464 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);} 1465 #endif 1466 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1467 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr); 1468 } 1469 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1470 } 1471 } 1472 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1473 } 1474 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1475 } 1476 } 1477 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1478 } 1479 PetscFunctionReturn(0); 1480 } 1481 1482 #undef __FUNCT__ 1483 #define __FUNCT__ "PetscSectionSymLabelSetLabel" 1484 /*@ 1485 PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries 1486 1487 Logically collective on sym 1488 1489 Input parameters: 1490 + sym - the section symmetries 1491 - label - the DMLabel describing the types of points 1492 1493 Level: developer: 1494 1495 .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms() 1496 @*/ 1497 PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label) 1498 { 1499 PetscSectionSym_Label *sl; 1500 PetscErrorCode ierr; 1501 1502 PetscFunctionBegin; 1503 PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1); 1504 sl = (PetscSectionSym_Label *) sym->data; 1505 if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);} 1506 if (label) { 1507 label->refct++; 1508 sl->label = label; 1509 ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr); 1510 ierr = PetscCalloc5(sl->numStrata+1,&sl->modes,sl->numStrata+1,&sl->sizes,sl->numStrata+1,&sl->perms,sl->numStrata+1,&sl->rots,sl->numStrata+1,&sl->minMaxOrients);CHKERRQ(ierr); 1511 } 1512 PetscFunctionReturn(0); 1513 } 1514 1515 #undef __FUNCT__ 1516 #define __FUNCT__ "PetscSectionSymLabelSetStratum" 1517 /*@C 1518 PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum 1519 1520 Logically collective on PetscSectionSym 1521 1522 InputParameters: 1523 + sys - the section symmetries 1524 . stratum - the stratum value in the label that we are assigning symmetries for 1525 . size - the number of dofs for points in the stratum of the label 1526 . minOrient - the smallest orientation for a point in this stratum 1527 . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient)) 1528 . mode - how sym should copy the perms and rots arrays 1529 . perms - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation. A NULL permutation is the identity 1530 + rots - NULL if there are no rotations, or (maxOrient - minOrient) sets of rotations, one for each orientation. A NULL set of orientations is the identity 1531 1532 Level: developer 1533 1534 .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel() 1535 @*/ 1536 PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots) 1537 { 1538 PetscInt i, j, k; 1539 PetscSectionSym_Label *sl; 1540 PetscErrorCode ierr; 1541 1542 PetscFunctionBegin; 1543 PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1); 1544 sl = (PetscSectionSym_Label *) sym->data; 1545 if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet"); 1546 for (i = 0; i <= sl->numStrata; i++) { 1547 PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue; 1548 1549 if (stratum == value) break; 1550 } 1551 if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,sl->label->name); 1552 sl->sizes[i] = size; 1553 sl->modes[i] = mode; 1554 sl->minMaxOrients[i][0] = minOrient; 1555 sl->minMaxOrients[i][1] = maxOrient; 1556 if (mode == PETSC_COPY_VALUES) { 1557 if (perms) { 1558 PetscInt **ownPerms; 1559 1560 ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr); 1561 for (j = 0; j < maxOrient-minOrient; j++) { 1562 if (perms[j]) { 1563 ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr); 1564 for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];} 1565 } 1566 } 1567 sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient]; 1568 } 1569 if (rots) { 1570 PetscScalar **ownRots; 1571 1572 ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr); 1573 for (j = 0; j < maxOrient-minOrient; j++) { 1574 if (rots[j]) { 1575 ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr); 1576 for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];} 1577 } 1578 } 1579 sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient]; 1580 } 1581 } else { 1582 sl->perms[i] = perms ? &perms[-minOrient] : NULL; 1583 sl->rots[i] = rots ? &rots[-minOrient] : NULL; 1584 } 1585 PetscFunctionReturn(0); 1586 } 1587 1588 #undef __FUNCT__ 1589 #define __FUNCT__ "PetscSectionSymGetPoints_Label" 1590 static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots) 1591 { 1592 PetscInt i, j, numStrata; 1593 PetscSectionSym_Label *sl; 1594 DMLabel label; 1595 PetscErrorCode ierr; 1596 1597 PetscFunctionBegin; 1598 sl = (PetscSectionSym_Label *) sym->data; 1599 numStrata = sl->numStrata; 1600 label = sl->label; 1601 for (i = 0; i < numPoints; i++) { 1602 PetscInt point = points[2*i]; 1603 PetscInt ornt = points[2*i+1]; 1604 1605 for (j = 0; j < numStrata; j++) { 1606 if (label->validIS[j]) { 1607 PetscInt k; 1608 1609 ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr); 1610 if (k >= 0) break; 1611 } else { 1612 PetscBool has; 1613 1614 PetscHashIHasKey(label->ht[j], point, has); 1615 if (has) break; 1616 } 1617 } 1618 if ((sl->minMaxOrients[j][1] > sl->minMaxOrients[j][0]) && (ornt < sl->minMaxOrients[j][0] || ornt >= sl->minMaxOrients[j][1])) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D orientation %D not in range [%D, %D) for stratum %D",point,ornt,sl->minMaxOrients[j][0],sl->minMaxOrients[j][1],j < numStrata ? label->stratumValues[j] : label->defaultValue); 1619 if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;} 1620 if (rots) {rots[i] = sl->rots[j] ? sl->rots[j][ornt] : NULL;} 1621 } 1622 PetscFunctionReturn(0); 1623 } 1624 1625 #undef __FUNCT__ 1626 #define __FUNCT__ "PetscSectionSymCreate_Label" 1627 PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym) 1628 { 1629 PetscSectionSym_Label *sl; 1630 PetscErrorCode ierr; 1631 1632 PetscFunctionBegin; 1633 ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr); 1634 sym->ops->getpoints = PetscSectionSymGetPoints_Label; 1635 sym->ops->view = PetscSectionSymView_Label; 1636 sym->ops->destroy = PetscSectionSymDestroy_Label; 1637 sym->data = (void *) sl; 1638 PetscFunctionReturn(0); 1639 } 1640 1641 #undef __FUNCT__ 1642 #define __FUNCT__ "PetscSectionSymCreateLabel" 1643 /*@ 1644 PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label 1645 1646 Collective 1647 1648 Input Parameters: 1649 + comm - the MPI communicator for the new symmetry 1650 - label - the label defining the strata 1651 1652 Output Parameters: 1653 . sym - the section symmetries 1654 1655 Level: developer 1656 1657 .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms() 1658 @*/ 1659 PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym) 1660 { 1661 PetscErrorCode ierr; 1662 1663 PetscFunctionBegin; 1664 ierr = DMInitializePackage();CHKERRQ(ierr); 1665 ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr); 1666 ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr); 1667 ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr); 1668 PetscFunctionReturn(0); 1669 } 1670