15fdea053SToby Isaac #include <petscdm.h> 2c58f1c22SToby Isaac #include <petsc/private/dmlabelimpl.h> /*I "petscdmlabel.h" I*/ 3c58f1c22SToby Isaac #include <petsc/private/isimpl.h> /*I "petscis.h" I*/ 4c58f1c22SToby Isaac #include <petscsf.h> 5c58f1c22SToby Isaac 6c58f1c22SToby Isaac /*@C 7c58f1c22SToby Isaac DMLabelCreate - Create a DMLabel object, which is a multimap 8c58f1c22SToby Isaac 9d67d17b1SMatthew G. Knepley Input parameters: 10d67d17b1SMatthew G. Knepley + comm - The communicator, usually PETSC_COMM_SELF 11d67d17b1SMatthew G. Knepley - 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 @*/ 20d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelCreate(MPI_Comm comm, const char name[], DMLabel *label) 21c58f1c22SToby Isaac { 22c58f1c22SToby Isaac PetscErrorCode ierr; 23c58f1c22SToby Isaac 24c58f1c22SToby Isaac PetscFunctionBegin; 25d67d17b1SMatthew G. Knepley PetscValidPointer(label,2); 26d67d17b1SMatthew G. Knepley ierr = DMInitializePackage();CHKERRQ(ierr); 27c58f1c22SToby Isaac 28d67d17b1SMatthew G. Knepley ierr = PetscHeaderCreate(*label,DMLABEL_CLASSID,"DMLabel","DMLabel","DM",comm,DMLabelDestroy,DMLabelView);CHKERRQ(ierr); 29d67d17b1SMatthew G. Knepley 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; 40d67d17b1SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *label, name);CHKERRQ(ierr); 41c58f1c22SToby Isaac PetscFunctionReturn(0); 42c58f1c22SToby Isaac } 43c58f1c22SToby Isaac 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 { 600c3c4a36SLisandro Dalcin PetscInt off = 0; 61ad8374ffSToby Isaac PetscInt *pointArray; 62c58f1c22SToby Isaac PetscErrorCode ierr; 63c58f1c22SToby Isaac 64ad8374ffSToby Isaac if (label->validIS[v]) return 0; 65c58f1c22SToby Isaac PetscFunctionBegin; 660c3c4a36SLisandro Dalcin if (v < 0 || v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeValid_Private\n", v); 67e8f14785SLisandro Dalcin ierr = PetscHSetIGetSize(label->ht[v], &label->stratumSizes[v]);CHKERRQ(ierr); 68c58f1c22SToby Isaac 69ad8374ffSToby Isaac ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr); 70e8f14785SLisandro Dalcin ierr = PetscHSetIGetElems(label->ht[v], &off, pointArray);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]); 72e8f14785SLisandro Dalcin PetscHSetIClear(label->ht[v]); 73ad8374ffSToby Isaac ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr); 74c58f1c22SToby Isaac if (label->bt) { 75c58f1c22SToby Isaac PetscInt p; 76c58f1c22SToby Isaac 77c58f1c22SToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) { 78ad8374ffSToby Isaac const PetscInt point = pointArray[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 } 84ad8374ffSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,label->stratumSizes[v],pointArray,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr); 85ad8374ffSToby Isaac ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr); 86ad8374ffSToby Isaac label->validIS[v] = PETSC_TRUE; 87d67d17b1SMatthew G. Knepley ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr); 88c58f1c22SToby Isaac PetscFunctionReturn(0); 89c58f1c22SToby Isaac } 90c58f1c22SToby Isaac 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 /* 117c58f1c22SToby Isaac DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format 118c58f1c22SToby Isaac 119c58f1c22SToby Isaac Input parameter: 120c58f1c22SToby Isaac + label - The DMLabel 121c58f1c22SToby Isaac - v - The stratum value 122c58f1c22SToby Isaac 123c58f1c22SToby Isaac Output parameter: 124c58f1c22SToby Isaac . label - The DMLabel with stratum in hash format 125c58f1c22SToby Isaac 126c58f1c22SToby Isaac Level: developer 127c58f1c22SToby Isaac 128c58f1c22SToby Isaac .seealso: DMLabelCreate() 129c58f1c22SToby Isaac */ 130c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v) 131c58f1c22SToby Isaac { 132c58f1c22SToby Isaac PetscInt p; 133ad8374ffSToby Isaac const PetscInt *points; 134c58f1c22SToby Isaac PetscErrorCode ierr; 135c58f1c22SToby Isaac 1360c3c4a36SLisandro Dalcin if (!label->validIS[v]) return 0; 137c58f1c22SToby Isaac PetscFunctionBegin; 1380c3c4a36SLisandro Dalcin if (v < 0 || v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeInvalid_Private\n", v); 139ad8374ffSToby Isaac if (label->points[v]) { 140ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 141e8f14785SLisandro Dalcin for (p = 0; p < label->stratumSizes[v]; ++p) { 142e8f14785SLisandro Dalcin ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr); 143e8f14785SLisandro Dalcin } 144ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 145ad8374ffSToby Isaac ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 146ad8374ffSToby Isaac } 147ad8374ffSToby Isaac label->validIS[v] = PETSC_FALSE; 148c58f1c22SToby Isaac PetscFunctionReturn(0); 149c58f1c22SToby Isaac } 150c58f1c22SToby Isaac 1510c3c4a36SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index) 1520c3c4a36SLisandro Dalcin { 1530c3c4a36SLisandro Dalcin PetscInt v; 1540c3c4a36SLisandro Dalcin PetscFunctionBegin; 155d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1560c3c4a36SLisandro Dalcin for (*index = -1, v = 0; v < label->numStrata; ++v) 1570c3c4a36SLisandro Dalcin if (label->stratumValues[v] == value) { *index = v; break; } 1580c3c4a36SLisandro Dalcin PetscFunctionReturn(0); 1590c3c4a36SLisandro Dalcin } 1600c3c4a36SLisandro Dalcin 1610c3c4a36SLisandro Dalcin static PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index) 162c58f1c22SToby Isaac { 163ad8374ffSToby Isaac PetscInt v, *tmpV, *tmpS; 164ad8374ffSToby Isaac IS *tmpP; 165e8f14785SLisandro Dalcin PetscHSetI *tmpH; 166c58f1c22SToby Isaac PetscBool *tmpB; 167c58f1c22SToby Isaac PetscErrorCode ierr; 168c58f1c22SToby Isaac 169c58f1c22SToby Isaac PetscFunctionBegin; 170d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 171c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpV);CHKERRQ(ierr); 172c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpS);CHKERRQ(ierr); 173c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpH);CHKERRQ(ierr); 174c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpP);CHKERRQ(ierr); 175c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpB);CHKERRQ(ierr); 176c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 177c58f1c22SToby Isaac tmpV[v] = label->stratumValues[v]; 178c58f1c22SToby Isaac tmpS[v] = label->stratumSizes[v]; 179c58f1c22SToby Isaac tmpH[v] = label->ht[v]; 180c58f1c22SToby Isaac tmpP[v] = label->points[v]; 181ad8374ffSToby Isaac tmpB[v] = label->validIS[v]; 182c58f1c22SToby Isaac } 183c58f1c22SToby Isaac tmpV[v] = value; 184c58f1c22SToby Isaac tmpS[v] = 0; 185e8f14785SLisandro Dalcin ierr = PetscHSetICreate(&tmpH[v]);CHKERRQ(ierr); 186ecb9a371SToby Isaac ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&tmpP[v]);CHKERRQ(ierr); 187c58f1c22SToby Isaac tmpB[v] = PETSC_TRUE; 1880c3c4a36SLisandro Dalcin 189c58f1c22SToby Isaac ++label->numStrata; 190c58f1c22SToby Isaac ierr = PetscFree(label->stratumValues);CHKERRQ(ierr); 191c58f1c22SToby Isaac ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr); 192c58f1c22SToby Isaac ierr = PetscFree(label->ht);CHKERRQ(ierr); 193c58f1c22SToby Isaac ierr = PetscFree(label->points);CHKERRQ(ierr); 194ad8374ffSToby Isaac ierr = PetscFree(label->validIS);CHKERRQ(ierr); 195c58f1c22SToby Isaac label->stratumValues = tmpV; 196c58f1c22SToby Isaac label->stratumSizes = tmpS; 197c58f1c22SToby Isaac label->ht = tmpH; 198c58f1c22SToby Isaac label->points = tmpP; 199ad8374ffSToby Isaac label->validIS = tmpB; 200c58f1c22SToby Isaac 2010c3c4a36SLisandro Dalcin *index = v; 2020c3c4a36SLisandro Dalcin PetscFunctionReturn(0); 2030c3c4a36SLisandro Dalcin } 2040c3c4a36SLisandro Dalcin 2050c3c4a36SLisandro Dalcin PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value) 2060c3c4a36SLisandro Dalcin { 2070c3c4a36SLisandro Dalcin PetscInt v; 2080c3c4a36SLisandro Dalcin PetscErrorCode ierr; 2090c3c4a36SLisandro Dalcin 2100c3c4a36SLisandro Dalcin PetscFunctionBegin; 211d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 2120c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 2130c3c4a36SLisandro Dalcin if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);} 214c58f1c22SToby Isaac PetscFunctionReturn(0); 215c58f1c22SToby Isaac } 216c58f1c22SToby Isaac 217c58f1c22SToby Isaac static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer) 218c58f1c22SToby Isaac { 219c58f1c22SToby Isaac PetscInt v; 220c58f1c22SToby Isaac PetscMPIInt rank; 221c58f1c22SToby Isaac PetscErrorCode ierr; 222c58f1c22SToby Isaac 223c58f1c22SToby Isaac PetscFunctionBegin; 224c58f1c22SToby Isaac ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr); 225c58f1c22SToby Isaac ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 226c58f1c22SToby Isaac if (label) { 227d67d17b1SMatthew G. Knepley const char *name; 228d67d17b1SMatthew G. Knepley 229d67d17b1SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr); 230d67d17b1SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", name);CHKERRQ(ierr); 231c58f1c22SToby Isaac if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, " Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);} 232c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 233c58f1c22SToby Isaac const PetscInt value = label->stratumValues[v]; 234ad8374ffSToby Isaac const PetscInt *points; 235c58f1c22SToby Isaac PetscInt p; 236c58f1c22SToby Isaac 237ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 238c58f1c22SToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) { 239ad8374ffSToby Isaac ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr); 240c58f1c22SToby Isaac } 241ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 242c58f1c22SToby Isaac } 243c58f1c22SToby Isaac } 244c58f1c22SToby Isaac ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 245c58f1c22SToby Isaac ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 246c58f1c22SToby Isaac PetscFunctionReturn(0); 247c58f1c22SToby Isaac } 248c58f1c22SToby Isaac 249c58f1c22SToby Isaac /*@C 250c58f1c22SToby Isaac DMLabelView - View the label 251c58f1c22SToby Isaac 252c58f1c22SToby Isaac Input Parameters: 253c58f1c22SToby Isaac + label - The DMLabel 254c58f1c22SToby Isaac - viewer - The PetscViewer 255c58f1c22SToby Isaac 256c58f1c22SToby Isaac Level: intermediate 257c58f1c22SToby Isaac 258c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy() 259c58f1c22SToby Isaac @*/ 260c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer) 261c58f1c22SToby Isaac { 262c58f1c22SToby Isaac PetscBool iascii; 263c58f1c22SToby Isaac PetscErrorCode ierr; 264c58f1c22SToby Isaac 265c58f1c22SToby Isaac PetscFunctionBegin; 266d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 267c58f1c22SToby Isaac PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 268c58f1c22SToby Isaac if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 269c58f1c22SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 270c58f1c22SToby Isaac if (iascii) { 271c58f1c22SToby Isaac ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr); 272c58f1c22SToby Isaac } 273c58f1c22SToby Isaac PetscFunctionReturn(0); 274c58f1c22SToby Isaac } 275c58f1c22SToby Isaac 27684f0b6dfSMatthew G. Knepley /*@ 277d67d17b1SMatthew G. Knepley DMLabelReset - Destroys internal data structures in a DMLabel 278d67d17b1SMatthew G. Knepley 279d67d17b1SMatthew G. Knepley Input Parameter: 280d67d17b1SMatthew G. Knepley . label - The DMLabel 281d67d17b1SMatthew G. Knepley 282d67d17b1SMatthew G. Knepley Level: beginner 283d67d17b1SMatthew G. Knepley 284d67d17b1SMatthew G. Knepley .seealso: DMLabelDestroy(), DMLabelCreate() 285d67d17b1SMatthew G. Knepley @*/ 286d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelReset(DMLabel label) 287d67d17b1SMatthew G. Knepley { 288d67d17b1SMatthew G. Knepley PetscInt v; 289d67d17b1SMatthew G. Knepley PetscErrorCode ierr; 290d67d17b1SMatthew G. Knepley 291d67d17b1SMatthew G. Knepley PetscFunctionBegin; 292d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 293d67d17b1SMatthew G. Knepley ierr = PetscFree(label->stratumValues);CHKERRQ(ierr); 294d67d17b1SMatthew G. Knepley ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr); 295d67d17b1SMatthew G. Knepley for (v = 0; v < label->numStrata; ++v) { 296d67d17b1SMatthew G. Knepley ierr = PetscHSetIDestroy(&label->ht[v]);CHKERRQ(ierr); 297d67d17b1SMatthew G. Knepley ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr); 298d67d17b1SMatthew G. Knepley } 299d67d17b1SMatthew G. Knepley ierr = PetscFree(label->ht);CHKERRQ(ierr); 300d67d17b1SMatthew G. Knepley ierr = PetscFree(label->points);CHKERRQ(ierr); 301d67d17b1SMatthew G. Knepley ierr = PetscFree(label->validIS);CHKERRQ(ierr); 302d67d17b1SMatthew G. Knepley ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 303d67d17b1SMatthew G. Knepley PetscFunctionReturn(0); 304d67d17b1SMatthew G. Knepley } 305d67d17b1SMatthew G. Knepley 306d67d17b1SMatthew G. Knepley /*@ 30784f0b6dfSMatthew G. Knepley DMLabelDestroy - Destroys a DMLabel 30884f0b6dfSMatthew G. Knepley 30984f0b6dfSMatthew G. Knepley Input Parameter: 31084f0b6dfSMatthew G. Knepley . label - The DMLabel 31184f0b6dfSMatthew G. Knepley 31284f0b6dfSMatthew G. Knepley Level: beginner 31384f0b6dfSMatthew G. Knepley 314d67d17b1SMatthew G. Knepley .seealso: DMLabelReset(), DMLabelCreate() 31584f0b6dfSMatthew G. Knepley @*/ 316c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label) 317c58f1c22SToby Isaac { 318c58f1c22SToby Isaac PetscErrorCode ierr; 319c58f1c22SToby Isaac 320c58f1c22SToby Isaac PetscFunctionBegin; 321d67d17b1SMatthew G. Knepley if (!*label) PetscFunctionReturn(0); 322d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific((*label),DMLABEL_CLASSID,1); 323d67d17b1SMatthew G. Knepley if (--((PetscObject)(*label))->refct > 0) { 324d67d17b1SMatthew G. Knepley *label = NULL; 325d67d17b1SMatthew G. Knepley PetscFunctionReturn(0); 3260c3c4a36SLisandro Dalcin } 327d67d17b1SMatthew G. Knepley ierr = DMLabelReset(*label);CHKERRQ(ierr); 328d67d17b1SMatthew G. Knepley ierr = PetscHeaderDestroy(label);CHKERRQ(ierr); 329c58f1c22SToby Isaac PetscFunctionReturn(0); 330c58f1c22SToby Isaac } 331c58f1c22SToby Isaac 33284f0b6dfSMatthew G. Knepley /*@ 33384f0b6dfSMatthew G. Knepley DMLabelDuplicate - Duplicates a DMLabel 33484f0b6dfSMatthew G. Knepley 33584f0b6dfSMatthew G. Knepley Input Parameter: 33684f0b6dfSMatthew G. Knepley . label - The DMLabel 33784f0b6dfSMatthew G. Knepley 33884f0b6dfSMatthew G. Knepley Output Parameter: 33984f0b6dfSMatthew G. Knepley . labelnew - location to put new vector 34084f0b6dfSMatthew G. Knepley 34184f0b6dfSMatthew G. Knepley Level: intermediate 34284f0b6dfSMatthew G. Knepley 34384f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelDestroy() 34484f0b6dfSMatthew G. Knepley @*/ 345c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew) 346c58f1c22SToby Isaac { 347d67d17b1SMatthew G. Knepley const char *name; 348ad8374ffSToby Isaac PetscInt v; 349c58f1c22SToby Isaac PetscErrorCode ierr; 350c58f1c22SToby Isaac 351c58f1c22SToby Isaac PetscFunctionBegin; 352d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 353c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 354d67d17b1SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr); 355d67d17b1SMatthew G. Knepley ierr = DMLabelCreate(PetscObjectComm((PetscObject) label), name, labelnew);CHKERRQ(ierr); 356c58f1c22SToby Isaac 357c58f1c22SToby Isaac (*labelnew)->numStrata = label->numStrata; 3585aa44df4SToby Isaac (*labelnew)->defaultValue = label->defaultValue; 359c58f1c22SToby Isaac if (label->numStrata) { 360c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr); 361c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr); 362c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr); 363c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr); 364ad8374ffSToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr); 365c58f1c22SToby Isaac /* Could eliminate unused space here */ 366c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 367e8f14785SLisandro Dalcin ierr = PetscHSetICreate(&(*labelnew)->ht[v]);CHKERRQ(ierr); 368ad8374ffSToby Isaac (*labelnew)->validIS[v] = PETSC_TRUE; 369c58f1c22SToby Isaac (*labelnew)->stratumValues[v] = label->stratumValues[v]; 370c58f1c22SToby Isaac (*labelnew)->stratumSizes[v] = label->stratumSizes[v]; 371ad8374ffSToby Isaac ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr); 372ad8374ffSToby Isaac (*labelnew)->points[v] = label->points[v]; 373c58f1c22SToby Isaac } 374c58f1c22SToby Isaac } 375c58f1c22SToby Isaac (*labelnew)->pStart = -1; 376c58f1c22SToby Isaac (*labelnew)->pEnd = -1; 377c58f1c22SToby Isaac (*labelnew)->bt = NULL; 378c58f1c22SToby Isaac PetscFunctionReturn(0); 379c58f1c22SToby Isaac } 380c58f1c22SToby Isaac 381c58f1c22SToby Isaac /* This can be hooked into SetValue(), ClearValue(), etc. for updating */ 382c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd) 383c58f1c22SToby Isaac { 384c58f1c22SToby Isaac PetscInt v; 385c58f1c22SToby Isaac PetscErrorCode ierr; 386c58f1c22SToby Isaac 387c58f1c22SToby Isaac PetscFunctionBegin; 388d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 3890c3c4a36SLisandro Dalcin ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr); 390c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 391c58f1c22SToby Isaac label->pStart = pStart; 392c58f1c22SToby Isaac label->pEnd = pEnd; 393c58f1c22SToby Isaac ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr); 394c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 395ad8374ffSToby Isaac const PetscInt *points; 396c58f1c22SToby Isaac PetscInt i; 397c58f1c22SToby Isaac 398ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 399c58f1c22SToby Isaac for (i = 0; i < label->stratumSizes[v]; ++i) { 400ad8374ffSToby Isaac const PetscInt point = points[i]; 401c58f1c22SToby Isaac 402c58f1c22SToby 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); 403c58f1c22SToby Isaac ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr); 404c58f1c22SToby Isaac } 405ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr); 406c58f1c22SToby Isaac } 407c58f1c22SToby Isaac PetscFunctionReturn(0); 408c58f1c22SToby Isaac } 409c58f1c22SToby Isaac 410c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label) 411c58f1c22SToby Isaac { 412c58f1c22SToby Isaac PetscErrorCode ierr; 413c58f1c22SToby Isaac 414c58f1c22SToby Isaac PetscFunctionBegin; 415d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 416c58f1c22SToby Isaac label->pStart = -1; 417c58f1c22SToby Isaac label->pEnd = -1; 4180c3c4a36SLisandro Dalcin ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 419c58f1c22SToby Isaac PetscFunctionReturn(0); 420c58f1c22SToby Isaac } 421c58f1c22SToby Isaac 422c58f1c22SToby Isaac /*@ 423c58f1c22SToby Isaac DMLabelHasValue - Determine whether a label assigns the value to any point 424c58f1c22SToby Isaac 425c58f1c22SToby Isaac Input Parameters: 426c58f1c22SToby Isaac + label - the DMLabel 427c58f1c22SToby Isaac - value - the value 428c58f1c22SToby Isaac 429c58f1c22SToby Isaac Output Parameter: 430c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point 431c58f1c22SToby Isaac 432c58f1c22SToby Isaac Level: developer 433c58f1c22SToby Isaac 434c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue() 435c58f1c22SToby Isaac @*/ 436c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains) 437c58f1c22SToby Isaac { 438c58f1c22SToby Isaac PetscInt v; 4390c3c4a36SLisandro Dalcin PetscErrorCode ierr; 440c58f1c22SToby Isaac 441c58f1c22SToby Isaac PetscFunctionBegin; 442d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 443c58f1c22SToby Isaac PetscValidPointer(contains, 3); 4440c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 4450c3c4a36SLisandro Dalcin *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE; 446c58f1c22SToby Isaac PetscFunctionReturn(0); 447c58f1c22SToby Isaac } 448c58f1c22SToby Isaac 449c58f1c22SToby Isaac /*@ 450c58f1c22SToby Isaac DMLabelHasPoint - Determine whether a label assigns a value to a point 451c58f1c22SToby Isaac 452c58f1c22SToby Isaac Input Parameters: 453c58f1c22SToby Isaac + label - the DMLabel 454c58f1c22SToby Isaac - point - the point 455c58f1c22SToby Isaac 456c58f1c22SToby Isaac Output Parameter: 457c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value 458c58f1c22SToby Isaac 459c58f1c22SToby Isaac Note: The user must call DMLabelCreateIndex() before this function. 460c58f1c22SToby Isaac 461c58f1c22SToby Isaac Level: developer 462c58f1c22SToby Isaac 463c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue() 464c58f1c22SToby Isaac @*/ 465c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains) 466c58f1c22SToby Isaac { 467c58f1c22SToby Isaac PetscErrorCode ierr; 468c58f1c22SToby Isaac 469c58f1c22SToby Isaac PetscFunctionBeginHot; 470d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 471c58f1c22SToby Isaac PetscValidPointer(contains, 3); 472c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 473c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG) 474c58f1c22SToby Isaac if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()"); 475c58f1c22SToby 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); 476c58f1c22SToby Isaac #endif 477c58f1c22SToby Isaac *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE; 478c58f1c22SToby Isaac PetscFunctionReturn(0); 479c58f1c22SToby Isaac } 480c58f1c22SToby Isaac 481c58f1c22SToby Isaac /*@ 482c58f1c22SToby Isaac DMLabelStratumHasPoint - Return true if the stratum contains a point 483c58f1c22SToby Isaac 484c58f1c22SToby Isaac Input Parameters: 485c58f1c22SToby Isaac + label - the DMLabel 486c58f1c22SToby Isaac . value - the stratum value 487c58f1c22SToby Isaac - point - the point 488c58f1c22SToby Isaac 489c58f1c22SToby Isaac Output Parameter: 490c58f1c22SToby Isaac . contains - true if the stratum contains the point 491c58f1c22SToby Isaac 492c58f1c22SToby Isaac Level: intermediate 493c58f1c22SToby Isaac 494c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue() 495c58f1c22SToby Isaac @*/ 496c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains) 497c58f1c22SToby Isaac { 498c58f1c22SToby Isaac PetscInt v; 499c58f1c22SToby Isaac PetscErrorCode ierr; 500c58f1c22SToby Isaac 5010c3c4a36SLisandro Dalcin PetscFunctionBeginHot; 502d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 503c58f1c22SToby Isaac PetscValidPointer(contains, 4); 504c58f1c22SToby Isaac *contains = PETSC_FALSE; 5050c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 5060c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 5070c3c4a36SLisandro Dalcin 508ad8374ffSToby Isaac if (label->validIS[v]) { 509c58f1c22SToby Isaac PetscInt i; 510c58f1c22SToby Isaac 511a2d74346SToby Isaac ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr); 5120c3c4a36SLisandro Dalcin if (i >= 0) *contains = PETSC_TRUE; 513c58f1c22SToby Isaac } else { 514c58f1c22SToby Isaac PetscBool has; 515c58f1c22SToby Isaac 516e8f14785SLisandro Dalcin PetscHSetIHas(label->ht[v], point, &has); 5170c3c4a36SLisandro Dalcin if (has) *contains = PETSC_TRUE; 518c58f1c22SToby Isaac } 519c58f1c22SToby Isaac PetscFunctionReturn(0); 520c58f1c22SToby Isaac } 521c58f1c22SToby Isaac 52284f0b6dfSMatthew G. Knepley /*@ 5235aa44df4SToby Isaac DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 5245aa44df4SToby Isaac When a label is created, it is initialized to -1. 5255aa44df4SToby Isaac 5265aa44df4SToby Isaac Input parameter: 5275aa44df4SToby Isaac . label - a DMLabel object 5285aa44df4SToby Isaac 5295aa44df4SToby Isaac Output parameter: 5305aa44df4SToby Isaac . defaultValue - the default value 5315aa44df4SToby Isaac 5325aa44df4SToby Isaac Level: beginner 5335aa44df4SToby Isaac 5345aa44df4SToby Isaac .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 53584f0b6dfSMatthew G. Knepley @*/ 5365aa44df4SToby Isaac PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue) 5375aa44df4SToby Isaac { 5385aa44df4SToby Isaac PetscFunctionBegin; 539d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 5405aa44df4SToby Isaac *defaultValue = label->defaultValue; 5415aa44df4SToby Isaac PetscFunctionReturn(0); 5425aa44df4SToby Isaac } 5435aa44df4SToby Isaac 54484f0b6dfSMatthew G. Knepley /*@ 5455aa44df4SToby Isaac DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 5465aa44df4SToby Isaac When a label is created, it is initialized to -1. 5475aa44df4SToby Isaac 5485aa44df4SToby Isaac Input parameter: 5495aa44df4SToby Isaac . label - a DMLabel object 5505aa44df4SToby Isaac 5515aa44df4SToby Isaac Output parameter: 5525aa44df4SToby Isaac . defaultValue - the default value 5535aa44df4SToby Isaac 5545aa44df4SToby Isaac Level: beginner 5555aa44df4SToby Isaac 5565aa44df4SToby Isaac .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 55784f0b6dfSMatthew G. Knepley @*/ 5585aa44df4SToby Isaac PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue) 5595aa44df4SToby Isaac { 5605aa44df4SToby Isaac PetscFunctionBegin; 561d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 5625aa44df4SToby Isaac label->defaultValue = defaultValue; 5635aa44df4SToby Isaac PetscFunctionReturn(0); 5645aa44df4SToby Isaac } 5655aa44df4SToby Isaac 566c58f1c22SToby Isaac /*@ 5675aa44df4SToby 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()) 568c58f1c22SToby Isaac 569c58f1c22SToby Isaac Input Parameters: 570c58f1c22SToby Isaac + label - the DMLabel 571c58f1c22SToby Isaac - point - the point 572c58f1c22SToby Isaac 573c58f1c22SToby Isaac Output Parameter: 5748e68afcfSMatthew G. Knepley . value - The point value, or the default value (-1 by default) 575c58f1c22SToby Isaac 576c58f1c22SToby Isaac Level: intermediate 577c58f1c22SToby Isaac 5785aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 579c58f1c22SToby Isaac @*/ 580c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value) 581c58f1c22SToby Isaac { 582c58f1c22SToby Isaac PetscInt v; 583c58f1c22SToby Isaac PetscErrorCode ierr; 584c58f1c22SToby Isaac 5850c3c4a36SLisandro Dalcin PetscFunctionBeginHot; 586d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 587c58f1c22SToby Isaac PetscValidPointer(value, 3); 5885aa44df4SToby Isaac *value = label->defaultValue; 589c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 590ad8374ffSToby Isaac if (label->validIS[v]) { 591c58f1c22SToby Isaac PetscInt i; 592c58f1c22SToby Isaac 593a2d74346SToby Isaac ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr); 594c58f1c22SToby Isaac if (i >= 0) { 595c58f1c22SToby Isaac *value = label->stratumValues[v]; 596c58f1c22SToby Isaac break; 597c58f1c22SToby Isaac } 598c58f1c22SToby Isaac } else { 599c58f1c22SToby Isaac PetscBool has; 600c58f1c22SToby Isaac 601e8f14785SLisandro Dalcin PetscHSetIHas(label->ht[v], point, &has); 602c58f1c22SToby Isaac if (has) { 603c58f1c22SToby Isaac *value = label->stratumValues[v]; 604c58f1c22SToby Isaac break; 605c58f1c22SToby Isaac } 606c58f1c22SToby Isaac } 607c58f1c22SToby Isaac } 608c58f1c22SToby Isaac PetscFunctionReturn(0); 609c58f1c22SToby Isaac } 610c58f1c22SToby Isaac 611c58f1c22SToby Isaac /*@ 612367003a6SStefano Zampini DMLabelSetValue - Set the value a label assigns to a point. If the value is the same as the label's default value (which is initially -1, and can be changed with DMLabelSetDefaultValue() to something different), then this function will do nothing. 613c58f1c22SToby Isaac 614c58f1c22SToby Isaac Input Parameters: 615c58f1c22SToby Isaac + label - the DMLabel 616c58f1c22SToby Isaac . point - the point 617c58f1c22SToby Isaac - value - The point value 618c58f1c22SToby Isaac 619c58f1c22SToby Isaac Level: intermediate 620c58f1c22SToby Isaac 6215aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 622c58f1c22SToby Isaac @*/ 623c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value) 624c58f1c22SToby Isaac { 625c58f1c22SToby Isaac PetscInt v; 626c58f1c22SToby Isaac PetscErrorCode ierr; 627c58f1c22SToby Isaac 628c58f1c22SToby Isaac PetscFunctionBegin; 629d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 6300c3c4a36SLisandro Dalcin /* Find label value, add new entry if needed */ 6315aa44df4SToby Isaac if (value == label->defaultValue) PetscFunctionReturn(0); 6320c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 6330c3c4a36SLisandro Dalcin if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);} 634c58f1c22SToby Isaac /* Set key */ 6350c3c4a36SLisandro Dalcin ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 636e8f14785SLisandro Dalcin ierr = PetscHSetIAdd(label->ht[v], point);CHKERRQ(ierr); 637c58f1c22SToby Isaac PetscFunctionReturn(0); 638c58f1c22SToby Isaac } 639c58f1c22SToby Isaac 640c58f1c22SToby Isaac /*@ 641c58f1c22SToby Isaac DMLabelClearValue - Clear the value a label assigns to a point 642c58f1c22SToby Isaac 643c58f1c22SToby Isaac Input Parameters: 644c58f1c22SToby Isaac + label - the DMLabel 645c58f1c22SToby Isaac . point - the point 646c58f1c22SToby Isaac - value - The point value 647c58f1c22SToby Isaac 648c58f1c22SToby Isaac Level: intermediate 649c58f1c22SToby Isaac 650c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue() 651c58f1c22SToby Isaac @*/ 652c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value) 653c58f1c22SToby Isaac { 654ad8374ffSToby Isaac PetscInt v; 655c58f1c22SToby Isaac PetscErrorCode ierr; 656c58f1c22SToby Isaac 657c58f1c22SToby Isaac PetscFunctionBegin; 658d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 659c58f1c22SToby Isaac /* Find label value */ 6600c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 6610c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 6620c3c4a36SLisandro Dalcin 663eeed21e7SToby Isaac if (label->bt) { 664eeed21e7SToby 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); 665eeed21e7SToby Isaac ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 666eeed21e7SToby Isaac } 6670c3c4a36SLisandro Dalcin 6680c3c4a36SLisandro Dalcin /* Delete key */ 6690c3c4a36SLisandro Dalcin ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 670e8f14785SLisandro Dalcin ierr = PetscHSetIDel(label->ht[v], point);CHKERRQ(ierr); 671c58f1c22SToby Isaac PetscFunctionReturn(0); 672c58f1c22SToby Isaac } 673c58f1c22SToby Isaac 674c58f1c22SToby Isaac /*@ 675c58f1c22SToby Isaac DMLabelInsertIS - Set all points in the IS to a value 676c58f1c22SToby Isaac 677c58f1c22SToby Isaac Input Parameters: 678c58f1c22SToby Isaac + label - the DMLabel 679c58f1c22SToby Isaac . is - the point IS 680c58f1c22SToby Isaac - value - The point value 681c58f1c22SToby Isaac 682c58f1c22SToby Isaac Level: intermediate 683c58f1c22SToby Isaac 684c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 685c58f1c22SToby Isaac @*/ 686c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value) 687c58f1c22SToby Isaac { 6880c3c4a36SLisandro Dalcin PetscInt v, n, p; 689c58f1c22SToby Isaac const PetscInt *points; 690c58f1c22SToby Isaac PetscErrorCode ierr; 691c58f1c22SToby Isaac 692c58f1c22SToby Isaac PetscFunctionBegin; 693d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 694c58f1c22SToby Isaac PetscValidHeaderSpecific(is, IS_CLASSID, 2); 6950c3c4a36SLisandro Dalcin /* Find label value, add new entry if needed */ 6960c3c4a36SLisandro Dalcin if (value == label->defaultValue) PetscFunctionReturn(0); 6970c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 6980c3c4a36SLisandro Dalcin if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);} 6990c3c4a36SLisandro Dalcin /* Set keys */ 7000c3c4a36SLisandro Dalcin ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 701c58f1c22SToby Isaac ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 702c58f1c22SToby Isaac ierr = ISGetIndices(is, &points);CHKERRQ(ierr); 703e8f14785SLisandro Dalcin for (p = 0; p < n; ++p) {ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);} 704c58f1c22SToby Isaac ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr); 705c58f1c22SToby Isaac PetscFunctionReturn(0); 706c58f1c22SToby Isaac } 707c58f1c22SToby Isaac 70884f0b6dfSMatthew G. Knepley /*@ 70984f0b6dfSMatthew G. Knepley DMLabelGetNumValues - Get the number of values that the DMLabel takes 71084f0b6dfSMatthew G. Knepley 71184f0b6dfSMatthew G. Knepley Input Parameter: 71284f0b6dfSMatthew G. Knepley . label - the DMLabel 71384f0b6dfSMatthew G. Knepley 71484f0b6dfSMatthew G. Knepley Output Paramater: 71584f0b6dfSMatthew G. Knepley . numValues - the number of values 71684f0b6dfSMatthew G. Knepley 71784f0b6dfSMatthew G. Knepley Level: intermediate 71884f0b6dfSMatthew G. Knepley 71984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 72084f0b6dfSMatthew G. Knepley @*/ 721c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues) 722c58f1c22SToby Isaac { 723c58f1c22SToby Isaac PetscFunctionBegin; 724d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 725c58f1c22SToby Isaac PetscValidPointer(numValues, 2); 726c58f1c22SToby Isaac *numValues = label->numStrata; 727c58f1c22SToby Isaac PetscFunctionReturn(0); 728c58f1c22SToby Isaac } 729c58f1c22SToby Isaac 73084f0b6dfSMatthew G. Knepley /*@ 73184f0b6dfSMatthew G. Knepley DMLabelGetValueIS - Get an IS of all values that the DMlabel takes 73284f0b6dfSMatthew G. Knepley 73384f0b6dfSMatthew G. Knepley Input Parameter: 73484f0b6dfSMatthew G. Knepley . label - the DMLabel 73584f0b6dfSMatthew G. Knepley 73684f0b6dfSMatthew G. Knepley Output Paramater: 73784f0b6dfSMatthew G. Knepley . is - the value IS 73884f0b6dfSMatthew G. Knepley 73984f0b6dfSMatthew G. Knepley Level: intermediate 74084f0b6dfSMatthew G. Knepley 74184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 74284f0b6dfSMatthew G. Knepley @*/ 743c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values) 744c58f1c22SToby Isaac { 745c58f1c22SToby Isaac PetscErrorCode ierr; 746c58f1c22SToby Isaac 747c58f1c22SToby Isaac PetscFunctionBegin; 748d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 749c58f1c22SToby Isaac PetscValidPointer(values, 2); 750c58f1c22SToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr); 751c58f1c22SToby Isaac PetscFunctionReturn(0); 752c58f1c22SToby Isaac } 753c58f1c22SToby Isaac 75484f0b6dfSMatthew G. Knepley /*@ 75584f0b6dfSMatthew G. Knepley DMLabelHasStratum - Determine whether points exist with the given value 75684f0b6dfSMatthew G. Knepley 75784f0b6dfSMatthew G. Knepley Input Parameters: 75884f0b6dfSMatthew G. Knepley + label - the DMLabel 75984f0b6dfSMatthew G. Knepley - value - the stratum value 76084f0b6dfSMatthew G. Knepley 76184f0b6dfSMatthew G. Knepley Output Paramater: 76284f0b6dfSMatthew G. Knepley . exists - Flag saying whether points exist 76384f0b6dfSMatthew G. Knepley 76484f0b6dfSMatthew G. Knepley Level: intermediate 76584f0b6dfSMatthew G. Knepley 76684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 76784f0b6dfSMatthew G. Knepley @*/ 768fada774cSMatthew G. Knepley PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists) 769fada774cSMatthew G. Knepley { 770fada774cSMatthew G. Knepley PetscInt v; 7710c3c4a36SLisandro Dalcin PetscErrorCode ierr; 772fada774cSMatthew G. Knepley 773fada774cSMatthew G. Knepley PetscFunctionBegin; 774d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 775fada774cSMatthew G. Knepley PetscValidPointer(exists, 3); 7760c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 7770c3c4a36SLisandro Dalcin *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE; 778fada774cSMatthew G. Knepley PetscFunctionReturn(0); 779fada774cSMatthew G. Knepley } 780fada774cSMatthew G. Knepley 78184f0b6dfSMatthew G. Knepley /*@ 78284f0b6dfSMatthew G. Knepley DMLabelGetStratumSize - Get the size of a stratum 78384f0b6dfSMatthew G. Knepley 78484f0b6dfSMatthew G. Knepley Input Parameters: 78584f0b6dfSMatthew G. Knepley + label - the DMLabel 78684f0b6dfSMatthew G. Knepley - value - the stratum value 78784f0b6dfSMatthew G. Knepley 78884f0b6dfSMatthew G. Knepley Output Paramater: 78984f0b6dfSMatthew G. Knepley . size - The number of points in the stratum 79084f0b6dfSMatthew G. Knepley 79184f0b6dfSMatthew G. Knepley Level: intermediate 79284f0b6dfSMatthew G. Knepley 79384f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 79484f0b6dfSMatthew G. Knepley @*/ 795c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size) 796c58f1c22SToby Isaac { 797c58f1c22SToby Isaac PetscInt v; 798c58f1c22SToby Isaac PetscErrorCode ierr; 799c58f1c22SToby Isaac 800c58f1c22SToby Isaac PetscFunctionBegin; 801d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 802c58f1c22SToby Isaac PetscValidPointer(size, 3); 803c58f1c22SToby Isaac *size = 0; 8040c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 8050c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 806c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 807c58f1c22SToby Isaac *size = label->stratumSizes[v]; 808c58f1c22SToby Isaac PetscFunctionReturn(0); 809c58f1c22SToby Isaac } 810c58f1c22SToby Isaac 81184f0b6dfSMatthew G. Knepley /*@ 81284f0b6dfSMatthew G. Knepley DMLabelGetStratumBounds - Get the largest and smallest point of a stratum 81384f0b6dfSMatthew G. Knepley 81484f0b6dfSMatthew G. Knepley Input Parameters: 81584f0b6dfSMatthew G. Knepley + label - the DMLabel 81684f0b6dfSMatthew G. Knepley - value - the stratum value 81784f0b6dfSMatthew G. Knepley 81884f0b6dfSMatthew G. Knepley Output Paramaters: 81984f0b6dfSMatthew G. Knepley + start - the smallest point in the stratum 82084f0b6dfSMatthew G. Knepley - end - the largest point in the stratum 82184f0b6dfSMatthew G. Knepley 82284f0b6dfSMatthew G. Knepley Level: intermediate 82384f0b6dfSMatthew G. Knepley 82484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 82584f0b6dfSMatthew G. Knepley @*/ 826c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end) 827c58f1c22SToby Isaac { 8280c3c4a36SLisandro Dalcin PetscInt v, min, max; 829c58f1c22SToby Isaac PetscErrorCode ierr; 830c58f1c22SToby Isaac 831c58f1c22SToby Isaac PetscFunctionBegin; 832d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 833c58f1c22SToby Isaac if (start) {PetscValidPointer(start, 3); *start = 0;} 834c58f1c22SToby Isaac if (end) {PetscValidPointer(end, 4); *end = 0;} 8350c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 8360c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 837c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 8380c3c4a36SLisandro Dalcin if (label->stratumSizes[v] <= 0) PetscFunctionReturn(0); 839d6cb179aSToby Isaac ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr); 840d6cb179aSToby Isaac if (start) *start = min; 841d6cb179aSToby Isaac if (end) *end = max+1; 842c58f1c22SToby Isaac PetscFunctionReturn(0); 843c58f1c22SToby Isaac } 844c58f1c22SToby Isaac 84584f0b6dfSMatthew G. Knepley /*@ 84684f0b6dfSMatthew G. Knepley DMLabelGetStratumIS - Get an IS with the stratum points 84784f0b6dfSMatthew G. Knepley 84884f0b6dfSMatthew G. Knepley Input Parameters: 84984f0b6dfSMatthew G. Knepley + label - the DMLabel 85084f0b6dfSMatthew G. Knepley - value - the stratum value 85184f0b6dfSMatthew G. Knepley 85284f0b6dfSMatthew G. Knepley Output Paramater: 85384f0b6dfSMatthew G. Knepley . points - The stratum points 85484f0b6dfSMatthew G. Knepley 85584f0b6dfSMatthew G. Knepley Level: intermediate 85684f0b6dfSMatthew G. Knepley 85784f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 85884f0b6dfSMatthew G. Knepley @*/ 859c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points) 860c58f1c22SToby Isaac { 861c58f1c22SToby Isaac PetscInt v; 862c58f1c22SToby Isaac PetscErrorCode ierr; 863c58f1c22SToby Isaac 864c58f1c22SToby Isaac PetscFunctionBegin; 865d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 866c58f1c22SToby Isaac PetscValidPointer(points, 3); 867c58f1c22SToby Isaac *points = NULL; 8680c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 8690c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 870c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 871ad8374ffSToby Isaac ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr); 872ad8374ffSToby Isaac *points = label->points[v]; 873c58f1c22SToby Isaac PetscFunctionReturn(0); 874c58f1c22SToby Isaac } 875c58f1c22SToby Isaac 87684f0b6dfSMatthew G. Knepley /*@ 87784f0b6dfSMatthew G. Knepley DMLabelSetStratumIS - Set the stratum points using an IS 87884f0b6dfSMatthew G. Knepley 87984f0b6dfSMatthew G. Knepley Input Parameters: 88084f0b6dfSMatthew G. Knepley + label - the DMLabel 88184f0b6dfSMatthew G. Knepley . value - the stratum value 88284f0b6dfSMatthew G. Knepley - points - The stratum points 88384f0b6dfSMatthew G. Knepley 88484f0b6dfSMatthew G. Knepley Level: intermediate 88584f0b6dfSMatthew G. Knepley 88684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 88784f0b6dfSMatthew G. Knepley @*/ 8884de306b1SToby Isaac PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is) 8894de306b1SToby Isaac { 8900c3c4a36SLisandro Dalcin PetscInt v; 8914de306b1SToby Isaac PetscErrorCode ierr; 8924de306b1SToby Isaac 8934de306b1SToby Isaac PetscFunctionBegin; 894d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 895d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(is, IS_CLASSID, 3); 8960c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 8970c3c4a36SLisandro Dalcin if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);} 8984de306b1SToby Isaac if (is == label->points[v]) PetscFunctionReturn(0); 8994de306b1SToby Isaac ierr = DMLabelClearStratum(label, value);CHKERRQ(ierr); 9004de306b1SToby Isaac ierr = ISGetLocalSize(is, &(label->stratumSizes[v]));CHKERRQ(ierr); 9014de306b1SToby Isaac ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 9024de306b1SToby Isaac ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 9030c3c4a36SLisandro Dalcin label->points[v] = is; 9040c3c4a36SLisandro Dalcin label->validIS[v] = PETSC_TRUE; 9054de306b1SToby Isaac if (label->bt) { 9064de306b1SToby Isaac const PetscInt *points; 9074de306b1SToby Isaac PetscInt p; 9084de306b1SToby Isaac 9094de306b1SToby Isaac ierr = ISGetIndices(is,&points);CHKERRQ(ierr); 9104de306b1SToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) { 9114de306b1SToby Isaac const PetscInt point = points[p]; 9124de306b1SToby Isaac 9134de306b1SToby 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); 9144de306b1SToby Isaac ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 9154de306b1SToby Isaac } 9164de306b1SToby Isaac } 9174de306b1SToby Isaac PetscFunctionReturn(0); 9184de306b1SToby Isaac } 9194de306b1SToby Isaac 92084f0b6dfSMatthew G. Knepley /*@ 92184f0b6dfSMatthew G. Knepley DMLabelClearStratum - Remove a stratum 9224de306b1SToby Isaac 92384f0b6dfSMatthew G. Knepley Input Parameters: 92484f0b6dfSMatthew G. Knepley + label - the DMLabel 92584f0b6dfSMatthew G. Knepley - value - the stratum value 92684f0b6dfSMatthew G. Knepley 92784f0b6dfSMatthew G. Knepley Level: intermediate 92884f0b6dfSMatthew G. Knepley 92984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 93084f0b6dfSMatthew G. Knepley @*/ 931c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value) 932c58f1c22SToby Isaac { 933c58f1c22SToby Isaac PetscInt v; 934c58f1c22SToby Isaac PetscErrorCode ierr; 935c58f1c22SToby Isaac 936c58f1c22SToby Isaac PetscFunctionBegin; 937d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 9380c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 9390c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 9404de306b1SToby Isaac if (label->validIS[v]) { 9414de306b1SToby Isaac if (label->bt) { 942c58f1c22SToby Isaac PetscInt i; 943ad8374ffSToby Isaac const PetscInt *points; 944c58f1c22SToby Isaac 945ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 946c58f1c22SToby Isaac for (i = 0; i < label->stratumSizes[v]; ++i) { 947ad8374ffSToby Isaac const PetscInt point = points[i]; 948c58f1c22SToby Isaac 949c58f1c22SToby 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); 950c58f1c22SToby Isaac ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 951c58f1c22SToby Isaac } 952ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr); 953c58f1c22SToby Isaac } 954c58f1c22SToby Isaac label->stratumSizes[v] = 0; 9550c3c4a36SLisandro Dalcin ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr); 9560c3c4a36SLisandro Dalcin ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr); 9570c3c4a36SLisandro Dalcin ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr); 958c58f1c22SToby Isaac } else { 959e8f14785SLisandro Dalcin PetscHSetIClear(label->ht[v]); 960c58f1c22SToby Isaac } 961c58f1c22SToby Isaac PetscFunctionReturn(0); 962c58f1c22SToby Isaac } 963c58f1c22SToby Isaac 96484f0b6dfSMatthew G. Knepley /*@ 96584f0b6dfSMatthew G. Knepley DMLabelFilter - Remove all points outside of [start, end) 96684f0b6dfSMatthew G. Knepley 96784f0b6dfSMatthew G. Knepley Input Parameters: 96884f0b6dfSMatthew G. Knepley + label - the DMLabel 96984f0b6dfSMatthew G. Knepley . start - the first point 97084f0b6dfSMatthew G. Knepley - end - the last point 97184f0b6dfSMatthew G. Knepley 97284f0b6dfSMatthew G. Knepley Level: intermediate 97384f0b6dfSMatthew G. Knepley 97484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 97584f0b6dfSMatthew G. Knepley @*/ 976c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end) 977c58f1c22SToby Isaac { 978c58f1c22SToby Isaac PetscInt v; 979c58f1c22SToby Isaac PetscErrorCode ierr; 980c58f1c22SToby Isaac 981c58f1c22SToby Isaac PetscFunctionBegin; 982d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 9830c3c4a36SLisandro Dalcin ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr); 984c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 985c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 986c58f1c22SToby Isaac PetscInt off, q; 987ad8374ffSToby Isaac const PetscInt *points; 988033405d5SLisandro Dalcin PetscInt numPointsNew = 0, *pointsNew = NULL; 989c58f1c22SToby Isaac 990ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 991033405d5SLisandro Dalcin for (q = 0; q < label->stratumSizes[v]; ++q) 992033405d5SLisandro Dalcin if (points[q] >= start && points[q] < end) 993033405d5SLisandro Dalcin numPointsNew++; 994033405d5SLisandro Dalcin ierr = PetscMalloc1(numPointsNew, &pointsNew);CHKERRQ(ierr); 995c58f1c22SToby Isaac for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) { 996033405d5SLisandro Dalcin if (points[q] >= start && points[q] < end) 997033405d5SLisandro Dalcin pointsNew[off++] = points[q]; 998ad8374ffSToby Isaac } 999ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 1000033405d5SLisandro Dalcin 1001033405d5SLisandro Dalcin label->stratumSizes[v] = numPointsNew; 1002033405d5SLisandro Dalcin ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr); 1003033405d5SLisandro Dalcin ierr = ISCreateGeneral(PETSC_COMM_SELF,numPointsNew, pointsNew, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr); 1004033405d5SLisandro Dalcin ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr); 1005c58f1c22SToby Isaac } 1006c58f1c22SToby Isaac ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr); 1007c58f1c22SToby Isaac PetscFunctionReturn(0); 1008c58f1c22SToby Isaac } 1009c58f1c22SToby Isaac 101084f0b6dfSMatthew G. Knepley /*@ 101184f0b6dfSMatthew G. Knepley DMLabelPermute - Create a new label with permuted points 101284f0b6dfSMatthew G. Knepley 101384f0b6dfSMatthew G. Knepley Input Parameters: 101484f0b6dfSMatthew G. Knepley + label - the DMLabel 101584f0b6dfSMatthew G. Knepley - permutation - the point permutation 101684f0b6dfSMatthew G. Knepley 101784f0b6dfSMatthew G. Knepley Output Parameter: 101884f0b6dfSMatthew G. Knepley . labelnew - the new label containing the permuted points 101984f0b6dfSMatthew G. Knepley 102084f0b6dfSMatthew G. Knepley Level: intermediate 102184f0b6dfSMatthew G. Knepley 102284f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 102384f0b6dfSMatthew G. Knepley @*/ 1024c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew) 1025c58f1c22SToby Isaac { 1026c58f1c22SToby Isaac const PetscInt *perm; 1027c58f1c22SToby Isaac PetscInt numValues, numPoints, v, q; 1028c58f1c22SToby Isaac PetscErrorCode ierr; 1029c58f1c22SToby Isaac 1030c58f1c22SToby Isaac PetscFunctionBegin; 1031d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1032d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(permutation, IS_CLASSID, 2); 1033c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 1034c58f1c22SToby Isaac ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr); 1035c58f1c22SToby Isaac ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr); 1036c58f1c22SToby Isaac ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr); 1037c58f1c22SToby Isaac ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr); 1038c58f1c22SToby Isaac for (v = 0; v < numValues; ++v) { 1039c58f1c22SToby Isaac const PetscInt size = (*labelNew)->stratumSizes[v]; 1040ad8374ffSToby Isaac const PetscInt *points; 1041ad8374ffSToby Isaac PetscInt *pointsNew; 1042c58f1c22SToby Isaac 1043ad8374ffSToby Isaac ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 1044ad8374ffSToby Isaac ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr); 1045c58f1c22SToby Isaac for (q = 0; q < size; ++q) { 1046ad8374ffSToby Isaac const PetscInt point = points[q]; 1047c58f1c22SToby Isaac 1048c58f1c22SToby 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); 1049ad8374ffSToby Isaac pointsNew[q] = perm[point]; 1050c58f1c22SToby Isaac } 1051ad8374ffSToby Isaac ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 1052ad8374ffSToby Isaac ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr); 1053ad8374ffSToby Isaac ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr); 1054fa8e8ae5SToby Isaac if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) { 1055fa8e8ae5SToby Isaac ierr = ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));CHKERRQ(ierr); 1056fa8e8ae5SToby Isaac ierr = PetscFree(pointsNew);CHKERRQ(ierr); 1057fa8e8ae5SToby Isaac } else { 1058ad8374ffSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr); 1059fa8e8ae5SToby Isaac } 1060ad8374ffSToby Isaac ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr); 1061c58f1c22SToby Isaac } 1062c58f1c22SToby Isaac ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr); 1063c58f1c22SToby Isaac if (label->bt) { 1064c58f1c22SToby Isaac ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 1065c58f1c22SToby Isaac ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr); 1066c58f1c22SToby Isaac } 1067c58f1c22SToby Isaac PetscFunctionReturn(0); 1068c58f1c22SToby Isaac } 1069c58f1c22SToby Isaac 107026c55118SMichael Lange PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata) 107126c55118SMichael Lange { 107226c55118SMichael Lange MPI_Comm comm; 107326c55118SMichael Lange PetscInt s, l, nroots, nleaves, dof, offset, size; 107426c55118SMichael Lange PetscInt *remoteOffsets, *rootStrata, *rootIdx; 107526c55118SMichael Lange PetscSection rootSection; 107626c55118SMichael Lange PetscSF labelSF; 107726c55118SMichael Lange PetscErrorCode ierr; 107826c55118SMichael Lange 107926c55118SMichael Lange PetscFunctionBegin; 108026c55118SMichael Lange if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 108126c55118SMichael Lange ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 108226c55118SMichael Lange /* Build a section of stratum values per point, generate the according SF 108326c55118SMichael Lange and distribute point-wise stratum values to leaves. */ 108426c55118SMichael Lange ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 108526c55118SMichael Lange ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 108626c55118SMichael Lange ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr); 108726c55118SMichael Lange if (label) { 108826c55118SMichael Lange for (s = 0; s < label->numStrata; ++s) { 1089ad8374ffSToby Isaac const PetscInt *points; 1090ad8374ffSToby Isaac 1091ad8374ffSToby Isaac ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 109226c55118SMichael Lange for (l = 0; l < label->stratumSizes[s]; l++) { 1093ad8374ffSToby Isaac ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr); 1094ad8374ffSToby Isaac ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr); 109526c55118SMichael Lange } 1096ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 109726c55118SMichael Lange } 109826c55118SMichael Lange } 109926c55118SMichael Lange ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 110026c55118SMichael Lange /* Create a point-wise array of stratum values */ 110126c55118SMichael Lange ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr); 110226c55118SMichael Lange ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr); 110326c55118SMichael Lange ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr); 110426c55118SMichael Lange if (label) { 110526c55118SMichael Lange for (s = 0; s < label->numStrata; ++s) { 1106ad8374ffSToby Isaac const PetscInt *points; 1107ad8374ffSToby Isaac 1108ad8374ffSToby Isaac ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 110926c55118SMichael Lange for (l = 0; l < label->stratumSizes[s]; l++) { 1110ad8374ffSToby Isaac const PetscInt p = points[l]; 111126c55118SMichael Lange ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr); 111226c55118SMichael Lange rootStrata[offset+rootIdx[p]++] = label->stratumValues[s]; 111326c55118SMichael Lange } 1114ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 111526c55118SMichael Lange } 111626c55118SMichael Lange } 111726c55118SMichael Lange /* Build SF that maps label points to remote processes */ 111826c55118SMichael Lange ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr); 111926c55118SMichael Lange ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr); 112026c55118SMichael Lange ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr); 112126c55118SMichael Lange ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 112226c55118SMichael Lange /* Send the strata for each point over the derived SF */ 112326c55118SMichael Lange ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr); 112426c55118SMichael Lange ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr); 112526c55118SMichael Lange ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 112626c55118SMichael Lange ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 112726c55118SMichael Lange /* Clean up */ 112826c55118SMichael Lange ierr = PetscFree(rootStrata);CHKERRQ(ierr); 112926c55118SMichael Lange ierr = PetscFree(rootIdx);CHKERRQ(ierr); 113026c55118SMichael Lange ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 113126c55118SMichael Lange ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr); 113226c55118SMichael Lange PetscFunctionReturn(0); 113326c55118SMichael Lange } 113426c55118SMichael Lange 113584f0b6dfSMatthew G. Knepley /*@ 113684f0b6dfSMatthew G. Knepley DMLabelDistribute - Create a new label pushed forward over the PetscSF 113784f0b6dfSMatthew G. Knepley 113884f0b6dfSMatthew G. Knepley Input Parameters: 113984f0b6dfSMatthew G. Knepley + label - the DMLabel 114084f0b6dfSMatthew G. Knepley - sf - the map from old to new distribution 114184f0b6dfSMatthew G. Knepley 114284f0b6dfSMatthew G. Knepley Output Parameter: 114384f0b6dfSMatthew G. Knepley . labelnew - the new redistributed label 114484f0b6dfSMatthew G. Knepley 114584f0b6dfSMatthew G. Knepley Level: intermediate 114684f0b6dfSMatthew G. Knepley 114784f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 114884f0b6dfSMatthew G. Knepley @*/ 1149c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew) 1150c58f1c22SToby Isaac { 1151c58f1c22SToby Isaac MPI_Comm comm; 115226c55118SMichael Lange PetscSection leafSection; 115326c55118SMichael Lange PetscInt p, pStart, pEnd, s, size, dof, offset, stratum; 115426c55118SMichael Lange PetscInt *leafStrata, *strataIdx; 1155ad8374ffSToby Isaac PetscInt **points; 1156d67d17b1SMatthew G. Knepley const char *lname = NULL; 1157c58f1c22SToby Isaac char *name; 1158c58f1c22SToby Isaac PetscInt nameSize; 1159e8f14785SLisandro Dalcin PetscHSetI stratumHash; 1160c58f1c22SToby Isaac size_t len = 0; 116126c55118SMichael Lange PetscMPIInt rank; 1162c58f1c22SToby Isaac PetscErrorCode ierr; 1163c58f1c22SToby Isaac 1164c58f1c22SToby Isaac PetscFunctionBegin; 1165d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 1166*f018e600SMatthew G. Knepley if (label) { 1167*f018e600SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1168*f018e600SMatthew G. Knepley ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 1169*f018e600SMatthew G. Knepley } 1170c58f1c22SToby Isaac ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 1171c58f1c22SToby Isaac ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1172c58f1c22SToby Isaac /* Bcast name */ 1173d67d17b1SMatthew G. Knepley if (!rank) { 1174d67d17b1SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr); 1175d67d17b1SMatthew G. Knepley ierr = PetscStrlen(lname, &len);CHKERRQ(ierr); 1176d67d17b1SMatthew G. Knepley } 1177c58f1c22SToby Isaac nameSize = len; 1178c58f1c22SToby Isaac ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1179c58f1c22SToby Isaac ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1180d67d17b1SMatthew G. Knepley if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);} 1181c58f1c22SToby Isaac ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1182d67d17b1SMatthew G. Knepley ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr); 1183c58f1c22SToby Isaac ierr = PetscFree(name);CHKERRQ(ierr); 118477d236dfSMichael Lange /* Bcast defaultValue */ 118577d236dfSMichael Lange if (!rank) (*labelNew)->defaultValue = label->defaultValue; 118677d236dfSMichael Lange ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 118726c55118SMichael Lange /* Distribute stratum values over the SF and get the point mapping on the receiver */ 118826c55118SMichael Lange ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr); 11895cbdf6fcSMichael Lange /* Determine received stratum values and initialise new label*/ 1190e8f14785SLisandro Dalcin ierr = PetscHSetICreate(&stratumHash);CHKERRQ(ierr); 119126c55118SMichael Lange ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr); 1192e8f14785SLisandro Dalcin for (p = 0; p < size; ++p) {ierr = PetscHSetIAdd(stratumHash, leafStrata[p]);CHKERRQ(ierr);} 1193e8f14785SLisandro Dalcin ierr = PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata);CHKERRQ(ierr); 1194ad8374ffSToby Isaac ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr); 1195ad8374ffSToby Isaac for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE; 11965cbdf6fcSMichael Lange ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr); 11975cbdf6fcSMichael Lange /* Turn leafStrata into indices rather than stratum values */ 11985cbdf6fcSMichael Lange offset = 0; 1199e8f14785SLisandro Dalcin ierr = PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr); 1200a205f1fdSToby Isaac ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr); 12015cbdf6fcSMichael Lange for (p = 0; p < size; ++p) { 1202231b9e6fSMatthew G. Knepley for (s = 0; s < (*labelNew)->numStrata; ++s) { 1203231b9e6fSMatthew G. Knepley if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;} 12045cbdf6fcSMichael Lange } 12055cbdf6fcSMichael Lange } 1206c58f1c22SToby Isaac /* Rebuild the point strata on the receiver */ 1207c58f1c22SToby Isaac ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr); 1208c58f1c22SToby Isaac ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1209c58f1c22SToby Isaac for (p=pStart; p<pEnd; p++) { 1210c58f1c22SToby Isaac ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1211c58f1c22SToby Isaac ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1212c58f1c22SToby Isaac for (s=0; s<dof; s++) { 1213c58f1c22SToby Isaac (*labelNew)->stratumSizes[leafStrata[offset+s]]++; 1214c58f1c22SToby Isaac } 1215c58f1c22SToby Isaac } 1216c58f1c22SToby Isaac ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr); 1217c58f1c22SToby Isaac ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr); 1218ad8374ffSToby Isaac ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr); 1219c58f1c22SToby Isaac for (s = 0; s < (*labelNew)->numStrata; ++s) { 1220e8f14785SLisandro Dalcin ierr = PetscHSetICreate(&(*labelNew)->ht[s]);CHKERRQ(ierr); 1221ad8374ffSToby Isaac ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr); 1222c58f1c22SToby Isaac } 1223c58f1c22SToby Isaac /* Insert points into new strata */ 1224c58f1c22SToby Isaac ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr); 1225c58f1c22SToby Isaac ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1226c58f1c22SToby Isaac for (p=pStart; p<pEnd; p++) { 1227c58f1c22SToby Isaac ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1228c58f1c22SToby Isaac ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1229c58f1c22SToby Isaac for (s=0; s<dof; s++) { 1230c58f1c22SToby Isaac stratum = leafStrata[offset+s]; 1231ad8374ffSToby Isaac points[stratum][strataIdx[stratum]++] = p; 1232c58f1c22SToby Isaac } 1233c58f1c22SToby Isaac } 1234ad8374ffSToby Isaac for (s = 0; s < (*labelNew)->numStrata; s++) { 1235ad8374ffSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr); 1236ad8374ffSToby Isaac ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr); 1237ad8374ffSToby Isaac } 1238ad8374ffSToby Isaac ierr = PetscFree(points);CHKERRQ(ierr); 1239e8f14785SLisandro Dalcin ierr = PetscHSetIDestroy(&stratumHash);CHKERRQ(ierr); 1240c58f1c22SToby Isaac ierr = PetscFree(leafStrata);CHKERRQ(ierr); 1241c58f1c22SToby Isaac ierr = PetscFree(strataIdx);CHKERRQ(ierr); 1242c58f1c22SToby Isaac ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 1243c58f1c22SToby Isaac PetscFunctionReturn(0); 1244c58f1c22SToby Isaac } 1245c58f1c22SToby Isaac 12467937d9ceSMichael Lange /*@ 12477937d9ceSMichael Lange DMLabelGather - Gather all label values from leafs into roots 12487937d9ceSMichael Lange 12497937d9ceSMichael Lange Input Parameters: 12507937d9ceSMichael Lange + label - the DMLabel 125184f0b6dfSMatthew G. Knepley - sf - the Star Forest point communication map 12527937d9ceSMichael Lange 125384f0b6dfSMatthew G. Knepley Output Parameters: 125484f0b6dfSMatthew G. Knepley . labelNew - the new DMLabel with localised leaf values 12557937d9ceSMichael Lange 12567937d9ceSMichael Lange Level: developer 12577937d9ceSMichael Lange 12587937d9ceSMichael Lange Note: This is the inverse operation to DMLabelDistribute. 12597937d9ceSMichael Lange 12607937d9ceSMichael Lange .seealso: DMLabelDistribute() 12617937d9ceSMichael Lange @*/ 12627937d9ceSMichael Lange PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew) 12637937d9ceSMichael Lange { 12647937d9ceSMichael Lange MPI_Comm comm; 12657937d9ceSMichael Lange PetscSection rootSection; 12667937d9ceSMichael Lange PetscSF sfLabel; 12677937d9ceSMichael Lange PetscSFNode *rootPoints, *leafPoints; 12687937d9ceSMichael Lange PetscInt p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset; 12697937d9ceSMichael Lange const PetscInt *rootDegree, *ilocal; 12707937d9ceSMichael Lange PetscInt *rootStrata; 1271d67d17b1SMatthew G. Knepley const char *lname; 12727937d9ceSMichael Lange char *name; 12737937d9ceSMichael Lange PetscInt nameSize; 12747937d9ceSMichael Lange size_t len = 0; 12759852e123SBarry Smith PetscMPIInt rank, size; 12767937d9ceSMichael Lange PetscErrorCode ierr; 12777937d9ceSMichael Lange 12787937d9ceSMichael Lange PetscFunctionBegin; 1279d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1280d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 12817937d9ceSMichael Lange ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 12827937d9ceSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 12839852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 12847937d9ceSMichael Lange /* Bcast name */ 1285d67d17b1SMatthew G. Knepley if (!rank) { 1286d67d17b1SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr); 1287d67d17b1SMatthew G. Knepley ierr = PetscStrlen(lname, &len);CHKERRQ(ierr); 1288d67d17b1SMatthew G. Knepley } 12897937d9ceSMichael Lange nameSize = len; 12907937d9ceSMichael Lange ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 12917937d9ceSMichael Lange ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1292d67d17b1SMatthew G. Knepley if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);} 12937937d9ceSMichael Lange ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1294d67d17b1SMatthew G. Knepley ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr); 12957937d9ceSMichael Lange ierr = PetscFree(name);CHKERRQ(ierr); 12967937d9ceSMichael Lange /* Gather rank/index pairs of leaves into local roots to build 12977937d9ceSMichael Lange an inverse, multi-rooted SF. Note that this ignores local leaf 12987937d9ceSMichael Lange indexing due to the use of the multiSF in PetscSFGather. */ 12997937d9ceSMichael Lange ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr); 1300dc53bc9bSMatthew G. Knepley ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr); 1301dc53bc9bSMatthew G. Knepley for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1; 13027937d9ceSMichael Lange for (p = 0; p < nleaves; p++) { 13038212dd46SStefano Zampini PetscInt ilp = ilocal ? ilocal[p] : p; 13048212dd46SStefano Zampini 13058212dd46SStefano Zampini leafPoints[ilp].index = ilp; 13068212dd46SStefano Zampini leafPoints[ilp].rank = rank; 13077937d9ceSMichael Lange } 13087937d9ceSMichael Lange ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr); 13097937d9ceSMichael Lange ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr); 13107937d9ceSMichael Lange for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p]; 13117937d9ceSMichael Lange ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr); 13127937d9ceSMichael Lange ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 13137937d9ceSMichael Lange ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 13147937d9ceSMichael Lange ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr); 13157937d9ceSMichael Lange ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 13167937d9ceSMichael Lange /* Migrate label over inverted SF to pull stratum values at leaves into roots. */ 13177937d9ceSMichael Lange ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr); 13187937d9ceSMichael Lange /* Rebuild the point strata on the receiver */ 13197937d9ceSMichael Lange for (p = 0, idx = 0; p < nroots; p++) { 13207937d9ceSMichael Lange for (d = 0; d < rootDegree[p]; d++) { 13217937d9ceSMichael Lange ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr); 13227937d9ceSMichael Lange ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr); 13237937d9ceSMichael Lange for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);} 13247937d9ceSMichael Lange } 13257937d9ceSMichael Lange idx += rootDegree[p]; 13267937d9ceSMichael Lange } 132777e0c0e7SMichael Lange ierr = PetscFree(leafPoints);CHKERRQ(ierr); 132877e0c0e7SMichael Lange ierr = PetscFree(rootStrata);CHKERRQ(ierr); 132977e0c0e7SMichael Lange ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 133077e0c0e7SMichael Lange ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr); 13317937d9ceSMichael Lange PetscFunctionReturn(0); 13327937d9ceSMichael Lange } 13337937d9ceSMichael Lange 133484f0b6dfSMatthew G. Knepley /*@ 133584f0b6dfSMatthew G. Knepley DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label 133684f0b6dfSMatthew G. Knepley 133784f0b6dfSMatthew G. Knepley Input Parameter: 133884f0b6dfSMatthew G. Knepley . label - the DMLabel 133984f0b6dfSMatthew G. Knepley 134084f0b6dfSMatthew G. Knepley Output Parameters: 134184f0b6dfSMatthew G. Knepley + section - the section giving offsets for each stratum 134284f0b6dfSMatthew G. Knepley - is - An IS containing all the label points 134384f0b6dfSMatthew G. Knepley 134484f0b6dfSMatthew G. Knepley Level: developer 134584f0b6dfSMatthew G. Knepley 134684f0b6dfSMatthew G. Knepley .seealso: DMLabelDistribute() 134784f0b6dfSMatthew G. Knepley @*/ 1348c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is) 1349c58f1c22SToby Isaac { 1350c58f1c22SToby Isaac IS vIS; 1351c58f1c22SToby Isaac const PetscInt *values; 1352c58f1c22SToby Isaac PetscInt *points; 1353c58f1c22SToby Isaac PetscInt nV, vS = 0, vE = 0, v, N; 1354c58f1c22SToby Isaac PetscErrorCode ierr; 1355c58f1c22SToby Isaac 1356c58f1c22SToby Isaac PetscFunctionBegin; 1357d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1358c58f1c22SToby Isaac ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr); 1359c58f1c22SToby Isaac ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr); 1360c58f1c22SToby Isaac ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr); 1361c58f1c22SToby Isaac if (nV) {vS = values[0]; vE = values[0]+1;} 1362c58f1c22SToby Isaac for (v = 1; v < nV; ++v) { 1363c58f1c22SToby Isaac vS = PetscMin(vS, values[v]); 1364c58f1c22SToby Isaac vE = PetscMax(vE, values[v]+1); 1365c58f1c22SToby Isaac } 1366c58f1c22SToby Isaac ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr); 1367c58f1c22SToby Isaac ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr); 1368c58f1c22SToby Isaac for (v = 0; v < nV; ++v) { 1369c58f1c22SToby Isaac PetscInt n; 1370c58f1c22SToby Isaac 1371c58f1c22SToby Isaac ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr); 1372c58f1c22SToby Isaac ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr); 1373c58f1c22SToby Isaac } 1374c58f1c22SToby Isaac ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 1375c58f1c22SToby Isaac ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr); 1376c58f1c22SToby Isaac ierr = PetscMalloc1(N, &points);CHKERRQ(ierr); 1377c58f1c22SToby Isaac for (v = 0; v < nV; ++v) { 1378c58f1c22SToby Isaac IS is; 1379c58f1c22SToby Isaac const PetscInt *spoints; 1380c58f1c22SToby Isaac PetscInt dof, off, p; 1381c58f1c22SToby Isaac 1382c58f1c22SToby Isaac ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr); 1383c58f1c22SToby Isaac ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr); 1384c58f1c22SToby Isaac ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr); 1385c58f1c22SToby Isaac ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr); 1386c58f1c22SToby Isaac for (p = 0; p < dof; ++p) points[off+p] = spoints[p]; 1387c58f1c22SToby Isaac ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr); 1388c58f1c22SToby Isaac ierr = ISDestroy(&is);CHKERRQ(ierr); 1389c58f1c22SToby Isaac } 1390c58f1c22SToby Isaac ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr); 1391c58f1c22SToby Isaac ierr = ISDestroy(&vIS);CHKERRQ(ierr); 1392c58f1c22SToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr); 1393c58f1c22SToby Isaac PetscFunctionReturn(0); 1394c58f1c22SToby Isaac } 1395c58f1c22SToby Isaac 139684f0b6dfSMatthew G. Knepley /*@ 1397c58f1c22SToby Isaac PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using 1398c58f1c22SToby Isaac the local section and an SF describing the section point overlap. 1399c58f1c22SToby Isaac 1400c58f1c22SToby Isaac Input Parameters: 1401c58f1c22SToby Isaac + s - The PetscSection for the local field layout 1402c58f1c22SToby Isaac . sf - The SF describing parallel layout of the section points 1403c58f1c22SToby Isaac . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs 1404c58f1c22SToby Isaac . label - The label specifying the points 1405c58f1c22SToby Isaac - labelValue - The label stratum specifying the points 1406c58f1c22SToby Isaac 1407c58f1c22SToby Isaac Output Parameter: 1408c58f1c22SToby Isaac . gsection - The PetscSection for the global field layout 1409c58f1c22SToby Isaac 1410c58f1c22SToby Isaac Note: This gives negative sizes and offsets to points not owned by this process 1411c58f1c22SToby Isaac 1412c58f1c22SToby Isaac Level: developer 1413c58f1c22SToby Isaac 1414c58f1c22SToby Isaac .seealso: PetscSectionCreate() 1415c58f1c22SToby Isaac @*/ 1416c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection) 1417c58f1c22SToby Isaac { 1418c58f1c22SToby Isaac PetscInt *neg = NULL, *tmpOff = NULL; 1419c58f1c22SToby Isaac PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots; 1420c58f1c22SToby Isaac PetscErrorCode ierr; 1421c58f1c22SToby Isaac 1422c58f1c22SToby Isaac PetscFunctionBegin; 1423d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1424d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 1425d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4); 1426c58f1c22SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr); 1427c58f1c22SToby Isaac ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 1428c58f1c22SToby Isaac ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr); 1429c58f1c22SToby Isaac ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 1430c58f1c22SToby Isaac if (nroots >= 0) { 1431c58f1c22SToby Isaac if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart); 1432c58f1c22SToby Isaac ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr); 1433c58f1c22SToby Isaac if (nroots > pEnd-pStart) { 1434c58f1c22SToby Isaac ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr); 1435c58f1c22SToby Isaac } else { 1436c58f1c22SToby Isaac tmpOff = &(*gsection)->atlasDof[-pStart]; 1437c58f1c22SToby Isaac } 1438c58f1c22SToby Isaac } 1439c58f1c22SToby Isaac /* Mark ghost points with negative dof */ 1440c58f1c22SToby Isaac for (p = pStart; p < pEnd; ++p) { 1441c58f1c22SToby Isaac PetscInt value; 1442c58f1c22SToby Isaac 1443c58f1c22SToby Isaac ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr); 1444c58f1c22SToby Isaac if (value != labelValue) continue; 1445c58f1c22SToby Isaac ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr); 1446c58f1c22SToby Isaac ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr); 1447c58f1c22SToby Isaac ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 1448c58f1c22SToby Isaac if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);} 1449c58f1c22SToby Isaac if (neg) neg[p] = -(dof+1); 1450c58f1c22SToby Isaac } 1451c58f1c22SToby Isaac ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr); 1452c58f1c22SToby Isaac if (nroots >= 0) { 1453c58f1c22SToby Isaac ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1454c58f1c22SToby Isaac ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1455c58f1c22SToby Isaac if (nroots > pEnd-pStart) { 1456c58f1c22SToby Isaac for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];} 1457c58f1c22SToby Isaac } 1458c58f1c22SToby Isaac } 1459c58f1c22SToby Isaac /* Calculate new sizes, get proccess offset, and calculate point offsets */ 1460c58f1c22SToby Isaac for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1461c58f1c22SToby Isaac cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0; 1462c58f1c22SToby Isaac (*gsection)->atlasOff[p] = off; 1463c58f1c22SToby Isaac off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0; 1464c58f1c22SToby Isaac } 1465c58f1c22SToby Isaac ierr = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr); 1466c58f1c22SToby Isaac globalOff -= off; 1467c58f1c22SToby Isaac for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1468c58f1c22SToby Isaac (*gsection)->atlasOff[p] += globalOff; 1469c58f1c22SToby Isaac if (neg) neg[p] = -((*gsection)->atlasOff[p]+1); 1470c58f1c22SToby Isaac } 1471c58f1c22SToby Isaac /* Put in negative offsets for ghost points */ 1472c58f1c22SToby Isaac if (nroots >= 0) { 1473c58f1c22SToby Isaac ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1474c58f1c22SToby Isaac ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1475c58f1c22SToby Isaac if (nroots > pEnd-pStart) { 1476c58f1c22SToby Isaac for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];} 1477c58f1c22SToby Isaac } 1478c58f1c22SToby Isaac } 1479c58f1c22SToby Isaac if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);} 1480c58f1c22SToby Isaac ierr = PetscFree(neg);CHKERRQ(ierr); 1481c58f1c22SToby Isaac PetscFunctionReturn(0); 1482c58f1c22SToby Isaac } 1483c58f1c22SToby Isaac 14845fdea053SToby Isaac typedef struct _n_PetscSectionSym_Label 14855fdea053SToby Isaac { 14865fdea053SToby Isaac DMLabel label; 14875fdea053SToby Isaac PetscCopyMode *modes; 14885fdea053SToby Isaac PetscInt *sizes; 14895fdea053SToby Isaac const PetscInt ***perms; 14905fdea053SToby Isaac const PetscScalar ***rots; 14915fdea053SToby Isaac PetscInt (*minMaxOrients)[2]; 14925fdea053SToby Isaac PetscInt numStrata; /* numStrata is only increasing, functions as a state */ 14935fdea053SToby Isaac } PetscSectionSym_Label; 14945fdea053SToby Isaac 14955fdea053SToby Isaac static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym) 14965fdea053SToby Isaac { 14975fdea053SToby Isaac PetscInt i, j; 14985fdea053SToby Isaac PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data; 14995fdea053SToby Isaac PetscErrorCode ierr; 15005fdea053SToby Isaac 15015fdea053SToby Isaac PetscFunctionBegin; 15025fdea053SToby Isaac for (i = 0; i <= sl->numStrata; i++) { 15035fdea053SToby Isaac if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) { 15045fdea053SToby Isaac for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 15055fdea053SToby Isaac if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);} 15065fdea053SToby Isaac if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);} 15075fdea053SToby Isaac } 15085fdea053SToby Isaac if (sl->perms[i]) { 15095fdea053SToby Isaac const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]]; 15105fdea053SToby Isaac 15115fdea053SToby Isaac ierr = PetscFree(perms);CHKERRQ(ierr); 15125fdea053SToby Isaac } 15135fdea053SToby Isaac if (sl->rots[i]) { 15145fdea053SToby Isaac const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]]; 15155fdea053SToby Isaac 15165fdea053SToby Isaac ierr = PetscFree(rots);CHKERRQ(ierr); 15175fdea053SToby Isaac } 15185fdea053SToby Isaac } 15195fdea053SToby Isaac } 15205fdea053SToby Isaac ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr); 15215fdea053SToby Isaac ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr); 15225fdea053SToby Isaac sl->numStrata = 0; 15235fdea053SToby Isaac PetscFunctionReturn(0); 15245fdea053SToby Isaac } 15255fdea053SToby Isaac 15265fdea053SToby Isaac static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym) 15275fdea053SToby Isaac { 15285fdea053SToby Isaac PetscErrorCode ierr; 15295fdea053SToby Isaac 15305fdea053SToby Isaac PetscFunctionBegin; 15315fdea053SToby Isaac ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr); 15325fdea053SToby Isaac ierr = PetscFree(sym->data);CHKERRQ(ierr); 15335fdea053SToby Isaac PetscFunctionReturn(0); 15345fdea053SToby Isaac } 15355fdea053SToby Isaac 15365fdea053SToby Isaac static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer) 15375fdea053SToby Isaac { 15385fdea053SToby Isaac PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data; 15395fdea053SToby Isaac PetscBool isAscii; 15405fdea053SToby Isaac DMLabel label = sl->label; 1541d67d17b1SMatthew G. Knepley const char *name; 15425fdea053SToby Isaac PetscErrorCode ierr; 15435fdea053SToby Isaac 15445fdea053SToby Isaac PetscFunctionBegin; 15455fdea053SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr); 15465fdea053SToby Isaac if (isAscii) { 15475fdea053SToby Isaac PetscInt i, j, k; 15485fdea053SToby Isaac PetscViewerFormat format; 15495fdea053SToby Isaac 15505fdea053SToby Isaac ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 15515fdea053SToby Isaac if (label) { 15525fdea053SToby Isaac ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 15535fdea053SToby Isaac if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 15545fdea053SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 15555fdea053SToby Isaac ierr = DMLabelView(label, viewer);CHKERRQ(ierr); 15565fdea053SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 15575fdea053SToby Isaac } else { 1558d67d17b1SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr); 1559d67d17b1SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer," Label '%s'\n",name);CHKERRQ(ierr); 15605fdea053SToby Isaac } 15615fdea053SToby Isaac } else { 15625fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr); 15635fdea053SToby Isaac } 15645fdea053SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 15655fdea053SToby Isaac for (i = 0; i <= sl->numStrata; i++) { 15665fdea053SToby Isaac PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue; 15675fdea053SToby Isaac 15685fdea053SToby Isaac if (!(sl->perms[i] || sl->rots[i])) { 15695fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr); 15705fdea053SToby Isaac } else { 15715fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr); 15725fdea053SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 15735fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr); 15745fdea053SToby Isaac if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 15755fdea053SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 15765fdea053SToby Isaac for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 15775fdea053SToby Isaac if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) { 15785fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr); 15795fdea053SToby Isaac } else { 15805fdea053SToby Isaac PetscInt tab; 15815fdea053SToby Isaac 15825fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr); 15835fdea053SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 15845fdea053SToby Isaac ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr); 15855fdea053SToby Isaac if (sl->perms[i] && sl->perms[i][j]) { 15865fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr); 15875fdea053SToby Isaac ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr); 15885fdea053SToby Isaac for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);} 15895fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 15905fdea053SToby Isaac ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr); 15915fdea053SToby Isaac } 15925fdea053SToby Isaac if (sl->rots[i] && sl->rots[i][j]) { 15935fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Rotations: ");CHKERRQ(ierr); 15945fdea053SToby Isaac ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr); 15955fdea053SToby Isaac #if defined(PETSC_USE_COMPLEX) 15965fdea053SToby Isaac for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f+i*%+f",PetscRealPart(sl->rots[i][j][k]),PetscImaginaryPart(sl->rots[i][j][k]));CHKERRQ(ierr);} 15975fdea053SToby Isaac #else 15985fdea053SToby Isaac for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);} 15995fdea053SToby Isaac #endif 16005fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 16015fdea053SToby Isaac ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr); 16025fdea053SToby Isaac } 16035fdea053SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 16045fdea053SToby Isaac } 16055fdea053SToby Isaac } 16065fdea053SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 16075fdea053SToby Isaac } 16085fdea053SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 16095fdea053SToby Isaac } 16105fdea053SToby Isaac } 16115fdea053SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 16125fdea053SToby Isaac } 16135fdea053SToby Isaac PetscFunctionReturn(0); 16145fdea053SToby Isaac } 16155fdea053SToby Isaac 16165fdea053SToby Isaac /*@ 16175fdea053SToby Isaac PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries 16185fdea053SToby Isaac 16195fdea053SToby Isaac Logically collective on sym 16205fdea053SToby Isaac 16215fdea053SToby Isaac Input parameters: 16225fdea053SToby Isaac + sym - the section symmetries 16235fdea053SToby Isaac - label - the DMLabel describing the types of points 16245fdea053SToby Isaac 16255fdea053SToby Isaac Level: developer: 16265fdea053SToby Isaac 16275fdea053SToby Isaac .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms() 16285fdea053SToby Isaac @*/ 16295fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label) 16305fdea053SToby Isaac { 16315fdea053SToby Isaac PetscSectionSym_Label *sl; 16325fdea053SToby Isaac PetscErrorCode ierr; 16335fdea053SToby Isaac 16345fdea053SToby Isaac PetscFunctionBegin; 16355fdea053SToby Isaac PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1); 16365fdea053SToby Isaac sl = (PetscSectionSym_Label *) sym->data; 16375fdea053SToby Isaac if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);} 16385fdea053SToby Isaac if (label) { 16395fdea053SToby Isaac sl->label = label; 1640d67d17b1SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) label);CHKERRQ(ierr); 16415fdea053SToby Isaac ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr); 16421a834cf9SToby Isaac ierr = PetscMalloc5(sl->numStrata+1,&sl->modes,sl->numStrata+1,&sl->sizes,sl->numStrata+1,&sl->perms,sl->numStrata+1,&sl->rots,sl->numStrata+1,&sl->minMaxOrients);CHKERRQ(ierr); 16431a834cf9SToby Isaac ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr); 16441a834cf9SToby Isaac ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr); 16451a834cf9SToby Isaac ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr); 16461a834cf9SToby Isaac ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr); 16471a834cf9SToby Isaac ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr); 16485fdea053SToby Isaac } 16495fdea053SToby Isaac PetscFunctionReturn(0); 16505fdea053SToby Isaac } 16515fdea053SToby Isaac 16525fdea053SToby Isaac /*@C 16535fdea053SToby Isaac PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum 16545fdea053SToby Isaac 16555fdea053SToby Isaac Logically collective on PetscSectionSym 16565fdea053SToby Isaac 16575fdea053SToby Isaac InputParameters: 16585fdea053SToby Isaac + sys - the section symmetries 16595fdea053SToby Isaac . stratum - the stratum value in the label that we are assigning symmetries for 16605fdea053SToby Isaac . size - the number of dofs for points in the stratum of the label 16615fdea053SToby Isaac . minOrient - the smallest orientation for a point in this stratum 16625fdea053SToby Isaac . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient)) 16635fdea053SToby Isaac . mode - how sym should copy the perms and rots arrays 16645fdea053SToby Isaac . perms - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation. A NULL permutation is the identity 16655fdea053SToby Isaac + rots - NULL if there are no rotations, or (maxOrient - minOrient) sets of rotations, one for each orientation. A NULL set of orientations is the identity 16665fdea053SToby Isaac 16675fdea053SToby Isaac Level: developer 16685fdea053SToby Isaac 16695fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel() 16705fdea053SToby Isaac @*/ 16715fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots) 16725fdea053SToby Isaac { 16735fdea053SToby Isaac PetscSectionSym_Label *sl; 1674d67d17b1SMatthew G. Knepley const char *name; 1675d67d17b1SMatthew G. Knepley PetscInt i, j, k; 16765fdea053SToby Isaac PetscErrorCode ierr; 16775fdea053SToby Isaac 16785fdea053SToby Isaac PetscFunctionBegin; 16795fdea053SToby Isaac PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1); 16805fdea053SToby Isaac sl = (PetscSectionSym_Label *) sym->data; 16815fdea053SToby Isaac if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet"); 16825fdea053SToby Isaac for (i = 0; i <= sl->numStrata; i++) { 16835fdea053SToby Isaac PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue; 16845fdea053SToby Isaac 16855fdea053SToby Isaac if (stratum == value) break; 16865fdea053SToby Isaac } 1687d67d17b1SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr); 1688d67d17b1SMatthew G. Knepley if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,name); 16895fdea053SToby Isaac sl->sizes[i] = size; 16905fdea053SToby Isaac sl->modes[i] = mode; 16915fdea053SToby Isaac sl->minMaxOrients[i][0] = minOrient; 16925fdea053SToby Isaac sl->minMaxOrients[i][1] = maxOrient; 16935fdea053SToby Isaac if (mode == PETSC_COPY_VALUES) { 16945fdea053SToby Isaac if (perms) { 16955fdea053SToby Isaac PetscInt **ownPerms; 16965fdea053SToby Isaac 16975fdea053SToby Isaac ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr); 16985fdea053SToby Isaac for (j = 0; j < maxOrient-minOrient; j++) { 16995fdea053SToby Isaac if (perms[j]) { 17005fdea053SToby Isaac ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr); 17015fdea053SToby Isaac for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];} 17025fdea053SToby Isaac } 17035fdea053SToby Isaac } 17045fdea053SToby Isaac sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient]; 17055fdea053SToby Isaac } 17065fdea053SToby Isaac if (rots) { 17075fdea053SToby Isaac PetscScalar **ownRots; 17085fdea053SToby Isaac 17095fdea053SToby Isaac ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr); 17105fdea053SToby Isaac for (j = 0; j < maxOrient-minOrient; j++) { 17115fdea053SToby Isaac if (rots[j]) { 17125fdea053SToby Isaac ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr); 17135fdea053SToby Isaac for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];} 17145fdea053SToby Isaac } 17155fdea053SToby Isaac } 17165fdea053SToby Isaac sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient]; 17175fdea053SToby Isaac } 17185fdea053SToby Isaac } else { 17195fdea053SToby Isaac sl->perms[i] = perms ? &perms[-minOrient] : NULL; 17205fdea053SToby Isaac sl->rots[i] = rots ? &rots[-minOrient] : NULL; 17215fdea053SToby Isaac } 17225fdea053SToby Isaac PetscFunctionReturn(0); 17235fdea053SToby Isaac } 17245fdea053SToby Isaac 17255fdea053SToby Isaac static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots) 17265fdea053SToby Isaac { 17275fdea053SToby Isaac PetscInt i, j, numStrata; 17285fdea053SToby Isaac PetscSectionSym_Label *sl; 17295fdea053SToby Isaac DMLabel label; 17305fdea053SToby Isaac PetscErrorCode ierr; 17315fdea053SToby Isaac 17325fdea053SToby Isaac PetscFunctionBegin; 17335fdea053SToby Isaac sl = (PetscSectionSym_Label *) sym->data; 17345fdea053SToby Isaac numStrata = sl->numStrata; 17355fdea053SToby Isaac label = sl->label; 17365fdea053SToby Isaac for (i = 0; i < numPoints; i++) { 17375fdea053SToby Isaac PetscInt point = points[2*i]; 17385fdea053SToby Isaac PetscInt ornt = points[2*i+1]; 17395fdea053SToby Isaac 17405fdea053SToby Isaac for (j = 0; j < numStrata; j++) { 17415fdea053SToby Isaac if (label->validIS[j]) { 17425fdea053SToby Isaac PetscInt k; 17435fdea053SToby Isaac 17445fdea053SToby Isaac ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr); 17455fdea053SToby Isaac if (k >= 0) break; 17465fdea053SToby Isaac } else { 17475fdea053SToby Isaac PetscBool has; 17485fdea053SToby Isaac 1749e8f14785SLisandro Dalcin PetscHSetIHas(label->ht[j], point, &has); 17505fdea053SToby Isaac if (has) break; 17515fdea053SToby Isaac } 17525fdea053SToby Isaac } 17535fdea053SToby Isaac if ((sl->minMaxOrients[j][1] > sl->minMaxOrients[j][0]) && (ornt < sl->minMaxOrients[j][0] || ornt >= sl->minMaxOrients[j][1])) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D orientation %D not in range [%D, %D) for stratum %D",point,ornt,sl->minMaxOrients[j][0],sl->minMaxOrients[j][1],j < numStrata ? label->stratumValues[j] : label->defaultValue); 17545fdea053SToby Isaac if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;} 17555fdea053SToby Isaac if (rots) {rots[i] = sl->rots[j] ? sl->rots[j][ornt] : NULL;} 17565fdea053SToby Isaac } 17575fdea053SToby Isaac PetscFunctionReturn(0); 17585fdea053SToby Isaac } 17595fdea053SToby Isaac 17605fdea053SToby Isaac PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym) 17615fdea053SToby Isaac { 17625fdea053SToby Isaac PetscSectionSym_Label *sl; 17635fdea053SToby Isaac PetscErrorCode ierr; 17645fdea053SToby Isaac 17655fdea053SToby Isaac PetscFunctionBegin; 17665fdea053SToby Isaac ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr); 17675fdea053SToby Isaac sym->ops->getpoints = PetscSectionSymGetPoints_Label; 17685fdea053SToby Isaac sym->ops->view = PetscSectionSymView_Label; 17695fdea053SToby Isaac sym->ops->destroy = PetscSectionSymDestroy_Label; 17705fdea053SToby Isaac sym->data = (void *) sl; 17715fdea053SToby Isaac PetscFunctionReturn(0); 17725fdea053SToby Isaac } 17735fdea053SToby Isaac 17745fdea053SToby Isaac /*@ 17755fdea053SToby Isaac PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label 17765fdea053SToby Isaac 17775fdea053SToby Isaac Collective 17785fdea053SToby Isaac 17795fdea053SToby Isaac Input Parameters: 17805fdea053SToby Isaac + comm - the MPI communicator for the new symmetry 17815fdea053SToby Isaac - label - the label defining the strata 17825fdea053SToby Isaac 17835fdea053SToby Isaac Output Parameters: 17845fdea053SToby Isaac . sym - the section symmetries 17855fdea053SToby Isaac 17865fdea053SToby Isaac Level: developer 17875fdea053SToby Isaac 17885fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms() 17895fdea053SToby Isaac @*/ 17905fdea053SToby Isaac PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym) 17915fdea053SToby Isaac { 17925fdea053SToby Isaac PetscErrorCode ierr; 17935fdea053SToby Isaac 17945fdea053SToby Isaac PetscFunctionBegin; 17955fdea053SToby Isaac ierr = DMInitializePackage();CHKERRQ(ierr); 17965fdea053SToby Isaac ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr); 17975fdea053SToby Isaac ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr); 17985fdea053SToby Isaac ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr); 17995fdea053SToby Isaac PetscFunctionReturn(0); 18005fdea053SToby Isaac } 1801