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