xref: /petsc/src/dm/label/dmlabel.c (revision 0e79e033d600ff7d649ae360e2604535d17718c2)
15fdea053SToby Isaac #include <petscdm.h>
2c58f1c22SToby Isaac #include <petsc/private/dmlabelimpl.h>   /*I      "petscdmlabel.h"   I*/
3c58f1c22SToby Isaac #include <petsc/private/isimpl.h>        /*I      "petscis.h"        I*/
4c58f1c22SToby Isaac #include <petscsf.h>
5c58f1c22SToby Isaac 
6c58f1c22SToby Isaac /*@C
7c58f1c22SToby Isaac   DMLabelCreate - Create a DMLabel object, which is a multimap
8c58f1c22SToby Isaac 
9d67d17b1SMatthew G. Knepley   Input parameters:
10d67d17b1SMatthew G. Knepley + comm - The communicator, usually PETSC_COMM_SELF
11d67d17b1SMatthew G. Knepley - name - The label name
12c58f1c22SToby Isaac 
13c58f1c22SToby Isaac   Output parameter:
14c58f1c22SToby Isaac . label - The DMLabel
15c58f1c22SToby Isaac 
16c58f1c22SToby Isaac   Level: beginner
17c58f1c22SToby Isaac 
18c58f1c22SToby Isaac .seealso: DMLabelDestroy()
19c58f1c22SToby Isaac @*/
20d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelCreate(MPI_Comm comm, const char name[], DMLabel *label)
21c58f1c22SToby Isaac {
22c58f1c22SToby Isaac   PetscErrorCode ierr;
23c58f1c22SToby Isaac 
24c58f1c22SToby Isaac   PetscFunctionBegin;
25d67d17b1SMatthew G. Knepley   PetscValidPointer(label,2);
26d67d17b1SMatthew G. Knepley   ierr = DMInitializePackage();CHKERRQ(ierr);
27c58f1c22SToby Isaac 
28d67d17b1SMatthew G. Knepley   ierr = PetscHeaderCreate(*label,DMLABEL_CLASSID,"DMLabel","DMLabel","DM",comm,DMLabelDestroy,DMLabelView);CHKERRQ(ierr);
29d67d17b1SMatthew G. Knepley 
30c58f1c22SToby Isaac   (*label)->numStrata      = 0;
315aa44df4SToby Isaac   (*label)->defaultValue   = -1;
32c58f1c22SToby Isaac   (*label)->stratumValues  = NULL;
33ad8374ffSToby Isaac   (*label)->validIS        = NULL;
34c58f1c22SToby Isaac   (*label)->stratumSizes   = NULL;
35c58f1c22SToby Isaac   (*label)->points         = NULL;
36c58f1c22SToby Isaac   (*label)->ht             = NULL;
37c58f1c22SToby Isaac   (*label)->pStart         = -1;
38c58f1c22SToby Isaac   (*label)->pEnd           = -1;
39c58f1c22SToby Isaac   (*label)->bt             = NULL;
40d67d17b1SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *label, name);CHKERRQ(ierr);
41c58f1c22SToby Isaac   PetscFunctionReturn(0);
42c58f1c22SToby Isaac }
43c58f1c22SToby Isaac 
44c58f1c22SToby Isaac /*
45c58f1c22SToby Isaac   DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format
46c58f1c22SToby Isaac 
47c58f1c22SToby Isaac   Input parameter:
48c58f1c22SToby Isaac + label - The DMLabel
49c58f1c22SToby Isaac - v - The stratum value
50c58f1c22SToby Isaac 
51c58f1c22SToby Isaac   Output parameter:
52c58f1c22SToby Isaac . label - The DMLabel with stratum in sorted list format
53c58f1c22SToby Isaac 
54c58f1c22SToby Isaac   Level: developer
55c58f1c22SToby Isaac 
56c58f1c22SToby Isaac .seealso: DMLabelCreate()
57c58f1c22SToby Isaac */
58c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v)
59c58f1c22SToby Isaac {
600c3c4a36SLisandro Dalcin   PetscInt       off = 0;
61ad8374ffSToby Isaac   PetscInt       *pointArray;
62c58f1c22SToby Isaac   PetscErrorCode ierr;
63c58f1c22SToby Isaac 
64ad8374ffSToby Isaac   if (label->validIS[v]) return 0;
65c58f1c22SToby Isaac   PetscFunctionBegin;
660c3c4a36SLisandro Dalcin   if (v < 0 || v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeValid_Private\n", v);
67e8f14785SLisandro Dalcin   ierr = PetscHSetIGetSize(label->ht[v], &label->stratumSizes[v]);CHKERRQ(ierr);
68c58f1c22SToby Isaac 
69ad8374ffSToby Isaac   ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr);
70e8f14785SLisandro Dalcin   ierr = PetscHSetIGetElems(label->ht[v], &off, pointArray);CHKERRQ(ierr);
71c58f1c22SToby Isaac   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]);
72e8f14785SLisandro Dalcin   PetscHSetIClear(label->ht[v]);
73ad8374ffSToby Isaac   ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr);
74c58f1c22SToby Isaac   if (label->bt) {
75c58f1c22SToby Isaac     PetscInt p;
76c58f1c22SToby Isaac 
77c58f1c22SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
78ad8374ffSToby Isaac       const PetscInt point = pointArray[p];
79c58f1c22SToby Isaac 
80c58f1c22SToby Isaac       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);
81c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
82c58f1c22SToby Isaac     }
83c58f1c22SToby Isaac   }
84ad8374ffSToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF,label->stratumSizes[v],pointArray,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
85ad8374ffSToby Isaac   ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
86ad8374ffSToby Isaac   label->validIS[v] = PETSC_TRUE;
87d67d17b1SMatthew G. Knepley   ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
88c58f1c22SToby Isaac   PetscFunctionReturn(0);
89c58f1c22SToby Isaac }
90c58f1c22SToby Isaac 
91c58f1c22SToby Isaac /*
92c58f1c22SToby Isaac   DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format
93c58f1c22SToby Isaac 
94c58f1c22SToby Isaac   Input parameter:
95c58f1c22SToby Isaac . label - The DMLabel
96c58f1c22SToby Isaac 
97c58f1c22SToby Isaac   Output parameter:
98c58f1c22SToby Isaac . label - The DMLabel with all strata in sorted list format
99c58f1c22SToby Isaac 
100c58f1c22SToby Isaac   Level: developer
101c58f1c22SToby Isaac 
102c58f1c22SToby Isaac .seealso: DMLabelCreate()
103c58f1c22SToby Isaac */
104c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label)
105c58f1c22SToby Isaac {
106c58f1c22SToby Isaac   PetscInt       v;
107c58f1c22SToby Isaac   PetscErrorCode ierr;
108c58f1c22SToby Isaac 
109c58f1c22SToby Isaac   PetscFunctionBegin;
110c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; v++){
111c58f1c22SToby Isaac     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
112c58f1c22SToby Isaac   }
113c58f1c22SToby Isaac   PetscFunctionReturn(0);
114c58f1c22SToby Isaac }
115c58f1c22SToby Isaac 
116c58f1c22SToby Isaac /*
117c58f1c22SToby Isaac   DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format
118c58f1c22SToby Isaac 
119c58f1c22SToby Isaac   Input parameter:
120c58f1c22SToby Isaac + label - The DMLabel
121c58f1c22SToby Isaac - v - The stratum value
122c58f1c22SToby Isaac 
123c58f1c22SToby Isaac   Output parameter:
124c58f1c22SToby Isaac . label - The DMLabel with stratum in hash format
125c58f1c22SToby Isaac 
126c58f1c22SToby Isaac   Level: developer
127c58f1c22SToby Isaac 
128c58f1c22SToby Isaac .seealso: DMLabelCreate()
129c58f1c22SToby Isaac */
130c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v)
131c58f1c22SToby Isaac {
132c58f1c22SToby Isaac   PetscInt       p;
133ad8374ffSToby Isaac   const PetscInt *points;
134c58f1c22SToby Isaac   PetscErrorCode ierr;
135c58f1c22SToby Isaac 
1360c3c4a36SLisandro Dalcin   if (!label->validIS[v]) return 0;
137c58f1c22SToby Isaac   PetscFunctionBegin;
1380c3c4a36SLisandro Dalcin   if (v < 0 || v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeInvalid_Private\n", v);
139ad8374ffSToby Isaac   if (label->points[v]) {
140ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
141e8f14785SLisandro Dalcin     for (p = 0; p < label->stratumSizes[v]; ++p) {
142e8f14785SLisandro Dalcin       ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);
143e8f14785SLisandro Dalcin     }
144ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
145ad8374ffSToby Isaac     ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
146ad8374ffSToby Isaac   }
147ad8374ffSToby Isaac   label->validIS[v] = PETSC_FALSE;
148c58f1c22SToby Isaac   PetscFunctionReturn(0);
149c58f1c22SToby Isaac }
150c58f1c22SToby Isaac 
1510c3c4a36SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index)
1520c3c4a36SLisandro Dalcin {
1530c3c4a36SLisandro Dalcin   PetscInt v;
154*0e79e033SBarry Smith 
1550c3c4a36SLisandro Dalcin   PetscFunctionBegin;
156*0e79e033SBarry Smith   *index = -1;
157d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
158*0e79e033SBarry Smith   for (v = 0; v < label->numStrata; ++v) {
159*0e79e033SBarry Smith     if (label->stratumValues[v] == value) {
160*0e79e033SBarry Smith       *index = v; break;
161*0e79e033SBarry Smith     }
162*0e79e033SBarry Smith   }
1630c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
1640c3c4a36SLisandro Dalcin }
1650c3c4a36SLisandro Dalcin 
1660c3c4a36SLisandro Dalcin static PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index)
167c58f1c22SToby Isaac {
168ad8374ffSToby Isaac   PetscInt       v, *tmpV, *tmpS;
169ad8374ffSToby Isaac   IS             *tmpP;
170e8f14785SLisandro Dalcin   PetscHSetI     *tmpH;
171c58f1c22SToby Isaac   PetscBool      *tmpB;
172c58f1c22SToby Isaac   PetscErrorCode ierr;
173c58f1c22SToby Isaac 
174c58f1c22SToby Isaac   PetscFunctionBegin;
175d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
176c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpV);CHKERRQ(ierr);
177c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpS);CHKERRQ(ierr);
178c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpH);CHKERRQ(ierr);
179c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpP);CHKERRQ(ierr);
180c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpB);CHKERRQ(ierr);
181c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
182c58f1c22SToby Isaac     tmpV[v] = label->stratumValues[v];
183c58f1c22SToby Isaac     tmpS[v] = label->stratumSizes[v];
184c58f1c22SToby Isaac     tmpH[v] = label->ht[v];
185c58f1c22SToby Isaac     tmpP[v] = label->points[v];
186ad8374ffSToby Isaac     tmpB[v] = label->validIS[v];
187c58f1c22SToby Isaac   }
188c58f1c22SToby Isaac   tmpV[v] = value;
189c58f1c22SToby Isaac   tmpS[v] = 0;
190e8f14785SLisandro Dalcin   ierr = PetscHSetICreate(&tmpH[v]);CHKERRQ(ierr);
191ecb9a371SToby Isaac   ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&tmpP[v]);CHKERRQ(ierr);
192c58f1c22SToby Isaac   tmpB[v] = PETSC_TRUE;
1930c3c4a36SLisandro Dalcin 
194c58f1c22SToby Isaac   ++label->numStrata;
195c58f1c22SToby Isaac   ierr = PetscFree(label->stratumValues);CHKERRQ(ierr);
196c58f1c22SToby Isaac   ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr);
197c58f1c22SToby Isaac   ierr = PetscFree(label->ht);CHKERRQ(ierr);
198c58f1c22SToby Isaac   ierr = PetscFree(label->points);CHKERRQ(ierr);
199ad8374ffSToby Isaac   ierr = PetscFree(label->validIS);CHKERRQ(ierr);
200c58f1c22SToby Isaac   label->stratumValues = tmpV;
201c58f1c22SToby Isaac   label->stratumSizes  = tmpS;
202c58f1c22SToby Isaac   label->ht            = tmpH;
203c58f1c22SToby Isaac   label->points        = tmpP;
204ad8374ffSToby Isaac   label->validIS       = tmpB;
205c58f1c22SToby Isaac 
2060c3c4a36SLisandro Dalcin   *index = v;
2070c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
2080c3c4a36SLisandro Dalcin }
2090c3c4a36SLisandro Dalcin 
2100c3c4a36SLisandro Dalcin PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
2110c3c4a36SLisandro Dalcin {
2120c3c4a36SLisandro Dalcin   PetscInt       v;
2130c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
2140c3c4a36SLisandro Dalcin 
2150c3c4a36SLisandro Dalcin   PetscFunctionBegin;
216d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
2170c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
2180c3c4a36SLisandro Dalcin   if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);}
219c58f1c22SToby Isaac   PetscFunctionReturn(0);
220c58f1c22SToby Isaac }
221c58f1c22SToby Isaac 
222c58f1c22SToby Isaac static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
223c58f1c22SToby Isaac {
224c58f1c22SToby Isaac   PetscInt       v;
225c58f1c22SToby Isaac   PetscMPIInt    rank;
226c58f1c22SToby Isaac   PetscErrorCode ierr;
227c58f1c22SToby Isaac 
228c58f1c22SToby Isaac   PetscFunctionBegin;
229c58f1c22SToby Isaac   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr);
230c58f1c22SToby Isaac   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
231c58f1c22SToby Isaac   if (label) {
232d67d17b1SMatthew G. Knepley     const char *name;
233d67d17b1SMatthew G. Knepley 
234d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
235d67d17b1SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", name);CHKERRQ(ierr);
236c58f1c22SToby Isaac     if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);}
237c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
238c58f1c22SToby Isaac       const PetscInt value = label->stratumValues[v];
239ad8374ffSToby Isaac       const PetscInt *points;
240c58f1c22SToby Isaac       PetscInt       p;
241c58f1c22SToby Isaac 
242ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
243c58f1c22SToby Isaac       for (p = 0; p < label->stratumSizes[v]; ++p) {
244ad8374ffSToby Isaac         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr);
245c58f1c22SToby Isaac       }
246ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
247c58f1c22SToby Isaac     }
248c58f1c22SToby Isaac   }
249c58f1c22SToby Isaac   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
250c58f1c22SToby Isaac   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
251c58f1c22SToby Isaac   PetscFunctionReturn(0);
252c58f1c22SToby Isaac }
253c58f1c22SToby Isaac 
254c58f1c22SToby Isaac /*@C
255c58f1c22SToby Isaac   DMLabelView - View the label
256c58f1c22SToby Isaac 
257c58f1c22SToby Isaac   Input Parameters:
258c58f1c22SToby Isaac + label - The DMLabel
259c58f1c22SToby Isaac - viewer - The PetscViewer
260c58f1c22SToby Isaac 
261c58f1c22SToby Isaac   Level: intermediate
262c58f1c22SToby Isaac 
263c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy()
264c58f1c22SToby Isaac @*/
265c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
266c58f1c22SToby Isaac {
267c58f1c22SToby Isaac   PetscBool      iascii;
268c58f1c22SToby Isaac   PetscErrorCode ierr;
269c58f1c22SToby Isaac 
270c58f1c22SToby Isaac   PetscFunctionBegin;
271d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
272c58f1c22SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
273c58f1c22SToby Isaac   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
274c58f1c22SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
275c58f1c22SToby Isaac   if (iascii) {
276c58f1c22SToby Isaac     ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr);
277c58f1c22SToby Isaac   }
278c58f1c22SToby Isaac   PetscFunctionReturn(0);
279c58f1c22SToby Isaac }
280c58f1c22SToby Isaac 
28184f0b6dfSMatthew G. Knepley /*@
282d67d17b1SMatthew G. Knepley   DMLabelReset - Destroys internal data structures in a DMLabel
283d67d17b1SMatthew G. Knepley 
284d67d17b1SMatthew G. Knepley   Input Parameter:
285d67d17b1SMatthew G. Knepley . label - The DMLabel
286d67d17b1SMatthew G. Knepley 
287d67d17b1SMatthew G. Knepley   Level: beginner
288d67d17b1SMatthew G. Knepley 
289d67d17b1SMatthew G. Knepley .seealso: DMLabelDestroy(), DMLabelCreate()
290d67d17b1SMatthew G. Knepley @*/
291d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelReset(DMLabel label)
292d67d17b1SMatthew G. Knepley {
293d67d17b1SMatthew G. Knepley   PetscInt       v;
294d67d17b1SMatthew G. Knepley   PetscErrorCode ierr;
295d67d17b1SMatthew G. Knepley 
296d67d17b1SMatthew G. Knepley   PetscFunctionBegin;
297d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
298d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->stratumValues);CHKERRQ(ierr);
299d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr);
300d67d17b1SMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
301d67d17b1SMatthew G. Knepley     ierr = PetscHSetIDestroy(&label->ht[v]);CHKERRQ(ierr);
302d67d17b1SMatthew G. Knepley     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
303d67d17b1SMatthew G. Knepley   }
304d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->ht);CHKERRQ(ierr);
305d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->points);CHKERRQ(ierr);
306d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->validIS);CHKERRQ(ierr);
307d67d17b1SMatthew G. Knepley   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
308d67d17b1SMatthew G. Knepley   PetscFunctionReturn(0);
309d67d17b1SMatthew G. Knepley }
310d67d17b1SMatthew G. Knepley 
311d67d17b1SMatthew G. Knepley /*@
31284f0b6dfSMatthew G. Knepley   DMLabelDestroy - Destroys a DMLabel
31384f0b6dfSMatthew G. Knepley 
31484f0b6dfSMatthew G. Knepley   Input Parameter:
31584f0b6dfSMatthew G. Knepley . label - The DMLabel
31684f0b6dfSMatthew G. Knepley 
31784f0b6dfSMatthew G. Knepley   Level: beginner
31884f0b6dfSMatthew G. Knepley 
319d67d17b1SMatthew G. Knepley .seealso: DMLabelReset(), DMLabelCreate()
32084f0b6dfSMatthew G. Knepley @*/
321c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label)
322c58f1c22SToby Isaac {
323c58f1c22SToby Isaac   PetscErrorCode ierr;
324c58f1c22SToby Isaac 
325c58f1c22SToby Isaac   PetscFunctionBegin;
326d67d17b1SMatthew G. Knepley   if (!*label) PetscFunctionReturn(0);
327d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific((*label),DMLABEL_CLASSID,1);
328d67d17b1SMatthew G. Knepley   if (--((PetscObject)(*label))->refct > 0) {
329d67d17b1SMatthew G. Knepley     *label = NULL;
330d67d17b1SMatthew G. Knepley     PetscFunctionReturn(0);
3310c3c4a36SLisandro Dalcin   }
332d67d17b1SMatthew G. Knepley   ierr = DMLabelReset(*label);CHKERRQ(ierr);
333d67d17b1SMatthew G. Knepley   ierr = PetscHeaderDestroy(label);CHKERRQ(ierr);
334c58f1c22SToby Isaac   PetscFunctionReturn(0);
335c58f1c22SToby Isaac }
336c58f1c22SToby Isaac 
33784f0b6dfSMatthew G. Knepley /*@
33884f0b6dfSMatthew G. Knepley   DMLabelDuplicate - Duplicates a DMLabel
33984f0b6dfSMatthew G. Knepley 
34084f0b6dfSMatthew G. Knepley   Input Parameter:
34184f0b6dfSMatthew G. Knepley . label - The DMLabel
34284f0b6dfSMatthew G. Knepley 
34384f0b6dfSMatthew G. Knepley   Output Parameter:
34484f0b6dfSMatthew G. Knepley . labelnew - location to put new vector
34584f0b6dfSMatthew G. Knepley 
34684f0b6dfSMatthew G. Knepley   Level: intermediate
34784f0b6dfSMatthew G. Knepley 
34884f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelDestroy()
34984f0b6dfSMatthew G. Knepley @*/
350c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
351c58f1c22SToby Isaac {
352d67d17b1SMatthew G. Knepley   const char    *name;
353ad8374ffSToby Isaac   PetscInt       v;
354c58f1c22SToby Isaac   PetscErrorCode ierr;
355c58f1c22SToby Isaac 
356c58f1c22SToby Isaac   PetscFunctionBegin;
357d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
358c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
359d67d17b1SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
360d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PetscObjectComm((PetscObject) label), name, labelnew);CHKERRQ(ierr);
361c58f1c22SToby Isaac 
362c58f1c22SToby Isaac   (*labelnew)->numStrata    = label->numStrata;
3635aa44df4SToby Isaac   (*labelnew)->defaultValue = label->defaultValue;
364c58f1c22SToby Isaac   if (label->numStrata) {
365c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr);
366c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr);
367c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr);
368c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr);
369ad8374ffSToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr);
370c58f1c22SToby Isaac     /* Could eliminate unused space here */
371c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
372e8f14785SLisandro Dalcin       ierr = PetscHSetICreate(&(*labelnew)->ht[v]);CHKERRQ(ierr);
373ad8374ffSToby Isaac       (*labelnew)->validIS[v]        = PETSC_TRUE;
374c58f1c22SToby Isaac       (*labelnew)->stratumValues[v]  = label->stratumValues[v];
375c58f1c22SToby Isaac       (*labelnew)->stratumSizes[v]   = label->stratumSizes[v];
376ad8374ffSToby Isaac       ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr);
377ad8374ffSToby Isaac       (*labelnew)->points[v]         = label->points[v];
378c58f1c22SToby Isaac     }
379c58f1c22SToby Isaac   }
380c58f1c22SToby Isaac   (*labelnew)->pStart = -1;
381c58f1c22SToby Isaac   (*labelnew)->pEnd   = -1;
382c58f1c22SToby Isaac   (*labelnew)->bt     = NULL;
383c58f1c22SToby Isaac   PetscFunctionReturn(0);
384c58f1c22SToby Isaac }
385c58f1c22SToby Isaac 
386c58f1c22SToby Isaac /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
387c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
388c58f1c22SToby Isaac {
389c58f1c22SToby Isaac   PetscInt       v;
390c58f1c22SToby Isaac   PetscErrorCode ierr;
391c58f1c22SToby Isaac 
392c58f1c22SToby Isaac   PetscFunctionBegin;
393d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
3940c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
395c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
396c58f1c22SToby Isaac   label->pStart = pStart;
397c58f1c22SToby Isaac   label->pEnd   = pEnd;
398c58f1c22SToby Isaac   ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr);
399c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
400ad8374ffSToby Isaac     const PetscInt *points;
401c58f1c22SToby Isaac     PetscInt       i;
402c58f1c22SToby Isaac 
403ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
404c58f1c22SToby Isaac     for (i = 0; i < label->stratumSizes[v]; ++i) {
405ad8374ffSToby Isaac       const PetscInt point = points[i];
406c58f1c22SToby Isaac 
407c58f1c22SToby Isaac       if ((point < pStart) || (point >= pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, pStart, pEnd);
408c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr);
409c58f1c22SToby Isaac     }
410ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
411c58f1c22SToby Isaac   }
412c58f1c22SToby Isaac   PetscFunctionReturn(0);
413c58f1c22SToby Isaac }
414c58f1c22SToby Isaac 
415c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label)
416c58f1c22SToby Isaac {
417c58f1c22SToby Isaac   PetscErrorCode ierr;
418c58f1c22SToby Isaac 
419c58f1c22SToby Isaac   PetscFunctionBegin;
420d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
421c58f1c22SToby Isaac   label->pStart = -1;
422c58f1c22SToby Isaac   label->pEnd   = -1;
4230c3c4a36SLisandro Dalcin   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
424c58f1c22SToby Isaac   PetscFunctionReturn(0);
425c58f1c22SToby Isaac }
426c58f1c22SToby Isaac 
427c58f1c22SToby Isaac /*@
428c58f1c22SToby Isaac   DMLabelHasValue - Determine whether a label assigns the value to any point
429c58f1c22SToby Isaac 
430c58f1c22SToby Isaac   Input Parameters:
431c58f1c22SToby Isaac + label - the DMLabel
432c58f1c22SToby Isaac - value - the value
433c58f1c22SToby Isaac 
434c58f1c22SToby Isaac   Output Parameter:
435c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point
436c58f1c22SToby Isaac 
437c58f1c22SToby Isaac   Level: developer
438c58f1c22SToby Isaac 
439c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
440c58f1c22SToby Isaac @*/
441c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
442c58f1c22SToby Isaac {
443c58f1c22SToby Isaac   PetscInt v;
4440c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
445c58f1c22SToby Isaac 
446c58f1c22SToby Isaac   PetscFunctionBegin;
447d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
448c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
4490c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
4500c3c4a36SLisandro Dalcin   *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE;
451c58f1c22SToby Isaac   PetscFunctionReturn(0);
452c58f1c22SToby Isaac }
453c58f1c22SToby Isaac 
454c58f1c22SToby Isaac /*@
455c58f1c22SToby Isaac   DMLabelHasPoint - Determine whether a label assigns a value to a point
456c58f1c22SToby Isaac 
457c58f1c22SToby Isaac   Input Parameters:
458c58f1c22SToby Isaac + label - the DMLabel
459c58f1c22SToby Isaac - point - the point
460c58f1c22SToby Isaac 
461c58f1c22SToby Isaac   Output Parameter:
462c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value
463c58f1c22SToby Isaac 
464c58f1c22SToby Isaac   Note: The user must call DMLabelCreateIndex() before this function.
465c58f1c22SToby Isaac 
466c58f1c22SToby Isaac   Level: developer
467c58f1c22SToby Isaac 
468c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
469c58f1c22SToby Isaac @*/
470c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
471c58f1c22SToby Isaac {
472c58f1c22SToby Isaac   PetscErrorCode ierr;
473c58f1c22SToby Isaac 
474c58f1c22SToby Isaac   PetscFunctionBeginHot;
475d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
476c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
477c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
478c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG)
479c58f1c22SToby Isaac   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
480c58f1c22SToby Isaac   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);
481c58f1c22SToby Isaac #endif
482c58f1c22SToby Isaac   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
483c58f1c22SToby Isaac   PetscFunctionReturn(0);
484c58f1c22SToby Isaac }
485c58f1c22SToby Isaac 
486c58f1c22SToby Isaac /*@
487c58f1c22SToby Isaac   DMLabelStratumHasPoint - Return true if the stratum contains a point
488c58f1c22SToby Isaac 
489c58f1c22SToby Isaac   Input Parameters:
490c58f1c22SToby Isaac + label - the DMLabel
491c58f1c22SToby Isaac . value - the stratum value
492c58f1c22SToby Isaac - point - the point
493c58f1c22SToby Isaac 
494c58f1c22SToby Isaac   Output Parameter:
495c58f1c22SToby Isaac . contains - true if the stratum contains the point
496c58f1c22SToby Isaac 
497c58f1c22SToby Isaac   Level: intermediate
498c58f1c22SToby Isaac 
499c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
500c58f1c22SToby Isaac @*/
501c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
502c58f1c22SToby Isaac {
503c58f1c22SToby Isaac   PetscInt       v;
504c58f1c22SToby Isaac   PetscErrorCode ierr;
505c58f1c22SToby Isaac 
5060c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
507d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
508c58f1c22SToby Isaac   PetscValidPointer(contains, 4);
509c58f1c22SToby Isaac   *contains = PETSC_FALSE;
5100c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
5110c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
5120c3c4a36SLisandro Dalcin 
513ad8374ffSToby Isaac   if (label->validIS[v]) {
514c58f1c22SToby Isaac     PetscInt i;
515c58f1c22SToby Isaac 
516a2d74346SToby Isaac     ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
5170c3c4a36SLisandro Dalcin     if (i >= 0) *contains = PETSC_TRUE;
518c58f1c22SToby Isaac   } else {
519c58f1c22SToby Isaac     PetscBool has;
520c58f1c22SToby Isaac 
521e8f14785SLisandro Dalcin     PetscHSetIHas(label->ht[v], point, &has);
5220c3c4a36SLisandro Dalcin     if (has) *contains = PETSC_TRUE;
523c58f1c22SToby Isaac   }
524c58f1c22SToby Isaac   PetscFunctionReturn(0);
525c58f1c22SToby Isaac }
526c58f1c22SToby Isaac 
52784f0b6dfSMatthew G. Knepley /*@
5285aa44df4SToby Isaac   DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
5295aa44df4SToby Isaac   When a label is created, it is initialized to -1.
5305aa44df4SToby Isaac 
5315aa44df4SToby Isaac   Input parameter:
5325aa44df4SToby Isaac . label - a DMLabel object
5335aa44df4SToby Isaac 
5345aa44df4SToby Isaac   Output parameter:
5355aa44df4SToby Isaac . defaultValue - the default value
5365aa44df4SToby Isaac 
5375aa44df4SToby Isaac   Level: beginner
5385aa44df4SToby Isaac 
5395aa44df4SToby Isaac .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
54084f0b6dfSMatthew G. Knepley @*/
5415aa44df4SToby Isaac PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
5425aa44df4SToby Isaac {
5435aa44df4SToby Isaac   PetscFunctionBegin;
544d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
5455aa44df4SToby Isaac   *defaultValue = label->defaultValue;
5465aa44df4SToby Isaac   PetscFunctionReturn(0);
5475aa44df4SToby Isaac }
5485aa44df4SToby Isaac 
54984f0b6dfSMatthew G. Knepley /*@
5505aa44df4SToby Isaac   DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
5515aa44df4SToby Isaac   When a label is created, it is initialized to -1.
5525aa44df4SToby Isaac 
5535aa44df4SToby Isaac   Input parameter:
5545aa44df4SToby Isaac . label - a DMLabel object
5555aa44df4SToby Isaac 
5565aa44df4SToby Isaac   Output parameter:
5575aa44df4SToby Isaac . defaultValue - the default value
5585aa44df4SToby Isaac 
5595aa44df4SToby Isaac   Level: beginner
5605aa44df4SToby Isaac 
5615aa44df4SToby Isaac .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
56284f0b6dfSMatthew G. Knepley @*/
5635aa44df4SToby Isaac PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
5645aa44df4SToby Isaac {
5655aa44df4SToby Isaac   PetscFunctionBegin;
566d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
5675aa44df4SToby Isaac   label->defaultValue = defaultValue;
5685aa44df4SToby Isaac   PetscFunctionReturn(0);
5695aa44df4SToby Isaac }
5705aa44df4SToby Isaac 
571c58f1c22SToby Isaac /*@
5725aa44df4SToby Isaac   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())
573c58f1c22SToby Isaac 
574c58f1c22SToby Isaac   Input Parameters:
575c58f1c22SToby Isaac + label - the DMLabel
576c58f1c22SToby Isaac - point - the point
577c58f1c22SToby Isaac 
578c58f1c22SToby Isaac   Output Parameter:
5798e68afcfSMatthew G. Knepley . value - The point value, or the default value (-1 by default)
580c58f1c22SToby Isaac 
581c58f1c22SToby Isaac   Level: intermediate
582c58f1c22SToby Isaac 
5835aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
584c58f1c22SToby Isaac @*/
585c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
586c58f1c22SToby Isaac {
587c58f1c22SToby Isaac   PetscInt       v;
588c58f1c22SToby Isaac   PetscErrorCode ierr;
589c58f1c22SToby Isaac 
5900c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
591d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
592c58f1c22SToby Isaac   PetscValidPointer(value, 3);
5935aa44df4SToby Isaac   *value = label->defaultValue;
594c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
595ad8374ffSToby Isaac     if (label->validIS[v]) {
596c58f1c22SToby Isaac       PetscInt i;
597c58f1c22SToby Isaac 
598a2d74346SToby Isaac       ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
599c58f1c22SToby Isaac       if (i >= 0) {
600c58f1c22SToby Isaac         *value = label->stratumValues[v];
601c58f1c22SToby Isaac         break;
602c58f1c22SToby Isaac       }
603c58f1c22SToby Isaac     } else {
604c58f1c22SToby Isaac       PetscBool has;
605c58f1c22SToby Isaac 
606e8f14785SLisandro Dalcin       PetscHSetIHas(label->ht[v], point, &has);
607c58f1c22SToby Isaac       if (has) {
608c58f1c22SToby Isaac         *value = label->stratumValues[v];
609c58f1c22SToby Isaac         break;
610c58f1c22SToby Isaac       }
611c58f1c22SToby Isaac     }
612c58f1c22SToby Isaac   }
613c58f1c22SToby Isaac   PetscFunctionReturn(0);
614c58f1c22SToby Isaac }
615c58f1c22SToby Isaac 
616c58f1c22SToby Isaac /*@
617367003a6SStefano Zampini   DMLabelSetValue - Set the value a label assigns to a point.  If the value is the same as the label's default value (which is initially -1, and can be changed with DMLabelSetDefaultValue() to something different), then this function will do nothing.
618c58f1c22SToby Isaac 
619c58f1c22SToby Isaac   Input Parameters:
620c58f1c22SToby Isaac + label - the DMLabel
621c58f1c22SToby Isaac . point - the point
622c58f1c22SToby Isaac - value - The point value
623c58f1c22SToby Isaac 
624c58f1c22SToby Isaac   Level: intermediate
625c58f1c22SToby Isaac 
6265aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
627c58f1c22SToby Isaac @*/
628c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
629c58f1c22SToby Isaac {
630c58f1c22SToby Isaac   PetscInt       v;
631c58f1c22SToby Isaac   PetscErrorCode ierr;
632c58f1c22SToby Isaac 
633c58f1c22SToby Isaac   PetscFunctionBegin;
634d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
6350c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
6365aa44df4SToby Isaac   if (value == label->defaultValue) PetscFunctionReturn(0);
6370c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
6380c3c4a36SLisandro Dalcin   if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);}
639c58f1c22SToby Isaac   /* Set key */
6400c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
641e8f14785SLisandro Dalcin   ierr = PetscHSetIAdd(label->ht[v], point);CHKERRQ(ierr);
642c58f1c22SToby Isaac   PetscFunctionReturn(0);
643c58f1c22SToby Isaac }
644c58f1c22SToby Isaac 
645c58f1c22SToby Isaac /*@
646c58f1c22SToby Isaac   DMLabelClearValue - Clear the value a label assigns to a point
647c58f1c22SToby Isaac 
648c58f1c22SToby Isaac   Input Parameters:
649c58f1c22SToby Isaac + label - the DMLabel
650c58f1c22SToby Isaac . point - the point
651c58f1c22SToby Isaac - value - The point value
652c58f1c22SToby Isaac 
653c58f1c22SToby Isaac   Level: intermediate
654c58f1c22SToby Isaac 
655c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
656c58f1c22SToby Isaac @*/
657c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
658c58f1c22SToby Isaac {
659ad8374ffSToby Isaac   PetscInt       v;
660c58f1c22SToby Isaac   PetscErrorCode ierr;
661c58f1c22SToby Isaac 
662c58f1c22SToby Isaac   PetscFunctionBegin;
663d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
664c58f1c22SToby Isaac   /* Find label value */
6650c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
6660c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
6670c3c4a36SLisandro Dalcin 
668eeed21e7SToby Isaac   if (label->bt) {
669eeed21e7SToby Isaac     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);
670eeed21e7SToby Isaac     ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
671eeed21e7SToby Isaac   }
6720c3c4a36SLisandro Dalcin 
6730c3c4a36SLisandro Dalcin   /* Delete key */
6740c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
675e8f14785SLisandro Dalcin   ierr = PetscHSetIDel(label->ht[v], point);CHKERRQ(ierr);
676c58f1c22SToby Isaac   PetscFunctionReturn(0);
677c58f1c22SToby Isaac }
678c58f1c22SToby Isaac 
679c58f1c22SToby Isaac /*@
680c58f1c22SToby Isaac   DMLabelInsertIS - Set all points in the IS to a value
681c58f1c22SToby Isaac 
682c58f1c22SToby Isaac   Input Parameters:
683c58f1c22SToby Isaac + label - the DMLabel
684c58f1c22SToby Isaac . is    - the point IS
685c58f1c22SToby Isaac - value - The point value
686c58f1c22SToby Isaac 
687c58f1c22SToby Isaac   Level: intermediate
688c58f1c22SToby Isaac 
689c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
690c58f1c22SToby Isaac @*/
691c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
692c58f1c22SToby Isaac {
6930c3c4a36SLisandro Dalcin   PetscInt        v, n, p;
694c58f1c22SToby Isaac   const PetscInt *points;
695c58f1c22SToby Isaac   PetscErrorCode  ierr;
696c58f1c22SToby Isaac 
697c58f1c22SToby Isaac   PetscFunctionBegin;
698d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
699c58f1c22SToby Isaac   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
7000c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
7010c3c4a36SLisandro Dalcin   if (value == label->defaultValue) PetscFunctionReturn(0);
7020c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
7030c3c4a36SLisandro Dalcin   if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);}
7040c3c4a36SLisandro Dalcin   /* Set keys */
7050c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
706c58f1c22SToby Isaac   ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr);
707c58f1c22SToby Isaac   ierr = ISGetIndices(is, &points);CHKERRQ(ierr);
708e8f14785SLisandro Dalcin   for (p = 0; p < n; ++p) {ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);}
709c58f1c22SToby Isaac   ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr);
710c58f1c22SToby Isaac   PetscFunctionReturn(0);
711c58f1c22SToby Isaac }
712c58f1c22SToby Isaac 
71384f0b6dfSMatthew G. Knepley /*@
71484f0b6dfSMatthew G. Knepley   DMLabelGetNumValues - Get the number of values that the DMLabel takes
71584f0b6dfSMatthew G. Knepley 
71684f0b6dfSMatthew G. Knepley   Input Parameter:
71784f0b6dfSMatthew G. Knepley . label - the DMLabel
71884f0b6dfSMatthew G. Knepley 
71984f0b6dfSMatthew G. Knepley   Output Paramater:
72084f0b6dfSMatthew G. Knepley . numValues - the number of values
72184f0b6dfSMatthew G. Knepley 
72284f0b6dfSMatthew G. Knepley   Level: intermediate
72384f0b6dfSMatthew G. Knepley 
72484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
72584f0b6dfSMatthew G. Knepley @*/
726c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
727c58f1c22SToby Isaac {
728c58f1c22SToby Isaac   PetscFunctionBegin;
729d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
730c58f1c22SToby Isaac   PetscValidPointer(numValues, 2);
731c58f1c22SToby Isaac   *numValues = label->numStrata;
732c58f1c22SToby Isaac   PetscFunctionReturn(0);
733c58f1c22SToby Isaac }
734c58f1c22SToby Isaac 
73584f0b6dfSMatthew G. Knepley /*@
73684f0b6dfSMatthew G. Knepley   DMLabelGetValueIS - Get an IS of all values that the DMlabel takes
73784f0b6dfSMatthew G. Knepley 
73884f0b6dfSMatthew G. Knepley   Input Parameter:
73984f0b6dfSMatthew G. Knepley . label - the DMLabel
74084f0b6dfSMatthew G. Knepley 
74184f0b6dfSMatthew G. Knepley   Output Paramater:
74284f0b6dfSMatthew G. Knepley . is    - the value IS
74384f0b6dfSMatthew G. Knepley 
74484f0b6dfSMatthew G. Knepley   Level: intermediate
74584f0b6dfSMatthew G. Knepley 
74684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
74784f0b6dfSMatthew G. Knepley @*/
748c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
749c58f1c22SToby Isaac {
750c58f1c22SToby Isaac   PetscErrorCode ierr;
751c58f1c22SToby Isaac 
752c58f1c22SToby Isaac   PetscFunctionBegin;
753d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
754c58f1c22SToby Isaac   PetscValidPointer(values, 2);
755c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr);
756c58f1c22SToby Isaac   PetscFunctionReturn(0);
757c58f1c22SToby Isaac }
758c58f1c22SToby Isaac 
75984f0b6dfSMatthew G. Knepley /*@
76084f0b6dfSMatthew G. Knepley   DMLabelHasStratum - Determine whether points exist with the given value
76184f0b6dfSMatthew G. Knepley 
76284f0b6dfSMatthew G. Knepley   Input Parameters:
76384f0b6dfSMatthew G. Knepley + label - the DMLabel
76484f0b6dfSMatthew G. Knepley - value - the stratum value
76584f0b6dfSMatthew G. Knepley 
76684f0b6dfSMatthew G. Knepley   Output Paramater:
76784f0b6dfSMatthew G. Knepley . exists - Flag saying whether points exist
76884f0b6dfSMatthew G. Knepley 
76984f0b6dfSMatthew G. Knepley   Level: intermediate
77084f0b6dfSMatthew G. Knepley 
77184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
77284f0b6dfSMatthew G. Knepley @*/
773fada774cSMatthew G. Knepley PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
774fada774cSMatthew G. Knepley {
775fada774cSMatthew G. Knepley   PetscInt       v;
7760c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
777fada774cSMatthew G. Knepley 
778fada774cSMatthew G. Knepley   PetscFunctionBegin;
779d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
780fada774cSMatthew G. Knepley   PetscValidPointer(exists, 3);
7810c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
7820c3c4a36SLisandro Dalcin   *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE;
783fada774cSMatthew G. Knepley   PetscFunctionReturn(0);
784fada774cSMatthew G. Knepley }
785fada774cSMatthew G. Knepley 
78684f0b6dfSMatthew G. Knepley /*@
78784f0b6dfSMatthew G. Knepley   DMLabelGetStratumSize - Get the size of a stratum
78884f0b6dfSMatthew G. Knepley 
78984f0b6dfSMatthew G. Knepley   Input Parameters:
79084f0b6dfSMatthew G. Knepley + label - the DMLabel
79184f0b6dfSMatthew G. Knepley - value - the stratum value
79284f0b6dfSMatthew G. Knepley 
79384f0b6dfSMatthew G. Knepley   Output Paramater:
79484f0b6dfSMatthew G. Knepley . size - The number of points in the stratum
79584f0b6dfSMatthew G. Knepley 
79684f0b6dfSMatthew G. Knepley   Level: intermediate
79784f0b6dfSMatthew G. Knepley 
79884f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
79984f0b6dfSMatthew G. Knepley @*/
800c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
801c58f1c22SToby Isaac {
802c58f1c22SToby Isaac   PetscInt       v;
803c58f1c22SToby Isaac   PetscErrorCode ierr;
804c58f1c22SToby Isaac 
805c58f1c22SToby Isaac   PetscFunctionBegin;
806d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
807c58f1c22SToby Isaac   PetscValidPointer(size, 3);
808c58f1c22SToby Isaac   *size = 0;
8090c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
8100c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
811c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
812c58f1c22SToby Isaac   *size = label->stratumSizes[v];
813c58f1c22SToby Isaac   PetscFunctionReturn(0);
814c58f1c22SToby Isaac }
815c58f1c22SToby Isaac 
81684f0b6dfSMatthew G. Knepley /*@
81784f0b6dfSMatthew G. Knepley   DMLabelGetStratumBounds - Get the largest and smallest point of a stratum
81884f0b6dfSMatthew G. Knepley 
81984f0b6dfSMatthew G. Knepley   Input Parameters:
82084f0b6dfSMatthew G. Knepley + label - the DMLabel
82184f0b6dfSMatthew G. Knepley - value - the stratum value
82284f0b6dfSMatthew G. Knepley 
82384f0b6dfSMatthew G. Knepley   Output Paramaters:
82484f0b6dfSMatthew G. Knepley + start - the smallest point in the stratum
82584f0b6dfSMatthew G. Knepley - end - the largest point in the stratum
82684f0b6dfSMatthew G. Knepley 
82784f0b6dfSMatthew G. Knepley   Level: intermediate
82884f0b6dfSMatthew G. Knepley 
82984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
83084f0b6dfSMatthew G. Knepley @*/
831c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
832c58f1c22SToby Isaac {
8330c3c4a36SLisandro Dalcin   PetscInt       v, min, max;
834c58f1c22SToby Isaac   PetscErrorCode ierr;
835c58f1c22SToby Isaac 
836c58f1c22SToby Isaac   PetscFunctionBegin;
837d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
838c58f1c22SToby Isaac   if (start) {PetscValidPointer(start, 3); *start = 0;}
839c58f1c22SToby Isaac   if (end)   {PetscValidPointer(end,   4); *end   = 0;}
8400c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
8410c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
842c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
8430c3c4a36SLisandro Dalcin   if (label->stratumSizes[v] <= 0) PetscFunctionReturn(0);
844d6cb179aSToby Isaac   ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr);
845d6cb179aSToby Isaac   if (start) *start = min;
846d6cb179aSToby Isaac   if (end)   *end   = max+1;
847c58f1c22SToby Isaac   PetscFunctionReturn(0);
848c58f1c22SToby Isaac }
849c58f1c22SToby Isaac 
85084f0b6dfSMatthew G. Knepley /*@
85184f0b6dfSMatthew G. Knepley   DMLabelGetStratumIS - Get an IS with the stratum points
85284f0b6dfSMatthew G. Knepley 
85384f0b6dfSMatthew G. Knepley   Input Parameters:
85484f0b6dfSMatthew G. Knepley + label - the DMLabel
85584f0b6dfSMatthew G. Knepley - value - the stratum value
85684f0b6dfSMatthew G. Knepley 
85784f0b6dfSMatthew G. Knepley   Output Paramater:
85884f0b6dfSMatthew G. Knepley . points - The stratum points
85984f0b6dfSMatthew G. Knepley 
86084f0b6dfSMatthew G. Knepley   Level: intermediate
86184f0b6dfSMatthew G. Knepley 
86284f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
86384f0b6dfSMatthew G. Knepley @*/
864c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
865c58f1c22SToby Isaac {
866c58f1c22SToby Isaac   PetscInt       v;
867c58f1c22SToby Isaac   PetscErrorCode ierr;
868c58f1c22SToby Isaac 
869c58f1c22SToby Isaac   PetscFunctionBegin;
870d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
871c58f1c22SToby Isaac   PetscValidPointer(points, 3);
872c58f1c22SToby Isaac   *points = NULL;
8730c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
8740c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
875c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
876ad8374ffSToby Isaac   ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr);
877ad8374ffSToby Isaac   *points = label->points[v];
878c58f1c22SToby Isaac   PetscFunctionReturn(0);
879c58f1c22SToby Isaac }
880c58f1c22SToby Isaac 
88184f0b6dfSMatthew G. Knepley /*@
88284f0b6dfSMatthew G. Knepley   DMLabelSetStratumIS - Set the stratum points using an IS
88384f0b6dfSMatthew G. Knepley 
88484f0b6dfSMatthew G. Knepley   Input Parameters:
88584f0b6dfSMatthew G. Knepley + label - the DMLabel
88684f0b6dfSMatthew G. Knepley . value - the stratum value
88784f0b6dfSMatthew G. Knepley - points - The stratum points
88884f0b6dfSMatthew G. Knepley 
88984f0b6dfSMatthew G. Knepley   Level: intermediate
89084f0b6dfSMatthew G. Knepley 
89184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
89284f0b6dfSMatthew G. Knepley @*/
8934de306b1SToby Isaac PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
8944de306b1SToby Isaac {
8950c3c4a36SLisandro Dalcin   PetscInt       v;
8964de306b1SToby Isaac   PetscErrorCode ierr;
8974de306b1SToby Isaac 
8984de306b1SToby Isaac   PetscFunctionBegin;
899d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
900d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
9010c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
9020c3c4a36SLisandro Dalcin   if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);}
9034de306b1SToby Isaac   if (is == label->points[v]) PetscFunctionReturn(0);
9044de306b1SToby Isaac   ierr = DMLabelClearStratum(label, value);CHKERRQ(ierr);
9054de306b1SToby Isaac   ierr = ISGetLocalSize(is, &(label->stratumSizes[v]));CHKERRQ(ierr);
9064de306b1SToby Isaac   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
9074de306b1SToby Isaac   ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
9080c3c4a36SLisandro Dalcin   label->points[v] = is;
9090c3c4a36SLisandro Dalcin   label->validIS[v] = PETSC_TRUE;
9104de306b1SToby Isaac   if (label->bt) {
9114de306b1SToby Isaac     const PetscInt *points;
9124de306b1SToby Isaac     PetscInt p;
9134de306b1SToby Isaac 
9144de306b1SToby Isaac     ierr = ISGetIndices(is,&points);CHKERRQ(ierr);
9154de306b1SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
9164de306b1SToby Isaac       const PetscInt point = points[p];
9174de306b1SToby Isaac 
9184de306b1SToby Isaac       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);
9194de306b1SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
9204de306b1SToby Isaac     }
9214de306b1SToby Isaac   }
9224de306b1SToby Isaac   PetscFunctionReturn(0);
9234de306b1SToby Isaac }
9244de306b1SToby Isaac 
92584f0b6dfSMatthew G. Knepley /*@
92684f0b6dfSMatthew G. Knepley   DMLabelClearStratum - Remove a stratum
9274de306b1SToby Isaac 
92884f0b6dfSMatthew G. Knepley   Input Parameters:
92984f0b6dfSMatthew G. Knepley + label - the DMLabel
93084f0b6dfSMatthew G. Knepley - value - the stratum value
93184f0b6dfSMatthew G. Knepley 
93284f0b6dfSMatthew G. Knepley   Level: intermediate
93384f0b6dfSMatthew G. Knepley 
93484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
93584f0b6dfSMatthew G. Knepley @*/
936c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
937c58f1c22SToby Isaac {
938c58f1c22SToby Isaac   PetscInt       v;
939c58f1c22SToby Isaac   PetscErrorCode ierr;
940c58f1c22SToby Isaac 
941c58f1c22SToby Isaac   PetscFunctionBegin;
942d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
9430c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
9440c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
9454de306b1SToby Isaac   if (label->validIS[v]) {
9464de306b1SToby Isaac     if (label->bt) {
947c58f1c22SToby Isaac       PetscInt       i;
948ad8374ffSToby Isaac       const PetscInt *points;
949c58f1c22SToby Isaac 
950ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
951c58f1c22SToby Isaac       for (i = 0; i < label->stratumSizes[v]; ++i) {
952ad8374ffSToby Isaac         const PetscInt point = points[i];
953c58f1c22SToby Isaac 
954c58f1c22SToby Isaac         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);
955c58f1c22SToby Isaac         ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
956c58f1c22SToby Isaac       }
957ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
958c58f1c22SToby Isaac     }
959c58f1c22SToby Isaac     label->stratumSizes[v] = 0;
9600c3c4a36SLisandro Dalcin     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
9610c3c4a36SLisandro Dalcin     ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr);
9620c3c4a36SLisandro Dalcin     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
963c58f1c22SToby Isaac   } else {
964e8f14785SLisandro Dalcin     PetscHSetIClear(label->ht[v]);
965c58f1c22SToby Isaac   }
966c58f1c22SToby Isaac   PetscFunctionReturn(0);
967c58f1c22SToby Isaac }
968c58f1c22SToby Isaac 
96984f0b6dfSMatthew G. Knepley /*@
97084f0b6dfSMatthew G. Knepley   DMLabelFilter - Remove all points outside of [start, end)
97184f0b6dfSMatthew G. Knepley 
97284f0b6dfSMatthew G. Knepley   Input Parameters:
97384f0b6dfSMatthew G. Knepley + label - the DMLabel
97484f0b6dfSMatthew G. Knepley . start - the first point
97584f0b6dfSMatthew G. Knepley - end - the last point
97684f0b6dfSMatthew G. Knepley 
97784f0b6dfSMatthew G. Knepley   Level: intermediate
97884f0b6dfSMatthew G. Knepley 
97984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
98084f0b6dfSMatthew G. Knepley @*/
981c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
982c58f1c22SToby Isaac {
983c58f1c22SToby Isaac   PetscInt       v;
984c58f1c22SToby Isaac   PetscErrorCode ierr;
985c58f1c22SToby Isaac 
986c58f1c22SToby Isaac   PetscFunctionBegin;
987d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
9880c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
989c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
990c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
991c58f1c22SToby Isaac     PetscInt off, q;
992ad8374ffSToby Isaac     const PetscInt *points;
993033405d5SLisandro Dalcin     PetscInt numPointsNew = 0, *pointsNew = NULL;
994c58f1c22SToby Isaac 
995ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
996033405d5SLisandro Dalcin     for (q = 0; q < label->stratumSizes[v]; ++q)
997033405d5SLisandro Dalcin       if (points[q] >= start && points[q] < end)
998033405d5SLisandro Dalcin         numPointsNew++;
999033405d5SLisandro Dalcin     ierr = PetscMalloc1(numPointsNew, &pointsNew);CHKERRQ(ierr);
1000c58f1c22SToby Isaac     for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) {
1001033405d5SLisandro Dalcin       if (points[q] >= start && points[q] < end)
1002033405d5SLisandro Dalcin         pointsNew[off++] = points[q];
1003ad8374ffSToby Isaac     }
1004ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
1005033405d5SLisandro Dalcin 
1006033405d5SLisandro Dalcin     label->stratumSizes[v] = numPointsNew;
1007033405d5SLisandro Dalcin     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
1008033405d5SLisandro Dalcin     ierr = ISCreateGeneral(PETSC_COMM_SELF,numPointsNew, pointsNew, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr);
1009033405d5SLisandro Dalcin     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
1010c58f1c22SToby Isaac   }
1011c58f1c22SToby Isaac   ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr);
1012c58f1c22SToby Isaac   PetscFunctionReturn(0);
1013c58f1c22SToby Isaac }
1014c58f1c22SToby Isaac 
101584f0b6dfSMatthew G. Knepley /*@
101684f0b6dfSMatthew G. Knepley   DMLabelPermute - Create a new label with permuted points
101784f0b6dfSMatthew G. Knepley 
101884f0b6dfSMatthew G. Knepley   Input Parameters:
101984f0b6dfSMatthew G. Knepley + label - the DMLabel
102084f0b6dfSMatthew G. Knepley - permutation - the point permutation
102184f0b6dfSMatthew G. Knepley 
102284f0b6dfSMatthew G. Knepley   Output Parameter:
102384f0b6dfSMatthew G. Knepley . labelnew - the new label containing the permuted points
102484f0b6dfSMatthew G. Knepley 
102584f0b6dfSMatthew G. Knepley   Level: intermediate
102684f0b6dfSMatthew G. Knepley 
102784f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
102884f0b6dfSMatthew G. Knepley @*/
1029c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
1030c58f1c22SToby Isaac {
1031c58f1c22SToby Isaac   const PetscInt *perm;
1032c58f1c22SToby Isaac   PetscInt        numValues, numPoints, v, q;
1033c58f1c22SToby Isaac   PetscErrorCode  ierr;
1034c58f1c22SToby Isaac 
1035c58f1c22SToby Isaac   PetscFunctionBegin;
1036d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1037d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(permutation, IS_CLASSID, 2);
1038c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1039c58f1c22SToby Isaac   ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr);
1040c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr);
1041c58f1c22SToby Isaac   ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr);
1042c58f1c22SToby Isaac   ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr);
1043c58f1c22SToby Isaac   for (v = 0; v < numValues; ++v) {
1044c58f1c22SToby Isaac     const PetscInt size   = (*labelNew)->stratumSizes[v];
1045ad8374ffSToby Isaac     const PetscInt *points;
1046ad8374ffSToby Isaac     PetscInt *pointsNew;
1047c58f1c22SToby Isaac 
1048ad8374ffSToby Isaac     ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1049ad8374ffSToby Isaac     ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr);
1050c58f1c22SToby Isaac     for (q = 0; q < size; ++q) {
1051ad8374ffSToby Isaac       const PetscInt point = points[q];
1052c58f1c22SToby Isaac 
1053c58f1c22SToby Isaac       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);
1054ad8374ffSToby Isaac       pointsNew[q] = perm[point];
1055c58f1c22SToby Isaac     }
1056ad8374ffSToby Isaac     ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1057ad8374ffSToby Isaac     ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr);
1058ad8374ffSToby Isaac     ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr);
1059fa8e8ae5SToby Isaac     if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) {
1060fa8e8ae5SToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));CHKERRQ(ierr);
1061fa8e8ae5SToby Isaac       ierr = PetscFree(pointsNew);CHKERRQ(ierr);
1062fa8e8ae5SToby Isaac     } else {
1063ad8374ffSToby Isaac       ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr);
1064fa8e8ae5SToby Isaac     }
1065ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr);
1066c58f1c22SToby Isaac   }
1067c58f1c22SToby Isaac   ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr);
1068c58f1c22SToby Isaac   if (label->bt) {
1069c58f1c22SToby Isaac     ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
1070c58f1c22SToby Isaac     ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
1071c58f1c22SToby Isaac   }
1072c58f1c22SToby Isaac   PetscFunctionReturn(0);
1073c58f1c22SToby Isaac }
1074c58f1c22SToby Isaac 
107526c55118SMichael Lange PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
107626c55118SMichael Lange {
107726c55118SMichael Lange   MPI_Comm       comm;
107826c55118SMichael Lange   PetscInt       s, l, nroots, nleaves, dof, offset, size;
107926c55118SMichael Lange   PetscInt      *remoteOffsets, *rootStrata, *rootIdx;
108026c55118SMichael Lange   PetscSection   rootSection;
108126c55118SMichael Lange   PetscSF        labelSF;
108226c55118SMichael Lange   PetscErrorCode ierr;
108326c55118SMichael Lange 
108426c55118SMichael Lange   PetscFunctionBegin;
108526c55118SMichael Lange   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
108626c55118SMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
108726c55118SMichael Lange   /* Build a section of stratum values per point, generate the according SF
108826c55118SMichael Lange      and distribute point-wise stratum values to leaves. */
108926c55118SMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
109026c55118SMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
109126c55118SMichael Lange   ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr);
109226c55118SMichael Lange   if (label) {
109326c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1094ad8374ffSToby Isaac       const PetscInt *points;
1095ad8374ffSToby Isaac 
1096ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
109726c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1098ad8374ffSToby Isaac         ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr);
1099ad8374ffSToby Isaac         ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr);
110026c55118SMichael Lange       }
1101ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
110226c55118SMichael Lange     }
110326c55118SMichael Lange   }
110426c55118SMichael Lange   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
110526c55118SMichael Lange   /* Create a point-wise array of stratum values */
110626c55118SMichael Lange   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
110726c55118SMichael Lange   ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr);
110826c55118SMichael Lange   ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr);
110926c55118SMichael Lange   if (label) {
111026c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1111ad8374ffSToby Isaac       const PetscInt *points;
1112ad8374ffSToby Isaac 
1113ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
111426c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1115ad8374ffSToby Isaac         const PetscInt p = points[l];
111626c55118SMichael Lange         ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr);
111726c55118SMichael Lange         rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
111826c55118SMichael Lange       }
1119ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
112026c55118SMichael Lange     }
112126c55118SMichael Lange   }
112226c55118SMichael Lange   /* Build SF that maps label points to remote processes */
112326c55118SMichael Lange   ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr);
112426c55118SMichael Lange   ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr);
112526c55118SMichael Lange   ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr);
112626c55118SMichael Lange   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
112726c55118SMichael Lange   /* Send the strata for each point over the derived SF */
112826c55118SMichael Lange   ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr);
112926c55118SMichael Lange   ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr);
113026c55118SMichael Lange   ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
113126c55118SMichael Lange   ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
113226c55118SMichael Lange   /* Clean up */
113326c55118SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
113426c55118SMichael Lange   ierr = PetscFree(rootIdx);CHKERRQ(ierr);
113526c55118SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
113626c55118SMichael Lange   ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr);
113726c55118SMichael Lange   PetscFunctionReturn(0);
113826c55118SMichael Lange }
113926c55118SMichael Lange 
114084f0b6dfSMatthew G. Knepley /*@
114184f0b6dfSMatthew G. Knepley   DMLabelDistribute - Create a new label pushed forward over the PetscSF
114284f0b6dfSMatthew G. Knepley 
114384f0b6dfSMatthew G. Knepley   Input Parameters:
114484f0b6dfSMatthew G. Knepley + label - the DMLabel
114584f0b6dfSMatthew G. Knepley - sf    - the map from old to new distribution
114684f0b6dfSMatthew G. Knepley 
114784f0b6dfSMatthew G. Knepley   Output Parameter:
114884f0b6dfSMatthew G. Knepley . labelnew - the new redistributed label
114984f0b6dfSMatthew G. Knepley 
115084f0b6dfSMatthew G. Knepley   Level: intermediate
115184f0b6dfSMatthew G. Knepley 
115284f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
115384f0b6dfSMatthew G. Knepley @*/
1154c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1155c58f1c22SToby Isaac {
1156c58f1c22SToby Isaac   MPI_Comm       comm;
115726c55118SMichael Lange   PetscSection   leafSection;
115826c55118SMichael Lange   PetscInt       p, pStart, pEnd, s, size, dof, offset, stratum;
115926c55118SMichael Lange   PetscInt      *leafStrata, *strataIdx;
1160ad8374ffSToby Isaac   PetscInt     **points;
1161d67d17b1SMatthew G. Knepley   const char    *lname = NULL;
1162c58f1c22SToby Isaac   char          *name;
1163c58f1c22SToby Isaac   PetscInt       nameSize;
1164e8f14785SLisandro Dalcin   PetscHSetI     stratumHash;
1165c58f1c22SToby Isaac   size_t         len = 0;
116626c55118SMichael Lange   PetscMPIInt    rank;
1167c58f1c22SToby Isaac   PetscErrorCode ierr;
1168c58f1c22SToby Isaac 
1169c58f1c22SToby Isaac   PetscFunctionBegin;
1170d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1171f018e600SMatthew G. Knepley   if (label) {
1172f018e600SMatthew G. Knepley     PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1173f018e600SMatthew G. Knepley     ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1174f018e600SMatthew G. Knepley   }
1175c58f1c22SToby Isaac   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1176c58f1c22SToby Isaac   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1177c58f1c22SToby Isaac   /* Bcast name */
1178d67d17b1SMatthew G. Knepley   if (!rank) {
1179d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1180d67d17b1SMatthew G. Knepley     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1181d67d17b1SMatthew G. Knepley   }
1182c58f1c22SToby Isaac   nameSize = len;
1183c58f1c22SToby Isaac   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1184c58f1c22SToby Isaac   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1185d67d17b1SMatthew G. Knepley   if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);}
1186c58f1c22SToby Isaac   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1187d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
1188c58f1c22SToby Isaac   ierr = PetscFree(name);CHKERRQ(ierr);
118977d236dfSMichael Lange   /* Bcast defaultValue */
119077d236dfSMichael Lange   if (!rank) (*labelNew)->defaultValue = label->defaultValue;
119177d236dfSMichael Lange   ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
119226c55118SMichael Lange   /* Distribute stratum values over the SF and get the point mapping on the receiver */
119326c55118SMichael Lange   ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr);
11945cbdf6fcSMichael Lange   /* Determine received stratum values and initialise new label*/
1195e8f14785SLisandro Dalcin   ierr = PetscHSetICreate(&stratumHash);CHKERRQ(ierr);
119626c55118SMichael Lange   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
1197e8f14785SLisandro Dalcin   for (p = 0; p < size; ++p) {ierr = PetscHSetIAdd(stratumHash, leafStrata[p]);CHKERRQ(ierr);}
1198e8f14785SLisandro Dalcin   ierr = PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata);CHKERRQ(ierr);
1199ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr);
1200ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
12015cbdf6fcSMichael Lange   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr);
12025cbdf6fcSMichael Lange   /* Turn leafStrata into indices rather than stratum values */
12035cbdf6fcSMichael Lange   offset = 0;
1204e8f14785SLisandro Dalcin   ierr = PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr);
1205a205f1fdSToby Isaac   ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr);
12065cbdf6fcSMichael Lange   for (p = 0; p < size; ++p) {
1207231b9e6fSMatthew G. Knepley     for (s = 0; s < (*labelNew)->numStrata; ++s) {
1208231b9e6fSMatthew G. Knepley       if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
12095cbdf6fcSMichael Lange     }
12105cbdf6fcSMichael Lange   }
1211c58f1c22SToby Isaac   /* Rebuild the point strata on the receiver */
1212c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr);
1213c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1214c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1215c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1216c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1217c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1218c58f1c22SToby Isaac       (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1219c58f1c22SToby Isaac     }
1220c58f1c22SToby Isaac   }
1221c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr);
1222c58f1c22SToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr);
1223ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr);
1224c58f1c22SToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1225e8f14785SLisandro Dalcin     ierr = PetscHSetICreate(&(*labelNew)->ht[s]);CHKERRQ(ierr);
1226ad8374ffSToby Isaac     ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr);
1227c58f1c22SToby Isaac   }
1228c58f1c22SToby Isaac   /* Insert points into new strata */
1229c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr);
1230c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1231c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1232c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1233c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1234c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1235c58f1c22SToby Isaac       stratum = leafStrata[offset+s];
1236ad8374ffSToby Isaac       points[stratum][strataIdx[stratum]++] = p;
1237c58f1c22SToby Isaac     }
1238c58f1c22SToby Isaac   }
1239ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; s++) {
1240ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr);
1241ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr);
1242ad8374ffSToby Isaac   }
1243ad8374ffSToby Isaac   ierr = PetscFree(points);CHKERRQ(ierr);
1244e8f14785SLisandro Dalcin   ierr = PetscHSetIDestroy(&stratumHash);CHKERRQ(ierr);
1245c58f1c22SToby Isaac   ierr = PetscFree(leafStrata);CHKERRQ(ierr);
1246c58f1c22SToby Isaac   ierr = PetscFree(strataIdx);CHKERRQ(ierr);
1247c58f1c22SToby Isaac   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
1248c58f1c22SToby Isaac   PetscFunctionReturn(0);
1249c58f1c22SToby Isaac }
1250c58f1c22SToby Isaac 
12517937d9ceSMichael Lange /*@
12527937d9ceSMichael Lange   DMLabelGather - Gather all label values from leafs into roots
12537937d9ceSMichael Lange 
12547937d9ceSMichael Lange   Input Parameters:
12557937d9ceSMichael Lange + label - the DMLabel
125684f0b6dfSMatthew G. Knepley - sf - the Star Forest point communication map
12577937d9ceSMichael Lange 
125884f0b6dfSMatthew G. Knepley   Output Parameters:
125984f0b6dfSMatthew G. Knepley . labelNew - the new DMLabel with localised leaf values
12607937d9ceSMichael Lange 
12617937d9ceSMichael Lange   Level: developer
12627937d9ceSMichael Lange 
12637937d9ceSMichael Lange   Note: This is the inverse operation to DMLabelDistribute.
12647937d9ceSMichael Lange 
12657937d9ceSMichael Lange .seealso: DMLabelDistribute()
12667937d9ceSMichael Lange @*/
12677937d9ceSMichael Lange PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
12687937d9ceSMichael Lange {
12697937d9ceSMichael Lange   MPI_Comm       comm;
12707937d9ceSMichael Lange   PetscSection   rootSection;
12717937d9ceSMichael Lange   PetscSF        sfLabel;
12727937d9ceSMichael Lange   PetscSFNode   *rootPoints, *leafPoints;
12737937d9ceSMichael Lange   PetscInt       p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
12747937d9ceSMichael Lange   const PetscInt *rootDegree, *ilocal;
12757937d9ceSMichael Lange   PetscInt       *rootStrata;
1276d67d17b1SMatthew G. Knepley   const char    *lname;
12777937d9ceSMichael Lange   char          *name;
12787937d9ceSMichael Lange   PetscInt       nameSize;
12797937d9ceSMichael Lange   size_t         len = 0;
12809852e123SBarry Smith   PetscMPIInt    rank, size;
12817937d9ceSMichael Lange   PetscErrorCode ierr;
12827937d9ceSMichael Lange 
12837937d9ceSMichael Lange   PetscFunctionBegin;
1284d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1285d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
12867937d9ceSMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
12877937d9ceSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
12889852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
12897937d9ceSMichael Lange   /* Bcast name */
1290d67d17b1SMatthew G. Knepley   if (!rank) {
1291d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1292d67d17b1SMatthew G. Knepley     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1293d67d17b1SMatthew G. Knepley   }
12947937d9ceSMichael Lange   nameSize = len;
12957937d9ceSMichael Lange   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
12967937d9ceSMichael Lange   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1297d67d17b1SMatthew G. Knepley   if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);}
12987937d9ceSMichael Lange   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1299d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
13007937d9ceSMichael Lange   ierr = PetscFree(name);CHKERRQ(ierr);
13017937d9ceSMichael Lange   /* Gather rank/index pairs of leaves into local roots to build
13027937d9ceSMichael Lange      an inverse, multi-rooted SF. Note that this ignores local leaf
13037937d9ceSMichael Lange      indexing due to the use of the multiSF in PetscSFGather. */
13047937d9ceSMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr);
1305dc53bc9bSMatthew G. Knepley   ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr);
1306dc53bc9bSMatthew G. Knepley   for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
13077937d9ceSMichael Lange   for (p = 0; p < nleaves; p++) {
13088212dd46SStefano Zampini     PetscInt ilp = ilocal ? ilocal[p] : p;
13098212dd46SStefano Zampini 
13108212dd46SStefano Zampini     leafPoints[ilp].index = ilp;
13118212dd46SStefano Zampini     leafPoints[ilp].rank  = rank;
13127937d9ceSMichael Lange   }
13137937d9ceSMichael Lange   ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr);
13147937d9ceSMichael Lange   ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr);
13157937d9ceSMichael Lange   for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
13167937d9ceSMichael Lange   ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr);
13177937d9ceSMichael Lange   ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
13187937d9ceSMichael Lange   ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
13197937d9ceSMichael Lange   ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr);
13207937d9ceSMichael Lange   ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
13217937d9ceSMichael Lange   /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
13227937d9ceSMichael Lange   ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr);
13237937d9ceSMichael Lange   /* Rebuild the point strata on the receiver */
13247937d9ceSMichael Lange   for (p = 0, idx = 0; p < nroots; p++) {
13257937d9ceSMichael Lange     for (d = 0; d < rootDegree[p]; d++) {
13267937d9ceSMichael Lange       ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr);
13277937d9ceSMichael Lange       ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr);
13287937d9ceSMichael Lange       for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);}
13297937d9ceSMichael Lange     }
13307937d9ceSMichael Lange     idx += rootDegree[p];
13317937d9ceSMichael Lange   }
133277e0c0e7SMichael Lange   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
133377e0c0e7SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
133477e0c0e7SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
133577e0c0e7SMichael Lange   ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr);
13367937d9ceSMichael Lange   PetscFunctionReturn(0);
13377937d9ceSMichael Lange }
13387937d9ceSMichael Lange 
133984f0b6dfSMatthew G. Knepley /*@
134084f0b6dfSMatthew G. Knepley   DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label
134184f0b6dfSMatthew G. Knepley 
134284f0b6dfSMatthew G. Knepley   Input Parameter:
134384f0b6dfSMatthew G. Knepley . label - the DMLabel
134484f0b6dfSMatthew G. Knepley 
134584f0b6dfSMatthew G. Knepley   Output Parameters:
134684f0b6dfSMatthew G. Knepley + section - the section giving offsets for each stratum
134784f0b6dfSMatthew G. Knepley - is - An IS containing all the label points
134884f0b6dfSMatthew G. Knepley 
134984f0b6dfSMatthew G. Knepley   Level: developer
135084f0b6dfSMatthew G. Knepley 
135184f0b6dfSMatthew G. Knepley .seealso: DMLabelDistribute()
135284f0b6dfSMatthew G. Knepley @*/
1353c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1354c58f1c22SToby Isaac {
1355c58f1c22SToby Isaac   IS              vIS;
1356c58f1c22SToby Isaac   const PetscInt *values;
1357c58f1c22SToby Isaac   PetscInt       *points;
1358c58f1c22SToby Isaac   PetscInt        nV, vS = 0, vE = 0, v, N;
1359c58f1c22SToby Isaac   PetscErrorCode  ierr;
1360c58f1c22SToby Isaac 
1361c58f1c22SToby Isaac   PetscFunctionBegin;
1362d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1363c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr);
1364c58f1c22SToby Isaac   ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
1365c58f1c22SToby Isaac   ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr);
1366c58f1c22SToby Isaac   if (nV) {vS = values[0]; vE = values[0]+1;}
1367c58f1c22SToby Isaac   for (v = 1; v < nV; ++v) {
1368c58f1c22SToby Isaac     vS = PetscMin(vS, values[v]);
1369c58f1c22SToby Isaac     vE = PetscMax(vE, values[v]+1);
1370c58f1c22SToby Isaac   }
1371c58f1c22SToby Isaac   ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr);
1372c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr);
1373c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1374c58f1c22SToby Isaac     PetscInt n;
1375c58f1c22SToby Isaac 
1376c58f1c22SToby Isaac     ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr);
1377c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr);
1378c58f1c22SToby Isaac   }
1379c58f1c22SToby Isaac   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
1380c58f1c22SToby Isaac   ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr);
1381c58f1c22SToby Isaac   ierr = PetscMalloc1(N, &points);CHKERRQ(ierr);
1382c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1383c58f1c22SToby Isaac     IS              is;
1384c58f1c22SToby Isaac     const PetscInt *spoints;
1385c58f1c22SToby Isaac     PetscInt        dof, off, p;
1386c58f1c22SToby Isaac 
1387c58f1c22SToby Isaac     ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr);
1388c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr);
1389c58f1c22SToby Isaac     ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr);
1390c58f1c22SToby Isaac     ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr);
1391c58f1c22SToby Isaac     for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1392c58f1c22SToby Isaac     ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr);
1393c58f1c22SToby Isaac     ierr = ISDestroy(&is);CHKERRQ(ierr);
1394c58f1c22SToby Isaac   }
1395c58f1c22SToby Isaac   ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr);
1396c58f1c22SToby Isaac   ierr = ISDestroy(&vIS);CHKERRQ(ierr);
1397c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1398c58f1c22SToby Isaac   PetscFunctionReturn(0);
1399c58f1c22SToby Isaac }
1400c58f1c22SToby Isaac 
140184f0b6dfSMatthew G. Knepley /*@
1402c58f1c22SToby Isaac   PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1403c58f1c22SToby Isaac   the local section and an SF describing the section point overlap.
1404c58f1c22SToby Isaac 
1405c58f1c22SToby Isaac   Input Parameters:
1406c58f1c22SToby Isaac   + s - The PetscSection for the local field layout
1407c58f1c22SToby Isaac   . sf - The SF describing parallel layout of the section points
1408c58f1c22SToby Isaac   . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1409c58f1c22SToby Isaac   . label - The label specifying the points
1410c58f1c22SToby Isaac   - labelValue - The label stratum specifying the points
1411c58f1c22SToby Isaac 
1412c58f1c22SToby Isaac   Output Parameter:
1413c58f1c22SToby Isaac   . gsection - The PetscSection for the global field layout
1414c58f1c22SToby Isaac 
1415c58f1c22SToby Isaac   Note: This gives negative sizes and offsets to points not owned by this process
1416c58f1c22SToby Isaac 
1417c58f1c22SToby Isaac   Level: developer
1418c58f1c22SToby Isaac 
1419c58f1c22SToby Isaac .seealso: PetscSectionCreate()
1420c58f1c22SToby Isaac @*/
1421c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1422c58f1c22SToby Isaac {
1423c58f1c22SToby Isaac   PetscInt      *neg = NULL, *tmpOff = NULL;
1424c58f1c22SToby Isaac   PetscInt       pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1425c58f1c22SToby Isaac   PetscErrorCode ierr;
1426c58f1c22SToby Isaac 
1427c58f1c22SToby Isaac   PetscFunctionBegin;
1428d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1429d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1430d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4);
1431c58f1c22SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr);
1432c58f1c22SToby Isaac   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1433c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr);
1434c58f1c22SToby Isaac   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1435c58f1c22SToby Isaac   if (nroots >= 0) {
1436c58f1c22SToby Isaac     if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1437c58f1c22SToby Isaac     ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr);
1438c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1439c58f1c22SToby Isaac       ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr);
1440c58f1c22SToby Isaac     } else {
1441c58f1c22SToby Isaac       tmpOff = &(*gsection)->atlasDof[-pStart];
1442c58f1c22SToby Isaac     }
1443c58f1c22SToby Isaac   }
1444c58f1c22SToby Isaac   /* Mark ghost points with negative dof */
1445c58f1c22SToby Isaac   for (p = pStart; p < pEnd; ++p) {
1446c58f1c22SToby Isaac     PetscInt value;
1447c58f1c22SToby Isaac 
1448c58f1c22SToby Isaac     ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr);
1449c58f1c22SToby Isaac     if (value != labelValue) continue;
1450c58f1c22SToby Isaac     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1451c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr);
1452c58f1c22SToby Isaac     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1453c58f1c22SToby Isaac     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);}
1454c58f1c22SToby Isaac     if (neg) neg[p] = -(dof+1);
1455c58f1c22SToby Isaac   }
1456c58f1c22SToby Isaac   ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr);
1457c58f1c22SToby Isaac   if (nroots >= 0) {
1458c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1459c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1460c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1461c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1462c58f1c22SToby Isaac     }
1463c58f1c22SToby Isaac   }
1464c58f1c22SToby Isaac   /* Calculate new sizes, get proccess offset, and calculate point offsets */
1465c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1466c58f1c22SToby Isaac     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1467c58f1c22SToby Isaac     (*gsection)->atlasOff[p] = off;
1468c58f1c22SToby Isaac     off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1469c58f1c22SToby Isaac   }
1470c58f1c22SToby Isaac   ierr       = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1471c58f1c22SToby Isaac   globalOff -= off;
1472c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1473c58f1c22SToby Isaac     (*gsection)->atlasOff[p] += globalOff;
1474c58f1c22SToby Isaac     if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1475c58f1c22SToby Isaac   }
1476c58f1c22SToby Isaac   /* Put in negative offsets for ghost points */
1477c58f1c22SToby Isaac   if (nroots >= 0) {
1478c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1479c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1480c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1481c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1482c58f1c22SToby Isaac     }
1483c58f1c22SToby Isaac   }
1484c58f1c22SToby Isaac   if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);}
1485c58f1c22SToby Isaac   ierr = PetscFree(neg);CHKERRQ(ierr);
1486c58f1c22SToby Isaac   PetscFunctionReturn(0);
1487c58f1c22SToby Isaac }
1488c58f1c22SToby Isaac 
14895fdea053SToby Isaac typedef struct _n_PetscSectionSym_Label
14905fdea053SToby Isaac {
14915fdea053SToby Isaac   DMLabel           label;
14925fdea053SToby Isaac   PetscCopyMode     *modes;
14935fdea053SToby Isaac   PetscInt          *sizes;
14945fdea053SToby Isaac   const PetscInt    ***perms;
14955fdea053SToby Isaac   const PetscScalar ***rots;
14965fdea053SToby Isaac   PetscInt          (*minMaxOrients)[2];
14975fdea053SToby Isaac   PetscInt          numStrata; /* numStrata is only increasing, functions as a state */
14985fdea053SToby Isaac } PetscSectionSym_Label;
14995fdea053SToby Isaac 
15005fdea053SToby Isaac static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
15015fdea053SToby Isaac {
15025fdea053SToby Isaac   PetscInt              i, j;
15035fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
15045fdea053SToby Isaac   PetscErrorCode        ierr;
15055fdea053SToby Isaac 
15065fdea053SToby Isaac   PetscFunctionBegin;
15075fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
15085fdea053SToby Isaac     if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
15095fdea053SToby Isaac       for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
15105fdea053SToby Isaac         if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);}
15115fdea053SToby Isaac         if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);}
15125fdea053SToby Isaac       }
15135fdea053SToby Isaac       if (sl->perms[i]) {
15145fdea053SToby Isaac         const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
15155fdea053SToby Isaac 
15165fdea053SToby Isaac         ierr = PetscFree(perms);CHKERRQ(ierr);
15175fdea053SToby Isaac       }
15185fdea053SToby Isaac       if (sl->rots[i]) {
15195fdea053SToby Isaac         const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
15205fdea053SToby Isaac 
15215fdea053SToby Isaac         ierr = PetscFree(rots);CHKERRQ(ierr);
15225fdea053SToby Isaac       }
15235fdea053SToby Isaac     }
15245fdea053SToby Isaac   }
15255fdea053SToby Isaac   ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr);
15265fdea053SToby Isaac   ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr);
15275fdea053SToby Isaac   sl->numStrata = 0;
15285fdea053SToby Isaac   PetscFunctionReturn(0);
15295fdea053SToby Isaac }
15305fdea053SToby Isaac 
15315fdea053SToby Isaac static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
15325fdea053SToby Isaac {
15335fdea053SToby Isaac   PetscErrorCode ierr;
15345fdea053SToby Isaac 
15355fdea053SToby Isaac   PetscFunctionBegin;
15365fdea053SToby Isaac   ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);
15375fdea053SToby Isaac   ierr = PetscFree(sym->data);CHKERRQ(ierr);
15385fdea053SToby Isaac   PetscFunctionReturn(0);
15395fdea053SToby Isaac }
15405fdea053SToby Isaac 
15415fdea053SToby Isaac static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
15425fdea053SToby Isaac {
15435fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
15445fdea053SToby Isaac   PetscBool             isAscii;
15455fdea053SToby Isaac   DMLabel               label = sl->label;
1546d67d17b1SMatthew G. Knepley   const char           *name;
15475fdea053SToby Isaac   PetscErrorCode        ierr;
15485fdea053SToby Isaac 
15495fdea053SToby Isaac   PetscFunctionBegin;
15505fdea053SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
15515fdea053SToby Isaac   if (isAscii) {
15525fdea053SToby Isaac     PetscInt          i, j, k;
15535fdea053SToby Isaac     PetscViewerFormat format;
15545fdea053SToby Isaac 
15555fdea053SToby Isaac     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
15565fdea053SToby Isaac     if (label) {
15575fdea053SToby Isaac       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
15585fdea053SToby Isaac       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
15595fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15605fdea053SToby Isaac         ierr = DMLabelView(label, viewer);CHKERRQ(ierr);
15615fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
15625fdea053SToby Isaac       } else {
1563d67d17b1SMatthew G. Knepley         ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1564d67d17b1SMatthew G. Knepley         ierr = PetscViewerASCIIPrintf(viewer,"  Label '%s'\n",name);CHKERRQ(ierr);
15655fdea053SToby Isaac       }
15665fdea053SToby Isaac     } else {
15675fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr);
15685fdea053SToby Isaac     }
15695fdea053SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15705fdea053SToby Isaac     for (i = 0; i <= sl->numStrata; i++) {
15715fdea053SToby Isaac       PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
15725fdea053SToby Isaac 
15735fdea053SToby Isaac       if (!(sl->perms[i] || sl->rots[i])) {
15745fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr);
15755fdea053SToby Isaac       } else {
15765fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr);
15775fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15785fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr);
15795fdea053SToby Isaac         if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
15805fdea053SToby Isaac           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15815fdea053SToby Isaac           for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
15825fdea053SToby Isaac             if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
15835fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr);
15845fdea053SToby Isaac             } else {
15855fdea053SToby Isaac               PetscInt tab;
15865fdea053SToby Isaac 
15875fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr);
15885fdea053SToby Isaac               ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15895fdea053SToby Isaac               ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr);
15905fdea053SToby Isaac               if (sl->perms[i] && sl->perms[i][j]) {
15915fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr);
15925fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
15935fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);}
15945fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
15955fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
15965fdea053SToby Isaac               }
15975fdea053SToby Isaac               if (sl->rots[i] && sl->rots[i][j]) {
15985fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Rotations:  ");CHKERRQ(ierr);
15995fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
16005fdea053SToby Isaac #if defined(PETSC_USE_COMPLEX)
16015fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f+i*%+f",PetscRealPart(sl->rots[i][j][k]),PetscImaginaryPart(sl->rots[i][j][k]));CHKERRQ(ierr);}
16025fdea053SToby Isaac #else
16035fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);}
16045fdea053SToby Isaac #endif
16055fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
16065fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
16075fdea053SToby Isaac               }
16085fdea053SToby Isaac               ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
16095fdea053SToby Isaac             }
16105fdea053SToby Isaac           }
16115fdea053SToby Isaac           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
16125fdea053SToby Isaac         }
16135fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
16145fdea053SToby Isaac       }
16155fdea053SToby Isaac     }
16165fdea053SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
16175fdea053SToby Isaac   }
16185fdea053SToby Isaac   PetscFunctionReturn(0);
16195fdea053SToby Isaac }
16205fdea053SToby Isaac 
16215fdea053SToby Isaac /*@
16225fdea053SToby Isaac   PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
16235fdea053SToby Isaac 
16245fdea053SToby Isaac   Logically collective on sym
16255fdea053SToby Isaac 
16265fdea053SToby Isaac   Input parameters:
16275fdea053SToby Isaac + sym - the section symmetries
16285fdea053SToby Isaac - label - the DMLabel describing the types of points
16295fdea053SToby Isaac 
16305fdea053SToby Isaac   Level: developer:
16315fdea053SToby Isaac 
16325fdea053SToby Isaac .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()
16335fdea053SToby Isaac @*/
16345fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
16355fdea053SToby Isaac {
16365fdea053SToby Isaac   PetscSectionSym_Label *sl;
16375fdea053SToby Isaac   PetscErrorCode        ierr;
16385fdea053SToby Isaac 
16395fdea053SToby Isaac   PetscFunctionBegin;
16405fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
16415fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
16425fdea053SToby Isaac   if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);}
16435fdea053SToby Isaac   if (label) {
16445fdea053SToby Isaac     sl->label = label;
1645d67d17b1SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) label);CHKERRQ(ierr);
16465fdea053SToby Isaac     ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr);
16471a834cf9SToby Isaac     ierr = PetscMalloc5(sl->numStrata+1,&sl->modes,sl->numStrata+1,&sl->sizes,sl->numStrata+1,&sl->perms,sl->numStrata+1,&sl->rots,sl->numStrata+1,&sl->minMaxOrients);CHKERRQ(ierr);
16481a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr);
16491a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr);
16501a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr);
16511a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr);
16521a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr);
16535fdea053SToby Isaac   }
16545fdea053SToby Isaac   PetscFunctionReturn(0);
16555fdea053SToby Isaac }
16565fdea053SToby Isaac 
16575fdea053SToby Isaac /*@C
16585fdea053SToby Isaac   PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
16595fdea053SToby Isaac 
16605fdea053SToby Isaac   Logically collective on PetscSectionSym
16615fdea053SToby Isaac 
16625fdea053SToby Isaac   InputParameters:
16635fdea053SToby Isaac + sys       - the section symmetries
16645fdea053SToby Isaac . stratum   - the stratum value in the label that we are assigning symmetries for
16655fdea053SToby Isaac . size      - the number of dofs for points in the stratum of the label
16665fdea053SToby Isaac . minOrient - the smallest orientation for a point in this stratum
16675fdea053SToby Isaac . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
16685fdea053SToby Isaac . mode      - how sym should copy the perms and rots arrays
16695fdea053SToby Isaac . perms     - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation.  A NULL permutation is the identity
16705fdea053SToby Isaac + 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
16715fdea053SToby Isaac 
16725fdea053SToby Isaac   Level: developer
16735fdea053SToby Isaac 
16745fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()
16755fdea053SToby Isaac @*/
16765fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
16775fdea053SToby Isaac {
16785fdea053SToby Isaac   PetscSectionSym_Label *sl;
1679d67d17b1SMatthew G. Knepley   const char            *name;
1680d67d17b1SMatthew G. Knepley   PetscInt               i, j, k;
16815fdea053SToby Isaac   PetscErrorCode         ierr;
16825fdea053SToby Isaac 
16835fdea053SToby Isaac   PetscFunctionBegin;
16845fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
16855fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
16865fdea053SToby Isaac   if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet");
16875fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
16885fdea053SToby Isaac     PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
16895fdea053SToby Isaac 
16905fdea053SToby Isaac     if (stratum == value) break;
16915fdea053SToby Isaac   }
1692d67d17b1SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1693d67d17b1SMatthew G. Knepley   if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,name);
16945fdea053SToby Isaac   sl->sizes[i] = size;
16955fdea053SToby Isaac   sl->modes[i] = mode;
16965fdea053SToby Isaac   sl->minMaxOrients[i][0] = minOrient;
16975fdea053SToby Isaac   sl->minMaxOrients[i][1] = maxOrient;
16985fdea053SToby Isaac   if (mode == PETSC_COPY_VALUES) {
16995fdea053SToby Isaac     if (perms) {
17005fdea053SToby Isaac       PetscInt    **ownPerms;
17015fdea053SToby Isaac 
17025fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr);
17035fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
17045fdea053SToby Isaac         if (perms[j]) {
17055fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr);
17065fdea053SToby Isaac           for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];}
17075fdea053SToby Isaac         }
17085fdea053SToby Isaac       }
17095fdea053SToby Isaac       sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient];
17105fdea053SToby Isaac     }
17115fdea053SToby Isaac     if (rots) {
17125fdea053SToby Isaac       PetscScalar **ownRots;
17135fdea053SToby Isaac 
17145fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr);
17155fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
17165fdea053SToby Isaac         if (rots[j]) {
17175fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr);
17185fdea053SToby Isaac           for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];}
17195fdea053SToby Isaac         }
17205fdea053SToby Isaac       }
17215fdea053SToby Isaac       sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient];
17225fdea053SToby Isaac     }
17235fdea053SToby Isaac   } else {
17245fdea053SToby Isaac     sl->perms[i] = perms ? &perms[-minOrient] : NULL;
17255fdea053SToby Isaac     sl->rots[i]  = rots ? &rots[-minOrient] : NULL;
17265fdea053SToby Isaac   }
17275fdea053SToby Isaac   PetscFunctionReturn(0);
17285fdea053SToby Isaac }
17295fdea053SToby Isaac 
17305fdea053SToby Isaac static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
17315fdea053SToby Isaac {
17325fdea053SToby Isaac   PetscInt              i, j, numStrata;
17335fdea053SToby Isaac   PetscSectionSym_Label *sl;
17345fdea053SToby Isaac   DMLabel               label;
17355fdea053SToby Isaac   PetscErrorCode        ierr;
17365fdea053SToby Isaac 
17375fdea053SToby Isaac   PetscFunctionBegin;
17385fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
17395fdea053SToby Isaac   numStrata = sl->numStrata;
17405fdea053SToby Isaac   label     = sl->label;
17415fdea053SToby Isaac   for (i = 0; i < numPoints; i++) {
17425fdea053SToby Isaac     PetscInt point = points[2*i];
17435fdea053SToby Isaac     PetscInt ornt  = points[2*i+1];
17445fdea053SToby Isaac 
17455fdea053SToby Isaac     for (j = 0; j < numStrata; j++) {
17465fdea053SToby Isaac       if (label->validIS[j]) {
17475fdea053SToby Isaac         PetscInt k;
17485fdea053SToby Isaac 
17495fdea053SToby Isaac         ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr);
17505fdea053SToby Isaac         if (k >= 0) break;
17515fdea053SToby Isaac       } else {
17525fdea053SToby Isaac         PetscBool has;
17535fdea053SToby Isaac 
1754e8f14785SLisandro Dalcin         PetscHSetIHas(label->ht[j], point, &has);
17555fdea053SToby Isaac         if (has) break;
17565fdea053SToby Isaac       }
17575fdea053SToby Isaac     }
17585fdea053SToby Isaac     if ((sl->minMaxOrients[j][1] > sl->minMaxOrients[j][0]) && (ornt < sl->minMaxOrients[j][0] || ornt >= sl->minMaxOrients[j][1])) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D orientation %D not in range [%D, %D) for stratum %D",point,ornt,sl->minMaxOrients[j][0],sl->minMaxOrients[j][1],j < numStrata ? label->stratumValues[j] : label->defaultValue);
17595fdea053SToby Isaac     if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;}
17605fdea053SToby Isaac     if (rots) {rots[i]  = sl->rots[j] ? sl->rots[j][ornt] : NULL;}
17615fdea053SToby Isaac   }
17625fdea053SToby Isaac   PetscFunctionReturn(0);
17635fdea053SToby Isaac }
17645fdea053SToby Isaac 
17655fdea053SToby Isaac PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
17665fdea053SToby Isaac {
17675fdea053SToby Isaac   PetscSectionSym_Label *sl;
17685fdea053SToby Isaac   PetscErrorCode        ierr;
17695fdea053SToby Isaac 
17705fdea053SToby Isaac   PetscFunctionBegin;
17715fdea053SToby Isaac   ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr);
17725fdea053SToby Isaac   sym->ops->getpoints = PetscSectionSymGetPoints_Label;
17735fdea053SToby Isaac   sym->ops->view      = PetscSectionSymView_Label;
17745fdea053SToby Isaac   sym->ops->destroy   = PetscSectionSymDestroy_Label;
17755fdea053SToby Isaac   sym->data           = (void *) sl;
17765fdea053SToby Isaac   PetscFunctionReturn(0);
17775fdea053SToby Isaac }
17785fdea053SToby Isaac 
17795fdea053SToby Isaac /*@
17805fdea053SToby Isaac   PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
17815fdea053SToby Isaac 
17825fdea053SToby Isaac   Collective
17835fdea053SToby Isaac 
17845fdea053SToby Isaac   Input Parameters:
17855fdea053SToby Isaac + comm - the MPI communicator for the new symmetry
17865fdea053SToby Isaac - label - the label defining the strata
17875fdea053SToby Isaac 
17885fdea053SToby Isaac   Output Parameters:
17895fdea053SToby Isaac . sym - the section symmetries
17905fdea053SToby Isaac 
17915fdea053SToby Isaac   Level: developer
17925fdea053SToby Isaac 
17935fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()
17945fdea053SToby Isaac @*/
17955fdea053SToby Isaac PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
17965fdea053SToby Isaac {
17975fdea053SToby Isaac   PetscErrorCode        ierr;
17985fdea053SToby Isaac 
17995fdea053SToby Isaac   PetscFunctionBegin;
18005fdea053SToby Isaac   ierr = DMInitializePackage();CHKERRQ(ierr);
18015fdea053SToby Isaac   ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr);
18025fdea053SToby Isaac   ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr);
18035fdea053SToby Isaac   ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr);
18045fdea053SToby Isaac   PetscFunctionReturn(0);
18055fdea053SToby Isaac }
1806