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