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