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