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