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 9c58f1c22SToby Isaac Input parameter: 10c58f1c22SToby Isaac . name - The label name 11c58f1c22SToby Isaac 12c58f1c22SToby Isaac Output parameter: 13c58f1c22SToby Isaac . label - The DMLabel 14c58f1c22SToby Isaac 15c58f1c22SToby Isaac Level: beginner 16c58f1c22SToby Isaac 17c58f1c22SToby Isaac .seealso: DMLabelDestroy() 18c58f1c22SToby Isaac @*/ 19c58f1c22SToby Isaac PetscErrorCode DMLabelCreate(const char name[], DMLabel *label) 20c58f1c22SToby Isaac { 21c58f1c22SToby Isaac PetscErrorCode ierr; 22c58f1c22SToby Isaac 23c58f1c22SToby Isaac PetscFunctionBegin; 24c58f1c22SToby Isaac ierr = PetscNew(label);CHKERRQ(ierr); 25c58f1c22SToby Isaac ierr = PetscStrallocpy(name, &(*label)->name);CHKERRQ(ierr); 26c58f1c22SToby Isaac 27c58f1c22SToby Isaac (*label)->refct = 1; 28c58f1c22SToby Isaac (*label)->state = -1; 29c58f1c22SToby Isaac (*label)->numStrata = 0; 305aa44df4SToby Isaac (*label)->defaultValue = -1; 31c58f1c22SToby Isaac (*label)->stratumValues = NULL; 32ad8374ffSToby Isaac (*label)->validIS = NULL; 33c58f1c22SToby Isaac (*label)->stratumSizes = NULL; 34c58f1c22SToby Isaac (*label)->points = NULL; 35c58f1c22SToby Isaac (*label)->ht = NULL; 36c58f1c22SToby Isaac (*label)->pStart = -1; 37c58f1c22SToby Isaac (*label)->pEnd = -1; 38c58f1c22SToby Isaac (*label)->bt = NULL; 39c58f1c22SToby Isaac PetscFunctionReturn(0); 40c58f1c22SToby Isaac } 41c58f1c22SToby Isaac 42c58f1c22SToby Isaac /* 43c58f1c22SToby Isaac DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format 44c58f1c22SToby Isaac 45c58f1c22SToby Isaac Input parameter: 46c58f1c22SToby Isaac + label - The DMLabel 47c58f1c22SToby Isaac - v - The stratum value 48c58f1c22SToby Isaac 49c58f1c22SToby Isaac Output parameter: 50c58f1c22SToby Isaac . label - The DMLabel with stratum in sorted list format 51c58f1c22SToby Isaac 52c58f1c22SToby Isaac Level: developer 53c58f1c22SToby Isaac 54c58f1c22SToby Isaac .seealso: DMLabelCreate() 55c58f1c22SToby Isaac */ 56c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v) 57c58f1c22SToby Isaac { 580c3c4a36SLisandro Dalcin PetscInt off = 0; 59ad8374ffSToby Isaac PetscInt *pointArray; 60c58f1c22SToby Isaac PetscErrorCode ierr; 61c58f1c22SToby Isaac 62ad8374ffSToby Isaac if (label->validIS[v]) return 0; 63c58f1c22SToby Isaac PetscFunctionBegin; 640c3c4a36SLisandro 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); 65*e8f14785SLisandro Dalcin ierr = PetscHSetIGetSize(label->ht[v], &label->stratumSizes[v]);CHKERRQ(ierr); 66c58f1c22SToby Isaac 67ad8374ffSToby Isaac ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr); 68*e8f14785SLisandro Dalcin ierr = PetscHSetIGetElems(label->ht[v], &off, pointArray);CHKERRQ(ierr); 69c58f1c22SToby 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]); 70*e8f14785SLisandro Dalcin PetscHSetIClear(label->ht[v]); 71ad8374ffSToby Isaac ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr); 72c58f1c22SToby Isaac if (label->bt) { 73c58f1c22SToby Isaac PetscInt p; 74c58f1c22SToby Isaac 75c58f1c22SToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) { 76ad8374ffSToby Isaac const PetscInt point = pointArray[p]; 77c58f1c22SToby Isaac 78c58f1c22SToby 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); 79c58f1c22SToby Isaac ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 80c58f1c22SToby Isaac } 81c58f1c22SToby Isaac } 82ad8374ffSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,label->stratumSizes[v],pointArray,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr); 83ad8374ffSToby Isaac ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr); 84ad8374ffSToby Isaac label->validIS[v] = PETSC_TRUE; 85c58f1c22SToby Isaac ++label->state; 86c58f1c22SToby Isaac PetscFunctionReturn(0); 87c58f1c22SToby Isaac } 88c58f1c22SToby Isaac 89c58f1c22SToby Isaac /* 90c58f1c22SToby Isaac DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format 91c58f1c22SToby Isaac 92c58f1c22SToby Isaac Input parameter: 93c58f1c22SToby Isaac . label - The DMLabel 94c58f1c22SToby Isaac 95c58f1c22SToby Isaac Output parameter: 96c58f1c22SToby Isaac . label - The DMLabel with all strata in sorted list format 97c58f1c22SToby Isaac 98c58f1c22SToby Isaac Level: developer 99c58f1c22SToby Isaac 100c58f1c22SToby Isaac .seealso: DMLabelCreate() 101c58f1c22SToby Isaac */ 102c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label) 103c58f1c22SToby Isaac { 104c58f1c22SToby Isaac PetscInt v; 105c58f1c22SToby Isaac PetscErrorCode ierr; 106c58f1c22SToby Isaac 107c58f1c22SToby Isaac PetscFunctionBegin; 108c58f1c22SToby Isaac for (v = 0; v < label->numStrata; v++){ 109c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 110c58f1c22SToby Isaac } 111c58f1c22SToby Isaac PetscFunctionReturn(0); 112c58f1c22SToby Isaac } 113c58f1c22SToby Isaac 114c58f1c22SToby Isaac /* 115c58f1c22SToby Isaac DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format 116c58f1c22SToby Isaac 117c58f1c22SToby Isaac Input parameter: 118c58f1c22SToby Isaac + label - The DMLabel 119c58f1c22SToby Isaac - v - The stratum value 120c58f1c22SToby Isaac 121c58f1c22SToby Isaac Output parameter: 122c58f1c22SToby Isaac . label - The DMLabel with stratum in hash format 123c58f1c22SToby Isaac 124c58f1c22SToby Isaac Level: developer 125c58f1c22SToby Isaac 126c58f1c22SToby Isaac .seealso: DMLabelCreate() 127c58f1c22SToby Isaac */ 128c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v) 129c58f1c22SToby Isaac { 130c58f1c22SToby Isaac PetscInt p; 131ad8374ffSToby Isaac const PetscInt *points; 132c58f1c22SToby Isaac PetscErrorCode ierr; 133c58f1c22SToby Isaac 1340c3c4a36SLisandro Dalcin if (!label->validIS[v]) return 0; 135c58f1c22SToby Isaac PetscFunctionBegin; 1360c3c4a36SLisandro 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); 137ad8374ffSToby Isaac if (label->points[v]) { 138ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 139*e8f14785SLisandro Dalcin for (p = 0; p < label->stratumSizes[v]; ++p) { 140*e8f14785SLisandro Dalcin ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr); 141*e8f14785SLisandro Dalcin } 142ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 143ad8374ffSToby Isaac ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 144ad8374ffSToby Isaac } 145ad8374ffSToby Isaac label->validIS[v] = PETSC_FALSE; 146c58f1c22SToby Isaac PetscFunctionReturn(0); 147c58f1c22SToby Isaac } 148c58f1c22SToby Isaac 149c58f1c22SToby Isaac PetscErrorCode DMLabelGetState(DMLabel label, PetscObjectState *state) 150c58f1c22SToby Isaac { 151c58f1c22SToby Isaac PetscFunctionBegin; 152c58f1c22SToby Isaac PetscValidPointer(state, 2); 153c58f1c22SToby Isaac *state = label->state; 154c58f1c22SToby Isaac PetscFunctionReturn(0); 155c58f1c22SToby Isaac } 156c58f1c22SToby Isaac 1570c3c4a36SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index) 1580c3c4a36SLisandro Dalcin { 1590c3c4a36SLisandro Dalcin PetscInt v; 1600c3c4a36SLisandro Dalcin PetscFunctionBegin; 1610c3c4a36SLisandro Dalcin for (*index = -1, v = 0; v < label->numStrata; ++v) 1620c3c4a36SLisandro Dalcin if (label->stratumValues[v] == value) { *index = v; break; } 1630c3c4a36SLisandro Dalcin PetscFunctionReturn(0); 1640c3c4a36SLisandro Dalcin } 1650c3c4a36SLisandro Dalcin 1660c3c4a36SLisandro Dalcin static PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index) 167c58f1c22SToby Isaac { 168ad8374ffSToby Isaac PetscInt v, *tmpV, *tmpS; 169ad8374ffSToby Isaac IS *tmpP; 170*e8f14785SLisandro Dalcin PetscHSetI *tmpH; 171c58f1c22SToby Isaac PetscBool *tmpB; 172c58f1c22SToby Isaac PetscErrorCode ierr; 173c58f1c22SToby Isaac 174c58f1c22SToby Isaac PetscFunctionBegin; 175c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpV);CHKERRQ(ierr); 176c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpS);CHKERRQ(ierr); 177c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpH);CHKERRQ(ierr); 178c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpP);CHKERRQ(ierr); 179c58f1c22SToby Isaac ierr = PetscMalloc1((label->numStrata+1), &tmpB);CHKERRQ(ierr); 180c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 181c58f1c22SToby Isaac tmpV[v] = label->stratumValues[v]; 182c58f1c22SToby Isaac tmpS[v] = label->stratumSizes[v]; 183c58f1c22SToby Isaac tmpH[v] = label->ht[v]; 184c58f1c22SToby Isaac tmpP[v] = label->points[v]; 185ad8374ffSToby Isaac tmpB[v] = label->validIS[v]; 186c58f1c22SToby Isaac } 187c58f1c22SToby Isaac tmpV[v] = value; 188c58f1c22SToby Isaac tmpS[v] = 0; 189*e8f14785SLisandro Dalcin ierr = PetscHSetICreate(&tmpH[v]);CHKERRQ(ierr); 190ecb9a371SToby Isaac ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&tmpP[v]);CHKERRQ(ierr); 191c58f1c22SToby Isaac tmpB[v] = PETSC_TRUE; 1920c3c4a36SLisandro Dalcin 193c58f1c22SToby Isaac ++label->numStrata; 194c58f1c22SToby Isaac ierr = PetscFree(label->stratumValues);CHKERRQ(ierr); 195c58f1c22SToby Isaac ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr); 196c58f1c22SToby Isaac ierr = PetscFree(label->ht);CHKERRQ(ierr); 197c58f1c22SToby Isaac ierr = PetscFree(label->points);CHKERRQ(ierr); 198ad8374ffSToby Isaac ierr = PetscFree(label->validIS);CHKERRQ(ierr); 199c58f1c22SToby Isaac label->stratumValues = tmpV; 200c58f1c22SToby Isaac label->stratumSizes = tmpS; 201c58f1c22SToby Isaac label->ht = tmpH; 202c58f1c22SToby Isaac label->points = tmpP; 203ad8374ffSToby Isaac label->validIS = tmpB; 204c58f1c22SToby Isaac 2050c3c4a36SLisandro Dalcin *index = v; 2060c3c4a36SLisandro Dalcin PetscFunctionReturn(0); 2070c3c4a36SLisandro Dalcin } 2080c3c4a36SLisandro Dalcin 2090c3c4a36SLisandro Dalcin PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value) 2100c3c4a36SLisandro Dalcin { 2110c3c4a36SLisandro Dalcin PetscInt v; 2120c3c4a36SLisandro Dalcin PetscErrorCode ierr; 2130c3c4a36SLisandro Dalcin 2140c3c4a36SLisandro Dalcin PetscFunctionBegin; 2150c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 2160c3c4a36SLisandro Dalcin if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);} 217c58f1c22SToby Isaac PetscFunctionReturn(0); 218c58f1c22SToby Isaac } 219c58f1c22SToby Isaac 220c58f1c22SToby Isaac /*@C 221c58f1c22SToby Isaac DMLabelGetName - Return the name of a DMLabel object 222c58f1c22SToby Isaac 223c58f1c22SToby Isaac Input parameter: 224c58f1c22SToby Isaac . label - The DMLabel 225c58f1c22SToby Isaac 226c58f1c22SToby Isaac Output parameter: 227c58f1c22SToby Isaac . name - The label name 228c58f1c22SToby Isaac 229c58f1c22SToby Isaac Level: beginner 230c58f1c22SToby Isaac 231c58f1c22SToby Isaac .seealso: DMLabelCreate() 232c58f1c22SToby Isaac @*/ 233c58f1c22SToby Isaac PetscErrorCode DMLabelGetName(DMLabel label, const char **name) 234c58f1c22SToby Isaac { 235c58f1c22SToby Isaac PetscFunctionBegin; 236c58f1c22SToby Isaac PetscValidPointer(name, 2); 237c58f1c22SToby Isaac *name = label->name; 238c58f1c22SToby Isaac PetscFunctionReturn(0); 239c58f1c22SToby Isaac } 240c58f1c22SToby Isaac 241c58f1c22SToby Isaac static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer) 242c58f1c22SToby Isaac { 243c58f1c22SToby Isaac PetscInt v; 244c58f1c22SToby Isaac PetscMPIInt rank; 245c58f1c22SToby Isaac PetscErrorCode ierr; 246c58f1c22SToby Isaac 247c58f1c22SToby Isaac PetscFunctionBegin; 248c58f1c22SToby Isaac ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr); 249c58f1c22SToby Isaac ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 250c58f1c22SToby Isaac if (label) { 251c58f1c22SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", label->name);CHKERRQ(ierr); 252c58f1c22SToby Isaac if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, " Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);} 253c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 254c58f1c22SToby Isaac const PetscInt value = label->stratumValues[v]; 255ad8374ffSToby Isaac const PetscInt *points; 256c58f1c22SToby Isaac PetscInt p; 257c58f1c22SToby Isaac 258ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr); 259c58f1c22SToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) { 260ad8374ffSToby Isaac ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr); 261c58f1c22SToby Isaac } 262ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 263c58f1c22SToby Isaac } 264c58f1c22SToby Isaac } 265c58f1c22SToby Isaac ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 266c58f1c22SToby Isaac ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 267c58f1c22SToby Isaac PetscFunctionReturn(0); 268c58f1c22SToby Isaac } 269c58f1c22SToby Isaac 270c58f1c22SToby Isaac /*@C 271c58f1c22SToby Isaac DMLabelView - View the label 272c58f1c22SToby Isaac 273c58f1c22SToby Isaac Input Parameters: 274c58f1c22SToby Isaac + label - The DMLabel 275c58f1c22SToby Isaac - viewer - The PetscViewer 276c58f1c22SToby Isaac 277c58f1c22SToby Isaac Level: intermediate 278c58f1c22SToby Isaac 279c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy() 280c58f1c22SToby Isaac @*/ 281c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer) 282c58f1c22SToby Isaac { 283c58f1c22SToby Isaac PetscBool iascii; 284c58f1c22SToby Isaac PetscErrorCode ierr; 285c58f1c22SToby Isaac 286c58f1c22SToby Isaac PetscFunctionBegin; 287c58f1c22SToby Isaac PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 288c58f1c22SToby Isaac if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 289c58f1c22SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 290c58f1c22SToby Isaac if (iascii) { 291c58f1c22SToby Isaac ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr); 292c58f1c22SToby Isaac } 293c58f1c22SToby Isaac PetscFunctionReturn(0); 294c58f1c22SToby Isaac } 295c58f1c22SToby Isaac 29684f0b6dfSMatthew G. Knepley /*@ 29784f0b6dfSMatthew G. Knepley DMLabelDestroy - Destroys a DMLabel 29884f0b6dfSMatthew G. Knepley 29984f0b6dfSMatthew G. Knepley Input Parameter: 30084f0b6dfSMatthew G. Knepley . label - The DMLabel 30184f0b6dfSMatthew G. Knepley 30284f0b6dfSMatthew G. Knepley Level: beginner 30384f0b6dfSMatthew G. Knepley 30484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate() 30584f0b6dfSMatthew G. Knepley @*/ 306c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label) 307c58f1c22SToby Isaac { 308c58f1c22SToby Isaac PetscInt v; 309c58f1c22SToby Isaac PetscErrorCode ierr; 310c58f1c22SToby Isaac 311c58f1c22SToby Isaac PetscFunctionBegin; 312c58f1c22SToby Isaac if (!(*label)) PetscFunctionReturn(0); 313c58f1c22SToby Isaac if (--(*label)->refct > 0) PetscFunctionReturn(0); 314c58f1c22SToby Isaac ierr = PetscFree((*label)->name);CHKERRQ(ierr); 315c58f1c22SToby Isaac ierr = PetscFree((*label)->stratumValues);CHKERRQ(ierr); 316c58f1c22SToby Isaac ierr = PetscFree((*label)->stratumSizes);CHKERRQ(ierr); 3170c3c4a36SLisandro Dalcin for (v = 0; v < (*label)->numStrata; ++v) { 318*e8f14785SLisandro Dalcin ierr = PetscHSetIDestroy(&(*label)->ht[v]);CHKERRQ(ierr); 319*e8f14785SLisandro Dalcin ierr = ISDestroy(&(*label)->points[v]);CHKERRQ(ierr); 3200c3c4a36SLisandro Dalcin } 3210c3c4a36SLisandro Dalcin ierr = PetscFree((*label)->ht);CHKERRQ(ierr); 322c58f1c22SToby Isaac ierr = PetscFree((*label)->points);CHKERRQ(ierr); 323ad8374ffSToby Isaac ierr = PetscFree((*label)->validIS);CHKERRQ(ierr); 324c58f1c22SToby Isaac ierr = PetscBTDestroy(&(*label)->bt);CHKERRQ(ierr); 325c58f1c22SToby Isaac ierr = PetscFree(*label);CHKERRQ(ierr); 326c58f1c22SToby Isaac PetscFunctionReturn(0); 327c58f1c22SToby Isaac } 328c58f1c22SToby Isaac 32984f0b6dfSMatthew G. Knepley /*@ 33084f0b6dfSMatthew G. Knepley DMLabelDuplicate - Duplicates a DMLabel 33184f0b6dfSMatthew G. Knepley 33284f0b6dfSMatthew G. Knepley Input Parameter: 33384f0b6dfSMatthew G. Knepley . label - The DMLabel 33484f0b6dfSMatthew G. Knepley 33584f0b6dfSMatthew G. Knepley Output Parameter: 33684f0b6dfSMatthew G. Knepley . labelnew - location to put new vector 33784f0b6dfSMatthew G. Knepley 33884f0b6dfSMatthew G. Knepley Level: intermediate 33984f0b6dfSMatthew G. Knepley 34084f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelDestroy() 34184f0b6dfSMatthew G. Knepley @*/ 342c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew) 343c58f1c22SToby Isaac { 344ad8374ffSToby Isaac PetscInt v; 345c58f1c22SToby Isaac PetscErrorCode ierr; 346c58f1c22SToby Isaac 347c58f1c22SToby Isaac PetscFunctionBegin; 348c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 349c58f1c22SToby Isaac ierr = PetscNew(labelnew);CHKERRQ(ierr); 350c58f1c22SToby Isaac ierr = PetscStrallocpy(label->name, &(*labelnew)->name);CHKERRQ(ierr); 351c58f1c22SToby Isaac 352c58f1c22SToby Isaac (*labelnew)->refct = 1; 353c58f1c22SToby Isaac (*labelnew)->numStrata = label->numStrata; 3545aa44df4SToby Isaac (*labelnew)->defaultValue = label->defaultValue; 355c58f1c22SToby Isaac if (label->numStrata) { 356c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr); 357c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr); 358c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr); 359c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr); 360ad8374ffSToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr); 361c58f1c22SToby Isaac /* Could eliminate unused space here */ 362c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 363*e8f14785SLisandro Dalcin ierr = PetscHSetICreate(&(*labelnew)->ht[v]);CHKERRQ(ierr); 364ad8374ffSToby Isaac (*labelnew)->validIS[v] = PETSC_TRUE; 365c58f1c22SToby Isaac (*labelnew)->stratumValues[v] = label->stratumValues[v]; 366c58f1c22SToby Isaac (*labelnew)->stratumSizes[v] = label->stratumSizes[v]; 367ad8374ffSToby Isaac ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr); 368ad8374ffSToby Isaac (*labelnew)->points[v] = label->points[v]; 369c58f1c22SToby Isaac } 370c58f1c22SToby Isaac } 371c58f1c22SToby Isaac (*labelnew)->pStart = -1; 372c58f1c22SToby Isaac (*labelnew)->pEnd = -1; 373c58f1c22SToby Isaac (*labelnew)->bt = NULL; 374c58f1c22SToby Isaac PetscFunctionReturn(0); 375c58f1c22SToby Isaac } 376c58f1c22SToby Isaac 377c58f1c22SToby Isaac /* This can be hooked into SetValue(), ClearValue(), etc. for updating */ 378c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd) 379c58f1c22SToby Isaac { 380c58f1c22SToby Isaac PetscInt v; 381c58f1c22SToby Isaac PetscErrorCode ierr; 382c58f1c22SToby Isaac 383c58f1c22SToby Isaac PetscFunctionBegin; 3840c3c4a36SLisandro Dalcin ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr); 385c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 386c58f1c22SToby Isaac label->pStart = pStart; 387c58f1c22SToby Isaac label->pEnd = pEnd; 388c58f1c22SToby Isaac ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr); 389c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 390ad8374ffSToby Isaac const PetscInt *points; 391c58f1c22SToby Isaac PetscInt i; 392c58f1c22SToby Isaac 393ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 394c58f1c22SToby Isaac for (i = 0; i < label->stratumSizes[v]; ++i) { 395ad8374ffSToby Isaac const PetscInt point = points[i]; 396c58f1c22SToby Isaac 397c58f1c22SToby 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); 398c58f1c22SToby Isaac ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr); 399c58f1c22SToby Isaac } 400ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr); 401c58f1c22SToby Isaac } 402c58f1c22SToby Isaac PetscFunctionReturn(0); 403c58f1c22SToby Isaac } 404c58f1c22SToby Isaac 405c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label) 406c58f1c22SToby Isaac { 407c58f1c22SToby Isaac PetscErrorCode ierr; 408c58f1c22SToby Isaac 409c58f1c22SToby Isaac PetscFunctionBegin; 410c58f1c22SToby Isaac label->pStart = -1; 411c58f1c22SToby Isaac label->pEnd = -1; 4120c3c4a36SLisandro Dalcin ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 413c58f1c22SToby Isaac PetscFunctionReturn(0); 414c58f1c22SToby Isaac } 415c58f1c22SToby Isaac 416c58f1c22SToby Isaac /*@ 417c58f1c22SToby Isaac DMLabelHasValue - Determine whether a label assigns the value to any point 418c58f1c22SToby Isaac 419c58f1c22SToby Isaac Input Parameters: 420c58f1c22SToby Isaac + label - the DMLabel 421c58f1c22SToby Isaac - value - the value 422c58f1c22SToby Isaac 423c58f1c22SToby Isaac Output Parameter: 424c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point 425c58f1c22SToby Isaac 426c58f1c22SToby Isaac Level: developer 427c58f1c22SToby Isaac 428c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue() 429c58f1c22SToby Isaac @*/ 430c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains) 431c58f1c22SToby Isaac { 432c58f1c22SToby Isaac PetscInt v; 4330c3c4a36SLisandro Dalcin PetscErrorCode ierr; 434c58f1c22SToby Isaac 435c58f1c22SToby Isaac PetscFunctionBegin; 436c58f1c22SToby Isaac PetscValidPointer(contains, 3); 4370c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 4380c3c4a36SLisandro Dalcin *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE; 439c58f1c22SToby Isaac PetscFunctionReturn(0); 440c58f1c22SToby Isaac } 441c58f1c22SToby Isaac 442c58f1c22SToby Isaac /*@ 443c58f1c22SToby Isaac DMLabelHasPoint - Determine whether a label assigns a value to a point 444c58f1c22SToby Isaac 445c58f1c22SToby Isaac Input Parameters: 446c58f1c22SToby Isaac + label - the DMLabel 447c58f1c22SToby Isaac - point - the point 448c58f1c22SToby Isaac 449c58f1c22SToby Isaac Output Parameter: 450c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value 451c58f1c22SToby Isaac 452c58f1c22SToby Isaac Note: The user must call DMLabelCreateIndex() before this function. 453c58f1c22SToby Isaac 454c58f1c22SToby Isaac Level: developer 455c58f1c22SToby Isaac 456c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue() 457c58f1c22SToby Isaac @*/ 458c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains) 459c58f1c22SToby Isaac { 460c58f1c22SToby Isaac PetscErrorCode ierr; 461c58f1c22SToby Isaac 462c58f1c22SToby Isaac PetscFunctionBeginHot; 463c58f1c22SToby Isaac PetscValidPointer(contains, 3); 464c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 465c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG) 466c58f1c22SToby Isaac if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()"); 467c58f1c22SToby 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); 468c58f1c22SToby Isaac #endif 469c58f1c22SToby Isaac *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE; 470c58f1c22SToby Isaac PetscFunctionReturn(0); 471c58f1c22SToby Isaac } 472c58f1c22SToby Isaac 473c58f1c22SToby Isaac /*@ 474c58f1c22SToby Isaac DMLabelStratumHasPoint - Return true if the stratum contains a point 475c58f1c22SToby Isaac 476c58f1c22SToby Isaac Input Parameters: 477c58f1c22SToby Isaac + label - the DMLabel 478c58f1c22SToby Isaac . value - the stratum value 479c58f1c22SToby Isaac - point - the point 480c58f1c22SToby Isaac 481c58f1c22SToby Isaac Output Parameter: 482c58f1c22SToby Isaac . contains - true if the stratum contains the point 483c58f1c22SToby Isaac 484c58f1c22SToby Isaac Level: intermediate 485c58f1c22SToby Isaac 486c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue() 487c58f1c22SToby Isaac @*/ 488c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains) 489c58f1c22SToby Isaac { 490c58f1c22SToby Isaac PetscInt v; 491c58f1c22SToby Isaac PetscErrorCode ierr; 492c58f1c22SToby Isaac 4930c3c4a36SLisandro Dalcin PetscFunctionBeginHot; 494c58f1c22SToby Isaac PetscValidPointer(contains, 4); 495c58f1c22SToby Isaac *contains = PETSC_FALSE; 4960c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 4970c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 4980c3c4a36SLisandro Dalcin 499ad8374ffSToby Isaac if (label->validIS[v]) { 500c58f1c22SToby Isaac PetscInt i; 501c58f1c22SToby Isaac 502a2d74346SToby Isaac ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr); 5030c3c4a36SLisandro Dalcin if (i >= 0) *contains = PETSC_TRUE; 504c58f1c22SToby Isaac } else { 505c58f1c22SToby Isaac PetscBool has; 506c58f1c22SToby Isaac 507*e8f14785SLisandro Dalcin PetscHSetIHas(label->ht[v], point, &has); 5080c3c4a36SLisandro Dalcin if (has) *contains = PETSC_TRUE; 509c58f1c22SToby Isaac } 510c58f1c22SToby Isaac PetscFunctionReturn(0); 511c58f1c22SToby Isaac } 512c58f1c22SToby Isaac 51384f0b6dfSMatthew G. Knepley /*@ 5145aa44df4SToby Isaac DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 5155aa44df4SToby Isaac When a label is created, it is initialized to -1. 5165aa44df4SToby Isaac 5175aa44df4SToby Isaac Input parameter: 5185aa44df4SToby Isaac . label - a DMLabel object 5195aa44df4SToby Isaac 5205aa44df4SToby Isaac Output parameter: 5215aa44df4SToby Isaac . defaultValue - the default value 5225aa44df4SToby Isaac 5235aa44df4SToby Isaac Level: beginner 5245aa44df4SToby Isaac 5255aa44df4SToby Isaac .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 52684f0b6dfSMatthew G. Knepley @*/ 5275aa44df4SToby Isaac PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue) 5285aa44df4SToby Isaac { 5295aa44df4SToby Isaac PetscFunctionBegin; 5305aa44df4SToby Isaac *defaultValue = label->defaultValue; 5315aa44df4SToby Isaac PetscFunctionReturn(0); 5325aa44df4SToby Isaac } 5335aa44df4SToby Isaac 53484f0b6dfSMatthew G. Knepley /*@ 5355aa44df4SToby Isaac DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 5365aa44df4SToby Isaac When a label is created, it is initialized to -1. 5375aa44df4SToby Isaac 5385aa44df4SToby Isaac Input parameter: 5395aa44df4SToby Isaac . label - a DMLabel object 5405aa44df4SToby Isaac 5415aa44df4SToby Isaac Output parameter: 5425aa44df4SToby Isaac . defaultValue - the default value 5435aa44df4SToby Isaac 5445aa44df4SToby Isaac Level: beginner 5455aa44df4SToby Isaac 5465aa44df4SToby Isaac .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 54784f0b6dfSMatthew G. Knepley @*/ 5485aa44df4SToby Isaac PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue) 5495aa44df4SToby Isaac { 5505aa44df4SToby Isaac PetscFunctionBegin; 5515aa44df4SToby Isaac label->defaultValue = defaultValue; 5525aa44df4SToby Isaac PetscFunctionReturn(0); 5535aa44df4SToby Isaac } 5545aa44df4SToby Isaac 555c58f1c22SToby Isaac /*@ 5565aa44df4SToby 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()) 557c58f1c22SToby Isaac 558c58f1c22SToby Isaac Input Parameters: 559c58f1c22SToby Isaac + label - the DMLabel 560c58f1c22SToby Isaac - point - the point 561c58f1c22SToby Isaac 562c58f1c22SToby Isaac Output Parameter: 563c58f1c22SToby Isaac . value - The point value, or -1 564c58f1c22SToby Isaac 565c58f1c22SToby Isaac Level: intermediate 566c58f1c22SToby Isaac 5675aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 568c58f1c22SToby Isaac @*/ 569c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value) 570c58f1c22SToby Isaac { 571c58f1c22SToby Isaac PetscInt v; 572c58f1c22SToby Isaac PetscErrorCode ierr; 573c58f1c22SToby Isaac 5740c3c4a36SLisandro Dalcin PetscFunctionBeginHot; 575c58f1c22SToby Isaac PetscValidPointer(value, 3); 5765aa44df4SToby Isaac *value = label->defaultValue; 577c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 578ad8374ffSToby Isaac if (label->validIS[v]) { 579c58f1c22SToby Isaac PetscInt i; 580c58f1c22SToby Isaac 581a2d74346SToby Isaac ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr); 582c58f1c22SToby Isaac if (i >= 0) { 583c58f1c22SToby Isaac *value = label->stratumValues[v]; 584c58f1c22SToby Isaac break; 585c58f1c22SToby Isaac } 586c58f1c22SToby Isaac } else { 587c58f1c22SToby Isaac PetscBool has; 588c58f1c22SToby Isaac 589*e8f14785SLisandro Dalcin PetscHSetIHas(label->ht[v], point, &has); 590c58f1c22SToby Isaac if (has) { 591c58f1c22SToby Isaac *value = label->stratumValues[v]; 592c58f1c22SToby Isaac break; 593c58f1c22SToby Isaac } 594c58f1c22SToby Isaac } 595c58f1c22SToby Isaac } 596c58f1c22SToby Isaac PetscFunctionReturn(0); 597c58f1c22SToby Isaac } 598c58f1c22SToby Isaac 599c58f1c22SToby Isaac /*@ 6005aa44df4SToby Isaac DMLabelSetValue - Set the value a label assigns to a point. If the value is the same as the label's default value (which is initially -1, and can be changed with DMLabelSetDefaultValue() to somethingg different), then this function will do nothing. 601c58f1c22SToby Isaac 602c58f1c22SToby Isaac Input Parameters: 603c58f1c22SToby Isaac + label - the DMLabel 604c58f1c22SToby Isaac . point - the point 605c58f1c22SToby Isaac - value - The point value 606c58f1c22SToby Isaac 607c58f1c22SToby Isaac Level: intermediate 608c58f1c22SToby Isaac 6095aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 610c58f1c22SToby Isaac @*/ 611c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value) 612c58f1c22SToby Isaac { 613c58f1c22SToby Isaac PetscInt v; 614c58f1c22SToby Isaac PetscErrorCode ierr; 615c58f1c22SToby Isaac 616c58f1c22SToby Isaac PetscFunctionBegin; 6170c3c4a36SLisandro Dalcin /* Find label value, add new entry if needed */ 6185aa44df4SToby Isaac if (value == label->defaultValue) PetscFunctionReturn(0); 6190c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 6200c3c4a36SLisandro Dalcin if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);} 621c58f1c22SToby Isaac /* Set key */ 6220c3c4a36SLisandro Dalcin ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 623*e8f14785SLisandro Dalcin ierr = PetscHSetIAdd(label->ht[v], point);CHKERRQ(ierr); 624c58f1c22SToby Isaac PetscFunctionReturn(0); 625c58f1c22SToby Isaac } 626c58f1c22SToby Isaac 627c58f1c22SToby Isaac /*@ 628c58f1c22SToby Isaac DMLabelClearValue - Clear the value a label assigns to a point 629c58f1c22SToby Isaac 630c58f1c22SToby Isaac Input Parameters: 631c58f1c22SToby Isaac + label - the DMLabel 632c58f1c22SToby Isaac . point - the point 633c58f1c22SToby Isaac - value - The point value 634c58f1c22SToby Isaac 635c58f1c22SToby Isaac Level: intermediate 636c58f1c22SToby Isaac 637c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue() 638c58f1c22SToby Isaac @*/ 639c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value) 640c58f1c22SToby Isaac { 641ad8374ffSToby Isaac PetscInt v; 642c58f1c22SToby Isaac PetscErrorCode ierr; 643c58f1c22SToby Isaac 644c58f1c22SToby Isaac PetscFunctionBegin; 645c58f1c22SToby Isaac /* Find label value */ 6460c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 6470c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 6480c3c4a36SLisandro Dalcin 649eeed21e7SToby Isaac if (label->bt) { 650eeed21e7SToby 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); 651eeed21e7SToby Isaac ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 652eeed21e7SToby Isaac } 6530c3c4a36SLisandro Dalcin 6540c3c4a36SLisandro Dalcin /* Delete key */ 6550c3c4a36SLisandro Dalcin ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 656*e8f14785SLisandro Dalcin ierr = PetscHSetIDel(label->ht[v], point);CHKERRQ(ierr); 657c58f1c22SToby Isaac PetscFunctionReturn(0); 658c58f1c22SToby Isaac } 659c58f1c22SToby Isaac 660c58f1c22SToby Isaac /*@ 661c58f1c22SToby Isaac DMLabelInsertIS - Set all points in the IS to a value 662c58f1c22SToby Isaac 663c58f1c22SToby Isaac Input Parameters: 664c58f1c22SToby Isaac + label - the DMLabel 665c58f1c22SToby Isaac . is - the point IS 666c58f1c22SToby Isaac - value - The point value 667c58f1c22SToby Isaac 668c58f1c22SToby Isaac Level: intermediate 669c58f1c22SToby Isaac 670c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 671c58f1c22SToby Isaac @*/ 672c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value) 673c58f1c22SToby Isaac { 6740c3c4a36SLisandro Dalcin PetscInt v, n, p; 675c58f1c22SToby Isaac const PetscInt *points; 676c58f1c22SToby Isaac PetscErrorCode ierr; 677c58f1c22SToby Isaac 678c58f1c22SToby Isaac PetscFunctionBegin; 679c58f1c22SToby Isaac PetscValidHeaderSpecific(is, IS_CLASSID, 2); 6800c3c4a36SLisandro Dalcin /* Find label value, add new entry if needed */ 6810c3c4a36SLisandro Dalcin if (value == label->defaultValue) PetscFunctionReturn(0); 6820c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 6830c3c4a36SLisandro Dalcin if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);} 6840c3c4a36SLisandro Dalcin /* Set keys */ 6850c3c4a36SLisandro Dalcin ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 686c58f1c22SToby Isaac ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 687c58f1c22SToby Isaac ierr = ISGetIndices(is, &points);CHKERRQ(ierr); 688*e8f14785SLisandro Dalcin for (p = 0; p < n; ++p) {ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);} 689c58f1c22SToby Isaac ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr); 690c58f1c22SToby Isaac PetscFunctionReturn(0); 691c58f1c22SToby Isaac } 692c58f1c22SToby Isaac 69384f0b6dfSMatthew G. Knepley /*@ 69484f0b6dfSMatthew G. Knepley DMLabelGetNumValues - Get the number of values that the DMLabel takes 69584f0b6dfSMatthew G. Knepley 69684f0b6dfSMatthew G. Knepley Input Parameter: 69784f0b6dfSMatthew G. Knepley . label - the DMLabel 69884f0b6dfSMatthew G. Knepley 69984f0b6dfSMatthew G. Knepley Output Paramater: 70084f0b6dfSMatthew G. Knepley . numValues - the number of values 70184f0b6dfSMatthew G. Knepley 70284f0b6dfSMatthew G. Knepley Level: intermediate 70384f0b6dfSMatthew G. Knepley 70484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 70584f0b6dfSMatthew G. Knepley @*/ 706c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues) 707c58f1c22SToby Isaac { 708c58f1c22SToby Isaac PetscFunctionBegin; 709c58f1c22SToby Isaac PetscValidPointer(numValues, 2); 710c58f1c22SToby Isaac *numValues = label->numStrata; 711c58f1c22SToby Isaac PetscFunctionReturn(0); 712c58f1c22SToby Isaac } 713c58f1c22SToby Isaac 71484f0b6dfSMatthew G. Knepley /*@ 71584f0b6dfSMatthew G. Knepley DMLabelGetValueIS - Get an IS of all values that the DMlabel takes 71684f0b6dfSMatthew G. Knepley 71784f0b6dfSMatthew G. Knepley Input Parameter: 71884f0b6dfSMatthew G. Knepley . label - the DMLabel 71984f0b6dfSMatthew G. Knepley 72084f0b6dfSMatthew G. Knepley Output Paramater: 72184f0b6dfSMatthew G. Knepley . is - the value IS 72284f0b6dfSMatthew G. Knepley 72384f0b6dfSMatthew G. Knepley Level: intermediate 72484f0b6dfSMatthew G. Knepley 72584f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 72684f0b6dfSMatthew G. Knepley @*/ 727c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values) 728c58f1c22SToby Isaac { 729c58f1c22SToby Isaac PetscErrorCode ierr; 730c58f1c22SToby Isaac 731c58f1c22SToby Isaac PetscFunctionBegin; 732c58f1c22SToby Isaac PetscValidPointer(values, 2); 733c58f1c22SToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr); 734c58f1c22SToby Isaac PetscFunctionReturn(0); 735c58f1c22SToby Isaac } 736c58f1c22SToby Isaac 73784f0b6dfSMatthew G. Knepley /*@ 73884f0b6dfSMatthew G. Knepley DMLabelHasStratum - Determine whether points exist with the given value 73984f0b6dfSMatthew G. Knepley 74084f0b6dfSMatthew G. Knepley Input Parameters: 74184f0b6dfSMatthew G. Knepley + label - the DMLabel 74284f0b6dfSMatthew G. Knepley - value - the stratum value 74384f0b6dfSMatthew G. Knepley 74484f0b6dfSMatthew G. Knepley Output Paramater: 74584f0b6dfSMatthew G. Knepley . exists - Flag saying whether points exist 74684f0b6dfSMatthew G. Knepley 74784f0b6dfSMatthew G. Knepley Level: intermediate 74884f0b6dfSMatthew G. Knepley 74984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 75084f0b6dfSMatthew G. Knepley @*/ 751fada774cSMatthew G. Knepley PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists) 752fada774cSMatthew G. Knepley { 753fada774cSMatthew G. Knepley PetscInt v; 7540c3c4a36SLisandro Dalcin PetscErrorCode ierr; 755fada774cSMatthew G. Knepley 756fada774cSMatthew G. Knepley PetscFunctionBegin; 757fada774cSMatthew G. Knepley PetscValidPointer(exists, 3); 7580c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 7590c3c4a36SLisandro Dalcin *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE; 760fada774cSMatthew G. Knepley PetscFunctionReturn(0); 761fada774cSMatthew G. Knepley } 762fada774cSMatthew G. Knepley 76384f0b6dfSMatthew G. Knepley /*@ 76484f0b6dfSMatthew G. Knepley DMLabelGetStratumSize - Get the size of a stratum 76584f0b6dfSMatthew G. Knepley 76684f0b6dfSMatthew G. Knepley Input Parameters: 76784f0b6dfSMatthew G. Knepley + label - the DMLabel 76884f0b6dfSMatthew G. Knepley - value - the stratum value 76984f0b6dfSMatthew G. Knepley 77084f0b6dfSMatthew G. Knepley Output Paramater: 77184f0b6dfSMatthew G. Knepley . size - The number of points in the stratum 77284f0b6dfSMatthew G. Knepley 77384f0b6dfSMatthew G. Knepley Level: intermediate 77484f0b6dfSMatthew G. Knepley 77584f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 77684f0b6dfSMatthew G. Knepley @*/ 777c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size) 778c58f1c22SToby Isaac { 779c58f1c22SToby Isaac PetscInt v; 780c58f1c22SToby Isaac PetscErrorCode ierr; 781c58f1c22SToby Isaac 782c58f1c22SToby Isaac PetscFunctionBegin; 783c58f1c22SToby Isaac PetscValidPointer(size, 3); 784c58f1c22SToby Isaac *size = 0; 7850c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 7860c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 787c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 788c58f1c22SToby Isaac *size = label->stratumSizes[v]; 789c58f1c22SToby Isaac PetscFunctionReturn(0); 790c58f1c22SToby Isaac } 791c58f1c22SToby Isaac 79284f0b6dfSMatthew G. Knepley /*@ 79384f0b6dfSMatthew G. Knepley DMLabelGetStratumBounds - Get the largest and smallest point of a stratum 79484f0b6dfSMatthew G. Knepley 79584f0b6dfSMatthew G. Knepley Input Parameters: 79684f0b6dfSMatthew G. Knepley + label - the DMLabel 79784f0b6dfSMatthew G. Knepley - value - the stratum value 79884f0b6dfSMatthew G. Knepley 79984f0b6dfSMatthew G. Knepley Output Paramaters: 80084f0b6dfSMatthew G. Knepley + start - the smallest point in the stratum 80184f0b6dfSMatthew G. Knepley - end - the largest point in the stratum 80284f0b6dfSMatthew G. Knepley 80384f0b6dfSMatthew G. Knepley Level: intermediate 80484f0b6dfSMatthew G. Knepley 80584f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 80684f0b6dfSMatthew G. Knepley @*/ 807c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end) 808c58f1c22SToby Isaac { 8090c3c4a36SLisandro Dalcin PetscInt v, min, max; 810c58f1c22SToby Isaac PetscErrorCode ierr; 811c58f1c22SToby Isaac 812c58f1c22SToby Isaac PetscFunctionBegin; 813c58f1c22SToby Isaac if (start) {PetscValidPointer(start, 3); *start = 0;} 814c58f1c22SToby Isaac if (end) {PetscValidPointer(end, 4); *end = 0;} 8150c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 8160c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 817c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 8180c3c4a36SLisandro Dalcin if (label->stratumSizes[v] <= 0) PetscFunctionReturn(0); 819d6cb179aSToby Isaac ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr); 820d6cb179aSToby Isaac if (start) *start = min; 821d6cb179aSToby Isaac if (end) *end = max+1; 822c58f1c22SToby Isaac PetscFunctionReturn(0); 823c58f1c22SToby Isaac } 824c58f1c22SToby Isaac 82584f0b6dfSMatthew G. Knepley /*@ 82684f0b6dfSMatthew G. Knepley DMLabelGetStratumIS - Get an IS with the stratum points 82784f0b6dfSMatthew G. Knepley 82884f0b6dfSMatthew G. Knepley Input Parameters: 82984f0b6dfSMatthew G. Knepley + label - the DMLabel 83084f0b6dfSMatthew G. Knepley - value - the stratum value 83184f0b6dfSMatthew G. Knepley 83284f0b6dfSMatthew G. Knepley Output Paramater: 83384f0b6dfSMatthew G. Knepley . points - The stratum points 83484f0b6dfSMatthew G. Knepley 83584f0b6dfSMatthew G. Knepley Level: intermediate 83684f0b6dfSMatthew G. Knepley 83784f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 83884f0b6dfSMatthew G. Knepley @*/ 839c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points) 840c58f1c22SToby Isaac { 841c58f1c22SToby Isaac PetscInt v; 842c58f1c22SToby Isaac PetscErrorCode ierr; 843c58f1c22SToby Isaac 844c58f1c22SToby Isaac PetscFunctionBegin; 845c58f1c22SToby Isaac PetscValidPointer(points, 3); 846c58f1c22SToby Isaac *points = NULL; 8470c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 8480c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 849c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 850ad8374ffSToby Isaac ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr); 851ad8374ffSToby Isaac *points = label->points[v]; 852c58f1c22SToby Isaac PetscFunctionReturn(0); 853c58f1c22SToby Isaac } 854c58f1c22SToby Isaac 85584f0b6dfSMatthew G. Knepley /*@ 85684f0b6dfSMatthew G. Knepley DMLabelSetStratumIS - Set the stratum points using an IS 85784f0b6dfSMatthew G. Knepley 85884f0b6dfSMatthew G. Knepley Input Parameters: 85984f0b6dfSMatthew G. Knepley + label - the DMLabel 86084f0b6dfSMatthew G. Knepley . value - the stratum value 86184f0b6dfSMatthew G. Knepley - points - The stratum points 86284f0b6dfSMatthew G. Knepley 86384f0b6dfSMatthew G. Knepley Level: intermediate 86484f0b6dfSMatthew G. Knepley 86584f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 86684f0b6dfSMatthew G. Knepley @*/ 8674de306b1SToby Isaac PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is) 8684de306b1SToby Isaac { 8690c3c4a36SLisandro Dalcin PetscInt v; 8704de306b1SToby Isaac PetscErrorCode ierr; 8714de306b1SToby Isaac 8724de306b1SToby Isaac PetscFunctionBegin; 8730c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 8740c3c4a36SLisandro Dalcin if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);} 8754de306b1SToby Isaac if (is == label->points[v]) PetscFunctionReturn(0); 8764de306b1SToby Isaac ierr = DMLabelClearStratum(label, value);CHKERRQ(ierr); 8774de306b1SToby Isaac ierr = ISGetLocalSize(is, &(label->stratumSizes[v]));CHKERRQ(ierr); 8784de306b1SToby Isaac ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 8794de306b1SToby Isaac ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 8800c3c4a36SLisandro Dalcin label->points[v] = is; 8810c3c4a36SLisandro Dalcin label->validIS[v] = PETSC_TRUE; 8824de306b1SToby Isaac if (label->bt) { 8834de306b1SToby Isaac const PetscInt *points; 8844de306b1SToby Isaac PetscInt p; 8854de306b1SToby Isaac 8864de306b1SToby Isaac ierr = ISGetIndices(is,&points);CHKERRQ(ierr); 8874de306b1SToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) { 8884de306b1SToby Isaac const PetscInt point = points[p]; 8894de306b1SToby Isaac 8904de306b1SToby 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); 8914de306b1SToby Isaac ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 8924de306b1SToby Isaac } 8934de306b1SToby Isaac } 8944de306b1SToby Isaac PetscFunctionReturn(0); 8954de306b1SToby Isaac } 8964de306b1SToby Isaac 89784f0b6dfSMatthew G. Knepley /*@ 89884f0b6dfSMatthew G. Knepley DMLabelClearStratum - Remove a stratum 8994de306b1SToby Isaac 90084f0b6dfSMatthew G. Knepley Input Parameters: 90184f0b6dfSMatthew G. Knepley + label - the DMLabel 90284f0b6dfSMatthew G. Knepley - value - the stratum value 90384f0b6dfSMatthew G. Knepley 90484f0b6dfSMatthew G. Knepley Level: intermediate 90584f0b6dfSMatthew G. Knepley 90684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 90784f0b6dfSMatthew G. Knepley @*/ 908c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value) 909c58f1c22SToby Isaac { 910c58f1c22SToby Isaac PetscInt v; 911c58f1c22SToby Isaac PetscErrorCode ierr; 912c58f1c22SToby Isaac 913c58f1c22SToby Isaac PetscFunctionBegin; 9140c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 9150c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 9164de306b1SToby Isaac if (label->validIS[v]) { 9174de306b1SToby Isaac if (label->bt) { 918c58f1c22SToby Isaac PetscInt i; 919ad8374ffSToby Isaac const PetscInt *points; 920c58f1c22SToby Isaac 921ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 922c58f1c22SToby Isaac for (i = 0; i < label->stratumSizes[v]; ++i) { 923ad8374ffSToby Isaac const PetscInt point = points[i]; 924c58f1c22SToby Isaac 925c58f1c22SToby 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); 926c58f1c22SToby Isaac ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 927c58f1c22SToby Isaac } 928ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr); 929c58f1c22SToby Isaac } 930c58f1c22SToby Isaac label->stratumSizes[v] = 0; 9310c3c4a36SLisandro Dalcin ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr); 9320c3c4a36SLisandro Dalcin ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr); 9330c3c4a36SLisandro Dalcin ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr); 934c58f1c22SToby Isaac } else { 935*e8f14785SLisandro Dalcin PetscHSetIClear(label->ht[v]); 936c58f1c22SToby Isaac } 937c58f1c22SToby Isaac PetscFunctionReturn(0); 938c58f1c22SToby Isaac } 939c58f1c22SToby Isaac 94084f0b6dfSMatthew G. Knepley /*@ 94184f0b6dfSMatthew G. Knepley DMLabelFilter - Remove all points outside of [start, end) 94284f0b6dfSMatthew G. Knepley 94384f0b6dfSMatthew G. Knepley Input Parameters: 94484f0b6dfSMatthew G. Knepley + label - the DMLabel 94584f0b6dfSMatthew G. Knepley . start - the first point 94684f0b6dfSMatthew G. Knepley - end - the last point 94784f0b6dfSMatthew G. Knepley 94884f0b6dfSMatthew G. Knepley Level: intermediate 94984f0b6dfSMatthew G. Knepley 95084f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 95184f0b6dfSMatthew G. Knepley @*/ 952c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end) 953c58f1c22SToby Isaac { 954c58f1c22SToby Isaac PetscInt v; 955c58f1c22SToby Isaac PetscErrorCode ierr; 956c58f1c22SToby Isaac 957c58f1c22SToby Isaac PetscFunctionBegin; 9580c3c4a36SLisandro Dalcin ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr); 959c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 960c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 961c58f1c22SToby Isaac PetscInt off, q; 962ad8374ffSToby Isaac const PetscInt *points; 963033405d5SLisandro Dalcin PetscInt numPointsNew = 0, *pointsNew = NULL; 964c58f1c22SToby Isaac 965ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 966033405d5SLisandro Dalcin for (q = 0; q < label->stratumSizes[v]; ++q) 967033405d5SLisandro Dalcin if (points[q] >= start && points[q] < end) 968033405d5SLisandro Dalcin numPointsNew++; 969033405d5SLisandro Dalcin ierr = PetscMalloc1(numPointsNew, &pointsNew);CHKERRQ(ierr); 970c58f1c22SToby Isaac for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) { 971033405d5SLisandro Dalcin if (points[q] >= start && points[q] < end) 972033405d5SLisandro Dalcin pointsNew[off++] = points[q]; 973ad8374ffSToby Isaac } 974ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 975033405d5SLisandro Dalcin 976033405d5SLisandro Dalcin label->stratumSizes[v] = numPointsNew; 977033405d5SLisandro Dalcin ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr); 978033405d5SLisandro Dalcin ierr = ISCreateGeneral(PETSC_COMM_SELF,numPointsNew, pointsNew, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr); 979033405d5SLisandro Dalcin ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr); 980c58f1c22SToby Isaac } 981c58f1c22SToby Isaac ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr); 982c58f1c22SToby Isaac PetscFunctionReturn(0); 983c58f1c22SToby Isaac } 984c58f1c22SToby Isaac 98584f0b6dfSMatthew G. Knepley /*@ 98684f0b6dfSMatthew G. Knepley DMLabelPermute - Create a new label with permuted points 98784f0b6dfSMatthew G. Knepley 98884f0b6dfSMatthew G. Knepley Input Parameters: 98984f0b6dfSMatthew G. Knepley + label - the DMLabel 99084f0b6dfSMatthew G. Knepley - permutation - the point permutation 99184f0b6dfSMatthew G. Knepley 99284f0b6dfSMatthew G. Knepley Output Parameter: 99384f0b6dfSMatthew G. Knepley . labelnew - the new label containing the permuted points 99484f0b6dfSMatthew G. Knepley 99584f0b6dfSMatthew G. Knepley Level: intermediate 99684f0b6dfSMatthew G. Knepley 99784f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 99884f0b6dfSMatthew G. Knepley @*/ 999c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew) 1000c58f1c22SToby Isaac { 1001c58f1c22SToby Isaac const PetscInt *perm; 1002c58f1c22SToby Isaac PetscInt numValues, numPoints, v, q; 1003c58f1c22SToby Isaac PetscErrorCode ierr; 1004c58f1c22SToby Isaac 1005c58f1c22SToby Isaac PetscFunctionBegin; 1006c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 1007c58f1c22SToby Isaac ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr); 1008c58f1c22SToby Isaac ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr); 1009c58f1c22SToby Isaac ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr); 1010c58f1c22SToby Isaac ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr); 1011c58f1c22SToby Isaac for (v = 0; v < numValues; ++v) { 1012c58f1c22SToby Isaac const PetscInt size = (*labelNew)->stratumSizes[v]; 1013ad8374ffSToby Isaac const PetscInt *points; 1014ad8374ffSToby Isaac PetscInt *pointsNew; 1015c58f1c22SToby Isaac 1016ad8374ffSToby Isaac ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 1017ad8374ffSToby Isaac ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr); 1018c58f1c22SToby Isaac for (q = 0; q < size; ++q) { 1019ad8374ffSToby Isaac const PetscInt point = points[q]; 1020c58f1c22SToby Isaac 1021c58f1c22SToby 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); 1022ad8374ffSToby Isaac pointsNew[q] = perm[point]; 1023c58f1c22SToby Isaac } 1024ad8374ffSToby Isaac ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 1025ad8374ffSToby Isaac ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr); 1026ad8374ffSToby Isaac ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr); 1027fa8e8ae5SToby Isaac if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) { 1028fa8e8ae5SToby Isaac ierr = ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));CHKERRQ(ierr); 1029fa8e8ae5SToby Isaac ierr = PetscFree(pointsNew);CHKERRQ(ierr); 1030fa8e8ae5SToby Isaac } else { 1031ad8374ffSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr); 1032fa8e8ae5SToby Isaac } 1033ad8374ffSToby Isaac ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr); 1034c58f1c22SToby Isaac } 1035c58f1c22SToby Isaac ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr); 1036c58f1c22SToby Isaac if (label->bt) { 1037c58f1c22SToby Isaac ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 1038c58f1c22SToby Isaac ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr); 1039c58f1c22SToby Isaac } 1040c58f1c22SToby Isaac PetscFunctionReturn(0); 1041c58f1c22SToby Isaac } 1042c58f1c22SToby Isaac 104326c55118SMichael Lange PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata) 104426c55118SMichael Lange { 104526c55118SMichael Lange MPI_Comm comm; 104626c55118SMichael Lange PetscInt s, l, nroots, nleaves, dof, offset, size; 104726c55118SMichael Lange PetscInt *remoteOffsets, *rootStrata, *rootIdx; 104826c55118SMichael Lange PetscSection rootSection; 104926c55118SMichael Lange PetscSF labelSF; 105026c55118SMichael Lange PetscErrorCode ierr; 105126c55118SMichael Lange 105226c55118SMichael Lange PetscFunctionBegin; 105326c55118SMichael Lange if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 105426c55118SMichael Lange ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 105526c55118SMichael Lange /* Build a section of stratum values per point, generate the according SF 105626c55118SMichael Lange and distribute point-wise stratum values to leaves. */ 105726c55118SMichael Lange ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 105826c55118SMichael Lange ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 105926c55118SMichael Lange ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr); 106026c55118SMichael Lange if (label) { 106126c55118SMichael Lange for (s = 0; s < label->numStrata; ++s) { 1062ad8374ffSToby Isaac const PetscInt *points; 1063ad8374ffSToby Isaac 1064ad8374ffSToby Isaac ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 106526c55118SMichael Lange for (l = 0; l < label->stratumSizes[s]; l++) { 1066ad8374ffSToby Isaac ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr); 1067ad8374ffSToby Isaac ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr); 106826c55118SMichael Lange } 1069ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 107026c55118SMichael Lange } 107126c55118SMichael Lange } 107226c55118SMichael Lange ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 107326c55118SMichael Lange /* Create a point-wise array of stratum values */ 107426c55118SMichael Lange ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr); 107526c55118SMichael Lange ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr); 107626c55118SMichael Lange ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr); 107726c55118SMichael Lange if (label) { 107826c55118SMichael Lange for (s = 0; s < label->numStrata; ++s) { 1079ad8374ffSToby Isaac const PetscInt *points; 1080ad8374ffSToby Isaac 1081ad8374ffSToby Isaac ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 108226c55118SMichael Lange for (l = 0; l < label->stratumSizes[s]; l++) { 1083ad8374ffSToby Isaac const PetscInt p = points[l]; 108426c55118SMichael Lange ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr); 108526c55118SMichael Lange rootStrata[offset+rootIdx[p]++] = label->stratumValues[s]; 108626c55118SMichael Lange } 1087ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 108826c55118SMichael Lange } 108926c55118SMichael Lange } 109026c55118SMichael Lange /* Build SF that maps label points to remote processes */ 109126c55118SMichael Lange ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr); 109226c55118SMichael Lange ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr); 109326c55118SMichael Lange ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr); 109426c55118SMichael Lange ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 109526c55118SMichael Lange /* Send the strata for each point over the derived SF */ 109626c55118SMichael Lange ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr); 109726c55118SMichael Lange ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr); 109826c55118SMichael Lange ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 109926c55118SMichael Lange ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 110026c55118SMichael Lange /* Clean up */ 110126c55118SMichael Lange ierr = PetscFree(rootStrata);CHKERRQ(ierr); 110226c55118SMichael Lange ierr = PetscFree(rootIdx);CHKERRQ(ierr); 110326c55118SMichael Lange ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 110426c55118SMichael Lange ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr); 110526c55118SMichael Lange PetscFunctionReturn(0); 110626c55118SMichael Lange } 110726c55118SMichael Lange 110884f0b6dfSMatthew G. Knepley /*@ 110984f0b6dfSMatthew G. Knepley DMLabelDistribute - Create a new label pushed forward over the PetscSF 111084f0b6dfSMatthew G. Knepley 111184f0b6dfSMatthew G. Knepley Input Parameters: 111284f0b6dfSMatthew G. Knepley + label - the DMLabel 111384f0b6dfSMatthew G. Knepley - sf - the map from old to new distribution 111484f0b6dfSMatthew G. Knepley 111584f0b6dfSMatthew G. Knepley Output Parameter: 111684f0b6dfSMatthew G. Knepley . labelnew - the new redistributed label 111784f0b6dfSMatthew G. Knepley 111884f0b6dfSMatthew G. Knepley Level: intermediate 111984f0b6dfSMatthew G. Knepley 112084f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 112184f0b6dfSMatthew G. Knepley @*/ 1122c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew) 1123c58f1c22SToby Isaac { 1124c58f1c22SToby Isaac MPI_Comm comm; 112526c55118SMichael Lange PetscSection leafSection; 112626c55118SMichael Lange PetscInt p, pStart, pEnd, s, size, dof, offset, stratum; 112726c55118SMichael Lange PetscInt *leafStrata, *strataIdx; 1128ad8374ffSToby Isaac PetscInt **points; 1129c58f1c22SToby Isaac char *name; 1130c58f1c22SToby Isaac PetscInt nameSize; 1131*e8f14785SLisandro Dalcin PetscHSetI stratumHash; 1132c58f1c22SToby Isaac size_t len = 0; 113326c55118SMichael Lange PetscMPIInt rank; 1134c58f1c22SToby Isaac PetscErrorCode ierr; 1135c58f1c22SToby Isaac 1136c58f1c22SToby Isaac PetscFunctionBegin; 1137c58f1c22SToby Isaac if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 1138c58f1c22SToby Isaac ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 1139c58f1c22SToby Isaac ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1140c58f1c22SToby Isaac /* Bcast name */ 1141c58f1c22SToby Isaac if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 1142c58f1c22SToby Isaac nameSize = len; 1143c58f1c22SToby Isaac ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1144c58f1c22SToby Isaac ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1145c58f1c22SToby Isaac if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 1146c58f1c22SToby Isaac ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1147c58f1c22SToby Isaac ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 1148c58f1c22SToby Isaac ierr = PetscFree(name);CHKERRQ(ierr); 114977d236dfSMichael Lange /* Bcast defaultValue */ 115077d236dfSMichael Lange if (!rank) (*labelNew)->defaultValue = label->defaultValue; 115177d236dfSMichael Lange ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 115226c55118SMichael Lange /* Distribute stratum values over the SF and get the point mapping on the receiver */ 115326c55118SMichael Lange ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr); 11545cbdf6fcSMichael Lange /* Determine received stratum values and initialise new label*/ 1155*e8f14785SLisandro Dalcin ierr = PetscHSetICreate(&stratumHash);CHKERRQ(ierr); 115626c55118SMichael Lange ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr); 1157*e8f14785SLisandro Dalcin for (p = 0; p < size; ++p) {ierr = PetscHSetIAdd(stratumHash, leafStrata[p]);CHKERRQ(ierr);} 1158*e8f14785SLisandro Dalcin ierr = PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata);CHKERRQ(ierr); 1159ad8374ffSToby Isaac ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr); 1160ad8374ffSToby Isaac for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE; 11615cbdf6fcSMichael Lange ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr); 11625cbdf6fcSMichael Lange /* Turn leafStrata into indices rather than stratum values */ 11635cbdf6fcSMichael Lange offset = 0; 1164*e8f14785SLisandro Dalcin ierr = PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr); 1165a205f1fdSToby Isaac ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr); 11665cbdf6fcSMichael Lange for (p = 0; p < size; ++p) { 1167231b9e6fSMatthew G. Knepley for (s = 0; s < (*labelNew)->numStrata; ++s) { 1168231b9e6fSMatthew G. Knepley if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;} 11695cbdf6fcSMichael Lange } 11705cbdf6fcSMichael Lange } 1171c58f1c22SToby Isaac /* Rebuild the point strata on the receiver */ 1172c58f1c22SToby Isaac ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr); 1173c58f1c22SToby Isaac ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1174c58f1c22SToby Isaac for (p=pStart; p<pEnd; p++) { 1175c58f1c22SToby Isaac ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1176c58f1c22SToby Isaac ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1177c58f1c22SToby Isaac for (s=0; s<dof; s++) { 1178c58f1c22SToby Isaac (*labelNew)->stratumSizes[leafStrata[offset+s]]++; 1179c58f1c22SToby Isaac } 1180c58f1c22SToby Isaac } 1181c58f1c22SToby Isaac ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr); 1182c58f1c22SToby Isaac ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr); 1183ad8374ffSToby Isaac ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr); 1184c58f1c22SToby Isaac for (s = 0; s < (*labelNew)->numStrata; ++s) { 1185*e8f14785SLisandro Dalcin ierr = PetscHSetICreate(&(*labelNew)->ht[s]);CHKERRQ(ierr); 1186ad8374ffSToby Isaac ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr); 1187c58f1c22SToby Isaac } 1188c58f1c22SToby Isaac /* Insert points into new strata */ 1189c58f1c22SToby Isaac ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr); 1190c58f1c22SToby Isaac ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1191c58f1c22SToby Isaac for (p=pStart; p<pEnd; p++) { 1192c58f1c22SToby Isaac ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1193c58f1c22SToby Isaac ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1194c58f1c22SToby Isaac for (s=0; s<dof; s++) { 1195c58f1c22SToby Isaac stratum = leafStrata[offset+s]; 1196ad8374ffSToby Isaac points[stratum][strataIdx[stratum]++] = p; 1197c58f1c22SToby Isaac } 1198c58f1c22SToby Isaac } 1199ad8374ffSToby Isaac for (s = 0; s < (*labelNew)->numStrata; s++) { 1200ad8374ffSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr); 1201ad8374ffSToby Isaac ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr); 1202ad8374ffSToby Isaac } 1203ad8374ffSToby Isaac ierr = PetscFree(points);CHKERRQ(ierr); 1204*e8f14785SLisandro Dalcin ierr = PetscHSetIDestroy(&stratumHash);CHKERRQ(ierr); 1205c58f1c22SToby Isaac ierr = PetscFree(leafStrata);CHKERRQ(ierr); 1206c58f1c22SToby Isaac ierr = PetscFree(strataIdx);CHKERRQ(ierr); 1207c58f1c22SToby Isaac ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 1208c58f1c22SToby Isaac PetscFunctionReturn(0); 1209c58f1c22SToby Isaac } 1210c58f1c22SToby Isaac 12117937d9ceSMichael Lange /*@ 12127937d9ceSMichael Lange DMLabelGather - Gather all label values from leafs into roots 12137937d9ceSMichael Lange 12147937d9ceSMichael Lange Input Parameters: 12157937d9ceSMichael Lange + label - the DMLabel 121684f0b6dfSMatthew G. Knepley - sf - the Star Forest point communication map 12177937d9ceSMichael Lange 121884f0b6dfSMatthew G. Knepley Output Parameters: 121984f0b6dfSMatthew G. Knepley . labelNew - the new DMLabel with localised leaf values 12207937d9ceSMichael Lange 12217937d9ceSMichael Lange Level: developer 12227937d9ceSMichael Lange 12237937d9ceSMichael Lange Note: This is the inverse operation to DMLabelDistribute. 12247937d9ceSMichael Lange 12257937d9ceSMichael Lange .seealso: DMLabelDistribute() 12267937d9ceSMichael Lange @*/ 12277937d9ceSMichael Lange PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew) 12287937d9ceSMichael Lange { 12297937d9ceSMichael Lange MPI_Comm comm; 12307937d9ceSMichael Lange PetscSection rootSection; 12317937d9ceSMichael Lange PetscSF sfLabel; 12327937d9ceSMichael Lange PetscSFNode *rootPoints, *leafPoints; 12337937d9ceSMichael Lange PetscInt p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset; 12347937d9ceSMichael Lange const PetscInt *rootDegree, *ilocal; 12357937d9ceSMichael Lange PetscInt *rootStrata; 12367937d9ceSMichael Lange char *name; 12377937d9ceSMichael Lange PetscInt nameSize; 12387937d9ceSMichael Lange size_t len = 0; 12399852e123SBarry Smith PetscMPIInt rank, size; 12407937d9ceSMichael Lange PetscErrorCode ierr; 12417937d9ceSMichael Lange 12427937d9ceSMichael Lange PetscFunctionBegin; 12437937d9ceSMichael Lange ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 12447937d9ceSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 12459852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 12467937d9ceSMichael Lange /* Bcast name */ 12477937d9ceSMichael Lange if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} 12487937d9ceSMichael Lange nameSize = len; 12497937d9ceSMichael Lange ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 12507937d9ceSMichael Lange ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 12517937d9ceSMichael Lange if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} 12527937d9ceSMichael Lange ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 12537937d9ceSMichael Lange ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); 12547937d9ceSMichael Lange ierr = PetscFree(name);CHKERRQ(ierr); 12557937d9ceSMichael Lange /* Gather rank/index pairs of leaves into local roots to build 12567937d9ceSMichael Lange an inverse, multi-rooted SF. Note that this ignores local leaf 12577937d9ceSMichael Lange indexing due to the use of the multiSF in PetscSFGather. */ 12587937d9ceSMichael Lange ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr); 1259dc53bc9bSMatthew G. Knepley ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr); 1260dc53bc9bSMatthew G. Knepley for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1; 12617937d9ceSMichael Lange for (p = 0; p < nleaves; p++) { 1262dc53bc9bSMatthew G. Knepley leafPoints[ilocal[p]].index = ilocal[p]; 1263dc53bc9bSMatthew G. Knepley leafPoints[ilocal[p]].rank = rank; 12647937d9ceSMichael Lange } 12657937d9ceSMichael Lange ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr); 12667937d9ceSMichael Lange ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr); 12677937d9ceSMichael Lange for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p]; 12687937d9ceSMichael Lange ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr); 12697937d9ceSMichael Lange ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 12707937d9ceSMichael Lange ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 12717937d9ceSMichael Lange ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr); 12727937d9ceSMichael Lange ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 12737937d9ceSMichael Lange /* Migrate label over inverted SF to pull stratum values at leaves into roots. */ 12747937d9ceSMichael Lange ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr); 12757937d9ceSMichael Lange /* Rebuild the point strata on the receiver */ 12767937d9ceSMichael Lange for (p = 0, idx = 0; p < nroots; p++) { 12777937d9ceSMichael Lange for (d = 0; d < rootDegree[p]; d++) { 12787937d9ceSMichael Lange ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr); 12797937d9ceSMichael Lange ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr); 12807937d9ceSMichael Lange for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);} 12817937d9ceSMichael Lange } 12827937d9ceSMichael Lange idx += rootDegree[p]; 12837937d9ceSMichael Lange } 128477e0c0e7SMichael Lange ierr = PetscFree(leafPoints);CHKERRQ(ierr); 128577e0c0e7SMichael Lange ierr = PetscFree(rootStrata);CHKERRQ(ierr); 128677e0c0e7SMichael Lange ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 128777e0c0e7SMichael Lange ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr); 12887937d9ceSMichael Lange PetscFunctionReturn(0); 12897937d9ceSMichael Lange } 12907937d9ceSMichael Lange 129184f0b6dfSMatthew G. Knepley /*@ 129284f0b6dfSMatthew G. Knepley DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label 129384f0b6dfSMatthew G. Knepley 129484f0b6dfSMatthew G. Knepley Input Parameter: 129584f0b6dfSMatthew G. Knepley . label - the DMLabel 129684f0b6dfSMatthew G. Knepley 129784f0b6dfSMatthew G. Knepley Output Parameters: 129884f0b6dfSMatthew G. Knepley + section - the section giving offsets for each stratum 129984f0b6dfSMatthew G. Knepley - is - An IS containing all the label points 130084f0b6dfSMatthew G. Knepley 130184f0b6dfSMatthew G. Knepley Level: developer 130284f0b6dfSMatthew G. Knepley 130384f0b6dfSMatthew G. Knepley .seealso: DMLabelDistribute() 130484f0b6dfSMatthew G. Knepley @*/ 1305c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is) 1306c58f1c22SToby Isaac { 1307c58f1c22SToby Isaac IS vIS; 1308c58f1c22SToby Isaac const PetscInt *values; 1309c58f1c22SToby Isaac PetscInt *points; 1310c58f1c22SToby Isaac PetscInt nV, vS = 0, vE = 0, v, N; 1311c58f1c22SToby Isaac PetscErrorCode ierr; 1312c58f1c22SToby Isaac 1313c58f1c22SToby Isaac PetscFunctionBegin; 1314c58f1c22SToby Isaac ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr); 1315c58f1c22SToby Isaac ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr); 1316c58f1c22SToby Isaac ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr); 1317c58f1c22SToby Isaac if (nV) {vS = values[0]; vE = values[0]+1;} 1318c58f1c22SToby Isaac for (v = 1; v < nV; ++v) { 1319c58f1c22SToby Isaac vS = PetscMin(vS, values[v]); 1320c58f1c22SToby Isaac vE = PetscMax(vE, values[v]+1); 1321c58f1c22SToby Isaac } 1322c58f1c22SToby Isaac ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr); 1323c58f1c22SToby Isaac ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr); 1324c58f1c22SToby Isaac for (v = 0; v < nV; ++v) { 1325c58f1c22SToby Isaac PetscInt n; 1326c58f1c22SToby Isaac 1327c58f1c22SToby Isaac ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr); 1328c58f1c22SToby Isaac ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr); 1329c58f1c22SToby Isaac } 1330c58f1c22SToby Isaac ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 1331c58f1c22SToby Isaac ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr); 1332c58f1c22SToby Isaac ierr = PetscMalloc1(N, &points);CHKERRQ(ierr); 1333c58f1c22SToby Isaac for (v = 0; v < nV; ++v) { 1334c58f1c22SToby Isaac IS is; 1335c58f1c22SToby Isaac const PetscInt *spoints; 1336c58f1c22SToby Isaac PetscInt dof, off, p; 1337c58f1c22SToby Isaac 1338c58f1c22SToby Isaac ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr); 1339c58f1c22SToby Isaac ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr); 1340c58f1c22SToby Isaac ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr); 1341c58f1c22SToby Isaac ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr); 1342c58f1c22SToby Isaac for (p = 0; p < dof; ++p) points[off+p] = spoints[p]; 1343c58f1c22SToby Isaac ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr); 1344c58f1c22SToby Isaac ierr = ISDestroy(&is);CHKERRQ(ierr); 1345c58f1c22SToby Isaac } 1346c58f1c22SToby Isaac ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr); 1347c58f1c22SToby Isaac ierr = ISDestroy(&vIS);CHKERRQ(ierr); 1348c58f1c22SToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr); 1349c58f1c22SToby Isaac PetscFunctionReturn(0); 1350c58f1c22SToby Isaac } 1351c58f1c22SToby Isaac 135284f0b6dfSMatthew G. Knepley /*@ 1353c58f1c22SToby Isaac PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using 1354c58f1c22SToby Isaac the local section and an SF describing the section point overlap. 1355c58f1c22SToby Isaac 1356c58f1c22SToby Isaac Input Parameters: 1357c58f1c22SToby Isaac + s - The PetscSection for the local field layout 1358c58f1c22SToby Isaac . sf - The SF describing parallel layout of the section points 1359c58f1c22SToby Isaac . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs 1360c58f1c22SToby Isaac . label - The label specifying the points 1361c58f1c22SToby Isaac - labelValue - The label stratum specifying the points 1362c58f1c22SToby Isaac 1363c58f1c22SToby Isaac Output Parameter: 1364c58f1c22SToby Isaac . gsection - The PetscSection for the global field layout 1365c58f1c22SToby Isaac 1366c58f1c22SToby Isaac Note: This gives negative sizes and offsets to points not owned by this process 1367c58f1c22SToby Isaac 1368c58f1c22SToby Isaac Level: developer 1369c58f1c22SToby Isaac 1370c58f1c22SToby Isaac .seealso: PetscSectionCreate() 1371c58f1c22SToby Isaac @*/ 1372c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection) 1373c58f1c22SToby Isaac { 1374c58f1c22SToby Isaac PetscInt *neg = NULL, *tmpOff = NULL; 1375c58f1c22SToby Isaac PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots; 1376c58f1c22SToby Isaac PetscErrorCode ierr; 1377c58f1c22SToby Isaac 1378c58f1c22SToby Isaac PetscFunctionBegin; 1379c58f1c22SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr); 1380c58f1c22SToby Isaac ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 1381c58f1c22SToby Isaac ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr); 1382c58f1c22SToby Isaac ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 1383c58f1c22SToby Isaac if (nroots >= 0) { 1384c58f1c22SToby Isaac if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart); 1385c58f1c22SToby Isaac ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr); 1386c58f1c22SToby Isaac if (nroots > pEnd-pStart) { 1387c58f1c22SToby Isaac ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr); 1388c58f1c22SToby Isaac } else { 1389c58f1c22SToby Isaac tmpOff = &(*gsection)->atlasDof[-pStart]; 1390c58f1c22SToby Isaac } 1391c58f1c22SToby Isaac } 1392c58f1c22SToby Isaac /* Mark ghost points with negative dof */ 1393c58f1c22SToby Isaac for (p = pStart; p < pEnd; ++p) { 1394c58f1c22SToby Isaac PetscInt value; 1395c58f1c22SToby Isaac 1396c58f1c22SToby Isaac ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr); 1397c58f1c22SToby Isaac if (value != labelValue) continue; 1398c58f1c22SToby Isaac ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr); 1399c58f1c22SToby Isaac ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr); 1400c58f1c22SToby Isaac ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 1401c58f1c22SToby Isaac if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);} 1402c58f1c22SToby Isaac if (neg) neg[p] = -(dof+1); 1403c58f1c22SToby Isaac } 1404c58f1c22SToby Isaac ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr); 1405c58f1c22SToby Isaac if (nroots >= 0) { 1406c58f1c22SToby Isaac ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1407c58f1c22SToby Isaac ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1408c58f1c22SToby Isaac if (nroots > pEnd-pStart) { 1409c58f1c22SToby Isaac for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];} 1410c58f1c22SToby Isaac } 1411c58f1c22SToby Isaac } 1412c58f1c22SToby Isaac /* Calculate new sizes, get proccess offset, and calculate point offsets */ 1413c58f1c22SToby Isaac for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1414c58f1c22SToby Isaac cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0; 1415c58f1c22SToby Isaac (*gsection)->atlasOff[p] = off; 1416c58f1c22SToby Isaac off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0; 1417c58f1c22SToby Isaac } 1418c58f1c22SToby Isaac ierr = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr); 1419c58f1c22SToby Isaac globalOff -= off; 1420c58f1c22SToby Isaac for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1421c58f1c22SToby Isaac (*gsection)->atlasOff[p] += globalOff; 1422c58f1c22SToby Isaac if (neg) neg[p] = -((*gsection)->atlasOff[p]+1); 1423c58f1c22SToby Isaac } 1424c58f1c22SToby Isaac /* Put in negative offsets for ghost points */ 1425c58f1c22SToby Isaac if (nroots >= 0) { 1426c58f1c22SToby Isaac ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1427c58f1c22SToby Isaac ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1428c58f1c22SToby Isaac if (nroots > pEnd-pStart) { 1429c58f1c22SToby Isaac for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];} 1430c58f1c22SToby Isaac } 1431c58f1c22SToby Isaac } 1432c58f1c22SToby Isaac if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);} 1433c58f1c22SToby Isaac ierr = PetscFree(neg);CHKERRQ(ierr); 1434c58f1c22SToby Isaac PetscFunctionReturn(0); 1435c58f1c22SToby Isaac } 1436c58f1c22SToby Isaac 14375fdea053SToby Isaac typedef struct _n_PetscSectionSym_Label 14385fdea053SToby Isaac { 14395fdea053SToby Isaac DMLabel label; 14405fdea053SToby Isaac PetscCopyMode *modes; 14415fdea053SToby Isaac PetscInt *sizes; 14425fdea053SToby Isaac const PetscInt ***perms; 14435fdea053SToby Isaac const PetscScalar ***rots; 14445fdea053SToby Isaac PetscInt (*minMaxOrients)[2]; 14455fdea053SToby Isaac PetscInt numStrata; /* numStrata is only increasing, functions as a state */ 14465fdea053SToby Isaac } PetscSectionSym_Label; 14475fdea053SToby Isaac 14485fdea053SToby Isaac static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym) 14495fdea053SToby Isaac { 14505fdea053SToby Isaac PetscInt i, j; 14515fdea053SToby Isaac PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data; 14525fdea053SToby Isaac PetscErrorCode ierr; 14535fdea053SToby Isaac 14545fdea053SToby Isaac PetscFunctionBegin; 14555fdea053SToby Isaac for (i = 0; i <= sl->numStrata; i++) { 14565fdea053SToby Isaac if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) { 14575fdea053SToby Isaac for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 14585fdea053SToby Isaac if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);} 14595fdea053SToby Isaac if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);} 14605fdea053SToby Isaac } 14615fdea053SToby Isaac if (sl->perms[i]) { 14625fdea053SToby Isaac const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]]; 14635fdea053SToby Isaac 14645fdea053SToby Isaac ierr = PetscFree(perms);CHKERRQ(ierr); 14655fdea053SToby Isaac } 14665fdea053SToby Isaac if (sl->rots[i]) { 14675fdea053SToby Isaac const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]]; 14685fdea053SToby Isaac 14695fdea053SToby Isaac ierr = PetscFree(rots);CHKERRQ(ierr); 14705fdea053SToby Isaac } 14715fdea053SToby Isaac } 14725fdea053SToby Isaac } 14735fdea053SToby Isaac ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr); 14745fdea053SToby Isaac ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr); 14755fdea053SToby Isaac sl->numStrata = 0; 14765fdea053SToby Isaac PetscFunctionReturn(0); 14775fdea053SToby Isaac } 14785fdea053SToby Isaac 14795fdea053SToby Isaac static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym) 14805fdea053SToby Isaac { 14815fdea053SToby Isaac PetscErrorCode ierr; 14825fdea053SToby Isaac 14835fdea053SToby Isaac PetscFunctionBegin; 14845fdea053SToby Isaac ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr); 14855fdea053SToby Isaac ierr = PetscFree(sym->data);CHKERRQ(ierr); 14865fdea053SToby Isaac PetscFunctionReturn(0); 14875fdea053SToby Isaac } 14885fdea053SToby Isaac 14895fdea053SToby Isaac static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer) 14905fdea053SToby Isaac { 14915fdea053SToby Isaac PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data; 14925fdea053SToby Isaac PetscBool isAscii; 14935fdea053SToby Isaac DMLabel label = sl->label; 14945fdea053SToby Isaac PetscErrorCode ierr; 14955fdea053SToby Isaac 14965fdea053SToby Isaac PetscFunctionBegin; 14975fdea053SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr); 14985fdea053SToby Isaac if (isAscii) { 14995fdea053SToby Isaac PetscInt i, j, k; 15005fdea053SToby Isaac PetscViewerFormat format; 15015fdea053SToby Isaac 15025fdea053SToby Isaac ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 15035fdea053SToby Isaac if (label) { 15045fdea053SToby Isaac ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 15055fdea053SToby Isaac if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 15065fdea053SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 15075fdea053SToby Isaac ierr = DMLabelView(label, viewer);CHKERRQ(ierr); 15085fdea053SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 15095fdea053SToby Isaac } else { 15105fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer," Label '%s'\n",sl->label->name);CHKERRQ(ierr); 15115fdea053SToby Isaac } 15125fdea053SToby Isaac } else { 15135fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr); 15145fdea053SToby Isaac } 15155fdea053SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 15165fdea053SToby Isaac for (i = 0; i <= sl->numStrata; i++) { 15175fdea053SToby Isaac PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue; 15185fdea053SToby Isaac 15195fdea053SToby Isaac if (!(sl->perms[i] || sl->rots[i])) { 15205fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr); 15215fdea053SToby Isaac } else { 15225fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr); 15235fdea053SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 15245fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr); 15255fdea053SToby Isaac if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 15265fdea053SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 15275fdea053SToby Isaac for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 15285fdea053SToby Isaac if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) { 15295fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr); 15305fdea053SToby Isaac } else { 15315fdea053SToby Isaac PetscInt tab; 15325fdea053SToby Isaac 15335fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr); 15345fdea053SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 15355fdea053SToby Isaac ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr); 15365fdea053SToby Isaac if (sl->perms[i] && sl->perms[i][j]) { 15375fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr); 15385fdea053SToby Isaac ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr); 15395fdea053SToby Isaac for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);} 15405fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 15415fdea053SToby Isaac ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr); 15425fdea053SToby Isaac } 15435fdea053SToby Isaac if (sl->rots[i] && sl->rots[i][j]) { 15445fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Rotations: ");CHKERRQ(ierr); 15455fdea053SToby Isaac ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr); 15465fdea053SToby Isaac #if defined(PETSC_USE_COMPLEX) 15475fdea053SToby 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);} 15485fdea053SToby Isaac #else 15495fdea053SToby Isaac for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);} 15505fdea053SToby Isaac #endif 15515fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 15525fdea053SToby Isaac ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr); 15535fdea053SToby Isaac } 15545fdea053SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 15555fdea053SToby Isaac } 15565fdea053SToby Isaac } 15575fdea053SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 15585fdea053SToby Isaac } 15595fdea053SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 15605fdea053SToby Isaac } 15615fdea053SToby Isaac } 15625fdea053SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 15635fdea053SToby Isaac } 15645fdea053SToby Isaac PetscFunctionReturn(0); 15655fdea053SToby Isaac } 15665fdea053SToby Isaac 15675fdea053SToby Isaac /*@ 15685fdea053SToby Isaac PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries 15695fdea053SToby Isaac 15705fdea053SToby Isaac Logically collective on sym 15715fdea053SToby Isaac 15725fdea053SToby Isaac Input parameters: 15735fdea053SToby Isaac + sym - the section symmetries 15745fdea053SToby Isaac - label - the DMLabel describing the types of points 15755fdea053SToby Isaac 15765fdea053SToby Isaac Level: developer: 15775fdea053SToby Isaac 15785fdea053SToby Isaac .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms() 15795fdea053SToby Isaac @*/ 15805fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label) 15815fdea053SToby Isaac { 15825fdea053SToby Isaac PetscSectionSym_Label *sl; 15835fdea053SToby Isaac PetscErrorCode ierr; 15845fdea053SToby Isaac 15855fdea053SToby Isaac PetscFunctionBegin; 15865fdea053SToby Isaac PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1); 15875fdea053SToby Isaac sl = (PetscSectionSym_Label *) sym->data; 15885fdea053SToby Isaac if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);} 15895fdea053SToby Isaac if (label) { 15905fdea053SToby Isaac label->refct++; 15915fdea053SToby Isaac sl->label = label; 15925fdea053SToby Isaac ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr); 15931a834cf9SToby 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); 15941a834cf9SToby Isaac ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr); 15951a834cf9SToby Isaac ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr); 15961a834cf9SToby Isaac ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr); 15971a834cf9SToby Isaac ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr); 15981a834cf9SToby Isaac ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr); 15995fdea053SToby Isaac } 16005fdea053SToby Isaac PetscFunctionReturn(0); 16015fdea053SToby Isaac } 16025fdea053SToby Isaac 16035fdea053SToby Isaac /*@C 16045fdea053SToby Isaac PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum 16055fdea053SToby Isaac 16065fdea053SToby Isaac Logically collective on PetscSectionSym 16075fdea053SToby Isaac 16085fdea053SToby Isaac InputParameters: 16095fdea053SToby Isaac + sys - the section symmetries 16105fdea053SToby Isaac . stratum - the stratum value in the label that we are assigning symmetries for 16115fdea053SToby Isaac . size - the number of dofs for points in the stratum of the label 16125fdea053SToby Isaac . minOrient - the smallest orientation for a point in this stratum 16135fdea053SToby Isaac . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient)) 16145fdea053SToby Isaac . mode - how sym should copy the perms and rots arrays 16155fdea053SToby Isaac . perms - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation. A NULL permutation is the identity 16165fdea053SToby 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 16175fdea053SToby Isaac 16185fdea053SToby Isaac Level: developer 16195fdea053SToby Isaac 16205fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel() 16215fdea053SToby Isaac @*/ 16225fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots) 16235fdea053SToby Isaac { 16245fdea053SToby Isaac PetscInt i, j, k; 16255fdea053SToby Isaac PetscSectionSym_Label *sl; 16265fdea053SToby Isaac PetscErrorCode ierr; 16275fdea053SToby Isaac 16285fdea053SToby Isaac PetscFunctionBegin; 16295fdea053SToby Isaac PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1); 16305fdea053SToby Isaac sl = (PetscSectionSym_Label *) sym->data; 16315fdea053SToby Isaac if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet"); 16325fdea053SToby Isaac for (i = 0; i <= sl->numStrata; i++) { 16335fdea053SToby Isaac PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue; 16345fdea053SToby Isaac 16355fdea053SToby Isaac if (stratum == value) break; 16365fdea053SToby Isaac } 16375fdea053SToby Isaac if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,sl->label->name); 16385fdea053SToby Isaac sl->sizes[i] = size; 16395fdea053SToby Isaac sl->modes[i] = mode; 16405fdea053SToby Isaac sl->minMaxOrients[i][0] = minOrient; 16415fdea053SToby Isaac sl->minMaxOrients[i][1] = maxOrient; 16425fdea053SToby Isaac if (mode == PETSC_COPY_VALUES) { 16435fdea053SToby Isaac if (perms) { 16445fdea053SToby Isaac PetscInt **ownPerms; 16455fdea053SToby Isaac 16465fdea053SToby Isaac ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr); 16475fdea053SToby Isaac for (j = 0; j < maxOrient-minOrient; j++) { 16485fdea053SToby Isaac if (perms[j]) { 16495fdea053SToby Isaac ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr); 16505fdea053SToby Isaac for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];} 16515fdea053SToby Isaac } 16525fdea053SToby Isaac } 16535fdea053SToby Isaac sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient]; 16545fdea053SToby Isaac } 16555fdea053SToby Isaac if (rots) { 16565fdea053SToby Isaac PetscScalar **ownRots; 16575fdea053SToby Isaac 16585fdea053SToby Isaac ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr); 16595fdea053SToby Isaac for (j = 0; j < maxOrient-minOrient; j++) { 16605fdea053SToby Isaac if (rots[j]) { 16615fdea053SToby Isaac ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr); 16625fdea053SToby Isaac for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];} 16635fdea053SToby Isaac } 16645fdea053SToby Isaac } 16655fdea053SToby Isaac sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient]; 16665fdea053SToby Isaac } 16675fdea053SToby Isaac } else { 16685fdea053SToby Isaac sl->perms[i] = perms ? &perms[-minOrient] : NULL; 16695fdea053SToby Isaac sl->rots[i] = rots ? &rots[-minOrient] : NULL; 16705fdea053SToby Isaac } 16715fdea053SToby Isaac PetscFunctionReturn(0); 16725fdea053SToby Isaac } 16735fdea053SToby Isaac 16745fdea053SToby Isaac static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots) 16755fdea053SToby Isaac { 16765fdea053SToby Isaac PetscInt i, j, numStrata; 16775fdea053SToby Isaac PetscSectionSym_Label *sl; 16785fdea053SToby Isaac DMLabel label; 16795fdea053SToby Isaac PetscErrorCode ierr; 16805fdea053SToby Isaac 16815fdea053SToby Isaac PetscFunctionBegin; 16825fdea053SToby Isaac sl = (PetscSectionSym_Label *) sym->data; 16835fdea053SToby Isaac numStrata = sl->numStrata; 16845fdea053SToby Isaac label = sl->label; 16855fdea053SToby Isaac for (i = 0; i < numPoints; i++) { 16865fdea053SToby Isaac PetscInt point = points[2*i]; 16875fdea053SToby Isaac PetscInt ornt = points[2*i+1]; 16885fdea053SToby Isaac 16895fdea053SToby Isaac for (j = 0; j < numStrata; j++) { 16905fdea053SToby Isaac if (label->validIS[j]) { 16915fdea053SToby Isaac PetscInt k; 16925fdea053SToby Isaac 16935fdea053SToby Isaac ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr); 16945fdea053SToby Isaac if (k >= 0) break; 16955fdea053SToby Isaac } else { 16965fdea053SToby Isaac PetscBool has; 16975fdea053SToby Isaac 1698*e8f14785SLisandro Dalcin PetscHSetIHas(label->ht[j], point, &has); 16995fdea053SToby Isaac if (has) break; 17005fdea053SToby Isaac } 17015fdea053SToby Isaac } 17025fdea053SToby 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); 17035fdea053SToby Isaac if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;} 17045fdea053SToby Isaac if (rots) {rots[i] = sl->rots[j] ? sl->rots[j][ornt] : NULL;} 17055fdea053SToby Isaac } 17065fdea053SToby Isaac PetscFunctionReturn(0); 17075fdea053SToby Isaac } 17085fdea053SToby Isaac 17095fdea053SToby Isaac PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym) 17105fdea053SToby Isaac { 17115fdea053SToby Isaac PetscSectionSym_Label *sl; 17125fdea053SToby Isaac PetscErrorCode ierr; 17135fdea053SToby Isaac 17145fdea053SToby Isaac PetscFunctionBegin; 17155fdea053SToby Isaac ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr); 17165fdea053SToby Isaac sym->ops->getpoints = PetscSectionSymGetPoints_Label; 17175fdea053SToby Isaac sym->ops->view = PetscSectionSymView_Label; 17185fdea053SToby Isaac sym->ops->destroy = PetscSectionSymDestroy_Label; 17195fdea053SToby Isaac sym->data = (void *) sl; 17205fdea053SToby Isaac PetscFunctionReturn(0); 17215fdea053SToby Isaac } 17225fdea053SToby Isaac 17235fdea053SToby Isaac /*@ 17245fdea053SToby Isaac PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label 17255fdea053SToby Isaac 17265fdea053SToby Isaac Collective 17275fdea053SToby Isaac 17285fdea053SToby Isaac Input Parameters: 17295fdea053SToby Isaac + comm - the MPI communicator for the new symmetry 17305fdea053SToby Isaac - label - the label defining the strata 17315fdea053SToby Isaac 17325fdea053SToby Isaac Output Parameters: 17335fdea053SToby Isaac . sym - the section symmetries 17345fdea053SToby Isaac 17355fdea053SToby Isaac Level: developer 17365fdea053SToby Isaac 17375fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms() 17385fdea053SToby Isaac @*/ 17395fdea053SToby Isaac PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym) 17405fdea053SToby Isaac { 17415fdea053SToby Isaac PetscErrorCode ierr; 17425fdea053SToby Isaac 17435fdea053SToby Isaac PetscFunctionBegin; 17445fdea053SToby Isaac ierr = DMInitializePackage();CHKERRQ(ierr); 17455fdea053SToby Isaac ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr); 17465fdea053SToby Isaac ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr); 17475fdea053SToby Isaac ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr); 17485fdea053SToby Isaac PetscFunctionReturn(0); 17495fdea053SToby Isaac } 1750