1c58f1c22SToby Isaac #include <petsc/private/dmlabelimpl.h> /*I "petscdmlabel.h" I*/ 2c58f1c22SToby Isaac #include <petsc/private/isimpl.h> /*I "petscis.h" I*/ 3c58f1c22SToby Isaac #include <petscsf.h> 4c58f1c22SToby Isaac 5c58f1c22SToby Isaac #undef __FUNCT__ 6c58f1c22SToby Isaac #define __FUNCT__ "DMLabelCreate" 7c58f1c22SToby Isaac /*@C 8c58f1c22SToby Isaac DMLabelCreate - Create a DMLabel object, which is a multimap 9c58f1c22SToby Isaac 10c58f1c22SToby Isaac Input parameter: 11c58f1c22SToby Isaac . name - The label name 12c58f1c22SToby Isaac 13c58f1c22SToby Isaac Output parameter: 14c58f1c22SToby Isaac . label - The DMLabel 15c58f1c22SToby Isaac 16c58f1c22SToby Isaac Level: beginner 17c58f1c22SToby Isaac 18c58f1c22SToby Isaac .seealso: DMLabelDestroy() 19c58f1c22SToby Isaac @*/ 20c58f1c22SToby Isaac PetscErrorCode DMLabelCreate(const char name[], DMLabel *label) 21c58f1c22SToby Isaac { 22c58f1c22SToby Isaac PetscErrorCode ierr; 23c58f1c22SToby Isaac 24c58f1c22SToby Isaac PetscFunctionBegin; 25c58f1c22SToby Isaac ierr = PetscNew(label);CHKERRQ(ierr); 26c58f1c22SToby Isaac ierr = PetscStrallocpy(name, &(*label)->name);CHKERRQ(ierr); 27c58f1c22SToby Isaac 28c58f1c22SToby Isaac (*label)->refct = 1; 29c58f1c22SToby Isaac (*label)->state = -1; 30c58f1c22SToby Isaac (*label)->numStrata = 0; 315aa44df4SToby Isaac (*label)->defaultValue = -1; 32c58f1c22SToby Isaac (*label)->stratumValues = NULL; 33ad8374ffSToby Isaac (*label)->validIS = NULL; 34c58f1c22SToby Isaac (*label)->stratumSizes = NULL; 35c58f1c22SToby Isaac (*label)->points = NULL; 36c58f1c22SToby Isaac (*label)->ht = NULL; 37c58f1c22SToby Isaac (*label)->pStart = -1; 38c58f1c22SToby Isaac (*label)->pEnd = -1; 39c58f1c22SToby Isaac (*label)->bt = NULL; 40c58f1c22SToby Isaac PetscFunctionReturn(0); 41c58f1c22SToby Isaac } 42c58f1c22SToby Isaac 43c58f1c22SToby Isaac #undef __FUNCT__ 44c58f1c22SToby Isaac #define __FUNCT__ "DMLabelMakeValid_Private" 45c58f1c22SToby Isaac /* 46c58f1c22SToby Isaac DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format 47c58f1c22SToby Isaac 48c58f1c22SToby Isaac Input parameter: 49c58f1c22SToby Isaac + label - The DMLabel 50c58f1c22SToby Isaac - v - The stratum value 51c58f1c22SToby Isaac 52c58f1c22SToby Isaac Output parameter: 53c58f1c22SToby Isaac . label - The DMLabel with stratum in sorted list format 54c58f1c22SToby Isaac 55c58f1c22SToby Isaac Level: developer 56c58f1c22SToby Isaac 57c58f1c22SToby Isaac .seealso: DMLabelCreate() 58c58f1c22SToby Isaac */ 59c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v) 60c58f1c22SToby Isaac { 61c58f1c22SToby Isaac PetscInt off; 62ad8374ffSToby Isaac PetscInt *pointArray; 63c58f1c22SToby Isaac PetscErrorCode ierr; 64c58f1c22SToby Isaac 65ad8374ffSToby Isaac if (label->validIS[v]) return 0; 66c58f1c22SToby Isaac if (v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeValid_Private\n", v); 67c58f1c22SToby Isaac PetscFunctionBegin; 68c58f1c22SToby Isaac PetscHashISize(label->ht[v], label->stratumSizes[v]); 69c58f1c22SToby Isaac 70ad8374ffSToby Isaac ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr); 71c58f1c22SToby Isaac off = 0; 72ad8374ffSToby Isaac ierr = PetscHashIGetKeys(label->ht[v], &off, pointArray);CHKERRQ(ierr); 73c58f1c22SToby Isaac if (off != label->stratumSizes[v]) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of contributed points %D from value %D should be %D", off, label->stratumValues[v], label->stratumSizes[v]); 74c58f1c22SToby Isaac PetscHashIClear(label->ht[v]); 75ad8374ffSToby Isaac ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr); 76c58f1c22SToby Isaac if (label->bt) { 77c58f1c22SToby Isaac PetscInt p; 78c58f1c22SToby Isaac 79c58f1c22SToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) { 80ad8374ffSToby Isaac const PetscInt point = pointArray[p]; 81c58f1c22SToby Isaac 82c58f1c22SToby Isaac 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); 83c58f1c22SToby Isaac ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 84c58f1c22SToby Isaac } 85c58f1c22SToby Isaac } 86ad8374ffSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,label->stratumSizes[v],pointArray,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr); 87ad8374ffSToby Isaac ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr); 88ad8374ffSToby Isaac label->validIS[v] = PETSC_TRUE; 89c58f1c22SToby Isaac ++label->state; 90c58f1c22SToby Isaac PetscFunctionReturn(0); 91c58f1c22SToby Isaac } 92c58f1c22SToby Isaac 93c58f1c22SToby Isaac #undef __FUNCT__ 94c58f1c22SToby Isaac #define __FUNCT__ "DMLabelMakeAllValid_Private" 95c58f1c22SToby Isaac /* 96c58f1c22SToby Isaac DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format 97c58f1c22SToby Isaac 98c58f1c22SToby Isaac Input parameter: 99c58f1c22SToby Isaac . label - The DMLabel 100c58f1c22SToby Isaac 101c58f1c22SToby Isaac Output parameter: 102c58f1c22SToby Isaac . label - The DMLabel with all strata in sorted list format 103c58f1c22SToby Isaac 104c58f1c22SToby Isaac Level: developer 105c58f1c22SToby Isaac 106c58f1c22SToby Isaac .seealso: DMLabelCreate() 107c58f1c22SToby Isaac */ 108c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label) 109c58f1c22SToby Isaac { 110c58f1c22SToby Isaac PetscInt v; 111c58f1c22SToby Isaac PetscErrorCode ierr; 112c58f1c22SToby Isaac 113c58f1c22SToby Isaac PetscFunctionBegin; 114c58f1c22SToby Isaac for (v = 0; v < label->numStrata; v++){ 115c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 116c58f1c22SToby Isaac } 117c58f1c22SToby Isaac PetscFunctionReturn(0); 118c58f1c22SToby Isaac } 119c58f1c22SToby Isaac 120c58f1c22SToby Isaac #undef __FUNCT__ 121c58f1c22SToby Isaac #define __FUNCT__ "DMLabelMakeInvalid_Private" 122c58f1c22SToby Isaac /* 123c58f1c22SToby Isaac DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format 124c58f1c22SToby Isaac 125c58f1c22SToby Isaac Input parameter: 126c58f1c22SToby Isaac + label - The DMLabel 127c58f1c22SToby Isaac - v - The stratum value 128c58f1c22SToby Isaac 129c58f1c22SToby Isaac Output parameter: 130c58f1c22SToby Isaac . label - The DMLabel with stratum in hash format 131c58f1c22SToby Isaac 132c58f1c22SToby Isaac Level: developer 133c58f1c22SToby Isaac 134c58f1c22SToby Isaac .seealso: DMLabelCreate() 135c58f1c22SToby Isaac */ 136c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v) 137c58f1c22SToby Isaac { 138c58f1c22SToby Isaac PETSC_UNUSED PetscHashIIter ret, iter; 139c58f1c22SToby Isaac PetscInt p; 140ad8374ffSToby Isaac const PetscInt *points; 141c58f1c22SToby Isaac PetscErrorCode ierr; 142c58f1c22SToby Isaac 143c58f1c22SToby Isaac PetscFunctionBegin; 144ad8374ffSToby Isaac if (!label->validIS[v]) PetscFunctionReturn(0); 145ad8374ffSToby Isaac if (label->points[v]) { 146ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 147ad8374ffSToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) PetscHashIPut(label->ht[v], points[p], ret, iter); 148ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 149ad8374ffSToby Isaac ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 150ad8374ffSToby Isaac } 151ad8374ffSToby Isaac label->validIS[v] = PETSC_FALSE; 152c58f1c22SToby Isaac PetscFunctionReturn(0); 153c58f1c22SToby Isaac } 154c58f1c22SToby Isaac 155c58f1c22SToby Isaac #undef __FUNCT__ 156c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetState" 157c58f1c22SToby Isaac PetscErrorCode DMLabelGetState(DMLabel label, PetscObjectState *state) 158c58f1c22SToby Isaac { 159c58f1c22SToby Isaac PetscFunctionBegin; 160c58f1c22SToby Isaac PetscValidPointer(state, 2); 161c58f1c22SToby Isaac *state = label->state; 162c58f1c22SToby Isaac PetscFunctionReturn(0); 163c58f1c22SToby Isaac } 164c58f1c22SToby Isaac 165c58f1c22SToby Isaac #undef __FUNCT__ 166c58f1c22SToby Isaac #define __FUNCT__ "DMLabelAddStratum" 167c58f1c22SToby Isaac PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value) 168c58f1c22SToby Isaac { 169ad8374ffSToby Isaac PetscInt v, *tmpV, *tmpS; 170ad8374ffSToby Isaac IS *tmpP; 171c58f1c22SToby Isaac PetscHashI *tmpH; 172c58f1c22SToby Isaac PetscBool *tmpB; 173c58f1c22SToby Isaac PetscErrorCode ierr; 174c58f1c22SToby Isaac 175c58f1c22SToby Isaac PetscFunctionBegin; 176c58f1c22SToby Isaac 177ad8374ffSToby Isaac for (v = 0; v < label->numStrata; v++) { 178ad8374ffSToby Isaac if (label->stratumValues[v] == value) PetscFunctionReturn(0); 179ad8374ffSToby Isaac } 180c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpV);CHKERRQ(ierr); 181c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpS);CHKERRQ(ierr); 182c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpH);CHKERRQ(ierr); 183c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpP);CHKERRQ(ierr); 184c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpB);CHKERRQ(ierr); 185c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 186c58f1c22SToby Isaac tmpV[v] = label->stratumValues[v]; 187c58f1c22SToby Isaac tmpS[v] = label->stratumSizes[v]; 188c58f1c22SToby Isaac tmpH[v] = label->ht[v]; 189c58f1c22SToby Isaac tmpP[v] = label->points[v]; 190ad8374ffSToby Isaac tmpB[v] = label->validIS[v]; 191c58f1c22SToby Isaac } 192c58f1c22SToby Isaac tmpV[v] = value; 193c58f1c22SToby Isaac tmpS[v] = 0; 194c58f1c22SToby Isaac PetscHashICreate(tmpH[v]); 195ad8374ffSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&tmpP[v]);CHKERRQ(ierr); 196c58f1c22SToby Isaac tmpB[v] = PETSC_TRUE; 197c58f1c22SToby Isaac ++label->numStrata; 198c58f1c22SToby Isaac ierr = PetscFree(label->stratumValues);CHKERRQ(ierr); 199c58f1c22SToby Isaac ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr); 200c58f1c22SToby Isaac ierr = PetscFree(label->ht);CHKERRQ(ierr); 201c58f1c22SToby Isaac ierr = PetscFree(label->points);CHKERRQ(ierr); 202ad8374ffSToby Isaac ierr = PetscFree(label->validIS);CHKERRQ(ierr); 203c58f1c22SToby Isaac label->stratumValues = tmpV; 204c58f1c22SToby Isaac label->stratumSizes = tmpS; 205c58f1c22SToby Isaac label->ht = tmpH; 206c58f1c22SToby Isaac label->points = tmpP; 207ad8374ffSToby Isaac label->validIS = tmpB; 208c58f1c22SToby Isaac 209c58f1c22SToby Isaac PetscFunctionReturn(0); 210c58f1c22SToby Isaac } 211c58f1c22SToby Isaac 212c58f1c22SToby Isaac #undef __FUNCT__ 213c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetName" 214c58f1c22SToby Isaac /*@C 215c58f1c22SToby Isaac DMLabelGetName - Return the name of a DMLabel object 216c58f1c22SToby Isaac 217c58f1c22SToby Isaac Input parameter: 218c58f1c22SToby Isaac . label - The DMLabel 219c58f1c22SToby Isaac 220c58f1c22SToby Isaac Output parameter: 221c58f1c22SToby Isaac . name - The label name 222c58f1c22SToby Isaac 223c58f1c22SToby Isaac Level: beginner 224c58f1c22SToby Isaac 225c58f1c22SToby Isaac .seealso: DMLabelCreate() 226c58f1c22SToby Isaac @*/ 227c58f1c22SToby Isaac PetscErrorCode DMLabelGetName(DMLabel label, const char **name) 228c58f1c22SToby Isaac { 229c58f1c22SToby Isaac PetscFunctionBegin; 230c58f1c22SToby Isaac PetscValidPointer(name, 2); 231c58f1c22SToby Isaac *name = label->name; 232c58f1c22SToby Isaac PetscFunctionReturn(0); 233c58f1c22SToby Isaac } 234c58f1c22SToby Isaac 235c58f1c22SToby Isaac #undef __FUNCT__ 236c58f1c22SToby Isaac #define __FUNCT__ "DMLabelView_Ascii" 237c58f1c22SToby Isaac static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer) 238c58f1c22SToby Isaac { 239c58f1c22SToby Isaac PetscInt v; 240c58f1c22SToby Isaac PetscMPIInt rank; 241c58f1c22SToby Isaac PetscErrorCode ierr; 242c58f1c22SToby Isaac 243c58f1c22SToby Isaac PetscFunctionBegin; 244c58f1c22SToby Isaac ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr); 245c58f1c22SToby Isaac ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 246c58f1c22SToby Isaac if (label) { 247c58f1c22SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", label->name);CHKERRQ(ierr); 248c58f1c22SToby Isaac if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, " Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);} 249c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 250c58f1c22SToby Isaac const PetscInt value = label->stratumValues[v]; 251ad8374ffSToby Isaac const PetscInt *points; 252c58f1c22SToby Isaac PetscInt p; 253c58f1c22SToby Isaac 254ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 255c58f1c22SToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) { 256ad8374ffSToby Isaac ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr); 257c58f1c22SToby Isaac } 258ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 259c58f1c22SToby Isaac } 260c58f1c22SToby Isaac } 261c58f1c22SToby Isaac ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 262c58f1c22SToby Isaac ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 263c58f1c22SToby Isaac PetscFunctionReturn(0); 264c58f1c22SToby Isaac } 265c58f1c22SToby Isaac 266c58f1c22SToby Isaac #undef __FUNCT__ 267c58f1c22SToby Isaac #define __FUNCT__ "DMLabelView" 268c58f1c22SToby Isaac /*@C 269c58f1c22SToby Isaac DMLabelView - View the label 270c58f1c22SToby Isaac 271c58f1c22SToby Isaac Input Parameters: 272c58f1c22SToby Isaac + label - The DMLabel 273c58f1c22SToby Isaac - viewer - The PetscViewer 274c58f1c22SToby Isaac 275c58f1c22SToby Isaac Level: intermediate 276c58f1c22SToby Isaac 277c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy() 278c58f1c22SToby Isaac @*/ 279c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer) 280c58f1c22SToby Isaac { 281c58f1c22SToby Isaac PetscBool iascii; 282c58f1c22SToby Isaac PetscErrorCode ierr; 283c58f1c22SToby Isaac 284c58f1c22SToby Isaac PetscFunctionBegin; 285c58f1c22SToby Isaac PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 286c58f1c22SToby Isaac if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 287c58f1c22SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 288c58f1c22SToby Isaac if (iascii) { 289c58f1c22SToby Isaac ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr); 290c58f1c22SToby Isaac } 291c58f1c22SToby Isaac PetscFunctionReturn(0); 292c58f1c22SToby Isaac } 293c58f1c22SToby Isaac 294c58f1c22SToby Isaac #undef __FUNCT__ 295c58f1c22SToby Isaac #define __FUNCT__ "DMLabelDestroy" 296c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label) 297c58f1c22SToby Isaac { 298c58f1c22SToby Isaac PetscInt v; 299c58f1c22SToby Isaac PetscErrorCode ierr; 300c58f1c22SToby Isaac 301c58f1c22SToby Isaac PetscFunctionBegin; 302c58f1c22SToby Isaac if (!(*label)) PetscFunctionReturn(0); 303c58f1c22SToby Isaac if (--(*label)->refct > 0) PetscFunctionReturn(0); 304c58f1c22SToby Isaac ierr = PetscFree((*label)->name);CHKERRQ(ierr); 305c58f1c22SToby Isaac ierr = PetscFree((*label)->stratumValues);CHKERRQ(ierr); 306c58f1c22SToby Isaac ierr = PetscFree((*label)->stratumSizes);CHKERRQ(ierr); 307ad8374ffSToby Isaac for (v = 0; v < (*label)->numStrata; ++v) {ierr = ISDestroy(&((*label)->points[v]));CHKERRQ(ierr);} 308c58f1c22SToby Isaac ierr = PetscFree((*label)->points);CHKERRQ(ierr); 309ad8374ffSToby Isaac ierr = PetscFree((*label)->validIS);CHKERRQ(ierr); 310c58f1c22SToby Isaac if ((*label)->ht) { 311c58f1c22SToby Isaac for (v = 0; v < (*label)->numStrata; ++v) {PetscHashIDestroy((*label)->ht[v]);} 312c58f1c22SToby Isaac ierr = PetscFree((*label)->ht);CHKERRQ(ierr); 313c58f1c22SToby Isaac } 314c58f1c22SToby Isaac ierr = PetscBTDestroy(&(*label)->bt);CHKERRQ(ierr); 315c58f1c22SToby Isaac ierr = PetscFree(*label);CHKERRQ(ierr); 316c58f1c22SToby Isaac PetscFunctionReturn(0); 317c58f1c22SToby Isaac } 318c58f1c22SToby Isaac 319c58f1c22SToby Isaac #undef __FUNCT__ 320c58f1c22SToby Isaac #define __FUNCT__ "DMLabelDuplicate" 321c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew) 322c58f1c22SToby Isaac { 323ad8374ffSToby Isaac PetscInt v; 324c58f1c22SToby Isaac PetscErrorCode ierr; 325c58f1c22SToby Isaac 326c58f1c22SToby Isaac PetscFunctionBegin; 327c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 328c58f1c22SToby Isaac ierr = PetscNew(labelnew);CHKERRQ(ierr); 329c58f1c22SToby Isaac ierr = PetscStrallocpy(label->name, &(*labelnew)->name);CHKERRQ(ierr); 330c58f1c22SToby Isaac 331c58f1c22SToby Isaac (*labelnew)->refct = 1; 332c58f1c22SToby Isaac (*labelnew)->numStrata = label->numStrata; 3335aa44df4SToby Isaac (*labelnew)->defaultValue = label->defaultValue; 334c58f1c22SToby Isaac if (label->numStrata) { 335c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr); 336c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr); 337c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr); 338c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr); 339ad8374ffSToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr); 340c58f1c22SToby Isaac /* Could eliminate unused space here */ 341c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 342c58f1c22SToby Isaac PetscHashICreate((*labelnew)->ht[v]); 343ad8374ffSToby Isaac (*labelnew)->validIS[v] = PETSC_TRUE; 344c58f1c22SToby Isaac (*labelnew)->stratumValues[v] = label->stratumValues[v]; 345c58f1c22SToby Isaac (*labelnew)->stratumSizes[v] = label->stratumSizes[v]; 346ad8374ffSToby Isaac ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr); 347ad8374ffSToby Isaac (*labelnew)->points[v] = label->points[v]; 348c58f1c22SToby Isaac } 349c58f1c22SToby Isaac } 350c58f1c22SToby Isaac (*labelnew)->pStart = -1; 351c58f1c22SToby Isaac (*labelnew)->pEnd = -1; 352c58f1c22SToby Isaac (*labelnew)->bt = NULL; 353c58f1c22SToby Isaac PetscFunctionReturn(0); 354c58f1c22SToby Isaac } 355c58f1c22SToby Isaac 356c58f1c22SToby Isaac #undef __FUNCT__ 357c58f1c22SToby Isaac #define __FUNCT__ "DMLabelCreateIndex" 358c58f1c22SToby Isaac /* This can be hooked into SetValue(), ClearValue(), etc. for updating */ 359c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd) 360c58f1c22SToby Isaac { 361c58f1c22SToby Isaac PetscInt v; 362c58f1c22SToby Isaac PetscErrorCode ierr; 363c58f1c22SToby Isaac 364c58f1c22SToby Isaac PetscFunctionBegin; 365c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 366c58f1c22SToby Isaac if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 367c58f1c22SToby Isaac label->pStart = pStart; 368c58f1c22SToby Isaac label->pEnd = pEnd; 369c58f1c22SToby Isaac ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr); 370c58f1c22SToby Isaac ierr = PetscBTMemzero(pEnd - pStart, label->bt);CHKERRQ(ierr); 371c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 372ad8374ffSToby Isaac const PetscInt *points; 373c58f1c22SToby Isaac PetscInt i; 374c58f1c22SToby Isaac 375ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 376c58f1c22SToby Isaac for (i = 0; i < label->stratumSizes[v]; ++i) { 377ad8374ffSToby Isaac const PetscInt point = points[i]; 378c58f1c22SToby Isaac 379c58f1c22SToby Isaac if ((point < pStart) || (point >= pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, pStart, pEnd); 380c58f1c22SToby Isaac ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr); 381c58f1c22SToby Isaac } 382ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 383c58f1c22SToby Isaac } 384c58f1c22SToby Isaac PetscFunctionReturn(0); 385c58f1c22SToby Isaac } 386c58f1c22SToby Isaac 387c58f1c22SToby Isaac #undef __FUNCT__ 388c58f1c22SToby Isaac #define __FUNCT__ "DMLabelDestroyIndex" 389c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label) 390c58f1c22SToby Isaac { 391c58f1c22SToby Isaac PetscErrorCode ierr; 392c58f1c22SToby Isaac 393c58f1c22SToby Isaac PetscFunctionBegin; 394c58f1c22SToby Isaac label->pStart = -1; 395c58f1c22SToby Isaac label->pEnd = -1; 396c58f1c22SToby Isaac if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 397c58f1c22SToby Isaac PetscFunctionReturn(0); 398c58f1c22SToby Isaac } 399c58f1c22SToby Isaac 400c58f1c22SToby Isaac #undef __FUNCT__ 401c58f1c22SToby Isaac #define __FUNCT__ "DMLabelHasValue" 402c58f1c22SToby Isaac /*@ 403c58f1c22SToby Isaac DMLabelHasValue - Determine whether a label assigns the value to any point 404c58f1c22SToby Isaac 405c58f1c22SToby Isaac Input Parameters: 406c58f1c22SToby Isaac + label - the DMLabel 407c58f1c22SToby Isaac - value - the value 408c58f1c22SToby Isaac 409c58f1c22SToby Isaac Output Parameter: 410c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point 411c58f1c22SToby Isaac 412c58f1c22SToby Isaac Level: developer 413c58f1c22SToby Isaac 414c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue() 415c58f1c22SToby Isaac @*/ 416c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains) 417c58f1c22SToby Isaac { 418c58f1c22SToby Isaac PetscInt v; 419c58f1c22SToby Isaac 420c58f1c22SToby Isaac PetscFunctionBegin; 421c58f1c22SToby Isaac PetscValidPointer(contains, 3); 422c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 423c58f1c22SToby Isaac if (value == label->stratumValues[v]) break; 424c58f1c22SToby Isaac } 425c58f1c22SToby Isaac *contains = (v < label->numStrata ? PETSC_TRUE : PETSC_FALSE); 426c58f1c22SToby Isaac PetscFunctionReturn(0); 427c58f1c22SToby Isaac } 428c58f1c22SToby Isaac 429c58f1c22SToby Isaac #undef __FUNCT__ 430c58f1c22SToby Isaac #define __FUNCT__ "DMLabelHasPoint" 431c58f1c22SToby Isaac /*@ 432c58f1c22SToby Isaac DMLabelHasPoint - Determine whether a label assigns a value to a point 433c58f1c22SToby Isaac 434c58f1c22SToby Isaac Input Parameters: 435c58f1c22SToby Isaac + label - the DMLabel 436c58f1c22SToby Isaac - point - the point 437c58f1c22SToby Isaac 438c58f1c22SToby Isaac Output Parameter: 439c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value 440c58f1c22SToby Isaac 441c58f1c22SToby Isaac Note: The user must call DMLabelCreateIndex() before this function. 442c58f1c22SToby Isaac 443c58f1c22SToby Isaac Level: developer 444c58f1c22SToby Isaac 445c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue() 446c58f1c22SToby Isaac @*/ 447c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains) 448c58f1c22SToby Isaac { 449c58f1c22SToby Isaac PetscErrorCode ierr; 450c58f1c22SToby Isaac 451c58f1c22SToby Isaac PetscFunctionBeginHot; 452c58f1c22SToby Isaac PetscValidPointer(contains, 3); 453c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 454c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG) 455c58f1c22SToby Isaac if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()"); 456c58f1c22SToby Isaac 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); 457c58f1c22SToby Isaac #endif 458c58f1c22SToby Isaac *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE; 459c58f1c22SToby Isaac PetscFunctionReturn(0); 460c58f1c22SToby Isaac } 461c58f1c22SToby Isaac 462c58f1c22SToby Isaac #undef __FUNCT__ 463c58f1c22SToby Isaac #define __FUNCT__ "DMLabelStratumHasPoint" 464c58f1c22SToby Isaac /*@ 465c58f1c22SToby Isaac DMLabelStratumHasPoint - Return true if the stratum contains a point 466c58f1c22SToby Isaac 467c58f1c22SToby Isaac Input Parameters: 468c58f1c22SToby Isaac + label - the DMLabel 469c58f1c22SToby Isaac . value - the stratum value 470c58f1c22SToby Isaac - point - the point 471c58f1c22SToby Isaac 472c58f1c22SToby Isaac Output Parameter: 473c58f1c22SToby Isaac . contains - true if the stratum contains the point 474c58f1c22SToby Isaac 475c58f1c22SToby Isaac Level: intermediate 476c58f1c22SToby Isaac 477c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue() 478c58f1c22SToby Isaac @*/ 479c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains) 480c58f1c22SToby Isaac { 481c58f1c22SToby Isaac PetscInt v; 482c58f1c22SToby Isaac PetscErrorCode ierr; 483c58f1c22SToby Isaac 484c58f1c22SToby Isaac PetscFunctionBegin; 485c58f1c22SToby Isaac PetscValidPointer(contains, 4); 486c58f1c22SToby Isaac *contains = PETSC_FALSE; 487c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 488c58f1c22SToby Isaac if (label->stratumValues[v] == value) { 489ad8374ffSToby Isaac if (label->validIS[v]) { 490c58f1c22SToby Isaac PetscInt i; 491c58f1c22SToby Isaac 492*a2d74346SToby Isaac ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr); 493c58f1c22SToby Isaac if (i >= 0) { 494c58f1c22SToby Isaac *contains = PETSC_TRUE; 495c58f1c22SToby Isaac break; 496c58f1c22SToby Isaac } 497c58f1c22SToby Isaac } else { 498c58f1c22SToby Isaac PetscBool has; 499c58f1c22SToby Isaac 500c58f1c22SToby Isaac PetscHashIHasKey(label->ht[v], point, has); 501c58f1c22SToby Isaac if (has) { 502c58f1c22SToby Isaac *contains = PETSC_TRUE; 503c58f1c22SToby Isaac break; 504c58f1c22SToby Isaac } 505c58f1c22SToby Isaac } 506c58f1c22SToby Isaac } 507c58f1c22SToby Isaac } 508c58f1c22SToby Isaac PetscFunctionReturn(0); 509c58f1c22SToby Isaac } 510c58f1c22SToby Isaac 511c58f1c22SToby Isaac #undef __FUNCT__ 5125aa44df4SToby Isaac #define __FUNCT__ "DMLabelGetDefaultValue" 5135aa44df4SToby Isaac /* 5145aa44df4SToby Isaac DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 5155aa44df4SToby Isaac When a label is created, it is initialized to -1. 5165aa44df4SToby Isaac 5175aa44df4SToby Isaac Input parameter: 5185aa44df4SToby Isaac . label - a DMLabel object 5195aa44df4SToby Isaac 5205aa44df4SToby Isaac Output parameter: 5215aa44df4SToby Isaac . defaultValue - the default value 5225aa44df4SToby Isaac 5235aa44df4SToby Isaac Level: beginner 5245aa44df4SToby Isaac 5255aa44df4SToby Isaac .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 5265aa44df4SToby Isaac */ 5275aa44df4SToby Isaac PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue) 5285aa44df4SToby Isaac { 5295aa44df4SToby Isaac PetscFunctionBegin; 5305aa44df4SToby Isaac *defaultValue = label->defaultValue; 5315aa44df4SToby Isaac PetscFunctionReturn(0); 5325aa44df4SToby Isaac } 5335aa44df4SToby Isaac 5345aa44df4SToby Isaac #undef __FUNCT__ 5355aa44df4SToby Isaac #define __FUNCT__ "DMLabelSetDefaultValue" 5365aa44df4SToby Isaac /* 5375aa44df4SToby Isaac DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 5385aa44df4SToby Isaac When a label is created, it is initialized to -1. 5395aa44df4SToby Isaac 5405aa44df4SToby Isaac Input parameter: 5415aa44df4SToby Isaac . label - a DMLabel object 5425aa44df4SToby Isaac 5435aa44df4SToby Isaac Output parameter: 5445aa44df4SToby Isaac . defaultValue - the default value 5455aa44df4SToby Isaac 5465aa44df4SToby Isaac Level: beginner 5475aa44df4SToby Isaac 5485aa44df4SToby Isaac .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 5495aa44df4SToby Isaac */ 5505aa44df4SToby Isaac PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue) 5515aa44df4SToby Isaac { 5525aa44df4SToby Isaac PetscFunctionBegin; 5535aa44df4SToby Isaac label->defaultValue = defaultValue; 5545aa44df4SToby Isaac PetscFunctionReturn(0); 5555aa44df4SToby Isaac } 5565aa44df4SToby Isaac 5575aa44df4SToby Isaac #undef __FUNCT__ 558c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetValue" 559c58f1c22SToby Isaac /*@ 5605aa44df4SToby Isaac 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()) 561c58f1c22SToby Isaac 562c58f1c22SToby Isaac Input Parameters: 563c58f1c22SToby Isaac + label - the DMLabel 564c58f1c22SToby Isaac - point - the point 565c58f1c22SToby Isaac 566c58f1c22SToby Isaac Output Parameter: 567c58f1c22SToby Isaac . value - The point value, or -1 568c58f1c22SToby Isaac 569c58f1c22SToby Isaac Level: intermediate 570c58f1c22SToby Isaac 5715aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 572c58f1c22SToby Isaac @*/ 573c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value) 574c58f1c22SToby Isaac { 575c58f1c22SToby Isaac PetscInt v; 576c58f1c22SToby Isaac PetscErrorCode ierr; 577c58f1c22SToby Isaac 578c58f1c22SToby Isaac PetscFunctionBegin; 579c58f1c22SToby Isaac PetscValidPointer(value, 3); 5805aa44df4SToby Isaac *value = label->defaultValue; 581c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 582ad8374ffSToby Isaac if (label->validIS[v]) { 583c58f1c22SToby Isaac PetscInt i; 584c58f1c22SToby Isaac 585*a2d74346SToby Isaac ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr); 586c58f1c22SToby Isaac if (i >= 0) { 587c58f1c22SToby Isaac *value = label->stratumValues[v]; 588c58f1c22SToby Isaac break; 589c58f1c22SToby Isaac } 590c58f1c22SToby Isaac } else { 591c58f1c22SToby Isaac PetscBool has; 592c58f1c22SToby Isaac 593c58f1c22SToby Isaac PetscHashIHasKey(label->ht[v], point, has); 594c58f1c22SToby Isaac if (has) { 595c58f1c22SToby Isaac *value = label->stratumValues[v]; 596c58f1c22SToby Isaac break; 597c58f1c22SToby Isaac } 598c58f1c22SToby Isaac } 599c58f1c22SToby Isaac } 600c58f1c22SToby Isaac PetscFunctionReturn(0); 601c58f1c22SToby Isaac } 602c58f1c22SToby Isaac 603c58f1c22SToby Isaac #undef __FUNCT__ 604c58f1c22SToby Isaac #define __FUNCT__ "DMLabelSetValue" 605c58f1c22SToby Isaac /*@ 6065aa44df4SToby Isaac 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 somethingg different), then this function will do nothing. 607c58f1c22SToby Isaac 608c58f1c22SToby Isaac Input Parameters: 609c58f1c22SToby Isaac + label - the DMLabel 610c58f1c22SToby Isaac . point - the point 611c58f1c22SToby Isaac - value - The point value 612c58f1c22SToby Isaac 613c58f1c22SToby Isaac Level: intermediate 614c58f1c22SToby Isaac 6155aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 616c58f1c22SToby Isaac @*/ 617c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value) 618c58f1c22SToby Isaac { 619c58f1c22SToby Isaac PETSC_UNUSED PetscHashIIter iter, ret; 620c58f1c22SToby Isaac PetscInt v; 621c58f1c22SToby Isaac PetscErrorCode ierr; 622c58f1c22SToby Isaac 623c58f1c22SToby Isaac PetscFunctionBegin; 624c58f1c22SToby Isaac /* Find, or add, label value */ 6255aa44df4SToby Isaac if (value == label->defaultValue) PetscFunctionReturn(0); 626c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 627c58f1c22SToby Isaac if (label->stratumValues[v] == value) break; 628c58f1c22SToby Isaac } 629c58f1c22SToby Isaac /* Create new table */ 630c58f1c22SToby Isaac if (v >= label->numStrata) {ierr = DMLabelAddStratum(label, value);CHKERRQ(ierr);} 631c58f1c22SToby Isaac ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 632c58f1c22SToby Isaac /* Set key */ 633c58f1c22SToby Isaac PetscHashIPut(label->ht[v], point, ret, iter); 634c58f1c22SToby Isaac PetscFunctionReturn(0); 635c58f1c22SToby Isaac } 636c58f1c22SToby Isaac 637c58f1c22SToby Isaac #undef __FUNCT__ 638c58f1c22SToby Isaac #define __FUNCT__ "DMLabelClearValue" 639c58f1c22SToby Isaac /*@ 640c58f1c22SToby Isaac DMLabelClearValue - Clear the value a label assigns to a point 641c58f1c22SToby Isaac 642c58f1c22SToby Isaac Input Parameters: 643c58f1c22SToby Isaac + label - the DMLabel 644c58f1c22SToby Isaac . point - the point 645c58f1c22SToby Isaac - value - The point value 646c58f1c22SToby Isaac 647c58f1c22SToby Isaac Level: intermediate 648c58f1c22SToby Isaac 649c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue() 650c58f1c22SToby Isaac @*/ 651c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value) 652c58f1c22SToby Isaac { 653ad8374ffSToby Isaac PetscInt v; 654c58f1c22SToby Isaac PetscErrorCode ierr; 655c58f1c22SToby Isaac 656c58f1c22SToby Isaac PetscFunctionBegin; 657c58f1c22SToby Isaac /* Find label value */ 658c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 659c58f1c22SToby Isaac if (label->stratumValues[v] == value) break; 660c58f1c22SToby Isaac } 661c58f1c22SToby Isaac if (v >= label->numStrata) PetscFunctionReturn(0); 662ad8374ffSToby Isaac if (label->validIS[v]) {ierr = DMLabelMakeInvalid_Private(label,v);CHKERRQ(ierr);} 663c58f1c22SToby Isaac ierr = PetscHashIDelKey(label->ht[v], point);CHKERRQ(ierr); 664c58f1c22SToby Isaac PetscFunctionReturn(0); 665c58f1c22SToby Isaac } 666c58f1c22SToby Isaac 667c58f1c22SToby Isaac #undef __FUNCT__ 668c58f1c22SToby Isaac #define __FUNCT__ "DMLabelInsertIS" 669c58f1c22SToby Isaac /*@ 670c58f1c22SToby Isaac DMLabelInsertIS - Set all points in the IS to a value 671c58f1c22SToby Isaac 672c58f1c22SToby Isaac Input Parameters: 673c58f1c22SToby Isaac + label - the DMLabel 674c58f1c22SToby Isaac . is - the point IS 675c58f1c22SToby Isaac - value - The point value 676c58f1c22SToby Isaac 677c58f1c22SToby Isaac Level: intermediate 678c58f1c22SToby Isaac 679c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 680c58f1c22SToby Isaac @*/ 681c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value) 682c58f1c22SToby Isaac { 683c58f1c22SToby Isaac const PetscInt *points; 684c58f1c22SToby Isaac PetscInt n, p; 685c58f1c22SToby Isaac PetscErrorCode ierr; 686c58f1c22SToby Isaac 687c58f1c22SToby Isaac PetscFunctionBegin; 688c58f1c22SToby Isaac PetscValidHeaderSpecific(is, IS_CLASSID, 2); 689c58f1c22SToby Isaac ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 690c58f1c22SToby Isaac ierr = ISGetIndices(is, &points);CHKERRQ(ierr); 691c58f1c22SToby Isaac for (p = 0; p < n; ++p) {ierr = DMLabelSetValue(label, points[p], value);CHKERRQ(ierr);} 692c58f1c22SToby Isaac ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr); 693c58f1c22SToby Isaac PetscFunctionReturn(0); 694c58f1c22SToby Isaac } 695c58f1c22SToby Isaac 696c58f1c22SToby Isaac #undef __FUNCT__ 697c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetNumValues" 698c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues) 699c58f1c22SToby Isaac { 700c58f1c22SToby Isaac PetscFunctionBegin; 701c58f1c22SToby Isaac PetscValidPointer(numValues, 2); 702c58f1c22SToby Isaac *numValues = label->numStrata; 703c58f1c22SToby Isaac PetscFunctionReturn(0); 704c58f1c22SToby Isaac } 705c58f1c22SToby Isaac 706c58f1c22SToby Isaac #undef __FUNCT__ 707c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetValueIS" 708c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values) 709c58f1c22SToby Isaac { 710c58f1c22SToby Isaac PetscErrorCode ierr; 711c58f1c22SToby Isaac 712c58f1c22SToby Isaac PetscFunctionBegin; 713c58f1c22SToby Isaac PetscValidPointer(values, 2); 714c58f1c22SToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr); 715c58f1c22SToby Isaac PetscFunctionReturn(0); 716c58f1c22SToby Isaac } 717c58f1c22SToby Isaac 718c58f1c22SToby Isaac #undef __FUNCT__ 719fada774cSMatthew G. Knepley #define __FUNCT__ "DMLabelHasStratum" 720fada774cSMatthew G. Knepley PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists) 721fada774cSMatthew G. Knepley { 722fada774cSMatthew G. Knepley PetscInt v; 723fada774cSMatthew G. Knepley 724fada774cSMatthew G. Knepley PetscFunctionBegin; 725fada774cSMatthew G. Knepley PetscValidPointer(exists, 3); 726fada774cSMatthew G. Knepley *exists = PETSC_FALSE; 727fada774cSMatthew G. Knepley for (v = 0; v < label->numStrata; ++v) { 728fada774cSMatthew G. Knepley if (label->stratumValues[v] == value) { 729fada774cSMatthew G. Knepley *exists = PETSC_TRUE; 730fada774cSMatthew G. Knepley break; 731fada774cSMatthew G. Knepley } 732fada774cSMatthew G. Knepley } 733fada774cSMatthew G. Knepley PetscFunctionReturn(0); 734fada774cSMatthew G. Knepley } 735fada774cSMatthew G. Knepley 736fada774cSMatthew G. Knepley #undef __FUNCT__ 737c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetStratumSize" 738c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size) 739c58f1c22SToby Isaac { 740c58f1c22SToby Isaac PetscInt v; 741c58f1c22SToby Isaac PetscErrorCode ierr; 742c58f1c22SToby Isaac 743c58f1c22SToby Isaac PetscFunctionBegin; 744c58f1c22SToby Isaac PetscValidPointer(size, 3); 745c58f1c22SToby Isaac *size = 0; 746c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 747c58f1c22SToby Isaac if (label->stratumValues[v] == value) { 748c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 749c58f1c22SToby Isaac *size = label->stratumSizes[v]; 750c58f1c22SToby Isaac break; 751c58f1c22SToby Isaac } 752c58f1c22SToby Isaac } 753c58f1c22SToby Isaac PetscFunctionReturn(0); 754c58f1c22SToby Isaac } 755c58f1c22SToby Isaac 756c58f1c22SToby Isaac #undef __FUNCT__ 757c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetStratumBounds" 758c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end) 759c58f1c22SToby Isaac { 760c58f1c22SToby Isaac PetscInt v; 761c58f1c22SToby Isaac PetscErrorCode ierr; 762c58f1c22SToby Isaac 763c58f1c22SToby Isaac PetscFunctionBegin; 764c58f1c22SToby Isaac if (start) {PetscValidPointer(start, 3); *start = 0;} 765c58f1c22SToby Isaac if (end) {PetscValidPointer(end, 4); *end = 0;} 766c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 767ad8374ffSToby Isaac const PetscInt *points; 768ad8374ffSToby Isaac 769c58f1c22SToby Isaac if (label->stratumValues[v] != value) continue; 770c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 771c58f1c22SToby Isaac if (label->stratumSizes[v] <= 0) break; 772ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 773ad8374ffSToby Isaac if (start) *start = points[0]; 774ad8374ffSToby Isaac if (end) *end = points[label->stratumSizes[v]-1]+1; 775ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 776c58f1c22SToby Isaac break; 777c58f1c22SToby Isaac } 778c58f1c22SToby Isaac PetscFunctionReturn(0); 779c58f1c22SToby Isaac } 780c58f1c22SToby Isaac 781c58f1c22SToby Isaac #undef __FUNCT__ 782c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetStratumIS" 783c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points) 784c58f1c22SToby Isaac { 785c58f1c22SToby Isaac PetscInt v; 786c58f1c22SToby Isaac PetscErrorCode ierr; 787c58f1c22SToby Isaac 788c58f1c22SToby Isaac PetscFunctionBegin; 789c58f1c22SToby Isaac PetscValidPointer(points, 3); 790c58f1c22SToby Isaac *points = NULL; 791c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 792c58f1c22SToby Isaac if (label->stratumValues[v] == value) { 793c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 794ad8374ffSToby Isaac if (label->validIS[v]) { 795ad8374ffSToby Isaac ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr); 796ad8374ffSToby Isaac *points = label->points[v]; 797c58f1c22SToby Isaac } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to implement this to speedup Stratify"); 798c58f1c22SToby Isaac break; 799c58f1c22SToby Isaac } 800c58f1c22SToby Isaac } 801c58f1c22SToby Isaac PetscFunctionReturn(0); 802c58f1c22SToby Isaac } 803c58f1c22SToby Isaac 804c58f1c22SToby Isaac #undef __FUNCT__ 8054de306b1SToby Isaac #define __FUNCT__ "DMLabelSetStratumIS" 8064de306b1SToby Isaac PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is) 8074de306b1SToby Isaac { 8084de306b1SToby Isaac PetscInt v, numStrata; 8094de306b1SToby Isaac PetscErrorCode ierr; 8104de306b1SToby Isaac 8114de306b1SToby Isaac PetscFunctionBegin; 8124de306b1SToby Isaac numStrata = label->numStrata; 8134de306b1SToby Isaac for (v = 0; v < numStrata; v++) { 8144de306b1SToby Isaac if (label->stratumValues[v] == value) break; 8154de306b1SToby Isaac } 8164de306b1SToby Isaac if (v >= numStrata) {ierr = DMLabelAddStratum(label,value);CHKERRQ(ierr);} 8174de306b1SToby Isaac if (is == label->points[v]) PetscFunctionReturn(0); 8184de306b1SToby Isaac ierr = DMLabelClearStratum(label,value);CHKERRQ(ierr); 8194de306b1SToby Isaac ierr = ISGetLocalSize(is,&(label->stratumSizes[v]));CHKERRQ(ierr); 8204de306b1SToby Isaac label->stratumValues[v] = value; 8214de306b1SToby Isaac label->validIS[v] = PETSC_TRUE; 8224de306b1SToby Isaac ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 8234de306b1SToby Isaac ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 8244de306b1SToby Isaac if (label->bt) { 8254de306b1SToby Isaac const PetscInt *points; 8264de306b1SToby Isaac PetscInt p; 8274de306b1SToby Isaac 8284de306b1SToby Isaac ierr = ISGetIndices(is,&points);CHKERRQ(ierr); 8294de306b1SToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) { 8304de306b1SToby Isaac const PetscInt point = points[p]; 8314de306b1SToby Isaac 8324de306b1SToby Isaac 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); 8334de306b1SToby Isaac ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 8344de306b1SToby Isaac } 8354de306b1SToby Isaac } 8364de306b1SToby Isaac label->points[v] = is; 8374de306b1SToby Isaac PetscFunctionReturn(0); 8384de306b1SToby Isaac } 8394de306b1SToby Isaac 8404de306b1SToby Isaac 8414de306b1SToby Isaac #undef __FUNCT__ 842c58f1c22SToby Isaac #define __FUNCT__ "DMLabelClearStratum" 843c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value) 844c58f1c22SToby Isaac { 845c58f1c22SToby Isaac PetscInt v; 846c58f1c22SToby Isaac PetscErrorCode ierr; 847c58f1c22SToby Isaac 848c58f1c22SToby Isaac PetscFunctionBegin; 849c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 850c58f1c22SToby Isaac if (label->stratumValues[v] == value) break; 851c58f1c22SToby Isaac } 852c58f1c22SToby Isaac if (v >= label->numStrata) PetscFunctionReturn(0); 8534de306b1SToby Isaac if (label->validIS[v]) { 8544de306b1SToby Isaac if (label->bt) { 855c58f1c22SToby Isaac PetscInt i; 856ad8374ffSToby Isaac const PetscInt *points; 857c58f1c22SToby Isaac 858ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 859c58f1c22SToby Isaac for (i = 0; i < label->stratumSizes[v]; ++i) { 860ad8374ffSToby Isaac const PetscInt point = points[i]; 861c58f1c22SToby Isaac 862c58f1c22SToby Isaac 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); 863c58f1c22SToby Isaac ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 864c58f1c22SToby Isaac } 865ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr); 866c58f1c22SToby Isaac } 867ad8374ffSToby Isaac ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 868c58f1c22SToby Isaac label->stratumSizes[v] = 0; 869ad8374ffSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr); 870ad8374ffSToby Isaac ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr); 871c58f1c22SToby Isaac } else { 872c58f1c22SToby Isaac PetscHashIClear(label->ht[v]); 873c58f1c22SToby Isaac } 874c58f1c22SToby Isaac PetscFunctionReturn(0); 875c58f1c22SToby Isaac } 876c58f1c22SToby Isaac 877c58f1c22SToby Isaac #undef __FUNCT__ 878c58f1c22SToby Isaac #define __FUNCT__ "DMLabelFilter" 879c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end) 880c58f1c22SToby Isaac { 881c58f1c22SToby Isaac PetscInt v; 882c58f1c22SToby Isaac PetscErrorCode ierr; 883c58f1c22SToby Isaac 884c58f1c22SToby Isaac PetscFunctionBegin; 885c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 886c58f1c22SToby Isaac label->pStart = start; 887c58f1c22SToby Isaac label->pEnd = end; 888c58f1c22SToby Isaac if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 889c58f1c22SToby Isaac /* Could squish offsets, but would only make sense if I reallocate the storage */ 890c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 891c58f1c22SToby Isaac PetscInt off, q; 892ad8374ffSToby Isaac const PetscInt *points; 893ad8374ffSToby Isaac PetscInt *pointsNew = NULL; 894c58f1c22SToby Isaac 895ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 896c58f1c22SToby Isaac for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) { 897ad8374ffSToby Isaac const PetscInt point = points[q]; 898c58f1c22SToby Isaac 899ad8374ffSToby Isaac if ((point < start) || (point >= end)) { 900ad8374ffSToby Isaac if (!pointsNew) { 901ad8374ffSToby Isaac ierr = PetscMalloc1(label->stratumSizes[v],&pointsNew);CHKERRQ(ierr); 902ad8374ffSToby Isaac ierr = PetscMemcpy(pointsNew,points,(size_t) off * sizeof(PetscInt));CHKERRQ(ierr); 903ad8374ffSToby Isaac } 904ad8374ffSToby Isaac continue; 905ad8374ffSToby Isaac } 906ad8374ffSToby Isaac if (pointsNew) { 907ad8374ffSToby Isaac pointsNew[off++] = point; 908ad8374ffSToby Isaac } 909ad8374ffSToby Isaac } 910ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 911ad8374ffSToby Isaac if (pointsNew) { 912ad8374ffSToby Isaac ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 913ad8374ffSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,off,pointsNew,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr); 914ad8374ffSToby Isaac ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr); 915c58f1c22SToby Isaac } 916c58f1c22SToby Isaac label->stratumSizes[v] = off; 917c58f1c22SToby Isaac } 918c58f1c22SToby Isaac ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr); 919c58f1c22SToby Isaac PetscFunctionReturn(0); 920c58f1c22SToby Isaac } 921c58f1c22SToby Isaac 922c58f1c22SToby Isaac #undef __FUNCT__ 923c58f1c22SToby Isaac #define __FUNCT__ "DMLabelPermute" 924c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew) 925c58f1c22SToby Isaac { 926c58f1c22SToby Isaac const PetscInt *perm; 927c58f1c22SToby Isaac PetscInt numValues, numPoints, v, q; 928c58f1c22SToby Isaac PetscErrorCode ierr; 929c58f1c22SToby Isaac 930c58f1c22SToby Isaac PetscFunctionBegin; 931c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 932c58f1c22SToby Isaac ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr); 933c58f1c22SToby Isaac ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr); 934c58f1c22SToby Isaac ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr); 935c58f1c22SToby Isaac ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr); 936c58f1c22SToby Isaac for (v = 0; v < numValues; ++v) { 937c58f1c22SToby Isaac const PetscInt size = (*labelNew)->stratumSizes[v]; 938ad8374ffSToby Isaac const PetscInt *points; 939ad8374ffSToby Isaac PetscInt *pointsNew; 940c58f1c22SToby Isaac 941ad8374ffSToby Isaac ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 942ad8374ffSToby Isaac ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr); 943c58f1c22SToby Isaac for (q = 0; q < size; ++q) { 944ad8374ffSToby Isaac const PetscInt point = points[q]; 945c58f1c22SToby Isaac 946c58f1c22SToby Isaac 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); 947ad8374ffSToby Isaac pointsNew[q] = perm[point]; 948c58f1c22SToby Isaac } 949ad8374ffSToby Isaac ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 950ad8374ffSToby Isaac ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr); 951ad8374ffSToby Isaac ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr); 952ad8374ffSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr); 953ad8374ffSToby Isaac ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr); 954c58f1c22SToby Isaac } 955c58f1c22SToby Isaac ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr); 956c58f1c22SToby Isaac if (label->bt) { 957c58f1c22SToby Isaac ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 958c58f1c22SToby Isaac ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr); 959c58f1c22SToby Isaac } 960c58f1c22SToby Isaac PetscFunctionReturn(0); 961c58f1c22SToby Isaac } 962c58f1c22SToby Isaac 963c58f1c22SToby Isaac #undef __FUNCT__ 96426c55118SMichael Lange #define __FUNCT__ "DMLabelDistribute_Internal" 96526c55118SMichael Lange PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata) 96626c55118SMichael Lange { 96726c55118SMichael Lange MPI_Comm comm; 96826c55118SMichael Lange PetscInt s, l, nroots, nleaves, dof, offset, size; 96926c55118SMichael Lange PetscInt *remoteOffsets, *rootStrata, *rootIdx; 97026c55118SMichael Lange PetscSection rootSection; 97126c55118SMichael Lange PetscSF labelSF; 97226c55118SMichael Lange PetscErrorCode ierr; 97326c55118SMichael Lange 97426c55118SMichael Lange PetscFunctionBegin; 97526c55118SMichael Lange if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 97626c55118SMichael Lange ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 97726c55118SMichael Lange /* Build a section of stratum values per point, generate the according SF 97826c55118SMichael Lange and distribute point-wise stratum values to leaves. */ 97926c55118SMichael Lange ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 98026c55118SMichael Lange ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 98126c55118SMichael Lange ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr); 98226c55118SMichael Lange if (label) { 98326c55118SMichael Lange for (s = 0; s < label->numStrata; ++s) { 984ad8374ffSToby Isaac const PetscInt *points; 985ad8374ffSToby Isaac 986ad8374ffSToby Isaac ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 98726c55118SMichael Lange for (l = 0; l < label->stratumSizes[s]; l++) { 988ad8374ffSToby Isaac ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr); 989ad8374ffSToby Isaac ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr); 99026c55118SMichael Lange } 991ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 99226c55118SMichael Lange } 99326c55118SMichael Lange } 99426c55118SMichael Lange ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 99526c55118SMichael Lange /* Create a point-wise array of stratum values */ 99626c55118SMichael Lange ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr); 99726c55118SMichael Lange ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr); 99826c55118SMichael Lange ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr); 99926c55118SMichael Lange if (label) { 100026c55118SMichael Lange for (s = 0; s < label->numStrata; ++s) { 1001ad8374ffSToby Isaac const PetscInt *points; 1002ad8374ffSToby Isaac 1003ad8374ffSToby Isaac ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 100426c55118SMichael Lange for (l = 0; l < label->stratumSizes[s]; l++) { 1005ad8374ffSToby Isaac const PetscInt p = points[l]; 100626c55118SMichael Lange ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr); 100726c55118SMichael Lange rootStrata[offset+rootIdx[p]++] = label->stratumValues[s]; 100826c55118SMichael Lange } 1009ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 101026c55118SMichael Lange } 101126c55118SMichael Lange } 101226c55118SMichael Lange /* Build SF that maps label points to remote processes */ 101326c55118SMichael Lange ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr); 101426c55118SMichael Lange ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr); 101526c55118SMichael Lange ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr); 101626c55118SMichael Lange ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 101726c55118SMichael Lange /* Send the strata for each point over the derived SF */ 101826c55118SMichael Lange ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr); 101926c55118SMichael Lange ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr); 102026c55118SMichael Lange ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 102126c55118SMichael Lange ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 102226c55118SMichael Lange /* Clean up */ 102326c55118SMichael Lange ierr = PetscFree(rootStrata);CHKERRQ(ierr); 102426c55118SMichael Lange ierr = PetscFree(rootIdx);CHKERRQ(ierr); 102526c55118SMichael Lange ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 102626c55118SMichael Lange ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr); 102726c55118SMichael Lange PetscFunctionReturn(0); 102826c55118SMichael Lange } 102926c55118SMichael Lange 103026c55118SMichael Lange #undef __FUNCT__ 1031c58f1c22SToby Isaac #define __FUNCT__ "DMLabelDistribute" 1032c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew) 1033c58f1c22SToby Isaac { 1034c58f1c22SToby Isaac MPI_Comm comm; 103526c55118SMichael Lange PetscSection leafSection; 103626c55118SMichael Lange PetscInt p, pStart, pEnd, s, size, dof, offset, stratum; 103726c55118SMichael Lange PetscInt *leafStrata, *strataIdx; 1038ad8374ffSToby Isaac PetscInt **points; 1039c58f1c22SToby Isaac char *name; 1040c58f1c22SToby Isaac PetscInt nameSize; 10415cbdf6fcSMichael Lange PetscHashI stratumHash; 10425cbdf6fcSMichael Lange PETSC_UNUSED PetscHashIIter ret, iter; 1043c58f1c22SToby Isaac size_t len = 0; 104426c55118SMichael Lange PetscMPIInt rank; 1045c58f1c22SToby Isaac PetscErrorCode ierr; 1046c58f1c22SToby Isaac 1047c58f1c22SToby Isaac PetscFunctionBegin; 1048c58f1c22SToby Isaac if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 1049c58f1c22SToby Isaac ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 1050c58f1c22SToby Isaac ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1051c58f1c22SToby Isaac /* Bcast name */ 1052c58f1c22SToby Isaac if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 1053c58f1c22SToby Isaac nameSize = len; 1054c58f1c22SToby Isaac ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1055c58f1c22SToby Isaac ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1056c58f1c22SToby Isaac if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 1057c58f1c22SToby Isaac ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1058c58f1c22SToby Isaac ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 1059c58f1c22SToby Isaac ierr = PetscFree(name);CHKERRQ(ierr); 106077d236dfSMichael Lange /* Bcast defaultValue */ 106177d236dfSMichael Lange if (!rank) (*labelNew)->defaultValue = label->defaultValue; 106277d236dfSMichael Lange ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 106326c55118SMichael Lange /* Distribute stratum values over the SF and get the point mapping on the receiver */ 106426c55118SMichael Lange ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr); 10655cbdf6fcSMichael Lange /* Determine received stratum values and initialise new label*/ 10665cbdf6fcSMichael Lange PetscHashICreate(stratumHash); 106726c55118SMichael Lange ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr); 10685cbdf6fcSMichael Lange for (p = 0; p < size; ++p) PetscHashIPut(stratumHash, leafStrata[p], ret, iter); 10695cbdf6fcSMichael Lange PetscHashISize(stratumHash, (*labelNew)->numStrata); 1070ad8374ffSToby Isaac ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr); 1071ad8374ffSToby Isaac for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE; 10725cbdf6fcSMichael Lange ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr); 10735cbdf6fcSMichael Lange /* Turn leafStrata into indices rather than stratum values */ 10745cbdf6fcSMichael Lange offset = 0; 10755cbdf6fcSMichael Lange ierr = PetscHashIGetKeys(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr); 10765cbdf6fcSMichael Lange for (p = 0; p < size; ++p) { 1077231b9e6fSMatthew G. Knepley for (s = 0; s < (*labelNew)->numStrata; ++s) { 1078231b9e6fSMatthew G. Knepley if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;} 10795cbdf6fcSMichael Lange } 10805cbdf6fcSMichael Lange } 1081c58f1c22SToby Isaac /* Rebuild the point strata on the receiver */ 1082c58f1c22SToby Isaac ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr); 1083c58f1c22SToby Isaac ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1084c58f1c22SToby Isaac for (p=pStart; p<pEnd; p++) { 1085c58f1c22SToby Isaac ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1086c58f1c22SToby Isaac ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1087c58f1c22SToby Isaac for (s=0; s<dof; s++) { 1088c58f1c22SToby Isaac (*labelNew)->stratumSizes[leafStrata[offset+s]]++; 1089c58f1c22SToby Isaac } 1090c58f1c22SToby Isaac } 1091c58f1c22SToby Isaac ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr); 1092c58f1c22SToby Isaac ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr); 1093ad8374ffSToby Isaac ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr); 1094c58f1c22SToby Isaac for (s = 0; s < (*labelNew)->numStrata; ++s) { 1095c58f1c22SToby Isaac PetscHashICreate((*labelNew)->ht[s]); 1096ad8374ffSToby Isaac ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr); 1097c58f1c22SToby Isaac } 1098c58f1c22SToby Isaac /* Insert points into new strata */ 1099c58f1c22SToby Isaac ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr); 1100c58f1c22SToby Isaac ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1101c58f1c22SToby Isaac for (p=pStart; p<pEnd; p++) { 1102c58f1c22SToby Isaac ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1103c58f1c22SToby Isaac ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1104c58f1c22SToby Isaac for (s=0; s<dof; s++) { 1105c58f1c22SToby Isaac stratum = leafStrata[offset+s]; 1106ad8374ffSToby Isaac points[stratum][strataIdx[stratum]++] = p; 1107c58f1c22SToby Isaac } 1108c58f1c22SToby Isaac } 1109ad8374ffSToby Isaac for (s = 0; s < (*labelNew)->numStrata; s++) { 1110ad8374ffSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr); 1111ad8374ffSToby Isaac ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr); 1112ad8374ffSToby Isaac } 1113ad8374ffSToby Isaac ierr = PetscFree(points);CHKERRQ(ierr); 11145cbdf6fcSMichael Lange PetscHashIDestroy(stratumHash); 1115c58f1c22SToby Isaac ierr = PetscFree(leafStrata);CHKERRQ(ierr); 1116c58f1c22SToby Isaac ierr = PetscFree(strataIdx);CHKERRQ(ierr); 1117c58f1c22SToby Isaac ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 1118c58f1c22SToby Isaac PetscFunctionReturn(0); 1119c58f1c22SToby Isaac } 1120c58f1c22SToby Isaac 1121c58f1c22SToby Isaac #undef __FUNCT__ 11227937d9ceSMichael Lange #define __FUNCT__ "DMLabelGather" 11237937d9ceSMichael Lange /*@ 11247937d9ceSMichael Lange DMLabelGather - Gather all label values from leafs into roots 11257937d9ceSMichael Lange 11267937d9ceSMichael Lange Input Parameters: 11277937d9ceSMichael Lange + label - the DMLabel 11287937d9ceSMichael Lange . point - the Star Forest point communication map 11297937d9ceSMichael Lange 11307937d9ceSMichael Lange Input Parameters: 11317937d9ceSMichael Lange + label - the new DMLabel with localised leaf values 11327937d9ceSMichael Lange 11337937d9ceSMichael Lange Level: developer 11347937d9ceSMichael Lange 11357937d9ceSMichael Lange Note: This is the inverse operation to DMLabelDistribute. 11367937d9ceSMichael Lange 11377937d9ceSMichael Lange .seealso: DMLabelDistribute() 11387937d9ceSMichael Lange @*/ 11397937d9ceSMichael Lange PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew) 11407937d9ceSMichael Lange { 11417937d9ceSMichael Lange MPI_Comm comm; 11427937d9ceSMichael Lange PetscSection rootSection; 11437937d9ceSMichael Lange PetscSF sfLabel; 11447937d9ceSMichael Lange PetscSFNode *rootPoints, *leafPoints; 11457937d9ceSMichael Lange PetscInt p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset; 11467937d9ceSMichael Lange const PetscInt *rootDegree, *ilocal; 11477937d9ceSMichael Lange PetscInt *rootStrata; 11487937d9ceSMichael Lange char *name; 11497937d9ceSMichael Lange PetscInt nameSize; 11507937d9ceSMichael Lange size_t len = 0; 11517937d9ceSMichael Lange PetscMPIInt rank, numProcs; 11527937d9ceSMichael Lange PetscErrorCode ierr; 11537937d9ceSMichael Lange 11547937d9ceSMichael Lange PetscFunctionBegin; 11557937d9ceSMichael Lange ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 11567937d9ceSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 11577937d9ceSMichael Lange ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 11587937d9ceSMichael Lange /* Bcast name */ 11597937d9ceSMichael Lange if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 11607937d9ceSMichael Lange nameSize = len; 11617937d9ceSMichael Lange ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 11627937d9ceSMichael Lange ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 11637937d9ceSMichael Lange if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 11647937d9ceSMichael Lange ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 11657937d9ceSMichael Lange ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 11667937d9ceSMichael Lange ierr = PetscFree(name);CHKERRQ(ierr); 11677937d9ceSMichael Lange /* Gather rank/index pairs of leaves into local roots to build 11687937d9ceSMichael Lange an inverse, multi-rooted SF. Note that this ignores local leaf 11697937d9ceSMichael Lange indexing due to the use of the multiSF in PetscSFGather. */ 11707937d9ceSMichael Lange ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr); 1171dc53bc9bSMatthew G. Knepley ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr); 1172dc53bc9bSMatthew G. Knepley for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1; 11737937d9ceSMichael Lange for (p = 0; p < nleaves; p++) { 1174dc53bc9bSMatthew G. Knepley leafPoints[ilocal[p]].index = ilocal[p]; 1175dc53bc9bSMatthew G. Knepley leafPoints[ilocal[p]].rank = rank; 11767937d9ceSMichael Lange } 11777937d9ceSMichael Lange ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr); 11787937d9ceSMichael Lange ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr); 11797937d9ceSMichael Lange for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p]; 11807937d9ceSMichael Lange ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr); 11817937d9ceSMichael Lange ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 11827937d9ceSMichael Lange ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 11837937d9ceSMichael Lange ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr); 11847937d9ceSMichael Lange ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 11857937d9ceSMichael Lange /* Migrate label over inverted SF to pull stratum values at leaves into roots. */ 11867937d9ceSMichael Lange ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr); 11877937d9ceSMichael Lange /* Rebuild the point strata on the receiver */ 11887937d9ceSMichael Lange for (p = 0, idx = 0; p < nroots; p++) { 11897937d9ceSMichael Lange for (d = 0; d < rootDegree[p]; d++) { 11907937d9ceSMichael Lange ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr); 11917937d9ceSMichael Lange ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr); 11927937d9ceSMichael Lange for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);} 11937937d9ceSMichael Lange } 11947937d9ceSMichael Lange idx += rootDegree[p]; 11957937d9ceSMichael Lange } 119677e0c0e7SMichael Lange ierr = PetscFree(leafPoints);CHKERRQ(ierr); 119777e0c0e7SMichael Lange ierr = PetscFree(rootStrata);CHKERRQ(ierr); 119877e0c0e7SMichael Lange ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 119977e0c0e7SMichael Lange ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr); 12007937d9ceSMichael Lange PetscFunctionReturn(0); 12017937d9ceSMichael Lange } 12027937d9ceSMichael Lange 12037937d9ceSMichael Lange #undef __FUNCT__ 1204c58f1c22SToby Isaac #define __FUNCT__ "DMLabelConvertToSection" 1205c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is) 1206c58f1c22SToby Isaac { 1207c58f1c22SToby Isaac IS vIS; 1208c58f1c22SToby Isaac const PetscInt *values; 1209c58f1c22SToby Isaac PetscInt *points; 1210c58f1c22SToby Isaac PetscInt nV, vS = 0, vE = 0, v, N; 1211c58f1c22SToby Isaac PetscErrorCode ierr; 1212c58f1c22SToby Isaac 1213c58f1c22SToby Isaac PetscFunctionBegin; 1214c58f1c22SToby Isaac ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr); 1215c58f1c22SToby Isaac ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr); 1216c58f1c22SToby Isaac ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr); 1217c58f1c22SToby Isaac if (nV) {vS = values[0]; vE = values[0]+1;} 1218c58f1c22SToby Isaac for (v = 1; v < nV; ++v) { 1219c58f1c22SToby Isaac vS = PetscMin(vS, values[v]); 1220c58f1c22SToby Isaac vE = PetscMax(vE, values[v]+1); 1221c58f1c22SToby Isaac } 1222c58f1c22SToby Isaac ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr); 1223c58f1c22SToby Isaac ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr); 1224c58f1c22SToby Isaac for (v = 0; v < nV; ++v) { 1225c58f1c22SToby Isaac PetscInt n; 1226c58f1c22SToby Isaac 1227c58f1c22SToby Isaac ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr); 1228c58f1c22SToby Isaac ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr); 1229c58f1c22SToby Isaac } 1230c58f1c22SToby Isaac ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 1231c58f1c22SToby Isaac ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr); 1232c58f1c22SToby Isaac ierr = PetscMalloc1(N, &points);CHKERRQ(ierr); 1233c58f1c22SToby Isaac for (v = 0; v < nV; ++v) { 1234c58f1c22SToby Isaac IS is; 1235c58f1c22SToby Isaac const PetscInt *spoints; 1236c58f1c22SToby Isaac PetscInt dof, off, p; 1237c58f1c22SToby Isaac 1238c58f1c22SToby Isaac ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr); 1239c58f1c22SToby Isaac ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr); 1240c58f1c22SToby Isaac ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr); 1241c58f1c22SToby Isaac ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr); 1242c58f1c22SToby Isaac for (p = 0; p < dof; ++p) points[off+p] = spoints[p]; 1243c58f1c22SToby Isaac ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr); 1244c58f1c22SToby Isaac ierr = ISDestroy(&is);CHKERRQ(ierr); 1245c58f1c22SToby Isaac } 1246c58f1c22SToby Isaac ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr); 1247c58f1c22SToby Isaac ierr = ISDestroy(&vIS);CHKERRQ(ierr); 1248c58f1c22SToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr); 1249c58f1c22SToby Isaac PetscFunctionReturn(0); 1250c58f1c22SToby Isaac } 1251c58f1c22SToby Isaac 1252c58f1c22SToby Isaac #undef __FUNCT__ 1253c58f1c22SToby Isaac #define __FUNCT__ "PetscSectionCreateGlobalSectionLabel" 1254c58f1c22SToby Isaac /*@C 1255c58f1c22SToby Isaac PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using 1256c58f1c22SToby Isaac the local section and an SF describing the section point overlap. 1257c58f1c22SToby Isaac 1258c58f1c22SToby Isaac Input Parameters: 1259c58f1c22SToby Isaac + s - The PetscSection for the local field layout 1260c58f1c22SToby Isaac . sf - The SF describing parallel layout of the section points 1261c58f1c22SToby Isaac . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs 1262c58f1c22SToby Isaac . label - The label specifying the points 1263c58f1c22SToby Isaac - labelValue - The label stratum specifying the points 1264c58f1c22SToby Isaac 1265c58f1c22SToby Isaac Output Parameter: 1266c58f1c22SToby Isaac . gsection - The PetscSection for the global field layout 1267c58f1c22SToby Isaac 1268c58f1c22SToby Isaac Note: This gives negative sizes and offsets to points not owned by this process 1269c58f1c22SToby Isaac 1270c58f1c22SToby Isaac Level: developer 1271c58f1c22SToby Isaac 1272c58f1c22SToby Isaac .seealso: PetscSectionCreate() 1273c58f1c22SToby Isaac @*/ 1274c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection) 1275c58f1c22SToby Isaac { 1276c58f1c22SToby Isaac PetscInt *neg = NULL, *tmpOff = NULL; 1277c58f1c22SToby Isaac PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots; 1278c58f1c22SToby Isaac PetscErrorCode ierr; 1279c58f1c22SToby Isaac 1280c58f1c22SToby Isaac PetscFunctionBegin; 1281c58f1c22SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr); 1282c58f1c22SToby Isaac ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 1283c58f1c22SToby Isaac ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr); 1284c58f1c22SToby Isaac ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 1285c58f1c22SToby Isaac if (nroots >= 0) { 1286c58f1c22SToby Isaac if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart); 1287c58f1c22SToby Isaac ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr); 1288c58f1c22SToby Isaac if (nroots > pEnd-pStart) { 1289c58f1c22SToby Isaac ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr); 1290c58f1c22SToby Isaac } else { 1291c58f1c22SToby Isaac tmpOff = &(*gsection)->atlasDof[-pStart]; 1292c58f1c22SToby Isaac } 1293c58f1c22SToby Isaac } 1294c58f1c22SToby Isaac /* Mark ghost points with negative dof */ 1295c58f1c22SToby Isaac for (p = pStart; p < pEnd; ++p) { 1296c58f1c22SToby Isaac PetscInt value; 1297c58f1c22SToby Isaac 1298c58f1c22SToby Isaac ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr); 1299c58f1c22SToby Isaac if (value != labelValue) continue; 1300c58f1c22SToby Isaac ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr); 1301c58f1c22SToby Isaac ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr); 1302c58f1c22SToby Isaac ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 1303c58f1c22SToby Isaac if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);} 1304c58f1c22SToby Isaac if (neg) neg[p] = -(dof+1); 1305c58f1c22SToby Isaac } 1306c58f1c22SToby Isaac ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr); 1307c58f1c22SToby Isaac if (nroots >= 0) { 1308c58f1c22SToby Isaac ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1309c58f1c22SToby Isaac ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1310c58f1c22SToby Isaac if (nroots > pEnd-pStart) { 1311c58f1c22SToby Isaac for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];} 1312c58f1c22SToby Isaac } 1313c58f1c22SToby Isaac } 1314c58f1c22SToby Isaac /* Calculate new sizes, get proccess offset, and calculate point offsets */ 1315c58f1c22SToby Isaac for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1316c58f1c22SToby Isaac cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0; 1317c58f1c22SToby Isaac (*gsection)->atlasOff[p] = off; 1318c58f1c22SToby Isaac off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0; 1319c58f1c22SToby Isaac } 1320c58f1c22SToby Isaac ierr = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr); 1321c58f1c22SToby Isaac globalOff -= off; 1322c58f1c22SToby Isaac for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1323c58f1c22SToby Isaac (*gsection)->atlasOff[p] += globalOff; 1324c58f1c22SToby Isaac if (neg) neg[p] = -((*gsection)->atlasOff[p]+1); 1325c58f1c22SToby Isaac } 1326c58f1c22SToby Isaac /* Put in negative offsets for ghost points */ 1327c58f1c22SToby Isaac if (nroots >= 0) { 1328c58f1c22SToby Isaac ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1329c58f1c22SToby Isaac ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1330c58f1c22SToby Isaac if (nroots > pEnd-pStart) { 1331c58f1c22SToby Isaac for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];} 1332c58f1c22SToby Isaac } 1333c58f1c22SToby Isaac } 1334c58f1c22SToby Isaac if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);} 1335c58f1c22SToby Isaac ierr = PetscFree(neg);CHKERRQ(ierr); 1336c58f1c22SToby Isaac PetscFunctionReturn(0); 1337c58f1c22SToby Isaac } 1338c58f1c22SToby Isaac 1339