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