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