xref: /petsc/src/dm/label/dmlabel.c (revision b9907514337d2f9686ada4b680cea2bfa7a6c371)
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;
40*b9907514SLisandro 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 {
61*b9907514SLisandro Dalcin   PetscInt       off = 0, *pointArray, p;
62c58f1c22SToby Isaac   PetscErrorCode ierr;
63c58f1c22SToby Isaac 
64*b9907514SLisandro 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);
70*b9907514SLisandro 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 
131*b9907514SLisandro 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 
146*b9907514SLisandro Dalcin #if !defined(DMLABEL_LOOKUP_THRESHOLD)
147*b9907514SLisandro Dalcin #define DMLABEL_LOOKUP_THRESHOLD 16
148*b9907514SLisandro Dalcin #endif
149*b9907514SLisandro Dalcin 
1500c3c4a36SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index)
1510c3c4a36SLisandro Dalcin {
1520c3c4a36SLisandro Dalcin   PetscInt       v;
153*b9907514SLisandro Dalcin   PetscErrorCode ierr;
1540e79e033SBarry Smith 
1550c3c4a36SLisandro Dalcin   PetscFunctionBegin;
1560e79e033SBarry Smith   *index = -1;
157*b9907514SLisandro Dalcin   if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) {
158*b9907514SLisandro Dalcin     for (v = 0; v < label->numStrata; ++v)
159*b9907514SLisandro Dalcin       if (label->stratumValues[v] == value) {*index = v; break;}
160*b9907514SLisandro Dalcin   } else {
161*b9907514SLisandro Dalcin     ierr = PetscHMapIGet(label->hmap, value, index);CHKERRQ(ierr);
162*b9907514SLisandro Dalcin #if defined(PETSC_USE_DEBUG)
163*b9907514SLisandro Dalcin     v = *index;
164*b9907514SLisandro Dalcin     if (v >= 0 && v >= label->numStrata) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent lookup using strata hash map");
165*b9907514SLisandro Dalcin     if (v >= 0 && label->stratumValues[v] != value) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent lookup using strata hash map");
166*b9907514SLisandro Dalcin #endif
1670e79e033SBarry Smith   }
1680c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
1690c3c4a36SLisandro Dalcin }
1700c3c4a36SLisandro Dalcin 
171*b9907514SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index)
172c58f1c22SToby Isaac {
173*b9907514SLisandro Dalcin   PetscInt       v;
174*b9907514SLisandro Dalcin   PetscInt      *tmpV;
175*b9907514SLisandro Dalcin   PetscInt      *tmpS;
176*b9907514SLisandro Dalcin   PetscHSetI    *tmpH, ht;
177*b9907514SLisandro Dalcin   IS            *tmpP, is;
178c58f1c22SToby Isaac   PetscBool     *tmpB;
179*b9907514SLisandro Dalcin   PetscHMapI     hmap = label->hmap;
180c58f1c22SToby Isaac   PetscErrorCode ierr;
181c58f1c22SToby Isaac 
182c58f1c22SToby Isaac   PetscFunctionBegin;
183*b9907514SLisandro Dalcin   v    = label->numStrata;
184*b9907514SLisandro Dalcin   tmpV = label->stratumValues;
185*b9907514SLisandro Dalcin   tmpS = label->stratumSizes;
186*b9907514SLisandro Dalcin   tmpH = label->ht;
187*b9907514SLisandro Dalcin   tmpP = label->points;
188*b9907514SLisandro Dalcin   tmpB = label->validIS;
189*b9907514SLisandro Dalcin   ierr = PetscRealloc((v+1)*sizeof(*tmpV), &tmpV);CHKERRQ(ierr);
190*b9907514SLisandro Dalcin   ierr = PetscRealloc((v+1)*sizeof(*tmpS), &tmpS);CHKERRQ(ierr);
191*b9907514SLisandro Dalcin   ierr = PetscRealloc((v+1)*sizeof(*tmpH), &tmpH);CHKERRQ(ierr);
192*b9907514SLisandro Dalcin   ierr = PetscRealloc((v+1)*sizeof(*tmpP), &tmpP);CHKERRQ(ierr);
193*b9907514SLisandro Dalcin   ierr = PetscRealloc((v+1)*sizeof(*tmpB), &tmpB);CHKERRQ(ierr);
194*b9907514SLisandro Dalcin   label->numStrata     = v+1;
195c58f1c22SToby Isaac   label->stratumValues = tmpV;
196c58f1c22SToby Isaac   label->stratumSizes  = tmpS;
197c58f1c22SToby Isaac   label->ht            = tmpH;
198c58f1c22SToby Isaac   label->points        = tmpP;
199ad8374ffSToby Isaac   label->validIS       = tmpB;
200*b9907514SLisandro Dalcin   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
201*b9907514SLisandro Dalcin   ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);CHKERRQ(ierr);
202*b9907514SLisandro Dalcin   ierr = PetscHMapISet(hmap, value, v);CHKERRQ(ierr);
203*b9907514SLisandro Dalcin   tmpV[v] = value;
204*b9907514SLisandro Dalcin   tmpS[v] = 0;
205*b9907514SLisandro Dalcin   tmpH[v] = ht;
206*b9907514SLisandro Dalcin   tmpP[v] = is;
207*b9907514SLisandro Dalcin   tmpB[v] = PETSC_TRUE;
2080c3c4a36SLisandro Dalcin   *index = v;
2090c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
2100c3c4a36SLisandro Dalcin }
2110c3c4a36SLisandro Dalcin 
212*b9907514SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupAddStratum(DMLabel label, PetscInt value, PetscInt *index)
213*b9907514SLisandro Dalcin {
214*b9907514SLisandro Dalcin   PetscErrorCode ierr;
215*b9907514SLisandro Dalcin   PetscFunctionBegin;
216*b9907514SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, index);CHKERRQ(ierr);
217*b9907514SLisandro Dalcin   if (*index < 0) {ierr = DMLabelNewStratum(label, value, index);CHKERRQ(ierr);}
218*b9907514SLisandro Dalcin   PetscFunctionReturn(0);
219*b9907514SLisandro Dalcin }
220*b9907514SLisandro Dalcin 
221*b9907514SLisandro Dalcin /*@
222*b9907514SLisandro Dalcin   DMLabelAddStratum - Adds a new stratum value in a DMLabel
223*b9907514SLisandro Dalcin 
224*b9907514SLisandro Dalcin   Input Parameter:
225*b9907514SLisandro Dalcin + label - The DMLabel
226*b9907514SLisandro Dalcin - value - The stratum value
227*b9907514SLisandro Dalcin 
228*b9907514SLisandro Dalcin   Level: beginner
229*b9907514SLisandro Dalcin 
230*b9907514SLisandro Dalcin .seealso:  DMLabelCreate(), DMLabelDestroy()
231*b9907514SLisandro Dalcin @*/
2320c3c4a36SLisandro Dalcin PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
2330c3c4a36SLisandro Dalcin {
2340c3c4a36SLisandro Dalcin   PetscInt       v;
2350c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
2360c3c4a36SLisandro Dalcin 
2370c3c4a36SLisandro Dalcin   PetscFunctionBegin;
238d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
239*b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
240*b9907514SLisandro Dalcin   PetscFunctionReturn(0);
241*b9907514SLisandro Dalcin }
242*b9907514SLisandro Dalcin 
243*b9907514SLisandro Dalcin /*@
244*b9907514SLisandro Dalcin   DMLabelAddStrata - Adds new stratum values in a DMLabel
245*b9907514SLisandro Dalcin 
246*b9907514SLisandro Dalcin   Input Parameter:
247*b9907514SLisandro Dalcin + label - The DMLabel
248*b9907514SLisandro Dalcin . numStrata - The number of stratum values
249*b9907514SLisandro Dalcin - stratumValues - The stratum values
250*b9907514SLisandro Dalcin 
251*b9907514SLisandro Dalcin   Level: beginner
252*b9907514SLisandro Dalcin 
253*b9907514SLisandro Dalcin .seealso:  DMLabelCreate(), DMLabelDestroy()
254*b9907514SLisandro Dalcin @*/
255*b9907514SLisandro Dalcin PetscErrorCode DMLabelAddStrata(DMLabel label, PetscInt numStrata, const PetscInt stratumValues[])
256*b9907514SLisandro Dalcin {
257*b9907514SLisandro Dalcin   PetscInt       *values, v;
258*b9907514SLisandro Dalcin   PetscErrorCode ierr;
259*b9907514SLisandro Dalcin 
260*b9907514SLisandro Dalcin   PetscFunctionBegin;
261*b9907514SLisandro Dalcin   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
262*b9907514SLisandro Dalcin   if (numStrata) PetscValidIntPointer(stratumValues, 3);
263*b9907514SLisandro Dalcin   ierr = PetscMalloc1(numStrata, &values);CHKERRQ(ierr);
264*b9907514SLisandro Dalcin   ierr = PetscMemcpy(values, stratumValues, numStrata*sizeof(PetscInt));CHKERRQ(ierr);
265*b9907514SLisandro Dalcin   ierr = PetscSortRemoveDupsInt(&numStrata, values);CHKERRQ(ierr);
266*b9907514SLisandro Dalcin   if (!label->numStrata) { /* Fast preallocation */
267*b9907514SLisandro Dalcin     PetscInt   *tmpV;
268*b9907514SLisandro Dalcin     PetscInt   *tmpS;
269*b9907514SLisandro Dalcin     PetscHSetI *tmpH, ht;
270*b9907514SLisandro Dalcin     IS         *tmpP, is;
271*b9907514SLisandro Dalcin     PetscBool  *tmpB;
272*b9907514SLisandro Dalcin     PetscHMapI  hmap = label->hmap;
273*b9907514SLisandro Dalcin 
274*b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpV);CHKERRQ(ierr);
275*b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpS);CHKERRQ(ierr);
276*b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpH);CHKERRQ(ierr);
277*b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpP);CHKERRQ(ierr);
278*b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpB);CHKERRQ(ierr);
279*b9907514SLisandro Dalcin     label->numStrata     = numStrata;
280*b9907514SLisandro Dalcin     label->stratumValues = tmpV;
281*b9907514SLisandro Dalcin     label->stratumSizes  = tmpS;
282*b9907514SLisandro Dalcin     label->ht            = tmpH;
283*b9907514SLisandro Dalcin     label->points        = tmpP;
284*b9907514SLisandro Dalcin     label->validIS       = tmpB;
285*b9907514SLisandro Dalcin     for (v = 0; v < numStrata; ++v) {
286*b9907514SLisandro Dalcin       ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
287*b9907514SLisandro Dalcin       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);CHKERRQ(ierr);
288*b9907514SLisandro Dalcin       ierr = PetscHMapISet(hmap, values[v], v);CHKERRQ(ierr);
289*b9907514SLisandro Dalcin       tmpV[v] = values[v];
290*b9907514SLisandro Dalcin       tmpS[v] = 0;
291*b9907514SLisandro Dalcin       tmpH[v] = ht;
292*b9907514SLisandro Dalcin       tmpP[v] = is;
293*b9907514SLisandro Dalcin       tmpB[v] = PETSC_TRUE;
294*b9907514SLisandro Dalcin     }
295*b9907514SLisandro Dalcin   } else {
296*b9907514SLisandro Dalcin     for (v = 0; v < numStrata; ++v) {
297*b9907514SLisandro Dalcin       ierr = DMLabelAddStratum(label, values[v]);CHKERRQ(ierr);
298*b9907514SLisandro Dalcin     }
299*b9907514SLisandro Dalcin   }
300*b9907514SLisandro Dalcin   ierr = PetscFree(values);CHKERRQ(ierr);
301*b9907514SLisandro Dalcin   PetscFunctionReturn(0);
302*b9907514SLisandro Dalcin }
303*b9907514SLisandro Dalcin 
304*b9907514SLisandro Dalcin /*@
305*b9907514SLisandro Dalcin   DMLabelAddStrataIS - Adds new stratum values in a DMLabel
306*b9907514SLisandro Dalcin 
307*b9907514SLisandro Dalcin   Input Parameter:
308*b9907514SLisandro Dalcin + label - The DMLabel
309*b9907514SLisandro Dalcin - valueIS - Index set with stratum values
310*b9907514SLisandro Dalcin 
311*b9907514SLisandro Dalcin   Level: beginner
312*b9907514SLisandro Dalcin 
313*b9907514SLisandro Dalcin .seealso:  DMLabelCreate(), DMLabelDestroy()
314*b9907514SLisandro Dalcin @*/
315*b9907514SLisandro Dalcin PetscErrorCode DMLabelAddStrataIS(DMLabel label, IS valueIS)
316*b9907514SLisandro Dalcin {
317*b9907514SLisandro Dalcin   PetscInt       numStrata;
318*b9907514SLisandro Dalcin   const PetscInt *stratumValues;
319*b9907514SLisandro Dalcin   PetscErrorCode ierr;
320*b9907514SLisandro Dalcin 
321*b9907514SLisandro Dalcin   PetscFunctionBegin;
322*b9907514SLisandro Dalcin   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
323*b9907514SLisandro Dalcin   PetscValidHeaderSpecific(valueIS, IS_CLASSID, 2);
324*b9907514SLisandro Dalcin   ierr = ISGetLocalSize(valueIS, &numStrata);CHKERRQ(ierr);
325*b9907514SLisandro Dalcin   ierr = ISGetIndices(valueIS, &stratumValues);CHKERRQ(ierr);
326*b9907514SLisandro Dalcin   ierr = DMLabelAddStrata(label, numStrata, stratumValues);CHKERRQ(ierr);
327c58f1c22SToby Isaac   PetscFunctionReturn(0);
328c58f1c22SToby Isaac }
329c58f1c22SToby Isaac 
330c58f1c22SToby Isaac static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
331c58f1c22SToby Isaac {
332c58f1c22SToby Isaac   PetscInt       v;
333c58f1c22SToby Isaac   PetscMPIInt    rank;
334c58f1c22SToby Isaac   PetscErrorCode ierr;
335c58f1c22SToby Isaac 
336c58f1c22SToby Isaac   PetscFunctionBegin;
337c58f1c22SToby Isaac   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr);
338c58f1c22SToby Isaac   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
339c58f1c22SToby Isaac   if (label) {
340d67d17b1SMatthew G. Knepley     const char *name;
341d67d17b1SMatthew G. Knepley 
342d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
343d67d17b1SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", name);CHKERRQ(ierr);
344c58f1c22SToby Isaac     if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);}
345c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
346c58f1c22SToby Isaac       const PetscInt value = label->stratumValues[v];
347ad8374ffSToby Isaac       const PetscInt *points;
348c58f1c22SToby Isaac       PetscInt       p;
349c58f1c22SToby Isaac 
350ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
351c58f1c22SToby Isaac       for (p = 0; p < label->stratumSizes[v]; ++p) {
352ad8374ffSToby Isaac         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr);
353c58f1c22SToby Isaac       }
354ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
355c58f1c22SToby Isaac     }
356c58f1c22SToby Isaac   }
357c58f1c22SToby Isaac   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
358c58f1c22SToby Isaac   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
359c58f1c22SToby Isaac   PetscFunctionReturn(0);
360c58f1c22SToby Isaac }
361c58f1c22SToby Isaac 
362c58f1c22SToby Isaac /*@C
363c58f1c22SToby Isaac   DMLabelView - View the label
364c58f1c22SToby Isaac 
365c58f1c22SToby Isaac   Input Parameters:
366c58f1c22SToby Isaac + label - The DMLabel
367c58f1c22SToby Isaac - viewer - The PetscViewer
368c58f1c22SToby Isaac 
369c58f1c22SToby Isaac   Level: intermediate
370c58f1c22SToby Isaac 
371c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy()
372c58f1c22SToby Isaac @*/
373c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
374c58f1c22SToby Isaac {
375c58f1c22SToby Isaac   PetscBool      iascii;
376c58f1c22SToby Isaac   PetscErrorCode ierr;
377c58f1c22SToby Isaac 
378c58f1c22SToby Isaac   PetscFunctionBegin;
379d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
380*b9907514SLisandro Dalcin   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)label), &viewer);CHKERRQ(ierr);}
381c58f1c22SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
382c58f1c22SToby Isaac   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
383c58f1c22SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
384c58f1c22SToby Isaac   if (iascii) {
385c58f1c22SToby Isaac     ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr);
386c58f1c22SToby Isaac   }
387c58f1c22SToby Isaac   PetscFunctionReturn(0);
388c58f1c22SToby Isaac }
389c58f1c22SToby Isaac 
39084f0b6dfSMatthew G. Knepley /*@
391d67d17b1SMatthew G. Knepley   DMLabelReset - Destroys internal data structures in a DMLabel
392d67d17b1SMatthew G. Knepley 
393d67d17b1SMatthew G. Knepley   Input Parameter:
394d67d17b1SMatthew G. Knepley . label - The DMLabel
395d67d17b1SMatthew G. Knepley 
396d67d17b1SMatthew G. Knepley   Level: beginner
397d67d17b1SMatthew G. Knepley 
398d67d17b1SMatthew G. Knepley .seealso: DMLabelDestroy(), DMLabelCreate()
399d67d17b1SMatthew G. Knepley @*/
400d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelReset(DMLabel label)
401d67d17b1SMatthew G. Knepley {
402d67d17b1SMatthew G. Knepley   PetscInt       v;
403d67d17b1SMatthew G. Knepley   PetscErrorCode ierr;
404d67d17b1SMatthew G. Knepley 
405d67d17b1SMatthew G. Knepley   PetscFunctionBegin;
406d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
407d67d17b1SMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
408d67d17b1SMatthew G. Knepley     ierr = PetscHSetIDestroy(&label->ht[v]);CHKERRQ(ierr);
409d67d17b1SMatthew G. Knepley     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
410d67d17b1SMatthew G. Knepley   }
411*b9907514SLisandro Dalcin   label->numStrata = 0;
412*b9907514SLisandro Dalcin   ierr = PetscFree(label->stratumValues);CHKERRQ(ierr);
413*b9907514SLisandro Dalcin   ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr);
414d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->ht);CHKERRQ(ierr);
415d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->points);CHKERRQ(ierr);
416d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->validIS);CHKERRQ(ierr);
417*b9907514SLisandro Dalcin   ierr = PetscHMapIReset(label->hmap);CHKERRQ(ierr);
418*b9907514SLisandro Dalcin   label->pStart = -1;
419*b9907514SLisandro Dalcin   label->pEnd   = -1;
420d67d17b1SMatthew G. Knepley   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
421d67d17b1SMatthew G. Knepley   PetscFunctionReturn(0);
422d67d17b1SMatthew G. Knepley }
423d67d17b1SMatthew G. Knepley 
424d67d17b1SMatthew G. Knepley /*@
42584f0b6dfSMatthew G. Knepley   DMLabelDestroy - Destroys a DMLabel
42684f0b6dfSMatthew G. Knepley 
42784f0b6dfSMatthew G. Knepley   Input Parameter:
42884f0b6dfSMatthew G. Knepley . label - The DMLabel
42984f0b6dfSMatthew G. Knepley 
43084f0b6dfSMatthew G. Knepley   Level: beginner
43184f0b6dfSMatthew G. Knepley 
432d67d17b1SMatthew G. Knepley .seealso: DMLabelReset(), DMLabelCreate()
43384f0b6dfSMatthew G. Knepley @*/
434c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label)
435c58f1c22SToby Isaac {
436c58f1c22SToby Isaac   PetscErrorCode ierr;
437c58f1c22SToby Isaac 
438c58f1c22SToby Isaac   PetscFunctionBegin;
439d67d17b1SMatthew G. Knepley   if (!*label) PetscFunctionReturn(0);
440d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific((*label),DMLABEL_CLASSID,1);
441*b9907514SLisandro Dalcin   if (--((PetscObject)(*label))->refct > 0) {*label = NULL; PetscFunctionReturn(0);}
442d67d17b1SMatthew G. Knepley   ierr = DMLabelReset(*label);CHKERRQ(ierr);
443*b9907514SLisandro Dalcin   ierr = PetscHMapIDestroy(&(*label)->hmap);CHKERRQ(ierr);
444d67d17b1SMatthew G. Knepley   ierr = PetscHeaderDestroy(label);CHKERRQ(ierr);
445c58f1c22SToby Isaac   PetscFunctionReturn(0);
446c58f1c22SToby Isaac }
447c58f1c22SToby Isaac 
44884f0b6dfSMatthew G. Knepley /*@
44984f0b6dfSMatthew G. Knepley   DMLabelDuplicate - Duplicates a DMLabel
45084f0b6dfSMatthew G. Knepley 
45184f0b6dfSMatthew G. Knepley   Input Parameter:
45284f0b6dfSMatthew G. Knepley . label - The DMLabel
45384f0b6dfSMatthew G. Knepley 
45484f0b6dfSMatthew G. Knepley   Output Parameter:
45584f0b6dfSMatthew G. Knepley . labelnew - location to put new vector
45684f0b6dfSMatthew G. Knepley 
45784f0b6dfSMatthew G. Knepley   Level: intermediate
45884f0b6dfSMatthew G. Knepley 
45984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelDestroy()
46084f0b6dfSMatthew G. Knepley @*/
461c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
462c58f1c22SToby Isaac {
463d67d17b1SMatthew G. Knepley   const char    *name;
464ad8374ffSToby Isaac   PetscInt       v;
465c58f1c22SToby Isaac   PetscErrorCode ierr;
466c58f1c22SToby Isaac 
467c58f1c22SToby Isaac   PetscFunctionBegin;
468d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
469c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
470d67d17b1SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
471d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PetscObjectComm((PetscObject) label), name, labelnew);CHKERRQ(ierr);
472c58f1c22SToby Isaac 
473c58f1c22SToby Isaac   (*labelnew)->numStrata    = label->numStrata;
4745aa44df4SToby Isaac   (*labelnew)->defaultValue = label->defaultValue;
475c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr);
476c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr);
477c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr);
478c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr);
479ad8374ffSToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr);
480c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
481e8f14785SLisandro Dalcin     ierr = PetscHSetICreate(&(*labelnew)->ht[v]);CHKERRQ(ierr);
482c58f1c22SToby Isaac     (*labelnew)->stratumValues[v]  = label->stratumValues[v];
483c58f1c22SToby Isaac     (*labelnew)->stratumSizes[v]   = label->stratumSizes[v];
484ad8374ffSToby Isaac     ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr);
485ad8374ffSToby Isaac     (*labelnew)->points[v]         = label->points[v];
486*b9907514SLisandro Dalcin     (*labelnew)->validIS[v]        = PETSC_TRUE;
487c58f1c22SToby Isaac   }
488*b9907514SLisandro Dalcin   ierr = PetscHMapIDuplicate(label->hmap,&(*labelnew)->hmap);CHKERRQ(ierr);
489c58f1c22SToby Isaac   (*labelnew)->pStart = -1;
490c58f1c22SToby Isaac   (*labelnew)->pEnd   = -1;
491c58f1c22SToby Isaac   (*labelnew)->bt     = NULL;
492c58f1c22SToby Isaac   PetscFunctionReturn(0);
493c58f1c22SToby Isaac }
494c58f1c22SToby Isaac 
495c6a43d28SMatthew G. Knepley /*@
496c6a43d28SMatthew G. Knepley   DMLabelComputeIndex - Create an index structure for membership determination, automatically determining the bounds
497c6a43d28SMatthew G. Knepley 
498c6a43d28SMatthew G. Knepley   Input Parameter:
499c6a43d28SMatthew G. Knepley . label  - The DMLabel
500c6a43d28SMatthew G. Knepley 
501c6a43d28SMatthew G. Knepley   Level: intermediate
502c6a43d28SMatthew G. Knepley 
503c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
504c6a43d28SMatthew G. Knepley @*/
505c6a43d28SMatthew G. Knepley PetscErrorCode DMLabelComputeIndex(DMLabel label)
506c6a43d28SMatthew G. Knepley {
507c6a43d28SMatthew G. Knepley   PetscInt       pStart = PETSC_MAX_INT, pEnd = -1, v;
508c6a43d28SMatthew G. Knepley   PetscErrorCode ierr;
509c6a43d28SMatthew G. Knepley 
510c6a43d28SMatthew G. Knepley   PetscFunctionBegin;
511c6a43d28SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
512c6a43d28SMatthew G. Knepley   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
513c6a43d28SMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
514c6a43d28SMatthew G. Knepley     const PetscInt *points;
515c6a43d28SMatthew G. Knepley     PetscInt       i;
516c6a43d28SMatthew G. Knepley 
517c6a43d28SMatthew G. Knepley     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
518c6a43d28SMatthew G. Knepley     for (i = 0; i < label->stratumSizes[v]; ++i) {
519c6a43d28SMatthew G. Knepley       const PetscInt point = points[i];
520c6a43d28SMatthew G. Knepley 
521c6a43d28SMatthew G. Knepley       pStart = PetscMin(point,   pStart);
522c6a43d28SMatthew G. Knepley       pEnd   = PetscMax(point+1, pEnd);
523c6a43d28SMatthew G. Knepley     }
524c6a43d28SMatthew G. Knepley     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
525c6a43d28SMatthew G. Knepley   }
526c6a43d28SMatthew G. Knepley   label->pStart = pStart == PETSC_MAX_INT ? -1 : pStart;
527c6a43d28SMatthew G. Knepley   label->pEnd   = pEnd;
528c6a43d28SMatthew G. Knepley   ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
529c6a43d28SMatthew G. Knepley   PetscFunctionReturn(0);
530c6a43d28SMatthew G. Knepley }
531c6a43d28SMatthew G. Knepley 
532c6a43d28SMatthew G. Knepley /*@
533c6a43d28SMatthew G. Knepley   DMLabelCreateIndex - Create an index structure for membership determination
534c6a43d28SMatthew G. Knepley 
535c6a43d28SMatthew G. Knepley   Input Parameters:
536c6a43d28SMatthew G. Knepley + label  - The DMLabel
537c6a43d28SMatthew G. Knepley . pStart - The smallest point
538c6a43d28SMatthew G. Knepley - pEnd   - The largest point + 1
539c6a43d28SMatthew G. Knepley 
540c6a43d28SMatthew G. Knepley   Level: intermediate
541c6a43d28SMatthew G. Knepley 
542c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelComputeIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
543c6a43d28SMatthew G. Knepley @*/
544c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
545c58f1c22SToby Isaac {
546c58f1c22SToby Isaac   PetscInt       v;
547c58f1c22SToby Isaac   PetscErrorCode ierr;
548c58f1c22SToby Isaac 
549c58f1c22SToby Isaac   PetscFunctionBegin;
550d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
5510c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
552c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
553c58f1c22SToby Isaac   label->pStart = pStart;
554c58f1c22SToby Isaac   label->pEnd   = pEnd;
555c6a43d28SMatthew G. Knepley   /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
556c58f1c22SToby Isaac   ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr);
557c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
558ad8374ffSToby Isaac     const PetscInt *points;
559c58f1c22SToby Isaac     PetscInt       i;
560c58f1c22SToby Isaac 
561ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
562c58f1c22SToby Isaac     for (i = 0; i < label->stratumSizes[v]; ++i) {
563ad8374ffSToby Isaac       const PetscInt point = points[i];
564c58f1c22SToby Isaac 
565c58f1c22SToby 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);
566c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr);
567c58f1c22SToby Isaac     }
568ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
569c58f1c22SToby Isaac   }
570c58f1c22SToby Isaac   PetscFunctionReturn(0);
571c58f1c22SToby Isaac }
572c58f1c22SToby Isaac 
573c6a43d28SMatthew G. Knepley /*@
574c6a43d28SMatthew G. Knepley   DMLabelDestroyIndex - Destroy the index structure
575c6a43d28SMatthew G. Knepley 
576c6a43d28SMatthew G. Knepley   Input Parameter:
577c6a43d28SMatthew G. Knepley . label - the DMLabel
578c6a43d28SMatthew G. Knepley 
579c6a43d28SMatthew G. Knepley   Level: intermediate
580c6a43d28SMatthew G. Knepley 
581c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
582c6a43d28SMatthew G. Knepley @*/
583c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label)
584c58f1c22SToby Isaac {
585c58f1c22SToby Isaac   PetscErrorCode ierr;
586c58f1c22SToby Isaac 
587c58f1c22SToby Isaac   PetscFunctionBegin;
588d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
589c58f1c22SToby Isaac   label->pStart = -1;
590c58f1c22SToby Isaac   label->pEnd   = -1;
5910c3c4a36SLisandro Dalcin   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
592c58f1c22SToby Isaac   PetscFunctionReturn(0);
593c58f1c22SToby Isaac }
594c58f1c22SToby Isaac 
595c58f1c22SToby Isaac /*@
596c6a43d28SMatthew G. Knepley   DMLabelGetBounds - Return the smallest and largest point in the label
597c6a43d28SMatthew G. Knepley 
598c6a43d28SMatthew G. Knepley   Input Parameter:
599c6a43d28SMatthew G. Knepley . label - the DMLabel
600c6a43d28SMatthew G. Knepley 
601c6a43d28SMatthew G. Knepley   Output Parameters:
602c6a43d28SMatthew G. Knepley + pStart - The smallest point
603c6a43d28SMatthew G. Knepley - pEnd   - The largest point + 1
604c6a43d28SMatthew G. Knepley 
605c6a43d28SMatthew G. Knepley   Note: This will compute an index for the label if one does not exist.
606c6a43d28SMatthew G. Knepley 
607c6a43d28SMatthew G. Knepley   Level: intermediate
608c6a43d28SMatthew G. Knepley 
609c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
610c6a43d28SMatthew G. Knepley @*/
611c6a43d28SMatthew G. Knepley PetscErrorCode DMLabelGetBounds(DMLabel label, PetscInt *pStart, PetscInt *pEnd)
612c6a43d28SMatthew G. Knepley {
613c6a43d28SMatthew G. Knepley   PetscErrorCode ierr;
614c6a43d28SMatthew G. Knepley 
615c6a43d28SMatthew G. Knepley   PetscFunctionBegin;
616c6a43d28SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
617c6a43d28SMatthew G. Knepley   if ((label->pStart == -1) && (label->pEnd == -1)) {ierr = DMLabelComputeIndex(label);CHKERRQ(ierr);}
618c6a43d28SMatthew G. Knepley   if (pStart) {
619c6a43d28SMatthew G. Knepley     PetscValidPointer(pStart, 2);
620c6a43d28SMatthew G. Knepley     *pStart = label->pStart;
621c6a43d28SMatthew G. Knepley   }
622c6a43d28SMatthew G. Knepley   if (pEnd) {
623c6a43d28SMatthew G. Knepley     PetscValidPointer(pEnd, 3);
624c6a43d28SMatthew G. Knepley     *pEnd = label->pEnd;
625c6a43d28SMatthew G. Knepley   }
626c6a43d28SMatthew G. Knepley   PetscFunctionReturn(0);
627c6a43d28SMatthew G. Knepley }
628c6a43d28SMatthew G. Knepley 
629c6a43d28SMatthew G. Knepley /*@
630c58f1c22SToby Isaac   DMLabelHasValue - Determine whether a label assigns the value to any point
631c58f1c22SToby Isaac 
632c58f1c22SToby Isaac   Input Parameters:
633c58f1c22SToby Isaac + label - the DMLabel
634c58f1c22SToby Isaac - value - the value
635c58f1c22SToby Isaac 
636c58f1c22SToby Isaac   Output Parameter:
637c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point
638c58f1c22SToby Isaac 
639c58f1c22SToby Isaac   Level: developer
640c58f1c22SToby Isaac 
641c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
642c58f1c22SToby Isaac @*/
643c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
644c58f1c22SToby Isaac {
645c58f1c22SToby Isaac   PetscInt v;
6460c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
647c58f1c22SToby Isaac 
648c58f1c22SToby Isaac   PetscFunctionBegin;
649d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
650c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
6510c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
6520c3c4a36SLisandro Dalcin   *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE;
653c58f1c22SToby Isaac   PetscFunctionReturn(0);
654c58f1c22SToby Isaac }
655c58f1c22SToby Isaac 
656c58f1c22SToby Isaac /*@
657c58f1c22SToby Isaac   DMLabelHasPoint - Determine whether a label assigns a value to a point
658c58f1c22SToby Isaac 
659c58f1c22SToby Isaac   Input Parameters:
660c58f1c22SToby Isaac + label - the DMLabel
661c58f1c22SToby Isaac - point - the point
662c58f1c22SToby Isaac 
663c58f1c22SToby Isaac   Output Parameter:
664c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value
665c58f1c22SToby Isaac 
666c58f1c22SToby Isaac   Note: The user must call DMLabelCreateIndex() before this function.
667c58f1c22SToby Isaac 
668c58f1c22SToby Isaac   Level: developer
669c58f1c22SToby Isaac 
670c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
671c58f1c22SToby Isaac @*/
672c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
673c58f1c22SToby Isaac {
674c58f1c22SToby Isaac   PetscErrorCode ierr;
675c58f1c22SToby Isaac 
676c58f1c22SToby Isaac   PetscFunctionBeginHot;
677d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
678c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
679c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
680c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG)
681c58f1c22SToby Isaac   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
682c58f1c22SToby 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);
683c58f1c22SToby Isaac #endif
684c58f1c22SToby Isaac   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
685c58f1c22SToby Isaac   PetscFunctionReturn(0);
686c58f1c22SToby Isaac }
687c58f1c22SToby Isaac 
688c58f1c22SToby Isaac /*@
689c58f1c22SToby Isaac   DMLabelStratumHasPoint - Return true if the stratum contains a point
690c58f1c22SToby Isaac 
691c58f1c22SToby Isaac   Input Parameters:
692c58f1c22SToby Isaac + label - the DMLabel
693c58f1c22SToby Isaac . value - the stratum value
694c58f1c22SToby Isaac - point - the point
695c58f1c22SToby Isaac 
696c58f1c22SToby Isaac   Output Parameter:
697c58f1c22SToby Isaac . contains - true if the stratum contains the point
698c58f1c22SToby Isaac 
699c58f1c22SToby Isaac   Level: intermediate
700c58f1c22SToby Isaac 
701c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
702c58f1c22SToby Isaac @*/
703c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
704c58f1c22SToby Isaac {
705c58f1c22SToby Isaac   PetscInt       v;
706c58f1c22SToby Isaac   PetscErrorCode ierr;
707c58f1c22SToby Isaac 
7080c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
709d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
710c58f1c22SToby Isaac   PetscValidPointer(contains, 4);
711c58f1c22SToby Isaac   *contains = PETSC_FALSE;
7120c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
7130c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
7140c3c4a36SLisandro Dalcin 
715ad8374ffSToby Isaac   if (label->validIS[v]) {
716c58f1c22SToby Isaac     PetscInt i;
717c58f1c22SToby Isaac 
718a2d74346SToby Isaac     ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
7190c3c4a36SLisandro Dalcin     if (i >= 0) *contains = PETSC_TRUE;
720c58f1c22SToby Isaac   } else {
721c58f1c22SToby Isaac     PetscBool has;
722c58f1c22SToby Isaac 
723*b9907514SLisandro Dalcin     ierr = PetscHSetIHas(label->ht[v], point, &has);CHKERRQ(ierr);
7240c3c4a36SLisandro Dalcin     if (has) *contains = PETSC_TRUE;
725c58f1c22SToby Isaac   }
726c58f1c22SToby Isaac   PetscFunctionReturn(0);
727c58f1c22SToby Isaac }
728c58f1c22SToby Isaac 
72984f0b6dfSMatthew G. Knepley /*@
7305aa44df4SToby Isaac   DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
7315aa44df4SToby Isaac   When a label is created, it is initialized to -1.
7325aa44df4SToby Isaac 
7335aa44df4SToby Isaac   Input parameter:
7345aa44df4SToby Isaac . label - a DMLabel object
7355aa44df4SToby Isaac 
7365aa44df4SToby Isaac   Output parameter:
7375aa44df4SToby Isaac . defaultValue - the default value
7385aa44df4SToby Isaac 
7395aa44df4SToby Isaac   Level: beginner
7405aa44df4SToby Isaac 
7415aa44df4SToby Isaac .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
74284f0b6dfSMatthew G. Knepley @*/
7435aa44df4SToby Isaac PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
7445aa44df4SToby Isaac {
7455aa44df4SToby Isaac   PetscFunctionBegin;
746d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
7475aa44df4SToby Isaac   *defaultValue = label->defaultValue;
7485aa44df4SToby Isaac   PetscFunctionReturn(0);
7495aa44df4SToby Isaac }
7505aa44df4SToby Isaac 
75184f0b6dfSMatthew G. Knepley /*@
7525aa44df4SToby Isaac   DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
7535aa44df4SToby Isaac   When a label is created, it is initialized to -1.
7545aa44df4SToby Isaac 
7555aa44df4SToby Isaac   Input parameter:
7565aa44df4SToby Isaac . label - a DMLabel object
7575aa44df4SToby Isaac 
7585aa44df4SToby Isaac   Output parameter:
7595aa44df4SToby Isaac . defaultValue - the default value
7605aa44df4SToby Isaac 
7615aa44df4SToby Isaac   Level: beginner
7625aa44df4SToby Isaac 
7635aa44df4SToby Isaac .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
76484f0b6dfSMatthew G. Knepley @*/
7655aa44df4SToby Isaac PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
7665aa44df4SToby Isaac {
7675aa44df4SToby Isaac   PetscFunctionBegin;
768d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
7695aa44df4SToby Isaac   label->defaultValue = defaultValue;
7705aa44df4SToby Isaac   PetscFunctionReturn(0);
7715aa44df4SToby Isaac }
7725aa44df4SToby Isaac 
773c58f1c22SToby Isaac /*@
7745aa44df4SToby 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())
775c58f1c22SToby Isaac 
776c58f1c22SToby Isaac   Input Parameters:
777c58f1c22SToby Isaac + label - the DMLabel
778c58f1c22SToby Isaac - point - the point
779c58f1c22SToby Isaac 
780c58f1c22SToby Isaac   Output Parameter:
7818e68afcfSMatthew G. Knepley . value - The point value, or the default value (-1 by default)
782c58f1c22SToby Isaac 
783c58f1c22SToby Isaac   Level: intermediate
784c58f1c22SToby Isaac 
7855aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
786c58f1c22SToby Isaac @*/
787c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
788c58f1c22SToby Isaac {
789c58f1c22SToby Isaac   PetscInt       v;
790c58f1c22SToby Isaac   PetscErrorCode ierr;
791c58f1c22SToby Isaac 
7920c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
793d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
794c58f1c22SToby Isaac   PetscValidPointer(value, 3);
7955aa44df4SToby Isaac   *value = label->defaultValue;
796c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
797ad8374ffSToby Isaac     if (label->validIS[v]) {
798c58f1c22SToby Isaac       PetscInt i;
799c58f1c22SToby Isaac 
800a2d74346SToby Isaac       ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
801c58f1c22SToby Isaac       if (i >= 0) {
802c58f1c22SToby Isaac         *value = label->stratumValues[v];
803c58f1c22SToby Isaac         break;
804c58f1c22SToby Isaac       }
805c58f1c22SToby Isaac     } else {
806c58f1c22SToby Isaac       PetscBool has;
807c58f1c22SToby Isaac 
808*b9907514SLisandro Dalcin       ierr = PetscHSetIHas(label->ht[v], point, &has);CHKERRQ(ierr);
809c58f1c22SToby Isaac       if (has) {
810c58f1c22SToby Isaac         *value = label->stratumValues[v];
811c58f1c22SToby Isaac         break;
812c58f1c22SToby Isaac       }
813c58f1c22SToby Isaac     }
814c58f1c22SToby Isaac   }
815c58f1c22SToby Isaac   PetscFunctionReturn(0);
816c58f1c22SToby Isaac }
817c58f1c22SToby Isaac 
818c58f1c22SToby Isaac /*@
819367003a6SStefano 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.
820c58f1c22SToby Isaac 
821c58f1c22SToby Isaac   Input Parameters:
822c58f1c22SToby Isaac + label - the DMLabel
823c58f1c22SToby Isaac . point - the point
824c58f1c22SToby Isaac - value - The point value
825c58f1c22SToby Isaac 
826c58f1c22SToby Isaac   Level: intermediate
827c58f1c22SToby Isaac 
8285aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
829c58f1c22SToby Isaac @*/
830c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
831c58f1c22SToby Isaac {
832c58f1c22SToby Isaac   PetscInt       v;
833c58f1c22SToby Isaac   PetscErrorCode ierr;
834c58f1c22SToby Isaac 
835c58f1c22SToby Isaac   PetscFunctionBegin;
836d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
8370c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
8385aa44df4SToby Isaac   if (value == label->defaultValue) PetscFunctionReturn(0);
839*b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
840c58f1c22SToby Isaac   /* Set key */
8410c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
842e8f14785SLisandro Dalcin   ierr = PetscHSetIAdd(label->ht[v], point);CHKERRQ(ierr);
843c58f1c22SToby Isaac   PetscFunctionReturn(0);
844c58f1c22SToby Isaac }
845c58f1c22SToby Isaac 
846c58f1c22SToby Isaac /*@
847c58f1c22SToby Isaac   DMLabelClearValue - Clear the value a label assigns to a point
848c58f1c22SToby Isaac 
849c58f1c22SToby Isaac   Input Parameters:
850c58f1c22SToby Isaac + label - the DMLabel
851c58f1c22SToby Isaac . point - the point
852c58f1c22SToby Isaac - value - The point value
853c58f1c22SToby Isaac 
854c58f1c22SToby Isaac   Level: intermediate
855c58f1c22SToby Isaac 
856c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
857c58f1c22SToby Isaac @*/
858c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
859c58f1c22SToby Isaac {
860ad8374ffSToby Isaac   PetscInt       v;
861c58f1c22SToby Isaac   PetscErrorCode ierr;
862c58f1c22SToby Isaac 
863c58f1c22SToby Isaac   PetscFunctionBegin;
864d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
865c58f1c22SToby Isaac   /* Find label value */
8660c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
8670c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
8680c3c4a36SLisandro Dalcin 
869eeed21e7SToby Isaac   if (label->bt) {
870eeed21e7SToby 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);
871eeed21e7SToby Isaac     ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
872eeed21e7SToby Isaac   }
8730c3c4a36SLisandro Dalcin 
8740c3c4a36SLisandro Dalcin   /* Delete key */
8750c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
876e8f14785SLisandro Dalcin   ierr = PetscHSetIDel(label->ht[v], point);CHKERRQ(ierr);
877c58f1c22SToby Isaac   PetscFunctionReturn(0);
878c58f1c22SToby Isaac }
879c58f1c22SToby Isaac 
880c58f1c22SToby Isaac /*@
881c58f1c22SToby Isaac   DMLabelInsertIS - Set all points in the IS to a value
882c58f1c22SToby Isaac 
883c58f1c22SToby Isaac   Input Parameters:
884c58f1c22SToby Isaac + label - the DMLabel
885c58f1c22SToby Isaac . is    - the point IS
886c58f1c22SToby Isaac - value - The point value
887c58f1c22SToby Isaac 
888c58f1c22SToby Isaac   Level: intermediate
889c58f1c22SToby Isaac 
890c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
891c58f1c22SToby Isaac @*/
892c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
893c58f1c22SToby Isaac {
8940c3c4a36SLisandro Dalcin   PetscInt        v, n, p;
895c58f1c22SToby Isaac   const PetscInt *points;
896c58f1c22SToby Isaac   PetscErrorCode  ierr;
897c58f1c22SToby Isaac 
898c58f1c22SToby Isaac   PetscFunctionBegin;
899d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
900c58f1c22SToby Isaac   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
9010c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
9020c3c4a36SLisandro Dalcin   if (value == label->defaultValue) PetscFunctionReturn(0);
903*b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
9040c3c4a36SLisandro Dalcin   /* Set keys */
9050c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
906c58f1c22SToby Isaac   ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr);
907c58f1c22SToby Isaac   ierr = ISGetIndices(is, &points);CHKERRQ(ierr);
908e8f14785SLisandro Dalcin   for (p = 0; p < n; ++p) {ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);}
909c58f1c22SToby Isaac   ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr);
910c58f1c22SToby Isaac   PetscFunctionReturn(0);
911c58f1c22SToby Isaac }
912c58f1c22SToby Isaac 
91384f0b6dfSMatthew G. Knepley /*@
91484f0b6dfSMatthew G. Knepley   DMLabelGetNumValues - Get the number of values that the DMLabel takes
91584f0b6dfSMatthew G. Knepley 
91684f0b6dfSMatthew G. Knepley   Input Parameter:
91784f0b6dfSMatthew G. Knepley . label - the DMLabel
91884f0b6dfSMatthew G. Knepley 
91984f0b6dfSMatthew G. Knepley   Output Paramater:
92084f0b6dfSMatthew G. Knepley . numValues - the number of values
92184f0b6dfSMatthew G. Knepley 
92284f0b6dfSMatthew G. Knepley   Level: intermediate
92384f0b6dfSMatthew G. Knepley 
92484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
92584f0b6dfSMatthew G. Knepley @*/
926c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
927c58f1c22SToby Isaac {
928c58f1c22SToby Isaac   PetscFunctionBegin;
929d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
930c58f1c22SToby Isaac   PetscValidPointer(numValues, 2);
931c58f1c22SToby Isaac   *numValues = label->numStrata;
932c58f1c22SToby Isaac   PetscFunctionReturn(0);
933c58f1c22SToby Isaac }
934c58f1c22SToby Isaac 
93584f0b6dfSMatthew G. Knepley /*@
93684f0b6dfSMatthew G. Knepley   DMLabelGetValueIS - Get an IS of all values that the DMlabel takes
93784f0b6dfSMatthew G. Knepley 
93884f0b6dfSMatthew G. Knepley   Input Parameter:
93984f0b6dfSMatthew G. Knepley . label - the DMLabel
94084f0b6dfSMatthew G. Knepley 
94184f0b6dfSMatthew G. Knepley   Output Paramater:
94284f0b6dfSMatthew G. Knepley . is    - the value IS
94384f0b6dfSMatthew G. Knepley 
94484f0b6dfSMatthew G. Knepley   Level: intermediate
94584f0b6dfSMatthew G. Knepley 
94684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
94784f0b6dfSMatthew G. Knepley @*/
948c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
949c58f1c22SToby Isaac {
950c58f1c22SToby Isaac   PetscErrorCode ierr;
951c58f1c22SToby Isaac 
952c58f1c22SToby Isaac   PetscFunctionBegin;
953d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
954c58f1c22SToby Isaac   PetscValidPointer(values, 2);
955c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr);
956c58f1c22SToby Isaac   PetscFunctionReturn(0);
957c58f1c22SToby Isaac }
958c58f1c22SToby Isaac 
95984f0b6dfSMatthew G. Knepley /*@
96084f0b6dfSMatthew G. Knepley   DMLabelHasStratum - Determine whether points exist with the given value
96184f0b6dfSMatthew G. Knepley 
96284f0b6dfSMatthew G. Knepley   Input Parameters:
96384f0b6dfSMatthew G. Knepley + label - the DMLabel
96484f0b6dfSMatthew G. Knepley - value - the stratum value
96584f0b6dfSMatthew G. Knepley 
96684f0b6dfSMatthew G. Knepley   Output Paramater:
96784f0b6dfSMatthew G. Knepley . exists - Flag saying whether points exist
96884f0b6dfSMatthew G. Knepley 
96984f0b6dfSMatthew G. Knepley   Level: intermediate
97084f0b6dfSMatthew G. Knepley 
97184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
97284f0b6dfSMatthew G. Knepley @*/
973fada774cSMatthew G. Knepley PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
974fada774cSMatthew G. Knepley {
975fada774cSMatthew G. Knepley   PetscInt       v;
9760c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
977fada774cSMatthew G. Knepley 
978fada774cSMatthew G. Knepley   PetscFunctionBegin;
979d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
980fada774cSMatthew G. Knepley   PetscValidPointer(exists, 3);
9810c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
9820c3c4a36SLisandro Dalcin   *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE;
983fada774cSMatthew G. Knepley   PetscFunctionReturn(0);
984fada774cSMatthew G. Knepley }
985fada774cSMatthew G. Knepley 
98684f0b6dfSMatthew G. Knepley /*@
98784f0b6dfSMatthew G. Knepley   DMLabelGetStratumSize - Get the size of a stratum
98884f0b6dfSMatthew G. Knepley 
98984f0b6dfSMatthew G. Knepley   Input Parameters:
99084f0b6dfSMatthew G. Knepley + label - the DMLabel
99184f0b6dfSMatthew G. Knepley - value - the stratum value
99284f0b6dfSMatthew G. Knepley 
99384f0b6dfSMatthew G. Knepley   Output Paramater:
99484f0b6dfSMatthew G. Knepley . size - The number of points in the stratum
99584f0b6dfSMatthew G. Knepley 
99684f0b6dfSMatthew G. Knepley   Level: intermediate
99784f0b6dfSMatthew G. Knepley 
99884f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
99984f0b6dfSMatthew G. Knepley @*/
1000c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
1001c58f1c22SToby Isaac {
1002c58f1c22SToby Isaac   PetscInt       v;
1003c58f1c22SToby Isaac   PetscErrorCode ierr;
1004c58f1c22SToby Isaac 
1005c58f1c22SToby Isaac   PetscFunctionBegin;
1006d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1007c58f1c22SToby Isaac   PetscValidPointer(size, 3);
1008c58f1c22SToby Isaac   *size = 0;
10090c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
10100c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
1011c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
1012c58f1c22SToby Isaac   *size = label->stratumSizes[v];
1013c58f1c22SToby Isaac   PetscFunctionReturn(0);
1014c58f1c22SToby Isaac }
1015c58f1c22SToby Isaac 
101684f0b6dfSMatthew G. Knepley /*@
101784f0b6dfSMatthew G. Knepley   DMLabelGetStratumBounds - Get the largest and smallest point of a stratum
101884f0b6dfSMatthew G. Knepley 
101984f0b6dfSMatthew G. Knepley   Input Parameters:
102084f0b6dfSMatthew G. Knepley + label - the DMLabel
102184f0b6dfSMatthew G. Knepley - value - the stratum value
102284f0b6dfSMatthew G. Knepley 
102384f0b6dfSMatthew G. Knepley   Output Paramaters:
102484f0b6dfSMatthew G. Knepley + start - the smallest point in the stratum
102584f0b6dfSMatthew G. Knepley - end - the largest point in the stratum
102684f0b6dfSMatthew G. Knepley 
102784f0b6dfSMatthew G. Knepley   Level: intermediate
102884f0b6dfSMatthew G. Knepley 
102984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
103084f0b6dfSMatthew G. Knepley @*/
1031c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
1032c58f1c22SToby Isaac {
10330c3c4a36SLisandro Dalcin   PetscInt       v, min, max;
1034c58f1c22SToby Isaac   PetscErrorCode ierr;
1035c58f1c22SToby Isaac 
1036c58f1c22SToby Isaac   PetscFunctionBegin;
1037d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1038c58f1c22SToby Isaac   if (start) {PetscValidPointer(start, 3); *start = 0;}
1039c58f1c22SToby Isaac   if (end)   {PetscValidPointer(end,   4); *end   = 0;}
10400c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
10410c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
1042c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
10430c3c4a36SLisandro Dalcin   if (label->stratumSizes[v] <= 0) PetscFunctionReturn(0);
1044d6cb179aSToby Isaac   ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr);
1045d6cb179aSToby Isaac   if (start) *start = min;
1046d6cb179aSToby Isaac   if (end)   *end   = max+1;
1047c58f1c22SToby Isaac   PetscFunctionReturn(0);
1048c58f1c22SToby Isaac }
1049c58f1c22SToby Isaac 
105084f0b6dfSMatthew G. Knepley /*@
105184f0b6dfSMatthew G. Knepley   DMLabelGetStratumIS - Get an IS with the stratum points
105284f0b6dfSMatthew G. Knepley 
105384f0b6dfSMatthew G. Knepley   Input Parameters:
105484f0b6dfSMatthew G. Knepley + label - the DMLabel
105584f0b6dfSMatthew G. Knepley - value - the stratum value
105684f0b6dfSMatthew G. Knepley 
105784f0b6dfSMatthew G. Knepley   Output Paramater:
105884f0b6dfSMatthew G. Knepley . points - The stratum points
105984f0b6dfSMatthew G. Knepley 
106084f0b6dfSMatthew G. Knepley   Level: intermediate
106184f0b6dfSMatthew G. Knepley 
106284f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
106384f0b6dfSMatthew G. Knepley @*/
1064c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
1065c58f1c22SToby Isaac {
1066c58f1c22SToby Isaac   PetscInt       v;
1067c58f1c22SToby Isaac   PetscErrorCode ierr;
1068c58f1c22SToby Isaac 
1069c58f1c22SToby Isaac   PetscFunctionBegin;
1070d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1071c58f1c22SToby Isaac   PetscValidPointer(points, 3);
1072c58f1c22SToby Isaac   *points = NULL;
10730c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
10740c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
1075c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
1076ad8374ffSToby Isaac   ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr);
1077ad8374ffSToby Isaac   *points = label->points[v];
1078c58f1c22SToby Isaac   PetscFunctionReturn(0);
1079c58f1c22SToby Isaac }
1080c58f1c22SToby Isaac 
108184f0b6dfSMatthew G. Knepley /*@
108284f0b6dfSMatthew G. Knepley   DMLabelSetStratumIS - Set the stratum points using an IS
108384f0b6dfSMatthew G. Knepley 
108484f0b6dfSMatthew G. Knepley   Input Parameters:
108584f0b6dfSMatthew G. Knepley + label - the DMLabel
108684f0b6dfSMatthew G. Knepley . value - the stratum value
108784f0b6dfSMatthew G. Knepley - points - The stratum points
108884f0b6dfSMatthew G. Knepley 
108984f0b6dfSMatthew G. Knepley   Level: intermediate
109084f0b6dfSMatthew G. Knepley 
109184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
109284f0b6dfSMatthew G. Knepley @*/
10934de306b1SToby Isaac PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
10944de306b1SToby Isaac {
10950c3c4a36SLisandro Dalcin   PetscInt       v;
10964de306b1SToby Isaac   PetscErrorCode ierr;
10974de306b1SToby Isaac 
10984de306b1SToby Isaac   PetscFunctionBegin;
1099d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1100d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
1101*b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
11024de306b1SToby Isaac   if (is == label->points[v]) PetscFunctionReturn(0);
11034de306b1SToby Isaac   ierr = DMLabelClearStratum(label, value);CHKERRQ(ierr);
11044de306b1SToby Isaac   ierr = ISGetLocalSize(is, &(label->stratumSizes[v]));CHKERRQ(ierr);
11054de306b1SToby Isaac   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
11064de306b1SToby Isaac   ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
11070c3c4a36SLisandro Dalcin   label->points[v] = is;
11080c3c4a36SLisandro Dalcin   label->validIS[v] = PETSC_TRUE;
11094de306b1SToby Isaac   if (label->bt) {
11104de306b1SToby Isaac     const PetscInt *points;
11114de306b1SToby Isaac     PetscInt p;
11124de306b1SToby Isaac 
11134de306b1SToby Isaac     ierr = ISGetIndices(is,&points);CHKERRQ(ierr);
11144de306b1SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
11154de306b1SToby Isaac       const PetscInt point = points[p];
11164de306b1SToby Isaac 
11174de306b1SToby 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);
11184de306b1SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
11194de306b1SToby Isaac     }
11204de306b1SToby Isaac   }
11214de306b1SToby Isaac   PetscFunctionReturn(0);
11224de306b1SToby Isaac }
11234de306b1SToby Isaac 
112484f0b6dfSMatthew G. Knepley /*@
112584f0b6dfSMatthew G. Knepley   DMLabelClearStratum - Remove a stratum
11264de306b1SToby Isaac 
112784f0b6dfSMatthew G. Knepley   Input Parameters:
112884f0b6dfSMatthew G. Knepley + label - the DMLabel
112984f0b6dfSMatthew G. Knepley - value - the stratum value
113084f0b6dfSMatthew G. Knepley 
113184f0b6dfSMatthew G. Knepley   Level: intermediate
113284f0b6dfSMatthew G. Knepley 
113384f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
113484f0b6dfSMatthew G. Knepley @*/
1135c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
1136c58f1c22SToby Isaac {
1137c58f1c22SToby Isaac   PetscInt       v;
1138c58f1c22SToby Isaac   PetscErrorCode ierr;
1139c58f1c22SToby Isaac 
1140c58f1c22SToby Isaac   PetscFunctionBegin;
1141d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
11420c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
11430c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
11444de306b1SToby Isaac   if (label->validIS[v]) {
11454de306b1SToby Isaac     if (label->bt) {
1146c58f1c22SToby Isaac       PetscInt       i;
1147ad8374ffSToby Isaac       const PetscInt *points;
1148c58f1c22SToby Isaac 
1149ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
1150c58f1c22SToby Isaac       for (i = 0; i < label->stratumSizes[v]; ++i) {
1151ad8374ffSToby Isaac         const PetscInt point = points[i];
1152c58f1c22SToby Isaac 
1153c58f1c22SToby 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);
1154c58f1c22SToby Isaac         ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
1155c58f1c22SToby Isaac       }
1156ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
1157c58f1c22SToby Isaac     }
1158c58f1c22SToby Isaac     label->stratumSizes[v] = 0;
11590c3c4a36SLisandro Dalcin     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
11600c3c4a36SLisandro Dalcin     ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr);
11610c3c4a36SLisandro Dalcin     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
1162c58f1c22SToby Isaac   } else {
1163*b9907514SLisandro Dalcin     ierr = PetscHSetIClear(label->ht[v]);CHKERRQ(ierr);
1164c58f1c22SToby Isaac   }
1165c58f1c22SToby Isaac   PetscFunctionReturn(0);
1166c58f1c22SToby Isaac }
1167c58f1c22SToby Isaac 
116884f0b6dfSMatthew G. Knepley /*@
116984f0b6dfSMatthew G. Knepley   DMLabelFilter - Remove all points outside of [start, end)
117084f0b6dfSMatthew G. Knepley 
117184f0b6dfSMatthew G. Knepley   Input Parameters:
117284f0b6dfSMatthew G. Knepley + label - the DMLabel
117384f0b6dfSMatthew G. Knepley . start - the first point
117484f0b6dfSMatthew G. Knepley - end - the last point
117584f0b6dfSMatthew G. Knepley 
117684f0b6dfSMatthew G. Knepley   Level: intermediate
117784f0b6dfSMatthew G. Knepley 
117884f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
117984f0b6dfSMatthew G. Knepley @*/
1180c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
1181c58f1c22SToby Isaac {
1182c58f1c22SToby Isaac   PetscInt       v;
1183c58f1c22SToby Isaac   PetscErrorCode ierr;
1184c58f1c22SToby Isaac 
1185c58f1c22SToby Isaac   PetscFunctionBegin;
1186d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
11870c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
1188c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1189c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
1190c58f1c22SToby Isaac     PetscInt off, q;
1191ad8374ffSToby Isaac     const PetscInt *points;
1192033405d5SLisandro Dalcin     PetscInt numPointsNew = 0, *pointsNew = NULL;
1193c58f1c22SToby Isaac 
1194ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
1195033405d5SLisandro Dalcin     for (q = 0; q < label->stratumSizes[v]; ++q)
1196033405d5SLisandro Dalcin       if (points[q] >= start && points[q] < end)
1197033405d5SLisandro Dalcin         numPointsNew++;
1198033405d5SLisandro Dalcin     ierr = PetscMalloc1(numPointsNew, &pointsNew);CHKERRQ(ierr);
1199c58f1c22SToby Isaac     for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) {
1200033405d5SLisandro Dalcin       if (points[q] >= start && points[q] < end)
1201033405d5SLisandro Dalcin         pointsNew[off++] = points[q];
1202ad8374ffSToby Isaac     }
1203ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
1204033405d5SLisandro Dalcin 
1205033405d5SLisandro Dalcin     label->stratumSizes[v] = numPointsNew;
1206033405d5SLisandro Dalcin     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
1207033405d5SLisandro Dalcin     ierr = ISCreateGeneral(PETSC_COMM_SELF,numPointsNew, pointsNew, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr);
1208033405d5SLisandro Dalcin     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
1209c58f1c22SToby Isaac   }
1210c58f1c22SToby Isaac   ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr);
1211c58f1c22SToby Isaac   PetscFunctionReturn(0);
1212c58f1c22SToby Isaac }
1213c58f1c22SToby Isaac 
121484f0b6dfSMatthew G. Knepley /*@
121584f0b6dfSMatthew G. Knepley   DMLabelPermute - Create a new label with permuted points
121684f0b6dfSMatthew G. Knepley 
121784f0b6dfSMatthew G. Knepley   Input Parameters:
121884f0b6dfSMatthew G. Knepley + label - the DMLabel
121984f0b6dfSMatthew G. Knepley - permutation - the point permutation
122084f0b6dfSMatthew G. Knepley 
122184f0b6dfSMatthew G. Knepley   Output Parameter:
122284f0b6dfSMatthew G. Knepley . labelnew - the new label containing the permuted points
122384f0b6dfSMatthew G. Knepley 
122484f0b6dfSMatthew G. Knepley   Level: intermediate
122584f0b6dfSMatthew G. Knepley 
122684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
122784f0b6dfSMatthew G. Knepley @*/
1228c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
1229c58f1c22SToby Isaac {
1230c58f1c22SToby Isaac   const PetscInt *perm;
1231c58f1c22SToby Isaac   PetscInt        numValues, numPoints, v, q;
1232c58f1c22SToby Isaac   PetscErrorCode  ierr;
1233c58f1c22SToby Isaac 
1234c58f1c22SToby Isaac   PetscFunctionBegin;
1235d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1236d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(permutation, IS_CLASSID, 2);
1237c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1238c58f1c22SToby Isaac   ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr);
1239c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr);
1240c58f1c22SToby Isaac   ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr);
1241c58f1c22SToby Isaac   ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr);
1242c58f1c22SToby Isaac   for (v = 0; v < numValues; ++v) {
1243c58f1c22SToby Isaac     const PetscInt size   = (*labelNew)->stratumSizes[v];
1244ad8374ffSToby Isaac     const PetscInt *points;
1245ad8374ffSToby Isaac     PetscInt *pointsNew;
1246c58f1c22SToby Isaac 
1247ad8374ffSToby Isaac     ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1248ad8374ffSToby Isaac     ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr);
1249c58f1c22SToby Isaac     for (q = 0; q < size; ++q) {
1250ad8374ffSToby Isaac       const PetscInt point = points[q];
1251c58f1c22SToby Isaac 
1252c58f1c22SToby 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);
1253ad8374ffSToby Isaac       pointsNew[q] = perm[point];
1254c58f1c22SToby Isaac     }
1255ad8374ffSToby Isaac     ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1256ad8374ffSToby Isaac     ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr);
1257ad8374ffSToby Isaac     ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr);
1258fa8e8ae5SToby Isaac     if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) {
1259fa8e8ae5SToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));CHKERRQ(ierr);
1260fa8e8ae5SToby Isaac       ierr = PetscFree(pointsNew);CHKERRQ(ierr);
1261fa8e8ae5SToby Isaac     } else {
1262ad8374ffSToby Isaac       ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr);
1263fa8e8ae5SToby Isaac     }
1264ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr);
1265c58f1c22SToby Isaac   }
1266c58f1c22SToby Isaac   ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr);
1267c58f1c22SToby Isaac   if (label->bt) {
1268c58f1c22SToby Isaac     ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
1269c58f1c22SToby Isaac     ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
1270c58f1c22SToby Isaac   }
1271c58f1c22SToby Isaac   PetscFunctionReturn(0);
1272c58f1c22SToby Isaac }
1273c58f1c22SToby Isaac 
127426c55118SMichael Lange PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
127526c55118SMichael Lange {
127626c55118SMichael Lange   MPI_Comm       comm;
127726c55118SMichael Lange   PetscInt       s, l, nroots, nleaves, dof, offset, size;
127826c55118SMichael Lange   PetscInt      *remoteOffsets, *rootStrata, *rootIdx;
127926c55118SMichael Lange   PetscSection   rootSection;
128026c55118SMichael Lange   PetscSF        labelSF;
128126c55118SMichael Lange   PetscErrorCode ierr;
128226c55118SMichael Lange 
128326c55118SMichael Lange   PetscFunctionBegin;
128426c55118SMichael Lange   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
128526c55118SMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
128626c55118SMichael Lange   /* Build a section of stratum values per point, generate the according SF
128726c55118SMichael Lange      and distribute point-wise stratum values to leaves. */
128826c55118SMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
128926c55118SMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
129026c55118SMichael Lange   ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr);
129126c55118SMichael Lange   if (label) {
129226c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1293ad8374ffSToby Isaac       const PetscInt *points;
1294ad8374ffSToby Isaac 
1295ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
129626c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1297ad8374ffSToby Isaac         ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr);
1298ad8374ffSToby Isaac         ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr);
129926c55118SMichael Lange       }
1300ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
130126c55118SMichael Lange     }
130226c55118SMichael Lange   }
130326c55118SMichael Lange   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
130426c55118SMichael Lange   /* Create a point-wise array of stratum values */
130526c55118SMichael Lange   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
130626c55118SMichael Lange   ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr);
130726c55118SMichael Lange   ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr);
130826c55118SMichael Lange   if (label) {
130926c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1310ad8374ffSToby Isaac       const PetscInt *points;
1311ad8374ffSToby Isaac 
1312ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
131326c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1314ad8374ffSToby Isaac         const PetscInt p = points[l];
131526c55118SMichael Lange         ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr);
131626c55118SMichael Lange         rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
131726c55118SMichael Lange       }
1318ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
131926c55118SMichael Lange     }
132026c55118SMichael Lange   }
132126c55118SMichael Lange   /* Build SF that maps label points to remote processes */
132226c55118SMichael Lange   ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr);
132326c55118SMichael Lange   ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr);
132426c55118SMichael Lange   ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr);
132526c55118SMichael Lange   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
132626c55118SMichael Lange   /* Send the strata for each point over the derived SF */
132726c55118SMichael Lange   ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr);
132826c55118SMichael Lange   ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr);
132926c55118SMichael Lange   ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
133026c55118SMichael Lange   ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
133126c55118SMichael Lange   /* Clean up */
133226c55118SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
133326c55118SMichael Lange   ierr = PetscFree(rootIdx);CHKERRQ(ierr);
133426c55118SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
133526c55118SMichael Lange   ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr);
133626c55118SMichael Lange   PetscFunctionReturn(0);
133726c55118SMichael Lange }
133826c55118SMichael Lange 
133984f0b6dfSMatthew G. Knepley /*@
134084f0b6dfSMatthew G. Knepley   DMLabelDistribute - Create a new label pushed forward over the PetscSF
134184f0b6dfSMatthew G. Knepley 
134284f0b6dfSMatthew G. Knepley   Input Parameters:
134384f0b6dfSMatthew G. Knepley + label - the DMLabel
134484f0b6dfSMatthew G. Knepley - sf    - the map from old to new distribution
134584f0b6dfSMatthew G. Knepley 
134684f0b6dfSMatthew G. Knepley   Output Parameter:
134784f0b6dfSMatthew G. Knepley . labelnew - the new redistributed label
134884f0b6dfSMatthew G. Knepley 
134984f0b6dfSMatthew G. Knepley   Level: intermediate
135084f0b6dfSMatthew G. Knepley 
135184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
135284f0b6dfSMatthew G. Knepley @*/
1353c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1354c58f1c22SToby Isaac {
1355c58f1c22SToby Isaac   MPI_Comm       comm;
135626c55118SMichael Lange   PetscSection   leafSection;
135726c55118SMichael Lange   PetscInt       p, pStart, pEnd, s, size, dof, offset, stratum;
135826c55118SMichael Lange   PetscInt      *leafStrata, *strataIdx;
1359ad8374ffSToby Isaac   PetscInt     **points;
1360d67d17b1SMatthew G. Knepley   const char    *lname = NULL;
1361c58f1c22SToby Isaac   char          *name;
1362c58f1c22SToby Isaac   PetscInt       nameSize;
1363e8f14785SLisandro Dalcin   PetscHSetI     stratumHash;
1364c58f1c22SToby Isaac   size_t         len = 0;
136526c55118SMichael Lange   PetscMPIInt    rank;
1366c58f1c22SToby Isaac   PetscErrorCode ierr;
1367c58f1c22SToby Isaac 
1368c58f1c22SToby Isaac   PetscFunctionBegin;
1369d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1370f018e600SMatthew G. Knepley   if (label) {
1371f018e600SMatthew G. Knepley     PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1372f018e600SMatthew G. Knepley     ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1373f018e600SMatthew G. Knepley   }
1374c58f1c22SToby Isaac   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1375c58f1c22SToby Isaac   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1376c58f1c22SToby Isaac   /* Bcast name */
1377d67d17b1SMatthew G. Knepley   if (!rank) {
1378d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1379d67d17b1SMatthew G. Knepley     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1380d67d17b1SMatthew G. Knepley   }
1381c58f1c22SToby Isaac   nameSize = len;
1382c58f1c22SToby Isaac   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1383c58f1c22SToby Isaac   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1384d67d17b1SMatthew G. Knepley   if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);}
1385c58f1c22SToby Isaac   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1386d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
1387c58f1c22SToby Isaac   ierr = PetscFree(name);CHKERRQ(ierr);
138877d236dfSMichael Lange   /* Bcast defaultValue */
138977d236dfSMichael Lange   if (!rank) (*labelNew)->defaultValue = label->defaultValue;
139077d236dfSMichael Lange   ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
139126c55118SMichael Lange   /* Distribute stratum values over the SF and get the point mapping on the receiver */
139226c55118SMichael Lange   ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr);
13935cbdf6fcSMichael Lange   /* Determine received stratum values and initialise new label*/
1394e8f14785SLisandro Dalcin   ierr = PetscHSetICreate(&stratumHash);CHKERRQ(ierr);
139526c55118SMichael Lange   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
1396e8f14785SLisandro Dalcin   for (p = 0; p < size; ++p) {ierr = PetscHSetIAdd(stratumHash, leafStrata[p]);CHKERRQ(ierr);}
1397e8f14785SLisandro Dalcin   ierr = PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata);CHKERRQ(ierr);
1398ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr);
1399ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
14005cbdf6fcSMichael Lange   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr);
14015cbdf6fcSMichael Lange   /* Turn leafStrata into indices rather than stratum values */
14025cbdf6fcSMichael Lange   offset = 0;
1403e8f14785SLisandro Dalcin   ierr = PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr);
1404a205f1fdSToby Isaac   ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr);
14055cbdf6fcSMichael Lange   for (p = 0; p < size; ++p) {
1406231b9e6fSMatthew G. Knepley     for (s = 0; s < (*labelNew)->numStrata; ++s) {
1407231b9e6fSMatthew G. Knepley       if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
14085cbdf6fcSMichael Lange     }
14095cbdf6fcSMichael Lange   }
1410c58f1c22SToby Isaac   /* Rebuild the point strata on the receiver */
1411c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr);
1412c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1413c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1414c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1415c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1416c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1417c58f1c22SToby Isaac       (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1418c58f1c22SToby Isaac     }
1419c58f1c22SToby Isaac   }
1420c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr);
1421c58f1c22SToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr);
1422ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr);
1423c58f1c22SToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1424e8f14785SLisandro Dalcin     ierr = PetscHSetICreate(&(*labelNew)->ht[s]);CHKERRQ(ierr);
1425ad8374ffSToby Isaac     ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr);
1426c58f1c22SToby Isaac   }
1427c58f1c22SToby Isaac   /* Insert points into new strata */
1428c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr);
1429c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1430c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1431c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1432c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1433c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1434c58f1c22SToby Isaac       stratum = leafStrata[offset+s];
1435ad8374ffSToby Isaac       points[stratum][strataIdx[stratum]++] = p;
1436c58f1c22SToby Isaac     }
1437c58f1c22SToby Isaac   }
1438ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; s++) {
1439ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr);
1440ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr);
1441ad8374ffSToby Isaac   }
1442ad8374ffSToby Isaac   ierr = PetscFree(points);CHKERRQ(ierr);
1443e8f14785SLisandro Dalcin   ierr = PetscHSetIDestroy(&stratumHash);CHKERRQ(ierr);
1444c58f1c22SToby Isaac   ierr = PetscFree(leafStrata);CHKERRQ(ierr);
1445c58f1c22SToby Isaac   ierr = PetscFree(strataIdx);CHKERRQ(ierr);
1446c58f1c22SToby Isaac   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
1447c58f1c22SToby Isaac   PetscFunctionReturn(0);
1448c58f1c22SToby Isaac }
1449c58f1c22SToby Isaac 
14507937d9ceSMichael Lange /*@
14517937d9ceSMichael Lange   DMLabelGather - Gather all label values from leafs into roots
14527937d9ceSMichael Lange 
14537937d9ceSMichael Lange   Input Parameters:
14547937d9ceSMichael Lange + label - the DMLabel
145584f0b6dfSMatthew G. Knepley - sf - the Star Forest point communication map
14567937d9ceSMichael Lange 
145784f0b6dfSMatthew G. Knepley   Output Parameters:
145884f0b6dfSMatthew G. Knepley . labelNew - the new DMLabel with localised leaf values
14597937d9ceSMichael Lange 
14607937d9ceSMichael Lange   Level: developer
14617937d9ceSMichael Lange 
14627937d9ceSMichael Lange   Note: This is the inverse operation to DMLabelDistribute.
14637937d9ceSMichael Lange 
14647937d9ceSMichael Lange .seealso: DMLabelDistribute()
14657937d9ceSMichael Lange @*/
14667937d9ceSMichael Lange PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
14677937d9ceSMichael Lange {
14687937d9ceSMichael Lange   MPI_Comm       comm;
14697937d9ceSMichael Lange   PetscSection   rootSection;
14707937d9ceSMichael Lange   PetscSF        sfLabel;
14717937d9ceSMichael Lange   PetscSFNode   *rootPoints, *leafPoints;
14727937d9ceSMichael Lange   PetscInt       p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
14737937d9ceSMichael Lange   const PetscInt *rootDegree, *ilocal;
14747937d9ceSMichael Lange   PetscInt       *rootStrata;
1475d67d17b1SMatthew G. Knepley   const char    *lname;
14767937d9ceSMichael Lange   char          *name;
14777937d9ceSMichael Lange   PetscInt       nameSize;
14787937d9ceSMichael Lange   size_t         len = 0;
14799852e123SBarry Smith   PetscMPIInt    rank, size;
14807937d9ceSMichael Lange   PetscErrorCode ierr;
14817937d9ceSMichael Lange 
14827937d9ceSMichael Lange   PetscFunctionBegin;
1483d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1484d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
14857937d9ceSMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
14867937d9ceSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
14879852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
14887937d9ceSMichael Lange   /* Bcast name */
1489d67d17b1SMatthew G. Knepley   if (!rank) {
1490d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1491d67d17b1SMatthew G. Knepley     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1492d67d17b1SMatthew G. Knepley   }
14937937d9ceSMichael Lange   nameSize = len;
14947937d9ceSMichael Lange   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
14957937d9ceSMichael Lange   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1496d67d17b1SMatthew G. Knepley   if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);}
14977937d9ceSMichael Lange   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1498d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
14997937d9ceSMichael Lange   ierr = PetscFree(name);CHKERRQ(ierr);
15007937d9ceSMichael Lange   /* Gather rank/index pairs of leaves into local roots to build
15017937d9ceSMichael Lange      an inverse, multi-rooted SF. Note that this ignores local leaf
15027937d9ceSMichael Lange      indexing due to the use of the multiSF in PetscSFGather. */
15037937d9ceSMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr);
1504dc53bc9bSMatthew G. Knepley   ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr);
1505dc53bc9bSMatthew G. Knepley   for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
15067937d9ceSMichael Lange   for (p = 0; p < nleaves; p++) {
15078212dd46SStefano Zampini     PetscInt ilp = ilocal ? ilocal[p] : p;
15088212dd46SStefano Zampini 
15098212dd46SStefano Zampini     leafPoints[ilp].index = ilp;
15108212dd46SStefano Zampini     leafPoints[ilp].rank  = rank;
15117937d9ceSMichael Lange   }
15127937d9ceSMichael Lange   ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr);
15137937d9ceSMichael Lange   ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr);
15147937d9ceSMichael Lange   for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
15157937d9ceSMichael Lange   ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr);
15167937d9ceSMichael Lange   ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
15177937d9ceSMichael Lange   ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
15187937d9ceSMichael Lange   ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr);
15197937d9ceSMichael Lange   ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
15207937d9ceSMichael Lange   /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
15217937d9ceSMichael Lange   ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr);
15227937d9ceSMichael Lange   /* Rebuild the point strata on the receiver */
15237937d9ceSMichael Lange   for (p = 0, idx = 0; p < nroots; p++) {
15247937d9ceSMichael Lange     for (d = 0; d < rootDegree[p]; d++) {
15257937d9ceSMichael Lange       ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr);
15267937d9ceSMichael Lange       ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr);
15277937d9ceSMichael Lange       for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);}
15287937d9ceSMichael Lange     }
15297937d9ceSMichael Lange     idx += rootDegree[p];
15307937d9ceSMichael Lange   }
153177e0c0e7SMichael Lange   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
153277e0c0e7SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
153377e0c0e7SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
153477e0c0e7SMichael Lange   ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr);
15357937d9ceSMichael Lange   PetscFunctionReturn(0);
15367937d9ceSMichael Lange }
15377937d9ceSMichael Lange 
153884f0b6dfSMatthew G. Knepley /*@
153984f0b6dfSMatthew G. Knepley   DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label
154084f0b6dfSMatthew G. Knepley 
154184f0b6dfSMatthew G. Knepley   Input Parameter:
154284f0b6dfSMatthew G. Knepley . label - the DMLabel
154384f0b6dfSMatthew G. Knepley 
154484f0b6dfSMatthew G. Knepley   Output Parameters:
154584f0b6dfSMatthew G. Knepley + section - the section giving offsets for each stratum
154684f0b6dfSMatthew G. Knepley - is - An IS containing all the label points
154784f0b6dfSMatthew G. Knepley 
154884f0b6dfSMatthew G. Knepley   Level: developer
154984f0b6dfSMatthew G. Knepley 
155084f0b6dfSMatthew G. Knepley .seealso: DMLabelDistribute()
155184f0b6dfSMatthew G. Knepley @*/
1552c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1553c58f1c22SToby Isaac {
1554c58f1c22SToby Isaac   IS              vIS;
1555c58f1c22SToby Isaac   const PetscInt *values;
1556c58f1c22SToby Isaac   PetscInt       *points;
1557c58f1c22SToby Isaac   PetscInt        nV, vS = 0, vE = 0, v, N;
1558c58f1c22SToby Isaac   PetscErrorCode  ierr;
1559c58f1c22SToby Isaac 
1560c58f1c22SToby Isaac   PetscFunctionBegin;
1561d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1562c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr);
1563c58f1c22SToby Isaac   ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
1564c58f1c22SToby Isaac   ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr);
1565c58f1c22SToby Isaac   if (nV) {vS = values[0]; vE = values[0]+1;}
1566c58f1c22SToby Isaac   for (v = 1; v < nV; ++v) {
1567c58f1c22SToby Isaac     vS = PetscMin(vS, values[v]);
1568c58f1c22SToby Isaac     vE = PetscMax(vE, values[v]+1);
1569c58f1c22SToby Isaac   }
1570c58f1c22SToby Isaac   ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr);
1571c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr);
1572c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1573c58f1c22SToby Isaac     PetscInt n;
1574c58f1c22SToby Isaac 
1575c58f1c22SToby Isaac     ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr);
1576c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr);
1577c58f1c22SToby Isaac   }
1578c58f1c22SToby Isaac   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
1579c58f1c22SToby Isaac   ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr);
1580c58f1c22SToby Isaac   ierr = PetscMalloc1(N, &points);CHKERRQ(ierr);
1581c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1582c58f1c22SToby Isaac     IS              is;
1583c58f1c22SToby Isaac     const PetscInt *spoints;
1584c58f1c22SToby Isaac     PetscInt        dof, off, p;
1585c58f1c22SToby Isaac 
1586c58f1c22SToby Isaac     ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr);
1587c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr);
1588c58f1c22SToby Isaac     ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr);
1589c58f1c22SToby Isaac     ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr);
1590c58f1c22SToby Isaac     for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1591c58f1c22SToby Isaac     ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr);
1592c58f1c22SToby Isaac     ierr = ISDestroy(&is);CHKERRQ(ierr);
1593c58f1c22SToby Isaac   }
1594c58f1c22SToby Isaac   ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr);
1595c58f1c22SToby Isaac   ierr = ISDestroy(&vIS);CHKERRQ(ierr);
1596c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1597c58f1c22SToby Isaac   PetscFunctionReturn(0);
1598c58f1c22SToby Isaac }
1599c58f1c22SToby Isaac 
160084f0b6dfSMatthew G. Knepley /*@
1601c58f1c22SToby Isaac   PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1602c58f1c22SToby Isaac   the local section and an SF describing the section point overlap.
1603c58f1c22SToby Isaac 
1604c58f1c22SToby Isaac   Input Parameters:
1605c58f1c22SToby Isaac   + s - The PetscSection for the local field layout
1606c58f1c22SToby Isaac   . sf - The SF describing parallel layout of the section points
1607c58f1c22SToby Isaac   . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1608c58f1c22SToby Isaac   . label - The label specifying the points
1609c58f1c22SToby Isaac   - labelValue - The label stratum specifying the points
1610c58f1c22SToby Isaac 
1611c58f1c22SToby Isaac   Output Parameter:
1612c58f1c22SToby Isaac   . gsection - The PetscSection for the global field layout
1613c58f1c22SToby Isaac 
1614c58f1c22SToby Isaac   Note: This gives negative sizes and offsets to points not owned by this process
1615c58f1c22SToby Isaac 
1616c58f1c22SToby Isaac   Level: developer
1617c58f1c22SToby Isaac 
1618c58f1c22SToby Isaac .seealso: PetscSectionCreate()
1619c58f1c22SToby Isaac @*/
1620c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1621c58f1c22SToby Isaac {
1622c58f1c22SToby Isaac   PetscInt      *neg = NULL, *tmpOff = NULL;
1623c58f1c22SToby Isaac   PetscInt       pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1624c58f1c22SToby Isaac   PetscErrorCode ierr;
1625c58f1c22SToby Isaac 
1626c58f1c22SToby Isaac   PetscFunctionBegin;
1627d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1628d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1629d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4);
1630c58f1c22SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr);
1631c58f1c22SToby Isaac   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1632c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr);
1633c58f1c22SToby Isaac   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1634c58f1c22SToby Isaac   if (nroots >= 0) {
1635c58f1c22SToby Isaac     if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1636c58f1c22SToby Isaac     ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr);
1637c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1638c58f1c22SToby Isaac       ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr);
1639c58f1c22SToby Isaac     } else {
1640c58f1c22SToby Isaac       tmpOff = &(*gsection)->atlasDof[-pStart];
1641c58f1c22SToby Isaac     }
1642c58f1c22SToby Isaac   }
1643c58f1c22SToby Isaac   /* Mark ghost points with negative dof */
1644c58f1c22SToby Isaac   for (p = pStart; p < pEnd; ++p) {
1645c58f1c22SToby Isaac     PetscInt value;
1646c58f1c22SToby Isaac 
1647c58f1c22SToby Isaac     ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr);
1648c58f1c22SToby Isaac     if (value != labelValue) continue;
1649c58f1c22SToby Isaac     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1650c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr);
1651c58f1c22SToby Isaac     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1652c58f1c22SToby Isaac     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);}
1653c58f1c22SToby Isaac     if (neg) neg[p] = -(dof+1);
1654c58f1c22SToby Isaac   }
1655c58f1c22SToby Isaac   ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr);
1656c58f1c22SToby Isaac   if (nroots >= 0) {
1657c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1658c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1659c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1660c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1661c58f1c22SToby Isaac     }
1662c58f1c22SToby Isaac   }
1663c58f1c22SToby Isaac   /* Calculate new sizes, get proccess offset, and calculate point offsets */
1664c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1665c58f1c22SToby Isaac     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1666c58f1c22SToby Isaac     (*gsection)->atlasOff[p] = off;
1667c58f1c22SToby Isaac     off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1668c58f1c22SToby Isaac   }
1669c58f1c22SToby Isaac   ierr       = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1670c58f1c22SToby Isaac   globalOff -= off;
1671c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1672c58f1c22SToby Isaac     (*gsection)->atlasOff[p] += globalOff;
1673c58f1c22SToby Isaac     if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1674c58f1c22SToby Isaac   }
1675c58f1c22SToby Isaac   /* Put in negative offsets for ghost points */
1676c58f1c22SToby Isaac   if (nroots >= 0) {
1677c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1678c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1679c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1680c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1681c58f1c22SToby Isaac     }
1682c58f1c22SToby Isaac   }
1683c58f1c22SToby Isaac   if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);}
1684c58f1c22SToby Isaac   ierr = PetscFree(neg);CHKERRQ(ierr);
1685c58f1c22SToby Isaac   PetscFunctionReturn(0);
1686c58f1c22SToby Isaac }
1687c58f1c22SToby Isaac 
16885fdea053SToby Isaac typedef struct _n_PetscSectionSym_Label
16895fdea053SToby Isaac {
16905fdea053SToby Isaac   DMLabel           label;
16915fdea053SToby Isaac   PetscCopyMode     *modes;
16925fdea053SToby Isaac   PetscInt          *sizes;
16935fdea053SToby Isaac   const PetscInt    ***perms;
16945fdea053SToby Isaac   const PetscScalar ***rots;
16955fdea053SToby Isaac   PetscInt          (*minMaxOrients)[2];
16965fdea053SToby Isaac   PetscInt          numStrata; /* numStrata is only increasing, functions as a state */
16975fdea053SToby Isaac } PetscSectionSym_Label;
16985fdea053SToby Isaac 
16995fdea053SToby Isaac static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
17005fdea053SToby Isaac {
17015fdea053SToby Isaac   PetscInt              i, j;
17025fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
17035fdea053SToby Isaac   PetscErrorCode        ierr;
17045fdea053SToby Isaac 
17055fdea053SToby Isaac   PetscFunctionBegin;
17065fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
17075fdea053SToby Isaac     if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
17085fdea053SToby Isaac       for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
17095fdea053SToby Isaac         if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);}
17105fdea053SToby Isaac         if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);}
17115fdea053SToby Isaac       }
17125fdea053SToby Isaac       if (sl->perms[i]) {
17135fdea053SToby Isaac         const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
17145fdea053SToby Isaac 
17155fdea053SToby Isaac         ierr = PetscFree(perms);CHKERRQ(ierr);
17165fdea053SToby Isaac       }
17175fdea053SToby Isaac       if (sl->rots[i]) {
17185fdea053SToby Isaac         const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
17195fdea053SToby Isaac 
17205fdea053SToby Isaac         ierr = PetscFree(rots);CHKERRQ(ierr);
17215fdea053SToby Isaac       }
17225fdea053SToby Isaac     }
17235fdea053SToby Isaac   }
17245fdea053SToby Isaac   ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr);
17255fdea053SToby Isaac   ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr);
17265fdea053SToby Isaac   sl->numStrata = 0;
17275fdea053SToby Isaac   PetscFunctionReturn(0);
17285fdea053SToby Isaac }
17295fdea053SToby Isaac 
17305fdea053SToby Isaac static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
17315fdea053SToby Isaac {
17325fdea053SToby Isaac   PetscErrorCode ierr;
17335fdea053SToby Isaac 
17345fdea053SToby Isaac   PetscFunctionBegin;
17355fdea053SToby Isaac   ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);
17365fdea053SToby Isaac   ierr = PetscFree(sym->data);CHKERRQ(ierr);
17375fdea053SToby Isaac   PetscFunctionReturn(0);
17385fdea053SToby Isaac }
17395fdea053SToby Isaac 
17405fdea053SToby Isaac static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
17415fdea053SToby Isaac {
17425fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
17435fdea053SToby Isaac   PetscBool             isAscii;
17445fdea053SToby Isaac   DMLabel               label = sl->label;
1745d67d17b1SMatthew G. Knepley   const char           *name;
17465fdea053SToby Isaac   PetscErrorCode        ierr;
17475fdea053SToby Isaac 
17485fdea053SToby Isaac   PetscFunctionBegin;
17495fdea053SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
17505fdea053SToby Isaac   if (isAscii) {
17515fdea053SToby Isaac     PetscInt          i, j, k;
17525fdea053SToby Isaac     PetscViewerFormat format;
17535fdea053SToby Isaac 
17545fdea053SToby Isaac     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
17555fdea053SToby Isaac     if (label) {
17565fdea053SToby Isaac       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
17575fdea053SToby Isaac       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
17585fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
17595fdea053SToby Isaac         ierr = DMLabelView(label, viewer);CHKERRQ(ierr);
17605fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
17615fdea053SToby Isaac       } else {
1762d67d17b1SMatthew G. Knepley         ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1763d67d17b1SMatthew G. Knepley         ierr = PetscViewerASCIIPrintf(viewer,"  Label '%s'\n",name);CHKERRQ(ierr);
17645fdea053SToby Isaac       }
17655fdea053SToby Isaac     } else {
17665fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr);
17675fdea053SToby Isaac     }
17685fdea053SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
17695fdea053SToby Isaac     for (i = 0; i <= sl->numStrata; i++) {
17705fdea053SToby Isaac       PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
17715fdea053SToby Isaac 
17725fdea053SToby Isaac       if (!(sl->perms[i] || sl->rots[i])) {
17735fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr);
17745fdea053SToby Isaac       } else {
17755fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr);
17765fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
17775fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr);
17785fdea053SToby Isaac         if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
17795fdea053SToby Isaac           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
17805fdea053SToby Isaac           for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
17815fdea053SToby Isaac             if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
17825fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr);
17835fdea053SToby Isaac             } else {
17845fdea053SToby Isaac               PetscInt tab;
17855fdea053SToby Isaac 
17865fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr);
17875fdea053SToby Isaac               ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
17885fdea053SToby Isaac               ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr);
17895fdea053SToby Isaac               if (sl->perms[i] && sl->perms[i][j]) {
17905fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr);
17915fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
17925fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);}
17935fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
17945fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
17955fdea053SToby Isaac               }
17965fdea053SToby Isaac               if (sl->rots[i] && sl->rots[i][j]) {
17975fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Rotations:  ");CHKERRQ(ierr);
17985fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
17995fdea053SToby Isaac #if defined(PETSC_USE_COMPLEX)
18005fdea053SToby 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);}
18015fdea053SToby Isaac #else
18025fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);}
18035fdea053SToby Isaac #endif
18045fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
18055fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
18065fdea053SToby Isaac               }
18075fdea053SToby Isaac               ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
18085fdea053SToby Isaac             }
18095fdea053SToby Isaac           }
18105fdea053SToby Isaac           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
18115fdea053SToby Isaac         }
18125fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
18135fdea053SToby Isaac       }
18145fdea053SToby Isaac     }
18155fdea053SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
18165fdea053SToby Isaac   }
18175fdea053SToby Isaac   PetscFunctionReturn(0);
18185fdea053SToby Isaac }
18195fdea053SToby Isaac 
18205fdea053SToby Isaac /*@
18215fdea053SToby Isaac   PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
18225fdea053SToby Isaac 
18235fdea053SToby Isaac   Logically collective on sym
18245fdea053SToby Isaac 
18255fdea053SToby Isaac   Input parameters:
18265fdea053SToby Isaac + sym - the section symmetries
18275fdea053SToby Isaac - label - the DMLabel describing the types of points
18285fdea053SToby Isaac 
18295fdea053SToby Isaac   Level: developer:
18305fdea053SToby Isaac 
18315fdea053SToby Isaac .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()
18325fdea053SToby Isaac @*/
18335fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
18345fdea053SToby Isaac {
18355fdea053SToby Isaac   PetscSectionSym_Label *sl;
18365fdea053SToby Isaac   PetscErrorCode        ierr;
18375fdea053SToby Isaac 
18385fdea053SToby Isaac   PetscFunctionBegin;
18395fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
18405fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
18415fdea053SToby Isaac   if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);}
18425fdea053SToby Isaac   if (label) {
18435fdea053SToby Isaac     sl->label = label;
1844d67d17b1SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) label);CHKERRQ(ierr);
18455fdea053SToby Isaac     ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr);
18461a834cf9SToby 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);
18471a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr);
18481a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr);
18491a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr);
18501a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr);
18511a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr);
18525fdea053SToby Isaac   }
18535fdea053SToby Isaac   PetscFunctionReturn(0);
18545fdea053SToby Isaac }
18555fdea053SToby Isaac 
18565fdea053SToby Isaac /*@C
18575fdea053SToby Isaac   PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
18585fdea053SToby Isaac 
18595fdea053SToby Isaac   Logically collective on PetscSectionSym
18605fdea053SToby Isaac 
18615fdea053SToby Isaac   InputParameters:
18625fdea053SToby Isaac + sys       - the section symmetries
18635fdea053SToby Isaac . stratum   - the stratum value in the label that we are assigning symmetries for
18645fdea053SToby Isaac . size      - the number of dofs for points in the stratum of the label
18655fdea053SToby Isaac . minOrient - the smallest orientation for a point in this stratum
18665fdea053SToby Isaac . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
18675fdea053SToby Isaac . mode      - how sym should copy the perms and rots arrays
18685fdea053SToby Isaac . perms     - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation.  A NULL permutation is the identity
18695fdea053SToby 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
18705fdea053SToby Isaac 
18715fdea053SToby Isaac   Level: developer
18725fdea053SToby Isaac 
18735fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()
18745fdea053SToby Isaac @*/
18755fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
18765fdea053SToby Isaac {
18775fdea053SToby Isaac   PetscSectionSym_Label *sl;
1878d67d17b1SMatthew G. Knepley   const char            *name;
1879d67d17b1SMatthew G. Knepley   PetscInt               i, j, k;
18805fdea053SToby Isaac   PetscErrorCode         ierr;
18815fdea053SToby Isaac 
18825fdea053SToby Isaac   PetscFunctionBegin;
18835fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
18845fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
18855fdea053SToby Isaac   if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet");
18865fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
18875fdea053SToby Isaac     PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
18885fdea053SToby Isaac 
18895fdea053SToby Isaac     if (stratum == value) break;
18905fdea053SToby Isaac   }
1891d67d17b1SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1892d67d17b1SMatthew G. Knepley   if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,name);
18935fdea053SToby Isaac   sl->sizes[i] = size;
18945fdea053SToby Isaac   sl->modes[i] = mode;
18955fdea053SToby Isaac   sl->minMaxOrients[i][0] = minOrient;
18965fdea053SToby Isaac   sl->minMaxOrients[i][1] = maxOrient;
18975fdea053SToby Isaac   if (mode == PETSC_COPY_VALUES) {
18985fdea053SToby Isaac     if (perms) {
18995fdea053SToby Isaac       PetscInt    **ownPerms;
19005fdea053SToby Isaac 
19015fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr);
19025fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
19035fdea053SToby Isaac         if (perms[j]) {
19045fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr);
19055fdea053SToby Isaac           for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];}
19065fdea053SToby Isaac         }
19075fdea053SToby Isaac       }
19085fdea053SToby Isaac       sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient];
19095fdea053SToby Isaac     }
19105fdea053SToby Isaac     if (rots) {
19115fdea053SToby Isaac       PetscScalar **ownRots;
19125fdea053SToby Isaac 
19135fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr);
19145fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
19155fdea053SToby Isaac         if (rots[j]) {
19165fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr);
19175fdea053SToby Isaac           for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];}
19185fdea053SToby Isaac         }
19195fdea053SToby Isaac       }
19205fdea053SToby Isaac       sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient];
19215fdea053SToby Isaac     }
19225fdea053SToby Isaac   } else {
19235fdea053SToby Isaac     sl->perms[i] = perms ? &perms[-minOrient] : NULL;
19245fdea053SToby Isaac     sl->rots[i]  = rots ? &rots[-minOrient] : NULL;
19255fdea053SToby Isaac   }
19265fdea053SToby Isaac   PetscFunctionReturn(0);
19275fdea053SToby Isaac }
19285fdea053SToby Isaac 
19295fdea053SToby Isaac static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
19305fdea053SToby Isaac {
19315fdea053SToby Isaac   PetscInt              i, j, numStrata;
19325fdea053SToby Isaac   PetscSectionSym_Label *sl;
19335fdea053SToby Isaac   DMLabel               label;
19345fdea053SToby Isaac   PetscErrorCode        ierr;
19355fdea053SToby Isaac 
19365fdea053SToby Isaac   PetscFunctionBegin;
19375fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
19385fdea053SToby Isaac   numStrata = sl->numStrata;
19395fdea053SToby Isaac   label     = sl->label;
19405fdea053SToby Isaac   for (i = 0; i < numPoints; i++) {
19415fdea053SToby Isaac     PetscInt point = points[2*i];
19425fdea053SToby Isaac     PetscInt ornt  = points[2*i+1];
19435fdea053SToby Isaac 
19445fdea053SToby Isaac     for (j = 0; j < numStrata; j++) {
19455fdea053SToby Isaac       if (label->validIS[j]) {
19465fdea053SToby Isaac         PetscInt k;
19475fdea053SToby Isaac 
19485fdea053SToby Isaac         ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr);
19495fdea053SToby Isaac         if (k >= 0) break;
19505fdea053SToby Isaac       } else {
19515fdea053SToby Isaac         PetscBool has;
19525fdea053SToby Isaac 
1953*b9907514SLisandro Dalcin         ierr = PetscHSetIHas(label->ht[j], point, &has);CHKERRQ(ierr);
19545fdea053SToby Isaac         if (has) break;
19555fdea053SToby Isaac       }
19565fdea053SToby Isaac     }
19575fdea053SToby 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);
19585fdea053SToby Isaac     if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;}
19595fdea053SToby Isaac     if (rots) {rots[i]  = sl->rots[j] ? sl->rots[j][ornt] : NULL;}
19605fdea053SToby Isaac   }
19615fdea053SToby Isaac   PetscFunctionReturn(0);
19625fdea053SToby Isaac }
19635fdea053SToby Isaac 
19645fdea053SToby Isaac PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
19655fdea053SToby Isaac {
19665fdea053SToby Isaac   PetscSectionSym_Label *sl;
19675fdea053SToby Isaac   PetscErrorCode        ierr;
19685fdea053SToby Isaac 
19695fdea053SToby Isaac   PetscFunctionBegin;
19705fdea053SToby Isaac   ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr);
19715fdea053SToby Isaac   sym->ops->getpoints = PetscSectionSymGetPoints_Label;
19725fdea053SToby Isaac   sym->ops->view      = PetscSectionSymView_Label;
19735fdea053SToby Isaac   sym->ops->destroy   = PetscSectionSymDestroy_Label;
19745fdea053SToby Isaac   sym->data           = (void *) sl;
19755fdea053SToby Isaac   PetscFunctionReturn(0);
19765fdea053SToby Isaac }
19775fdea053SToby Isaac 
19785fdea053SToby Isaac /*@
19795fdea053SToby Isaac   PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
19805fdea053SToby Isaac 
19815fdea053SToby Isaac   Collective
19825fdea053SToby Isaac 
19835fdea053SToby Isaac   Input Parameters:
19845fdea053SToby Isaac + comm - the MPI communicator for the new symmetry
19855fdea053SToby Isaac - label - the label defining the strata
19865fdea053SToby Isaac 
19875fdea053SToby Isaac   Output Parameters:
19885fdea053SToby Isaac . sym - the section symmetries
19895fdea053SToby Isaac 
19905fdea053SToby Isaac   Level: developer
19915fdea053SToby Isaac 
19925fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()
19935fdea053SToby Isaac @*/
19945fdea053SToby Isaac PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
19955fdea053SToby Isaac {
19965fdea053SToby Isaac   PetscErrorCode        ierr;
19975fdea053SToby Isaac 
19985fdea053SToby Isaac   PetscFunctionBegin;
19995fdea053SToby Isaac   ierr = DMInitializePackage();CHKERRQ(ierr);
20005fdea053SToby Isaac   ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr);
20015fdea053SToby Isaac   ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr);
20025fdea053SToby Isaac   ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr);
20035fdea053SToby Isaac   PetscFunctionReturn(0);
20045fdea053SToby Isaac }
2005