xref: /petsc/src/dm/label/dmlabel.c (revision ea844a1adba28e307898a6e94a1f422fd7d2db48)
15fdea053SToby Isaac #include <petscdm.h>
2c58f1c22SToby Isaac #include <petsc/private/dmlabelimpl.h>   /*I      "petscdmlabel.h"   I*/
3*ea844a1aSMatthew Knepley #include <petsc/private/sectionimpl.h>   /*I      "petscsection.h"   I*/
4c58f1c22SToby Isaac #include <petscsf.h>
5*ea844a1aSMatthew Knepley #include <petscsection.h>
6c58f1c22SToby Isaac 
7c58f1c22SToby Isaac /*@C
8c58f1c22SToby Isaac   DMLabelCreate - Create a DMLabel object, which is a multimap
9c58f1c22SToby Isaac 
105b5e7992SMatthew G. Knepley   Collective
115b5e7992SMatthew G. Knepley 
12d67d17b1SMatthew G. Knepley   Input parameters:
13d67d17b1SMatthew G. Knepley + comm - The communicator, usually PETSC_COMM_SELF
14d67d17b1SMatthew G. Knepley - name - The label name
15c58f1c22SToby Isaac 
16c58f1c22SToby Isaac   Output parameter:
17c58f1c22SToby Isaac . label - The DMLabel
18c58f1c22SToby Isaac 
19c58f1c22SToby Isaac   Level: beginner
20c58f1c22SToby Isaac 
21c58f1c22SToby Isaac .seealso: DMLabelDestroy()
22c58f1c22SToby Isaac @*/
23d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelCreate(MPI_Comm comm, const char name[], DMLabel *label)
24c58f1c22SToby Isaac {
25c58f1c22SToby Isaac   PetscErrorCode ierr;
26c58f1c22SToby Isaac 
27c58f1c22SToby Isaac   PetscFunctionBegin;
28d67d17b1SMatthew G. Knepley   PetscValidPointer(label,2);
29d67d17b1SMatthew G. Knepley   ierr = DMInitializePackage();CHKERRQ(ierr);
30c58f1c22SToby Isaac 
31d67d17b1SMatthew G. Knepley   ierr = PetscHeaderCreate(*label,DMLABEL_CLASSID,"DMLabel","DMLabel","DM",comm,DMLabelDestroy,DMLabelView);CHKERRQ(ierr);
32d67d17b1SMatthew G. Knepley 
33c58f1c22SToby Isaac   (*label)->numStrata      = 0;
345aa44df4SToby Isaac   (*label)->defaultValue   = -1;
35c58f1c22SToby Isaac   (*label)->stratumValues  = NULL;
36ad8374ffSToby Isaac   (*label)->validIS        = NULL;
37c58f1c22SToby Isaac   (*label)->stratumSizes   = NULL;
38c58f1c22SToby Isaac   (*label)->points         = NULL;
39c58f1c22SToby Isaac   (*label)->ht             = NULL;
40c58f1c22SToby Isaac   (*label)->pStart         = -1;
41c58f1c22SToby Isaac   (*label)->pEnd           = -1;
42c58f1c22SToby Isaac   (*label)->bt             = NULL;
43b9907514SLisandro Dalcin   ierr = PetscHMapICreate(&(*label)->hmap);CHKERRQ(ierr);
44d67d17b1SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *label, name);CHKERRQ(ierr);
45c58f1c22SToby Isaac   PetscFunctionReturn(0);
46c58f1c22SToby Isaac }
47c58f1c22SToby Isaac 
48c58f1c22SToby Isaac /*
49c58f1c22SToby Isaac   DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format
50c58f1c22SToby Isaac 
515b5e7992SMatthew G. Knepley   Not collective
525b5e7992SMatthew G. Knepley 
53c58f1c22SToby Isaac   Input parameter:
54c58f1c22SToby Isaac + label - The DMLabel
55c58f1c22SToby Isaac - v - The stratum value
56c58f1c22SToby Isaac 
57c58f1c22SToby Isaac   Output parameter:
58c58f1c22SToby Isaac . label - The DMLabel with stratum in sorted list format
59c58f1c22SToby Isaac 
60c58f1c22SToby Isaac   Level: developer
61c58f1c22SToby Isaac 
62c58f1c22SToby Isaac .seealso: DMLabelCreate()
63c58f1c22SToby Isaac */
64c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v)
65c58f1c22SToby Isaac {
66277ea44aSLisandro Dalcin   IS             is;
67b9907514SLisandro Dalcin   PetscInt       off = 0, *pointArray, p;
68c58f1c22SToby Isaac   PetscErrorCode ierr;
69c58f1c22SToby Isaac 
70b9907514SLisandro Dalcin   if (PetscLikely(v >= 0 && v < label->numStrata) && label->validIS[v]) return 0;
71c58f1c22SToby Isaac   PetscFunctionBegin;
720c3c4a36SLisandro 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);
73e8f14785SLisandro Dalcin   ierr = PetscHSetIGetSize(label->ht[v], &label->stratumSizes[v]);CHKERRQ(ierr);
74ad8374ffSToby Isaac   ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr);
75e8f14785SLisandro Dalcin   ierr = PetscHSetIGetElems(label->ht[v], &off, pointArray);CHKERRQ(ierr);
76b9907514SLisandro Dalcin   ierr = PetscHSetIClear(label->ht[v]);CHKERRQ(ierr);
77ad8374ffSToby Isaac   ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr);
78c58f1c22SToby Isaac   if (label->bt) {
79c58f1c22SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
80ad8374ffSToby Isaac       const PetscInt point = pointArray[p];
81c58f1c22SToby 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);
82c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
83c58f1c22SToby Isaac     }
84c58f1c22SToby Isaac   }
85277ea44aSLisandro Dalcin   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->stratumSizes[v], pointArray, PETSC_OWN_POINTER, &is);CHKERRQ(ierr);
86277ea44aSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject) is, "indices");CHKERRQ(ierr);
87277ea44aSLisandro Dalcin   label->points[v]  = is;
88ad8374ffSToby Isaac   label->validIS[v] = PETSC_TRUE;
89d67d17b1SMatthew G. Knepley   ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
90c58f1c22SToby Isaac   PetscFunctionReturn(0);
91c58f1c22SToby Isaac }
92c58f1c22SToby Isaac 
93c58f1c22SToby Isaac /*
94c58f1c22SToby Isaac   DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format
95c58f1c22SToby Isaac 
965b5e7992SMatthew G. Knepley   Not collective
975b5e7992SMatthew G. Knepley 
98c58f1c22SToby Isaac   Input parameter:
99c58f1c22SToby Isaac . label - The DMLabel
100c58f1c22SToby Isaac 
101c58f1c22SToby Isaac   Output parameter:
102c58f1c22SToby Isaac . label - The DMLabel with all strata in sorted list format
103c58f1c22SToby Isaac 
104c58f1c22SToby Isaac   Level: developer
105c58f1c22SToby Isaac 
106c58f1c22SToby Isaac .seealso: DMLabelCreate()
107c58f1c22SToby Isaac */
108c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label)
109c58f1c22SToby Isaac {
110c58f1c22SToby Isaac   PetscInt       v;
111c58f1c22SToby Isaac   PetscErrorCode ierr;
112c58f1c22SToby Isaac 
113c58f1c22SToby Isaac   PetscFunctionBegin;
114c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; v++){
115c58f1c22SToby Isaac     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
116c58f1c22SToby Isaac   }
117c58f1c22SToby Isaac   PetscFunctionReturn(0);
118c58f1c22SToby Isaac }
119c58f1c22SToby Isaac 
120c58f1c22SToby Isaac /*
121c58f1c22SToby Isaac   DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format
122c58f1c22SToby Isaac 
1235b5e7992SMatthew G. Knepley   Not collective
1245b5e7992SMatthew G. Knepley 
125c58f1c22SToby Isaac   Input parameter:
126c58f1c22SToby Isaac + label - The DMLabel
127c58f1c22SToby Isaac - v - The stratum value
128c58f1c22SToby Isaac 
129c58f1c22SToby Isaac   Output parameter:
130c58f1c22SToby Isaac . label - The DMLabel with stratum in hash format
131c58f1c22SToby Isaac 
132c58f1c22SToby Isaac   Level: developer
133c58f1c22SToby Isaac 
134c58f1c22SToby Isaac .seealso: DMLabelCreate()
135c58f1c22SToby Isaac */
136c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v)
137c58f1c22SToby Isaac {
138c58f1c22SToby Isaac   PetscInt       p;
139ad8374ffSToby Isaac   const PetscInt *points;
140c58f1c22SToby Isaac   PetscErrorCode ierr;
141c58f1c22SToby Isaac 
142b9907514SLisandro Dalcin   if (PetscLikely(v >= 0 && v < label->numStrata) && !label->validIS[v]) return 0;
143c58f1c22SToby Isaac   PetscFunctionBegin;
1440c3c4a36SLisandro 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);
145ad8374ffSToby Isaac   if (label->points[v]) {
146ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
147e8f14785SLisandro Dalcin     for (p = 0; p < label->stratumSizes[v]; ++p) {
148e8f14785SLisandro Dalcin       ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);
149e8f14785SLisandro Dalcin     }
150ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
151ad8374ffSToby Isaac     ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
152ad8374ffSToby Isaac   }
153ad8374ffSToby Isaac   label->validIS[v] = PETSC_FALSE;
154c58f1c22SToby Isaac   PetscFunctionReturn(0);
155c58f1c22SToby Isaac }
156c58f1c22SToby Isaac 
157b9907514SLisandro Dalcin #if !defined(DMLABEL_LOOKUP_THRESHOLD)
158b9907514SLisandro Dalcin #define DMLABEL_LOOKUP_THRESHOLD 16
159b9907514SLisandro Dalcin #endif
160b9907514SLisandro Dalcin 
1610c3c4a36SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index)
1620c3c4a36SLisandro Dalcin {
1630c3c4a36SLisandro Dalcin   PetscInt       v;
164b9907514SLisandro Dalcin   PetscErrorCode ierr;
1650e79e033SBarry Smith 
1660c3c4a36SLisandro Dalcin   PetscFunctionBegin;
1670e79e033SBarry Smith   *index = -1;
168b9907514SLisandro Dalcin   if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) {
169b9907514SLisandro Dalcin     for (v = 0; v < label->numStrata; ++v)
170b9907514SLisandro Dalcin       if (label->stratumValues[v] == value) {*index = v; break;}
171b9907514SLisandro Dalcin   } else {
172b9907514SLisandro Dalcin     ierr = PetscHMapIGet(label->hmap, value, index);CHKERRQ(ierr);
1730e79e033SBarry Smith   }
17490e9b2aeSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
17590e9b2aeSLisandro Dalcin   { /* Check strata hash map consistency */
17690e9b2aeSLisandro Dalcin     PetscInt len, loc = -1;
17790e9b2aeSLisandro Dalcin     ierr = PetscHMapIGetSize(label->hmap, &len);CHKERRQ(ierr);
17890e9b2aeSLisandro Dalcin     if (len != label->numStrata) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map size");
17990e9b2aeSLisandro Dalcin     if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) {
18090e9b2aeSLisandro Dalcin       ierr = PetscHMapIGet(label->hmap, value, &loc);CHKERRQ(ierr);
18190e9b2aeSLisandro Dalcin     } else {
18290e9b2aeSLisandro Dalcin       for (v = 0; v < label->numStrata; ++v)
18390e9b2aeSLisandro Dalcin         if (label->stratumValues[v] == value) {loc = v; break;}
18490e9b2aeSLisandro Dalcin     }
18590e9b2aeSLisandro Dalcin     if (loc != *index) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map lookup");
18690e9b2aeSLisandro Dalcin   }
18790e9b2aeSLisandro Dalcin #endif
1880c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
1890c3c4a36SLisandro Dalcin }
1900c3c4a36SLisandro Dalcin 
191b9907514SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index)
192c58f1c22SToby Isaac {
193b9907514SLisandro Dalcin   PetscInt       v;
194b9907514SLisandro Dalcin   PetscInt      *tmpV;
195b9907514SLisandro Dalcin   PetscInt      *tmpS;
196b9907514SLisandro Dalcin   PetscHSetI    *tmpH, ht;
197b9907514SLisandro Dalcin   IS            *tmpP, is;
198c58f1c22SToby Isaac   PetscBool     *tmpB;
199b9907514SLisandro Dalcin   PetscHMapI     hmap = label->hmap;
200c58f1c22SToby Isaac   PetscErrorCode ierr;
201c58f1c22SToby Isaac 
202c58f1c22SToby Isaac   PetscFunctionBegin;
203b9907514SLisandro Dalcin   v    = label->numStrata;
204b9907514SLisandro Dalcin   tmpV = label->stratumValues;
205b9907514SLisandro Dalcin   tmpS = label->stratumSizes;
206b9907514SLisandro Dalcin   tmpH = label->ht;
207b9907514SLisandro Dalcin   tmpP = label->points;
208b9907514SLisandro Dalcin   tmpB = label->validIS;
2098e1f8cf0SLisandro Dalcin   { /* TODO: PetscRealloc() is broken, use malloc+memcpy+free  */
2108e1f8cf0SLisandro Dalcin     PetscInt   *oldV = tmpV;
2118e1f8cf0SLisandro Dalcin     PetscInt   *oldS = tmpS;
2128e1f8cf0SLisandro Dalcin     PetscHSetI *oldH = tmpH;
2138e1f8cf0SLisandro Dalcin     IS         *oldP = tmpP;
2148e1f8cf0SLisandro Dalcin     PetscBool  *oldB = tmpB;
2158e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpV), &tmpV);CHKERRQ(ierr);
2168e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpS), &tmpS);CHKERRQ(ierr);
2178e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpH), &tmpH);CHKERRQ(ierr);
2188e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpP), &tmpP);CHKERRQ(ierr);
2198e1f8cf0SLisandro Dalcin     ierr = PetscMalloc((v+1)*sizeof(*tmpB), &tmpB);CHKERRQ(ierr);
220580bdb30SBarry Smith     ierr = PetscArraycpy(tmpV, oldV, v);CHKERRQ(ierr);
221580bdb30SBarry Smith     ierr = PetscArraycpy(tmpS, oldS, v);CHKERRQ(ierr);
222580bdb30SBarry Smith     ierr = PetscArraycpy(tmpH, oldH, v);CHKERRQ(ierr);
223580bdb30SBarry Smith     ierr = PetscArraycpy(tmpP, oldP, v);CHKERRQ(ierr);
224580bdb30SBarry Smith     ierr = PetscArraycpy(tmpB, oldB, v);CHKERRQ(ierr);
2258e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldV);CHKERRQ(ierr);
2268e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldS);CHKERRQ(ierr);
2278e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldH);CHKERRQ(ierr);
2288e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldP);CHKERRQ(ierr);
2298e1f8cf0SLisandro Dalcin     ierr = PetscFree(oldB);CHKERRQ(ierr);
2308e1f8cf0SLisandro Dalcin   }
231b9907514SLisandro Dalcin   label->numStrata     = v+1;
232c58f1c22SToby Isaac   label->stratumValues = tmpV;
233c58f1c22SToby Isaac   label->stratumSizes  = tmpS;
234c58f1c22SToby Isaac   label->ht            = tmpH;
235c58f1c22SToby Isaac   label->points        = tmpP;
236ad8374ffSToby Isaac   label->validIS       = tmpB;
237b9907514SLisandro Dalcin   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
238b9907514SLisandro Dalcin   ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);CHKERRQ(ierr);
239b9907514SLisandro Dalcin   ierr = PetscHMapISet(hmap, value, v);CHKERRQ(ierr);
240b9907514SLisandro Dalcin   tmpV[v] = value;
241b9907514SLisandro Dalcin   tmpS[v] = 0;
242b9907514SLisandro Dalcin   tmpH[v] = ht;
243b9907514SLisandro Dalcin   tmpP[v] = is;
244b9907514SLisandro Dalcin   tmpB[v] = PETSC_TRUE;
245277ea44aSLisandro Dalcin   ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
2460c3c4a36SLisandro Dalcin   *index = v;
2470c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
2480c3c4a36SLisandro Dalcin }
2490c3c4a36SLisandro Dalcin 
250b9907514SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupAddStratum(DMLabel label, PetscInt value, PetscInt *index)
251b9907514SLisandro Dalcin {
252b9907514SLisandro Dalcin   PetscErrorCode ierr;
253b9907514SLisandro Dalcin   PetscFunctionBegin;
254b9907514SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, index);CHKERRQ(ierr);
255b9907514SLisandro Dalcin   if (*index < 0) {ierr = DMLabelNewStratum(label, value, index);CHKERRQ(ierr);}
256b9907514SLisandro Dalcin   PetscFunctionReturn(0);
257b9907514SLisandro Dalcin }
258b9907514SLisandro Dalcin 
259b9907514SLisandro Dalcin /*@
260b9907514SLisandro Dalcin   DMLabelAddStratum - Adds a new stratum value in a DMLabel
261b9907514SLisandro Dalcin 
262b9907514SLisandro Dalcin   Input Parameter:
263b9907514SLisandro Dalcin + label - The DMLabel
264b9907514SLisandro Dalcin - value - The stratum value
265b9907514SLisandro Dalcin 
266b9907514SLisandro Dalcin   Level: beginner
267b9907514SLisandro Dalcin 
268b9907514SLisandro Dalcin .seealso:  DMLabelCreate(), DMLabelDestroy()
269b9907514SLisandro Dalcin @*/
2700c3c4a36SLisandro Dalcin PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
2710c3c4a36SLisandro Dalcin {
2720c3c4a36SLisandro Dalcin   PetscInt       v;
2730c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
2740c3c4a36SLisandro Dalcin 
2750c3c4a36SLisandro Dalcin   PetscFunctionBegin;
276d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
277b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
278b9907514SLisandro Dalcin   PetscFunctionReturn(0);
279b9907514SLisandro Dalcin }
280b9907514SLisandro Dalcin 
281b9907514SLisandro Dalcin /*@
282b9907514SLisandro Dalcin   DMLabelAddStrata - Adds new stratum values in a DMLabel
283b9907514SLisandro Dalcin 
2845b5e7992SMatthew G. Knepley   Not collective
2855b5e7992SMatthew G. Knepley 
286b9907514SLisandro Dalcin   Input Parameter:
287b9907514SLisandro Dalcin + label - The DMLabel
288b9907514SLisandro Dalcin . numStrata - The number of stratum values
289b9907514SLisandro Dalcin - stratumValues - The stratum values
290b9907514SLisandro Dalcin 
291b9907514SLisandro Dalcin   Level: beginner
292b9907514SLisandro Dalcin 
293b9907514SLisandro Dalcin .seealso:  DMLabelCreate(), DMLabelDestroy()
294b9907514SLisandro Dalcin @*/
295b9907514SLisandro Dalcin PetscErrorCode DMLabelAddStrata(DMLabel label, PetscInt numStrata, const PetscInt stratumValues[])
296b9907514SLisandro Dalcin {
297b9907514SLisandro Dalcin   PetscInt       *values, v;
298b9907514SLisandro Dalcin   PetscErrorCode ierr;
299b9907514SLisandro Dalcin 
300b9907514SLisandro Dalcin   PetscFunctionBegin;
301b9907514SLisandro Dalcin   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
302b9907514SLisandro Dalcin   if (numStrata) PetscValidIntPointer(stratumValues, 3);
303b9907514SLisandro Dalcin   ierr = PetscMalloc1(numStrata, &values);CHKERRQ(ierr);
304580bdb30SBarry Smith   ierr = PetscArraycpy(values, stratumValues, numStrata);CHKERRQ(ierr);
305b9907514SLisandro Dalcin   ierr = PetscSortRemoveDupsInt(&numStrata, values);CHKERRQ(ierr);
306b9907514SLisandro Dalcin   if (!label->numStrata) { /* Fast preallocation */
307b9907514SLisandro Dalcin     PetscInt   *tmpV;
308b9907514SLisandro Dalcin     PetscInt   *tmpS;
309b9907514SLisandro Dalcin     PetscHSetI *tmpH, ht;
310b9907514SLisandro Dalcin     IS         *tmpP, is;
311b9907514SLisandro Dalcin     PetscBool  *tmpB;
312b9907514SLisandro Dalcin     PetscHMapI  hmap = label->hmap;
313b9907514SLisandro Dalcin 
314b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpV);CHKERRQ(ierr);
315b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpS);CHKERRQ(ierr);
316b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpH);CHKERRQ(ierr);
317b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpP);CHKERRQ(ierr);
318b9907514SLisandro Dalcin     ierr = PetscMalloc1(numStrata, &tmpB);CHKERRQ(ierr);
319b9907514SLisandro Dalcin     label->numStrata     = numStrata;
320b9907514SLisandro Dalcin     label->stratumValues = tmpV;
321b9907514SLisandro Dalcin     label->stratumSizes  = tmpS;
322b9907514SLisandro Dalcin     label->ht            = tmpH;
323b9907514SLisandro Dalcin     label->points        = tmpP;
324b9907514SLisandro Dalcin     label->validIS       = tmpB;
325b9907514SLisandro Dalcin     for (v = 0; v < numStrata; ++v) {
326b9907514SLisandro Dalcin       ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
327b9907514SLisandro Dalcin       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);CHKERRQ(ierr);
328b9907514SLisandro Dalcin       ierr = PetscHMapISet(hmap, values[v], v);CHKERRQ(ierr);
329b9907514SLisandro Dalcin       tmpV[v] = values[v];
330b9907514SLisandro Dalcin       tmpS[v] = 0;
331b9907514SLisandro Dalcin       tmpH[v] = ht;
332b9907514SLisandro Dalcin       tmpP[v] = is;
333b9907514SLisandro Dalcin       tmpB[v] = PETSC_TRUE;
334b9907514SLisandro Dalcin     }
335277ea44aSLisandro Dalcin     ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
336b9907514SLisandro Dalcin   } else {
337b9907514SLisandro Dalcin     for (v = 0; v < numStrata; ++v) {
338b9907514SLisandro Dalcin       ierr = DMLabelAddStratum(label, values[v]);CHKERRQ(ierr);
339b9907514SLisandro Dalcin     }
340b9907514SLisandro Dalcin   }
341b9907514SLisandro Dalcin   ierr = PetscFree(values);CHKERRQ(ierr);
342b9907514SLisandro Dalcin   PetscFunctionReturn(0);
343b9907514SLisandro Dalcin }
344b9907514SLisandro Dalcin 
345b9907514SLisandro Dalcin /*@
346b9907514SLisandro Dalcin   DMLabelAddStrataIS - Adds new stratum values in a DMLabel
347b9907514SLisandro Dalcin 
3485b5e7992SMatthew G. Knepley   Not collective
3495b5e7992SMatthew G. Knepley 
350b9907514SLisandro Dalcin   Input Parameter:
351b9907514SLisandro Dalcin + label - The DMLabel
352b9907514SLisandro Dalcin - valueIS - Index set with stratum values
353b9907514SLisandro Dalcin 
354b9907514SLisandro Dalcin   Level: beginner
355b9907514SLisandro Dalcin 
356b9907514SLisandro Dalcin .seealso:  DMLabelCreate(), DMLabelDestroy()
357b9907514SLisandro Dalcin @*/
358b9907514SLisandro Dalcin PetscErrorCode DMLabelAddStrataIS(DMLabel label, IS valueIS)
359b9907514SLisandro Dalcin {
360b9907514SLisandro Dalcin   PetscInt       numStrata;
361b9907514SLisandro Dalcin   const PetscInt *stratumValues;
362b9907514SLisandro Dalcin   PetscErrorCode ierr;
363b9907514SLisandro Dalcin 
364b9907514SLisandro Dalcin   PetscFunctionBegin;
365b9907514SLisandro Dalcin   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
366b9907514SLisandro Dalcin   PetscValidHeaderSpecific(valueIS, IS_CLASSID, 2);
367b9907514SLisandro Dalcin   ierr = ISGetLocalSize(valueIS, &numStrata);CHKERRQ(ierr);
368b9907514SLisandro Dalcin   ierr = ISGetIndices(valueIS, &stratumValues);CHKERRQ(ierr);
369b9907514SLisandro Dalcin   ierr = DMLabelAddStrata(label, numStrata, stratumValues);CHKERRQ(ierr);
370c58f1c22SToby Isaac   PetscFunctionReturn(0);
371c58f1c22SToby Isaac }
372c58f1c22SToby Isaac 
373c58f1c22SToby Isaac static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
374c58f1c22SToby Isaac {
375c58f1c22SToby Isaac   PetscInt       v;
376c58f1c22SToby Isaac   PetscMPIInt    rank;
377c58f1c22SToby Isaac   PetscErrorCode ierr;
378c58f1c22SToby Isaac 
379c58f1c22SToby Isaac   PetscFunctionBegin;
380c58f1c22SToby Isaac   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr);
381c58f1c22SToby Isaac   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
382c58f1c22SToby Isaac   if (label) {
383d67d17b1SMatthew G. Knepley     const char *name;
384d67d17b1SMatthew G. Knepley 
385d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
386d67d17b1SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", name);CHKERRQ(ierr);
387c58f1c22SToby Isaac     if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);}
388c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
389c58f1c22SToby Isaac       const PetscInt value = label->stratumValues[v];
390ad8374ffSToby Isaac       const PetscInt *points;
391c58f1c22SToby Isaac       PetscInt       p;
392c58f1c22SToby Isaac 
393ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
394c58f1c22SToby Isaac       for (p = 0; p < label->stratumSizes[v]; ++p) {
395ad8374ffSToby Isaac         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr);
396c58f1c22SToby Isaac       }
397ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
398c58f1c22SToby Isaac     }
399c58f1c22SToby Isaac   }
400c58f1c22SToby Isaac   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
401c58f1c22SToby Isaac   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
402c58f1c22SToby Isaac   PetscFunctionReturn(0);
403c58f1c22SToby Isaac }
404c58f1c22SToby Isaac 
405c58f1c22SToby Isaac /*@C
406c58f1c22SToby Isaac   DMLabelView - View the label
407c58f1c22SToby Isaac 
4085b5e7992SMatthew G. Knepley   Collective on viewer
4095b5e7992SMatthew G. Knepley 
410c58f1c22SToby Isaac   Input Parameters:
411c58f1c22SToby Isaac + label - The DMLabel
412c58f1c22SToby Isaac - viewer - The PetscViewer
413c58f1c22SToby Isaac 
414c58f1c22SToby Isaac   Level: intermediate
415c58f1c22SToby Isaac 
416c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy()
417c58f1c22SToby Isaac @*/
418c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
419c58f1c22SToby Isaac {
420c58f1c22SToby Isaac   PetscBool      iascii;
421c58f1c22SToby Isaac   PetscErrorCode ierr;
422c58f1c22SToby Isaac 
423c58f1c22SToby Isaac   PetscFunctionBegin;
424d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
425b9907514SLisandro Dalcin   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)label), &viewer);CHKERRQ(ierr);}
426c58f1c22SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
427c58f1c22SToby Isaac   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
428c58f1c22SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
429c58f1c22SToby Isaac   if (iascii) {
430c58f1c22SToby Isaac     ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr);
431c58f1c22SToby Isaac   }
432c58f1c22SToby Isaac   PetscFunctionReturn(0);
433c58f1c22SToby Isaac }
434c58f1c22SToby Isaac 
43584f0b6dfSMatthew G. Knepley /*@
436d67d17b1SMatthew G. Knepley   DMLabelReset - Destroys internal data structures in a DMLabel
437d67d17b1SMatthew G. Knepley 
4385b5e7992SMatthew G. Knepley   Not collective
4395b5e7992SMatthew G. Knepley 
440d67d17b1SMatthew G. Knepley   Input Parameter:
441d67d17b1SMatthew G. Knepley . label - The DMLabel
442d67d17b1SMatthew G. Knepley 
443d67d17b1SMatthew G. Knepley   Level: beginner
444d67d17b1SMatthew G. Knepley 
445d67d17b1SMatthew G. Knepley .seealso: DMLabelDestroy(), DMLabelCreate()
446d67d17b1SMatthew G. Knepley @*/
447d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelReset(DMLabel label)
448d67d17b1SMatthew G. Knepley {
449d67d17b1SMatthew G. Knepley   PetscInt       v;
450d67d17b1SMatthew G. Knepley   PetscErrorCode ierr;
451d67d17b1SMatthew G. Knepley 
452d67d17b1SMatthew G. Knepley   PetscFunctionBegin;
453d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
454d67d17b1SMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
455d67d17b1SMatthew G. Knepley     ierr = PetscHSetIDestroy(&label->ht[v]);CHKERRQ(ierr);
456d67d17b1SMatthew G. Knepley     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
457d67d17b1SMatthew G. Knepley   }
458b9907514SLisandro Dalcin   label->numStrata = 0;
459b9907514SLisandro Dalcin   ierr = PetscFree(label->stratumValues);CHKERRQ(ierr);
460b9907514SLisandro Dalcin   ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr);
461d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->ht);CHKERRQ(ierr);
462d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->points);CHKERRQ(ierr);
463d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->validIS);CHKERRQ(ierr);
464b9907514SLisandro Dalcin   ierr = PetscHMapIReset(label->hmap);CHKERRQ(ierr);
465b9907514SLisandro Dalcin   label->pStart = -1;
466b9907514SLisandro Dalcin   label->pEnd   = -1;
467d67d17b1SMatthew G. Knepley   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
468d67d17b1SMatthew G. Knepley   PetscFunctionReturn(0);
469d67d17b1SMatthew G. Knepley }
470d67d17b1SMatthew G. Knepley 
471d67d17b1SMatthew G. Knepley /*@
47284f0b6dfSMatthew G. Knepley   DMLabelDestroy - Destroys a DMLabel
47384f0b6dfSMatthew G. Knepley 
4745b5e7992SMatthew G. Knepley   Collective on label
4755b5e7992SMatthew G. Knepley 
47684f0b6dfSMatthew G. Knepley   Input Parameter:
47784f0b6dfSMatthew G. Knepley . label - The DMLabel
47884f0b6dfSMatthew G. Knepley 
47984f0b6dfSMatthew G. Knepley   Level: beginner
48084f0b6dfSMatthew G. Knepley 
481d67d17b1SMatthew G. Knepley .seealso: DMLabelReset(), DMLabelCreate()
48284f0b6dfSMatthew G. Knepley @*/
483c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label)
484c58f1c22SToby Isaac {
485c58f1c22SToby Isaac   PetscErrorCode ierr;
486c58f1c22SToby Isaac 
487c58f1c22SToby Isaac   PetscFunctionBegin;
488d67d17b1SMatthew G. Knepley   if (!*label) PetscFunctionReturn(0);
489d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific((*label),DMLABEL_CLASSID,1);
490b9907514SLisandro Dalcin   if (--((PetscObject)(*label))->refct > 0) {*label = NULL; PetscFunctionReturn(0);}
491d67d17b1SMatthew G. Knepley   ierr = DMLabelReset(*label);CHKERRQ(ierr);
492b9907514SLisandro Dalcin   ierr = PetscHMapIDestroy(&(*label)->hmap);CHKERRQ(ierr);
493d67d17b1SMatthew G. Knepley   ierr = PetscHeaderDestroy(label);CHKERRQ(ierr);
494c58f1c22SToby Isaac   PetscFunctionReturn(0);
495c58f1c22SToby Isaac }
496c58f1c22SToby Isaac 
49784f0b6dfSMatthew G. Knepley /*@
49884f0b6dfSMatthew G. Knepley   DMLabelDuplicate - Duplicates a DMLabel
49984f0b6dfSMatthew G. Knepley 
5005b5e7992SMatthew G. Knepley   Collective on label
5015b5e7992SMatthew G. Knepley 
50284f0b6dfSMatthew G. Knepley   Input Parameter:
50384f0b6dfSMatthew G. Knepley . label - The DMLabel
50484f0b6dfSMatthew G. Knepley 
50584f0b6dfSMatthew G. Knepley   Output Parameter:
50684f0b6dfSMatthew G. Knepley . labelnew - location to put new vector
50784f0b6dfSMatthew G. Knepley 
50884f0b6dfSMatthew G. Knepley   Level: intermediate
50984f0b6dfSMatthew G. Knepley 
51084f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelDestroy()
51184f0b6dfSMatthew G. Knepley @*/
512c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
513c58f1c22SToby Isaac {
514d67d17b1SMatthew G. Knepley   const char    *name;
515ad8374ffSToby Isaac   PetscInt       v;
516c58f1c22SToby Isaac   PetscErrorCode ierr;
517c58f1c22SToby Isaac 
518c58f1c22SToby Isaac   PetscFunctionBegin;
519d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
520c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
521d67d17b1SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
522d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PetscObjectComm((PetscObject) label), name, labelnew);CHKERRQ(ierr);
523c58f1c22SToby Isaac 
524c58f1c22SToby Isaac   (*labelnew)->numStrata    = label->numStrata;
5255aa44df4SToby Isaac   (*labelnew)->defaultValue = label->defaultValue;
526c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr);
527c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr);
528c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr);
529c58f1c22SToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr);
530ad8374ffSToby Isaac   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr);
531c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
532e8f14785SLisandro Dalcin     ierr = PetscHSetICreate(&(*labelnew)->ht[v]);CHKERRQ(ierr);
533c58f1c22SToby Isaac     (*labelnew)->stratumValues[v]  = label->stratumValues[v];
534c58f1c22SToby Isaac     (*labelnew)->stratumSizes[v]   = label->stratumSizes[v];
535ad8374ffSToby Isaac     ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr);
536ad8374ffSToby Isaac     (*labelnew)->points[v]         = label->points[v];
537b9907514SLisandro Dalcin     (*labelnew)->validIS[v]        = PETSC_TRUE;
538c58f1c22SToby Isaac   }
539f14fe40dSLisandro Dalcin   ierr = PetscHMapIDestroy(&(*labelnew)->hmap);CHKERRQ(ierr);
540b9907514SLisandro Dalcin   ierr = PetscHMapIDuplicate(label->hmap,&(*labelnew)->hmap);CHKERRQ(ierr);
541c58f1c22SToby Isaac   (*labelnew)->pStart = -1;
542c58f1c22SToby Isaac   (*labelnew)->pEnd   = -1;
543c58f1c22SToby Isaac   (*labelnew)->bt     = NULL;
544c58f1c22SToby Isaac   PetscFunctionReturn(0);
545c58f1c22SToby Isaac }
546c58f1c22SToby Isaac 
547c6a43d28SMatthew G. Knepley /*@
548c6a43d28SMatthew G. Knepley   DMLabelComputeIndex - Create an index structure for membership determination, automatically determining the bounds
549c6a43d28SMatthew G. Knepley 
5505b5e7992SMatthew G. Knepley   Not collective
5515b5e7992SMatthew G. Knepley 
552c6a43d28SMatthew G. Knepley   Input Parameter:
553c6a43d28SMatthew G. Knepley . label  - The DMLabel
554c6a43d28SMatthew G. Knepley 
555c6a43d28SMatthew G. Knepley   Level: intermediate
556c6a43d28SMatthew G. Knepley 
557c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
558c6a43d28SMatthew G. Knepley @*/
559c6a43d28SMatthew G. Knepley PetscErrorCode DMLabelComputeIndex(DMLabel label)
560c6a43d28SMatthew G. Knepley {
561c6a43d28SMatthew G. Knepley   PetscInt       pStart = PETSC_MAX_INT, pEnd = -1, v;
562c6a43d28SMatthew G. Knepley   PetscErrorCode ierr;
563c6a43d28SMatthew G. Knepley 
564c6a43d28SMatthew G. Knepley   PetscFunctionBegin;
565c6a43d28SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
566c6a43d28SMatthew G. Knepley   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
567c6a43d28SMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
568c6a43d28SMatthew G. Knepley     const PetscInt *points;
569c6a43d28SMatthew G. Knepley     PetscInt       i;
570c6a43d28SMatthew G. Knepley 
571c6a43d28SMatthew G. Knepley     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
572c6a43d28SMatthew G. Knepley     for (i = 0; i < label->stratumSizes[v]; ++i) {
573c6a43d28SMatthew G. Knepley       const PetscInt point = points[i];
574c6a43d28SMatthew G. Knepley 
575c6a43d28SMatthew G. Knepley       pStart = PetscMin(point,   pStart);
576c6a43d28SMatthew G. Knepley       pEnd   = PetscMax(point+1, pEnd);
577c6a43d28SMatthew G. Knepley     }
578c6a43d28SMatthew G. Knepley     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
579c6a43d28SMatthew G. Knepley   }
580c6a43d28SMatthew G. Knepley   label->pStart = pStart == PETSC_MAX_INT ? -1 : pStart;
581c6a43d28SMatthew G. Knepley   label->pEnd   = pEnd;
582c6a43d28SMatthew G. Knepley   ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
583c6a43d28SMatthew G. Knepley   PetscFunctionReturn(0);
584c6a43d28SMatthew G. Knepley }
585c6a43d28SMatthew G. Knepley 
586c6a43d28SMatthew G. Knepley /*@
587c6a43d28SMatthew G. Knepley   DMLabelCreateIndex - Create an index structure for membership determination
588c6a43d28SMatthew G. Knepley 
5895b5e7992SMatthew G. Knepley   Not collective
5905b5e7992SMatthew G. Knepley 
591c6a43d28SMatthew G. Knepley   Input Parameters:
592c6a43d28SMatthew G. Knepley + label  - The DMLabel
593c6a43d28SMatthew G. Knepley . pStart - The smallest point
594c6a43d28SMatthew G. Knepley - pEnd   - The largest point + 1
595c6a43d28SMatthew G. Knepley 
596c6a43d28SMatthew G. Knepley   Level: intermediate
597c6a43d28SMatthew G. Knepley 
598c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelComputeIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
599c6a43d28SMatthew G. Knepley @*/
600c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
601c58f1c22SToby Isaac {
602c58f1c22SToby Isaac   PetscInt       v;
603c58f1c22SToby Isaac   PetscErrorCode ierr;
604c58f1c22SToby Isaac 
605c58f1c22SToby Isaac   PetscFunctionBegin;
606d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
6070c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
608c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
609c58f1c22SToby Isaac   label->pStart = pStart;
610c58f1c22SToby Isaac   label->pEnd   = pEnd;
611c6a43d28SMatthew G. Knepley   /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
612c58f1c22SToby Isaac   ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr);
613c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
614ad8374ffSToby Isaac     const PetscInt *points;
615c58f1c22SToby Isaac     PetscInt       i;
616c58f1c22SToby Isaac 
617ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
618c58f1c22SToby Isaac     for (i = 0; i < label->stratumSizes[v]; ++i) {
619ad8374ffSToby Isaac       const PetscInt point = points[i];
620c58f1c22SToby Isaac 
621c58f1c22SToby 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);
622c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr);
623c58f1c22SToby Isaac     }
624ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
625c58f1c22SToby Isaac   }
626c58f1c22SToby Isaac   PetscFunctionReturn(0);
627c58f1c22SToby Isaac }
628c58f1c22SToby Isaac 
629c6a43d28SMatthew G. Knepley /*@
630c6a43d28SMatthew G. Knepley   DMLabelDestroyIndex - Destroy the index structure
631c6a43d28SMatthew G. Knepley 
6325b5e7992SMatthew G. Knepley   Not collective
6335b5e7992SMatthew G. Knepley 
634c6a43d28SMatthew G. Knepley   Input Parameter:
635c6a43d28SMatthew G. Knepley . label - the DMLabel
636c6a43d28SMatthew G. Knepley 
637c6a43d28SMatthew G. Knepley   Level: intermediate
638c6a43d28SMatthew G. Knepley 
639c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
640c6a43d28SMatthew G. Knepley @*/
641c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label)
642c58f1c22SToby Isaac {
643c58f1c22SToby Isaac   PetscErrorCode ierr;
644c58f1c22SToby Isaac 
645c58f1c22SToby Isaac   PetscFunctionBegin;
646d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
647c58f1c22SToby Isaac   label->pStart = -1;
648c58f1c22SToby Isaac   label->pEnd   = -1;
6490c3c4a36SLisandro Dalcin   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
650c58f1c22SToby Isaac   PetscFunctionReturn(0);
651c58f1c22SToby Isaac }
652c58f1c22SToby Isaac 
653c58f1c22SToby Isaac /*@
654c6a43d28SMatthew G. Knepley   DMLabelGetBounds - Return the smallest and largest point in the label
655c6a43d28SMatthew G. Knepley 
6565b5e7992SMatthew G. Knepley   Not collective
6575b5e7992SMatthew G. Knepley 
658c6a43d28SMatthew G. Knepley   Input Parameter:
659c6a43d28SMatthew G. Knepley . label - the DMLabel
660c6a43d28SMatthew G. Knepley 
661c6a43d28SMatthew G. Knepley   Output Parameters:
662c6a43d28SMatthew G. Knepley + pStart - The smallest point
663c6a43d28SMatthew G. Knepley - pEnd   - The largest point + 1
664c6a43d28SMatthew G. Knepley 
665c6a43d28SMatthew G. Knepley   Note: This will compute an index for the label if one does not exist.
666c6a43d28SMatthew G. Knepley 
667c6a43d28SMatthew G. Knepley   Level: intermediate
668c6a43d28SMatthew G. Knepley 
669c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
670c6a43d28SMatthew G. Knepley @*/
671c6a43d28SMatthew G. Knepley PetscErrorCode DMLabelGetBounds(DMLabel label, PetscInt *pStart, PetscInt *pEnd)
672c6a43d28SMatthew G. Knepley {
673c6a43d28SMatthew G. Knepley   PetscErrorCode ierr;
674c6a43d28SMatthew G. Knepley 
675c6a43d28SMatthew G. Knepley   PetscFunctionBegin;
676c6a43d28SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
677c6a43d28SMatthew G. Knepley   if ((label->pStart == -1) && (label->pEnd == -1)) {ierr = DMLabelComputeIndex(label);CHKERRQ(ierr);}
678c6a43d28SMatthew G. Knepley   if (pStart) {
679534a8f05SLisandro Dalcin     PetscValidIntPointer(pStart, 2);
680c6a43d28SMatthew G. Knepley     *pStart = label->pStart;
681c6a43d28SMatthew G. Knepley   }
682c6a43d28SMatthew G. Knepley   if (pEnd) {
683534a8f05SLisandro Dalcin     PetscValidIntPointer(pEnd, 3);
684c6a43d28SMatthew G. Knepley     *pEnd = label->pEnd;
685c6a43d28SMatthew G. Knepley   }
686c6a43d28SMatthew G. Knepley   PetscFunctionReturn(0);
687c6a43d28SMatthew G. Knepley }
688c6a43d28SMatthew G. Knepley 
689c6a43d28SMatthew G. Knepley /*@
690c58f1c22SToby Isaac   DMLabelHasValue - Determine whether a label assigns the value to any point
691c58f1c22SToby Isaac 
6925b5e7992SMatthew G. Knepley   Not collective
6935b5e7992SMatthew G. Knepley 
694c58f1c22SToby Isaac   Input Parameters:
695c58f1c22SToby Isaac + label - the DMLabel
696c58f1c22SToby Isaac - value - the value
697c58f1c22SToby Isaac 
698c58f1c22SToby Isaac   Output Parameter:
699c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point
700c58f1c22SToby Isaac 
701c58f1c22SToby Isaac   Level: developer
702c58f1c22SToby Isaac 
703c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
704c58f1c22SToby Isaac @*/
705c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
706c58f1c22SToby Isaac {
707c58f1c22SToby Isaac   PetscInt v;
7080c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
709c58f1c22SToby Isaac 
710c58f1c22SToby Isaac   PetscFunctionBegin;
711d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
712534a8f05SLisandro Dalcin   PetscValidBoolPointer(contains, 3);
7130c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
7140c3c4a36SLisandro Dalcin   *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE;
715c58f1c22SToby Isaac   PetscFunctionReturn(0);
716c58f1c22SToby Isaac }
717c58f1c22SToby Isaac 
718c58f1c22SToby Isaac /*@
719c58f1c22SToby Isaac   DMLabelHasPoint - Determine whether a label assigns a value to a point
720c58f1c22SToby Isaac 
7215b5e7992SMatthew G. Knepley   Not collective
7225b5e7992SMatthew G. Knepley 
723c58f1c22SToby Isaac   Input Parameters:
724c58f1c22SToby Isaac + label - the DMLabel
725c58f1c22SToby Isaac - point - the point
726c58f1c22SToby Isaac 
727c58f1c22SToby Isaac   Output Parameter:
728c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value
729c58f1c22SToby Isaac 
730c58f1c22SToby Isaac   Note: The user must call DMLabelCreateIndex() before this function.
731c58f1c22SToby Isaac 
732c58f1c22SToby Isaac   Level: developer
733c58f1c22SToby Isaac 
734c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
735c58f1c22SToby Isaac @*/
736c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
737c58f1c22SToby Isaac {
738c58f1c22SToby Isaac   PetscErrorCode ierr;
739c58f1c22SToby Isaac 
740c58f1c22SToby Isaac   PetscFunctionBeginHot;
741d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
742534a8f05SLisandro Dalcin   PetscValidBoolPointer(contains, 3);
743c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
744c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG)
745c58f1c22SToby Isaac   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
746c58f1c22SToby 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);
747c58f1c22SToby Isaac #endif
748c58f1c22SToby Isaac   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
749c58f1c22SToby Isaac   PetscFunctionReturn(0);
750c58f1c22SToby Isaac }
751c58f1c22SToby Isaac 
752c58f1c22SToby Isaac /*@
753c58f1c22SToby Isaac   DMLabelStratumHasPoint - Return true if the stratum contains a point
754c58f1c22SToby Isaac 
7555b5e7992SMatthew G. Knepley   Not collective
7565b5e7992SMatthew G. Knepley 
757c58f1c22SToby Isaac   Input Parameters:
758c58f1c22SToby Isaac + label - the DMLabel
759c58f1c22SToby Isaac . value - the stratum value
760c58f1c22SToby Isaac - point - the point
761c58f1c22SToby Isaac 
762c58f1c22SToby Isaac   Output Parameter:
763c58f1c22SToby Isaac . contains - true if the stratum contains the point
764c58f1c22SToby Isaac 
765c58f1c22SToby Isaac   Level: intermediate
766c58f1c22SToby Isaac 
767c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
768c58f1c22SToby Isaac @*/
769c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
770c58f1c22SToby Isaac {
771c58f1c22SToby Isaac   PetscInt       v;
772c58f1c22SToby Isaac   PetscErrorCode ierr;
773c58f1c22SToby Isaac 
7740c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
775d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
776534a8f05SLisandro Dalcin   PetscValidBoolPointer(contains, 4);
777c58f1c22SToby Isaac   *contains = PETSC_FALSE;
7780c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
7790c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
7800c3c4a36SLisandro Dalcin 
781ad8374ffSToby Isaac   if (label->validIS[v]) {
782c58f1c22SToby Isaac     PetscInt i;
783c58f1c22SToby Isaac 
784a2d74346SToby Isaac     ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
7850c3c4a36SLisandro Dalcin     if (i >= 0) *contains = PETSC_TRUE;
786c58f1c22SToby Isaac   } else {
787c58f1c22SToby Isaac     PetscBool has;
788c58f1c22SToby Isaac 
789b9907514SLisandro Dalcin     ierr = PetscHSetIHas(label->ht[v], point, &has);CHKERRQ(ierr);
7900c3c4a36SLisandro Dalcin     if (has) *contains = PETSC_TRUE;
791c58f1c22SToby Isaac   }
792c58f1c22SToby Isaac   PetscFunctionReturn(0);
793c58f1c22SToby Isaac }
794c58f1c22SToby Isaac 
79584f0b6dfSMatthew G. Knepley /*@
7965aa44df4SToby Isaac   DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
7975aa44df4SToby Isaac   When a label is created, it is initialized to -1.
7985aa44df4SToby Isaac 
7995b5e7992SMatthew G. Knepley   Not collective
8005b5e7992SMatthew G. Knepley 
8015aa44df4SToby Isaac   Input parameter:
8025aa44df4SToby Isaac . label - a DMLabel object
8035aa44df4SToby Isaac 
8045aa44df4SToby Isaac   Output parameter:
8055aa44df4SToby Isaac . defaultValue - the default value
8065aa44df4SToby Isaac 
8075aa44df4SToby Isaac   Level: beginner
8085aa44df4SToby Isaac 
8095aa44df4SToby Isaac .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
81084f0b6dfSMatthew G. Knepley @*/
8115aa44df4SToby Isaac PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
8125aa44df4SToby Isaac {
8135aa44df4SToby Isaac   PetscFunctionBegin;
814d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
8155aa44df4SToby Isaac   *defaultValue = label->defaultValue;
8165aa44df4SToby Isaac   PetscFunctionReturn(0);
8175aa44df4SToby Isaac }
8185aa44df4SToby Isaac 
81984f0b6dfSMatthew G. Knepley /*@
8205aa44df4SToby Isaac   DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
8215aa44df4SToby Isaac   When a label is created, it is initialized to -1.
8225aa44df4SToby Isaac 
8235b5e7992SMatthew G. Knepley   Not collective
8245b5e7992SMatthew G. Knepley 
8255aa44df4SToby Isaac   Input parameter:
8265aa44df4SToby Isaac . label - a DMLabel object
8275aa44df4SToby Isaac 
8285aa44df4SToby Isaac   Output parameter:
8295aa44df4SToby Isaac . defaultValue - the default value
8305aa44df4SToby Isaac 
8315aa44df4SToby Isaac   Level: beginner
8325aa44df4SToby Isaac 
8335aa44df4SToby Isaac .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
83484f0b6dfSMatthew G. Knepley @*/
8355aa44df4SToby Isaac PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
8365aa44df4SToby Isaac {
8375aa44df4SToby Isaac   PetscFunctionBegin;
838d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
8395aa44df4SToby Isaac   label->defaultValue = defaultValue;
8405aa44df4SToby Isaac   PetscFunctionReturn(0);
8415aa44df4SToby Isaac }
8425aa44df4SToby Isaac 
843c58f1c22SToby Isaac /*@
8445aa44df4SToby 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())
845c58f1c22SToby Isaac 
8465b5e7992SMatthew G. Knepley   Not collective
8475b5e7992SMatthew G. Knepley 
848c58f1c22SToby Isaac   Input Parameters:
849c58f1c22SToby Isaac + label - the DMLabel
850c58f1c22SToby Isaac - point - the point
851c58f1c22SToby Isaac 
852c58f1c22SToby Isaac   Output Parameter:
8538e68afcfSMatthew G. Knepley . value - The point value, or the default value (-1 by default)
854c58f1c22SToby Isaac 
855c58f1c22SToby Isaac   Level: intermediate
856c58f1c22SToby Isaac 
8575aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
858c58f1c22SToby Isaac @*/
859c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
860c58f1c22SToby Isaac {
861c58f1c22SToby Isaac   PetscInt       v;
862c58f1c22SToby Isaac   PetscErrorCode ierr;
863c58f1c22SToby Isaac 
8640c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
865d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
866c58f1c22SToby Isaac   PetscValidPointer(value, 3);
8675aa44df4SToby Isaac   *value = label->defaultValue;
868c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
869ad8374ffSToby Isaac     if (label->validIS[v]) {
870c58f1c22SToby Isaac       PetscInt i;
871c58f1c22SToby Isaac 
872a2d74346SToby Isaac       ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
873c58f1c22SToby Isaac       if (i >= 0) {
874c58f1c22SToby Isaac         *value = label->stratumValues[v];
875c58f1c22SToby Isaac         break;
876c58f1c22SToby Isaac       }
877c58f1c22SToby Isaac     } else {
878c58f1c22SToby Isaac       PetscBool has;
879c58f1c22SToby Isaac 
880b9907514SLisandro Dalcin       ierr = PetscHSetIHas(label->ht[v], point, &has);CHKERRQ(ierr);
881c58f1c22SToby Isaac       if (has) {
882c58f1c22SToby Isaac         *value = label->stratumValues[v];
883c58f1c22SToby Isaac         break;
884c58f1c22SToby Isaac       }
885c58f1c22SToby Isaac     }
886c58f1c22SToby Isaac   }
887c58f1c22SToby Isaac   PetscFunctionReturn(0);
888c58f1c22SToby Isaac }
889c58f1c22SToby Isaac 
890c58f1c22SToby Isaac /*@
891367003a6SStefano 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.
892c58f1c22SToby Isaac 
8935b5e7992SMatthew G. Knepley   Not collective
8945b5e7992SMatthew G. Knepley 
895c58f1c22SToby Isaac   Input Parameters:
896c58f1c22SToby Isaac + label - the DMLabel
897c58f1c22SToby Isaac . point - the point
898c58f1c22SToby Isaac - value - The point value
899c58f1c22SToby Isaac 
900c58f1c22SToby Isaac   Level: intermediate
901c58f1c22SToby Isaac 
9025aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
903c58f1c22SToby Isaac @*/
904c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
905c58f1c22SToby Isaac {
906c58f1c22SToby Isaac   PetscInt       v;
907c58f1c22SToby Isaac   PetscErrorCode ierr;
908c58f1c22SToby Isaac 
909c58f1c22SToby Isaac   PetscFunctionBegin;
910d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
9110c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
9125aa44df4SToby Isaac   if (value == label->defaultValue) PetscFunctionReturn(0);
913b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
914c58f1c22SToby Isaac   /* Set key */
9150c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
916e8f14785SLisandro Dalcin   ierr = PetscHSetIAdd(label->ht[v], point);CHKERRQ(ierr);
917c58f1c22SToby Isaac   PetscFunctionReturn(0);
918c58f1c22SToby Isaac }
919c58f1c22SToby Isaac 
920c58f1c22SToby Isaac /*@
921c58f1c22SToby Isaac   DMLabelClearValue - Clear the value a label assigns to a point
922c58f1c22SToby Isaac 
9235b5e7992SMatthew G. Knepley   Not collective
9245b5e7992SMatthew G. Knepley 
925c58f1c22SToby Isaac   Input Parameters:
926c58f1c22SToby Isaac + label - the DMLabel
927c58f1c22SToby Isaac . point - the point
928c58f1c22SToby Isaac - value - The point value
929c58f1c22SToby Isaac 
930c58f1c22SToby Isaac   Level: intermediate
931c58f1c22SToby Isaac 
932c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
933c58f1c22SToby Isaac @*/
934c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
935c58f1c22SToby Isaac {
936ad8374ffSToby Isaac   PetscInt       v;
937c58f1c22SToby Isaac   PetscErrorCode ierr;
938c58f1c22SToby Isaac 
939c58f1c22SToby Isaac   PetscFunctionBegin;
940d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
941c58f1c22SToby Isaac   /* Find label value */
9420c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
9430c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
9440c3c4a36SLisandro Dalcin 
945eeed21e7SToby Isaac   if (label->bt) {
946eeed21e7SToby 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);
947eeed21e7SToby Isaac     ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
948eeed21e7SToby Isaac   }
9490c3c4a36SLisandro Dalcin 
9500c3c4a36SLisandro Dalcin   /* Delete key */
9510c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
952e8f14785SLisandro Dalcin   ierr = PetscHSetIDel(label->ht[v], point);CHKERRQ(ierr);
953c58f1c22SToby Isaac   PetscFunctionReturn(0);
954c58f1c22SToby Isaac }
955c58f1c22SToby Isaac 
956c58f1c22SToby Isaac /*@
957c58f1c22SToby Isaac   DMLabelInsertIS - Set all points in the IS to a value
958c58f1c22SToby Isaac 
9595b5e7992SMatthew G. Knepley   Not collective
9605b5e7992SMatthew G. Knepley 
961c58f1c22SToby Isaac   Input Parameters:
962c58f1c22SToby Isaac + label - the DMLabel
963c58f1c22SToby Isaac . is    - the point IS
964c58f1c22SToby Isaac - value - The point value
965c58f1c22SToby Isaac 
966c58f1c22SToby Isaac   Level: intermediate
967c58f1c22SToby Isaac 
968c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
969c58f1c22SToby Isaac @*/
970c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
971c58f1c22SToby Isaac {
9720c3c4a36SLisandro Dalcin   PetscInt        v, n, p;
973c58f1c22SToby Isaac   const PetscInt *points;
974c58f1c22SToby Isaac   PetscErrorCode  ierr;
975c58f1c22SToby Isaac 
976c58f1c22SToby Isaac   PetscFunctionBegin;
977d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
978c58f1c22SToby Isaac   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
9790c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
9800c3c4a36SLisandro Dalcin   if (value == label->defaultValue) PetscFunctionReturn(0);
981b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
9820c3c4a36SLisandro Dalcin   /* Set keys */
9830c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
984c58f1c22SToby Isaac   ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr);
985c58f1c22SToby Isaac   ierr = ISGetIndices(is, &points);CHKERRQ(ierr);
986e8f14785SLisandro Dalcin   for (p = 0; p < n; ++p) {ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);}
987c58f1c22SToby Isaac   ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr);
988c58f1c22SToby Isaac   PetscFunctionReturn(0);
989c58f1c22SToby Isaac }
990c58f1c22SToby Isaac 
99184f0b6dfSMatthew G. Knepley /*@
99284f0b6dfSMatthew G. Knepley   DMLabelGetNumValues - Get the number of values that the DMLabel takes
99384f0b6dfSMatthew G. Knepley 
9945b5e7992SMatthew G. Knepley   Not collective
9955b5e7992SMatthew G. Knepley 
99684f0b6dfSMatthew G. Knepley   Input Parameter:
99784f0b6dfSMatthew G. Knepley . label - the DMLabel
99884f0b6dfSMatthew G. Knepley 
99984f0b6dfSMatthew G. Knepley   Output Paramater:
100084f0b6dfSMatthew G. Knepley . numValues - the number of values
100184f0b6dfSMatthew G. Knepley 
100284f0b6dfSMatthew G. Knepley   Level: intermediate
100384f0b6dfSMatthew G. Knepley 
100484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
100584f0b6dfSMatthew G. Knepley @*/
1006c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
1007c58f1c22SToby Isaac {
1008c58f1c22SToby Isaac   PetscFunctionBegin;
1009d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1010c58f1c22SToby Isaac   PetscValidPointer(numValues, 2);
1011c58f1c22SToby Isaac   *numValues = label->numStrata;
1012c58f1c22SToby Isaac   PetscFunctionReturn(0);
1013c58f1c22SToby Isaac }
1014c58f1c22SToby Isaac 
101584f0b6dfSMatthew G. Knepley /*@
101684f0b6dfSMatthew G. Knepley   DMLabelGetValueIS - Get an IS of all values that the DMlabel takes
101784f0b6dfSMatthew G. Knepley 
10185b5e7992SMatthew G. Knepley   Not collective
10195b5e7992SMatthew G. Knepley 
102084f0b6dfSMatthew G. Knepley   Input Parameter:
102184f0b6dfSMatthew G. Knepley . label - the DMLabel
102284f0b6dfSMatthew G. Knepley 
102384f0b6dfSMatthew G. Knepley   Output Paramater:
102484f0b6dfSMatthew G. Knepley . is    - the value IS
102584f0b6dfSMatthew G. Knepley 
102684f0b6dfSMatthew G. Knepley   Level: intermediate
102784f0b6dfSMatthew G. Knepley 
102884f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
102984f0b6dfSMatthew G. Knepley @*/
1030c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
1031c58f1c22SToby Isaac {
1032c58f1c22SToby Isaac   PetscErrorCode ierr;
1033c58f1c22SToby Isaac 
1034c58f1c22SToby Isaac   PetscFunctionBegin;
1035d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1036c58f1c22SToby Isaac   PetscValidPointer(values, 2);
1037c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr);
1038c58f1c22SToby Isaac   PetscFunctionReturn(0);
1039c58f1c22SToby Isaac }
1040c58f1c22SToby Isaac 
104184f0b6dfSMatthew G. Knepley /*@
104284f0b6dfSMatthew G. Knepley   DMLabelHasStratum - Determine whether points exist with the given value
104384f0b6dfSMatthew G. Knepley 
10445b5e7992SMatthew G. Knepley   Not collective
10455b5e7992SMatthew G. Knepley 
104684f0b6dfSMatthew G. Knepley   Input Parameters:
104784f0b6dfSMatthew G. Knepley + label - the DMLabel
104884f0b6dfSMatthew G. Knepley - value - the stratum value
104984f0b6dfSMatthew G. Knepley 
105084f0b6dfSMatthew G. Knepley   Output Paramater:
105184f0b6dfSMatthew G. Knepley . exists - Flag saying whether points exist
105284f0b6dfSMatthew G. Knepley 
105384f0b6dfSMatthew G. Knepley   Level: intermediate
105484f0b6dfSMatthew G. Knepley 
105584f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
105684f0b6dfSMatthew G. Knepley @*/
1057fada774cSMatthew G. Knepley PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
1058fada774cSMatthew G. Knepley {
1059fada774cSMatthew G. Knepley   PetscInt       v;
10600c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
1061fada774cSMatthew G. Knepley 
1062fada774cSMatthew G. Knepley   PetscFunctionBegin;
1063d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1064fada774cSMatthew G. Knepley   PetscValidPointer(exists, 3);
10650c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
10660c3c4a36SLisandro Dalcin   *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE;
1067fada774cSMatthew G. Knepley   PetscFunctionReturn(0);
1068fada774cSMatthew G. Knepley }
1069fada774cSMatthew G. Knepley 
107084f0b6dfSMatthew G. Knepley /*@
107184f0b6dfSMatthew G. Knepley   DMLabelGetStratumSize - Get the size of a stratum
107284f0b6dfSMatthew G. Knepley 
10735b5e7992SMatthew G. Knepley   Not collective
10745b5e7992SMatthew G. Knepley 
107584f0b6dfSMatthew G. Knepley   Input Parameters:
107684f0b6dfSMatthew G. Knepley + label - the DMLabel
107784f0b6dfSMatthew G. Knepley - value - the stratum value
107884f0b6dfSMatthew G. Knepley 
107984f0b6dfSMatthew G. Knepley   Output Paramater:
108084f0b6dfSMatthew G. Knepley . size - The number of points in the stratum
108184f0b6dfSMatthew G. Knepley 
108284f0b6dfSMatthew G. Knepley   Level: intermediate
108384f0b6dfSMatthew G. Knepley 
108484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
108584f0b6dfSMatthew G. Knepley @*/
1086c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
1087c58f1c22SToby Isaac {
1088c58f1c22SToby Isaac   PetscInt       v;
1089c58f1c22SToby Isaac   PetscErrorCode ierr;
1090c58f1c22SToby Isaac 
1091c58f1c22SToby Isaac   PetscFunctionBegin;
1092d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1093c58f1c22SToby Isaac   PetscValidPointer(size, 3);
1094c58f1c22SToby Isaac   *size = 0;
10950c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
10960c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
1097c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
1098c58f1c22SToby Isaac   *size = label->stratumSizes[v];
1099c58f1c22SToby Isaac   PetscFunctionReturn(0);
1100c58f1c22SToby Isaac }
1101c58f1c22SToby Isaac 
110284f0b6dfSMatthew G. Knepley /*@
110384f0b6dfSMatthew G. Knepley   DMLabelGetStratumBounds - Get the largest and smallest point of a stratum
110484f0b6dfSMatthew G. Knepley 
11055b5e7992SMatthew G. Knepley   Not collective
11065b5e7992SMatthew G. Knepley 
110784f0b6dfSMatthew G. Knepley   Input Parameters:
110884f0b6dfSMatthew G. Knepley + label - the DMLabel
110984f0b6dfSMatthew G. Knepley - value - the stratum value
111084f0b6dfSMatthew G. Knepley 
111184f0b6dfSMatthew G. Knepley   Output Paramaters:
111284f0b6dfSMatthew G. Knepley + start - the smallest point in the stratum
111384f0b6dfSMatthew G. Knepley - end - the largest point in the stratum
111484f0b6dfSMatthew G. Knepley 
111584f0b6dfSMatthew G. Knepley   Level: intermediate
111684f0b6dfSMatthew G. Knepley 
111784f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
111884f0b6dfSMatthew G. Knepley @*/
1119c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
1120c58f1c22SToby Isaac {
11210c3c4a36SLisandro Dalcin   PetscInt       v, min, max;
1122c58f1c22SToby Isaac   PetscErrorCode ierr;
1123c58f1c22SToby Isaac 
1124c58f1c22SToby Isaac   PetscFunctionBegin;
1125d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1126c58f1c22SToby Isaac   if (start) {PetscValidPointer(start, 3); *start = 0;}
1127c58f1c22SToby Isaac   if (end)   {PetscValidPointer(end,   4); *end   = 0;}
11280c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
11290c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
1130c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
11310c3c4a36SLisandro Dalcin   if (label->stratumSizes[v] <= 0) PetscFunctionReturn(0);
1132d6cb179aSToby Isaac   ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr);
1133d6cb179aSToby Isaac   if (start) *start = min;
1134d6cb179aSToby Isaac   if (end)   *end   = max+1;
1135c58f1c22SToby Isaac   PetscFunctionReturn(0);
1136c58f1c22SToby Isaac }
1137c58f1c22SToby Isaac 
113884f0b6dfSMatthew G. Knepley /*@
113984f0b6dfSMatthew G. Knepley   DMLabelGetStratumIS - Get an IS with the stratum points
114084f0b6dfSMatthew G. Knepley 
11415b5e7992SMatthew G. Knepley   Not collective
11425b5e7992SMatthew G. Knepley 
114384f0b6dfSMatthew G. Knepley   Input Parameters:
114484f0b6dfSMatthew G. Knepley + label - the DMLabel
114584f0b6dfSMatthew G. Knepley - value - the stratum value
114684f0b6dfSMatthew G. Knepley 
114784f0b6dfSMatthew G. Knepley   Output Paramater:
114884f0b6dfSMatthew G. Knepley . points - The stratum points
114984f0b6dfSMatthew G. Knepley 
115084f0b6dfSMatthew G. Knepley   Level: intermediate
115184f0b6dfSMatthew G. Knepley 
115284f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
115384f0b6dfSMatthew G. Knepley @*/
1154c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
1155c58f1c22SToby Isaac {
1156c58f1c22SToby Isaac   PetscInt       v;
1157c58f1c22SToby Isaac   PetscErrorCode ierr;
1158c58f1c22SToby Isaac 
1159c58f1c22SToby Isaac   PetscFunctionBegin;
1160d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1161c58f1c22SToby Isaac   PetscValidPointer(points, 3);
1162c58f1c22SToby Isaac   *points = NULL;
11630c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
11640c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
1165c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
1166ad8374ffSToby Isaac   ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr);
1167ad8374ffSToby Isaac   *points = label->points[v];
1168c58f1c22SToby Isaac   PetscFunctionReturn(0);
1169c58f1c22SToby Isaac }
1170c58f1c22SToby Isaac 
117184f0b6dfSMatthew G. Knepley /*@
117284f0b6dfSMatthew G. Knepley   DMLabelSetStratumIS - Set the stratum points using an IS
117384f0b6dfSMatthew G. Knepley 
11745b5e7992SMatthew G. Knepley   Not collective
11755b5e7992SMatthew G. Knepley 
117684f0b6dfSMatthew G. Knepley   Input Parameters:
117784f0b6dfSMatthew G. Knepley + label - the DMLabel
117884f0b6dfSMatthew G. Knepley . value - the stratum value
117984f0b6dfSMatthew G. Knepley - points - The stratum points
118084f0b6dfSMatthew G. Knepley 
118184f0b6dfSMatthew G. Knepley   Level: intermediate
118284f0b6dfSMatthew G. Knepley 
118384f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
118484f0b6dfSMatthew G. Knepley @*/
11854de306b1SToby Isaac PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
11864de306b1SToby Isaac {
11870c3c4a36SLisandro Dalcin   PetscInt       v;
11884de306b1SToby Isaac   PetscErrorCode ierr;
11894de306b1SToby Isaac 
11904de306b1SToby Isaac   PetscFunctionBegin;
1191d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1192d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
1193b9907514SLisandro Dalcin   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
11944de306b1SToby Isaac   if (is == label->points[v]) PetscFunctionReturn(0);
11954de306b1SToby Isaac   ierr = DMLabelClearStratum(label, value);CHKERRQ(ierr);
11964de306b1SToby Isaac   ierr = ISGetLocalSize(is, &(label->stratumSizes[v]));CHKERRQ(ierr);
11974de306b1SToby Isaac   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
11984de306b1SToby Isaac   ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
11990c3c4a36SLisandro Dalcin   label->points[v]  = is;
12000c3c4a36SLisandro Dalcin   label->validIS[v] = PETSC_TRUE;
1201277ea44aSLisandro Dalcin   ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
12024de306b1SToby Isaac   if (label->bt) {
12034de306b1SToby Isaac     const PetscInt *points;
12044de306b1SToby Isaac     PetscInt p;
12054de306b1SToby Isaac 
12064de306b1SToby Isaac     ierr = ISGetIndices(is,&points);CHKERRQ(ierr);
12074de306b1SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
12084de306b1SToby Isaac       const PetscInt point = points[p];
12094de306b1SToby Isaac 
12104de306b1SToby 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);
12114de306b1SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
12124de306b1SToby Isaac     }
12134de306b1SToby Isaac   }
12144de306b1SToby Isaac   PetscFunctionReturn(0);
12154de306b1SToby Isaac }
12164de306b1SToby Isaac 
121784f0b6dfSMatthew G. Knepley /*@
121884f0b6dfSMatthew G. Knepley   DMLabelClearStratum - Remove a stratum
12194de306b1SToby Isaac 
12205b5e7992SMatthew G. Knepley   Not collective
12215b5e7992SMatthew G. Knepley 
122284f0b6dfSMatthew G. Knepley   Input Parameters:
122384f0b6dfSMatthew G. Knepley + label - the DMLabel
122484f0b6dfSMatthew G. Knepley - value - the stratum value
122584f0b6dfSMatthew G. Knepley 
122684f0b6dfSMatthew G. Knepley   Level: intermediate
122784f0b6dfSMatthew G. Knepley 
122884f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
122984f0b6dfSMatthew G. Knepley @*/
1230c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
1231c58f1c22SToby Isaac {
1232c58f1c22SToby Isaac   PetscInt       v;
1233c58f1c22SToby Isaac   PetscErrorCode ierr;
1234c58f1c22SToby Isaac 
1235c58f1c22SToby Isaac   PetscFunctionBegin;
1236d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
12370c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
12380c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
12394de306b1SToby Isaac   if (label->validIS[v]) {
12404de306b1SToby Isaac     if (label->bt) {
1241c58f1c22SToby Isaac       PetscInt       i;
1242ad8374ffSToby Isaac       const PetscInt *points;
1243c58f1c22SToby Isaac 
1244ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
1245c58f1c22SToby Isaac       for (i = 0; i < label->stratumSizes[v]; ++i) {
1246ad8374ffSToby Isaac         const PetscInt point = points[i];
1247c58f1c22SToby Isaac 
1248c58f1c22SToby 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);
1249c58f1c22SToby Isaac         ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
1250c58f1c22SToby Isaac       }
1251ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
1252c58f1c22SToby Isaac     }
1253c58f1c22SToby Isaac     label->stratumSizes[v] = 0;
12540c3c4a36SLisandro Dalcin     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
1255277ea44aSLisandro Dalcin     ierr = ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &label->points[v]);CHKERRQ(ierr);
12560c3c4a36SLisandro Dalcin     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
1257277ea44aSLisandro Dalcin     ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
1258c58f1c22SToby Isaac   } else {
1259b9907514SLisandro Dalcin     ierr = PetscHSetIClear(label->ht[v]);CHKERRQ(ierr);
1260c58f1c22SToby Isaac   }
1261c58f1c22SToby Isaac   PetscFunctionReturn(0);
1262c58f1c22SToby Isaac }
1263c58f1c22SToby Isaac 
126484f0b6dfSMatthew G. Knepley /*@
126584f0b6dfSMatthew G. Knepley   DMLabelFilter - Remove all points outside of [start, end)
126684f0b6dfSMatthew G. Knepley 
12675b5e7992SMatthew G. Knepley   Not collective
12685b5e7992SMatthew G. Knepley 
126984f0b6dfSMatthew G. Knepley   Input Parameters:
127084f0b6dfSMatthew G. Knepley + label - the DMLabel
127122cd2cdaSVaclav Hapla . start - the first point kept
127222cd2cdaSVaclav Hapla - end - one more than the last point kept
127384f0b6dfSMatthew G. Knepley 
127484f0b6dfSMatthew G. Knepley   Level: intermediate
127584f0b6dfSMatthew G. Knepley 
127684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
127784f0b6dfSMatthew G. Knepley @*/
1278c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
1279c58f1c22SToby Isaac {
1280c58f1c22SToby Isaac   PetscInt       v;
1281c58f1c22SToby Isaac   PetscErrorCode ierr;
1282c58f1c22SToby Isaac 
1283c58f1c22SToby Isaac   PetscFunctionBegin;
1284d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
12850c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
1286c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1287c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
12889106b633SVaclav Hapla     ierr = ISGeneralFilter(label->points[v], start, end);CHKERRQ(ierr);
1289c58f1c22SToby Isaac   }
1290c58f1c22SToby Isaac   ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr);
1291c58f1c22SToby Isaac   PetscFunctionReturn(0);
1292c58f1c22SToby Isaac }
1293c58f1c22SToby Isaac 
129484f0b6dfSMatthew G. Knepley /*@
129584f0b6dfSMatthew G. Knepley   DMLabelPermute - Create a new label with permuted points
129684f0b6dfSMatthew G. Knepley 
12975b5e7992SMatthew G. Knepley   Not collective
12985b5e7992SMatthew G. Knepley 
129984f0b6dfSMatthew G. Knepley   Input Parameters:
130084f0b6dfSMatthew G. Knepley + label - the DMLabel
130184f0b6dfSMatthew G. Knepley - permutation - the point permutation
130284f0b6dfSMatthew G. Knepley 
130384f0b6dfSMatthew G. Knepley   Output Parameter:
130484f0b6dfSMatthew G. Knepley . labelnew - the new label containing the permuted points
130584f0b6dfSMatthew G. Knepley 
130684f0b6dfSMatthew G. Knepley   Level: intermediate
130784f0b6dfSMatthew G. Knepley 
130884f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
130984f0b6dfSMatthew G. Knepley @*/
1310c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
1311c58f1c22SToby Isaac {
1312c58f1c22SToby Isaac   const PetscInt *perm;
1313c58f1c22SToby Isaac   PetscInt        numValues, numPoints, v, q;
1314c58f1c22SToby Isaac   PetscErrorCode  ierr;
1315c58f1c22SToby Isaac 
1316c58f1c22SToby Isaac   PetscFunctionBegin;
1317d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1318d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(permutation, IS_CLASSID, 2);
1319c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1320c58f1c22SToby Isaac   ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr);
1321c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr);
1322c58f1c22SToby Isaac   ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr);
1323c58f1c22SToby Isaac   ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr);
1324c58f1c22SToby Isaac   for (v = 0; v < numValues; ++v) {
1325c58f1c22SToby Isaac     const PetscInt size   = (*labelNew)->stratumSizes[v];
1326ad8374ffSToby Isaac     const PetscInt *points;
1327ad8374ffSToby Isaac     PetscInt *pointsNew;
1328c58f1c22SToby Isaac 
1329ad8374ffSToby Isaac     ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1330ad8374ffSToby Isaac     ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr);
1331c58f1c22SToby Isaac     for (q = 0; q < size; ++q) {
1332ad8374ffSToby Isaac       const PetscInt point = points[q];
1333c58f1c22SToby Isaac 
1334c58f1c22SToby 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);
1335ad8374ffSToby Isaac       pointsNew[q] = perm[point];
1336c58f1c22SToby Isaac     }
1337ad8374ffSToby Isaac     ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1338ad8374ffSToby Isaac     ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr);
1339ad8374ffSToby Isaac     ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr);
1340fa8e8ae5SToby Isaac     if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) {
1341fa8e8ae5SToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));CHKERRQ(ierr);
1342fa8e8ae5SToby Isaac       ierr = PetscFree(pointsNew);CHKERRQ(ierr);
1343fa8e8ae5SToby Isaac     } else {
1344ad8374ffSToby Isaac       ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr);
1345fa8e8ae5SToby Isaac     }
1346ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr);
1347c58f1c22SToby Isaac   }
1348c58f1c22SToby Isaac   ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr);
1349c58f1c22SToby Isaac   if (label->bt) {
1350c58f1c22SToby Isaac     ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
1351c58f1c22SToby Isaac     ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
1352c58f1c22SToby Isaac   }
1353c58f1c22SToby Isaac   PetscFunctionReturn(0);
1354c58f1c22SToby Isaac }
1355c58f1c22SToby Isaac 
135626c55118SMichael Lange PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
135726c55118SMichael Lange {
135826c55118SMichael Lange   MPI_Comm       comm;
135926c55118SMichael Lange   PetscInt       s, l, nroots, nleaves, dof, offset, size;
136026c55118SMichael Lange   PetscInt      *remoteOffsets, *rootStrata, *rootIdx;
136126c55118SMichael Lange   PetscSection   rootSection;
136226c55118SMichael Lange   PetscSF        labelSF;
136326c55118SMichael Lange   PetscErrorCode ierr;
136426c55118SMichael Lange 
136526c55118SMichael Lange   PetscFunctionBegin;
136626c55118SMichael Lange   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
136726c55118SMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
136826c55118SMichael Lange   /* Build a section of stratum values per point, generate the according SF
136926c55118SMichael Lange      and distribute point-wise stratum values to leaves. */
137026c55118SMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
137126c55118SMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
137226c55118SMichael Lange   ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr);
137326c55118SMichael Lange   if (label) {
137426c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1375ad8374ffSToby Isaac       const PetscInt *points;
1376ad8374ffSToby Isaac 
1377ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
137826c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1379ad8374ffSToby Isaac         ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr);
1380ad8374ffSToby Isaac         ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr);
138126c55118SMichael Lange       }
1382ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
138326c55118SMichael Lange     }
138426c55118SMichael Lange   }
138526c55118SMichael Lange   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
138626c55118SMichael Lange   /* Create a point-wise array of stratum values */
138726c55118SMichael Lange   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
138826c55118SMichael Lange   ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr);
138926c55118SMichael Lange   ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr);
139026c55118SMichael Lange   if (label) {
139126c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1392ad8374ffSToby Isaac       const PetscInt *points;
1393ad8374ffSToby Isaac 
1394ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
139526c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1396ad8374ffSToby Isaac         const PetscInt p = points[l];
139726c55118SMichael Lange         ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr);
139826c55118SMichael Lange         rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
139926c55118SMichael Lange       }
1400ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
140126c55118SMichael Lange     }
140226c55118SMichael Lange   }
140326c55118SMichael Lange   /* Build SF that maps label points to remote processes */
140426c55118SMichael Lange   ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr);
140526c55118SMichael Lange   ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr);
140626c55118SMichael Lange   ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr);
140726c55118SMichael Lange   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
140826c55118SMichael Lange   /* Send the strata for each point over the derived SF */
140926c55118SMichael Lange   ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr);
141026c55118SMichael Lange   ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr);
141126c55118SMichael Lange   ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
141226c55118SMichael Lange   ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
141326c55118SMichael Lange   /* Clean up */
141426c55118SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
141526c55118SMichael Lange   ierr = PetscFree(rootIdx);CHKERRQ(ierr);
141626c55118SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
141726c55118SMichael Lange   ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr);
141826c55118SMichael Lange   PetscFunctionReturn(0);
141926c55118SMichael Lange }
142026c55118SMichael Lange 
142184f0b6dfSMatthew G. Knepley /*@
142284f0b6dfSMatthew G. Knepley   DMLabelDistribute - Create a new label pushed forward over the PetscSF
142384f0b6dfSMatthew G. Knepley 
14245b5e7992SMatthew G. Knepley   Collective on sf
14255b5e7992SMatthew G. Knepley 
142684f0b6dfSMatthew G. Knepley   Input Parameters:
142784f0b6dfSMatthew G. Knepley + label - the DMLabel
142884f0b6dfSMatthew G. Knepley - sf    - the map from old to new distribution
142984f0b6dfSMatthew G. Knepley 
143084f0b6dfSMatthew G. Knepley   Output Parameter:
143184f0b6dfSMatthew G. Knepley . labelnew - the new redistributed label
143284f0b6dfSMatthew G. Knepley 
143384f0b6dfSMatthew G. Knepley   Level: intermediate
143484f0b6dfSMatthew G. Knepley 
143584f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
143684f0b6dfSMatthew G. Knepley @*/
1437c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1438c58f1c22SToby Isaac {
1439c58f1c22SToby Isaac   MPI_Comm       comm;
144026c55118SMichael Lange   PetscSection   leafSection;
144126c55118SMichael Lange   PetscInt       p, pStart, pEnd, s, size, dof, offset, stratum;
144226c55118SMichael Lange   PetscInt      *leafStrata, *strataIdx;
1443ad8374ffSToby Isaac   PetscInt     **points;
1444d67d17b1SMatthew G. Knepley   const char    *lname = NULL;
1445c58f1c22SToby Isaac   char          *name;
1446c58f1c22SToby Isaac   PetscInt       nameSize;
1447e8f14785SLisandro Dalcin   PetscHSetI     stratumHash;
1448c58f1c22SToby Isaac   size_t         len = 0;
144926c55118SMichael Lange   PetscMPIInt    rank;
1450c58f1c22SToby Isaac   PetscErrorCode ierr;
1451c58f1c22SToby Isaac 
1452c58f1c22SToby Isaac   PetscFunctionBegin;
1453d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1454f018e600SMatthew G. Knepley   if (label) {
1455f018e600SMatthew G. Knepley     PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1456f018e600SMatthew G. Knepley     ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1457f018e600SMatthew G. Knepley   }
1458c58f1c22SToby Isaac   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1459c58f1c22SToby Isaac   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1460c58f1c22SToby Isaac   /* Bcast name */
1461d67d17b1SMatthew G. Knepley   if (!rank) {
1462d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1463d67d17b1SMatthew G. Knepley     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1464d67d17b1SMatthew G. Knepley   }
1465c58f1c22SToby Isaac   nameSize = len;
1466c58f1c22SToby Isaac   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1467c58f1c22SToby Isaac   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1468580bdb30SBarry Smith   if (!rank) {ierr = PetscArraycpy(name, lname, nameSize+1);CHKERRQ(ierr);}
1469c58f1c22SToby Isaac   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1470d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
1471c58f1c22SToby Isaac   ierr = PetscFree(name);CHKERRQ(ierr);
147277d236dfSMichael Lange   /* Bcast defaultValue */
147377d236dfSMichael Lange   if (!rank) (*labelNew)->defaultValue = label->defaultValue;
147477d236dfSMichael Lange   ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
147526c55118SMichael Lange   /* Distribute stratum values over the SF and get the point mapping on the receiver */
147626c55118SMichael Lange   ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr);
14775cbdf6fcSMichael Lange   /* Determine received stratum values and initialise new label*/
1478e8f14785SLisandro Dalcin   ierr = PetscHSetICreate(&stratumHash);CHKERRQ(ierr);
147926c55118SMichael Lange   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
1480e8f14785SLisandro Dalcin   for (p = 0; p < size; ++p) {ierr = PetscHSetIAdd(stratumHash, leafStrata[p]);CHKERRQ(ierr);}
1481e8f14785SLisandro Dalcin   ierr = PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata);CHKERRQ(ierr);
1482ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr);
1483ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
14845cbdf6fcSMichael Lange   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr);
14855cbdf6fcSMichael Lange   /* Turn leafStrata into indices rather than stratum values */
14865cbdf6fcSMichael Lange   offset = 0;
1487e8f14785SLisandro Dalcin   ierr = PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr);
1488a205f1fdSToby Isaac   ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr);
148990e9b2aeSLisandro Dalcin   for (s = 0; s < (*labelNew)->numStrata; ++s) {
149090e9b2aeSLisandro Dalcin     ierr = PetscHMapISet((*labelNew)->hmap, (*labelNew)->stratumValues[s], s);CHKERRQ(ierr);
149190e9b2aeSLisandro Dalcin   }
14925cbdf6fcSMichael Lange   for (p = 0; p < size; ++p) {
1493231b9e6fSMatthew G. Knepley     for (s = 0; s < (*labelNew)->numStrata; ++s) {
1494231b9e6fSMatthew G. Knepley       if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
14955cbdf6fcSMichael Lange     }
14965cbdf6fcSMichael Lange   }
1497c58f1c22SToby Isaac   /* Rebuild the point strata on the receiver */
1498c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr);
1499c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1500c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1501c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1502c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1503c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1504c58f1c22SToby Isaac       (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1505c58f1c22SToby Isaac     }
1506c58f1c22SToby Isaac   }
1507c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr);
1508c58f1c22SToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr);
1509ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr);
1510c58f1c22SToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1511e8f14785SLisandro Dalcin     ierr = PetscHSetICreate(&(*labelNew)->ht[s]);CHKERRQ(ierr);
1512ad8374ffSToby Isaac     ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr);
1513c58f1c22SToby Isaac   }
1514c58f1c22SToby Isaac   /* Insert points into new strata */
1515c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr);
1516c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1517c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1518c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1519c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1520c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1521c58f1c22SToby Isaac       stratum = leafStrata[offset+s];
1522ad8374ffSToby Isaac       points[stratum][strataIdx[stratum]++] = p;
1523c58f1c22SToby Isaac     }
1524c58f1c22SToby Isaac   }
1525ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; s++) {
1526ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr);
1527ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr);
1528ad8374ffSToby Isaac   }
1529ad8374ffSToby Isaac   ierr = PetscFree(points);CHKERRQ(ierr);
1530e8f14785SLisandro Dalcin   ierr = PetscHSetIDestroy(&stratumHash);CHKERRQ(ierr);
1531c58f1c22SToby Isaac   ierr = PetscFree(leafStrata);CHKERRQ(ierr);
1532c58f1c22SToby Isaac   ierr = PetscFree(strataIdx);CHKERRQ(ierr);
1533c58f1c22SToby Isaac   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
1534c58f1c22SToby Isaac   PetscFunctionReturn(0);
1535c58f1c22SToby Isaac }
1536c58f1c22SToby Isaac 
15377937d9ceSMichael Lange /*@
15387937d9ceSMichael Lange   DMLabelGather - Gather all label values from leafs into roots
15397937d9ceSMichael Lange 
15405b5e7992SMatthew G. Knepley   Collective on sf
15415b5e7992SMatthew G. Knepley 
15427937d9ceSMichael Lange   Input Parameters:
15437937d9ceSMichael Lange + label - the DMLabel
154484f0b6dfSMatthew G. Knepley - sf - the Star Forest point communication map
15457937d9ceSMichael Lange 
154684f0b6dfSMatthew G. Knepley   Output Parameters:
154784f0b6dfSMatthew G. Knepley . labelNew - the new DMLabel with localised leaf values
15487937d9ceSMichael Lange 
15497937d9ceSMichael Lange   Level: developer
15507937d9ceSMichael Lange 
15517937d9ceSMichael Lange   Note: This is the inverse operation to DMLabelDistribute.
15527937d9ceSMichael Lange 
15537937d9ceSMichael Lange .seealso: DMLabelDistribute()
15547937d9ceSMichael Lange @*/
15557937d9ceSMichael Lange PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
15567937d9ceSMichael Lange {
15577937d9ceSMichael Lange   MPI_Comm       comm;
15587937d9ceSMichael Lange   PetscSection   rootSection;
15597937d9ceSMichael Lange   PetscSF        sfLabel;
15607937d9ceSMichael Lange   PetscSFNode   *rootPoints, *leafPoints;
15617937d9ceSMichael Lange   PetscInt       p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
15627937d9ceSMichael Lange   const PetscInt *rootDegree, *ilocal;
15637937d9ceSMichael Lange   PetscInt       *rootStrata;
1564d67d17b1SMatthew G. Knepley   const char    *lname;
15657937d9ceSMichael Lange   char          *name;
15667937d9ceSMichael Lange   PetscInt       nameSize;
15677937d9ceSMichael Lange   size_t         len = 0;
15689852e123SBarry Smith   PetscMPIInt    rank, size;
15697937d9ceSMichael Lange   PetscErrorCode ierr;
15707937d9ceSMichael Lange 
15717937d9ceSMichael Lange   PetscFunctionBegin;
1572d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1573d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
15747937d9ceSMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
15757937d9ceSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
15769852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
15777937d9ceSMichael Lange   /* Bcast name */
1578d67d17b1SMatthew G. Knepley   if (!rank) {
1579d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1580d67d17b1SMatthew G. Knepley     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1581d67d17b1SMatthew G. Knepley   }
15827937d9ceSMichael Lange   nameSize = len;
15837937d9ceSMichael Lange   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
15847937d9ceSMichael Lange   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1585580bdb30SBarry Smith   if (!rank) {ierr = PetscArraycpy(name, lname, nameSize+1);CHKERRQ(ierr);}
15867937d9ceSMichael Lange   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1587d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
15887937d9ceSMichael Lange   ierr = PetscFree(name);CHKERRQ(ierr);
15897937d9ceSMichael Lange   /* Gather rank/index pairs of leaves into local roots to build
15907937d9ceSMichael Lange      an inverse, multi-rooted SF. Note that this ignores local leaf
15917937d9ceSMichael Lange      indexing due to the use of the multiSF in PetscSFGather. */
15927937d9ceSMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr);
1593dc53bc9bSMatthew G. Knepley   ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr);
1594dc53bc9bSMatthew G. Knepley   for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
15957937d9ceSMichael Lange   for (p = 0; p < nleaves; p++) {
15968212dd46SStefano Zampini     PetscInt ilp = ilocal ? ilocal[p] : p;
15978212dd46SStefano Zampini 
15988212dd46SStefano Zampini     leafPoints[ilp].index = ilp;
15998212dd46SStefano Zampini     leafPoints[ilp].rank  = rank;
16007937d9ceSMichael Lange   }
16017937d9ceSMichael Lange   ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr);
16027937d9ceSMichael Lange   ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr);
16037937d9ceSMichael Lange   for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
16047937d9ceSMichael Lange   ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr);
16057937d9ceSMichael Lange   ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
16067937d9ceSMichael Lange   ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
16077937d9ceSMichael Lange   ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr);
16087937d9ceSMichael Lange   ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
16097937d9ceSMichael Lange   /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
16107937d9ceSMichael Lange   ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr);
16117937d9ceSMichael Lange   /* Rebuild the point strata on the receiver */
16127937d9ceSMichael Lange   for (p = 0, idx = 0; p < nroots; p++) {
16137937d9ceSMichael Lange     for (d = 0; d < rootDegree[p]; d++) {
16147937d9ceSMichael Lange       ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr);
16157937d9ceSMichael Lange       ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr);
16167937d9ceSMichael Lange       for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);}
16177937d9ceSMichael Lange     }
16187937d9ceSMichael Lange     idx += rootDegree[p];
16197937d9ceSMichael Lange   }
162077e0c0e7SMichael Lange   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
162177e0c0e7SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
162277e0c0e7SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
162377e0c0e7SMichael Lange   ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr);
16247937d9ceSMichael Lange   PetscFunctionReturn(0);
16257937d9ceSMichael Lange }
16267937d9ceSMichael Lange 
162784f0b6dfSMatthew G. Knepley /*@
162884f0b6dfSMatthew G. Knepley   DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label
162984f0b6dfSMatthew G. Knepley 
16305b5e7992SMatthew G. Knepley   Not collective
16315b5e7992SMatthew G. Knepley 
163284f0b6dfSMatthew G. Knepley   Input Parameter:
163384f0b6dfSMatthew G. Knepley . label - the DMLabel
163484f0b6dfSMatthew G. Knepley 
163584f0b6dfSMatthew G. Knepley   Output Parameters:
163684f0b6dfSMatthew G. Knepley + section - the section giving offsets for each stratum
163784f0b6dfSMatthew G. Knepley - is - An IS containing all the label points
163884f0b6dfSMatthew G. Knepley 
163984f0b6dfSMatthew G. Knepley   Level: developer
164084f0b6dfSMatthew G. Knepley 
164184f0b6dfSMatthew G. Knepley .seealso: DMLabelDistribute()
164284f0b6dfSMatthew G. Knepley @*/
1643c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1644c58f1c22SToby Isaac {
1645c58f1c22SToby Isaac   IS              vIS;
1646c58f1c22SToby Isaac   const PetscInt *values;
1647c58f1c22SToby Isaac   PetscInt       *points;
1648c58f1c22SToby Isaac   PetscInt        nV, vS = 0, vE = 0, v, N;
1649c58f1c22SToby Isaac   PetscErrorCode  ierr;
1650c58f1c22SToby Isaac 
1651c58f1c22SToby Isaac   PetscFunctionBegin;
1652d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1653c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr);
1654c58f1c22SToby Isaac   ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
1655c58f1c22SToby Isaac   ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr);
1656c58f1c22SToby Isaac   if (nV) {vS = values[0]; vE = values[0]+1;}
1657c58f1c22SToby Isaac   for (v = 1; v < nV; ++v) {
1658c58f1c22SToby Isaac     vS = PetscMin(vS, values[v]);
1659c58f1c22SToby Isaac     vE = PetscMax(vE, values[v]+1);
1660c58f1c22SToby Isaac   }
1661c58f1c22SToby Isaac   ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr);
1662c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr);
1663c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1664c58f1c22SToby Isaac     PetscInt n;
1665c58f1c22SToby Isaac 
1666c58f1c22SToby Isaac     ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr);
1667c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr);
1668c58f1c22SToby Isaac   }
1669c58f1c22SToby Isaac   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
1670c58f1c22SToby Isaac   ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr);
1671c58f1c22SToby Isaac   ierr = PetscMalloc1(N, &points);CHKERRQ(ierr);
1672c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1673c58f1c22SToby Isaac     IS              is;
1674c58f1c22SToby Isaac     const PetscInt *spoints;
1675c58f1c22SToby Isaac     PetscInt        dof, off, p;
1676c58f1c22SToby Isaac 
1677c58f1c22SToby Isaac     ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr);
1678c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr);
1679c58f1c22SToby Isaac     ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr);
1680c58f1c22SToby Isaac     ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr);
1681c58f1c22SToby Isaac     for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1682c58f1c22SToby Isaac     ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr);
1683c58f1c22SToby Isaac     ierr = ISDestroy(&is);CHKERRQ(ierr);
1684c58f1c22SToby Isaac   }
1685c58f1c22SToby Isaac   ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr);
1686c58f1c22SToby Isaac   ierr = ISDestroy(&vIS);CHKERRQ(ierr);
1687c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1688c58f1c22SToby Isaac   PetscFunctionReturn(0);
1689c58f1c22SToby Isaac }
1690c58f1c22SToby Isaac 
169184f0b6dfSMatthew G. Knepley /*@
1692c58f1c22SToby Isaac   PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1693c58f1c22SToby Isaac   the local section and an SF describing the section point overlap.
1694c58f1c22SToby Isaac 
16955b5e7992SMatthew G. Knepley   Collective on sf
16965b5e7992SMatthew G. Knepley 
1697c58f1c22SToby Isaac   Input Parameters:
1698c58f1c22SToby Isaac   + s - The PetscSection for the local field layout
1699c58f1c22SToby Isaac   . sf - The SF describing parallel layout of the section points
1700c58f1c22SToby Isaac   . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1701c58f1c22SToby Isaac   . label - The label specifying the points
1702c58f1c22SToby Isaac   - labelValue - The label stratum specifying the points
1703c58f1c22SToby Isaac 
1704c58f1c22SToby Isaac   Output Parameter:
1705c58f1c22SToby Isaac   . gsection - The PetscSection for the global field layout
1706c58f1c22SToby Isaac 
1707c58f1c22SToby Isaac   Note: This gives negative sizes and offsets to points not owned by this process
1708c58f1c22SToby Isaac 
1709c58f1c22SToby Isaac   Level: developer
1710c58f1c22SToby Isaac 
1711c58f1c22SToby Isaac .seealso: PetscSectionCreate()
1712c58f1c22SToby Isaac @*/
1713c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1714c58f1c22SToby Isaac {
1715c58f1c22SToby Isaac   PetscInt      *neg = NULL, *tmpOff = NULL;
1716c58f1c22SToby Isaac   PetscInt       pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1717c58f1c22SToby Isaac   PetscErrorCode ierr;
1718c58f1c22SToby Isaac 
1719c58f1c22SToby Isaac   PetscFunctionBegin;
1720d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1721d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1722d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4);
1723c58f1c22SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr);
1724c58f1c22SToby Isaac   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1725c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr);
1726c58f1c22SToby Isaac   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1727c58f1c22SToby Isaac   if (nroots >= 0) {
1728c58f1c22SToby Isaac     if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1729c58f1c22SToby Isaac     ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr);
1730c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1731c58f1c22SToby Isaac       ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr);
1732c58f1c22SToby Isaac     } else {
1733c58f1c22SToby Isaac       tmpOff = &(*gsection)->atlasDof[-pStart];
1734c58f1c22SToby Isaac     }
1735c58f1c22SToby Isaac   }
1736c58f1c22SToby Isaac   /* Mark ghost points with negative dof */
1737c58f1c22SToby Isaac   for (p = pStart; p < pEnd; ++p) {
1738c58f1c22SToby Isaac     PetscInt value;
1739c58f1c22SToby Isaac 
1740c58f1c22SToby Isaac     ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr);
1741c58f1c22SToby Isaac     if (value != labelValue) continue;
1742c58f1c22SToby Isaac     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1743c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr);
1744c58f1c22SToby Isaac     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1745c58f1c22SToby Isaac     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);}
1746c58f1c22SToby Isaac     if (neg) neg[p] = -(dof+1);
1747c58f1c22SToby Isaac   }
1748c58f1c22SToby Isaac   ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr);
1749c58f1c22SToby Isaac   if (nroots >= 0) {
1750c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1751c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1752c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1753c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1754c58f1c22SToby Isaac     }
1755c58f1c22SToby Isaac   }
1756c58f1c22SToby Isaac   /* Calculate new sizes, get proccess offset, and calculate point offsets */
1757c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1758c58f1c22SToby Isaac     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1759c58f1c22SToby Isaac     (*gsection)->atlasOff[p] = off;
1760c58f1c22SToby Isaac     off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1761c58f1c22SToby Isaac   }
1762c58f1c22SToby Isaac   ierr       = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1763c58f1c22SToby Isaac   globalOff -= off;
1764c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1765c58f1c22SToby Isaac     (*gsection)->atlasOff[p] += globalOff;
1766c58f1c22SToby Isaac     if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1767c58f1c22SToby Isaac   }
1768c58f1c22SToby Isaac   /* Put in negative offsets for ghost points */
1769c58f1c22SToby Isaac   if (nroots >= 0) {
1770c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1771c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1772c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1773c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1774c58f1c22SToby Isaac     }
1775c58f1c22SToby Isaac   }
1776c58f1c22SToby Isaac   if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);}
1777c58f1c22SToby Isaac   ierr = PetscFree(neg);CHKERRQ(ierr);
1778c58f1c22SToby Isaac   PetscFunctionReturn(0);
1779c58f1c22SToby Isaac }
1780c58f1c22SToby Isaac 
17815fdea053SToby Isaac typedef struct _n_PetscSectionSym_Label
17825fdea053SToby Isaac {
17835fdea053SToby Isaac   DMLabel           label;
17845fdea053SToby Isaac   PetscCopyMode     *modes;
17855fdea053SToby Isaac   PetscInt          *sizes;
17865fdea053SToby Isaac   const PetscInt    ***perms;
17875fdea053SToby Isaac   const PetscScalar ***rots;
17885fdea053SToby Isaac   PetscInt          (*minMaxOrients)[2];
17895fdea053SToby Isaac   PetscInt          numStrata; /* numStrata is only increasing, functions as a state */
17905fdea053SToby Isaac } PetscSectionSym_Label;
17915fdea053SToby Isaac 
17925fdea053SToby Isaac static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
17935fdea053SToby Isaac {
17945fdea053SToby Isaac   PetscInt              i, j;
17955fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
17965fdea053SToby Isaac   PetscErrorCode        ierr;
17975fdea053SToby Isaac 
17985fdea053SToby Isaac   PetscFunctionBegin;
17995fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
18005fdea053SToby Isaac     if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
18015fdea053SToby Isaac       for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
18025fdea053SToby Isaac         if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);}
18035fdea053SToby Isaac         if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);}
18045fdea053SToby Isaac       }
18055fdea053SToby Isaac       if (sl->perms[i]) {
18065fdea053SToby Isaac         const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
18075fdea053SToby Isaac 
18085fdea053SToby Isaac         ierr = PetscFree(perms);CHKERRQ(ierr);
18095fdea053SToby Isaac       }
18105fdea053SToby Isaac       if (sl->rots[i]) {
18115fdea053SToby Isaac         const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
18125fdea053SToby Isaac 
18135fdea053SToby Isaac         ierr = PetscFree(rots);CHKERRQ(ierr);
18145fdea053SToby Isaac       }
18155fdea053SToby Isaac     }
18165fdea053SToby Isaac   }
18175fdea053SToby Isaac   ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr);
18185fdea053SToby Isaac   ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr);
18195fdea053SToby Isaac   sl->numStrata = 0;
18205fdea053SToby Isaac   PetscFunctionReturn(0);
18215fdea053SToby Isaac }
18225fdea053SToby Isaac 
18235fdea053SToby Isaac static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
18245fdea053SToby Isaac {
18255fdea053SToby Isaac   PetscErrorCode ierr;
18265fdea053SToby Isaac 
18275fdea053SToby Isaac   PetscFunctionBegin;
18285fdea053SToby Isaac   ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);
18295fdea053SToby Isaac   ierr = PetscFree(sym->data);CHKERRQ(ierr);
18305fdea053SToby Isaac   PetscFunctionReturn(0);
18315fdea053SToby Isaac }
18325fdea053SToby Isaac 
18335fdea053SToby Isaac static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
18345fdea053SToby Isaac {
18355fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
18365fdea053SToby Isaac   PetscBool             isAscii;
18375fdea053SToby Isaac   DMLabel               label = sl->label;
1838d67d17b1SMatthew G. Knepley   const char           *name;
18395fdea053SToby Isaac   PetscErrorCode        ierr;
18405fdea053SToby Isaac 
18415fdea053SToby Isaac   PetscFunctionBegin;
18425fdea053SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
18435fdea053SToby Isaac   if (isAscii) {
18445fdea053SToby Isaac     PetscInt          i, j, k;
18455fdea053SToby Isaac     PetscViewerFormat format;
18465fdea053SToby Isaac 
18475fdea053SToby Isaac     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
18485fdea053SToby Isaac     if (label) {
18495fdea053SToby Isaac       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
18505fdea053SToby Isaac       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
18515fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18525fdea053SToby Isaac         ierr = DMLabelView(label, viewer);CHKERRQ(ierr);
18535fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
18545fdea053SToby Isaac       } else {
1855d67d17b1SMatthew G. Knepley         ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1856d67d17b1SMatthew G. Knepley         ierr = PetscViewerASCIIPrintf(viewer,"  Label '%s'\n",name);CHKERRQ(ierr);
18575fdea053SToby Isaac       }
18585fdea053SToby Isaac     } else {
18595fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr);
18605fdea053SToby Isaac     }
18615fdea053SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18625fdea053SToby Isaac     for (i = 0; i <= sl->numStrata; i++) {
18635fdea053SToby Isaac       PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
18645fdea053SToby Isaac 
18655fdea053SToby Isaac       if (!(sl->perms[i] || sl->rots[i])) {
18665fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr);
18675fdea053SToby Isaac       } else {
18685fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr);
18695fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18705fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr);
18715fdea053SToby Isaac         if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
18725fdea053SToby Isaac           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18735fdea053SToby Isaac           for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
18745fdea053SToby Isaac             if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
18755fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr);
18765fdea053SToby Isaac             } else {
18775fdea053SToby Isaac               PetscInt tab;
18785fdea053SToby Isaac 
18795fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr);
18805fdea053SToby Isaac               ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18815fdea053SToby Isaac               ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr);
18825fdea053SToby Isaac               if (sl->perms[i] && sl->perms[i][j]) {
18835fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr);
18845fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
18855fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);}
18865fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
18875fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
18885fdea053SToby Isaac               }
18895fdea053SToby Isaac               if (sl->rots[i] && sl->rots[i][j]) {
18905fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Rotations:  ");CHKERRQ(ierr);
18915fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
18925fdea053SToby Isaac #if defined(PETSC_USE_COMPLEX)
18935fdea053SToby 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);}
18945fdea053SToby Isaac #else
18955fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);}
18965fdea053SToby Isaac #endif
18975fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
18985fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
18995fdea053SToby Isaac               }
19005fdea053SToby Isaac               ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
19015fdea053SToby Isaac             }
19025fdea053SToby Isaac           }
19035fdea053SToby Isaac           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
19045fdea053SToby Isaac         }
19055fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
19065fdea053SToby Isaac       }
19075fdea053SToby Isaac     }
19085fdea053SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
19095fdea053SToby Isaac   }
19105fdea053SToby Isaac   PetscFunctionReturn(0);
19115fdea053SToby Isaac }
19125fdea053SToby Isaac 
19135fdea053SToby Isaac /*@
19145fdea053SToby Isaac   PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
19155fdea053SToby Isaac 
19165fdea053SToby Isaac   Logically collective on sym
19175fdea053SToby Isaac 
19185fdea053SToby Isaac   Input parameters:
19195fdea053SToby Isaac + sym - the section symmetries
19205fdea053SToby Isaac - label - the DMLabel describing the types of points
19215fdea053SToby Isaac 
19225fdea053SToby Isaac   Level: developer:
19235fdea053SToby Isaac 
19245fdea053SToby Isaac .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()
19255fdea053SToby Isaac @*/
19265fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
19275fdea053SToby Isaac {
19285fdea053SToby Isaac   PetscSectionSym_Label *sl;
19295fdea053SToby Isaac   PetscErrorCode        ierr;
19305fdea053SToby Isaac 
19315fdea053SToby Isaac   PetscFunctionBegin;
19325fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
19335fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
19345fdea053SToby Isaac   if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);}
19355fdea053SToby Isaac   if (label) {
19365fdea053SToby Isaac     sl->label = label;
1937d67d17b1SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) label);CHKERRQ(ierr);
19385fdea053SToby Isaac     ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr);
19391a834cf9SToby 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);
19401a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr);
19411a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr);
19421a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr);
19431a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr);
19441a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr);
19455fdea053SToby Isaac   }
19465fdea053SToby Isaac   PetscFunctionReturn(0);
19475fdea053SToby Isaac }
19485fdea053SToby Isaac 
19495fdea053SToby Isaac /*@C
19505fdea053SToby Isaac   PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
19515fdea053SToby Isaac 
19525b5e7992SMatthew G. Knepley   Logically collective on sym
19535fdea053SToby Isaac 
19545fdea053SToby Isaac   InputParameters:
19555b5e7992SMatthew G. Knepley + sym       - the section symmetries
19565fdea053SToby Isaac . stratum   - the stratum value in the label that we are assigning symmetries for
19575fdea053SToby Isaac . size      - the number of dofs for points in the stratum of the label
19585fdea053SToby Isaac . minOrient - the smallest orientation for a point in this stratum
19595fdea053SToby Isaac . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
19605fdea053SToby Isaac . mode      - how sym should copy the perms and rots arrays
19615fdea053SToby Isaac . perms     - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation.  A NULL permutation is the identity
1962a2b725a8SWilliam Gropp - rots      - NULL if there are no rotations, or (maxOrient - minOrient) sets of rotations, one for each orientation.  A NULL set of orientations is the identity
19635fdea053SToby Isaac 
19645fdea053SToby Isaac   Level: developer
19655fdea053SToby Isaac 
19665fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()
19675fdea053SToby Isaac @*/
19685fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
19695fdea053SToby Isaac {
19705fdea053SToby Isaac   PetscSectionSym_Label *sl;
1971d67d17b1SMatthew G. Knepley   const char            *name;
1972d67d17b1SMatthew G. Knepley   PetscInt               i, j, k;
19735fdea053SToby Isaac   PetscErrorCode         ierr;
19745fdea053SToby Isaac 
19755fdea053SToby Isaac   PetscFunctionBegin;
19765fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
19775fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
19785fdea053SToby Isaac   if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet");
19795fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
19805fdea053SToby Isaac     PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
19815fdea053SToby Isaac 
19825fdea053SToby Isaac     if (stratum == value) break;
19835fdea053SToby Isaac   }
1984d67d17b1SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1985d67d17b1SMatthew G. Knepley   if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,name);
19865fdea053SToby Isaac   sl->sizes[i] = size;
19875fdea053SToby Isaac   sl->modes[i] = mode;
19885fdea053SToby Isaac   sl->minMaxOrients[i][0] = minOrient;
19895fdea053SToby Isaac   sl->minMaxOrients[i][1] = maxOrient;
19905fdea053SToby Isaac   if (mode == PETSC_COPY_VALUES) {
19915fdea053SToby Isaac     if (perms) {
19925fdea053SToby Isaac       PetscInt    **ownPerms;
19935fdea053SToby Isaac 
19945fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr);
19955fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
19965fdea053SToby Isaac         if (perms[j]) {
19975fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr);
19985fdea053SToby Isaac           for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];}
19995fdea053SToby Isaac         }
20005fdea053SToby Isaac       }
20015fdea053SToby Isaac       sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient];
20025fdea053SToby Isaac     }
20035fdea053SToby Isaac     if (rots) {
20045fdea053SToby Isaac       PetscScalar **ownRots;
20055fdea053SToby Isaac 
20065fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr);
20075fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
20085fdea053SToby Isaac         if (rots[j]) {
20095fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr);
20105fdea053SToby Isaac           for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];}
20115fdea053SToby Isaac         }
20125fdea053SToby Isaac       }
20135fdea053SToby Isaac       sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient];
20145fdea053SToby Isaac     }
20155fdea053SToby Isaac   } else {
20165fdea053SToby Isaac     sl->perms[i] = perms ? &perms[-minOrient] : NULL;
20175fdea053SToby Isaac     sl->rots[i]  = rots ? &rots[-minOrient] : NULL;
20185fdea053SToby Isaac   }
20195fdea053SToby Isaac   PetscFunctionReturn(0);
20205fdea053SToby Isaac }
20215fdea053SToby Isaac 
20225fdea053SToby Isaac static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
20235fdea053SToby Isaac {
20245fdea053SToby Isaac   PetscInt              i, j, numStrata;
20255fdea053SToby Isaac   PetscSectionSym_Label *sl;
20265fdea053SToby Isaac   DMLabel               label;
20275fdea053SToby Isaac   PetscErrorCode        ierr;
20285fdea053SToby Isaac 
20295fdea053SToby Isaac   PetscFunctionBegin;
20305fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
20315fdea053SToby Isaac   numStrata = sl->numStrata;
20325fdea053SToby Isaac   label     = sl->label;
20335fdea053SToby Isaac   for (i = 0; i < numPoints; i++) {
20345fdea053SToby Isaac     PetscInt point = points[2*i];
20355fdea053SToby Isaac     PetscInt ornt  = points[2*i+1];
20365fdea053SToby Isaac 
20375fdea053SToby Isaac     for (j = 0; j < numStrata; j++) {
20385fdea053SToby Isaac       if (label->validIS[j]) {
20395fdea053SToby Isaac         PetscInt k;
20405fdea053SToby Isaac 
20415fdea053SToby Isaac         ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr);
20425fdea053SToby Isaac         if (k >= 0) break;
20435fdea053SToby Isaac       } else {
20445fdea053SToby Isaac         PetscBool has;
20455fdea053SToby Isaac 
2046b9907514SLisandro Dalcin         ierr = PetscHSetIHas(label->ht[j], point, &has);CHKERRQ(ierr);
20475fdea053SToby Isaac         if (has) break;
20485fdea053SToby Isaac       }
20495fdea053SToby Isaac     }
20505fdea053SToby 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);
20515fdea053SToby Isaac     if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;}
20525fdea053SToby Isaac     if (rots) {rots[i]  = sl->rots[j] ? sl->rots[j][ornt] : NULL;}
20535fdea053SToby Isaac   }
20545fdea053SToby Isaac   PetscFunctionReturn(0);
20555fdea053SToby Isaac }
20565fdea053SToby Isaac 
20575fdea053SToby Isaac PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
20585fdea053SToby Isaac {
20595fdea053SToby Isaac   PetscSectionSym_Label *sl;
20605fdea053SToby Isaac   PetscErrorCode        ierr;
20615fdea053SToby Isaac 
20625fdea053SToby Isaac   PetscFunctionBegin;
20635fdea053SToby Isaac   ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr);
20645fdea053SToby Isaac   sym->ops->getpoints = PetscSectionSymGetPoints_Label;
20655fdea053SToby Isaac   sym->ops->view      = PetscSectionSymView_Label;
20665fdea053SToby Isaac   sym->ops->destroy   = PetscSectionSymDestroy_Label;
20675fdea053SToby Isaac   sym->data           = (void *) sl;
20685fdea053SToby Isaac   PetscFunctionReturn(0);
20695fdea053SToby Isaac }
20705fdea053SToby Isaac 
20715fdea053SToby Isaac /*@
20725fdea053SToby Isaac   PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
20735fdea053SToby Isaac 
20745fdea053SToby Isaac   Collective
20755fdea053SToby Isaac 
20765fdea053SToby Isaac   Input Parameters:
20775fdea053SToby Isaac + comm - the MPI communicator for the new symmetry
20785fdea053SToby Isaac - label - the label defining the strata
20795fdea053SToby Isaac 
20805fdea053SToby Isaac   Output Parameters:
20815fdea053SToby Isaac . sym - the section symmetries
20825fdea053SToby Isaac 
20835fdea053SToby Isaac   Level: developer
20845fdea053SToby Isaac 
20855fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()
20865fdea053SToby Isaac @*/
20875fdea053SToby Isaac PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
20885fdea053SToby Isaac {
20895fdea053SToby Isaac   PetscErrorCode        ierr;
20905fdea053SToby Isaac 
20915fdea053SToby Isaac   PetscFunctionBegin;
20925fdea053SToby Isaac   ierr = DMInitializePackage();CHKERRQ(ierr);
20935fdea053SToby Isaac   ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr);
20945fdea053SToby Isaac   ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr);
20955fdea053SToby Isaac   ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr);
20965fdea053SToby Isaac   PetscFunctionReturn(0);
20975fdea053SToby Isaac }
2098