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