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