xref: /petsc/src/dm/label/dmlabel.c (revision 8fcddce65efd55a8fe3f87d4c08c15577ce4cbef)
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 /*@
387   DMLabelComputeIndex - Create an index structure for membership determination, automatically determining the bounds
388 
389   Input Parameter:
390 . label  - The DMLabel
391 
392   Level: intermediate
393 
394 .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
395 @*/
396 PetscErrorCode DMLabelComputeIndex(DMLabel label)
397 {
398   PetscInt       pStart = PETSC_MAX_INT, pEnd = -1, v;
399   PetscErrorCode ierr;
400 
401   PetscFunctionBegin;
402   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
403   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
404   for (v = 0; v < label->numStrata; ++v) {
405     const PetscInt *points;
406     PetscInt       i;
407 
408     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
409     for (i = 0; i < label->stratumSizes[v]; ++i) {
410       const PetscInt point = points[i];
411 
412       pStart = PetscMin(point,   pStart);
413       pEnd   = PetscMax(point+1, pEnd);
414     }
415     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
416   }
417   label->pStart = pStart == PETSC_MAX_INT ? -1 : pStart;
418   label->pEnd   = pEnd;
419   ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
420   PetscFunctionReturn(0);
421 }
422 
423 /*@
424   DMLabelCreateIndex - Create an index structure for membership determination
425 
426   Input Parameters:
427 + label  - The DMLabel
428 . pStart - The smallest point
429 - pEnd   - The largest point + 1
430 
431   Level: intermediate
432 
433 .seealso: DMLabelHasPoint(), DMLabelComputeIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
434 @*/
435 PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
436 {
437   PetscInt       v;
438   PetscErrorCode ierr;
439 
440   PetscFunctionBegin;
441   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
442   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
443   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
444   label->pStart = pStart;
445   label->pEnd   = pEnd;
446   /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
447   ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr);
448   for (v = 0; v < label->numStrata; ++v) {
449     const PetscInt *points;
450     PetscInt       i;
451 
452     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
453     for (i = 0; i < label->stratumSizes[v]; ++i) {
454       const PetscInt point = points[i];
455 
456       if ((point < pStart) || (point >= pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, pStart, pEnd);
457       ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr);
458     }
459     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
460   }
461   PetscFunctionReturn(0);
462 }
463 
464 /*@
465   DMLabelDestroyIndex - Destroy the index structure
466 
467   Input Parameter:
468 . label - the DMLabel
469 
470   Level: intermediate
471 
472 .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
473 @*/
474 PetscErrorCode DMLabelDestroyIndex(DMLabel label)
475 {
476   PetscErrorCode ierr;
477 
478   PetscFunctionBegin;
479   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
480   label->pStart = -1;
481   label->pEnd   = -1;
482   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
483   PetscFunctionReturn(0);
484 }
485 
486 /*@
487   DMLabelGetBounds - Return the smallest and largest point in the label
488 
489   Input Parameter:
490 . label - the DMLabel
491 
492   Output Parameters:
493 + pStart - The smallest point
494 - pEnd   - The largest point + 1
495 
496   Note: This will compute an index for the label if one does not exist.
497 
498   Level: intermediate
499 
500 .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
501 @*/
502 PetscErrorCode DMLabelGetBounds(DMLabel label, PetscInt *pStart, PetscInt *pEnd)
503 {
504   PetscErrorCode ierr;
505 
506   PetscFunctionBegin;
507   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
508   if ((label->pStart == -1) && (label->pEnd == -1)) {ierr = DMLabelComputeIndex(label);CHKERRQ(ierr);}
509   if (pStart) {
510     PetscValidPointer(pStart, 2);
511     *pStart = label->pStart;
512   }
513   if (pEnd) {
514     PetscValidPointer(pEnd, 3);
515     *pEnd = label->pEnd;
516   }
517   PetscFunctionReturn(0);
518 }
519 
520 /*@
521   DMLabelHasValue - Determine whether a label assigns the value to any point
522 
523   Input Parameters:
524 + label - the DMLabel
525 - value - the value
526 
527   Output Parameter:
528 . contains - Flag indicating whether the label maps this value to any point
529 
530   Level: developer
531 
532 .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
533 @*/
534 PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
535 {
536   PetscInt v;
537   PetscErrorCode ierr;
538 
539   PetscFunctionBegin;
540   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
541   PetscValidPointer(contains, 3);
542   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
543   *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE;
544   PetscFunctionReturn(0);
545 }
546 
547 /*@
548   DMLabelHasPoint - Determine whether a label assigns a value to a point
549 
550   Input Parameters:
551 + label - the DMLabel
552 - point - the point
553 
554   Output Parameter:
555 . contains - Flag indicating whether the label maps this point to a value
556 
557   Note: The user must call DMLabelCreateIndex() before this function.
558 
559   Level: developer
560 
561 .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
562 @*/
563 PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
564 {
565   PetscErrorCode ierr;
566 
567   PetscFunctionBeginHot;
568   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
569   PetscValidPointer(contains, 3);
570   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
571 #if defined(PETSC_USE_DEBUG)
572   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
573   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);
574 #endif
575   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
576   PetscFunctionReturn(0);
577 }
578 
579 /*@
580   DMLabelStratumHasPoint - Return true if the stratum contains a point
581 
582   Input Parameters:
583 + label - the DMLabel
584 . value - the stratum value
585 - point - the point
586 
587   Output Parameter:
588 . contains - true if the stratum contains the point
589 
590   Level: intermediate
591 
592 .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
593 @*/
594 PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
595 {
596   PetscInt       v;
597   PetscErrorCode ierr;
598 
599   PetscFunctionBeginHot;
600   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
601   PetscValidPointer(contains, 4);
602   *contains = PETSC_FALSE;
603   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
604   if (v < 0) PetscFunctionReturn(0);
605 
606   if (label->validIS[v]) {
607     PetscInt i;
608 
609     ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
610     if (i >= 0) *contains = PETSC_TRUE;
611   } else {
612     PetscBool has;
613 
614     PetscHSetIHas(label->ht[v], point, &has);
615     if (has) *contains = PETSC_TRUE;
616   }
617   PetscFunctionReturn(0);
618 }
619 
620 /*@
621   DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
622   When a label is created, it is initialized to -1.
623 
624   Input parameter:
625 . label - a DMLabel object
626 
627   Output parameter:
628 . defaultValue - the default value
629 
630   Level: beginner
631 
632 .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
633 @*/
634 PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
635 {
636   PetscFunctionBegin;
637   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
638   *defaultValue = label->defaultValue;
639   PetscFunctionReturn(0);
640 }
641 
642 /*@
643   DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
644   When a label is created, it is initialized to -1.
645 
646   Input parameter:
647 . label - a DMLabel object
648 
649   Output parameter:
650 . defaultValue - the default value
651 
652   Level: beginner
653 
654 .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
655 @*/
656 PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
657 {
658   PetscFunctionBegin;
659   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
660   label->defaultValue = defaultValue;
661   PetscFunctionReturn(0);
662 }
663 
664 /*@
665   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())
666 
667   Input Parameters:
668 + label - the DMLabel
669 - point - the point
670 
671   Output Parameter:
672 . value - The point value, or the default value (-1 by default)
673 
674   Level: intermediate
675 
676 .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
677 @*/
678 PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
679 {
680   PetscInt       v;
681   PetscErrorCode ierr;
682 
683   PetscFunctionBeginHot;
684   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
685   PetscValidPointer(value, 3);
686   *value = label->defaultValue;
687   for (v = 0; v < label->numStrata; ++v) {
688     if (label->validIS[v]) {
689       PetscInt i;
690 
691       ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
692       if (i >= 0) {
693         *value = label->stratumValues[v];
694         break;
695       }
696     } else {
697       PetscBool has;
698 
699       PetscHSetIHas(label->ht[v], point, &has);
700       if (has) {
701         *value = label->stratumValues[v];
702         break;
703       }
704     }
705   }
706   PetscFunctionReturn(0);
707 }
708 
709 /*@
710   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.
711 
712   Input Parameters:
713 + label - the DMLabel
714 . point - the point
715 - value - The point value
716 
717   Level: intermediate
718 
719 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
720 @*/
721 PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
722 {
723   PetscInt       v;
724   PetscErrorCode ierr;
725 
726   PetscFunctionBegin;
727   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
728   /* Find label value, add new entry if needed */
729   if (value == label->defaultValue) PetscFunctionReturn(0);
730   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
731   if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);}
732   /* Set key */
733   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
734   ierr = PetscHSetIAdd(label->ht[v], point);CHKERRQ(ierr);
735   PetscFunctionReturn(0);
736 }
737 
738 /*@
739   DMLabelClearValue - Clear the value a label assigns to a point
740 
741   Input Parameters:
742 + label - the DMLabel
743 . point - the point
744 - value - The point value
745 
746   Level: intermediate
747 
748 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
749 @*/
750 PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
751 {
752   PetscInt       v;
753   PetscErrorCode ierr;
754 
755   PetscFunctionBegin;
756   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
757   /* Find label value */
758   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
759   if (v < 0) PetscFunctionReturn(0);
760 
761   if (label->bt) {
762     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);
763     ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
764   }
765 
766   /* Delete key */
767   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
768   ierr = PetscHSetIDel(label->ht[v], point);CHKERRQ(ierr);
769   PetscFunctionReturn(0);
770 }
771 
772 /*@
773   DMLabelInsertIS - Set all points in the IS to a value
774 
775   Input Parameters:
776 + label - the DMLabel
777 . is    - the point IS
778 - value - The point value
779 
780   Level: intermediate
781 
782 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
783 @*/
784 PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
785 {
786   PetscInt        v, n, p;
787   const PetscInt *points;
788   PetscErrorCode  ierr;
789 
790   PetscFunctionBegin;
791   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
792   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
793   /* Find label value, add new entry if needed */
794   if (value == label->defaultValue) PetscFunctionReturn(0);
795   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
796   if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);}
797   /* Set keys */
798   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
799   ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr);
800   ierr = ISGetIndices(is, &points);CHKERRQ(ierr);
801   for (p = 0; p < n; ++p) {ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);}
802   ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr);
803   PetscFunctionReturn(0);
804 }
805 
806 /*@
807   DMLabelGetNumValues - Get the number of values that the DMLabel takes
808 
809   Input Parameter:
810 . label - the DMLabel
811 
812   Output Paramater:
813 . numValues - the number of values
814 
815   Level: intermediate
816 
817 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
818 @*/
819 PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
820 {
821   PetscFunctionBegin;
822   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
823   PetscValidPointer(numValues, 2);
824   *numValues = label->numStrata;
825   PetscFunctionReturn(0);
826 }
827 
828 /*@
829   DMLabelGetValueIS - Get an IS of all values that the DMlabel takes
830 
831   Input Parameter:
832 . label - the DMLabel
833 
834   Output Paramater:
835 . is    - the value IS
836 
837   Level: intermediate
838 
839 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
840 @*/
841 PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
842 {
843   PetscErrorCode ierr;
844 
845   PetscFunctionBegin;
846   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
847   PetscValidPointer(values, 2);
848   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr);
849   PetscFunctionReturn(0);
850 }
851 
852 /*@
853   DMLabelHasStratum - Determine whether points exist with the given value
854 
855   Input Parameters:
856 + label - the DMLabel
857 - value - the stratum value
858 
859   Output Paramater:
860 . exists - Flag saying whether points exist
861 
862   Level: intermediate
863 
864 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
865 @*/
866 PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
867 {
868   PetscInt       v;
869   PetscErrorCode ierr;
870 
871   PetscFunctionBegin;
872   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
873   PetscValidPointer(exists, 3);
874   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
875   *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE;
876   PetscFunctionReturn(0);
877 }
878 
879 /*@
880   DMLabelGetStratumSize - Get the size of a stratum
881 
882   Input Parameters:
883 + label - the DMLabel
884 - value - the stratum value
885 
886   Output Paramater:
887 . size - The number of points in the stratum
888 
889   Level: intermediate
890 
891 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
892 @*/
893 PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
894 {
895   PetscInt       v;
896   PetscErrorCode ierr;
897 
898   PetscFunctionBegin;
899   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
900   PetscValidPointer(size, 3);
901   *size = 0;
902   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
903   if (v < 0) PetscFunctionReturn(0);
904   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
905   *size = label->stratumSizes[v];
906   PetscFunctionReturn(0);
907 }
908 
909 /*@
910   DMLabelGetStratumBounds - Get the largest and smallest point of a stratum
911 
912   Input Parameters:
913 + label - the DMLabel
914 - value - the stratum value
915 
916   Output Paramaters:
917 + start - the smallest point in the stratum
918 - end - the largest point in the stratum
919 
920   Level: intermediate
921 
922 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
923 @*/
924 PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
925 {
926   PetscInt       v, min, max;
927   PetscErrorCode ierr;
928 
929   PetscFunctionBegin;
930   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
931   if (start) {PetscValidPointer(start, 3); *start = 0;}
932   if (end)   {PetscValidPointer(end,   4); *end   = 0;}
933   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
934   if (v < 0) PetscFunctionReturn(0);
935   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
936   if (label->stratumSizes[v] <= 0) PetscFunctionReturn(0);
937   ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr);
938   if (start) *start = min;
939   if (end)   *end   = max+1;
940   PetscFunctionReturn(0);
941 }
942 
943 /*@
944   DMLabelGetStratumIS - Get an IS with the stratum points
945 
946   Input Parameters:
947 + label - the DMLabel
948 - value - the stratum value
949 
950   Output Paramater:
951 . points - The stratum points
952 
953   Level: intermediate
954 
955 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
956 @*/
957 PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
958 {
959   PetscInt       v;
960   PetscErrorCode ierr;
961 
962   PetscFunctionBegin;
963   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
964   PetscValidPointer(points, 3);
965   *points = NULL;
966   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
967   if (v < 0) PetscFunctionReturn(0);
968   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
969   ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr);
970   *points = label->points[v];
971   PetscFunctionReturn(0);
972 }
973 
974 /*@
975   DMLabelSetStratumIS - Set the stratum points using an IS
976 
977   Input Parameters:
978 + label - the DMLabel
979 . value - the stratum value
980 - points - The stratum points
981 
982   Level: intermediate
983 
984 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
985 @*/
986 PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
987 {
988   PetscInt       v;
989   PetscErrorCode ierr;
990 
991   PetscFunctionBegin;
992   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
993   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
994   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
995   if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);}
996   if (is == label->points[v]) PetscFunctionReturn(0);
997   ierr = DMLabelClearStratum(label, value);CHKERRQ(ierr);
998   ierr = ISGetLocalSize(is, &(label->stratumSizes[v]));CHKERRQ(ierr);
999   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1000   ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
1001   label->points[v] = is;
1002   label->validIS[v] = PETSC_TRUE;
1003   if (label->bt) {
1004     const PetscInt *points;
1005     PetscInt p;
1006 
1007     ierr = ISGetIndices(is,&points);CHKERRQ(ierr);
1008     for (p = 0; p < label->stratumSizes[v]; ++p) {
1009       const PetscInt point = points[p];
1010 
1011       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);
1012       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
1013     }
1014   }
1015   PetscFunctionReturn(0);
1016 }
1017 
1018 /*@
1019   DMLabelClearStratum - Remove a stratum
1020 
1021   Input Parameters:
1022 + label - the DMLabel
1023 - value - the stratum value
1024 
1025   Level: intermediate
1026 
1027 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1028 @*/
1029 PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
1030 {
1031   PetscInt       v;
1032   PetscErrorCode ierr;
1033 
1034   PetscFunctionBegin;
1035   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1036   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
1037   if (v < 0) PetscFunctionReturn(0);
1038   if (label->validIS[v]) {
1039     if (label->bt) {
1040       PetscInt       i;
1041       const PetscInt *points;
1042 
1043       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
1044       for (i = 0; i < label->stratumSizes[v]; ++i) {
1045         const PetscInt point = points[i];
1046 
1047         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);
1048         ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
1049       }
1050       ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
1051     }
1052     label->stratumSizes[v] = 0;
1053     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
1054     ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr);
1055     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
1056   } else {
1057     PetscHSetIClear(label->ht[v]);
1058   }
1059   PetscFunctionReturn(0);
1060 }
1061 
1062 /*@
1063   DMLabelFilter - Remove all points outside of [start, end)
1064 
1065   Input Parameters:
1066 + label - the DMLabel
1067 . start - the first point
1068 - end - the last point
1069 
1070   Level: intermediate
1071 
1072 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1073 @*/
1074 PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
1075 {
1076   PetscInt       v;
1077   PetscErrorCode ierr;
1078 
1079   PetscFunctionBegin;
1080   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1081   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
1082   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1083   for (v = 0; v < label->numStrata; ++v) {
1084     PetscInt off, q;
1085     const PetscInt *points;
1086     PetscInt numPointsNew = 0, *pointsNew = NULL;
1087 
1088     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
1089     for (q = 0; q < label->stratumSizes[v]; ++q)
1090       if (points[q] >= start && points[q] < end)
1091         numPointsNew++;
1092     ierr = PetscMalloc1(numPointsNew, &pointsNew);CHKERRQ(ierr);
1093     for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) {
1094       if (points[q] >= start && points[q] < end)
1095         pointsNew[off++] = points[q];
1096     }
1097     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
1098 
1099     label->stratumSizes[v] = numPointsNew;
1100     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
1101     ierr = ISCreateGeneral(PETSC_COMM_SELF,numPointsNew, pointsNew, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr);
1102     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
1103   }
1104   ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr);
1105   PetscFunctionReturn(0);
1106 }
1107 
1108 /*@
1109   DMLabelPermute - Create a new label with permuted points
1110 
1111   Input Parameters:
1112 + label - the DMLabel
1113 - permutation - the point permutation
1114 
1115   Output Parameter:
1116 . labelnew - the new label containing the permuted points
1117 
1118   Level: intermediate
1119 
1120 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1121 @*/
1122 PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
1123 {
1124   const PetscInt *perm;
1125   PetscInt        numValues, numPoints, v, q;
1126   PetscErrorCode  ierr;
1127 
1128   PetscFunctionBegin;
1129   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1130   PetscValidHeaderSpecific(permutation, IS_CLASSID, 2);
1131   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1132   ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr);
1133   ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr);
1134   ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr);
1135   ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr);
1136   for (v = 0; v < numValues; ++v) {
1137     const PetscInt size   = (*labelNew)->stratumSizes[v];
1138     const PetscInt *points;
1139     PetscInt *pointsNew;
1140 
1141     ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1142     ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr);
1143     for (q = 0; q < size; ++q) {
1144       const PetscInt point = points[q];
1145 
1146       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);
1147       pointsNew[q] = perm[point];
1148     }
1149     ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1150     ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr);
1151     ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr);
1152     if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) {
1153       ierr = ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));CHKERRQ(ierr);
1154       ierr = PetscFree(pointsNew);CHKERRQ(ierr);
1155     } else {
1156       ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr);
1157     }
1158     ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr);
1159   }
1160   ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr);
1161   if (label->bt) {
1162     ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
1163     ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
1164   }
1165   PetscFunctionReturn(0);
1166 }
1167 
1168 PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
1169 {
1170   MPI_Comm       comm;
1171   PetscInt       s, l, nroots, nleaves, dof, offset, size;
1172   PetscInt      *remoteOffsets, *rootStrata, *rootIdx;
1173   PetscSection   rootSection;
1174   PetscSF        labelSF;
1175   PetscErrorCode ierr;
1176 
1177   PetscFunctionBegin;
1178   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
1179   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1180   /* Build a section of stratum values per point, generate the according SF
1181      and distribute point-wise stratum values to leaves. */
1182   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
1183   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
1184   ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr);
1185   if (label) {
1186     for (s = 0; s < label->numStrata; ++s) {
1187       const PetscInt *points;
1188 
1189       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
1190       for (l = 0; l < label->stratumSizes[s]; l++) {
1191         ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr);
1192         ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr);
1193       }
1194       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
1195     }
1196   }
1197   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
1198   /* Create a point-wise array of stratum values */
1199   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
1200   ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr);
1201   ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr);
1202   if (label) {
1203     for (s = 0; s < label->numStrata; ++s) {
1204       const PetscInt *points;
1205 
1206       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
1207       for (l = 0; l < label->stratumSizes[s]; l++) {
1208         const PetscInt p = points[l];
1209         ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr);
1210         rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
1211       }
1212       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
1213     }
1214   }
1215   /* Build SF that maps label points to remote processes */
1216   ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr);
1217   ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr);
1218   ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr);
1219   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
1220   /* Send the strata for each point over the derived SF */
1221   ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr);
1222   ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr);
1223   ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
1224   ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
1225   /* Clean up */
1226   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
1227   ierr = PetscFree(rootIdx);CHKERRQ(ierr);
1228   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
1229   ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr);
1230   PetscFunctionReturn(0);
1231 }
1232 
1233 /*@
1234   DMLabelDistribute - Create a new label pushed forward over the PetscSF
1235 
1236   Input Parameters:
1237 + label - the DMLabel
1238 - sf    - the map from old to new distribution
1239 
1240   Output Parameter:
1241 . labelnew - the new redistributed label
1242 
1243   Level: intermediate
1244 
1245 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
1246 @*/
1247 PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1248 {
1249   MPI_Comm       comm;
1250   PetscSection   leafSection;
1251   PetscInt       p, pStart, pEnd, s, size, dof, offset, stratum;
1252   PetscInt      *leafStrata, *strataIdx;
1253   PetscInt     **points;
1254   const char    *lname = NULL;
1255   char          *name;
1256   PetscInt       nameSize;
1257   PetscHSetI     stratumHash;
1258   size_t         len = 0;
1259   PetscMPIInt    rank;
1260   PetscErrorCode ierr;
1261 
1262   PetscFunctionBegin;
1263   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1264   if (label) {
1265     PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1266     ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1267   }
1268   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1269   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1270   /* Bcast name */
1271   if (!rank) {
1272     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1273     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1274   }
1275   nameSize = len;
1276   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1277   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1278   if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);}
1279   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1280   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
1281   ierr = PetscFree(name);CHKERRQ(ierr);
1282   /* Bcast defaultValue */
1283   if (!rank) (*labelNew)->defaultValue = label->defaultValue;
1284   ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1285   /* Distribute stratum values over the SF and get the point mapping on the receiver */
1286   ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr);
1287   /* Determine received stratum values and initialise new label*/
1288   ierr = PetscHSetICreate(&stratumHash);CHKERRQ(ierr);
1289   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
1290   for (p = 0; p < size; ++p) {ierr = PetscHSetIAdd(stratumHash, leafStrata[p]);CHKERRQ(ierr);}
1291   ierr = PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata);CHKERRQ(ierr);
1292   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr);
1293   for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
1294   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr);
1295   /* Turn leafStrata into indices rather than stratum values */
1296   offset = 0;
1297   ierr = PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr);
1298   ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr);
1299   for (p = 0; p < size; ++p) {
1300     for (s = 0; s < (*labelNew)->numStrata; ++s) {
1301       if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
1302     }
1303   }
1304   /* Rebuild the point strata on the receiver */
1305   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr);
1306   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1307   for (p=pStart; p<pEnd; p++) {
1308     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1309     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1310     for (s=0; s<dof; s++) {
1311       (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1312     }
1313   }
1314   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr);
1315   ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr);
1316   ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr);
1317   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1318     ierr = PetscHSetICreate(&(*labelNew)->ht[s]);CHKERRQ(ierr);
1319     ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr);
1320   }
1321   /* Insert points into new strata */
1322   ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr);
1323   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1324   for (p=pStart; p<pEnd; p++) {
1325     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1326     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1327     for (s=0; s<dof; s++) {
1328       stratum = leafStrata[offset+s];
1329       points[stratum][strataIdx[stratum]++] = p;
1330     }
1331   }
1332   for (s = 0; s < (*labelNew)->numStrata; s++) {
1333     ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr);
1334     ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr);
1335   }
1336   ierr = PetscFree(points);CHKERRQ(ierr);
1337   ierr = PetscHSetIDestroy(&stratumHash);CHKERRQ(ierr);
1338   ierr = PetscFree(leafStrata);CHKERRQ(ierr);
1339   ierr = PetscFree(strataIdx);CHKERRQ(ierr);
1340   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
1341   PetscFunctionReturn(0);
1342 }
1343 
1344 /*@
1345   DMLabelGather - Gather all label values from leafs into roots
1346 
1347   Input Parameters:
1348 + label - the DMLabel
1349 - sf - the Star Forest point communication map
1350 
1351   Output Parameters:
1352 . labelNew - the new DMLabel with localised leaf values
1353 
1354   Level: developer
1355 
1356   Note: This is the inverse operation to DMLabelDistribute.
1357 
1358 .seealso: DMLabelDistribute()
1359 @*/
1360 PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
1361 {
1362   MPI_Comm       comm;
1363   PetscSection   rootSection;
1364   PetscSF        sfLabel;
1365   PetscSFNode   *rootPoints, *leafPoints;
1366   PetscInt       p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
1367   const PetscInt *rootDegree, *ilocal;
1368   PetscInt       *rootStrata;
1369   const char    *lname;
1370   char          *name;
1371   PetscInt       nameSize;
1372   size_t         len = 0;
1373   PetscMPIInt    rank, size;
1374   PetscErrorCode ierr;
1375 
1376   PetscFunctionBegin;
1377   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1378   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1379   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1380   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1381   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1382   /* Bcast name */
1383   if (!rank) {
1384     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1385     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1386   }
1387   nameSize = len;
1388   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1389   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1390   if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);}
1391   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1392   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
1393   ierr = PetscFree(name);CHKERRQ(ierr);
1394   /* Gather rank/index pairs of leaves into local roots to build
1395      an inverse, multi-rooted SF. Note that this ignores local leaf
1396      indexing due to the use of the multiSF in PetscSFGather. */
1397   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr);
1398   ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr);
1399   for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
1400   for (p = 0; p < nleaves; p++) {
1401     PetscInt ilp = ilocal ? ilocal[p] : p;
1402 
1403     leafPoints[ilp].index = ilp;
1404     leafPoints[ilp].rank  = rank;
1405   }
1406   ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr);
1407   ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr);
1408   for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
1409   ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr);
1410   ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
1411   ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
1412   ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr);
1413   ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
1414   /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
1415   ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr);
1416   /* Rebuild the point strata on the receiver */
1417   for (p = 0, idx = 0; p < nroots; p++) {
1418     for (d = 0; d < rootDegree[p]; d++) {
1419       ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr);
1420       ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr);
1421       for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);}
1422     }
1423     idx += rootDegree[p];
1424   }
1425   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
1426   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
1427   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
1428   ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr);
1429   PetscFunctionReturn(0);
1430 }
1431 
1432 /*@
1433   DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label
1434 
1435   Input Parameter:
1436 . label - the DMLabel
1437 
1438   Output Parameters:
1439 + section - the section giving offsets for each stratum
1440 - is - An IS containing all the label points
1441 
1442   Level: developer
1443 
1444 .seealso: DMLabelDistribute()
1445 @*/
1446 PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1447 {
1448   IS              vIS;
1449   const PetscInt *values;
1450   PetscInt       *points;
1451   PetscInt        nV, vS = 0, vE = 0, v, N;
1452   PetscErrorCode  ierr;
1453 
1454   PetscFunctionBegin;
1455   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1456   ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr);
1457   ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
1458   ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr);
1459   if (nV) {vS = values[0]; vE = values[0]+1;}
1460   for (v = 1; v < nV; ++v) {
1461     vS = PetscMin(vS, values[v]);
1462     vE = PetscMax(vE, values[v]+1);
1463   }
1464   ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr);
1465   ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr);
1466   for (v = 0; v < nV; ++v) {
1467     PetscInt n;
1468 
1469     ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr);
1470     ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr);
1471   }
1472   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
1473   ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr);
1474   ierr = PetscMalloc1(N, &points);CHKERRQ(ierr);
1475   for (v = 0; v < nV; ++v) {
1476     IS              is;
1477     const PetscInt *spoints;
1478     PetscInt        dof, off, p;
1479 
1480     ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr);
1481     ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr);
1482     ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr);
1483     ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr);
1484     for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1485     ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr);
1486     ierr = ISDestroy(&is);CHKERRQ(ierr);
1487   }
1488   ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr);
1489   ierr = ISDestroy(&vIS);CHKERRQ(ierr);
1490   ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1491   PetscFunctionReturn(0);
1492 }
1493 
1494 /*@
1495   PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1496   the local section and an SF describing the section point overlap.
1497 
1498   Input Parameters:
1499   + s - The PetscSection for the local field layout
1500   . sf - The SF describing parallel layout of the section points
1501   . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1502   . label - The label specifying the points
1503   - labelValue - The label stratum specifying the points
1504 
1505   Output Parameter:
1506   . gsection - The PetscSection for the global field layout
1507 
1508   Note: This gives negative sizes and offsets to points not owned by this process
1509 
1510   Level: developer
1511 
1512 .seealso: PetscSectionCreate()
1513 @*/
1514 PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1515 {
1516   PetscInt      *neg = NULL, *tmpOff = NULL;
1517   PetscInt       pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1518   PetscErrorCode ierr;
1519 
1520   PetscFunctionBegin;
1521   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1522   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1523   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4);
1524   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr);
1525   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1526   ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr);
1527   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1528   if (nroots >= 0) {
1529     if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1530     ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr);
1531     if (nroots > pEnd-pStart) {
1532       ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr);
1533     } else {
1534       tmpOff = &(*gsection)->atlasDof[-pStart];
1535     }
1536   }
1537   /* Mark ghost points with negative dof */
1538   for (p = pStart; p < pEnd; ++p) {
1539     PetscInt value;
1540 
1541     ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr);
1542     if (value != labelValue) continue;
1543     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1544     ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr);
1545     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1546     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);}
1547     if (neg) neg[p] = -(dof+1);
1548   }
1549   ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr);
1550   if (nroots >= 0) {
1551     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1552     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1553     if (nroots > pEnd-pStart) {
1554       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1555     }
1556   }
1557   /* Calculate new sizes, get proccess offset, and calculate point offsets */
1558   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1559     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1560     (*gsection)->atlasOff[p] = off;
1561     off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1562   }
1563   ierr       = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1564   globalOff -= off;
1565   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1566     (*gsection)->atlasOff[p] += globalOff;
1567     if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1568   }
1569   /* Put in negative offsets for ghost points */
1570   if (nroots >= 0) {
1571     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1572     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1573     if (nroots > pEnd-pStart) {
1574       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1575     }
1576   }
1577   if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);}
1578   ierr = PetscFree(neg);CHKERRQ(ierr);
1579   PetscFunctionReturn(0);
1580 }
1581 
1582 typedef struct _n_PetscSectionSym_Label
1583 {
1584   DMLabel           label;
1585   PetscCopyMode     *modes;
1586   PetscInt          *sizes;
1587   const PetscInt    ***perms;
1588   const PetscScalar ***rots;
1589   PetscInt          (*minMaxOrients)[2];
1590   PetscInt          numStrata; /* numStrata is only increasing, functions as a state */
1591 } PetscSectionSym_Label;
1592 
1593 static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
1594 {
1595   PetscInt              i, j;
1596   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
1597   PetscErrorCode        ierr;
1598 
1599   PetscFunctionBegin;
1600   for (i = 0; i <= sl->numStrata; i++) {
1601     if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
1602       for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
1603         if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);}
1604         if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);}
1605       }
1606       if (sl->perms[i]) {
1607         const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
1608 
1609         ierr = PetscFree(perms);CHKERRQ(ierr);
1610       }
1611       if (sl->rots[i]) {
1612         const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
1613 
1614         ierr = PetscFree(rots);CHKERRQ(ierr);
1615       }
1616     }
1617   }
1618   ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr);
1619   ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr);
1620   sl->numStrata = 0;
1621   PetscFunctionReturn(0);
1622 }
1623 
1624 static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
1625 {
1626   PetscErrorCode ierr;
1627 
1628   PetscFunctionBegin;
1629   ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);
1630   ierr = PetscFree(sym->data);CHKERRQ(ierr);
1631   PetscFunctionReturn(0);
1632 }
1633 
1634 static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
1635 {
1636   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
1637   PetscBool             isAscii;
1638   DMLabel               label = sl->label;
1639   const char           *name;
1640   PetscErrorCode        ierr;
1641 
1642   PetscFunctionBegin;
1643   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
1644   if (isAscii) {
1645     PetscInt          i, j, k;
1646     PetscViewerFormat format;
1647 
1648     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1649     if (label) {
1650       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1651       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1652         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1653         ierr = DMLabelView(label, viewer);CHKERRQ(ierr);
1654         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1655       } else {
1656         ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1657         ierr = PetscViewerASCIIPrintf(viewer,"  Label '%s'\n",name);CHKERRQ(ierr);
1658       }
1659     } else {
1660       ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr);
1661     }
1662     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1663     for (i = 0; i <= sl->numStrata; i++) {
1664       PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
1665 
1666       if (!(sl->perms[i] || sl->rots[i])) {
1667         ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr);
1668       } else {
1669       ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr);
1670         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1671         ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr);
1672         if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1673           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1674           for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
1675             if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
1676               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr);
1677             } else {
1678               PetscInt tab;
1679 
1680               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr);
1681               ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1682               ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr);
1683               if (sl->perms[i] && sl->perms[i][j]) {
1684                 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr);
1685                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
1686                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);}
1687                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1688                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
1689               }
1690               if (sl->rots[i] && sl->rots[i][j]) {
1691                 ierr = PetscViewerASCIIPrintf(viewer,"Rotations:  ");CHKERRQ(ierr);
1692                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
1693 #if defined(PETSC_USE_COMPLEX)
1694                 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);}
1695 #else
1696                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);}
1697 #endif
1698                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1699                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
1700               }
1701               ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1702             }
1703           }
1704           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1705         }
1706         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1707       }
1708     }
1709     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1710   }
1711   PetscFunctionReturn(0);
1712 }
1713 
1714 /*@
1715   PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
1716 
1717   Logically collective on sym
1718 
1719   Input parameters:
1720 + sym - the section symmetries
1721 - label - the DMLabel describing the types of points
1722 
1723   Level: developer:
1724 
1725 .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()
1726 @*/
1727 PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
1728 {
1729   PetscSectionSym_Label *sl;
1730   PetscErrorCode        ierr;
1731 
1732   PetscFunctionBegin;
1733   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
1734   sl = (PetscSectionSym_Label *) sym->data;
1735   if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);}
1736   if (label) {
1737     sl->label = label;
1738     ierr = PetscObjectReference((PetscObject) label);CHKERRQ(ierr);
1739     ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr);
1740     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);
1741     ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr);
1742     ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr);
1743     ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr);
1744     ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr);
1745     ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr);
1746   }
1747   PetscFunctionReturn(0);
1748 }
1749 
1750 /*@C
1751   PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
1752 
1753   Logically collective on PetscSectionSym
1754 
1755   InputParameters:
1756 + sys       - the section symmetries
1757 . stratum   - the stratum value in the label that we are assigning symmetries for
1758 . size      - the number of dofs for points in the stratum of the label
1759 . minOrient - the smallest orientation for a point in this stratum
1760 . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
1761 . mode      - how sym should copy the perms and rots arrays
1762 . perms     - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation.  A NULL permutation is the identity
1763 + 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
1764 
1765   Level: developer
1766 
1767 .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()
1768 @*/
1769 PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
1770 {
1771   PetscSectionSym_Label *sl;
1772   const char            *name;
1773   PetscInt               i, j, k;
1774   PetscErrorCode         ierr;
1775 
1776   PetscFunctionBegin;
1777   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
1778   sl = (PetscSectionSym_Label *) sym->data;
1779   if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet");
1780   for (i = 0; i <= sl->numStrata; i++) {
1781     PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
1782 
1783     if (stratum == value) break;
1784   }
1785   ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1786   if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,name);
1787   sl->sizes[i] = size;
1788   sl->modes[i] = mode;
1789   sl->minMaxOrients[i][0] = minOrient;
1790   sl->minMaxOrients[i][1] = maxOrient;
1791   if (mode == PETSC_COPY_VALUES) {
1792     if (perms) {
1793       PetscInt    **ownPerms;
1794 
1795       ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr);
1796       for (j = 0; j < maxOrient-minOrient; j++) {
1797         if (perms[j]) {
1798           ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr);
1799           for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];}
1800         }
1801       }
1802       sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient];
1803     }
1804     if (rots) {
1805       PetscScalar **ownRots;
1806 
1807       ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr);
1808       for (j = 0; j < maxOrient-minOrient; j++) {
1809         if (rots[j]) {
1810           ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr);
1811           for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];}
1812         }
1813       }
1814       sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient];
1815     }
1816   } else {
1817     sl->perms[i] = perms ? &perms[-minOrient] : NULL;
1818     sl->rots[i]  = rots ? &rots[-minOrient] : NULL;
1819   }
1820   PetscFunctionReturn(0);
1821 }
1822 
1823 static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
1824 {
1825   PetscInt              i, j, numStrata;
1826   PetscSectionSym_Label *sl;
1827   DMLabel               label;
1828   PetscErrorCode        ierr;
1829 
1830   PetscFunctionBegin;
1831   sl = (PetscSectionSym_Label *) sym->data;
1832   numStrata = sl->numStrata;
1833   label     = sl->label;
1834   for (i = 0; i < numPoints; i++) {
1835     PetscInt point = points[2*i];
1836     PetscInt ornt  = points[2*i+1];
1837 
1838     for (j = 0; j < numStrata; j++) {
1839       if (label->validIS[j]) {
1840         PetscInt k;
1841 
1842         ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr);
1843         if (k >= 0) break;
1844       } else {
1845         PetscBool has;
1846 
1847         PetscHSetIHas(label->ht[j], point, &has);
1848         if (has) break;
1849       }
1850     }
1851     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);
1852     if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;}
1853     if (rots) {rots[i]  = sl->rots[j] ? sl->rots[j][ornt] : NULL;}
1854   }
1855   PetscFunctionReturn(0);
1856 }
1857 
1858 PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
1859 {
1860   PetscSectionSym_Label *sl;
1861   PetscErrorCode        ierr;
1862 
1863   PetscFunctionBegin;
1864   ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr);
1865   sym->ops->getpoints = PetscSectionSymGetPoints_Label;
1866   sym->ops->view      = PetscSectionSymView_Label;
1867   sym->ops->destroy   = PetscSectionSymDestroy_Label;
1868   sym->data           = (void *) sl;
1869   PetscFunctionReturn(0);
1870 }
1871 
1872 /*@
1873   PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
1874 
1875   Collective
1876 
1877   Input Parameters:
1878 + comm - the MPI communicator for the new symmetry
1879 - label - the label defining the strata
1880 
1881   Output Parameters:
1882 . sym - the section symmetries
1883 
1884   Level: developer
1885 
1886 .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()
1887 @*/
1888 PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
1889 {
1890   PetscErrorCode        ierr;
1891 
1892   PetscFunctionBegin;
1893   ierr = DMInitializePackage();CHKERRQ(ierr);
1894   ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr);
1895   ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr);
1896   ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr);
1897   PetscFunctionReturn(0);
1898 }
1899