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