xref: /petsc/src/dm/label/dmlabel.c (revision 71f1c950b7bbc7c3d25e91463aa371a2ce720a45)
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 PetscFunctionList DMLabelList              = NULL;
8 PetscBool         DMLabelRegisterAllCalled = PETSC_FALSE;
9 
10 /*@
11   DMLabelCreate - Create a `DMLabel` object, which is a multimap
12 
13   Collective
14 
15   Input Parameters:
16 + comm - The communicator, usually `PETSC_COMM_SELF`
17 - name - The label name
18 
19   Output Parameter:
20 . label - The `DMLabel`
21 
22   Level: beginner
23 
24   Notes:
25   The label name is actually usually the `PetscObject` name.
26   One can get/set it with `PetscObjectGetName()`/`PetscObjectSetName()`.
27 
28 .seealso: `DMLabel`, `DM`, `DMLabelDestroy()`
29 @*/
30 PetscErrorCode DMLabelCreate(MPI_Comm comm, const char name[], DMLabel *label)
31 {
32   PetscFunctionBegin;
33   PetscAssertPointer(label, 3);
34   PetscCall(DMInitializePackage());
35 
36   PetscCall(PetscHeaderCreate(*label, DMLABEL_CLASSID, "DMLabel", "DMLabel", "DM", comm, DMLabelDestroy, DMLabelView));
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   PetscCall(PetscHMapICreate(&(*label)->hmap));
48   PetscCall(PetscObjectSetName((PetscObject)*label, name));
49   PetscCall(DMLabelSetType(*label, DMLABELCONCRETE));
50   PetscFunctionReturn(PETSC_SUCCESS);
51 }
52 
53 /*@
54   DMLabelSetUp - SetUp a `DMLabel` object
55 
56   Collective
57 
58   Input Parameters:
59 . label - The `DMLabel`
60 
61   Level: intermediate
62 
63 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelDestroy()`
64 @*/
65 PetscErrorCode DMLabelSetUp(DMLabel label)
66 {
67   PetscFunctionBegin;
68   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
69   PetscTryTypeMethod(label, setup);
70   PetscFunctionReturn(PETSC_SUCCESS);
71 }
72 
73 /*
74   DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format
75 
76   Not collective
77 
78   Input parameter:
79 + label - The `DMLabel`
80 - v - The stratum value
81 
82   Output parameter:
83 . label - The `DMLabel` with stratum in sorted list format
84 
85   Level: developer
86 
87 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`
88 */
89 static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v)
90 {
91   IS       is;
92   PetscInt off = 0, *pointArray, p;
93 
94   PetscFunctionBegin;
95   if ((PetscLikely(v >= 0 && v < label->numStrata) && label->validIS[v]) || label->readonly) PetscFunctionReturn(PETSC_SUCCESS);
96   PetscCheck(v >= 0 && v < label->numStrata, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %" PetscInt_FMT " in DMLabelMakeValid_Private", v);
97   PetscCall(PetscHSetIGetSize(label->ht[v], &label->stratumSizes[v]));
98   PetscCall(PetscMalloc1(label->stratumSizes[v], &pointArray));
99   PetscCall(PetscHSetIGetElems(label->ht[v], &off, pointArray));
100   PetscCall(PetscHSetIClear(label->ht[v]));
101   PetscCall(PetscSortInt(label->stratumSizes[v], pointArray));
102   if (label->bt) {
103     for (p = 0; p < label->stratumSizes[v]; ++p) {
104       const PetscInt point = pointArray[p];
105       PetscCheck(!(point < label->pStart) && !(point >= label->pEnd), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %" PetscInt_FMT " is not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, label->pStart, label->pEnd);
106       PetscCall(PetscBTSet(label->bt, point - label->pStart));
107     }
108   }
109   if (label->stratumSizes[v] > 0 && pointArray[label->stratumSizes[v] - 1] == pointArray[0] + label->stratumSizes[v] - 1) {
110     PetscCall(ISCreateStride(PETSC_COMM_SELF, label->stratumSizes[v], pointArray[0], 1, &is));
111     PetscCall(PetscFree(pointArray));
112   } else {
113     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, label->stratumSizes[v], pointArray, PETSC_OWN_POINTER, &is));
114   }
115   PetscCall(ISSetInfo(is, IS_SORTED, IS_LOCAL, PETSC_TRUE, PETSC_TRUE));
116   PetscCall(PetscObjectSetName((PetscObject)is, "indices"));
117   label->points[v]  = is;
118   label->validIS[v] = PETSC_TRUE;
119   PetscCall(PetscObjectStateIncrease((PetscObject)label));
120   PetscFunctionReturn(PETSC_SUCCESS);
121 }
122 
123 /*
124   DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format
125 
126   Not Collective
127 
128   Input parameter:
129 . label - The `DMLabel`
130 
131   Output parameter:
132 . label - The `DMLabel` with all strata in sorted list format
133 
134   Level: developer
135 
136 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`
137 */
138 static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label)
139 {
140   PetscInt v;
141 
142   PetscFunctionBegin;
143   for (v = 0; v < label->numStrata; v++) PetscCall(DMLabelMakeValid_Private(label, v));
144   PetscFunctionReturn(PETSC_SUCCESS);
145 }
146 
147 /*
148   DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format
149 
150   Not Collective
151 
152   Input parameter:
153 + label - The `DMLabel`
154 - v - The stratum value
155 
156   Output parameter:
157 . label - The `DMLabel` with stratum in hash format
158 
159   Level: developer
160 
161 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`
162 */
163 static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v)
164 {
165   PetscInt        p;
166   const PetscInt *points;
167 
168   PetscFunctionBegin;
169   if ((PetscLikely(v >= 0 && v < label->numStrata) && !label->validIS[v]) || label->readonly) PetscFunctionReturn(PETSC_SUCCESS);
170   PetscCheck(v >= 0 && v < label->numStrata, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %" PetscInt_FMT " in DMLabelMakeInvalid_Private", v);
171   if (label->points[v]) {
172     PetscCall(ISGetIndices(label->points[v], &points));
173     for (p = 0; p < label->stratumSizes[v]; ++p) PetscCall(PetscHSetIAdd(label->ht[v], points[p]));
174     PetscCall(ISRestoreIndices(label->points[v], &points));
175     PetscCall(ISDestroy(&label->points[v]));
176   }
177   label->validIS[v] = PETSC_FALSE;
178   PetscFunctionReturn(PETSC_SUCCESS);
179 }
180 
181 PetscErrorCode DMLabelMakeAllInvalid_Internal(DMLabel label)
182 {
183   PetscInt v;
184 
185   PetscFunctionBegin;
186   for (v = 0; v < label->numStrata; v++) PetscCall(DMLabelMakeInvalid_Private(label, v));
187   PetscFunctionReturn(PETSC_SUCCESS);
188 }
189 
190 #if !defined(DMLABEL_LOOKUP_THRESHOLD)
191   #define DMLABEL_LOOKUP_THRESHOLD 16
192 #endif
193 
194 PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index)
195 {
196   PetscInt v;
197 
198   PetscFunctionBegin;
199   *index = -1;
200   if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD || label->readonly) {
201     for (v = 0; v < label->numStrata; ++v)
202       if (label->stratumValues[v] == value) {
203         *index = v;
204         break;
205       }
206   } else {
207     PetscCall(PetscHMapIGet(label->hmap, value, index));
208   }
209   if (PetscDefined(USE_DEBUG) && !label->readonly) { /* Check strata hash map consistency */
210     PetscInt len, loc = -1;
211     PetscCall(PetscHMapIGetSize(label->hmap, &len));
212     PetscCheck(len == label->numStrata, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map size");
213     if (label->numStrata <= DMLABEL_LOOKUP_THRESHOLD) {
214       PetscCall(PetscHMapIGet(label->hmap, value, &loc));
215     } else {
216       for (v = 0; v < label->numStrata; ++v)
217         if (label->stratumValues[v] == value) {
218           loc = v;
219           break;
220         }
221     }
222     PetscCheck(loc == *index, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent strata hash map lookup");
223   }
224   PetscFunctionReturn(PETSC_SUCCESS);
225 }
226 
227 static inline PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index)
228 {
229   PetscInt    v;
230   PetscInt   *tmpV;
231   PetscInt   *tmpS;
232   PetscHSetI *tmpH, ht;
233   IS         *tmpP, is;
234   PetscBool  *tmpB;
235   PetscHMapI  hmap = label->hmap;
236 
237   PetscFunctionBegin;
238   v    = label->numStrata;
239   tmpV = label->stratumValues;
240   tmpS = label->stratumSizes;
241   tmpH = label->ht;
242   tmpP = label->points;
243   tmpB = label->validIS;
244   { /* TODO: PetscRealloc() is broken, use malloc+memcpy+free  */
245     PetscInt   *oldV = tmpV;
246     PetscInt   *oldS = tmpS;
247     PetscHSetI *oldH = tmpH;
248     IS         *oldP = tmpP;
249     PetscBool  *oldB = tmpB;
250     PetscCall(PetscMalloc((v + 1) * sizeof(*tmpV), &tmpV));
251     PetscCall(PetscMalloc((v + 1) * sizeof(*tmpS), &tmpS));
252     PetscCall(PetscCalloc((v + 1) * sizeof(*tmpH), &tmpH));
253     PetscCall(PetscCalloc((v + 1) * sizeof(*tmpP), &tmpP));
254     PetscCall(PetscMalloc((v + 1) * sizeof(*tmpB), &tmpB));
255     PetscCall(PetscArraycpy(tmpV, oldV, v));
256     PetscCall(PetscArraycpy(tmpS, oldS, v));
257     PetscCall(PetscArraycpy(tmpH, oldH, v));
258     PetscCall(PetscArraycpy(tmpP, oldP, v));
259     PetscCall(PetscArraycpy(tmpB, oldB, v));
260     PetscCall(PetscFree(oldV));
261     PetscCall(PetscFree(oldS));
262     PetscCall(PetscFree(oldH));
263     PetscCall(PetscFree(oldP));
264     PetscCall(PetscFree(oldB));
265   }
266   label->numStrata     = v + 1;
267   label->stratumValues = tmpV;
268   label->stratumSizes  = tmpS;
269   label->ht            = tmpH;
270   label->points        = tmpP;
271   label->validIS       = tmpB;
272   PetscCall(PetscHSetICreate(&ht));
273   PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &is));
274   PetscCall(PetscHMapISet(hmap, value, v));
275   tmpV[v] = value;
276   tmpS[v] = 0;
277   tmpH[v] = ht;
278   tmpP[v] = is;
279   tmpB[v] = PETSC_TRUE;
280   PetscCall(PetscObjectStateIncrease((PetscObject)label));
281   *index = v;
282   PetscFunctionReturn(PETSC_SUCCESS);
283 }
284 
285 static inline PetscErrorCode DMLabelLookupAddStratum(DMLabel label, PetscInt value, PetscInt *index)
286 {
287   PetscFunctionBegin;
288   PetscCall(DMLabelLookupStratum(label, value, index));
289   if (*index < 0) PetscCall(DMLabelNewStratum(label, value, index));
290   PetscFunctionReturn(PETSC_SUCCESS);
291 }
292 
293 PetscErrorCode DMLabelGetStratumSize_Private(DMLabel label, PetscInt v, PetscInt *size)
294 {
295   PetscFunctionBegin;
296   *size = 0;
297   if (v < 0) PetscFunctionReturn(PETSC_SUCCESS);
298   if (label->readonly || label->validIS[v]) {
299     *size = label->stratumSizes[v];
300   } else {
301     PetscCall(PetscHSetIGetSize(label->ht[v], size));
302   }
303   PetscFunctionReturn(PETSC_SUCCESS);
304 }
305 
306 /*@
307   DMLabelAddStratum - Adds a new stratum value in a `DMLabel`
308 
309   Input Parameters:
310 + label - The `DMLabel`
311 - value - The stratum value
312 
313   Level: beginner
314 
315 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelDestroy()`
316 @*/
317 PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
318 {
319   PetscInt v;
320 
321   PetscFunctionBegin;
322   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
323   PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered");
324   PetscCall(DMLabelLookupAddStratum(label, value, &v));
325   PetscFunctionReturn(PETSC_SUCCESS);
326 }
327 
328 /*@
329   DMLabelAddStrata - Adds new stratum values in a `DMLabel`
330 
331   Not Collective
332 
333   Input Parameters:
334 + label         - The `DMLabel`
335 . numStrata     - The number of stratum values
336 - stratumValues - The stratum values
337 
338   Level: beginner
339 
340 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelDestroy()`
341 @*/
342 PetscErrorCode DMLabelAddStrata(DMLabel label, PetscInt numStrata, const PetscInt stratumValues[])
343 {
344   PetscInt *values, v;
345 
346   PetscFunctionBegin;
347   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
348   if (numStrata) PetscAssertPointer(stratumValues, 3);
349   PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered");
350   PetscCall(PetscMalloc1(numStrata, &values));
351   PetscCall(PetscArraycpy(values, stratumValues, numStrata));
352   PetscCall(PetscSortRemoveDupsInt(&numStrata, values));
353   if (!label->numStrata) { /* Fast preallocation */
354     PetscInt   *tmpV;
355     PetscInt   *tmpS;
356     PetscHSetI *tmpH, ht;
357     IS         *tmpP, is;
358     PetscBool  *tmpB;
359     PetscHMapI  hmap = label->hmap;
360 
361     PetscCall(PetscMalloc1(numStrata, &tmpV));
362     PetscCall(PetscMalloc1(numStrata, &tmpS));
363     PetscCall(PetscCalloc1(numStrata, &tmpH));
364     PetscCall(PetscCalloc1(numStrata, &tmpP));
365     PetscCall(PetscMalloc1(numStrata, &tmpB));
366     label->numStrata     = numStrata;
367     label->stratumValues = tmpV;
368     label->stratumSizes  = tmpS;
369     label->ht            = tmpH;
370     label->points        = tmpP;
371     label->validIS       = tmpB;
372     for (v = 0; v < numStrata; ++v) {
373       PetscCall(PetscHSetICreate(&ht));
374       PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &is));
375       PetscCall(PetscHMapISet(hmap, values[v], v));
376       tmpV[v] = values[v];
377       tmpS[v] = 0;
378       tmpH[v] = ht;
379       tmpP[v] = is;
380       tmpB[v] = PETSC_TRUE;
381     }
382     PetscCall(PetscObjectStateIncrease((PetscObject)label));
383   } else {
384     for (v = 0; v < numStrata; ++v) PetscCall(DMLabelAddStratum(label, values[v]));
385   }
386   PetscCall(PetscFree(values));
387   PetscFunctionReturn(PETSC_SUCCESS);
388 }
389 
390 /*@
391   DMLabelAddStrataIS - Adds new stratum values in a `DMLabel`
392 
393   Not Collective
394 
395   Input Parameters:
396 + label   - The `DMLabel`
397 - valueIS - Index set with stratum values
398 
399   Level: beginner
400 
401 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelDestroy()`
402 @*/
403 PetscErrorCode DMLabelAddStrataIS(DMLabel label, IS valueIS)
404 {
405   PetscInt        numStrata;
406   const PetscInt *stratumValues;
407 
408   PetscFunctionBegin;
409   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
410   PetscValidHeaderSpecific(valueIS, IS_CLASSID, 2);
411   PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered");
412   PetscCall(ISGetLocalSize(valueIS, &numStrata));
413   PetscCall(ISGetIndices(valueIS, &stratumValues));
414   PetscCall(DMLabelAddStrata(label, numStrata, stratumValues));
415   PetscFunctionReturn(PETSC_SUCCESS);
416 }
417 
418 static PetscErrorCode DMLabelView_Concrete_Ascii(DMLabel label, PetscViewer viewer)
419 {
420   PetscInt    v;
421   PetscMPIInt rank;
422 
423   PetscFunctionBegin;
424   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
425   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
426   if (label) {
427     const char *name;
428 
429     PetscCall(PetscObjectGetName((PetscObject)label, &name));
430     PetscCall(PetscViewerASCIIPrintf(viewer, "Label '%s':\n", name));
431     if (label->bt) PetscCall(PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%" PetscInt_FMT ", %" PetscInt_FMT ")\n", label->pStart, label->pEnd));
432     for (v = 0; v < label->numStrata; ++v) {
433       const PetscInt  value = label->stratumValues[v];
434       const PetscInt *points;
435       PetscInt        p;
436 
437       PetscCall(ISGetIndices(label->points[v], &points));
438       for (p = 0; p < label->stratumSizes[v]; ++p) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %" PetscInt_FMT " (%" PetscInt_FMT ")\n", rank, points[p], value));
439       PetscCall(ISRestoreIndices(label->points[v], &points));
440     }
441   }
442   PetscCall(PetscViewerFlush(viewer));
443   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
444   PetscFunctionReturn(PETSC_SUCCESS);
445 }
446 
447 static PetscErrorCode DMLabelView_Concrete(DMLabel label, PetscViewer viewer)
448 {
449   PetscBool isascii;
450 
451   PetscFunctionBegin;
452   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
453   if (isascii) PetscCall(DMLabelView_Concrete_Ascii(label, viewer));
454   PetscFunctionReturn(PETSC_SUCCESS);
455 }
456 
457 /*@
458   DMLabelView - View the label
459 
460   Collective
461 
462   Input Parameters:
463 + label  - The `DMLabel`
464 - viewer - The `PetscViewer`
465 
466   Level: intermediate
467 
468 .seealso: `DMLabel`, `PetscViewer`, `DM`, `DMLabelCreate()`, `DMLabelDestroy()`
469 @*/
470 PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
471 {
472   PetscFunctionBegin;
473   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
474   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)label), &viewer));
475   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
476   PetscCall(DMLabelMakeAllValid_Private(label));
477   PetscUseTypeMethod(label, view, viewer);
478   PetscFunctionReturn(PETSC_SUCCESS);
479 }
480 
481 /*@
482   DMLabelViewFromOptions - View a `DMLabel` in a particular way based on a request in the options database
483 
484   Collective
485 
486   Input Parameters:
487 + label - the `DMLabel` object
488 . obj   - optional object that provides the prefix for the options database (if `NULL` then the prefix in `obj` is used)
489 - name  - option string that is used to activate viewing
490 
491   Level: intermediate
492 
493   Note:
494   See `PetscObjectViewFromOptions()` for a list of values that can be provided in the options database to determine how the `DMLabel` is viewed
495 
496 .seealso: [](ch_dmbase), `DMLabel`, `DMLabelView()`, `PetscObjectViewFromOptions()`, `DMLabelCreate()`
497 @*/
498 PetscErrorCode DMLabelViewFromOptions(DMLabel label, PeOp PetscObject obj, const char name[])
499 {
500   PetscFunctionBegin;
501   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
502   PetscCall(PetscObjectViewFromOptions((PetscObject)label, obj, name));
503   PetscFunctionReturn(PETSC_SUCCESS);
504 }
505 
506 /*@
507   DMLabelReset - Destroys internal data structures in a `DMLabel`
508 
509   Not Collective
510 
511   Input Parameter:
512 . label - The `DMLabel`
513 
514   Level: beginner
515 
516 .seealso: `DMLabel`, `DM`, `DMLabelDestroy()`, `DMLabelCreate()`
517 @*/
518 PetscErrorCode DMLabelReset(DMLabel label)
519 {
520   PetscInt v;
521 
522   PetscFunctionBegin;
523   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
524   for (v = 0; v < label->numStrata; ++v) {
525     if (label->ht[v]) PetscCall(PetscHSetIDestroy(&label->ht[v]));
526     PetscCall(ISDestroy(&label->points[v]));
527   }
528   label->numStrata = 0;
529   PetscCall(PetscFree(label->stratumValues));
530   PetscCall(PetscFree(label->stratumSizes));
531   PetscCall(PetscFree(label->ht));
532   PetscCall(PetscFree(label->points));
533   PetscCall(PetscFree(label->validIS));
534   PetscCall(PetscHMapIReset(label->hmap));
535   label->pStart = -1;
536   label->pEnd   = -1;
537   PetscCall(PetscBTDestroy(&label->bt));
538   PetscFunctionReturn(PETSC_SUCCESS);
539 }
540 
541 /*@
542   DMLabelDestroy - Destroys a `DMLabel`
543 
544   Collective
545 
546   Input Parameter:
547 . label - The `DMLabel`
548 
549   Level: beginner
550 
551 .seealso: `DMLabel`, `DM`, `DMLabelReset()`, `DMLabelCreate()`
552 @*/
553 PetscErrorCode DMLabelDestroy(DMLabel *label)
554 {
555   PetscFunctionBegin;
556   if (!*label) PetscFunctionReturn(PETSC_SUCCESS);
557   PetscValidHeaderSpecific(*label, DMLABEL_CLASSID, 1);
558   if (--((PetscObject)*label)->refct > 0) {
559     *label = NULL;
560     PetscFunctionReturn(PETSC_SUCCESS);
561   }
562   PetscCall(DMLabelReset(*label));
563   PetscCall(PetscHMapIDestroy(&(*label)->hmap));
564   PetscCall(PetscHeaderDestroy(label));
565   PetscFunctionReturn(PETSC_SUCCESS);
566 }
567 
568 static PetscErrorCode DMLabelDuplicate_Concrete(DMLabel label, DMLabel *labelnew)
569 {
570   PetscFunctionBegin;
571   for (PetscInt v = 0; v < label->numStrata; ++v) {
572     PetscCall(PetscHSetICreate(&(*labelnew)->ht[v]));
573     PetscCall(PetscObjectReference((PetscObject)label->points[v]));
574     (*labelnew)->points[v] = label->points[v];
575   }
576   PetscCall(PetscHMapIDestroy(&(*labelnew)->hmap));
577   PetscCall(PetscHMapIDuplicate(label->hmap, &(*labelnew)->hmap));
578   PetscFunctionReturn(PETSC_SUCCESS);
579 }
580 
581 /*@
582   DMLabelDuplicate - Duplicates a `DMLabel`
583 
584   Collective
585 
586   Input Parameter:
587 . label - The `DMLabel`
588 
589   Output Parameter:
590 . labelnew - new label
591 
592   Level: intermediate
593 
594 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelDestroy()`
595 @*/
596 PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
597 {
598   const char *name;
599 
600   PetscFunctionBegin;
601   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
602   PetscCall(DMLabelMakeAllValid_Private(label));
603   PetscCall(PetscObjectGetName((PetscObject)label, &name));
604   PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)label), name, labelnew));
605 
606   (*labelnew)->numStrata    = label->numStrata;
607   (*labelnew)->defaultValue = label->defaultValue;
608   (*labelnew)->readonly     = label->readonly;
609   PetscCall(PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues));
610   PetscCall(PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes));
611   PetscCall(PetscCalloc1(label->numStrata, &(*labelnew)->ht));
612   PetscCall(PetscCalloc1(label->numStrata, &(*labelnew)->points));
613   PetscCall(PetscMalloc1(label->numStrata, &(*labelnew)->validIS));
614   for (PetscInt v = 0; v < label->numStrata; ++v) {
615     (*labelnew)->stratumValues[v] = label->stratumValues[v];
616     (*labelnew)->stratumSizes[v]  = label->stratumSizes[v];
617     (*labelnew)->validIS[v]       = PETSC_TRUE;
618   }
619   (*labelnew)->pStart = -1;
620   (*labelnew)->pEnd   = -1;
621   (*labelnew)->bt     = NULL;
622   PetscUseTypeMethod(label, duplicate, labelnew);
623   PetscFunctionReturn(PETSC_SUCCESS);
624 }
625 
626 /*@C
627   DMLabelCompare - Compare two `DMLabel` objects
628 
629   Collective; No Fortran Support
630 
631   Input Parameters:
632 + comm - Comm over which to compare labels
633 . l0   - First `DMLabel`
634 - l1   - Second `DMLabel`
635 
636   Output Parameters:
637 + equal   - (Optional) Flag whether the two labels are equal
638 - message - (Optional) Message describing the difference
639 
640   Level: intermediate
641 
642   Notes:
643   The output flag equal is the same on all processes.
644   If it is passed as `NULL` and difference is found, an error is thrown on all processes.
645   Make sure to pass `NULL` on all processes.
646 
647   The output message is set independently on each rank.
648   It is set to `NULL` if no difference was found on the current rank. It must be freed by user.
649   If message is passed as `NULL` and difference is found, the difference description is printed to stderr in synchronized manner.
650   Make sure to pass `NULL` on all processes.
651 
652   For the comparison, we ignore the order of stratum values, and strata with no points.
653 
654   The communicator needs to be specified because currently `DMLabel` can live on `PETSC_COMM_SELF` even if the underlying `DM` is parallel.
655 
656   Developer Note:
657   Fortran stub cannot be generated automatically because `message` must be freed with `PetscFree()`
658 
659 .seealso: `DMLabel`, `DM`, `DMCompareLabels()`, `DMLabelGetNumValues()`, `DMLabelGetDefaultValue()`, `DMLabelGetNonEmptyStratumValuesIS()`, `DMLabelGetStratumIS()`
660 @*/
661 PetscErrorCode DMLabelCompare(MPI_Comm comm, DMLabel l0, DMLabel l1, PetscBool *equal, char **message)
662 {
663   const char *name0, *name1;
664   char        msg[PETSC_MAX_PATH_LEN] = "";
665   PetscBool   eq;
666   PetscMPIInt rank;
667 
668   PetscFunctionBegin;
669   PetscValidHeaderSpecific(l0, DMLABEL_CLASSID, 2);
670   PetscValidHeaderSpecific(l1, DMLABEL_CLASSID, 3);
671   if (equal) PetscAssertPointer(equal, 4);
672   if (message) PetscAssertPointer(message, 5);
673   PetscCallMPI(MPI_Comm_rank(comm, &rank));
674   PetscCall(PetscObjectGetName((PetscObject)l0, &name0));
675   PetscCall(PetscObjectGetName((PetscObject)l1, &name1));
676   {
677     PetscInt v0, v1;
678 
679     PetscCall(DMLabelGetDefaultValue(l0, &v0));
680     PetscCall(DMLabelGetDefaultValue(l1, &v1));
681     eq = (PetscBool)(v0 == v1);
682     if (!eq) PetscCall(PetscSNPrintf(msg, sizeof(msg), "Default value of DMLabel l0 \"%s\" = %" PetscInt_FMT " != %" PetscInt_FMT " = Default value of DMLabel l1 \"%s\"", name0, v0, v1, name1));
683     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPI_C_BOOL, MPI_LAND, comm));
684     if (!eq) goto finish;
685   }
686   {
687     IS is0, is1;
688 
689     PetscCall(DMLabelGetNonEmptyStratumValuesIS(l0, &is0));
690     PetscCall(DMLabelGetNonEmptyStratumValuesIS(l1, &is1));
691     PetscCall(ISEqual(is0, is1, &eq));
692     PetscCall(ISDestroy(&is0));
693     PetscCall(ISDestroy(&is1));
694     if (!eq) PetscCall(PetscSNPrintf(msg, sizeof(msg), "Stratum values in DMLabel l0 \"%s\" are different than in DMLabel l1 \"%s\"", name0, name1));
695     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPI_C_BOOL, MPI_LAND, comm));
696     if (!eq) goto finish;
697   }
698   {
699     PetscInt i, nValues;
700 
701     PetscCall(DMLabelGetNumValues(l0, &nValues));
702     for (i = 0; i < nValues; i++) {
703       const PetscInt v = l0->stratumValues[i];
704       PetscInt       n;
705       IS             is0, is1;
706 
707       PetscCall(DMLabelGetStratumSize_Private(l0, i, &n));
708       if (!n) continue;
709       PetscCall(DMLabelGetStratumIS(l0, v, &is0));
710       PetscCall(DMLabelGetStratumIS(l1, v, &is1));
711       PetscCall(ISEqualUnsorted(is0, is1, &eq));
712       PetscCall(ISDestroy(&is0));
713       PetscCall(ISDestroy(&is1));
714       if (!eq) {
715         PetscCall(PetscSNPrintf(msg, sizeof(msg), "Stratum #%" PetscInt_FMT " with value %" PetscInt_FMT " contains different points in DMLabel l0 \"%s\" and DMLabel l1 \"%s\"", i, v, name0, name1));
716         break;
717       }
718     }
719     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &eq, 1, MPI_C_BOOL, MPI_LAND, comm));
720   }
721 finish:
722   /* If message output arg not set, print to stderr */
723   if (message) {
724     *message = NULL;
725     if (msg[0]) PetscCall(PetscStrallocpy(msg, message));
726   } else {
727     if (msg[0]) PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] %s\n", rank, msg));
728     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR));
729   }
730   /* If same output arg not ser and labels are not equal, throw error */
731   if (equal) *equal = eq;
732   else PetscCheck(eq, comm, PETSC_ERR_ARG_INCOMP, "DMLabels l0 \"%s\" and l1 \"%s\" are not equal", name0, name1);
733   PetscFunctionReturn(PETSC_SUCCESS);
734 }
735 
736 /*@
737   DMLabelComputeIndex - Create an index structure for membership determination, automatically determining the bounds
738 
739   Not Collective
740 
741   Input Parameter:
742 . label - The `DMLabel`
743 
744   Level: intermediate
745 
746 .seealso: `DMLabel`, `DM`, `DMLabelHasPoint()`, `DMLabelCreateIndex()`, `DMLabelDestroyIndex()`, `DMLabelGetValue()`, `DMLabelSetValue()`
747 @*/
748 PetscErrorCode DMLabelComputeIndex(DMLabel label)
749 {
750   PetscInt pStart = PETSC_INT_MAX, pEnd = -1, v;
751 
752   PetscFunctionBegin;
753   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
754   PetscCall(DMLabelMakeAllValid_Private(label));
755   for (v = 0; v < label->numStrata; ++v) {
756     const PetscInt *points;
757     PetscInt        i;
758 
759     PetscCall(ISGetIndices(label->points[v], &points));
760     for (i = 0; i < label->stratumSizes[v]; ++i) {
761       const PetscInt point = points[i];
762 
763       pStart = PetscMin(point, pStart);
764       pEnd   = PetscMax(point + 1, pEnd);
765     }
766     PetscCall(ISRestoreIndices(label->points[v], &points));
767   }
768   label->pStart = pStart == PETSC_INT_MAX ? -1 : pStart;
769   label->pEnd   = pEnd;
770   PetscCall(DMLabelCreateIndex(label, label->pStart, label->pEnd));
771   PetscFunctionReturn(PETSC_SUCCESS);
772 }
773 
774 /*@
775   DMLabelCreateIndex - Create an index structure for membership determination
776 
777   Not Collective
778 
779   Input Parameters:
780 + label  - The `DMLabel`
781 . pStart - The smallest point
782 - pEnd   - The largest point + 1
783 
784   Level: intermediate
785 
786 .seealso: `DMLabel`, `DM`, `DMLabelHasPoint()`, `DMLabelComputeIndex()`, `DMLabelDestroyIndex()`, `DMLabelGetValue()`, `DMLabelSetValue()`
787 @*/
788 PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
789 {
790   PetscInt v;
791 
792   PetscFunctionBegin;
793   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
794   PetscCall(DMLabelDestroyIndex(label));
795   PetscCall(DMLabelMakeAllValid_Private(label));
796   label->pStart = pStart;
797   label->pEnd   = pEnd;
798   /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
799   PetscCall(PetscBTCreate(pEnd - pStart, &label->bt));
800   for (v = 0; v < label->numStrata; ++v) {
801     IS              pointIS;
802     const PetscInt *points;
803     PetscInt        i;
804 
805     PetscUseTypeMethod(label, getstratumis, v, &pointIS);
806     PetscCall(ISGetIndices(pointIS, &points));
807     for (i = 0; i < label->stratumSizes[v]; ++i) {
808       const PetscInt point = points[i];
809 
810       PetscCheck(!(point < pStart) && !(point >= pEnd), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %" PetscInt_FMT " in stratum %" PetscInt_FMT " is not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, label->stratumValues[v], pStart, pEnd);
811       PetscCall(PetscBTSet(label->bt, point - pStart));
812     }
813     PetscCall(ISRestoreIndices(label->points[v], &points));
814     PetscCall(ISDestroy(&pointIS));
815   }
816   PetscFunctionReturn(PETSC_SUCCESS);
817 }
818 
819 /*@
820   DMLabelDestroyIndex - Destroy the index structure
821 
822   Not Collective
823 
824   Input Parameter:
825 . label - the `DMLabel`
826 
827   Level: intermediate
828 
829 .seealso: `DMLabel`, `DM`, `DMLabelHasPoint()`, `DMLabelCreateIndex()`, `DMLabelGetValue()`, `DMLabelSetValue()`
830 @*/
831 PetscErrorCode DMLabelDestroyIndex(DMLabel label)
832 {
833   PetscFunctionBegin;
834   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
835   label->pStart = -1;
836   label->pEnd   = -1;
837   PetscCall(PetscBTDestroy(&label->bt));
838   PetscFunctionReturn(PETSC_SUCCESS);
839 }
840 
841 /*@
842   DMLabelGetBounds - Return the smallest and largest point in the label
843 
844   Not Collective
845 
846   Input Parameter:
847 . label - the `DMLabel`
848 
849   Output Parameters:
850 + pStart - The smallest point
851 - pEnd   - The largest point + 1
852 
853   Level: intermediate
854 
855   Note:
856   This will compute an index for the label if one does not exist.
857 
858 .seealso: `DMLabel`, `DM`, `DMLabelHasPoint()`, `DMLabelCreateIndex()`, `DMLabelGetValue()`, `DMLabelSetValue()`
859 @*/
860 PetscErrorCode DMLabelGetBounds(DMLabel label, PetscInt *pStart, PetscInt *pEnd)
861 {
862   PetscFunctionBegin;
863   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
864   if ((label->pStart == -1) && (label->pEnd == -1)) PetscCall(DMLabelComputeIndex(label));
865   if (pStart) {
866     PetscAssertPointer(pStart, 2);
867     *pStart = label->pStart;
868   }
869   if (pEnd) {
870     PetscAssertPointer(pEnd, 3);
871     *pEnd = label->pEnd;
872   }
873   PetscFunctionReturn(PETSC_SUCCESS);
874 }
875 
876 /*@
877   DMLabelHasValue - Determine whether a label assigns the value to any point
878 
879   Not Collective
880 
881   Input Parameters:
882 + label - the `DMLabel`
883 - value - the value
884 
885   Output Parameter:
886 . contains - Flag indicating whether the label maps this value to any point
887 
888   Level: developer
889 
890 .seealso: `DMLabel`, `DM`, `DMLabelHasPoint()`, `DMLabelGetValue()`, `DMLabelSetValue()`
891 @*/
892 PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
893 {
894   PetscInt v;
895 
896   PetscFunctionBegin;
897   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
898   PetscAssertPointer(contains, 3);
899   PetscCall(DMLabelLookupStratum(label, value, &v));
900   *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE;
901   PetscFunctionReturn(PETSC_SUCCESS);
902 }
903 
904 /*@
905   DMLabelHasPoint - Determine whether a label assigns a value to a point
906 
907   Not Collective
908 
909   Input Parameters:
910 + label - the `DMLabel`
911 - point - the point
912 
913   Output Parameter:
914 . contains - Flag indicating whether the label maps this point to a value
915 
916   Level: developer
917 
918   Note:
919   The user must call `DMLabelCreateIndex()` before this function.
920 
921 .seealso: `DMLabel`, `DM`, `DMLabelCreateIndex()`, `DMLabelGetValue()`, `DMLabelSetValue()`
922 @*/
923 PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
924 {
925   PetscFunctionBeginHot;
926   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
927   PetscAssertPointer(contains, 3);
928   PetscCall(DMLabelMakeAllValid_Private(label));
929   if (PetscDefined(USE_DEBUG)) {
930     PetscCheck(label->bt, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
931     PetscCheck(!(point < label->pStart) && !(point >= label->pEnd), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %" PetscInt_FMT " is not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, label->pStart, label->pEnd);
932   }
933   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
934   PetscFunctionReturn(PETSC_SUCCESS);
935 }
936 
937 /*@
938   DMLabelStratumHasPoint - Return true if the stratum contains a point
939 
940   Not Collective
941 
942   Input Parameters:
943 + label - the `DMLabel`
944 . value - the stratum value
945 - point - the point
946 
947   Output Parameter:
948 . contains - true if the stratum contains the point
949 
950   Level: intermediate
951 
952 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelSetValue()`, `DMLabelClearValue()`
953 @*/
954 PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
955 {
956   PetscFunctionBeginHot;
957   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
958   PetscAssertPointer(contains, 4);
959   if (value == label->defaultValue) {
960     PetscInt pointVal;
961 
962     PetscCall(DMLabelGetValue(label, point, &pointVal));
963     *contains = (PetscBool)(pointVal == value);
964   } else {
965     PetscInt v;
966 
967     PetscCall(DMLabelLookupStratum(label, value, &v));
968     if (v >= 0) {
969       if (label->validIS[v] || label->readonly) {
970         IS       is;
971         PetscInt i;
972 
973         PetscUseTypeMethod(label, getstratumis, v, &is);
974         PetscCall(ISLocate(is, point, &i));
975         PetscCall(ISDestroy(&is));
976         *contains = (PetscBool)(i >= 0);
977       } else {
978         PetscCall(PetscHSetIHas(label->ht[v], point, contains));
979       }
980     } else { // value is not present
981       *contains = PETSC_FALSE;
982     }
983   }
984   PetscFunctionReturn(PETSC_SUCCESS);
985 }
986 
987 /*@
988   DMLabelGetDefaultValue - Get the default value returned by `DMLabelGetValue()` if a point has not been explicitly given a value.
989   When a label is created, it is initialized to -1.
990 
991   Not Collective
992 
993   Input Parameter:
994 . label - a `DMLabel` object
995 
996   Output Parameter:
997 . defaultValue - the default value
998 
999   Level: beginner
1000 
1001 .seealso: `DMLabel`, `DM`, `DMLabelSetDefaultValue()`, `DMLabelGetValue()`, `DMLabelSetValue()`
1002 @*/
1003 PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
1004 {
1005   PetscFunctionBegin;
1006   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1007   *defaultValue = label->defaultValue;
1008   PetscFunctionReturn(PETSC_SUCCESS);
1009 }
1010 
1011 /*@
1012   DMLabelSetDefaultValue - Set the default value returned by `DMLabelGetValue()` if a point has not been explicitly given a value.
1013   When a label is created, it is initialized to -1.
1014 
1015   Not Collective
1016 
1017   Input Parameter:
1018 . label - a `DMLabel` object
1019 
1020   Output Parameter:
1021 . defaultValue - the default value
1022 
1023   Level: beginner
1024 
1025 .seealso: `DMLabel`, `DM`, `DMLabelGetDefaultValue()`, `DMLabelGetValue()`, `DMLabelSetValue()`
1026 @*/
1027 PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
1028 {
1029   PetscFunctionBegin;
1030   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1031   label->defaultValue = defaultValue;
1032   PetscFunctionReturn(PETSC_SUCCESS);
1033 }
1034 
1035 /*@
1036   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
1037   `DMLabelSetDefaultValue()`)
1038 
1039   Not Collective
1040 
1041   Input Parameters:
1042 + label - the `DMLabel`
1043 - point - the point
1044 
1045   Output Parameter:
1046 . value - The point value, or the default value (-1 by default)
1047 
1048   Level: intermediate
1049 
1050   Note:
1051   A label may assign multiple values to a point.  No guarantees are made about which value is returned in that case.
1052   Use `DMLabelStratumHasPoint()` to check for inclusion in a specific value stratum.
1053 
1054 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelSetValue()`, `DMLabelClearValue()`, `DMLabelGetDefaultValue()`, `DMLabelSetDefaultValue()`
1055 @*/
1056 PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
1057 {
1058   PetscInt v;
1059 
1060   PetscFunctionBeginHot;
1061   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1062   PetscAssertPointer(value, 3);
1063   *value = label->defaultValue;
1064   for (v = 0; v < label->numStrata; ++v) {
1065     if (label->validIS[v] || label->readonly) {
1066       IS       is;
1067       PetscInt i;
1068 
1069       PetscUseTypeMethod(label, getstratumis, v, &is);
1070       PetscCall(ISLocate(label->points[v], point, &i));
1071       PetscCall(ISDestroy(&is));
1072       if (i >= 0) {
1073         *value = label->stratumValues[v];
1074         break;
1075       }
1076     } else {
1077       PetscBool has;
1078 
1079       PetscCall(PetscHSetIHas(label->ht[v], point, &has));
1080       if (has) {
1081         *value = label->stratumValues[v];
1082         break;
1083       }
1084     }
1085   }
1086   PetscFunctionReturn(PETSC_SUCCESS);
1087 }
1088 
1089 /*@
1090   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
1091   be changed with `DMLabelSetDefaultValue()` to something different), then this function will do nothing.
1092 
1093   Not Collective
1094 
1095   Input Parameters:
1096 + label - the `DMLabel`
1097 . point - the point
1098 - value - The point value
1099 
1100   Level: intermediate
1101 
1102 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelClearValue()`, `DMLabelGetDefaultValue()`, `DMLabelSetDefaultValue()`
1103 @*/
1104 PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
1105 {
1106   PetscInt v;
1107 
1108   PetscFunctionBegin;
1109   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1110   /* Find label value, add new entry if needed */
1111   if (value == label->defaultValue) PetscFunctionReturn(PETSC_SUCCESS);
1112   PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered");
1113   PetscCall(DMLabelLookupAddStratum(label, value, &v));
1114   /* Set key */
1115   PetscCall(DMLabelMakeInvalid_Private(label, v));
1116   PetscCall(PetscHSetIAdd(label->ht[v], point));
1117   PetscFunctionReturn(PETSC_SUCCESS);
1118 }
1119 
1120 /*@
1121   DMLabelClearValue - Clear the value a label assigns to a point
1122 
1123   Not Collective
1124 
1125   Input Parameters:
1126 + label - the `DMLabel`
1127 . point - the point
1128 - value - The point value
1129 
1130   Level: intermediate
1131 
1132 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`
1133 @*/
1134 PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
1135 {
1136   PetscInt v;
1137 
1138   PetscFunctionBegin;
1139   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1140   PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered");
1141   /* Find label value */
1142   PetscCall(DMLabelLookupStratum(label, value, &v));
1143   if (v < 0) PetscFunctionReturn(PETSC_SUCCESS);
1144 
1145   if (label->bt) {
1146     PetscCheck(!(point < label->pStart) && !(point >= label->pEnd), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %" PetscInt_FMT " is not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, label->pStart, label->pEnd);
1147     PetscCall(PetscBTClear(label->bt, point - label->pStart));
1148   }
1149 
1150   /* Delete key */
1151   PetscCall(DMLabelMakeInvalid_Private(label, v));
1152   PetscCall(PetscHSetIDel(label->ht[v], point));
1153   PetscFunctionReturn(PETSC_SUCCESS);
1154 }
1155 
1156 /*@
1157   DMLabelInsertIS - Set all points in the `IS` to a value
1158 
1159   Not Collective
1160 
1161   Input Parameters:
1162 + label - the `DMLabel`
1163 . is    - the point `IS`
1164 - value - The point value
1165 
1166   Level: intermediate
1167 
1168 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1169 @*/
1170 PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
1171 {
1172   PetscInt        v, n, p;
1173   const PetscInt *points;
1174 
1175   PetscFunctionBegin;
1176   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1177   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
1178   /* Find label value, add new entry if needed */
1179   if (value == label->defaultValue) PetscFunctionReturn(PETSC_SUCCESS);
1180   PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered");
1181   PetscCall(DMLabelLookupAddStratum(label, value, &v));
1182   /* Set keys */
1183   PetscCall(DMLabelMakeInvalid_Private(label, v));
1184   PetscCall(ISGetLocalSize(is, &n));
1185   PetscCall(ISGetIndices(is, &points));
1186   for (p = 0; p < n; ++p) PetscCall(PetscHSetIAdd(label->ht[v], points[p]));
1187   PetscCall(ISRestoreIndices(is, &points));
1188   PetscFunctionReturn(PETSC_SUCCESS);
1189 }
1190 
1191 /*@
1192   DMLabelGetNumValues - Get the number of values that the `DMLabel` takes
1193 
1194   Not Collective
1195 
1196   Input Parameter:
1197 . label - the `DMLabel`
1198 
1199   Output Parameter:
1200 . numValues - the number of values
1201 
1202   Level: intermediate
1203 
1204 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1205 @*/
1206 PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
1207 {
1208   PetscFunctionBegin;
1209   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1210   PetscAssertPointer(numValues, 2);
1211   *numValues = label->numStrata;
1212   PetscFunctionReturn(PETSC_SUCCESS);
1213 }
1214 
1215 /*@
1216   DMLabelGetValueIS - Get an `IS` of all values that the `DMlabel` takes
1217 
1218   Not Collective
1219 
1220   Input Parameter:
1221 . label - the `DMLabel`
1222 
1223   Output Parameter:
1224 . values - the value `IS`
1225 
1226   Level: intermediate
1227 
1228   Notes:
1229   The output `IS` should be destroyed when no longer needed.
1230   Strata which are allocated but empty [`DMLabelGetStratumSize()` yields 0] are counted.
1231   If you need to count only nonempty strata, use `DMLabelGetNonEmptyStratumValuesIS()`.
1232 
1233 .seealso: `DMLabel`, `DM`, `DMLabelGetNonEmptyStratumValuesIS()`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1234 @*/
1235 PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
1236 {
1237   PetscFunctionBegin;
1238   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1239   PetscAssertPointer(values, 2);
1240   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values));
1241   PetscFunctionReturn(PETSC_SUCCESS);
1242 }
1243 
1244 /*@
1245   DMLabelGetValueBounds - Return the smallest and largest value in the label
1246 
1247   Not Collective
1248 
1249   Input Parameter:
1250 . label - the `DMLabel`
1251 
1252   Output Parameters:
1253 + minValue - The smallest value
1254 - maxValue - The largest value
1255 
1256   Level: intermediate
1257 
1258 .seealso: `DMLabel`, `DM`, `DMLabelGetBounds()`, `DMLabelGetValue()`, `DMLabelSetValue()`
1259 @*/
1260 PetscErrorCode DMLabelGetValueBounds(DMLabel label, PetscInt *minValue, PetscInt *maxValue)
1261 {
1262   PetscInt min = PETSC_INT_MAX, max = PETSC_INT_MIN;
1263 
1264   PetscFunctionBegin;
1265   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1266   for (PetscInt v = 0; v < label->numStrata; ++v) {
1267     min = PetscMin(min, label->stratumValues[v]);
1268     max = PetscMax(max, label->stratumValues[v]);
1269   }
1270   if (minValue) {
1271     PetscAssertPointer(minValue, 2);
1272     *minValue = min;
1273   }
1274   if (maxValue) {
1275     PetscAssertPointer(maxValue, 3);
1276     *maxValue = max;
1277   }
1278   PetscFunctionReturn(PETSC_SUCCESS);
1279 }
1280 
1281 /*@
1282   DMLabelGetNonEmptyStratumValuesIS - Get an `IS` of all values that the `DMlabel` takes
1283 
1284   Not Collective
1285 
1286   Input Parameter:
1287 . label - the `DMLabel`
1288 
1289   Output Parameter:
1290 . values - the value `IS`
1291 
1292   Level: intermediate
1293 
1294   Notes:
1295   The output `IS` should be destroyed when no longer needed.
1296   This is similar to `DMLabelGetValueIS()` but counts only nonempty strata.
1297 
1298 .seealso: `DMLabel`, `DM`, `DMLabelGetValueIS()`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1299 @*/
1300 PetscErrorCode DMLabelGetNonEmptyStratumValuesIS(DMLabel label, IS *values)
1301 {
1302   PetscInt  i, j;
1303   PetscInt *valuesArr;
1304 
1305   PetscFunctionBegin;
1306   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1307   PetscAssertPointer(values, 2);
1308   PetscCall(PetscMalloc1(label->numStrata, &valuesArr));
1309   for (i = 0, j = 0; i < label->numStrata; i++) {
1310     PetscInt n;
1311 
1312     PetscCall(DMLabelGetStratumSize_Private(label, i, &n));
1313     if (n) valuesArr[j++] = label->stratumValues[i];
1314   }
1315   if (j == label->numStrata) {
1316     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values));
1317   } else {
1318     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, j, valuesArr, PETSC_COPY_VALUES, values));
1319   }
1320   PetscCall(PetscFree(valuesArr));
1321   PetscFunctionReturn(PETSC_SUCCESS);
1322 }
1323 
1324 /*@
1325   DMLabelGetValueIndex - Get the index of a given value in the list of values for the `DMlabel`, or -1 if it is not present
1326 
1327   Not Collective
1328 
1329   Input Parameters:
1330 + label - the `DMLabel`
1331 - value - the value
1332 
1333   Output Parameter:
1334 . index - the index of value in the list of values
1335 
1336   Level: intermediate
1337 
1338 .seealso: `DMLabel`, `DM`, `DMLabelGetValueIS()`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1339 @*/
1340 PetscErrorCode DMLabelGetValueIndex(DMLabel label, PetscInt value, PetscInt *index)
1341 {
1342   PetscInt v;
1343 
1344   PetscFunctionBegin;
1345   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1346   PetscAssertPointer(index, 3);
1347   /* Do not assume they are sorted */
1348   for (v = 0; v < label->numStrata; ++v)
1349     if (label->stratumValues[v] == value) break;
1350   if (v >= label->numStrata) *index = -1;
1351   else *index = v;
1352   PetscFunctionReturn(PETSC_SUCCESS);
1353 }
1354 
1355 /*@
1356   DMLabelHasStratum - Determine whether points exist with the given value
1357 
1358   Not Collective
1359 
1360   Input Parameters:
1361 + label - the `DMLabel`
1362 - value - the stratum value
1363 
1364   Output Parameter:
1365 . exists - Flag saying whether points exist
1366 
1367   Level: intermediate
1368 
1369 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1370 @*/
1371 PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
1372 {
1373   PetscInt v;
1374 
1375   PetscFunctionBegin;
1376   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1377   PetscAssertPointer(exists, 3);
1378   PetscCall(DMLabelLookupStratum(label, value, &v));
1379   *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE;
1380   PetscFunctionReturn(PETSC_SUCCESS);
1381 }
1382 
1383 /*@
1384   DMLabelGetStratumSize - Get the size of a stratum
1385 
1386   Not Collective
1387 
1388   Input Parameters:
1389 + label - the `DMLabel`
1390 - value - the stratum value
1391 
1392   Output Parameter:
1393 . size - The number of points in the stratum
1394 
1395   Level: intermediate
1396 
1397 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1398 @*/
1399 PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
1400 {
1401   PetscInt v;
1402 
1403   PetscFunctionBegin;
1404   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1405   PetscAssertPointer(size, 3);
1406   PetscCall(DMLabelLookupStratum(label, value, &v));
1407   PetscCall(DMLabelGetStratumSize_Private(label, v, size));
1408   PetscFunctionReturn(PETSC_SUCCESS);
1409 }
1410 
1411 /*@
1412   DMLabelGetStratumBounds - Get the largest and smallest point of a stratum
1413 
1414   Not Collective
1415 
1416   Input Parameters:
1417 + label - the `DMLabel`
1418 - value - the stratum value
1419 
1420   Output Parameters:
1421 + start - the smallest point in the stratum
1422 - end   - the largest point in the stratum
1423 
1424   Level: intermediate
1425 
1426 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1427 @*/
1428 PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
1429 {
1430   IS       is;
1431   PetscInt v, min, max;
1432 
1433   PetscFunctionBegin;
1434   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1435   if (start) {
1436     PetscAssertPointer(start, 3);
1437     *start = -1;
1438   }
1439   if (end) {
1440     PetscAssertPointer(end, 4);
1441     *end = -1;
1442   }
1443   PetscCall(DMLabelLookupStratum(label, value, &v));
1444   if (v < 0) PetscFunctionReturn(PETSC_SUCCESS);
1445   PetscCall(DMLabelMakeValid_Private(label, v));
1446   if (label->stratumSizes[v] <= 0) PetscFunctionReturn(PETSC_SUCCESS);
1447   PetscUseTypeMethod(label, getstratumis, v, &is);
1448   PetscCall(ISGetMinMax(is, &min, &max));
1449   PetscCall(ISDestroy(&is));
1450   if (start) *start = min;
1451   if (end) *end = max + 1;
1452   PetscFunctionReturn(PETSC_SUCCESS);
1453 }
1454 
1455 static PetscErrorCode DMLabelGetStratumIS_Concrete(DMLabel label, PetscInt v, IS *pointIS)
1456 {
1457   PetscFunctionBegin;
1458   PetscCall(PetscObjectReference((PetscObject)label->points[v]));
1459   *pointIS = label->points[v];
1460   PetscFunctionReturn(PETSC_SUCCESS);
1461 }
1462 
1463 /*@
1464   DMLabelGetStratumIS - Get an `IS` with the stratum points
1465 
1466   Not Collective
1467 
1468   Input Parameters:
1469 + label - the `DMLabel`
1470 - value - the stratum value
1471 
1472   Output Parameter:
1473 . points - The stratum points
1474 
1475   Level: intermediate
1476 
1477   Notes:
1478   The output `IS` should be destroyed when no longer needed.
1479   Returns `NULL` if the stratum is empty.
1480 
1481 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1482 @*/
1483 PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
1484 {
1485   PetscInt v;
1486 
1487   PetscFunctionBegin;
1488   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1489   PetscAssertPointer(points, 3);
1490   *points = NULL;
1491   PetscCall(DMLabelLookupStratum(label, value, &v));
1492   if (v < 0) PetscFunctionReturn(PETSC_SUCCESS);
1493   PetscCall(DMLabelMakeValid_Private(label, v));
1494   PetscUseTypeMethod(label, getstratumis, v, points);
1495   PetscFunctionReturn(PETSC_SUCCESS);
1496 }
1497 
1498 /*@
1499   DMLabelSetStratumIS - Set the stratum points using an `IS`
1500 
1501   Not Collective
1502 
1503   Input Parameters:
1504 + label - the `DMLabel`
1505 . value - the stratum value
1506 - is    - The stratum points
1507 
1508   Level: intermediate
1509 
1510 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1511 @*/
1512 PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
1513 {
1514   PetscInt v;
1515 
1516   PetscFunctionBegin;
1517   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1518   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
1519   PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered");
1520   PetscCall(DMLabelLookupAddStratum(label, value, &v));
1521   if (is == label->points[v]) PetscFunctionReturn(PETSC_SUCCESS);
1522   PetscCall(DMLabelClearStratum(label, value));
1523   PetscCall(ISGetLocalSize(is, &label->stratumSizes[v]));
1524   PetscCall(PetscObjectReference((PetscObject)is));
1525   PetscCall(ISDestroy(&label->points[v]));
1526   label->points[v]  = is;
1527   label->validIS[v] = PETSC_TRUE;
1528   PetscCall(PetscObjectStateIncrease((PetscObject)label));
1529   if (label->bt) {
1530     const PetscInt *points;
1531     PetscInt        p;
1532 
1533     PetscCall(ISGetIndices(is, &points));
1534     for (p = 0; p < label->stratumSizes[v]; ++p) {
1535       const PetscInt point = points[p];
1536 
1537       PetscCheck(!(point < label->pStart) && !(point >= label->pEnd), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %" PetscInt_FMT " is not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, label->pStart, label->pEnd);
1538       PetscCall(PetscBTSet(label->bt, point - label->pStart));
1539     }
1540   }
1541   PetscFunctionReturn(PETSC_SUCCESS);
1542 }
1543 
1544 /*@
1545   DMLabelClearStratum - Remove a stratum
1546 
1547   Not Collective
1548 
1549   Input Parameters:
1550 + label - the `DMLabel`
1551 - value - the stratum value
1552 
1553   Level: intermediate
1554 
1555 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1556 @*/
1557 PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
1558 {
1559   PetscInt v;
1560 
1561   PetscFunctionBegin;
1562   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1563   PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered");
1564   PetscCall(DMLabelLookupStratum(label, value, &v));
1565   if (v < 0) PetscFunctionReturn(PETSC_SUCCESS);
1566   if (label->validIS[v]) {
1567     if (label->bt) {
1568       PetscInt        i;
1569       const PetscInt *points;
1570 
1571       PetscCall(ISGetIndices(label->points[v], &points));
1572       for (i = 0; i < label->stratumSizes[v]; ++i) {
1573         const PetscInt point = points[i];
1574 
1575         PetscCheck(!(point < label->pStart) && !(point >= label->pEnd), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %" PetscInt_FMT " is not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, label->pStart, label->pEnd);
1576         PetscCall(PetscBTClear(label->bt, point - label->pStart));
1577       }
1578       PetscCall(ISRestoreIndices(label->points[v], &points));
1579     }
1580     label->stratumSizes[v] = 0;
1581     PetscCall(ISDestroy(&label->points[v]));
1582     PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &label->points[v]));
1583     PetscCall(PetscObjectSetName((PetscObject)label->points[v], "indices"));
1584     PetscCall(PetscObjectStateIncrease((PetscObject)label));
1585   } else {
1586     PetscCall(PetscHSetIClear(label->ht[v]));
1587   }
1588   PetscFunctionReturn(PETSC_SUCCESS);
1589 }
1590 
1591 /*@
1592   DMLabelSetStratumBounds - Efficiently give a contiguous set of points a given label value
1593 
1594   Not Collective
1595 
1596   Input Parameters:
1597 + label  - The `DMLabel`
1598 . value  - The label value for all points
1599 . pStart - The first point
1600 - pEnd   - A point beyond all marked points
1601 
1602   Level: intermediate
1603 
1604   Note:
1605   The marks points are [`pStart`, `pEnd`), and only the bounds are stored.
1606 
1607 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelSetStratumIS()`, `DMLabelGetStratumIS()`
1608 @*/
1609 PetscErrorCode DMLabelSetStratumBounds(DMLabel label, PetscInt value, PetscInt pStart, PetscInt pEnd)
1610 {
1611   IS pIS;
1612 
1613   PetscFunctionBegin;
1614   PetscCall(ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &pIS));
1615   PetscCall(DMLabelSetStratumIS(label, value, pIS));
1616   PetscCall(ISDestroy(&pIS));
1617   PetscFunctionReturn(PETSC_SUCCESS);
1618 }
1619 
1620 /*@
1621   DMLabelGetStratumPointIndex - Get the index of a point in a given stratum
1622 
1623   Not Collective
1624 
1625   Input Parameters:
1626 + label - The `DMLabel`
1627 . value - The label value
1628 - p     - A point with this value
1629 
1630   Output Parameter:
1631 . index - The index of this point in the stratum, or -1 if the point is not in the stratum or the stratum does not exist
1632 
1633   Level: intermediate
1634 
1635 .seealso: `DMLabel`, `DM`, `DMLabelGetValueIndex()`, `DMLabelGetStratumIS()`, `DMLabelCreate()`
1636 @*/
1637 PetscErrorCode DMLabelGetStratumPointIndex(DMLabel label, PetscInt value, PetscInt p, PetscInt *index)
1638 {
1639   IS       pointIS;
1640   PetscInt v;
1641 
1642   PetscFunctionBegin;
1643   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1644   PetscAssertPointer(index, 4);
1645   *index = -1;
1646   PetscCall(DMLabelLookupStratum(label, value, &v));
1647   if (v < 0) PetscFunctionReturn(PETSC_SUCCESS);
1648   PetscCall(DMLabelMakeValid_Private(label, v));
1649   PetscUseTypeMethod(label, getstratumis, v, &pointIS);
1650   PetscCall(ISLocate(pointIS, p, index));
1651   PetscCall(ISDestroy(&pointIS));
1652   PetscFunctionReturn(PETSC_SUCCESS);
1653 }
1654 
1655 /*@
1656   DMLabelFilter - Remove all points outside of [`start`, `end`)
1657 
1658   Not Collective
1659 
1660   Input Parameters:
1661 + label - the `DMLabel`
1662 . start - the first point kept
1663 - end   - one more than the last point kept
1664 
1665   Level: intermediate
1666 
1667 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1668 @*/
1669 PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
1670 {
1671   PetscInt v;
1672 
1673   PetscFunctionBegin;
1674   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1675   PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered");
1676   PetscCall(DMLabelDestroyIndex(label));
1677   PetscCall(DMLabelMakeAllValid_Private(label));
1678   for (v = 0; v < label->numStrata; ++v) {
1679     PetscCall(ISGeneralFilter(label->points[v], start, end));
1680     PetscCall(ISGetLocalSize(label->points[v], &label->stratumSizes[v]));
1681   }
1682   PetscCall(DMLabelCreateIndex(label, start, end));
1683   PetscFunctionReturn(PETSC_SUCCESS);
1684 }
1685 
1686 /*@
1687   DMLabelPermute - Create a new label with permuted points
1688 
1689   Not Collective
1690 
1691   Input Parameters:
1692 + label       - the `DMLabel`
1693 - permutation - the point permutation
1694 
1695   Output Parameter:
1696 . labelNew - the new label containing the permuted points
1697 
1698   Level: intermediate
1699 
1700 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1701 @*/
1702 PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
1703 {
1704   const PetscInt *perm;
1705   PetscInt        numValues, numPoints, v, q;
1706 
1707   PetscFunctionBegin;
1708   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1709   PetscValidHeaderSpecific(permutation, IS_CLASSID, 2);
1710   PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered");
1711   PetscCall(DMLabelMakeAllValid_Private(label));
1712   PetscCall(DMLabelDuplicate(label, labelNew));
1713   PetscCall(DMLabelGetNumValues(*labelNew, &numValues));
1714   PetscCall(ISGetLocalSize(permutation, &numPoints));
1715   PetscCall(ISGetIndices(permutation, &perm));
1716   for (v = 0; v < numValues; ++v) {
1717     const PetscInt  size = (*labelNew)->stratumSizes[v];
1718     const PetscInt *points;
1719     PetscInt       *pointsNew;
1720 
1721     PetscCall(ISGetIndices((*labelNew)->points[v], &points));
1722     PetscCall(PetscCalloc1(size, &pointsNew));
1723     for (q = 0; q < size; ++q) {
1724       const PetscInt point = points[q];
1725 
1726       PetscCheck(!(point < 0) && !(point >= numPoints), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ") for the remapping", point, numPoints);
1727       pointsNew[q] = perm[point];
1728     }
1729     PetscCall(ISRestoreIndices((*labelNew)->points[v], &points));
1730     PetscCall(PetscSortInt(size, pointsNew));
1731     PetscCall(ISDestroy(&(*labelNew)->points[v]));
1732     if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) {
1733       PetscCall(ISCreateStride(PETSC_COMM_SELF, size, pointsNew[0], 1, &((*labelNew)->points[v])));
1734       PetscCall(PetscFree(pointsNew));
1735     } else {
1736       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size, pointsNew, PETSC_OWN_POINTER, &((*labelNew)->points[v])));
1737     }
1738     PetscCall(PetscObjectSetName((PetscObject)((*labelNew)->points[v]), "indices"));
1739   }
1740   PetscCall(ISRestoreIndices(permutation, &perm));
1741   if (label->bt) {
1742     PetscCall(PetscBTDestroy(&label->bt));
1743     PetscCall(DMLabelCreateIndex(label, label->pStart, label->pEnd));
1744   }
1745   PetscFunctionReturn(PETSC_SUCCESS);
1746 }
1747 
1748 /*@
1749   DMLabelPermuteValues - Permute the values in a label
1750 
1751   Not collective
1752 
1753   Input Parameters:
1754 + label       - the `DMLabel`
1755 - permutation - the value permutation, permutation[old value] = new value
1756 
1757   Output Parameter:
1758 . label - the `DMLabel` now with permuted values
1759 
1760   Note:
1761   The modification is done in-place
1762 
1763   Level: intermediate
1764 
1765 .seealso: `DMLabelRewriteValues()`, `DMLabel`, `DM`, `DMLabelPermute()`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1766 @*/
1767 PetscErrorCode DMLabelPermuteValues(DMLabel label, IS permutation)
1768 {
1769   PetscInt Nv, Np;
1770 
1771   PetscFunctionBegin;
1772   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1773   PetscValidHeaderSpecific(permutation, IS_CLASSID, 2);
1774   PetscCall(DMLabelGetNumValues(label, &Nv));
1775   PetscCall(ISGetLocalSize(permutation, &Np));
1776   PetscCheck(Np == Nv, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_SIZ, "Permutation has size %" PetscInt_FMT " != %" PetscInt_FMT " number of label values", Np, Nv);
1777   if (PetscDefined(USE_DEBUG)) {
1778     PetscBool flg;
1779     PetscCall(ISGetInfo(permutation, IS_PERMUTATION, IS_LOCAL, PETSC_TRUE, &flg));
1780     PetscCheck(flg, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "IS is not a permutation");
1781   }
1782   PetscCall(DMLabelRewriteValues(label, permutation));
1783   PetscFunctionReturn(PETSC_SUCCESS);
1784 }
1785 
1786 /*@
1787   DMLabelRewriteValues - Permute the values in a label, but some may be omitted
1788 
1789   Not collective
1790 
1791   Input Parameters:
1792 + label       - the `DMLabel`
1793 - permutation - the value permutation, permutation[old value] = new value, but some maybe omitted
1794 
1795   Output Parameter:
1796 . label - the `DMLabel` now with permuted values
1797 
1798   Note:
1799   The modification is done in-place
1800 
1801   Level: intermediate
1802 
1803 .seealso: `DMLabelPermuteValues()`, `DMLabel`, `DM`, `DMLabelPermute()`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1804 @*/
1805 PetscErrorCode DMLabelRewriteValues(DMLabel label, IS permutation)
1806 {
1807   const PetscInt *perm;
1808   PetscInt        Nv, Np;
1809 
1810   PetscFunctionBegin;
1811   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1812   PetscValidHeaderSpecific(permutation, IS_CLASSID, 2);
1813   PetscCall(DMLabelMakeAllValid_Private(label));
1814   PetscCall(DMLabelGetNumValues(label, &Nv));
1815   PetscCall(ISGetLocalSize(permutation, &Np));
1816   PetscCheck(Np >= Nv, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_SIZ, "Permutation has size %" PetscInt_FMT " < %" PetscInt_FMT " number of label values", Np, Nv);
1817   PetscCall(ISGetIndices(permutation, &perm));
1818   for (PetscInt v = 0; v < Nv; ++v) label->stratumValues[v] = perm[label->stratumValues[v]];
1819   PetscCall(ISRestoreIndices(permutation, &perm));
1820   PetscFunctionReturn(PETSC_SUCCESS);
1821 }
1822 
1823 static PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
1824 {
1825   MPI_Comm     comm;
1826   PetscInt     s, l, nroots, nleaves, offset, size;
1827   PetscInt    *remoteOffsets, *rootStrata, *rootIdx;
1828   PetscSection rootSection;
1829   PetscSF      labelSF;
1830 
1831   PetscFunctionBegin;
1832   if (label) PetscCall(DMLabelMakeAllValid_Private(label));
1833   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
1834   /* Build a section of stratum values per point, generate the according SF
1835      and distribute point-wise stratum values to leaves. */
1836   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL));
1837   PetscCall(PetscSectionCreate(comm, &rootSection));
1838   PetscCall(PetscSectionSetChart(rootSection, 0, nroots));
1839   if (label) {
1840     for (s = 0; s < label->numStrata; ++s) {
1841       const PetscInt *points;
1842 
1843       PetscCall(ISGetIndices(label->points[s], &points));
1844       for (l = 0; l < label->stratumSizes[s]; l++) PetscCall(PetscSectionAddDof(rootSection, points[l], 1));
1845       PetscCall(ISRestoreIndices(label->points[s], &points));
1846     }
1847   }
1848   PetscCall(PetscSectionSetUp(rootSection));
1849   /* Create a point-wise array of stratum values */
1850   PetscCall(PetscSectionGetStorageSize(rootSection, &size));
1851   PetscCall(PetscMalloc1(size, &rootStrata));
1852   PetscCall(PetscCalloc1(nroots, &rootIdx));
1853   if (label) {
1854     for (s = 0; s < label->numStrata; ++s) {
1855       const PetscInt *points;
1856 
1857       PetscCall(ISGetIndices(label->points[s], &points));
1858       for (l = 0; l < label->stratumSizes[s]; l++) {
1859         const PetscInt p = points[l];
1860         PetscCall(PetscSectionGetOffset(rootSection, p, &offset));
1861         rootStrata[offset + rootIdx[p]++] = label->stratumValues[s];
1862       }
1863       PetscCall(ISRestoreIndices(label->points[s], &points));
1864     }
1865   }
1866   /* Build SF that maps label points to remote processes */
1867   PetscCall(PetscSectionCreate(comm, leafSection));
1868   PetscCall(PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection));
1869   PetscCall(PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF));
1870   PetscCall(PetscFree(remoteOffsets));
1871   /* Send the strata for each point over the derived SF */
1872   PetscCall(PetscSectionGetStorageSize(*leafSection, &size));
1873   PetscCall(PetscMalloc1(size, leafStrata));
1874   PetscCall(PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata, MPI_REPLACE));
1875   PetscCall(PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata, MPI_REPLACE));
1876   /* Clean up */
1877   PetscCall(PetscFree(rootStrata));
1878   PetscCall(PetscFree(rootIdx));
1879   PetscCall(PetscSectionDestroy(&rootSection));
1880   PetscCall(PetscSFDestroy(&labelSF));
1881   PetscFunctionReturn(PETSC_SUCCESS);
1882 }
1883 
1884 /*@
1885   DMLabelDistribute - Create a new label pushed forward over the `PetscSF`
1886 
1887   Collective
1888 
1889   Input Parameters:
1890 + label - the `DMLabel`
1891 - sf    - the map from old to new distribution
1892 
1893   Output Parameter:
1894 . labelNew - the new redistributed label
1895 
1896   Level: intermediate
1897 
1898 .seealso: `DMLabel`, `DM`, `DMLabelCreate()`, `DMLabelGetValue()`, `DMLabelSetValue()`, `DMLabelClearValue()`
1899 @*/
1900 PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1901 {
1902   MPI_Comm     comm;
1903   PetscSection leafSection;
1904   PetscInt     p, pStart, pEnd, s, size, dof, offset, stratum;
1905   PetscInt    *leafStrata, *strataIdx;
1906   PetscInt   **points;
1907   const char  *lname = NULL;
1908   char        *name;
1909   PetscMPIInt  nameSize;
1910   PetscHSetI   stratumHash;
1911   size_t       len = 0;
1912   PetscMPIInt  rank;
1913 
1914   PetscFunctionBegin;
1915   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1916   if (label) {
1917     PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1918     PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered");
1919     PetscCall(DMLabelMakeAllValid_Private(label));
1920   }
1921   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
1922   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1923   /* Bcast name */
1924   if (rank == 0) {
1925     PetscCall(PetscObjectGetName((PetscObject)label, &lname));
1926     PetscCall(PetscStrlen(lname, &len));
1927   }
1928   PetscCall(PetscMPIIntCast(len, &nameSize));
1929   PetscCallMPI(MPI_Bcast(&nameSize, 1, MPI_INT, 0, comm));
1930   PetscCall(PetscMalloc1(nameSize + 1, &name));
1931   if (rank == 0) PetscCall(PetscArraycpy(name, lname, nameSize + 1));
1932   PetscCallMPI(MPI_Bcast(name, nameSize + 1, MPI_CHAR, 0, comm));
1933   PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, labelNew));
1934   PetscCall(PetscFree(name));
1935   /* Bcast defaultValue */
1936   if (rank == 0) (*labelNew)->defaultValue = label->defaultValue;
1937   PetscCallMPI(MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm));
1938   /* Distribute stratum values over the SF and get the point mapping on the receiver */
1939   PetscCall(DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata));
1940   /* Determine received stratum values and initialise new label*/
1941   PetscCall(PetscHSetICreate(&stratumHash));
1942   PetscCall(PetscSectionGetStorageSize(leafSection, &size));
1943   for (p = 0; p < size; ++p) PetscCall(PetscHSetIAdd(stratumHash, leafStrata[p]));
1944   PetscCall(PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata));
1945   PetscCall(PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS));
1946   for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
1947   PetscCall(PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues));
1948   /* Turn leafStrata into indices rather than stratum values */
1949   offset = 0;
1950   PetscCall(PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues));
1951   PetscCall(PetscSortInt((*labelNew)->numStrata, (*labelNew)->stratumValues));
1952   for (s = 0; s < (*labelNew)->numStrata; ++s) PetscCall(PetscHMapISet((*labelNew)->hmap, (*labelNew)->stratumValues[s], s));
1953   for (p = 0; p < size; ++p) {
1954     for (s = 0; s < (*labelNew)->numStrata; ++s) {
1955       if (leafStrata[p] == (*labelNew)->stratumValues[s]) {
1956         leafStrata[p] = s;
1957         break;
1958       }
1959     }
1960   }
1961   /* Rebuild the point strata on the receiver */
1962   PetscCall(PetscCalloc1((*labelNew)->numStrata, &(*labelNew)->stratumSizes));
1963   PetscCall(PetscSectionGetChart(leafSection, &pStart, &pEnd));
1964   for (p = pStart; p < pEnd; p++) {
1965     PetscCall(PetscSectionGetDof(leafSection, p, &dof));
1966     PetscCall(PetscSectionGetOffset(leafSection, p, &offset));
1967     for (s = 0; s < dof; s++) (*labelNew)->stratumSizes[leafStrata[offset + s]]++;
1968   }
1969   PetscCall(PetscCalloc1((*labelNew)->numStrata, &(*labelNew)->ht));
1970   PetscCall(PetscCalloc1((*labelNew)->numStrata, &(*labelNew)->points));
1971   PetscCall(PetscCalloc1((*labelNew)->numStrata, &points));
1972   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1973     PetscCall(PetscHSetICreate(&(*labelNew)->ht[s]));
1974     PetscCall(PetscMalloc1((*labelNew)->stratumSizes[s], &points[s]));
1975   }
1976   /* Insert points into new strata */
1977   PetscCall(PetscCalloc1((*labelNew)->numStrata, &strataIdx));
1978   PetscCall(PetscSectionGetChart(leafSection, &pStart, &pEnd));
1979   for (p = pStart; p < pEnd; p++) {
1980     PetscCall(PetscSectionGetDof(leafSection, p, &dof));
1981     PetscCall(PetscSectionGetOffset(leafSection, p, &offset));
1982     for (s = 0; s < dof; s++) {
1983       stratum                               = leafStrata[offset + s];
1984       points[stratum][strataIdx[stratum]++] = p;
1985     }
1986   }
1987   for (s = 0; s < (*labelNew)->numStrata; s++) {
1988     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, (*labelNew)->stratumSizes[s], &points[s][0], PETSC_OWN_POINTER, &((*labelNew)->points[s])));
1989     PetscCall(PetscObjectSetName((PetscObject)((*labelNew)->points[s]), "indices"));
1990   }
1991   PetscCall(PetscFree(points));
1992   PetscCall(PetscHSetIDestroy(&stratumHash));
1993   PetscCall(PetscFree(leafStrata));
1994   PetscCall(PetscFree(strataIdx));
1995   PetscCall(PetscSectionDestroy(&leafSection));
1996   PetscFunctionReturn(PETSC_SUCCESS);
1997 }
1998 
1999 /*@
2000   DMLabelGather - Gather all label values from leafs into roots
2001 
2002   Collective
2003 
2004   Input Parameters:
2005 + label - the `DMLabel`
2006 - sf    - the `PetscSF` communication map
2007 
2008   Output Parameter:
2009 . labelNew - the new `DMLabel` with localised leaf values
2010 
2011   Level: developer
2012 
2013   Note:
2014   This is the inverse operation to `DMLabelDistribute()`.
2015 
2016 .seealso: `DMLabel`, `DM`, `DMLabelDistribute()`
2017 @*/
2018 PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
2019 {
2020   MPI_Comm        comm;
2021   PetscSection    rootSection;
2022   PetscSF         sfLabel;
2023   PetscSFNode    *rootPoints, *leafPoints;
2024   PetscInt        p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
2025   const PetscInt *rootDegree, *ilocal;
2026   PetscInt       *rootStrata;
2027   const char     *lname;
2028   char           *name;
2029   PetscMPIInt     nameSize;
2030   size_t          len = 0;
2031   PetscMPIInt     rank, size;
2032 
2033   PetscFunctionBegin;
2034   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
2035   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
2036   PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered");
2037   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
2038   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2039   PetscCallMPI(MPI_Comm_size(comm, &size));
2040   /* Bcast name */
2041   if (rank == 0) {
2042     PetscCall(PetscObjectGetName((PetscObject)label, &lname));
2043     PetscCall(PetscStrlen(lname, &len));
2044   }
2045   PetscCall(PetscMPIIntCast(len, &nameSize));
2046   PetscCallMPI(MPI_Bcast(&nameSize, 1, MPI_INT, 0, comm));
2047   PetscCall(PetscMalloc1(nameSize + 1, &name));
2048   if (rank == 0) PetscCall(PetscArraycpy(name, lname, nameSize + 1));
2049   PetscCallMPI(MPI_Bcast(name, nameSize + 1, MPI_CHAR, 0, comm));
2050   PetscCall(DMLabelCreate(PETSC_COMM_SELF, name, labelNew));
2051   PetscCall(PetscFree(name));
2052   /* Gather rank/index pairs of leaves into local roots to build
2053      an inverse, multi-rooted SF. Note that this ignores local leaf
2054      indexing due to the use of the multiSF in PetscSFGather. */
2055   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL));
2056   PetscCall(PetscMalloc1(nroots, &leafPoints));
2057   for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
2058   for (p = 0; p < nleaves; p++) {
2059     PetscInt ilp = ilocal ? ilocal[p] : p;
2060 
2061     leafPoints[ilp].index = ilp;
2062     leafPoints[ilp].rank  = rank;
2063   }
2064   PetscCall(PetscSFComputeDegreeBegin(sf, &rootDegree));
2065   PetscCall(PetscSFComputeDegreeEnd(sf, &rootDegree));
2066   for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
2067   PetscCall(PetscMalloc1(nmultiroots, &rootPoints));
2068   PetscCall(PetscSFGatherBegin(sf, MPIU_SF_NODE, leafPoints, rootPoints));
2069   PetscCall(PetscSFGatherEnd(sf, MPIU_SF_NODE, leafPoints, rootPoints));
2070   PetscCall(PetscSFCreate(comm, &sfLabel));
2071   PetscCall(PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER));
2072   /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
2073   PetscCall(DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata));
2074   /* Rebuild the point strata on the receiver */
2075   for (p = 0, idx = 0; p < nroots; p++) {
2076     for (d = 0; d < rootDegree[p]; d++) {
2077       PetscCall(PetscSectionGetDof(rootSection, idx + d, &dof));
2078       PetscCall(PetscSectionGetOffset(rootSection, idx + d, &offset));
2079       for (s = 0; s < dof; s++) PetscCall(DMLabelSetValue(*labelNew, p, rootStrata[offset + s]));
2080     }
2081     idx += rootDegree[p];
2082   }
2083   PetscCall(PetscFree(leafPoints));
2084   PetscCall(PetscFree(rootStrata));
2085   PetscCall(PetscSectionDestroy(&rootSection));
2086   PetscCall(PetscSFDestroy(&sfLabel));
2087   PetscFunctionReturn(PETSC_SUCCESS);
2088 }
2089 
2090 static PetscErrorCode DMLabelPropagateInit_Internal(DMLabel label, PetscSF pointSF, PetscInt valArray[])
2091 {
2092   const PetscInt *degree;
2093   const PetscInt *points;
2094   PetscInt        Nr, r, Nl, l, val, defVal;
2095 
2096   PetscFunctionBegin;
2097   PetscCall(DMLabelGetDefaultValue(label, &defVal));
2098   /* Add in leaves */
2099   PetscCall(PetscSFGetGraph(pointSF, &Nr, &Nl, &points, NULL));
2100   for (l = 0; l < Nl; ++l) {
2101     PetscCall(DMLabelGetValue(label, points[l], &val));
2102     if (val != defVal) valArray[points[l]] = val;
2103   }
2104   /* Add in shared roots */
2105   PetscCall(PetscSFComputeDegreeBegin(pointSF, &degree));
2106   PetscCall(PetscSFComputeDegreeEnd(pointSF, &degree));
2107   for (r = 0; r < Nr; ++r) {
2108     if (degree[r]) {
2109       PetscCall(DMLabelGetValue(label, r, &val));
2110       if (val != defVal) valArray[r] = val;
2111     }
2112   }
2113   PetscFunctionReturn(PETSC_SUCCESS);
2114 }
2115 
2116 static PetscErrorCode DMLabelPropagateFini_Internal(DMLabel label, PetscSF pointSF, PetscInt valArray[], PetscErrorCode (*markPoint)(DMLabel, PetscInt, PetscInt, void *), void *ctx)
2117 {
2118   const PetscInt *degree;
2119   const PetscInt *points;
2120   PetscInt        Nr, r, Nl, l, val, defVal;
2121 
2122   PetscFunctionBegin;
2123   PetscCall(DMLabelGetDefaultValue(label, &defVal));
2124   /* Read out leaves */
2125   PetscCall(PetscSFGetGraph(pointSF, &Nr, &Nl, &points, NULL));
2126   for (l = 0; l < Nl; ++l) {
2127     const PetscInt p    = points[l];
2128     const PetscInt cval = valArray[p];
2129 
2130     if (cval != defVal) {
2131       PetscCall(DMLabelGetValue(label, p, &val));
2132       if (val == defVal) {
2133         PetscCall(DMLabelSetValue(label, p, cval));
2134         if (markPoint) PetscCall((*markPoint)(label, p, cval, ctx));
2135       }
2136     }
2137   }
2138   /* Read out shared roots */
2139   PetscCall(PetscSFComputeDegreeBegin(pointSF, &degree));
2140   PetscCall(PetscSFComputeDegreeEnd(pointSF, &degree));
2141   for (r = 0; r < Nr; ++r) {
2142     if (degree[r]) {
2143       const PetscInt cval = valArray[r];
2144 
2145       if (cval != defVal) {
2146         PetscCall(DMLabelGetValue(label, r, &val));
2147         if (val == defVal) {
2148           PetscCall(DMLabelSetValue(label, r, cval));
2149           if (markPoint) PetscCall((*markPoint)(label, r, cval, ctx));
2150         }
2151       }
2152     }
2153   }
2154   PetscFunctionReturn(PETSC_SUCCESS);
2155 }
2156 
2157 /*@
2158   DMLabelPropagateBegin - Setup a cycle of label propagation
2159 
2160   Collective
2161 
2162   Input Parameters:
2163 + label - The `DMLabel` to propagate across processes
2164 - sf    - The `PetscSF` describing parallel layout of the label points
2165 
2166   Level: intermediate
2167 
2168 .seealso: `DMLabel`, `DM`, `DMLabelPropagateEnd()`, `DMLabelPropagatePush()`
2169 @*/
2170 PetscErrorCode DMLabelPropagateBegin(DMLabel label, PetscSF sf)
2171 {
2172   PetscInt    Nr, r, defVal;
2173   PetscMPIInt size;
2174 
2175   PetscFunctionBegin;
2176   PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered");
2177   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sf), &size));
2178   if (size > 1) {
2179     PetscCall(DMLabelGetDefaultValue(label, &defVal));
2180     PetscCall(PetscSFGetGraph(sf, &Nr, NULL, NULL, NULL));
2181     if (Nr >= 0) PetscCall(PetscMalloc1(Nr, &label->propArray));
2182     for (r = 0; r < Nr; ++r) label->propArray[r] = defVal;
2183   }
2184   PetscFunctionReturn(PETSC_SUCCESS);
2185 }
2186 
2187 /*@
2188   DMLabelPropagateEnd - Tear down a cycle of label propagation
2189 
2190   Collective
2191 
2192   Input Parameters:
2193 + label   - The `DMLabel` to propagate across processes
2194 - pointSF - The `PetscSF` describing parallel layout of the label points
2195 
2196   Level: intermediate
2197 
2198 .seealso: `DMLabel`, `DM`, `DMLabelPropagateBegin()`, `DMLabelPropagatePush()`
2199 @*/
2200 PetscErrorCode DMLabelPropagateEnd(DMLabel label, PetscSF pointSF)
2201 {
2202   PetscFunctionBegin;
2203   PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered");
2204   PetscCall(PetscFree(label->propArray));
2205   label->propArray = NULL;
2206   PetscFunctionReturn(PETSC_SUCCESS);
2207 }
2208 
2209 /*@C
2210   DMLabelPropagatePush - Tear down a cycle of label propagation
2211 
2212   Collective
2213 
2214   Input Parameters:
2215 + label     - The `DMLabel` to propagate across processes
2216 . pointSF   - The `PetscSF` describing parallel layout of the label points
2217 . markPoint - An optional callback that is called when a point is marked, or `NULL`
2218 - ctx       - An optional user context for the callback, or `NULL`
2219 
2220   Calling sequence of `markPoint`:
2221 + label - The `DMLabel`
2222 . p     - The point being marked
2223 . val   - The label value for `p`
2224 - ctx   - An optional user context
2225 
2226   Level: intermediate
2227 
2228 .seealso: `DMLabel`, `DM`, `DMLabelPropagateBegin()`, `DMLabelPropagateEnd()`
2229 @*/
2230 PetscErrorCode DMLabelPropagatePush(DMLabel label, PetscSF pointSF, PetscErrorCode (*markPoint)(DMLabel label, PetscInt p, PetscInt val, void *ctx), void *ctx)
2231 {
2232   PetscInt   *valArray = label->propArray, Nr;
2233   PetscMPIInt size;
2234 
2235   PetscFunctionBegin;
2236   PetscCheck(!label->readonly, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_WRONG, "Read-only labels cannot be altered");
2237   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pointSF), &size));
2238   PetscCall(PetscSFGetGraph(pointSF, &Nr, NULL, NULL, NULL));
2239   if (size > 1 && Nr >= 0) {
2240     /* Communicate marked edges
2241        The current implementation allocates an array the size of the number of root. We put the label values into the
2242        array, and then call PetscSFReduce()+PetscSFBcast() to make the marks consistent.
2243 
2244        TODO: We could use in-place communication with a different SF
2245        We use MPI_SUM for the Reduce, and check the result against the rootdegree. If sum >= rootdegree+1, then the edge has
2246        already been marked. If not, it might have been handled on the process in this round, but we add it anyway.
2247 
2248        In order to update the queue with the new edges from the label communication, we use BcastAnOp(MPI_SUM), so that new
2249        values will have 1+0=1 and old values will have 1+1=2. Loop over these, resetting the values to 1, and adding any new
2250        edge to the queue.
2251     */
2252     PetscCall(DMLabelPropagateInit_Internal(label, pointSF, valArray));
2253     PetscCall(PetscSFReduceBegin(pointSF, MPIU_INT, valArray, valArray, MPI_MAX));
2254     PetscCall(PetscSFReduceEnd(pointSF, MPIU_INT, valArray, valArray, MPI_MAX));
2255     PetscCall(PetscSFBcastBegin(pointSF, MPIU_INT, valArray, valArray, MPI_REPLACE));
2256     PetscCall(PetscSFBcastEnd(pointSF, MPIU_INT, valArray, valArray, MPI_REPLACE));
2257     PetscCall(DMLabelPropagateFini_Internal(label, pointSF, valArray, markPoint, ctx));
2258   }
2259   PetscFunctionReturn(PETSC_SUCCESS);
2260 }
2261 
2262 /*@
2263   DMLabelConvertToSection - Make a `PetscSection`/`IS` pair that encodes the label
2264 
2265   Not Collective
2266 
2267   Input Parameter:
2268 . label - the `DMLabel`
2269 
2270   Output Parameters:
2271 + section - the section giving offsets for each stratum
2272 - is      - An `IS` containing all the label points
2273 
2274   Level: developer
2275 
2276 .seealso: `DMLabel`, `DM`, `DMLabelDistribute()`
2277 @*/
2278 PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
2279 {
2280   IS              vIS;
2281   const PetscInt *values;
2282   PetscInt       *points;
2283   PetscInt        nV, vS = 0, vE = 0, v, N;
2284 
2285   PetscFunctionBegin;
2286   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
2287   PetscCall(DMLabelGetNumValues(label, &nV));
2288   PetscCall(DMLabelGetValueIS(label, &vIS));
2289   PetscCall(ISGetIndices(vIS, &values));
2290   if (nV) {
2291     vS = values[0];
2292     vE = values[0] + 1;
2293   }
2294   for (v = 1; v < nV; ++v) {
2295     vS = PetscMin(vS, values[v]);
2296     vE = PetscMax(vE, values[v] + 1);
2297   }
2298   PetscCall(PetscSectionCreate(PETSC_COMM_SELF, section));
2299   PetscCall(PetscSectionSetChart(*section, vS, vE));
2300   for (v = 0; v < nV; ++v) {
2301     PetscInt n;
2302 
2303     PetscCall(DMLabelGetStratumSize(label, values[v], &n));
2304     PetscCall(PetscSectionSetDof(*section, values[v], n));
2305   }
2306   PetscCall(PetscSectionSetUp(*section));
2307   PetscCall(PetscSectionGetStorageSize(*section, &N));
2308   PetscCall(PetscMalloc1(N, &points));
2309   for (v = 0; v < nV; ++v) {
2310     IS              is;
2311     const PetscInt *spoints;
2312     PetscInt        dof, off, p;
2313 
2314     PetscCall(PetscSectionGetDof(*section, values[v], &dof));
2315     PetscCall(PetscSectionGetOffset(*section, values[v], &off));
2316     PetscCall(DMLabelGetStratumIS(label, values[v], &is));
2317     PetscCall(ISGetIndices(is, &spoints));
2318     for (p = 0; p < dof; ++p) points[off + p] = spoints[p];
2319     PetscCall(ISRestoreIndices(is, &spoints));
2320     PetscCall(ISDestroy(&is));
2321   }
2322   PetscCall(ISRestoreIndices(vIS, &values));
2323   PetscCall(ISDestroy(&vIS));
2324   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is));
2325   PetscFunctionReturn(PETSC_SUCCESS);
2326 }
2327 
2328 /*@C
2329   DMLabelRegister - Adds a new label component implementation
2330 
2331   Not Collective
2332 
2333   Input Parameters:
2334 + name        - The name of a new user-defined creation routine
2335 - create_func - The creation routine itself
2336 
2337   Notes:
2338   `DMLabelRegister()` may be called multiple times to add several user-defined labels
2339 
2340   Example Usage:
2341 .vb
2342   DMLabelRegister("my_label", MyLabelCreate);
2343 .ve
2344 
2345   Then, your label type can be chosen with the procedural interface via
2346 .vb
2347   DMLabelCreate(MPI_Comm, DMLabel *);
2348   DMLabelSetType(DMLabel, "my_label");
2349 .ve
2350   or at runtime via the option
2351 .vb
2352   -dm_label_type my_label
2353 .ve
2354 
2355   Level: advanced
2356 
2357 .seealso: `DMLabel`, `DM`, `DMLabelType`, `DMLabelRegisterAll()`, `DMLabelRegisterDestroy()`
2358 @*/
2359 PetscErrorCode DMLabelRegister(const char name[], PetscErrorCode (*create_func)(DMLabel))
2360 {
2361   PetscFunctionBegin;
2362   PetscCall(DMInitializePackage());
2363   PetscCall(PetscFunctionListAdd(&DMLabelList, name, create_func));
2364   PetscFunctionReturn(PETSC_SUCCESS);
2365 }
2366 
2367 PETSC_EXTERN PetscErrorCode DMLabelCreate_Concrete(DMLabel);
2368 PETSC_EXTERN PetscErrorCode DMLabelCreate_Ephemeral(DMLabel);
2369 
2370 /*@C
2371   DMLabelRegisterAll - Registers all of the `DMLabel` implementations in the `DM` package.
2372 
2373   Not Collective
2374 
2375   Level: advanced
2376 
2377 .seealso: `DMLabel`, `DM`, `DMRegisterAll()`, `DMLabelRegisterDestroy()`
2378 @*/
2379 PetscErrorCode DMLabelRegisterAll(void)
2380 {
2381   PetscFunctionBegin;
2382   if (DMLabelRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
2383   DMLabelRegisterAllCalled = PETSC_TRUE;
2384 
2385   PetscCall(DMLabelRegister(DMLABELCONCRETE, DMLabelCreate_Concrete));
2386   PetscCall(DMLabelRegister(DMLABELEPHEMERAL, DMLabelCreate_Ephemeral));
2387   PetscFunctionReturn(PETSC_SUCCESS);
2388 }
2389 
2390 /*@C
2391   DMLabelRegisterDestroy - This function destroys the `DMLabel` registry. It is called from `PetscFinalize()`.
2392 
2393   Level: developer
2394 
2395 .seealso: `DMLabel`, `DM`, `PetscInitialize()`
2396 @*/
2397 PetscErrorCode DMLabelRegisterDestroy(void)
2398 {
2399   PetscFunctionBegin;
2400   PetscCall(PetscFunctionListDestroy(&DMLabelList));
2401   DMLabelRegisterAllCalled = PETSC_FALSE;
2402   PetscFunctionReturn(PETSC_SUCCESS);
2403 }
2404 
2405 /*@
2406   DMLabelSetType - Sets the particular implementation for a label.
2407 
2408   Collective
2409 
2410   Input Parameters:
2411 + label  - The label
2412 - method - The name of the label type
2413 
2414   Options Database Key:
2415 . -dm_label_type <type> - Sets the label type; use -help for a list of available types or see `DMLabelType`
2416 
2417   Level: intermediate
2418 
2419 .seealso: `DMLabel`, `DM`, `DMLabelGetType()`, `DMLabelCreate()`
2420 @*/
2421 PetscErrorCode DMLabelSetType(DMLabel label, DMLabelType method)
2422 {
2423   PetscErrorCode (*r)(DMLabel);
2424   PetscBool match;
2425 
2426   PetscFunctionBegin;
2427   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
2428   PetscCall(PetscObjectTypeCompare((PetscObject)label, method, &match));
2429   if (match) PetscFunctionReturn(PETSC_SUCCESS);
2430 
2431   PetscCall(DMLabelRegisterAll());
2432   PetscCall(PetscFunctionListFind(DMLabelList, method, &r));
2433   PetscCheck(r, PetscObjectComm((PetscObject)label), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DMLabel type: %s", method);
2434 
2435   PetscTryTypeMethod(label, destroy);
2436   PetscCall(PetscMemzero(label->ops, sizeof(*label->ops)));
2437   PetscCall(PetscObjectChangeTypeName((PetscObject)label, method));
2438   PetscCall((*r)(label));
2439   PetscFunctionReturn(PETSC_SUCCESS);
2440 }
2441 
2442 /*@
2443   DMLabelGetType - Gets the type name (as a string) from the label.
2444 
2445   Not Collective
2446 
2447   Input Parameter:
2448 . label - The `DMLabel`
2449 
2450   Output Parameter:
2451 . type - The `DMLabel` type name
2452 
2453   Level: intermediate
2454 
2455 .seealso: `DMLabel`, `DM`, `DMLabelSetType()`, `DMLabelCreate()`
2456 @*/
2457 PetscErrorCode DMLabelGetType(DMLabel label, DMLabelType *type)
2458 {
2459   PetscFunctionBegin;
2460   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
2461   PetscAssertPointer(type, 2);
2462   PetscCall(DMLabelRegisterAll());
2463   *type = ((PetscObject)label)->type_name;
2464   PetscFunctionReturn(PETSC_SUCCESS);
2465 }
2466 
2467 static PetscErrorCode DMLabelInitialize_Concrete(DMLabel label)
2468 {
2469   PetscFunctionBegin;
2470   label->ops->view         = DMLabelView_Concrete;
2471   label->ops->setup        = NULL;
2472   label->ops->duplicate    = DMLabelDuplicate_Concrete;
2473   label->ops->getstratumis = DMLabelGetStratumIS_Concrete;
2474   PetscFunctionReturn(PETSC_SUCCESS);
2475 }
2476 
2477 PETSC_EXTERN PetscErrorCode DMLabelCreate_Concrete(DMLabel label)
2478 {
2479   PetscFunctionBegin;
2480   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
2481   PetscCall(DMLabelInitialize_Concrete(label));
2482   PetscFunctionReturn(PETSC_SUCCESS);
2483 }
2484 
2485 /*@
2486   PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
2487   the local section and an `PetscSF` describing the section point overlap.
2488 
2489   Collective
2490 
2491   Input Parameters:
2492 + s                  - The `PetscSection` for the local field layout
2493 . sf                 - The `PetscSF` describing parallel layout of the section points
2494 . includeConstraints - By default this is `PETSC_FALSE`, meaning that the global field vector will not possess constrained dofs
2495 . label              - The label specifying the points
2496 - labelValue         - The label stratum specifying the points
2497 
2498   Output Parameter:
2499 . gsection - The `PetscSection` for the global field layout
2500 
2501   Level: developer
2502 
2503   Note:
2504   This gives negative sizes and offsets to points not owned by this process
2505 
2506 .seealso: `DMLabel`, `DM`, `PetscSectionCreate()`
2507 @*/
2508 PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
2509 {
2510   PetscInt *neg = NULL, *tmpOff = NULL;
2511   PetscInt  pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
2512 
2513   PetscFunctionBegin;
2514   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
2515   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
2516   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4);
2517   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)s), gsection));
2518   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
2519   PetscCall(PetscSectionSetChart(*gsection, pStart, pEnd));
2520   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
2521   if (nroots >= 0) {
2522     PetscCheck(nroots >= pEnd - pStart, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %" PetscInt_FMT " < %" PetscInt_FMT " section size", nroots, pEnd - pStart);
2523     PetscCall(PetscCalloc1(nroots, &neg));
2524     if (nroots > pEnd - pStart) {
2525       PetscCall(PetscCalloc1(nroots, &tmpOff));
2526     } else {
2527       tmpOff = &(*gsection)->atlasDof[-pStart];
2528     }
2529   }
2530   /* Mark ghost points with negative dof */
2531   for (p = pStart; p < pEnd; ++p) {
2532     PetscInt value;
2533 
2534     PetscCall(DMLabelGetValue(label, p, &value));
2535     if (value != labelValue) continue;
2536     PetscCall(PetscSectionGetDof(s, p, &dof));
2537     PetscCall(PetscSectionSetDof(*gsection, p, dof));
2538     PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
2539     if (!includeConstraints && cdof > 0) PetscCall(PetscSectionSetConstraintDof(*gsection, p, cdof));
2540     if (neg) neg[p] = -(dof + 1);
2541   }
2542   PetscCall(PetscSectionSetUpBC(*gsection));
2543   if (nroots >= 0) {
2544     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE));
2545     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE));
2546     if (nroots > pEnd - pStart) {
2547       for (p = pStart; p < pEnd; ++p) {
2548         if (tmpOff[p] < 0) (*gsection)->atlasDof[p - pStart] = tmpOff[p];
2549       }
2550     }
2551   }
2552   /* Calculate new sizes, get process offset, and calculate point offsets */
2553   for (p = 0, off = 0; p < pEnd - pStart; ++p) {
2554     cdof                     = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
2555     (*gsection)->atlasOff[p] = off;
2556     off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p] - cdof : 0;
2557   }
2558   PetscCallMPI(MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)s)));
2559   globalOff -= off;
2560   for (p = 0, off = 0; p < pEnd - pStart; ++p) {
2561     (*gsection)->atlasOff[p] += globalOff;
2562     if (neg) neg[p] = -((*gsection)->atlasOff[p] + 1);
2563   }
2564   /* Put in negative offsets for ghost points */
2565   if (nroots >= 0) {
2566     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE));
2567     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE));
2568     if (nroots > pEnd - pStart) {
2569       for (p = pStart; p < pEnd; ++p) {
2570         if (tmpOff[p] < 0) (*gsection)->atlasOff[p - pStart] = tmpOff[p];
2571       }
2572     }
2573   }
2574   if (nroots >= 0 && nroots > pEnd - pStart) PetscCall(PetscFree(tmpOff));
2575   PetscCall(PetscFree(neg));
2576   PetscFunctionReturn(PETSC_SUCCESS);
2577 }
2578 
2579 typedef struct _n_PetscSectionSym_Label {
2580   DMLabel              label;
2581   PetscCopyMode       *modes;
2582   PetscInt            *sizes;
2583   const PetscInt    ***perms;
2584   const PetscScalar ***rots;
2585   PetscInt (*minMaxOrients)[2];
2586   PetscInt numStrata; /* numStrata is only increasing, functions as a state */
2587 } PetscSectionSym_Label;
2588 
2589 static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
2590 {
2591   PetscInt               i, j;
2592   PetscSectionSym_Label *sl = (PetscSectionSym_Label *)sym->data;
2593 
2594   PetscFunctionBegin;
2595   for (i = 0; i <= sl->numStrata; i++) {
2596     if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
2597       for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
2598         if (sl->perms[i]) PetscCall(PetscFree(sl->perms[i][j]));
2599         if (sl->rots[i]) PetscCall(PetscFree(sl->rots[i][j]));
2600       }
2601       if (sl->perms[i]) {
2602         const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
2603 
2604         PetscCall(PetscFree(perms));
2605       }
2606       if (sl->rots[i]) {
2607         const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
2608 
2609         PetscCall(PetscFree(rots));
2610       }
2611     }
2612   }
2613   PetscCall(PetscFree5(sl->modes, sl->sizes, sl->perms, sl->rots, sl->minMaxOrients));
2614   PetscCall(DMLabelDestroy(&sl->label));
2615   sl->numStrata = 0;
2616   PetscFunctionReturn(PETSC_SUCCESS);
2617 }
2618 
2619 static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
2620 {
2621   PetscFunctionBegin;
2622   PetscCall(PetscSectionSymLabelReset(sym));
2623   PetscCall(PetscFree(sym->data));
2624   PetscFunctionReturn(PETSC_SUCCESS);
2625 }
2626 
2627 static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
2628 {
2629   PetscSectionSym_Label *sl = (PetscSectionSym_Label *)sym->data;
2630   PetscBool              isAscii;
2631   DMLabel                label = sl->label;
2632   const char            *name;
2633 
2634   PetscFunctionBegin;
2635   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isAscii));
2636   if (isAscii) {
2637     PetscInt          i, j, k;
2638     PetscViewerFormat format;
2639 
2640     PetscCall(PetscViewerGetFormat(viewer, &format));
2641     if (label) {
2642       PetscCall(PetscViewerGetFormat(viewer, &format));
2643       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
2644         PetscCall(PetscViewerASCIIPushTab(viewer));
2645         PetscCall(DMLabelView(label, viewer));
2646         PetscCall(PetscViewerASCIIPopTab(viewer));
2647       } else {
2648         PetscCall(PetscObjectGetName((PetscObject)sl->label, &name));
2649         PetscCall(PetscViewerASCIIPrintf(viewer, "  Label '%s'\n", name));
2650       }
2651     } else {
2652       PetscCall(PetscViewerASCIIPrintf(viewer, "No label given\n"));
2653     }
2654     PetscCall(PetscViewerASCIIPushTab(viewer));
2655     for (i = 0; i <= sl->numStrata; i++) {
2656       PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
2657 
2658       if (!(sl->perms[i] || sl->rots[i])) {
2659         PetscCall(PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %" PetscInt_FMT " (%" PetscInt_FMT " dofs per point): no symmetries\n", value, sl->sizes[i]));
2660       } else {
2661         PetscCall(PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %" PetscInt_FMT " (%" PetscInt_FMT " dofs per point):\n", value, sl->sizes[i]));
2662         PetscCall(PetscViewerASCIIPushTab(viewer));
2663         PetscCall(PetscViewerASCIIPrintf(viewer, "Orientation range: [%" PetscInt_FMT ", %" PetscInt_FMT ")\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]));
2664         if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
2665           PetscCall(PetscViewerASCIIPushTab(viewer));
2666           for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
2667             if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
2668               PetscCall(PetscViewerASCIIPrintf(viewer, "Orientation %" PetscInt_FMT ": identity\n", j));
2669             } else {
2670               PetscInt tab;
2671 
2672               PetscCall(PetscViewerASCIIPrintf(viewer, "Orientation %" PetscInt_FMT ":\n", j));
2673               PetscCall(PetscViewerASCIIPushTab(viewer));
2674               PetscCall(PetscViewerASCIIGetTab(viewer, &tab));
2675               if (sl->perms[i] && sl->perms[i][j]) {
2676                 PetscCall(PetscViewerASCIIPrintf(viewer, "Permutation:"));
2677                 PetscCall(PetscViewerASCIISetTab(viewer, 0));
2678                 for (k = 0; k < sl->sizes[i]; k++) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, sl->perms[i][j][k]));
2679                 PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
2680                 PetscCall(PetscViewerASCIISetTab(viewer, tab));
2681               }
2682               if (sl->rots[i] && sl->rots[i][j]) {
2683                 PetscCall(PetscViewerASCIIPrintf(viewer, "Rotations:  "));
2684                 PetscCall(PetscViewerASCIISetTab(viewer, 0));
2685 #if defined(PETSC_USE_COMPLEX)
2686                 for (k = 0; k < sl->sizes[i]; k++) PetscCall(PetscViewerASCIIPrintf(viewer, " %+g+i*%+g", (double)PetscRealPart(sl->rots[i][j][k]), (double)PetscImaginaryPart(sl->rots[i][j][k])));
2687 #else
2688                 for (k = 0; k < sl->sizes[i]; k++) PetscCall(PetscViewerASCIIPrintf(viewer, " %+g", (double)sl->rots[i][j][k]));
2689 #endif
2690                 PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
2691                 PetscCall(PetscViewerASCIISetTab(viewer, tab));
2692               }
2693               PetscCall(PetscViewerASCIIPopTab(viewer));
2694             }
2695           }
2696           PetscCall(PetscViewerASCIIPopTab(viewer));
2697         }
2698         PetscCall(PetscViewerASCIIPopTab(viewer));
2699       }
2700     }
2701     PetscCall(PetscViewerASCIIPopTab(viewer));
2702   }
2703   PetscFunctionReturn(PETSC_SUCCESS);
2704 }
2705 
2706 /*@
2707   PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
2708 
2709   Logically
2710 
2711   Input Parameters:
2712 + sym   - the section symmetries
2713 - label - the `DMLabel` describing the types of points
2714 
2715   Level: developer:
2716 
2717 .seealso: `DMLabel`, `DM`, `PetscSectionSymLabelSetStratum()`, `PetscSectionSymCreateLabel()`, `PetscSectionGetPointSyms()`
2718 @*/
2719 PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
2720 {
2721   PetscSectionSym_Label *sl;
2722 
2723   PetscFunctionBegin;
2724   PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1);
2725   sl = (PetscSectionSym_Label *)sym->data;
2726   if (sl->label && sl->label != label) PetscCall(PetscSectionSymLabelReset(sym));
2727   if (label) {
2728     sl->label = label;
2729     PetscCall(PetscObjectReference((PetscObject)label));
2730     PetscCall(DMLabelGetNumValues(label, &sl->numStrata));
2731     PetscCall(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));
2732     PetscCall(PetscMemzero((void *)sl->modes, (sl->numStrata + 1) * sizeof(PetscCopyMode)));
2733     PetscCall(PetscMemzero((void *)sl->sizes, (sl->numStrata + 1) * sizeof(PetscInt)));
2734     PetscCall(PetscMemzero((void *)sl->perms, (sl->numStrata + 1) * sizeof(const PetscInt **)));
2735     PetscCall(PetscMemzero((void *)sl->rots, (sl->numStrata + 1) * sizeof(const PetscScalar **)));
2736     PetscCall(PetscMemzero((void *)sl->minMaxOrients, (sl->numStrata + 1) * sizeof(PetscInt[2])));
2737   }
2738   PetscFunctionReturn(PETSC_SUCCESS);
2739 }
2740 
2741 /*@C
2742   PetscSectionSymLabelGetStratum - get the symmetries for the orientations of a stratum
2743 
2744   Logically Collective
2745 
2746   Input Parameters:
2747 + sym     - the section symmetries
2748 - stratum - the stratum value in the label that we are assigning symmetries for
2749 
2750   Output Parameters:
2751 + size      - the number of dofs for points in the `stratum` of the label
2752 . minOrient - the smallest orientation for a point in this `stratum`
2753 . maxOrient - one greater than the largest orientation for a ppoint in this `stratum` (i.e., orientations are in the range [`minOrient`, `maxOrient`))
2754 . perms     - `NULL` if there are no permutations, or (`maxOrient` - `minOrient`) permutations, one for each orientation.  A `NULL` permutation is the identity
2755 - 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
2756 
2757   Level: developer
2758 
2759 .seealso: `DMLabel`, `DM`, `PetscSectionSymLabelSetStratum()`, `PetscSectionSymCreate()`, `PetscSectionSetSym()`, `PetscSectionGetPointSyms()`, `PetscSectionSymCreateLabel()`
2760 @*/
2761 PetscErrorCode PetscSectionSymLabelGetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt *size, PetscInt *minOrient, PetscInt *maxOrient, const PetscInt ***perms, const PetscScalar ***rots)
2762 {
2763   PetscSectionSym_Label *sl;
2764   const char            *name;
2765   PetscInt               i;
2766 
2767   PetscFunctionBegin;
2768   PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1);
2769   sl = (PetscSectionSym_Label *)sym->data;
2770   PetscCheck(sl->label, PetscObjectComm((PetscObject)sym), PETSC_ERR_ARG_WRONGSTATE, "No label set yet");
2771   for (i = 0; i <= sl->numStrata; i++) {
2772     PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
2773 
2774     if (stratum == value) break;
2775   }
2776   PetscCall(PetscObjectGetName((PetscObject)sl->label, &name));
2777   PetscCheck(i <= sl->numStrata, PetscObjectComm((PetscObject)sym), PETSC_ERR_ARG_OUTOFRANGE, "Stratum %" PetscInt_FMT " not found in label %s", stratum, name);
2778   if (size) {
2779     PetscAssertPointer(size, 3);
2780     *size = sl->sizes[i];
2781   }
2782   if (minOrient) {
2783     PetscAssertPointer(minOrient, 4);
2784     *minOrient = sl->minMaxOrients[i][0];
2785   }
2786   if (maxOrient) {
2787     PetscAssertPointer(maxOrient, 5);
2788     *maxOrient = sl->minMaxOrients[i][1];
2789   }
2790   if (perms) {
2791     PetscAssertPointer(perms, 6);
2792     *perms = PetscSafePointerPlusOffset(sl->perms[i], sl->minMaxOrients[i][0]);
2793   }
2794   if (rots) {
2795     PetscAssertPointer(rots, 7);
2796     *rots = PetscSafePointerPlusOffset(sl->rots[i], sl->minMaxOrients[i][0]);
2797   }
2798   PetscFunctionReturn(PETSC_SUCCESS);
2799 }
2800 
2801 /*@C
2802   PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
2803 
2804   Logically
2805 
2806   Input Parameters:
2807 + sym       - the section symmetries
2808 . stratum   - the stratum value in the label that we are assigning symmetries for
2809 . size      - the number of dofs for points in the `stratum` of the label
2810 . minOrient - the smallest orientation for a point in this `stratum`
2811 . maxOrient - one greater than the largest orientation for a point in this `stratum` (i.e., orientations are in the range [`minOrient`, `maxOrient`))
2812 . mode      - how `sym` should copy the `perms` and `rots` arrays
2813 . perms     - `NULL` if there are no permutations, or (`maxOrient` - `minOrient`) permutations, one for each orientation.  A `NULL` permutation is the identity
2814 - 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
2815 
2816   Level: developer
2817 
2818 .seealso: `DMLabel`, `DM`, `PetscSectionSymLabelGetStratum()`, `PetscSectionSymCreate()`, `PetscSectionSetSym()`, `PetscSectionGetPointSyms()`, `PetscSectionSymCreateLabel()`
2819 @*/
2820 PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
2821 {
2822   PetscSectionSym_Label *sl;
2823   const char            *name;
2824   PetscInt               i, j, k;
2825 
2826   PetscFunctionBegin;
2827   PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1);
2828   sl = (PetscSectionSym_Label *)sym->data;
2829   PetscCheck(sl->label, PetscObjectComm((PetscObject)sym), PETSC_ERR_ARG_WRONGSTATE, "No label set yet");
2830   for (i = 0; i <= sl->numStrata; i++) {
2831     PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
2832 
2833     if (stratum == value) break;
2834   }
2835   PetscCall(PetscObjectGetName((PetscObject)sl->label, &name));
2836   PetscCheck(i <= sl->numStrata, PetscObjectComm((PetscObject)sym), PETSC_ERR_ARG_OUTOFRANGE, "Stratum %" PetscInt_FMT " not found in label %s", stratum, name);
2837   sl->sizes[i]            = size;
2838   sl->modes[i]            = mode;
2839   sl->minMaxOrients[i][0] = minOrient;
2840   sl->minMaxOrients[i][1] = maxOrient;
2841   if (mode == PETSC_COPY_VALUES) {
2842     if (perms) {
2843       PetscInt **ownPerms;
2844 
2845       PetscCall(PetscCalloc1(maxOrient - minOrient, &ownPerms));
2846       for (j = 0; j < maxOrient - minOrient; j++) {
2847         if (perms[j]) {
2848           PetscCall(PetscMalloc1(size, &ownPerms[j]));
2849           for (k = 0; k < size; k++) ownPerms[j][k] = perms[j][k];
2850         }
2851       }
2852       sl->perms[i] = (const PetscInt **)&ownPerms[-minOrient];
2853     }
2854     if (rots) {
2855       PetscScalar **ownRots;
2856 
2857       PetscCall(PetscCalloc1(maxOrient - minOrient, &ownRots));
2858       for (j = 0; j < maxOrient - minOrient; j++) {
2859         if (rots[j]) {
2860           PetscCall(PetscMalloc1(size, &ownRots[j]));
2861           for (k = 0; k < size; k++) ownRots[j][k] = rots[j][k];
2862         }
2863       }
2864       sl->rots[i] = (const PetscScalar **)&ownRots[-minOrient];
2865     }
2866   } else {
2867     sl->perms[i] = PetscSafePointerPlusOffset(perms, -minOrient);
2868     sl->rots[i]  = PetscSafePointerPlusOffset(rots, -minOrient);
2869   }
2870   PetscFunctionReturn(PETSC_SUCCESS);
2871 }
2872 
2873 static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
2874 {
2875   PetscInt               i, j, numStrata;
2876   PetscSectionSym_Label *sl;
2877   DMLabel                label;
2878 
2879   PetscFunctionBegin;
2880   sl        = (PetscSectionSym_Label *)sym->data;
2881   numStrata = sl->numStrata;
2882   label     = sl->label;
2883   for (i = 0; i < numPoints; i++) {
2884     PetscInt point = points[2 * i];
2885     PetscInt ornt  = points[2 * i + 1];
2886 
2887     for (j = 0; j < numStrata; j++) {
2888       if (label->validIS[j]) {
2889         PetscInt k;
2890 
2891         PetscCall(ISLocate(label->points[j], point, &k));
2892         if (k >= 0) break;
2893       } else {
2894         PetscBool has;
2895 
2896         PetscCall(PetscHSetIHas(label->ht[j], point, &has));
2897         if (has) break;
2898       }
2899     }
2900     PetscCheck(!(sl->minMaxOrients[j][1] > sl->minMaxOrients[j][0]) || !(ornt < sl->minMaxOrients[j][0] || ornt >= sl->minMaxOrients[j][1]), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "point %" PetscInt_FMT " orientation %" PetscInt_FMT " not in range [%" PetscInt_FMT ", %" PetscInt_FMT ") for stratum %" PetscInt_FMT, point, ornt, sl->minMaxOrients[j][0], sl->minMaxOrients[j][1],
2901                j < numStrata ? label->stratumValues[j] : label->defaultValue);
2902     if (perms) perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;
2903     if (rots) rots[i] = sl->rots[j] ? sl->rots[j][ornt] : NULL;
2904   }
2905   PetscFunctionReturn(PETSC_SUCCESS);
2906 }
2907 
2908 static PetscErrorCode PetscSectionSymCopy_Label(PetscSectionSym sym, PetscSectionSym nsym)
2909 {
2910   PetscSectionSym_Label *sl = (PetscSectionSym_Label *)nsym->data;
2911   IS                     valIS;
2912   const PetscInt        *values;
2913   PetscInt               Nv, v;
2914 
2915   PetscFunctionBegin;
2916   PetscCall(DMLabelGetNumValues(sl->label, &Nv));
2917   PetscCall(DMLabelGetValueIS(sl->label, &valIS));
2918   PetscCall(ISGetIndices(valIS, &values));
2919   for (v = 0; v < Nv; ++v) {
2920     const PetscInt      val = values[v];
2921     PetscInt            size, minOrient, maxOrient;
2922     const PetscInt    **perms;
2923     const PetscScalar **rots;
2924 
2925     PetscCall(PetscSectionSymLabelGetStratum(sym, val, &size, &minOrient, &maxOrient, &perms, &rots));
2926     PetscCall(PetscSectionSymLabelSetStratum(nsym, val, size, minOrient, maxOrient, PETSC_COPY_VALUES, perms, rots));
2927   }
2928   PetscCall(ISDestroy(&valIS));
2929   PetscFunctionReturn(PETSC_SUCCESS);
2930 }
2931 
2932 static PetscErrorCode PetscSectionSymDistribute_Label(PetscSectionSym sym, PetscSF migrationSF, PetscSectionSym *dsym)
2933 {
2934   PetscSectionSym_Label *sl = (PetscSectionSym_Label *)sym->data;
2935   DMLabel                dlabel;
2936 
2937   PetscFunctionBegin;
2938   PetscCall(DMLabelDistribute(sl->label, migrationSF, &dlabel));
2939   PetscCall(PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)sym), dlabel, dsym));
2940   PetscCall(DMLabelDestroy(&dlabel));
2941   PetscCall(PetscSectionSymCopy(sym, *dsym));
2942   PetscFunctionReturn(PETSC_SUCCESS);
2943 }
2944 
2945 PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
2946 {
2947   PetscSectionSym_Label *sl;
2948 
2949   PetscFunctionBegin;
2950   PetscCall(PetscNew(&sl));
2951   sym->ops->getpoints  = PetscSectionSymGetPoints_Label;
2952   sym->ops->distribute = PetscSectionSymDistribute_Label;
2953   sym->ops->copy       = PetscSectionSymCopy_Label;
2954   sym->ops->view       = PetscSectionSymView_Label;
2955   sym->ops->destroy    = PetscSectionSymDestroy_Label;
2956   sym->data            = (void *)sl;
2957   PetscFunctionReturn(PETSC_SUCCESS);
2958 }
2959 
2960 /*@
2961   PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
2962 
2963   Collective
2964 
2965   Input Parameters:
2966 + comm  - the MPI communicator for the new symmetry
2967 - label - the label defining the strata
2968 
2969   Output Parameter:
2970 . sym - the section symmetries
2971 
2972   Level: developer
2973 
2974 .seealso: `DMLabel`, `DM`, `PetscSectionSymCreate()`, `PetscSectionSetSym()`, `PetscSectionGetSym()`, `PetscSectionSymLabelSetStratum()`, `PetscSectionGetPointSyms()`
2975 @*/
2976 PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
2977 {
2978   PetscFunctionBegin;
2979   PetscCall(DMInitializePackage());
2980   PetscCall(PetscSectionSymCreate(comm, sym));
2981   PetscCall(PetscSectionSymSetType(*sym, PETSCSECTIONSYMLABEL));
2982   PetscCall(PetscSectionSymLabelSetLabel(*sym, label));
2983   PetscFunctionReturn(PETSC_SUCCESS);
2984 }
2985