xref: /petsc/src/dm/label/dmlabel.c (revision fc8a9adeb7fcdc98711d755fa2dc544ddccf0f3e)
1 #include <petscdm.h>
2 #include <petsc/private/dmlabelimpl.h>   /*I      "petscdmlabel.h"   I*/
3 #include <petsc/private/isimpl.h>        /*I      "petscis.h"        I*/
4 #include <petscsf.h>
5 
6 /*@C
7   DMLabelCreate - Create a DMLabel object, which is a multimap
8 
9   Input parameters:
10 + comm - The communicator, usually PETSC_COMM_SELF
11 - name - The label name
12 
13   Output parameter:
14 . label - The DMLabel
15 
16   Level: beginner
17 
18 .seealso: DMLabelDestroy()
19 @*/
20 PetscErrorCode DMLabelCreate(MPI_Comm comm, const char name[], DMLabel *label)
21 {
22   PetscErrorCode ierr;
23 
24   PetscFunctionBegin;
25   PetscValidPointer(label,2);
26   ierr = DMInitializePackage();CHKERRQ(ierr);
27 
28   ierr = PetscHeaderCreate(*label,DMLABEL_CLASSID,"DMLabel","DMLabel","DM",comm,DMLabelDestroy,DMLabelView);CHKERRQ(ierr);
29 
30   (*label)->numStrata      = 0;
31   (*label)->defaultValue   = -1;
32   (*label)->stratumValues  = NULL;
33   (*label)->validIS        = NULL;
34   (*label)->stratumSizes   = NULL;
35   (*label)->points         = NULL;
36   (*label)->ht             = NULL;
37   (*label)->pStart         = -1;
38   (*label)->pEnd           = -1;
39   (*label)->bt             = NULL;
40   ierr = PetscHMapICreate(&(*label)->hmap);CHKERRQ(ierr);
41   ierr = PetscObjectSetName((PetscObject) *label, name);CHKERRQ(ierr);
42   PetscFunctionReturn(0);
43 }
44 
45 /*
46   DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format
47 
48   Input parameter:
49 + label - The DMLabel
50 - v - The stratum value
51 
52   Output parameter:
53 . label - The DMLabel with stratum in sorted list format
54 
55   Level: developer
56 
57 .seealso: DMLabelCreate()
58 */
59 static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v)
60 {
61   IS             is;
62   PetscInt       off = 0, *pointArray, p;
63   PetscErrorCode ierr;
64 
65   if (PetscLikely(v >= 0 && v < label->numStrata) && label->validIS[v]) return 0;
66   PetscFunctionBegin;
67   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);
68   ierr = PetscHSetIGetSize(label->ht[v], &label->stratumSizes[v]);CHKERRQ(ierr);
69   ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr);
70   ierr = PetscHSetIGetElems(label->ht[v], &off, pointArray);CHKERRQ(ierr);
71   ierr = PetscHSetIClear(label->ht[v]);CHKERRQ(ierr);
72   ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr);
73   if (label->bt) {
74     for (p = 0; p < label->stratumSizes[v]; ++p) {
75       const PetscInt point = pointArray[p];
76       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);
77       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
78     }
79   }
80   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->stratumSizes[v], pointArray, PETSC_OWN_POINTER, &is);CHKERRQ(ierr);
81   ierr = PetscObjectSetName((PetscObject) is, "indices");CHKERRQ(ierr);
82   label->points[v]  = is;
83   label->validIS[v] = PETSC_TRUE;
84   ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
85   PetscFunctionReturn(0);
86 }
87 
88 /*
89   DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format
90 
91   Input parameter:
92 . label - The DMLabel
93 
94   Output parameter:
95 . label - The DMLabel with all strata in sorted list format
96 
97   Level: developer
98 
99 .seealso: DMLabelCreate()
100 */
101 static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label)
102 {
103   PetscInt       v;
104   PetscErrorCode ierr;
105 
106   PetscFunctionBegin;
107   for (v = 0; v < label->numStrata; v++){
108     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
109   }
110   PetscFunctionReturn(0);
111 }
112 
113 /*
114   DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format
115 
116   Input parameter:
117 + label - The DMLabel
118 - v - The stratum value
119 
120   Output parameter:
121 . label - The DMLabel with stratum in hash format
122 
123   Level: developer
124 
125 .seealso: DMLabelCreate()
126 */
127 static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v)
128 {
129   PetscInt       p;
130   const PetscInt *points;
131   PetscErrorCode ierr;
132 
133   if (PetscLikely(v >= 0 && v < label->numStrata) && !label->validIS[v]) return 0;
134   PetscFunctionBegin;
135   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);
136   if (label->points[v]) {
137     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
138     for (p = 0; p < label->stratumSizes[v]; ++p) {
139       ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);
140     }
141     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
142     ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
143   }
144   label->validIS[v] = PETSC_FALSE;
145   PetscFunctionReturn(0);
146 }
147 
148 #if !defined(DMLABEL_LOOKUP_THRESHOLD)
149 #define DMLABEL_LOOKUP_THRESHOLD 16
150 #endif
151 
152 PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index)
153 {
154   PetscInt       v;
155   PetscErrorCode ierr;
156 
157   PetscFunctionBegin;
158   *index = -1;
159   if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) {
160     for (v = 0; v < label->numStrata; ++v)
161       if (label->stratumValues[v] == value) {*index = v; break;}
162   } else {
163     ierr = PetscHMapIGet(label->hmap, value, index);CHKERRQ(ierr);
164   }
165 #if defined(PETSC_USE_DEBUG)
166   { /* Check strata hash map consistency */
167     PetscInt len, loc = -1;
168     ierr = PetscHMapIGetSize(label->hmap, &len);CHKERRQ(ierr);
169     if (len != label->numStrata) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map size");
170     if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) {
171       ierr = PetscHMapIGet(label->hmap, value, &loc);CHKERRQ(ierr);
172     } else {
173       for (v = 0; v < label->numStrata; ++v)
174         if (label->stratumValues[v] == value) {loc = v; break;}
175     }
176     if (loc != *index) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map lookup");
177   }
178 #endif
179   PetscFunctionReturn(0);
180 }
181 
182 PETSC_STATIC_INLINE PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index)
183 {
184   PetscInt       v;
185   PetscInt      *tmpV;
186   PetscInt      *tmpS;
187   PetscHSetI    *tmpH, ht;
188   IS            *tmpP, is;
189   PetscBool     *tmpB;
190   PetscHMapI     hmap = label->hmap;
191   PetscErrorCode ierr;
192 
193   PetscFunctionBegin;
194   v    = label->numStrata;
195   tmpV = label->stratumValues;
196   tmpS = label->stratumSizes;
197   tmpH = label->ht;
198   tmpP = label->points;
199   tmpB = label->validIS;
200   { /* TODO: PetscRealloc() is broken, use malloc+memcpy+free  */
201     PetscInt   *oldV = tmpV;
202     PetscInt   *oldS = tmpS;
203     PetscHSetI *oldH = tmpH;
204     IS         *oldP = tmpP;
205     PetscBool  *oldB = tmpB;
206     ierr = PetscMalloc((v+1)*sizeof(*tmpV), &tmpV);CHKERRQ(ierr);
207     ierr = PetscMalloc((v+1)*sizeof(*tmpS), &tmpS);CHKERRQ(ierr);
208     ierr = PetscMalloc((v+1)*sizeof(*tmpH), &tmpH);CHKERRQ(ierr);
209     ierr = PetscMalloc((v+1)*sizeof(*tmpP), &tmpP);CHKERRQ(ierr);
210     ierr = PetscMalloc((v+1)*sizeof(*tmpB), &tmpB);CHKERRQ(ierr);
211     ierr = PetscMemcpy(tmpV, oldV, v*sizeof(*tmpV));CHKERRQ(ierr);
212     ierr = PetscMemcpy(tmpS, oldS, v*sizeof(*tmpS));CHKERRQ(ierr);
213     ierr = PetscMemcpy(tmpH, oldH, v*sizeof(*tmpH));CHKERRQ(ierr);
214     ierr = PetscMemcpy(tmpP, oldP, v*sizeof(*tmpP));CHKERRQ(ierr);
215     ierr = PetscMemcpy(tmpB, oldB, v*sizeof(*tmpB));CHKERRQ(ierr);
216     ierr = PetscFree(oldV);CHKERRQ(ierr);
217     ierr = PetscFree(oldS);CHKERRQ(ierr);
218     ierr = PetscFree(oldH);CHKERRQ(ierr);
219     ierr = PetscFree(oldP);CHKERRQ(ierr);
220     ierr = PetscFree(oldB);CHKERRQ(ierr);
221   }
222   label->numStrata     = v+1;
223   label->stratumValues = tmpV;
224   label->stratumSizes  = tmpS;
225   label->ht            = tmpH;
226   label->points        = tmpP;
227   label->validIS       = tmpB;
228   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
229   ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);CHKERRQ(ierr);
230   ierr = PetscHMapISet(hmap, value, v);CHKERRQ(ierr);
231   tmpV[v] = value;
232   tmpS[v] = 0;
233   tmpH[v] = ht;
234   tmpP[v] = is;
235   tmpB[v] = PETSC_TRUE;
236   ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
237   *index = v;
238   PetscFunctionReturn(0);
239 }
240 
241 PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupAddStratum(DMLabel label, PetscInt value, PetscInt *index)
242 {
243   PetscErrorCode ierr;
244   PetscFunctionBegin;
245   ierr = DMLabelLookupStratum(label, value, index);CHKERRQ(ierr);
246   if (*index < 0) {ierr = DMLabelNewStratum(label, value, index);CHKERRQ(ierr);}
247   PetscFunctionReturn(0);
248 }
249 
250 /*@
251   DMLabelAddStratum - Adds a new stratum value in a DMLabel
252 
253   Input Parameter:
254 + label - The DMLabel
255 - value - The stratum value
256 
257   Level: beginner
258 
259 .seealso:  DMLabelCreate(), DMLabelDestroy()
260 @*/
261 PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
262 {
263   PetscInt       v;
264   PetscErrorCode ierr;
265 
266   PetscFunctionBegin;
267   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
268   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
269   PetscFunctionReturn(0);
270 }
271 
272 /*@
273   DMLabelAddStrata - Adds new stratum values in a DMLabel
274 
275   Input Parameter:
276 + label - The DMLabel
277 . numStrata - The number of stratum values
278 - stratumValues - The stratum values
279 
280   Level: beginner
281 
282 .seealso:  DMLabelCreate(), DMLabelDestroy()
283 @*/
284 PetscErrorCode DMLabelAddStrata(DMLabel label, PetscInt numStrata, const PetscInt stratumValues[])
285 {
286   PetscInt       *values, v;
287   PetscErrorCode ierr;
288 
289   PetscFunctionBegin;
290   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
291   if (numStrata) PetscValidIntPointer(stratumValues, 3);
292   ierr = PetscMalloc1(numStrata, &values);CHKERRQ(ierr);
293   ierr = PetscMemcpy(values, stratumValues, numStrata*sizeof(PetscInt));CHKERRQ(ierr);
294   ierr = PetscSortRemoveDupsInt(&numStrata, values);CHKERRQ(ierr);
295   if (!label->numStrata) { /* Fast preallocation */
296     PetscInt   *tmpV;
297     PetscInt   *tmpS;
298     PetscHSetI *tmpH, ht;
299     IS         *tmpP, is;
300     PetscBool  *tmpB;
301     PetscHMapI  hmap = label->hmap;
302 
303     ierr = PetscMalloc1(numStrata, &tmpV);CHKERRQ(ierr);
304     ierr = PetscMalloc1(numStrata, &tmpS);CHKERRQ(ierr);
305     ierr = PetscMalloc1(numStrata, &tmpH);CHKERRQ(ierr);
306     ierr = PetscMalloc1(numStrata, &tmpP);CHKERRQ(ierr);
307     ierr = PetscMalloc1(numStrata, &tmpB);CHKERRQ(ierr);
308     label->numStrata     = numStrata;
309     label->stratumValues = tmpV;
310     label->stratumSizes  = tmpS;
311     label->ht            = tmpH;
312     label->points        = tmpP;
313     label->validIS       = tmpB;
314     for (v = 0; v < numStrata; ++v) {
315       ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
316       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is);CHKERRQ(ierr);
317       ierr = PetscHMapISet(hmap, values[v], v);CHKERRQ(ierr);
318       tmpV[v] = values[v];
319       tmpS[v] = 0;
320       tmpH[v] = ht;
321       tmpP[v] = is;
322       tmpB[v] = PETSC_TRUE;
323     }
324     ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
325   } else {
326     for (v = 0; v < numStrata; ++v) {
327       ierr = DMLabelAddStratum(label, values[v]);CHKERRQ(ierr);
328     }
329   }
330   ierr = PetscFree(values);CHKERRQ(ierr);
331   PetscFunctionReturn(0);
332 }
333 
334 /*@
335   DMLabelAddStrataIS - Adds new stratum values in a DMLabel
336 
337   Input Parameter:
338 + label - The DMLabel
339 - valueIS - Index set with stratum values
340 
341   Level: beginner
342 
343 .seealso:  DMLabelCreate(), DMLabelDestroy()
344 @*/
345 PetscErrorCode DMLabelAddStrataIS(DMLabel label, IS valueIS)
346 {
347   PetscInt       numStrata;
348   const PetscInt *stratumValues;
349   PetscErrorCode ierr;
350 
351   PetscFunctionBegin;
352   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
353   PetscValidHeaderSpecific(valueIS, IS_CLASSID, 2);
354   ierr = ISGetLocalSize(valueIS, &numStrata);CHKERRQ(ierr);
355   ierr = ISGetIndices(valueIS, &stratumValues);CHKERRQ(ierr);
356   ierr = DMLabelAddStrata(label, numStrata, stratumValues);CHKERRQ(ierr);
357   PetscFunctionReturn(0);
358 }
359 
360 static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
361 {
362   PetscInt       v;
363   PetscMPIInt    rank;
364   PetscErrorCode ierr;
365 
366   PetscFunctionBegin;
367   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr);
368   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
369   if (label) {
370     const char *name;
371 
372     ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
373     ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", name);CHKERRQ(ierr);
374     if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);}
375     for (v = 0; v < label->numStrata; ++v) {
376       const PetscInt value = label->stratumValues[v];
377       const PetscInt *points;
378       PetscInt       p;
379 
380       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
381       for (p = 0; p < label->stratumSizes[v]; ++p) {
382         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr);
383       }
384       ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
385     }
386   }
387   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
388   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
389   PetscFunctionReturn(0);
390 }
391 
392 /*@C
393   DMLabelView - View the label
394 
395   Input Parameters:
396 + label - The DMLabel
397 - viewer - The PetscViewer
398 
399   Level: intermediate
400 
401 .seealso: DMLabelCreate(), DMLabelDestroy()
402 @*/
403 PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
404 {
405   PetscBool      iascii;
406   PetscErrorCode ierr;
407 
408   PetscFunctionBegin;
409   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
410   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)label), &viewer);CHKERRQ(ierr);}
411   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
412   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
413   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
414   if (iascii) {
415     ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr);
416   }
417   PetscFunctionReturn(0);
418 }
419 
420 /*@
421   DMLabelReset - Destroys internal data structures in a DMLabel
422 
423   Input Parameter:
424 . label - The DMLabel
425 
426   Level: beginner
427 
428 .seealso: DMLabelDestroy(), DMLabelCreate()
429 @*/
430 PetscErrorCode DMLabelReset(DMLabel label)
431 {
432   PetscInt       v;
433   PetscErrorCode ierr;
434 
435   PetscFunctionBegin;
436   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
437   for (v = 0; v < label->numStrata; ++v) {
438     ierr = PetscHSetIDestroy(&label->ht[v]);CHKERRQ(ierr);
439     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
440   }
441   label->numStrata = 0;
442   ierr = PetscFree(label->stratumValues);CHKERRQ(ierr);
443   ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr);
444   ierr = PetscFree(label->ht);CHKERRQ(ierr);
445   ierr = PetscFree(label->points);CHKERRQ(ierr);
446   ierr = PetscFree(label->validIS);CHKERRQ(ierr);
447   ierr = PetscHMapIReset(label->hmap);CHKERRQ(ierr);
448   label->pStart = -1;
449   label->pEnd   = -1;
450   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
451   PetscFunctionReturn(0);
452 }
453 
454 /*@
455   DMLabelDestroy - Destroys a DMLabel
456 
457   Input Parameter:
458 . label - The DMLabel
459 
460   Level: beginner
461 
462 .seealso: DMLabelReset(), DMLabelCreate()
463 @*/
464 PetscErrorCode DMLabelDestroy(DMLabel *label)
465 {
466   PetscErrorCode ierr;
467 
468   PetscFunctionBegin;
469   if (!*label) PetscFunctionReturn(0);
470   PetscValidHeaderSpecific((*label),DMLABEL_CLASSID,1);
471   if (--((PetscObject)(*label))->refct > 0) {*label = NULL; PetscFunctionReturn(0);}
472   ierr = DMLabelReset(*label);CHKERRQ(ierr);
473   ierr = PetscHMapIDestroy(&(*label)->hmap);CHKERRQ(ierr);
474   ierr = PetscHeaderDestroy(label);CHKERRQ(ierr);
475   PetscFunctionReturn(0);
476 }
477 
478 /*@
479   DMLabelDuplicate - Duplicates a DMLabel
480 
481   Input Parameter:
482 . label - The DMLabel
483 
484   Output Parameter:
485 . labelnew - location to put new vector
486 
487   Level: intermediate
488 
489 .seealso: DMLabelCreate(), DMLabelDestroy()
490 @*/
491 PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
492 {
493   const char    *name;
494   PetscInt       v;
495   PetscErrorCode ierr;
496 
497   PetscFunctionBegin;
498   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
499   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
500   ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
501   ierr = DMLabelCreate(PetscObjectComm((PetscObject) label), name, labelnew);CHKERRQ(ierr);
502 
503   (*labelnew)->numStrata    = label->numStrata;
504   (*labelnew)->defaultValue = label->defaultValue;
505   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr);
506   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr);
507   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr);
508   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr);
509   ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr);
510   for (v = 0; v < label->numStrata; ++v) {
511     ierr = PetscHSetICreate(&(*labelnew)->ht[v]);CHKERRQ(ierr);
512     (*labelnew)->stratumValues[v]  = label->stratumValues[v];
513     (*labelnew)->stratumSizes[v]   = label->stratumSizes[v];
514     ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr);
515     (*labelnew)->points[v]         = label->points[v];
516     (*labelnew)->validIS[v]        = PETSC_TRUE;
517   }
518   ierr = PetscHMapIDestroy(&(*labelnew)->hmap);CHKERRQ(ierr);
519   ierr = PetscHMapIDuplicate(label->hmap,&(*labelnew)->hmap);CHKERRQ(ierr);
520   (*labelnew)->pStart = -1;
521   (*labelnew)->pEnd   = -1;
522   (*labelnew)->bt     = NULL;
523   PetscFunctionReturn(0);
524 }
525 
526 /*@
527   DMLabelComputeIndex - Create an index structure for membership determination, automatically determining the bounds
528 
529   Input Parameter:
530 . label  - The DMLabel
531 
532   Level: intermediate
533 
534 .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
535 @*/
536 PetscErrorCode DMLabelComputeIndex(DMLabel label)
537 {
538   PetscInt       pStart = PETSC_MAX_INT, pEnd = -1, v;
539   PetscErrorCode ierr;
540 
541   PetscFunctionBegin;
542   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
543   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
544   for (v = 0; v < label->numStrata; ++v) {
545     const PetscInt *points;
546     PetscInt       i;
547 
548     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
549     for (i = 0; i < label->stratumSizes[v]; ++i) {
550       const PetscInt point = points[i];
551 
552       pStart = PetscMin(point,   pStart);
553       pEnd   = PetscMax(point+1, pEnd);
554     }
555     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
556   }
557   label->pStart = pStart == PETSC_MAX_INT ? -1 : pStart;
558   label->pEnd   = pEnd;
559   ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
560   PetscFunctionReturn(0);
561 }
562 
563 /*@
564   DMLabelCreateIndex - Create an index structure for membership determination
565 
566   Input Parameters:
567 + label  - The DMLabel
568 . pStart - The smallest point
569 - pEnd   - The largest point + 1
570 
571   Level: intermediate
572 
573 .seealso: DMLabelHasPoint(), DMLabelComputeIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
574 @*/
575 PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
576 {
577   PetscInt       v;
578   PetscErrorCode ierr;
579 
580   PetscFunctionBegin;
581   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
582   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
583   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
584   label->pStart = pStart;
585   label->pEnd   = pEnd;
586   /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
587   ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr);
588   for (v = 0; v < label->numStrata; ++v) {
589     const PetscInt *points;
590     PetscInt       i;
591 
592     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
593     for (i = 0; i < label->stratumSizes[v]; ++i) {
594       const PetscInt point = points[i];
595 
596       if ((point < pStart) || (point >= pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, pStart, pEnd);
597       ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr);
598     }
599     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
600   }
601   PetscFunctionReturn(0);
602 }
603 
604 /*@
605   DMLabelDestroyIndex - Destroy the index structure
606 
607   Input Parameter:
608 . label - the DMLabel
609 
610   Level: intermediate
611 
612 .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
613 @*/
614 PetscErrorCode DMLabelDestroyIndex(DMLabel label)
615 {
616   PetscErrorCode ierr;
617 
618   PetscFunctionBegin;
619   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
620   label->pStart = -1;
621   label->pEnd   = -1;
622   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
623   PetscFunctionReturn(0);
624 }
625 
626 /*@
627   DMLabelGetBounds - Return the smallest and largest point in the label
628 
629   Input Parameter:
630 . label - the DMLabel
631 
632   Output Parameters:
633 + pStart - The smallest point
634 - pEnd   - The largest point + 1
635 
636   Note: This will compute an index for the label if one does not exist.
637 
638   Level: intermediate
639 
640 .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
641 @*/
642 PetscErrorCode DMLabelGetBounds(DMLabel label, PetscInt *pStart, PetscInt *pEnd)
643 {
644   PetscErrorCode ierr;
645 
646   PetscFunctionBegin;
647   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
648   if ((label->pStart == -1) && (label->pEnd == -1)) {ierr = DMLabelComputeIndex(label);CHKERRQ(ierr);}
649   if (pStart) {
650     PetscValidPointer(pStart, 2);
651     *pStart = label->pStart;
652   }
653   if (pEnd) {
654     PetscValidPointer(pEnd, 3);
655     *pEnd = label->pEnd;
656   }
657   PetscFunctionReturn(0);
658 }
659 
660 /*@
661   DMLabelHasValue - Determine whether a label assigns the value to any point
662 
663   Input Parameters:
664 + label - the DMLabel
665 - value - the value
666 
667   Output Parameter:
668 . contains - Flag indicating whether the label maps this value to any point
669 
670   Level: developer
671 
672 .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
673 @*/
674 PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
675 {
676   PetscInt v;
677   PetscErrorCode ierr;
678 
679   PetscFunctionBegin;
680   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
681   PetscValidPointer(contains, 3);
682   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
683   *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE;
684   PetscFunctionReturn(0);
685 }
686 
687 /*@
688   DMLabelHasPoint - Determine whether a label assigns a value to a point
689 
690   Input Parameters:
691 + label - the DMLabel
692 - point - the point
693 
694   Output Parameter:
695 . contains - Flag indicating whether the label maps this point to a value
696 
697   Note: The user must call DMLabelCreateIndex() before this function.
698 
699   Level: developer
700 
701 .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
702 @*/
703 PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
704 {
705   PetscErrorCode ierr;
706 
707   PetscFunctionBeginHot;
708   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
709   PetscValidPointer(contains, 3);
710   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
711 #if defined(PETSC_USE_DEBUG)
712   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
713   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);
714 #endif
715   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
716   PetscFunctionReturn(0);
717 }
718 
719 /*@
720   DMLabelStratumHasPoint - Return true if the stratum contains a point
721 
722   Input Parameters:
723 + label - the DMLabel
724 . value - the stratum value
725 - point - the point
726 
727   Output Parameter:
728 . contains - true if the stratum contains the point
729 
730   Level: intermediate
731 
732 .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
733 @*/
734 PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
735 {
736   PetscInt       v;
737   PetscErrorCode ierr;
738 
739   PetscFunctionBeginHot;
740   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
741   PetscValidPointer(contains, 4);
742   *contains = PETSC_FALSE;
743   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
744   if (v < 0) PetscFunctionReturn(0);
745 
746   if (label->validIS[v]) {
747     PetscInt i;
748 
749     ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
750     if (i >= 0) *contains = PETSC_TRUE;
751   } else {
752     PetscBool has;
753 
754     ierr = PetscHSetIHas(label->ht[v], point, &has);CHKERRQ(ierr);
755     if (has) *contains = PETSC_TRUE;
756   }
757   PetscFunctionReturn(0);
758 }
759 
760 /*@
761   DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
762   When a label is created, it is initialized to -1.
763 
764   Input parameter:
765 . label - a DMLabel object
766 
767   Output parameter:
768 . defaultValue - the default value
769 
770   Level: beginner
771 
772 .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
773 @*/
774 PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
775 {
776   PetscFunctionBegin;
777   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
778   *defaultValue = label->defaultValue;
779   PetscFunctionReturn(0);
780 }
781 
782 /*@
783   DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
784   When a label is created, it is initialized to -1.
785 
786   Input parameter:
787 . label - a DMLabel object
788 
789   Output parameter:
790 . defaultValue - the default value
791 
792   Level: beginner
793 
794 .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
795 @*/
796 PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
797 {
798   PetscFunctionBegin;
799   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
800   label->defaultValue = defaultValue;
801   PetscFunctionReturn(0);
802 }
803 
804 /*@
805   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())
806 
807   Input Parameters:
808 + label - the DMLabel
809 - point - the point
810 
811   Output Parameter:
812 . value - The point value, or the default value (-1 by default)
813 
814   Level: intermediate
815 
816 .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
817 @*/
818 PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
819 {
820   PetscInt       v;
821   PetscErrorCode ierr;
822 
823   PetscFunctionBeginHot;
824   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
825   PetscValidPointer(value, 3);
826   *value = label->defaultValue;
827   for (v = 0; v < label->numStrata; ++v) {
828     if (label->validIS[v]) {
829       PetscInt i;
830 
831       ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
832       if (i >= 0) {
833         *value = label->stratumValues[v];
834         break;
835       }
836     } else {
837       PetscBool has;
838 
839       ierr = PetscHSetIHas(label->ht[v], point, &has);CHKERRQ(ierr);
840       if (has) {
841         *value = label->stratumValues[v];
842         break;
843       }
844     }
845   }
846   PetscFunctionReturn(0);
847 }
848 
849 /*@
850   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.
851 
852   Input Parameters:
853 + label - the DMLabel
854 . point - the point
855 - value - The point value
856 
857   Level: intermediate
858 
859 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
860 @*/
861 PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
862 {
863   PetscInt       v;
864   PetscErrorCode ierr;
865 
866   PetscFunctionBegin;
867   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
868   /* Find label value, add new entry if needed */
869   if (value == label->defaultValue) PetscFunctionReturn(0);
870   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
871   /* Set key */
872   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
873   ierr = PetscHSetIAdd(label->ht[v], point);CHKERRQ(ierr);
874   PetscFunctionReturn(0);
875 }
876 
877 /*@
878   DMLabelClearValue - Clear the value a label assigns to a point
879 
880   Input Parameters:
881 + label - the DMLabel
882 . point - the point
883 - value - The point value
884 
885   Level: intermediate
886 
887 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
888 @*/
889 PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
890 {
891   PetscInt       v;
892   PetscErrorCode ierr;
893 
894   PetscFunctionBegin;
895   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
896   /* Find label value */
897   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
898   if (v < 0) PetscFunctionReturn(0);
899 
900   if (label->bt) {
901     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);
902     ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
903   }
904 
905   /* Delete key */
906   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
907   ierr = PetscHSetIDel(label->ht[v], point);CHKERRQ(ierr);
908   PetscFunctionReturn(0);
909 }
910 
911 /*@
912   DMLabelInsertIS - Set all points in the IS to a value
913 
914   Input Parameters:
915 + label - the DMLabel
916 . is    - the point IS
917 - value - The point value
918 
919   Level: intermediate
920 
921 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
922 @*/
923 PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
924 {
925   PetscInt        v, n, p;
926   const PetscInt *points;
927   PetscErrorCode  ierr;
928 
929   PetscFunctionBegin;
930   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
931   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
932   /* Find label value, add new entry if needed */
933   if (value == label->defaultValue) PetscFunctionReturn(0);
934   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
935   /* Set keys */
936   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
937   ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr);
938   ierr = ISGetIndices(is, &points);CHKERRQ(ierr);
939   for (p = 0; p < n; ++p) {ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);}
940   ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr);
941   PetscFunctionReturn(0);
942 }
943 
944 /*@
945   DMLabelGetNumValues - Get the number of values that the DMLabel takes
946 
947   Input Parameter:
948 . label - the DMLabel
949 
950   Output Paramater:
951 . numValues - the number of values
952 
953   Level: intermediate
954 
955 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
956 @*/
957 PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
958 {
959   PetscFunctionBegin;
960   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
961   PetscValidPointer(numValues, 2);
962   *numValues = label->numStrata;
963   PetscFunctionReturn(0);
964 }
965 
966 /*@
967   DMLabelGetValueIS - Get an IS of all values that the DMlabel takes
968 
969   Input Parameter:
970 . label - the DMLabel
971 
972   Output Paramater:
973 . is    - the value IS
974 
975   Level: intermediate
976 
977 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
978 @*/
979 PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
980 {
981   PetscErrorCode ierr;
982 
983   PetscFunctionBegin;
984   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
985   PetscValidPointer(values, 2);
986   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr);
987   PetscFunctionReturn(0);
988 }
989 
990 /*@
991   DMLabelHasStratum - Determine whether points exist with the given value
992 
993   Input Parameters:
994 + label - the DMLabel
995 - value - the stratum value
996 
997   Output Paramater:
998 . exists - Flag saying whether points exist
999 
1000   Level: intermediate
1001 
1002 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1003 @*/
1004 PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
1005 {
1006   PetscInt       v;
1007   PetscErrorCode ierr;
1008 
1009   PetscFunctionBegin;
1010   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1011   PetscValidPointer(exists, 3);
1012   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
1013   *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE;
1014   PetscFunctionReturn(0);
1015 }
1016 
1017 /*@
1018   DMLabelGetStratumSize - Get the size of a stratum
1019 
1020   Input Parameters:
1021 + label - the DMLabel
1022 - value - the stratum value
1023 
1024   Output Paramater:
1025 . size - The number of points in the stratum
1026 
1027   Level: intermediate
1028 
1029 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1030 @*/
1031 PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
1032 {
1033   PetscInt       v;
1034   PetscErrorCode ierr;
1035 
1036   PetscFunctionBegin;
1037   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1038   PetscValidPointer(size, 3);
1039   *size = 0;
1040   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
1041   if (v < 0) PetscFunctionReturn(0);
1042   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
1043   *size = label->stratumSizes[v];
1044   PetscFunctionReturn(0);
1045 }
1046 
1047 /*@
1048   DMLabelGetStratumBounds - Get the largest and smallest point of a stratum
1049 
1050   Input Parameters:
1051 + label - the DMLabel
1052 - value - the stratum value
1053 
1054   Output Paramaters:
1055 + start - the smallest point in the stratum
1056 - end - the largest point in the stratum
1057 
1058   Level: intermediate
1059 
1060 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1061 @*/
1062 PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
1063 {
1064   PetscInt       v, min, max;
1065   PetscErrorCode ierr;
1066 
1067   PetscFunctionBegin;
1068   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1069   if (start) {PetscValidPointer(start, 3); *start = 0;}
1070   if (end)   {PetscValidPointer(end,   4); *end   = 0;}
1071   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
1072   if (v < 0) PetscFunctionReturn(0);
1073   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
1074   if (label->stratumSizes[v] <= 0) PetscFunctionReturn(0);
1075   ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr);
1076   if (start) *start = min;
1077   if (end)   *end   = max+1;
1078   PetscFunctionReturn(0);
1079 }
1080 
1081 /*@
1082   DMLabelGetStratumIS - Get an IS with the stratum points
1083 
1084   Input Parameters:
1085 + label - the DMLabel
1086 - value - the stratum value
1087 
1088   Output Paramater:
1089 . points - The stratum points
1090 
1091   Level: intermediate
1092 
1093 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1094 @*/
1095 PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
1096 {
1097   PetscInt       v;
1098   PetscErrorCode ierr;
1099 
1100   PetscFunctionBegin;
1101   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1102   PetscValidPointer(points, 3);
1103   *points = NULL;
1104   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
1105   if (v < 0) PetscFunctionReturn(0);
1106   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
1107   ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr);
1108   *points = label->points[v];
1109   PetscFunctionReturn(0);
1110 }
1111 
1112 /*@
1113   DMLabelSetStratumIS - Set the stratum points using an IS
1114 
1115   Input Parameters:
1116 + label - the DMLabel
1117 . value - the stratum value
1118 - points - The stratum points
1119 
1120   Level: intermediate
1121 
1122 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1123 @*/
1124 PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
1125 {
1126   PetscInt       v;
1127   PetscErrorCode ierr;
1128 
1129   PetscFunctionBegin;
1130   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1131   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
1132   ierr = DMLabelLookupAddStratum(label, value, &v);CHKERRQ(ierr);
1133   if (is == label->points[v]) PetscFunctionReturn(0);
1134   ierr = DMLabelClearStratum(label, value);CHKERRQ(ierr);
1135   ierr = ISGetLocalSize(is, &(label->stratumSizes[v]));CHKERRQ(ierr);
1136   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1137   ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
1138   label->points[v]  = is;
1139   label->validIS[v] = PETSC_TRUE;
1140   ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
1141   if (label->bt) {
1142     const PetscInt *points;
1143     PetscInt p;
1144 
1145     ierr = ISGetIndices(is,&points);CHKERRQ(ierr);
1146     for (p = 0; p < label->stratumSizes[v]; ++p) {
1147       const PetscInt point = points[p];
1148 
1149       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);
1150       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
1151     }
1152   }
1153   PetscFunctionReturn(0);
1154 }
1155 
1156 /*@
1157   DMLabelClearStratum - Remove a stratum
1158 
1159   Input Parameters:
1160 + label - the DMLabel
1161 - value - the stratum value
1162 
1163   Level: intermediate
1164 
1165 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1166 @*/
1167 PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
1168 {
1169   PetscInt       v;
1170   PetscErrorCode ierr;
1171 
1172   PetscFunctionBegin;
1173   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1174   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
1175   if (v < 0) PetscFunctionReturn(0);
1176   if (label->validIS[v]) {
1177     if (label->bt) {
1178       PetscInt       i;
1179       const PetscInt *points;
1180 
1181       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
1182       for (i = 0; i < label->stratumSizes[v]; ++i) {
1183         const PetscInt point = points[i];
1184 
1185         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);
1186         ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
1187       }
1188       ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
1189     }
1190     label->stratumSizes[v] = 0;
1191     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
1192     ierr = ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &label->points[v]);CHKERRQ(ierr);
1193     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
1194     ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
1195   } else {
1196     ierr = PetscHSetIClear(label->ht[v]);CHKERRQ(ierr);
1197   }
1198   PetscFunctionReturn(0);
1199 }
1200 
1201 /*@
1202   DMLabelFilter - Remove all points outside of [start, end)
1203 
1204   Input Parameters:
1205 + label - the DMLabel
1206 . start - the first point
1207 - end - the last point
1208 
1209   Level: intermediate
1210 
1211 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1212 @*/
1213 PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
1214 {
1215   PetscInt       v;
1216   PetscErrorCode ierr;
1217 
1218   PetscFunctionBegin;
1219   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1220   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
1221   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1222   for (v = 0; v < label->numStrata; ++v) {
1223     PetscInt off, q;
1224     const PetscInt *points;
1225     PetscInt numPointsNew = 0, *pointsNew = NULL;
1226 
1227     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
1228     for (q = 0; q < label->stratumSizes[v]; ++q)
1229       if (points[q] >= start && points[q] < end)
1230         numPointsNew++;
1231     ierr = PetscMalloc1(numPointsNew, &pointsNew);CHKERRQ(ierr);
1232     for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) {
1233       if (points[q] >= start && points[q] < end)
1234         pointsNew[off++] = points[q];
1235     }
1236     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
1237 
1238     label->stratumSizes[v] = numPointsNew;
1239     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
1240     ierr = ISCreateGeneral(PETSC_COMM_SELF,numPointsNew, pointsNew, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr);
1241     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
1242   }
1243   ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr);
1244   PetscFunctionReturn(0);
1245 }
1246 
1247 /*@
1248   DMLabelPermute - Create a new label with permuted points
1249 
1250   Input Parameters:
1251 + label - the DMLabel
1252 - permutation - the point permutation
1253 
1254   Output Parameter:
1255 . labelnew - the new label containing the permuted points
1256 
1257   Level: intermediate
1258 
1259 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1260 @*/
1261 PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
1262 {
1263   const PetscInt *perm;
1264   PetscInt        numValues, numPoints, v, q;
1265   PetscErrorCode  ierr;
1266 
1267   PetscFunctionBegin;
1268   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1269   PetscValidHeaderSpecific(permutation, IS_CLASSID, 2);
1270   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1271   ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr);
1272   ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr);
1273   ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr);
1274   ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr);
1275   for (v = 0; v < numValues; ++v) {
1276     const PetscInt size   = (*labelNew)->stratumSizes[v];
1277     const PetscInt *points;
1278     PetscInt *pointsNew;
1279 
1280     ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1281     ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr);
1282     for (q = 0; q < size; ++q) {
1283       const PetscInt point = points[q];
1284 
1285       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);
1286       pointsNew[q] = perm[point];
1287     }
1288     ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1289     ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr);
1290     ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr);
1291     if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) {
1292       ierr = ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));CHKERRQ(ierr);
1293       ierr = PetscFree(pointsNew);CHKERRQ(ierr);
1294     } else {
1295       ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr);
1296     }
1297     ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr);
1298   }
1299   ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr);
1300   if (label->bt) {
1301     ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
1302     ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
1303   }
1304   PetscFunctionReturn(0);
1305 }
1306 
1307 PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
1308 {
1309   MPI_Comm       comm;
1310   PetscInt       s, l, nroots, nleaves, dof, offset, size;
1311   PetscInt      *remoteOffsets, *rootStrata, *rootIdx;
1312   PetscSection   rootSection;
1313   PetscSF        labelSF;
1314   PetscErrorCode ierr;
1315 
1316   PetscFunctionBegin;
1317   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
1318   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1319   /* Build a section of stratum values per point, generate the according SF
1320      and distribute point-wise stratum values to leaves. */
1321   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
1322   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
1323   ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr);
1324   if (label) {
1325     for (s = 0; s < label->numStrata; ++s) {
1326       const PetscInt *points;
1327 
1328       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
1329       for (l = 0; l < label->stratumSizes[s]; l++) {
1330         ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr);
1331         ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr);
1332       }
1333       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
1334     }
1335   }
1336   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
1337   /* Create a point-wise array of stratum values */
1338   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
1339   ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr);
1340   ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr);
1341   if (label) {
1342     for (s = 0; s < label->numStrata; ++s) {
1343       const PetscInt *points;
1344 
1345       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
1346       for (l = 0; l < label->stratumSizes[s]; l++) {
1347         const PetscInt p = points[l];
1348         ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr);
1349         rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
1350       }
1351       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
1352     }
1353   }
1354   /* Build SF that maps label points to remote processes */
1355   ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr);
1356   ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr);
1357   ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr);
1358   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
1359   /* Send the strata for each point over the derived SF */
1360   ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr);
1361   ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr);
1362   ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
1363   ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
1364   /* Clean up */
1365   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
1366   ierr = PetscFree(rootIdx);CHKERRQ(ierr);
1367   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
1368   ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr);
1369   PetscFunctionReturn(0);
1370 }
1371 
1372 /*@
1373   DMLabelDistribute - Create a new label pushed forward over the PetscSF
1374 
1375   Input Parameters:
1376 + label - the DMLabel
1377 - sf    - the map from old to new distribution
1378 
1379   Output Parameter:
1380 . labelnew - the new redistributed label
1381 
1382   Level: intermediate
1383 
1384 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1385 @*/
1386 PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1387 {
1388   MPI_Comm       comm;
1389   PetscSection   leafSection;
1390   PetscInt       p, pStart, pEnd, s, size, dof, offset, stratum;
1391   PetscInt      *leafStrata, *strataIdx;
1392   PetscInt     **points;
1393   const char    *lname = NULL;
1394   char          *name;
1395   PetscInt       nameSize;
1396   PetscHSetI     stratumHash;
1397   size_t         len = 0;
1398   PetscMPIInt    rank;
1399   PetscErrorCode ierr;
1400 
1401   PetscFunctionBegin;
1402   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1403   if (label) {
1404     PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1405     ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1406   }
1407   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1408   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1409   /* Bcast name */
1410   if (!rank) {
1411     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1412     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1413   }
1414   nameSize = len;
1415   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1416   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1417   if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);}
1418   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1419   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
1420   ierr = PetscFree(name);CHKERRQ(ierr);
1421   /* Bcast defaultValue */
1422   if (!rank) (*labelNew)->defaultValue = label->defaultValue;
1423   ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1424   /* Distribute stratum values over the SF and get the point mapping on the receiver */
1425   ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr);
1426   /* Determine received stratum values and initialise new label*/
1427   ierr = PetscHSetICreate(&stratumHash);CHKERRQ(ierr);
1428   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
1429   for (p = 0; p < size; ++p) {ierr = PetscHSetIAdd(stratumHash, leafStrata[p]);CHKERRQ(ierr);}
1430   ierr = PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata);CHKERRQ(ierr);
1431   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr);
1432   for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
1433   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr);
1434   /* Turn leafStrata into indices rather than stratum values */
1435   offset = 0;
1436   ierr = PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr);
1437   ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr);
1438   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1439     ierr = PetscHMapISet((*labelNew)->hmap, (*labelNew)->stratumValues[s], s);CHKERRQ(ierr);
1440   }
1441   for (p = 0; p < size; ++p) {
1442     for (s = 0; s < (*labelNew)->numStrata; ++s) {
1443       if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
1444     }
1445   }
1446   /* Rebuild the point strata on the receiver */
1447   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr);
1448   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1449   for (p=pStart; p<pEnd; p++) {
1450     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1451     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1452     for (s=0; s<dof; s++) {
1453       (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1454     }
1455   }
1456   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr);
1457   ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr);
1458   ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr);
1459   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1460     ierr = PetscHSetICreate(&(*labelNew)->ht[s]);CHKERRQ(ierr);
1461     ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr);
1462   }
1463   /* Insert points into new strata */
1464   ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr);
1465   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1466   for (p=pStart; p<pEnd; p++) {
1467     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1468     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1469     for (s=0; s<dof; s++) {
1470       stratum = leafStrata[offset+s];
1471       points[stratum][strataIdx[stratum]++] = p;
1472     }
1473   }
1474   for (s = 0; s < (*labelNew)->numStrata; s++) {
1475     ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr);
1476     ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr);
1477   }
1478   ierr = PetscFree(points);CHKERRQ(ierr);
1479   ierr = PetscHSetIDestroy(&stratumHash);CHKERRQ(ierr);
1480   ierr = PetscFree(leafStrata);CHKERRQ(ierr);
1481   ierr = PetscFree(strataIdx);CHKERRQ(ierr);
1482   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
1483   PetscFunctionReturn(0);
1484 }
1485 
1486 /*@
1487   DMLabelGather - Gather all label values from leafs into roots
1488 
1489   Input Parameters:
1490 + label - the DMLabel
1491 - sf - the Star Forest point communication map
1492 
1493   Output Parameters:
1494 . labelNew - the new DMLabel with localised leaf values
1495 
1496   Level: developer
1497 
1498   Note: This is the inverse operation to DMLabelDistribute.
1499 
1500 .seealso: DMLabelDistribute()
1501 @*/
1502 PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
1503 {
1504   MPI_Comm       comm;
1505   PetscSection   rootSection;
1506   PetscSF        sfLabel;
1507   PetscSFNode   *rootPoints, *leafPoints;
1508   PetscInt       p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
1509   const PetscInt *rootDegree, *ilocal;
1510   PetscInt       *rootStrata;
1511   const char    *lname;
1512   char          *name;
1513   PetscInt       nameSize;
1514   size_t         len = 0;
1515   PetscMPIInt    rank, size;
1516   PetscErrorCode ierr;
1517 
1518   PetscFunctionBegin;
1519   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1520   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1521   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1522   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1523   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1524   /* Bcast name */
1525   if (!rank) {
1526     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1527     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1528   }
1529   nameSize = len;
1530   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1531   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1532   if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);}
1533   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1534   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
1535   ierr = PetscFree(name);CHKERRQ(ierr);
1536   /* Gather rank/index pairs of leaves into local roots to build
1537      an inverse, multi-rooted SF. Note that this ignores local leaf
1538      indexing due to the use of the multiSF in PetscSFGather. */
1539   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr);
1540   ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr);
1541   for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
1542   for (p = 0; p < nleaves; p++) {
1543     PetscInt ilp = ilocal ? ilocal[p] : p;
1544 
1545     leafPoints[ilp].index = ilp;
1546     leafPoints[ilp].rank  = rank;
1547   }
1548   ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr);
1549   ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr);
1550   for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
1551   ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr);
1552   ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
1553   ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
1554   ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr);
1555   ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
1556   /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
1557   ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr);
1558   /* Rebuild the point strata on the receiver */
1559   for (p = 0, idx = 0; p < nroots; p++) {
1560     for (d = 0; d < rootDegree[p]; d++) {
1561       ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr);
1562       ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr);
1563       for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);}
1564     }
1565     idx += rootDegree[p];
1566   }
1567   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
1568   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
1569   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
1570   ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr);
1571   PetscFunctionReturn(0);
1572 }
1573 
1574 /*@
1575   DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label
1576 
1577   Input Parameter:
1578 . label - the DMLabel
1579 
1580   Output Parameters:
1581 + section - the section giving offsets for each stratum
1582 - is - An IS containing all the label points
1583 
1584   Level: developer
1585 
1586 .seealso: DMLabelDistribute()
1587 @*/
1588 PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1589 {
1590   IS              vIS;
1591   const PetscInt *values;
1592   PetscInt       *points;
1593   PetscInt        nV, vS = 0, vE = 0, v, N;
1594   PetscErrorCode  ierr;
1595 
1596   PetscFunctionBegin;
1597   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1598   ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr);
1599   ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
1600   ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr);
1601   if (nV) {vS = values[0]; vE = values[0]+1;}
1602   for (v = 1; v < nV; ++v) {
1603     vS = PetscMin(vS, values[v]);
1604     vE = PetscMax(vE, values[v]+1);
1605   }
1606   ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr);
1607   ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr);
1608   for (v = 0; v < nV; ++v) {
1609     PetscInt n;
1610 
1611     ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr);
1612     ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr);
1613   }
1614   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
1615   ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr);
1616   ierr = PetscMalloc1(N, &points);CHKERRQ(ierr);
1617   for (v = 0; v < nV; ++v) {
1618     IS              is;
1619     const PetscInt *spoints;
1620     PetscInt        dof, off, p;
1621 
1622     ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr);
1623     ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr);
1624     ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr);
1625     ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr);
1626     for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1627     ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr);
1628     ierr = ISDestroy(&is);CHKERRQ(ierr);
1629   }
1630   ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr);
1631   ierr = ISDestroy(&vIS);CHKERRQ(ierr);
1632   ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1633   PetscFunctionReturn(0);
1634 }
1635 
1636 /*@
1637   PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1638   the local section and an SF describing the section point overlap.
1639 
1640   Input Parameters:
1641   + s - The PetscSection for the local field layout
1642   . sf - The SF describing parallel layout of the section points
1643   . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1644   . label - The label specifying the points
1645   - labelValue - The label stratum specifying the points
1646 
1647   Output Parameter:
1648   . gsection - The PetscSection for the global field layout
1649 
1650   Note: This gives negative sizes and offsets to points not owned by this process
1651 
1652   Level: developer
1653 
1654 .seealso: PetscSectionCreate()
1655 @*/
1656 PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1657 {
1658   PetscInt      *neg = NULL, *tmpOff = NULL;
1659   PetscInt       pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1660   PetscErrorCode ierr;
1661 
1662   PetscFunctionBegin;
1663   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1664   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1665   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4);
1666   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr);
1667   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1668   ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr);
1669   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1670   if (nroots >= 0) {
1671     if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1672     ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr);
1673     if (nroots > pEnd-pStart) {
1674       ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr);
1675     } else {
1676       tmpOff = &(*gsection)->atlasDof[-pStart];
1677     }
1678   }
1679   /* Mark ghost points with negative dof */
1680   for (p = pStart; p < pEnd; ++p) {
1681     PetscInt value;
1682 
1683     ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr);
1684     if (value != labelValue) continue;
1685     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1686     ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr);
1687     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1688     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);}
1689     if (neg) neg[p] = -(dof+1);
1690   }
1691   ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr);
1692   if (nroots >= 0) {
1693     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1694     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1695     if (nroots > pEnd-pStart) {
1696       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1697     }
1698   }
1699   /* Calculate new sizes, get proccess offset, and calculate point offsets */
1700   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1701     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1702     (*gsection)->atlasOff[p] = off;
1703     off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1704   }
1705   ierr       = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1706   globalOff -= off;
1707   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1708     (*gsection)->atlasOff[p] += globalOff;
1709     if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1710   }
1711   /* Put in negative offsets for ghost points */
1712   if (nroots >= 0) {
1713     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1714     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1715     if (nroots > pEnd-pStart) {
1716       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1717     }
1718   }
1719   if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);}
1720   ierr = PetscFree(neg);CHKERRQ(ierr);
1721   PetscFunctionReturn(0);
1722 }
1723 
1724 typedef struct _n_PetscSectionSym_Label
1725 {
1726   DMLabel           label;
1727   PetscCopyMode     *modes;
1728   PetscInt          *sizes;
1729   const PetscInt    ***perms;
1730   const PetscScalar ***rots;
1731   PetscInt          (*minMaxOrients)[2];
1732   PetscInt          numStrata; /* numStrata is only increasing, functions as a state */
1733 } PetscSectionSym_Label;
1734 
1735 static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
1736 {
1737   PetscInt              i, j;
1738   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
1739   PetscErrorCode        ierr;
1740 
1741   PetscFunctionBegin;
1742   for (i = 0; i <= sl->numStrata; i++) {
1743     if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
1744       for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
1745         if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);}
1746         if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);}
1747       }
1748       if (sl->perms[i]) {
1749         const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
1750 
1751         ierr = PetscFree(perms);CHKERRQ(ierr);
1752       }
1753       if (sl->rots[i]) {
1754         const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
1755 
1756         ierr = PetscFree(rots);CHKERRQ(ierr);
1757       }
1758     }
1759   }
1760   ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr);
1761   ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr);
1762   sl->numStrata = 0;
1763   PetscFunctionReturn(0);
1764 }
1765 
1766 static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
1767 {
1768   PetscErrorCode ierr;
1769 
1770   PetscFunctionBegin;
1771   ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);
1772   ierr = PetscFree(sym->data);CHKERRQ(ierr);
1773   PetscFunctionReturn(0);
1774 }
1775 
1776 static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
1777 {
1778   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
1779   PetscBool             isAscii;
1780   DMLabel               label = sl->label;
1781   const char           *name;
1782   PetscErrorCode        ierr;
1783 
1784   PetscFunctionBegin;
1785   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
1786   if (isAscii) {
1787     PetscInt          i, j, k;
1788     PetscViewerFormat format;
1789 
1790     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1791     if (label) {
1792       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1793       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1794         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1795         ierr = DMLabelView(label, viewer);CHKERRQ(ierr);
1796         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1797       } else {
1798         ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1799         ierr = PetscViewerASCIIPrintf(viewer,"  Label '%s'\n",name);CHKERRQ(ierr);
1800       }
1801     } else {
1802       ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr);
1803     }
1804     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1805     for (i = 0; i <= sl->numStrata; i++) {
1806       PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
1807 
1808       if (!(sl->perms[i] || sl->rots[i])) {
1809         ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr);
1810       } else {
1811       ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr);
1812         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1813         ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr);
1814         if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1815           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1816           for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
1817             if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
1818               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr);
1819             } else {
1820               PetscInt tab;
1821 
1822               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr);
1823               ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1824               ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr);
1825               if (sl->perms[i] && sl->perms[i][j]) {
1826                 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr);
1827                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
1828                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);}
1829                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1830                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
1831               }
1832               if (sl->rots[i] && sl->rots[i][j]) {
1833                 ierr = PetscViewerASCIIPrintf(viewer,"Rotations:  ");CHKERRQ(ierr);
1834                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
1835 #if defined(PETSC_USE_COMPLEX)
1836                 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);}
1837 #else
1838                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);}
1839 #endif
1840                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1841                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
1842               }
1843               ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1844             }
1845           }
1846           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1847         }
1848         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1849       }
1850     }
1851     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1852   }
1853   PetscFunctionReturn(0);
1854 }
1855 
1856 /*@
1857   PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
1858 
1859   Logically collective on sym
1860 
1861   Input parameters:
1862 + sym - the section symmetries
1863 - label - the DMLabel describing the types of points
1864 
1865   Level: developer:
1866 
1867 .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()
1868 @*/
1869 PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
1870 {
1871   PetscSectionSym_Label *sl;
1872   PetscErrorCode        ierr;
1873 
1874   PetscFunctionBegin;
1875   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
1876   sl = (PetscSectionSym_Label *) sym->data;
1877   if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);}
1878   if (label) {
1879     sl->label = label;
1880     ierr = PetscObjectReference((PetscObject) label);CHKERRQ(ierr);
1881     ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr);
1882     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);
1883     ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr);
1884     ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr);
1885     ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr);
1886     ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr);
1887     ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr);
1888   }
1889   PetscFunctionReturn(0);
1890 }
1891 
1892 /*@C
1893   PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
1894 
1895   Logically collective on PetscSectionSym
1896 
1897   InputParameters:
1898 + sys       - the section symmetries
1899 . stratum   - the stratum value in the label that we are assigning symmetries for
1900 . size      - the number of dofs for points in the stratum of the label
1901 . minOrient - the smallest orientation for a point in this stratum
1902 . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
1903 . mode      - how sym should copy the perms and rots arrays
1904 . perms     - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation.  A NULL permutation is the identity
1905 + 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
1906 
1907   Level: developer
1908 
1909 .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()
1910 @*/
1911 PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
1912 {
1913   PetscSectionSym_Label *sl;
1914   const char            *name;
1915   PetscInt               i, j, k;
1916   PetscErrorCode         ierr;
1917 
1918   PetscFunctionBegin;
1919   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
1920   sl = (PetscSectionSym_Label *) sym->data;
1921   if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet");
1922   for (i = 0; i <= sl->numStrata; i++) {
1923     PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
1924 
1925     if (stratum == value) break;
1926   }
1927   ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1928   if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,name);
1929   sl->sizes[i] = size;
1930   sl->modes[i] = mode;
1931   sl->minMaxOrients[i][0] = minOrient;
1932   sl->minMaxOrients[i][1] = maxOrient;
1933   if (mode == PETSC_COPY_VALUES) {
1934     if (perms) {
1935       PetscInt    **ownPerms;
1936 
1937       ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr);
1938       for (j = 0; j < maxOrient-minOrient; j++) {
1939         if (perms[j]) {
1940           ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr);
1941           for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];}
1942         }
1943       }
1944       sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient];
1945     }
1946     if (rots) {
1947       PetscScalar **ownRots;
1948 
1949       ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr);
1950       for (j = 0; j < maxOrient-minOrient; j++) {
1951         if (rots[j]) {
1952           ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr);
1953           for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];}
1954         }
1955       }
1956       sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient];
1957     }
1958   } else {
1959     sl->perms[i] = perms ? &perms[-minOrient] : NULL;
1960     sl->rots[i]  = rots ? &rots[-minOrient] : NULL;
1961   }
1962   PetscFunctionReturn(0);
1963 }
1964 
1965 static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
1966 {
1967   PetscInt              i, j, numStrata;
1968   PetscSectionSym_Label *sl;
1969   DMLabel               label;
1970   PetscErrorCode        ierr;
1971 
1972   PetscFunctionBegin;
1973   sl = (PetscSectionSym_Label *) sym->data;
1974   numStrata = sl->numStrata;
1975   label     = sl->label;
1976   for (i = 0; i < numPoints; i++) {
1977     PetscInt point = points[2*i];
1978     PetscInt ornt  = points[2*i+1];
1979 
1980     for (j = 0; j < numStrata; j++) {
1981       if (label->validIS[j]) {
1982         PetscInt k;
1983 
1984         ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr);
1985         if (k >= 0) break;
1986       } else {
1987         PetscBool has;
1988 
1989         ierr = PetscHSetIHas(label->ht[j], point, &has);CHKERRQ(ierr);
1990         if (has) break;
1991       }
1992     }
1993     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);
1994     if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;}
1995     if (rots) {rots[i]  = sl->rots[j] ? sl->rots[j][ornt] : NULL;}
1996   }
1997   PetscFunctionReturn(0);
1998 }
1999 
2000 PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
2001 {
2002   PetscSectionSym_Label *sl;
2003   PetscErrorCode        ierr;
2004 
2005   PetscFunctionBegin;
2006   ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr);
2007   sym->ops->getpoints = PetscSectionSymGetPoints_Label;
2008   sym->ops->view      = PetscSectionSymView_Label;
2009   sym->ops->destroy   = PetscSectionSymDestroy_Label;
2010   sym->data           = (void *) sl;
2011   PetscFunctionReturn(0);
2012 }
2013 
2014 /*@
2015   PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
2016 
2017   Collective
2018 
2019   Input Parameters:
2020 + comm - the MPI communicator for the new symmetry
2021 - label - the label defining the strata
2022 
2023   Output Parameters:
2024 . sym - the section symmetries
2025 
2026   Level: developer
2027 
2028 .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()
2029 @*/
2030 PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
2031 {
2032   PetscErrorCode        ierr;
2033 
2034   PetscFunctionBegin;
2035   ierr = DMInitializePackage();CHKERRQ(ierr);
2036   ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr);
2037   ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr);
2038   ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr);
2039   PetscFunctionReturn(0);
2040 }
2041