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