xref: /petsc/src/dm/label/dmlabel.c (revision 277ea44a122bf0a8226b47f1967981bcacd44b51)
15fdea053SToby Isaac #include <petscdm.h>
2c58f1c22SToby Isaac #include <petsc/private/dmlabelimpl.h>   /*I      "petscdmlabel.h"   I*/
3c58f1c22SToby Isaac #include <petsc/private/isimpl.h>        /*I      "petscis.h"        I*/
4c58f1c22SToby Isaac #include <petscsf.h>
5c58f1c22SToby Isaac 
6c58f1c22SToby Isaac /*@C
7c58f1c22SToby Isaac   DMLabelCreate - Create a DMLabel object, which is a multimap
8c58f1c22SToby Isaac 
9d67d17b1SMatthew G. Knepley   Input parameters:
10d67d17b1SMatthew G. Knepley + comm - The communicator, usually PETSC_COMM_SELF
11d67d17b1SMatthew G. Knepley - name - The label name
12c58f1c22SToby Isaac 
13c58f1c22SToby Isaac   Output parameter:
14c58f1c22SToby Isaac . label - The DMLabel
15c58f1c22SToby Isaac 
16c58f1c22SToby Isaac   Level: beginner
17c58f1c22SToby Isaac 
18c58f1c22SToby Isaac .seealso: DMLabelDestroy()
19c58f1c22SToby Isaac @*/
20d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelCreate(MPI_Comm comm, const char name[], DMLabel *label)
21c58f1c22SToby Isaac {
22c58f1c22SToby Isaac   PetscErrorCode ierr;
23c58f1c22SToby Isaac 
24c58f1c22SToby Isaac   PetscFunctionBegin;
25d67d17b1SMatthew G. Knepley   PetscValidPointer(label,2);
26d67d17b1SMatthew G. Knepley   ierr = DMInitializePackage();CHKERRQ(ierr);
27c58f1c22SToby Isaac 
28d67d17b1SMatthew G. Knepley   ierr = PetscHeaderCreate(*label,DMLABEL_CLASSID,"DMLabel","DMLabel","DM",comm,DMLabelDestroy,DMLabelView);CHKERRQ(ierr);
29d67d17b1SMatthew G. Knepley 
30c58f1c22SToby Isaac   (*label)->numStrata      = 0;
315aa44df4SToby Isaac   (*label)->defaultValue   = -1;
32c58f1c22SToby Isaac   (*label)->stratumValues  = NULL;
33ad8374ffSToby Isaac   (*label)->validIS        = NULL;
34c58f1c22SToby Isaac   (*label)->stratumSizes   = NULL;
35c58f1c22SToby Isaac   (*label)->points         = NULL;
36c58f1c22SToby Isaac   (*label)->ht             = NULL;
37c58f1c22SToby Isaac   (*label)->pStart         = -1;
38c58f1c22SToby Isaac   (*label)->pEnd           = -1;
39c58f1c22SToby Isaac   (*label)->bt             = NULL;
40b9907514SLisandro Dalcin   ierr = PetscHMapICreate(&(*label)->hmap);CHKERRQ(ierr);
41d67d17b1SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *label, name);CHKERRQ(ierr);
42c58f1c22SToby Isaac   PetscFunctionReturn(0);
43c58f1c22SToby Isaac }
44c58f1c22SToby Isaac 
45c58f1c22SToby Isaac /*
46c58f1c22SToby Isaac   DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format
47c58f1c22SToby Isaac 
48c58f1c22SToby Isaac   Input parameter:
49c58f1c22SToby Isaac + label - The DMLabel
50c58f1c22SToby Isaac - v - The stratum value
51c58f1c22SToby Isaac 
52c58f1c22SToby Isaac   Output parameter:
53c58f1c22SToby Isaac . label - The DMLabel with stratum in sorted list format
54c58f1c22SToby Isaac 
55c58f1c22SToby Isaac   Level: developer
56c58f1c22SToby Isaac 
57c58f1c22SToby Isaac .seealso: DMLabelCreate()
58c58f1c22SToby Isaac */
59c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v)
60c58f1c22SToby Isaac {
61*277ea44aSLisandro Dalcin   IS             is;
62b9907514SLisandro Dalcin   PetscInt       off = 0, *pointArray, p;
63c58f1c22SToby Isaac   PetscErrorCode ierr;
64c58f1c22SToby Isaac 
65b9907514SLisandro Dalcin   if (PetscLikely(v >= 0 && v < label->numStrata) && label->validIS[v]) return 0;
66c58f1c22SToby Isaac   PetscFunctionBegin;
670c3c4a36SLisandro 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);
68e8f14785SLisandro Dalcin   ierr = PetscHSetIGetSize(label->ht[v], &label->stratumSizes[v]);CHKERRQ(ierr);
69ad8374ffSToby Isaac   ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr);
70e8f14785SLisandro Dalcin   ierr = PetscHSetIGetElems(label->ht[v], &off, pointArray);CHKERRQ(ierr);
71b9907514SLisandro Dalcin   ierr = PetscHSetIClear(label->ht[v]);CHKERRQ(ierr);
72ad8374ffSToby Isaac   ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr);
73c58f1c22SToby Isaac   if (label->bt) {
74c58f1c22SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
75ad8374ffSToby Isaac       const PetscInt point = pointArray[p];
76c58f1c22SToby 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);
77c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
78c58f1c22SToby Isaac     }
79c58f1c22SToby Isaac   }
80*277ea44aSLisandro Dalcin   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->stratumSizes[v], pointArray, PETSC_OWN_POINTER, &is);CHKERRQ(ierr);
81*277ea44aSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject) is, "indices");CHKERRQ(ierr);
82*277ea44aSLisandro Dalcin   label->points[v]  = is;
83ad8374ffSToby Isaac   label->validIS[v] = PETSC_TRUE;
84d67d17b1SMatthew G. Knepley   ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
85c58f1c22SToby Isaac   PetscFunctionReturn(0);
86c58f1c22SToby Isaac }
87c58f1c22SToby Isaac 
88c58f1c22SToby Isaac /*
89c58f1c22SToby Isaac   DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format
90c58f1c22SToby Isaac 
91c58f1c22SToby Isaac   Input parameter:
92c58f1c22SToby Isaac . label - The DMLabel
93c58f1c22SToby Isaac 
94c58f1c22SToby Isaac   Output parameter:
95c58f1c22SToby Isaac . label - The DMLabel with all strata in sorted list format
96c58f1c22SToby Isaac 
97c58f1c22SToby Isaac   Level: developer
98c58f1c22SToby Isaac 
99c58f1c22SToby Isaac .seealso: DMLabelCreate()
100c58f1c22SToby Isaac */
101c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label)
102c58f1c22SToby Isaac {
103c58f1c22SToby Isaac   PetscInt       v;
104c58f1c22SToby Isaac   PetscErrorCode ierr;
105c58f1c22SToby Isaac 
106c58f1c22SToby Isaac   PetscFunctionBegin;
107c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; v++){
108c58f1c22SToby Isaac     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
109c58f1c22SToby Isaac   }
110c58f1c22SToby Isaac   PetscFunctionReturn(0);
111c58f1c22SToby Isaac }
112c58f1c22SToby Isaac 
113c58f1c22SToby Isaac /*
114c58f1c22SToby Isaac   DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format
115c58f1c22SToby Isaac 
116c58f1c22SToby Isaac   Input parameter:
117c58f1c22SToby Isaac + label - The DMLabel
118c58f1c22SToby Isaac - v - The stratum value
119c58f1c22SToby Isaac 
120c58f1c22SToby Isaac   Output parameter:
121c58f1c22SToby Isaac . label - The DMLabel with stratum in hash format
122c58f1c22SToby Isaac 
123c58f1c22SToby Isaac   Level: developer
124c58f1c22SToby Isaac 
125c58f1c22SToby Isaac .seealso: DMLabelCreate()
126c58f1c22SToby Isaac */
127c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v)
128c58f1c22SToby Isaac {
129c58f1c22SToby Isaac   PetscInt       p;
130ad8374ffSToby Isaac   const PetscInt *points;
131c58f1c22SToby Isaac   PetscErrorCode ierr;
132c58f1c22SToby Isaac 
133b9907514SLisandro Dalcin   if (PetscLikely(v >= 0 && v < label->numStrata) && !label->validIS[v]) return 0;
134c58f1c22SToby Isaac   PetscFunctionBegin;
1350c3c4a36SLisandro 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);
136ad8374ffSToby Isaac   if (label->points[v]) {
137ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
138e8f14785SLisandro Dalcin     for (p = 0; p < label->stratumSizes[v]; ++p) {
139e8f14785SLisandro Dalcin       ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);
140e8f14785SLisandro Dalcin     }
141ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
142ad8374ffSToby Isaac     ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
143ad8374ffSToby Isaac   }
144ad8374ffSToby Isaac   label->validIS[v] = PETSC_FALSE;
145c58f1c22SToby Isaac   PetscFunctionReturn(0);
146c58f1c22SToby Isaac }
147c58f1c22SToby Isaac 
148b9907514SLisandro Dalcin #if !defined(DMLABEL_LOOKUP_THRESHOLD)
149b9907514SLisandro Dalcin #define DMLABEL_LOOKUP_THRESHOLD 16
150b9907514SLisandro Dalcin #endif
151b9907514SLisandro Dalcin 
1520c3c4a36SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index)
1530c3c4a36SLisandro Dalcin {
1540c3c4a36SLisandro Dalcin   PetscInt       v;
155b9907514SLisandro Dalcin   PetscErrorCode ierr;
1560e79e033SBarry Smith 
1570c3c4a36SLisandro Dalcin   PetscFunctionBegin;
1580e79e033SBarry Smith   *index = -1;
159b9907514SLisandro Dalcin   if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) {
160b9907514SLisandro Dalcin     for (v = 0; v < label->numStrata; ++v)
161b9907514SLisandro Dalcin       if (label->stratumValues[v] == value) {*index = v; break;}
162b9907514SLisandro Dalcin   } else {
163b9907514SLisandro Dalcin     ierr = PetscHMapIGet(label->hmap, value, index);CHKERRQ(ierr);
1640e79e033SBarry Smith   }
16590e9b2aeSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
16690e9b2aeSLisandro Dalcin   { /* Check strata hash map consistency */
16790e9b2aeSLisandro Dalcin     PetscInt len, loc = -1;
16890e9b2aeSLisandro Dalcin     ierr = PetscHMapIGetSize(label->hmap, &len);CHKERRQ(ierr);
16990e9b2aeSLisandro Dalcin     if (len != label->numStrata) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map size");
17090e9b2aeSLisandro Dalcin     if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) {
17190e9b2aeSLisandro Dalcin       ierr = PetscHMapIGet(label->hmap, value, &loc);CHKERRQ(ierr);
17290e9b2aeSLisandro Dalcin     } else {
17390e9b2aeSLisandro Dalcin       for (v = 0; v < label->numStrata; ++v)
17490e9b2aeSLisandro Dalcin         if (label->stratumValues[v] == value) {loc = v; break;}
17590e9b2aeSLisandro Dalcin     }
17690e9b2aeSLisandro Dalcin     if (loc != *index) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map lookup");
17790e9b2aeSLisandro Dalcin   }
17890e9b2aeSLisandro Dalcin #endif
1790c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
1800c3c4a36SLisandro Dalcin }
1810c3c4a36SLisandro Dalcin 
182b9907514SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index)
183c58f1c22SToby Isaac {
184b9907514SLisandro Dalcin   PetscInt       v;
185b9907514SLisandro Dalcin   PetscInt      *tmpV;
186b9907514SLisandro Dalcin   PetscInt      *tmpS;
187b9907514SLisandro Dalcin   PetscHSetI    *tmpH, ht;
188b9907514SLisandro Dalcin   IS            *tmpP, is;
189c58f1c22SToby Isaac   PetscBool     *tmpB;
190b9907514SLisandro Dalcin   PetscHMapI     hmap = label->hmap;
191c58f1c22SToby Isaac   PetscErrorCode ierr;
192c58f1c22SToby Isaac 
193c58f1c22SToby Isaac   PetscFunctionBegin;
194b9907514SLisandro Dalcin   v    = label->numStrata;
195b9907514SLisandro Dalcin   tmpV = label->stratumValues;
196b9907514SLisandro Dalcin   tmpS = label->stratumSizes;
197b9907514SLisandro Dalcin   tmpH = label->ht;
198b9907514SLisandro Dalcin   tmpP = label->points;
199b9907514SLisandro Dalcin   tmpB = label->validIS;
2008e1f8cf0SLisandro Dalcin   { /* TODO: PetscRealloc() is broken, use malloc+memcpy+free  */
2018e1f8cf0SLisandro Dalcin     PetscInt   *oldV = tmpV;
2028e1f8cf0SLisandro Dalcin     PetscInt   *oldS = tmpS;
2038e1f8cf0SLisandro Dalcin     PetscHSetI *oldH = tmpH;
2048e1f8cf0SLisandro Dalcin     IS         *oldP = tmpP;
2058e1f8cf0SLisandro Dalcin     PetscBool  *oldB = tmpB;
2068e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpV), &tmpV);CHKERRQ(ierr);
2078e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpS), &tmpS);CHKERRQ(ierr);
2088e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpH), &tmpH);CHKERRQ(ierr);
2098e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpP), &tmpP);CHKERRQ(ierr);
2108e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpB), &tmpB);CHKERRQ(ierr);
2118e1f8cf0SLisandro Dalcin     ierr = PetscMemcpy(tmpV, oldV, v*sizeof(*tmpV));CHKERRQ(ierr);
2128e1f8cf0SLisandro Dalcin     ierr = PetscMemcpy(tmpS, oldS, v*sizeof(*tmpS));CHKERRQ(ierr);
2138e1f8cf0SLisandro Dalcin     ierr = PetscMemcpy(tmpH, oldH, v*sizeof(*tmpH));CHKERRQ(ierr);
2148e1f8cf0SLisandro Dalcin     ierr = PetscMemcpy(tmpP, oldP, v*sizeof(*tmpP));CHKERRQ(ierr);
2158e1f8cf0SLisandro Dalcin     ierr = PetscMemcpy(tmpB, oldB, v*sizeof(*tmpB));CHKERRQ(ierr);
2168e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldV);CHKERRQ(ierr);
2178e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldS);CHKERRQ(ierr);
2188e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldH);CHKERRQ(ierr);
2198e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldP);CHKERRQ(ierr);
2208e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldB);CHKERRQ(ierr);
2218e1f8cf0SLisandro Dalcin   }
222b9907514SLisandro Dalcin   label->numStrata     = v+1;
223c58f1c22SToby Isaac   label->stratumValues = tmpV;
224c58f1c22SToby Isaac   label->stratumSizes  = tmpS;
225c58f1c22SToby Isaac   label->ht            = tmpH;
226c58f1c22SToby Isaac   label->points        = tmpP;
227ad8374ffSToby Isaac   label->validIS       = tmpB;
228b9907514SLisandro Dalcin   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
229b9907514SLisandro Dalcin   ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);CHKERRQ(ierr);
230b9907514SLisandro Dalcin   ierr = PetscHMapISet(hmap, value, v);CHKERRQ(ierr);
231b9907514SLisandro Dalcin   tmpV[v] = value;
232b9907514SLisandro Dalcin   tmpS[v] = 0;
233b9907514SLisandro Dalcin   tmpH[v] = ht;
234b9907514SLisandro Dalcin   tmpP[v] = is;
235b9907514SLisandro Dalcin   tmpB[v] = PETSC_TRUE;
236*277ea44aSLisandro Dalcin   ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
2370c3c4a36SLisandro Dalcin   *index = v;
2380c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
2390c3c4a36SLisandro Dalcin }
2400c3c4a36SLisandro Dalcin 
241b9907514SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupAddStratum(DMLabel label, PetscInt value, PetscInt *index)
242b9907514SLisandro Dalcin {
243b9907514SLisandro Dalcin   PetscErrorCode ierr;
244b9907514SLisandro Dalcin   PetscFunctionBegin;
245b9907514SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, index);CHKERRQ(ierr);
246b9907514SLisandro Dalcin   if (*index < 0) {ierr = DMLabelNewStratum(label, value, index);CHKERRQ(ierr);}
247b9907514SLisandro Dalcin   PetscFunctionReturn(0);
248b9907514SLisandro Dalcin }
249b9907514SLisandro Dalcin 
250b9907514SLisandro Dalcin /*@
251b9907514SLisandro Dalcin   DMLabelAddStratum - Adds a new stratum value in a DMLabel
252b9907514SLisandro Dalcin 
253b9907514SLisandro Dalcin   Input Parameter:
254b9907514SLisandro Dalcin + label - The DMLabel
255b9907514SLisandro Dalcin - value - The stratum value
256b9907514SLisandro Dalcin 
257b9907514SLisandro Dalcin   Level: beginner
258b9907514SLisandro Dalcin 
259b9907514SLisandro Dalcin .seealso:  DMLabelCreate(), DMLabelDestroy()
260b9907514SLisandro Dalcin @*/
2610c3c4a36SLisandro Dalcin PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
2620c3c4a36SLisandro Dalcin {
2630c3c4a36SLisandro Dalcin   PetscInt       v;
2640c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
2650c3c4a36SLisandro Dalcin 
2660c3c4a36SLisandro Dalcin   PetscFunctionBegin;
267d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
268b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
269b9907514SLisandro Dalcin   PetscFunctionReturn(0);
270b9907514SLisandro Dalcin }
271b9907514SLisandro Dalcin 
272b9907514SLisandro Dalcin /*@
273b9907514SLisandro Dalcin   DMLabelAddStrata - Adds new stratum values in a DMLabel
274b9907514SLisandro Dalcin 
275b9907514SLisandro Dalcin   Input Parameter:
276b9907514SLisandro Dalcin + label - The DMLabel
277b9907514SLisandro Dalcin . numStrata - The number of stratum values
278b9907514SLisandro Dalcin - stratumValues - The stratum values
279b9907514SLisandro Dalcin 
280b9907514SLisandro Dalcin   Level: beginner
281b9907514SLisandro Dalcin 
282b9907514SLisandro Dalcin .seealso:  DMLabelCreate(), DMLabelDestroy()
283b9907514SLisandro Dalcin @*/
284b9907514SLisandro Dalcin PetscErrorCode DMLabelAddStrata(DMLabel label, PetscInt numStrata, const PetscInt stratumValues[])
285b9907514SLisandro Dalcin {
286b9907514SLisandro Dalcin   PetscInt       *values, v;
287b9907514SLisandro Dalcin   PetscErrorCode ierr;
288b9907514SLisandro Dalcin 
289b9907514SLisandro Dalcin   PetscFunctionBegin;
290b9907514SLisandro Dalcin   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
291b9907514SLisandro Dalcin   if (numStrata) PetscValidIntPointer(stratumValues, 3);
292b9907514SLisandro Dalcin   ierr = PetscMalloc1(numStrata, &values);CHKERRQ(ierr);
293b9907514SLisandro Dalcin   ierr = PetscMemcpy(values, stratumValues, numStrata*sizeof(PetscInt));CHKERRQ(ierr);
294b9907514SLisandro Dalcin   ierr = PetscSortRemoveDupsInt(&numStrata, values);CHKERRQ(ierr);
295b9907514SLisandro Dalcin   if (!label->numStrata) { /* Fast preallocation */
296b9907514SLisandro Dalcin     PetscInt   *tmpV;
297b9907514SLisandro Dalcin     PetscInt   *tmpS;
298b9907514SLisandro Dalcin     PetscHSetI *tmpH, ht;
299b9907514SLisandro Dalcin     IS         *tmpP, is;
300b9907514SLisandro Dalcin     PetscBool  *tmpB;
301b9907514SLisandro Dalcin     PetscHMapI  hmap = label->hmap;
302b9907514SLisandro Dalcin 
303b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpV);CHKERRQ(ierr);
304b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpS);CHKERRQ(ierr);
305b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpH);CHKERRQ(ierr);
306b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpP);CHKERRQ(ierr);
307b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpB);CHKERRQ(ierr);
308b9907514SLisandro Dalcin     label->numStrata     = numStrata;
309b9907514SLisandro Dalcin     label->stratumValues = tmpV;
310b9907514SLisandro Dalcin     label->stratumSizes  = tmpS;
311b9907514SLisandro Dalcin     label->ht            = tmpH;
312b9907514SLisandro Dalcin     label->points        = tmpP;
313b9907514SLisandro Dalcin     label->validIS       = tmpB;
314b9907514SLisandro Dalcin     for (v = 0; v < numStrata; ++v) {
315b9907514SLisandro Dalcin       ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
316b9907514SLisandro Dalcin       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);CHKERRQ(ierr);
317b9907514SLisandro Dalcin       ierr = PetscHMapISet(hmap, values[v], v);CHKERRQ(ierr);
318b9907514SLisandro Dalcin       tmpV[v] = values[v];
319b9907514SLisandro Dalcin       tmpS[v] = 0;
320b9907514SLisandro Dalcin       tmpH[v] = ht;
321b9907514SLisandro Dalcin       tmpP[v] = is;
322b9907514SLisandro Dalcin       tmpB[v] = PETSC_TRUE;
323b9907514SLisandro Dalcin     }
324*277ea44aSLisandro Dalcin     ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
325b9907514SLisandro Dalcin   } else {
326b9907514SLisandro Dalcin     for (v = 0; v < numStrata; ++v) {
327b9907514SLisandro Dalcin       ierr = DMLabelAddStratum(label, values[v]);CHKERRQ(ierr);
328b9907514SLisandro Dalcin     }
329b9907514SLisandro Dalcin   }
330b9907514SLisandro Dalcin   ierr = PetscFree(values);CHKERRQ(ierr);
331b9907514SLisandro Dalcin   PetscFunctionReturn(0);
332b9907514SLisandro Dalcin }
333b9907514SLisandro Dalcin 
334b9907514SLisandro Dalcin /*@
335b9907514SLisandro Dalcin   DMLabelAddStrataIS - Adds new stratum values in a DMLabel
336b9907514SLisandro Dalcin 
337b9907514SLisandro Dalcin   Input Parameter:
338b9907514SLisandro Dalcin + label - The DMLabel
339b9907514SLisandro Dalcin - valueIS - Index set with stratum values
340b9907514SLisandro Dalcin 
341b9907514SLisandro Dalcin   Level: beginner
342b9907514SLisandro Dalcin 
343b9907514SLisandro Dalcin .seealso:  DMLabelCreate(), DMLabelDestroy()
344b9907514SLisandro Dalcin @*/
345b9907514SLisandro Dalcin PetscErrorCode DMLabelAddStrataIS(DMLabel label, IS valueIS)
346b9907514SLisandro Dalcin {
347b9907514SLisandro Dalcin   PetscInt       numStrata;
348b9907514SLisandro Dalcin   const PetscInt *stratumValues;
349b9907514SLisandro Dalcin   PetscErrorCode ierr;
350b9907514SLisandro Dalcin 
351b9907514SLisandro Dalcin   PetscFunctionBegin;
352b9907514SLisandro Dalcin   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
353b9907514SLisandro Dalcin   PetscValidHeaderSpecific(valueIS, IS_CLASSID, 2);
354b9907514SLisandro Dalcin   ierr = ISGetLocalSize(valueIS, &numStrata);CHKERRQ(ierr);
355b9907514SLisandro Dalcin   ierr = ISGetIndices(valueIS, &stratumValues);CHKERRQ(ierr);
356b9907514SLisandro Dalcin   ierr = DMLabelAddStrata(label, numStrata, stratumValues);CHKERRQ(ierr);
357c58f1c22SToby Isaac   PetscFunctionReturn(0);
358c58f1c22SToby Isaac }
359c58f1c22SToby Isaac 
360c58f1c22SToby Isaac static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
361c58f1c22SToby Isaac {
362c58f1c22SToby Isaac   PetscInt       v;
363c58f1c22SToby Isaac   PetscMPIInt    rank;
364c58f1c22SToby Isaac   PetscErrorCode ierr;
365c58f1c22SToby Isaac 
366c58f1c22SToby Isaac   PetscFunctionBegin;
367c58f1c22SToby Isaac   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr);
368c58f1c22SToby Isaac   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
369c58f1c22SToby Isaac   if (label) {
370d67d17b1SMatthew G. Knepley     const char *name;
371d67d17b1SMatthew G. Knepley 
372d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
373d67d17b1SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", name);CHKERRQ(ierr);
374c58f1c22SToby Isaac     if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);}
375c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
376c58f1c22SToby Isaac       const PetscInt value = label->stratumValues[v];
377ad8374ffSToby Isaac       const PetscInt *points;
378c58f1c22SToby Isaac       PetscInt       p;
379c58f1c22SToby Isaac 
380ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
381c58f1c22SToby Isaac       for (p = 0; p < label->stratumSizes[v]; ++p) {
382ad8374ffSToby Isaac         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr);
383c58f1c22SToby Isaac       }
384ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
385c58f1c22SToby Isaac     }
386c58f1c22SToby Isaac   }
387c58f1c22SToby Isaac   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
388c58f1c22SToby Isaac   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
389c58f1c22SToby Isaac   PetscFunctionReturn(0);
390c58f1c22SToby Isaac }
391c58f1c22SToby Isaac 
392c58f1c22SToby Isaac /*@C
393c58f1c22SToby Isaac   DMLabelView - View the label
394c58f1c22SToby Isaac 
395c58f1c22SToby Isaac   Input Parameters:
396c58f1c22SToby Isaac + label - The DMLabel
397c58f1c22SToby Isaac - viewer - The PetscViewer
398c58f1c22SToby Isaac 
399c58f1c22SToby Isaac   Level: intermediate
400c58f1c22SToby Isaac 
401c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy()
402c58f1c22SToby Isaac @*/
403c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
404c58f1c22SToby Isaac {
405c58f1c22SToby Isaac   PetscBool      iascii;
406c58f1c22SToby Isaac   PetscErrorCode ierr;
407c58f1c22SToby Isaac 
408c58f1c22SToby Isaac   PetscFunctionBegin;
409d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
410b9907514SLisandro Dalcin   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)label), &viewer);CHKERRQ(ierr);}
411c58f1c22SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
412c58f1c22SToby Isaac   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
413c58f1c22SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
414c58f1c22SToby Isaac   if (iascii) {
415c58f1c22SToby Isaac     ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr);
416c58f1c22SToby Isaac   }
417c58f1c22SToby Isaac   PetscFunctionReturn(0);
418c58f1c22SToby Isaac }
419c58f1c22SToby Isaac 
42084f0b6dfSMatthew G. Knepley /*@
421d67d17b1SMatthew G. Knepley   DMLabelReset - Destroys internal data structures in a DMLabel
422d67d17b1SMatthew G. Knepley 
423d67d17b1SMatthew G. Knepley   Input Parameter:
424d67d17b1SMatthew G. Knepley . label - The DMLabel
425d67d17b1SMatthew G. Knepley 
426d67d17b1SMatthew G. Knepley   Level: beginner
427d67d17b1SMatthew G. Knepley 
428d67d17b1SMatthew G. Knepley .seealso: DMLabelDestroy(), DMLabelCreate()
429d67d17b1SMatthew G. Knepley @*/
430d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelReset(DMLabel label)
431d67d17b1SMatthew G. Knepley {
432d67d17b1SMatthew G. Knepley   PetscInt       v;
433d67d17b1SMatthew G. Knepley   PetscErrorCode ierr;
434d67d17b1SMatthew G. Knepley 
435d67d17b1SMatthew G. Knepley   PetscFunctionBegin;
436d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
437d67d17b1SMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
438d67d17b1SMatthew G. Knepley     ierr = PetscHSetIDestroy(&label->ht[v]);CHKERRQ(ierr);
439d67d17b1SMatthew G. Knepley     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
440d67d17b1SMatthew G. Knepley   }
441b9907514SLisandro Dalcin   label->numStrata = 0;
442b9907514SLisandro Dalcin   ierr = PetscFree(label->stratumValues);CHKERRQ(ierr);
443b9907514SLisandro Dalcin   ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr);
444d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->ht);CHKERRQ(ierr);
445d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->points);CHKERRQ(ierr);
446d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->validIS);CHKERRQ(ierr);
447b9907514SLisandro Dalcin   ierr = PetscHMapIReset(label->hmap);CHKERRQ(ierr);
448b9907514SLisandro Dalcin   label->pStart = -1;
449b9907514SLisandro Dalcin   label->pEnd   = -1;
450d67d17b1SMatthew G. Knepley   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
451d67d17b1SMatthew G. Knepley   PetscFunctionReturn(0);
452d67d17b1SMatthew G. Knepley }
453d67d17b1SMatthew G. Knepley 
454d67d17b1SMatthew G. Knepley /*@
45584f0b6dfSMatthew G. Knepley   DMLabelDestroy - Destroys a DMLabel
45684f0b6dfSMatthew G. Knepley 
45784f0b6dfSMatthew G. Knepley   Input Parameter:
45884f0b6dfSMatthew G. Knepley . label - The DMLabel
45984f0b6dfSMatthew G. Knepley 
46084f0b6dfSMatthew G. Knepley   Level: beginner
46184f0b6dfSMatthew G. Knepley 
462d67d17b1SMatthew G. Knepley .seealso: DMLabelReset(), DMLabelCreate()
46384f0b6dfSMatthew G. Knepley @*/
464c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label)
465c58f1c22SToby Isaac {
466c58f1c22SToby Isaac   PetscErrorCode ierr;
467c58f1c22SToby Isaac 
468c58f1c22SToby Isaac   PetscFunctionBegin;
469d67d17b1SMatthew G. Knepley   if (!*label) PetscFunctionReturn(0);
470d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific((*label),DMLABEL_CLASSID,1);
471b9907514SLisandro Dalcin   if (--((PetscObject)(*label))->refct > 0) {*label = NULL; PetscFunctionReturn(0);}
472d67d17b1SMatthew G. Knepley   ierr = DMLabelReset(*label);CHKERRQ(ierr);
473b9907514SLisandro Dalcin   ierr = PetscHMapIDestroy(&(*label)->hmap);CHKERRQ(ierr);
474d67d17b1SMatthew G. Knepley   ierr = PetscHeaderDestroy(label);CHKERRQ(ierr);
475c58f1c22SToby Isaac   PetscFunctionReturn(0);
476c58f1c22SToby Isaac }
477c58f1c22SToby Isaac 
47884f0b6dfSMatthew G. Knepley /*@
47984f0b6dfSMatthew G. Knepley   DMLabelDuplicate - Duplicates a DMLabel
48084f0b6dfSMatthew G. Knepley 
48184f0b6dfSMatthew G. Knepley   Input Parameter:
48284f0b6dfSMatthew G. Knepley . label - The DMLabel
48384f0b6dfSMatthew G. Knepley 
48484f0b6dfSMatthew G. Knepley   Output Parameter:
48584f0b6dfSMatthew G. Knepley . labelnew - location to put new vector
48684f0b6dfSMatthew G. Knepley 
48784f0b6dfSMatthew G. Knepley   Level: intermediate
48884f0b6dfSMatthew G. Knepley 
48984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelDestroy()
49084f0b6dfSMatthew G. Knepley @*/
491c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
492c58f1c22SToby Isaac {
493d67d17b1SMatthew G. Knepley   const char    *name;
494ad8374ffSToby Isaac   PetscInt       v;
495c58f1c22SToby Isaac   PetscErrorCode ierr;
496c58f1c22SToby Isaac 
497c58f1c22SToby Isaac   PetscFunctionBegin;
498d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
499c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
500d67d17b1SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
501d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PetscObjectComm((PetscObject) label), name, labelnew);CHKERRQ(ierr);
502c58f1c22SToby Isaac 
503c58f1c22SToby Isaac   (*labelnew)->numStrata    = label->numStrata;
5045aa44df4SToby Isaac   (*labelnew)->defaultValue = label->defaultValue;
505c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr);
506c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr);
507c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr);
508c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr);
509ad8374ffSToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr);
510c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
511e8f14785SLisandro Dalcin     ierr = PetscHSetICreate(&(*labelnew)->ht[v]);CHKERRQ(ierr);
512c58f1c22SToby Isaac     (*labelnew)->stratumValues[v]  = label->stratumValues[v];
513c58f1c22SToby Isaac     (*labelnew)->stratumSizes[v]   = label->stratumSizes[v];
514ad8374ffSToby Isaac     ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr);
515ad8374ffSToby Isaac     (*labelnew)->points[v]         = label->points[v];
516b9907514SLisandro Dalcin     (*labelnew)->validIS[v]        = PETSC_TRUE;
517c58f1c22SToby Isaac   }
518f14fe40dSLisandro Dalcin   ierr = PetscHMapIDestroy(&(*labelnew)->hmap);CHKERRQ(ierr);
519b9907514SLisandro Dalcin   ierr = PetscHMapIDuplicate(label->hmap,&(*labelnew)->hmap);CHKERRQ(ierr);
520c58f1c22SToby Isaac   (*labelnew)->pStart = -1;
521c58f1c22SToby Isaac   (*labelnew)->pEnd   = -1;
522c58f1c22SToby Isaac   (*labelnew)->bt     = NULL;
523c58f1c22SToby Isaac   PetscFunctionReturn(0);
524c58f1c22SToby Isaac }
525c58f1c22SToby Isaac 
526c6a43d28SMatthew G. Knepley /*@
527c6a43d28SMatthew G. Knepley   DMLabelComputeIndex - Create an index structure for membership determination, automatically determining the bounds
528c6a43d28SMatthew G. Knepley 
529c6a43d28SMatthew G. Knepley   Input Parameter:
530c6a43d28SMatthew G. Knepley . label  - The DMLabel
531c6a43d28SMatthew G. Knepley 
532c6a43d28SMatthew G. Knepley   Level: intermediate
533c6a43d28SMatthew G. Knepley 
534c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
535c6a43d28SMatthew G. Knepley @*/
536c6a43d28SMatthew G. Knepley PetscErrorCode DMLabelComputeIndex(DMLabel label)
537c6a43d28SMatthew G. Knepley {
538c6a43d28SMatthew G. Knepley   PetscInt       pStart = PETSC_MAX_INT, pEnd = -1, v;
539c6a43d28SMatthew G. Knepley   PetscErrorCode ierr;
540c6a43d28SMatthew G. Knepley 
541c6a43d28SMatthew G. Knepley   PetscFunctionBegin;
542c6a43d28SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
543c6a43d28SMatthew G. Knepley   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
544c6a43d28SMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
545c6a43d28SMatthew G. Knepley     const PetscInt *points;
546c6a43d28SMatthew G. Knepley     PetscInt       i;
547c6a43d28SMatthew G. Knepley 
548c6a43d28SMatthew G. Knepley     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
549c6a43d28SMatthew G. Knepley     for (i = 0; i < label->stratumSizes[v]; ++i) {
550c6a43d28SMatthew G. Knepley       const PetscInt point = points[i];
551c6a43d28SMatthew G. Knepley 
552c6a43d28SMatthew G. Knepley       pStart = PetscMin(point,   pStart);
553c6a43d28SMatthew G. Knepley       pEnd   = PetscMax(point+1, pEnd);
554c6a43d28SMatthew G. Knepley     }
555c6a43d28SMatthew G. Knepley     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
556c6a43d28SMatthew G. Knepley   }
557c6a43d28SMatthew G. Knepley   label->pStart = pStart == PETSC_MAX_INT ? -1 : pStart;
558c6a43d28SMatthew G. Knepley   label->pEnd   = pEnd;
559c6a43d28SMatthew G. Knepley   ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
560c6a43d28SMatthew G. Knepley   PetscFunctionReturn(0);
561c6a43d28SMatthew G. Knepley }
562c6a43d28SMatthew G. Knepley 
563c6a43d28SMatthew G. Knepley /*@
564c6a43d28SMatthew G. Knepley   DMLabelCreateIndex - Create an index structure for membership determination
565c6a43d28SMatthew G. Knepley 
566c6a43d28SMatthew G. Knepley   Input Parameters:
567c6a43d28SMatthew G. Knepley + label  - The DMLabel
568c6a43d28SMatthew G. Knepley . pStart - The smallest point
569c6a43d28SMatthew G. Knepley - pEnd   - The largest point + 1
570c6a43d28SMatthew G. Knepley 
571c6a43d28SMatthew G. Knepley   Level: intermediate
572c6a43d28SMatthew G. Knepley 
573c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelComputeIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
574c6a43d28SMatthew G. Knepley @*/
575c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
576c58f1c22SToby Isaac {
577c58f1c22SToby Isaac   PetscInt       v;
578c58f1c22SToby Isaac   PetscErrorCode ierr;
579c58f1c22SToby Isaac 
580c58f1c22SToby Isaac   PetscFunctionBegin;
581d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
5820c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
583c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
584c58f1c22SToby Isaac   label->pStart = pStart;
585c58f1c22SToby Isaac   label->pEnd   = pEnd;
586c6a43d28SMatthew G. Knepley   /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
587c58f1c22SToby Isaac   ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr);
588c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
589ad8374ffSToby Isaac     const PetscInt *points;
590c58f1c22SToby Isaac     PetscInt       i;
591c58f1c22SToby Isaac 
592ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
593c58f1c22SToby Isaac     for (i = 0; i < label->stratumSizes[v]; ++i) {
594ad8374ffSToby Isaac       const PetscInt point = points[i];
595c58f1c22SToby Isaac 
596c58f1c22SToby 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);
597c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr);
598c58f1c22SToby Isaac     }
599ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
600c58f1c22SToby Isaac   }
601c58f1c22SToby Isaac   PetscFunctionReturn(0);
602c58f1c22SToby Isaac }
603c58f1c22SToby Isaac 
604c6a43d28SMatthew G. Knepley /*@
605c6a43d28SMatthew G. Knepley   DMLabelDestroyIndex - Destroy the index structure
606c6a43d28SMatthew G. Knepley 
607c6a43d28SMatthew G. Knepley   Input Parameter:
608c6a43d28SMatthew G. Knepley . label - the DMLabel
609c6a43d28SMatthew G. Knepley 
610c6a43d28SMatthew G. Knepley   Level: intermediate
611c6a43d28SMatthew G. Knepley 
612c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
613c6a43d28SMatthew G. Knepley @*/
614c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label)
615c58f1c22SToby Isaac {
616c58f1c22SToby Isaac   PetscErrorCode ierr;
617c58f1c22SToby Isaac 
618c58f1c22SToby Isaac   PetscFunctionBegin;
619d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
620c58f1c22SToby Isaac   label->pStart = -1;
621c58f1c22SToby Isaac   label->pEnd   = -1;
6220c3c4a36SLisandro Dalcin   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
623c58f1c22SToby Isaac   PetscFunctionReturn(0);
624c58f1c22SToby Isaac }
625c58f1c22SToby Isaac 
626c58f1c22SToby Isaac /*@
627c6a43d28SMatthew G. Knepley   DMLabelGetBounds - Return the smallest and largest point in the label
628c6a43d28SMatthew G. Knepley 
629c6a43d28SMatthew G. Knepley   Input Parameter:
630c6a43d28SMatthew G. Knepley . label - the DMLabel
631c6a43d28SMatthew G. Knepley 
632c6a43d28SMatthew G. Knepley   Output Parameters:
633c6a43d28SMatthew G. Knepley + pStart - The smallest point
634c6a43d28SMatthew G. Knepley - pEnd   - The largest point + 1
635c6a43d28SMatthew G. Knepley 
636c6a43d28SMatthew G. Knepley   Note: This will compute an index for the label if one does not exist.
637c6a43d28SMatthew G. Knepley 
638c6a43d28SMatthew G. Knepley   Level: intermediate
639c6a43d28SMatthew G. Knepley 
640c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
641c6a43d28SMatthew G. Knepley @*/
642c6a43d28SMatthew G. Knepley PetscErrorCode DMLabelGetBounds(DMLabel label, PetscInt *pStart, PetscInt *pEnd)
643c6a43d28SMatthew G. Knepley {
644c6a43d28SMatthew G. Knepley   PetscErrorCode ierr;
645c6a43d28SMatthew G. Knepley 
646c6a43d28SMatthew G. Knepley   PetscFunctionBegin;
647c6a43d28SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
648c6a43d28SMatthew G. Knepley   if ((label->pStart == -1) && (label->pEnd == -1)) {ierr = DMLabelComputeIndex(label);CHKERRQ(ierr);}
649c6a43d28SMatthew G. Knepley   if (pStart) {
650c6a43d28SMatthew G. Knepley     PetscValidPointer(pStart, 2);
651c6a43d28SMatthew G. Knepley     *pStart = label->pStart;
652c6a43d28SMatthew G. Knepley   }
653c6a43d28SMatthew G. Knepley   if (pEnd) {
654c6a43d28SMatthew G. Knepley     PetscValidPointer(pEnd, 3);
655c6a43d28SMatthew G. Knepley     *pEnd = label->pEnd;
656c6a43d28SMatthew G. Knepley   }
657c6a43d28SMatthew G. Knepley   PetscFunctionReturn(0);
658c6a43d28SMatthew G. Knepley }
659c6a43d28SMatthew G. Knepley 
660c6a43d28SMatthew G. Knepley /*@
661c58f1c22SToby Isaac   DMLabelHasValue - Determine whether a label assigns the value to any point
662c58f1c22SToby Isaac 
663c58f1c22SToby Isaac   Input Parameters:
664c58f1c22SToby Isaac + label - the DMLabel
665c58f1c22SToby Isaac - value - the value
666c58f1c22SToby Isaac 
667c58f1c22SToby Isaac   Output Parameter:
668c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point
669c58f1c22SToby Isaac 
670c58f1c22SToby Isaac   Level: developer
671c58f1c22SToby Isaac 
672c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
673c58f1c22SToby Isaac @*/
674c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
675c58f1c22SToby Isaac {
676c58f1c22SToby Isaac   PetscInt v;
6770c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
678c58f1c22SToby Isaac 
679c58f1c22SToby Isaac   PetscFunctionBegin;
680d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
681c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
6820c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
6830c3c4a36SLisandro Dalcin   *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE;
684c58f1c22SToby Isaac   PetscFunctionReturn(0);
685c58f1c22SToby Isaac }
686c58f1c22SToby Isaac 
687c58f1c22SToby Isaac /*@
688c58f1c22SToby Isaac   DMLabelHasPoint - Determine whether a label assigns a value to a point
689c58f1c22SToby Isaac 
690c58f1c22SToby Isaac   Input Parameters:
691c58f1c22SToby Isaac + label - the DMLabel
692c58f1c22SToby Isaac - point - the point
693c58f1c22SToby Isaac 
694c58f1c22SToby Isaac   Output Parameter:
695c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value
696c58f1c22SToby Isaac 
697c58f1c22SToby Isaac   Note: The user must call DMLabelCreateIndex() before this function.
698c58f1c22SToby Isaac 
699c58f1c22SToby Isaac   Level: developer
700c58f1c22SToby Isaac 
701c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
702c58f1c22SToby Isaac @*/
703c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
704c58f1c22SToby Isaac {
705c58f1c22SToby Isaac   PetscErrorCode ierr;
706c58f1c22SToby Isaac 
707c58f1c22SToby Isaac   PetscFunctionBeginHot;
708d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
709c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
710c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
711c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG)
712c58f1c22SToby Isaac   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
713c58f1c22SToby 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);
714c58f1c22SToby Isaac #endif
715c58f1c22SToby Isaac   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
716c58f1c22SToby Isaac   PetscFunctionReturn(0);
717c58f1c22SToby Isaac }
718c58f1c22SToby Isaac 
719c58f1c22SToby Isaac /*@
720c58f1c22SToby Isaac   DMLabelStratumHasPoint - Return true if the stratum contains a point
721c58f1c22SToby Isaac 
722c58f1c22SToby Isaac   Input Parameters:
723c58f1c22SToby Isaac + label - the DMLabel
724c58f1c22SToby Isaac . value - the stratum value
725c58f1c22SToby Isaac - point - the point
726c58f1c22SToby Isaac 
727c58f1c22SToby Isaac   Output Parameter:
728c58f1c22SToby Isaac . contains - true if the stratum contains the point
729c58f1c22SToby Isaac 
730c58f1c22SToby Isaac   Level: intermediate
731c58f1c22SToby Isaac 
732c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
733c58f1c22SToby Isaac @*/
734c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
735c58f1c22SToby Isaac {
736c58f1c22SToby Isaac   PetscInt       v;
737c58f1c22SToby Isaac   PetscErrorCode ierr;
738c58f1c22SToby Isaac 
7390c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
740d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
741c58f1c22SToby Isaac   PetscValidPointer(contains, 4);
742c58f1c22SToby Isaac   *contains = PETSC_FALSE;
7430c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
7440c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
7450c3c4a36SLisandro Dalcin 
746ad8374ffSToby Isaac   if (label->validIS[v]) {
747c58f1c22SToby Isaac     PetscInt i;
748c58f1c22SToby Isaac 
749a2d74346SToby Isaac     ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
7500c3c4a36SLisandro Dalcin     if (i >= 0) *contains = PETSC_TRUE;
751c58f1c22SToby Isaac   } else {
752c58f1c22SToby Isaac     PetscBool has;
753c58f1c22SToby Isaac 
754b9907514SLisandro Dalcin     ierr = PetscHSetIHas(label->ht[v], point, &has);CHKERRQ(ierr);
7550c3c4a36SLisandro Dalcin     if (has) *contains = PETSC_TRUE;
756c58f1c22SToby Isaac   }
757c58f1c22SToby Isaac   PetscFunctionReturn(0);
758c58f1c22SToby Isaac }
759c58f1c22SToby Isaac 
76084f0b6dfSMatthew G. Knepley /*@
7615aa44df4SToby Isaac   DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
7625aa44df4SToby Isaac   When a label is created, it is initialized to -1.
7635aa44df4SToby Isaac 
7645aa44df4SToby Isaac   Input parameter:
7655aa44df4SToby Isaac . label - a DMLabel object
7665aa44df4SToby Isaac 
7675aa44df4SToby Isaac   Output parameter:
7685aa44df4SToby Isaac . defaultValue - the default value
7695aa44df4SToby Isaac 
7705aa44df4SToby Isaac   Level: beginner
7715aa44df4SToby Isaac 
7725aa44df4SToby Isaac .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
77384f0b6dfSMatthew G. Knepley @*/
7745aa44df4SToby Isaac PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
7755aa44df4SToby Isaac {
7765aa44df4SToby Isaac   PetscFunctionBegin;
777d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
7785aa44df4SToby Isaac   *defaultValue = label->defaultValue;
7795aa44df4SToby Isaac   PetscFunctionReturn(0);
7805aa44df4SToby Isaac }
7815aa44df4SToby Isaac 
78284f0b6dfSMatthew G. Knepley /*@
7835aa44df4SToby Isaac   DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
7845aa44df4SToby Isaac   When a label is created, it is initialized to -1.
7855aa44df4SToby Isaac 
7865aa44df4SToby Isaac   Input parameter:
7875aa44df4SToby Isaac . label - a DMLabel object
7885aa44df4SToby Isaac 
7895aa44df4SToby Isaac   Output parameter:
7905aa44df4SToby Isaac . defaultValue - the default value
7915aa44df4SToby Isaac 
7925aa44df4SToby Isaac   Level: beginner
7935aa44df4SToby Isaac 
7945aa44df4SToby Isaac .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
79584f0b6dfSMatthew G. Knepley @*/
7965aa44df4SToby Isaac PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
7975aa44df4SToby Isaac {
7985aa44df4SToby Isaac   PetscFunctionBegin;
799d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
8005aa44df4SToby Isaac   label->defaultValue = defaultValue;
8015aa44df4SToby Isaac   PetscFunctionReturn(0);
8025aa44df4SToby Isaac }
8035aa44df4SToby Isaac 
804c58f1c22SToby Isaac /*@
8055aa44df4SToby 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())
806c58f1c22SToby Isaac 
807c58f1c22SToby Isaac   Input Parameters:
808c58f1c22SToby Isaac + label - the DMLabel
809c58f1c22SToby Isaac - point - the point
810c58f1c22SToby Isaac 
811c58f1c22SToby Isaac   Output Parameter:
8128e68afcfSMatthew G. Knepley . value - The point value, or the default value (-1 by default)
813c58f1c22SToby Isaac 
814c58f1c22SToby Isaac   Level: intermediate
815c58f1c22SToby Isaac 
8165aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
817c58f1c22SToby Isaac @*/
818c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
819c58f1c22SToby Isaac {
820c58f1c22SToby Isaac   PetscInt       v;
821c58f1c22SToby Isaac   PetscErrorCode ierr;
822c58f1c22SToby Isaac 
8230c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
824d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
825c58f1c22SToby Isaac   PetscValidPointer(value, 3);
8265aa44df4SToby Isaac   *value = label->defaultValue;
827c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
828ad8374ffSToby Isaac     if (label->validIS[v]) {
829c58f1c22SToby Isaac       PetscInt i;
830c58f1c22SToby Isaac 
831a2d74346SToby Isaac       ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
832c58f1c22SToby Isaac       if (i >= 0) {
833c58f1c22SToby Isaac         *value = label->stratumValues[v];
834c58f1c22SToby Isaac         break;
835c58f1c22SToby Isaac       }
836c58f1c22SToby Isaac     } else {
837c58f1c22SToby Isaac       PetscBool has;
838c58f1c22SToby Isaac 
839b9907514SLisandro Dalcin       ierr = PetscHSetIHas(label->ht[v], point, &has);CHKERRQ(ierr);
840c58f1c22SToby Isaac       if (has) {
841c58f1c22SToby Isaac         *value = label->stratumValues[v];
842c58f1c22SToby Isaac         break;
843c58f1c22SToby Isaac       }
844c58f1c22SToby Isaac     }
845c58f1c22SToby Isaac   }
846c58f1c22SToby Isaac   PetscFunctionReturn(0);
847c58f1c22SToby Isaac }
848c58f1c22SToby Isaac 
849c58f1c22SToby Isaac /*@
850367003a6SStefano 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.
851c58f1c22SToby Isaac 
852c58f1c22SToby Isaac   Input Parameters:
853c58f1c22SToby Isaac + label - the DMLabel
854c58f1c22SToby Isaac . point - the point
855c58f1c22SToby Isaac - value - The point value
856c58f1c22SToby Isaac 
857c58f1c22SToby Isaac   Level: intermediate
858c58f1c22SToby Isaac 
8595aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
860c58f1c22SToby Isaac @*/
861c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
862c58f1c22SToby Isaac {
863c58f1c22SToby Isaac   PetscInt       v;
864c58f1c22SToby Isaac   PetscErrorCode ierr;
865c58f1c22SToby Isaac 
866c58f1c22SToby Isaac   PetscFunctionBegin;
867d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
8680c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
8695aa44df4SToby Isaac   if (value == label->defaultValue) PetscFunctionReturn(0);
870b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
871c58f1c22SToby Isaac   /* Set key */
8720c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
873e8f14785SLisandro Dalcin   ierr = PetscHSetIAdd(label->ht[v], point);CHKERRQ(ierr);
874c58f1c22SToby Isaac   PetscFunctionReturn(0);
875c58f1c22SToby Isaac }
876c58f1c22SToby Isaac 
877c58f1c22SToby Isaac /*@
878c58f1c22SToby Isaac   DMLabelClearValue - Clear the value a label assigns to a point
879c58f1c22SToby Isaac 
880c58f1c22SToby Isaac   Input Parameters:
881c58f1c22SToby Isaac + label - the DMLabel
882c58f1c22SToby Isaac . point - the point
883c58f1c22SToby Isaac - value - The point value
884c58f1c22SToby Isaac 
885c58f1c22SToby Isaac   Level: intermediate
886c58f1c22SToby Isaac 
887c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
888c58f1c22SToby Isaac @*/
889c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
890c58f1c22SToby Isaac {
891ad8374ffSToby Isaac   PetscInt       v;
892c58f1c22SToby Isaac   PetscErrorCode ierr;
893c58f1c22SToby Isaac 
894c58f1c22SToby Isaac   PetscFunctionBegin;
895d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
896c58f1c22SToby Isaac   /* Find label value */
8970c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
8980c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
8990c3c4a36SLisandro Dalcin 
900eeed21e7SToby Isaac   if (label->bt) {
901eeed21e7SToby 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);
902eeed21e7SToby Isaac     ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
903eeed21e7SToby Isaac   }
9040c3c4a36SLisandro Dalcin 
9050c3c4a36SLisandro Dalcin   /* Delete key */
9060c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
907e8f14785SLisandro Dalcin   ierr = PetscHSetIDel(label->ht[v], point);CHKERRQ(ierr);
908c58f1c22SToby Isaac   PetscFunctionReturn(0);
909c58f1c22SToby Isaac }
910c58f1c22SToby Isaac 
911c58f1c22SToby Isaac /*@
912c58f1c22SToby Isaac   DMLabelInsertIS - Set all points in the IS to a value
913c58f1c22SToby Isaac 
914c58f1c22SToby Isaac   Input Parameters:
915c58f1c22SToby Isaac + label - the DMLabel
916c58f1c22SToby Isaac . is    - the point IS
917c58f1c22SToby Isaac - value - The point value
918c58f1c22SToby Isaac 
919c58f1c22SToby Isaac   Level: intermediate
920c58f1c22SToby Isaac 
921c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
922c58f1c22SToby Isaac @*/
923c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
924c58f1c22SToby Isaac {
9250c3c4a36SLisandro Dalcin   PetscInt        v, n, p;
926c58f1c22SToby Isaac   const PetscInt *points;
927c58f1c22SToby Isaac   PetscErrorCode  ierr;
928c58f1c22SToby Isaac 
929c58f1c22SToby Isaac   PetscFunctionBegin;
930d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
931c58f1c22SToby Isaac   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
9320c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
9330c3c4a36SLisandro Dalcin   if (value == label->defaultValue) PetscFunctionReturn(0);
934b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
9350c3c4a36SLisandro Dalcin   /* Set keys */
9360c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
937c58f1c22SToby Isaac   ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr);
938c58f1c22SToby Isaac   ierr = ISGetIndices(is, &points);CHKERRQ(ierr);
939e8f14785SLisandro Dalcin   for (p = 0; p < n; ++p) {ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);}
940c58f1c22SToby Isaac   ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr);
941c58f1c22SToby Isaac   PetscFunctionReturn(0);
942c58f1c22SToby Isaac }
943c58f1c22SToby Isaac 
94484f0b6dfSMatthew G. Knepley /*@
94584f0b6dfSMatthew G. Knepley   DMLabelGetNumValues - Get the number of values that the DMLabel takes
94684f0b6dfSMatthew G. Knepley 
94784f0b6dfSMatthew G. Knepley   Input Parameter:
94884f0b6dfSMatthew G. Knepley . label - the DMLabel
94984f0b6dfSMatthew G. Knepley 
95084f0b6dfSMatthew G. Knepley   Output Paramater:
95184f0b6dfSMatthew G. Knepley . numValues - the number of values
95284f0b6dfSMatthew G. Knepley 
95384f0b6dfSMatthew G. Knepley   Level: intermediate
95484f0b6dfSMatthew G. Knepley 
95584f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
95684f0b6dfSMatthew G. Knepley @*/
957c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
958c58f1c22SToby Isaac {
959c58f1c22SToby Isaac   PetscFunctionBegin;
960d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
961c58f1c22SToby Isaac   PetscValidPointer(numValues, 2);
962c58f1c22SToby Isaac   *numValues = label->numStrata;
963c58f1c22SToby Isaac   PetscFunctionReturn(0);
964c58f1c22SToby Isaac }
965c58f1c22SToby Isaac 
96684f0b6dfSMatthew G. Knepley /*@
96784f0b6dfSMatthew G. Knepley   DMLabelGetValueIS - Get an IS of all values that the DMlabel takes
96884f0b6dfSMatthew G. Knepley 
96984f0b6dfSMatthew G. Knepley   Input Parameter:
97084f0b6dfSMatthew G. Knepley . label - the DMLabel
97184f0b6dfSMatthew G. Knepley 
97284f0b6dfSMatthew G. Knepley   Output Paramater:
97384f0b6dfSMatthew G. Knepley . is    - the value IS
97484f0b6dfSMatthew G. Knepley 
97584f0b6dfSMatthew G. Knepley   Level: intermediate
97684f0b6dfSMatthew G. Knepley 
97784f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
97884f0b6dfSMatthew G. Knepley @*/
979c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
980c58f1c22SToby Isaac {
981c58f1c22SToby Isaac   PetscErrorCode ierr;
982c58f1c22SToby Isaac 
983c58f1c22SToby Isaac   PetscFunctionBegin;
984d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
985c58f1c22SToby Isaac   PetscValidPointer(values, 2);
986c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr);
987c58f1c22SToby Isaac   PetscFunctionReturn(0);
988c58f1c22SToby Isaac }
989c58f1c22SToby Isaac 
99084f0b6dfSMatthew G. Knepley /*@
99184f0b6dfSMatthew G. Knepley   DMLabelHasStratum - Determine whether points exist with the given value
99284f0b6dfSMatthew G. Knepley 
99384f0b6dfSMatthew G. Knepley   Input Parameters:
99484f0b6dfSMatthew G. Knepley + label - the DMLabel
99584f0b6dfSMatthew G. Knepley - value - the stratum value
99684f0b6dfSMatthew G. Knepley 
99784f0b6dfSMatthew G. Knepley   Output Paramater:
99884f0b6dfSMatthew G. Knepley . exists - Flag saying whether points exist
99984f0b6dfSMatthew G. Knepley 
100084f0b6dfSMatthew G. Knepley   Level: intermediate
100184f0b6dfSMatthew G. Knepley 
100284f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
100384f0b6dfSMatthew G. Knepley @*/
1004fada774cSMatthew G. Knepley PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
1005fada774cSMatthew G. Knepley {
1006fada774cSMatthew G. Knepley   PetscInt       v;
10070c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
1008fada774cSMatthew G. Knepley 
1009fada774cSMatthew G. Knepley   PetscFunctionBegin;
1010d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1011fada774cSMatthew G. Knepley   PetscValidPointer(exists, 3);
10120c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
10130c3c4a36SLisandro Dalcin   *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE;
1014fada774cSMatthew G. Knepley   PetscFunctionReturn(0);
1015fada774cSMatthew G. Knepley }
1016fada774cSMatthew G. Knepley 
101784f0b6dfSMatthew G. Knepley /*@
101884f0b6dfSMatthew G. Knepley   DMLabelGetStratumSize - Get the size of a stratum
101984f0b6dfSMatthew G. Knepley 
102084f0b6dfSMatthew G. Knepley   Input Parameters:
102184f0b6dfSMatthew G. Knepley + label - the DMLabel
102284f0b6dfSMatthew G. Knepley - value - the stratum value
102384f0b6dfSMatthew G. Knepley 
102484f0b6dfSMatthew G. Knepley   Output Paramater:
102584f0b6dfSMatthew G. Knepley . size - The number of points 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 DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
1032c58f1c22SToby Isaac {
1033c58f1c22SToby Isaac   PetscInt       v;
1034c58f1c22SToby Isaac   PetscErrorCode ierr;
1035c58f1c22SToby Isaac 
1036c58f1c22SToby Isaac   PetscFunctionBegin;
1037d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1038c58f1c22SToby Isaac   PetscValidPointer(size, 3);
1039c58f1c22SToby Isaac   *size = 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);
1043c58f1c22SToby Isaac   *size = label->stratumSizes[v];
1044c58f1c22SToby Isaac   PetscFunctionReturn(0);
1045c58f1c22SToby Isaac }
1046c58f1c22SToby Isaac 
104784f0b6dfSMatthew G. Knepley /*@
104884f0b6dfSMatthew G. Knepley   DMLabelGetStratumBounds - Get the largest and smallest point of a stratum
104984f0b6dfSMatthew G. Knepley 
105084f0b6dfSMatthew G. Knepley   Input Parameters:
105184f0b6dfSMatthew G. Knepley + label - the DMLabel
105284f0b6dfSMatthew G. Knepley - value - the stratum value
105384f0b6dfSMatthew G. Knepley 
105484f0b6dfSMatthew G. Knepley   Output Paramaters:
105584f0b6dfSMatthew G. Knepley + start - the smallest point in the stratum
105684f0b6dfSMatthew G. Knepley - end - the largest point in the stratum
105784f0b6dfSMatthew G. Knepley 
105884f0b6dfSMatthew G. Knepley   Level: intermediate
105984f0b6dfSMatthew G. Knepley 
106084f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
106184f0b6dfSMatthew G. Knepley @*/
1062c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
1063c58f1c22SToby Isaac {
10640c3c4a36SLisandro Dalcin   PetscInt       v, min, max;
1065c58f1c22SToby Isaac   PetscErrorCode ierr;
1066c58f1c22SToby Isaac 
1067c58f1c22SToby Isaac   PetscFunctionBegin;
1068d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1069c58f1c22SToby Isaac   if (start) {PetscValidPointer(start, 3); *start = 0;}
1070c58f1c22SToby Isaac   if (end)   {PetscValidPointer(end,   4); *end   = 0;}
10710c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
10720c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
1073c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
10740c3c4a36SLisandro Dalcin   if (label->stratumSizes[v] <= 0) PetscFunctionReturn(0);
1075d6cb179aSToby Isaac   ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr);
1076d6cb179aSToby Isaac   if (start) *start = min;
1077d6cb179aSToby Isaac   if (end)   *end   = max+1;
1078c58f1c22SToby Isaac   PetscFunctionReturn(0);
1079c58f1c22SToby Isaac }
1080c58f1c22SToby Isaac 
108184f0b6dfSMatthew G. Knepley /*@
108284f0b6dfSMatthew G. Knepley   DMLabelGetStratumIS - Get an IS with the stratum points
108384f0b6dfSMatthew G. Knepley 
108484f0b6dfSMatthew G. Knepley   Input Parameters:
108584f0b6dfSMatthew G. Knepley + label - the DMLabel
108684f0b6dfSMatthew G. Knepley - value - the stratum value
108784f0b6dfSMatthew G. Knepley 
108884f0b6dfSMatthew G. Knepley   Output Paramater:
108984f0b6dfSMatthew G. Knepley . points - The stratum points
109084f0b6dfSMatthew G. Knepley 
109184f0b6dfSMatthew G. Knepley   Level: intermediate
109284f0b6dfSMatthew G. Knepley 
109384f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
109484f0b6dfSMatthew G. Knepley @*/
1095c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
1096c58f1c22SToby Isaac {
1097c58f1c22SToby Isaac   PetscInt       v;
1098c58f1c22SToby Isaac   PetscErrorCode ierr;
1099c58f1c22SToby Isaac 
1100c58f1c22SToby Isaac   PetscFunctionBegin;
1101d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1102c58f1c22SToby Isaac   PetscValidPointer(points, 3);
1103c58f1c22SToby Isaac   *points = NULL;
11040c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
11050c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
1106c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
1107ad8374ffSToby Isaac   ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr);
1108ad8374ffSToby Isaac   *points = label->points[v];
1109c58f1c22SToby Isaac   PetscFunctionReturn(0);
1110c58f1c22SToby Isaac }
1111c58f1c22SToby Isaac 
111284f0b6dfSMatthew G. Knepley /*@
111384f0b6dfSMatthew G. Knepley   DMLabelSetStratumIS - Set the stratum points using an IS
111484f0b6dfSMatthew G. Knepley 
111584f0b6dfSMatthew G. Knepley   Input Parameters:
111684f0b6dfSMatthew G. Knepley + label - the DMLabel
111784f0b6dfSMatthew G. Knepley . value - the stratum value
111884f0b6dfSMatthew G. Knepley - points - The stratum points
111984f0b6dfSMatthew G. Knepley 
112084f0b6dfSMatthew G. Knepley   Level: intermediate
112184f0b6dfSMatthew G. Knepley 
112284f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
112384f0b6dfSMatthew G. Knepley @*/
11244de306b1SToby Isaac PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
11254de306b1SToby Isaac {
11260c3c4a36SLisandro Dalcin   PetscInt       v;
11274de306b1SToby Isaac   PetscErrorCode ierr;
11284de306b1SToby Isaac 
11294de306b1SToby Isaac   PetscFunctionBegin;
1130d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1131d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
1132b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
11334de306b1SToby Isaac   if (is == label->points[v]) PetscFunctionReturn(0);
11344de306b1SToby Isaac   ierr = DMLabelClearStratum(label, value);CHKERRQ(ierr);
11354de306b1SToby Isaac   ierr = ISGetLocalSize(is, &(label->stratumSizes[v]));CHKERRQ(ierr);
11364de306b1SToby Isaac   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
11374de306b1SToby Isaac   ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
11380c3c4a36SLisandro Dalcin   label->points[v]  = is;
11390c3c4a36SLisandro Dalcin   label->validIS[v] = PETSC_TRUE;
1140*277ea44aSLisandro Dalcin   ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
11414de306b1SToby Isaac   if (label->bt) {
11424de306b1SToby Isaac     const PetscInt *points;
11434de306b1SToby Isaac     PetscInt p;
11444de306b1SToby Isaac 
11454de306b1SToby Isaac     ierr = ISGetIndices(is,&points);CHKERRQ(ierr);
11464de306b1SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
11474de306b1SToby Isaac       const PetscInt point = points[p];
11484de306b1SToby Isaac 
11494de306b1SToby 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);
11504de306b1SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
11514de306b1SToby Isaac     }
11524de306b1SToby Isaac   }
11534de306b1SToby Isaac   PetscFunctionReturn(0);
11544de306b1SToby Isaac }
11554de306b1SToby Isaac 
115684f0b6dfSMatthew G. Knepley /*@
115784f0b6dfSMatthew G. Knepley   DMLabelClearStratum - Remove a stratum
11584de306b1SToby Isaac 
115984f0b6dfSMatthew G. Knepley   Input Parameters:
116084f0b6dfSMatthew G. Knepley + label - the DMLabel
116184f0b6dfSMatthew G. Knepley - value - the stratum value
116284f0b6dfSMatthew G. Knepley 
116384f0b6dfSMatthew G. Knepley   Level: intermediate
116484f0b6dfSMatthew G. Knepley 
116584f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
116684f0b6dfSMatthew G. Knepley @*/
1167c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
1168c58f1c22SToby Isaac {
1169c58f1c22SToby Isaac   PetscInt       v;
1170c58f1c22SToby Isaac   PetscErrorCode ierr;
1171c58f1c22SToby Isaac 
1172c58f1c22SToby Isaac   PetscFunctionBegin;
1173d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
11740c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
11750c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
11764de306b1SToby Isaac   if (label->validIS[v]) {
11774de306b1SToby Isaac     if (label->bt) {
1178c58f1c22SToby Isaac       PetscInt       i;
1179ad8374ffSToby Isaac       const PetscInt *points;
1180c58f1c22SToby Isaac 
1181ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
1182c58f1c22SToby Isaac       for (i = 0; i < label->stratumSizes[v]; ++i) {
1183ad8374ffSToby Isaac         const PetscInt point = points[i];
1184c58f1c22SToby Isaac 
1185c58f1c22SToby 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);
1186c58f1c22SToby Isaac         ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
1187c58f1c22SToby Isaac       }
1188ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
1189c58f1c22SToby Isaac     }
1190c58f1c22SToby Isaac     label->stratumSizes[v] = 0;
11910c3c4a36SLisandro Dalcin     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
1192*277ea44aSLisandro Dalcin     ierr = ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &label->points[v]);CHKERRQ(ierr);
11930c3c4a36SLisandro Dalcin     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
1194*277ea44aSLisandro Dalcin     ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
1195c58f1c22SToby Isaac   } else {
1196b9907514SLisandro Dalcin     ierr = PetscHSetIClear(label->ht[v]);CHKERRQ(ierr);
1197c58f1c22SToby Isaac   }
1198c58f1c22SToby Isaac   PetscFunctionReturn(0);
1199c58f1c22SToby Isaac }
1200c58f1c22SToby Isaac 
120184f0b6dfSMatthew G. Knepley /*@
120284f0b6dfSMatthew G. Knepley   DMLabelFilter - Remove all points outside of [start, end)
120384f0b6dfSMatthew G. Knepley 
120484f0b6dfSMatthew G. Knepley   Input Parameters:
120584f0b6dfSMatthew G. Knepley + label - the DMLabel
120684f0b6dfSMatthew G. Knepley . start - the first point
120784f0b6dfSMatthew G. Knepley - end - the last point
120884f0b6dfSMatthew G. Knepley 
120984f0b6dfSMatthew G. Knepley   Level: intermediate
121084f0b6dfSMatthew G. Knepley 
121184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
121284f0b6dfSMatthew G. Knepley @*/
1213c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
1214c58f1c22SToby Isaac {
1215c58f1c22SToby Isaac   PetscInt       v;
1216c58f1c22SToby Isaac   PetscErrorCode ierr;
1217c58f1c22SToby Isaac 
1218c58f1c22SToby Isaac   PetscFunctionBegin;
1219d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
12200c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
1221c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1222c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
1223c58f1c22SToby Isaac     PetscInt off, q;
1224ad8374ffSToby Isaac     const PetscInt *points;
1225033405d5SLisandro Dalcin     PetscInt numPointsNew = 0, *pointsNew = NULL;
1226c58f1c22SToby Isaac 
1227ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
1228033405d5SLisandro Dalcin     for (q = 0; q < label->stratumSizes[v]; ++q)
1229033405d5SLisandro Dalcin       if (points[q] >= start && points[q] < end)
1230033405d5SLisandro Dalcin         numPointsNew++;
1231033405d5SLisandro Dalcin     ierr = PetscMalloc1(numPointsNew, &pointsNew);CHKERRQ(ierr);
1232c58f1c22SToby Isaac     for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) {
1233033405d5SLisandro Dalcin       if (points[q] >= start && points[q] < end)
1234033405d5SLisandro Dalcin         pointsNew[off++] = points[q];
1235ad8374ffSToby Isaac     }
1236ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
1237033405d5SLisandro Dalcin 
1238033405d5SLisandro Dalcin     label->stratumSizes[v] = numPointsNew;
1239033405d5SLisandro Dalcin     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
1240033405d5SLisandro Dalcin     ierr = ISCreateGeneral(PETSC_COMM_SELF,numPointsNew, pointsNew, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr);
1241033405d5SLisandro Dalcin     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
1242c58f1c22SToby Isaac   }
1243c58f1c22SToby Isaac   ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr);
1244c58f1c22SToby Isaac   PetscFunctionReturn(0);
1245c58f1c22SToby Isaac }
1246c58f1c22SToby Isaac 
124784f0b6dfSMatthew G. Knepley /*@
124884f0b6dfSMatthew G. Knepley   DMLabelPermute - Create a new label with permuted points
124984f0b6dfSMatthew G. Knepley 
125084f0b6dfSMatthew G. Knepley   Input Parameters:
125184f0b6dfSMatthew G. Knepley + label - the DMLabel
125284f0b6dfSMatthew G. Knepley - permutation - the point permutation
125384f0b6dfSMatthew G. Knepley 
125484f0b6dfSMatthew G. Knepley   Output Parameter:
125584f0b6dfSMatthew G. Knepley . labelnew - the new label containing the permuted points
125684f0b6dfSMatthew G. Knepley 
125784f0b6dfSMatthew G. Knepley   Level: intermediate
125884f0b6dfSMatthew G. Knepley 
125984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
126084f0b6dfSMatthew G. Knepley @*/
1261c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
1262c58f1c22SToby Isaac {
1263c58f1c22SToby Isaac   const PetscInt *perm;
1264c58f1c22SToby Isaac   PetscInt        numValues, numPoints, v, q;
1265c58f1c22SToby Isaac   PetscErrorCode  ierr;
1266c58f1c22SToby Isaac 
1267c58f1c22SToby Isaac   PetscFunctionBegin;
1268d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1269d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(permutation, IS_CLASSID, 2);
1270c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1271c58f1c22SToby Isaac   ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr);
1272c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr);
1273c58f1c22SToby Isaac   ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr);
1274c58f1c22SToby Isaac   ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr);
1275c58f1c22SToby Isaac   for (v = 0; v < numValues; ++v) {
1276c58f1c22SToby Isaac     const PetscInt size   = (*labelNew)->stratumSizes[v];
1277ad8374ffSToby Isaac     const PetscInt *points;
1278ad8374ffSToby Isaac     PetscInt *pointsNew;
1279c58f1c22SToby Isaac 
1280ad8374ffSToby Isaac     ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1281ad8374ffSToby Isaac     ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr);
1282c58f1c22SToby Isaac     for (q = 0; q < size; ++q) {
1283ad8374ffSToby Isaac       const PetscInt point = points[q];
1284c58f1c22SToby Isaac 
1285c58f1c22SToby 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);
1286ad8374ffSToby Isaac       pointsNew[q] = perm[point];
1287c58f1c22SToby Isaac     }
1288ad8374ffSToby Isaac     ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1289ad8374ffSToby Isaac     ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr);
1290ad8374ffSToby Isaac     ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr);
1291fa8e8ae5SToby Isaac     if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) {
1292fa8e8ae5SToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));CHKERRQ(ierr);
1293fa8e8ae5SToby Isaac       ierr = PetscFree(pointsNew);CHKERRQ(ierr);
1294fa8e8ae5SToby Isaac     } else {
1295ad8374ffSToby Isaac       ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr);
1296fa8e8ae5SToby Isaac     }
1297ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr);
1298c58f1c22SToby Isaac   }
1299c58f1c22SToby Isaac   ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr);
1300c58f1c22SToby Isaac   if (label->bt) {
1301c58f1c22SToby Isaac     ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
1302c58f1c22SToby Isaac     ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
1303c58f1c22SToby Isaac   }
1304c58f1c22SToby Isaac   PetscFunctionReturn(0);
1305c58f1c22SToby Isaac }
1306c58f1c22SToby Isaac 
130726c55118SMichael Lange PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
130826c55118SMichael Lange {
130926c55118SMichael Lange   MPI_Comm       comm;
131026c55118SMichael Lange   PetscInt       s, l, nroots, nleaves, dof, offset, size;
131126c55118SMichael Lange   PetscInt      *remoteOffsets, *rootStrata, *rootIdx;
131226c55118SMichael Lange   PetscSection   rootSection;
131326c55118SMichael Lange   PetscSF        labelSF;
131426c55118SMichael Lange   PetscErrorCode ierr;
131526c55118SMichael Lange 
131626c55118SMichael Lange   PetscFunctionBegin;
131726c55118SMichael Lange   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
131826c55118SMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
131926c55118SMichael Lange   /* Build a section of stratum values per point, generate the according SF
132026c55118SMichael Lange      and distribute point-wise stratum values to leaves. */
132126c55118SMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
132226c55118SMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
132326c55118SMichael Lange   ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr);
132426c55118SMichael Lange   if (label) {
132526c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1326ad8374ffSToby Isaac       const PetscInt *points;
1327ad8374ffSToby Isaac 
1328ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
132926c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1330ad8374ffSToby Isaac         ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr);
1331ad8374ffSToby Isaac         ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr);
133226c55118SMichael Lange       }
1333ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
133426c55118SMichael Lange     }
133526c55118SMichael Lange   }
133626c55118SMichael Lange   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
133726c55118SMichael Lange   /* Create a point-wise array of stratum values */
133826c55118SMichael Lange   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
133926c55118SMichael Lange   ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr);
134026c55118SMichael Lange   ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr);
134126c55118SMichael Lange   if (label) {
134226c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1343ad8374ffSToby Isaac       const PetscInt *points;
1344ad8374ffSToby Isaac 
1345ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
134626c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1347ad8374ffSToby Isaac         const PetscInt p = points[l];
134826c55118SMichael Lange         ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr);
134926c55118SMichael Lange         rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
135026c55118SMichael Lange       }
1351ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
135226c55118SMichael Lange     }
135326c55118SMichael Lange   }
135426c55118SMichael Lange   /* Build SF that maps label points to remote processes */
135526c55118SMichael Lange   ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr);
135626c55118SMichael Lange   ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr);
135726c55118SMichael Lange   ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr);
135826c55118SMichael Lange   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
135926c55118SMichael Lange   /* Send the strata for each point over the derived SF */
136026c55118SMichael Lange   ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr);
136126c55118SMichael Lange   ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr);
136226c55118SMichael Lange   ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
136326c55118SMichael Lange   ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
136426c55118SMichael Lange   /* Clean up */
136526c55118SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
136626c55118SMichael Lange   ierr = PetscFree(rootIdx);CHKERRQ(ierr);
136726c55118SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
136826c55118SMichael Lange   ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr);
136926c55118SMichael Lange   PetscFunctionReturn(0);
137026c55118SMichael Lange }
137126c55118SMichael Lange 
137284f0b6dfSMatthew G. Knepley /*@
137384f0b6dfSMatthew G. Knepley   DMLabelDistribute - Create a new label pushed forward over the PetscSF
137484f0b6dfSMatthew G. Knepley 
137584f0b6dfSMatthew G. Knepley   Input Parameters:
137684f0b6dfSMatthew G. Knepley + label - the DMLabel
137784f0b6dfSMatthew G. Knepley - sf    - the map from old to new distribution
137884f0b6dfSMatthew G. Knepley 
137984f0b6dfSMatthew G. Knepley   Output Parameter:
138084f0b6dfSMatthew G. Knepley . labelnew - the new redistributed label
138184f0b6dfSMatthew G. Knepley 
138284f0b6dfSMatthew G. Knepley   Level: intermediate
138384f0b6dfSMatthew G. Knepley 
138484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
138584f0b6dfSMatthew G. Knepley @*/
1386c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1387c58f1c22SToby Isaac {
1388c58f1c22SToby Isaac   MPI_Comm       comm;
138926c55118SMichael Lange   PetscSection   leafSection;
139026c55118SMichael Lange   PetscInt       p, pStart, pEnd, s, size, dof, offset, stratum;
139126c55118SMichael Lange   PetscInt      *leafStrata, *strataIdx;
1392ad8374ffSToby Isaac   PetscInt     **points;
1393d67d17b1SMatthew G. Knepley   const char    *lname = NULL;
1394c58f1c22SToby Isaac   char          *name;
1395c58f1c22SToby Isaac   PetscInt       nameSize;
1396e8f14785SLisandro Dalcin   PetscHSetI     stratumHash;
1397c58f1c22SToby Isaac   size_t         len = 0;
139826c55118SMichael Lange   PetscMPIInt    rank;
1399c58f1c22SToby Isaac   PetscErrorCode ierr;
1400c58f1c22SToby Isaac 
1401c58f1c22SToby Isaac   PetscFunctionBegin;
1402d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1403f018e600SMatthew G. Knepley   if (label) {
1404f018e600SMatthew G. Knepley     PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1405f018e600SMatthew G. Knepley     ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1406f018e600SMatthew G. Knepley   }
1407c58f1c22SToby Isaac   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1408c58f1c22SToby Isaac   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1409c58f1c22SToby Isaac   /* Bcast name */
1410d67d17b1SMatthew G. Knepley   if (!rank) {
1411d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1412d67d17b1SMatthew G. Knepley     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1413d67d17b1SMatthew G. Knepley   }
1414c58f1c22SToby Isaac   nameSize = len;
1415c58f1c22SToby Isaac   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1416c58f1c22SToby Isaac   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1417d67d17b1SMatthew G. Knepley   if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);}
1418c58f1c22SToby Isaac   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1419d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
1420c58f1c22SToby Isaac   ierr = PetscFree(name);CHKERRQ(ierr);
142177d236dfSMichael Lange   /* Bcast defaultValue */
142277d236dfSMichael Lange   if (!rank) (*labelNew)->defaultValue = label->defaultValue;
142377d236dfSMichael Lange   ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
142426c55118SMichael Lange   /* Distribute stratum values over the SF and get the point mapping on the receiver */
142526c55118SMichael Lange   ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr);
14265cbdf6fcSMichael Lange   /* Determine received stratum values and initialise new label*/
1427e8f14785SLisandro Dalcin   ierr = PetscHSetICreate(&stratumHash);CHKERRQ(ierr);
142826c55118SMichael Lange   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
1429e8f14785SLisandro Dalcin   for (p = 0; p < size; ++p) {ierr = PetscHSetIAdd(stratumHash, leafStrata[p]);CHKERRQ(ierr);}
1430e8f14785SLisandro Dalcin   ierr = PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata);CHKERRQ(ierr);
1431ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr);
1432ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
14335cbdf6fcSMichael Lange   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr);
14345cbdf6fcSMichael Lange   /* Turn leafStrata into indices rather than stratum values */
14355cbdf6fcSMichael Lange   offset = 0;
1436e8f14785SLisandro Dalcin   ierr = PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr);
1437a205f1fdSToby Isaac   ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr);
143890e9b2aeSLisandro Dalcin   for (s = 0; s < (*labelNew)->numStrata; ++s) {
143990e9b2aeSLisandro Dalcin     ierr = PetscHMapISet((*labelNew)->hmap, (*labelNew)->stratumValues[s], s);CHKERRQ(ierr);
144090e9b2aeSLisandro Dalcin   }
14415cbdf6fcSMichael Lange   for (p = 0; p < size; ++p) {
1442231b9e6fSMatthew G. Knepley     for (s = 0; s < (*labelNew)->numStrata; ++s) {
1443231b9e6fSMatthew G. Knepley       if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
14445cbdf6fcSMichael Lange     }
14455cbdf6fcSMichael Lange   }
1446c58f1c22SToby Isaac   /* Rebuild the point strata on the receiver */
1447c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr);
1448c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1449c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1450c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1451c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1452c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1453c58f1c22SToby Isaac       (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1454c58f1c22SToby Isaac     }
1455c58f1c22SToby Isaac   }
1456c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr);
1457c58f1c22SToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr);
1458ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr);
1459c58f1c22SToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1460e8f14785SLisandro Dalcin     ierr = PetscHSetICreate(&(*labelNew)->ht[s]);CHKERRQ(ierr);
1461ad8374ffSToby Isaac     ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr);
1462c58f1c22SToby Isaac   }
1463c58f1c22SToby Isaac   /* Insert points into new strata */
1464c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr);
1465c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1466c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1467c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1468c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1469c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1470c58f1c22SToby Isaac       stratum = leafStrata[offset+s];
1471ad8374ffSToby Isaac       points[stratum][strataIdx[stratum]++] = p;
1472c58f1c22SToby Isaac     }
1473c58f1c22SToby Isaac   }
1474ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; s++) {
1475ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr);
1476ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr);
1477ad8374ffSToby Isaac   }
1478ad8374ffSToby Isaac   ierr = PetscFree(points);CHKERRQ(ierr);
1479e8f14785SLisandro Dalcin   ierr = PetscHSetIDestroy(&stratumHash);CHKERRQ(ierr);
1480c58f1c22SToby Isaac   ierr = PetscFree(leafStrata);CHKERRQ(ierr);
1481c58f1c22SToby Isaac   ierr = PetscFree(strataIdx);CHKERRQ(ierr);
1482c58f1c22SToby Isaac   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
1483c58f1c22SToby Isaac   PetscFunctionReturn(0);
1484c58f1c22SToby Isaac }
1485c58f1c22SToby Isaac 
14867937d9ceSMichael Lange /*@
14877937d9ceSMichael Lange   DMLabelGather - Gather all label values from leafs into roots
14887937d9ceSMichael Lange 
14897937d9ceSMichael Lange   Input Parameters:
14907937d9ceSMichael Lange + label - the DMLabel
149184f0b6dfSMatthew G. Knepley - sf - the Star Forest point communication map
14927937d9ceSMichael Lange 
149384f0b6dfSMatthew G. Knepley   Output Parameters:
149484f0b6dfSMatthew G. Knepley . labelNew - the new DMLabel with localised leaf values
14957937d9ceSMichael Lange 
14967937d9ceSMichael Lange   Level: developer
14977937d9ceSMichael Lange 
14987937d9ceSMichael Lange   Note: This is the inverse operation to DMLabelDistribute.
14997937d9ceSMichael Lange 
15007937d9ceSMichael Lange .seealso: DMLabelDistribute()
15017937d9ceSMichael Lange @*/
15027937d9ceSMichael Lange PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
15037937d9ceSMichael Lange {
15047937d9ceSMichael Lange   MPI_Comm       comm;
15057937d9ceSMichael Lange   PetscSection   rootSection;
15067937d9ceSMichael Lange   PetscSF        sfLabel;
15077937d9ceSMichael Lange   PetscSFNode   *rootPoints, *leafPoints;
15087937d9ceSMichael Lange   PetscInt       p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
15097937d9ceSMichael Lange   const PetscInt *rootDegree, *ilocal;
15107937d9ceSMichael Lange   PetscInt       *rootStrata;
1511d67d17b1SMatthew G. Knepley   const char    *lname;
15127937d9ceSMichael Lange   char          *name;
15137937d9ceSMichael Lange   PetscInt       nameSize;
15147937d9ceSMichael Lange   size_t         len = 0;
15159852e123SBarry Smith   PetscMPIInt    rank, size;
15167937d9ceSMichael Lange   PetscErrorCode ierr;
15177937d9ceSMichael Lange 
15187937d9ceSMichael Lange   PetscFunctionBegin;
1519d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1520d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
15217937d9ceSMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
15227937d9ceSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
15239852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
15247937d9ceSMichael Lange   /* Bcast name */
1525d67d17b1SMatthew G. Knepley   if (!rank) {
1526d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1527d67d17b1SMatthew G. Knepley     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1528d67d17b1SMatthew G. Knepley   }
15297937d9ceSMichael Lange   nameSize = len;
15307937d9ceSMichael Lange   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
15317937d9ceSMichael Lange   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1532d67d17b1SMatthew G. Knepley   if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);}
15337937d9ceSMichael Lange   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1534d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
15357937d9ceSMichael Lange   ierr = PetscFree(name);CHKERRQ(ierr);
15367937d9ceSMichael Lange   /* Gather rank/index pairs of leaves into local roots to build
15377937d9ceSMichael Lange      an inverse, multi-rooted SF. Note that this ignores local leaf
15387937d9ceSMichael Lange      indexing due to the use of the multiSF in PetscSFGather. */
15397937d9ceSMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr);
1540dc53bc9bSMatthew G. Knepley   ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr);
1541dc53bc9bSMatthew G. Knepley   for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
15427937d9ceSMichael Lange   for (p = 0; p < nleaves; p++) {
15438212dd46SStefano Zampini     PetscInt ilp = ilocal ? ilocal[p] : p;
15448212dd46SStefano Zampini 
15458212dd46SStefano Zampini     leafPoints[ilp].index = ilp;
15468212dd46SStefano Zampini     leafPoints[ilp].rank  = rank;
15477937d9ceSMichael Lange   }
15487937d9ceSMichael Lange   ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr);
15497937d9ceSMichael Lange   ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr);
15507937d9ceSMichael Lange   for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
15517937d9ceSMichael Lange   ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr);
15527937d9ceSMichael Lange   ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
15537937d9ceSMichael Lange   ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
15547937d9ceSMichael Lange   ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr);
15557937d9ceSMichael Lange   ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
15567937d9ceSMichael Lange   /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
15577937d9ceSMichael Lange   ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr);
15587937d9ceSMichael Lange   /* Rebuild the point strata on the receiver */
15597937d9ceSMichael Lange   for (p = 0, idx = 0; p < nroots; p++) {
15607937d9ceSMichael Lange     for (d = 0; d < rootDegree[p]; d++) {
15617937d9ceSMichael Lange       ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr);
15627937d9ceSMichael Lange       ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr);
15637937d9ceSMichael Lange       for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);}
15647937d9ceSMichael Lange     }
15657937d9ceSMichael Lange     idx += rootDegree[p];
15667937d9ceSMichael Lange   }
156777e0c0e7SMichael Lange   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
156877e0c0e7SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
156977e0c0e7SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
157077e0c0e7SMichael Lange   ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr);
15717937d9ceSMichael Lange   PetscFunctionReturn(0);
15727937d9ceSMichael Lange }
15737937d9ceSMichael Lange 
157484f0b6dfSMatthew G. Knepley /*@
157584f0b6dfSMatthew G. Knepley   DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label
157684f0b6dfSMatthew G. Knepley 
157784f0b6dfSMatthew G. Knepley   Input Parameter:
157884f0b6dfSMatthew G. Knepley . label - the DMLabel
157984f0b6dfSMatthew G. Knepley 
158084f0b6dfSMatthew G. Knepley   Output Parameters:
158184f0b6dfSMatthew G. Knepley + section - the section giving offsets for each stratum
158284f0b6dfSMatthew G. Knepley - is - An IS containing all the label points
158384f0b6dfSMatthew G. Knepley 
158484f0b6dfSMatthew G. Knepley   Level: developer
158584f0b6dfSMatthew G. Knepley 
158684f0b6dfSMatthew G. Knepley .seealso: DMLabelDistribute()
158784f0b6dfSMatthew G. Knepley @*/
1588c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1589c58f1c22SToby Isaac {
1590c58f1c22SToby Isaac   IS              vIS;
1591c58f1c22SToby Isaac   const PetscInt *values;
1592c58f1c22SToby Isaac   PetscInt       *points;
1593c58f1c22SToby Isaac   PetscInt        nV, vS = 0, vE = 0, v, N;
1594c58f1c22SToby Isaac   PetscErrorCode  ierr;
1595c58f1c22SToby Isaac 
1596c58f1c22SToby Isaac   PetscFunctionBegin;
1597d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1598c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr);
1599c58f1c22SToby Isaac   ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
1600c58f1c22SToby Isaac   ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr);
1601c58f1c22SToby Isaac   if (nV) {vS = values[0]; vE = values[0]+1;}
1602c58f1c22SToby Isaac   for (v = 1; v < nV; ++v) {
1603c58f1c22SToby Isaac     vS = PetscMin(vS, values[v]);
1604c58f1c22SToby Isaac     vE = PetscMax(vE, values[v]+1);
1605c58f1c22SToby Isaac   }
1606c58f1c22SToby Isaac   ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr);
1607c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr);
1608c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1609c58f1c22SToby Isaac     PetscInt n;
1610c58f1c22SToby Isaac 
1611c58f1c22SToby Isaac     ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr);
1612c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr);
1613c58f1c22SToby Isaac   }
1614c58f1c22SToby Isaac   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
1615c58f1c22SToby Isaac   ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr);
1616c58f1c22SToby Isaac   ierr = PetscMalloc1(N, &points);CHKERRQ(ierr);
1617c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1618c58f1c22SToby Isaac     IS              is;
1619c58f1c22SToby Isaac     const PetscInt *spoints;
1620c58f1c22SToby Isaac     PetscInt        dof, off, p;
1621c58f1c22SToby Isaac 
1622c58f1c22SToby Isaac     ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr);
1623c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr);
1624c58f1c22SToby Isaac     ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr);
1625c58f1c22SToby Isaac     ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr);
1626c58f1c22SToby Isaac     for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1627c58f1c22SToby Isaac     ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr);
1628c58f1c22SToby Isaac     ierr = ISDestroy(&is);CHKERRQ(ierr);
1629c58f1c22SToby Isaac   }
1630c58f1c22SToby Isaac   ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr);
1631c58f1c22SToby Isaac   ierr = ISDestroy(&vIS);CHKERRQ(ierr);
1632c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1633c58f1c22SToby Isaac   PetscFunctionReturn(0);
1634c58f1c22SToby Isaac }
1635c58f1c22SToby Isaac 
163684f0b6dfSMatthew G. Knepley /*@
1637c58f1c22SToby Isaac   PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1638c58f1c22SToby Isaac   the local section and an SF describing the section point overlap.
1639c58f1c22SToby Isaac 
1640c58f1c22SToby Isaac   Input Parameters:
1641c58f1c22SToby Isaac   + s - The PetscSection for the local field layout
1642c58f1c22SToby Isaac   . sf - The SF describing parallel layout of the section points
1643c58f1c22SToby Isaac   . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1644c58f1c22SToby Isaac   . label - The label specifying the points
1645c58f1c22SToby Isaac   - labelValue - The label stratum specifying the points
1646c58f1c22SToby Isaac 
1647c58f1c22SToby Isaac   Output Parameter:
1648c58f1c22SToby Isaac   . gsection - The PetscSection for the global field layout
1649c58f1c22SToby Isaac 
1650c58f1c22SToby Isaac   Note: This gives negative sizes and offsets to points not owned by this process
1651c58f1c22SToby Isaac 
1652c58f1c22SToby Isaac   Level: developer
1653c58f1c22SToby Isaac 
1654c58f1c22SToby Isaac .seealso: PetscSectionCreate()
1655c58f1c22SToby Isaac @*/
1656c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1657c58f1c22SToby Isaac {
1658c58f1c22SToby Isaac   PetscInt      *neg = NULL, *tmpOff = NULL;
1659c58f1c22SToby Isaac   PetscInt       pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1660c58f1c22SToby Isaac   PetscErrorCode ierr;
1661c58f1c22SToby Isaac 
1662c58f1c22SToby Isaac   PetscFunctionBegin;
1663d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1664d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1665d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4);
1666c58f1c22SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr);
1667c58f1c22SToby Isaac   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1668c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr);
1669c58f1c22SToby Isaac   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1670c58f1c22SToby Isaac   if (nroots >= 0) {
1671c58f1c22SToby Isaac     if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1672c58f1c22SToby Isaac     ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr);
1673c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1674c58f1c22SToby Isaac       ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr);
1675c58f1c22SToby Isaac     } else {
1676c58f1c22SToby Isaac       tmpOff = &(*gsection)->atlasDof[-pStart];
1677c58f1c22SToby Isaac     }
1678c58f1c22SToby Isaac   }
1679c58f1c22SToby Isaac   /* Mark ghost points with negative dof */
1680c58f1c22SToby Isaac   for (p = pStart; p < pEnd; ++p) {
1681c58f1c22SToby Isaac     PetscInt value;
1682c58f1c22SToby Isaac 
1683c58f1c22SToby Isaac     ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr);
1684c58f1c22SToby Isaac     if (value != labelValue) continue;
1685c58f1c22SToby Isaac     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1686c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr);
1687c58f1c22SToby Isaac     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1688c58f1c22SToby Isaac     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);}
1689c58f1c22SToby Isaac     if (neg) neg[p] = -(dof+1);
1690c58f1c22SToby Isaac   }
1691c58f1c22SToby Isaac   ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr);
1692c58f1c22SToby Isaac   if (nroots >= 0) {
1693c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1694c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1695c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1696c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1697c58f1c22SToby Isaac     }
1698c58f1c22SToby Isaac   }
1699c58f1c22SToby Isaac   /* Calculate new sizes, get proccess offset, and calculate point offsets */
1700c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1701c58f1c22SToby Isaac     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1702c58f1c22SToby Isaac     (*gsection)->atlasOff[p] = off;
1703c58f1c22SToby Isaac     off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1704c58f1c22SToby Isaac   }
1705c58f1c22SToby Isaac   ierr       = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1706c58f1c22SToby Isaac   globalOff -= off;
1707c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1708c58f1c22SToby Isaac     (*gsection)->atlasOff[p] += globalOff;
1709c58f1c22SToby Isaac     if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1710c58f1c22SToby Isaac   }
1711c58f1c22SToby Isaac   /* Put in negative offsets for ghost points */
1712c58f1c22SToby Isaac   if (nroots >= 0) {
1713c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1714c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1715c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1716c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1717c58f1c22SToby Isaac     }
1718c58f1c22SToby Isaac   }
1719c58f1c22SToby Isaac   if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);}
1720c58f1c22SToby Isaac   ierr = PetscFree(neg);CHKERRQ(ierr);
1721c58f1c22SToby Isaac   PetscFunctionReturn(0);
1722c58f1c22SToby Isaac }
1723c58f1c22SToby Isaac 
17245fdea053SToby Isaac typedef struct _n_PetscSectionSym_Label
17255fdea053SToby Isaac {
17265fdea053SToby Isaac   DMLabel           label;
17275fdea053SToby Isaac   PetscCopyMode     *modes;
17285fdea053SToby Isaac   PetscInt          *sizes;
17295fdea053SToby Isaac   const PetscInt    ***perms;
17305fdea053SToby Isaac   const PetscScalar ***rots;
17315fdea053SToby Isaac   PetscInt          (*minMaxOrients)[2];
17325fdea053SToby Isaac   PetscInt          numStrata; /* numStrata is only increasing, functions as a state */
17335fdea053SToby Isaac } PetscSectionSym_Label;
17345fdea053SToby Isaac 
17355fdea053SToby Isaac static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
17365fdea053SToby Isaac {
17375fdea053SToby Isaac   PetscInt              i, j;
17385fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
17395fdea053SToby Isaac   PetscErrorCode        ierr;
17405fdea053SToby Isaac 
17415fdea053SToby Isaac   PetscFunctionBegin;
17425fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
17435fdea053SToby Isaac     if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
17445fdea053SToby Isaac       for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
17455fdea053SToby Isaac         if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);}
17465fdea053SToby Isaac         if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);}
17475fdea053SToby Isaac       }
17485fdea053SToby Isaac       if (sl->perms[i]) {
17495fdea053SToby Isaac         const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
17505fdea053SToby Isaac 
17515fdea053SToby Isaac         ierr = PetscFree(perms);CHKERRQ(ierr);
17525fdea053SToby Isaac       }
17535fdea053SToby Isaac       if (sl->rots[i]) {
17545fdea053SToby Isaac         const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
17555fdea053SToby Isaac 
17565fdea053SToby Isaac         ierr = PetscFree(rots);CHKERRQ(ierr);
17575fdea053SToby Isaac       }
17585fdea053SToby Isaac     }
17595fdea053SToby Isaac   }
17605fdea053SToby Isaac   ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr);
17615fdea053SToby Isaac   ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr);
17625fdea053SToby Isaac   sl->numStrata = 0;
17635fdea053SToby Isaac   PetscFunctionReturn(0);
17645fdea053SToby Isaac }
17655fdea053SToby Isaac 
17665fdea053SToby Isaac static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
17675fdea053SToby Isaac {
17685fdea053SToby Isaac   PetscErrorCode ierr;
17695fdea053SToby Isaac 
17705fdea053SToby Isaac   PetscFunctionBegin;
17715fdea053SToby Isaac   ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);
17725fdea053SToby Isaac   ierr = PetscFree(sym->data);CHKERRQ(ierr);
17735fdea053SToby Isaac   PetscFunctionReturn(0);
17745fdea053SToby Isaac }
17755fdea053SToby Isaac 
17765fdea053SToby Isaac static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
17775fdea053SToby Isaac {
17785fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
17795fdea053SToby Isaac   PetscBool             isAscii;
17805fdea053SToby Isaac   DMLabel               label = sl->label;
1781d67d17b1SMatthew G. Knepley   const char           *name;
17825fdea053SToby Isaac   PetscErrorCode        ierr;
17835fdea053SToby Isaac 
17845fdea053SToby Isaac   PetscFunctionBegin;
17855fdea053SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
17865fdea053SToby Isaac   if (isAscii) {
17875fdea053SToby Isaac     PetscInt          i, j, k;
17885fdea053SToby Isaac     PetscViewerFormat format;
17895fdea053SToby Isaac 
17905fdea053SToby Isaac     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
17915fdea053SToby Isaac     if (label) {
17925fdea053SToby Isaac       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
17935fdea053SToby Isaac       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
17945fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
17955fdea053SToby Isaac         ierr = DMLabelView(label, viewer);CHKERRQ(ierr);
17965fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
17975fdea053SToby Isaac       } else {
1798d67d17b1SMatthew G. Knepley         ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1799d67d17b1SMatthew G. Knepley         ierr = PetscViewerASCIIPrintf(viewer,"  Label '%s'\n",name);CHKERRQ(ierr);
18005fdea053SToby Isaac       }
18015fdea053SToby Isaac     } else {
18025fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr);
18035fdea053SToby Isaac     }
18045fdea053SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18055fdea053SToby Isaac     for (i = 0; i <= sl->numStrata; i++) {
18065fdea053SToby Isaac       PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
18075fdea053SToby Isaac 
18085fdea053SToby Isaac       if (!(sl->perms[i] || sl->rots[i])) {
18095fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr);
18105fdea053SToby Isaac       } else {
18115fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr);
18125fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18135fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr);
18145fdea053SToby Isaac         if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
18155fdea053SToby Isaac           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18165fdea053SToby Isaac           for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
18175fdea053SToby Isaac             if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
18185fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr);
18195fdea053SToby Isaac             } else {
18205fdea053SToby Isaac               PetscInt tab;
18215fdea053SToby Isaac 
18225fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr);
18235fdea053SToby Isaac               ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18245fdea053SToby Isaac               ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr);
18255fdea053SToby Isaac               if (sl->perms[i] && sl->perms[i][j]) {
18265fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr);
18275fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
18285fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);}
18295fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
18305fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
18315fdea053SToby Isaac               }
18325fdea053SToby Isaac               if (sl->rots[i] && sl->rots[i][j]) {
18335fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Rotations:  ");CHKERRQ(ierr);
18345fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
18355fdea053SToby Isaac #if defined(PETSC_USE_COMPLEX)
18365fdea053SToby 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);}
18375fdea053SToby Isaac #else
18385fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);}
18395fdea053SToby Isaac #endif
18405fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
18415fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
18425fdea053SToby Isaac               }
18435fdea053SToby Isaac               ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
18445fdea053SToby Isaac             }
18455fdea053SToby Isaac           }
18465fdea053SToby Isaac           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
18475fdea053SToby Isaac         }
18485fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
18495fdea053SToby Isaac       }
18505fdea053SToby Isaac     }
18515fdea053SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
18525fdea053SToby Isaac   }
18535fdea053SToby Isaac   PetscFunctionReturn(0);
18545fdea053SToby Isaac }
18555fdea053SToby Isaac 
18565fdea053SToby Isaac /*@
18575fdea053SToby Isaac   PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
18585fdea053SToby Isaac 
18595fdea053SToby Isaac   Logically collective on sym
18605fdea053SToby Isaac 
18615fdea053SToby Isaac   Input parameters:
18625fdea053SToby Isaac + sym - the section symmetries
18635fdea053SToby Isaac - label - the DMLabel describing the types of points
18645fdea053SToby Isaac 
18655fdea053SToby Isaac   Level: developer:
18665fdea053SToby Isaac 
18675fdea053SToby Isaac .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()
18685fdea053SToby Isaac @*/
18695fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
18705fdea053SToby Isaac {
18715fdea053SToby Isaac   PetscSectionSym_Label *sl;
18725fdea053SToby Isaac   PetscErrorCode        ierr;
18735fdea053SToby Isaac 
18745fdea053SToby Isaac   PetscFunctionBegin;
18755fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
18765fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
18775fdea053SToby Isaac   if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);}
18785fdea053SToby Isaac   if (label) {
18795fdea053SToby Isaac     sl->label = label;
1880d67d17b1SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) label);CHKERRQ(ierr);
18815fdea053SToby Isaac     ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr);
18821a834cf9SToby 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);
18831a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr);
18841a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr);
18851a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr);
18861a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr);
18871a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr);
18885fdea053SToby Isaac   }
18895fdea053SToby Isaac   PetscFunctionReturn(0);
18905fdea053SToby Isaac }
18915fdea053SToby Isaac 
18925fdea053SToby Isaac /*@C
18935fdea053SToby Isaac   PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
18945fdea053SToby Isaac 
18955fdea053SToby Isaac   Logically collective on PetscSectionSym
18965fdea053SToby Isaac 
18975fdea053SToby Isaac   InputParameters:
18985fdea053SToby Isaac + sys       - the section symmetries
18995fdea053SToby Isaac . stratum   - the stratum value in the label that we are assigning symmetries for
19005fdea053SToby Isaac . size      - the number of dofs for points in the stratum of the label
19015fdea053SToby Isaac . minOrient - the smallest orientation for a point in this stratum
19025fdea053SToby Isaac . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
19035fdea053SToby Isaac . mode      - how sym should copy the perms and rots arrays
19045fdea053SToby Isaac . perms     - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation.  A NULL permutation is the identity
19055fdea053SToby 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
19065fdea053SToby Isaac 
19075fdea053SToby Isaac   Level: developer
19085fdea053SToby Isaac 
19095fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()
19105fdea053SToby Isaac @*/
19115fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
19125fdea053SToby Isaac {
19135fdea053SToby Isaac   PetscSectionSym_Label *sl;
1914d67d17b1SMatthew G. Knepley   const char            *name;
1915d67d17b1SMatthew G. Knepley   PetscInt               i, j, k;
19165fdea053SToby Isaac   PetscErrorCode         ierr;
19175fdea053SToby Isaac 
19185fdea053SToby Isaac   PetscFunctionBegin;
19195fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
19205fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
19215fdea053SToby Isaac   if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet");
19225fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
19235fdea053SToby Isaac     PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
19245fdea053SToby Isaac 
19255fdea053SToby Isaac     if (stratum == value) break;
19265fdea053SToby Isaac   }
1927d67d17b1SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1928d67d17b1SMatthew G. Knepley   if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,name);
19295fdea053SToby Isaac   sl->sizes[i] = size;
19305fdea053SToby Isaac   sl->modes[i] = mode;
19315fdea053SToby Isaac   sl->minMaxOrients[i][0] = minOrient;
19325fdea053SToby Isaac   sl->minMaxOrients[i][1] = maxOrient;
19335fdea053SToby Isaac   if (mode == PETSC_COPY_VALUES) {
19345fdea053SToby Isaac     if (perms) {
19355fdea053SToby Isaac       PetscInt    **ownPerms;
19365fdea053SToby Isaac 
19375fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr);
19385fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
19395fdea053SToby Isaac         if (perms[j]) {
19405fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr);
19415fdea053SToby Isaac           for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];}
19425fdea053SToby Isaac         }
19435fdea053SToby Isaac       }
19445fdea053SToby Isaac       sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient];
19455fdea053SToby Isaac     }
19465fdea053SToby Isaac     if (rots) {
19475fdea053SToby Isaac       PetscScalar **ownRots;
19485fdea053SToby Isaac 
19495fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr);
19505fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
19515fdea053SToby Isaac         if (rots[j]) {
19525fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr);
19535fdea053SToby Isaac           for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];}
19545fdea053SToby Isaac         }
19555fdea053SToby Isaac       }
19565fdea053SToby Isaac       sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient];
19575fdea053SToby Isaac     }
19585fdea053SToby Isaac   } else {
19595fdea053SToby Isaac     sl->perms[i] = perms ? &perms[-minOrient] : NULL;
19605fdea053SToby Isaac     sl->rots[i]  = rots ? &rots[-minOrient] : NULL;
19615fdea053SToby Isaac   }
19625fdea053SToby Isaac   PetscFunctionReturn(0);
19635fdea053SToby Isaac }
19645fdea053SToby Isaac 
19655fdea053SToby Isaac static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
19665fdea053SToby Isaac {
19675fdea053SToby Isaac   PetscInt              i, j, numStrata;
19685fdea053SToby Isaac   PetscSectionSym_Label *sl;
19695fdea053SToby Isaac   DMLabel               label;
19705fdea053SToby Isaac   PetscErrorCode        ierr;
19715fdea053SToby Isaac 
19725fdea053SToby Isaac   PetscFunctionBegin;
19735fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
19745fdea053SToby Isaac   numStrata = sl->numStrata;
19755fdea053SToby Isaac   label     = sl->label;
19765fdea053SToby Isaac   for (i = 0; i < numPoints; i++) {
19775fdea053SToby Isaac     PetscInt point = points[2*i];
19785fdea053SToby Isaac     PetscInt ornt  = points[2*i+1];
19795fdea053SToby Isaac 
19805fdea053SToby Isaac     for (j = 0; j < numStrata; j++) {
19815fdea053SToby Isaac       if (label->validIS[j]) {
19825fdea053SToby Isaac         PetscInt k;
19835fdea053SToby Isaac 
19845fdea053SToby Isaac         ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr);
19855fdea053SToby Isaac         if (k >= 0) break;
19865fdea053SToby Isaac       } else {
19875fdea053SToby Isaac         PetscBool has;
19885fdea053SToby Isaac 
1989b9907514SLisandro Dalcin         ierr = PetscHSetIHas(label->ht[j], point, &has);CHKERRQ(ierr);
19905fdea053SToby Isaac         if (has) break;
19915fdea053SToby Isaac       }
19925fdea053SToby Isaac     }
19935fdea053SToby 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);
19945fdea053SToby Isaac     if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;}
19955fdea053SToby Isaac     if (rots) {rots[i]  = sl->rots[j] ? sl->rots[j][ornt] : NULL;}
19965fdea053SToby Isaac   }
19975fdea053SToby Isaac   PetscFunctionReturn(0);
19985fdea053SToby Isaac }
19995fdea053SToby Isaac 
20005fdea053SToby Isaac PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
20015fdea053SToby Isaac {
20025fdea053SToby Isaac   PetscSectionSym_Label *sl;
20035fdea053SToby Isaac   PetscErrorCode        ierr;
20045fdea053SToby Isaac 
20055fdea053SToby Isaac   PetscFunctionBegin;
20065fdea053SToby Isaac   ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr);
20075fdea053SToby Isaac   sym->ops->getpoints = PetscSectionSymGetPoints_Label;
20085fdea053SToby Isaac   sym->ops->view      = PetscSectionSymView_Label;
20095fdea053SToby Isaac   sym->ops->destroy   = PetscSectionSymDestroy_Label;
20105fdea053SToby Isaac   sym->data           = (void *) sl;
20115fdea053SToby Isaac   PetscFunctionReturn(0);
20125fdea053SToby Isaac }
20135fdea053SToby Isaac 
20145fdea053SToby Isaac /*@
20155fdea053SToby Isaac   PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
20165fdea053SToby Isaac 
20175fdea053SToby Isaac   Collective
20185fdea053SToby Isaac 
20195fdea053SToby Isaac   Input Parameters:
20205fdea053SToby Isaac + comm - the MPI communicator for the new symmetry
20215fdea053SToby Isaac - label - the label defining the strata
20225fdea053SToby Isaac 
20235fdea053SToby Isaac   Output Parameters:
20245fdea053SToby Isaac . sym - the section symmetries
20255fdea053SToby Isaac 
20265fdea053SToby Isaac   Level: developer
20275fdea053SToby Isaac 
20285fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()
20295fdea053SToby Isaac @*/
20305fdea053SToby Isaac PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
20315fdea053SToby Isaac {
20325fdea053SToby Isaac   PetscErrorCode        ierr;
20335fdea053SToby Isaac 
20345fdea053SToby Isaac   PetscFunctionBegin;
20355fdea053SToby Isaac   ierr = DMInitializePackage();CHKERRQ(ierr);
20365fdea053SToby Isaac   ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr);
20375fdea053SToby Isaac   ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr);
20385fdea053SToby Isaac   ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr);
20395fdea053SToby Isaac   PetscFunctionReturn(0);
20405fdea053SToby Isaac }
2041