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