xref: /petsc/src/dm/label/dmlabel.c (revision fb2dfec1ce55457db93eeb2a8a1569a386563980)
1 #include <petsc/private/dmlabelimpl.h>   /*I      "petscdmlabel.h"   I*/
2 #include <petsc/private/isimpl.h>        /*I      "petscis.h"        I*/
3 #include <petscsf.h>
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "DMLabelCreate"
7 /*@C
8   DMLabelCreate - Create a DMLabel object, which is a multimap
9 
10   Input parameter:
11 . name - The label name
12 
13   Output parameter:
14 . label - The DMLabel
15 
16   Level: beginner
17 
18 .seealso: DMLabelDestroy()
19 @*/
20 PetscErrorCode DMLabelCreate(const char name[], DMLabel *label)
21 {
22   PetscErrorCode ierr;
23 
24   PetscFunctionBegin;
25   ierr = PetscNew(label);CHKERRQ(ierr);
26   ierr = PetscStrallocpy(name, &(*label)->name);CHKERRQ(ierr);
27 
28   (*label)->refct          = 1;
29   (*label)->state          = -1;
30   (*label)->numStrata      = 0;
31   (*label)->defaultValue   = -1;
32   (*label)->stratumValues  = NULL;
33   (*label)->validIS        = NULL;
34   (*label)->stratumSizes   = NULL;
35   (*label)->points         = NULL;
36   (*label)->ht             = NULL;
37   (*label)->pStart         = -1;
38   (*label)->pEnd           = -1;
39   (*label)->bt             = NULL;
40   PetscFunctionReturn(0);
41 }
42 
43 #undef __FUNCT__
44 #define __FUNCT__ "DMLabelMakeValid_Private"
45 /*
46   DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format
47 
48   Input parameter:
49 + label - The DMLabel
50 - v - The stratum value
51 
52   Output parameter:
53 . label - The DMLabel with stratum in sorted list format
54 
55   Level: developer
56 
57 .seealso: DMLabelCreate()
58 */
59 static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v)
60 {
61   PetscInt       off;
62   PetscInt       *pointArray;
63   PetscErrorCode ierr;
64 
65   if (label->validIS[v]) return 0;
66   if (v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeValid_Private\n", v);
67   PetscFunctionBegin;
68   PetscHashISize(label->ht[v], label->stratumSizes[v]);
69 
70   ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr);
71   off  = 0;
72   ierr = PetscHashIGetKeys(label->ht[v], &off, pointArray);CHKERRQ(ierr);
73   if (off != label->stratumSizes[v]) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of contributed points %D from value %D should be %D", off, label->stratumValues[v], label->stratumSizes[v]);
74   PetscHashIClear(label->ht[v]);
75   ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr);
76   if (label->bt) {
77     PetscInt p;
78 
79     for (p = 0; p < label->stratumSizes[v]; ++p) {
80       const PetscInt point = pointArray[p];
81 
82       if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, label->pStart, label->pEnd);
83       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
84     }
85   }
86   ierr = ISCreateGeneral(PETSC_COMM_SELF,label->stratumSizes[v],pointArray,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
87   ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
88   label->validIS[v] = PETSC_TRUE;
89   ++label->state;
90   PetscFunctionReturn(0);
91 }
92 
93 #undef __FUNCT__
94 #define __FUNCT__ "DMLabelMakeAllValid_Private"
95 /*
96   DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format
97 
98   Input parameter:
99 . label - The DMLabel
100 
101   Output parameter:
102 . label - The DMLabel with all strata in sorted list format
103 
104   Level: developer
105 
106 .seealso: DMLabelCreate()
107 */
108 static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label)
109 {
110   PetscInt       v;
111   PetscErrorCode ierr;
112 
113   PetscFunctionBegin;
114   for (v = 0; v < label->numStrata; v++){
115     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
116   }
117   PetscFunctionReturn(0);
118 }
119 
120 #undef __FUNCT__
121 #define __FUNCT__ "DMLabelMakeInvalid_Private"
122 /*
123   DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format
124 
125   Input parameter:
126 + label - The DMLabel
127 - v - The stratum value
128 
129   Output parameter:
130 . label - The DMLabel with stratum in hash format
131 
132   Level: developer
133 
134 .seealso: DMLabelCreate()
135 */
136 static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v)
137 {
138   PETSC_UNUSED PetscHashIIter ret, iter;
139   PetscInt                    p;
140   const PetscInt              *points;
141   PetscErrorCode ierr;
142 
143   PetscFunctionBegin;
144   if (!label->validIS[v]) PetscFunctionReturn(0);
145   if (label->points[v]) {
146     ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
147     for (p = 0; p < label->stratumSizes[v]; ++p) PetscHashIPut(label->ht[v], points[p], ret, iter);
148     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
149     ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
150   }
151   label->validIS[v] = PETSC_FALSE;
152   PetscFunctionReturn(0);
153 }
154 
155 #undef __FUNCT__
156 #define __FUNCT__ "DMLabelGetState"
157 PetscErrorCode DMLabelGetState(DMLabel label, PetscObjectState *state)
158 {
159   PetscFunctionBegin;
160   PetscValidPointer(state, 2);
161   *state = label->state;
162   PetscFunctionReturn(0);
163 }
164 
165 #undef __FUNCT__
166 #define __FUNCT__ "DMLabelAddStratum"
167 PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
168 {
169   PetscInt    v, *tmpV, *tmpS;
170   IS         *tmpP;
171   PetscHashI *tmpH;
172   PetscBool  *tmpB;
173   PetscErrorCode ierr;
174 
175   PetscFunctionBegin;
176 
177   for (v = 0; v < label->numStrata; v++) {
178     if (label->stratumValues[v] == value) PetscFunctionReturn(0);
179   }
180   ierr = PetscMalloc1((label->numStrata+1), &tmpV);CHKERRQ(ierr);
181   ierr = PetscMalloc1((label->numStrata+1), &tmpS);CHKERRQ(ierr);
182   ierr = PetscMalloc1((label->numStrata+1), &tmpH);CHKERRQ(ierr);
183   ierr = PetscMalloc1((label->numStrata+1), &tmpP);CHKERRQ(ierr);
184   ierr = PetscMalloc1((label->numStrata+1), &tmpB);CHKERRQ(ierr);
185   for (v = 0; v < label->numStrata; ++v) {
186     tmpV[v] = label->stratumValues[v];
187     tmpS[v] = label->stratumSizes[v];
188     tmpH[v] = label->ht[v];
189     tmpP[v] = label->points[v];
190     tmpB[v] = label->validIS[v];
191   }
192   tmpV[v] = value;
193   tmpS[v] = 0;
194   PetscHashICreate(tmpH[v]);
195   ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&tmpP[v]);CHKERRQ(ierr);
196   tmpB[v] = PETSC_TRUE;
197   ++label->numStrata;
198   ierr = PetscFree(label->stratumValues);CHKERRQ(ierr);
199   ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr);
200   ierr = PetscFree(label->ht);CHKERRQ(ierr);
201   ierr = PetscFree(label->points);CHKERRQ(ierr);
202   ierr = PetscFree(label->validIS);CHKERRQ(ierr);
203   label->stratumValues = tmpV;
204   label->stratumSizes  = tmpS;
205   label->ht            = tmpH;
206   label->points        = tmpP;
207   label->validIS       = tmpB;
208 
209   PetscFunctionReturn(0);
210 }
211 
212 #undef __FUNCT__
213 #define __FUNCT__ "DMLabelGetName"
214 /*@C
215   DMLabelGetName - Return the name of a DMLabel object
216 
217   Input parameter:
218 . label - The DMLabel
219 
220   Output parameter:
221 . name - The label name
222 
223   Level: beginner
224 
225 .seealso: DMLabelCreate()
226 @*/
227 PetscErrorCode DMLabelGetName(DMLabel label, const char **name)
228 {
229   PetscFunctionBegin;
230   PetscValidPointer(name, 2);
231   *name = label->name;
232   PetscFunctionReturn(0);
233 }
234 
235 #undef __FUNCT__
236 #define __FUNCT__ "DMLabelView_Ascii"
237 static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
238 {
239   PetscInt       v;
240   PetscMPIInt    rank;
241   PetscErrorCode ierr;
242 
243   PetscFunctionBegin;
244   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr);
245   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
246   if (label) {
247     ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", label->name);CHKERRQ(ierr);
248     if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);}
249     for (v = 0; v < label->numStrata; ++v) {
250       const PetscInt value = label->stratumValues[v];
251       const PetscInt *points;
252       PetscInt       p;
253 
254       ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
255       for (p = 0; p < label->stratumSizes[v]; ++p) {
256         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr);
257       }
258       ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
259     }
260   }
261   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
262   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
263   PetscFunctionReturn(0);
264 }
265 
266 #undef __FUNCT__
267 #define __FUNCT__ "DMLabelView"
268 /*@C
269   DMLabelView - View the label
270 
271   Input Parameters:
272 + label - The DMLabel
273 - viewer - The PetscViewer
274 
275   Level: intermediate
276 
277 .seealso: DMLabelCreate(), DMLabelDestroy()
278 @*/
279 PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
280 {
281   PetscBool      iascii;
282   PetscErrorCode ierr;
283 
284   PetscFunctionBegin;
285   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
286   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
287   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
288   if (iascii) {
289     ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr);
290   }
291   PetscFunctionReturn(0);
292 }
293 
294 #undef __FUNCT__
295 #define __FUNCT__ "DMLabelDestroy"
296 PetscErrorCode DMLabelDestroy(DMLabel *label)
297 {
298   PetscInt       v;
299   PetscErrorCode ierr;
300 
301   PetscFunctionBegin;
302   if (!(*label)) PetscFunctionReturn(0);
303   if (--(*label)->refct > 0) PetscFunctionReturn(0);
304   ierr = PetscFree((*label)->name);CHKERRQ(ierr);
305   ierr = PetscFree((*label)->stratumValues);CHKERRQ(ierr);
306   ierr = PetscFree((*label)->stratumSizes);CHKERRQ(ierr);
307   for (v = 0; v < (*label)->numStrata; ++v) {ierr = ISDestroy(&((*label)->points[v]));CHKERRQ(ierr);}
308   ierr = PetscFree((*label)->points);CHKERRQ(ierr);
309   ierr = PetscFree((*label)->validIS);CHKERRQ(ierr);
310   if ((*label)->ht) {
311     for (v = 0; v < (*label)->numStrata; ++v) {PetscHashIDestroy((*label)->ht[v]);}
312     ierr = PetscFree((*label)->ht);CHKERRQ(ierr);
313   }
314   ierr = PetscBTDestroy(&(*label)->bt);CHKERRQ(ierr);
315   ierr = PetscFree(*label);CHKERRQ(ierr);
316   PetscFunctionReturn(0);
317 }
318 
319 #undef __FUNCT__
320 #define __FUNCT__ "DMLabelDuplicate"
321 PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
322 {
323   PetscInt       v;
324   PetscErrorCode ierr;
325 
326   PetscFunctionBegin;
327   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
328   ierr = PetscNew(labelnew);CHKERRQ(ierr);
329   ierr = PetscStrallocpy(label->name, &(*labelnew)->name);CHKERRQ(ierr);
330 
331   (*labelnew)->refct        = 1;
332   (*labelnew)->numStrata    = label->numStrata;
333   (*labelnew)->defaultValue = label->defaultValue;
334   if (label->numStrata) {
335     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr);
336     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr);
337     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr);
338     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr);
339     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr);
340     /* Could eliminate unused space here */
341     for (v = 0; v < label->numStrata; ++v) {
342       PetscHashICreate((*labelnew)->ht[v]);
343       (*labelnew)->validIS[v]        = PETSC_TRUE;
344       (*labelnew)->stratumValues[v]  = label->stratumValues[v];
345       (*labelnew)->stratumSizes[v]   = label->stratumSizes[v];
346       ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr);
347       (*labelnew)->points[v]         = label->points[v];
348     }
349   }
350   (*labelnew)->pStart = -1;
351   (*labelnew)->pEnd   = -1;
352   (*labelnew)->bt     = NULL;
353   PetscFunctionReturn(0);
354 }
355 
356 #undef __FUNCT__
357 #define __FUNCT__ "DMLabelCreateIndex"
358 /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
359 PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
360 {
361   PetscInt       v;
362   PetscErrorCode ierr;
363 
364   PetscFunctionBegin;
365   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
366   if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);}
367   label->pStart = pStart;
368   label->pEnd   = pEnd;
369   ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr);
370   ierr = PetscBTMemzero(pEnd - pStart, label->bt);CHKERRQ(ierr);
371   for (v = 0; v < label->numStrata; ++v) {
372     const PetscInt *points;
373     PetscInt       i;
374 
375     ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
376     for (i = 0; i < label->stratumSizes[v]; ++i) {
377       const PetscInt point = points[i];
378 
379       if ((point < pStart) || (point >= pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, pStart, pEnd);
380       ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr);
381     }
382     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
383   }
384   PetscFunctionReturn(0);
385 }
386 
387 #undef __FUNCT__
388 #define __FUNCT__ "DMLabelDestroyIndex"
389 PetscErrorCode DMLabelDestroyIndex(DMLabel label)
390 {
391   PetscErrorCode ierr;
392 
393   PetscFunctionBegin;
394   label->pStart = -1;
395   label->pEnd   = -1;
396   if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);}
397   PetscFunctionReturn(0);
398 }
399 
400 #undef __FUNCT__
401 #define __FUNCT__ "DMLabelHasValue"
402 /*@
403   DMLabelHasValue - Determine whether a label assigns the value to any point
404 
405   Input Parameters:
406 + label - the DMLabel
407 - value - the value
408 
409   Output Parameter:
410 . contains - Flag indicating whether the label maps this value to any point
411 
412   Level: developer
413 
414 .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
415 @*/
416 PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
417 {
418   PetscInt v;
419 
420   PetscFunctionBegin;
421   PetscValidPointer(contains, 3);
422   for (v = 0; v < label->numStrata; ++v) {
423     if (value == label->stratumValues[v]) break;
424   }
425   *contains = (v < label->numStrata ? PETSC_TRUE : PETSC_FALSE);
426   PetscFunctionReturn(0);
427 }
428 
429 #undef __FUNCT__
430 #define __FUNCT__ "DMLabelHasPoint"
431 /*@
432   DMLabelHasPoint - Determine whether a label assigns a value to a point
433 
434   Input Parameters:
435 + label - the DMLabel
436 - point - the point
437 
438   Output Parameter:
439 . contains - Flag indicating whether the label maps this point to a value
440 
441   Note: The user must call DMLabelCreateIndex() before this function.
442 
443   Level: developer
444 
445 .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
446 @*/
447 PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
448 {
449   PetscErrorCode ierr;
450 
451   PetscFunctionBeginHot;
452   PetscValidPointer(contains, 3);
453   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
454 #if defined(PETSC_USE_DEBUG)
455   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
456   if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, label->pStart, label->pEnd);
457 #endif
458   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
459   PetscFunctionReturn(0);
460 }
461 
462 #undef __FUNCT__
463 #define __FUNCT__ "DMLabelStratumHasPoint"
464 /*@
465   DMLabelStratumHasPoint - Return true if the stratum contains a point
466 
467   Input Parameters:
468 + label - the DMLabel
469 . value - the stratum value
470 - point - the point
471 
472   Output Parameter:
473 . contains - true if the stratum contains the point
474 
475   Level: intermediate
476 
477 .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
478 @*/
479 PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
480 {
481   PetscInt       v;
482   PetscErrorCode ierr;
483 
484   PetscFunctionBegin;
485   PetscValidPointer(contains, 4);
486   *contains = PETSC_FALSE;
487   for (v = 0; v < label->numStrata; ++v) {
488     if (label->stratumValues[v] == value) {
489       if (label->validIS[v]) {
490         PetscInt i;
491 
492         ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr);
493         if (i >= 0) {
494           *contains = PETSC_TRUE;
495           break;
496         }
497       } else {
498         PetscBool has;
499 
500         PetscHashIHasKey(label->ht[v], point, has);
501         if (has) {
502           *contains = PETSC_TRUE;
503           break;
504         }
505       }
506     }
507   }
508   PetscFunctionReturn(0);
509 }
510 
511 #undef __FUNCT__
512 #define __FUNCT__ "DMLabelGetDefaultValue"
513 /*
514   DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
515   When a label is created, it is initialized to -1.
516 
517   Input parameter:
518 . label - a DMLabel object
519 
520   Output parameter:
521 . defaultValue - the default value
522 
523   Level: beginner
524 
525 .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
526 */
527 PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
528 {
529   PetscFunctionBegin;
530   *defaultValue = label->defaultValue;
531   PetscFunctionReturn(0);
532 }
533 
534 #undef __FUNCT__
535 #define __FUNCT__ "DMLabelSetDefaultValue"
536 /*
537   DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
538   When a label is created, it is initialized to -1.
539 
540   Input parameter:
541 . label - a DMLabel object
542 
543   Output parameter:
544 . defaultValue - the default value
545 
546   Level: beginner
547 
548 .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
549 */
550 PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
551 {
552   PetscFunctionBegin;
553   label->defaultValue = defaultValue;
554   PetscFunctionReturn(0);
555 }
556 
557 #undef __FUNCT__
558 #define __FUNCT__ "DMLabelGetValue"
559 /*@
560   DMLabelGetValue - Return the value a label assigns to a point, or the label's default value (which is initially -1, and can be changed with DMLabelSetDefaultValue())
561 
562   Input Parameters:
563 + label - the DMLabel
564 - point - the point
565 
566   Output Parameter:
567 . value - The point value, or -1
568 
569   Level: intermediate
570 
571 .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
572 @*/
573 PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
574 {
575   PetscInt       v;
576   PetscErrorCode ierr;
577 
578   PetscFunctionBegin;
579   PetscValidPointer(value, 3);
580   *value = label->defaultValue;
581   for (v = 0; v < label->numStrata; ++v) {
582     if (label->validIS[v]) {
583       PetscInt i;
584 
585       ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr);
586       if (i >= 0) {
587         *value = label->stratumValues[v];
588         break;
589       }
590     } else {
591       PetscBool has;
592 
593       PetscHashIHasKey(label->ht[v], point, has);
594       if (has) {
595         *value = label->stratumValues[v];
596         break;
597       }
598     }
599   }
600   PetscFunctionReturn(0);
601 }
602 
603 #undef __FUNCT__
604 #define __FUNCT__ "DMLabelSetValue"
605 /*@
606   DMLabelSetValue - Set the value a label assigns to a point.  If the value is the same as the label's default value (which is initially -1, and can be changed with DMLabelSetDefaultValue() to somethingg different), then this function will do nothing.
607 
608   Input Parameters:
609 + label - the DMLabel
610 . point - the point
611 - value - The point value
612 
613   Level: intermediate
614 
615 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
616 @*/
617 PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
618 {
619   PETSC_UNUSED PetscHashIIter iter, ret;
620   PetscInt                    v;
621   PetscErrorCode              ierr;
622 
623   PetscFunctionBegin;
624   /* Find, or add, label value */
625   if (value == label->defaultValue) PetscFunctionReturn(0);
626   for (v = 0; v < label->numStrata; ++v) {
627     if (label->stratumValues[v] == value) break;
628   }
629   /* Create new table */
630   if (v >= label->numStrata) {ierr = DMLabelAddStratum(label, value);CHKERRQ(ierr);}
631   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
632   /* Set key */
633   PetscHashIPut(label->ht[v], point, ret, iter);
634   PetscFunctionReturn(0);
635 }
636 
637 #undef __FUNCT__
638 #define __FUNCT__ "DMLabelClearValue"
639 /*@
640   DMLabelClearValue - Clear the value a label assigns to a point
641 
642   Input Parameters:
643 + label - the DMLabel
644 . point - the point
645 - value - The point value
646 
647   Level: intermediate
648 
649 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
650 @*/
651 PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
652 {
653   PetscInt       v;
654   PetscErrorCode ierr;
655 
656   PetscFunctionBegin;
657   /* Find label value */
658   for (v = 0; v < label->numStrata; ++v) {
659     if (label->stratumValues[v] == value) break;
660   }
661   if (v >= label->numStrata) PetscFunctionReturn(0);
662   if (label->validIS[v]) {
663     ierr = DMLabelMakeInvalid_Private(label,v);CHKERRQ(ierr);
664   }
665   if (label->bt) {
666     if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, label->pStart, label->pEnd);
667     ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
668   }
669   ierr = PetscHashIDelKey(label->ht[v], point);CHKERRQ(ierr);
670   PetscFunctionReturn(0);
671 }
672 
673 #undef __FUNCT__
674 #define __FUNCT__ "DMLabelInsertIS"
675 /*@
676   DMLabelInsertIS - Set all points in the IS to a value
677 
678   Input Parameters:
679 + label - the DMLabel
680 . is    - the point IS
681 - value - The point value
682 
683   Level: intermediate
684 
685 .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
686 @*/
687 PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
688 {
689   const PetscInt *points;
690   PetscInt        n, p;
691   PetscErrorCode  ierr;
692 
693   PetscFunctionBegin;
694   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
695   ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr);
696   ierr = ISGetIndices(is, &points);CHKERRQ(ierr);
697   for (p = 0; p < n; ++p) {ierr = DMLabelSetValue(label, points[p], value);CHKERRQ(ierr);}
698   ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr);
699   PetscFunctionReturn(0);
700 }
701 
702 #undef __FUNCT__
703 #define __FUNCT__ "DMLabelGetNumValues"
704 PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
705 {
706   PetscFunctionBegin;
707   PetscValidPointer(numValues, 2);
708   *numValues = label->numStrata;
709   PetscFunctionReturn(0);
710 }
711 
712 #undef __FUNCT__
713 #define __FUNCT__ "DMLabelGetValueIS"
714 PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
715 {
716   PetscErrorCode ierr;
717 
718   PetscFunctionBegin;
719   PetscValidPointer(values, 2);
720   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr);
721   PetscFunctionReturn(0);
722 }
723 
724 #undef __FUNCT__
725 #define __FUNCT__ "DMLabelHasStratum"
726 PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
727 {
728   PetscInt v;
729 
730   PetscFunctionBegin;
731   PetscValidPointer(exists, 3);
732   *exists = PETSC_FALSE;
733   for (v = 0; v < label->numStrata; ++v) {
734     if (label->stratumValues[v] == value) {
735       *exists = PETSC_TRUE;
736       break;
737     }
738   }
739   PetscFunctionReturn(0);
740 }
741 
742 #undef __FUNCT__
743 #define __FUNCT__ "DMLabelGetStratumSize"
744 PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
745 {
746   PetscInt       v;
747   PetscErrorCode ierr;
748 
749   PetscFunctionBegin;
750   PetscValidPointer(size, 3);
751   *size = 0;
752   for (v = 0; v < label->numStrata; ++v) {
753     if (label->stratumValues[v] == value) {
754       ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
755       *size = label->stratumSizes[v];
756       break;
757     }
758   }
759   PetscFunctionReturn(0);
760 }
761 
762 #undef __FUNCT__
763 #define __FUNCT__ "DMLabelGetStratumBounds"
764 PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
765 {
766   PetscInt       v;
767   PetscErrorCode ierr;
768 
769   PetscFunctionBegin;
770   if (start) {PetscValidPointer(start, 3); *start = 0;}
771   if (end)   {PetscValidPointer(end,   4); *end   = 0;}
772   for (v = 0; v < label->numStrata; ++v) {
773     PetscInt min, max;
774 
775     if (label->stratumValues[v] != value) continue;
776     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
777     if (label->stratumSizes[v]  <= 0)     break;
778     ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr);
779     if (start) *start = min;
780     if (end)   *end   = max+1;
781     break;
782   }
783   PetscFunctionReturn(0);
784 }
785 
786 #undef __FUNCT__
787 #define __FUNCT__ "DMLabelGetStratumIS"
788 PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
789 {
790   PetscInt       v;
791   PetscErrorCode ierr;
792 
793   PetscFunctionBegin;
794   PetscValidPointer(points, 3);
795   *points = NULL;
796   for (v = 0; v < label->numStrata; ++v) {
797     if (label->stratumValues[v] == value) {
798       ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
799       if (label->validIS[v]) {
800         ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr);
801         *points = label->points[v];
802       } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to implement this to speedup Stratify");
803       break;
804     }
805   }
806   PetscFunctionReturn(0);
807 }
808 
809 #undef __FUNCT__
810 #define __FUNCT__ "DMLabelSetStratumIS"
811 PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
812 {
813   PetscInt       v, numStrata;
814   PetscErrorCode ierr;
815 
816   PetscFunctionBegin;
817   numStrata = label->numStrata;
818   for (v = 0; v < numStrata; v++) {
819     if (label->stratumValues[v] == value) break;
820   }
821   if (v >= numStrata) {ierr = DMLabelAddStratum(label,value);CHKERRQ(ierr);}
822   if (is == label->points[v]) PetscFunctionReturn(0);
823   ierr = DMLabelClearStratum(label,value);CHKERRQ(ierr);
824   ierr = ISGetLocalSize(is,&(label->stratumSizes[v]));CHKERRQ(ierr);
825   label->stratumValues[v] = value;
826   label->validIS[v] = PETSC_TRUE;
827   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
828   ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
829   if (label->bt) {
830     const PetscInt *points;
831     PetscInt p;
832 
833     ierr = ISGetIndices(is,&points);CHKERRQ(ierr);
834     for (p = 0; p < label->stratumSizes[v]; ++p) {
835       const PetscInt point = points[p];
836 
837       if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, label->pStart, label->pEnd);
838       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
839     }
840   }
841   label->points[v] = is;
842   PetscFunctionReturn(0);
843 }
844 
845 
846 #undef __FUNCT__
847 #define __FUNCT__ "DMLabelClearStratum"
848 PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
849 {
850   PetscInt       v;
851   PetscErrorCode ierr;
852 
853   PetscFunctionBegin;
854   for (v = 0; v < label->numStrata; ++v) {
855     if (label->stratumValues[v] == value) break;
856   }
857   if (v >= label->numStrata) PetscFunctionReturn(0);
858   if (label->validIS[v]) {
859     if (label->bt) {
860       PetscInt       i;
861       const PetscInt *points;
862 
863       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
864       for (i = 0; i < label->stratumSizes[v]; ++i) {
865         const PetscInt point = points[i];
866 
867         if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, label->pStart, label->pEnd);
868         ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
869       }
870       ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
871     }
872     ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
873     label->stratumSizes[v] = 0;
874     ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
875     ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
876   } else {
877     PetscHashIClear(label->ht[v]);
878   }
879   PetscFunctionReturn(0);
880 }
881 
882 #undef __FUNCT__
883 #define __FUNCT__ "DMLabelFilter"
884 PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
885 {
886   PetscInt       v;
887   PetscErrorCode ierr;
888 
889   PetscFunctionBegin;
890   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
891   label->pStart = start;
892   label->pEnd   = end;
893   if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);}
894   /* Could squish offsets, but would only make sense if I reallocate the storage */
895   for (v = 0; v < label->numStrata; ++v) {
896     PetscInt off, q;
897     const PetscInt *points;
898     PetscInt *pointsNew = NULL;
899 
900     ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
901     for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) {
902       const PetscInt point = points[q];
903 
904       if ((point < start) || (point >= end)) {
905         if (!pointsNew) {
906           ierr = PetscMalloc1(label->stratumSizes[v],&pointsNew);CHKERRQ(ierr);
907           ierr = PetscMemcpy(pointsNew,points,(size_t) off * sizeof(PetscInt));CHKERRQ(ierr);
908         }
909         continue;
910       }
911       if (pointsNew) {
912         pointsNew[off++] = point;
913       }
914     }
915     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
916     if (pointsNew) {
917       ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
918       ierr = ISCreateGeneral(PETSC_COMM_SELF,off,pointsNew,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
919       ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
920     }
921     label->stratumSizes[v] = off;
922   }
923   ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr);
924   PetscFunctionReturn(0);
925 }
926 
927 #undef __FUNCT__
928 #define __FUNCT__ "DMLabelPermute"
929 PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
930 {
931   const PetscInt *perm;
932   PetscInt        numValues, numPoints, v, q;
933   PetscErrorCode  ierr;
934 
935   PetscFunctionBegin;
936   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
937   ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr);
938   ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr);
939   ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr);
940   ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr);
941   for (v = 0; v < numValues; ++v) {
942     const PetscInt size   = (*labelNew)->stratumSizes[v];
943     const PetscInt *points;
944     PetscInt *pointsNew;
945 
946     ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
947     ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr);
948     for (q = 0; q < size; ++q) {
949       const PetscInt point = points[q];
950 
951       if ((point < 0) || (point >= numPoints)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [0, %D) for the remapping", point, numPoints);
952       pointsNew[q] = perm[point];
953     }
954     ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
955     ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr);
956     ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr);
957     ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr);
958     ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr);
959   }
960   ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr);
961   if (label->bt) {
962     ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
963     ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
964   }
965   PetscFunctionReturn(0);
966 }
967 
968 #undef __FUNCT__
969 #define __FUNCT__ "DMLabelDistribute_Internal"
970 PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
971 {
972   MPI_Comm       comm;
973   PetscInt       s, l, nroots, nleaves, dof, offset, size;
974   PetscInt      *remoteOffsets, *rootStrata, *rootIdx;
975   PetscSection   rootSection;
976   PetscSF        labelSF;
977   PetscErrorCode ierr;
978 
979   PetscFunctionBegin;
980   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
981   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
982   /* Build a section of stratum values per point, generate the according SF
983      and distribute point-wise stratum values to leaves. */
984   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
985   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
986   ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr);
987   if (label) {
988     for (s = 0; s < label->numStrata; ++s) {
989       const PetscInt *points;
990 
991       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
992       for (l = 0; l < label->stratumSizes[s]; l++) {
993         ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr);
994         ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr);
995       }
996       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
997     }
998   }
999   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
1000   /* Create a point-wise array of stratum values */
1001   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
1002   ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr);
1003   ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr);
1004   if (label) {
1005     for (s = 0; s < label->numStrata; ++s) {
1006       const PetscInt *points;
1007 
1008       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
1009       for (l = 0; l < label->stratumSizes[s]; l++) {
1010         const PetscInt p = points[l];
1011         ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr);
1012         rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
1013       }
1014       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
1015     }
1016   }
1017   /* Build SF that maps label points to remote processes */
1018   ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr);
1019   ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr);
1020   ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr);
1021   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
1022   /* Send the strata for each point over the derived SF */
1023   ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr);
1024   ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr);
1025   ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
1026   ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
1027   /* Clean up */
1028   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
1029   ierr = PetscFree(rootIdx);CHKERRQ(ierr);
1030   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
1031   ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr);
1032   PetscFunctionReturn(0);
1033 }
1034 
1035 #undef __FUNCT__
1036 #define __FUNCT__ "DMLabelDistribute"
1037 PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1038 {
1039   MPI_Comm       comm;
1040   PetscSection   leafSection;
1041   PetscInt       p, pStart, pEnd, s, size, dof, offset, stratum;
1042   PetscInt      *leafStrata, *strataIdx;
1043   PetscInt     **points;
1044   char          *name;
1045   PetscInt       nameSize;
1046   PetscHashI     stratumHash;
1047   PETSC_UNUSED   PetscHashIIter ret, iter;
1048   size_t         len = 0;
1049   PetscMPIInt    rank;
1050   PetscErrorCode ierr;
1051 
1052   PetscFunctionBegin;
1053   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
1054   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1055   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1056   /* Bcast name */
1057   if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);}
1058   nameSize = len;
1059   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1060   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1061   if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);}
1062   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1063   ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr);
1064   ierr = PetscFree(name);CHKERRQ(ierr);
1065   /* Bcast defaultValue */
1066   if (!rank) (*labelNew)->defaultValue = label->defaultValue;
1067   ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1068   /* Distribute stratum values over the SF and get the point mapping on the receiver */
1069   ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr);
1070   /* Determine received stratum values and initialise new label*/
1071   PetscHashICreate(stratumHash);
1072   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
1073   for (p = 0; p < size; ++p) PetscHashIPut(stratumHash, leafStrata[p], ret, iter);
1074   PetscHashISize(stratumHash, (*labelNew)->numStrata);
1075   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr);
1076   for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
1077   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr);
1078   /* Turn leafStrata into indices rather than stratum values */
1079   offset = 0;
1080   ierr = PetscHashIGetKeys(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr);
1081   for (p = 0; p < size; ++p) {
1082     for (s = 0; s < (*labelNew)->numStrata; ++s) {
1083       if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
1084     }
1085   }
1086   /* Rebuild the point strata on the receiver */
1087   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr);
1088   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1089   for (p=pStart; p<pEnd; p++) {
1090     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1091     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1092     for (s=0; s<dof; s++) {
1093       (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1094     }
1095   }
1096   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr);
1097   ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr);
1098   ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr);
1099   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1100     PetscHashICreate((*labelNew)->ht[s]);
1101     ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr);
1102   }
1103   /* Insert points into new strata */
1104   ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr);
1105   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1106   for (p=pStart; p<pEnd; p++) {
1107     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1108     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1109     for (s=0; s<dof; s++) {
1110       stratum = leafStrata[offset+s];
1111       points[stratum][strataIdx[stratum]++] = p;
1112     }
1113   }
1114   for (s = 0; s < (*labelNew)->numStrata; s++) {
1115     ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr);
1116     ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr);
1117   }
1118   ierr = PetscFree(points);CHKERRQ(ierr);
1119   PetscHashIDestroy(stratumHash);
1120   ierr = PetscFree(leafStrata);CHKERRQ(ierr);
1121   ierr = PetscFree(strataIdx);CHKERRQ(ierr);
1122   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
1123   PetscFunctionReturn(0);
1124 }
1125 
1126 #undef __FUNCT__
1127 #define __FUNCT__ "DMLabelGather"
1128 /*@
1129   DMLabelGather - Gather all label values from leafs into roots
1130 
1131   Input Parameters:
1132 + label - the DMLabel
1133 . point - the Star Forest point communication map
1134 
1135   Input Parameters:
1136 + label - the new DMLabel with localised leaf values
1137 
1138   Level: developer
1139 
1140   Note: This is the inverse operation to DMLabelDistribute.
1141 
1142 .seealso: DMLabelDistribute()
1143 @*/
1144 PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
1145 {
1146   MPI_Comm       comm;
1147   PetscSection   rootSection;
1148   PetscSF        sfLabel;
1149   PetscSFNode   *rootPoints, *leafPoints;
1150   PetscInt       p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
1151   const PetscInt *rootDegree, *ilocal;
1152   PetscInt       *rootStrata;
1153   char          *name;
1154   PetscInt       nameSize;
1155   size_t         len = 0;
1156   PetscMPIInt    rank, numProcs;
1157   PetscErrorCode ierr;
1158 
1159   PetscFunctionBegin;
1160   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1161   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1162   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
1163   /* Bcast name */
1164   if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);}
1165   nameSize = len;
1166   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1167   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1168   if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);}
1169   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1170   ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr);
1171   ierr = PetscFree(name);CHKERRQ(ierr);
1172   /* Gather rank/index pairs of leaves into local roots to build
1173      an inverse, multi-rooted SF. Note that this ignores local leaf
1174      indexing due to the use of the multiSF in PetscSFGather. */
1175   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr);
1176   ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr);
1177   for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
1178   for (p = 0; p < nleaves; p++) {
1179     leafPoints[ilocal[p]].index = ilocal[p];
1180     leafPoints[ilocal[p]].rank  = rank;
1181   }
1182   ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr);
1183   ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr);
1184   for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
1185   ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr);
1186   ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
1187   ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
1188   ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr);
1189   ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
1190   /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
1191   ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr);
1192   /* Rebuild the point strata on the receiver */
1193   for (p = 0, idx = 0; p < nroots; p++) {
1194     for (d = 0; d < rootDegree[p]; d++) {
1195       ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr);
1196       ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr);
1197       for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);}
1198     }
1199     idx += rootDegree[p];
1200   }
1201   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
1202   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
1203   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
1204   ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr);
1205   PetscFunctionReturn(0);
1206 }
1207 
1208 #undef __FUNCT__
1209 #define __FUNCT__ "DMLabelConvertToSection"
1210 PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1211 {
1212   IS              vIS;
1213   const PetscInt *values;
1214   PetscInt       *points;
1215   PetscInt        nV, vS = 0, vE = 0, v, N;
1216   PetscErrorCode  ierr;
1217 
1218   PetscFunctionBegin;
1219   ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr);
1220   ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
1221   ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr);
1222   if (nV) {vS = values[0]; vE = values[0]+1;}
1223   for (v = 1; v < nV; ++v) {
1224     vS = PetscMin(vS, values[v]);
1225     vE = PetscMax(vE, values[v]+1);
1226   }
1227   ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr);
1228   ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr);
1229   for (v = 0; v < nV; ++v) {
1230     PetscInt n;
1231 
1232     ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr);
1233     ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr);
1234   }
1235   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
1236   ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr);
1237   ierr = PetscMalloc1(N, &points);CHKERRQ(ierr);
1238   for (v = 0; v < nV; ++v) {
1239     IS              is;
1240     const PetscInt *spoints;
1241     PetscInt        dof, off, p;
1242 
1243     ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr);
1244     ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr);
1245     ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr);
1246     ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr);
1247     for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1248     ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr);
1249     ierr = ISDestroy(&is);CHKERRQ(ierr);
1250   }
1251   ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr);
1252   ierr = ISDestroy(&vIS);CHKERRQ(ierr);
1253   ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1254   PetscFunctionReturn(0);
1255 }
1256 
1257 #undef __FUNCT__
1258 #define __FUNCT__ "PetscSectionCreateGlobalSectionLabel"
1259 /*@C
1260   PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1261   the local section and an SF describing the section point overlap.
1262 
1263   Input Parameters:
1264   + s - The PetscSection for the local field layout
1265   . sf - The SF describing parallel layout of the section points
1266   . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1267   . label - The label specifying the points
1268   - labelValue - The label stratum specifying the points
1269 
1270   Output Parameter:
1271   . gsection - The PetscSection for the global field layout
1272 
1273   Note: This gives negative sizes and offsets to points not owned by this process
1274 
1275   Level: developer
1276 
1277 .seealso: PetscSectionCreate()
1278 @*/
1279 PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1280 {
1281   PetscInt      *neg = NULL, *tmpOff = NULL;
1282   PetscInt       pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1283   PetscErrorCode ierr;
1284 
1285   PetscFunctionBegin;
1286   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr);
1287   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1288   ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr);
1289   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1290   if (nroots >= 0) {
1291     if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1292     ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr);
1293     if (nroots > pEnd-pStart) {
1294       ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr);
1295     } else {
1296       tmpOff = &(*gsection)->atlasDof[-pStart];
1297     }
1298   }
1299   /* Mark ghost points with negative dof */
1300   for (p = pStart; p < pEnd; ++p) {
1301     PetscInt value;
1302 
1303     ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr);
1304     if (value != labelValue) continue;
1305     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1306     ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr);
1307     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1308     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);}
1309     if (neg) neg[p] = -(dof+1);
1310   }
1311   ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr);
1312   if (nroots >= 0) {
1313     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1314     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1315     if (nroots > pEnd-pStart) {
1316       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1317     }
1318   }
1319   /* Calculate new sizes, get proccess offset, and calculate point offsets */
1320   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1321     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1322     (*gsection)->atlasOff[p] = off;
1323     off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1324   }
1325   ierr       = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1326   globalOff -= off;
1327   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1328     (*gsection)->atlasOff[p] += globalOff;
1329     if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1330   }
1331   /* Put in negative offsets for ghost points */
1332   if (nroots >= 0) {
1333     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1334     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1335     if (nroots > pEnd-pStart) {
1336       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1337     }
1338   }
1339   if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);}
1340   ierr = PetscFree(neg);CHKERRQ(ierr);
1341   PetscFunctionReturn(0);
1342 }
1343 
1344