xref: /petsc/src/dm/label/dmlabel.c (revision f14fe40db710f1db1d8b4863e660105d93d93376)
15fdea053SToby Isaac #include <petscdm.h>
2c58f1c22SToby Isaac #include <petsc/private/dmlabelimpl.h>   /*I      "petscdmlabel.h"   I*/
3c58f1c22SToby Isaac #include <petsc/private/isimpl.h>        /*I      "petscis.h"        I*/
4c58f1c22SToby Isaac #include <petscsf.h>
5c58f1c22SToby Isaac 
6c58f1c22SToby Isaac /*@C
7c58f1c22SToby Isaac   DMLabelCreate - Create a DMLabel object, which is a multimap
8c58f1c22SToby Isaac 
9d67d17b1SMatthew G. Knepley   Input parameters:
10d67d17b1SMatthew G. Knepley + comm - The communicator, usually PETSC_COMM_SELF
11d67d17b1SMatthew G. Knepley - name - The label name
12c58f1c22SToby Isaac 
13c58f1c22SToby Isaac   Output parameter:
14c58f1c22SToby Isaac . label - The DMLabel
15c58f1c22SToby Isaac 
16c58f1c22SToby Isaac   Level: beginner
17c58f1c22SToby Isaac 
18c58f1c22SToby Isaac .seealso: DMLabelDestroy()
19c58f1c22SToby Isaac @*/
20d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelCreate(MPI_Comm comm, const char name[], DMLabel *label)
21c58f1c22SToby Isaac {
22c58f1c22SToby Isaac   PetscErrorCode ierr;
23c58f1c22SToby Isaac 
24c58f1c22SToby Isaac   PetscFunctionBegin;
25d67d17b1SMatthew G. Knepley   PetscValidPointer(label,2);
26d67d17b1SMatthew G. Knepley   ierr = DMInitializePackage();CHKERRQ(ierr);
27c58f1c22SToby Isaac 
28d67d17b1SMatthew G. Knepley   ierr = PetscHeaderCreate(*label,DMLABEL_CLASSID,"DMLabel","DMLabel","DM",comm,DMLabelDestroy,DMLabelView);CHKERRQ(ierr);
29d67d17b1SMatthew G. Knepley 
30c58f1c22SToby Isaac   (*label)->numStrata      = 0;
315aa44df4SToby Isaac   (*label)->defaultValue   = -1;
32c58f1c22SToby Isaac   (*label)->stratumValues  = NULL;
33ad8374ffSToby Isaac   (*label)->validIS        = NULL;
34c58f1c22SToby Isaac   (*label)->stratumSizes   = NULL;
35c58f1c22SToby Isaac   (*label)->points         = NULL;
36c58f1c22SToby Isaac   (*label)->ht             = NULL;
37c58f1c22SToby Isaac   (*label)->pStart         = -1;
38c58f1c22SToby Isaac   (*label)->pEnd           = -1;
39c58f1c22SToby Isaac   (*label)->bt             = NULL;
40b9907514SLisandro Dalcin   ierr = PetscHMapICreate(&(*label)->hmap);CHKERRQ(ierr);
41d67d17b1SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *label, name);CHKERRQ(ierr);
42c58f1c22SToby Isaac   PetscFunctionReturn(0);
43c58f1c22SToby Isaac }
44c58f1c22SToby Isaac 
45c58f1c22SToby Isaac /*
46c58f1c22SToby Isaac   DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format
47c58f1c22SToby Isaac 
48c58f1c22SToby Isaac   Input parameter:
49c58f1c22SToby Isaac + label - The DMLabel
50c58f1c22SToby Isaac - v - The stratum value
51c58f1c22SToby Isaac 
52c58f1c22SToby Isaac   Output parameter:
53c58f1c22SToby Isaac . label - The DMLabel with stratum in sorted list format
54c58f1c22SToby Isaac 
55c58f1c22SToby Isaac   Level: developer
56c58f1c22SToby Isaac 
57c58f1c22SToby Isaac .seealso: DMLabelCreate()
58c58f1c22SToby Isaac */
59c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v)
60c58f1c22SToby Isaac {
61b9907514SLisandro Dalcin   PetscInt       off = 0, *pointArray, p;
62c58f1c22SToby Isaac   PetscErrorCode ierr;
63c58f1c22SToby Isaac 
64b9907514SLisandro Dalcin   if (PetscLikely(v >= 0 && v < label->numStrata) && label->validIS[v]) return 0;
65c58f1c22SToby Isaac   PetscFunctionBegin;
660c3c4a36SLisandro Dalcin   if (v < 0 || v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeValid_Private\n", v);
67e8f14785SLisandro Dalcin   ierr = PetscHSetIGetSize(label->ht[v], &label->stratumSizes[v]);CHKERRQ(ierr);
68ad8374ffSToby Isaac   ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr);
69e8f14785SLisandro Dalcin   ierr = PetscHSetIGetElems(label->ht[v], &off, pointArray);CHKERRQ(ierr);
70b9907514SLisandro Dalcin   ierr = PetscHSetIClear(label->ht[v]);CHKERRQ(ierr);
71ad8374ffSToby Isaac   ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr);
72c58f1c22SToby Isaac   if (label->bt) {
73c58f1c22SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
74ad8374ffSToby Isaac       const PetscInt point = pointArray[p];
75c58f1c22SToby 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);
76c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
77c58f1c22SToby Isaac     }
78c58f1c22SToby Isaac   }
79ad8374ffSToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF,label->stratumSizes[v],pointArray,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
80ad8374ffSToby Isaac   ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
81ad8374ffSToby Isaac   label->validIS[v] = PETSC_TRUE;
82d67d17b1SMatthew G. Knepley   ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
83c58f1c22SToby Isaac   PetscFunctionReturn(0);
84c58f1c22SToby Isaac }
85c58f1c22SToby Isaac 
86c58f1c22SToby Isaac /*
87c58f1c22SToby Isaac   DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format
88c58f1c22SToby Isaac 
89c58f1c22SToby Isaac   Input parameter:
90c58f1c22SToby Isaac . label - The DMLabel
91c58f1c22SToby Isaac 
92c58f1c22SToby Isaac   Output parameter:
93c58f1c22SToby Isaac . label - The DMLabel with all strata in sorted list format
94c58f1c22SToby Isaac 
95c58f1c22SToby Isaac   Level: developer
96c58f1c22SToby Isaac 
97c58f1c22SToby Isaac .seealso: DMLabelCreate()
98c58f1c22SToby Isaac */
99c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label)
100c58f1c22SToby Isaac {
101c58f1c22SToby Isaac   PetscInt       v;
102c58f1c22SToby Isaac   PetscErrorCode ierr;
103c58f1c22SToby Isaac 
104c58f1c22SToby Isaac   PetscFunctionBegin;
105c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; v++){
106c58f1c22SToby Isaac     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
107c58f1c22SToby Isaac   }
108c58f1c22SToby Isaac   PetscFunctionReturn(0);
109c58f1c22SToby Isaac }
110c58f1c22SToby Isaac 
111c58f1c22SToby Isaac /*
112c58f1c22SToby Isaac   DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format
113c58f1c22SToby Isaac 
114c58f1c22SToby Isaac   Input parameter:
115c58f1c22SToby Isaac + label - The DMLabel
116c58f1c22SToby Isaac - v - The stratum value
117c58f1c22SToby Isaac 
118c58f1c22SToby Isaac   Output parameter:
119c58f1c22SToby Isaac . label - The DMLabel with stratum in hash format
120c58f1c22SToby Isaac 
121c58f1c22SToby Isaac   Level: developer
122c58f1c22SToby Isaac 
123c58f1c22SToby Isaac .seealso: DMLabelCreate()
124c58f1c22SToby Isaac */
125c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v)
126c58f1c22SToby Isaac {
127c58f1c22SToby Isaac   PetscInt       p;
128ad8374ffSToby Isaac   const PetscInt *points;
129c58f1c22SToby Isaac   PetscErrorCode ierr;
130c58f1c22SToby Isaac 
131b9907514SLisandro Dalcin   if (PetscLikely(v >= 0 && v < label->numStrata) && !label->validIS[v]) return 0;
132c58f1c22SToby Isaac   PetscFunctionBegin;
1330c3c4a36SLisandro 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);
134ad8374ffSToby Isaac   if (label->points[v]) {
135ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
136e8f14785SLisandro Dalcin     for (p = 0; p < label->stratumSizes[v]; ++p) {
137e8f14785SLisandro Dalcin       ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);
138e8f14785SLisandro Dalcin     }
139ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
140ad8374ffSToby Isaac     ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
141ad8374ffSToby Isaac   }
142ad8374ffSToby Isaac   label->validIS[v] = PETSC_FALSE;
143c58f1c22SToby Isaac   PetscFunctionReturn(0);
144c58f1c22SToby Isaac }
145c58f1c22SToby Isaac 
146b9907514SLisandro Dalcin #if !defined(DMLABEL_LOOKUP_THRESHOLD)
147b9907514SLisandro Dalcin #define DMLABEL_LOOKUP_THRESHOLD 16
148b9907514SLisandro Dalcin #endif
149b9907514SLisandro Dalcin 
1500c3c4a36SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index)
1510c3c4a36SLisandro Dalcin {
1520c3c4a36SLisandro Dalcin   PetscInt       v;
153b9907514SLisandro Dalcin   PetscErrorCode ierr;
1540e79e033SBarry Smith 
1550c3c4a36SLisandro Dalcin   PetscFunctionBegin;
1560e79e033SBarry Smith   *index = -1;
157b9907514SLisandro Dalcin   if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) {
158b9907514SLisandro Dalcin     for (v = 0; v < label->numStrata; ++v)
159b9907514SLisandro Dalcin       if (label->stratumValues[v] == value) {*index = v; break;}
160b9907514SLisandro Dalcin   } else {
161b9907514SLisandro Dalcin     ierr = PetscHMapIGet(label->hmap, value, index);CHKERRQ(ierr);
162b9907514SLisandro Dalcin #if defined(PETSC_USE_DEBUG)
163b9907514SLisandro Dalcin     v = *index;
164b9907514SLisandro Dalcin     if (v >= 0 && v >= label->numStrata) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent lookup using strata hash map");
165b9907514SLisandro Dalcin     if (v >= 0 && label->stratumValues[v] != value) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent lookup using strata hash map");
166b9907514SLisandro Dalcin #endif
1670e79e033SBarry Smith   }
1680c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
1690c3c4a36SLisandro Dalcin }
1700c3c4a36SLisandro Dalcin 
171b9907514SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index)
172c58f1c22SToby Isaac {
173b9907514SLisandro Dalcin   PetscInt       v;
174b9907514SLisandro Dalcin   PetscInt      *tmpV;
175b9907514SLisandro Dalcin   PetscInt      *tmpS;
176b9907514SLisandro Dalcin   PetscHSetI    *tmpH, ht;
177b9907514SLisandro Dalcin   IS            *tmpP, is;
178c58f1c22SToby Isaac   PetscBool     *tmpB;
179b9907514SLisandro Dalcin   PetscHMapI     hmap = label->hmap;
180c58f1c22SToby Isaac   PetscErrorCode ierr;
181c58f1c22SToby Isaac 
182c58f1c22SToby Isaac   PetscFunctionBegin;
183b9907514SLisandro Dalcin   v    = label->numStrata;
184b9907514SLisandro Dalcin   tmpV = label->stratumValues;
185b9907514SLisandro Dalcin   tmpS = label->stratumSizes;
186b9907514SLisandro Dalcin   tmpH = label->ht;
187b9907514SLisandro Dalcin   tmpP = label->points;
188b9907514SLisandro Dalcin   tmpB = label->validIS;
1898e1f8cf0SLisandro Dalcin   { /* TODO: PetscRealloc() is broken, use malloc+memcpy+free  */
1908e1f8cf0SLisandro Dalcin     PetscInt   *oldV = tmpV;
1918e1f8cf0SLisandro Dalcin     PetscInt   *oldS = tmpS;
1928e1f8cf0SLisandro Dalcin     PetscHSetI *oldH = tmpH;
1938e1f8cf0SLisandro Dalcin     IS         *oldP = tmpP;
1948e1f8cf0SLisandro Dalcin     PetscBool  *oldB = tmpB;
1958e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpV), &tmpV);CHKERRQ(ierr);
1968e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpS), &tmpS);CHKERRQ(ierr);
1978e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpH), &tmpH);CHKERRQ(ierr);
1988e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpP), &tmpP);CHKERRQ(ierr);
1998e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpB), &tmpB);CHKERRQ(ierr);
2008e1f8cf0SLisandro Dalcin     ierr = PetscMemcpy(tmpV, oldV, v*sizeof(*tmpV));CHKERRQ(ierr);
2018e1f8cf0SLisandro Dalcin     ierr = PetscMemcpy(tmpS, oldS, v*sizeof(*tmpS));CHKERRQ(ierr);
2028e1f8cf0SLisandro Dalcin     ierr = PetscMemcpy(tmpH, oldH, v*sizeof(*tmpH));CHKERRQ(ierr);
2038e1f8cf0SLisandro Dalcin     ierr = PetscMemcpy(tmpP, oldP, v*sizeof(*tmpP));CHKERRQ(ierr);
2048e1f8cf0SLisandro Dalcin     ierr = PetscMemcpy(tmpB, oldB, v*sizeof(*tmpB));CHKERRQ(ierr);
2058e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldV);CHKERRQ(ierr);
2068e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldS);CHKERRQ(ierr);
2078e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldH);CHKERRQ(ierr);
2088e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldP);CHKERRQ(ierr);
2098e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldB);CHKERRQ(ierr);
2108e1f8cf0SLisandro Dalcin   }
211b9907514SLisandro Dalcin   label->numStrata     = v+1;
212c58f1c22SToby Isaac   label->stratumValues = tmpV;
213c58f1c22SToby Isaac   label->stratumSizes  = tmpS;
214c58f1c22SToby Isaac   label->ht            = tmpH;
215c58f1c22SToby Isaac   label->points        = tmpP;
216ad8374ffSToby Isaac   label->validIS       = tmpB;
217b9907514SLisandro Dalcin   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
218b9907514SLisandro Dalcin   ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);CHKERRQ(ierr);
219b9907514SLisandro Dalcin   ierr = PetscHMapISet(hmap, value, v);CHKERRQ(ierr);
220b9907514SLisandro Dalcin   tmpV[v] = value;
221b9907514SLisandro Dalcin   tmpS[v] = 0;
222b9907514SLisandro Dalcin   tmpH[v] = ht;
223b9907514SLisandro Dalcin   tmpP[v] = is;
224b9907514SLisandro Dalcin   tmpB[v] = PETSC_TRUE;
2250c3c4a36SLisandro Dalcin   *index = v;
2260c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
2270c3c4a36SLisandro Dalcin }
2280c3c4a36SLisandro Dalcin 
229b9907514SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupAddStratum(DMLabel label, PetscInt value, PetscInt *index)
230b9907514SLisandro Dalcin {
231b9907514SLisandro Dalcin   PetscErrorCode ierr;
232b9907514SLisandro Dalcin   PetscFunctionBegin;
233b9907514SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, index);CHKERRQ(ierr);
234b9907514SLisandro Dalcin   if (*index < 0) {ierr = DMLabelNewStratum(label, value, index);CHKERRQ(ierr);}
235b9907514SLisandro Dalcin   PetscFunctionReturn(0);
236b9907514SLisandro Dalcin }
237b9907514SLisandro Dalcin 
238b9907514SLisandro Dalcin /*@
239b9907514SLisandro Dalcin   DMLabelAddStratum - Adds a new stratum value in a DMLabel
240b9907514SLisandro Dalcin 
241b9907514SLisandro Dalcin   Input Parameter:
242b9907514SLisandro Dalcin + label - The DMLabel
243b9907514SLisandro Dalcin - value - The stratum value
244b9907514SLisandro Dalcin 
245b9907514SLisandro Dalcin   Level: beginner
246b9907514SLisandro Dalcin 
247b9907514SLisandro Dalcin .seealso:  DMLabelCreate(), DMLabelDestroy()
248b9907514SLisandro Dalcin @*/
2490c3c4a36SLisandro Dalcin PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
2500c3c4a36SLisandro Dalcin {
2510c3c4a36SLisandro Dalcin   PetscInt       v;
2520c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
2530c3c4a36SLisandro Dalcin 
2540c3c4a36SLisandro Dalcin   PetscFunctionBegin;
255d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
256b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
257b9907514SLisandro Dalcin   PetscFunctionReturn(0);
258b9907514SLisandro Dalcin }
259b9907514SLisandro Dalcin 
260b9907514SLisandro Dalcin /*@
261b9907514SLisandro Dalcin   DMLabelAddStrata - Adds new stratum values in a DMLabel
262b9907514SLisandro Dalcin 
263b9907514SLisandro Dalcin   Input Parameter:
264b9907514SLisandro Dalcin + label - The DMLabel
265b9907514SLisandro Dalcin . numStrata - The number of stratum values
266b9907514SLisandro Dalcin - stratumValues - The stratum values
267b9907514SLisandro Dalcin 
268b9907514SLisandro Dalcin   Level: beginner
269b9907514SLisandro Dalcin 
270b9907514SLisandro Dalcin .seealso:  DMLabelCreate(), DMLabelDestroy()
271b9907514SLisandro Dalcin @*/
272b9907514SLisandro Dalcin PetscErrorCode DMLabelAddStrata(DMLabel label, PetscInt numStrata, const PetscInt stratumValues[])
273b9907514SLisandro Dalcin {
274b9907514SLisandro Dalcin   PetscInt       *values, v;
275b9907514SLisandro Dalcin   PetscErrorCode ierr;
276b9907514SLisandro Dalcin 
277b9907514SLisandro Dalcin   PetscFunctionBegin;
278b9907514SLisandro Dalcin   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
279b9907514SLisandro Dalcin   if (numStrata) PetscValidIntPointer(stratumValues, 3);
280b9907514SLisandro Dalcin   ierr = PetscMalloc1(numStrata, &values);CHKERRQ(ierr);
281b9907514SLisandro Dalcin   ierr = PetscMemcpy(values, stratumValues, numStrata*sizeof(PetscInt));CHKERRQ(ierr);
282b9907514SLisandro Dalcin   ierr = PetscSortRemoveDupsInt(&numStrata, values);CHKERRQ(ierr);
283b9907514SLisandro Dalcin   if (!label->numStrata) { /* Fast preallocation */
284b9907514SLisandro Dalcin     PetscInt   *tmpV;
285b9907514SLisandro Dalcin     PetscInt   *tmpS;
286b9907514SLisandro Dalcin     PetscHSetI *tmpH, ht;
287b9907514SLisandro Dalcin     IS         *tmpP, is;
288b9907514SLisandro Dalcin     PetscBool  *tmpB;
289b9907514SLisandro Dalcin     PetscHMapI  hmap = label->hmap;
290b9907514SLisandro Dalcin 
291b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpV);CHKERRQ(ierr);
292b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpS);CHKERRQ(ierr);
293b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpH);CHKERRQ(ierr);
294b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpP);CHKERRQ(ierr);
295b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpB);CHKERRQ(ierr);
296b9907514SLisandro Dalcin     label->numStrata     = numStrata;
297b9907514SLisandro Dalcin     label->stratumValues = tmpV;
298b9907514SLisandro Dalcin     label->stratumSizes  = tmpS;
299b9907514SLisandro Dalcin     label->ht            = tmpH;
300b9907514SLisandro Dalcin     label->points        = tmpP;
301b9907514SLisandro Dalcin     label->validIS       = tmpB;
302b9907514SLisandro Dalcin     for (v = 0; v < numStrata; ++v) {
303b9907514SLisandro Dalcin       ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
304b9907514SLisandro Dalcin       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);CHKERRQ(ierr);
305b9907514SLisandro Dalcin       ierr = PetscHMapISet(hmap, values[v], v);CHKERRQ(ierr);
306b9907514SLisandro Dalcin       tmpV[v] = values[v];
307b9907514SLisandro Dalcin       tmpS[v] = 0;
308b9907514SLisandro Dalcin       tmpH[v] = ht;
309b9907514SLisandro Dalcin       tmpP[v] = is;
310b9907514SLisandro Dalcin       tmpB[v] = PETSC_TRUE;
311b9907514SLisandro Dalcin     }
312b9907514SLisandro Dalcin   } else {
313b9907514SLisandro Dalcin     for (v = 0; v < numStrata; ++v) {
314b9907514SLisandro Dalcin       ierr = DMLabelAddStratum(label, values[v]);CHKERRQ(ierr);
315b9907514SLisandro Dalcin     }
316b9907514SLisandro Dalcin   }
317b9907514SLisandro Dalcin   ierr = PetscFree(values);CHKERRQ(ierr);
318b9907514SLisandro Dalcin   PetscFunctionReturn(0);
319b9907514SLisandro Dalcin }
320b9907514SLisandro Dalcin 
321b9907514SLisandro Dalcin /*@
322b9907514SLisandro Dalcin   DMLabelAddStrataIS - Adds new stratum values in a DMLabel
323b9907514SLisandro Dalcin 
324b9907514SLisandro Dalcin   Input Parameter:
325b9907514SLisandro Dalcin + label - The DMLabel
326b9907514SLisandro Dalcin - valueIS - Index set with stratum values
327b9907514SLisandro Dalcin 
328b9907514SLisandro Dalcin   Level: beginner
329b9907514SLisandro Dalcin 
330b9907514SLisandro Dalcin .seealso:  DMLabelCreate(), DMLabelDestroy()
331b9907514SLisandro Dalcin @*/
332b9907514SLisandro Dalcin PetscErrorCode DMLabelAddStrataIS(DMLabel label, IS valueIS)
333b9907514SLisandro Dalcin {
334b9907514SLisandro Dalcin   PetscInt       numStrata;
335b9907514SLisandro Dalcin   const PetscInt *stratumValues;
336b9907514SLisandro Dalcin   PetscErrorCode ierr;
337b9907514SLisandro Dalcin 
338b9907514SLisandro Dalcin   PetscFunctionBegin;
339b9907514SLisandro Dalcin   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
340b9907514SLisandro Dalcin   PetscValidHeaderSpecific(valueIS, IS_CLASSID, 2);
341b9907514SLisandro Dalcin   ierr = ISGetLocalSize(valueIS, &numStrata);CHKERRQ(ierr);
342b9907514SLisandro Dalcin   ierr = ISGetIndices(valueIS, &stratumValues);CHKERRQ(ierr);
343b9907514SLisandro Dalcin   ierr = DMLabelAddStrata(label, numStrata, stratumValues);CHKERRQ(ierr);
344c58f1c22SToby Isaac   PetscFunctionReturn(0);
345c58f1c22SToby Isaac }
346c58f1c22SToby Isaac 
347c58f1c22SToby Isaac static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
348c58f1c22SToby Isaac {
349c58f1c22SToby Isaac   PetscInt       v;
350c58f1c22SToby Isaac   PetscMPIInt    rank;
351c58f1c22SToby Isaac   PetscErrorCode ierr;
352c58f1c22SToby Isaac 
353c58f1c22SToby Isaac   PetscFunctionBegin;
354c58f1c22SToby Isaac   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr);
355c58f1c22SToby Isaac   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
356c58f1c22SToby Isaac   if (label) {
357d67d17b1SMatthew G. Knepley     const char *name;
358d67d17b1SMatthew G. Knepley 
359d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
360d67d17b1SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", name);CHKERRQ(ierr);
361c58f1c22SToby Isaac     if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);}
362c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
363c58f1c22SToby Isaac       const PetscInt value = label->stratumValues[v];
364ad8374ffSToby Isaac       const PetscInt *points;
365c58f1c22SToby Isaac       PetscInt       p;
366c58f1c22SToby Isaac 
367ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
368c58f1c22SToby Isaac       for (p = 0; p < label->stratumSizes[v]; ++p) {
369ad8374ffSToby Isaac         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr);
370c58f1c22SToby Isaac       }
371ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
372c58f1c22SToby Isaac     }
373c58f1c22SToby Isaac   }
374c58f1c22SToby Isaac   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
375c58f1c22SToby Isaac   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
376c58f1c22SToby Isaac   PetscFunctionReturn(0);
377c58f1c22SToby Isaac }
378c58f1c22SToby Isaac 
379c58f1c22SToby Isaac /*@C
380c58f1c22SToby Isaac   DMLabelView - View the label
381c58f1c22SToby Isaac 
382c58f1c22SToby Isaac   Input Parameters:
383c58f1c22SToby Isaac + label - The DMLabel
384c58f1c22SToby Isaac - viewer - The PetscViewer
385c58f1c22SToby Isaac 
386c58f1c22SToby Isaac   Level: intermediate
387c58f1c22SToby Isaac 
388c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy()
389c58f1c22SToby Isaac @*/
390c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
391c58f1c22SToby Isaac {
392c58f1c22SToby Isaac   PetscBool      iascii;
393c58f1c22SToby Isaac   PetscErrorCode ierr;
394c58f1c22SToby Isaac 
395c58f1c22SToby Isaac   PetscFunctionBegin;
396d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
397b9907514SLisandro Dalcin   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)label), &viewer);CHKERRQ(ierr);}
398c58f1c22SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
399c58f1c22SToby Isaac   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
400c58f1c22SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
401c58f1c22SToby Isaac   if (iascii) {
402c58f1c22SToby Isaac     ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr);
403c58f1c22SToby Isaac   }
404c58f1c22SToby Isaac   PetscFunctionReturn(0);
405c58f1c22SToby Isaac }
406c58f1c22SToby Isaac 
40784f0b6dfSMatthew G. Knepley /*@
408d67d17b1SMatthew G. Knepley   DMLabelReset - Destroys internal data structures in a DMLabel
409d67d17b1SMatthew G. Knepley 
410d67d17b1SMatthew G. Knepley   Input Parameter:
411d67d17b1SMatthew G. Knepley . label - The DMLabel
412d67d17b1SMatthew G. Knepley 
413d67d17b1SMatthew G. Knepley   Level: beginner
414d67d17b1SMatthew G. Knepley 
415d67d17b1SMatthew G. Knepley .seealso: DMLabelDestroy(), DMLabelCreate()
416d67d17b1SMatthew G. Knepley @*/
417d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelReset(DMLabel label)
418d67d17b1SMatthew G. Knepley {
419d67d17b1SMatthew G. Knepley   PetscInt       v;
420d67d17b1SMatthew G. Knepley   PetscErrorCode ierr;
421d67d17b1SMatthew G. Knepley 
422d67d17b1SMatthew G. Knepley   PetscFunctionBegin;
423d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
424d67d17b1SMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
425d67d17b1SMatthew G. Knepley     ierr = PetscHSetIDestroy(&label->ht[v]);CHKERRQ(ierr);
426d67d17b1SMatthew G. Knepley     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
427d67d17b1SMatthew G. Knepley   }
428b9907514SLisandro Dalcin   label->numStrata = 0;
429b9907514SLisandro Dalcin   ierr = PetscFree(label->stratumValues);CHKERRQ(ierr);
430b9907514SLisandro Dalcin   ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr);
431d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->ht);CHKERRQ(ierr);
432d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->points);CHKERRQ(ierr);
433d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->validIS);CHKERRQ(ierr);
434b9907514SLisandro Dalcin   ierr = PetscHMapIReset(label->hmap);CHKERRQ(ierr);
435b9907514SLisandro Dalcin   label->pStart = -1;
436b9907514SLisandro Dalcin   label->pEnd   = -1;
437d67d17b1SMatthew G. Knepley   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
438d67d17b1SMatthew G. Knepley   PetscFunctionReturn(0);
439d67d17b1SMatthew G. Knepley }
440d67d17b1SMatthew G. Knepley 
441d67d17b1SMatthew G. Knepley /*@
44284f0b6dfSMatthew G. Knepley   DMLabelDestroy - Destroys a DMLabel
44384f0b6dfSMatthew G. Knepley 
44484f0b6dfSMatthew G. Knepley   Input Parameter:
44584f0b6dfSMatthew G. Knepley . label - The DMLabel
44684f0b6dfSMatthew G. Knepley 
44784f0b6dfSMatthew G. Knepley   Level: beginner
44884f0b6dfSMatthew G. Knepley 
449d67d17b1SMatthew G. Knepley .seealso: DMLabelReset(), DMLabelCreate()
45084f0b6dfSMatthew G. Knepley @*/
451c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label)
452c58f1c22SToby Isaac {
453c58f1c22SToby Isaac   PetscErrorCode ierr;
454c58f1c22SToby Isaac 
455c58f1c22SToby Isaac   PetscFunctionBegin;
456d67d17b1SMatthew G. Knepley   if (!*label) PetscFunctionReturn(0);
457d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific((*label),DMLABEL_CLASSID,1);
458b9907514SLisandro Dalcin   if (--((PetscObject)(*label))->refct > 0) {*label = NULL; PetscFunctionReturn(0);}
459d67d17b1SMatthew G. Knepley   ierr = DMLabelReset(*label);CHKERRQ(ierr);
460b9907514SLisandro Dalcin   ierr = PetscHMapIDestroy(&(*label)->hmap);CHKERRQ(ierr);
461d67d17b1SMatthew G. Knepley   ierr = PetscHeaderDestroy(label);CHKERRQ(ierr);
462c58f1c22SToby Isaac   PetscFunctionReturn(0);
463c58f1c22SToby Isaac }
464c58f1c22SToby Isaac 
46584f0b6dfSMatthew G. Knepley /*@
46684f0b6dfSMatthew G. Knepley   DMLabelDuplicate - Duplicates a DMLabel
46784f0b6dfSMatthew G. Knepley 
46884f0b6dfSMatthew G. Knepley   Input Parameter:
46984f0b6dfSMatthew G. Knepley . label - The DMLabel
47084f0b6dfSMatthew G. Knepley 
47184f0b6dfSMatthew G. Knepley   Output Parameter:
47284f0b6dfSMatthew G. Knepley . labelnew - location to put new vector
47384f0b6dfSMatthew G. Knepley 
47484f0b6dfSMatthew G. Knepley   Level: intermediate
47584f0b6dfSMatthew G. Knepley 
47684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelDestroy()
47784f0b6dfSMatthew G. Knepley @*/
478c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
479c58f1c22SToby Isaac {
480d67d17b1SMatthew G. Knepley   const char    *name;
481ad8374ffSToby Isaac   PetscInt       v;
482c58f1c22SToby Isaac   PetscErrorCode ierr;
483c58f1c22SToby Isaac 
484c58f1c22SToby Isaac   PetscFunctionBegin;
485d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
486c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
487d67d17b1SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
488d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PetscObjectComm((PetscObject) label), name, labelnew);CHKERRQ(ierr);
489c58f1c22SToby Isaac 
490c58f1c22SToby Isaac   (*labelnew)->numStrata    = label->numStrata;
4915aa44df4SToby Isaac   (*labelnew)->defaultValue = label->defaultValue;
492c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr);
493c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr);
494c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr);
495c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr);
496ad8374ffSToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr);
497c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
498e8f14785SLisandro Dalcin     ierr = PetscHSetICreate(&(*labelnew)->ht[v]);CHKERRQ(ierr);
499c58f1c22SToby Isaac     (*labelnew)->stratumValues[v]  = label->stratumValues[v];
500c58f1c22SToby Isaac     (*labelnew)->stratumSizes[v]   = label->stratumSizes[v];
501ad8374ffSToby Isaac     ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr);
502ad8374ffSToby Isaac     (*labelnew)->points[v]         = label->points[v];
503b9907514SLisandro Dalcin     (*labelnew)->validIS[v]        = PETSC_TRUE;
504c58f1c22SToby Isaac   }
505*f14fe40dSLisandro Dalcin   ierr = PetscHMapIDestroy(&(*labelnew)->hmap);CHKERRQ(ierr);
506b9907514SLisandro Dalcin   ierr = PetscHMapIDuplicate(label->hmap,&(*labelnew)->hmap);CHKERRQ(ierr);
507c58f1c22SToby Isaac   (*labelnew)->pStart = -1;
508c58f1c22SToby Isaac   (*labelnew)->pEnd   = -1;
509c58f1c22SToby Isaac   (*labelnew)->bt     = NULL;
510c58f1c22SToby Isaac   PetscFunctionReturn(0);
511c58f1c22SToby Isaac }
512c58f1c22SToby Isaac 
513c6a43d28SMatthew G. Knepley /*@
514c6a43d28SMatthew G. Knepley   DMLabelComputeIndex - Create an index structure for membership determination, automatically determining the bounds
515c6a43d28SMatthew G. Knepley 
516c6a43d28SMatthew G. Knepley   Input Parameter:
517c6a43d28SMatthew G. Knepley . label  - The DMLabel
518c6a43d28SMatthew G. Knepley 
519c6a43d28SMatthew G. Knepley   Level: intermediate
520c6a43d28SMatthew G. Knepley 
521c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
522c6a43d28SMatthew G. Knepley @*/
523c6a43d28SMatthew G. Knepley PetscErrorCode DMLabelComputeIndex(DMLabel label)
524c6a43d28SMatthew G. Knepley {
525c6a43d28SMatthew G. Knepley   PetscInt       pStart = PETSC_MAX_INT, pEnd = -1, v;
526c6a43d28SMatthew G. Knepley   PetscErrorCode ierr;
527c6a43d28SMatthew G. Knepley 
528c6a43d28SMatthew G. Knepley   PetscFunctionBegin;
529c6a43d28SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
530c6a43d28SMatthew G. Knepley   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
531c6a43d28SMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
532c6a43d28SMatthew G. Knepley     const PetscInt *points;
533c6a43d28SMatthew G. Knepley     PetscInt       i;
534c6a43d28SMatthew G. Knepley 
535c6a43d28SMatthew G. Knepley     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
536c6a43d28SMatthew G. Knepley     for (i = 0; i < label->stratumSizes[v]; ++i) {
537c6a43d28SMatthew G. Knepley       const PetscInt point = points[i];
538c6a43d28SMatthew G. Knepley 
539c6a43d28SMatthew G. Knepley       pStart = PetscMin(point,   pStart);
540c6a43d28SMatthew G. Knepley       pEnd   = PetscMax(point+1, pEnd);
541c6a43d28SMatthew G. Knepley     }
542c6a43d28SMatthew G. Knepley     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
543c6a43d28SMatthew G. Knepley   }
544c6a43d28SMatthew G. Knepley   label->pStart = pStart == PETSC_MAX_INT ? -1 : pStart;
545c6a43d28SMatthew G. Knepley   label->pEnd   = pEnd;
546c6a43d28SMatthew G. Knepley   ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
547c6a43d28SMatthew G. Knepley   PetscFunctionReturn(0);
548c6a43d28SMatthew G. Knepley }
549c6a43d28SMatthew G. Knepley 
550c6a43d28SMatthew G. Knepley /*@
551c6a43d28SMatthew G. Knepley   DMLabelCreateIndex - Create an index structure for membership determination
552c6a43d28SMatthew G. Knepley 
553c6a43d28SMatthew G. Knepley   Input Parameters:
554c6a43d28SMatthew G. Knepley + label  - The DMLabel
555c6a43d28SMatthew G. Knepley . pStart - The smallest point
556c6a43d28SMatthew G. Knepley - pEnd   - The largest point + 1
557c6a43d28SMatthew G. Knepley 
558c6a43d28SMatthew G. Knepley   Level: intermediate
559c6a43d28SMatthew G. Knepley 
560c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelComputeIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
561c6a43d28SMatthew G. Knepley @*/
562c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
563c58f1c22SToby Isaac {
564c58f1c22SToby Isaac   PetscInt       v;
565c58f1c22SToby Isaac   PetscErrorCode ierr;
566c58f1c22SToby Isaac 
567c58f1c22SToby Isaac   PetscFunctionBegin;
568d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
5690c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
570c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
571c58f1c22SToby Isaac   label->pStart = pStart;
572c58f1c22SToby Isaac   label->pEnd   = pEnd;
573c6a43d28SMatthew G. Knepley   /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
574c58f1c22SToby Isaac   ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr);
575c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
576ad8374ffSToby Isaac     const PetscInt *points;
577c58f1c22SToby Isaac     PetscInt       i;
578c58f1c22SToby Isaac 
579ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
580c58f1c22SToby Isaac     for (i = 0; i < label->stratumSizes[v]; ++i) {
581ad8374ffSToby Isaac       const PetscInt point = points[i];
582c58f1c22SToby Isaac 
583c58f1c22SToby 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);
584c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr);
585c58f1c22SToby Isaac     }
586ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
587c58f1c22SToby Isaac   }
588c58f1c22SToby Isaac   PetscFunctionReturn(0);
589c58f1c22SToby Isaac }
590c58f1c22SToby Isaac 
591c6a43d28SMatthew G. Knepley /*@
592c6a43d28SMatthew G. Knepley   DMLabelDestroyIndex - Destroy the index structure
593c6a43d28SMatthew G. Knepley 
594c6a43d28SMatthew G. Knepley   Input Parameter:
595c6a43d28SMatthew G. Knepley . label - the DMLabel
596c6a43d28SMatthew G. Knepley 
597c6a43d28SMatthew G. Knepley   Level: intermediate
598c6a43d28SMatthew G. Knepley 
599c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
600c6a43d28SMatthew G. Knepley @*/
601c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label)
602c58f1c22SToby Isaac {
603c58f1c22SToby Isaac   PetscErrorCode ierr;
604c58f1c22SToby Isaac 
605c58f1c22SToby Isaac   PetscFunctionBegin;
606d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
607c58f1c22SToby Isaac   label->pStart = -1;
608c58f1c22SToby Isaac   label->pEnd   = -1;
6090c3c4a36SLisandro Dalcin   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
610c58f1c22SToby Isaac   PetscFunctionReturn(0);
611c58f1c22SToby Isaac }
612c58f1c22SToby Isaac 
613c58f1c22SToby Isaac /*@
614c6a43d28SMatthew G. Knepley   DMLabelGetBounds - Return the smallest and largest point in the label
615c6a43d28SMatthew G. Knepley 
616c6a43d28SMatthew G. Knepley   Input Parameter:
617c6a43d28SMatthew G. Knepley . label - the DMLabel
618c6a43d28SMatthew G. Knepley 
619c6a43d28SMatthew G. Knepley   Output Parameters:
620c6a43d28SMatthew G. Knepley + pStart - The smallest point
621c6a43d28SMatthew G. Knepley - pEnd   - The largest point + 1
622c6a43d28SMatthew G. Knepley 
623c6a43d28SMatthew G. Knepley   Note: This will compute an index for the label if one does not exist.
624c6a43d28SMatthew G. Knepley 
625c6a43d28SMatthew G. Knepley   Level: intermediate
626c6a43d28SMatthew G. Knepley 
627c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
628c6a43d28SMatthew G. Knepley @*/
629c6a43d28SMatthew G. Knepley PetscErrorCode DMLabelGetBounds(DMLabel label, PetscInt *pStart, PetscInt *pEnd)
630c6a43d28SMatthew G. Knepley {
631c6a43d28SMatthew G. Knepley   PetscErrorCode ierr;
632c6a43d28SMatthew G. Knepley 
633c6a43d28SMatthew G. Knepley   PetscFunctionBegin;
634c6a43d28SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
635c6a43d28SMatthew G. Knepley   if ((label->pStart == -1) && (label->pEnd == -1)) {ierr = DMLabelComputeIndex(label);CHKERRQ(ierr);}
636c6a43d28SMatthew G. Knepley   if (pStart) {
637c6a43d28SMatthew G. Knepley     PetscValidPointer(pStart, 2);
638c6a43d28SMatthew G. Knepley     *pStart = label->pStart;
639c6a43d28SMatthew G. Knepley   }
640c6a43d28SMatthew G. Knepley   if (pEnd) {
641c6a43d28SMatthew G. Knepley     PetscValidPointer(pEnd, 3);
642c6a43d28SMatthew G. Knepley     *pEnd = label->pEnd;
643c6a43d28SMatthew G. Knepley   }
644c6a43d28SMatthew G. Knepley   PetscFunctionReturn(0);
645c6a43d28SMatthew G. Knepley }
646c6a43d28SMatthew G. Knepley 
647c6a43d28SMatthew G. Knepley /*@
648c58f1c22SToby Isaac   DMLabelHasValue - Determine whether a label assigns the value to any point
649c58f1c22SToby Isaac 
650c58f1c22SToby Isaac   Input Parameters:
651c58f1c22SToby Isaac + label - the DMLabel
652c58f1c22SToby Isaac - value - the value
653c58f1c22SToby Isaac 
654c58f1c22SToby Isaac   Output Parameter:
655c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point
656c58f1c22SToby Isaac 
657c58f1c22SToby Isaac   Level: developer
658c58f1c22SToby Isaac 
659c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
660c58f1c22SToby Isaac @*/
661c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
662c58f1c22SToby Isaac {
663c58f1c22SToby Isaac   PetscInt v;
6640c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
665c58f1c22SToby Isaac 
666c58f1c22SToby Isaac   PetscFunctionBegin;
667d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
668c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
6690c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
6700c3c4a36SLisandro Dalcin   *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE;
671c58f1c22SToby Isaac   PetscFunctionReturn(0);
672c58f1c22SToby Isaac }
673c58f1c22SToby Isaac 
674c58f1c22SToby Isaac /*@
675c58f1c22SToby Isaac   DMLabelHasPoint - Determine whether a label assigns a value to a point
676c58f1c22SToby Isaac 
677c58f1c22SToby Isaac   Input Parameters:
678c58f1c22SToby Isaac + label - the DMLabel
679c58f1c22SToby Isaac - point - the point
680c58f1c22SToby Isaac 
681c58f1c22SToby Isaac   Output Parameter:
682c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value
683c58f1c22SToby Isaac 
684c58f1c22SToby Isaac   Note: The user must call DMLabelCreateIndex() before this function.
685c58f1c22SToby Isaac 
686c58f1c22SToby Isaac   Level: developer
687c58f1c22SToby Isaac 
688c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
689c58f1c22SToby Isaac @*/
690c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
691c58f1c22SToby Isaac {
692c58f1c22SToby Isaac   PetscErrorCode ierr;
693c58f1c22SToby Isaac 
694c58f1c22SToby Isaac   PetscFunctionBeginHot;
695d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
696c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
697c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
698c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG)
699c58f1c22SToby Isaac   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
700c58f1c22SToby 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);
701c58f1c22SToby Isaac #endif
702c58f1c22SToby Isaac   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
703c58f1c22SToby Isaac   PetscFunctionReturn(0);
704c58f1c22SToby Isaac }
705c58f1c22SToby Isaac 
706c58f1c22SToby Isaac /*@
707c58f1c22SToby Isaac   DMLabelStratumHasPoint - Return true if the stratum contains a point
708c58f1c22SToby Isaac 
709c58f1c22SToby Isaac   Input Parameters:
710c58f1c22SToby Isaac + label - the DMLabel
711c58f1c22SToby Isaac . value - the stratum value
712c58f1c22SToby Isaac - point - the point
713c58f1c22SToby Isaac 
714c58f1c22SToby Isaac   Output Parameter:
715c58f1c22SToby Isaac . contains - true if the stratum contains the point
716c58f1c22SToby Isaac 
717c58f1c22SToby Isaac   Level: intermediate
718c58f1c22SToby Isaac 
719c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
720c58f1c22SToby Isaac @*/
721c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
722c58f1c22SToby Isaac {
723c58f1c22SToby Isaac   PetscInt       v;
724c58f1c22SToby Isaac   PetscErrorCode ierr;
725c58f1c22SToby Isaac 
7260c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
727d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
728c58f1c22SToby Isaac   PetscValidPointer(contains, 4);
729c58f1c22SToby Isaac   *contains = PETSC_FALSE;
7300c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
7310c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
7320c3c4a36SLisandro Dalcin 
733ad8374ffSToby Isaac   if (label->validIS[v]) {
734c58f1c22SToby Isaac     PetscInt i;
735c58f1c22SToby Isaac 
736a2d74346SToby Isaac     ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
7370c3c4a36SLisandro Dalcin     if (i >= 0) *contains = PETSC_TRUE;
738c58f1c22SToby Isaac   } else {
739c58f1c22SToby Isaac     PetscBool has;
740c58f1c22SToby Isaac 
741b9907514SLisandro Dalcin     ierr = PetscHSetIHas(label->ht[v], point, &has);CHKERRQ(ierr);
7420c3c4a36SLisandro Dalcin     if (has) *contains = PETSC_TRUE;
743c58f1c22SToby Isaac   }
744c58f1c22SToby Isaac   PetscFunctionReturn(0);
745c58f1c22SToby Isaac }
746c58f1c22SToby Isaac 
74784f0b6dfSMatthew G. Knepley /*@
7485aa44df4SToby Isaac   DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
7495aa44df4SToby Isaac   When a label is created, it is initialized to -1.
7505aa44df4SToby Isaac 
7515aa44df4SToby Isaac   Input parameter:
7525aa44df4SToby Isaac . label - a DMLabel object
7535aa44df4SToby Isaac 
7545aa44df4SToby Isaac   Output parameter:
7555aa44df4SToby Isaac . defaultValue - the default value
7565aa44df4SToby Isaac 
7575aa44df4SToby Isaac   Level: beginner
7585aa44df4SToby Isaac 
7595aa44df4SToby Isaac .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
76084f0b6dfSMatthew G. Knepley @*/
7615aa44df4SToby Isaac PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
7625aa44df4SToby Isaac {
7635aa44df4SToby Isaac   PetscFunctionBegin;
764d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
7655aa44df4SToby Isaac   *defaultValue = label->defaultValue;
7665aa44df4SToby Isaac   PetscFunctionReturn(0);
7675aa44df4SToby Isaac }
7685aa44df4SToby Isaac 
76984f0b6dfSMatthew G. Knepley /*@
7705aa44df4SToby Isaac   DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
7715aa44df4SToby Isaac   When a label is created, it is initialized to -1.
7725aa44df4SToby Isaac 
7735aa44df4SToby Isaac   Input parameter:
7745aa44df4SToby Isaac . label - a DMLabel object
7755aa44df4SToby Isaac 
7765aa44df4SToby Isaac   Output parameter:
7775aa44df4SToby Isaac . defaultValue - the default value
7785aa44df4SToby Isaac 
7795aa44df4SToby Isaac   Level: beginner
7805aa44df4SToby Isaac 
7815aa44df4SToby Isaac .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
78284f0b6dfSMatthew G. Knepley @*/
7835aa44df4SToby Isaac PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
7845aa44df4SToby Isaac {
7855aa44df4SToby Isaac   PetscFunctionBegin;
786d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
7875aa44df4SToby Isaac   label->defaultValue = defaultValue;
7885aa44df4SToby Isaac   PetscFunctionReturn(0);
7895aa44df4SToby Isaac }
7905aa44df4SToby Isaac 
791c58f1c22SToby Isaac /*@
7925aa44df4SToby 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())
793c58f1c22SToby Isaac 
794c58f1c22SToby Isaac   Input Parameters:
795c58f1c22SToby Isaac + label - the DMLabel
796c58f1c22SToby Isaac - point - the point
797c58f1c22SToby Isaac 
798c58f1c22SToby Isaac   Output Parameter:
7998e68afcfSMatthew G. Knepley . value - The point value, or the default value (-1 by default)
800c58f1c22SToby Isaac 
801c58f1c22SToby Isaac   Level: intermediate
802c58f1c22SToby Isaac 
8035aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
804c58f1c22SToby Isaac @*/
805c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
806c58f1c22SToby Isaac {
807c58f1c22SToby Isaac   PetscInt       v;
808c58f1c22SToby Isaac   PetscErrorCode ierr;
809c58f1c22SToby Isaac 
8100c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
811d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
812c58f1c22SToby Isaac   PetscValidPointer(value, 3);
8135aa44df4SToby Isaac   *value = label->defaultValue;
814c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
815ad8374ffSToby Isaac     if (label->validIS[v]) {
816c58f1c22SToby Isaac       PetscInt i;
817c58f1c22SToby Isaac 
818a2d74346SToby Isaac       ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
819c58f1c22SToby Isaac       if (i >= 0) {
820c58f1c22SToby Isaac         *value = label->stratumValues[v];
821c58f1c22SToby Isaac         break;
822c58f1c22SToby Isaac       }
823c58f1c22SToby Isaac     } else {
824c58f1c22SToby Isaac       PetscBool has;
825c58f1c22SToby Isaac 
826b9907514SLisandro Dalcin       ierr = PetscHSetIHas(label->ht[v], point, &has);CHKERRQ(ierr);
827c58f1c22SToby Isaac       if (has) {
828c58f1c22SToby Isaac         *value = label->stratumValues[v];
829c58f1c22SToby Isaac         break;
830c58f1c22SToby Isaac       }
831c58f1c22SToby Isaac     }
832c58f1c22SToby Isaac   }
833c58f1c22SToby Isaac   PetscFunctionReturn(0);
834c58f1c22SToby Isaac }
835c58f1c22SToby Isaac 
836c58f1c22SToby Isaac /*@
837367003a6SStefano 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.
838c58f1c22SToby Isaac 
839c58f1c22SToby Isaac   Input Parameters:
840c58f1c22SToby Isaac + label - the DMLabel
841c58f1c22SToby Isaac . point - the point
842c58f1c22SToby Isaac - value - The point value
843c58f1c22SToby Isaac 
844c58f1c22SToby Isaac   Level: intermediate
845c58f1c22SToby Isaac 
8465aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
847c58f1c22SToby Isaac @*/
848c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
849c58f1c22SToby Isaac {
850c58f1c22SToby Isaac   PetscInt       v;
851c58f1c22SToby Isaac   PetscErrorCode ierr;
852c58f1c22SToby Isaac 
853c58f1c22SToby Isaac   PetscFunctionBegin;
854d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
8550c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
8565aa44df4SToby Isaac   if (value == label->defaultValue) PetscFunctionReturn(0);
857b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
858c58f1c22SToby Isaac   /* Set key */
8590c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
860e8f14785SLisandro Dalcin   ierr = PetscHSetIAdd(label->ht[v], point);CHKERRQ(ierr);
861c58f1c22SToby Isaac   PetscFunctionReturn(0);
862c58f1c22SToby Isaac }
863c58f1c22SToby Isaac 
864c58f1c22SToby Isaac /*@
865c58f1c22SToby Isaac   DMLabelClearValue - Clear the value a label assigns to a point
866c58f1c22SToby Isaac 
867c58f1c22SToby Isaac   Input Parameters:
868c58f1c22SToby Isaac + label - the DMLabel
869c58f1c22SToby Isaac . point - the point
870c58f1c22SToby Isaac - value - The point value
871c58f1c22SToby Isaac 
872c58f1c22SToby Isaac   Level: intermediate
873c58f1c22SToby Isaac 
874c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
875c58f1c22SToby Isaac @*/
876c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
877c58f1c22SToby Isaac {
878ad8374ffSToby Isaac   PetscInt       v;
879c58f1c22SToby Isaac   PetscErrorCode ierr;
880c58f1c22SToby Isaac 
881c58f1c22SToby Isaac   PetscFunctionBegin;
882d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
883c58f1c22SToby Isaac   /* Find label value */
8840c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
8850c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
8860c3c4a36SLisandro Dalcin 
887eeed21e7SToby Isaac   if (label->bt) {
888eeed21e7SToby 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);
889eeed21e7SToby Isaac     ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
890eeed21e7SToby Isaac   }
8910c3c4a36SLisandro Dalcin 
8920c3c4a36SLisandro Dalcin   /* Delete key */
8930c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
894e8f14785SLisandro Dalcin   ierr = PetscHSetIDel(label->ht[v], point);CHKERRQ(ierr);
895c58f1c22SToby Isaac   PetscFunctionReturn(0);
896c58f1c22SToby Isaac }
897c58f1c22SToby Isaac 
898c58f1c22SToby Isaac /*@
899c58f1c22SToby Isaac   DMLabelInsertIS - Set all points in the IS to a value
900c58f1c22SToby Isaac 
901c58f1c22SToby Isaac   Input Parameters:
902c58f1c22SToby Isaac + label - the DMLabel
903c58f1c22SToby Isaac . is    - the point IS
904c58f1c22SToby Isaac - value - The point value
905c58f1c22SToby Isaac 
906c58f1c22SToby Isaac   Level: intermediate
907c58f1c22SToby Isaac 
908c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
909c58f1c22SToby Isaac @*/
910c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
911c58f1c22SToby Isaac {
9120c3c4a36SLisandro Dalcin   PetscInt        v, n, p;
913c58f1c22SToby Isaac   const PetscInt *points;
914c58f1c22SToby Isaac   PetscErrorCode  ierr;
915c58f1c22SToby Isaac 
916c58f1c22SToby Isaac   PetscFunctionBegin;
917d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
918c58f1c22SToby Isaac   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
9190c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
9200c3c4a36SLisandro Dalcin   if (value == label->defaultValue) PetscFunctionReturn(0);
921b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
9220c3c4a36SLisandro Dalcin   /* Set keys */
9230c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
924c58f1c22SToby Isaac   ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr);
925c58f1c22SToby Isaac   ierr = ISGetIndices(is, &points);CHKERRQ(ierr);
926e8f14785SLisandro Dalcin   for (p = 0; p < n; ++p) {ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);}
927c58f1c22SToby Isaac   ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr);
928c58f1c22SToby Isaac   PetscFunctionReturn(0);
929c58f1c22SToby Isaac }
930c58f1c22SToby Isaac 
93184f0b6dfSMatthew G. Knepley /*@
93284f0b6dfSMatthew G. Knepley   DMLabelGetNumValues - Get the number of values that the DMLabel takes
93384f0b6dfSMatthew G. Knepley 
93484f0b6dfSMatthew G. Knepley   Input Parameter:
93584f0b6dfSMatthew G. Knepley . label - the DMLabel
93684f0b6dfSMatthew G. Knepley 
93784f0b6dfSMatthew G. Knepley   Output Paramater:
93884f0b6dfSMatthew G. Knepley . numValues - the number of values
93984f0b6dfSMatthew G. Knepley 
94084f0b6dfSMatthew G. Knepley   Level: intermediate
94184f0b6dfSMatthew G. Knepley 
94284f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
94384f0b6dfSMatthew G. Knepley @*/
944c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
945c58f1c22SToby Isaac {
946c58f1c22SToby Isaac   PetscFunctionBegin;
947d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
948c58f1c22SToby Isaac   PetscValidPointer(numValues, 2);
949c58f1c22SToby Isaac   *numValues = label->numStrata;
950c58f1c22SToby Isaac   PetscFunctionReturn(0);
951c58f1c22SToby Isaac }
952c58f1c22SToby Isaac 
95384f0b6dfSMatthew G. Knepley /*@
95484f0b6dfSMatthew G. Knepley   DMLabelGetValueIS - Get an IS of all values that the DMlabel takes
95584f0b6dfSMatthew G. Knepley 
95684f0b6dfSMatthew G. Knepley   Input Parameter:
95784f0b6dfSMatthew G. Knepley . label - the DMLabel
95884f0b6dfSMatthew G. Knepley 
95984f0b6dfSMatthew G. Knepley   Output Paramater:
96084f0b6dfSMatthew G. Knepley . is    - the value IS
96184f0b6dfSMatthew G. Knepley 
96284f0b6dfSMatthew G. Knepley   Level: intermediate
96384f0b6dfSMatthew G. Knepley 
96484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
96584f0b6dfSMatthew G. Knepley @*/
966c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
967c58f1c22SToby Isaac {
968c58f1c22SToby Isaac   PetscErrorCode ierr;
969c58f1c22SToby Isaac 
970c58f1c22SToby Isaac   PetscFunctionBegin;
971d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
972c58f1c22SToby Isaac   PetscValidPointer(values, 2);
973c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr);
974c58f1c22SToby Isaac   PetscFunctionReturn(0);
975c58f1c22SToby Isaac }
976c58f1c22SToby Isaac 
97784f0b6dfSMatthew G. Knepley /*@
97884f0b6dfSMatthew G. Knepley   DMLabelHasStratum - Determine whether points exist with the given value
97984f0b6dfSMatthew G. Knepley 
98084f0b6dfSMatthew G. Knepley   Input Parameters:
98184f0b6dfSMatthew G. Knepley + label - the DMLabel
98284f0b6dfSMatthew G. Knepley - value - the stratum value
98384f0b6dfSMatthew G. Knepley 
98484f0b6dfSMatthew G. Knepley   Output Paramater:
98584f0b6dfSMatthew G. Knepley . exists - Flag saying whether points exist
98684f0b6dfSMatthew G. Knepley 
98784f0b6dfSMatthew G. Knepley   Level: intermediate
98884f0b6dfSMatthew G. Knepley 
98984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
99084f0b6dfSMatthew G. Knepley @*/
991fada774cSMatthew G. Knepley PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
992fada774cSMatthew G. Knepley {
993fada774cSMatthew G. Knepley   PetscInt       v;
9940c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
995fada774cSMatthew G. Knepley 
996fada774cSMatthew G. Knepley   PetscFunctionBegin;
997d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
998fada774cSMatthew G. Knepley   PetscValidPointer(exists, 3);
9990c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
10000c3c4a36SLisandro Dalcin   *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE;
1001fada774cSMatthew G. Knepley   PetscFunctionReturn(0);
1002fada774cSMatthew G. Knepley }
1003fada774cSMatthew G. Knepley 
100484f0b6dfSMatthew G. Knepley /*@
100584f0b6dfSMatthew G. Knepley   DMLabelGetStratumSize - Get the size of a stratum
100684f0b6dfSMatthew G. Knepley 
100784f0b6dfSMatthew G. Knepley   Input Parameters:
100884f0b6dfSMatthew G. Knepley + label - the DMLabel
100984f0b6dfSMatthew G. Knepley - value - the stratum value
101084f0b6dfSMatthew G. Knepley 
101184f0b6dfSMatthew G. Knepley   Output Paramater:
101284f0b6dfSMatthew G. Knepley . size - The number of points in the stratum
101384f0b6dfSMatthew G. Knepley 
101484f0b6dfSMatthew G. Knepley   Level: intermediate
101584f0b6dfSMatthew G. Knepley 
101684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
101784f0b6dfSMatthew G. Knepley @*/
1018c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
1019c58f1c22SToby Isaac {
1020c58f1c22SToby Isaac   PetscInt       v;
1021c58f1c22SToby Isaac   PetscErrorCode ierr;
1022c58f1c22SToby Isaac 
1023c58f1c22SToby Isaac   PetscFunctionBegin;
1024d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1025c58f1c22SToby Isaac   PetscValidPointer(size, 3);
1026c58f1c22SToby Isaac   *size = 0;
10270c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
10280c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
1029c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
1030c58f1c22SToby Isaac   *size = label->stratumSizes[v];
1031c58f1c22SToby Isaac   PetscFunctionReturn(0);
1032c58f1c22SToby Isaac }
1033c58f1c22SToby Isaac 
103484f0b6dfSMatthew G. Knepley /*@
103584f0b6dfSMatthew G. Knepley   DMLabelGetStratumBounds - Get the largest and smallest point of a stratum
103684f0b6dfSMatthew G. Knepley 
103784f0b6dfSMatthew G. Knepley   Input Parameters:
103884f0b6dfSMatthew G. Knepley + label - the DMLabel
103984f0b6dfSMatthew G. Knepley - value - the stratum value
104084f0b6dfSMatthew G. Knepley 
104184f0b6dfSMatthew G. Knepley   Output Paramaters:
104284f0b6dfSMatthew G. Knepley + start - the smallest point in the stratum
104384f0b6dfSMatthew G. Knepley - end - the largest point in the stratum
104484f0b6dfSMatthew G. Knepley 
104584f0b6dfSMatthew G. Knepley   Level: intermediate
104684f0b6dfSMatthew G. Knepley 
104784f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
104884f0b6dfSMatthew G. Knepley @*/
1049c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
1050c58f1c22SToby Isaac {
10510c3c4a36SLisandro Dalcin   PetscInt       v, min, max;
1052c58f1c22SToby Isaac   PetscErrorCode ierr;
1053c58f1c22SToby Isaac 
1054c58f1c22SToby Isaac   PetscFunctionBegin;
1055d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1056c58f1c22SToby Isaac   if (start) {PetscValidPointer(start, 3); *start = 0;}
1057c58f1c22SToby Isaac   if (end)   {PetscValidPointer(end,   4); *end   = 0;}
10580c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
10590c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
1060c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
10610c3c4a36SLisandro Dalcin   if (label->stratumSizes[v] <= 0) PetscFunctionReturn(0);
1062d6cb179aSToby Isaac   ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr);
1063d6cb179aSToby Isaac   if (start) *start = min;
1064d6cb179aSToby Isaac   if (end)   *end   = max+1;
1065c58f1c22SToby Isaac   PetscFunctionReturn(0);
1066c58f1c22SToby Isaac }
1067c58f1c22SToby Isaac 
106884f0b6dfSMatthew G. Knepley /*@
106984f0b6dfSMatthew G. Knepley   DMLabelGetStratumIS - Get an IS with the stratum points
107084f0b6dfSMatthew G. Knepley 
107184f0b6dfSMatthew G. Knepley   Input Parameters:
107284f0b6dfSMatthew G. Knepley + label - the DMLabel
107384f0b6dfSMatthew G. Knepley - value - the stratum value
107484f0b6dfSMatthew G. Knepley 
107584f0b6dfSMatthew G. Knepley   Output Paramater:
107684f0b6dfSMatthew G. Knepley . points - The stratum points
107784f0b6dfSMatthew G. Knepley 
107884f0b6dfSMatthew G. Knepley   Level: intermediate
107984f0b6dfSMatthew G. Knepley 
108084f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
108184f0b6dfSMatthew G. Knepley @*/
1082c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
1083c58f1c22SToby Isaac {
1084c58f1c22SToby Isaac   PetscInt       v;
1085c58f1c22SToby Isaac   PetscErrorCode ierr;
1086c58f1c22SToby Isaac 
1087c58f1c22SToby Isaac   PetscFunctionBegin;
1088d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1089c58f1c22SToby Isaac   PetscValidPointer(points, 3);
1090c58f1c22SToby Isaac   *points = NULL;
10910c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
10920c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
1093c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
1094ad8374ffSToby Isaac   ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr);
1095ad8374ffSToby Isaac   *points = label->points[v];
1096c58f1c22SToby Isaac   PetscFunctionReturn(0);
1097c58f1c22SToby Isaac }
1098c58f1c22SToby Isaac 
109984f0b6dfSMatthew G. Knepley /*@
110084f0b6dfSMatthew G. Knepley   DMLabelSetStratumIS - Set the stratum points using an IS
110184f0b6dfSMatthew G. Knepley 
110284f0b6dfSMatthew G. Knepley   Input Parameters:
110384f0b6dfSMatthew G. Knepley + label - the DMLabel
110484f0b6dfSMatthew G. Knepley . value - the stratum value
110584f0b6dfSMatthew G. Knepley - points - The stratum points
110684f0b6dfSMatthew G. Knepley 
110784f0b6dfSMatthew G. Knepley   Level: intermediate
110884f0b6dfSMatthew G. Knepley 
110984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
111084f0b6dfSMatthew G. Knepley @*/
11114de306b1SToby Isaac PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
11124de306b1SToby Isaac {
11130c3c4a36SLisandro Dalcin   PetscInt       v;
11144de306b1SToby Isaac   PetscErrorCode ierr;
11154de306b1SToby Isaac 
11164de306b1SToby Isaac   PetscFunctionBegin;
1117d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1118d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
1119b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
11204de306b1SToby Isaac   if (is == label->points[v]) PetscFunctionReturn(0);
11214de306b1SToby Isaac   ierr = DMLabelClearStratum(label, value);CHKERRQ(ierr);
11224de306b1SToby Isaac   ierr = ISGetLocalSize(is, &(label->stratumSizes[v]));CHKERRQ(ierr);
11234de306b1SToby Isaac   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
11244de306b1SToby Isaac   ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
11250c3c4a36SLisandro Dalcin   label->points[v] = is;
11260c3c4a36SLisandro Dalcin   label->validIS[v] = PETSC_TRUE;
11274de306b1SToby Isaac   if (label->bt) {
11284de306b1SToby Isaac     const PetscInt *points;
11294de306b1SToby Isaac     PetscInt p;
11304de306b1SToby Isaac 
11314de306b1SToby Isaac     ierr = ISGetIndices(is,&points);CHKERRQ(ierr);
11324de306b1SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
11334de306b1SToby Isaac       const PetscInt point = points[p];
11344de306b1SToby Isaac 
11354de306b1SToby 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);
11364de306b1SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
11374de306b1SToby Isaac     }
11384de306b1SToby Isaac   }
11394de306b1SToby Isaac   PetscFunctionReturn(0);
11404de306b1SToby Isaac }
11414de306b1SToby Isaac 
114284f0b6dfSMatthew G. Knepley /*@
114384f0b6dfSMatthew G. Knepley   DMLabelClearStratum - Remove a stratum
11444de306b1SToby Isaac 
114584f0b6dfSMatthew G. Knepley   Input Parameters:
114684f0b6dfSMatthew G. Knepley + label - the DMLabel
114784f0b6dfSMatthew G. Knepley - value - the stratum value
114884f0b6dfSMatthew G. Knepley 
114984f0b6dfSMatthew G. Knepley   Level: intermediate
115084f0b6dfSMatthew G. Knepley 
115184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
115284f0b6dfSMatthew G. Knepley @*/
1153c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
1154c58f1c22SToby Isaac {
1155c58f1c22SToby Isaac   PetscInt       v;
1156c58f1c22SToby Isaac   PetscErrorCode ierr;
1157c58f1c22SToby Isaac 
1158c58f1c22SToby Isaac   PetscFunctionBegin;
1159d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
11600c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
11610c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
11624de306b1SToby Isaac   if (label->validIS[v]) {
11634de306b1SToby Isaac     if (label->bt) {
1164c58f1c22SToby Isaac       PetscInt       i;
1165ad8374ffSToby Isaac       const PetscInt *points;
1166c58f1c22SToby Isaac 
1167ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
1168c58f1c22SToby Isaac       for (i = 0; i < label->stratumSizes[v]; ++i) {
1169ad8374ffSToby Isaac         const PetscInt point = points[i];
1170c58f1c22SToby Isaac 
1171c58f1c22SToby 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);
1172c58f1c22SToby Isaac         ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
1173c58f1c22SToby Isaac       }
1174ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
1175c58f1c22SToby Isaac     }
1176c58f1c22SToby Isaac     label->stratumSizes[v] = 0;
11770c3c4a36SLisandro Dalcin     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
11780c3c4a36SLisandro Dalcin     ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr);
11790c3c4a36SLisandro Dalcin     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
1180c58f1c22SToby Isaac   } else {
1181b9907514SLisandro Dalcin     ierr = PetscHSetIClear(label->ht[v]);CHKERRQ(ierr);
1182c58f1c22SToby Isaac   }
1183c58f1c22SToby Isaac   PetscFunctionReturn(0);
1184c58f1c22SToby Isaac }
1185c58f1c22SToby Isaac 
118684f0b6dfSMatthew G. Knepley /*@
118784f0b6dfSMatthew G. Knepley   DMLabelFilter - Remove all points outside of [start, end)
118884f0b6dfSMatthew G. Knepley 
118984f0b6dfSMatthew G. Knepley   Input Parameters:
119084f0b6dfSMatthew G. Knepley + label - the DMLabel
119184f0b6dfSMatthew G. Knepley . start - the first point
119284f0b6dfSMatthew G. Knepley - end - the last point
119384f0b6dfSMatthew G. Knepley 
119484f0b6dfSMatthew G. Knepley   Level: intermediate
119584f0b6dfSMatthew G. Knepley 
119684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
119784f0b6dfSMatthew G. Knepley @*/
1198c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
1199c58f1c22SToby Isaac {
1200c58f1c22SToby Isaac   PetscInt       v;
1201c58f1c22SToby Isaac   PetscErrorCode ierr;
1202c58f1c22SToby Isaac 
1203c58f1c22SToby Isaac   PetscFunctionBegin;
1204d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
12050c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
1206c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1207c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
1208c58f1c22SToby Isaac     PetscInt off, q;
1209ad8374ffSToby Isaac     const PetscInt *points;
1210033405d5SLisandro Dalcin     PetscInt numPointsNew = 0, *pointsNew = NULL;
1211c58f1c22SToby Isaac 
1212ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
1213033405d5SLisandro Dalcin     for (q = 0; q < label->stratumSizes[v]; ++q)
1214033405d5SLisandro Dalcin       if (points[q] >= start && points[q] < end)
1215033405d5SLisandro Dalcin         numPointsNew++;
1216033405d5SLisandro Dalcin     ierr = PetscMalloc1(numPointsNew, &pointsNew);CHKERRQ(ierr);
1217c58f1c22SToby Isaac     for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) {
1218033405d5SLisandro Dalcin       if (points[q] >= start && points[q] < end)
1219033405d5SLisandro Dalcin         pointsNew[off++] = points[q];
1220ad8374ffSToby Isaac     }
1221ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
1222033405d5SLisandro Dalcin 
1223033405d5SLisandro Dalcin     label->stratumSizes[v] = numPointsNew;
1224033405d5SLisandro Dalcin     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
1225033405d5SLisandro Dalcin     ierr = ISCreateGeneral(PETSC_COMM_SELF,numPointsNew, pointsNew, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr);
1226033405d5SLisandro Dalcin     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
1227c58f1c22SToby Isaac   }
1228c58f1c22SToby Isaac   ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr);
1229c58f1c22SToby Isaac   PetscFunctionReturn(0);
1230c58f1c22SToby Isaac }
1231c58f1c22SToby Isaac 
123284f0b6dfSMatthew G. Knepley /*@
123384f0b6dfSMatthew G. Knepley   DMLabelPermute - Create a new label with permuted points
123484f0b6dfSMatthew G. Knepley 
123584f0b6dfSMatthew G. Knepley   Input Parameters:
123684f0b6dfSMatthew G. Knepley + label - the DMLabel
123784f0b6dfSMatthew G. Knepley - permutation - the point permutation
123884f0b6dfSMatthew G. Knepley 
123984f0b6dfSMatthew G. Knepley   Output Parameter:
124084f0b6dfSMatthew G. Knepley . labelnew - the new label containing the permuted points
124184f0b6dfSMatthew G. Knepley 
124284f0b6dfSMatthew G. Knepley   Level: intermediate
124384f0b6dfSMatthew G. Knepley 
124484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
124584f0b6dfSMatthew G. Knepley @*/
1246c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
1247c58f1c22SToby Isaac {
1248c58f1c22SToby Isaac   const PetscInt *perm;
1249c58f1c22SToby Isaac   PetscInt        numValues, numPoints, v, q;
1250c58f1c22SToby Isaac   PetscErrorCode  ierr;
1251c58f1c22SToby Isaac 
1252c58f1c22SToby Isaac   PetscFunctionBegin;
1253d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1254d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(permutation, IS_CLASSID, 2);
1255c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1256c58f1c22SToby Isaac   ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr);
1257c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr);
1258c58f1c22SToby Isaac   ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr);
1259c58f1c22SToby Isaac   ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr);
1260c58f1c22SToby Isaac   for (v = 0; v < numValues; ++v) {
1261c58f1c22SToby Isaac     const PetscInt size   = (*labelNew)->stratumSizes[v];
1262ad8374ffSToby Isaac     const PetscInt *points;
1263ad8374ffSToby Isaac     PetscInt *pointsNew;
1264c58f1c22SToby Isaac 
1265ad8374ffSToby Isaac     ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1266ad8374ffSToby Isaac     ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr);
1267c58f1c22SToby Isaac     for (q = 0; q < size; ++q) {
1268ad8374ffSToby Isaac       const PetscInt point = points[q];
1269c58f1c22SToby Isaac 
1270c58f1c22SToby 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);
1271ad8374ffSToby Isaac       pointsNew[q] = perm[point];
1272c58f1c22SToby Isaac     }
1273ad8374ffSToby Isaac     ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1274ad8374ffSToby Isaac     ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr);
1275ad8374ffSToby Isaac     ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr);
1276fa8e8ae5SToby Isaac     if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) {
1277fa8e8ae5SToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));CHKERRQ(ierr);
1278fa8e8ae5SToby Isaac       ierr = PetscFree(pointsNew);CHKERRQ(ierr);
1279fa8e8ae5SToby Isaac     } else {
1280ad8374ffSToby Isaac       ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr);
1281fa8e8ae5SToby Isaac     }
1282ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr);
1283c58f1c22SToby Isaac   }
1284c58f1c22SToby Isaac   ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr);
1285c58f1c22SToby Isaac   if (label->bt) {
1286c58f1c22SToby Isaac     ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
1287c58f1c22SToby Isaac     ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
1288c58f1c22SToby Isaac   }
1289c58f1c22SToby Isaac   PetscFunctionReturn(0);
1290c58f1c22SToby Isaac }
1291c58f1c22SToby Isaac 
129226c55118SMichael Lange PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
129326c55118SMichael Lange {
129426c55118SMichael Lange   MPI_Comm       comm;
129526c55118SMichael Lange   PetscInt       s, l, nroots, nleaves, dof, offset, size;
129626c55118SMichael Lange   PetscInt      *remoteOffsets, *rootStrata, *rootIdx;
129726c55118SMichael Lange   PetscSection   rootSection;
129826c55118SMichael Lange   PetscSF        labelSF;
129926c55118SMichael Lange   PetscErrorCode ierr;
130026c55118SMichael Lange 
130126c55118SMichael Lange   PetscFunctionBegin;
130226c55118SMichael Lange   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
130326c55118SMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
130426c55118SMichael Lange   /* Build a section of stratum values per point, generate the according SF
130526c55118SMichael Lange      and distribute point-wise stratum values to leaves. */
130626c55118SMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
130726c55118SMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
130826c55118SMichael Lange   ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr);
130926c55118SMichael Lange   if (label) {
131026c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1311ad8374ffSToby Isaac       const PetscInt *points;
1312ad8374ffSToby Isaac 
1313ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
131426c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1315ad8374ffSToby Isaac         ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr);
1316ad8374ffSToby Isaac         ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr);
131726c55118SMichael Lange       }
1318ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
131926c55118SMichael Lange     }
132026c55118SMichael Lange   }
132126c55118SMichael Lange   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
132226c55118SMichael Lange   /* Create a point-wise array of stratum values */
132326c55118SMichael Lange   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
132426c55118SMichael Lange   ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr);
132526c55118SMichael Lange   ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr);
132626c55118SMichael Lange   if (label) {
132726c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1328ad8374ffSToby Isaac       const PetscInt *points;
1329ad8374ffSToby Isaac 
1330ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
133126c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1332ad8374ffSToby Isaac         const PetscInt p = points[l];
133326c55118SMichael Lange         ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr);
133426c55118SMichael Lange         rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
133526c55118SMichael Lange       }
1336ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
133726c55118SMichael Lange     }
133826c55118SMichael Lange   }
133926c55118SMichael Lange   /* Build SF that maps label points to remote processes */
134026c55118SMichael Lange   ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr);
134126c55118SMichael Lange   ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr);
134226c55118SMichael Lange   ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr);
134326c55118SMichael Lange   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
134426c55118SMichael Lange   /* Send the strata for each point over the derived SF */
134526c55118SMichael Lange   ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr);
134626c55118SMichael Lange   ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr);
134726c55118SMichael Lange   ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
134826c55118SMichael Lange   ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
134926c55118SMichael Lange   /* Clean up */
135026c55118SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
135126c55118SMichael Lange   ierr = PetscFree(rootIdx);CHKERRQ(ierr);
135226c55118SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
135326c55118SMichael Lange   ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr);
135426c55118SMichael Lange   PetscFunctionReturn(0);
135526c55118SMichael Lange }
135626c55118SMichael Lange 
135784f0b6dfSMatthew G. Knepley /*@
135884f0b6dfSMatthew G. Knepley   DMLabelDistribute - Create a new label pushed forward over the PetscSF
135984f0b6dfSMatthew G. Knepley 
136084f0b6dfSMatthew G. Knepley   Input Parameters:
136184f0b6dfSMatthew G. Knepley + label - the DMLabel
136284f0b6dfSMatthew G. Knepley - sf    - the map from old to new distribution
136384f0b6dfSMatthew G. Knepley 
136484f0b6dfSMatthew G. Knepley   Output Parameter:
136584f0b6dfSMatthew G. Knepley . labelnew - the new redistributed label
136684f0b6dfSMatthew G. Knepley 
136784f0b6dfSMatthew G. Knepley   Level: intermediate
136884f0b6dfSMatthew G. Knepley 
136984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
137084f0b6dfSMatthew G. Knepley @*/
1371c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1372c58f1c22SToby Isaac {
1373c58f1c22SToby Isaac   MPI_Comm       comm;
137426c55118SMichael Lange   PetscSection   leafSection;
137526c55118SMichael Lange   PetscInt       p, pStart, pEnd, s, size, dof, offset, stratum;
137626c55118SMichael Lange   PetscInt      *leafStrata, *strataIdx;
1377ad8374ffSToby Isaac   PetscInt     **points;
1378d67d17b1SMatthew G. Knepley   const char    *lname = NULL;
1379c58f1c22SToby Isaac   char          *name;
1380c58f1c22SToby Isaac   PetscInt       nameSize;
1381e8f14785SLisandro Dalcin   PetscHSetI     stratumHash;
1382c58f1c22SToby Isaac   size_t         len = 0;
138326c55118SMichael Lange   PetscMPIInt    rank;
1384c58f1c22SToby Isaac   PetscErrorCode ierr;
1385c58f1c22SToby Isaac 
1386c58f1c22SToby Isaac   PetscFunctionBegin;
1387d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1388f018e600SMatthew G. Knepley   if (label) {
1389f018e600SMatthew G. Knepley     PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1390f018e600SMatthew G. Knepley     ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1391f018e600SMatthew G. Knepley   }
1392c58f1c22SToby Isaac   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1393c58f1c22SToby Isaac   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1394c58f1c22SToby Isaac   /* Bcast name */
1395d67d17b1SMatthew G. Knepley   if (!rank) {
1396d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1397d67d17b1SMatthew G. Knepley     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1398d67d17b1SMatthew G. Knepley   }
1399c58f1c22SToby Isaac   nameSize = len;
1400c58f1c22SToby Isaac   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1401c58f1c22SToby Isaac   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1402d67d17b1SMatthew G. Knepley   if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);}
1403c58f1c22SToby Isaac   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1404d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
1405c58f1c22SToby Isaac   ierr = PetscFree(name);CHKERRQ(ierr);
140677d236dfSMichael Lange   /* Bcast defaultValue */
140777d236dfSMichael Lange   if (!rank) (*labelNew)->defaultValue = label->defaultValue;
140877d236dfSMichael Lange   ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
140926c55118SMichael Lange   /* Distribute stratum values over the SF and get the point mapping on the receiver */
141026c55118SMichael Lange   ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr);
14115cbdf6fcSMichael Lange   /* Determine received stratum values and initialise new label*/
1412e8f14785SLisandro Dalcin   ierr = PetscHSetICreate(&stratumHash);CHKERRQ(ierr);
141326c55118SMichael Lange   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
1414e8f14785SLisandro Dalcin   for (p = 0; p < size; ++p) {ierr = PetscHSetIAdd(stratumHash, leafStrata[p]);CHKERRQ(ierr);}
1415e8f14785SLisandro Dalcin   ierr = PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata);CHKERRQ(ierr);
1416ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr);
1417ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
14185cbdf6fcSMichael Lange   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr);
14195cbdf6fcSMichael Lange   /* Turn leafStrata into indices rather than stratum values */
14205cbdf6fcSMichael Lange   offset = 0;
1421e8f14785SLisandro Dalcin   ierr = PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr);
1422a205f1fdSToby Isaac   ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr);
14235cbdf6fcSMichael Lange   for (p = 0; p < size; ++p) {
1424231b9e6fSMatthew G. Knepley     for (s = 0; s < (*labelNew)->numStrata; ++s) {
1425231b9e6fSMatthew G. Knepley       if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
14265cbdf6fcSMichael Lange     }
14275cbdf6fcSMichael Lange   }
1428c58f1c22SToby Isaac   /* Rebuild the point strata on the receiver */
1429c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr);
1430c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1431c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1432c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1433c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1434c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1435c58f1c22SToby Isaac       (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1436c58f1c22SToby Isaac     }
1437c58f1c22SToby Isaac   }
1438c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr);
1439c58f1c22SToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr);
1440ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr);
1441c58f1c22SToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1442e8f14785SLisandro Dalcin     ierr = PetscHSetICreate(&(*labelNew)->ht[s]);CHKERRQ(ierr);
1443ad8374ffSToby Isaac     ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr);
1444c58f1c22SToby Isaac   }
1445c58f1c22SToby Isaac   /* Insert points into new strata */
1446c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr);
1447c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1448c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1449c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1450c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1451c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1452c58f1c22SToby Isaac       stratum = leafStrata[offset+s];
1453ad8374ffSToby Isaac       points[stratum][strataIdx[stratum]++] = p;
1454c58f1c22SToby Isaac     }
1455c58f1c22SToby Isaac   }
1456ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; s++) {
1457ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr);
1458ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr);
1459ad8374ffSToby Isaac   }
1460ad8374ffSToby Isaac   ierr = PetscFree(points);CHKERRQ(ierr);
1461e8f14785SLisandro Dalcin   ierr = PetscHSetIDestroy(&stratumHash);CHKERRQ(ierr);
1462c58f1c22SToby Isaac   ierr = PetscFree(leafStrata);CHKERRQ(ierr);
1463c58f1c22SToby Isaac   ierr = PetscFree(strataIdx);CHKERRQ(ierr);
1464c58f1c22SToby Isaac   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
1465c58f1c22SToby Isaac   PetscFunctionReturn(0);
1466c58f1c22SToby Isaac }
1467c58f1c22SToby Isaac 
14687937d9ceSMichael Lange /*@
14697937d9ceSMichael Lange   DMLabelGather - Gather all label values from leafs into roots
14707937d9ceSMichael Lange 
14717937d9ceSMichael Lange   Input Parameters:
14727937d9ceSMichael Lange + label - the DMLabel
147384f0b6dfSMatthew G. Knepley - sf - the Star Forest point communication map
14747937d9ceSMichael Lange 
147584f0b6dfSMatthew G. Knepley   Output Parameters:
147684f0b6dfSMatthew G. Knepley . labelNew - the new DMLabel with localised leaf values
14777937d9ceSMichael Lange 
14787937d9ceSMichael Lange   Level: developer
14797937d9ceSMichael Lange 
14807937d9ceSMichael Lange   Note: This is the inverse operation to DMLabelDistribute.
14817937d9ceSMichael Lange 
14827937d9ceSMichael Lange .seealso: DMLabelDistribute()
14837937d9ceSMichael Lange @*/
14847937d9ceSMichael Lange PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
14857937d9ceSMichael Lange {
14867937d9ceSMichael Lange   MPI_Comm       comm;
14877937d9ceSMichael Lange   PetscSection   rootSection;
14887937d9ceSMichael Lange   PetscSF        sfLabel;
14897937d9ceSMichael Lange   PetscSFNode   *rootPoints, *leafPoints;
14907937d9ceSMichael Lange   PetscInt       p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
14917937d9ceSMichael Lange   const PetscInt *rootDegree, *ilocal;
14927937d9ceSMichael Lange   PetscInt       *rootStrata;
1493d67d17b1SMatthew G. Knepley   const char    *lname;
14947937d9ceSMichael Lange   char          *name;
14957937d9ceSMichael Lange   PetscInt       nameSize;
14967937d9ceSMichael Lange   size_t         len = 0;
14979852e123SBarry Smith   PetscMPIInt    rank, size;
14987937d9ceSMichael Lange   PetscErrorCode ierr;
14997937d9ceSMichael Lange 
15007937d9ceSMichael Lange   PetscFunctionBegin;
1501d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1502d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
15037937d9ceSMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
15047937d9ceSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
15059852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
15067937d9ceSMichael Lange   /* Bcast name */
1507d67d17b1SMatthew G. Knepley   if (!rank) {
1508d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1509d67d17b1SMatthew G. Knepley     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1510d67d17b1SMatthew G. Knepley   }
15117937d9ceSMichael Lange   nameSize = len;
15127937d9ceSMichael Lange   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
15137937d9ceSMichael Lange   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1514d67d17b1SMatthew G. Knepley   if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);}
15157937d9ceSMichael Lange   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1516d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
15177937d9ceSMichael Lange   ierr = PetscFree(name);CHKERRQ(ierr);
15187937d9ceSMichael Lange   /* Gather rank/index pairs of leaves into local roots to build
15197937d9ceSMichael Lange      an inverse, multi-rooted SF. Note that this ignores local leaf
15207937d9ceSMichael Lange      indexing due to the use of the multiSF in PetscSFGather. */
15217937d9ceSMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr);
1522dc53bc9bSMatthew G. Knepley   ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr);
1523dc53bc9bSMatthew G. Knepley   for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
15247937d9ceSMichael Lange   for (p = 0; p < nleaves; p++) {
15258212dd46SStefano Zampini     PetscInt ilp = ilocal ? ilocal[p] : p;
15268212dd46SStefano Zampini 
15278212dd46SStefano Zampini     leafPoints[ilp].index = ilp;
15288212dd46SStefano Zampini     leafPoints[ilp].rank  = rank;
15297937d9ceSMichael Lange   }
15307937d9ceSMichael Lange   ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr);
15317937d9ceSMichael Lange   ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr);
15327937d9ceSMichael Lange   for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
15337937d9ceSMichael Lange   ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr);
15347937d9ceSMichael Lange   ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
15357937d9ceSMichael Lange   ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
15367937d9ceSMichael Lange   ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr);
15377937d9ceSMichael Lange   ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
15387937d9ceSMichael Lange   /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
15397937d9ceSMichael Lange   ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr);
15407937d9ceSMichael Lange   /* Rebuild the point strata on the receiver */
15417937d9ceSMichael Lange   for (p = 0, idx = 0; p < nroots; p++) {
15427937d9ceSMichael Lange     for (d = 0; d < rootDegree[p]; d++) {
15437937d9ceSMichael Lange       ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr);
15447937d9ceSMichael Lange       ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr);
15457937d9ceSMichael Lange       for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);}
15467937d9ceSMichael Lange     }
15477937d9ceSMichael Lange     idx += rootDegree[p];
15487937d9ceSMichael Lange   }
154977e0c0e7SMichael Lange   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
155077e0c0e7SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
155177e0c0e7SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
155277e0c0e7SMichael Lange   ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr);
15537937d9ceSMichael Lange   PetscFunctionReturn(0);
15547937d9ceSMichael Lange }
15557937d9ceSMichael Lange 
155684f0b6dfSMatthew G. Knepley /*@
155784f0b6dfSMatthew G. Knepley   DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label
155884f0b6dfSMatthew G. Knepley 
155984f0b6dfSMatthew G. Knepley   Input Parameter:
156084f0b6dfSMatthew G. Knepley . label - the DMLabel
156184f0b6dfSMatthew G. Knepley 
156284f0b6dfSMatthew G. Knepley   Output Parameters:
156384f0b6dfSMatthew G. Knepley + section - the section giving offsets for each stratum
156484f0b6dfSMatthew G. Knepley - is - An IS containing all the label points
156584f0b6dfSMatthew G. Knepley 
156684f0b6dfSMatthew G. Knepley   Level: developer
156784f0b6dfSMatthew G. Knepley 
156884f0b6dfSMatthew G. Knepley .seealso: DMLabelDistribute()
156984f0b6dfSMatthew G. Knepley @*/
1570c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1571c58f1c22SToby Isaac {
1572c58f1c22SToby Isaac   IS              vIS;
1573c58f1c22SToby Isaac   const PetscInt *values;
1574c58f1c22SToby Isaac   PetscInt       *points;
1575c58f1c22SToby Isaac   PetscInt        nV, vS = 0, vE = 0, v, N;
1576c58f1c22SToby Isaac   PetscErrorCode  ierr;
1577c58f1c22SToby Isaac 
1578c58f1c22SToby Isaac   PetscFunctionBegin;
1579d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1580c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr);
1581c58f1c22SToby Isaac   ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
1582c58f1c22SToby Isaac   ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr);
1583c58f1c22SToby Isaac   if (nV) {vS = values[0]; vE = values[0]+1;}
1584c58f1c22SToby Isaac   for (v = 1; v < nV; ++v) {
1585c58f1c22SToby Isaac     vS = PetscMin(vS, values[v]);
1586c58f1c22SToby Isaac     vE = PetscMax(vE, values[v]+1);
1587c58f1c22SToby Isaac   }
1588c58f1c22SToby Isaac   ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr);
1589c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr);
1590c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1591c58f1c22SToby Isaac     PetscInt n;
1592c58f1c22SToby Isaac 
1593c58f1c22SToby Isaac     ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr);
1594c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr);
1595c58f1c22SToby Isaac   }
1596c58f1c22SToby Isaac   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
1597c58f1c22SToby Isaac   ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr);
1598c58f1c22SToby Isaac   ierr = PetscMalloc1(N, &points);CHKERRQ(ierr);
1599c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1600c58f1c22SToby Isaac     IS              is;
1601c58f1c22SToby Isaac     const PetscInt *spoints;
1602c58f1c22SToby Isaac     PetscInt        dof, off, p;
1603c58f1c22SToby Isaac 
1604c58f1c22SToby Isaac     ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr);
1605c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr);
1606c58f1c22SToby Isaac     ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr);
1607c58f1c22SToby Isaac     ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr);
1608c58f1c22SToby Isaac     for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1609c58f1c22SToby Isaac     ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr);
1610c58f1c22SToby Isaac     ierr = ISDestroy(&is);CHKERRQ(ierr);
1611c58f1c22SToby Isaac   }
1612c58f1c22SToby Isaac   ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr);
1613c58f1c22SToby Isaac   ierr = ISDestroy(&vIS);CHKERRQ(ierr);
1614c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1615c58f1c22SToby Isaac   PetscFunctionReturn(0);
1616c58f1c22SToby Isaac }
1617c58f1c22SToby Isaac 
161884f0b6dfSMatthew G. Knepley /*@
1619c58f1c22SToby Isaac   PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1620c58f1c22SToby Isaac   the local section and an SF describing the section point overlap.
1621c58f1c22SToby Isaac 
1622c58f1c22SToby Isaac   Input Parameters:
1623c58f1c22SToby Isaac   + s - The PetscSection for the local field layout
1624c58f1c22SToby Isaac   . sf - The SF describing parallel layout of the section points
1625c58f1c22SToby Isaac   . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1626c58f1c22SToby Isaac   . label - The label specifying the points
1627c58f1c22SToby Isaac   - labelValue - The label stratum specifying the points
1628c58f1c22SToby Isaac 
1629c58f1c22SToby Isaac   Output Parameter:
1630c58f1c22SToby Isaac   . gsection - The PetscSection for the global field layout
1631c58f1c22SToby Isaac 
1632c58f1c22SToby Isaac   Note: This gives negative sizes and offsets to points not owned by this process
1633c58f1c22SToby Isaac 
1634c58f1c22SToby Isaac   Level: developer
1635c58f1c22SToby Isaac 
1636c58f1c22SToby Isaac .seealso: PetscSectionCreate()
1637c58f1c22SToby Isaac @*/
1638c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1639c58f1c22SToby Isaac {
1640c58f1c22SToby Isaac   PetscInt      *neg = NULL, *tmpOff = NULL;
1641c58f1c22SToby Isaac   PetscInt       pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1642c58f1c22SToby Isaac   PetscErrorCode ierr;
1643c58f1c22SToby Isaac 
1644c58f1c22SToby Isaac   PetscFunctionBegin;
1645d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1646d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1647d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4);
1648c58f1c22SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr);
1649c58f1c22SToby Isaac   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1650c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr);
1651c58f1c22SToby Isaac   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1652c58f1c22SToby Isaac   if (nroots >= 0) {
1653c58f1c22SToby Isaac     if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1654c58f1c22SToby Isaac     ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr);
1655c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1656c58f1c22SToby Isaac       ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr);
1657c58f1c22SToby Isaac     } else {
1658c58f1c22SToby Isaac       tmpOff = &(*gsection)->atlasDof[-pStart];
1659c58f1c22SToby Isaac     }
1660c58f1c22SToby Isaac   }
1661c58f1c22SToby Isaac   /* Mark ghost points with negative dof */
1662c58f1c22SToby Isaac   for (p = pStart; p < pEnd; ++p) {
1663c58f1c22SToby Isaac     PetscInt value;
1664c58f1c22SToby Isaac 
1665c58f1c22SToby Isaac     ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr);
1666c58f1c22SToby Isaac     if (value != labelValue) continue;
1667c58f1c22SToby Isaac     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1668c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr);
1669c58f1c22SToby Isaac     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1670c58f1c22SToby Isaac     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);}
1671c58f1c22SToby Isaac     if (neg) neg[p] = -(dof+1);
1672c58f1c22SToby Isaac   }
1673c58f1c22SToby Isaac   ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr);
1674c58f1c22SToby Isaac   if (nroots >= 0) {
1675c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1676c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1677c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1678c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1679c58f1c22SToby Isaac     }
1680c58f1c22SToby Isaac   }
1681c58f1c22SToby Isaac   /* Calculate new sizes, get proccess offset, and calculate point offsets */
1682c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1683c58f1c22SToby Isaac     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1684c58f1c22SToby Isaac     (*gsection)->atlasOff[p] = off;
1685c58f1c22SToby Isaac     off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1686c58f1c22SToby Isaac   }
1687c58f1c22SToby Isaac   ierr       = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1688c58f1c22SToby Isaac   globalOff -= off;
1689c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1690c58f1c22SToby Isaac     (*gsection)->atlasOff[p] += globalOff;
1691c58f1c22SToby Isaac     if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1692c58f1c22SToby Isaac   }
1693c58f1c22SToby Isaac   /* Put in negative offsets for ghost points */
1694c58f1c22SToby Isaac   if (nroots >= 0) {
1695c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1696c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1697c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1698c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1699c58f1c22SToby Isaac     }
1700c58f1c22SToby Isaac   }
1701c58f1c22SToby Isaac   if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);}
1702c58f1c22SToby Isaac   ierr = PetscFree(neg);CHKERRQ(ierr);
1703c58f1c22SToby Isaac   PetscFunctionReturn(0);
1704c58f1c22SToby Isaac }
1705c58f1c22SToby Isaac 
17065fdea053SToby Isaac typedef struct _n_PetscSectionSym_Label
17075fdea053SToby Isaac {
17085fdea053SToby Isaac   DMLabel           label;
17095fdea053SToby Isaac   PetscCopyMode     *modes;
17105fdea053SToby Isaac   PetscInt          *sizes;
17115fdea053SToby Isaac   const PetscInt    ***perms;
17125fdea053SToby Isaac   const PetscScalar ***rots;
17135fdea053SToby Isaac   PetscInt          (*minMaxOrients)[2];
17145fdea053SToby Isaac   PetscInt          numStrata; /* numStrata is only increasing, functions as a state */
17155fdea053SToby Isaac } PetscSectionSym_Label;
17165fdea053SToby Isaac 
17175fdea053SToby Isaac static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
17185fdea053SToby Isaac {
17195fdea053SToby Isaac   PetscInt              i, j;
17205fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
17215fdea053SToby Isaac   PetscErrorCode        ierr;
17225fdea053SToby Isaac 
17235fdea053SToby Isaac   PetscFunctionBegin;
17245fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
17255fdea053SToby Isaac     if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
17265fdea053SToby Isaac       for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
17275fdea053SToby Isaac         if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);}
17285fdea053SToby Isaac         if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);}
17295fdea053SToby Isaac       }
17305fdea053SToby Isaac       if (sl->perms[i]) {
17315fdea053SToby Isaac         const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
17325fdea053SToby Isaac 
17335fdea053SToby Isaac         ierr = PetscFree(perms);CHKERRQ(ierr);
17345fdea053SToby Isaac       }
17355fdea053SToby Isaac       if (sl->rots[i]) {
17365fdea053SToby Isaac         const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
17375fdea053SToby Isaac 
17385fdea053SToby Isaac         ierr = PetscFree(rots);CHKERRQ(ierr);
17395fdea053SToby Isaac       }
17405fdea053SToby Isaac     }
17415fdea053SToby Isaac   }
17425fdea053SToby Isaac   ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr);
17435fdea053SToby Isaac   ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr);
17445fdea053SToby Isaac   sl->numStrata = 0;
17455fdea053SToby Isaac   PetscFunctionReturn(0);
17465fdea053SToby Isaac }
17475fdea053SToby Isaac 
17485fdea053SToby Isaac static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
17495fdea053SToby Isaac {
17505fdea053SToby Isaac   PetscErrorCode ierr;
17515fdea053SToby Isaac 
17525fdea053SToby Isaac   PetscFunctionBegin;
17535fdea053SToby Isaac   ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);
17545fdea053SToby Isaac   ierr = PetscFree(sym->data);CHKERRQ(ierr);
17555fdea053SToby Isaac   PetscFunctionReturn(0);
17565fdea053SToby Isaac }
17575fdea053SToby Isaac 
17585fdea053SToby Isaac static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
17595fdea053SToby Isaac {
17605fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
17615fdea053SToby Isaac   PetscBool             isAscii;
17625fdea053SToby Isaac   DMLabel               label = sl->label;
1763d67d17b1SMatthew G. Knepley   const char           *name;
17645fdea053SToby Isaac   PetscErrorCode        ierr;
17655fdea053SToby Isaac 
17665fdea053SToby Isaac   PetscFunctionBegin;
17675fdea053SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
17685fdea053SToby Isaac   if (isAscii) {
17695fdea053SToby Isaac     PetscInt          i, j, k;
17705fdea053SToby Isaac     PetscViewerFormat format;
17715fdea053SToby Isaac 
17725fdea053SToby Isaac     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
17735fdea053SToby Isaac     if (label) {
17745fdea053SToby Isaac       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
17755fdea053SToby Isaac       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
17765fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
17775fdea053SToby Isaac         ierr = DMLabelView(label, viewer);CHKERRQ(ierr);
17785fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
17795fdea053SToby Isaac       } else {
1780d67d17b1SMatthew G. Knepley         ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1781d67d17b1SMatthew G. Knepley         ierr = PetscViewerASCIIPrintf(viewer,"  Label '%s'\n",name);CHKERRQ(ierr);
17825fdea053SToby Isaac       }
17835fdea053SToby Isaac     } else {
17845fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr);
17855fdea053SToby Isaac     }
17865fdea053SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
17875fdea053SToby Isaac     for (i = 0; i <= sl->numStrata; i++) {
17885fdea053SToby Isaac       PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
17895fdea053SToby Isaac 
17905fdea053SToby Isaac       if (!(sl->perms[i] || sl->rots[i])) {
17915fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr);
17925fdea053SToby Isaac       } else {
17935fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr);
17945fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
17955fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr);
17965fdea053SToby Isaac         if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
17975fdea053SToby Isaac           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
17985fdea053SToby Isaac           for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
17995fdea053SToby Isaac             if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
18005fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr);
18015fdea053SToby Isaac             } else {
18025fdea053SToby Isaac               PetscInt tab;
18035fdea053SToby Isaac 
18045fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr);
18055fdea053SToby Isaac               ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18065fdea053SToby Isaac               ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr);
18075fdea053SToby Isaac               if (sl->perms[i] && sl->perms[i][j]) {
18085fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr);
18095fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
18105fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);}
18115fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
18125fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
18135fdea053SToby Isaac               }
18145fdea053SToby Isaac               if (sl->rots[i] && sl->rots[i][j]) {
18155fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Rotations:  ");CHKERRQ(ierr);
18165fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
18175fdea053SToby Isaac #if defined(PETSC_USE_COMPLEX)
18185fdea053SToby 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);}
18195fdea053SToby Isaac #else
18205fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);}
18215fdea053SToby Isaac #endif
18225fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
18235fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
18245fdea053SToby Isaac               }
18255fdea053SToby Isaac               ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
18265fdea053SToby Isaac             }
18275fdea053SToby Isaac           }
18285fdea053SToby Isaac           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
18295fdea053SToby Isaac         }
18305fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
18315fdea053SToby Isaac       }
18325fdea053SToby Isaac     }
18335fdea053SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
18345fdea053SToby Isaac   }
18355fdea053SToby Isaac   PetscFunctionReturn(0);
18365fdea053SToby Isaac }
18375fdea053SToby Isaac 
18385fdea053SToby Isaac /*@
18395fdea053SToby Isaac   PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
18405fdea053SToby Isaac 
18415fdea053SToby Isaac   Logically collective on sym
18425fdea053SToby Isaac 
18435fdea053SToby Isaac   Input parameters:
18445fdea053SToby Isaac + sym - the section symmetries
18455fdea053SToby Isaac - label - the DMLabel describing the types of points
18465fdea053SToby Isaac 
18475fdea053SToby Isaac   Level: developer:
18485fdea053SToby Isaac 
18495fdea053SToby Isaac .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()
18505fdea053SToby Isaac @*/
18515fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
18525fdea053SToby Isaac {
18535fdea053SToby Isaac   PetscSectionSym_Label *sl;
18545fdea053SToby Isaac   PetscErrorCode        ierr;
18555fdea053SToby Isaac 
18565fdea053SToby Isaac   PetscFunctionBegin;
18575fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
18585fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
18595fdea053SToby Isaac   if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);}
18605fdea053SToby Isaac   if (label) {
18615fdea053SToby Isaac     sl->label = label;
1862d67d17b1SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) label);CHKERRQ(ierr);
18635fdea053SToby Isaac     ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr);
18641a834cf9SToby 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);
18651a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr);
18661a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr);
18671a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr);
18681a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr);
18691a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr);
18705fdea053SToby Isaac   }
18715fdea053SToby Isaac   PetscFunctionReturn(0);
18725fdea053SToby Isaac }
18735fdea053SToby Isaac 
18745fdea053SToby Isaac /*@C
18755fdea053SToby Isaac   PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
18765fdea053SToby Isaac 
18775fdea053SToby Isaac   Logically collective on PetscSectionSym
18785fdea053SToby Isaac 
18795fdea053SToby Isaac   InputParameters:
18805fdea053SToby Isaac + sys       - the section symmetries
18815fdea053SToby Isaac . stratum   - the stratum value in the label that we are assigning symmetries for
18825fdea053SToby Isaac . size      - the number of dofs for points in the stratum of the label
18835fdea053SToby Isaac . minOrient - the smallest orientation for a point in this stratum
18845fdea053SToby Isaac . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
18855fdea053SToby Isaac . mode      - how sym should copy the perms and rots arrays
18865fdea053SToby Isaac . perms     - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation.  A NULL permutation is the identity
18875fdea053SToby Isaac + rots      - NULL if there are no rotations, or (maxOrient - minOrient) sets of rotations, one for each orientation.  A NULL set of orientations is the identity
18885fdea053SToby Isaac 
18895fdea053SToby Isaac   Level: developer
18905fdea053SToby Isaac 
18915fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()
18925fdea053SToby Isaac @*/
18935fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
18945fdea053SToby Isaac {
18955fdea053SToby Isaac   PetscSectionSym_Label *sl;
1896d67d17b1SMatthew G. Knepley   const char            *name;
1897d67d17b1SMatthew G. Knepley   PetscInt               i, j, k;
18985fdea053SToby Isaac   PetscErrorCode         ierr;
18995fdea053SToby Isaac 
19005fdea053SToby Isaac   PetscFunctionBegin;
19015fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
19025fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
19035fdea053SToby Isaac   if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet");
19045fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
19055fdea053SToby Isaac     PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
19065fdea053SToby Isaac 
19075fdea053SToby Isaac     if (stratum == value) break;
19085fdea053SToby Isaac   }
1909d67d17b1SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1910d67d17b1SMatthew G. Knepley   if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,name);
19115fdea053SToby Isaac   sl->sizes[i] = size;
19125fdea053SToby Isaac   sl->modes[i] = mode;
19135fdea053SToby Isaac   sl->minMaxOrients[i][0] = minOrient;
19145fdea053SToby Isaac   sl->minMaxOrients[i][1] = maxOrient;
19155fdea053SToby Isaac   if (mode == PETSC_COPY_VALUES) {
19165fdea053SToby Isaac     if (perms) {
19175fdea053SToby Isaac       PetscInt    **ownPerms;
19185fdea053SToby Isaac 
19195fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr);
19205fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
19215fdea053SToby Isaac         if (perms[j]) {
19225fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr);
19235fdea053SToby Isaac           for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];}
19245fdea053SToby Isaac         }
19255fdea053SToby Isaac       }
19265fdea053SToby Isaac       sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient];
19275fdea053SToby Isaac     }
19285fdea053SToby Isaac     if (rots) {
19295fdea053SToby Isaac       PetscScalar **ownRots;
19305fdea053SToby Isaac 
19315fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr);
19325fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
19335fdea053SToby Isaac         if (rots[j]) {
19345fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr);
19355fdea053SToby Isaac           for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];}
19365fdea053SToby Isaac         }
19375fdea053SToby Isaac       }
19385fdea053SToby Isaac       sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient];
19395fdea053SToby Isaac     }
19405fdea053SToby Isaac   } else {
19415fdea053SToby Isaac     sl->perms[i] = perms ? &perms[-minOrient] : NULL;
19425fdea053SToby Isaac     sl->rots[i]  = rots ? &rots[-minOrient] : NULL;
19435fdea053SToby Isaac   }
19445fdea053SToby Isaac   PetscFunctionReturn(0);
19455fdea053SToby Isaac }
19465fdea053SToby Isaac 
19475fdea053SToby Isaac static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
19485fdea053SToby Isaac {
19495fdea053SToby Isaac   PetscInt              i, j, numStrata;
19505fdea053SToby Isaac   PetscSectionSym_Label *sl;
19515fdea053SToby Isaac   DMLabel               label;
19525fdea053SToby Isaac   PetscErrorCode        ierr;
19535fdea053SToby Isaac 
19545fdea053SToby Isaac   PetscFunctionBegin;
19555fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
19565fdea053SToby Isaac   numStrata = sl->numStrata;
19575fdea053SToby Isaac   label     = sl->label;
19585fdea053SToby Isaac   for (i = 0; i < numPoints; i++) {
19595fdea053SToby Isaac     PetscInt point = points[2*i];
19605fdea053SToby Isaac     PetscInt ornt  = points[2*i+1];
19615fdea053SToby Isaac 
19625fdea053SToby Isaac     for (j = 0; j < numStrata; j++) {
19635fdea053SToby Isaac       if (label->validIS[j]) {
19645fdea053SToby Isaac         PetscInt k;
19655fdea053SToby Isaac 
19665fdea053SToby Isaac         ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr);
19675fdea053SToby Isaac         if (k >= 0) break;
19685fdea053SToby Isaac       } else {
19695fdea053SToby Isaac         PetscBool has;
19705fdea053SToby Isaac 
1971b9907514SLisandro Dalcin         ierr = PetscHSetIHas(label->ht[j], point, &has);CHKERRQ(ierr);
19725fdea053SToby Isaac         if (has) break;
19735fdea053SToby Isaac       }
19745fdea053SToby Isaac     }
19755fdea053SToby 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);
19765fdea053SToby Isaac     if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;}
19775fdea053SToby Isaac     if (rots) {rots[i]  = sl->rots[j] ? sl->rots[j][ornt] : NULL;}
19785fdea053SToby Isaac   }
19795fdea053SToby Isaac   PetscFunctionReturn(0);
19805fdea053SToby Isaac }
19815fdea053SToby Isaac 
19825fdea053SToby Isaac PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
19835fdea053SToby Isaac {
19845fdea053SToby Isaac   PetscSectionSym_Label *sl;
19855fdea053SToby Isaac   PetscErrorCode        ierr;
19865fdea053SToby Isaac 
19875fdea053SToby Isaac   PetscFunctionBegin;
19885fdea053SToby Isaac   ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr);
19895fdea053SToby Isaac   sym->ops->getpoints = PetscSectionSymGetPoints_Label;
19905fdea053SToby Isaac   sym->ops->view      = PetscSectionSymView_Label;
19915fdea053SToby Isaac   sym->ops->destroy   = PetscSectionSymDestroy_Label;
19925fdea053SToby Isaac   sym->data           = (void *) sl;
19935fdea053SToby Isaac   PetscFunctionReturn(0);
19945fdea053SToby Isaac }
19955fdea053SToby Isaac 
19965fdea053SToby Isaac /*@
19975fdea053SToby Isaac   PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
19985fdea053SToby Isaac 
19995fdea053SToby Isaac   Collective
20005fdea053SToby Isaac 
20015fdea053SToby Isaac   Input Parameters:
20025fdea053SToby Isaac + comm - the MPI communicator for the new symmetry
20035fdea053SToby Isaac - label - the label defining the strata
20045fdea053SToby Isaac 
20055fdea053SToby Isaac   Output Parameters:
20065fdea053SToby Isaac . sym - the section symmetries
20075fdea053SToby Isaac 
20085fdea053SToby Isaac   Level: developer
20095fdea053SToby Isaac 
20105fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()
20115fdea053SToby Isaac @*/
20125fdea053SToby Isaac PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
20135fdea053SToby Isaac {
20145fdea053SToby Isaac   PetscErrorCode        ierr;
20155fdea053SToby Isaac 
20165fdea053SToby Isaac   PetscFunctionBegin;
20175fdea053SToby Isaac   ierr = DMInitializePackage();CHKERRQ(ierr);
20185fdea053SToby Isaac   ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr);
20195fdea053SToby Isaac   ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr);
20205fdea053SToby Isaac   ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr);
20215fdea053SToby Isaac   PetscFunctionReturn(0);
20225fdea053SToby Isaac }
2023