xref: /petsc/src/dm/label/dmlabel.c (revision a205f1fd2f49e1bbb35bc9f87225a617b4893275)
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 
9c58f1c22SToby Isaac   Input parameter:
10c58f1c22SToby Isaac . name - The label name
11c58f1c22SToby Isaac 
12c58f1c22SToby Isaac   Output parameter:
13c58f1c22SToby Isaac . label - The DMLabel
14c58f1c22SToby Isaac 
15c58f1c22SToby Isaac   Level: beginner
16c58f1c22SToby Isaac 
17c58f1c22SToby Isaac .seealso: DMLabelDestroy()
18c58f1c22SToby Isaac @*/
19c58f1c22SToby Isaac PetscErrorCode DMLabelCreate(const char name[], DMLabel *label)
20c58f1c22SToby Isaac {
21c58f1c22SToby Isaac   PetscErrorCode ierr;
22c58f1c22SToby Isaac 
23c58f1c22SToby Isaac   PetscFunctionBegin;
24c58f1c22SToby Isaac   ierr = PetscNew(label);CHKERRQ(ierr);
25c58f1c22SToby Isaac   ierr = PetscStrallocpy(name, &(*label)->name);CHKERRQ(ierr);
26c58f1c22SToby Isaac 
27c58f1c22SToby Isaac   (*label)->refct          = 1;
28c58f1c22SToby Isaac   (*label)->state          = -1;
29c58f1c22SToby Isaac   (*label)->numStrata      = 0;
305aa44df4SToby Isaac   (*label)->defaultValue   = -1;
31c58f1c22SToby Isaac   (*label)->stratumValues  = NULL;
32ad8374ffSToby Isaac   (*label)->validIS        = NULL;
33c58f1c22SToby Isaac   (*label)->stratumSizes   = NULL;
34c58f1c22SToby Isaac   (*label)->points         = NULL;
35c58f1c22SToby Isaac   (*label)->ht             = NULL;
36c58f1c22SToby Isaac   (*label)->pStart         = -1;
37c58f1c22SToby Isaac   (*label)->pEnd           = -1;
38c58f1c22SToby Isaac   (*label)->bt             = NULL;
39c58f1c22SToby Isaac   PetscFunctionReturn(0);
40c58f1c22SToby Isaac }
41c58f1c22SToby Isaac 
42c58f1c22SToby Isaac /*
43c58f1c22SToby Isaac   DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format
44c58f1c22SToby Isaac 
45c58f1c22SToby Isaac   Input parameter:
46c58f1c22SToby Isaac + label - The DMLabel
47c58f1c22SToby Isaac - v - The stratum value
48c58f1c22SToby Isaac 
49c58f1c22SToby Isaac   Output parameter:
50c58f1c22SToby Isaac . label - The DMLabel with stratum in sorted list format
51c58f1c22SToby Isaac 
52c58f1c22SToby Isaac   Level: developer
53c58f1c22SToby Isaac 
54c58f1c22SToby Isaac .seealso: DMLabelCreate()
55c58f1c22SToby Isaac */
56c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v)
57c58f1c22SToby Isaac {
58c58f1c22SToby Isaac   PetscInt       off;
59ad8374ffSToby Isaac   PetscInt       *pointArray;
60c58f1c22SToby Isaac   PetscErrorCode ierr;
61c58f1c22SToby Isaac 
62ad8374ffSToby Isaac   if (label->validIS[v]) return 0;
63c58f1c22SToby Isaac   if (v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeValid_Private\n", v);
64c58f1c22SToby Isaac   PetscFunctionBegin;
65c58f1c22SToby Isaac   PetscHashISize(label->ht[v], label->stratumSizes[v]);
66c58f1c22SToby Isaac 
67ad8374ffSToby Isaac   ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr);
68c58f1c22SToby Isaac   off  = 0;
69ad8374ffSToby Isaac   ierr = PetscHashIGetKeys(label->ht[v], &off, pointArray);CHKERRQ(ierr);
70c58f1c22SToby 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]);
71c58f1c22SToby Isaac   PetscHashIClear(label->ht[v]);
72ad8374ffSToby Isaac   ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr);
73c58f1c22SToby Isaac   if (label->bt) {
74c58f1c22SToby Isaac     PetscInt p;
75c58f1c22SToby Isaac 
76c58f1c22SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
77ad8374ffSToby Isaac       const PetscInt point = pointArray[p];
78c58f1c22SToby Isaac 
79c58f1c22SToby 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);
80c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
81c58f1c22SToby Isaac     }
82c58f1c22SToby Isaac   }
83ad8374ffSToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF,label->stratumSizes[v],pointArray,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
84ad8374ffSToby Isaac   ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
85ad8374ffSToby Isaac   label->validIS[v] = PETSC_TRUE;
86c58f1c22SToby Isaac   ++label->state;
87c58f1c22SToby Isaac   PetscFunctionReturn(0);
88c58f1c22SToby Isaac }
89c58f1c22SToby Isaac 
90c58f1c22SToby Isaac /*
91c58f1c22SToby Isaac   DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format
92c58f1c22SToby Isaac 
93c58f1c22SToby Isaac   Input parameter:
94c58f1c22SToby Isaac . label - The DMLabel
95c58f1c22SToby Isaac 
96c58f1c22SToby Isaac   Output parameter:
97c58f1c22SToby Isaac . label - The DMLabel with all strata in sorted list format
98c58f1c22SToby Isaac 
99c58f1c22SToby Isaac   Level: developer
100c58f1c22SToby Isaac 
101c58f1c22SToby Isaac .seealso: DMLabelCreate()
102c58f1c22SToby Isaac */
103c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label)
104c58f1c22SToby Isaac {
105c58f1c22SToby Isaac   PetscInt       v;
106c58f1c22SToby Isaac   PetscErrorCode ierr;
107c58f1c22SToby Isaac 
108c58f1c22SToby Isaac   PetscFunctionBegin;
109c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; v++){
110c58f1c22SToby Isaac     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
111c58f1c22SToby Isaac   }
112c58f1c22SToby Isaac   PetscFunctionReturn(0);
113c58f1c22SToby Isaac }
114c58f1c22SToby Isaac 
115c58f1c22SToby Isaac /*
116c58f1c22SToby Isaac   DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format
117c58f1c22SToby Isaac 
118c58f1c22SToby Isaac   Input parameter:
119c58f1c22SToby Isaac + label - The DMLabel
120c58f1c22SToby Isaac - v - The stratum value
121c58f1c22SToby Isaac 
122c58f1c22SToby Isaac   Output parameter:
123c58f1c22SToby Isaac . label - The DMLabel with stratum in hash format
124c58f1c22SToby Isaac 
125c58f1c22SToby Isaac   Level: developer
126c58f1c22SToby Isaac 
127c58f1c22SToby Isaac .seealso: DMLabelCreate()
128c58f1c22SToby Isaac */
129c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v)
130c58f1c22SToby Isaac {
131c58f1c22SToby Isaac   PETSC_UNUSED PetscHashIIter ret, iter;
132c58f1c22SToby Isaac   PetscInt                    p;
133ad8374ffSToby Isaac   const PetscInt              *points;
134c58f1c22SToby Isaac   PetscErrorCode ierr;
135c58f1c22SToby Isaac 
136c58f1c22SToby Isaac   PetscFunctionBegin;
137ad8374ffSToby Isaac   if (!label->validIS[v]) PetscFunctionReturn(0);
138ad8374ffSToby Isaac   if (label->points[v]) {
139ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
140ad8374ffSToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) PetscHashIPut(label->ht[v], points[p], ret, iter);
141ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
142ad8374ffSToby Isaac     ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
143ad8374ffSToby Isaac   }
144ad8374ffSToby Isaac   label->validIS[v] = PETSC_FALSE;
145c58f1c22SToby Isaac   PetscFunctionReturn(0);
146c58f1c22SToby Isaac }
147c58f1c22SToby Isaac 
148c58f1c22SToby Isaac PetscErrorCode DMLabelGetState(DMLabel label, PetscObjectState *state)
149c58f1c22SToby Isaac {
150c58f1c22SToby Isaac   PetscFunctionBegin;
151c58f1c22SToby Isaac   PetscValidPointer(state, 2);
152c58f1c22SToby Isaac   *state = label->state;
153c58f1c22SToby Isaac   PetscFunctionReturn(0);
154c58f1c22SToby Isaac }
155c58f1c22SToby Isaac 
156c58f1c22SToby Isaac PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
157c58f1c22SToby Isaac {
158ad8374ffSToby Isaac   PetscInt    v, *tmpV, *tmpS;
159ad8374ffSToby Isaac   IS         *tmpP;
160c58f1c22SToby Isaac   PetscHashI *tmpH;
161c58f1c22SToby Isaac   PetscBool  *tmpB;
162c58f1c22SToby Isaac   PetscErrorCode ierr;
163c58f1c22SToby Isaac 
164c58f1c22SToby Isaac   PetscFunctionBegin;
165c58f1c22SToby Isaac 
166ad8374ffSToby Isaac   for (v = 0; v < label->numStrata; v++) {
167ad8374ffSToby Isaac     if (label->stratumValues[v] == value) PetscFunctionReturn(0);
168ad8374ffSToby Isaac   }
169c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpV);CHKERRQ(ierr);
170c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpS);CHKERRQ(ierr);
171c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpH);CHKERRQ(ierr);
172c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpP);CHKERRQ(ierr);
173c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpB);CHKERRQ(ierr);
174c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
175c58f1c22SToby Isaac     tmpV[v] = label->stratumValues[v];
176c58f1c22SToby Isaac     tmpS[v] = label->stratumSizes[v];
177c58f1c22SToby Isaac     tmpH[v] = label->ht[v];
178c58f1c22SToby Isaac     tmpP[v] = label->points[v];
179ad8374ffSToby Isaac     tmpB[v] = label->validIS[v];
180c58f1c22SToby Isaac   }
181c58f1c22SToby Isaac   tmpV[v] = value;
182c58f1c22SToby Isaac   tmpS[v] = 0;
183c58f1c22SToby Isaac   PetscHashICreate(tmpH[v]);
184ad8374ffSToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&tmpP[v]);CHKERRQ(ierr);
185c58f1c22SToby Isaac   tmpB[v] = PETSC_TRUE;
186c58f1c22SToby Isaac   ++label->numStrata;
187c58f1c22SToby Isaac   ierr = PetscFree(label->stratumValues);CHKERRQ(ierr);
188c58f1c22SToby Isaac   ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr);
189c58f1c22SToby Isaac   ierr = PetscFree(label->ht);CHKERRQ(ierr);
190c58f1c22SToby Isaac   ierr = PetscFree(label->points);CHKERRQ(ierr);
191ad8374ffSToby Isaac   ierr = PetscFree(label->validIS);CHKERRQ(ierr);
192c58f1c22SToby Isaac   label->stratumValues = tmpV;
193c58f1c22SToby Isaac   label->stratumSizes  = tmpS;
194c58f1c22SToby Isaac   label->ht            = tmpH;
195c58f1c22SToby Isaac   label->points        = tmpP;
196ad8374ffSToby Isaac   label->validIS       = tmpB;
197c58f1c22SToby Isaac 
198c58f1c22SToby Isaac   PetscFunctionReturn(0);
199c58f1c22SToby Isaac }
200c58f1c22SToby Isaac 
201c58f1c22SToby Isaac /*@C
202c58f1c22SToby Isaac   DMLabelGetName - Return the name of a DMLabel object
203c58f1c22SToby Isaac 
204c58f1c22SToby Isaac   Input parameter:
205c58f1c22SToby Isaac . label - The DMLabel
206c58f1c22SToby Isaac 
207c58f1c22SToby Isaac   Output parameter:
208c58f1c22SToby Isaac . name - The label name
209c58f1c22SToby Isaac 
210c58f1c22SToby Isaac   Level: beginner
211c58f1c22SToby Isaac 
212c58f1c22SToby Isaac .seealso: DMLabelCreate()
213c58f1c22SToby Isaac @*/
214c58f1c22SToby Isaac PetscErrorCode DMLabelGetName(DMLabel label, const char **name)
215c58f1c22SToby Isaac {
216c58f1c22SToby Isaac   PetscFunctionBegin;
217c58f1c22SToby Isaac   PetscValidPointer(name, 2);
218c58f1c22SToby Isaac   *name = label->name;
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) {
232c58f1c22SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", label->name);CHKERRQ(ierr);
233c58f1c22SToby Isaac     if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);}
234c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
235c58f1c22SToby Isaac       const PetscInt value = label->stratumValues[v];
236ad8374ffSToby Isaac       const PetscInt *points;
237c58f1c22SToby Isaac       PetscInt       p;
238c58f1c22SToby Isaac 
239ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
240c58f1c22SToby Isaac       for (p = 0; p < label->stratumSizes[v]; ++p) {
241ad8374ffSToby Isaac         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr);
242c58f1c22SToby Isaac       }
243ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
244c58f1c22SToby Isaac     }
245c58f1c22SToby Isaac   }
246c58f1c22SToby Isaac   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
247c58f1c22SToby Isaac   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
248c58f1c22SToby Isaac   PetscFunctionReturn(0);
249c58f1c22SToby Isaac }
250c58f1c22SToby Isaac 
251c58f1c22SToby Isaac /*@C
252c58f1c22SToby Isaac   DMLabelView - View the label
253c58f1c22SToby Isaac 
254c58f1c22SToby Isaac   Input Parameters:
255c58f1c22SToby Isaac + label - The DMLabel
256c58f1c22SToby Isaac - viewer - The PetscViewer
257c58f1c22SToby Isaac 
258c58f1c22SToby Isaac   Level: intermediate
259c58f1c22SToby Isaac 
260c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy()
261c58f1c22SToby Isaac @*/
262c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
263c58f1c22SToby Isaac {
264c58f1c22SToby Isaac   PetscBool      iascii;
265c58f1c22SToby Isaac   PetscErrorCode ierr;
266c58f1c22SToby Isaac 
267c58f1c22SToby Isaac   PetscFunctionBegin;
268c58f1c22SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
269c58f1c22SToby Isaac   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
270c58f1c22SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
271c58f1c22SToby Isaac   if (iascii) {
272c58f1c22SToby Isaac     ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr);
273c58f1c22SToby Isaac   }
274c58f1c22SToby Isaac   PetscFunctionReturn(0);
275c58f1c22SToby Isaac }
276c58f1c22SToby Isaac 
277c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label)
278c58f1c22SToby Isaac {
279c58f1c22SToby Isaac   PetscInt       v;
280c58f1c22SToby Isaac   PetscErrorCode ierr;
281c58f1c22SToby Isaac 
282c58f1c22SToby Isaac   PetscFunctionBegin;
283c58f1c22SToby Isaac   if (!(*label)) PetscFunctionReturn(0);
284c58f1c22SToby Isaac   if (--(*label)->refct > 0) PetscFunctionReturn(0);
285c58f1c22SToby Isaac   ierr = PetscFree((*label)->name);CHKERRQ(ierr);
286c58f1c22SToby Isaac   ierr = PetscFree((*label)->stratumValues);CHKERRQ(ierr);
287c58f1c22SToby Isaac   ierr = PetscFree((*label)->stratumSizes);CHKERRQ(ierr);
288ad8374ffSToby Isaac   for (v = 0; v < (*label)->numStrata; ++v) {ierr = ISDestroy(&((*label)->points[v]));CHKERRQ(ierr);}
289c58f1c22SToby Isaac   ierr = PetscFree((*label)->points);CHKERRQ(ierr);
290ad8374ffSToby Isaac   ierr = PetscFree((*label)->validIS);CHKERRQ(ierr);
291c58f1c22SToby Isaac   if ((*label)->ht) {
292c58f1c22SToby Isaac     for (v = 0; v < (*label)->numStrata; ++v) {PetscHashIDestroy((*label)->ht[v]);}
293c58f1c22SToby Isaac     ierr = PetscFree((*label)->ht);CHKERRQ(ierr);
294c58f1c22SToby Isaac   }
295c58f1c22SToby Isaac   ierr = PetscBTDestroy(&(*label)->bt);CHKERRQ(ierr);
296c58f1c22SToby Isaac   ierr = PetscFree(*label);CHKERRQ(ierr);
297c58f1c22SToby Isaac   PetscFunctionReturn(0);
298c58f1c22SToby Isaac }
299c58f1c22SToby Isaac 
300c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
301c58f1c22SToby Isaac {
302ad8374ffSToby Isaac   PetscInt       v;
303c58f1c22SToby Isaac   PetscErrorCode ierr;
304c58f1c22SToby Isaac 
305c58f1c22SToby Isaac   PetscFunctionBegin;
306c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
307c58f1c22SToby Isaac   ierr = PetscNew(labelnew);CHKERRQ(ierr);
308c58f1c22SToby Isaac   ierr = PetscStrallocpy(label->name, &(*labelnew)->name);CHKERRQ(ierr);
309c58f1c22SToby Isaac 
310c58f1c22SToby Isaac   (*labelnew)->refct        = 1;
311c58f1c22SToby Isaac   (*labelnew)->numStrata    = label->numStrata;
3125aa44df4SToby Isaac   (*labelnew)->defaultValue = label->defaultValue;
313c58f1c22SToby Isaac   if (label->numStrata) {
314c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr);
315c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr);
316c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr);
317c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr);
318ad8374ffSToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr);
319c58f1c22SToby Isaac     /* Could eliminate unused space here */
320c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
321c58f1c22SToby Isaac       PetscHashICreate((*labelnew)->ht[v]);
322ad8374ffSToby Isaac       (*labelnew)->validIS[v]        = PETSC_TRUE;
323c58f1c22SToby Isaac       (*labelnew)->stratumValues[v]  = label->stratumValues[v];
324c58f1c22SToby Isaac       (*labelnew)->stratumSizes[v]   = label->stratumSizes[v];
325ad8374ffSToby Isaac       ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr);
326ad8374ffSToby Isaac       (*labelnew)->points[v]         = label->points[v];
327c58f1c22SToby Isaac     }
328c58f1c22SToby Isaac   }
329c58f1c22SToby Isaac   (*labelnew)->pStart = -1;
330c58f1c22SToby Isaac   (*labelnew)->pEnd   = -1;
331c58f1c22SToby Isaac   (*labelnew)->bt     = NULL;
332c58f1c22SToby Isaac   PetscFunctionReturn(0);
333c58f1c22SToby Isaac }
334c58f1c22SToby Isaac 
335c58f1c22SToby Isaac /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
336c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
337c58f1c22SToby Isaac {
338c58f1c22SToby Isaac   PetscInt       v;
339c58f1c22SToby Isaac   PetscErrorCode ierr;
340c58f1c22SToby Isaac 
341c58f1c22SToby Isaac   PetscFunctionBegin;
342c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
343c58f1c22SToby Isaac   if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);}
344c58f1c22SToby Isaac   label->pStart = pStart;
345c58f1c22SToby Isaac   label->pEnd   = pEnd;
346c58f1c22SToby Isaac   ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr);
347c58f1c22SToby Isaac   ierr = PetscBTMemzero(pEnd - pStart, label->bt);CHKERRQ(ierr);
348c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
349ad8374ffSToby Isaac     const PetscInt *points;
350c58f1c22SToby Isaac     PetscInt       i;
351c58f1c22SToby Isaac 
352ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
353c58f1c22SToby Isaac     for (i = 0; i < label->stratumSizes[v]; ++i) {
354ad8374ffSToby Isaac       const PetscInt point = points[i];
355c58f1c22SToby Isaac 
356c58f1c22SToby 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);
357c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr);
358c58f1c22SToby Isaac     }
359ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
360c58f1c22SToby Isaac   }
361c58f1c22SToby Isaac   PetscFunctionReturn(0);
362c58f1c22SToby Isaac }
363c58f1c22SToby Isaac 
364c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label)
365c58f1c22SToby Isaac {
366c58f1c22SToby Isaac   PetscErrorCode ierr;
367c58f1c22SToby Isaac 
368c58f1c22SToby Isaac   PetscFunctionBegin;
369c58f1c22SToby Isaac   label->pStart = -1;
370c58f1c22SToby Isaac   label->pEnd   = -1;
371c58f1c22SToby Isaac   if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);}
372c58f1c22SToby Isaac   PetscFunctionReturn(0);
373c58f1c22SToby Isaac }
374c58f1c22SToby Isaac 
375c58f1c22SToby Isaac /*@
376c58f1c22SToby Isaac   DMLabelHasValue - Determine whether a label assigns the value to any point
377c58f1c22SToby Isaac 
378c58f1c22SToby Isaac   Input Parameters:
379c58f1c22SToby Isaac + label - the DMLabel
380c58f1c22SToby Isaac - value - the value
381c58f1c22SToby Isaac 
382c58f1c22SToby Isaac   Output Parameter:
383c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point
384c58f1c22SToby Isaac 
385c58f1c22SToby Isaac   Level: developer
386c58f1c22SToby Isaac 
387c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
388c58f1c22SToby Isaac @*/
389c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
390c58f1c22SToby Isaac {
391c58f1c22SToby Isaac   PetscInt v;
392c58f1c22SToby Isaac 
393c58f1c22SToby Isaac   PetscFunctionBegin;
394c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
395c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
396c58f1c22SToby Isaac     if (value == label->stratumValues[v]) break;
397c58f1c22SToby Isaac   }
398c58f1c22SToby Isaac   *contains = (v < label->numStrata ? PETSC_TRUE : PETSC_FALSE);
399c58f1c22SToby Isaac   PetscFunctionReturn(0);
400c58f1c22SToby Isaac }
401c58f1c22SToby Isaac 
402c58f1c22SToby Isaac /*@
403c58f1c22SToby Isaac   DMLabelHasPoint - Determine whether a label assigns a value to a point
404c58f1c22SToby Isaac 
405c58f1c22SToby Isaac   Input Parameters:
406c58f1c22SToby Isaac + label - the DMLabel
407c58f1c22SToby Isaac - point - the point
408c58f1c22SToby Isaac 
409c58f1c22SToby Isaac   Output Parameter:
410c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value
411c58f1c22SToby Isaac 
412c58f1c22SToby Isaac   Note: The user must call DMLabelCreateIndex() before this function.
413c58f1c22SToby Isaac 
414c58f1c22SToby Isaac   Level: developer
415c58f1c22SToby Isaac 
416c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
417c58f1c22SToby Isaac @*/
418c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
419c58f1c22SToby Isaac {
420c58f1c22SToby Isaac   PetscErrorCode ierr;
421c58f1c22SToby Isaac 
422c58f1c22SToby Isaac   PetscFunctionBeginHot;
423c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
424c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
425c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG)
426c58f1c22SToby Isaac   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
427c58f1c22SToby 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);
428c58f1c22SToby Isaac #endif
429c58f1c22SToby Isaac   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
430c58f1c22SToby Isaac   PetscFunctionReturn(0);
431c58f1c22SToby Isaac }
432c58f1c22SToby Isaac 
433c58f1c22SToby Isaac /*@
434c58f1c22SToby Isaac   DMLabelStratumHasPoint - Return true if the stratum contains a point
435c58f1c22SToby Isaac 
436c58f1c22SToby Isaac   Input Parameters:
437c58f1c22SToby Isaac + label - the DMLabel
438c58f1c22SToby Isaac . value - the stratum value
439c58f1c22SToby Isaac - point - the point
440c58f1c22SToby Isaac 
441c58f1c22SToby Isaac   Output Parameter:
442c58f1c22SToby Isaac . contains - true if the stratum contains the point
443c58f1c22SToby Isaac 
444c58f1c22SToby Isaac   Level: intermediate
445c58f1c22SToby Isaac 
446c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
447c58f1c22SToby Isaac @*/
448c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
449c58f1c22SToby Isaac {
450c58f1c22SToby Isaac   PetscInt       v;
451c58f1c22SToby Isaac   PetscErrorCode ierr;
452c58f1c22SToby Isaac 
453c58f1c22SToby Isaac   PetscFunctionBegin;
454c58f1c22SToby Isaac   PetscValidPointer(contains, 4);
455c58f1c22SToby Isaac   *contains = PETSC_FALSE;
456c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
457c58f1c22SToby Isaac     if (label->stratumValues[v] == value) {
458ad8374ffSToby Isaac       if (label->validIS[v]) {
459c58f1c22SToby Isaac         PetscInt i;
460c58f1c22SToby Isaac 
461a2d74346SToby Isaac         ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr);
462c58f1c22SToby Isaac         if (i >= 0) {
463c58f1c22SToby Isaac           *contains = PETSC_TRUE;
464c58f1c22SToby Isaac           break;
465c58f1c22SToby Isaac         }
466c58f1c22SToby Isaac       } else {
467c58f1c22SToby Isaac         PetscBool has;
468c58f1c22SToby Isaac 
469c58f1c22SToby Isaac         PetscHashIHasKey(label->ht[v], point, has);
470c58f1c22SToby Isaac         if (has) {
471c58f1c22SToby Isaac           *contains = PETSC_TRUE;
472c58f1c22SToby Isaac           break;
473c58f1c22SToby Isaac         }
474c58f1c22SToby Isaac       }
475c58f1c22SToby Isaac     }
476c58f1c22SToby Isaac   }
477c58f1c22SToby Isaac   PetscFunctionReturn(0);
478c58f1c22SToby Isaac }
479c58f1c22SToby Isaac 
4805aa44df4SToby Isaac /*
4815aa44df4SToby Isaac   DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
4825aa44df4SToby Isaac   When a label is created, it is initialized to -1.
4835aa44df4SToby Isaac 
4845aa44df4SToby Isaac   Input parameter:
4855aa44df4SToby Isaac . label - a DMLabel object
4865aa44df4SToby Isaac 
4875aa44df4SToby Isaac   Output parameter:
4885aa44df4SToby Isaac . defaultValue - the default value
4895aa44df4SToby Isaac 
4905aa44df4SToby Isaac   Level: beginner
4915aa44df4SToby Isaac 
4925aa44df4SToby Isaac .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
4935aa44df4SToby Isaac */
4945aa44df4SToby Isaac PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
4955aa44df4SToby Isaac {
4965aa44df4SToby Isaac   PetscFunctionBegin;
4975aa44df4SToby Isaac   *defaultValue = label->defaultValue;
4985aa44df4SToby Isaac   PetscFunctionReturn(0);
4995aa44df4SToby Isaac }
5005aa44df4SToby Isaac 
5015aa44df4SToby Isaac /*
5025aa44df4SToby Isaac   DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
5035aa44df4SToby Isaac   When a label is created, it is initialized to -1.
5045aa44df4SToby Isaac 
5055aa44df4SToby Isaac   Input parameter:
5065aa44df4SToby Isaac . label - a DMLabel object
5075aa44df4SToby Isaac 
5085aa44df4SToby Isaac   Output parameter:
5095aa44df4SToby Isaac . defaultValue - the default value
5105aa44df4SToby Isaac 
5115aa44df4SToby Isaac   Level: beginner
5125aa44df4SToby Isaac 
5135aa44df4SToby Isaac .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
5145aa44df4SToby Isaac */
5155aa44df4SToby Isaac PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
5165aa44df4SToby Isaac {
5175aa44df4SToby Isaac   PetscFunctionBegin;
5185aa44df4SToby Isaac   label->defaultValue = defaultValue;
5195aa44df4SToby Isaac   PetscFunctionReturn(0);
5205aa44df4SToby Isaac }
5215aa44df4SToby Isaac 
522c58f1c22SToby Isaac /*@
5235aa44df4SToby 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())
524c58f1c22SToby Isaac 
525c58f1c22SToby Isaac   Input Parameters:
526c58f1c22SToby Isaac + label - the DMLabel
527c58f1c22SToby Isaac - point - the point
528c58f1c22SToby Isaac 
529c58f1c22SToby Isaac   Output Parameter:
530c58f1c22SToby Isaac . value - The point value, or -1
531c58f1c22SToby Isaac 
532c58f1c22SToby Isaac   Level: intermediate
533c58f1c22SToby Isaac 
5345aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
535c58f1c22SToby Isaac @*/
536c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
537c58f1c22SToby Isaac {
538c58f1c22SToby Isaac   PetscInt       v;
539c58f1c22SToby Isaac   PetscErrorCode ierr;
540c58f1c22SToby Isaac 
541c58f1c22SToby Isaac   PetscFunctionBegin;
542c58f1c22SToby Isaac   PetscValidPointer(value, 3);
5435aa44df4SToby Isaac   *value = label->defaultValue;
544c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
545ad8374ffSToby Isaac     if (label->validIS[v]) {
546c58f1c22SToby Isaac       PetscInt i;
547c58f1c22SToby Isaac 
548a2d74346SToby Isaac       ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr);
549c58f1c22SToby Isaac       if (i >= 0) {
550c58f1c22SToby Isaac         *value = label->stratumValues[v];
551c58f1c22SToby Isaac         break;
552c58f1c22SToby Isaac       }
553c58f1c22SToby Isaac     } else {
554c58f1c22SToby Isaac       PetscBool has;
555c58f1c22SToby Isaac 
556c58f1c22SToby Isaac       PetscHashIHasKey(label->ht[v], point, has);
557c58f1c22SToby Isaac       if (has) {
558c58f1c22SToby Isaac         *value = label->stratumValues[v];
559c58f1c22SToby Isaac         break;
560c58f1c22SToby Isaac       }
561c58f1c22SToby Isaac     }
562c58f1c22SToby Isaac   }
563c58f1c22SToby Isaac   PetscFunctionReturn(0);
564c58f1c22SToby Isaac }
565c58f1c22SToby Isaac 
566c58f1c22SToby Isaac /*@
5675aa44df4SToby Isaac   DMLabelSetValue - Set the value a label assigns to a point.  If the value is the same as the label's default value (which is initially -1, and can be changed with DMLabelSetDefaultValue() to somethingg different), then this function will do nothing.
568c58f1c22SToby Isaac 
569c58f1c22SToby Isaac   Input Parameters:
570c58f1c22SToby Isaac + label - the DMLabel
571c58f1c22SToby Isaac . point - the point
572c58f1c22SToby Isaac - value - The point value
573c58f1c22SToby Isaac 
574c58f1c22SToby Isaac   Level: intermediate
575c58f1c22SToby Isaac 
5765aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
577c58f1c22SToby Isaac @*/
578c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
579c58f1c22SToby Isaac {
580c58f1c22SToby Isaac   PETSC_UNUSED PetscHashIIter iter, ret;
581c58f1c22SToby Isaac   PetscInt                    v;
582c58f1c22SToby Isaac   PetscErrorCode              ierr;
583c58f1c22SToby Isaac 
584c58f1c22SToby Isaac   PetscFunctionBegin;
585c58f1c22SToby Isaac   /* Find, or add, label value */
5865aa44df4SToby Isaac   if (value == label->defaultValue) PetscFunctionReturn(0);
587c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
588c58f1c22SToby Isaac     if (label->stratumValues[v] == value) break;
589c58f1c22SToby Isaac   }
590c58f1c22SToby Isaac   /* Create new table */
591c58f1c22SToby Isaac   if (v >= label->numStrata) {ierr = DMLabelAddStratum(label, value);CHKERRQ(ierr);}
592c58f1c22SToby Isaac   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
593c58f1c22SToby Isaac   /* Set key */
594c58f1c22SToby Isaac   PetscHashIPut(label->ht[v], point, ret, iter);
595c58f1c22SToby Isaac   PetscFunctionReturn(0);
596c58f1c22SToby Isaac }
597c58f1c22SToby Isaac 
598c58f1c22SToby Isaac /*@
599c58f1c22SToby Isaac   DMLabelClearValue - Clear the value a label assigns to a point
600c58f1c22SToby Isaac 
601c58f1c22SToby Isaac   Input Parameters:
602c58f1c22SToby Isaac + label - the DMLabel
603c58f1c22SToby Isaac . point - the point
604c58f1c22SToby Isaac - value - The point value
605c58f1c22SToby Isaac 
606c58f1c22SToby Isaac   Level: intermediate
607c58f1c22SToby Isaac 
608c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
609c58f1c22SToby Isaac @*/
610c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
611c58f1c22SToby Isaac {
612ad8374ffSToby Isaac   PetscInt       v;
613c58f1c22SToby Isaac   PetscErrorCode ierr;
614c58f1c22SToby Isaac 
615c58f1c22SToby Isaac   PetscFunctionBegin;
616c58f1c22SToby Isaac   /* Find label value */
617c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
618c58f1c22SToby Isaac     if (label->stratumValues[v] == value) break;
619c58f1c22SToby Isaac   }
620c58f1c22SToby Isaac   if (v >= label->numStrata) PetscFunctionReturn(0);
621eeed21e7SToby Isaac   if (label->validIS[v]) {
622eeed21e7SToby Isaac     ierr = DMLabelMakeInvalid_Private(label,v);CHKERRQ(ierr);
623eeed21e7SToby Isaac   }
624eeed21e7SToby Isaac   if (label->bt) {
625eeed21e7SToby 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);
626eeed21e7SToby Isaac     ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
627eeed21e7SToby Isaac   }
628c58f1c22SToby Isaac   ierr = PetscHashIDelKey(label->ht[v], point);CHKERRQ(ierr);
629c58f1c22SToby Isaac   PetscFunctionReturn(0);
630c58f1c22SToby Isaac }
631c58f1c22SToby Isaac 
632c58f1c22SToby Isaac /*@
633c58f1c22SToby Isaac   DMLabelInsertIS - Set all points in the IS to a value
634c58f1c22SToby Isaac 
635c58f1c22SToby Isaac   Input Parameters:
636c58f1c22SToby Isaac + label - the DMLabel
637c58f1c22SToby Isaac . is    - the point IS
638c58f1c22SToby Isaac - value - The point value
639c58f1c22SToby Isaac 
640c58f1c22SToby Isaac   Level: intermediate
641c58f1c22SToby Isaac 
642c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
643c58f1c22SToby Isaac @*/
644c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
645c58f1c22SToby Isaac {
646c58f1c22SToby Isaac   const PetscInt *points;
647c58f1c22SToby Isaac   PetscInt        n, p;
648c58f1c22SToby Isaac   PetscErrorCode  ierr;
649c58f1c22SToby Isaac 
650c58f1c22SToby Isaac   PetscFunctionBegin;
651c58f1c22SToby Isaac   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
652c58f1c22SToby Isaac   ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr);
653c58f1c22SToby Isaac   ierr = ISGetIndices(is, &points);CHKERRQ(ierr);
654c58f1c22SToby Isaac   for (p = 0; p < n; ++p) {ierr = DMLabelSetValue(label, points[p], value);CHKERRQ(ierr);}
655c58f1c22SToby Isaac   ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr);
656c58f1c22SToby Isaac   PetscFunctionReturn(0);
657c58f1c22SToby Isaac }
658c58f1c22SToby Isaac 
659c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
660c58f1c22SToby Isaac {
661c58f1c22SToby Isaac   PetscFunctionBegin;
662c58f1c22SToby Isaac   PetscValidPointer(numValues, 2);
663c58f1c22SToby Isaac   *numValues = label->numStrata;
664c58f1c22SToby Isaac   PetscFunctionReturn(0);
665c58f1c22SToby Isaac }
666c58f1c22SToby Isaac 
667c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
668c58f1c22SToby Isaac {
669c58f1c22SToby Isaac   PetscErrorCode ierr;
670c58f1c22SToby Isaac 
671c58f1c22SToby Isaac   PetscFunctionBegin;
672c58f1c22SToby Isaac   PetscValidPointer(values, 2);
673c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr);
674c58f1c22SToby Isaac   PetscFunctionReturn(0);
675c58f1c22SToby Isaac }
676c58f1c22SToby Isaac 
677fada774cSMatthew G. Knepley PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
678fada774cSMatthew G. Knepley {
679fada774cSMatthew G. Knepley   PetscInt v;
680fada774cSMatthew G. Knepley 
681fada774cSMatthew G. Knepley   PetscFunctionBegin;
682fada774cSMatthew G. Knepley   PetscValidPointer(exists, 3);
683fada774cSMatthew G. Knepley   *exists = PETSC_FALSE;
684fada774cSMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
685fada774cSMatthew G. Knepley     if (label->stratumValues[v] == value) {
686fada774cSMatthew G. Knepley       *exists = PETSC_TRUE;
687fada774cSMatthew G. Knepley       break;
688fada774cSMatthew G. Knepley     }
689fada774cSMatthew G. Knepley   }
690fada774cSMatthew G. Knepley   PetscFunctionReturn(0);
691fada774cSMatthew G. Knepley }
692fada774cSMatthew G. Knepley 
693c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
694c58f1c22SToby Isaac {
695c58f1c22SToby Isaac   PetscInt       v;
696c58f1c22SToby Isaac   PetscErrorCode ierr;
697c58f1c22SToby Isaac 
698c58f1c22SToby Isaac   PetscFunctionBegin;
699c58f1c22SToby Isaac   PetscValidPointer(size, 3);
700c58f1c22SToby Isaac   *size = 0;
701c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
702c58f1c22SToby Isaac     if (label->stratumValues[v] == value) {
703c58f1c22SToby Isaac       ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
704c58f1c22SToby Isaac       *size = label->stratumSizes[v];
705c58f1c22SToby Isaac       break;
706c58f1c22SToby Isaac     }
707c58f1c22SToby Isaac   }
708c58f1c22SToby Isaac   PetscFunctionReturn(0);
709c58f1c22SToby Isaac }
710c58f1c22SToby Isaac 
711c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
712c58f1c22SToby Isaac {
713c58f1c22SToby Isaac   PetscInt       v;
714c58f1c22SToby Isaac   PetscErrorCode ierr;
715c58f1c22SToby Isaac 
716c58f1c22SToby Isaac   PetscFunctionBegin;
717c58f1c22SToby Isaac   if (start) {PetscValidPointer(start, 3); *start = 0;}
718c58f1c22SToby Isaac   if (end)   {PetscValidPointer(end,   4); *end   = 0;}
719c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
720d6cb179aSToby Isaac     PetscInt min, max;
721ad8374ffSToby Isaac 
722c58f1c22SToby Isaac     if (label->stratumValues[v] != value) continue;
723c58f1c22SToby Isaac     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
724c58f1c22SToby Isaac     if (label->stratumSizes[v]  <= 0)     break;
725d6cb179aSToby Isaac     ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr);
726d6cb179aSToby Isaac     if (start) *start = min;
727d6cb179aSToby Isaac     if (end)   *end   = max+1;
728c58f1c22SToby Isaac     break;
729c58f1c22SToby Isaac   }
730c58f1c22SToby Isaac   PetscFunctionReturn(0);
731c58f1c22SToby Isaac }
732c58f1c22SToby Isaac 
733c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
734c58f1c22SToby Isaac {
735c58f1c22SToby Isaac   PetscInt       v;
736c58f1c22SToby Isaac   PetscErrorCode ierr;
737c58f1c22SToby Isaac 
738c58f1c22SToby Isaac   PetscFunctionBegin;
739c58f1c22SToby Isaac   PetscValidPointer(points, 3);
740c58f1c22SToby Isaac   *points = NULL;
741c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
742c58f1c22SToby Isaac     if (label->stratumValues[v] == value) {
743c58f1c22SToby Isaac       ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
744ad8374ffSToby Isaac       if (label->validIS[v]) {
745ad8374ffSToby Isaac         ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr);
746ad8374ffSToby Isaac         *points = label->points[v];
747c58f1c22SToby Isaac       } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to implement this to speedup Stratify");
748c58f1c22SToby Isaac       break;
749c58f1c22SToby Isaac     }
750c58f1c22SToby Isaac   }
751c58f1c22SToby Isaac   PetscFunctionReturn(0);
752c58f1c22SToby Isaac }
753c58f1c22SToby Isaac 
7544de306b1SToby Isaac PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
7554de306b1SToby Isaac {
7564de306b1SToby Isaac   PetscInt       v, numStrata;
7574de306b1SToby Isaac   PetscErrorCode ierr;
7584de306b1SToby Isaac 
7594de306b1SToby Isaac   PetscFunctionBegin;
7604de306b1SToby Isaac   numStrata = label->numStrata;
7614de306b1SToby Isaac   for (v = 0; v < numStrata; v++) {
7624de306b1SToby Isaac     if (label->stratumValues[v] == value) break;
7634de306b1SToby Isaac   }
7644de306b1SToby Isaac   if (v >= numStrata) {ierr = DMLabelAddStratum(label,value);CHKERRQ(ierr);}
7654de306b1SToby Isaac   if (is == label->points[v]) PetscFunctionReturn(0);
7664de306b1SToby Isaac   ierr = DMLabelClearStratum(label,value);CHKERRQ(ierr);
7674de306b1SToby Isaac   ierr = ISGetLocalSize(is,&(label->stratumSizes[v]));CHKERRQ(ierr);
7684de306b1SToby Isaac   label->stratumValues[v] = value;
7694de306b1SToby Isaac   label->validIS[v] = PETSC_TRUE;
7704de306b1SToby Isaac   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
7714de306b1SToby Isaac   ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
7724de306b1SToby Isaac   if (label->bt) {
7734de306b1SToby Isaac     const PetscInt *points;
7744de306b1SToby Isaac     PetscInt p;
7754de306b1SToby Isaac 
7764de306b1SToby Isaac     ierr = ISGetIndices(is,&points);CHKERRQ(ierr);
7774de306b1SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
7784de306b1SToby Isaac       const PetscInt point = points[p];
7794de306b1SToby Isaac 
7804de306b1SToby 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);
7814de306b1SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
7824de306b1SToby Isaac     }
7834de306b1SToby Isaac   }
7844de306b1SToby Isaac   label->points[v] = is;
7854de306b1SToby Isaac   PetscFunctionReturn(0);
7864de306b1SToby Isaac }
7874de306b1SToby Isaac 
7884de306b1SToby Isaac 
789c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
790c58f1c22SToby Isaac {
791c58f1c22SToby Isaac   PetscInt       v;
792c58f1c22SToby Isaac   PetscErrorCode ierr;
793c58f1c22SToby Isaac 
794c58f1c22SToby Isaac   PetscFunctionBegin;
795c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
796c58f1c22SToby Isaac     if (label->stratumValues[v] == value) break;
797c58f1c22SToby Isaac   }
798c58f1c22SToby Isaac   if (v >= label->numStrata) PetscFunctionReturn(0);
7994de306b1SToby Isaac   if (label->validIS[v]) {
8004de306b1SToby Isaac     if (label->bt) {
801c58f1c22SToby Isaac       PetscInt       i;
802ad8374ffSToby Isaac       const PetscInt *points;
803c58f1c22SToby Isaac 
804ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
805c58f1c22SToby Isaac       for (i = 0; i < label->stratumSizes[v]; ++i) {
806ad8374ffSToby Isaac         const PetscInt point = points[i];
807c58f1c22SToby Isaac 
808c58f1c22SToby 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);
809c58f1c22SToby Isaac         ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
810c58f1c22SToby Isaac       }
811ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
812c58f1c22SToby Isaac     }
813ad8374ffSToby Isaac     ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
814c58f1c22SToby Isaac     label->stratumSizes[v] = 0;
815ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
816ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
817c58f1c22SToby Isaac   } else {
818c58f1c22SToby Isaac     PetscHashIClear(label->ht[v]);
819c58f1c22SToby Isaac   }
820c58f1c22SToby Isaac   PetscFunctionReturn(0);
821c58f1c22SToby Isaac }
822c58f1c22SToby Isaac 
823c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
824c58f1c22SToby Isaac {
825c58f1c22SToby Isaac   PetscInt       v;
826c58f1c22SToby Isaac   PetscErrorCode ierr;
827c58f1c22SToby Isaac 
828c58f1c22SToby Isaac   PetscFunctionBegin;
829c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
830c58f1c22SToby Isaac   label->pStart = start;
831c58f1c22SToby Isaac   label->pEnd   = end;
832c58f1c22SToby Isaac   if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);}
833c58f1c22SToby Isaac   /* Could squish offsets, but would only make sense if I reallocate the storage */
834c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
835c58f1c22SToby Isaac     PetscInt off, q;
836ad8374ffSToby Isaac     const PetscInt *points;
837ad8374ffSToby Isaac     PetscInt *pointsNew = NULL;
838c58f1c22SToby Isaac 
839ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
840c58f1c22SToby Isaac     for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) {
841ad8374ffSToby Isaac       const PetscInt point = points[q];
842c58f1c22SToby Isaac 
843ad8374ffSToby Isaac       if ((point < start) || (point >= end)) {
844ad8374ffSToby Isaac         if (!pointsNew) {
845ad8374ffSToby Isaac           ierr = PetscMalloc1(label->stratumSizes[v],&pointsNew);CHKERRQ(ierr);
846ad8374ffSToby Isaac           ierr = PetscMemcpy(pointsNew,points,(size_t) off * sizeof(PetscInt));CHKERRQ(ierr);
847ad8374ffSToby Isaac         }
848ad8374ffSToby Isaac         continue;
849ad8374ffSToby Isaac       }
850ad8374ffSToby Isaac       if (pointsNew) {
851ad8374ffSToby Isaac         pointsNew[off++] = point;
852ad8374ffSToby Isaac       }
853ad8374ffSToby Isaac     }
854ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
855ad8374ffSToby Isaac     if (pointsNew) {
856ad8374ffSToby Isaac       ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
857ad8374ffSToby Isaac       ierr = ISCreateGeneral(PETSC_COMM_SELF,off,pointsNew,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
858ad8374ffSToby Isaac       ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
859c58f1c22SToby Isaac     }
860c58f1c22SToby Isaac     label->stratumSizes[v] = off;
861c58f1c22SToby Isaac   }
862c58f1c22SToby Isaac   ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr);
863c58f1c22SToby Isaac   PetscFunctionReturn(0);
864c58f1c22SToby Isaac }
865c58f1c22SToby Isaac 
866c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
867c58f1c22SToby Isaac {
868c58f1c22SToby Isaac   const PetscInt *perm;
869c58f1c22SToby Isaac   PetscInt        numValues, numPoints, v, q;
870c58f1c22SToby Isaac   PetscErrorCode  ierr;
871c58f1c22SToby Isaac 
872c58f1c22SToby Isaac   PetscFunctionBegin;
873c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
874c58f1c22SToby Isaac   ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr);
875c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr);
876c58f1c22SToby Isaac   ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr);
877c58f1c22SToby Isaac   ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr);
878c58f1c22SToby Isaac   for (v = 0; v < numValues; ++v) {
879c58f1c22SToby Isaac     const PetscInt size   = (*labelNew)->stratumSizes[v];
880ad8374ffSToby Isaac     const PetscInt *points;
881ad8374ffSToby Isaac     PetscInt *pointsNew;
882c58f1c22SToby Isaac 
883ad8374ffSToby Isaac     ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
884ad8374ffSToby Isaac     ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr);
885c58f1c22SToby Isaac     for (q = 0; q < size; ++q) {
886ad8374ffSToby Isaac       const PetscInt point = points[q];
887c58f1c22SToby Isaac 
888c58f1c22SToby 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);
889ad8374ffSToby Isaac       pointsNew[q] = perm[point];
890c58f1c22SToby Isaac     }
891ad8374ffSToby Isaac     ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
892ad8374ffSToby Isaac     ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr);
893ad8374ffSToby Isaac     ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr);
894ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr);
895ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr);
896c58f1c22SToby Isaac   }
897c58f1c22SToby Isaac   ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr);
898c58f1c22SToby Isaac   if (label->bt) {
899c58f1c22SToby Isaac     ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
900c58f1c22SToby Isaac     ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
901c58f1c22SToby Isaac   }
902c58f1c22SToby Isaac   PetscFunctionReturn(0);
903c58f1c22SToby Isaac }
904c58f1c22SToby Isaac 
90526c55118SMichael Lange PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
90626c55118SMichael Lange {
90726c55118SMichael Lange   MPI_Comm       comm;
90826c55118SMichael Lange   PetscInt       s, l, nroots, nleaves, dof, offset, size;
90926c55118SMichael Lange   PetscInt      *remoteOffsets, *rootStrata, *rootIdx;
91026c55118SMichael Lange   PetscSection   rootSection;
91126c55118SMichael Lange   PetscSF        labelSF;
91226c55118SMichael Lange   PetscErrorCode ierr;
91326c55118SMichael Lange 
91426c55118SMichael Lange   PetscFunctionBegin;
91526c55118SMichael Lange   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
91626c55118SMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
91726c55118SMichael Lange   /* Build a section of stratum values per point, generate the according SF
91826c55118SMichael Lange      and distribute point-wise stratum values to leaves. */
91926c55118SMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
92026c55118SMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
92126c55118SMichael Lange   ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr);
92226c55118SMichael Lange   if (label) {
92326c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
924ad8374ffSToby Isaac       const PetscInt *points;
925ad8374ffSToby Isaac 
926ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
92726c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
928ad8374ffSToby Isaac         ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr);
929ad8374ffSToby Isaac         ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr);
93026c55118SMichael Lange       }
931ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
93226c55118SMichael Lange     }
93326c55118SMichael Lange   }
93426c55118SMichael Lange   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
93526c55118SMichael Lange   /* Create a point-wise array of stratum values */
93626c55118SMichael Lange   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
93726c55118SMichael Lange   ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr);
93826c55118SMichael Lange   ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr);
93926c55118SMichael Lange   if (label) {
94026c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
941ad8374ffSToby Isaac       const PetscInt *points;
942ad8374ffSToby Isaac 
943ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
94426c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
945ad8374ffSToby Isaac         const PetscInt p = points[l];
94626c55118SMichael Lange         ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr);
94726c55118SMichael Lange         rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
94826c55118SMichael Lange       }
949ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
95026c55118SMichael Lange     }
95126c55118SMichael Lange   }
95226c55118SMichael Lange   /* Build SF that maps label points to remote processes */
95326c55118SMichael Lange   ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr);
95426c55118SMichael Lange   ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr);
95526c55118SMichael Lange   ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr);
95626c55118SMichael Lange   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
95726c55118SMichael Lange   /* Send the strata for each point over the derived SF */
95826c55118SMichael Lange   ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr);
95926c55118SMichael Lange   ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr);
96026c55118SMichael Lange   ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
96126c55118SMichael Lange   ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
96226c55118SMichael Lange   /* Clean up */
96326c55118SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
96426c55118SMichael Lange   ierr = PetscFree(rootIdx);CHKERRQ(ierr);
96526c55118SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
96626c55118SMichael Lange   ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr);
96726c55118SMichael Lange   PetscFunctionReturn(0);
96826c55118SMichael Lange }
96926c55118SMichael Lange 
970c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
971c58f1c22SToby Isaac {
972c58f1c22SToby Isaac   MPI_Comm       comm;
97326c55118SMichael Lange   PetscSection   leafSection;
97426c55118SMichael Lange   PetscInt       p, pStart, pEnd, s, size, dof, offset, stratum;
97526c55118SMichael Lange   PetscInt      *leafStrata, *strataIdx;
976ad8374ffSToby Isaac   PetscInt     **points;
977c58f1c22SToby Isaac   char          *name;
978c58f1c22SToby Isaac   PetscInt       nameSize;
9795cbdf6fcSMichael Lange   PetscHashI     stratumHash;
9805cbdf6fcSMichael Lange   PETSC_UNUSED   PetscHashIIter ret, iter;
981c58f1c22SToby Isaac   size_t         len = 0;
98226c55118SMichael Lange   PetscMPIInt    rank;
983c58f1c22SToby Isaac   PetscErrorCode ierr;
984c58f1c22SToby Isaac 
985c58f1c22SToby Isaac   PetscFunctionBegin;
986c58f1c22SToby Isaac   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
987c58f1c22SToby Isaac   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
988c58f1c22SToby Isaac   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
989c58f1c22SToby Isaac   /* Bcast name */
990c58f1c22SToby Isaac   if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);}
991c58f1c22SToby Isaac   nameSize = len;
992c58f1c22SToby Isaac   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
993c58f1c22SToby Isaac   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
994c58f1c22SToby Isaac   if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);}
995c58f1c22SToby Isaac   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
996c58f1c22SToby Isaac   ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr);
997c58f1c22SToby Isaac   ierr = PetscFree(name);CHKERRQ(ierr);
99877d236dfSMichael Lange   /* Bcast defaultValue */
99977d236dfSMichael Lange   if (!rank) (*labelNew)->defaultValue = label->defaultValue;
100077d236dfSMichael Lange   ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
100126c55118SMichael Lange   /* Distribute stratum values over the SF and get the point mapping on the receiver */
100226c55118SMichael Lange   ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr);
10035cbdf6fcSMichael Lange   /* Determine received stratum values and initialise new label*/
10045cbdf6fcSMichael Lange   PetscHashICreate(stratumHash);
100526c55118SMichael Lange   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
10065cbdf6fcSMichael Lange   for (p = 0; p < size; ++p) PetscHashIPut(stratumHash, leafStrata[p], ret, iter);
10075cbdf6fcSMichael Lange   PetscHashISize(stratumHash, (*labelNew)->numStrata);
1008ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr);
1009ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
10105cbdf6fcSMichael Lange   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr);
10115cbdf6fcSMichael Lange   /* Turn leafStrata into indices rather than stratum values */
10125cbdf6fcSMichael Lange   offset = 0;
10135cbdf6fcSMichael Lange   ierr = PetscHashIGetKeys(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr);
1014*a205f1fdSToby Isaac   ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr);
10155cbdf6fcSMichael Lange   for (p = 0; p < size; ++p) {
1016231b9e6fSMatthew G. Knepley     for (s = 0; s < (*labelNew)->numStrata; ++s) {
1017231b9e6fSMatthew G. Knepley       if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
10185cbdf6fcSMichael Lange     }
10195cbdf6fcSMichael Lange   }
1020c58f1c22SToby Isaac   /* Rebuild the point strata on the receiver */
1021c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr);
1022c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1023c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1024c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1025c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1026c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1027c58f1c22SToby Isaac       (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1028c58f1c22SToby Isaac     }
1029c58f1c22SToby Isaac   }
1030c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr);
1031c58f1c22SToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr);
1032ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr);
1033c58f1c22SToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1034c58f1c22SToby Isaac     PetscHashICreate((*labelNew)->ht[s]);
1035ad8374ffSToby Isaac     ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr);
1036c58f1c22SToby Isaac   }
1037c58f1c22SToby Isaac   /* Insert points into new strata */
1038c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr);
1039c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1040c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1041c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1042c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1043c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1044c58f1c22SToby Isaac       stratum = leafStrata[offset+s];
1045ad8374ffSToby Isaac       points[stratum][strataIdx[stratum]++] = p;
1046c58f1c22SToby Isaac     }
1047c58f1c22SToby Isaac   }
1048ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; s++) {
1049ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr);
1050ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr);
1051ad8374ffSToby Isaac   }
1052ad8374ffSToby Isaac   ierr = PetscFree(points);CHKERRQ(ierr);
10535cbdf6fcSMichael Lange   PetscHashIDestroy(stratumHash);
1054c58f1c22SToby Isaac   ierr = PetscFree(leafStrata);CHKERRQ(ierr);
1055c58f1c22SToby Isaac   ierr = PetscFree(strataIdx);CHKERRQ(ierr);
1056c58f1c22SToby Isaac   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
1057c58f1c22SToby Isaac   PetscFunctionReturn(0);
1058c58f1c22SToby Isaac }
1059c58f1c22SToby Isaac 
10607937d9ceSMichael Lange /*@
10617937d9ceSMichael Lange   DMLabelGather - Gather all label values from leafs into roots
10627937d9ceSMichael Lange 
10637937d9ceSMichael Lange   Input Parameters:
10647937d9ceSMichael Lange + label - the DMLabel
10657937d9ceSMichael Lange . point - the Star Forest point communication map
10667937d9ceSMichael Lange 
10677937d9ceSMichael Lange   Input Parameters:
10687937d9ceSMichael Lange + label - the new DMLabel with localised leaf values
10697937d9ceSMichael Lange 
10707937d9ceSMichael Lange   Level: developer
10717937d9ceSMichael Lange 
10727937d9ceSMichael Lange   Note: This is the inverse operation to DMLabelDistribute.
10737937d9ceSMichael Lange 
10747937d9ceSMichael Lange .seealso: DMLabelDistribute()
10757937d9ceSMichael Lange @*/
10767937d9ceSMichael Lange PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
10777937d9ceSMichael Lange {
10787937d9ceSMichael Lange   MPI_Comm       comm;
10797937d9ceSMichael Lange   PetscSection   rootSection;
10807937d9ceSMichael Lange   PetscSF        sfLabel;
10817937d9ceSMichael Lange   PetscSFNode   *rootPoints, *leafPoints;
10827937d9ceSMichael Lange   PetscInt       p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
10837937d9ceSMichael Lange   const PetscInt *rootDegree, *ilocal;
10847937d9ceSMichael Lange   PetscInt       *rootStrata;
10857937d9ceSMichael Lange   char          *name;
10867937d9ceSMichael Lange   PetscInt       nameSize;
10877937d9ceSMichael Lange   size_t         len = 0;
10889852e123SBarry Smith   PetscMPIInt    rank, size;
10897937d9ceSMichael Lange   PetscErrorCode ierr;
10907937d9ceSMichael Lange 
10917937d9ceSMichael Lange   PetscFunctionBegin;
10927937d9ceSMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
10937937d9ceSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
10949852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
10957937d9ceSMichael Lange   /* Bcast name */
10967937d9ceSMichael Lange   if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);}
10977937d9ceSMichael Lange   nameSize = len;
10987937d9ceSMichael Lange   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
10997937d9ceSMichael Lange   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
11007937d9ceSMichael Lange   if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);}
11017937d9ceSMichael Lange   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
11027937d9ceSMichael Lange   ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr);
11037937d9ceSMichael Lange   ierr = PetscFree(name);CHKERRQ(ierr);
11047937d9ceSMichael Lange   /* Gather rank/index pairs of leaves into local roots to build
11057937d9ceSMichael Lange      an inverse, multi-rooted SF. Note that this ignores local leaf
11067937d9ceSMichael Lange      indexing due to the use of the multiSF in PetscSFGather. */
11077937d9ceSMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr);
1108dc53bc9bSMatthew G. Knepley   ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr);
1109dc53bc9bSMatthew G. Knepley   for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
11107937d9ceSMichael Lange   for (p = 0; p < nleaves; p++) {
1111dc53bc9bSMatthew G. Knepley     leafPoints[ilocal[p]].index = ilocal[p];
1112dc53bc9bSMatthew G. Knepley     leafPoints[ilocal[p]].rank  = rank;
11137937d9ceSMichael Lange   }
11147937d9ceSMichael Lange   ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr);
11157937d9ceSMichael Lange   ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr);
11167937d9ceSMichael Lange   for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
11177937d9ceSMichael Lange   ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr);
11187937d9ceSMichael Lange   ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
11197937d9ceSMichael Lange   ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
11207937d9ceSMichael Lange   ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr);
11217937d9ceSMichael Lange   ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
11227937d9ceSMichael Lange   /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
11237937d9ceSMichael Lange   ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr);
11247937d9ceSMichael Lange   /* Rebuild the point strata on the receiver */
11257937d9ceSMichael Lange   for (p = 0, idx = 0; p < nroots; p++) {
11267937d9ceSMichael Lange     for (d = 0; d < rootDegree[p]; d++) {
11277937d9ceSMichael Lange       ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr);
11287937d9ceSMichael Lange       ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr);
11297937d9ceSMichael Lange       for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);}
11307937d9ceSMichael Lange     }
11317937d9ceSMichael Lange     idx += rootDegree[p];
11327937d9ceSMichael Lange   }
113377e0c0e7SMichael Lange   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
113477e0c0e7SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
113577e0c0e7SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
113677e0c0e7SMichael Lange   ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr);
11377937d9ceSMichael Lange   PetscFunctionReturn(0);
11387937d9ceSMichael Lange }
11397937d9ceSMichael Lange 
1140c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1141c58f1c22SToby Isaac {
1142c58f1c22SToby Isaac   IS              vIS;
1143c58f1c22SToby Isaac   const PetscInt *values;
1144c58f1c22SToby Isaac   PetscInt       *points;
1145c58f1c22SToby Isaac   PetscInt        nV, vS = 0, vE = 0, v, N;
1146c58f1c22SToby Isaac   PetscErrorCode  ierr;
1147c58f1c22SToby Isaac 
1148c58f1c22SToby Isaac   PetscFunctionBegin;
1149c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr);
1150c58f1c22SToby Isaac   ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
1151c58f1c22SToby Isaac   ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr);
1152c58f1c22SToby Isaac   if (nV) {vS = values[0]; vE = values[0]+1;}
1153c58f1c22SToby Isaac   for (v = 1; v < nV; ++v) {
1154c58f1c22SToby Isaac     vS = PetscMin(vS, values[v]);
1155c58f1c22SToby Isaac     vE = PetscMax(vE, values[v]+1);
1156c58f1c22SToby Isaac   }
1157c58f1c22SToby Isaac   ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr);
1158c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr);
1159c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1160c58f1c22SToby Isaac     PetscInt n;
1161c58f1c22SToby Isaac 
1162c58f1c22SToby Isaac     ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr);
1163c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr);
1164c58f1c22SToby Isaac   }
1165c58f1c22SToby Isaac   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
1166c58f1c22SToby Isaac   ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr);
1167c58f1c22SToby Isaac   ierr = PetscMalloc1(N, &points);CHKERRQ(ierr);
1168c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1169c58f1c22SToby Isaac     IS              is;
1170c58f1c22SToby Isaac     const PetscInt *spoints;
1171c58f1c22SToby Isaac     PetscInt        dof, off, p;
1172c58f1c22SToby Isaac 
1173c58f1c22SToby Isaac     ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr);
1174c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr);
1175c58f1c22SToby Isaac     ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr);
1176c58f1c22SToby Isaac     ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr);
1177c58f1c22SToby Isaac     for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1178c58f1c22SToby Isaac     ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr);
1179c58f1c22SToby Isaac     ierr = ISDestroy(&is);CHKERRQ(ierr);
1180c58f1c22SToby Isaac   }
1181c58f1c22SToby Isaac   ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr);
1182c58f1c22SToby Isaac   ierr = ISDestroy(&vIS);CHKERRQ(ierr);
1183c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1184c58f1c22SToby Isaac   PetscFunctionReturn(0);
1185c58f1c22SToby Isaac }
1186c58f1c22SToby Isaac 
1187c58f1c22SToby Isaac /*@C
1188c58f1c22SToby Isaac   PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1189c58f1c22SToby Isaac   the local section and an SF describing the section point overlap.
1190c58f1c22SToby Isaac 
1191c58f1c22SToby Isaac   Input Parameters:
1192c58f1c22SToby Isaac   + s - The PetscSection for the local field layout
1193c58f1c22SToby Isaac   . sf - The SF describing parallel layout of the section points
1194c58f1c22SToby Isaac   . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1195c58f1c22SToby Isaac   . label - The label specifying the points
1196c58f1c22SToby Isaac   - labelValue - The label stratum specifying the points
1197c58f1c22SToby Isaac 
1198c58f1c22SToby Isaac   Output Parameter:
1199c58f1c22SToby Isaac   . gsection - The PetscSection for the global field layout
1200c58f1c22SToby Isaac 
1201c58f1c22SToby Isaac   Note: This gives negative sizes and offsets to points not owned by this process
1202c58f1c22SToby Isaac 
1203c58f1c22SToby Isaac   Level: developer
1204c58f1c22SToby Isaac 
1205c58f1c22SToby Isaac .seealso: PetscSectionCreate()
1206c58f1c22SToby Isaac @*/
1207c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1208c58f1c22SToby Isaac {
1209c58f1c22SToby Isaac   PetscInt      *neg = NULL, *tmpOff = NULL;
1210c58f1c22SToby Isaac   PetscInt       pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1211c58f1c22SToby Isaac   PetscErrorCode ierr;
1212c58f1c22SToby Isaac 
1213c58f1c22SToby Isaac   PetscFunctionBegin;
1214c58f1c22SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr);
1215c58f1c22SToby Isaac   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1216c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr);
1217c58f1c22SToby Isaac   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1218c58f1c22SToby Isaac   if (nroots >= 0) {
1219c58f1c22SToby Isaac     if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1220c58f1c22SToby Isaac     ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr);
1221c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1222c58f1c22SToby Isaac       ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr);
1223c58f1c22SToby Isaac     } else {
1224c58f1c22SToby Isaac       tmpOff = &(*gsection)->atlasDof[-pStart];
1225c58f1c22SToby Isaac     }
1226c58f1c22SToby Isaac   }
1227c58f1c22SToby Isaac   /* Mark ghost points with negative dof */
1228c58f1c22SToby Isaac   for (p = pStart; p < pEnd; ++p) {
1229c58f1c22SToby Isaac     PetscInt value;
1230c58f1c22SToby Isaac 
1231c58f1c22SToby Isaac     ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr);
1232c58f1c22SToby Isaac     if (value != labelValue) continue;
1233c58f1c22SToby Isaac     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1234c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr);
1235c58f1c22SToby Isaac     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1236c58f1c22SToby Isaac     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);}
1237c58f1c22SToby Isaac     if (neg) neg[p] = -(dof+1);
1238c58f1c22SToby Isaac   }
1239c58f1c22SToby Isaac   ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr);
1240c58f1c22SToby Isaac   if (nroots >= 0) {
1241c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1242c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1243c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1244c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1245c58f1c22SToby Isaac     }
1246c58f1c22SToby Isaac   }
1247c58f1c22SToby Isaac   /* Calculate new sizes, get proccess offset, and calculate point offsets */
1248c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1249c58f1c22SToby Isaac     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1250c58f1c22SToby Isaac     (*gsection)->atlasOff[p] = off;
1251c58f1c22SToby Isaac     off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1252c58f1c22SToby Isaac   }
1253c58f1c22SToby Isaac   ierr       = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1254c58f1c22SToby Isaac   globalOff -= off;
1255c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1256c58f1c22SToby Isaac     (*gsection)->atlasOff[p] += globalOff;
1257c58f1c22SToby Isaac     if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1258c58f1c22SToby Isaac   }
1259c58f1c22SToby Isaac   /* Put in negative offsets for ghost points */
1260c58f1c22SToby Isaac   if (nroots >= 0) {
1261c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1262c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1263c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1264c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1265c58f1c22SToby Isaac     }
1266c58f1c22SToby Isaac   }
1267c58f1c22SToby Isaac   if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);}
1268c58f1c22SToby Isaac   ierr = PetscFree(neg);CHKERRQ(ierr);
1269c58f1c22SToby Isaac   PetscFunctionReturn(0);
1270c58f1c22SToby Isaac }
1271c58f1c22SToby Isaac 
12725fdea053SToby Isaac typedef struct _n_PetscSectionSym_Label
12735fdea053SToby Isaac {
12745fdea053SToby Isaac   DMLabel           label;
12755fdea053SToby Isaac   PetscCopyMode     *modes;
12765fdea053SToby Isaac   PetscInt          *sizes;
12775fdea053SToby Isaac   const PetscInt    ***perms;
12785fdea053SToby Isaac   const PetscScalar ***rots;
12795fdea053SToby Isaac   PetscInt          (*minMaxOrients)[2];
12805fdea053SToby Isaac   PetscInt          numStrata; /* numStrata is only increasing, functions as a state */
12815fdea053SToby Isaac } PetscSectionSym_Label;
12825fdea053SToby Isaac 
12835fdea053SToby Isaac static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
12845fdea053SToby Isaac {
12855fdea053SToby Isaac   PetscInt              i, j;
12865fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
12875fdea053SToby Isaac   PetscErrorCode        ierr;
12885fdea053SToby Isaac 
12895fdea053SToby Isaac   PetscFunctionBegin;
12905fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
12915fdea053SToby Isaac     if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
12925fdea053SToby Isaac       for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
12935fdea053SToby Isaac         if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);}
12945fdea053SToby Isaac         if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);}
12955fdea053SToby Isaac       }
12965fdea053SToby Isaac       if (sl->perms[i]) {
12975fdea053SToby Isaac         const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
12985fdea053SToby Isaac 
12995fdea053SToby Isaac         ierr = PetscFree(perms);CHKERRQ(ierr);
13005fdea053SToby Isaac       }
13015fdea053SToby Isaac       if (sl->rots[i]) {
13025fdea053SToby Isaac         const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
13035fdea053SToby Isaac 
13045fdea053SToby Isaac         ierr = PetscFree(rots);CHKERRQ(ierr);
13055fdea053SToby Isaac       }
13065fdea053SToby Isaac     }
13075fdea053SToby Isaac   }
13085fdea053SToby Isaac   ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr);
13095fdea053SToby Isaac   ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr);
13105fdea053SToby Isaac   sl->numStrata = 0;
13115fdea053SToby Isaac   PetscFunctionReturn(0);
13125fdea053SToby Isaac }
13135fdea053SToby Isaac 
13145fdea053SToby Isaac static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
13155fdea053SToby Isaac {
13165fdea053SToby Isaac   PetscErrorCode ierr;
13175fdea053SToby Isaac 
13185fdea053SToby Isaac   PetscFunctionBegin;
13195fdea053SToby Isaac   ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);
13205fdea053SToby Isaac   ierr = PetscFree(sym->data);CHKERRQ(ierr);
13215fdea053SToby Isaac   PetscFunctionReturn(0);
13225fdea053SToby Isaac }
13235fdea053SToby Isaac 
13245fdea053SToby Isaac static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
13255fdea053SToby Isaac {
13265fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
13275fdea053SToby Isaac   PetscBool             isAscii;
13285fdea053SToby Isaac   DMLabel               label = sl->label;
13295fdea053SToby Isaac   PetscErrorCode        ierr;
13305fdea053SToby Isaac 
13315fdea053SToby Isaac   PetscFunctionBegin;
13325fdea053SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
13335fdea053SToby Isaac   if (isAscii) {
13345fdea053SToby Isaac     PetscInt          i, j, k;
13355fdea053SToby Isaac     PetscViewerFormat format;
13365fdea053SToby Isaac 
13375fdea053SToby Isaac     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
13385fdea053SToby Isaac     if (label) {
13395fdea053SToby Isaac       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
13405fdea053SToby Isaac       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
13415fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
13425fdea053SToby Isaac         ierr = DMLabelView(label, viewer);CHKERRQ(ierr);
13435fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
13445fdea053SToby Isaac       } else {
13455fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer,"  Label '%s'\n",sl->label->name);CHKERRQ(ierr);
13465fdea053SToby Isaac       }
13475fdea053SToby Isaac     } else {
13485fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr);
13495fdea053SToby Isaac     }
13505fdea053SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
13515fdea053SToby Isaac     for (i = 0; i <= sl->numStrata; i++) {
13525fdea053SToby Isaac       PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
13535fdea053SToby Isaac 
13545fdea053SToby Isaac       if (!(sl->perms[i] || sl->rots[i])) {
13555fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr);
13565fdea053SToby Isaac       } else {
13575fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr);
13585fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
13595fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr);
13605fdea053SToby Isaac         if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
13615fdea053SToby Isaac           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
13625fdea053SToby Isaac           for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
13635fdea053SToby Isaac             if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
13645fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr);
13655fdea053SToby Isaac             } else {
13665fdea053SToby Isaac               PetscInt tab;
13675fdea053SToby Isaac 
13685fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr);
13695fdea053SToby Isaac               ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
13705fdea053SToby Isaac               ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr);
13715fdea053SToby Isaac               if (sl->perms[i] && sl->perms[i][j]) {
13725fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr);
13735fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
13745fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);}
13755fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
13765fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
13775fdea053SToby Isaac               }
13785fdea053SToby Isaac               if (sl->rots[i] && sl->rots[i][j]) {
13795fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Rotations:  ");CHKERRQ(ierr);
13805fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
13815fdea053SToby Isaac #if defined(PETSC_USE_COMPLEX)
13825fdea053SToby 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);}
13835fdea053SToby Isaac #else
13845fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);}
13855fdea053SToby Isaac #endif
13865fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
13875fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
13885fdea053SToby Isaac               }
13895fdea053SToby Isaac               ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
13905fdea053SToby Isaac             }
13915fdea053SToby Isaac           }
13925fdea053SToby Isaac           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
13935fdea053SToby Isaac         }
13945fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
13955fdea053SToby Isaac       }
13965fdea053SToby Isaac     }
13975fdea053SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
13985fdea053SToby Isaac   }
13995fdea053SToby Isaac   PetscFunctionReturn(0);
14005fdea053SToby Isaac }
14015fdea053SToby Isaac 
14025fdea053SToby Isaac /*@
14035fdea053SToby Isaac   PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
14045fdea053SToby Isaac 
14055fdea053SToby Isaac   Logically collective on sym
14065fdea053SToby Isaac 
14075fdea053SToby Isaac   Input parameters:
14085fdea053SToby Isaac + sym - the section symmetries
14095fdea053SToby Isaac - label - the DMLabel describing the types of points
14105fdea053SToby Isaac 
14115fdea053SToby Isaac   Level: developer:
14125fdea053SToby Isaac 
14135fdea053SToby Isaac .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()
14145fdea053SToby Isaac @*/
14155fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
14165fdea053SToby Isaac {
14175fdea053SToby Isaac   PetscSectionSym_Label *sl;
14185fdea053SToby Isaac   PetscErrorCode ierr;
14195fdea053SToby Isaac 
14205fdea053SToby Isaac   PetscFunctionBegin;
14215fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
14225fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
14235fdea053SToby Isaac   if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);}
14245fdea053SToby Isaac   if (label) {
14255fdea053SToby Isaac     label->refct++;
14265fdea053SToby Isaac     sl->label = label;
14275fdea053SToby Isaac     ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr);
14281a834cf9SToby 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);
14291a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr);
14301a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr);
14311a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr);
14321a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr);
14331a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr);
14345fdea053SToby Isaac   }
14355fdea053SToby Isaac   PetscFunctionReturn(0);
14365fdea053SToby Isaac }
14375fdea053SToby Isaac 
14385fdea053SToby Isaac /*@C
14395fdea053SToby Isaac   PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
14405fdea053SToby Isaac 
14415fdea053SToby Isaac   Logically collective on PetscSectionSym
14425fdea053SToby Isaac 
14435fdea053SToby Isaac   InputParameters:
14445fdea053SToby Isaac + sys       - the section symmetries
14455fdea053SToby Isaac . stratum   - the stratum value in the label that we are assigning symmetries for
14465fdea053SToby Isaac . size      - the number of dofs for points in the stratum of the label
14475fdea053SToby Isaac . minOrient - the smallest orientation for a point in this stratum
14485fdea053SToby Isaac . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
14495fdea053SToby Isaac . mode      - how sym should copy the perms and rots arrays
14505fdea053SToby Isaac . perms     - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation.  A NULL permutation is the identity
14515fdea053SToby 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
14525fdea053SToby Isaac 
14535fdea053SToby Isaac   Level: developer
14545fdea053SToby Isaac 
14555fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()
14565fdea053SToby Isaac @*/
14575fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
14585fdea053SToby Isaac {
14595fdea053SToby Isaac   PetscInt       i, j, k;
14605fdea053SToby Isaac   PetscSectionSym_Label *sl;
14615fdea053SToby Isaac   PetscErrorCode ierr;
14625fdea053SToby Isaac 
14635fdea053SToby Isaac   PetscFunctionBegin;
14645fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
14655fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
14665fdea053SToby Isaac   if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet");
14675fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
14685fdea053SToby Isaac     PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
14695fdea053SToby Isaac 
14705fdea053SToby Isaac     if (stratum == value) break;
14715fdea053SToby Isaac   }
14725fdea053SToby Isaac   if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,sl->label->name);
14735fdea053SToby Isaac   sl->sizes[i] = size;
14745fdea053SToby Isaac   sl->modes[i] = mode;
14755fdea053SToby Isaac   sl->minMaxOrients[i][0] = minOrient;
14765fdea053SToby Isaac   sl->minMaxOrients[i][1] = maxOrient;
14775fdea053SToby Isaac   if (mode == PETSC_COPY_VALUES) {
14785fdea053SToby Isaac     if (perms) {
14795fdea053SToby Isaac       PetscInt    **ownPerms;
14805fdea053SToby Isaac 
14815fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr);
14825fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
14835fdea053SToby Isaac         if (perms[j]) {
14845fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr);
14855fdea053SToby Isaac           for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];}
14865fdea053SToby Isaac         }
14875fdea053SToby Isaac       }
14885fdea053SToby Isaac       sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient];
14895fdea053SToby Isaac     }
14905fdea053SToby Isaac     if (rots) {
14915fdea053SToby Isaac       PetscScalar **ownRots;
14925fdea053SToby Isaac 
14935fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr);
14945fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
14955fdea053SToby Isaac         if (rots[j]) {
14965fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr);
14975fdea053SToby Isaac           for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];}
14985fdea053SToby Isaac         }
14995fdea053SToby Isaac       }
15005fdea053SToby Isaac       sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient];
15015fdea053SToby Isaac     }
15025fdea053SToby Isaac   } else {
15035fdea053SToby Isaac     sl->perms[i] = perms ? &perms[-minOrient] : NULL;
15045fdea053SToby Isaac     sl->rots[i]  = rots ? &rots[-minOrient] : NULL;
15055fdea053SToby Isaac   }
15065fdea053SToby Isaac   PetscFunctionReturn(0);
15075fdea053SToby Isaac }
15085fdea053SToby Isaac 
15095fdea053SToby Isaac static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
15105fdea053SToby Isaac {
15115fdea053SToby Isaac   PetscInt              i, j, numStrata;
15125fdea053SToby Isaac   PetscSectionSym_Label *sl;
15135fdea053SToby Isaac   DMLabel               label;
15145fdea053SToby Isaac   PetscErrorCode        ierr;
15155fdea053SToby Isaac 
15165fdea053SToby Isaac   PetscFunctionBegin;
15175fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
15185fdea053SToby Isaac   numStrata = sl->numStrata;
15195fdea053SToby Isaac   label     = sl->label;
15205fdea053SToby Isaac   for (i = 0; i < numPoints; i++) {
15215fdea053SToby Isaac     PetscInt point = points[2*i];
15225fdea053SToby Isaac     PetscInt ornt  = points[2*i+1];
15235fdea053SToby Isaac 
15245fdea053SToby Isaac     for (j = 0; j < numStrata; j++) {
15255fdea053SToby Isaac       if (label->validIS[j]) {
15265fdea053SToby Isaac         PetscInt k;
15275fdea053SToby Isaac 
15285fdea053SToby Isaac         ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr);
15295fdea053SToby Isaac         if (k >= 0) break;
15305fdea053SToby Isaac       } else {
15315fdea053SToby Isaac         PetscBool has;
15325fdea053SToby Isaac 
15335fdea053SToby Isaac         PetscHashIHasKey(label->ht[j], point, has);
15345fdea053SToby Isaac         if (has) break;
15355fdea053SToby Isaac       }
15365fdea053SToby Isaac     }
15375fdea053SToby 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);
15385fdea053SToby Isaac     if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;}
15395fdea053SToby Isaac     if (rots) {rots[i]  = sl->rots[j] ? sl->rots[j][ornt] : NULL;}
15405fdea053SToby Isaac   }
15415fdea053SToby Isaac   PetscFunctionReturn(0);
15425fdea053SToby Isaac }
15435fdea053SToby Isaac 
15445fdea053SToby Isaac PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
15455fdea053SToby Isaac {
15465fdea053SToby Isaac   PetscSectionSym_Label *sl;
15475fdea053SToby Isaac   PetscErrorCode        ierr;
15485fdea053SToby Isaac 
15495fdea053SToby Isaac   PetscFunctionBegin;
15505fdea053SToby Isaac   ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr);
15515fdea053SToby Isaac   sym->ops->getpoints = PetscSectionSymGetPoints_Label;
15525fdea053SToby Isaac   sym->ops->view      = PetscSectionSymView_Label;
15535fdea053SToby Isaac   sym->ops->destroy   = PetscSectionSymDestroy_Label;
15545fdea053SToby Isaac   sym->data           = (void *) sl;
15555fdea053SToby Isaac   PetscFunctionReturn(0);
15565fdea053SToby Isaac }
15575fdea053SToby Isaac 
15585fdea053SToby Isaac /*@
15595fdea053SToby Isaac   PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
15605fdea053SToby Isaac 
15615fdea053SToby Isaac   Collective
15625fdea053SToby Isaac 
15635fdea053SToby Isaac   Input Parameters:
15645fdea053SToby Isaac + comm - the MPI communicator for the new symmetry
15655fdea053SToby Isaac - label - the label defining the strata
15665fdea053SToby Isaac 
15675fdea053SToby Isaac   Output Parameters:
15685fdea053SToby Isaac . sym - the section symmetries
15695fdea053SToby Isaac 
15705fdea053SToby Isaac   Level: developer
15715fdea053SToby Isaac 
15725fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()
15735fdea053SToby Isaac @*/
15745fdea053SToby Isaac PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
15755fdea053SToby Isaac {
15765fdea053SToby Isaac   PetscErrorCode        ierr;
15775fdea053SToby Isaac 
15785fdea053SToby Isaac   PetscFunctionBegin;
15795fdea053SToby Isaac   ierr = DMInitializePackage();CHKERRQ(ierr);
15805fdea053SToby Isaac   ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr);
15815fdea053SToby Isaac   ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr);
15825fdea053SToby Isaac   ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr);
15835fdea053SToby Isaac   PetscFunctionReturn(0);
15845fdea053SToby Isaac }
1585