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