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 95b5e7992SMatthew G. Knepley Collective 105b5e7992SMatthew G. Knepley 11d67d17b1SMatthew G. Knepley Input parameters: 12d67d17b1SMatthew G. Knepley + comm - The communicator, usually PETSC_COMM_SELF 13d67d17b1SMatthew G. Knepley - name - The label name 14c58f1c22SToby Isaac 15c58f1c22SToby Isaac Output parameter: 16c58f1c22SToby Isaac . label - The DMLabel 17c58f1c22SToby Isaac 18c58f1c22SToby Isaac Level: beginner 19c58f1c22SToby Isaac 20c58f1c22SToby Isaac .seealso: DMLabelDestroy() 21c58f1c22SToby Isaac @*/ 22d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelCreate(MPI_Comm comm, const char name[], DMLabel *label) 23c58f1c22SToby Isaac { 24c58f1c22SToby Isaac PetscErrorCode ierr; 25c58f1c22SToby Isaac 26c58f1c22SToby Isaac PetscFunctionBegin; 27d67d17b1SMatthew G. Knepley PetscValidPointer(label,2); 28d67d17b1SMatthew G. Knepley ierr = DMInitializePackage();CHKERRQ(ierr); 29c58f1c22SToby Isaac 30d67d17b1SMatthew G. Knepley ierr = PetscHeaderCreate(*label,DMLABEL_CLASSID,"DMLabel","DMLabel","DM",comm,DMLabelDestroy,DMLabelView);CHKERRQ(ierr); 31d67d17b1SMatthew G. Knepley 32c58f1c22SToby Isaac (*label)->numStrata = 0; 335aa44df4SToby Isaac (*label)->defaultValue = -1; 34c58f1c22SToby Isaac (*label)->stratumValues = NULL; 35ad8374ffSToby Isaac (*label)->validIS = NULL; 36c58f1c22SToby Isaac (*label)->stratumSizes = NULL; 37c58f1c22SToby Isaac (*label)->points = NULL; 38c58f1c22SToby Isaac (*label)->ht = NULL; 39c58f1c22SToby Isaac (*label)->pStart = -1; 40c58f1c22SToby Isaac (*label)->pEnd = -1; 41c58f1c22SToby Isaac (*label)->bt = NULL; 42b9907514SLisandro Dalcin ierr = PetscHMapICreate(&(*label)->hmap);CHKERRQ(ierr); 43d67d17b1SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *label, name);CHKERRQ(ierr); 44c58f1c22SToby Isaac PetscFunctionReturn(0); 45c58f1c22SToby Isaac } 46c58f1c22SToby Isaac 47c58f1c22SToby Isaac /* 48c58f1c22SToby Isaac DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format 49c58f1c22SToby Isaac 505b5e7992SMatthew G. Knepley Not collective 515b5e7992SMatthew G. Knepley 52c58f1c22SToby Isaac Input parameter: 53c58f1c22SToby Isaac + label - The DMLabel 54c58f1c22SToby Isaac - v - The stratum value 55c58f1c22SToby Isaac 56c58f1c22SToby Isaac Output parameter: 57c58f1c22SToby Isaac . label - The DMLabel with stratum in sorted list format 58c58f1c22SToby Isaac 59c58f1c22SToby Isaac Level: developer 60c58f1c22SToby Isaac 61c58f1c22SToby Isaac .seealso: DMLabelCreate() 62c58f1c22SToby Isaac */ 63c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v) 64c58f1c22SToby Isaac { 65277ea44aSLisandro Dalcin IS is; 66b9907514SLisandro Dalcin PetscInt off = 0, *pointArray, p; 67c58f1c22SToby Isaac PetscErrorCode ierr; 68c58f1c22SToby Isaac 69b9907514SLisandro Dalcin if (PetscLikely(v >= 0 && v < label->numStrata) && label->validIS[v]) return 0; 70c58f1c22SToby Isaac PetscFunctionBegin; 710c3c4a36SLisandro 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); 72e8f14785SLisandro Dalcin ierr = PetscHSetIGetSize(label->ht[v], &label->stratumSizes[v]);CHKERRQ(ierr); 73ad8374ffSToby Isaac ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr); 74e8f14785SLisandro Dalcin ierr = PetscHSetIGetElems(label->ht[v], &off, pointArray);CHKERRQ(ierr); 75b9907514SLisandro Dalcin ierr = PetscHSetIClear(label->ht[v]);CHKERRQ(ierr); 76ad8374ffSToby Isaac ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr); 77c58f1c22SToby Isaac if (label->bt) { 78c58f1c22SToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) { 79ad8374ffSToby Isaac const PetscInt point = pointArray[p]; 80c58f1c22SToby Isaac if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, label->pStart, label->pEnd); 81c58f1c22SToby Isaac ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 82c58f1c22SToby Isaac } 83c58f1c22SToby Isaac } 84277ea44aSLisandro Dalcin ierr = ISCreateGeneral(PETSC_COMM_SELF, label->stratumSizes[v], pointArray, PETSC_OWN_POINTER, &is);CHKERRQ(ierr); 85277ea44aSLisandro Dalcin ierr = PetscObjectSetName((PetscObject) is, "indices");CHKERRQ(ierr); 86277ea44aSLisandro Dalcin label->points[v] = is; 87ad8374ffSToby Isaac label->validIS[v] = PETSC_TRUE; 88d67d17b1SMatthew G. Knepley ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr); 89c58f1c22SToby Isaac PetscFunctionReturn(0); 90c58f1c22SToby Isaac } 91c58f1c22SToby Isaac 92c58f1c22SToby Isaac /* 93c58f1c22SToby Isaac DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format 94c58f1c22SToby Isaac 955b5e7992SMatthew G. Knepley Not collective 965b5e7992SMatthew G. Knepley 97c58f1c22SToby Isaac Input parameter: 98c58f1c22SToby Isaac . label - The DMLabel 99c58f1c22SToby Isaac 100c58f1c22SToby Isaac Output parameter: 101c58f1c22SToby Isaac . label - The DMLabel with all strata in sorted list format 102c58f1c22SToby Isaac 103c58f1c22SToby Isaac Level: developer 104c58f1c22SToby Isaac 105c58f1c22SToby Isaac .seealso: DMLabelCreate() 106c58f1c22SToby Isaac */ 107c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label) 108c58f1c22SToby Isaac { 109c58f1c22SToby Isaac PetscInt v; 110c58f1c22SToby Isaac PetscErrorCode ierr; 111c58f1c22SToby Isaac 112c58f1c22SToby Isaac PetscFunctionBegin; 113c58f1c22SToby Isaac for (v = 0; v < label->numStrata; v++){ 114c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 115c58f1c22SToby Isaac } 116c58f1c22SToby Isaac PetscFunctionReturn(0); 117c58f1c22SToby Isaac } 118c58f1c22SToby Isaac 119c58f1c22SToby Isaac /* 120c58f1c22SToby Isaac DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format 121c58f1c22SToby Isaac 1225b5e7992SMatthew G. Knepley Not collective 1235b5e7992SMatthew G. Knepley 124c58f1c22SToby Isaac Input parameter: 125c58f1c22SToby Isaac + label - The DMLabel 126c58f1c22SToby Isaac - v - The stratum value 127c58f1c22SToby Isaac 128c58f1c22SToby Isaac Output parameter: 129c58f1c22SToby Isaac . label - The DMLabel with stratum in hash format 130c58f1c22SToby Isaac 131c58f1c22SToby Isaac Level: developer 132c58f1c22SToby Isaac 133c58f1c22SToby Isaac .seealso: DMLabelCreate() 134c58f1c22SToby Isaac */ 135c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v) 136c58f1c22SToby Isaac { 137c58f1c22SToby Isaac PetscInt p; 138ad8374ffSToby Isaac const PetscInt *points; 139c58f1c22SToby Isaac PetscErrorCode ierr; 140c58f1c22SToby Isaac 141b9907514SLisandro Dalcin if (PetscLikely(v >= 0 && v < label->numStrata) && !label->validIS[v]) return 0; 142c58f1c22SToby Isaac PetscFunctionBegin; 1430c3c4a36SLisandro 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); 144ad8374ffSToby Isaac if (label->points[v]) { 145ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 146e8f14785SLisandro Dalcin for (p = 0; p < label->stratumSizes[v]; ++p) { 147e8f14785SLisandro Dalcin ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr); 148e8f14785SLisandro Dalcin } 149ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 150ad8374ffSToby Isaac ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 151ad8374ffSToby Isaac } 152ad8374ffSToby Isaac label->validIS[v] = PETSC_FALSE; 153c58f1c22SToby Isaac PetscFunctionReturn(0); 154c58f1c22SToby Isaac } 155c58f1c22SToby Isaac 156b9907514SLisandro Dalcin #if !defined(DMLABEL_LOOKUP_THRESHOLD) 157b9907514SLisandro Dalcin #define DMLABEL_LOOKUP_THRESHOLD 16 158b9907514SLisandro Dalcin #endif 159b9907514SLisandro Dalcin 1600c3c4a36SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index) 1610c3c4a36SLisandro Dalcin { 1620c3c4a36SLisandro Dalcin PetscInt v; 163b9907514SLisandro Dalcin PetscErrorCode ierr; 1640e79e033SBarry Smith 1650c3c4a36SLisandro Dalcin PetscFunctionBegin; 1660e79e033SBarry Smith *index = -1; 167b9907514SLisandro Dalcin if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) { 168b9907514SLisandro Dalcin for (v = 0; v < label->numStrata; ++v) 169b9907514SLisandro Dalcin if (label->stratumValues[v] == value) {*index = v; break;} 170b9907514SLisandro Dalcin } else { 171b9907514SLisandro Dalcin ierr = PetscHMapIGet(label->hmap, value, index);CHKERRQ(ierr); 1720e79e033SBarry Smith } 17390e9b2aeSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 17490e9b2aeSLisandro Dalcin { /* Check strata hash map consistency */ 17590e9b2aeSLisandro Dalcin PetscInt len, loc = -1; 17690e9b2aeSLisandro Dalcin ierr = PetscHMapIGetSize(label->hmap, &len);CHKERRQ(ierr); 17790e9b2aeSLisandro Dalcin if (len != label->numStrata) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map size"); 17890e9b2aeSLisandro Dalcin if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) { 17990e9b2aeSLisandro Dalcin ierr = PetscHMapIGet(label->hmap, value, &loc);CHKERRQ(ierr); 18090e9b2aeSLisandro Dalcin } else { 18190e9b2aeSLisandro Dalcin for (v = 0; v < label->numStrata; ++v) 18290e9b2aeSLisandro Dalcin if (label->stratumValues[v] == value) {loc = v; break;} 18390e9b2aeSLisandro Dalcin } 18490e9b2aeSLisandro Dalcin if (loc != *index) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map lookup"); 18590e9b2aeSLisandro Dalcin } 18690e9b2aeSLisandro Dalcin #endif 1870c3c4a36SLisandro Dalcin PetscFunctionReturn(0); 1880c3c4a36SLisandro Dalcin } 1890c3c4a36SLisandro Dalcin 190b9907514SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index) 191c58f1c22SToby Isaac { 192b9907514SLisandro Dalcin PetscInt v; 193b9907514SLisandro Dalcin PetscInt *tmpV; 194b9907514SLisandro Dalcin PetscInt *tmpS; 195b9907514SLisandro Dalcin PetscHSetI *tmpH, ht; 196b9907514SLisandro Dalcin IS *tmpP, is; 197c58f1c22SToby Isaac PetscBool *tmpB; 198b9907514SLisandro Dalcin PetscHMapI hmap = label->hmap; 199c58f1c22SToby Isaac PetscErrorCode ierr; 200c58f1c22SToby Isaac 201c58f1c22SToby Isaac PetscFunctionBegin; 202b9907514SLisandro Dalcin v = label->numStrata; 203b9907514SLisandro Dalcin tmpV = label->stratumValues; 204b9907514SLisandro Dalcin tmpS = label->stratumSizes; 205b9907514SLisandro Dalcin tmpH = label->ht; 206b9907514SLisandro Dalcin tmpP = label->points; 207b9907514SLisandro Dalcin tmpB = label->validIS; 2088e1f8cf0SLisandro Dalcin { /* TODO: PetscRealloc() is broken, use malloc+memcpy+free */ 2098e1f8cf0SLisandro Dalcin PetscInt *oldV = tmpV; 2108e1f8cf0SLisandro Dalcin PetscInt *oldS = tmpS; 2118e1f8cf0SLisandro Dalcin PetscHSetI *oldH = tmpH; 2128e1f8cf0SLisandro Dalcin IS *oldP = tmpP; 2138e1f8cf0SLisandro Dalcin PetscBool *oldB = tmpB; 2148e1f8cf0SLisandro Dalcin ierr = PetscMalloc((v+1)*sizeof(*tmpV), &tmpV);CHKERRQ(ierr); 2158e1f8cf0SLisandro Dalcin ierr = PetscMalloc((v+1)*sizeof(*tmpS), &tmpS);CHKERRQ(ierr); 2168e1f8cf0SLisandro Dalcin ierr = PetscMalloc((v+1)*sizeof(*tmpH), &tmpH);CHKERRQ(ierr); 2178e1f8cf0SLisandro Dalcin ierr = PetscMalloc((v+1)*sizeof(*tmpP), &tmpP);CHKERRQ(ierr); 2188e1f8cf0SLisandro Dalcin ierr = PetscMalloc((v+1)*sizeof(*tmpB), &tmpB);CHKERRQ(ierr); 219580bdb30SBarry Smith ierr = PetscArraycpy(tmpV, oldV, v);CHKERRQ(ierr); 220580bdb30SBarry Smith ierr = PetscArraycpy(tmpS, oldS, v);CHKERRQ(ierr); 221580bdb30SBarry Smith ierr = PetscArraycpy(tmpH, oldH, v);CHKERRQ(ierr); 222580bdb30SBarry Smith ierr = PetscArraycpy(tmpP, oldP, v);CHKERRQ(ierr); 223580bdb30SBarry Smith ierr = PetscArraycpy(tmpB, oldB, v);CHKERRQ(ierr); 2248e1f8cf0SLisandro Dalcin ierr = PetscFree(oldV);CHKERRQ(ierr); 2258e1f8cf0SLisandro Dalcin ierr = PetscFree(oldS);CHKERRQ(ierr); 2268e1f8cf0SLisandro Dalcin ierr = PetscFree(oldH);CHKERRQ(ierr); 2278e1f8cf0SLisandro Dalcin ierr = PetscFree(oldP);CHKERRQ(ierr); 2288e1f8cf0SLisandro Dalcin ierr = PetscFree(oldB);CHKERRQ(ierr); 2298e1f8cf0SLisandro Dalcin } 230b9907514SLisandro Dalcin label->numStrata = v+1; 231c58f1c22SToby Isaac label->stratumValues = tmpV; 232c58f1c22SToby Isaac label->stratumSizes = tmpS; 233c58f1c22SToby Isaac label->ht = tmpH; 234c58f1c22SToby Isaac label->points = tmpP; 235ad8374ffSToby Isaac label->validIS = tmpB; 236b9907514SLisandro Dalcin ierr = PetscHSetICreate(&ht);CHKERRQ(ierr); 237b9907514SLisandro Dalcin ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);CHKERRQ(ierr); 238b9907514SLisandro Dalcin ierr = PetscHMapISet(hmap, value, v);CHKERRQ(ierr); 239b9907514SLisandro Dalcin tmpV[v] = value; 240b9907514SLisandro Dalcin tmpS[v] = 0; 241b9907514SLisandro Dalcin tmpH[v] = ht; 242b9907514SLisandro Dalcin tmpP[v] = is; 243b9907514SLisandro Dalcin tmpB[v] = PETSC_TRUE; 244277ea44aSLisandro Dalcin ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr); 2450c3c4a36SLisandro Dalcin *index = v; 2460c3c4a36SLisandro Dalcin PetscFunctionReturn(0); 2470c3c4a36SLisandro Dalcin } 2480c3c4a36SLisandro Dalcin 249b9907514SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupAddStratum(DMLabel label, PetscInt value, PetscInt *index) 250b9907514SLisandro Dalcin { 251b9907514SLisandro Dalcin PetscErrorCode ierr; 252b9907514SLisandro Dalcin PetscFunctionBegin; 253b9907514SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, index);CHKERRQ(ierr); 254b9907514SLisandro Dalcin if (*index < 0) {ierr = DMLabelNewStratum(label, value, index);CHKERRQ(ierr);} 255b9907514SLisandro Dalcin PetscFunctionReturn(0); 256b9907514SLisandro Dalcin } 257b9907514SLisandro Dalcin 258b9907514SLisandro Dalcin /*@ 259b9907514SLisandro Dalcin DMLabelAddStratum - Adds a new stratum value in a DMLabel 260b9907514SLisandro Dalcin 261b9907514SLisandro Dalcin Input Parameter: 262b9907514SLisandro Dalcin + label - The DMLabel 263b9907514SLisandro Dalcin - value - The stratum value 264b9907514SLisandro Dalcin 265b9907514SLisandro Dalcin Level: beginner 266b9907514SLisandro Dalcin 267b9907514SLisandro Dalcin .seealso: DMLabelCreate(), DMLabelDestroy() 268b9907514SLisandro Dalcin @*/ 2690c3c4a36SLisandro Dalcin PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value) 2700c3c4a36SLisandro Dalcin { 2710c3c4a36SLisandro Dalcin PetscInt v; 2720c3c4a36SLisandro Dalcin PetscErrorCode ierr; 2730c3c4a36SLisandro Dalcin 2740c3c4a36SLisandro Dalcin PetscFunctionBegin; 275d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 276b9907514SLisandro Dalcin ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr); 277b9907514SLisandro Dalcin PetscFunctionReturn(0); 278b9907514SLisandro Dalcin } 279b9907514SLisandro Dalcin 280b9907514SLisandro Dalcin /*@ 281b9907514SLisandro Dalcin DMLabelAddStrata - Adds new stratum values in a DMLabel 282b9907514SLisandro Dalcin 2835b5e7992SMatthew G. Knepley Not collective 2845b5e7992SMatthew G. Knepley 285b9907514SLisandro Dalcin Input Parameter: 286b9907514SLisandro Dalcin + label - The DMLabel 287b9907514SLisandro Dalcin . numStrata - The number of stratum values 288b9907514SLisandro Dalcin - stratumValues - The stratum values 289b9907514SLisandro Dalcin 290b9907514SLisandro Dalcin Level: beginner 291b9907514SLisandro Dalcin 292b9907514SLisandro Dalcin .seealso: DMLabelCreate(), DMLabelDestroy() 293b9907514SLisandro Dalcin @*/ 294b9907514SLisandro Dalcin PetscErrorCode DMLabelAddStrata(DMLabel label, PetscInt numStrata, const PetscInt stratumValues[]) 295b9907514SLisandro Dalcin { 296b9907514SLisandro Dalcin PetscInt *values, v; 297b9907514SLisandro Dalcin PetscErrorCode ierr; 298b9907514SLisandro Dalcin 299b9907514SLisandro Dalcin PetscFunctionBegin; 300b9907514SLisandro Dalcin PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 301b9907514SLisandro Dalcin if (numStrata) PetscValidIntPointer(stratumValues, 3); 302b9907514SLisandro Dalcin ierr = PetscMalloc1(numStrata, &values);CHKERRQ(ierr); 303580bdb30SBarry Smith ierr = PetscArraycpy(values, stratumValues, numStrata);CHKERRQ(ierr); 304b9907514SLisandro Dalcin ierr = PetscSortRemoveDupsInt(&numStrata, values);CHKERRQ(ierr); 305b9907514SLisandro Dalcin if (!label->numStrata) { /* Fast preallocation */ 306b9907514SLisandro Dalcin PetscInt *tmpV; 307b9907514SLisandro Dalcin PetscInt *tmpS; 308b9907514SLisandro Dalcin PetscHSetI *tmpH, ht; 309b9907514SLisandro Dalcin IS *tmpP, is; 310b9907514SLisandro Dalcin PetscBool *tmpB; 311b9907514SLisandro Dalcin PetscHMapI hmap = label->hmap; 312b9907514SLisandro Dalcin 313b9907514SLisandro Dalcin ierr = PetscMalloc1(numStrata, &tmpV);CHKERRQ(ierr); 314b9907514SLisandro Dalcin ierr = PetscMalloc1(numStrata, &tmpS);CHKERRQ(ierr); 315b9907514SLisandro Dalcin ierr = PetscMalloc1(numStrata, &tmpH);CHKERRQ(ierr); 316b9907514SLisandro Dalcin ierr = PetscMalloc1(numStrata, &tmpP);CHKERRQ(ierr); 317b9907514SLisandro Dalcin ierr = PetscMalloc1(numStrata, &tmpB);CHKERRQ(ierr); 318b9907514SLisandro Dalcin label->numStrata = numStrata; 319b9907514SLisandro Dalcin label->stratumValues = tmpV; 320b9907514SLisandro Dalcin label->stratumSizes = tmpS; 321b9907514SLisandro Dalcin label->ht = tmpH; 322b9907514SLisandro Dalcin label->points = tmpP; 323b9907514SLisandro Dalcin label->validIS = tmpB; 324b9907514SLisandro Dalcin for (v = 0; v < numStrata; ++v) { 325b9907514SLisandro Dalcin ierr = PetscHSetICreate(&ht);CHKERRQ(ierr); 326b9907514SLisandro Dalcin ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);CHKERRQ(ierr); 327b9907514SLisandro Dalcin ierr = PetscHMapISet(hmap, values[v], v);CHKERRQ(ierr); 328b9907514SLisandro Dalcin tmpV[v] = values[v]; 329b9907514SLisandro Dalcin tmpS[v] = 0; 330b9907514SLisandro Dalcin tmpH[v] = ht; 331b9907514SLisandro Dalcin tmpP[v] = is; 332b9907514SLisandro Dalcin tmpB[v] = PETSC_TRUE; 333b9907514SLisandro Dalcin } 334277ea44aSLisandro Dalcin ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr); 335b9907514SLisandro Dalcin } else { 336b9907514SLisandro Dalcin for (v = 0; v < numStrata; ++v) { 337b9907514SLisandro Dalcin ierr = DMLabelAddStratum(label, values[v]);CHKERRQ(ierr); 338b9907514SLisandro Dalcin } 339b9907514SLisandro Dalcin } 340b9907514SLisandro Dalcin ierr = PetscFree(values);CHKERRQ(ierr); 341b9907514SLisandro Dalcin PetscFunctionReturn(0); 342b9907514SLisandro Dalcin } 343b9907514SLisandro Dalcin 344b9907514SLisandro Dalcin /*@ 345b9907514SLisandro Dalcin DMLabelAddStrataIS - Adds new stratum values in a DMLabel 346b9907514SLisandro Dalcin 3475b5e7992SMatthew G. Knepley Not collective 3485b5e7992SMatthew G. Knepley 349b9907514SLisandro Dalcin Input Parameter: 350b9907514SLisandro Dalcin + label - The DMLabel 351b9907514SLisandro Dalcin - valueIS - Index set with stratum values 352b9907514SLisandro Dalcin 353b9907514SLisandro Dalcin Level: beginner 354b9907514SLisandro Dalcin 355b9907514SLisandro Dalcin .seealso: DMLabelCreate(), DMLabelDestroy() 356b9907514SLisandro Dalcin @*/ 357b9907514SLisandro Dalcin PetscErrorCode DMLabelAddStrataIS(DMLabel label, IS valueIS) 358b9907514SLisandro Dalcin { 359b9907514SLisandro Dalcin PetscInt numStrata; 360b9907514SLisandro Dalcin const PetscInt *stratumValues; 361b9907514SLisandro Dalcin PetscErrorCode ierr; 362b9907514SLisandro Dalcin 363b9907514SLisandro Dalcin PetscFunctionBegin; 364b9907514SLisandro Dalcin PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 365b9907514SLisandro Dalcin PetscValidHeaderSpecific(valueIS, IS_CLASSID, 2); 366b9907514SLisandro Dalcin ierr = ISGetLocalSize(valueIS, &numStrata);CHKERRQ(ierr); 367b9907514SLisandro Dalcin ierr = ISGetIndices(valueIS, &stratumValues);CHKERRQ(ierr); 368b9907514SLisandro Dalcin ierr = DMLabelAddStrata(label, numStrata, stratumValues);CHKERRQ(ierr); 369c58f1c22SToby Isaac PetscFunctionReturn(0); 370c58f1c22SToby Isaac } 371c58f1c22SToby Isaac 372c58f1c22SToby Isaac static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer) 373c58f1c22SToby Isaac { 374c58f1c22SToby Isaac PetscInt v; 375c58f1c22SToby Isaac PetscMPIInt rank; 376c58f1c22SToby Isaac PetscErrorCode ierr; 377c58f1c22SToby Isaac 378c58f1c22SToby Isaac PetscFunctionBegin; 379c58f1c22SToby Isaac ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr); 380c58f1c22SToby Isaac ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 381c58f1c22SToby Isaac if (label) { 382d67d17b1SMatthew G. Knepley const char *name; 383d67d17b1SMatthew G. Knepley 384d67d17b1SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr); 385d67d17b1SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", name);CHKERRQ(ierr); 386c58f1c22SToby Isaac if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, " Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);} 387c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 388c58f1c22SToby Isaac const PetscInt value = label->stratumValues[v]; 389ad8374ffSToby Isaac const PetscInt *points; 390c58f1c22SToby Isaac PetscInt p; 391c58f1c22SToby Isaac 392ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 393c58f1c22SToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) { 394ad8374ffSToby Isaac ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr); 395c58f1c22SToby Isaac } 396ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr); 397c58f1c22SToby Isaac } 398c58f1c22SToby Isaac } 399c58f1c22SToby Isaac ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 400c58f1c22SToby Isaac ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 401c58f1c22SToby Isaac PetscFunctionReturn(0); 402c58f1c22SToby Isaac } 403c58f1c22SToby Isaac 404c58f1c22SToby Isaac /*@C 405c58f1c22SToby Isaac DMLabelView - View the label 406c58f1c22SToby Isaac 4075b5e7992SMatthew G. Knepley Collective on viewer 4085b5e7992SMatthew G. Knepley 409c58f1c22SToby Isaac Input Parameters: 410c58f1c22SToby Isaac + label - The DMLabel 411c58f1c22SToby Isaac - viewer - The PetscViewer 412c58f1c22SToby Isaac 413c58f1c22SToby Isaac Level: intermediate 414c58f1c22SToby Isaac 415c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy() 416c58f1c22SToby Isaac @*/ 417c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer) 418c58f1c22SToby Isaac { 419c58f1c22SToby Isaac PetscBool iascii; 420c58f1c22SToby Isaac PetscErrorCode ierr; 421c58f1c22SToby Isaac 422c58f1c22SToby Isaac PetscFunctionBegin; 423d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 424b9907514SLisandro Dalcin if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)label), &viewer);CHKERRQ(ierr);} 425c58f1c22SToby Isaac PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 426c58f1c22SToby Isaac if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 427c58f1c22SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 428c58f1c22SToby Isaac if (iascii) { 429c58f1c22SToby Isaac ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr); 430c58f1c22SToby Isaac } 431c58f1c22SToby Isaac PetscFunctionReturn(0); 432c58f1c22SToby Isaac } 433c58f1c22SToby Isaac 43484f0b6dfSMatthew G. Knepley /*@ 435d67d17b1SMatthew G. Knepley DMLabelReset - Destroys internal data structures in a DMLabel 436d67d17b1SMatthew G. Knepley 4375b5e7992SMatthew G. Knepley Not collective 4385b5e7992SMatthew G. Knepley 439d67d17b1SMatthew G. Knepley Input Parameter: 440d67d17b1SMatthew G. Knepley . label - The DMLabel 441d67d17b1SMatthew G. Knepley 442d67d17b1SMatthew G. Knepley Level: beginner 443d67d17b1SMatthew G. Knepley 444d67d17b1SMatthew G. Knepley .seealso: DMLabelDestroy(), DMLabelCreate() 445d67d17b1SMatthew G. Knepley @*/ 446d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelReset(DMLabel label) 447d67d17b1SMatthew G. Knepley { 448d67d17b1SMatthew G. Knepley PetscInt v; 449d67d17b1SMatthew G. Knepley PetscErrorCode ierr; 450d67d17b1SMatthew G. Knepley 451d67d17b1SMatthew G. Knepley PetscFunctionBegin; 452d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 453d67d17b1SMatthew G. Knepley for (v = 0; v < label->numStrata; ++v) { 454d67d17b1SMatthew G. Knepley ierr = PetscHSetIDestroy(&label->ht[v]);CHKERRQ(ierr); 455d67d17b1SMatthew G. Knepley ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr); 456d67d17b1SMatthew G. Knepley } 457b9907514SLisandro Dalcin label->numStrata = 0; 458b9907514SLisandro Dalcin ierr = PetscFree(label->stratumValues);CHKERRQ(ierr); 459b9907514SLisandro Dalcin ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr); 460d67d17b1SMatthew G. Knepley ierr = PetscFree(label->ht);CHKERRQ(ierr); 461d67d17b1SMatthew G. Knepley ierr = PetscFree(label->points);CHKERRQ(ierr); 462d67d17b1SMatthew G. Knepley ierr = PetscFree(label->validIS);CHKERRQ(ierr); 463b9907514SLisandro Dalcin ierr = PetscHMapIReset(label->hmap);CHKERRQ(ierr); 464b9907514SLisandro Dalcin label->pStart = -1; 465b9907514SLisandro Dalcin label->pEnd = -1; 466d67d17b1SMatthew G. Knepley ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 467d67d17b1SMatthew G. Knepley PetscFunctionReturn(0); 468d67d17b1SMatthew G. Knepley } 469d67d17b1SMatthew G. Knepley 470d67d17b1SMatthew G. Knepley /*@ 47184f0b6dfSMatthew G. Knepley DMLabelDestroy - Destroys a DMLabel 47284f0b6dfSMatthew G. Knepley 4735b5e7992SMatthew G. Knepley Collective on label 4745b5e7992SMatthew G. Knepley 47584f0b6dfSMatthew G. Knepley Input Parameter: 47684f0b6dfSMatthew G. Knepley . label - The DMLabel 47784f0b6dfSMatthew G. Knepley 47884f0b6dfSMatthew G. Knepley Level: beginner 47984f0b6dfSMatthew G. Knepley 480d67d17b1SMatthew G. Knepley .seealso: DMLabelReset(), DMLabelCreate() 48184f0b6dfSMatthew G. Knepley @*/ 482c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label) 483c58f1c22SToby Isaac { 484c58f1c22SToby Isaac PetscErrorCode ierr; 485c58f1c22SToby Isaac 486c58f1c22SToby Isaac PetscFunctionBegin; 487d67d17b1SMatthew G. Knepley if (!*label) PetscFunctionReturn(0); 488d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific((*label),DMLABEL_CLASSID,1); 489b9907514SLisandro Dalcin if (--((PetscObject)(*label))->refct > 0) {*label = NULL; PetscFunctionReturn(0);} 490d67d17b1SMatthew G. Knepley ierr = DMLabelReset(*label);CHKERRQ(ierr); 491b9907514SLisandro Dalcin ierr = PetscHMapIDestroy(&(*label)->hmap);CHKERRQ(ierr); 492d67d17b1SMatthew G. Knepley ierr = PetscHeaderDestroy(label);CHKERRQ(ierr); 493c58f1c22SToby Isaac PetscFunctionReturn(0); 494c58f1c22SToby Isaac } 495c58f1c22SToby Isaac 49684f0b6dfSMatthew G. Knepley /*@ 49784f0b6dfSMatthew G. Knepley DMLabelDuplicate - Duplicates a DMLabel 49884f0b6dfSMatthew G. Knepley 4995b5e7992SMatthew G. Knepley Collective on label 5005b5e7992SMatthew G. Knepley 50184f0b6dfSMatthew G. Knepley Input Parameter: 50284f0b6dfSMatthew G. Knepley . label - The DMLabel 50384f0b6dfSMatthew G. Knepley 50484f0b6dfSMatthew G. Knepley Output Parameter: 50584f0b6dfSMatthew G. Knepley . labelnew - location to put new vector 50684f0b6dfSMatthew G. Knepley 50784f0b6dfSMatthew G. Knepley Level: intermediate 50884f0b6dfSMatthew G. Knepley 50984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelDestroy() 51084f0b6dfSMatthew G. Knepley @*/ 511c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew) 512c58f1c22SToby Isaac { 513d67d17b1SMatthew G. Knepley const char *name; 514ad8374ffSToby Isaac PetscInt v; 515c58f1c22SToby Isaac PetscErrorCode ierr; 516c58f1c22SToby Isaac 517c58f1c22SToby Isaac PetscFunctionBegin; 518d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 519c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 520d67d17b1SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr); 521d67d17b1SMatthew G. Knepley ierr = DMLabelCreate(PetscObjectComm((PetscObject) label), name, labelnew);CHKERRQ(ierr); 522c58f1c22SToby Isaac 523c58f1c22SToby Isaac (*labelnew)->numStrata = label->numStrata; 5245aa44df4SToby Isaac (*labelnew)->defaultValue = label->defaultValue; 525c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr); 526c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr); 527c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr); 528c58f1c22SToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr); 529ad8374ffSToby Isaac ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr); 530c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 531e8f14785SLisandro Dalcin ierr = PetscHSetICreate(&(*labelnew)->ht[v]);CHKERRQ(ierr); 532c58f1c22SToby Isaac (*labelnew)->stratumValues[v] = label->stratumValues[v]; 533c58f1c22SToby Isaac (*labelnew)->stratumSizes[v] = label->stratumSizes[v]; 534ad8374ffSToby Isaac ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr); 535ad8374ffSToby Isaac (*labelnew)->points[v] = label->points[v]; 536b9907514SLisandro Dalcin (*labelnew)->validIS[v] = PETSC_TRUE; 537c58f1c22SToby Isaac } 538f14fe40dSLisandro Dalcin ierr = PetscHMapIDestroy(&(*labelnew)->hmap);CHKERRQ(ierr); 539b9907514SLisandro Dalcin ierr = PetscHMapIDuplicate(label->hmap,&(*labelnew)->hmap);CHKERRQ(ierr); 540c58f1c22SToby Isaac (*labelnew)->pStart = -1; 541c58f1c22SToby Isaac (*labelnew)->pEnd = -1; 542c58f1c22SToby Isaac (*labelnew)->bt = NULL; 543c58f1c22SToby Isaac PetscFunctionReturn(0); 544c58f1c22SToby Isaac } 545c58f1c22SToby Isaac 546c6a43d28SMatthew G. Knepley /*@ 547c6a43d28SMatthew G. Knepley DMLabelComputeIndex - Create an index structure for membership determination, automatically determining the bounds 548c6a43d28SMatthew G. Knepley 5495b5e7992SMatthew G. Knepley Not collective 5505b5e7992SMatthew G. Knepley 551c6a43d28SMatthew G. Knepley Input Parameter: 552c6a43d28SMatthew G. Knepley . label - The DMLabel 553c6a43d28SMatthew G. Knepley 554c6a43d28SMatthew G. Knepley Level: intermediate 555c6a43d28SMatthew G. Knepley 556c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue() 557c6a43d28SMatthew G. Knepley @*/ 558c6a43d28SMatthew G. Knepley PetscErrorCode DMLabelComputeIndex(DMLabel label) 559c6a43d28SMatthew G. Knepley { 560c6a43d28SMatthew G. Knepley PetscInt pStart = PETSC_MAX_INT, pEnd = -1, v; 561c6a43d28SMatthew G. Knepley PetscErrorCode ierr; 562c6a43d28SMatthew G. Knepley 563c6a43d28SMatthew G. Knepley PetscFunctionBegin; 564c6a43d28SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 565c6a43d28SMatthew G. Knepley ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 566c6a43d28SMatthew G. Knepley for (v = 0; v < label->numStrata; ++v) { 567c6a43d28SMatthew G. Knepley const PetscInt *points; 568c6a43d28SMatthew G. Knepley PetscInt i; 569c6a43d28SMatthew G. Knepley 570c6a43d28SMatthew G. Knepley ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 571c6a43d28SMatthew G. Knepley for (i = 0; i < label->stratumSizes[v]; ++i) { 572c6a43d28SMatthew G. Knepley const PetscInt point = points[i]; 573c6a43d28SMatthew G. Knepley 574c6a43d28SMatthew G. Knepley pStart = PetscMin(point, pStart); 575c6a43d28SMatthew G. Knepley pEnd = PetscMax(point+1, pEnd); 576c6a43d28SMatthew G. Knepley } 577c6a43d28SMatthew G. Knepley ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr); 578c6a43d28SMatthew G. Knepley } 579c6a43d28SMatthew G. Knepley label->pStart = pStart == PETSC_MAX_INT ? -1 : pStart; 580c6a43d28SMatthew G. Knepley label->pEnd = pEnd; 581c6a43d28SMatthew G. Knepley ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr); 582c6a43d28SMatthew G. Knepley PetscFunctionReturn(0); 583c6a43d28SMatthew G. Knepley } 584c6a43d28SMatthew G. Knepley 585c6a43d28SMatthew G. Knepley /*@ 586c6a43d28SMatthew G. Knepley DMLabelCreateIndex - Create an index structure for membership determination 587c6a43d28SMatthew G. Knepley 5885b5e7992SMatthew G. Knepley Not collective 5895b5e7992SMatthew G. Knepley 590c6a43d28SMatthew G. Knepley Input Parameters: 591c6a43d28SMatthew G. Knepley + label - The DMLabel 592c6a43d28SMatthew G. Knepley . pStart - The smallest point 593c6a43d28SMatthew G. Knepley - pEnd - The largest point + 1 594c6a43d28SMatthew G. Knepley 595c6a43d28SMatthew G. Knepley Level: intermediate 596c6a43d28SMatthew G. Knepley 597c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelComputeIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue() 598c6a43d28SMatthew G. Knepley @*/ 599c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd) 600c58f1c22SToby Isaac { 601c58f1c22SToby Isaac PetscInt v; 602c58f1c22SToby Isaac PetscErrorCode ierr; 603c58f1c22SToby Isaac 604c58f1c22SToby Isaac PetscFunctionBegin; 605d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 6060c3c4a36SLisandro Dalcin ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr); 607c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 608c58f1c22SToby Isaac label->pStart = pStart; 609c58f1c22SToby Isaac label->pEnd = pEnd; 610c6a43d28SMatthew G. Knepley /* This can be hooked into SetValue(), ClearValue(), etc. for updating */ 611c58f1c22SToby Isaac ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr); 612c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 613ad8374ffSToby Isaac const PetscInt *points; 614c58f1c22SToby Isaac PetscInt i; 615c58f1c22SToby Isaac 616ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 617c58f1c22SToby Isaac for (i = 0; i < label->stratumSizes[v]; ++i) { 618ad8374ffSToby Isaac const PetscInt point = points[i]; 619c58f1c22SToby Isaac 620c58f1c22SToby 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); 621c58f1c22SToby Isaac ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr); 622c58f1c22SToby Isaac } 623ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr); 624c58f1c22SToby Isaac } 625c58f1c22SToby Isaac PetscFunctionReturn(0); 626c58f1c22SToby Isaac } 627c58f1c22SToby Isaac 628c6a43d28SMatthew G. Knepley /*@ 629c6a43d28SMatthew G. Knepley DMLabelDestroyIndex - Destroy the index structure 630c6a43d28SMatthew G. Knepley 6315b5e7992SMatthew G. Knepley Not collective 6325b5e7992SMatthew G. Knepley 633c6a43d28SMatthew G. Knepley Input Parameter: 634c6a43d28SMatthew G. Knepley . label - the DMLabel 635c6a43d28SMatthew G. Knepley 636c6a43d28SMatthew G. Knepley Level: intermediate 637c6a43d28SMatthew G. Knepley 638c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue() 639c6a43d28SMatthew G. Knepley @*/ 640c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label) 641c58f1c22SToby Isaac { 642c58f1c22SToby Isaac PetscErrorCode ierr; 643c58f1c22SToby Isaac 644c58f1c22SToby Isaac PetscFunctionBegin; 645d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 646c58f1c22SToby Isaac label->pStart = -1; 647c58f1c22SToby Isaac label->pEnd = -1; 6480c3c4a36SLisandro Dalcin ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 649c58f1c22SToby Isaac PetscFunctionReturn(0); 650c58f1c22SToby Isaac } 651c58f1c22SToby Isaac 652c58f1c22SToby Isaac /*@ 653c6a43d28SMatthew G. Knepley DMLabelGetBounds - Return the smallest and largest point in the label 654c6a43d28SMatthew G. Knepley 6555b5e7992SMatthew G. Knepley Not collective 6565b5e7992SMatthew G. Knepley 657c6a43d28SMatthew G. Knepley Input Parameter: 658c6a43d28SMatthew G. Knepley . label - the DMLabel 659c6a43d28SMatthew G. Knepley 660c6a43d28SMatthew G. Knepley Output Parameters: 661c6a43d28SMatthew G. Knepley + pStart - The smallest point 662c6a43d28SMatthew G. Knepley - pEnd - The largest point + 1 663c6a43d28SMatthew G. Knepley 664c6a43d28SMatthew G. Knepley Note: This will compute an index for the label if one does not exist. 665c6a43d28SMatthew G. Knepley 666c6a43d28SMatthew G. Knepley Level: intermediate 667c6a43d28SMatthew G. Knepley 668c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue() 669c6a43d28SMatthew G. Knepley @*/ 670c6a43d28SMatthew G. Knepley PetscErrorCode DMLabelGetBounds(DMLabel label, PetscInt *pStart, PetscInt *pEnd) 671c6a43d28SMatthew G. Knepley { 672c6a43d28SMatthew G. Knepley PetscErrorCode ierr; 673c6a43d28SMatthew G. Knepley 674c6a43d28SMatthew G. Knepley PetscFunctionBegin; 675c6a43d28SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 676c6a43d28SMatthew G. Knepley if ((label->pStart == -1) && (label->pEnd == -1)) {ierr = DMLabelComputeIndex(label);CHKERRQ(ierr);} 677c6a43d28SMatthew G. Knepley if (pStart) { 678534a8f05SLisandro Dalcin PetscValidIntPointer(pStart, 2); 679c6a43d28SMatthew G. Knepley *pStart = label->pStart; 680c6a43d28SMatthew G. Knepley } 681c6a43d28SMatthew G. Knepley if (pEnd) { 682534a8f05SLisandro Dalcin PetscValidIntPointer(pEnd, 3); 683c6a43d28SMatthew G. Knepley *pEnd = label->pEnd; 684c6a43d28SMatthew G. Knepley } 685c6a43d28SMatthew G. Knepley PetscFunctionReturn(0); 686c6a43d28SMatthew G. Knepley } 687c6a43d28SMatthew G. Knepley 688c6a43d28SMatthew G. Knepley /*@ 689c58f1c22SToby Isaac DMLabelHasValue - Determine whether a label assigns the value to any point 690c58f1c22SToby Isaac 6915b5e7992SMatthew G. Knepley Not collective 6925b5e7992SMatthew G. Knepley 693c58f1c22SToby Isaac Input Parameters: 694c58f1c22SToby Isaac + label - the DMLabel 695c58f1c22SToby Isaac - value - the value 696c58f1c22SToby Isaac 697c58f1c22SToby Isaac Output Parameter: 698c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point 699c58f1c22SToby Isaac 700c58f1c22SToby Isaac Level: developer 701c58f1c22SToby Isaac 702c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue() 703c58f1c22SToby Isaac @*/ 704c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains) 705c58f1c22SToby Isaac { 706c58f1c22SToby Isaac PetscInt v; 7070c3c4a36SLisandro Dalcin PetscErrorCode ierr; 708c58f1c22SToby Isaac 709c58f1c22SToby Isaac PetscFunctionBegin; 710d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 711534a8f05SLisandro Dalcin PetscValidBoolPointer(contains, 3); 7120c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 7130c3c4a36SLisandro Dalcin *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE; 714c58f1c22SToby Isaac PetscFunctionReturn(0); 715c58f1c22SToby Isaac } 716c58f1c22SToby Isaac 717c58f1c22SToby Isaac /*@ 718c58f1c22SToby Isaac DMLabelHasPoint - Determine whether a label assigns a value to a point 719c58f1c22SToby Isaac 7205b5e7992SMatthew G. Knepley Not collective 7215b5e7992SMatthew G. Knepley 722c58f1c22SToby Isaac Input Parameters: 723c58f1c22SToby Isaac + label - the DMLabel 724c58f1c22SToby Isaac - point - the point 725c58f1c22SToby Isaac 726c58f1c22SToby Isaac Output Parameter: 727c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value 728c58f1c22SToby Isaac 729c58f1c22SToby Isaac Note: The user must call DMLabelCreateIndex() before this function. 730c58f1c22SToby Isaac 731c58f1c22SToby Isaac Level: developer 732c58f1c22SToby Isaac 733c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue() 734c58f1c22SToby Isaac @*/ 735c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains) 736c58f1c22SToby Isaac { 737c58f1c22SToby Isaac PetscErrorCode ierr; 738c58f1c22SToby Isaac 739c58f1c22SToby Isaac PetscFunctionBeginHot; 740d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 741534a8f05SLisandro Dalcin PetscValidBoolPointer(contains, 3); 742c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 743c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG) 744c58f1c22SToby Isaac if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()"); 745c58f1c22SToby 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); 746c58f1c22SToby Isaac #endif 747c58f1c22SToby Isaac *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE; 748c58f1c22SToby Isaac PetscFunctionReturn(0); 749c58f1c22SToby Isaac } 750c58f1c22SToby Isaac 751c58f1c22SToby Isaac /*@ 752c58f1c22SToby Isaac DMLabelStratumHasPoint - Return true if the stratum contains a point 753c58f1c22SToby Isaac 7545b5e7992SMatthew G. Knepley Not collective 7555b5e7992SMatthew G. Knepley 756c58f1c22SToby Isaac Input Parameters: 757c58f1c22SToby Isaac + label - the DMLabel 758c58f1c22SToby Isaac . value - the stratum value 759c58f1c22SToby Isaac - point - the point 760c58f1c22SToby Isaac 761c58f1c22SToby Isaac Output Parameter: 762c58f1c22SToby Isaac . contains - true if the stratum contains the point 763c58f1c22SToby Isaac 764c58f1c22SToby Isaac Level: intermediate 765c58f1c22SToby Isaac 766c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue() 767c58f1c22SToby Isaac @*/ 768c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains) 769c58f1c22SToby Isaac { 770c58f1c22SToby Isaac PetscInt v; 771c58f1c22SToby Isaac PetscErrorCode ierr; 772c58f1c22SToby Isaac 7730c3c4a36SLisandro Dalcin PetscFunctionBeginHot; 774d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 775534a8f05SLisandro Dalcin PetscValidBoolPointer(contains, 4); 776c58f1c22SToby Isaac *contains = PETSC_FALSE; 7770c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 7780c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 7790c3c4a36SLisandro Dalcin 780ad8374ffSToby Isaac if (label->validIS[v]) { 781c58f1c22SToby Isaac PetscInt i; 782c58f1c22SToby Isaac 783a2d74346SToby Isaac ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr); 7840c3c4a36SLisandro Dalcin if (i >= 0) *contains = PETSC_TRUE; 785c58f1c22SToby Isaac } else { 786c58f1c22SToby Isaac PetscBool has; 787c58f1c22SToby Isaac 788b9907514SLisandro Dalcin ierr = PetscHSetIHas(label->ht[v], point, &has);CHKERRQ(ierr); 7890c3c4a36SLisandro Dalcin if (has) *contains = PETSC_TRUE; 790c58f1c22SToby Isaac } 791c58f1c22SToby Isaac PetscFunctionReturn(0); 792c58f1c22SToby Isaac } 793c58f1c22SToby Isaac 79484f0b6dfSMatthew G. Knepley /*@ 7955aa44df4SToby Isaac DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 7965aa44df4SToby Isaac When a label is created, it is initialized to -1. 7975aa44df4SToby Isaac 7985b5e7992SMatthew G. Knepley Not collective 7995b5e7992SMatthew G. Knepley 8005aa44df4SToby Isaac Input parameter: 8015aa44df4SToby Isaac . label - a DMLabel object 8025aa44df4SToby Isaac 8035aa44df4SToby Isaac Output parameter: 8045aa44df4SToby Isaac . defaultValue - the default value 8055aa44df4SToby Isaac 8065aa44df4SToby Isaac Level: beginner 8075aa44df4SToby Isaac 8085aa44df4SToby Isaac .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 80984f0b6dfSMatthew G. Knepley @*/ 8105aa44df4SToby Isaac PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue) 8115aa44df4SToby Isaac { 8125aa44df4SToby Isaac PetscFunctionBegin; 813d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 8145aa44df4SToby Isaac *defaultValue = label->defaultValue; 8155aa44df4SToby Isaac PetscFunctionReturn(0); 8165aa44df4SToby Isaac } 8175aa44df4SToby Isaac 81884f0b6dfSMatthew G. Knepley /*@ 8195aa44df4SToby Isaac DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value. 8205aa44df4SToby Isaac When a label is created, it is initialized to -1. 8215aa44df4SToby Isaac 8225b5e7992SMatthew G. Knepley Not collective 8235b5e7992SMatthew G. Knepley 8245aa44df4SToby Isaac Input parameter: 8255aa44df4SToby Isaac . label - a DMLabel object 8265aa44df4SToby Isaac 8275aa44df4SToby Isaac Output parameter: 8285aa44df4SToby Isaac . defaultValue - the default value 8295aa44df4SToby Isaac 8305aa44df4SToby Isaac Level: beginner 8315aa44df4SToby Isaac 8325aa44df4SToby Isaac .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue() 83384f0b6dfSMatthew G. Knepley @*/ 8345aa44df4SToby Isaac PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue) 8355aa44df4SToby Isaac { 8365aa44df4SToby Isaac PetscFunctionBegin; 837d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 8385aa44df4SToby Isaac label->defaultValue = defaultValue; 8395aa44df4SToby Isaac PetscFunctionReturn(0); 8405aa44df4SToby Isaac } 8415aa44df4SToby Isaac 842c58f1c22SToby Isaac /*@ 8435aa44df4SToby 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()) 844c58f1c22SToby Isaac 8455b5e7992SMatthew G. Knepley Not collective 8465b5e7992SMatthew G. Knepley 847c58f1c22SToby Isaac Input Parameters: 848c58f1c22SToby Isaac + label - the DMLabel 849c58f1c22SToby Isaac - point - the point 850c58f1c22SToby Isaac 851c58f1c22SToby Isaac Output Parameter: 8528e68afcfSMatthew G. Knepley . value - The point value, or the default value (-1 by default) 853c58f1c22SToby Isaac 854c58f1c22SToby Isaac Level: intermediate 855c58f1c22SToby Isaac 8565aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 857c58f1c22SToby Isaac @*/ 858c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value) 859c58f1c22SToby Isaac { 860c58f1c22SToby Isaac PetscInt v; 861c58f1c22SToby Isaac PetscErrorCode ierr; 862c58f1c22SToby Isaac 8630c3c4a36SLisandro Dalcin PetscFunctionBeginHot; 864d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 865c58f1c22SToby Isaac PetscValidPointer(value, 3); 8665aa44df4SToby Isaac *value = label->defaultValue; 867c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 868ad8374ffSToby Isaac if (label->validIS[v]) { 869c58f1c22SToby Isaac PetscInt i; 870c58f1c22SToby Isaac 871a2d74346SToby Isaac ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr); 872c58f1c22SToby Isaac if (i >= 0) { 873c58f1c22SToby Isaac *value = label->stratumValues[v]; 874c58f1c22SToby Isaac break; 875c58f1c22SToby Isaac } 876c58f1c22SToby Isaac } else { 877c58f1c22SToby Isaac PetscBool has; 878c58f1c22SToby Isaac 879b9907514SLisandro Dalcin ierr = PetscHSetIHas(label->ht[v], point, &has);CHKERRQ(ierr); 880c58f1c22SToby Isaac if (has) { 881c58f1c22SToby Isaac *value = label->stratumValues[v]; 882c58f1c22SToby Isaac break; 883c58f1c22SToby Isaac } 884c58f1c22SToby Isaac } 885c58f1c22SToby Isaac } 886c58f1c22SToby Isaac PetscFunctionReturn(0); 887c58f1c22SToby Isaac } 888c58f1c22SToby Isaac 889c58f1c22SToby Isaac /*@ 890367003a6SStefano Zampini DMLabelSetValue - Set the value a label assigns to a point. If the value is the same as the label's default value (which is initially -1, and can be changed with DMLabelSetDefaultValue() to something different), then this function will do nothing. 891c58f1c22SToby Isaac 8925b5e7992SMatthew G. Knepley Not collective 8935b5e7992SMatthew G. Knepley 894c58f1c22SToby Isaac Input Parameters: 895c58f1c22SToby Isaac + label - the DMLabel 896c58f1c22SToby Isaac . point - the point 897c58f1c22SToby Isaac - value - The point value 898c58f1c22SToby Isaac 899c58f1c22SToby Isaac Level: intermediate 900c58f1c22SToby Isaac 9015aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue() 902c58f1c22SToby Isaac @*/ 903c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value) 904c58f1c22SToby Isaac { 905c58f1c22SToby Isaac PetscInt v; 906c58f1c22SToby Isaac PetscErrorCode ierr; 907c58f1c22SToby Isaac 908c58f1c22SToby Isaac PetscFunctionBegin; 909d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 9100c3c4a36SLisandro Dalcin /* Find label value, add new entry if needed */ 9115aa44df4SToby Isaac if (value == label->defaultValue) PetscFunctionReturn(0); 912b9907514SLisandro Dalcin ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr); 913c58f1c22SToby Isaac /* Set key */ 9140c3c4a36SLisandro Dalcin ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 915e8f14785SLisandro Dalcin ierr = PetscHSetIAdd(label->ht[v], point);CHKERRQ(ierr); 916c58f1c22SToby Isaac PetscFunctionReturn(0); 917c58f1c22SToby Isaac } 918c58f1c22SToby Isaac 919c58f1c22SToby Isaac /*@ 920c58f1c22SToby Isaac DMLabelClearValue - Clear the value a label assigns to a point 921c58f1c22SToby Isaac 9225b5e7992SMatthew G. Knepley Not collective 9235b5e7992SMatthew G. Knepley 924c58f1c22SToby Isaac Input Parameters: 925c58f1c22SToby Isaac + label - the DMLabel 926c58f1c22SToby Isaac . point - the point 927c58f1c22SToby Isaac - value - The point value 928c58f1c22SToby Isaac 929c58f1c22SToby Isaac Level: intermediate 930c58f1c22SToby Isaac 931c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue() 932c58f1c22SToby Isaac @*/ 933c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value) 934c58f1c22SToby Isaac { 935ad8374ffSToby Isaac PetscInt v; 936c58f1c22SToby Isaac PetscErrorCode ierr; 937c58f1c22SToby Isaac 938c58f1c22SToby Isaac PetscFunctionBegin; 939d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 940c58f1c22SToby Isaac /* Find label value */ 9410c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 9420c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 9430c3c4a36SLisandro Dalcin 944eeed21e7SToby Isaac if (label->bt) { 945eeed21e7SToby 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); 946eeed21e7SToby Isaac ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 947eeed21e7SToby Isaac } 9480c3c4a36SLisandro Dalcin 9490c3c4a36SLisandro Dalcin /* Delete key */ 9500c3c4a36SLisandro Dalcin ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 951e8f14785SLisandro Dalcin ierr = PetscHSetIDel(label->ht[v], point);CHKERRQ(ierr); 952c58f1c22SToby Isaac PetscFunctionReturn(0); 953c58f1c22SToby Isaac } 954c58f1c22SToby Isaac 955c58f1c22SToby Isaac /*@ 956c58f1c22SToby Isaac DMLabelInsertIS - Set all points in the IS to a value 957c58f1c22SToby Isaac 9585b5e7992SMatthew G. Knepley Not collective 9595b5e7992SMatthew G. Knepley 960c58f1c22SToby Isaac Input Parameters: 961c58f1c22SToby Isaac + label - the DMLabel 962c58f1c22SToby Isaac . is - the point IS 963c58f1c22SToby Isaac - value - The point value 964c58f1c22SToby Isaac 965c58f1c22SToby Isaac Level: intermediate 966c58f1c22SToby Isaac 967c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 968c58f1c22SToby Isaac @*/ 969c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value) 970c58f1c22SToby Isaac { 9710c3c4a36SLisandro Dalcin PetscInt v, n, p; 972c58f1c22SToby Isaac const PetscInt *points; 973c58f1c22SToby Isaac PetscErrorCode ierr; 974c58f1c22SToby Isaac 975c58f1c22SToby Isaac PetscFunctionBegin; 976d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 977c58f1c22SToby Isaac PetscValidHeaderSpecific(is, IS_CLASSID, 2); 9780c3c4a36SLisandro Dalcin /* Find label value, add new entry if needed */ 9790c3c4a36SLisandro Dalcin if (value == label->defaultValue) PetscFunctionReturn(0); 980b9907514SLisandro Dalcin ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr); 9810c3c4a36SLisandro Dalcin /* Set keys */ 9820c3c4a36SLisandro Dalcin ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr); 983c58f1c22SToby Isaac ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 984c58f1c22SToby Isaac ierr = ISGetIndices(is, &points);CHKERRQ(ierr); 985e8f14785SLisandro Dalcin for (p = 0; p < n; ++p) {ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);} 986c58f1c22SToby Isaac ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr); 987c58f1c22SToby Isaac PetscFunctionReturn(0); 988c58f1c22SToby Isaac } 989c58f1c22SToby Isaac 99084f0b6dfSMatthew G. Knepley /*@ 99184f0b6dfSMatthew G. Knepley DMLabelGetNumValues - Get the number of values that the DMLabel takes 99284f0b6dfSMatthew G. Knepley 9935b5e7992SMatthew G. Knepley Not collective 9945b5e7992SMatthew G. Knepley 99584f0b6dfSMatthew G. Knepley Input Parameter: 99684f0b6dfSMatthew G. Knepley . label - the DMLabel 99784f0b6dfSMatthew G. Knepley 99884f0b6dfSMatthew G. Knepley Output Paramater: 99984f0b6dfSMatthew G. Knepley . numValues - the number of values 100084f0b6dfSMatthew G. Knepley 100184f0b6dfSMatthew G. Knepley Level: intermediate 100284f0b6dfSMatthew G. Knepley 100384f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 100484f0b6dfSMatthew G. Knepley @*/ 1005c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues) 1006c58f1c22SToby Isaac { 1007c58f1c22SToby Isaac PetscFunctionBegin; 1008d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1009c58f1c22SToby Isaac PetscValidPointer(numValues, 2); 1010c58f1c22SToby Isaac *numValues = label->numStrata; 1011c58f1c22SToby Isaac PetscFunctionReturn(0); 1012c58f1c22SToby Isaac } 1013c58f1c22SToby Isaac 101484f0b6dfSMatthew G. Knepley /*@ 101584f0b6dfSMatthew G. Knepley DMLabelGetValueIS - Get an IS of all values that the DMlabel takes 101684f0b6dfSMatthew G. Knepley 10175b5e7992SMatthew G. Knepley Not collective 10185b5e7992SMatthew G. Knepley 101984f0b6dfSMatthew G. Knepley Input Parameter: 102084f0b6dfSMatthew G. Knepley . label - the DMLabel 102184f0b6dfSMatthew G. Knepley 102284f0b6dfSMatthew G. Knepley Output Paramater: 102384f0b6dfSMatthew G. Knepley . is - the value IS 102484f0b6dfSMatthew G. Knepley 102584f0b6dfSMatthew G. Knepley Level: intermediate 102684f0b6dfSMatthew G. Knepley 102784f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 102884f0b6dfSMatthew G. Knepley @*/ 1029c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values) 1030c58f1c22SToby Isaac { 1031c58f1c22SToby Isaac PetscErrorCode ierr; 1032c58f1c22SToby Isaac 1033c58f1c22SToby Isaac PetscFunctionBegin; 1034d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1035c58f1c22SToby Isaac PetscValidPointer(values, 2); 1036c58f1c22SToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr); 1037c58f1c22SToby Isaac PetscFunctionReturn(0); 1038c58f1c22SToby Isaac } 1039c58f1c22SToby Isaac 104084f0b6dfSMatthew G. Knepley /*@ 104184f0b6dfSMatthew G. Knepley DMLabelHasStratum - Determine whether points exist with the given value 104284f0b6dfSMatthew G. Knepley 10435b5e7992SMatthew G. Knepley Not collective 10445b5e7992SMatthew G. Knepley 104584f0b6dfSMatthew G. Knepley Input Parameters: 104684f0b6dfSMatthew G. Knepley + label - the DMLabel 104784f0b6dfSMatthew G. Knepley - value - the stratum value 104884f0b6dfSMatthew G. Knepley 104984f0b6dfSMatthew G. Knepley Output Paramater: 105084f0b6dfSMatthew G. Knepley . exists - Flag saying whether points exist 105184f0b6dfSMatthew G. Knepley 105284f0b6dfSMatthew G. Knepley Level: intermediate 105384f0b6dfSMatthew G. Knepley 105484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 105584f0b6dfSMatthew G. Knepley @*/ 1056fada774cSMatthew G. Knepley PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists) 1057fada774cSMatthew G. Knepley { 1058fada774cSMatthew G. Knepley PetscInt v; 10590c3c4a36SLisandro Dalcin PetscErrorCode ierr; 1060fada774cSMatthew G. Knepley 1061fada774cSMatthew G. Knepley PetscFunctionBegin; 1062d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1063fada774cSMatthew G. Knepley PetscValidPointer(exists, 3); 10640c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 10650c3c4a36SLisandro Dalcin *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE; 1066fada774cSMatthew G. Knepley PetscFunctionReturn(0); 1067fada774cSMatthew G. Knepley } 1068fada774cSMatthew G. Knepley 106984f0b6dfSMatthew G. Knepley /*@ 107084f0b6dfSMatthew G. Knepley DMLabelGetStratumSize - Get the size of a stratum 107184f0b6dfSMatthew G. Knepley 10725b5e7992SMatthew G. Knepley Not collective 10735b5e7992SMatthew G. Knepley 107484f0b6dfSMatthew G. Knepley Input Parameters: 107584f0b6dfSMatthew G. Knepley + label - the DMLabel 107684f0b6dfSMatthew G. Knepley - value - the stratum value 107784f0b6dfSMatthew G. Knepley 107884f0b6dfSMatthew G. Knepley Output Paramater: 107984f0b6dfSMatthew G. Knepley . size - The number of points in the stratum 108084f0b6dfSMatthew G. Knepley 108184f0b6dfSMatthew G. Knepley Level: intermediate 108284f0b6dfSMatthew G. Knepley 108384f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 108484f0b6dfSMatthew G. Knepley @*/ 1085c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size) 1086c58f1c22SToby Isaac { 1087c58f1c22SToby Isaac PetscInt v; 1088c58f1c22SToby Isaac PetscErrorCode ierr; 1089c58f1c22SToby Isaac 1090c58f1c22SToby Isaac PetscFunctionBegin; 1091d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1092c58f1c22SToby Isaac PetscValidPointer(size, 3); 1093c58f1c22SToby Isaac *size = 0; 10940c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 10950c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 1096c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 1097c58f1c22SToby Isaac *size = label->stratumSizes[v]; 1098c58f1c22SToby Isaac PetscFunctionReturn(0); 1099c58f1c22SToby Isaac } 1100c58f1c22SToby Isaac 110184f0b6dfSMatthew G. Knepley /*@ 110284f0b6dfSMatthew G. Knepley DMLabelGetStratumBounds - Get the largest and smallest point of a stratum 110384f0b6dfSMatthew G. Knepley 11045b5e7992SMatthew G. Knepley Not collective 11055b5e7992SMatthew G. Knepley 110684f0b6dfSMatthew G. Knepley Input Parameters: 110784f0b6dfSMatthew G. Knepley + label - the DMLabel 110884f0b6dfSMatthew G. Knepley - value - the stratum value 110984f0b6dfSMatthew G. Knepley 111084f0b6dfSMatthew G. Knepley Output Paramaters: 111184f0b6dfSMatthew G. Knepley + start - the smallest point in the stratum 111284f0b6dfSMatthew G. Knepley - end - the largest point in the stratum 111384f0b6dfSMatthew G. Knepley 111484f0b6dfSMatthew G. Knepley Level: intermediate 111584f0b6dfSMatthew G. Knepley 111684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 111784f0b6dfSMatthew G. Knepley @*/ 1118c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end) 1119c58f1c22SToby Isaac { 11200c3c4a36SLisandro Dalcin PetscInt v, min, max; 1121c58f1c22SToby Isaac PetscErrorCode ierr; 1122c58f1c22SToby Isaac 1123c58f1c22SToby Isaac PetscFunctionBegin; 1124d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1125c58f1c22SToby Isaac if (start) {PetscValidPointer(start, 3); *start = 0;} 1126c58f1c22SToby Isaac if (end) {PetscValidPointer(end, 4); *end = 0;} 11270c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 11280c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 1129c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 11300c3c4a36SLisandro Dalcin if (label->stratumSizes[v] <= 0) PetscFunctionReturn(0); 1131d6cb179aSToby Isaac ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr); 1132d6cb179aSToby Isaac if (start) *start = min; 1133d6cb179aSToby Isaac if (end) *end = max+1; 1134c58f1c22SToby Isaac PetscFunctionReturn(0); 1135c58f1c22SToby Isaac } 1136c58f1c22SToby Isaac 113784f0b6dfSMatthew G. Knepley /*@ 113884f0b6dfSMatthew G. Knepley DMLabelGetStratumIS - Get an IS with the stratum points 113984f0b6dfSMatthew G. Knepley 11405b5e7992SMatthew G. Knepley Not collective 11415b5e7992SMatthew G. Knepley 114284f0b6dfSMatthew G. Knepley Input Parameters: 114384f0b6dfSMatthew G. Knepley + label - the DMLabel 114484f0b6dfSMatthew G. Knepley - value - the stratum value 114584f0b6dfSMatthew G. Knepley 114684f0b6dfSMatthew G. Knepley Output Paramater: 114784f0b6dfSMatthew G. Knepley . points - The stratum points 114884f0b6dfSMatthew G. Knepley 114984f0b6dfSMatthew G. Knepley Level: intermediate 115084f0b6dfSMatthew G. Knepley 1151*593ce467SVaclav Hapla Notes: 1152*593ce467SVaclav Hapla The output IS should be destroyed when no longer needed. 1153*593ce467SVaclav Hapla Returns NULL if the stratum is empty. 1154*593ce467SVaclav Hapla 115584f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 115684f0b6dfSMatthew G. Knepley @*/ 1157c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points) 1158c58f1c22SToby Isaac { 1159c58f1c22SToby Isaac PetscInt v; 1160c58f1c22SToby Isaac PetscErrorCode ierr; 1161c58f1c22SToby Isaac 1162c58f1c22SToby Isaac PetscFunctionBegin; 1163d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1164c58f1c22SToby Isaac PetscValidPointer(points, 3); 1165c58f1c22SToby Isaac *points = NULL; 11660c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 11670c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 1168c58f1c22SToby Isaac ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr); 1169ad8374ffSToby Isaac ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr); 1170ad8374ffSToby Isaac *points = label->points[v]; 1171c58f1c22SToby Isaac PetscFunctionReturn(0); 1172c58f1c22SToby Isaac } 1173c58f1c22SToby Isaac 117484f0b6dfSMatthew G. Knepley /*@ 117584f0b6dfSMatthew G. Knepley DMLabelSetStratumIS - Set the stratum points using an IS 117684f0b6dfSMatthew G. Knepley 11775b5e7992SMatthew G. Knepley Not collective 11785b5e7992SMatthew G. Knepley 117984f0b6dfSMatthew G. Knepley Input Parameters: 118084f0b6dfSMatthew G. Knepley + label - the DMLabel 118184f0b6dfSMatthew G. Knepley . value - the stratum value 118284f0b6dfSMatthew G. Knepley - points - The stratum points 118384f0b6dfSMatthew G. Knepley 118484f0b6dfSMatthew G. Knepley Level: intermediate 118584f0b6dfSMatthew G. Knepley 118684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 118784f0b6dfSMatthew G. Knepley @*/ 11884de306b1SToby Isaac PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is) 11894de306b1SToby Isaac { 11900c3c4a36SLisandro Dalcin PetscInt v; 11914de306b1SToby Isaac PetscErrorCode ierr; 11924de306b1SToby Isaac 11934de306b1SToby Isaac PetscFunctionBegin; 1194d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1195d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(is, IS_CLASSID, 3); 1196b9907514SLisandro Dalcin ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr); 11974de306b1SToby Isaac if (is == label->points[v]) PetscFunctionReturn(0); 11984de306b1SToby Isaac ierr = DMLabelClearStratum(label, value);CHKERRQ(ierr); 11994de306b1SToby Isaac ierr = ISGetLocalSize(is, &(label->stratumSizes[v]));CHKERRQ(ierr); 12004de306b1SToby Isaac ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 12014de306b1SToby Isaac ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr); 12020c3c4a36SLisandro Dalcin label->points[v] = is; 12030c3c4a36SLisandro Dalcin label->validIS[v] = PETSC_TRUE; 1204277ea44aSLisandro Dalcin ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr); 12054de306b1SToby Isaac if (label->bt) { 12064de306b1SToby Isaac const PetscInt *points; 12074de306b1SToby Isaac PetscInt p; 12084de306b1SToby Isaac 12094de306b1SToby Isaac ierr = ISGetIndices(is,&points);CHKERRQ(ierr); 12104de306b1SToby Isaac for (p = 0; p < label->stratumSizes[v]; ++p) { 12114de306b1SToby Isaac const PetscInt point = points[p]; 12124de306b1SToby Isaac 12134de306b1SToby 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); 12144de306b1SToby Isaac ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); 12154de306b1SToby Isaac } 12164de306b1SToby Isaac } 12174de306b1SToby Isaac PetscFunctionReturn(0); 12184de306b1SToby Isaac } 12194de306b1SToby Isaac 122084f0b6dfSMatthew G. Knepley /*@ 122184f0b6dfSMatthew G. Knepley DMLabelClearStratum - Remove a stratum 12224de306b1SToby Isaac 12235b5e7992SMatthew G. Knepley Not collective 12245b5e7992SMatthew G. Knepley 122584f0b6dfSMatthew G. Knepley Input Parameters: 122684f0b6dfSMatthew G. Knepley + label - the DMLabel 122784f0b6dfSMatthew G. Knepley - value - the stratum value 122884f0b6dfSMatthew G. Knepley 122984f0b6dfSMatthew G. Knepley Level: intermediate 123084f0b6dfSMatthew G. Knepley 123184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 123284f0b6dfSMatthew G. Knepley @*/ 1233c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value) 1234c58f1c22SToby Isaac { 1235c58f1c22SToby Isaac PetscInt v; 1236c58f1c22SToby Isaac PetscErrorCode ierr; 1237c58f1c22SToby Isaac 1238c58f1c22SToby Isaac PetscFunctionBegin; 1239d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 12400c3c4a36SLisandro Dalcin ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr); 12410c3c4a36SLisandro Dalcin if (v < 0) PetscFunctionReturn(0); 12424de306b1SToby Isaac if (label->validIS[v]) { 12434de306b1SToby Isaac if (label->bt) { 1244c58f1c22SToby Isaac PetscInt i; 1245ad8374ffSToby Isaac const PetscInt *points; 1246c58f1c22SToby Isaac 1247ad8374ffSToby Isaac ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr); 1248c58f1c22SToby Isaac for (i = 0; i < label->stratumSizes[v]; ++i) { 1249ad8374ffSToby Isaac const PetscInt point = points[i]; 1250c58f1c22SToby Isaac 1251c58f1c22SToby 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); 1252c58f1c22SToby Isaac ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr); 1253c58f1c22SToby Isaac } 1254ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr); 1255c58f1c22SToby Isaac } 1256c58f1c22SToby Isaac label->stratumSizes[v] = 0; 12570c3c4a36SLisandro Dalcin ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr); 1258277ea44aSLisandro Dalcin ierr = ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &label->points[v]);CHKERRQ(ierr); 12590c3c4a36SLisandro Dalcin ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr); 1260277ea44aSLisandro Dalcin ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr); 1261c58f1c22SToby Isaac } else { 1262b9907514SLisandro Dalcin ierr = PetscHSetIClear(label->ht[v]);CHKERRQ(ierr); 1263c58f1c22SToby Isaac } 1264c58f1c22SToby Isaac PetscFunctionReturn(0); 1265c58f1c22SToby Isaac } 1266c58f1c22SToby Isaac 126784f0b6dfSMatthew G. Knepley /*@ 126884f0b6dfSMatthew G. Knepley DMLabelFilter - Remove all points outside of [start, end) 126984f0b6dfSMatthew G. Knepley 12705b5e7992SMatthew G. Knepley Not collective 12715b5e7992SMatthew G. Knepley 127284f0b6dfSMatthew G. Knepley Input Parameters: 127384f0b6dfSMatthew G. Knepley + label - the DMLabel 127422cd2cdaSVaclav Hapla . start - the first point kept 127522cd2cdaSVaclav Hapla - end - one more than the last point kept 127684f0b6dfSMatthew G. Knepley 127784f0b6dfSMatthew G. Knepley Level: intermediate 127884f0b6dfSMatthew G. Knepley 127984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 128084f0b6dfSMatthew G. Knepley @*/ 1281c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end) 1282c58f1c22SToby Isaac { 1283c58f1c22SToby Isaac PetscInt v; 1284c58f1c22SToby Isaac PetscErrorCode ierr; 1285c58f1c22SToby Isaac 1286c58f1c22SToby Isaac PetscFunctionBegin; 1287d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 12880c3c4a36SLisandro Dalcin ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr); 1289c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 1290c58f1c22SToby Isaac for (v = 0; v < label->numStrata; ++v) { 12919106b633SVaclav Hapla ierr = ISGeneralFilter(label->points[v], start, end);CHKERRQ(ierr); 1292c58f1c22SToby Isaac } 1293c58f1c22SToby Isaac ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr); 1294c58f1c22SToby Isaac PetscFunctionReturn(0); 1295c58f1c22SToby Isaac } 1296c58f1c22SToby Isaac 129784f0b6dfSMatthew G. Knepley /*@ 129884f0b6dfSMatthew G. Knepley DMLabelPermute - Create a new label with permuted points 129984f0b6dfSMatthew G. Knepley 13005b5e7992SMatthew G. Knepley Not collective 13015b5e7992SMatthew G. Knepley 130284f0b6dfSMatthew G. Knepley Input Parameters: 130384f0b6dfSMatthew G. Knepley + label - the DMLabel 130484f0b6dfSMatthew G. Knepley - permutation - the point permutation 130584f0b6dfSMatthew G. Knepley 130684f0b6dfSMatthew G. Knepley Output Parameter: 130784f0b6dfSMatthew G. Knepley . labelnew - the new label containing the permuted points 130884f0b6dfSMatthew G. Knepley 130984f0b6dfSMatthew G. Knepley Level: intermediate 131084f0b6dfSMatthew G. Knepley 131184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 131284f0b6dfSMatthew G. Knepley @*/ 1313c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew) 1314c58f1c22SToby Isaac { 1315c58f1c22SToby Isaac const PetscInt *perm; 1316c58f1c22SToby Isaac PetscInt numValues, numPoints, v, q; 1317c58f1c22SToby Isaac PetscErrorCode ierr; 1318c58f1c22SToby Isaac 1319c58f1c22SToby Isaac PetscFunctionBegin; 1320d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1321d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(permutation, IS_CLASSID, 2); 1322c58f1c22SToby Isaac ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 1323c58f1c22SToby Isaac ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr); 1324c58f1c22SToby Isaac ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr); 1325c58f1c22SToby Isaac ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr); 1326c58f1c22SToby Isaac ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr); 1327c58f1c22SToby Isaac for (v = 0; v < numValues; ++v) { 1328c58f1c22SToby Isaac const PetscInt size = (*labelNew)->stratumSizes[v]; 1329ad8374ffSToby Isaac const PetscInt *points; 1330ad8374ffSToby Isaac PetscInt *pointsNew; 1331c58f1c22SToby Isaac 1332ad8374ffSToby Isaac ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 1333ad8374ffSToby Isaac ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr); 1334c58f1c22SToby Isaac for (q = 0; q < size; ++q) { 1335ad8374ffSToby Isaac const PetscInt point = points[q]; 1336c58f1c22SToby Isaac 1337c58f1c22SToby 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); 1338ad8374ffSToby Isaac pointsNew[q] = perm[point]; 1339c58f1c22SToby Isaac } 1340ad8374ffSToby Isaac ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr); 1341ad8374ffSToby Isaac ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr); 1342ad8374ffSToby Isaac ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr); 1343fa8e8ae5SToby Isaac if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) { 1344fa8e8ae5SToby Isaac ierr = ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));CHKERRQ(ierr); 1345fa8e8ae5SToby Isaac ierr = PetscFree(pointsNew);CHKERRQ(ierr); 1346fa8e8ae5SToby Isaac } else { 1347ad8374ffSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr); 1348fa8e8ae5SToby Isaac } 1349ad8374ffSToby Isaac ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr); 1350c58f1c22SToby Isaac } 1351c58f1c22SToby Isaac ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr); 1352c58f1c22SToby Isaac if (label->bt) { 1353c58f1c22SToby Isaac ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr); 1354c58f1c22SToby Isaac ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr); 1355c58f1c22SToby Isaac } 1356c58f1c22SToby Isaac PetscFunctionReturn(0); 1357c58f1c22SToby Isaac } 1358c58f1c22SToby Isaac 135926c55118SMichael Lange PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata) 136026c55118SMichael Lange { 136126c55118SMichael Lange MPI_Comm comm; 136226c55118SMichael Lange PetscInt s, l, nroots, nleaves, dof, offset, size; 136326c55118SMichael Lange PetscInt *remoteOffsets, *rootStrata, *rootIdx; 136426c55118SMichael Lange PetscSection rootSection; 136526c55118SMichael Lange PetscSF labelSF; 136626c55118SMichael Lange PetscErrorCode ierr; 136726c55118SMichael Lange 136826c55118SMichael Lange PetscFunctionBegin; 136926c55118SMichael Lange if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} 137026c55118SMichael Lange ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 137126c55118SMichael Lange /* Build a section of stratum values per point, generate the according SF 137226c55118SMichael Lange and distribute point-wise stratum values to leaves. */ 137326c55118SMichael Lange ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 137426c55118SMichael Lange ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 137526c55118SMichael Lange ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr); 137626c55118SMichael Lange if (label) { 137726c55118SMichael Lange for (s = 0; s < label->numStrata; ++s) { 1378ad8374ffSToby Isaac const PetscInt *points; 1379ad8374ffSToby Isaac 1380ad8374ffSToby Isaac ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 138126c55118SMichael Lange for (l = 0; l < label->stratumSizes[s]; l++) { 1382ad8374ffSToby Isaac ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr); 1383ad8374ffSToby Isaac ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr); 138426c55118SMichael Lange } 1385ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 138626c55118SMichael Lange } 138726c55118SMichael Lange } 138826c55118SMichael Lange ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 138926c55118SMichael Lange /* Create a point-wise array of stratum values */ 139026c55118SMichael Lange ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr); 139126c55118SMichael Lange ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr); 139226c55118SMichael Lange ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr); 139326c55118SMichael Lange if (label) { 139426c55118SMichael Lange for (s = 0; s < label->numStrata; ++s) { 1395ad8374ffSToby Isaac const PetscInt *points; 1396ad8374ffSToby Isaac 1397ad8374ffSToby Isaac ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr); 139826c55118SMichael Lange for (l = 0; l < label->stratumSizes[s]; l++) { 1399ad8374ffSToby Isaac const PetscInt p = points[l]; 140026c55118SMichael Lange ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr); 140126c55118SMichael Lange rootStrata[offset+rootIdx[p]++] = label->stratumValues[s]; 140226c55118SMichael Lange } 1403ad8374ffSToby Isaac ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr); 140426c55118SMichael Lange } 140526c55118SMichael Lange } 140626c55118SMichael Lange /* Build SF that maps label points to remote processes */ 140726c55118SMichael Lange ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr); 140826c55118SMichael Lange ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr); 140926c55118SMichael Lange ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr); 141026c55118SMichael Lange ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 141126c55118SMichael Lange /* Send the strata for each point over the derived SF */ 141226c55118SMichael Lange ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr); 141326c55118SMichael Lange ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr); 141426c55118SMichael Lange ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 141526c55118SMichael Lange ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); 141626c55118SMichael Lange /* Clean up */ 141726c55118SMichael Lange ierr = PetscFree(rootStrata);CHKERRQ(ierr); 141826c55118SMichael Lange ierr = PetscFree(rootIdx);CHKERRQ(ierr); 141926c55118SMichael Lange ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 142026c55118SMichael Lange ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr); 142126c55118SMichael Lange PetscFunctionReturn(0); 142226c55118SMichael Lange } 142326c55118SMichael Lange 142484f0b6dfSMatthew G. Knepley /*@ 142584f0b6dfSMatthew G. Knepley DMLabelDistribute - Create a new label pushed forward over the PetscSF 142684f0b6dfSMatthew G. Knepley 14275b5e7992SMatthew G. Knepley Collective on sf 14285b5e7992SMatthew G. Knepley 142984f0b6dfSMatthew G. Knepley Input Parameters: 143084f0b6dfSMatthew G. Knepley + label - the DMLabel 143184f0b6dfSMatthew G. Knepley - sf - the map from old to new distribution 143284f0b6dfSMatthew G. Knepley 143384f0b6dfSMatthew G. Knepley Output Parameter: 143484f0b6dfSMatthew G. Knepley . labelnew - the new redistributed label 143584f0b6dfSMatthew G. Knepley 143684f0b6dfSMatthew G. Knepley Level: intermediate 143784f0b6dfSMatthew G. Knepley 143884f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue() 143984f0b6dfSMatthew G. Knepley @*/ 1440c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew) 1441c58f1c22SToby Isaac { 1442c58f1c22SToby Isaac MPI_Comm comm; 144326c55118SMichael Lange PetscSection leafSection; 144426c55118SMichael Lange PetscInt p, pStart, pEnd, s, size, dof, offset, stratum; 144526c55118SMichael Lange PetscInt *leafStrata, *strataIdx; 1446ad8374ffSToby Isaac PetscInt **points; 1447d67d17b1SMatthew G. Knepley const char *lname = NULL; 1448c58f1c22SToby Isaac char *name; 1449c58f1c22SToby Isaac PetscInt nameSize; 1450e8f14785SLisandro Dalcin PetscHSetI stratumHash; 1451c58f1c22SToby Isaac size_t len = 0; 145226c55118SMichael Lange PetscMPIInt rank; 1453c58f1c22SToby Isaac PetscErrorCode ierr; 1454c58f1c22SToby Isaac 1455c58f1c22SToby Isaac PetscFunctionBegin; 1456d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 1457f018e600SMatthew G. Knepley if (label) { 1458f018e600SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1459f018e600SMatthew G. Knepley ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); 1460f018e600SMatthew G. Knepley } 1461c58f1c22SToby Isaac ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 1462c58f1c22SToby Isaac ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1463c58f1c22SToby Isaac /* Bcast name */ 1464d67d17b1SMatthew G. Knepley if (!rank) { 1465d67d17b1SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr); 1466d67d17b1SMatthew G. Knepley ierr = PetscStrlen(lname, &len);CHKERRQ(ierr); 1467d67d17b1SMatthew G. Knepley } 1468c58f1c22SToby Isaac nameSize = len; 1469c58f1c22SToby Isaac ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1470c58f1c22SToby Isaac ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1471580bdb30SBarry Smith if (!rank) {ierr = PetscArraycpy(name, lname, nameSize+1);CHKERRQ(ierr);} 1472c58f1c22SToby Isaac ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1473d67d17b1SMatthew G. Knepley ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr); 1474c58f1c22SToby Isaac ierr = PetscFree(name);CHKERRQ(ierr); 147577d236dfSMichael Lange /* Bcast defaultValue */ 147677d236dfSMichael Lange if (!rank) (*labelNew)->defaultValue = label->defaultValue; 147777d236dfSMichael Lange ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 147826c55118SMichael Lange /* Distribute stratum values over the SF and get the point mapping on the receiver */ 147926c55118SMichael Lange ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr); 14805cbdf6fcSMichael Lange /* Determine received stratum values and initialise new label*/ 1481e8f14785SLisandro Dalcin ierr = PetscHSetICreate(&stratumHash);CHKERRQ(ierr); 148226c55118SMichael Lange ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr); 1483e8f14785SLisandro Dalcin for (p = 0; p < size; ++p) {ierr = PetscHSetIAdd(stratumHash, leafStrata[p]);CHKERRQ(ierr);} 1484e8f14785SLisandro Dalcin ierr = PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata);CHKERRQ(ierr); 1485ad8374ffSToby Isaac ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr); 1486ad8374ffSToby Isaac for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE; 14875cbdf6fcSMichael Lange ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr); 14885cbdf6fcSMichael Lange /* Turn leafStrata into indices rather than stratum values */ 14895cbdf6fcSMichael Lange offset = 0; 1490e8f14785SLisandro Dalcin ierr = PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr); 1491a205f1fdSToby Isaac ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr); 149290e9b2aeSLisandro Dalcin for (s = 0; s < (*labelNew)->numStrata; ++s) { 149390e9b2aeSLisandro Dalcin ierr = PetscHMapISet((*labelNew)->hmap, (*labelNew)->stratumValues[s], s);CHKERRQ(ierr); 149490e9b2aeSLisandro Dalcin } 14955cbdf6fcSMichael Lange for (p = 0; p < size; ++p) { 1496231b9e6fSMatthew G. Knepley for (s = 0; s < (*labelNew)->numStrata; ++s) { 1497231b9e6fSMatthew G. Knepley if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;} 14985cbdf6fcSMichael Lange } 14995cbdf6fcSMichael Lange } 1500c58f1c22SToby Isaac /* Rebuild the point strata on the receiver */ 1501c58f1c22SToby Isaac ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr); 1502c58f1c22SToby Isaac ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1503c58f1c22SToby Isaac for (p=pStart; p<pEnd; p++) { 1504c58f1c22SToby Isaac ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1505c58f1c22SToby Isaac ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1506c58f1c22SToby Isaac for (s=0; s<dof; s++) { 1507c58f1c22SToby Isaac (*labelNew)->stratumSizes[leafStrata[offset+s]]++; 1508c58f1c22SToby Isaac } 1509c58f1c22SToby Isaac } 1510c58f1c22SToby Isaac ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr); 1511c58f1c22SToby Isaac ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr); 1512ad8374ffSToby Isaac ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr); 1513c58f1c22SToby Isaac for (s = 0; s < (*labelNew)->numStrata; ++s) { 1514e8f14785SLisandro Dalcin ierr = PetscHSetICreate(&(*labelNew)->ht[s]);CHKERRQ(ierr); 1515ad8374ffSToby Isaac ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr); 1516c58f1c22SToby Isaac } 1517c58f1c22SToby Isaac /* Insert points into new strata */ 1518c58f1c22SToby Isaac ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr); 1519c58f1c22SToby Isaac ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); 1520c58f1c22SToby Isaac for (p=pStart; p<pEnd; p++) { 1521c58f1c22SToby Isaac ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); 1522c58f1c22SToby Isaac ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); 1523c58f1c22SToby Isaac for (s=0; s<dof; s++) { 1524c58f1c22SToby Isaac stratum = leafStrata[offset+s]; 1525ad8374ffSToby Isaac points[stratum][strataIdx[stratum]++] = p; 1526c58f1c22SToby Isaac } 1527c58f1c22SToby Isaac } 1528ad8374ffSToby Isaac for (s = 0; s < (*labelNew)->numStrata; s++) { 1529ad8374ffSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr); 1530ad8374ffSToby Isaac ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr); 1531ad8374ffSToby Isaac } 1532ad8374ffSToby Isaac ierr = PetscFree(points);CHKERRQ(ierr); 1533e8f14785SLisandro Dalcin ierr = PetscHSetIDestroy(&stratumHash);CHKERRQ(ierr); 1534c58f1c22SToby Isaac ierr = PetscFree(leafStrata);CHKERRQ(ierr); 1535c58f1c22SToby Isaac ierr = PetscFree(strataIdx);CHKERRQ(ierr); 1536c58f1c22SToby Isaac ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 1537c58f1c22SToby Isaac PetscFunctionReturn(0); 1538c58f1c22SToby Isaac } 1539c58f1c22SToby Isaac 15407937d9ceSMichael Lange /*@ 15417937d9ceSMichael Lange DMLabelGather - Gather all label values from leafs into roots 15427937d9ceSMichael Lange 15435b5e7992SMatthew G. Knepley Collective on sf 15445b5e7992SMatthew G. Knepley 15457937d9ceSMichael Lange Input Parameters: 15467937d9ceSMichael Lange + label - the DMLabel 154784f0b6dfSMatthew G. Knepley - sf - the Star Forest point communication map 15487937d9ceSMichael Lange 154984f0b6dfSMatthew G. Knepley Output Parameters: 155084f0b6dfSMatthew G. Knepley . labelNew - the new DMLabel with localised leaf values 15517937d9ceSMichael Lange 15527937d9ceSMichael Lange Level: developer 15537937d9ceSMichael Lange 15547937d9ceSMichael Lange Note: This is the inverse operation to DMLabelDistribute. 15557937d9ceSMichael Lange 15567937d9ceSMichael Lange .seealso: DMLabelDistribute() 15577937d9ceSMichael Lange @*/ 15587937d9ceSMichael Lange PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew) 15597937d9ceSMichael Lange { 15607937d9ceSMichael Lange MPI_Comm comm; 15617937d9ceSMichael Lange PetscSection rootSection; 15627937d9ceSMichael Lange PetscSF sfLabel; 15637937d9ceSMichael Lange PetscSFNode *rootPoints, *leafPoints; 15647937d9ceSMichael Lange PetscInt p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset; 15657937d9ceSMichael Lange const PetscInt *rootDegree, *ilocal; 15667937d9ceSMichael Lange PetscInt *rootStrata; 1567d67d17b1SMatthew G. Knepley const char *lname; 15687937d9ceSMichael Lange char *name; 15697937d9ceSMichael Lange PetscInt nameSize; 15707937d9ceSMichael Lange size_t len = 0; 15719852e123SBarry Smith PetscMPIInt rank, size; 15727937d9ceSMichael Lange PetscErrorCode ierr; 15737937d9ceSMichael Lange 15747937d9ceSMichael Lange PetscFunctionBegin; 1575d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1576d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 15777937d9ceSMichael Lange ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); 15787937d9ceSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 15799852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 15807937d9ceSMichael Lange /* Bcast name */ 1581d67d17b1SMatthew G. Knepley if (!rank) { 1582d67d17b1SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr); 1583d67d17b1SMatthew G. Knepley ierr = PetscStrlen(lname, &len);CHKERRQ(ierr); 1584d67d17b1SMatthew G. Knepley } 15857937d9ceSMichael Lange nameSize = len; 15867937d9ceSMichael Lange ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 15877937d9ceSMichael Lange ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); 1588580bdb30SBarry Smith if (!rank) {ierr = PetscArraycpy(name, lname, nameSize+1);CHKERRQ(ierr);} 15897937d9ceSMichael Lange ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1590d67d17b1SMatthew G. Knepley ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr); 15917937d9ceSMichael Lange ierr = PetscFree(name);CHKERRQ(ierr); 15927937d9ceSMichael Lange /* Gather rank/index pairs of leaves into local roots to build 15937937d9ceSMichael Lange an inverse, multi-rooted SF. Note that this ignores local leaf 15947937d9ceSMichael Lange indexing due to the use of the multiSF in PetscSFGather. */ 15957937d9ceSMichael Lange ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr); 1596dc53bc9bSMatthew G. Knepley ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr); 1597dc53bc9bSMatthew G. Knepley for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1; 15987937d9ceSMichael Lange for (p = 0; p < nleaves; p++) { 15998212dd46SStefano Zampini PetscInt ilp = ilocal ? ilocal[p] : p; 16008212dd46SStefano Zampini 16018212dd46SStefano Zampini leafPoints[ilp].index = ilp; 16028212dd46SStefano Zampini leafPoints[ilp].rank = rank; 16037937d9ceSMichael Lange } 16047937d9ceSMichael Lange ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr); 16057937d9ceSMichael Lange ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr); 16067937d9ceSMichael Lange for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p]; 16077937d9ceSMichael Lange ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr); 16087937d9ceSMichael Lange ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 16097937d9ceSMichael Lange ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr); 16107937d9ceSMichael Lange ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr); 16117937d9ceSMichael Lange ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 16127937d9ceSMichael Lange /* Migrate label over inverted SF to pull stratum values at leaves into roots. */ 16137937d9ceSMichael Lange ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr); 16147937d9ceSMichael Lange /* Rebuild the point strata on the receiver */ 16157937d9ceSMichael Lange for (p = 0, idx = 0; p < nroots; p++) { 16167937d9ceSMichael Lange for (d = 0; d < rootDegree[p]; d++) { 16177937d9ceSMichael Lange ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr); 16187937d9ceSMichael Lange ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr); 16197937d9ceSMichael Lange for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);} 16207937d9ceSMichael Lange } 16217937d9ceSMichael Lange idx += rootDegree[p]; 16227937d9ceSMichael Lange } 162377e0c0e7SMichael Lange ierr = PetscFree(leafPoints);CHKERRQ(ierr); 162477e0c0e7SMichael Lange ierr = PetscFree(rootStrata);CHKERRQ(ierr); 162577e0c0e7SMichael Lange ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 162677e0c0e7SMichael Lange ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr); 16277937d9ceSMichael Lange PetscFunctionReturn(0); 16287937d9ceSMichael Lange } 16297937d9ceSMichael Lange 163084f0b6dfSMatthew G. Knepley /*@ 163184f0b6dfSMatthew G. Knepley DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label 163284f0b6dfSMatthew G. Knepley 16335b5e7992SMatthew G. Knepley Not collective 16345b5e7992SMatthew G. Knepley 163584f0b6dfSMatthew G. Knepley Input Parameter: 163684f0b6dfSMatthew G. Knepley . label - the DMLabel 163784f0b6dfSMatthew G. Knepley 163884f0b6dfSMatthew G. Knepley Output Parameters: 163984f0b6dfSMatthew G. Knepley + section - the section giving offsets for each stratum 164084f0b6dfSMatthew G. Knepley - is - An IS containing all the label points 164184f0b6dfSMatthew G. Knepley 164284f0b6dfSMatthew G. Knepley Level: developer 164384f0b6dfSMatthew G. Knepley 164484f0b6dfSMatthew G. Knepley .seealso: DMLabelDistribute() 164584f0b6dfSMatthew G. Knepley @*/ 1646c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is) 1647c58f1c22SToby Isaac { 1648c58f1c22SToby Isaac IS vIS; 1649c58f1c22SToby Isaac const PetscInt *values; 1650c58f1c22SToby Isaac PetscInt *points; 1651c58f1c22SToby Isaac PetscInt nV, vS = 0, vE = 0, v, N; 1652c58f1c22SToby Isaac PetscErrorCode ierr; 1653c58f1c22SToby Isaac 1654c58f1c22SToby Isaac PetscFunctionBegin; 1655d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1); 1656c58f1c22SToby Isaac ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr); 1657c58f1c22SToby Isaac ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr); 1658c58f1c22SToby Isaac ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr); 1659c58f1c22SToby Isaac if (nV) {vS = values[0]; vE = values[0]+1;} 1660c58f1c22SToby Isaac for (v = 1; v < nV; ++v) { 1661c58f1c22SToby Isaac vS = PetscMin(vS, values[v]); 1662c58f1c22SToby Isaac vE = PetscMax(vE, values[v]+1); 1663c58f1c22SToby Isaac } 1664c58f1c22SToby Isaac ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr); 1665c58f1c22SToby Isaac ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr); 1666c58f1c22SToby Isaac for (v = 0; v < nV; ++v) { 1667c58f1c22SToby Isaac PetscInt n; 1668c58f1c22SToby Isaac 1669c58f1c22SToby Isaac ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr); 1670c58f1c22SToby Isaac ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr); 1671c58f1c22SToby Isaac } 1672c58f1c22SToby Isaac ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 1673c58f1c22SToby Isaac ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr); 1674c58f1c22SToby Isaac ierr = PetscMalloc1(N, &points);CHKERRQ(ierr); 1675c58f1c22SToby Isaac for (v = 0; v < nV; ++v) { 1676c58f1c22SToby Isaac IS is; 1677c58f1c22SToby Isaac const PetscInt *spoints; 1678c58f1c22SToby Isaac PetscInt dof, off, p; 1679c58f1c22SToby Isaac 1680c58f1c22SToby Isaac ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr); 1681c58f1c22SToby Isaac ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr); 1682c58f1c22SToby Isaac ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr); 1683c58f1c22SToby Isaac ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr); 1684c58f1c22SToby Isaac for (p = 0; p < dof; ++p) points[off+p] = spoints[p]; 1685c58f1c22SToby Isaac ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr); 1686c58f1c22SToby Isaac ierr = ISDestroy(&is);CHKERRQ(ierr); 1687c58f1c22SToby Isaac } 1688c58f1c22SToby Isaac ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr); 1689c58f1c22SToby Isaac ierr = ISDestroy(&vIS);CHKERRQ(ierr); 1690c58f1c22SToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr); 1691c58f1c22SToby Isaac PetscFunctionReturn(0); 1692c58f1c22SToby Isaac } 1693c58f1c22SToby Isaac 169484f0b6dfSMatthew G. Knepley /*@ 1695c58f1c22SToby Isaac PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using 1696c58f1c22SToby Isaac the local section and an SF describing the section point overlap. 1697c58f1c22SToby Isaac 16985b5e7992SMatthew G. Knepley Collective on sf 16995b5e7992SMatthew G. Knepley 1700c58f1c22SToby Isaac Input Parameters: 1701c58f1c22SToby Isaac + s - The PetscSection for the local field layout 1702c58f1c22SToby Isaac . sf - The SF describing parallel layout of the section points 1703c58f1c22SToby Isaac . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs 1704c58f1c22SToby Isaac . label - The label specifying the points 1705c58f1c22SToby Isaac - labelValue - The label stratum specifying the points 1706c58f1c22SToby Isaac 1707c58f1c22SToby Isaac Output Parameter: 1708c58f1c22SToby Isaac . gsection - The PetscSection for the global field layout 1709c58f1c22SToby Isaac 1710c58f1c22SToby Isaac Note: This gives negative sizes and offsets to points not owned by this process 1711c58f1c22SToby Isaac 1712c58f1c22SToby Isaac Level: developer 1713c58f1c22SToby Isaac 1714c58f1c22SToby Isaac .seealso: PetscSectionCreate() 1715c58f1c22SToby Isaac @*/ 1716c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection) 1717c58f1c22SToby Isaac { 1718c58f1c22SToby Isaac PetscInt *neg = NULL, *tmpOff = NULL; 1719c58f1c22SToby Isaac PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots; 1720c58f1c22SToby Isaac PetscErrorCode ierr; 1721c58f1c22SToby Isaac 1722c58f1c22SToby Isaac PetscFunctionBegin; 1723d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1724d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 1725d67d17b1SMatthew G. Knepley PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4); 1726c58f1c22SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr); 1727c58f1c22SToby Isaac ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 1728c58f1c22SToby Isaac ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr); 1729c58f1c22SToby Isaac ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 1730c58f1c22SToby Isaac if (nroots >= 0) { 1731c58f1c22SToby Isaac if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart); 1732c58f1c22SToby Isaac ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr); 1733c58f1c22SToby Isaac if (nroots > pEnd-pStart) { 1734c58f1c22SToby Isaac ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr); 1735c58f1c22SToby Isaac } else { 1736c58f1c22SToby Isaac tmpOff = &(*gsection)->atlasDof[-pStart]; 1737c58f1c22SToby Isaac } 1738c58f1c22SToby Isaac } 1739c58f1c22SToby Isaac /* Mark ghost points with negative dof */ 1740c58f1c22SToby Isaac for (p = pStart; p < pEnd; ++p) { 1741c58f1c22SToby Isaac PetscInt value; 1742c58f1c22SToby Isaac 1743c58f1c22SToby Isaac ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr); 1744c58f1c22SToby Isaac if (value != labelValue) continue; 1745c58f1c22SToby Isaac ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr); 1746c58f1c22SToby Isaac ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr); 1747c58f1c22SToby Isaac ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 1748c58f1c22SToby Isaac if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);} 1749c58f1c22SToby Isaac if (neg) neg[p] = -(dof+1); 1750c58f1c22SToby Isaac } 1751c58f1c22SToby Isaac ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr); 1752c58f1c22SToby Isaac if (nroots >= 0) { 1753c58f1c22SToby Isaac ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1754c58f1c22SToby Isaac ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1755c58f1c22SToby Isaac if (nroots > pEnd-pStart) { 1756c58f1c22SToby Isaac for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];} 1757c58f1c22SToby Isaac } 1758c58f1c22SToby Isaac } 1759c58f1c22SToby Isaac /* Calculate new sizes, get proccess offset, and calculate point offsets */ 1760c58f1c22SToby Isaac for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1761c58f1c22SToby Isaac cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0; 1762c58f1c22SToby Isaac (*gsection)->atlasOff[p] = off; 1763c58f1c22SToby Isaac off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0; 1764c58f1c22SToby Isaac } 1765c58f1c22SToby Isaac ierr = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr); 1766c58f1c22SToby Isaac globalOff -= off; 1767c58f1c22SToby Isaac for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1768c58f1c22SToby Isaac (*gsection)->atlasOff[p] += globalOff; 1769c58f1c22SToby Isaac if (neg) neg[p] = -((*gsection)->atlasOff[p]+1); 1770c58f1c22SToby Isaac } 1771c58f1c22SToby Isaac /* Put in negative offsets for ghost points */ 1772c58f1c22SToby Isaac if (nroots >= 0) { 1773c58f1c22SToby Isaac ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1774c58f1c22SToby Isaac ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr); 1775c58f1c22SToby Isaac if (nroots > pEnd-pStart) { 1776c58f1c22SToby Isaac for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];} 1777c58f1c22SToby Isaac } 1778c58f1c22SToby Isaac } 1779c58f1c22SToby Isaac if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);} 1780c58f1c22SToby Isaac ierr = PetscFree(neg);CHKERRQ(ierr); 1781c58f1c22SToby Isaac PetscFunctionReturn(0); 1782c58f1c22SToby Isaac } 1783c58f1c22SToby Isaac 17845fdea053SToby Isaac typedef struct _n_PetscSectionSym_Label 17855fdea053SToby Isaac { 17865fdea053SToby Isaac DMLabel label; 17875fdea053SToby Isaac PetscCopyMode *modes; 17885fdea053SToby Isaac PetscInt *sizes; 17895fdea053SToby Isaac const PetscInt ***perms; 17905fdea053SToby Isaac const PetscScalar ***rots; 17915fdea053SToby Isaac PetscInt (*minMaxOrients)[2]; 17925fdea053SToby Isaac PetscInt numStrata; /* numStrata is only increasing, functions as a state */ 17935fdea053SToby Isaac } PetscSectionSym_Label; 17945fdea053SToby Isaac 17955fdea053SToby Isaac static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym) 17965fdea053SToby Isaac { 17975fdea053SToby Isaac PetscInt i, j; 17985fdea053SToby Isaac PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data; 17995fdea053SToby Isaac PetscErrorCode ierr; 18005fdea053SToby Isaac 18015fdea053SToby Isaac PetscFunctionBegin; 18025fdea053SToby Isaac for (i = 0; i <= sl->numStrata; i++) { 18035fdea053SToby Isaac if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) { 18045fdea053SToby Isaac for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 18055fdea053SToby Isaac if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);} 18065fdea053SToby Isaac if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);} 18075fdea053SToby Isaac } 18085fdea053SToby Isaac if (sl->perms[i]) { 18095fdea053SToby Isaac const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]]; 18105fdea053SToby Isaac 18115fdea053SToby Isaac ierr = PetscFree(perms);CHKERRQ(ierr); 18125fdea053SToby Isaac } 18135fdea053SToby Isaac if (sl->rots[i]) { 18145fdea053SToby Isaac const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]]; 18155fdea053SToby Isaac 18165fdea053SToby Isaac ierr = PetscFree(rots);CHKERRQ(ierr); 18175fdea053SToby Isaac } 18185fdea053SToby Isaac } 18195fdea053SToby Isaac } 18205fdea053SToby Isaac ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr); 18215fdea053SToby Isaac ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr); 18225fdea053SToby Isaac sl->numStrata = 0; 18235fdea053SToby Isaac PetscFunctionReturn(0); 18245fdea053SToby Isaac } 18255fdea053SToby Isaac 18265fdea053SToby Isaac static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym) 18275fdea053SToby Isaac { 18285fdea053SToby Isaac PetscErrorCode ierr; 18295fdea053SToby Isaac 18305fdea053SToby Isaac PetscFunctionBegin; 18315fdea053SToby Isaac ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr); 18325fdea053SToby Isaac ierr = PetscFree(sym->data);CHKERRQ(ierr); 18335fdea053SToby Isaac PetscFunctionReturn(0); 18345fdea053SToby Isaac } 18355fdea053SToby Isaac 18365fdea053SToby Isaac static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer) 18375fdea053SToby Isaac { 18385fdea053SToby Isaac PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data; 18395fdea053SToby Isaac PetscBool isAscii; 18405fdea053SToby Isaac DMLabel label = sl->label; 1841d67d17b1SMatthew G. Knepley const char *name; 18425fdea053SToby Isaac PetscErrorCode ierr; 18435fdea053SToby Isaac 18445fdea053SToby Isaac PetscFunctionBegin; 18455fdea053SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr); 18465fdea053SToby Isaac if (isAscii) { 18475fdea053SToby Isaac PetscInt i, j, k; 18485fdea053SToby Isaac PetscViewerFormat format; 18495fdea053SToby Isaac 18505fdea053SToby Isaac ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 18515fdea053SToby Isaac if (label) { 18525fdea053SToby Isaac ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 18535fdea053SToby Isaac if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 18545fdea053SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 18555fdea053SToby Isaac ierr = DMLabelView(label, viewer);CHKERRQ(ierr); 18565fdea053SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 18575fdea053SToby Isaac } else { 1858d67d17b1SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr); 1859d67d17b1SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer," Label '%s'\n",name);CHKERRQ(ierr); 18605fdea053SToby Isaac } 18615fdea053SToby Isaac } else { 18625fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr); 18635fdea053SToby Isaac } 18645fdea053SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 18655fdea053SToby Isaac for (i = 0; i <= sl->numStrata; i++) { 18665fdea053SToby Isaac PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue; 18675fdea053SToby Isaac 18685fdea053SToby Isaac if (!(sl->perms[i] || sl->rots[i])) { 18695fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr); 18705fdea053SToby Isaac } else { 18715fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr); 18725fdea053SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 18735fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr); 18745fdea053SToby Isaac if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 18755fdea053SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 18765fdea053SToby Isaac for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) { 18775fdea053SToby Isaac if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) { 18785fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr); 18795fdea053SToby Isaac } else { 18805fdea053SToby Isaac PetscInt tab; 18815fdea053SToby Isaac 18825fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr); 18835fdea053SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 18845fdea053SToby Isaac ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr); 18855fdea053SToby Isaac if (sl->perms[i] && sl->perms[i][j]) { 18865fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr); 18875fdea053SToby Isaac ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr); 18885fdea053SToby Isaac for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);} 18895fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 18905fdea053SToby Isaac ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr); 18915fdea053SToby Isaac } 18925fdea053SToby Isaac if (sl->rots[i] && sl->rots[i][j]) { 18935fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Rotations: ");CHKERRQ(ierr); 18945fdea053SToby Isaac ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr); 18955fdea053SToby Isaac #if defined(PETSC_USE_COMPLEX) 18965fdea053SToby 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);} 18975fdea053SToby Isaac #else 18985fdea053SToby Isaac for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);} 18995fdea053SToby Isaac #endif 19005fdea053SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 19015fdea053SToby Isaac ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr); 19025fdea053SToby Isaac } 19035fdea053SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 19045fdea053SToby Isaac } 19055fdea053SToby Isaac } 19065fdea053SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 19075fdea053SToby Isaac } 19085fdea053SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 19095fdea053SToby Isaac } 19105fdea053SToby Isaac } 19115fdea053SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 19125fdea053SToby Isaac } 19135fdea053SToby Isaac PetscFunctionReturn(0); 19145fdea053SToby Isaac } 19155fdea053SToby Isaac 19165fdea053SToby Isaac /*@ 19175fdea053SToby Isaac PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries 19185fdea053SToby Isaac 19195fdea053SToby Isaac Logically collective on sym 19205fdea053SToby Isaac 19215fdea053SToby Isaac Input parameters: 19225fdea053SToby Isaac + sym - the section symmetries 19235fdea053SToby Isaac - label - the DMLabel describing the types of points 19245fdea053SToby Isaac 19255fdea053SToby Isaac Level: developer: 19265fdea053SToby Isaac 19275fdea053SToby Isaac .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms() 19285fdea053SToby Isaac @*/ 19295fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label) 19305fdea053SToby Isaac { 19315fdea053SToby Isaac PetscSectionSym_Label *sl; 19325fdea053SToby Isaac PetscErrorCode ierr; 19335fdea053SToby Isaac 19345fdea053SToby Isaac PetscFunctionBegin; 19355fdea053SToby Isaac PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1); 19365fdea053SToby Isaac sl = (PetscSectionSym_Label *) sym->data; 19375fdea053SToby Isaac if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);} 19385fdea053SToby Isaac if (label) { 19395fdea053SToby Isaac sl->label = label; 1940d67d17b1SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) label);CHKERRQ(ierr); 19415fdea053SToby Isaac ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr); 19421a834cf9SToby 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); 19431a834cf9SToby Isaac ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr); 19441a834cf9SToby Isaac ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr); 19451a834cf9SToby Isaac ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr); 19461a834cf9SToby Isaac ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr); 19471a834cf9SToby Isaac ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr); 19485fdea053SToby Isaac } 19495fdea053SToby Isaac PetscFunctionReturn(0); 19505fdea053SToby Isaac } 19515fdea053SToby Isaac 19525fdea053SToby Isaac /*@C 19535fdea053SToby Isaac PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum 19545fdea053SToby Isaac 19555b5e7992SMatthew G. Knepley Logically collective on sym 19565fdea053SToby Isaac 19575fdea053SToby Isaac InputParameters: 19585b5e7992SMatthew G. Knepley + sym - the section symmetries 19595fdea053SToby Isaac . stratum - the stratum value in the label that we are assigning symmetries for 19605fdea053SToby Isaac . size - the number of dofs for points in the stratum of the label 19615fdea053SToby Isaac . minOrient - the smallest orientation for a point in this stratum 19625fdea053SToby Isaac . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient)) 19635fdea053SToby Isaac . mode - how sym should copy the perms and rots arrays 19645fdea053SToby Isaac . perms - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation. A NULL permutation is the identity 1965a2b725a8SWilliam Gropp - 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 19665fdea053SToby Isaac 19675fdea053SToby Isaac Level: developer 19685fdea053SToby Isaac 19695fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel() 19705fdea053SToby Isaac @*/ 19715fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots) 19725fdea053SToby Isaac { 19735fdea053SToby Isaac PetscSectionSym_Label *sl; 1974d67d17b1SMatthew G. Knepley const char *name; 1975d67d17b1SMatthew G. Knepley PetscInt i, j, k; 19765fdea053SToby Isaac PetscErrorCode ierr; 19775fdea053SToby Isaac 19785fdea053SToby Isaac PetscFunctionBegin; 19795fdea053SToby Isaac PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1); 19805fdea053SToby Isaac sl = (PetscSectionSym_Label *) sym->data; 19815fdea053SToby Isaac if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet"); 19825fdea053SToby Isaac for (i = 0; i <= sl->numStrata; i++) { 19835fdea053SToby Isaac PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue; 19845fdea053SToby Isaac 19855fdea053SToby Isaac if (stratum == value) break; 19865fdea053SToby Isaac } 1987d67d17b1SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr); 1988d67d17b1SMatthew G. Knepley if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,name); 19895fdea053SToby Isaac sl->sizes[i] = size; 19905fdea053SToby Isaac sl->modes[i] = mode; 19915fdea053SToby Isaac sl->minMaxOrients[i][0] = minOrient; 19925fdea053SToby Isaac sl->minMaxOrients[i][1] = maxOrient; 19935fdea053SToby Isaac if (mode == PETSC_COPY_VALUES) { 19945fdea053SToby Isaac if (perms) { 19955fdea053SToby Isaac PetscInt **ownPerms; 19965fdea053SToby Isaac 19975fdea053SToby Isaac ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr); 19985fdea053SToby Isaac for (j = 0; j < maxOrient-minOrient; j++) { 19995fdea053SToby Isaac if (perms[j]) { 20005fdea053SToby Isaac ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr); 20015fdea053SToby Isaac for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];} 20025fdea053SToby Isaac } 20035fdea053SToby Isaac } 20045fdea053SToby Isaac sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient]; 20055fdea053SToby Isaac } 20065fdea053SToby Isaac if (rots) { 20075fdea053SToby Isaac PetscScalar **ownRots; 20085fdea053SToby Isaac 20095fdea053SToby Isaac ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr); 20105fdea053SToby Isaac for (j = 0; j < maxOrient-minOrient; j++) { 20115fdea053SToby Isaac if (rots[j]) { 20125fdea053SToby Isaac ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr); 20135fdea053SToby Isaac for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];} 20145fdea053SToby Isaac } 20155fdea053SToby Isaac } 20165fdea053SToby Isaac sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient]; 20175fdea053SToby Isaac } 20185fdea053SToby Isaac } else { 20195fdea053SToby Isaac sl->perms[i] = perms ? &perms[-minOrient] : NULL; 20205fdea053SToby Isaac sl->rots[i] = rots ? &rots[-minOrient] : NULL; 20215fdea053SToby Isaac } 20225fdea053SToby Isaac PetscFunctionReturn(0); 20235fdea053SToby Isaac } 20245fdea053SToby Isaac 20255fdea053SToby Isaac static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots) 20265fdea053SToby Isaac { 20275fdea053SToby Isaac PetscInt i, j, numStrata; 20285fdea053SToby Isaac PetscSectionSym_Label *sl; 20295fdea053SToby Isaac DMLabel label; 20305fdea053SToby Isaac PetscErrorCode ierr; 20315fdea053SToby Isaac 20325fdea053SToby Isaac PetscFunctionBegin; 20335fdea053SToby Isaac sl = (PetscSectionSym_Label *) sym->data; 20345fdea053SToby Isaac numStrata = sl->numStrata; 20355fdea053SToby Isaac label = sl->label; 20365fdea053SToby Isaac for (i = 0; i < numPoints; i++) { 20375fdea053SToby Isaac PetscInt point = points[2*i]; 20385fdea053SToby Isaac PetscInt ornt = points[2*i+1]; 20395fdea053SToby Isaac 20405fdea053SToby Isaac for (j = 0; j < numStrata; j++) { 20415fdea053SToby Isaac if (label->validIS[j]) { 20425fdea053SToby Isaac PetscInt k; 20435fdea053SToby Isaac 20445fdea053SToby Isaac ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr); 20455fdea053SToby Isaac if (k >= 0) break; 20465fdea053SToby Isaac } else { 20475fdea053SToby Isaac PetscBool has; 20485fdea053SToby Isaac 2049b9907514SLisandro Dalcin ierr = PetscHSetIHas(label->ht[j], point, &has);CHKERRQ(ierr); 20505fdea053SToby Isaac if (has) break; 20515fdea053SToby Isaac } 20525fdea053SToby Isaac } 20535fdea053SToby 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); 20545fdea053SToby Isaac if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;} 20555fdea053SToby Isaac if (rots) {rots[i] = sl->rots[j] ? sl->rots[j][ornt] : NULL;} 20565fdea053SToby Isaac } 20575fdea053SToby Isaac PetscFunctionReturn(0); 20585fdea053SToby Isaac } 20595fdea053SToby Isaac 20605fdea053SToby Isaac PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym) 20615fdea053SToby Isaac { 20625fdea053SToby Isaac PetscSectionSym_Label *sl; 20635fdea053SToby Isaac PetscErrorCode ierr; 20645fdea053SToby Isaac 20655fdea053SToby Isaac PetscFunctionBegin; 20665fdea053SToby Isaac ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr); 20675fdea053SToby Isaac sym->ops->getpoints = PetscSectionSymGetPoints_Label; 20685fdea053SToby Isaac sym->ops->view = PetscSectionSymView_Label; 20695fdea053SToby Isaac sym->ops->destroy = PetscSectionSymDestroy_Label; 20705fdea053SToby Isaac sym->data = (void *) sl; 20715fdea053SToby Isaac PetscFunctionReturn(0); 20725fdea053SToby Isaac } 20735fdea053SToby Isaac 20745fdea053SToby Isaac /*@ 20755fdea053SToby Isaac PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label 20765fdea053SToby Isaac 20775fdea053SToby Isaac Collective 20785fdea053SToby Isaac 20795fdea053SToby Isaac Input Parameters: 20805fdea053SToby Isaac + comm - the MPI communicator for the new symmetry 20815fdea053SToby Isaac - label - the label defining the strata 20825fdea053SToby Isaac 20835fdea053SToby Isaac Output Parameters: 20845fdea053SToby Isaac . sym - the section symmetries 20855fdea053SToby Isaac 20865fdea053SToby Isaac Level: developer 20875fdea053SToby Isaac 20885fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms() 20895fdea053SToby Isaac @*/ 20905fdea053SToby Isaac PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym) 20915fdea053SToby Isaac { 20925fdea053SToby Isaac PetscErrorCode ierr; 20935fdea053SToby Isaac 20945fdea053SToby Isaac PetscFunctionBegin; 20955fdea053SToby Isaac ierr = DMInitializePackage();CHKERRQ(ierr); 20965fdea053SToby Isaac ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr); 20975fdea053SToby Isaac ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr); 20985fdea053SToby Isaac ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr); 20995fdea053SToby Isaac PetscFunctionReturn(0); 21005fdea053SToby Isaac } 2101