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