xref: /petsc/src/dm/label/dmlabel.c (revision ba2698f1808e61e430b749a0784ecbf435cb3659)
15fdea053SToby Isaac #include <petscdm.h>
2c58f1c22SToby Isaac #include <petsc/private/dmlabelimpl.h>   /*I      "petscdmlabel.h"   I*/
3ea844a1aSMatthew Knepley #include <petsc/private/sectionimpl.h>   /*I      "petscsection.h"   I*/
4c58f1c22SToby Isaac #include <petscsf.h>
5ea844a1aSMatthew Knepley #include <petscsection.h>
6c58f1c22SToby Isaac 
7c58f1c22SToby Isaac /*@C
8c58f1c22SToby Isaac   DMLabelCreate - Create a DMLabel object, which is a multimap
9c58f1c22SToby Isaac 
105b5e7992SMatthew G. Knepley   Collective
115b5e7992SMatthew G. Knepley 
12d67d17b1SMatthew G. Knepley   Input parameters:
13d67d17b1SMatthew G. Knepley + comm - The communicator, usually PETSC_COMM_SELF
14d67d17b1SMatthew G. Knepley - name - The label name
15c58f1c22SToby Isaac 
16c58f1c22SToby Isaac   Output parameter:
17c58f1c22SToby Isaac . label - The DMLabel
18c58f1c22SToby Isaac 
19c58f1c22SToby Isaac   Level: beginner
20c58f1c22SToby Isaac 
2105ab7a9fSVaclav Hapla   Notes:
2205ab7a9fSVaclav Hapla   The label name is actually usual PetscObject name.
2305ab7a9fSVaclav Hapla   One can get/set it with PetscObjectGetName()/PetscObjectSetName().
2405ab7a9fSVaclav Hapla 
25c58f1c22SToby Isaac .seealso: DMLabelDestroy()
26c58f1c22SToby Isaac @*/
27d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelCreate(MPI_Comm comm, const char name[], DMLabel *label)
28c58f1c22SToby Isaac {
29c58f1c22SToby Isaac   PetscErrorCode ierr;
30c58f1c22SToby Isaac 
31c58f1c22SToby Isaac   PetscFunctionBegin;
32d67d17b1SMatthew G. Knepley   PetscValidPointer(label,2);
33d67d17b1SMatthew G. Knepley   ierr = DMInitializePackage();CHKERRQ(ierr);
34c58f1c22SToby Isaac 
35d67d17b1SMatthew G. Knepley   ierr = PetscHeaderCreate(*label,DMLABEL_CLASSID,"DMLabel","DMLabel","DM",comm,DMLabelDestroy,DMLabelView);CHKERRQ(ierr);
36d67d17b1SMatthew G. Knepley 
37c58f1c22SToby Isaac   (*label)->numStrata      = 0;
385aa44df4SToby Isaac   (*label)->defaultValue   = -1;
39c58f1c22SToby Isaac   (*label)->stratumValues  = NULL;
40ad8374ffSToby Isaac   (*label)->validIS        = NULL;
41c58f1c22SToby Isaac   (*label)->stratumSizes   = NULL;
42c58f1c22SToby Isaac   (*label)->points         = NULL;
43c58f1c22SToby Isaac   (*label)->ht             = NULL;
44c58f1c22SToby Isaac   (*label)->pStart         = -1;
45c58f1c22SToby Isaac   (*label)->pEnd           = -1;
46c58f1c22SToby Isaac   (*label)->bt             = NULL;
47b9907514SLisandro Dalcin   ierr = PetscHMapICreate(&(*label)->hmap);CHKERRQ(ierr);
48d67d17b1SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *label, name);CHKERRQ(ierr);
49c58f1c22SToby Isaac   PetscFunctionReturn(0);
50c58f1c22SToby Isaac }
51c58f1c22SToby Isaac 
52c58f1c22SToby Isaac /*
53c58f1c22SToby Isaac   DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format
54c58f1c22SToby Isaac 
555b5e7992SMatthew G. Knepley   Not collective
565b5e7992SMatthew G. Knepley 
57c58f1c22SToby Isaac   Input parameter:
58c58f1c22SToby Isaac + label - The DMLabel
59c58f1c22SToby Isaac - v - The stratum value
60c58f1c22SToby Isaac 
61c58f1c22SToby Isaac   Output parameter:
62c58f1c22SToby Isaac . label - The DMLabel with stratum in sorted list format
63c58f1c22SToby Isaac 
64c58f1c22SToby Isaac   Level: developer
65c58f1c22SToby Isaac 
66c58f1c22SToby Isaac .seealso: DMLabelCreate()
67c58f1c22SToby Isaac */
68c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v)
69c58f1c22SToby Isaac {
70277ea44aSLisandro Dalcin   IS             is;
71b9907514SLisandro Dalcin   PetscInt       off = 0, *pointArray, p;
72c58f1c22SToby Isaac   PetscErrorCode ierr;
73c58f1c22SToby Isaac 
74b9907514SLisandro Dalcin   if (PetscLikely(v >= 0 && v < label->numStrata) && label->validIS[v]) return 0;
75c58f1c22SToby Isaac   PetscFunctionBegin;
760c3c4a36SLisandro 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);
77e8f14785SLisandro Dalcin   ierr = PetscHSetIGetSize(label->ht[v], &label->stratumSizes[v]);CHKERRQ(ierr);
78ad8374ffSToby Isaac   ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr);
79e8f14785SLisandro Dalcin   ierr = PetscHSetIGetElems(label->ht[v], &off, pointArray);CHKERRQ(ierr);
80b9907514SLisandro Dalcin   ierr = PetscHSetIClear(label->ht[v]);CHKERRQ(ierr);
81ad8374ffSToby Isaac   ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr);
82c58f1c22SToby Isaac   if (label->bt) {
83c58f1c22SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
84ad8374ffSToby Isaac       const PetscInt point = pointArray[p];
85c58f1c22SToby 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);
86c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
87c58f1c22SToby Isaac     }
88c58f1c22SToby Isaac   }
89*ba2698f1SMatthew G. Knepley   if (label->stratumSizes[v] > 0 && pointArray[label->stratumSizes[v]-1] == pointArray[0] + label->stratumSizes[v]-1) {
90*ba2698f1SMatthew G. Knepley     ierr = ISCreateStride(PETSC_COMM_SELF, label->stratumSizes[v], pointArray[0], 1, &is);CHKERRQ(ierr);
91*ba2698f1SMatthew G. Knepley     ierr = PetscFree(pointArray);CHKERRQ(ierr);
92*ba2698f1SMatthew G. Knepley   } else {
93277ea44aSLisandro Dalcin     ierr = ISCreateGeneral(PETSC_COMM_SELF, label->stratumSizes[v], pointArray, PETSC_OWN_POINTER, &is);CHKERRQ(ierr);
94*ba2698f1SMatthew G. Knepley   }
95277ea44aSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject) is, "indices");CHKERRQ(ierr);
96277ea44aSLisandro Dalcin   label->points[v]  = is;
97ad8374ffSToby Isaac   label->validIS[v] = PETSC_TRUE;
98d67d17b1SMatthew G. Knepley   ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
99c58f1c22SToby Isaac   PetscFunctionReturn(0);
100c58f1c22SToby Isaac }
101c58f1c22SToby Isaac 
102c58f1c22SToby Isaac /*
103c58f1c22SToby Isaac   DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format
104c58f1c22SToby Isaac 
1055b5e7992SMatthew G. Knepley   Not collective
1065b5e7992SMatthew G. Knepley 
107c58f1c22SToby Isaac   Input parameter:
108c58f1c22SToby Isaac . label - The DMLabel
109c58f1c22SToby Isaac 
110c58f1c22SToby Isaac   Output parameter:
111c58f1c22SToby Isaac . label - The DMLabel with all strata in sorted list format
112c58f1c22SToby Isaac 
113c58f1c22SToby Isaac   Level: developer
114c58f1c22SToby Isaac 
115c58f1c22SToby Isaac .seealso: DMLabelCreate()
116c58f1c22SToby Isaac */
117c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label)
118c58f1c22SToby Isaac {
119c58f1c22SToby Isaac   PetscInt       v;
120c58f1c22SToby Isaac   PetscErrorCode ierr;
121c58f1c22SToby Isaac 
122c58f1c22SToby Isaac   PetscFunctionBegin;
123c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; v++){
124c58f1c22SToby Isaac     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
125c58f1c22SToby Isaac   }
126c58f1c22SToby Isaac   PetscFunctionReturn(0);
127c58f1c22SToby Isaac }
128c58f1c22SToby Isaac 
129c58f1c22SToby Isaac /*
130c58f1c22SToby Isaac   DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format
131c58f1c22SToby Isaac 
1325b5e7992SMatthew G. Knepley   Not collective
1335b5e7992SMatthew G. Knepley 
134c58f1c22SToby Isaac   Input parameter:
135c58f1c22SToby Isaac + label - The DMLabel
136c58f1c22SToby Isaac - v - The stratum value
137c58f1c22SToby Isaac 
138c58f1c22SToby Isaac   Output parameter:
139c58f1c22SToby Isaac . label - The DMLabel with stratum in hash format
140c58f1c22SToby Isaac 
141c58f1c22SToby Isaac   Level: developer
142c58f1c22SToby Isaac 
143c58f1c22SToby Isaac .seealso: DMLabelCreate()
144c58f1c22SToby Isaac */
145c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v)
146c58f1c22SToby Isaac {
147c58f1c22SToby Isaac   PetscInt       p;
148ad8374ffSToby Isaac   const PetscInt *points;
149c58f1c22SToby Isaac   PetscErrorCode ierr;
150c58f1c22SToby Isaac 
151b9907514SLisandro Dalcin   if (PetscLikely(v >= 0 && v < label->numStrata) && !label->validIS[v]) return 0;
152c58f1c22SToby Isaac   PetscFunctionBegin;
1530c3c4a36SLisandro 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);
154ad8374ffSToby Isaac   if (label->points[v]) {
155ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
156e8f14785SLisandro Dalcin     for (p = 0; p < label->stratumSizes[v]; ++p) {
157e8f14785SLisandro Dalcin       ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);
158e8f14785SLisandro Dalcin     }
159ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
160ad8374ffSToby Isaac     ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
161ad8374ffSToby Isaac   }
162ad8374ffSToby Isaac   label->validIS[v] = PETSC_FALSE;
163c58f1c22SToby Isaac   PetscFunctionReturn(0);
164c58f1c22SToby Isaac }
165c58f1c22SToby Isaac 
166b9907514SLisandro Dalcin #if !defined(DMLABEL_LOOKUP_THRESHOLD)
167b9907514SLisandro Dalcin #define DMLABEL_LOOKUP_THRESHOLD 16
168b9907514SLisandro Dalcin #endif
169b9907514SLisandro Dalcin 
1700c3c4a36SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index)
1710c3c4a36SLisandro Dalcin {
1720c3c4a36SLisandro Dalcin   PetscInt       v;
173b9907514SLisandro Dalcin   PetscErrorCode ierr;
1740e79e033SBarry Smith 
1750c3c4a36SLisandro Dalcin   PetscFunctionBegin;
1760e79e033SBarry Smith   *index = -1;
177b9907514SLisandro Dalcin   if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) {
178b9907514SLisandro Dalcin     for (v = 0; v < label->numStrata; ++v)
179b9907514SLisandro Dalcin       if (label->stratumValues[v] == value) {*index = v; break;}
180b9907514SLisandro Dalcin   } else {
181b9907514SLisandro Dalcin     ierr = PetscHMapIGet(label->hmap, value, index);CHKERRQ(ierr);
1820e79e033SBarry Smith   }
18390e9b2aeSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
18490e9b2aeSLisandro Dalcin   { /* Check strata hash map consistency */
18590e9b2aeSLisandro Dalcin     PetscInt len, loc = -1;
18690e9b2aeSLisandro Dalcin     ierr = PetscHMapIGetSize(label->hmap, &len);CHKERRQ(ierr);
18790e9b2aeSLisandro Dalcin     if (len != label->numStrata) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map size");
18890e9b2aeSLisandro Dalcin     if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) {
18990e9b2aeSLisandro Dalcin       ierr = PetscHMapIGet(label->hmap, value, &loc);CHKERRQ(ierr);
19090e9b2aeSLisandro Dalcin     } else {
19190e9b2aeSLisandro Dalcin       for (v = 0; v < label->numStrata; ++v)
19290e9b2aeSLisandro Dalcin         if (label->stratumValues[v] == value) {loc = v; break;}
19390e9b2aeSLisandro Dalcin     }
19490e9b2aeSLisandro Dalcin     if (loc != *index) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map lookup");
19590e9b2aeSLisandro Dalcin   }
19690e9b2aeSLisandro Dalcin #endif
1970c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
1980c3c4a36SLisandro Dalcin }
1990c3c4a36SLisandro Dalcin 
200b9907514SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index)
201c58f1c22SToby Isaac {
202b9907514SLisandro Dalcin   PetscInt       v;
203b9907514SLisandro Dalcin   PetscInt      *tmpV;
204b9907514SLisandro Dalcin   PetscInt      *tmpS;
205b9907514SLisandro Dalcin   PetscHSetI    *tmpH, ht;
206b9907514SLisandro Dalcin   IS            *tmpP, is;
207c58f1c22SToby Isaac   PetscBool     *tmpB;
208b9907514SLisandro Dalcin   PetscHMapI     hmap = label->hmap;
209c58f1c22SToby Isaac   PetscErrorCode ierr;
210c58f1c22SToby Isaac 
211c58f1c22SToby Isaac   PetscFunctionBegin;
212b9907514SLisandro Dalcin   v    = label->numStrata;
213b9907514SLisandro Dalcin   tmpV = label->stratumValues;
214b9907514SLisandro Dalcin   tmpS = label->stratumSizes;
215b9907514SLisandro Dalcin   tmpH = label->ht;
216b9907514SLisandro Dalcin   tmpP = label->points;
217b9907514SLisandro Dalcin   tmpB = label->validIS;
2188e1f8cf0SLisandro Dalcin   { /* TODO: PetscRealloc() is broken, use malloc+memcpy+free  */
2198e1f8cf0SLisandro Dalcin     PetscInt   *oldV = tmpV;
2208e1f8cf0SLisandro Dalcin     PetscInt   *oldS = tmpS;
2218e1f8cf0SLisandro Dalcin     PetscHSetI *oldH = tmpH;
2228e1f8cf0SLisandro Dalcin     IS         *oldP = tmpP;
2238e1f8cf0SLisandro Dalcin     PetscBool  *oldB = tmpB;
2248e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpV), &tmpV);CHKERRQ(ierr);
2258e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpS), &tmpS);CHKERRQ(ierr);
2268e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpH), &tmpH);CHKERRQ(ierr);
2278e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpP), &tmpP);CHKERRQ(ierr);
2288e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpB), &tmpB);CHKERRQ(ierr);
229580bdb30SBarry Smith     ierr = PetscArraycpy(tmpV, oldV, v);CHKERRQ(ierr);
230580bdb30SBarry Smith     ierr = PetscArraycpy(tmpS, oldS, v);CHKERRQ(ierr);
231580bdb30SBarry Smith     ierr = PetscArraycpy(tmpH, oldH, v);CHKERRQ(ierr);
232580bdb30SBarry Smith     ierr = PetscArraycpy(tmpP, oldP, v);CHKERRQ(ierr);
233580bdb30SBarry Smith     ierr = PetscArraycpy(tmpB, oldB, v);CHKERRQ(ierr);
2348e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldV);CHKERRQ(ierr);
2358e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldS);CHKERRQ(ierr);
2368e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldH);CHKERRQ(ierr);
2378e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldP);CHKERRQ(ierr);
2388e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldB);CHKERRQ(ierr);
2398e1f8cf0SLisandro Dalcin   }
240b9907514SLisandro Dalcin   label->numStrata     = v+1;
241c58f1c22SToby Isaac   label->stratumValues = tmpV;
242c58f1c22SToby Isaac   label->stratumSizes  = tmpS;
243c58f1c22SToby Isaac   label->ht            = tmpH;
244c58f1c22SToby Isaac   label->points        = tmpP;
245ad8374ffSToby Isaac   label->validIS       = tmpB;
246b9907514SLisandro Dalcin   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
247b9907514SLisandro Dalcin   ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);CHKERRQ(ierr);
248b9907514SLisandro Dalcin   ierr = PetscHMapISet(hmap, value, v);CHKERRQ(ierr);
249b9907514SLisandro Dalcin   tmpV[v] = value;
250b9907514SLisandro Dalcin   tmpS[v] = 0;
251b9907514SLisandro Dalcin   tmpH[v] = ht;
252b9907514SLisandro Dalcin   tmpP[v] = is;
253b9907514SLisandro Dalcin   tmpB[v] = PETSC_TRUE;
254277ea44aSLisandro Dalcin   ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
2550c3c4a36SLisandro Dalcin   *index = v;
2560c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
2570c3c4a36SLisandro Dalcin }
2580c3c4a36SLisandro Dalcin 
259b9907514SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupAddStratum(DMLabel label, PetscInt value, PetscInt *index)
260b9907514SLisandro Dalcin {
261b9907514SLisandro Dalcin   PetscErrorCode ierr;
262b9907514SLisandro Dalcin   PetscFunctionBegin;
263b9907514SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, index);CHKERRQ(ierr);
264b9907514SLisandro Dalcin   if (*index < 0) {ierr = DMLabelNewStratum(label, value, index);CHKERRQ(ierr);}
265b9907514SLisandro Dalcin   PetscFunctionReturn(0);
266b9907514SLisandro Dalcin }
267b9907514SLisandro Dalcin 
268b9907514SLisandro Dalcin /*@
269b9907514SLisandro Dalcin   DMLabelAddStratum - Adds a new stratum value in a DMLabel
270b9907514SLisandro Dalcin 
271b9907514SLisandro Dalcin   Input Parameter:
272b9907514SLisandro Dalcin + label - The DMLabel
273b9907514SLisandro Dalcin - value - The stratum value
274b9907514SLisandro Dalcin 
275b9907514SLisandro Dalcin   Level: beginner
276b9907514SLisandro Dalcin 
277b9907514SLisandro Dalcin .seealso:  DMLabelCreate(), DMLabelDestroy()
278b9907514SLisandro Dalcin @*/
2790c3c4a36SLisandro Dalcin PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
2800c3c4a36SLisandro Dalcin {
2810c3c4a36SLisandro Dalcin   PetscInt       v;
2820c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
2830c3c4a36SLisandro Dalcin 
2840c3c4a36SLisandro Dalcin   PetscFunctionBegin;
285d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
286b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
287b9907514SLisandro Dalcin   PetscFunctionReturn(0);
288b9907514SLisandro Dalcin }
289b9907514SLisandro Dalcin 
290b9907514SLisandro Dalcin /*@
291b9907514SLisandro Dalcin   DMLabelAddStrata - Adds new stratum values in a DMLabel
292b9907514SLisandro Dalcin 
2935b5e7992SMatthew G. Knepley   Not collective
2945b5e7992SMatthew G. Knepley 
295b9907514SLisandro Dalcin   Input Parameter:
296b9907514SLisandro Dalcin + label - The DMLabel
297b9907514SLisandro Dalcin . numStrata - The number of stratum values
298b9907514SLisandro Dalcin - stratumValues - The stratum values
299b9907514SLisandro Dalcin 
300b9907514SLisandro Dalcin   Level: beginner
301b9907514SLisandro Dalcin 
302b9907514SLisandro Dalcin .seealso:  DMLabelCreate(), DMLabelDestroy()
303b9907514SLisandro Dalcin @*/
304b9907514SLisandro Dalcin PetscErrorCode DMLabelAddStrata(DMLabel label, PetscInt numStrata, const PetscInt stratumValues[])
305b9907514SLisandro Dalcin {
306b9907514SLisandro Dalcin   PetscInt       *values, v;
307b9907514SLisandro Dalcin   PetscErrorCode ierr;
308b9907514SLisandro Dalcin 
309b9907514SLisandro Dalcin   PetscFunctionBegin;
310b9907514SLisandro Dalcin   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
311b9907514SLisandro Dalcin   if (numStrata) PetscValidIntPointer(stratumValues, 3);
312b9907514SLisandro Dalcin   ierr = PetscMalloc1(numStrata, &values);CHKERRQ(ierr);
313580bdb30SBarry Smith   ierr = PetscArraycpy(values, stratumValues, numStrata);CHKERRQ(ierr);
314b9907514SLisandro Dalcin   ierr = PetscSortRemoveDupsInt(&numStrata, values);CHKERRQ(ierr);
315b9907514SLisandro Dalcin   if (!label->numStrata) { /* Fast preallocation */
316b9907514SLisandro Dalcin     PetscInt   *tmpV;
317b9907514SLisandro Dalcin     PetscInt   *tmpS;
318b9907514SLisandro Dalcin     PetscHSetI *tmpH, ht;
319b9907514SLisandro Dalcin     IS         *tmpP, is;
320b9907514SLisandro Dalcin     PetscBool  *tmpB;
321b9907514SLisandro Dalcin     PetscHMapI  hmap = label->hmap;
322b9907514SLisandro Dalcin 
323b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpV);CHKERRQ(ierr);
324b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpS);CHKERRQ(ierr);
325b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpH);CHKERRQ(ierr);
326b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpP);CHKERRQ(ierr);
327b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpB);CHKERRQ(ierr);
328b9907514SLisandro Dalcin     label->numStrata     = numStrata;
329b9907514SLisandro Dalcin     label->stratumValues = tmpV;
330b9907514SLisandro Dalcin     label->stratumSizes  = tmpS;
331b9907514SLisandro Dalcin     label->ht            = tmpH;
332b9907514SLisandro Dalcin     label->points        = tmpP;
333b9907514SLisandro Dalcin     label->validIS       = tmpB;
334b9907514SLisandro Dalcin     for (v = 0; v < numStrata; ++v) {
335b9907514SLisandro Dalcin       ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
336b9907514SLisandro Dalcin       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);CHKERRQ(ierr);
337b9907514SLisandro Dalcin       ierr = PetscHMapISet(hmap, values[v], v);CHKERRQ(ierr);
338b9907514SLisandro Dalcin       tmpV[v] = values[v];
339b9907514SLisandro Dalcin       tmpS[v] = 0;
340b9907514SLisandro Dalcin       tmpH[v] = ht;
341b9907514SLisandro Dalcin       tmpP[v] = is;
342b9907514SLisandro Dalcin       tmpB[v] = PETSC_TRUE;
343b9907514SLisandro Dalcin     }
344277ea44aSLisandro Dalcin     ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
345b9907514SLisandro Dalcin   } else {
346b9907514SLisandro Dalcin     for (v = 0; v < numStrata; ++v) {
347b9907514SLisandro Dalcin       ierr = DMLabelAddStratum(label, values[v]);CHKERRQ(ierr);
348b9907514SLisandro Dalcin     }
349b9907514SLisandro Dalcin   }
350b9907514SLisandro Dalcin   ierr = PetscFree(values);CHKERRQ(ierr);
351b9907514SLisandro Dalcin   PetscFunctionReturn(0);
352b9907514SLisandro Dalcin }
353b9907514SLisandro Dalcin 
354b9907514SLisandro Dalcin /*@
355b9907514SLisandro Dalcin   DMLabelAddStrataIS - Adds new stratum values in a DMLabel
356b9907514SLisandro Dalcin 
3575b5e7992SMatthew G. Knepley   Not collective
3585b5e7992SMatthew G. Knepley 
359b9907514SLisandro Dalcin   Input Parameter:
360b9907514SLisandro Dalcin + label - The DMLabel
361b9907514SLisandro Dalcin - valueIS - Index set with stratum values
362b9907514SLisandro Dalcin 
363b9907514SLisandro Dalcin   Level: beginner
364b9907514SLisandro Dalcin 
365b9907514SLisandro Dalcin .seealso:  DMLabelCreate(), DMLabelDestroy()
366b9907514SLisandro Dalcin @*/
367b9907514SLisandro Dalcin PetscErrorCode DMLabelAddStrataIS(DMLabel label, IS valueIS)
368b9907514SLisandro Dalcin {
369b9907514SLisandro Dalcin   PetscInt       numStrata;
370b9907514SLisandro Dalcin   const PetscInt *stratumValues;
371b9907514SLisandro Dalcin   PetscErrorCode ierr;
372b9907514SLisandro Dalcin 
373b9907514SLisandro Dalcin   PetscFunctionBegin;
374b9907514SLisandro Dalcin   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
375b9907514SLisandro Dalcin   PetscValidHeaderSpecific(valueIS, IS_CLASSID, 2);
376b9907514SLisandro Dalcin   ierr = ISGetLocalSize(valueIS, &numStrata);CHKERRQ(ierr);
377b9907514SLisandro Dalcin   ierr = ISGetIndices(valueIS, &stratumValues);CHKERRQ(ierr);
378b9907514SLisandro Dalcin   ierr = DMLabelAddStrata(label, numStrata, stratumValues);CHKERRQ(ierr);
379c58f1c22SToby Isaac   PetscFunctionReturn(0);
380c58f1c22SToby Isaac }
381c58f1c22SToby Isaac 
382c58f1c22SToby Isaac static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
383c58f1c22SToby Isaac {
384c58f1c22SToby Isaac   PetscInt       v;
385c58f1c22SToby Isaac   PetscMPIInt    rank;
386c58f1c22SToby Isaac   PetscErrorCode ierr;
387c58f1c22SToby Isaac 
388c58f1c22SToby Isaac   PetscFunctionBegin;
389c58f1c22SToby Isaac   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr);
390c58f1c22SToby Isaac   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
391c58f1c22SToby Isaac   if (label) {
392d67d17b1SMatthew G. Knepley     const char *name;
393d67d17b1SMatthew G. Knepley 
394d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
395d67d17b1SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", name);CHKERRQ(ierr);
396c58f1c22SToby Isaac     if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);}
397c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
398c58f1c22SToby Isaac       const PetscInt value = label->stratumValues[v];
399ad8374ffSToby Isaac       const PetscInt *points;
400c58f1c22SToby Isaac       PetscInt       p;
401c58f1c22SToby Isaac 
402ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
403c58f1c22SToby Isaac       for (p = 0; p < label->stratumSizes[v]; ++p) {
404ad8374ffSToby Isaac         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr);
405c58f1c22SToby Isaac       }
406ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
407c58f1c22SToby Isaac     }
408c58f1c22SToby Isaac   }
409c58f1c22SToby Isaac   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
410c58f1c22SToby Isaac   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
411c58f1c22SToby Isaac   PetscFunctionReturn(0);
412c58f1c22SToby Isaac }
413c58f1c22SToby Isaac 
414c58f1c22SToby Isaac /*@C
415c58f1c22SToby Isaac   DMLabelView - View the label
416c58f1c22SToby Isaac 
4175b5e7992SMatthew G. Knepley   Collective on viewer
4185b5e7992SMatthew G. Knepley 
419c58f1c22SToby Isaac   Input Parameters:
420c58f1c22SToby Isaac + label - The DMLabel
421c58f1c22SToby Isaac - viewer - The PetscViewer
422c58f1c22SToby Isaac 
423c58f1c22SToby Isaac   Level: intermediate
424c58f1c22SToby Isaac 
425c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy()
426c58f1c22SToby Isaac @*/
427c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
428c58f1c22SToby Isaac {
429c58f1c22SToby Isaac   PetscBool      iascii;
430c58f1c22SToby Isaac   PetscErrorCode ierr;
431c58f1c22SToby Isaac 
432c58f1c22SToby Isaac   PetscFunctionBegin;
433d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
434b9907514SLisandro Dalcin   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)label), &viewer);CHKERRQ(ierr);}
435c58f1c22SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
436c58f1c22SToby Isaac   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
437c58f1c22SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
438c58f1c22SToby Isaac   if (iascii) {
439c58f1c22SToby Isaac     ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr);
440c58f1c22SToby Isaac   }
441c58f1c22SToby Isaac   PetscFunctionReturn(0);
442c58f1c22SToby Isaac }
443c58f1c22SToby Isaac 
44484f0b6dfSMatthew G. Knepley /*@
445d67d17b1SMatthew G. Knepley   DMLabelReset - Destroys internal data structures in a DMLabel
446d67d17b1SMatthew G. Knepley 
4475b5e7992SMatthew G. Knepley   Not collective
4485b5e7992SMatthew G. Knepley 
449d67d17b1SMatthew G. Knepley   Input Parameter:
450d67d17b1SMatthew G. Knepley . label - The DMLabel
451d67d17b1SMatthew G. Knepley 
452d67d17b1SMatthew G. Knepley   Level: beginner
453d67d17b1SMatthew G. Knepley 
454d67d17b1SMatthew G. Knepley .seealso: DMLabelDestroy(), DMLabelCreate()
455d67d17b1SMatthew G. Knepley @*/
456d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelReset(DMLabel label)
457d67d17b1SMatthew G. Knepley {
458d67d17b1SMatthew G. Knepley   PetscInt       v;
459d67d17b1SMatthew G. Knepley   PetscErrorCode ierr;
460d67d17b1SMatthew G. Knepley 
461d67d17b1SMatthew G. Knepley   PetscFunctionBegin;
462d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
463d67d17b1SMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
464d67d17b1SMatthew G. Knepley     ierr = PetscHSetIDestroy(&label->ht[v]);CHKERRQ(ierr);
465d67d17b1SMatthew G. Knepley     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
466d67d17b1SMatthew G. Knepley   }
467b9907514SLisandro Dalcin   label->numStrata = 0;
468b9907514SLisandro Dalcin   ierr = PetscFree(label->stratumValues);CHKERRQ(ierr);
469b9907514SLisandro Dalcin   ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr);
470d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->ht);CHKERRQ(ierr);
471d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->points);CHKERRQ(ierr);
472d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->validIS);CHKERRQ(ierr);
473b9907514SLisandro Dalcin   ierr = PetscHMapIReset(label->hmap);CHKERRQ(ierr);
474b9907514SLisandro Dalcin   label->pStart = -1;
475b9907514SLisandro Dalcin   label->pEnd   = -1;
476d67d17b1SMatthew G. Knepley   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
477d67d17b1SMatthew G. Knepley   PetscFunctionReturn(0);
478d67d17b1SMatthew G. Knepley }
479d67d17b1SMatthew G. Knepley 
480d67d17b1SMatthew G. Knepley /*@
48184f0b6dfSMatthew G. Knepley   DMLabelDestroy - Destroys a DMLabel
48284f0b6dfSMatthew G. Knepley 
4835b5e7992SMatthew G. Knepley   Collective on label
4845b5e7992SMatthew G. Knepley 
48584f0b6dfSMatthew G. Knepley   Input Parameter:
48684f0b6dfSMatthew G. Knepley . label - The DMLabel
48784f0b6dfSMatthew G. Knepley 
48884f0b6dfSMatthew G. Knepley   Level: beginner
48984f0b6dfSMatthew G. Knepley 
490d67d17b1SMatthew G. Knepley .seealso: DMLabelReset(), DMLabelCreate()
49184f0b6dfSMatthew G. Knepley @*/
492c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label)
493c58f1c22SToby Isaac {
494c58f1c22SToby Isaac   PetscErrorCode ierr;
495c58f1c22SToby Isaac 
496c58f1c22SToby Isaac   PetscFunctionBegin;
497d67d17b1SMatthew G. Knepley   if (!*label) PetscFunctionReturn(0);
498d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific((*label),DMLABEL_CLASSID,1);
499b9907514SLisandro Dalcin   if (--((PetscObject)(*label))->refct > 0) {*label = NULL; PetscFunctionReturn(0);}
500d67d17b1SMatthew G. Knepley   ierr = DMLabelReset(*label);CHKERRQ(ierr);
501b9907514SLisandro Dalcin   ierr = PetscHMapIDestroy(&(*label)->hmap);CHKERRQ(ierr);
502d67d17b1SMatthew G. Knepley   ierr = PetscHeaderDestroy(label);CHKERRQ(ierr);
503c58f1c22SToby Isaac   PetscFunctionReturn(0);
504c58f1c22SToby Isaac }
505c58f1c22SToby Isaac 
50684f0b6dfSMatthew G. Knepley /*@
50784f0b6dfSMatthew G. Knepley   DMLabelDuplicate - Duplicates a DMLabel
50884f0b6dfSMatthew G. Knepley 
5095b5e7992SMatthew G. Knepley   Collective on label
5105b5e7992SMatthew G. Knepley 
51184f0b6dfSMatthew G. Knepley   Input Parameter:
51284f0b6dfSMatthew G. Knepley . label - The DMLabel
51384f0b6dfSMatthew G. Knepley 
51484f0b6dfSMatthew G. Knepley   Output Parameter:
51584f0b6dfSMatthew G. Knepley . labelnew - location to put new vector
51684f0b6dfSMatthew G. Knepley 
51784f0b6dfSMatthew G. Knepley   Level: intermediate
51884f0b6dfSMatthew G. Knepley 
51984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelDestroy()
52084f0b6dfSMatthew G. Knepley @*/
521c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
522c58f1c22SToby Isaac {
523d67d17b1SMatthew G. Knepley   const char    *name;
524ad8374ffSToby Isaac   PetscInt       v;
525c58f1c22SToby Isaac   PetscErrorCode ierr;
526c58f1c22SToby Isaac 
527c58f1c22SToby Isaac   PetscFunctionBegin;
528d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
529c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
530d67d17b1SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
531d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PetscObjectComm((PetscObject) label), name, labelnew);CHKERRQ(ierr);
532c58f1c22SToby Isaac 
533c58f1c22SToby Isaac   (*labelnew)->numStrata    = label->numStrata;
5345aa44df4SToby Isaac   (*labelnew)->defaultValue = label->defaultValue;
535c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr);
536c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr);
537c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr);
538c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr);
539ad8374ffSToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr);
540c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
541e8f14785SLisandro Dalcin     ierr = PetscHSetICreate(&(*labelnew)->ht[v]);CHKERRQ(ierr);
542c58f1c22SToby Isaac     (*labelnew)->stratumValues[v]  = label->stratumValues[v];
543c58f1c22SToby Isaac     (*labelnew)->stratumSizes[v]   = label->stratumSizes[v];
544ad8374ffSToby Isaac     ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr);
545ad8374ffSToby Isaac     (*labelnew)->points[v]         = label->points[v];
546b9907514SLisandro Dalcin     (*labelnew)->validIS[v]        = PETSC_TRUE;
547c58f1c22SToby Isaac   }
548f14fe40dSLisandro Dalcin   ierr = PetscHMapIDestroy(&(*labelnew)->hmap);CHKERRQ(ierr);
549b9907514SLisandro Dalcin   ierr = PetscHMapIDuplicate(label->hmap,&(*labelnew)->hmap);CHKERRQ(ierr);
550c58f1c22SToby Isaac   (*labelnew)->pStart = -1;
551c58f1c22SToby Isaac   (*labelnew)->pEnd   = -1;
552c58f1c22SToby Isaac   (*labelnew)->bt     = NULL;
553c58f1c22SToby Isaac   PetscFunctionReturn(0);
554c58f1c22SToby Isaac }
555c58f1c22SToby Isaac 
556c6a43d28SMatthew G. Knepley /*@
557c6a43d28SMatthew G. Knepley   DMLabelComputeIndex - Create an index structure for membership determination, automatically determining the bounds
558c6a43d28SMatthew G. Knepley 
5595b5e7992SMatthew G. Knepley   Not collective
5605b5e7992SMatthew G. Knepley 
561c6a43d28SMatthew G. Knepley   Input Parameter:
562c6a43d28SMatthew G. Knepley . label  - The DMLabel
563c6a43d28SMatthew G. Knepley 
564c6a43d28SMatthew G. Knepley   Level: intermediate
565c6a43d28SMatthew G. Knepley 
566c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
567c6a43d28SMatthew G. Knepley @*/
568c6a43d28SMatthew G. Knepley PetscErrorCode DMLabelComputeIndex(DMLabel label)
569c6a43d28SMatthew G. Knepley {
570c6a43d28SMatthew G. Knepley   PetscInt       pStart = PETSC_MAX_INT, pEnd = -1, v;
571c6a43d28SMatthew G. Knepley   PetscErrorCode ierr;
572c6a43d28SMatthew G. Knepley 
573c6a43d28SMatthew G. Knepley   PetscFunctionBegin;
574c6a43d28SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
575c6a43d28SMatthew G. Knepley   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
576c6a43d28SMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
577c6a43d28SMatthew G. Knepley     const PetscInt *points;
578c6a43d28SMatthew G. Knepley     PetscInt       i;
579c6a43d28SMatthew G. Knepley 
580c6a43d28SMatthew G. Knepley     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
581c6a43d28SMatthew G. Knepley     for (i = 0; i < label->stratumSizes[v]; ++i) {
582c6a43d28SMatthew G. Knepley       const PetscInt point = points[i];
583c6a43d28SMatthew G. Knepley 
584c6a43d28SMatthew G. Knepley       pStart = PetscMin(point,   pStart);
585c6a43d28SMatthew G. Knepley       pEnd   = PetscMax(point+1, pEnd);
586c6a43d28SMatthew G. Knepley     }
587c6a43d28SMatthew G. Knepley     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
588c6a43d28SMatthew G. Knepley   }
589c6a43d28SMatthew G. Knepley   label->pStart = pStart == PETSC_MAX_INT ? -1 : pStart;
590c6a43d28SMatthew G. Knepley   label->pEnd   = pEnd;
591c6a43d28SMatthew G. Knepley   ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
592c6a43d28SMatthew G. Knepley   PetscFunctionReturn(0);
593c6a43d28SMatthew G. Knepley }
594c6a43d28SMatthew G. Knepley 
595c6a43d28SMatthew G. Knepley /*@
596c6a43d28SMatthew G. Knepley   DMLabelCreateIndex - Create an index structure for membership determination
597c6a43d28SMatthew G. Knepley 
5985b5e7992SMatthew G. Knepley   Not collective
5995b5e7992SMatthew G. Knepley 
600c6a43d28SMatthew G. Knepley   Input Parameters:
601c6a43d28SMatthew G. Knepley + label  - The DMLabel
602c6a43d28SMatthew G. Knepley . pStart - The smallest point
603c6a43d28SMatthew G. Knepley - pEnd   - The largest point + 1
604c6a43d28SMatthew G. Knepley 
605c6a43d28SMatthew G. Knepley   Level: intermediate
606c6a43d28SMatthew G. Knepley 
607c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelComputeIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
608c6a43d28SMatthew G. Knepley @*/
609c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
610c58f1c22SToby Isaac {
611c58f1c22SToby Isaac   PetscInt       v;
612c58f1c22SToby Isaac   PetscErrorCode ierr;
613c58f1c22SToby Isaac 
614c58f1c22SToby Isaac   PetscFunctionBegin;
615d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
6160c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
617c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
618c58f1c22SToby Isaac   label->pStart = pStart;
619c58f1c22SToby Isaac   label->pEnd   = pEnd;
620c6a43d28SMatthew G. Knepley   /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
621c58f1c22SToby Isaac   ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr);
622c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
623ad8374ffSToby Isaac     const PetscInt *points;
624c58f1c22SToby Isaac     PetscInt       i;
625c58f1c22SToby Isaac 
626ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
627c58f1c22SToby Isaac     for (i = 0; i < label->stratumSizes[v]; ++i) {
628ad8374ffSToby Isaac       const PetscInt point = points[i];
629c58f1c22SToby Isaac 
630c58f1c22SToby 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);
631c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr);
632c58f1c22SToby Isaac     }
633ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
634c58f1c22SToby Isaac   }
635c58f1c22SToby Isaac   PetscFunctionReturn(0);
636c58f1c22SToby Isaac }
637c58f1c22SToby Isaac 
638c6a43d28SMatthew G. Knepley /*@
639c6a43d28SMatthew G. Knepley   DMLabelDestroyIndex - Destroy the index structure
640c6a43d28SMatthew G. Knepley 
6415b5e7992SMatthew G. Knepley   Not collective
6425b5e7992SMatthew G. Knepley 
643c6a43d28SMatthew G. Knepley   Input Parameter:
644c6a43d28SMatthew G. Knepley . label - the DMLabel
645c6a43d28SMatthew G. Knepley 
646c6a43d28SMatthew G. Knepley   Level: intermediate
647c6a43d28SMatthew G. Knepley 
648c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
649c6a43d28SMatthew G. Knepley @*/
650c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label)
651c58f1c22SToby Isaac {
652c58f1c22SToby Isaac   PetscErrorCode ierr;
653c58f1c22SToby Isaac 
654c58f1c22SToby Isaac   PetscFunctionBegin;
655d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
656c58f1c22SToby Isaac   label->pStart = -1;
657c58f1c22SToby Isaac   label->pEnd   = -1;
6580c3c4a36SLisandro Dalcin   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
659c58f1c22SToby Isaac   PetscFunctionReturn(0);
660c58f1c22SToby Isaac }
661c58f1c22SToby Isaac 
662c58f1c22SToby Isaac /*@
663c6a43d28SMatthew G. Knepley   DMLabelGetBounds - Return the smallest and largest point in the label
664c6a43d28SMatthew G. Knepley 
6655b5e7992SMatthew G. Knepley   Not collective
6665b5e7992SMatthew G. Knepley 
667c6a43d28SMatthew G. Knepley   Input Parameter:
668c6a43d28SMatthew G. Knepley . label - the DMLabel
669c6a43d28SMatthew G. Knepley 
670c6a43d28SMatthew G. Knepley   Output Parameters:
671c6a43d28SMatthew G. Knepley + pStart - The smallest point
672c6a43d28SMatthew G. Knepley - pEnd   - The largest point + 1
673c6a43d28SMatthew G. Knepley 
674c6a43d28SMatthew G. Knepley   Note: This will compute an index for the label if one does not exist.
675c6a43d28SMatthew G. Knepley 
676c6a43d28SMatthew G. Knepley   Level: intermediate
677c6a43d28SMatthew G. Knepley 
678c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
679c6a43d28SMatthew G. Knepley @*/
680c6a43d28SMatthew G. Knepley PetscErrorCode DMLabelGetBounds(DMLabel label, PetscInt *pStart, PetscInt *pEnd)
681c6a43d28SMatthew G. Knepley {
682c6a43d28SMatthew G. Knepley   PetscErrorCode ierr;
683c6a43d28SMatthew G. Knepley 
684c6a43d28SMatthew G. Knepley   PetscFunctionBegin;
685c6a43d28SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
686c6a43d28SMatthew G. Knepley   if ((label->pStart == -1) && (label->pEnd == -1)) {ierr = DMLabelComputeIndex(label);CHKERRQ(ierr);}
687c6a43d28SMatthew G. Knepley   if (pStart) {
688534a8f05SLisandro Dalcin     PetscValidIntPointer(pStart, 2);
689c6a43d28SMatthew G. Knepley     *pStart = label->pStart;
690c6a43d28SMatthew G. Knepley   }
691c6a43d28SMatthew G. Knepley   if (pEnd) {
692534a8f05SLisandro Dalcin     PetscValidIntPointer(pEnd, 3);
693c6a43d28SMatthew G. Knepley     *pEnd = label->pEnd;
694c6a43d28SMatthew G. Knepley   }
695c6a43d28SMatthew G. Knepley   PetscFunctionReturn(0);
696c6a43d28SMatthew G. Knepley }
697c6a43d28SMatthew G. Knepley 
698c6a43d28SMatthew G. Knepley /*@
699c58f1c22SToby Isaac   DMLabelHasValue - Determine whether a label assigns the value to any point
700c58f1c22SToby Isaac 
7015b5e7992SMatthew G. Knepley   Not collective
7025b5e7992SMatthew G. Knepley 
703c58f1c22SToby Isaac   Input Parameters:
704c58f1c22SToby Isaac + label - the DMLabel
705c58f1c22SToby Isaac - value - the value
706c58f1c22SToby Isaac 
707c58f1c22SToby Isaac   Output Parameter:
708c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point
709c58f1c22SToby Isaac 
710c58f1c22SToby Isaac   Level: developer
711c58f1c22SToby Isaac 
712c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
713c58f1c22SToby Isaac @*/
714c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
715c58f1c22SToby Isaac {
716c58f1c22SToby Isaac   PetscInt v;
7170c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
718c58f1c22SToby Isaac 
719c58f1c22SToby Isaac   PetscFunctionBegin;
720d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
721534a8f05SLisandro Dalcin   PetscValidBoolPointer(contains, 3);
7220c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
7230c3c4a36SLisandro Dalcin   *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE;
724c58f1c22SToby Isaac   PetscFunctionReturn(0);
725c58f1c22SToby Isaac }
726c58f1c22SToby Isaac 
727c58f1c22SToby Isaac /*@
728c58f1c22SToby Isaac   DMLabelHasPoint - Determine whether a label assigns a value to a point
729c58f1c22SToby Isaac 
7305b5e7992SMatthew G. Knepley   Not collective
7315b5e7992SMatthew G. Knepley 
732c58f1c22SToby Isaac   Input Parameters:
733c58f1c22SToby Isaac + label - the DMLabel
734c58f1c22SToby Isaac - point - the point
735c58f1c22SToby Isaac 
736c58f1c22SToby Isaac   Output Parameter:
737c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value
738c58f1c22SToby Isaac 
739c58f1c22SToby Isaac   Note: The user must call DMLabelCreateIndex() before this function.
740c58f1c22SToby Isaac 
741c58f1c22SToby Isaac   Level: developer
742c58f1c22SToby Isaac 
743c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
744c58f1c22SToby Isaac @*/
745c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
746c58f1c22SToby Isaac {
747c58f1c22SToby Isaac   PetscErrorCode ierr;
748c58f1c22SToby Isaac 
749c58f1c22SToby Isaac   PetscFunctionBeginHot;
750d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
751534a8f05SLisandro Dalcin   PetscValidBoolPointer(contains, 3);
752c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
753c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG)
754c58f1c22SToby Isaac   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
755c58f1c22SToby 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);
756c58f1c22SToby Isaac #endif
757c58f1c22SToby Isaac   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
758c58f1c22SToby Isaac   PetscFunctionReturn(0);
759c58f1c22SToby Isaac }
760c58f1c22SToby Isaac 
761c58f1c22SToby Isaac /*@
762c58f1c22SToby Isaac   DMLabelStratumHasPoint - Return true if the stratum contains a point
763c58f1c22SToby Isaac 
7645b5e7992SMatthew G. Knepley   Not collective
7655b5e7992SMatthew G. Knepley 
766c58f1c22SToby Isaac   Input Parameters:
767c58f1c22SToby Isaac + label - the DMLabel
768c58f1c22SToby Isaac . value - the stratum value
769c58f1c22SToby Isaac - point - the point
770c58f1c22SToby Isaac 
771c58f1c22SToby Isaac   Output Parameter:
772c58f1c22SToby Isaac . contains - true if the stratum contains the point
773c58f1c22SToby Isaac 
774c58f1c22SToby Isaac   Level: intermediate
775c58f1c22SToby Isaac 
776c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
777c58f1c22SToby Isaac @*/
778c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
779c58f1c22SToby Isaac {
780c58f1c22SToby Isaac   PetscInt       v;
781c58f1c22SToby Isaac   PetscErrorCode ierr;
782c58f1c22SToby Isaac 
7830c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
784d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
785534a8f05SLisandro Dalcin   PetscValidBoolPointer(contains, 4);
786c58f1c22SToby Isaac   *contains = PETSC_FALSE;
7870c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
7880c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
7890c3c4a36SLisandro Dalcin 
790ad8374ffSToby Isaac   if (label->validIS[v]) {
791c58f1c22SToby Isaac     PetscInt i;
792c58f1c22SToby Isaac 
793a2d74346SToby Isaac     ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
7940c3c4a36SLisandro Dalcin     if (i >= 0) *contains = PETSC_TRUE;
795c58f1c22SToby Isaac   } else {
796c58f1c22SToby Isaac     PetscBool has;
797c58f1c22SToby Isaac 
798b9907514SLisandro Dalcin     ierr = PetscHSetIHas(label->ht[v], point, &has);CHKERRQ(ierr);
7990c3c4a36SLisandro Dalcin     if (has) *contains = PETSC_TRUE;
800c58f1c22SToby Isaac   }
801c58f1c22SToby Isaac   PetscFunctionReturn(0);
802c58f1c22SToby Isaac }
803c58f1c22SToby Isaac 
80484f0b6dfSMatthew G. Knepley /*@
8055aa44df4SToby Isaac   DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
8065aa44df4SToby Isaac   When a label is created, it is initialized to -1.
8075aa44df4SToby Isaac 
8085b5e7992SMatthew G. Knepley   Not collective
8095b5e7992SMatthew G. Knepley 
8105aa44df4SToby Isaac   Input parameter:
8115aa44df4SToby Isaac . label - a DMLabel object
8125aa44df4SToby Isaac 
8135aa44df4SToby Isaac   Output parameter:
8145aa44df4SToby Isaac . defaultValue - the default value
8155aa44df4SToby Isaac 
8165aa44df4SToby Isaac   Level: beginner
8175aa44df4SToby Isaac 
8185aa44df4SToby Isaac .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
81984f0b6dfSMatthew G. Knepley @*/
8205aa44df4SToby Isaac PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
8215aa44df4SToby Isaac {
8225aa44df4SToby Isaac   PetscFunctionBegin;
823d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
8245aa44df4SToby Isaac   *defaultValue = label->defaultValue;
8255aa44df4SToby Isaac   PetscFunctionReturn(0);
8265aa44df4SToby Isaac }
8275aa44df4SToby Isaac 
82884f0b6dfSMatthew G. Knepley /*@
8295aa44df4SToby Isaac   DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
8305aa44df4SToby Isaac   When a label is created, it is initialized to -1.
8315aa44df4SToby Isaac 
8325b5e7992SMatthew G. Knepley   Not collective
8335b5e7992SMatthew G. Knepley 
8345aa44df4SToby Isaac   Input parameter:
8355aa44df4SToby Isaac . label - a DMLabel object
8365aa44df4SToby Isaac 
8375aa44df4SToby Isaac   Output parameter:
8385aa44df4SToby Isaac . defaultValue - the default value
8395aa44df4SToby Isaac 
8405aa44df4SToby Isaac   Level: beginner
8415aa44df4SToby Isaac 
8425aa44df4SToby Isaac .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
84384f0b6dfSMatthew G. Knepley @*/
8445aa44df4SToby Isaac PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
8455aa44df4SToby Isaac {
8465aa44df4SToby Isaac   PetscFunctionBegin;
847d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
8485aa44df4SToby Isaac   label->defaultValue = defaultValue;
8495aa44df4SToby Isaac   PetscFunctionReturn(0);
8505aa44df4SToby Isaac }
8515aa44df4SToby Isaac 
852c58f1c22SToby Isaac /*@
8535aa44df4SToby 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())
854c58f1c22SToby Isaac 
8555b5e7992SMatthew G. Knepley   Not collective
8565b5e7992SMatthew G. Knepley 
857c58f1c22SToby Isaac   Input Parameters:
858c58f1c22SToby Isaac + label - the DMLabel
859c58f1c22SToby Isaac - point - the point
860c58f1c22SToby Isaac 
861c58f1c22SToby Isaac   Output Parameter:
8628e68afcfSMatthew G. Knepley . value - The point value, or the default value (-1 by default)
863c58f1c22SToby Isaac 
864c58f1c22SToby Isaac   Level: intermediate
865c58f1c22SToby Isaac 
8665aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
867c58f1c22SToby Isaac @*/
868c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
869c58f1c22SToby Isaac {
870c58f1c22SToby Isaac   PetscInt       v;
871c58f1c22SToby Isaac   PetscErrorCode ierr;
872c58f1c22SToby Isaac 
8730c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
874d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
875c58f1c22SToby Isaac   PetscValidPointer(value, 3);
8765aa44df4SToby Isaac   *value = label->defaultValue;
877c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
878ad8374ffSToby Isaac     if (label->validIS[v]) {
879c58f1c22SToby Isaac       PetscInt i;
880c58f1c22SToby Isaac 
881a2d74346SToby Isaac       ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
882c58f1c22SToby Isaac       if (i >= 0) {
883c58f1c22SToby Isaac         *value = label->stratumValues[v];
884c58f1c22SToby Isaac         break;
885c58f1c22SToby Isaac       }
886c58f1c22SToby Isaac     } else {
887c58f1c22SToby Isaac       PetscBool has;
888c58f1c22SToby Isaac 
889b9907514SLisandro Dalcin       ierr = PetscHSetIHas(label->ht[v], point, &has);CHKERRQ(ierr);
890c58f1c22SToby Isaac       if (has) {
891c58f1c22SToby Isaac         *value = label->stratumValues[v];
892c58f1c22SToby Isaac         break;
893c58f1c22SToby Isaac       }
894c58f1c22SToby Isaac     }
895c58f1c22SToby Isaac   }
896c58f1c22SToby Isaac   PetscFunctionReturn(0);
897c58f1c22SToby Isaac }
898c58f1c22SToby Isaac 
899c58f1c22SToby Isaac /*@
900367003a6SStefano 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.
901c58f1c22SToby Isaac 
9025b5e7992SMatthew G. Knepley   Not collective
9035b5e7992SMatthew G. Knepley 
904c58f1c22SToby Isaac   Input Parameters:
905c58f1c22SToby Isaac + label - the DMLabel
906c58f1c22SToby Isaac . point - the point
907c58f1c22SToby Isaac - value - The point value
908c58f1c22SToby Isaac 
909c58f1c22SToby Isaac   Level: intermediate
910c58f1c22SToby Isaac 
9115aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
912c58f1c22SToby Isaac @*/
913c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
914c58f1c22SToby Isaac {
915c58f1c22SToby Isaac   PetscInt       v;
916c58f1c22SToby Isaac   PetscErrorCode ierr;
917c58f1c22SToby Isaac 
918c58f1c22SToby Isaac   PetscFunctionBegin;
919d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
9200c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
9215aa44df4SToby Isaac   if (value == label->defaultValue) PetscFunctionReturn(0);
922b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
923c58f1c22SToby Isaac   /* Set key */
9240c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
925e8f14785SLisandro Dalcin   ierr = PetscHSetIAdd(label->ht[v], point);CHKERRQ(ierr);
926c58f1c22SToby Isaac   PetscFunctionReturn(0);
927c58f1c22SToby Isaac }
928c58f1c22SToby Isaac 
929c58f1c22SToby Isaac /*@
930c58f1c22SToby Isaac   DMLabelClearValue - Clear the value a label assigns to a point
931c58f1c22SToby Isaac 
9325b5e7992SMatthew G. Knepley   Not collective
9335b5e7992SMatthew G. Knepley 
934c58f1c22SToby Isaac   Input Parameters:
935c58f1c22SToby Isaac + label - the DMLabel
936c58f1c22SToby Isaac . point - the point
937c58f1c22SToby Isaac - value - The point value
938c58f1c22SToby Isaac 
939c58f1c22SToby Isaac   Level: intermediate
940c58f1c22SToby Isaac 
941c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
942c58f1c22SToby Isaac @*/
943c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
944c58f1c22SToby Isaac {
945ad8374ffSToby Isaac   PetscInt       v;
946c58f1c22SToby Isaac   PetscErrorCode ierr;
947c58f1c22SToby Isaac 
948c58f1c22SToby Isaac   PetscFunctionBegin;
949d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
950c58f1c22SToby Isaac   /* Find label value */
9510c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
9520c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
9530c3c4a36SLisandro Dalcin 
954eeed21e7SToby Isaac   if (label->bt) {
955eeed21e7SToby 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);
956eeed21e7SToby Isaac     ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
957eeed21e7SToby Isaac   }
9580c3c4a36SLisandro Dalcin 
9590c3c4a36SLisandro Dalcin   /* Delete key */
9600c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
961e8f14785SLisandro Dalcin   ierr = PetscHSetIDel(label->ht[v], point);CHKERRQ(ierr);
962c58f1c22SToby Isaac   PetscFunctionReturn(0);
963c58f1c22SToby Isaac }
964c58f1c22SToby Isaac 
965c58f1c22SToby Isaac /*@
966c58f1c22SToby Isaac   DMLabelInsertIS - Set all points in the IS to a value
967c58f1c22SToby Isaac 
9685b5e7992SMatthew G. Knepley   Not collective
9695b5e7992SMatthew G. Knepley 
970c58f1c22SToby Isaac   Input Parameters:
971c58f1c22SToby Isaac + label - the DMLabel
972c58f1c22SToby Isaac . is    - the point IS
973c58f1c22SToby Isaac - value - The point value
974c58f1c22SToby Isaac 
975c58f1c22SToby Isaac   Level: intermediate
976c58f1c22SToby Isaac 
977c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
978c58f1c22SToby Isaac @*/
979c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
980c58f1c22SToby Isaac {
9810c3c4a36SLisandro Dalcin   PetscInt        v, n, p;
982c58f1c22SToby Isaac   const PetscInt *points;
983c58f1c22SToby Isaac   PetscErrorCode  ierr;
984c58f1c22SToby Isaac 
985c58f1c22SToby Isaac   PetscFunctionBegin;
986d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
987c58f1c22SToby Isaac   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
9880c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
9890c3c4a36SLisandro Dalcin   if (value == label->defaultValue) PetscFunctionReturn(0);
990b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
9910c3c4a36SLisandro Dalcin   /* Set keys */
9920c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
993c58f1c22SToby Isaac   ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr);
994c58f1c22SToby Isaac   ierr = ISGetIndices(is, &points);CHKERRQ(ierr);
995e8f14785SLisandro Dalcin   for (p = 0; p < n; ++p) {ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);}
996c58f1c22SToby Isaac   ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr);
997c58f1c22SToby Isaac   PetscFunctionReturn(0);
998c58f1c22SToby Isaac }
999c58f1c22SToby Isaac 
100084f0b6dfSMatthew G. Knepley /*@
100184f0b6dfSMatthew G. Knepley   DMLabelGetNumValues - Get the number of values that the DMLabel takes
100284f0b6dfSMatthew G. Knepley 
10035b5e7992SMatthew G. Knepley   Not collective
10045b5e7992SMatthew G. Knepley 
100584f0b6dfSMatthew G. Knepley   Input Parameter:
100684f0b6dfSMatthew G. Knepley . label - the DMLabel
100784f0b6dfSMatthew G. Knepley 
100884f0b6dfSMatthew G. Knepley   Output Paramater:
100984f0b6dfSMatthew G. Knepley . numValues - the number of values
101084f0b6dfSMatthew G. Knepley 
101184f0b6dfSMatthew G. Knepley   Level: intermediate
101284f0b6dfSMatthew G. Knepley 
101384f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
101484f0b6dfSMatthew G. Knepley @*/
1015c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
1016c58f1c22SToby Isaac {
1017c58f1c22SToby Isaac   PetscFunctionBegin;
1018d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1019c58f1c22SToby Isaac   PetscValidPointer(numValues, 2);
1020c58f1c22SToby Isaac   *numValues = label->numStrata;
1021c58f1c22SToby Isaac   PetscFunctionReturn(0);
1022c58f1c22SToby Isaac }
1023c58f1c22SToby Isaac 
102484f0b6dfSMatthew G. Knepley /*@
102584f0b6dfSMatthew G. Knepley   DMLabelGetValueIS - Get an IS of all values that the DMlabel takes
102684f0b6dfSMatthew G. Knepley 
10275b5e7992SMatthew G. Knepley   Not collective
10285b5e7992SMatthew G. Knepley 
102984f0b6dfSMatthew G. Knepley   Input Parameter:
103084f0b6dfSMatthew G. Knepley . label - the DMLabel
103184f0b6dfSMatthew G. Knepley 
103284f0b6dfSMatthew G. Knepley   Output Paramater:
103384f0b6dfSMatthew G. Knepley . is    - the value IS
103484f0b6dfSMatthew G. Knepley 
103584f0b6dfSMatthew G. Knepley   Level: intermediate
103684f0b6dfSMatthew G. Knepley 
103784f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
103884f0b6dfSMatthew G. Knepley @*/
1039c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
1040c58f1c22SToby Isaac {
1041c58f1c22SToby Isaac   PetscErrorCode ierr;
1042c58f1c22SToby Isaac 
1043c58f1c22SToby Isaac   PetscFunctionBegin;
1044d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1045c58f1c22SToby Isaac   PetscValidPointer(values, 2);
1046c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr);
1047c58f1c22SToby Isaac   PetscFunctionReturn(0);
1048c58f1c22SToby Isaac }
1049c58f1c22SToby Isaac 
105084f0b6dfSMatthew G. Knepley /*@
105184f0b6dfSMatthew G. Knepley   DMLabelHasStratum - Determine whether points exist with the given value
105284f0b6dfSMatthew G. Knepley 
10535b5e7992SMatthew G. Knepley   Not collective
10545b5e7992SMatthew G. Knepley 
105584f0b6dfSMatthew G. Knepley   Input Parameters:
105684f0b6dfSMatthew G. Knepley + label - the DMLabel
105784f0b6dfSMatthew G. Knepley - value - the stratum value
105884f0b6dfSMatthew G. Knepley 
105984f0b6dfSMatthew G. Knepley   Output Paramater:
106084f0b6dfSMatthew G. Knepley . exists - Flag saying whether points exist
106184f0b6dfSMatthew G. Knepley 
106284f0b6dfSMatthew G. Knepley   Level: intermediate
106384f0b6dfSMatthew G. Knepley 
106484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
106584f0b6dfSMatthew G. Knepley @*/
1066fada774cSMatthew G. Knepley PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
1067fada774cSMatthew G. Knepley {
1068fada774cSMatthew G. Knepley   PetscInt       v;
10690c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
1070fada774cSMatthew G. Knepley 
1071fada774cSMatthew G. Knepley   PetscFunctionBegin;
1072d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1073fada774cSMatthew G. Knepley   PetscValidPointer(exists, 3);
10740c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
10750c3c4a36SLisandro Dalcin   *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE;
1076fada774cSMatthew G. Knepley   PetscFunctionReturn(0);
1077fada774cSMatthew G. Knepley }
1078fada774cSMatthew G. Knepley 
107984f0b6dfSMatthew G. Knepley /*@
108084f0b6dfSMatthew G. Knepley   DMLabelGetStratumSize - Get the size of a stratum
108184f0b6dfSMatthew G. Knepley 
10825b5e7992SMatthew G. Knepley   Not collective
10835b5e7992SMatthew G. Knepley 
108484f0b6dfSMatthew G. Knepley   Input Parameters:
108584f0b6dfSMatthew G. Knepley + label - the DMLabel
108684f0b6dfSMatthew G. Knepley - value - the stratum value
108784f0b6dfSMatthew G. Knepley 
108884f0b6dfSMatthew G. Knepley   Output Paramater:
108984f0b6dfSMatthew G. Knepley . size - The number of points in the stratum
109084f0b6dfSMatthew G. Knepley 
109184f0b6dfSMatthew G. Knepley   Level: intermediate
109284f0b6dfSMatthew G. Knepley 
109384f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
109484f0b6dfSMatthew G. Knepley @*/
1095c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
1096c58f1c22SToby Isaac {
1097c58f1c22SToby Isaac   PetscInt       v;
1098c58f1c22SToby Isaac   PetscErrorCode ierr;
1099c58f1c22SToby Isaac 
1100c58f1c22SToby Isaac   PetscFunctionBegin;
1101d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1102c58f1c22SToby Isaac   PetscValidPointer(size, 3);
1103c58f1c22SToby Isaac   *size = 0;
11040c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
11050c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
1106c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
1107c58f1c22SToby Isaac   *size = label->stratumSizes[v];
1108c58f1c22SToby Isaac   PetscFunctionReturn(0);
1109c58f1c22SToby Isaac }
1110c58f1c22SToby Isaac 
111184f0b6dfSMatthew G. Knepley /*@
111284f0b6dfSMatthew G. Knepley   DMLabelGetStratumBounds - Get the largest and smallest point of a stratum
111384f0b6dfSMatthew G. Knepley 
11145b5e7992SMatthew G. Knepley   Not collective
11155b5e7992SMatthew G. Knepley 
111684f0b6dfSMatthew G. Knepley   Input Parameters:
111784f0b6dfSMatthew G. Knepley + label - the DMLabel
111884f0b6dfSMatthew G. Knepley - value - the stratum value
111984f0b6dfSMatthew G. Knepley 
112084f0b6dfSMatthew G. Knepley   Output Paramaters:
112184f0b6dfSMatthew G. Knepley + start - the smallest point in the stratum
112284f0b6dfSMatthew G. Knepley - end - the largest point in the stratum
112384f0b6dfSMatthew G. Knepley 
112484f0b6dfSMatthew G. Knepley   Level: intermediate
112584f0b6dfSMatthew G. Knepley 
112684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
112784f0b6dfSMatthew G. Knepley @*/
1128c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
1129c58f1c22SToby Isaac {
11300c3c4a36SLisandro Dalcin   PetscInt       v, min, max;
1131c58f1c22SToby Isaac   PetscErrorCode ierr;
1132c58f1c22SToby Isaac 
1133c58f1c22SToby Isaac   PetscFunctionBegin;
1134d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1135c58f1c22SToby Isaac   if (start) {PetscValidPointer(start, 3); *start = 0;}
1136c58f1c22SToby Isaac   if (end)   {PetscValidPointer(end,   4); *end   = 0;}
11370c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
11380c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
1139c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
11400c3c4a36SLisandro Dalcin   if (label->stratumSizes[v] <= 0) PetscFunctionReturn(0);
1141d6cb179aSToby Isaac   ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr);
1142d6cb179aSToby Isaac   if (start) *start = min;
1143d6cb179aSToby Isaac   if (end)   *end   = max+1;
1144c58f1c22SToby Isaac   PetscFunctionReturn(0);
1145c58f1c22SToby Isaac }
1146c58f1c22SToby Isaac 
114784f0b6dfSMatthew G. Knepley /*@
114884f0b6dfSMatthew G. Knepley   DMLabelGetStratumIS - Get an IS with the stratum points
114984f0b6dfSMatthew G. Knepley 
11505b5e7992SMatthew G. Knepley   Not collective
11515b5e7992SMatthew G. Knepley 
115284f0b6dfSMatthew G. Knepley   Input Parameters:
115384f0b6dfSMatthew G. Knepley + label - the DMLabel
115484f0b6dfSMatthew G. Knepley - value - the stratum value
115584f0b6dfSMatthew G. Knepley 
115684f0b6dfSMatthew G. Knepley   Output Paramater:
115784f0b6dfSMatthew G. Knepley . points - The stratum points
115884f0b6dfSMatthew G. Knepley 
115984f0b6dfSMatthew G. Knepley   Level: intermediate
116084f0b6dfSMatthew G. Knepley 
1161593ce467SVaclav Hapla   Notes:
1162593ce467SVaclav Hapla   The output IS should be destroyed when no longer needed.
1163593ce467SVaclav Hapla   Returns NULL if the stratum is empty.
1164593ce467SVaclav Hapla 
116584f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
116684f0b6dfSMatthew G. Knepley @*/
1167c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
1168c58f1c22SToby Isaac {
1169c58f1c22SToby Isaac   PetscInt       v;
1170c58f1c22SToby Isaac   PetscErrorCode ierr;
1171c58f1c22SToby Isaac 
1172c58f1c22SToby Isaac   PetscFunctionBegin;
1173d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1174c58f1c22SToby Isaac   PetscValidPointer(points, 3);
1175c58f1c22SToby Isaac   *points = NULL;
11760c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
11770c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
1178c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
1179ad8374ffSToby Isaac   ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr);
1180ad8374ffSToby Isaac   *points = label->points[v];
1181c58f1c22SToby Isaac   PetscFunctionReturn(0);
1182c58f1c22SToby Isaac }
1183c58f1c22SToby Isaac 
118484f0b6dfSMatthew G. Knepley /*@
118584f0b6dfSMatthew G. Knepley   DMLabelSetStratumIS - Set the stratum points using an IS
118684f0b6dfSMatthew G. Knepley 
11875b5e7992SMatthew G. Knepley   Not collective
11885b5e7992SMatthew G. Knepley 
118984f0b6dfSMatthew G. Knepley   Input Parameters:
119084f0b6dfSMatthew G. Knepley + label - the DMLabel
119184f0b6dfSMatthew G. Knepley . value - the stratum value
119284f0b6dfSMatthew G. Knepley - points - The stratum points
119384f0b6dfSMatthew G. Knepley 
119484f0b6dfSMatthew G. Knepley   Level: intermediate
119584f0b6dfSMatthew G. Knepley 
119684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
119784f0b6dfSMatthew G. Knepley @*/
11984de306b1SToby Isaac PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
11994de306b1SToby Isaac {
12000c3c4a36SLisandro Dalcin   PetscInt       v;
12014de306b1SToby Isaac   PetscErrorCode ierr;
12024de306b1SToby Isaac 
12034de306b1SToby Isaac   PetscFunctionBegin;
1204d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1205d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
1206b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
12074de306b1SToby Isaac   if (is == label->points[v]) PetscFunctionReturn(0);
12084de306b1SToby Isaac   ierr = DMLabelClearStratum(label, value);CHKERRQ(ierr);
12094de306b1SToby Isaac   ierr = ISGetLocalSize(is, &(label->stratumSizes[v]));CHKERRQ(ierr);
12104de306b1SToby Isaac   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
12114de306b1SToby Isaac   ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
12120c3c4a36SLisandro Dalcin   label->points[v]  = is;
12130c3c4a36SLisandro Dalcin   label->validIS[v] = PETSC_TRUE;
1214277ea44aSLisandro Dalcin   ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
12154de306b1SToby Isaac   if (label->bt) {
12164de306b1SToby Isaac     const PetscInt *points;
12174de306b1SToby Isaac     PetscInt p;
12184de306b1SToby Isaac 
12194de306b1SToby Isaac     ierr = ISGetIndices(is,&points);CHKERRQ(ierr);
12204de306b1SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
12214de306b1SToby Isaac       const PetscInt point = points[p];
12224de306b1SToby Isaac 
12234de306b1SToby 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);
12244de306b1SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
12254de306b1SToby Isaac     }
12264de306b1SToby Isaac   }
12274de306b1SToby Isaac   PetscFunctionReturn(0);
12284de306b1SToby Isaac }
12294de306b1SToby Isaac 
123084f0b6dfSMatthew G. Knepley /*@
123184f0b6dfSMatthew G. Knepley   DMLabelClearStratum - Remove a stratum
12324de306b1SToby Isaac 
12335b5e7992SMatthew G. Knepley   Not collective
12345b5e7992SMatthew G. Knepley 
123584f0b6dfSMatthew G. Knepley   Input Parameters:
123684f0b6dfSMatthew G. Knepley + label - the DMLabel
123784f0b6dfSMatthew G. Knepley - value - the stratum value
123884f0b6dfSMatthew G. Knepley 
123984f0b6dfSMatthew G. Knepley   Level: intermediate
124084f0b6dfSMatthew G. Knepley 
124184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
124284f0b6dfSMatthew G. Knepley @*/
1243c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
1244c58f1c22SToby Isaac {
1245c58f1c22SToby Isaac   PetscInt       v;
1246c58f1c22SToby Isaac   PetscErrorCode ierr;
1247c58f1c22SToby Isaac 
1248c58f1c22SToby Isaac   PetscFunctionBegin;
1249d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
12500c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
12510c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
12524de306b1SToby Isaac   if (label->validIS[v]) {
12534de306b1SToby Isaac     if (label->bt) {
1254c58f1c22SToby Isaac       PetscInt       i;
1255ad8374ffSToby Isaac       const PetscInt *points;
1256c58f1c22SToby Isaac 
1257ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
1258c58f1c22SToby Isaac       for (i = 0; i < label->stratumSizes[v]; ++i) {
1259ad8374ffSToby Isaac         const PetscInt point = points[i];
1260c58f1c22SToby Isaac 
1261c58f1c22SToby 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);
1262c58f1c22SToby Isaac         ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
1263c58f1c22SToby Isaac       }
1264ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
1265c58f1c22SToby Isaac     }
1266c58f1c22SToby Isaac     label->stratumSizes[v] = 0;
12670c3c4a36SLisandro Dalcin     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
1268277ea44aSLisandro Dalcin     ierr = ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &label->points[v]);CHKERRQ(ierr);
12690c3c4a36SLisandro Dalcin     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
1270277ea44aSLisandro Dalcin     ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
1271c58f1c22SToby Isaac   } else {
1272b9907514SLisandro Dalcin     ierr = PetscHSetIClear(label->ht[v]);CHKERRQ(ierr);
1273c58f1c22SToby Isaac   }
1274c58f1c22SToby Isaac   PetscFunctionReturn(0);
1275c58f1c22SToby Isaac }
1276c58f1c22SToby Isaac 
127784f0b6dfSMatthew G. Knepley /*@
127884f0b6dfSMatthew G. Knepley   DMLabelFilter - Remove all points outside of [start, end)
127984f0b6dfSMatthew G. Knepley 
12805b5e7992SMatthew G. Knepley   Not collective
12815b5e7992SMatthew G. Knepley 
128284f0b6dfSMatthew G. Knepley   Input Parameters:
128384f0b6dfSMatthew G. Knepley + label - the DMLabel
128422cd2cdaSVaclav Hapla . start - the first point kept
128522cd2cdaSVaclav Hapla - end - one more than the last point kept
128684f0b6dfSMatthew G. Knepley 
128784f0b6dfSMatthew G. Knepley   Level: intermediate
128884f0b6dfSMatthew G. Knepley 
128984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
129084f0b6dfSMatthew G. Knepley @*/
1291c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
1292c58f1c22SToby Isaac {
1293c58f1c22SToby Isaac   PetscInt       v;
1294c58f1c22SToby Isaac   PetscErrorCode ierr;
1295c58f1c22SToby Isaac 
1296c58f1c22SToby Isaac   PetscFunctionBegin;
1297d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
12980c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
1299c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1300c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
13019106b633SVaclav Hapla     ierr = ISGeneralFilter(label->points[v], start, end);CHKERRQ(ierr);
1302c58f1c22SToby Isaac   }
1303c58f1c22SToby Isaac   ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr);
1304c58f1c22SToby Isaac   PetscFunctionReturn(0);
1305c58f1c22SToby Isaac }
1306c58f1c22SToby Isaac 
130784f0b6dfSMatthew G. Knepley /*@
130884f0b6dfSMatthew G. Knepley   DMLabelPermute - Create a new label with permuted points
130984f0b6dfSMatthew G. Knepley 
13105b5e7992SMatthew G. Knepley   Not collective
13115b5e7992SMatthew G. Knepley 
131284f0b6dfSMatthew G. Knepley   Input Parameters:
131384f0b6dfSMatthew G. Knepley + label - the DMLabel
131484f0b6dfSMatthew G. Knepley - permutation - the point permutation
131584f0b6dfSMatthew G. Knepley 
131684f0b6dfSMatthew G. Knepley   Output Parameter:
131784f0b6dfSMatthew G. Knepley . labelnew - the new label containing the permuted points
131884f0b6dfSMatthew G. Knepley 
131984f0b6dfSMatthew G. Knepley   Level: intermediate
132084f0b6dfSMatthew G. Knepley 
132184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
132284f0b6dfSMatthew G. Knepley @*/
1323c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
1324c58f1c22SToby Isaac {
1325c58f1c22SToby Isaac   const PetscInt *perm;
1326c58f1c22SToby Isaac   PetscInt        numValues, numPoints, v, q;
1327c58f1c22SToby Isaac   PetscErrorCode  ierr;
1328c58f1c22SToby Isaac 
1329c58f1c22SToby Isaac   PetscFunctionBegin;
1330d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1331d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(permutation, IS_CLASSID, 2);
1332c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1333c58f1c22SToby Isaac   ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr);
1334c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr);
1335c58f1c22SToby Isaac   ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr);
1336c58f1c22SToby Isaac   ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr);
1337c58f1c22SToby Isaac   for (v = 0; v < numValues; ++v) {
1338c58f1c22SToby Isaac     const PetscInt size   = (*labelNew)->stratumSizes[v];
1339ad8374ffSToby Isaac     const PetscInt *points;
1340ad8374ffSToby Isaac     PetscInt *pointsNew;
1341c58f1c22SToby Isaac 
1342ad8374ffSToby Isaac     ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1343ad8374ffSToby Isaac     ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr);
1344c58f1c22SToby Isaac     for (q = 0; q < size; ++q) {
1345ad8374ffSToby Isaac       const PetscInt point = points[q];
1346c58f1c22SToby Isaac 
1347c58f1c22SToby 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);
1348ad8374ffSToby Isaac       pointsNew[q] = perm[point];
1349c58f1c22SToby Isaac     }
1350ad8374ffSToby Isaac     ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1351ad8374ffSToby Isaac     ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr);
1352ad8374ffSToby Isaac     ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr);
1353fa8e8ae5SToby Isaac     if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) {
1354fa8e8ae5SToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));CHKERRQ(ierr);
1355fa8e8ae5SToby Isaac       ierr = PetscFree(pointsNew);CHKERRQ(ierr);
1356fa8e8ae5SToby Isaac     } else {
1357ad8374ffSToby Isaac       ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr);
1358fa8e8ae5SToby Isaac     }
1359ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr);
1360c58f1c22SToby Isaac   }
1361c58f1c22SToby Isaac   ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr);
1362c58f1c22SToby Isaac   if (label->bt) {
1363c58f1c22SToby Isaac     ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
1364c58f1c22SToby Isaac     ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
1365c58f1c22SToby Isaac   }
1366c58f1c22SToby Isaac   PetscFunctionReturn(0);
1367c58f1c22SToby Isaac }
1368c58f1c22SToby Isaac 
136926c55118SMichael Lange PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
137026c55118SMichael Lange {
137126c55118SMichael Lange   MPI_Comm       comm;
137226c55118SMichael Lange   PetscInt       s, l, nroots, nleaves, dof, offset, size;
137326c55118SMichael Lange   PetscInt      *remoteOffsets, *rootStrata, *rootIdx;
137426c55118SMichael Lange   PetscSection   rootSection;
137526c55118SMichael Lange   PetscSF        labelSF;
137626c55118SMichael Lange   PetscErrorCode ierr;
137726c55118SMichael Lange 
137826c55118SMichael Lange   PetscFunctionBegin;
137926c55118SMichael Lange   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
138026c55118SMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
138126c55118SMichael Lange   /* Build a section of stratum values per point, generate the according SF
138226c55118SMichael Lange      and distribute point-wise stratum values to leaves. */
138326c55118SMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
138426c55118SMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
138526c55118SMichael Lange   ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr);
138626c55118SMichael Lange   if (label) {
138726c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1388ad8374ffSToby Isaac       const PetscInt *points;
1389ad8374ffSToby Isaac 
1390ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
139126c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1392ad8374ffSToby Isaac         ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr);
1393ad8374ffSToby Isaac         ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr);
139426c55118SMichael Lange       }
1395ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
139626c55118SMichael Lange     }
139726c55118SMichael Lange   }
139826c55118SMichael Lange   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
139926c55118SMichael Lange   /* Create a point-wise array of stratum values */
140026c55118SMichael Lange   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
140126c55118SMichael Lange   ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr);
140226c55118SMichael Lange   ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr);
140326c55118SMichael Lange   if (label) {
140426c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1405ad8374ffSToby Isaac       const PetscInt *points;
1406ad8374ffSToby Isaac 
1407ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
140826c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1409ad8374ffSToby Isaac         const PetscInt p = points[l];
141026c55118SMichael Lange         ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr);
141126c55118SMichael Lange         rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
141226c55118SMichael Lange       }
1413ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
141426c55118SMichael Lange     }
141526c55118SMichael Lange   }
141626c55118SMichael Lange   /* Build SF that maps label points to remote processes */
141726c55118SMichael Lange   ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr);
141826c55118SMichael Lange   ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr);
141926c55118SMichael Lange   ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr);
142026c55118SMichael Lange   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
142126c55118SMichael Lange   /* Send the strata for each point over the derived SF */
142226c55118SMichael Lange   ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr);
142326c55118SMichael Lange   ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr);
142426c55118SMichael Lange   ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
142526c55118SMichael Lange   ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
142626c55118SMichael Lange   /* Clean up */
142726c55118SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
142826c55118SMichael Lange   ierr = PetscFree(rootIdx);CHKERRQ(ierr);
142926c55118SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
143026c55118SMichael Lange   ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr);
143126c55118SMichael Lange   PetscFunctionReturn(0);
143226c55118SMichael Lange }
143326c55118SMichael Lange 
143484f0b6dfSMatthew G. Knepley /*@
143584f0b6dfSMatthew G. Knepley   DMLabelDistribute - Create a new label pushed forward over the PetscSF
143684f0b6dfSMatthew G. Knepley 
14375b5e7992SMatthew G. Knepley   Collective on sf
14385b5e7992SMatthew G. Knepley 
143984f0b6dfSMatthew G. Knepley   Input Parameters:
144084f0b6dfSMatthew G. Knepley + label - the DMLabel
144184f0b6dfSMatthew G. Knepley - sf    - the map from old to new distribution
144284f0b6dfSMatthew G. Knepley 
144384f0b6dfSMatthew G. Knepley   Output Parameter:
144484f0b6dfSMatthew G. Knepley . labelnew - the new redistributed label
144584f0b6dfSMatthew G. Knepley 
144684f0b6dfSMatthew G. Knepley   Level: intermediate
144784f0b6dfSMatthew G. Knepley 
144884f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
144984f0b6dfSMatthew G. Knepley @*/
1450c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1451c58f1c22SToby Isaac {
1452c58f1c22SToby Isaac   MPI_Comm       comm;
145326c55118SMichael Lange   PetscSection   leafSection;
145426c55118SMichael Lange   PetscInt       p, pStart, pEnd, s, size, dof, offset, stratum;
145526c55118SMichael Lange   PetscInt      *leafStrata, *strataIdx;
1456ad8374ffSToby Isaac   PetscInt     **points;
1457d67d17b1SMatthew G. Knepley   const char    *lname = NULL;
1458c58f1c22SToby Isaac   char          *name;
1459c58f1c22SToby Isaac   PetscInt       nameSize;
1460e8f14785SLisandro Dalcin   PetscHSetI     stratumHash;
1461c58f1c22SToby Isaac   size_t         len = 0;
146226c55118SMichael Lange   PetscMPIInt    rank;
1463c58f1c22SToby Isaac   PetscErrorCode ierr;
1464c58f1c22SToby Isaac 
1465c58f1c22SToby Isaac   PetscFunctionBegin;
1466d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1467f018e600SMatthew G. Knepley   if (label) {
1468f018e600SMatthew G. Knepley     PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1469f018e600SMatthew G. Knepley     ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1470f018e600SMatthew G. Knepley   }
1471c58f1c22SToby Isaac   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1472c58f1c22SToby Isaac   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1473c58f1c22SToby Isaac   /* Bcast name */
1474d67d17b1SMatthew G. Knepley   if (!rank) {
1475d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1476d67d17b1SMatthew G. Knepley     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1477d67d17b1SMatthew G. Knepley   }
1478c58f1c22SToby Isaac   nameSize = len;
1479c58f1c22SToby Isaac   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1480c58f1c22SToby Isaac   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1481580bdb30SBarry Smith   if (!rank) {ierr = PetscArraycpy(name, lname, nameSize+1);CHKERRQ(ierr);}
1482c58f1c22SToby Isaac   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1483d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
1484c58f1c22SToby Isaac   ierr = PetscFree(name);CHKERRQ(ierr);
148577d236dfSMichael Lange   /* Bcast defaultValue */
148677d236dfSMichael Lange   if (!rank) (*labelNew)->defaultValue = label->defaultValue;
148777d236dfSMichael Lange   ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
148826c55118SMichael Lange   /* Distribute stratum values over the SF and get the point mapping on the receiver */
148926c55118SMichael Lange   ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr);
14905cbdf6fcSMichael Lange   /* Determine received stratum values and initialise new label*/
1491e8f14785SLisandro Dalcin   ierr = PetscHSetICreate(&stratumHash);CHKERRQ(ierr);
149226c55118SMichael Lange   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
1493e8f14785SLisandro Dalcin   for (p = 0; p < size; ++p) {ierr = PetscHSetIAdd(stratumHash, leafStrata[p]);CHKERRQ(ierr);}
1494e8f14785SLisandro Dalcin   ierr = PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata);CHKERRQ(ierr);
1495ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr);
1496ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
14975cbdf6fcSMichael Lange   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr);
14985cbdf6fcSMichael Lange   /* Turn leafStrata into indices rather than stratum values */
14995cbdf6fcSMichael Lange   offset = 0;
1500e8f14785SLisandro Dalcin   ierr = PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr);
1501a205f1fdSToby Isaac   ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr);
150290e9b2aeSLisandro Dalcin   for (s = 0; s < (*labelNew)->numStrata; ++s) {
150390e9b2aeSLisandro Dalcin     ierr = PetscHMapISet((*labelNew)->hmap, (*labelNew)->stratumValues[s], s);CHKERRQ(ierr);
150490e9b2aeSLisandro Dalcin   }
15055cbdf6fcSMichael Lange   for (p = 0; p < size; ++p) {
1506231b9e6fSMatthew G. Knepley     for (s = 0; s < (*labelNew)->numStrata; ++s) {
1507231b9e6fSMatthew G. Knepley       if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
15085cbdf6fcSMichael Lange     }
15095cbdf6fcSMichael Lange   }
1510c58f1c22SToby Isaac   /* Rebuild the point strata on the receiver */
1511c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr);
1512c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1513c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1514c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1515c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1516c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1517c58f1c22SToby Isaac       (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1518c58f1c22SToby Isaac     }
1519c58f1c22SToby Isaac   }
1520c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr);
1521c58f1c22SToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr);
1522ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr);
1523c58f1c22SToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1524e8f14785SLisandro Dalcin     ierr = PetscHSetICreate(&(*labelNew)->ht[s]);CHKERRQ(ierr);
1525ad8374ffSToby Isaac     ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr);
1526c58f1c22SToby Isaac   }
1527c58f1c22SToby Isaac   /* Insert points into new strata */
1528c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr);
1529c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1530c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1531c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1532c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1533c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1534c58f1c22SToby Isaac       stratum = leafStrata[offset+s];
1535ad8374ffSToby Isaac       points[stratum][strataIdx[stratum]++] = p;
1536c58f1c22SToby Isaac     }
1537c58f1c22SToby Isaac   }
1538ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; s++) {
1539ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr);
1540ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr);
1541ad8374ffSToby Isaac   }
1542ad8374ffSToby Isaac   ierr = PetscFree(points);CHKERRQ(ierr);
1543e8f14785SLisandro Dalcin   ierr = PetscHSetIDestroy(&stratumHash);CHKERRQ(ierr);
1544c58f1c22SToby Isaac   ierr = PetscFree(leafStrata);CHKERRQ(ierr);
1545c58f1c22SToby Isaac   ierr = PetscFree(strataIdx);CHKERRQ(ierr);
1546c58f1c22SToby Isaac   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
1547c58f1c22SToby Isaac   PetscFunctionReturn(0);
1548c58f1c22SToby Isaac }
1549c58f1c22SToby Isaac 
15507937d9ceSMichael Lange /*@
15517937d9ceSMichael Lange   DMLabelGather - Gather all label values from leafs into roots
15527937d9ceSMichael Lange 
15535b5e7992SMatthew G. Knepley   Collective on sf
15545b5e7992SMatthew G. Knepley 
15557937d9ceSMichael Lange   Input Parameters:
15567937d9ceSMichael Lange + label - the DMLabel
155784f0b6dfSMatthew G. Knepley - sf - the Star Forest point communication map
15587937d9ceSMichael Lange 
155984f0b6dfSMatthew G. Knepley   Output Parameters:
156084f0b6dfSMatthew G. Knepley . labelNew - the new DMLabel with localised leaf values
15617937d9ceSMichael Lange 
15627937d9ceSMichael Lange   Level: developer
15637937d9ceSMichael Lange 
15647937d9ceSMichael Lange   Note: This is the inverse operation to DMLabelDistribute.
15657937d9ceSMichael Lange 
15667937d9ceSMichael Lange .seealso: DMLabelDistribute()
15677937d9ceSMichael Lange @*/
15687937d9ceSMichael Lange PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
15697937d9ceSMichael Lange {
15707937d9ceSMichael Lange   MPI_Comm       comm;
15717937d9ceSMichael Lange   PetscSection   rootSection;
15727937d9ceSMichael Lange   PetscSF        sfLabel;
15737937d9ceSMichael Lange   PetscSFNode   *rootPoints, *leafPoints;
15747937d9ceSMichael Lange   PetscInt       p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
15757937d9ceSMichael Lange   const PetscInt *rootDegree, *ilocal;
15767937d9ceSMichael Lange   PetscInt       *rootStrata;
1577d67d17b1SMatthew G. Knepley   const char    *lname;
15787937d9ceSMichael Lange   char          *name;
15797937d9ceSMichael Lange   PetscInt       nameSize;
15807937d9ceSMichael Lange   size_t         len = 0;
15819852e123SBarry Smith   PetscMPIInt    rank, size;
15827937d9ceSMichael Lange   PetscErrorCode ierr;
15837937d9ceSMichael Lange 
15847937d9ceSMichael Lange   PetscFunctionBegin;
1585d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1586d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
15877937d9ceSMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
15887937d9ceSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
15899852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
15907937d9ceSMichael Lange   /* Bcast name */
1591d67d17b1SMatthew G. Knepley   if (!rank) {
1592d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1593d67d17b1SMatthew G. Knepley     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1594d67d17b1SMatthew G. Knepley   }
15957937d9ceSMichael Lange   nameSize = len;
15967937d9ceSMichael Lange   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
15977937d9ceSMichael Lange   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1598580bdb30SBarry Smith   if (!rank) {ierr = PetscArraycpy(name, lname, nameSize+1);CHKERRQ(ierr);}
15997937d9ceSMichael Lange   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1600d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
16017937d9ceSMichael Lange   ierr = PetscFree(name);CHKERRQ(ierr);
16027937d9ceSMichael Lange   /* Gather rank/index pairs of leaves into local roots to build
16037937d9ceSMichael Lange      an inverse, multi-rooted SF. Note that this ignores local leaf
16047937d9ceSMichael Lange      indexing due to the use of the multiSF in PetscSFGather. */
16057937d9ceSMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr);
1606dc53bc9bSMatthew G. Knepley   ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr);
1607dc53bc9bSMatthew G. Knepley   for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
16087937d9ceSMichael Lange   for (p = 0; p < nleaves; p++) {
16098212dd46SStefano Zampini     PetscInt ilp = ilocal ? ilocal[p] : p;
16108212dd46SStefano Zampini 
16118212dd46SStefano Zampini     leafPoints[ilp].index = ilp;
16128212dd46SStefano Zampini     leafPoints[ilp].rank  = rank;
16137937d9ceSMichael Lange   }
16147937d9ceSMichael Lange   ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr);
16157937d9ceSMichael Lange   ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr);
16167937d9ceSMichael Lange   for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
16177937d9ceSMichael Lange   ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr);
16187937d9ceSMichael Lange   ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
16197937d9ceSMichael Lange   ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
16207937d9ceSMichael Lange   ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr);
16217937d9ceSMichael Lange   ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
16227937d9ceSMichael Lange   /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
16237937d9ceSMichael Lange   ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr);
16247937d9ceSMichael Lange   /* Rebuild the point strata on the receiver */
16257937d9ceSMichael Lange   for (p = 0, idx = 0; p < nroots; p++) {
16267937d9ceSMichael Lange     for (d = 0; d < rootDegree[p]; d++) {
16277937d9ceSMichael Lange       ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr);
16287937d9ceSMichael Lange       ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr);
16297937d9ceSMichael Lange       for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);}
16307937d9ceSMichael Lange     }
16317937d9ceSMichael Lange     idx += rootDegree[p];
16327937d9ceSMichael Lange   }
163377e0c0e7SMichael Lange   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
163477e0c0e7SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
163577e0c0e7SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
163677e0c0e7SMichael Lange   ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr);
16377937d9ceSMichael Lange   PetscFunctionReturn(0);
16387937d9ceSMichael Lange }
16397937d9ceSMichael Lange 
164084f0b6dfSMatthew G. Knepley /*@
164184f0b6dfSMatthew G. Knepley   DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label
164284f0b6dfSMatthew G. Knepley 
16435b5e7992SMatthew G. Knepley   Not collective
16445b5e7992SMatthew G. Knepley 
164584f0b6dfSMatthew G. Knepley   Input Parameter:
164684f0b6dfSMatthew G. Knepley . label - the DMLabel
164784f0b6dfSMatthew G. Knepley 
164884f0b6dfSMatthew G. Knepley   Output Parameters:
164984f0b6dfSMatthew G. Knepley + section - the section giving offsets for each stratum
165084f0b6dfSMatthew G. Knepley - is - An IS containing all the label points
165184f0b6dfSMatthew G. Knepley 
165284f0b6dfSMatthew G. Knepley   Level: developer
165384f0b6dfSMatthew G. Knepley 
165484f0b6dfSMatthew G. Knepley .seealso: DMLabelDistribute()
165584f0b6dfSMatthew G. Knepley @*/
1656c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1657c58f1c22SToby Isaac {
1658c58f1c22SToby Isaac   IS              vIS;
1659c58f1c22SToby Isaac   const PetscInt *values;
1660c58f1c22SToby Isaac   PetscInt       *points;
1661c58f1c22SToby Isaac   PetscInt        nV, vS = 0, vE = 0, v, N;
1662c58f1c22SToby Isaac   PetscErrorCode  ierr;
1663c58f1c22SToby Isaac 
1664c58f1c22SToby Isaac   PetscFunctionBegin;
1665d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1666c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr);
1667c58f1c22SToby Isaac   ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
1668c58f1c22SToby Isaac   ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr);
1669c58f1c22SToby Isaac   if (nV) {vS = values[0]; vE = values[0]+1;}
1670c58f1c22SToby Isaac   for (v = 1; v < nV; ++v) {
1671c58f1c22SToby Isaac     vS = PetscMin(vS, values[v]);
1672c58f1c22SToby Isaac     vE = PetscMax(vE, values[v]+1);
1673c58f1c22SToby Isaac   }
1674c58f1c22SToby Isaac   ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr);
1675c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr);
1676c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1677c58f1c22SToby Isaac     PetscInt n;
1678c58f1c22SToby Isaac 
1679c58f1c22SToby Isaac     ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr);
1680c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr);
1681c58f1c22SToby Isaac   }
1682c58f1c22SToby Isaac   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
1683c58f1c22SToby Isaac   ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr);
1684c58f1c22SToby Isaac   ierr = PetscMalloc1(N, &points);CHKERRQ(ierr);
1685c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1686c58f1c22SToby Isaac     IS              is;
1687c58f1c22SToby Isaac     const PetscInt *spoints;
1688c58f1c22SToby Isaac     PetscInt        dof, off, p;
1689c58f1c22SToby Isaac 
1690c58f1c22SToby Isaac     ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr);
1691c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr);
1692c58f1c22SToby Isaac     ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr);
1693c58f1c22SToby Isaac     ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr);
1694c58f1c22SToby Isaac     for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1695c58f1c22SToby Isaac     ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr);
1696c58f1c22SToby Isaac     ierr = ISDestroy(&is);CHKERRQ(ierr);
1697c58f1c22SToby Isaac   }
1698c58f1c22SToby Isaac   ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr);
1699c58f1c22SToby Isaac   ierr = ISDestroy(&vIS);CHKERRQ(ierr);
1700c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1701c58f1c22SToby Isaac   PetscFunctionReturn(0);
1702c58f1c22SToby Isaac }
1703c58f1c22SToby Isaac 
170484f0b6dfSMatthew G. Knepley /*@
1705c58f1c22SToby Isaac   PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1706c58f1c22SToby Isaac   the local section and an SF describing the section point overlap.
1707c58f1c22SToby Isaac 
17085b5e7992SMatthew G. Knepley   Collective on sf
17095b5e7992SMatthew G. Knepley 
1710c58f1c22SToby Isaac   Input Parameters:
1711c58f1c22SToby Isaac   + s - The PetscSection for the local field layout
1712c58f1c22SToby Isaac   . sf - The SF describing parallel layout of the section points
1713c58f1c22SToby Isaac   . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1714c58f1c22SToby Isaac   . label - The label specifying the points
1715c58f1c22SToby Isaac   - labelValue - The label stratum specifying the points
1716c58f1c22SToby Isaac 
1717c58f1c22SToby Isaac   Output Parameter:
1718c58f1c22SToby Isaac   . gsection - The PetscSection for the global field layout
1719c58f1c22SToby Isaac 
1720c58f1c22SToby Isaac   Note: This gives negative sizes and offsets to points not owned by this process
1721c58f1c22SToby Isaac 
1722c58f1c22SToby Isaac   Level: developer
1723c58f1c22SToby Isaac 
1724c58f1c22SToby Isaac .seealso: PetscSectionCreate()
1725c58f1c22SToby Isaac @*/
1726c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1727c58f1c22SToby Isaac {
1728c58f1c22SToby Isaac   PetscInt      *neg = NULL, *tmpOff = NULL;
1729c58f1c22SToby Isaac   PetscInt       pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1730c58f1c22SToby Isaac   PetscErrorCode ierr;
1731c58f1c22SToby Isaac 
1732c58f1c22SToby Isaac   PetscFunctionBegin;
1733d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1734d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1735d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4);
1736c58f1c22SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr);
1737c58f1c22SToby Isaac   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1738c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr);
1739c58f1c22SToby Isaac   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1740c58f1c22SToby Isaac   if (nroots >= 0) {
1741c58f1c22SToby Isaac     if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1742c58f1c22SToby Isaac     ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr);
1743c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1744c58f1c22SToby Isaac       ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr);
1745c58f1c22SToby Isaac     } else {
1746c58f1c22SToby Isaac       tmpOff = &(*gsection)->atlasDof[-pStart];
1747c58f1c22SToby Isaac     }
1748c58f1c22SToby Isaac   }
1749c58f1c22SToby Isaac   /* Mark ghost points with negative dof */
1750c58f1c22SToby Isaac   for (p = pStart; p < pEnd; ++p) {
1751c58f1c22SToby Isaac     PetscInt value;
1752c58f1c22SToby Isaac 
1753c58f1c22SToby Isaac     ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr);
1754c58f1c22SToby Isaac     if (value != labelValue) continue;
1755c58f1c22SToby Isaac     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1756c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr);
1757c58f1c22SToby Isaac     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1758c58f1c22SToby Isaac     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);}
1759c58f1c22SToby Isaac     if (neg) neg[p] = -(dof+1);
1760c58f1c22SToby Isaac   }
1761c58f1c22SToby Isaac   ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr);
1762c58f1c22SToby Isaac   if (nroots >= 0) {
1763c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1764c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1765c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1766c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1767c58f1c22SToby Isaac     }
1768c58f1c22SToby Isaac   }
1769c58f1c22SToby Isaac   /* Calculate new sizes, get proccess offset, and calculate point offsets */
1770c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1771c58f1c22SToby Isaac     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1772c58f1c22SToby Isaac     (*gsection)->atlasOff[p] = off;
1773c58f1c22SToby Isaac     off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1774c58f1c22SToby Isaac   }
1775c58f1c22SToby Isaac   ierr       = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1776c58f1c22SToby Isaac   globalOff -= off;
1777c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1778c58f1c22SToby Isaac     (*gsection)->atlasOff[p] += globalOff;
1779c58f1c22SToby Isaac     if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1780c58f1c22SToby Isaac   }
1781c58f1c22SToby Isaac   /* Put in negative offsets for ghost points */
1782c58f1c22SToby Isaac   if (nroots >= 0) {
1783c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1784c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1785c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1786c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1787c58f1c22SToby Isaac     }
1788c58f1c22SToby Isaac   }
1789c58f1c22SToby Isaac   if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);}
1790c58f1c22SToby Isaac   ierr = PetscFree(neg);CHKERRQ(ierr);
1791c58f1c22SToby Isaac   PetscFunctionReturn(0);
1792c58f1c22SToby Isaac }
1793c58f1c22SToby Isaac 
17945fdea053SToby Isaac typedef struct _n_PetscSectionSym_Label
17955fdea053SToby Isaac {
17965fdea053SToby Isaac   DMLabel           label;
17975fdea053SToby Isaac   PetscCopyMode     *modes;
17985fdea053SToby Isaac   PetscInt          *sizes;
17995fdea053SToby Isaac   const PetscInt    ***perms;
18005fdea053SToby Isaac   const PetscScalar ***rots;
18015fdea053SToby Isaac   PetscInt          (*minMaxOrients)[2];
18025fdea053SToby Isaac   PetscInt          numStrata; /* numStrata is only increasing, functions as a state */
18035fdea053SToby Isaac } PetscSectionSym_Label;
18045fdea053SToby Isaac 
18055fdea053SToby Isaac static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
18065fdea053SToby Isaac {
18075fdea053SToby Isaac   PetscInt              i, j;
18085fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
18095fdea053SToby Isaac   PetscErrorCode        ierr;
18105fdea053SToby Isaac 
18115fdea053SToby Isaac   PetscFunctionBegin;
18125fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
18135fdea053SToby Isaac     if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
18145fdea053SToby Isaac       for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
18155fdea053SToby Isaac         if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);}
18165fdea053SToby Isaac         if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);}
18175fdea053SToby Isaac       }
18185fdea053SToby Isaac       if (sl->perms[i]) {
18195fdea053SToby Isaac         const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
18205fdea053SToby Isaac 
18215fdea053SToby Isaac         ierr = PetscFree(perms);CHKERRQ(ierr);
18225fdea053SToby Isaac       }
18235fdea053SToby Isaac       if (sl->rots[i]) {
18245fdea053SToby Isaac         const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
18255fdea053SToby Isaac 
18265fdea053SToby Isaac         ierr = PetscFree(rots);CHKERRQ(ierr);
18275fdea053SToby Isaac       }
18285fdea053SToby Isaac     }
18295fdea053SToby Isaac   }
18305fdea053SToby Isaac   ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr);
18315fdea053SToby Isaac   ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr);
18325fdea053SToby Isaac   sl->numStrata = 0;
18335fdea053SToby Isaac   PetscFunctionReturn(0);
18345fdea053SToby Isaac }
18355fdea053SToby Isaac 
18365fdea053SToby Isaac static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
18375fdea053SToby Isaac {
18385fdea053SToby Isaac   PetscErrorCode ierr;
18395fdea053SToby Isaac 
18405fdea053SToby Isaac   PetscFunctionBegin;
18415fdea053SToby Isaac   ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);
18425fdea053SToby Isaac   ierr = PetscFree(sym->data);CHKERRQ(ierr);
18435fdea053SToby Isaac   PetscFunctionReturn(0);
18445fdea053SToby Isaac }
18455fdea053SToby Isaac 
18465fdea053SToby Isaac static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
18475fdea053SToby Isaac {
18485fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
18495fdea053SToby Isaac   PetscBool             isAscii;
18505fdea053SToby Isaac   DMLabel               label = sl->label;
1851d67d17b1SMatthew G. Knepley   const char           *name;
18525fdea053SToby Isaac   PetscErrorCode        ierr;
18535fdea053SToby Isaac 
18545fdea053SToby Isaac   PetscFunctionBegin;
18555fdea053SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
18565fdea053SToby Isaac   if (isAscii) {
18575fdea053SToby Isaac     PetscInt          i, j, k;
18585fdea053SToby Isaac     PetscViewerFormat format;
18595fdea053SToby Isaac 
18605fdea053SToby Isaac     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
18615fdea053SToby Isaac     if (label) {
18625fdea053SToby Isaac       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
18635fdea053SToby Isaac       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
18645fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18655fdea053SToby Isaac         ierr = DMLabelView(label, viewer);CHKERRQ(ierr);
18665fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
18675fdea053SToby Isaac       } else {
1868d67d17b1SMatthew G. Knepley         ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1869d67d17b1SMatthew G. Knepley         ierr = PetscViewerASCIIPrintf(viewer,"  Label '%s'\n",name);CHKERRQ(ierr);
18705fdea053SToby Isaac       }
18715fdea053SToby Isaac     } else {
18725fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr);
18735fdea053SToby Isaac     }
18745fdea053SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18755fdea053SToby Isaac     for (i = 0; i <= sl->numStrata; i++) {
18765fdea053SToby Isaac       PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
18775fdea053SToby Isaac 
18785fdea053SToby Isaac       if (!(sl->perms[i] || sl->rots[i])) {
18795fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr);
18805fdea053SToby Isaac       } else {
18815fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr);
18825fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18835fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr);
18845fdea053SToby Isaac         if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
18855fdea053SToby Isaac           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18865fdea053SToby Isaac           for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
18875fdea053SToby Isaac             if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
18885fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr);
18895fdea053SToby Isaac             } else {
18905fdea053SToby Isaac               PetscInt tab;
18915fdea053SToby Isaac 
18925fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr);
18935fdea053SToby Isaac               ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18945fdea053SToby Isaac               ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr);
18955fdea053SToby Isaac               if (sl->perms[i] && sl->perms[i][j]) {
18965fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr);
18975fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
18985fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);}
18995fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
19005fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
19015fdea053SToby Isaac               }
19025fdea053SToby Isaac               if (sl->rots[i] && sl->rots[i][j]) {
19035fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Rotations:  ");CHKERRQ(ierr);
19045fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
19055fdea053SToby Isaac #if defined(PETSC_USE_COMPLEX)
19065fdea053SToby 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);}
19075fdea053SToby Isaac #else
19085fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);}
19095fdea053SToby Isaac #endif
19105fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
19115fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
19125fdea053SToby Isaac               }
19135fdea053SToby Isaac               ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
19145fdea053SToby Isaac             }
19155fdea053SToby Isaac           }
19165fdea053SToby Isaac           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
19175fdea053SToby Isaac         }
19185fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
19195fdea053SToby Isaac       }
19205fdea053SToby Isaac     }
19215fdea053SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
19225fdea053SToby Isaac   }
19235fdea053SToby Isaac   PetscFunctionReturn(0);
19245fdea053SToby Isaac }
19255fdea053SToby Isaac 
19265fdea053SToby Isaac /*@
19275fdea053SToby Isaac   PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
19285fdea053SToby Isaac 
19295fdea053SToby Isaac   Logically collective on sym
19305fdea053SToby Isaac 
19315fdea053SToby Isaac   Input parameters:
19325fdea053SToby Isaac + sym - the section symmetries
19335fdea053SToby Isaac - label - the DMLabel describing the types of points
19345fdea053SToby Isaac 
19355fdea053SToby Isaac   Level: developer:
19365fdea053SToby Isaac 
19375fdea053SToby Isaac .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()
19385fdea053SToby Isaac @*/
19395fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
19405fdea053SToby Isaac {
19415fdea053SToby Isaac   PetscSectionSym_Label *sl;
19425fdea053SToby Isaac   PetscErrorCode        ierr;
19435fdea053SToby Isaac 
19445fdea053SToby Isaac   PetscFunctionBegin;
19455fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
19465fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
19475fdea053SToby Isaac   if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);}
19485fdea053SToby Isaac   if (label) {
19495fdea053SToby Isaac     sl->label = label;
1950d67d17b1SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) label);CHKERRQ(ierr);
19515fdea053SToby Isaac     ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr);
19521a834cf9SToby 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);
19531a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr);
19541a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr);
19551a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr);
19561a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr);
19571a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr);
19585fdea053SToby Isaac   }
19595fdea053SToby Isaac   PetscFunctionReturn(0);
19605fdea053SToby Isaac }
19615fdea053SToby Isaac 
19625fdea053SToby Isaac /*@C
19635fdea053SToby Isaac   PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
19645fdea053SToby Isaac 
19655b5e7992SMatthew G. Knepley   Logically collective on sym
19665fdea053SToby Isaac 
19675fdea053SToby Isaac   InputParameters:
19685b5e7992SMatthew G. Knepley + sym       - the section symmetries
19695fdea053SToby Isaac . stratum   - the stratum value in the label that we are assigning symmetries for
19705fdea053SToby Isaac . size      - the number of dofs for points in the stratum of the label
19715fdea053SToby Isaac . minOrient - the smallest orientation for a point in this stratum
19725fdea053SToby Isaac . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
19735fdea053SToby Isaac . mode      - how sym should copy the perms and rots arrays
19745fdea053SToby Isaac . perms     - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation.  A NULL permutation is the identity
1975a2b725a8SWilliam 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
19765fdea053SToby Isaac 
19775fdea053SToby Isaac   Level: developer
19785fdea053SToby Isaac 
19795fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()
19805fdea053SToby Isaac @*/
19815fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
19825fdea053SToby Isaac {
19835fdea053SToby Isaac   PetscSectionSym_Label *sl;
1984d67d17b1SMatthew G. Knepley   const char            *name;
1985d67d17b1SMatthew G. Knepley   PetscInt               i, j, k;
19865fdea053SToby Isaac   PetscErrorCode         ierr;
19875fdea053SToby Isaac 
19885fdea053SToby Isaac   PetscFunctionBegin;
19895fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
19905fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
19915fdea053SToby Isaac   if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet");
19925fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
19935fdea053SToby Isaac     PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
19945fdea053SToby Isaac 
19955fdea053SToby Isaac     if (stratum == value) break;
19965fdea053SToby Isaac   }
1997d67d17b1SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1998d67d17b1SMatthew G. Knepley   if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,name);
19995fdea053SToby Isaac   sl->sizes[i] = size;
20005fdea053SToby Isaac   sl->modes[i] = mode;
20015fdea053SToby Isaac   sl->minMaxOrients[i][0] = minOrient;
20025fdea053SToby Isaac   sl->minMaxOrients[i][1] = maxOrient;
20035fdea053SToby Isaac   if (mode == PETSC_COPY_VALUES) {
20045fdea053SToby Isaac     if (perms) {
20055fdea053SToby Isaac       PetscInt    **ownPerms;
20065fdea053SToby Isaac 
20075fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr);
20085fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
20095fdea053SToby Isaac         if (perms[j]) {
20105fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr);
20115fdea053SToby Isaac           for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];}
20125fdea053SToby Isaac         }
20135fdea053SToby Isaac       }
20145fdea053SToby Isaac       sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient];
20155fdea053SToby Isaac     }
20165fdea053SToby Isaac     if (rots) {
20175fdea053SToby Isaac       PetscScalar **ownRots;
20185fdea053SToby Isaac 
20195fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr);
20205fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
20215fdea053SToby Isaac         if (rots[j]) {
20225fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr);
20235fdea053SToby Isaac           for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];}
20245fdea053SToby Isaac         }
20255fdea053SToby Isaac       }
20265fdea053SToby Isaac       sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient];
20275fdea053SToby Isaac     }
20285fdea053SToby Isaac   } else {
20295fdea053SToby Isaac     sl->perms[i] = perms ? &perms[-minOrient] : NULL;
20305fdea053SToby Isaac     sl->rots[i]  = rots ? &rots[-minOrient] : NULL;
20315fdea053SToby Isaac   }
20325fdea053SToby Isaac   PetscFunctionReturn(0);
20335fdea053SToby Isaac }
20345fdea053SToby Isaac 
20355fdea053SToby Isaac static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
20365fdea053SToby Isaac {
20375fdea053SToby Isaac   PetscInt              i, j, numStrata;
20385fdea053SToby Isaac   PetscSectionSym_Label *sl;
20395fdea053SToby Isaac   DMLabel               label;
20405fdea053SToby Isaac   PetscErrorCode        ierr;
20415fdea053SToby Isaac 
20425fdea053SToby Isaac   PetscFunctionBegin;
20435fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
20445fdea053SToby Isaac   numStrata = sl->numStrata;
20455fdea053SToby Isaac   label     = sl->label;
20465fdea053SToby Isaac   for (i = 0; i < numPoints; i++) {
20475fdea053SToby Isaac     PetscInt point = points[2*i];
20485fdea053SToby Isaac     PetscInt ornt  = points[2*i+1];
20495fdea053SToby Isaac 
20505fdea053SToby Isaac     for (j = 0; j < numStrata; j++) {
20515fdea053SToby Isaac       if (label->validIS[j]) {
20525fdea053SToby Isaac         PetscInt k;
20535fdea053SToby Isaac 
20545fdea053SToby Isaac         ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr);
20555fdea053SToby Isaac         if (k >= 0) break;
20565fdea053SToby Isaac       } else {
20575fdea053SToby Isaac         PetscBool has;
20585fdea053SToby Isaac 
2059b9907514SLisandro Dalcin         ierr = PetscHSetIHas(label->ht[j], point, &has);CHKERRQ(ierr);
20605fdea053SToby Isaac         if (has) break;
20615fdea053SToby Isaac       }
20625fdea053SToby Isaac     }
20635fdea053SToby 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);
20645fdea053SToby Isaac     if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;}
20655fdea053SToby Isaac     if (rots) {rots[i]  = sl->rots[j] ? sl->rots[j][ornt] : NULL;}
20665fdea053SToby Isaac   }
20675fdea053SToby Isaac   PetscFunctionReturn(0);
20685fdea053SToby Isaac }
20695fdea053SToby Isaac 
20705fdea053SToby Isaac PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
20715fdea053SToby Isaac {
20725fdea053SToby Isaac   PetscSectionSym_Label *sl;
20735fdea053SToby Isaac   PetscErrorCode        ierr;
20745fdea053SToby Isaac 
20755fdea053SToby Isaac   PetscFunctionBegin;
20765fdea053SToby Isaac   ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr);
20775fdea053SToby Isaac   sym->ops->getpoints = PetscSectionSymGetPoints_Label;
20785fdea053SToby Isaac   sym->ops->view      = PetscSectionSymView_Label;
20795fdea053SToby Isaac   sym->ops->destroy   = PetscSectionSymDestroy_Label;
20805fdea053SToby Isaac   sym->data           = (void *) sl;
20815fdea053SToby Isaac   PetscFunctionReturn(0);
20825fdea053SToby Isaac }
20835fdea053SToby Isaac 
20845fdea053SToby Isaac /*@
20855fdea053SToby Isaac   PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
20865fdea053SToby Isaac 
20875fdea053SToby Isaac   Collective
20885fdea053SToby Isaac 
20895fdea053SToby Isaac   Input Parameters:
20905fdea053SToby Isaac + comm - the MPI communicator for the new symmetry
20915fdea053SToby Isaac - label - the label defining the strata
20925fdea053SToby Isaac 
20935fdea053SToby Isaac   Output Parameters:
20945fdea053SToby Isaac . sym - the section symmetries
20955fdea053SToby Isaac 
20965fdea053SToby Isaac   Level: developer
20975fdea053SToby Isaac 
20985fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()
20995fdea053SToby Isaac @*/
21005fdea053SToby Isaac PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
21015fdea053SToby Isaac {
21025fdea053SToby Isaac   PetscErrorCode        ierr;
21035fdea053SToby Isaac 
21045fdea053SToby Isaac   PetscFunctionBegin;
21055fdea053SToby Isaac   ierr = DMInitializePackage();CHKERRQ(ierr);
21065fdea053SToby Isaac   ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr);
21075fdea053SToby Isaac   ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr);
21085fdea053SToby Isaac   ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr);
21095fdea053SToby Isaac   PetscFunctionReturn(0);
21105fdea053SToby Isaac }
2111