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