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; 31c58f1c22SToby Isaac (*label)->stratumValues = NULL; 32c58f1c22SToby Isaac (*label)->arrayValid = NULL; 33c58f1c22SToby Isaac (*label)->stratumSizes = NULL; 34c58f1c22SToby Isaac (*label)->points = NULL; 35c58f1c22SToby Isaac (*label)->ht = NULL; 36c58f1c22SToby Isaac (*label)->pStart = -1; 37c58f1c22SToby Isaac (*label)->pEnd = -1; 38c58f1c22SToby Isaac (*label)->bt = NULL; 39c58f1c22SToby Isaac PetscFunctionReturn(0); 40c58f1c22SToby Isaac } 41c58f1c22SToby Isaac 42c58f1c22SToby Isaac #undef __FUNCT__ 43c58f1c22SToby Isaac #define __FUNCT__ "DMLabelMakeValid_Private" 44c58f1c22SToby Isaac /* 45c58f1c22SToby Isaac DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format 46c58f1c22SToby Isaac 47c58f1c22SToby Isaac Input parameter: 48c58f1c22SToby Isaac + label - The DMLabel 49c58f1c22SToby Isaac - v - The stratum value 50c58f1c22SToby Isaac 51c58f1c22SToby Isaac Output parameter: 52c58f1c22SToby Isaac . label - The DMLabel with stratum in sorted list format 53c58f1c22SToby Isaac 54c58f1c22SToby Isaac Level: developer 55c58f1c22SToby Isaac 56c58f1c22SToby Isaac .seealso: DMLabelCreate() 57c58f1c22SToby Isaac */ 58c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v) 59c58f1c22SToby Isaac { 60c58f1c22SToby Isaac PetscInt off; 61c58f1c22SToby Isaac PetscErrorCode ierr; 62c58f1c22SToby Isaac 63c58f1c22SToby Isaac if (label->arrayValid[v]) return 0; 64c58f1c22SToby Isaac if (v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeValid_Private\n", v); 65c58f1c22SToby Isaac PetscFunctionBegin; 66c58f1c22SToby Isaac PetscHashISize(label->ht[v], label->stratumSizes[v]); 67c58f1c22SToby Isaac 68c58f1c22SToby Isaac ierr = PetscMalloc1(label->stratumSizes[v], &label->points[v]);CHKERRQ(ierr); 69c58f1c22SToby Isaac off = 0; 70c58f1c22SToby Isaac ierr = PetscHashIGetKeys(label->ht[v], &off, &(label->points[v][0]));CHKERRQ(ierr); 71c58f1c22SToby 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]); 72c58f1c22SToby Isaac PetscHashIClear(label->ht[v]); 73c58f1c22SToby Isaac ierr = PetscSortInt(label->stratumSizes[v], label->points[v]);CHKERRQ(ierr); 74c58f1c22SToby Isaac if (label->bt) { 75c58f1c22SToby Isaac PetscInt p; 76c58f1c22SToby Isaac 77c58f1c22SToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) { 78c58f1c22SToby Isaac const PetscInt point = label->points[v][p]; 79c58f1c22SToby Isaac 80c58f1c22SToby 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); 81c58f1c22SToby Isaac ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 82c58f1c22SToby Isaac } 83c58f1c22SToby Isaac } 84c58f1c22SToby Isaac label->arrayValid[v] = PETSC_TRUE; 85c58f1c22SToby Isaac ++label->state; 86c58f1c22SToby Isaac PetscFunctionReturn(0); 87c58f1c22SToby Isaac } 88c58f1c22SToby Isaac 89c58f1c22SToby Isaac #undef __FUNCT__ 90c58f1c22SToby Isaac #define __FUNCT__ "DMLabelMakeAllValid_Private" 91c58f1c22SToby Isaac /* 92c58f1c22SToby Isaac DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format 93c58f1c22SToby Isaac 94c58f1c22SToby Isaac Input parameter: 95c58f1c22SToby Isaac . label - The DMLabel 96c58f1c22SToby Isaac 97c58f1c22SToby Isaac Output parameter: 98c58f1c22SToby Isaac . label - The DMLabel with all strata in sorted list format 99c58f1c22SToby Isaac 100c58f1c22SToby Isaac Level: developer 101c58f1c22SToby Isaac 102c58f1c22SToby Isaac .seealso: DMLabelCreate() 103c58f1c22SToby Isaac */ 104c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label) 105c58f1c22SToby Isaac { 106c58f1c22SToby Isaac PetscInt v; 107c58f1c22SToby Isaac PetscErrorCode ierr; 108c58f1c22SToby Isaac 109c58f1c22SToby Isaac PetscFunctionBegin; 110c58f1c22SToby Isaac for (v = 0; v < label->numStrata; v++){ 111c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 112c58f1c22SToby Isaac } 113c58f1c22SToby Isaac PetscFunctionReturn(0); 114c58f1c22SToby Isaac } 115c58f1c22SToby Isaac 116c58f1c22SToby Isaac #undef __FUNCT__ 117c58f1c22SToby Isaac #define __FUNCT__ "DMLabelMakeInvalid_Private" 118c58f1c22SToby Isaac /* 119c58f1c22SToby Isaac DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format 120c58f1c22SToby Isaac 121c58f1c22SToby Isaac Input parameter: 122c58f1c22SToby Isaac + label - The DMLabel 123c58f1c22SToby Isaac - v - The stratum value 124c58f1c22SToby Isaac 125c58f1c22SToby Isaac Output parameter: 126c58f1c22SToby Isaac . label - The DMLabel with stratum in hash format 127c58f1c22SToby Isaac 128c58f1c22SToby Isaac Level: developer 129c58f1c22SToby Isaac 130c58f1c22SToby Isaac .seealso: DMLabelCreate() 131c58f1c22SToby Isaac */ 132c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v) 133c58f1c22SToby Isaac { 134c58f1c22SToby Isaac PETSC_UNUSED PetscHashIIter ret, iter; 135c58f1c22SToby Isaac PetscInt p; 136c58f1c22SToby Isaac PetscErrorCode ierr; 137c58f1c22SToby Isaac 138c58f1c22SToby Isaac PetscFunctionBegin; 139c58f1c22SToby Isaac if (!label->arrayValid[v]) PetscFunctionReturn(0); 140c58f1c22SToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) PetscHashIPut(label->ht[v], label->points[v][p], ret, iter); 141c58f1c22SToby Isaac ierr = PetscFree(label->points[v]);CHKERRQ(ierr); 142c58f1c22SToby Isaac label->arrayValid[v] = PETSC_FALSE; 143c58f1c22SToby Isaac PetscFunctionReturn(0); 144c58f1c22SToby Isaac } 145c58f1c22SToby Isaac 146c58f1c22SToby Isaac #undef __FUNCT__ 147c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetState" 148c58f1c22SToby Isaac PetscErrorCode DMLabelGetState(DMLabel label, PetscObjectState *state) 149c58f1c22SToby Isaac { 150c58f1c22SToby Isaac PetscFunctionBegin; 151c58f1c22SToby Isaac PetscValidPointer(state, 2); 152c58f1c22SToby Isaac *state = label->state; 153c58f1c22SToby Isaac PetscFunctionReturn(0); 154c58f1c22SToby Isaac } 155c58f1c22SToby Isaac 156c58f1c22SToby Isaac #undef __FUNCT__ 157c58f1c22SToby Isaac #define __FUNCT__ "DMLabelAddStratum" 158c58f1c22SToby Isaac PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value) 159c58f1c22SToby Isaac { 160c58f1c22SToby Isaac PetscInt v, *tmpV, *tmpS, **tmpP; 161c58f1c22SToby Isaac PetscHashI *tmpH; 162c58f1c22SToby Isaac PetscBool *tmpB; 163c58f1c22SToby Isaac PetscErrorCode ierr; 164c58f1c22SToby Isaac 165c58f1c22SToby Isaac PetscFunctionBegin; 166c58f1c22SToby Isaac 167c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpV);CHKERRQ(ierr); 168c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpS);CHKERRQ(ierr); 169c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpH);CHKERRQ(ierr); 170c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpP);CHKERRQ(ierr); 171c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpB);CHKERRQ(ierr); 172c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 173c58f1c22SToby Isaac tmpV[v] = label->stratumValues[v]; 174c58f1c22SToby Isaac tmpS[v] = label->stratumSizes[v]; 175c58f1c22SToby Isaac tmpH[v] = label->ht[v]; 176c58f1c22SToby Isaac tmpP[v] = label->points[v]; 177c58f1c22SToby Isaac tmpB[v] = label->arrayValid[v]; 178c58f1c22SToby Isaac } 179c58f1c22SToby Isaac tmpV[v] = value; 180c58f1c22SToby Isaac tmpS[v] = 0; 181c58f1c22SToby Isaac PetscHashICreate(tmpH[v]); 182c58f1c22SToby Isaac tmpP[v] = NULL; 183c58f1c22SToby Isaac tmpB[v] = PETSC_TRUE; 184c58f1c22SToby Isaac ++label->numStrata; 185c58f1c22SToby Isaac ierr = PetscFree(label->stratumValues);CHKERRQ(ierr); 186c58f1c22SToby Isaac ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr); 187c58f1c22SToby Isaac ierr = PetscFree(label->ht);CHKERRQ(ierr); 188c58f1c22SToby Isaac ierr = PetscFree(label->points);CHKERRQ(ierr); 189c58f1c22SToby Isaac ierr = PetscFree(label->arrayValid);CHKERRQ(ierr); 190c58f1c22SToby Isaac label->stratumValues = tmpV; 191c58f1c22SToby Isaac label->stratumSizes = tmpS; 192c58f1c22SToby Isaac label->ht = tmpH; 193c58f1c22SToby Isaac label->points = tmpP; 194c58f1c22SToby Isaac label->arrayValid = tmpB; 195c58f1c22SToby Isaac 196c58f1c22SToby Isaac PetscFunctionReturn(0); 197c58f1c22SToby Isaac } 198c58f1c22SToby Isaac 199c58f1c22SToby Isaac #undef __FUNCT__ 200c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetName" 201c58f1c22SToby Isaac /*@C 202c58f1c22SToby Isaac DMLabelGetName - Return the name of a DMLabel object 203c58f1c22SToby Isaac 204c58f1c22SToby Isaac Input parameter: 205c58f1c22SToby Isaac . label - The DMLabel 206c58f1c22SToby Isaac 207c58f1c22SToby Isaac Output parameter: 208c58f1c22SToby Isaac . name - The label name 209c58f1c22SToby Isaac 210c58f1c22SToby Isaac Level: beginner 211c58f1c22SToby Isaac 212c58f1c22SToby Isaac .seealso: DMLabelCreate() 213c58f1c22SToby Isaac @*/ 214c58f1c22SToby Isaac PetscErrorCode DMLabelGetName(DMLabel label, const char **name) 215c58f1c22SToby Isaac { 216c58f1c22SToby Isaac PetscFunctionBegin; 217c58f1c22SToby Isaac PetscValidPointer(name, 2); 218c58f1c22SToby Isaac *name = label->name; 219c58f1c22SToby Isaac PetscFunctionReturn(0); 220c58f1c22SToby Isaac } 221c58f1c22SToby Isaac 222c58f1c22SToby Isaac #undef __FUNCT__ 223c58f1c22SToby Isaac #define __FUNCT__ "DMLabelView_Ascii" 224c58f1c22SToby Isaac static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer) 225c58f1c22SToby Isaac { 226c58f1c22SToby Isaac PetscInt v; 227c58f1c22SToby Isaac PetscMPIInt rank; 228c58f1c22SToby Isaac PetscErrorCode ierr; 229c58f1c22SToby Isaac 230c58f1c22SToby Isaac PetscFunctionBegin; 231c58f1c22SToby Isaac ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr); 232c58f1c22SToby Isaac ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 233c58f1c22SToby Isaac if (label) { 234c58f1c22SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", label->name);CHKERRQ(ierr); 235c58f1c22SToby Isaac if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, " Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);} 236c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 237c58f1c22SToby Isaac const PetscInt value = label->stratumValues[v]; 238c58f1c22SToby Isaac PetscInt p; 239c58f1c22SToby Isaac 240c58f1c22SToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) { 241c58f1c22SToby Isaac ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, label->points[v][p], value);CHKERRQ(ierr); 242c58f1c22SToby Isaac } 243c58f1c22SToby Isaac } 244c58f1c22SToby Isaac } 245c58f1c22SToby Isaac ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 246c58f1c22SToby Isaac ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 247c58f1c22SToby Isaac PetscFunctionReturn(0); 248c58f1c22SToby Isaac } 249c58f1c22SToby Isaac 250c58f1c22SToby Isaac #undef __FUNCT__ 251c58f1c22SToby Isaac #define __FUNCT__ "DMLabelView" 252c58f1c22SToby Isaac /*@C 253c58f1c22SToby Isaac DMLabelView - View the label 254c58f1c22SToby Isaac 255c58f1c22SToby Isaac Input Parameters: 256c58f1c22SToby Isaac + label - The DMLabel 257c58f1c22SToby Isaac - viewer - The PetscViewer 258c58f1c22SToby Isaac 259c58f1c22SToby Isaac Level: intermediate 260c58f1c22SToby Isaac 261c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy() 262c58f1c22SToby Isaac @*/ 263c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer) 264c58f1c22SToby Isaac { 265c58f1c22SToby Isaac PetscBool iascii; 266c58f1c22SToby Isaac PetscErrorCode ierr; 267c58f1c22SToby Isaac 268c58f1c22SToby Isaac PetscFunctionBegin; 269c58f1c22SToby Isaac PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 270c58f1c22SToby Isaac if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 271c58f1c22SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 272c58f1c22SToby Isaac if (iascii) { 273c58f1c22SToby Isaac ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr); 274c58f1c22SToby Isaac } 275c58f1c22SToby Isaac PetscFunctionReturn(0); 276c58f1c22SToby Isaac } 277c58f1c22SToby Isaac 278c58f1c22SToby Isaac #undef __FUNCT__ 279c58f1c22SToby Isaac #define __FUNCT__ "DMLabelDestroy" 280c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label) 281c58f1c22SToby Isaac { 282c58f1c22SToby Isaac PetscInt v; 283c58f1c22SToby Isaac PetscErrorCode ierr; 284c58f1c22SToby Isaac 285c58f1c22SToby Isaac PetscFunctionBegin; 286c58f1c22SToby Isaac if (!(*label)) PetscFunctionReturn(0); 287c58f1c22SToby Isaac if (--(*label)->refct > 0) PetscFunctionReturn(0); 288c58f1c22SToby Isaac ierr = PetscFree((*label)->name);CHKERRQ(ierr); 289c58f1c22SToby Isaac ierr = PetscFree((*label)->stratumValues);CHKERRQ(ierr); 290c58f1c22SToby Isaac ierr = PetscFree((*label)->stratumSizes);CHKERRQ(ierr); 291c58f1c22SToby Isaac for (v = 0; v < (*label)->numStrata; ++v) {ierr = PetscFree((*label)->points[v]);CHKERRQ(ierr);} 292c58f1c22SToby Isaac ierr = PetscFree((*label)->points);CHKERRQ(ierr); 293c58f1c22SToby Isaac ierr = PetscFree((*label)->arrayValid);CHKERRQ(ierr); 294c58f1c22SToby Isaac if ((*label)->ht) { 295c58f1c22SToby Isaac for (v = 0; v < (*label)->numStrata; ++v) {PetscHashIDestroy((*label)->ht[v]);} 296c58f1c22SToby Isaac ierr = PetscFree((*label)->ht);CHKERRQ(ierr); 297c58f1c22SToby Isaac } 298c58f1c22SToby Isaac ierr = PetscBTDestroy(&(*label)->bt);CHKERRQ(ierr); 299c58f1c22SToby Isaac ierr = PetscFree(*label);CHKERRQ(ierr); 300c58f1c22SToby Isaac PetscFunctionReturn(0); 301c58f1c22SToby Isaac } 302c58f1c22SToby Isaac 303c58f1c22SToby Isaac #undef __FUNCT__ 304c58f1c22SToby Isaac #define __FUNCT__ "DMLabelDuplicate" 305c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew) 306c58f1c22SToby Isaac { 307c58f1c22SToby Isaac PetscInt v, q; 308c58f1c22SToby Isaac PetscErrorCode ierr; 309c58f1c22SToby Isaac 310c58f1c22SToby Isaac PetscFunctionBegin; 311c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 312c58f1c22SToby Isaac ierr = PetscNew(labelnew);CHKERRQ(ierr); 313c58f1c22SToby Isaac ierr = PetscStrallocpy(label->name, &(*labelnew)->name);CHKERRQ(ierr); 314c58f1c22SToby Isaac 315c58f1c22SToby Isaac (*labelnew)->refct = 1; 316c58f1c22SToby Isaac (*labelnew)->numStrata = label->numStrata; 317c58f1c22SToby Isaac if (label->numStrata) { 318c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr); 319c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr); 320c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr); 321c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr); 322c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->arrayValid);CHKERRQ(ierr); 323c58f1c22SToby Isaac /* Could eliminate unused space here */ 324c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 325c58f1c22SToby Isaac ierr = PetscMalloc1(label->stratumSizes[v], &(*labelnew)->points[v]);CHKERRQ(ierr); 326c58f1c22SToby Isaac PetscHashICreate((*labelnew)->ht[v]); 327c58f1c22SToby Isaac (*labelnew)->arrayValid[v] = PETSC_TRUE; 328c58f1c22SToby Isaac (*labelnew)->stratumValues[v] = label->stratumValues[v]; 329c58f1c22SToby Isaac (*labelnew)->stratumSizes[v] = label->stratumSizes[v]; 330c58f1c22SToby Isaac for (q = 0; q < label->stratumSizes[v]; ++q) { 331c58f1c22SToby Isaac (*labelnew)->points[v][q] = label->points[v][q]; 332c58f1c22SToby Isaac } 333c58f1c22SToby Isaac } 334c58f1c22SToby Isaac } 335c58f1c22SToby Isaac (*labelnew)->pStart = -1; 336c58f1c22SToby Isaac (*labelnew)->pEnd = -1; 337c58f1c22SToby Isaac (*labelnew)->bt = NULL; 338c58f1c22SToby Isaac PetscFunctionReturn(0); 339c58f1c22SToby Isaac } 340c58f1c22SToby Isaac 341c58f1c22SToby Isaac #undef __FUNCT__ 342c58f1c22SToby Isaac #define __FUNCT__ "DMLabelCreateIndex" 343c58f1c22SToby Isaac /* This can be hooked into SetValue(), ClearValue(), etc. for updating */ 344c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd) 345c58f1c22SToby Isaac { 346c58f1c22SToby Isaac PetscInt v; 347c58f1c22SToby Isaac PetscErrorCode ierr; 348c58f1c22SToby Isaac 349c58f1c22SToby Isaac PetscFunctionBegin; 350c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 351c58f1c22SToby Isaac if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 352c58f1c22SToby Isaac label->pStart = pStart; 353c58f1c22SToby Isaac label->pEnd = pEnd; 354c58f1c22SToby Isaac ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr); 355c58f1c22SToby Isaac ierr = PetscBTMemzero(pEnd - pStart, label->bt);CHKERRQ(ierr); 356c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 357c58f1c22SToby Isaac PetscInt i; 358c58f1c22SToby Isaac 359c58f1c22SToby Isaac for (i = 0; i < label->stratumSizes[v]; ++i) { 360c58f1c22SToby Isaac const PetscInt point = label->points[v][i]; 361c58f1c22SToby Isaac 362c58f1c22SToby 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); 363c58f1c22SToby Isaac ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr); 364c58f1c22SToby Isaac } 365c58f1c22SToby Isaac } 366c58f1c22SToby Isaac PetscFunctionReturn(0); 367c58f1c22SToby Isaac } 368c58f1c22SToby Isaac 369c58f1c22SToby Isaac #undef __FUNCT__ 370c58f1c22SToby Isaac #define __FUNCT__ "DMLabelDestroyIndex" 371c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label) 372c58f1c22SToby Isaac { 373c58f1c22SToby Isaac PetscErrorCode ierr; 374c58f1c22SToby Isaac 375c58f1c22SToby Isaac PetscFunctionBegin; 376c58f1c22SToby Isaac label->pStart = -1; 377c58f1c22SToby Isaac label->pEnd = -1; 378c58f1c22SToby Isaac if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 379c58f1c22SToby Isaac PetscFunctionReturn(0); 380c58f1c22SToby Isaac } 381c58f1c22SToby Isaac 382c58f1c22SToby Isaac #undef __FUNCT__ 383c58f1c22SToby Isaac #define __FUNCT__ "DMLabelHasValue" 384c58f1c22SToby Isaac /*@ 385c58f1c22SToby Isaac DMLabelHasValue - Determine whether a label assigns the value to any point 386c58f1c22SToby Isaac 387c58f1c22SToby Isaac Input Parameters: 388c58f1c22SToby Isaac + label - the DMLabel 389c58f1c22SToby Isaac - value - the value 390c58f1c22SToby Isaac 391c58f1c22SToby Isaac Output Parameter: 392c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point 393c58f1c22SToby Isaac 394c58f1c22SToby Isaac Level: developer 395c58f1c22SToby Isaac 396c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue() 397c58f1c22SToby Isaac @*/ 398c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains) 399c58f1c22SToby Isaac { 400c58f1c22SToby Isaac PetscInt v; 401c58f1c22SToby Isaac 402c58f1c22SToby Isaac PetscFunctionBegin; 403c58f1c22SToby Isaac PetscValidPointer(contains, 3); 404c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 405c58f1c22SToby Isaac if (value == label->stratumValues[v]) break; 406c58f1c22SToby Isaac } 407c58f1c22SToby Isaac *contains = (v < label->numStrata ? PETSC_TRUE : PETSC_FALSE); 408c58f1c22SToby Isaac PetscFunctionReturn(0); 409c58f1c22SToby Isaac } 410c58f1c22SToby Isaac 411c58f1c22SToby Isaac #undef __FUNCT__ 412c58f1c22SToby Isaac #define __FUNCT__ "DMLabelHasPoint" 413c58f1c22SToby Isaac /*@ 414c58f1c22SToby Isaac DMLabelHasPoint - Determine whether a label assigns a value to a point 415c58f1c22SToby Isaac 416c58f1c22SToby Isaac Input Parameters: 417c58f1c22SToby Isaac + label - the DMLabel 418c58f1c22SToby Isaac - point - the point 419c58f1c22SToby Isaac 420c58f1c22SToby Isaac Output Parameter: 421c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value 422c58f1c22SToby Isaac 423c58f1c22SToby Isaac Note: The user must call DMLabelCreateIndex() before this function. 424c58f1c22SToby Isaac 425c58f1c22SToby Isaac Level: developer 426c58f1c22SToby Isaac 427c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue() 428c58f1c22SToby Isaac @*/ 429c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains) 430c58f1c22SToby Isaac { 431c58f1c22SToby Isaac PetscErrorCode ierr; 432c58f1c22SToby Isaac 433c58f1c22SToby Isaac PetscFunctionBeginHot; 434c58f1c22SToby Isaac PetscValidPointer(contains, 3); 435c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 436c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG) 437c58f1c22SToby Isaac if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()"); 438c58f1c22SToby 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); 439c58f1c22SToby Isaac #endif 440c58f1c22SToby Isaac *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE; 441c58f1c22SToby Isaac PetscFunctionReturn(0); 442c58f1c22SToby Isaac } 443c58f1c22SToby Isaac 444c58f1c22SToby Isaac #undef __FUNCT__ 445c58f1c22SToby Isaac #define __FUNCT__ "DMLabelStratumHasPoint" 446c58f1c22SToby Isaac /*@ 447c58f1c22SToby Isaac DMLabelStratumHasPoint - Return true if the stratum contains a point 448c58f1c22SToby Isaac 449c58f1c22SToby Isaac Input Parameters: 450c58f1c22SToby Isaac + label - the DMLabel 451c58f1c22SToby Isaac . value - the stratum value 452c58f1c22SToby Isaac - point - the point 453c58f1c22SToby Isaac 454c58f1c22SToby Isaac Output Parameter: 455c58f1c22SToby Isaac . contains - true if the stratum contains the point 456c58f1c22SToby Isaac 457c58f1c22SToby Isaac Level: intermediate 458c58f1c22SToby Isaac 459c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue() 460c58f1c22SToby Isaac @*/ 461c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains) 462c58f1c22SToby Isaac { 463c58f1c22SToby Isaac PetscInt v; 464c58f1c22SToby Isaac PetscErrorCode ierr; 465c58f1c22SToby Isaac 466c58f1c22SToby Isaac PetscFunctionBegin; 467c58f1c22SToby Isaac PetscValidPointer(contains, 4); 468c58f1c22SToby Isaac *contains = PETSC_FALSE; 469c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 470c58f1c22SToby Isaac if (label->stratumValues[v] == value) { 471c58f1c22SToby Isaac if (label->arrayValid[v]) { 472c58f1c22SToby Isaac PetscInt i; 473c58f1c22SToby Isaac 474c58f1c22SToby Isaac ierr = PetscFindInt(point, label->stratumSizes[v], &label->points[v][0], &i);CHKERRQ(ierr); 475c58f1c22SToby Isaac if (i >= 0) { 476c58f1c22SToby Isaac *contains = PETSC_TRUE; 477c58f1c22SToby Isaac break; 478c58f1c22SToby Isaac } 479c58f1c22SToby Isaac } else { 480c58f1c22SToby Isaac PetscBool has; 481c58f1c22SToby Isaac 482c58f1c22SToby Isaac PetscHashIHasKey(label->ht[v], point, has); 483c58f1c22SToby Isaac if (has) { 484c58f1c22SToby Isaac *contains = PETSC_TRUE; 485c58f1c22SToby Isaac break; 486c58f1c22SToby Isaac } 487c58f1c22SToby Isaac } 488c58f1c22SToby Isaac } 489c58f1c22SToby Isaac } 490c58f1c22SToby Isaac PetscFunctionReturn(0); 491c58f1c22SToby Isaac } 492c58f1c22SToby Isaac 493c58f1c22SToby Isaac #undef __FUNCT__ 494c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetValue" 495c58f1c22SToby Isaac /*@ 496c58f1c22SToby Isaac DMLabelGetValue - Return the value a label assigns to a point, or -1 497c58f1c22SToby Isaac 498c58f1c22SToby Isaac Input Parameters: 499c58f1c22SToby Isaac + label - the DMLabel 500c58f1c22SToby Isaac - point - the point 501c58f1c22SToby Isaac 502c58f1c22SToby Isaac Output Parameter: 503c58f1c22SToby Isaac . value - The point value, or -1 504c58f1c22SToby Isaac 505c58f1c22SToby Isaac Level: intermediate 506c58f1c22SToby Isaac 507c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue() 508c58f1c22SToby Isaac @*/ 509c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value) 510c58f1c22SToby Isaac { 511c58f1c22SToby Isaac PetscInt v; 512c58f1c22SToby Isaac PetscErrorCode ierr; 513c58f1c22SToby Isaac 514c58f1c22SToby Isaac PetscFunctionBegin; 515c58f1c22SToby Isaac PetscValidPointer(value, 3); 516c58f1c22SToby Isaac *value = -1; 517c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 518c58f1c22SToby Isaac if (label->arrayValid[v]) { 519c58f1c22SToby Isaac PetscInt i; 520c58f1c22SToby Isaac 521c58f1c22SToby Isaac ierr = PetscFindInt(point, label->stratumSizes[v], &label->points[v][0], &i);CHKERRQ(ierr); 522c58f1c22SToby Isaac if (i >= 0) { 523c58f1c22SToby Isaac *value = label->stratumValues[v]; 524c58f1c22SToby Isaac break; 525c58f1c22SToby Isaac } 526c58f1c22SToby Isaac } else { 527c58f1c22SToby Isaac PetscBool has; 528c58f1c22SToby Isaac 529c58f1c22SToby Isaac PetscHashIHasKey(label->ht[v], point, has); 530c58f1c22SToby Isaac if (has) { 531c58f1c22SToby Isaac *value = label->stratumValues[v]; 532c58f1c22SToby Isaac break; 533c58f1c22SToby Isaac } 534c58f1c22SToby Isaac } 535c58f1c22SToby Isaac } 536c58f1c22SToby Isaac PetscFunctionReturn(0); 537c58f1c22SToby Isaac } 538c58f1c22SToby Isaac 539c58f1c22SToby Isaac #undef __FUNCT__ 540c58f1c22SToby Isaac #define __FUNCT__ "DMLabelSetValue" 541c58f1c22SToby Isaac /*@ 542c58f1c22SToby Isaac DMLabelSetValue - Set the value a label assigns to a point 543c58f1c22SToby Isaac 544c58f1c22SToby Isaac Input Parameters: 545c58f1c22SToby Isaac + label - the DMLabel 546c58f1c22SToby Isaac . point - the point 547c58f1c22SToby Isaac - value - The point value 548c58f1c22SToby Isaac 549c58f1c22SToby Isaac Level: intermediate 550c58f1c22SToby Isaac 551c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue() 552c58f1c22SToby Isaac @*/ 553c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value) 554c58f1c22SToby Isaac { 555c58f1c22SToby Isaac PETSC_UNUSED PetscHashIIter iter, ret; 556c58f1c22SToby Isaac PetscInt v; 557c58f1c22SToby Isaac PetscErrorCode ierr; 558c58f1c22SToby Isaac 559c58f1c22SToby Isaac PetscFunctionBegin; 560c58f1c22SToby Isaac /* Find, or add, label value */ 561c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 562c58f1c22SToby Isaac if (label->stratumValues[v] == value) break; 563c58f1c22SToby Isaac } 564c58f1c22SToby Isaac /* Create new table */ 565c58f1c22SToby Isaac if (v >= label->numStrata) {ierr = DMLabelAddStratum(label, value);CHKERRQ(ierr);} 566c58f1c22SToby Isaac ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 567c58f1c22SToby Isaac /* Set key */ 568c58f1c22SToby Isaac PetscHashIPut(label->ht[v], point, ret, iter); 569c58f1c22SToby Isaac PetscFunctionReturn(0); 570c58f1c22SToby Isaac } 571c58f1c22SToby Isaac 572c58f1c22SToby Isaac #undef __FUNCT__ 573c58f1c22SToby Isaac #define __FUNCT__ "DMLabelClearValue" 574c58f1c22SToby Isaac /*@ 575c58f1c22SToby Isaac DMLabelClearValue - Clear the value a label assigns to a point 576c58f1c22SToby Isaac 577c58f1c22SToby Isaac Input Parameters: 578c58f1c22SToby Isaac + label - the DMLabel 579c58f1c22SToby Isaac . point - the point 580c58f1c22SToby Isaac - value - The point value 581c58f1c22SToby Isaac 582c58f1c22SToby Isaac Level: intermediate 583c58f1c22SToby Isaac 584c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue() 585c58f1c22SToby Isaac @*/ 586c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value) 587c58f1c22SToby Isaac { 588c58f1c22SToby Isaac PetscInt v, p; 589c58f1c22SToby Isaac PetscErrorCode ierr; 590c58f1c22SToby Isaac 591c58f1c22SToby Isaac PetscFunctionBegin; 592c58f1c22SToby Isaac /* Find label value */ 593c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 594c58f1c22SToby Isaac if (label->stratumValues[v] == value) break; 595c58f1c22SToby Isaac } 596c58f1c22SToby Isaac if (v >= label->numStrata) PetscFunctionReturn(0); 597c58f1c22SToby Isaac if (label->arrayValid[v]) { 598c58f1c22SToby Isaac /* Check whether point exists */ 599c58f1c22SToby Isaac ierr = PetscFindInt(point, label->stratumSizes[v], &label->points[v][0], &p);CHKERRQ(ierr); 600c58f1c22SToby Isaac if (p >= 0) { 601c58f1c22SToby Isaac ierr = PetscMemmove(&label->points[v][p], &label->points[v][p+1], (label->stratumSizes[v]-p-1) * sizeof(PetscInt));CHKERRQ(ierr); 602c58f1c22SToby Isaac --label->stratumSizes[v]; 603c58f1c22SToby Isaac if (label->bt) { 604c58f1c22SToby 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); 605c58f1c22SToby Isaac ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 606c58f1c22SToby Isaac } 607c58f1c22SToby Isaac } 608c58f1c22SToby Isaac } else { 609c58f1c22SToby Isaac ierr = PetscHashIDelKey(label->ht[v], point);CHKERRQ(ierr); 610c58f1c22SToby Isaac } 611c58f1c22SToby Isaac PetscFunctionReturn(0); 612c58f1c22SToby Isaac } 613c58f1c22SToby Isaac 614c58f1c22SToby Isaac #undef __FUNCT__ 615c58f1c22SToby Isaac #define __FUNCT__ "DMLabelInsertIS" 616c58f1c22SToby Isaac /*@ 617c58f1c22SToby Isaac DMLabelInsertIS - Set all points in the IS to a value 618c58f1c22SToby Isaac 619c58f1c22SToby Isaac Input Parameters: 620c58f1c22SToby Isaac + label - the DMLabel 621c58f1c22SToby Isaac . is - the point IS 622c58f1c22SToby Isaac - value - The point value 623c58f1c22SToby Isaac 624c58f1c22SToby Isaac Level: intermediate 625c58f1c22SToby Isaac 626c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 627c58f1c22SToby Isaac @*/ 628c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value) 629c58f1c22SToby Isaac { 630c58f1c22SToby Isaac const PetscInt *points; 631c58f1c22SToby Isaac PetscInt n, p; 632c58f1c22SToby Isaac PetscErrorCode ierr; 633c58f1c22SToby Isaac 634c58f1c22SToby Isaac PetscFunctionBegin; 635c58f1c22SToby Isaac PetscValidHeaderSpecific(is, IS_CLASSID, 2); 636c58f1c22SToby Isaac ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 637c58f1c22SToby Isaac ierr = ISGetIndices(is, &points);CHKERRQ(ierr); 638c58f1c22SToby Isaac for (p = 0; p < n; ++p) {ierr = DMLabelSetValue(label, points[p], value);CHKERRQ(ierr);} 639c58f1c22SToby Isaac ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr); 640c58f1c22SToby Isaac PetscFunctionReturn(0); 641c58f1c22SToby Isaac } 642c58f1c22SToby Isaac 643c58f1c22SToby Isaac #undef __FUNCT__ 644c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetNumValues" 645c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues) 646c58f1c22SToby Isaac { 647c58f1c22SToby Isaac PetscFunctionBegin; 648c58f1c22SToby Isaac PetscValidPointer(numValues, 2); 649c58f1c22SToby Isaac *numValues = label->numStrata; 650c58f1c22SToby Isaac PetscFunctionReturn(0); 651c58f1c22SToby Isaac } 652c58f1c22SToby Isaac 653c58f1c22SToby Isaac #undef __FUNCT__ 654c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetValueIS" 655c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values) 656c58f1c22SToby Isaac { 657c58f1c22SToby Isaac PetscErrorCode ierr; 658c58f1c22SToby Isaac 659c58f1c22SToby Isaac PetscFunctionBegin; 660c58f1c22SToby Isaac PetscValidPointer(values, 2); 661c58f1c22SToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr); 662c58f1c22SToby Isaac PetscFunctionReturn(0); 663c58f1c22SToby Isaac } 664c58f1c22SToby Isaac 665c58f1c22SToby Isaac #undef __FUNCT__ 666c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetStratumSize" 667c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size) 668c58f1c22SToby Isaac { 669c58f1c22SToby Isaac PetscInt v; 670c58f1c22SToby Isaac PetscErrorCode ierr; 671c58f1c22SToby Isaac 672c58f1c22SToby Isaac PetscFunctionBegin; 673c58f1c22SToby Isaac PetscValidPointer(size, 3); 674c58f1c22SToby Isaac *size = 0; 675c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 676c58f1c22SToby Isaac if (label->stratumValues[v] == value) { 677c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 678c58f1c22SToby Isaac *size = label->stratumSizes[v]; 679c58f1c22SToby Isaac break; 680c58f1c22SToby Isaac } 681c58f1c22SToby Isaac } 682c58f1c22SToby Isaac PetscFunctionReturn(0); 683c58f1c22SToby Isaac } 684c58f1c22SToby Isaac 685c58f1c22SToby Isaac #undef __FUNCT__ 686c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetStratumBounds" 687c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end) 688c58f1c22SToby Isaac { 689c58f1c22SToby Isaac PetscInt v; 690c58f1c22SToby Isaac PetscErrorCode ierr; 691c58f1c22SToby Isaac 692c58f1c22SToby Isaac PetscFunctionBegin; 693c58f1c22SToby Isaac if (start) {PetscValidPointer(start, 3); *start = 0;} 694c58f1c22SToby Isaac if (end) {PetscValidPointer(end, 4); *end = 0;} 695c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 696c58f1c22SToby Isaac if (label->stratumValues[v] != value) continue; 697c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 698c58f1c22SToby Isaac if (label->stratumSizes[v] <= 0) break; 699c58f1c22SToby Isaac if (start) *start = label->points[v][0]; 700c58f1c22SToby Isaac if (end) *end = label->points[v][label->stratumSizes[v]-1]+1; 701c58f1c22SToby Isaac break; 702c58f1c22SToby Isaac } 703c58f1c22SToby Isaac PetscFunctionReturn(0); 704c58f1c22SToby Isaac } 705c58f1c22SToby Isaac 706c58f1c22SToby Isaac #undef __FUNCT__ 707c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetStratumIS" 708c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points) 709c58f1c22SToby Isaac { 710c58f1c22SToby Isaac PetscInt v; 711c58f1c22SToby Isaac PetscErrorCode ierr; 712c58f1c22SToby Isaac 713c58f1c22SToby Isaac PetscFunctionBegin; 714c58f1c22SToby Isaac PetscValidPointer(points, 3); 715c58f1c22SToby Isaac *points = NULL; 716c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 717c58f1c22SToby Isaac if (label->stratumValues[v] == value) { 718c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 719c58f1c22SToby Isaac if (label->arrayValid[v]) { 720c58f1c22SToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF, label->stratumSizes[v], &label->points[v][0], PETSC_COPY_VALUES, points);CHKERRQ(ierr); 721c58f1c22SToby Isaac ierr = PetscObjectSetName((PetscObject) *points, "indices");CHKERRQ(ierr); 722c58f1c22SToby Isaac } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to implement this to speedup Stratify"); 723c58f1c22SToby Isaac break; 724c58f1c22SToby Isaac } 725c58f1c22SToby Isaac } 726c58f1c22SToby Isaac PetscFunctionReturn(0); 727c58f1c22SToby Isaac } 728c58f1c22SToby Isaac 729c58f1c22SToby Isaac #undef __FUNCT__ 730c58f1c22SToby Isaac #define __FUNCT__ "DMLabelClearStratum" 731c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value) 732c58f1c22SToby Isaac { 733c58f1c22SToby Isaac PetscInt v; 734c58f1c22SToby Isaac PetscErrorCode ierr; 735c58f1c22SToby Isaac 736c58f1c22SToby Isaac PetscFunctionBegin; 737c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 738c58f1c22SToby Isaac if (label->stratumValues[v] == value) break; 739c58f1c22SToby Isaac } 740c58f1c22SToby Isaac if (v >= label->numStrata) PetscFunctionReturn(0); 741c58f1c22SToby Isaac if (label->bt) { 742c58f1c22SToby Isaac PetscInt i; 743c58f1c22SToby Isaac 744c58f1c22SToby Isaac for (i = 0; i < label->stratumSizes[v]; ++i) { 745c58f1c22SToby Isaac const PetscInt point = label->points[v][i]; 746c58f1c22SToby Isaac 747c58f1c22SToby 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); 748c58f1c22SToby Isaac ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 749c58f1c22SToby Isaac } 750c58f1c22SToby Isaac } 751c58f1c22SToby Isaac if (label->arrayValid[v]) { 752c58f1c22SToby Isaac label->stratumSizes[v] = 0; 753c58f1c22SToby Isaac } else { 754c58f1c22SToby Isaac PetscHashIClear(label->ht[v]); 755c58f1c22SToby Isaac } 756c58f1c22SToby Isaac PetscFunctionReturn(0); 757c58f1c22SToby Isaac } 758c58f1c22SToby Isaac 759c58f1c22SToby Isaac #undef __FUNCT__ 760c58f1c22SToby Isaac #define __FUNCT__ "DMLabelFilter" 761c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end) 762c58f1c22SToby Isaac { 763c58f1c22SToby Isaac PetscInt v; 764c58f1c22SToby Isaac PetscErrorCode ierr; 765c58f1c22SToby Isaac 766c58f1c22SToby Isaac PetscFunctionBegin; 767c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 768c58f1c22SToby Isaac label->pStart = start; 769c58f1c22SToby Isaac label->pEnd = end; 770c58f1c22SToby Isaac if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);} 771c58f1c22SToby Isaac /* Could squish offsets, but would only make sense if I reallocate the storage */ 772c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 773c58f1c22SToby Isaac PetscInt off, q; 774c58f1c22SToby Isaac 775c58f1c22SToby Isaac for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) { 776c58f1c22SToby Isaac const PetscInt point = label->points[v][q]; 777c58f1c22SToby Isaac 778c58f1c22SToby Isaac if ((point < start) || (point >= end)) continue; 779c58f1c22SToby Isaac label->points[v][off++] = point; 780c58f1c22SToby Isaac } 781c58f1c22SToby Isaac label->stratumSizes[v] = off; 782c58f1c22SToby Isaac } 783c58f1c22SToby Isaac ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr); 784c58f1c22SToby Isaac PetscFunctionReturn(0); 785c58f1c22SToby Isaac } 786c58f1c22SToby Isaac 787c58f1c22SToby Isaac #undef __FUNCT__ 788c58f1c22SToby Isaac #define __FUNCT__ "DMLabelPermute" 789c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew) 790c58f1c22SToby Isaac { 791c58f1c22SToby Isaac const PetscInt *perm; 792c58f1c22SToby Isaac PetscInt numValues, numPoints, v, q; 793c58f1c22SToby Isaac PetscErrorCode ierr; 794c58f1c22SToby Isaac 795c58f1c22SToby Isaac PetscFunctionBegin; 796c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 797c58f1c22SToby Isaac ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr); 798c58f1c22SToby Isaac ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr); 799c58f1c22SToby Isaac ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr); 800c58f1c22SToby Isaac ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr); 801c58f1c22SToby Isaac for (v = 0; v < numValues; ++v) { 802c58f1c22SToby Isaac const PetscInt size = (*labelNew)->stratumSizes[v]; 803c58f1c22SToby Isaac 804c58f1c22SToby Isaac for (q = 0; q < size; ++q) { 805c58f1c22SToby Isaac const PetscInt point = (*labelNew)->points[v][q]; 806c58f1c22SToby Isaac 807c58f1c22SToby 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); 808c58f1c22SToby Isaac (*labelNew)->points[v][q] = perm[point]; 809c58f1c22SToby Isaac } 810c58f1c22SToby Isaac ierr = PetscSortInt(size, &(*labelNew)->points[v][0]);CHKERRQ(ierr); 811c58f1c22SToby Isaac } 812c58f1c22SToby Isaac ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr); 813c58f1c22SToby Isaac if (label->bt) { 814c58f1c22SToby Isaac ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 815c58f1c22SToby Isaac ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr); 816c58f1c22SToby Isaac } 817c58f1c22SToby Isaac PetscFunctionReturn(0); 818c58f1c22SToby Isaac } 819c58f1c22SToby Isaac 820c58f1c22SToby Isaac #undef __FUNCT__ 821c58f1c22SToby Isaac #define __FUNCT__ "DMLabelDistribute" 822c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew) 823c58f1c22SToby Isaac { 824c58f1c22SToby Isaac MPI_Comm comm; 825c58f1c22SToby Isaac PetscSection rootSection, leafSection; 826c58f1c22SToby Isaac PetscSF labelSF; 827c58f1c22SToby Isaac PetscInt p, pStart, pEnd, l, lStart, lEnd, s, nroots, nleaves, size, dof, offset, stratum; 828c58f1c22SToby Isaac PetscInt *remoteOffsets, *rootStrata, *rootIdx, *leafStrata, *strataIdx; 829c58f1c22SToby Isaac char *name; 830c58f1c22SToby Isaac PetscInt nameSize; 831*5cbdf6fcSMichael Lange PetscHashI stratumHash; 832*5cbdf6fcSMichael Lange PETSC_UNUSED PetscHashIIter ret, iter; 833c58f1c22SToby Isaac size_t len = 0; 834c58f1c22SToby Isaac PetscMPIInt rank, numProcs; 835c58f1c22SToby Isaac PetscErrorCode ierr; 836c58f1c22SToby Isaac 837c58f1c22SToby Isaac PetscFunctionBegin; 838c58f1c22SToby Isaac if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 839c58f1c22SToby Isaac ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 840c58f1c22SToby Isaac ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 841c58f1c22SToby Isaac ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 842c58f1c22SToby Isaac /* Bcast name */ 843c58f1c22SToby Isaac if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 844c58f1c22SToby Isaac nameSize = len; 845c58f1c22SToby Isaac ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 846c58f1c22SToby Isaac ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 847c58f1c22SToby Isaac if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 848c58f1c22SToby Isaac ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 849c58f1c22SToby Isaac ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 850c58f1c22SToby Isaac ierr = PetscFree(name);CHKERRQ(ierr); 851c58f1c22SToby Isaac /* Build a section detailing strata-per-point, distribute and build SF 852c58f1c22SToby Isaac from that and then send our points. */ 853c58f1c22SToby Isaac ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 854c58f1c22SToby Isaac ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 855c58f1c22SToby Isaac ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr); 856c58f1c22SToby Isaac if (label) { 857c58f1c22SToby Isaac for (s = 0; s < label->numStrata; ++s) { 858c58f1c22SToby Isaac lStart = 0; 859c58f1c22SToby Isaac lEnd = label->stratumSizes[s]; 860c58f1c22SToby Isaac for (l=lStart; l<lEnd; l++) { 861c58f1c22SToby Isaac ierr = PetscSectionGetDof(rootSection, label->points[s][l], &dof);CHKERRQ(ierr); 862c58f1c22SToby Isaac ierr = PetscSectionSetDof(rootSection, label->points[s][l], dof+1);CHKERRQ(ierr); 863c58f1c22SToby Isaac } 864c58f1c22SToby Isaac } 865c58f1c22SToby Isaac } 866c58f1c22SToby Isaac ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 867c58f1c22SToby Isaac 868c58f1c22SToby Isaac /* Create a point-wise array of point strata */ 869c58f1c22SToby Isaac ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr); 870c58f1c22SToby Isaac ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr); 871c58f1c22SToby Isaac ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr); 872c58f1c22SToby Isaac if (label) { 873c58f1c22SToby Isaac for (s = 0; s < label->numStrata; ++s) { 874c58f1c22SToby Isaac lStart = 0; 875c58f1c22SToby Isaac lEnd = label->stratumSizes[s]; 876c58f1c22SToby Isaac for (l=lStart; l<lEnd; l++) { 877c58f1c22SToby Isaac p = label->points[s][l]; 878c58f1c22SToby Isaac ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr); 879*5cbdf6fcSMichael Lange /* Send stratum values for each point */ 880*5cbdf6fcSMichael Lange rootStrata[offset+rootIdx[p]++] = label->stratumValues[s]; 881c58f1c22SToby Isaac } 882c58f1c22SToby Isaac } 883c58f1c22SToby Isaac } 884c58f1c22SToby Isaac 885c58f1c22SToby Isaac /* Build SF that maps label points to remote processes */ 886c58f1c22SToby Isaac ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr); 887c58f1c22SToby Isaac ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, leafSection);CHKERRQ(ierr); 888c58f1c22SToby Isaac ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, leafSection, &labelSF);CHKERRQ(ierr); 889d2a29c01SMatthew G. Knepley ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 890c58f1c22SToby Isaac 891c58f1c22SToby Isaac /* Send the strata for each point over the derived SF */ 892c58f1c22SToby Isaac ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr); 893c58f1c22SToby Isaac ierr = PetscMalloc1(size, &leafStrata);CHKERRQ(ierr); 894c58f1c22SToby Isaac ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, leafStrata);CHKERRQ(ierr); 895c58f1c22SToby Isaac ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, leafStrata);CHKERRQ(ierr); 896*5cbdf6fcSMichael Lange /* Determine received stratum values and initialise new label*/ 897*5cbdf6fcSMichael Lange PetscHashICreate(stratumHash); 898*5cbdf6fcSMichael Lange for (p = 0; p < size; ++p) PetscHashIPut(stratumHash, leafStrata[p], ret, iter); 899*5cbdf6fcSMichael Lange PetscHashISize(stratumHash, (*labelNew)->numStrata); 900*5cbdf6fcSMichael Lange ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->arrayValid);CHKERRQ(ierr); 901*5cbdf6fcSMichael Lange for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->arrayValid[s] = PETSC_TRUE; 902*5cbdf6fcSMichael Lange ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr); 903*5cbdf6fcSMichael Lange /* Turn leafStrata into indices rather than stratum values */ 904*5cbdf6fcSMichael Lange offset = 0; 905*5cbdf6fcSMichael Lange ierr = PetscHashIGetKeys(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr); 906*5cbdf6fcSMichael Lange for (s = 0; s < (*labelNew)->numStrata; ++s) { 907*5cbdf6fcSMichael Lange for (p = 0; p < size; ++p) { 908*5cbdf6fcSMichael Lange if (leafStrata[p] == (*labelNew)->stratumValues[s]) leafStrata[p] = s; 909*5cbdf6fcSMichael Lange } 910*5cbdf6fcSMichael Lange } 911c58f1c22SToby Isaac /* Rebuild the point strata on the receiver */ 912c58f1c22SToby Isaac ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr); 913c58f1c22SToby Isaac ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 914c58f1c22SToby Isaac for (p=pStart; p<pEnd; p++) { 915c58f1c22SToby Isaac ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 916c58f1c22SToby Isaac ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 917c58f1c22SToby Isaac for (s=0; s<dof; s++) { 918c58f1c22SToby Isaac (*labelNew)->stratumSizes[leafStrata[offset+s]]++; 919c58f1c22SToby Isaac } 920c58f1c22SToby Isaac } 921c58f1c22SToby Isaac ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr); 922c58f1c22SToby Isaac ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr); 923c58f1c22SToby Isaac for (s = 0; s < (*labelNew)->numStrata; ++s) { 924c58f1c22SToby Isaac PetscHashICreate((*labelNew)->ht[s]); 925c58f1c22SToby Isaac ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(*labelNew)->points[s]);CHKERRQ(ierr); 926c58f1c22SToby Isaac } 927c58f1c22SToby Isaac 928c58f1c22SToby Isaac /* Insert points into new strata */ 929c58f1c22SToby Isaac ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr); 930c58f1c22SToby Isaac ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 931c58f1c22SToby Isaac for (p=pStart; p<pEnd; p++) { 932c58f1c22SToby Isaac ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 933c58f1c22SToby Isaac ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 934c58f1c22SToby Isaac for (s=0; s<dof; s++) { 935c58f1c22SToby Isaac stratum = leafStrata[offset+s]; 936c58f1c22SToby Isaac (*labelNew)->points[stratum][strataIdx[stratum]++] = p; 937c58f1c22SToby Isaac } 938c58f1c22SToby Isaac } 939*5cbdf6fcSMichael Lange PetscHashIDestroy(stratumHash); 940c58f1c22SToby Isaac ierr = PetscFree(rootStrata);CHKERRQ(ierr); 941c58f1c22SToby Isaac ierr = PetscFree(leafStrata);CHKERRQ(ierr); 942c58f1c22SToby Isaac ierr = PetscFree(rootIdx);CHKERRQ(ierr); 943c58f1c22SToby Isaac ierr = PetscFree(strataIdx);CHKERRQ(ierr); 944c58f1c22SToby Isaac ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 945c58f1c22SToby Isaac ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 946c58f1c22SToby Isaac ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr); 947c58f1c22SToby Isaac PetscFunctionReturn(0); 948c58f1c22SToby Isaac } 949c58f1c22SToby Isaac 950c58f1c22SToby Isaac #undef __FUNCT__ 951c58f1c22SToby Isaac #define __FUNCT__ "DMLabelConvertToSection" 952c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is) 953c58f1c22SToby Isaac { 954c58f1c22SToby Isaac IS vIS; 955c58f1c22SToby Isaac const PetscInt *values; 956c58f1c22SToby Isaac PetscInt *points; 957c58f1c22SToby Isaac PetscInt nV, vS = 0, vE = 0, v, N; 958c58f1c22SToby Isaac PetscErrorCode ierr; 959c58f1c22SToby Isaac 960c58f1c22SToby Isaac PetscFunctionBegin; 961c58f1c22SToby Isaac ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr); 962c58f1c22SToby Isaac ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr); 963c58f1c22SToby Isaac ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr); 964c58f1c22SToby Isaac if (nV) {vS = values[0]; vE = values[0]+1;} 965c58f1c22SToby Isaac for (v = 1; v < nV; ++v) { 966c58f1c22SToby Isaac vS = PetscMin(vS, values[v]); 967c58f1c22SToby Isaac vE = PetscMax(vE, values[v]+1); 968c58f1c22SToby Isaac } 969c58f1c22SToby Isaac ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr); 970c58f1c22SToby Isaac ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr); 971c58f1c22SToby Isaac for (v = 0; v < nV; ++v) { 972c58f1c22SToby Isaac PetscInt n; 973c58f1c22SToby Isaac 974c58f1c22SToby Isaac ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr); 975c58f1c22SToby Isaac ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr); 976c58f1c22SToby Isaac } 977c58f1c22SToby Isaac ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 978c58f1c22SToby Isaac ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr); 979c58f1c22SToby Isaac ierr = PetscMalloc1(N, &points);CHKERRQ(ierr); 980c58f1c22SToby Isaac for (v = 0; v < nV; ++v) { 981c58f1c22SToby Isaac IS is; 982c58f1c22SToby Isaac const PetscInt *spoints; 983c58f1c22SToby Isaac PetscInt dof, off, p; 984c58f1c22SToby Isaac 985c58f1c22SToby Isaac ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr); 986c58f1c22SToby Isaac ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr); 987c58f1c22SToby Isaac ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr); 988c58f1c22SToby Isaac ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr); 989c58f1c22SToby Isaac for (p = 0; p < dof; ++p) points[off+p] = spoints[p]; 990c58f1c22SToby Isaac ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr); 991c58f1c22SToby Isaac ierr = ISDestroy(&is);CHKERRQ(ierr); 992c58f1c22SToby Isaac } 993c58f1c22SToby Isaac ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr); 994c58f1c22SToby Isaac ierr = ISDestroy(&vIS);CHKERRQ(ierr); 995c58f1c22SToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr); 996c58f1c22SToby Isaac PetscFunctionReturn(0); 997c58f1c22SToby Isaac } 998c58f1c22SToby Isaac 999c58f1c22SToby Isaac #undef __FUNCT__ 1000c58f1c22SToby Isaac #define __FUNCT__ "PetscSectionCreateGlobalSectionLabel" 1001c58f1c22SToby Isaac /*@C 1002c58f1c22SToby Isaac PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using 1003c58f1c22SToby Isaac the local section and an SF describing the section point overlap. 1004c58f1c22SToby Isaac 1005c58f1c22SToby Isaac Input Parameters: 1006c58f1c22SToby Isaac + s - The PetscSection for the local field layout 1007c58f1c22SToby Isaac . sf - The SF describing parallel layout of the section points 1008c58f1c22SToby Isaac . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs 1009c58f1c22SToby Isaac . label - The label specifying the points 1010c58f1c22SToby Isaac - labelValue - The label stratum specifying the points 1011c58f1c22SToby Isaac 1012c58f1c22SToby Isaac Output Parameter: 1013c58f1c22SToby Isaac . gsection - The PetscSection for the global field layout 1014c58f1c22SToby Isaac 1015c58f1c22SToby Isaac Note: This gives negative sizes and offsets to points not owned by this process 1016c58f1c22SToby Isaac 1017c58f1c22SToby Isaac Level: developer 1018c58f1c22SToby Isaac 1019c58f1c22SToby Isaac .seealso: PetscSectionCreate() 1020c58f1c22SToby Isaac @*/ 1021c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection) 1022c58f1c22SToby Isaac { 1023c58f1c22SToby Isaac PetscInt *neg = NULL, *tmpOff = NULL; 1024c58f1c22SToby Isaac PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots; 1025c58f1c22SToby Isaac PetscErrorCode ierr; 1026c58f1c22SToby Isaac 1027c58f1c22SToby Isaac PetscFunctionBegin; 1028c58f1c22SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr); 1029c58f1c22SToby Isaac ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 1030c58f1c22SToby Isaac ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr); 1031c58f1c22SToby Isaac ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 1032c58f1c22SToby Isaac if (nroots >= 0) { 1033c58f1c22SToby Isaac if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart); 1034c58f1c22SToby Isaac ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr); 1035c58f1c22SToby Isaac if (nroots > pEnd-pStart) { 1036c58f1c22SToby Isaac ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr); 1037c58f1c22SToby Isaac } else { 1038c58f1c22SToby Isaac tmpOff = &(*gsection)->atlasDof[-pStart]; 1039c58f1c22SToby Isaac } 1040c58f1c22SToby Isaac } 1041c58f1c22SToby Isaac /* Mark ghost points with negative dof */ 1042c58f1c22SToby Isaac for (p = pStart; p < pEnd; ++p) { 1043c58f1c22SToby Isaac PetscInt value; 1044c58f1c22SToby Isaac 1045c58f1c22SToby Isaac ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr); 1046c58f1c22SToby Isaac if (value != labelValue) continue; 1047c58f1c22SToby Isaac ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr); 1048c58f1c22SToby Isaac ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr); 1049c58f1c22SToby Isaac ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 1050c58f1c22SToby Isaac if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);} 1051c58f1c22SToby Isaac if (neg) neg[p] = -(dof+1); 1052c58f1c22SToby Isaac } 1053c58f1c22SToby Isaac ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr); 1054c58f1c22SToby Isaac if (nroots >= 0) { 1055c58f1c22SToby Isaac ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1056c58f1c22SToby Isaac ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1057c58f1c22SToby Isaac if (nroots > pEnd-pStart) { 1058c58f1c22SToby Isaac for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];} 1059c58f1c22SToby Isaac } 1060c58f1c22SToby Isaac } 1061c58f1c22SToby Isaac /* Calculate new sizes, get proccess offset, and calculate point offsets */ 1062c58f1c22SToby Isaac for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1063c58f1c22SToby Isaac cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0; 1064c58f1c22SToby Isaac (*gsection)->atlasOff[p] = off; 1065c58f1c22SToby Isaac off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0; 1066c58f1c22SToby Isaac } 1067c58f1c22SToby Isaac ierr = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr); 1068c58f1c22SToby Isaac globalOff -= off; 1069c58f1c22SToby Isaac for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1070c58f1c22SToby Isaac (*gsection)->atlasOff[p] += globalOff; 1071c58f1c22SToby Isaac if (neg) neg[p] = -((*gsection)->atlasOff[p]+1); 1072c58f1c22SToby Isaac } 1073c58f1c22SToby Isaac /* Put in negative offsets for ghost points */ 1074c58f1c22SToby Isaac if (nroots >= 0) { 1075c58f1c22SToby Isaac ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1076c58f1c22SToby Isaac ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1077c58f1c22SToby Isaac if (nroots > pEnd-pStart) { 1078c58f1c22SToby Isaac for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];} 1079c58f1c22SToby Isaac } 1080c58f1c22SToby Isaac } 1081c58f1c22SToby Isaac if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);} 1082c58f1c22SToby Isaac ierr = PetscFree(neg);CHKERRQ(ierr); 1083c58f1c22SToby Isaac PetscFunctionReturn(0); 1084c58f1c22SToby Isaac } 1085c58f1c22SToby Isaac 1086