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