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