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