xref: /petsc/src/dm/label/dmlabel.c (revision 1a834cf98066b46fe06eb0a17d5d9fc7108bb4a0)
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 #undef __FUNCT__
7c58f1c22SToby Isaac #define __FUNCT__ "DMLabelCreate"
8c58f1c22SToby Isaac /*@C
9c58f1c22SToby Isaac   DMLabelCreate - Create a DMLabel object, which is a multimap
10c58f1c22SToby Isaac 
11c58f1c22SToby Isaac   Input parameter:
12c58f1c22SToby Isaac . name - The label name
13c58f1c22SToby Isaac 
14c58f1c22SToby Isaac   Output parameter:
15c58f1c22SToby Isaac . label - The DMLabel
16c58f1c22SToby Isaac 
17c58f1c22SToby Isaac   Level: beginner
18c58f1c22SToby Isaac 
19c58f1c22SToby Isaac .seealso: DMLabelDestroy()
20c58f1c22SToby Isaac @*/
21c58f1c22SToby Isaac PetscErrorCode DMLabelCreate(const char name[], DMLabel *label)
22c58f1c22SToby Isaac {
23c58f1c22SToby Isaac   PetscErrorCode ierr;
24c58f1c22SToby Isaac 
25c58f1c22SToby Isaac   PetscFunctionBegin;
26c58f1c22SToby Isaac   ierr = PetscNew(label);CHKERRQ(ierr);
27c58f1c22SToby Isaac   ierr = PetscStrallocpy(name, &(*label)->name);CHKERRQ(ierr);
28c58f1c22SToby Isaac 
29c58f1c22SToby Isaac   (*label)->refct          = 1;
30c58f1c22SToby Isaac   (*label)->state          = -1;
31c58f1c22SToby Isaac   (*label)->numStrata      = 0;
325aa44df4SToby Isaac   (*label)->defaultValue   = -1;
33c58f1c22SToby Isaac   (*label)->stratumValues  = NULL;
34ad8374ffSToby Isaac   (*label)->validIS        = NULL;
35c58f1c22SToby Isaac   (*label)->stratumSizes   = NULL;
36c58f1c22SToby Isaac   (*label)->points         = NULL;
37c58f1c22SToby Isaac   (*label)->ht             = NULL;
38c58f1c22SToby Isaac   (*label)->pStart         = -1;
39c58f1c22SToby Isaac   (*label)->pEnd           = -1;
40c58f1c22SToby Isaac   (*label)->bt             = NULL;
41c58f1c22SToby Isaac   PetscFunctionReturn(0);
42c58f1c22SToby Isaac }
43c58f1c22SToby Isaac 
44c58f1c22SToby Isaac #undef __FUNCT__
45c58f1c22SToby Isaac #define __FUNCT__ "DMLabelMakeValid_Private"
46c58f1c22SToby Isaac /*
47c58f1c22SToby Isaac   DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format
48c58f1c22SToby Isaac 
49c58f1c22SToby Isaac   Input parameter:
50c58f1c22SToby Isaac + label - The DMLabel
51c58f1c22SToby Isaac - v - The stratum value
52c58f1c22SToby Isaac 
53c58f1c22SToby Isaac   Output parameter:
54c58f1c22SToby Isaac . label - The DMLabel with stratum in sorted list format
55c58f1c22SToby Isaac 
56c58f1c22SToby Isaac   Level: developer
57c58f1c22SToby Isaac 
58c58f1c22SToby Isaac .seealso: DMLabelCreate()
59c58f1c22SToby Isaac */
60c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v)
61c58f1c22SToby Isaac {
62c58f1c22SToby Isaac   PetscInt       off;
63ad8374ffSToby Isaac   PetscInt       *pointArray;
64c58f1c22SToby Isaac   PetscErrorCode ierr;
65c58f1c22SToby Isaac 
66ad8374ffSToby Isaac   if (label->validIS[v]) return 0;
67c58f1c22SToby Isaac   if (v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeValid_Private\n", v);
68c58f1c22SToby Isaac   PetscFunctionBegin;
69c58f1c22SToby Isaac   PetscHashISize(label->ht[v], label->stratumSizes[v]);
70c58f1c22SToby Isaac 
71ad8374ffSToby Isaac   ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr);
72c58f1c22SToby Isaac   off  = 0;
73ad8374ffSToby Isaac   ierr = PetscHashIGetKeys(label->ht[v], &off, pointArray);CHKERRQ(ierr);
74c58f1c22SToby 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]);
75c58f1c22SToby Isaac   PetscHashIClear(label->ht[v]);
76ad8374ffSToby Isaac   ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr);
77c58f1c22SToby Isaac   if (label->bt) {
78c58f1c22SToby Isaac     PetscInt p;
79c58f1c22SToby Isaac 
80c58f1c22SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
81ad8374ffSToby Isaac       const PetscInt point = pointArray[p];
82c58f1c22SToby Isaac 
83c58f1c22SToby 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);
84c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
85c58f1c22SToby Isaac     }
86c58f1c22SToby Isaac   }
87ad8374ffSToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF,label->stratumSizes[v],pointArray,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
88ad8374ffSToby Isaac   ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
89ad8374ffSToby Isaac   label->validIS[v] = PETSC_TRUE;
90c58f1c22SToby Isaac   ++label->state;
91c58f1c22SToby Isaac   PetscFunctionReturn(0);
92c58f1c22SToby Isaac }
93c58f1c22SToby Isaac 
94c58f1c22SToby Isaac #undef __FUNCT__
95c58f1c22SToby Isaac #define __FUNCT__ "DMLabelMakeAllValid_Private"
96c58f1c22SToby Isaac /*
97c58f1c22SToby Isaac   DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format
98c58f1c22SToby Isaac 
99c58f1c22SToby Isaac   Input parameter:
100c58f1c22SToby Isaac . label - The DMLabel
101c58f1c22SToby Isaac 
102c58f1c22SToby Isaac   Output parameter:
103c58f1c22SToby Isaac . label - The DMLabel with all strata in sorted list format
104c58f1c22SToby Isaac 
105c58f1c22SToby Isaac   Level: developer
106c58f1c22SToby Isaac 
107c58f1c22SToby Isaac .seealso: DMLabelCreate()
108c58f1c22SToby Isaac */
109c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label)
110c58f1c22SToby Isaac {
111c58f1c22SToby Isaac   PetscInt       v;
112c58f1c22SToby Isaac   PetscErrorCode ierr;
113c58f1c22SToby Isaac 
114c58f1c22SToby Isaac   PetscFunctionBegin;
115c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; v++){
116c58f1c22SToby Isaac     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
117c58f1c22SToby Isaac   }
118c58f1c22SToby Isaac   PetscFunctionReturn(0);
119c58f1c22SToby Isaac }
120c58f1c22SToby Isaac 
121c58f1c22SToby Isaac #undef __FUNCT__
122c58f1c22SToby Isaac #define __FUNCT__ "DMLabelMakeInvalid_Private"
123c58f1c22SToby Isaac /*
124c58f1c22SToby Isaac   DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format
125c58f1c22SToby Isaac 
126c58f1c22SToby Isaac   Input parameter:
127c58f1c22SToby Isaac + label - The DMLabel
128c58f1c22SToby Isaac - v - The stratum value
129c58f1c22SToby Isaac 
130c58f1c22SToby Isaac   Output parameter:
131c58f1c22SToby Isaac . label - The DMLabel with stratum in hash format
132c58f1c22SToby Isaac 
133c58f1c22SToby Isaac   Level: developer
134c58f1c22SToby Isaac 
135c58f1c22SToby Isaac .seealso: DMLabelCreate()
136c58f1c22SToby Isaac */
137c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v)
138c58f1c22SToby Isaac {
139c58f1c22SToby Isaac   PETSC_UNUSED PetscHashIIter ret, iter;
140c58f1c22SToby Isaac   PetscInt                    p;
141ad8374ffSToby Isaac   const PetscInt              *points;
142c58f1c22SToby Isaac   PetscErrorCode ierr;
143c58f1c22SToby Isaac 
144c58f1c22SToby Isaac   PetscFunctionBegin;
145ad8374ffSToby Isaac   if (!label->validIS[v]) PetscFunctionReturn(0);
146ad8374ffSToby Isaac   if (label->points[v]) {
147ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
148ad8374ffSToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) PetscHashIPut(label->ht[v], points[p], ret, iter);
149ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
150ad8374ffSToby Isaac     ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
151ad8374ffSToby Isaac   }
152ad8374ffSToby Isaac   label->validIS[v] = PETSC_FALSE;
153c58f1c22SToby Isaac   PetscFunctionReturn(0);
154c58f1c22SToby Isaac }
155c58f1c22SToby Isaac 
156c58f1c22SToby Isaac #undef __FUNCT__
157c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetState"
158c58f1c22SToby Isaac PetscErrorCode DMLabelGetState(DMLabel label, PetscObjectState *state)
159c58f1c22SToby Isaac {
160c58f1c22SToby Isaac   PetscFunctionBegin;
161c58f1c22SToby Isaac   PetscValidPointer(state, 2);
162c58f1c22SToby Isaac   *state = label->state;
163c58f1c22SToby Isaac   PetscFunctionReturn(0);
164c58f1c22SToby Isaac }
165c58f1c22SToby Isaac 
166c58f1c22SToby Isaac #undef __FUNCT__
167c58f1c22SToby Isaac #define __FUNCT__ "DMLabelAddStratum"
168c58f1c22SToby Isaac PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
169c58f1c22SToby Isaac {
170ad8374ffSToby Isaac   PetscInt    v, *tmpV, *tmpS;
171ad8374ffSToby Isaac   IS         *tmpP;
172c58f1c22SToby Isaac   PetscHashI *tmpH;
173c58f1c22SToby Isaac   PetscBool  *tmpB;
174c58f1c22SToby Isaac   PetscErrorCode ierr;
175c58f1c22SToby Isaac 
176c58f1c22SToby Isaac   PetscFunctionBegin;
177c58f1c22SToby Isaac 
178ad8374ffSToby Isaac   for (v = 0; v < label->numStrata; v++) {
179ad8374ffSToby Isaac     if (label->stratumValues[v] == value) PetscFunctionReturn(0);
180ad8374ffSToby Isaac   }
181c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpV);CHKERRQ(ierr);
182c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpS);CHKERRQ(ierr);
183c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpH);CHKERRQ(ierr);
184c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpP);CHKERRQ(ierr);
185c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpB);CHKERRQ(ierr);
186c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
187c58f1c22SToby Isaac     tmpV[v] = label->stratumValues[v];
188c58f1c22SToby Isaac     tmpS[v] = label->stratumSizes[v];
189c58f1c22SToby Isaac     tmpH[v] = label->ht[v];
190c58f1c22SToby Isaac     tmpP[v] = label->points[v];
191ad8374ffSToby Isaac     tmpB[v] = label->validIS[v];
192c58f1c22SToby Isaac   }
193c58f1c22SToby Isaac   tmpV[v] = value;
194c58f1c22SToby Isaac   tmpS[v] = 0;
195c58f1c22SToby Isaac   PetscHashICreate(tmpH[v]);
196ad8374ffSToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&tmpP[v]);CHKERRQ(ierr);
197c58f1c22SToby Isaac   tmpB[v] = PETSC_TRUE;
198c58f1c22SToby Isaac   ++label->numStrata;
199c58f1c22SToby Isaac   ierr = PetscFree(label->stratumValues);CHKERRQ(ierr);
200c58f1c22SToby Isaac   ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr);
201c58f1c22SToby Isaac   ierr = PetscFree(label->ht);CHKERRQ(ierr);
202c58f1c22SToby Isaac   ierr = PetscFree(label->points);CHKERRQ(ierr);
203ad8374ffSToby Isaac   ierr = PetscFree(label->validIS);CHKERRQ(ierr);
204c58f1c22SToby Isaac   label->stratumValues = tmpV;
205c58f1c22SToby Isaac   label->stratumSizes  = tmpS;
206c58f1c22SToby Isaac   label->ht            = tmpH;
207c58f1c22SToby Isaac   label->points        = tmpP;
208ad8374ffSToby Isaac   label->validIS       = tmpB;
209c58f1c22SToby Isaac 
210c58f1c22SToby Isaac   PetscFunctionReturn(0);
211c58f1c22SToby Isaac }
212c58f1c22SToby Isaac 
213c58f1c22SToby Isaac #undef __FUNCT__
214c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetName"
215c58f1c22SToby Isaac /*@C
216c58f1c22SToby Isaac   DMLabelGetName - Return the name of a DMLabel object
217c58f1c22SToby Isaac 
218c58f1c22SToby Isaac   Input parameter:
219c58f1c22SToby Isaac . label - The DMLabel
220c58f1c22SToby Isaac 
221c58f1c22SToby Isaac   Output parameter:
222c58f1c22SToby Isaac . name - The label name
223c58f1c22SToby Isaac 
224c58f1c22SToby Isaac   Level: beginner
225c58f1c22SToby Isaac 
226c58f1c22SToby Isaac .seealso: DMLabelCreate()
227c58f1c22SToby Isaac @*/
228c58f1c22SToby Isaac PetscErrorCode DMLabelGetName(DMLabel label, const char **name)
229c58f1c22SToby Isaac {
230c58f1c22SToby Isaac   PetscFunctionBegin;
231c58f1c22SToby Isaac   PetscValidPointer(name, 2);
232c58f1c22SToby Isaac   *name = label->name;
233c58f1c22SToby Isaac   PetscFunctionReturn(0);
234c58f1c22SToby Isaac }
235c58f1c22SToby Isaac 
236c58f1c22SToby Isaac #undef __FUNCT__
237c58f1c22SToby Isaac #define __FUNCT__ "DMLabelView_Ascii"
238c58f1c22SToby Isaac static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
239c58f1c22SToby Isaac {
240c58f1c22SToby Isaac   PetscInt       v;
241c58f1c22SToby Isaac   PetscMPIInt    rank;
242c58f1c22SToby Isaac   PetscErrorCode ierr;
243c58f1c22SToby Isaac 
244c58f1c22SToby Isaac   PetscFunctionBegin;
245c58f1c22SToby Isaac   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr);
246c58f1c22SToby Isaac   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
247c58f1c22SToby Isaac   if (label) {
248c58f1c22SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", label->name);CHKERRQ(ierr);
249c58f1c22SToby Isaac     if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);}
250c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
251c58f1c22SToby Isaac       const PetscInt value = label->stratumValues[v];
252ad8374ffSToby Isaac       const PetscInt *points;
253c58f1c22SToby Isaac       PetscInt       p;
254c58f1c22SToby Isaac 
255ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
256c58f1c22SToby Isaac       for (p = 0; p < label->stratumSizes[v]; ++p) {
257ad8374ffSToby Isaac         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr);
258c58f1c22SToby Isaac       }
259ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
260c58f1c22SToby Isaac     }
261c58f1c22SToby Isaac   }
262c58f1c22SToby Isaac   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
263c58f1c22SToby Isaac   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
264c58f1c22SToby Isaac   PetscFunctionReturn(0);
265c58f1c22SToby Isaac }
266c58f1c22SToby Isaac 
267c58f1c22SToby Isaac #undef __FUNCT__
268c58f1c22SToby Isaac #define __FUNCT__ "DMLabelView"
269c58f1c22SToby Isaac /*@C
270c58f1c22SToby Isaac   DMLabelView - View the label
271c58f1c22SToby Isaac 
272c58f1c22SToby Isaac   Input Parameters:
273c58f1c22SToby Isaac + label - The DMLabel
274c58f1c22SToby Isaac - viewer - The PetscViewer
275c58f1c22SToby Isaac 
276c58f1c22SToby Isaac   Level: intermediate
277c58f1c22SToby Isaac 
278c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy()
279c58f1c22SToby Isaac @*/
280c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
281c58f1c22SToby Isaac {
282c58f1c22SToby Isaac   PetscBool      iascii;
283c58f1c22SToby Isaac   PetscErrorCode ierr;
284c58f1c22SToby Isaac 
285c58f1c22SToby Isaac   PetscFunctionBegin;
286c58f1c22SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
287c58f1c22SToby Isaac   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
288c58f1c22SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
289c58f1c22SToby Isaac   if (iascii) {
290c58f1c22SToby Isaac     ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr);
291c58f1c22SToby Isaac   }
292c58f1c22SToby Isaac   PetscFunctionReturn(0);
293c58f1c22SToby Isaac }
294c58f1c22SToby Isaac 
295c58f1c22SToby Isaac #undef __FUNCT__
296c58f1c22SToby Isaac #define __FUNCT__ "DMLabelDestroy"
297c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label)
298c58f1c22SToby Isaac {
299c58f1c22SToby Isaac   PetscInt       v;
300c58f1c22SToby Isaac   PetscErrorCode ierr;
301c58f1c22SToby Isaac 
302c58f1c22SToby Isaac   PetscFunctionBegin;
303c58f1c22SToby Isaac   if (!(*label)) PetscFunctionReturn(0);
304c58f1c22SToby Isaac   if (--(*label)->refct > 0) PetscFunctionReturn(0);
305c58f1c22SToby Isaac   ierr = PetscFree((*label)->name);CHKERRQ(ierr);
306c58f1c22SToby Isaac   ierr = PetscFree((*label)->stratumValues);CHKERRQ(ierr);
307c58f1c22SToby Isaac   ierr = PetscFree((*label)->stratumSizes);CHKERRQ(ierr);
308ad8374ffSToby Isaac   for (v = 0; v < (*label)->numStrata; ++v) {ierr = ISDestroy(&((*label)->points[v]));CHKERRQ(ierr);}
309c58f1c22SToby Isaac   ierr = PetscFree((*label)->points);CHKERRQ(ierr);
310ad8374ffSToby Isaac   ierr = PetscFree((*label)->validIS);CHKERRQ(ierr);
311c58f1c22SToby Isaac   if ((*label)->ht) {
312c58f1c22SToby Isaac     for (v = 0; v < (*label)->numStrata; ++v) {PetscHashIDestroy((*label)->ht[v]);}
313c58f1c22SToby Isaac     ierr = PetscFree((*label)->ht);CHKERRQ(ierr);
314c58f1c22SToby Isaac   }
315c58f1c22SToby Isaac   ierr = PetscBTDestroy(&(*label)->bt);CHKERRQ(ierr);
316c58f1c22SToby Isaac   ierr = PetscFree(*label);CHKERRQ(ierr);
317c58f1c22SToby Isaac   PetscFunctionReturn(0);
318c58f1c22SToby Isaac }
319c58f1c22SToby Isaac 
320c58f1c22SToby Isaac #undef __FUNCT__
321c58f1c22SToby Isaac #define __FUNCT__ "DMLabelDuplicate"
322c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
323c58f1c22SToby Isaac {
324ad8374ffSToby Isaac   PetscInt       v;
325c58f1c22SToby Isaac   PetscErrorCode ierr;
326c58f1c22SToby Isaac 
327c58f1c22SToby Isaac   PetscFunctionBegin;
328c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
329c58f1c22SToby Isaac   ierr = PetscNew(labelnew);CHKERRQ(ierr);
330c58f1c22SToby Isaac   ierr = PetscStrallocpy(label->name, &(*labelnew)->name);CHKERRQ(ierr);
331c58f1c22SToby Isaac 
332c58f1c22SToby Isaac   (*labelnew)->refct        = 1;
333c58f1c22SToby Isaac   (*labelnew)->numStrata    = label->numStrata;
3345aa44df4SToby Isaac   (*labelnew)->defaultValue = label->defaultValue;
335c58f1c22SToby Isaac   if (label->numStrata) {
336c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr);
337c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr);
338c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr);
339c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr);
340ad8374ffSToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr);
341c58f1c22SToby Isaac     /* Could eliminate unused space here */
342c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
343c58f1c22SToby Isaac       PetscHashICreate((*labelnew)->ht[v]);
344ad8374ffSToby Isaac       (*labelnew)->validIS[v]        = PETSC_TRUE;
345c58f1c22SToby Isaac       (*labelnew)->stratumValues[v]  = label->stratumValues[v];
346c58f1c22SToby Isaac       (*labelnew)->stratumSizes[v]   = label->stratumSizes[v];
347ad8374ffSToby Isaac       ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr);
348ad8374ffSToby Isaac       (*labelnew)->points[v]         = label->points[v];
349c58f1c22SToby Isaac     }
350c58f1c22SToby Isaac   }
351c58f1c22SToby Isaac   (*labelnew)->pStart = -1;
352c58f1c22SToby Isaac   (*labelnew)->pEnd   = -1;
353c58f1c22SToby Isaac   (*labelnew)->bt     = NULL;
354c58f1c22SToby Isaac   PetscFunctionReturn(0);
355c58f1c22SToby Isaac }
356c58f1c22SToby Isaac 
357c58f1c22SToby Isaac #undef __FUNCT__
358c58f1c22SToby Isaac #define __FUNCT__ "DMLabelCreateIndex"
359c58f1c22SToby Isaac /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
360c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
361c58f1c22SToby Isaac {
362c58f1c22SToby Isaac   PetscInt       v;
363c58f1c22SToby Isaac   PetscErrorCode ierr;
364c58f1c22SToby Isaac 
365c58f1c22SToby Isaac   PetscFunctionBegin;
366c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
367c58f1c22SToby Isaac   if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);}
368c58f1c22SToby Isaac   label->pStart = pStart;
369c58f1c22SToby Isaac   label->pEnd   = pEnd;
370c58f1c22SToby Isaac   ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr);
371c58f1c22SToby Isaac   ierr = PetscBTMemzero(pEnd - pStart, label->bt);CHKERRQ(ierr);
372c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
373ad8374ffSToby Isaac     const PetscInt *points;
374c58f1c22SToby Isaac     PetscInt       i;
375c58f1c22SToby Isaac 
376ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
377c58f1c22SToby Isaac     for (i = 0; i < label->stratumSizes[v]; ++i) {
378ad8374ffSToby Isaac       const PetscInt point = points[i];
379c58f1c22SToby Isaac 
380c58f1c22SToby 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);
381c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr);
382c58f1c22SToby Isaac     }
383ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
384c58f1c22SToby Isaac   }
385c58f1c22SToby Isaac   PetscFunctionReturn(0);
386c58f1c22SToby Isaac }
387c58f1c22SToby Isaac 
388c58f1c22SToby Isaac #undef __FUNCT__
389c58f1c22SToby Isaac #define __FUNCT__ "DMLabelDestroyIndex"
390c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label)
391c58f1c22SToby Isaac {
392c58f1c22SToby Isaac   PetscErrorCode ierr;
393c58f1c22SToby Isaac 
394c58f1c22SToby Isaac   PetscFunctionBegin;
395c58f1c22SToby Isaac   label->pStart = -1;
396c58f1c22SToby Isaac   label->pEnd   = -1;
397c58f1c22SToby Isaac   if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);}
398c58f1c22SToby Isaac   PetscFunctionReturn(0);
399c58f1c22SToby Isaac }
400c58f1c22SToby Isaac 
401c58f1c22SToby Isaac #undef __FUNCT__
402c58f1c22SToby Isaac #define __FUNCT__ "DMLabelHasValue"
403c58f1c22SToby Isaac /*@
404c58f1c22SToby Isaac   DMLabelHasValue - Determine whether a label assigns the value to any point
405c58f1c22SToby Isaac 
406c58f1c22SToby Isaac   Input Parameters:
407c58f1c22SToby Isaac + label - the DMLabel
408c58f1c22SToby Isaac - value - the value
409c58f1c22SToby Isaac 
410c58f1c22SToby Isaac   Output Parameter:
411c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point
412c58f1c22SToby Isaac 
413c58f1c22SToby Isaac   Level: developer
414c58f1c22SToby Isaac 
415c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
416c58f1c22SToby Isaac @*/
417c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
418c58f1c22SToby Isaac {
419c58f1c22SToby Isaac   PetscInt v;
420c58f1c22SToby Isaac 
421c58f1c22SToby Isaac   PetscFunctionBegin;
422c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
423c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
424c58f1c22SToby Isaac     if (value == label->stratumValues[v]) break;
425c58f1c22SToby Isaac   }
426c58f1c22SToby Isaac   *contains = (v < label->numStrata ? PETSC_TRUE : PETSC_FALSE);
427c58f1c22SToby Isaac   PetscFunctionReturn(0);
428c58f1c22SToby Isaac }
429c58f1c22SToby Isaac 
430c58f1c22SToby Isaac #undef __FUNCT__
431c58f1c22SToby Isaac #define __FUNCT__ "DMLabelHasPoint"
432c58f1c22SToby Isaac /*@
433c58f1c22SToby Isaac   DMLabelHasPoint - Determine whether a label assigns a value to a point
434c58f1c22SToby Isaac 
435c58f1c22SToby Isaac   Input Parameters:
436c58f1c22SToby Isaac + label - the DMLabel
437c58f1c22SToby Isaac - point - the point
438c58f1c22SToby Isaac 
439c58f1c22SToby Isaac   Output Parameter:
440c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value
441c58f1c22SToby Isaac 
442c58f1c22SToby Isaac   Note: The user must call DMLabelCreateIndex() before this function.
443c58f1c22SToby Isaac 
444c58f1c22SToby Isaac   Level: developer
445c58f1c22SToby Isaac 
446c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
447c58f1c22SToby Isaac @*/
448c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
449c58f1c22SToby Isaac {
450c58f1c22SToby Isaac   PetscErrorCode ierr;
451c58f1c22SToby Isaac 
452c58f1c22SToby Isaac   PetscFunctionBeginHot;
453c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
454c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
455c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG)
456c58f1c22SToby Isaac   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
457c58f1c22SToby 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);
458c58f1c22SToby Isaac #endif
459c58f1c22SToby Isaac   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
460c58f1c22SToby Isaac   PetscFunctionReturn(0);
461c58f1c22SToby Isaac }
462c58f1c22SToby Isaac 
463c58f1c22SToby Isaac #undef __FUNCT__
464c58f1c22SToby Isaac #define __FUNCT__ "DMLabelStratumHasPoint"
465c58f1c22SToby Isaac /*@
466c58f1c22SToby Isaac   DMLabelStratumHasPoint - Return true if the stratum contains a point
467c58f1c22SToby Isaac 
468c58f1c22SToby Isaac   Input Parameters:
469c58f1c22SToby Isaac + label - the DMLabel
470c58f1c22SToby Isaac . value - the stratum value
471c58f1c22SToby Isaac - point - the point
472c58f1c22SToby Isaac 
473c58f1c22SToby Isaac   Output Parameter:
474c58f1c22SToby Isaac . contains - true if the stratum contains the point
475c58f1c22SToby Isaac 
476c58f1c22SToby Isaac   Level: intermediate
477c58f1c22SToby Isaac 
478c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
479c58f1c22SToby Isaac @*/
480c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
481c58f1c22SToby Isaac {
482c58f1c22SToby Isaac   PetscInt       v;
483c58f1c22SToby Isaac   PetscErrorCode ierr;
484c58f1c22SToby Isaac 
485c58f1c22SToby Isaac   PetscFunctionBegin;
486c58f1c22SToby Isaac   PetscValidPointer(contains, 4);
487c58f1c22SToby Isaac   *contains = PETSC_FALSE;
488c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
489c58f1c22SToby Isaac     if (label->stratumValues[v] == value) {
490ad8374ffSToby Isaac       if (label->validIS[v]) {
491c58f1c22SToby Isaac         PetscInt i;
492c58f1c22SToby Isaac 
493a2d74346SToby Isaac         ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr);
494c58f1c22SToby Isaac         if (i >= 0) {
495c58f1c22SToby Isaac           *contains = PETSC_TRUE;
496c58f1c22SToby Isaac           break;
497c58f1c22SToby Isaac         }
498c58f1c22SToby Isaac       } else {
499c58f1c22SToby Isaac         PetscBool has;
500c58f1c22SToby Isaac 
501c58f1c22SToby Isaac         PetscHashIHasKey(label->ht[v], point, has);
502c58f1c22SToby Isaac         if (has) {
503c58f1c22SToby Isaac           *contains = PETSC_TRUE;
504c58f1c22SToby Isaac           break;
505c58f1c22SToby Isaac         }
506c58f1c22SToby Isaac       }
507c58f1c22SToby Isaac     }
508c58f1c22SToby Isaac   }
509c58f1c22SToby Isaac   PetscFunctionReturn(0);
510c58f1c22SToby Isaac }
511c58f1c22SToby Isaac 
512c58f1c22SToby Isaac #undef __FUNCT__
5135aa44df4SToby Isaac #define __FUNCT__ "DMLabelGetDefaultValue"
5145aa44df4SToby Isaac /*
5155aa44df4SToby Isaac   DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
5165aa44df4SToby Isaac   When a label is created, it is initialized to -1.
5175aa44df4SToby Isaac 
5185aa44df4SToby Isaac   Input parameter:
5195aa44df4SToby Isaac . label - a DMLabel object
5205aa44df4SToby Isaac 
5215aa44df4SToby Isaac   Output parameter:
5225aa44df4SToby Isaac . defaultValue - the default value
5235aa44df4SToby Isaac 
5245aa44df4SToby Isaac   Level: beginner
5255aa44df4SToby Isaac 
5265aa44df4SToby Isaac .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
5275aa44df4SToby Isaac */
5285aa44df4SToby Isaac PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
5295aa44df4SToby Isaac {
5305aa44df4SToby Isaac   PetscFunctionBegin;
5315aa44df4SToby Isaac   *defaultValue = label->defaultValue;
5325aa44df4SToby Isaac   PetscFunctionReturn(0);
5335aa44df4SToby Isaac }
5345aa44df4SToby Isaac 
5355aa44df4SToby Isaac #undef __FUNCT__
5365aa44df4SToby Isaac #define __FUNCT__ "DMLabelSetDefaultValue"
5375aa44df4SToby Isaac /*
5385aa44df4SToby Isaac   DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
5395aa44df4SToby Isaac   When a label is created, it is initialized to -1.
5405aa44df4SToby Isaac 
5415aa44df4SToby Isaac   Input parameter:
5425aa44df4SToby Isaac . label - a DMLabel object
5435aa44df4SToby Isaac 
5445aa44df4SToby Isaac   Output parameter:
5455aa44df4SToby Isaac . defaultValue - the default value
5465aa44df4SToby Isaac 
5475aa44df4SToby Isaac   Level: beginner
5485aa44df4SToby Isaac 
5495aa44df4SToby Isaac .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
5505aa44df4SToby Isaac */
5515aa44df4SToby Isaac PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
5525aa44df4SToby Isaac {
5535aa44df4SToby Isaac   PetscFunctionBegin;
5545aa44df4SToby Isaac   label->defaultValue = defaultValue;
5555aa44df4SToby Isaac   PetscFunctionReturn(0);
5565aa44df4SToby Isaac }
5575aa44df4SToby Isaac 
5585aa44df4SToby Isaac #undef __FUNCT__
559c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetValue"
560c58f1c22SToby Isaac /*@
5615aa44df4SToby 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())
562c58f1c22SToby Isaac 
563c58f1c22SToby Isaac   Input Parameters:
564c58f1c22SToby Isaac + label - the DMLabel
565c58f1c22SToby Isaac - point - the point
566c58f1c22SToby Isaac 
567c58f1c22SToby Isaac   Output Parameter:
568c58f1c22SToby Isaac . value - The point value, or -1
569c58f1c22SToby Isaac 
570c58f1c22SToby Isaac   Level: intermediate
571c58f1c22SToby Isaac 
5725aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
573c58f1c22SToby Isaac @*/
574c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
575c58f1c22SToby Isaac {
576c58f1c22SToby Isaac   PetscInt       v;
577c58f1c22SToby Isaac   PetscErrorCode ierr;
578c58f1c22SToby Isaac 
579c58f1c22SToby Isaac   PetscFunctionBegin;
580c58f1c22SToby Isaac   PetscValidPointer(value, 3);
5815aa44df4SToby Isaac   *value = label->defaultValue;
582c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
583ad8374ffSToby Isaac     if (label->validIS[v]) {
584c58f1c22SToby Isaac       PetscInt i;
585c58f1c22SToby Isaac 
586a2d74346SToby Isaac       ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr);
587c58f1c22SToby Isaac       if (i >= 0) {
588c58f1c22SToby Isaac         *value = label->stratumValues[v];
589c58f1c22SToby Isaac         break;
590c58f1c22SToby Isaac       }
591c58f1c22SToby Isaac     } else {
592c58f1c22SToby Isaac       PetscBool has;
593c58f1c22SToby Isaac 
594c58f1c22SToby Isaac       PetscHashIHasKey(label->ht[v], point, has);
595c58f1c22SToby Isaac       if (has) {
596c58f1c22SToby Isaac         *value = label->stratumValues[v];
597c58f1c22SToby Isaac         break;
598c58f1c22SToby Isaac       }
599c58f1c22SToby Isaac     }
600c58f1c22SToby Isaac   }
601c58f1c22SToby Isaac   PetscFunctionReturn(0);
602c58f1c22SToby Isaac }
603c58f1c22SToby Isaac 
604c58f1c22SToby Isaac #undef __FUNCT__
605c58f1c22SToby Isaac #define __FUNCT__ "DMLabelSetValue"
606c58f1c22SToby Isaac /*@
6075aa44df4SToby 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.
608c58f1c22SToby Isaac 
609c58f1c22SToby Isaac   Input Parameters:
610c58f1c22SToby Isaac + label - the DMLabel
611c58f1c22SToby Isaac . point - the point
612c58f1c22SToby Isaac - value - The point value
613c58f1c22SToby Isaac 
614c58f1c22SToby Isaac   Level: intermediate
615c58f1c22SToby Isaac 
6165aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
617c58f1c22SToby Isaac @*/
618c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
619c58f1c22SToby Isaac {
620c58f1c22SToby Isaac   PETSC_UNUSED PetscHashIIter iter, ret;
621c58f1c22SToby Isaac   PetscInt                    v;
622c58f1c22SToby Isaac   PetscErrorCode              ierr;
623c58f1c22SToby Isaac 
624c58f1c22SToby Isaac   PetscFunctionBegin;
625c58f1c22SToby Isaac   /* Find, or add, label value */
6265aa44df4SToby Isaac   if (value == label->defaultValue) PetscFunctionReturn(0);
627c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
628c58f1c22SToby Isaac     if (label->stratumValues[v] == value) break;
629c58f1c22SToby Isaac   }
630c58f1c22SToby Isaac   /* Create new table */
631c58f1c22SToby Isaac   if (v >= label->numStrata) {ierr = DMLabelAddStratum(label, value);CHKERRQ(ierr);}
632c58f1c22SToby Isaac   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
633c58f1c22SToby Isaac   /* Set key */
634c58f1c22SToby Isaac   PetscHashIPut(label->ht[v], point, ret, iter);
635c58f1c22SToby Isaac   PetscFunctionReturn(0);
636c58f1c22SToby Isaac }
637c58f1c22SToby Isaac 
638c58f1c22SToby Isaac #undef __FUNCT__
639c58f1c22SToby Isaac #define __FUNCT__ "DMLabelClearValue"
640c58f1c22SToby Isaac /*@
641c58f1c22SToby Isaac   DMLabelClearValue - Clear the value a label assigns to a point
642c58f1c22SToby Isaac 
643c58f1c22SToby Isaac   Input Parameters:
644c58f1c22SToby Isaac + label - the DMLabel
645c58f1c22SToby Isaac . point - the point
646c58f1c22SToby Isaac - value - The point value
647c58f1c22SToby Isaac 
648c58f1c22SToby Isaac   Level: intermediate
649c58f1c22SToby Isaac 
650c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
651c58f1c22SToby Isaac @*/
652c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
653c58f1c22SToby Isaac {
654ad8374ffSToby Isaac   PetscInt       v;
655c58f1c22SToby Isaac   PetscErrorCode ierr;
656c58f1c22SToby Isaac 
657c58f1c22SToby Isaac   PetscFunctionBegin;
658c58f1c22SToby Isaac   /* Find label value */
659c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
660c58f1c22SToby Isaac     if (label->stratumValues[v] == value) break;
661c58f1c22SToby Isaac   }
662c58f1c22SToby Isaac   if (v >= label->numStrata) PetscFunctionReturn(0);
663eeed21e7SToby Isaac   if (label->validIS[v]) {
664eeed21e7SToby Isaac     ierr = DMLabelMakeInvalid_Private(label,v);CHKERRQ(ierr);
665eeed21e7SToby Isaac   }
666eeed21e7SToby Isaac   if (label->bt) {
667eeed21e7SToby 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);
668eeed21e7SToby Isaac     ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
669eeed21e7SToby Isaac   }
670c58f1c22SToby Isaac   ierr = PetscHashIDelKey(label->ht[v], point);CHKERRQ(ierr);
671c58f1c22SToby Isaac   PetscFunctionReturn(0);
672c58f1c22SToby Isaac }
673c58f1c22SToby Isaac 
674c58f1c22SToby Isaac #undef __FUNCT__
675c58f1c22SToby Isaac #define __FUNCT__ "DMLabelInsertIS"
676c58f1c22SToby Isaac /*@
677c58f1c22SToby Isaac   DMLabelInsertIS - Set all points in the IS to a value
678c58f1c22SToby Isaac 
679c58f1c22SToby Isaac   Input Parameters:
680c58f1c22SToby Isaac + label - the DMLabel
681c58f1c22SToby Isaac . is    - the point IS
682c58f1c22SToby Isaac - value - The point value
683c58f1c22SToby Isaac 
684c58f1c22SToby Isaac   Level: intermediate
685c58f1c22SToby Isaac 
686c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
687c58f1c22SToby Isaac @*/
688c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
689c58f1c22SToby Isaac {
690c58f1c22SToby Isaac   const PetscInt *points;
691c58f1c22SToby Isaac   PetscInt        n, p;
692c58f1c22SToby Isaac   PetscErrorCode  ierr;
693c58f1c22SToby Isaac 
694c58f1c22SToby Isaac   PetscFunctionBegin;
695c58f1c22SToby Isaac   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
696c58f1c22SToby Isaac   ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr);
697c58f1c22SToby Isaac   ierr = ISGetIndices(is, &points);CHKERRQ(ierr);
698c58f1c22SToby Isaac   for (p = 0; p < n; ++p) {ierr = DMLabelSetValue(label, points[p], value);CHKERRQ(ierr);}
699c58f1c22SToby Isaac   ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr);
700c58f1c22SToby Isaac   PetscFunctionReturn(0);
701c58f1c22SToby Isaac }
702c58f1c22SToby Isaac 
703c58f1c22SToby Isaac #undef __FUNCT__
704c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetNumValues"
705c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
706c58f1c22SToby Isaac {
707c58f1c22SToby Isaac   PetscFunctionBegin;
708c58f1c22SToby Isaac   PetscValidPointer(numValues, 2);
709c58f1c22SToby Isaac   *numValues = label->numStrata;
710c58f1c22SToby Isaac   PetscFunctionReturn(0);
711c58f1c22SToby Isaac }
712c58f1c22SToby Isaac 
713c58f1c22SToby Isaac #undef __FUNCT__
714c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetValueIS"
715c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
716c58f1c22SToby Isaac {
717c58f1c22SToby Isaac   PetscErrorCode ierr;
718c58f1c22SToby Isaac 
719c58f1c22SToby Isaac   PetscFunctionBegin;
720c58f1c22SToby Isaac   PetscValidPointer(values, 2);
721c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr);
722c58f1c22SToby Isaac   PetscFunctionReturn(0);
723c58f1c22SToby Isaac }
724c58f1c22SToby Isaac 
725c58f1c22SToby Isaac #undef __FUNCT__
726fada774cSMatthew G. Knepley #define __FUNCT__ "DMLabelHasStratum"
727fada774cSMatthew G. Knepley PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
728fada774cSMatthew G. Knepley {
729fada774cSMatthew G. Knepley   PetscInt v;
730fada774cSMatthew G. Knepley 
731fada774cSMatthew G. Knepley   PetscFunctionBegin;
732fada774cSMatthew G. Knepley   PetscValidPointer(exists, 3);
733fada774cSMatthew G. Knepley   *exists = PETSC_FALSE;
734fada774cSMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
735fada774cSMatthew G. Knepley     if (label->stratumValues[v] == value) {
736fada774cSMatthew G. Knepley       *exists = PETSC_TRUE;
737fada774cSMatthew G. Knepley       break;
738fada774cSMatthew G. Knepley     }
739fada774cSMatthew G. Knepley   }
740fada774cSMatthew G. Knepley   PetscFunctionReturn(0);
741fada774cSMatthew G. Knepley }
742fada774cSMatthew G. Knepley 
743fada774cSMatthew G. Knepley #undef __FUNCT__
744c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetStratumSize"
745c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
746c58f1c22SToby Isaac {
747c58f1c22SToby Isaac   PetscInt       v;
748c58f1c22SToby Isaac   PetscErrorCode ierr;
749c58f1c22SToby Isaac 
750c58f1c22SToby Isaac   PetscFunctionBegin;
751c58f1c22SToby Isaac   PetscValidPointer(size, 3);
752c58f1c22SToby Isaac   *size = 0;
753c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
754c58f1c22SToby Isaac     if (label->stratumValues[v] == value) {
755c58f1c22SToby Isaac       ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
756c58f1c22SToby Isaac       *size = label->stratumSizes[v];
757c58f1c22SToby Isaac       break;
758c58f1c22SToby Isaac     }
759c58f1c22SToby Isaac   }
760c58f1c22SToby Isaac   PetscFunctionReturn(0);
761c58f1c22SToby Isaac }
762c58f1c22SToby Isaac 
763c58f1c22SToby Isaac #undef __FUNCT__
764c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetStratumBounds"
765c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
766c58f1c22SToby Isaac {
767c58f1c22SToby Isaac   PetscInt       v;
768c58f1c22SToby Isaac   PetscErrorCode ierr;
769c58f1c22SToby Isaac 
770c58f1c22SToby Isaac   PetscFunctionBegin;
771c58f1c22SToby Isaac   if (start) {PetscValidPointer(start, 3); *start = 0;}
772c58f1c22SToby Isaac   if (end)   {PetscValidPointer(end,   4); *end   = 0;}
773c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
774ad8374ffSToby Isaac     const PetscInt *points;
775ad8374ffSToby Isaac 
776c58f1c22SToby Isaac     if (label->stratumValues[v] != value) continue;
777c58f1c22SToby Isaac     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
778c58f1c22SToby Isaac     if (label->stratumSizes[v]  <= 0)     break;
779ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
780ad8374ffSToby Isaac     if (start) *start = points[0];
781ad8374ffSToby Isaac     if (end)   *end   = points[label->stratumSizes[v]-1]+1;
782ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
783c58f1c22SToby Isaac     break;
784c58f1c22SToby Isaac   }
785c58f1c22SToby Isaac   PetscFunctionReturn(0);
786c58f1c22SToby Isaac }
787c58f1c22SToby Isaac 
788c58f1c22SToby Isaac #undef __FUNCT__
789c58f1c22SToby Isaac #define __FUNCT__ "DMLabelGetStratumIS"
790c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
791c58f1c22SToby Isaac {
792c58f1c22SToby Isaac   PetscInt       v;
793c58f1c22SToby Isaac   PetscErrorCode ierr;
794c58f1c22SToby Isaac 
795c58f1c22SToby Isaac   PetscFunctionBegin;
796c58f1c22SToby Isaac   PetscValidPointer(points, 3);
797c58f1c22SToby Isaac   *points = NULL;
798c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
799c58f1c22SToby Isaac     if (label->stratumValues[v] == value) {
800c58f1c22SToby Isaac       ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
801ad8374ffSToby Isaac       if (label->validIS[v]) {
802ad8374ffSToby Isaac         ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr);
803ad8374ffSToby Isaac         *points = label->points[v];
804c58f1c22SToby Isaac       } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to implement this to speedup Stratify");
805c58f1c22SToby Isaac       break;
806c58f1c22SToby Isaac     }
807c58f1c22SToby Isaac   }
808c58f1c22SToby Isaac   PetscFunctionReturn(0);
809c58f1c22SToby Isaac }
810c58f1c22SToby Isaac 
811c58f1c22SToby Isaac #undef __FUNCT__
8124de306b1SToby Isaac #define __FUNCT__ "DMLabelSetStratumIS"
8134de306b1SToby Isaac PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
8144de306b1SToby Isaac {
8154de306b1SToby Isaac   PetscInt       v, numStrata;
8164de306b1SToby Isaac   PetscErrorCode ierr;
8174de306b1SToby Isaac 
8184de306b1SToby Isaac   PetscFunctionBegin;
8194de306b1SToby Isaac   numStrata = label->numStrata;
8204de306b1SToby Isaac   for (v = 0; v < numStrata; v++) {
8214de306b1SToby Isaac     if (label->stratumValues[v] == value) break;
8224de306b1SToby Isaac   }
8234de306b1SToby Isaac   if (v >= numStrata) {ierr = DMLabelAddStratum(label,value);CHKERRQ(ierr);}
8244de306b1SToby Isaac   if (is == label->points[v]) PetscFunctionReturn(0);
8254de306b1SToby Isaac   ierr = DMLabelClearStratum(label,value);CHKERRQ(ierr);
8264de306b1SToby Isaac   ierr = ISGetLocalSize(is,&(label->stratumSizes[v]));CHKERRQ(ierr);
8274de306b1SToby Isaac   label->stratumValues[v] = value;
8284de306b1SToby Isaac   label->validIS[v] = PETSC_TRUE;
8294de306b1SToby Isaac   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
8304de306b1SToby Isaac   ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
8314de306b1SToby Isaac   if (label->bt) {
8324de306b1SToby Isaac     const PetscInt *points;
8334de306b1SToby Isaac     PetscInt p;
8344de306b1SToby Isaac 
8354de306b1SToby Isaac     ierr = ISGetIndices(is,&points);CHKERRQ(ierr);
8364de306b1SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
8374de306b1SToby Isaac       const PetscInt point = points[p];
8384de306b1SToby Isaac 
8394de306b1SToby 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);
8404de306b1SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
8414de306b1SToby Isaac     }
8424de306b1SToby Isaac   }
8434de306b1SToby Isaac   label->points[v] = is;
8444de306b1SToby Isaac   PetscFunctionReturn(0);
8454de306b1SToby Isaac }
8464de306b1SToby Isaac 
8474de306b1SToby Isaac 
8484de306b1SToby Isaac #undef __FUNCT__
849c58f1c22SToby Isaac #define __FUNCT__ "DMLabelClearStratum"
850c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
851c58f1c22SToby Isaac {
852c58f1c22SToby Isaac   PetscInt       v;
853c58f1c22SToby Isaac   PetscErrorCode ierr;
854c58f1c22SToby Isaac 
855c58f1c22SToby Isaac   PetscFunctionBegin;
856c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
857c58f1c22SToby Isaac     if (label->stratumValues[v] == value) break;
858c58f1c22SToby Isaac   }
859c58f1c22SToby Isaac   if (v >= label->numStrata) PetscFunctionReturn(0);
8604de306b1SToby Isaac   if (label->validIS[v]) {
8614de306b1SToby Isaac     if (label->bt) {
862c58f1c22SToby Isaac       PetscInt       i;
863ad8374ffSToby Isaac       const PetscInt *points;
864c58f1c22SToby Isaac 
865ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
866c58f1c22SToby Isaac       for (i = 0; i < label->stratumSizes[v]; ++i) {
867ad8374ffSToby Isaac         const PetscInt point = points[i];
868c58f1c22SToby Isaac 
869c58f1c22SToby 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);
870c58f1c22SToby Isaac         ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
871c58f1c22SToby Isaac       }
872ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
873c58f1c22SToby Isaac     }
874ad8374ffSToby Isaac     ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
875c58f1c22SToby Isaac     label->stratumSizes[v] = 0;
876ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
877ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
878c58f1c22SToby Isaac   } else {
879c58f1c22SToby Isaac     PetscHashIClear(label->ht[v]);
880c58f1c22SToby Isaac   }
881c58f1c22SToby Isaac   PetscFunctionReturn(0);
882c58f1c22SToby Isaac }
883c58f1c22SToby Isaac 
884c58f1c22SToby Isaac #undef __FUNCT__
885c58f1c22SToby Isaac #define __FUNCT__ "DMLabelFilter"
886c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
887c58f1c22SToby Isaac {
888c58f1c22SToby Isaac   PetscInt       v;
889c58f1c22SToby Isaac   PetscErrorCode ierr;
890c58f1c22SToby Isaac 
891c58f1c22SToby Isaac   PetscFunctionBegin;
892c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
893c58f1c22SToby Isaac   label->pStart = start;
894c58f1c22SToby Isaac   label->pEnd   = end;
895c58f1c22SToby Isaac   if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);}
896c58f1c22SToby Isaac   /* Could squish offsets, but would only make sense if I reallocate the storage */
897c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
898c58f1c22SToby Isaac     PetscInt off, q;
899ad8374ffSToby Isaac     const PetscInt *points;
900ad8374ffSToby Isaac     PetscInt *pointsNew = NULL;
901c58f1c22SToby Isaac 
902ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
903c58f1c22SToby Isaac     for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) {
904ad8374ffSToby Isaac       const PetscInt point = points[q];
905c58f1c22SToby Isaac 
906ad8374ffSToby Isaac       if ((point < start) || (point >= end)) {
907ad8374ffSToby Isaac         if (!pointsNew) {
908ad8374ffSToby Isaac           ierr = PetscMalloc1(label->stratumSizes[v],&pointsNew);CHKERRQ(ierr);
909ad8374ffSToby Isaac           ierr = PetscMemcpy(pointsNew,points,(size_t) off * sizeof(PetscInt));CHKERRQ(ierr);
910ad8374ffSToby Isaac         }
911ad8374ffSToby Isaac         continue;
912ad8374ffSToby Isaac       }
913ad8374ffSToby Isaac       if (pointsNew) {
914ad8374ffSToby Isaac         pointsNew[off++] = point;
915ad8374ffSToby Isaac       }
916ad8374ffSToby Isaac     }
917ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
918ad8374ffSToby Isaac     if (pointsNew) {
919ad8374ffSToby Isaac       ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
920ad8374ffSToby Isaac       ierr = ISCreateGeneral(PETSC_COMM_SELF,off,pointsNew,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
921ad8374ffSToby Isaac       ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
922c58f1c22SToby Isaac     }
923c58f1c22SToby Isaac     label->stratumSizes[v] = off;
924c58f1c22SToby Isaac   }
925c58f1c22SToby Isaac   ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr);
926c58f1c22SToby Isaac   PetscFunctionReturn(0);
927c58f1c22SToby Isaac }
928c58f1c22SToby Isaac 
929c58f1c22SToby Isaac #undef __FUNCT__
930c58f1c22SToby Isaac #define __FUNCT__ "DMLabelPermute"
931c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
932c58f1c22SToby Isaac {
933c58f1c22SToby Isaac   const PetscInt *perm;
934c58f1c22SToby Isaac   PetscInt        numValues, numPoints, v, q;
935c58f1c22SToby Isaac   PetscErrorCode  ierr;
936c58f1c22SToby Isaac 
937c58f1c22SToby Isaac   PetscFunctionBegin;
938c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
939c58f1c22SToby Isaac   ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr);
940c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr);
941c58f1c22SToby Isaac   ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr);
942c58f1c22SToby Isaac   ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr);
943c58f1c22SToby Isaac   for (v = 0; v < numValues; ++v) {
944c58f1c22SToby Isaac     const PetscInt size   = (*labelNew)->stratumSizes[v];
945ad8374ffSToby Isaac     const PetscInt *points;
946ad8374ffSToby Isaac     PetscInt *pointsNew;
947c58f1c22SToby Isaac 
948ad8374ffSToby Isaac     ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
949ad8374ffSToby Isaac     ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr);
950c58f1c22SToby Isaac     for (q = 0; q < size; ++q) {
951ad8374ffSToby Isaac       const PetscInt point = points[q];
952c58f1c22SToby Isaac 
953c58f1c22SToby 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);
954ad8374ffSToby Isaac       pointsNew[q] = perm[point];
955c58f1c22SToby Isaac     }
956ad8374ffSToby Isaac     ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
957ad8374ffSToby Isaac     ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr);
958ad8374ffSToby Isaac     ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr);
959ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr);
960ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr);
961c58f1c22SToby Isaac   }
962c58f1c22SToby Isaac   ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr);
963c58f1c22SToby Isaac   if (label->bt) {
964c58f1c22SToby Isaac     ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
965c58f1c22SToby Isaac     ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
966c58f1c22SToby Isaac   }
967c58f1c22SToby Isaac   PetscFunctionReturn(0);
968c58f1c22SToby Isaac }
969c58f1c22SToby Isaac 
970c58f1c22SToby Isaac #undef __FUNCT__
97126c55118SMichael Lange #define __FUNCT__ "DMLabelDistribute_Internal"
97226c55118SMichael Lange PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
97326c55118SMichael Lange {
97426c55118SMichael Lange   MPI_Comm       comm;
97526c55118SMichael Lange   PetscInt       s, l, nroots, nleaves, dof, offset, size;
97626c55118SMichael Lange   PetscInt      *remoteOffsets, *rootStrata, *rootIdx;
97726c55118SMichael Lange   PetscSection   rootSection;
97826c55118SMichael Lange   PetscSF        labelSF;
97926c55118SMichael Lange   PetscErrorCode ierr;
98026c55118SMichael Lange 
98126c55118SMichael Lange   PetscFunctionBegin;
98226c55118SMichael Lange   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
98326c55118SMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
98426c55118SMichael Lange   /* Build a section of stratum values per point, generate the according SF
98526c55118SMichael Lange      and distribute point-wise stratum values to leaves. */
98626c55118SMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
98726c55118SMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
98826c55118SMichael Lange   ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr);
98926c55118SMichael Lange   if (label) {
99026c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
991ad8374ffSToby Isaac       const PetscInt *points;
992ad8374ffSToby Isaac 
993ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
99426c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
995ad8374ffSToby Isaac         ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr);
996ad8374ffSToby Isaac         ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr);
99726c55118SMichael Lange       }
998ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
99926c55118SMichael Lange     }
100026c55118SMichael Lange   }
100126c55118SMichael Lange   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
100226c55118SMichael Lange   /* Create a point-wise array of stratum values */
100326c55118SMichael Lange   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
100426c55118SMichael Lange   ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr);
100526c55118SMichael Lange   ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr);
100626c55118SMichael Lange   if (label) {
100726c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1008ad8374ffSToby Isaac       const PetscInt *points;
1009ad8374ffSToby Isaac 
1010ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
101126c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1012ad8374ffSToby Isaac         const PetscInt p = points[l];
101326c55118SMichael Lange         ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr);
101426c55118SMichael Lange         rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
101526c55118SMichael Lange       }
1016ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
101726c55118SMichael Lange     }
101826c55118SMichael Lange   }
101926c55118SMichael Lange   /* Build SF that maps label points to remote processes */
102026c55118SMichael Lange   ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr);
102126c55118SMichael Lange   ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr);
102226c55118SMichael Lange   ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr);
102326c55118SMichael Lange   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
102426c55118SMichael Lange   /* Send the strata for each point over the derived SF */
102526c55118SMichael Lange   ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr);
102626c55118SMichael Lange   ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr);
102726c55118SMichael Lange   ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
102826c55118SMichael Lange   ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
102926c55118SMichael Lange   /* Clean up */
103026c55118SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
103126c55118SMichael Lange   ierr = PetscFree(rootIdx);CHKERRQ(ierr);
103226c55118SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
103326c55118SMichael Lange   ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr);
103426c55118SMichael Lange   PetscFunctionReturn(0);
103526c55118SMichael Lange }
103626c55118SMichael Lange 
103726c55118SMichael Lange #undef __FUNCT__
1038c58f1c22SToby Isaac #define __FUNCT__ "DMLabelDistribute"
1039c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1040c58f1c22SToby Isaac {
1041c58f1c22SToby Isaac   MPI_Comm       comm;
104226c55118SMichael Lange   PetscSection   leafSection;
104326c55118SMichael Lange   PetscInt       p, pStart, pEnd, s, size, dof, offset, stratum;
104426c55118SMichael Lange   PetscInt      *leafStrata, *strataIdx;
1045ad8374ffSToby Isaac   PetscInt     **points;
1046c58f1c22SToby Isaac   char          *name;
1047c58f1c22SToby Isaac   PetscInt       nameSize;
10485cbdf6fcSMichael Lange   PetscHashI     stratumHash;
10495cbdf6fcSMichael Lange   PETSC_UNUSED   PetscHashIIter ret, iter;
1050c58f1c22SToby Isaac   size_t         len = 0;
105126c55118SMichael Lange   PetscMPIInt    rank;
1052c58f1c22SToby Isaac   PetscErrorCode ierr;
1053c58f1c22SToby Isaac 
1054c58f1c22SToby Isaac   PetscFunctionBegin;
1055c58f1c22SToby Isaac   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
1056c58f1c22SToby Isaac   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1057c58f1c22SToby Isaac   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1058c58f1c22SToby Isaac   /* Bcast name */
1059c58f1c22SToby Isaac   if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);}
1060c58f1c22SToby Isaac   nameSize = len;
1061c58f1c22SToby Isaac   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1062c58f1c22SToby Isaac   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1063c58f1c22SToby Isaac   if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);}
1064c58f1c22SToby Isaac   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1065c58f1c22SToby Isaac   ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr);
1066c58f1c22SToby Isaac   ierr = PetscFree(name);CHKERRQ(ierr);
106777d236dfSMichael Lange   /* Bcast defaultValue */
106877d236dfSMichael Lange   if (!rank) (*labelNew)->defaultValue = label->defaultValue;
106977d236dfSMichael Lange   ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
107026c55118SMichael Lange   /* Distribute stratum values over the SF and get the point mapping on the receiver */
107126c55118SMichael Lange   ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr);
10725cbdf6fcSMichael Lange   /* Determine received stratum values and initialise new label*/
10735cbdf6fcSMichael Lange   PetscHashICreate(stratumHash);
107426c55118SMichael Lange   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
10755cbdf6fcSMichael Lange   for (p = 0; p < size; ++p) PetscHashIPut(stratumHash, leafStrata[p], ret, iter);
10765cbdf6fcSMichael Lange   PetscHashISize(stratumHash, (*labelNew)->numStrata);
1077ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr);
1078ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
10795cbdf6fcSMichael Lange   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr);
10805cbdf6fcSMichael Lange   /* Turn leafStrata into indices rather than stratum values */
10815cbdf6fcSMichael Lange   offset = 0;
10825cbdf6fcSMichael Lange   ierr = PetscHashIGetKeys(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr);
10835cbdf6fcSMichael Lange   for (p = 0; p < size; ++p) {
1084231b9e6fSMatthew G. Knepley     for (s = 0; s < (*labelNew)->numStrata; ++s) {
1085231b9e6fSMatthew G. Knepley       if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
10865cbdf6fcSMichael Lange     }
10875cbdf6fcSMichael Lange   }
1088c58f1c22SToby Isaac   /* Rebuild the point strata on the receiver */
1089c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr);
1090c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1091c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1092c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1093c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1094c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1095c58f1c22SToby Isaac       (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1096c58f1c22SToby Isaac     }
1097c58f1c22SToby Isaac   }
1098c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr);
1099c58f1c22SToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr);
1100ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr);
1101c58f1c22SToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1102c58f1c22SToby Isaac     PetscHashICreate((*labelNew)->ht[s]);
1103ad8374ffSToby Isaac     ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr);
1104c58f1c22SToby Isaac   }
1105c58f1c22SToby Isaac   /* Insert points into new strata */
1106c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr);
1107c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1108c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1109c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1110c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1111c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1112c58f1c22SToby Isaac       stratum = leafStrata[offset+s];
1113ad8374ffSToby Isaac       points[stratum][strataIdx[stratum]++] = p;
1114c58f1c22SToby Isaac     }
1115c58f1c22SToby Isaac   }
1116ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; s++) {
1117ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr);
1118ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr);
1119ad8374ffSToby Isaac   }
1120ad8374ffSToby Isaac   ierr = PetscFree(points);CHKERRQ(ierr);
11215cbdf6fcSMichael Lange   PetscHashIDestroy(stratumHash);
1122c58f1c22SToby Isaac   ierr = PetscFree(leafStrata);CHKERRQ(ierr);
1123c58f1c22SToby Isaac   ierr = PetscFree(strataIdx);CHKERRQ(ierr);
1124c58f1c22SToby Isaac   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
1125c58f1c22SToby Isaac   PetscFunctionReturn(0);
1126c58f1c22SToby Isaac }
1127c58f1c22SToby Isaac 
1128c58f1c22SToby Isaac #undef __FUNCT__
11297937d9ceSMichael Lange #define __FUNCT__ "DMLabelGather"
11307937d9ceSMichael Lange /*@
11317937d9ceSMichael Lange   DMLabelGather - Gather all label values from leafs into roots
11327937d9ceSMichael Lange 
11337937d9ceSMichael Lange   Input Parameters:
11347937d9ceSMichael Lange + label - the DMLabel
11357937d9ceSMichael Lange . point - the Star Forest point communication map
11367937d9ceSMichael Lange 
11377937d9ceSMichael Lange   Input Parameters:
11387937d9ceSMichael Lange + label - the new DMLabel with localised leaf values
11397937d9ceSMichael Lange 
11407937d9ceSMichael Lange   Level: developer
11417937d9ceSMichael Lange 
11427937d9ceSMichael Lange   Note: This is the inverse operation to DMLabelDistribute.
11437937d9ceSMichael Lange 
11447937d9ceSMichael Lange .seealso: DMLabelDistribute()
11457937d9ceSMichael Lange @*/
11467937d9ceSMichael Lange PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
11477937d9ceSMichael Lange {
11487937d9ceSMichael Lange   MPI_Comm       comm;
11497937d9ceSMichael Lange   PetscSection   rootSection;
11507937d9ceSMichael Lange   PetscSF        sfLabel;
11517937d9ceSMichael Lange   PetscSFNode   *rootPoints, *leafPoints;
11527937d9ceSMichael Lange   PetscInt       p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
11537937d9ceSMichael Lange   const PetscInt *rootDegree, *ilocal;
11547937d9ceSMichael Lange   PetscInt       *rootStrata;
11557937d9ceSMichael Lange   char          *name;
11567937d9ceSMichael Lange   PetscInt       nameSize;
11577937d9ceSMichael Lange   size_t         len = 0;
11587937d9ceSMichael Lange   PetscMPIInt    rank, numProcs;
11597937d9ceSMichael Lange   PetscErrorCode ierr;
11607937d9ceSMichael Lange 
11617937d9ceSMichael Lange   PetscFunctionBegin;
11627937d9ceSMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
11637937d9ceSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
11647937d9ceSMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
11657937d9ceSMichael Lange   /* Bcast name */
11667937d9ceSMichael Lange   if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);}
11677937d9ceSMichael Lange   nameSize = len;
11687937d9ceSMichael Lange   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
11697937d9ceSMichael Lange   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
11707937d9ceSMichael Lange   if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);}
11717937d9ceSMichael Lange   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
11727937d9ceSMichael Lange   ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr);
11737937d9ceSMichael Lange   ierr = PetscFree(name);CHKERRQ(ierr);
11747937d9ceSMichael Lange   /* Gather rank/index pairs of leaves into local roots to build
11757937d9ceSMichael Lange      an inverse, multi-rooted SF. Note that this ignores local leaf
11767937d9ceSMichael Lange      indexing due to the use of the multiSF in PetscSFGather. */
11777937d9ceSMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr);
1178dc53bc9bSMatthew G. Knepley   ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr);
1179dc53bc9bSMatthew G. Knepley   for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
11807937d9ceSMichael Lange   for (p = 0; p < nleaves; p++) {
1181dc53bc9bSMatthew G. Knepley     leafPoints[ilocal[p]].index = ilocal[p];
1182dc53bc9bSMatthew G. Knepley     leafPoints[ilocal[p]].rank  = rank;
11837937d9ceSMichael Lange   }
11847937d9ceSMichael Lange   ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr);
11857937d9ceSMichael Lange   ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr);
11867937d9ceSMichael Lange   for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
11877937d9ceSMichael Lange   ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr);
11887937d9ceSMichael Lange   ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
11897937d9ceSMichael Lange   ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
11907937d9ceSMichael Lange   ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr);
11917937d9ceSMichael Lange   ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
11927937d9ceSMichael Lange   /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
11937937d9ceSMichael Lange   ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr);
11947937d9ceSMichael Lange   /* Rebuild the point strata on the receiver */
11957937d9ceSMichael Lange   for (p = 0, idx = 0; p < nroots; p++) {
11967937d9ceSMichael Lange     for (d = 0; d < rootDegree[p]; d++) {
11977937d9ceSMichael Lange       ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr);
11987937d9ceSMichael Lange       ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr);
11997937d9ceSMichael Lange       for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);}
12007937d9ceSMichael Lange     }
12017937d9ceSMichael Lange     idx += rootDegree[p];
12027937d9ceSMichael Lange   }
120377e0c0e7SMichael Lange   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
120477e0c0e7SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
120577e0c0e7SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
120677e0c0e7SMichael Lange   ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr);
12077937d9ceSMichael Lange   PetscFunctionReturn(0);
12087937d9ceSMichael Lange }
12097937d9ceSMichael Lange 
12107937d9ceSMichael Lange #undef __FUNCT__
1211c58f1c22SToby Isaac #define __FUNCT__ "DMLabelConvertToSection"
1212c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1213c58f1c22SToby Isaac {
1214c58f1c22SToby Isaac   IS              vIS;
1215c58f1c22SToby Isaac   const PetscInt *values;
1216c58f1c22SToby Isaac   PetscInt       *points;
1217c58f1c22SToby Isaac   PetscInt        nV, vS = 0, vE = 0, v, N;
1218c58f1c22SToby Isaac   PetscErrorCode  ierr;
1219c58f1c22SToby Isaac 
1220c58f1c22SToby Isaac   PetscFunctionBegin;
1221c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr);
1222c58f1c22SToby Isaac   ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
1223c58f1c22SToby Isaac   ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr);
1224c58f1c22SToby Isaac   if (nV) {vS = values[0]; vE = values[0]+1;}
1225c58f1c22SToby Isaac   for (v = 1; v < nV; ++v) {
1226c58f1c22SToby Isaac     vS = PetscMin(vS, values[v]);
1227c58f1c22SToby Isaac     vE = PetscMax(vE, values[v]+1);
1228c58f1c22SToby Isaac   }
1229c58f1c22SToby Isaac   ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr);
1230c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr);
1231c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1232c58f1c22SToby Isaac     PetscInt n;
1233c58f1c22SToby Isaac 
1234c58f1c22SToby Isaac     ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr);
1235c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr);
1236c58f1c22SToby Isaac   }
1237c58f1c22SToby Isaac   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
1238c58f1c22SToby Isaac   ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr);
1239c58f1c22SToby Isaac   ierr = PetscMalloc1(N, &points);CHKERRQ(ierr);
1240c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1241c58f1c22SToby Isaac     IS              is;
1242c58f1c22SToby Isaac     const PetscInt *spoints;
1243c58f1c22SToby Isaac     PetscInt        dof, off, p;
1244c58f1c22SToby Isaac 
1245c58f1c22SToby Isaac     ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr);
1246c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr);
1247c58f1c22SToby Isaac     ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr);
1248c58f1c22SToby Isaac     ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr);
1249c58f1c22SToby Isaac     for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1250c58f1c22SToby Isaac     ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr);
1251c58f1c22SToby Isaac     ierr = ISDestroy(&is);CHKERRQ(ierr);
1252c58f1c22SToby Isaac   }
1253c58f1c22SToby Isaac   ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr);
1254c58f1c22SToby Isaac   ierr = ISDestroy(&vIS);CHKERRQ(ierr);
1255c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1256c58f1c22SToby Isaac   PetscFunctionReturn(0);
1257c58f1c22SToby Isaac }
1258c58f1c22SToby Isaac 
1259c58f1c22SToby Isaac #undef __FUNCT__
1260c58f1c22SToby Isaac #define __FUNCT__ "PetscSectionCreateGlobalSectionLabel"
1261c58f1c22SToby Isaac /*@C
1262c58f1c22SToby Isaac   PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1263c58f1c22SToby Isaac   the local section and an SF describing the section point overlap.
1264c58f1c22SToby Isaac 
1265c58f1c22SToby Isaac   Input Parameters:
1266c58f1c22SToby Isaac   + s - The PetscSection for the local field layout
1267c58f1c22SToby Isaac   . sf - The SF describing parallel layout of the section points
1268c58f1c22SToby Isaac   . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1269c58f1c22SToby Isaac   . label - The label specifying the points
1270c58f1c22SToby Isaac   - labelValue - The label stratum specifying the points
1271c58f1c22SToby Isaac 
1272c58f1c22SToby Isaac   Output Parameter:
1273c58f1c22SToby Isaac   . gsection - The PetscSection for the global field layout
1274c58f1c22SToby Isaac 
1275c58f1c22SToby Isaac   Note: This gives negative sizes and offsets to points not owned by this process
1276c58f1c22SToby Isaac 
1277c58f1c22SToby Isaac   Level: developer
1278c58f1c22SToby Isaac 
1279c58f1c22SToby Isaac .seealso: PetscSectionCreate()
1280c58f1c22SToby Isaac @*/
1281c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1282c58f1c22SToby Isaac {
1283c58f1c22SToby Isaac   PetscInt      *neg = NULL, *tmpOff = NULL;
1284c58f1c22SToby Isaac   PetscInt       pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1285c58f1c22SToby Isaac   PetscErrorCode ierr;
1286c58f1c22SToby Isaac 
1287c58f1c22SToby Isaac   PetscFunctionBegin;
1288c58f1c22SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr);
1289c58f1c22SToby Isaac   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1290c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr);
1291c58f1c22SToby Isaac   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1292c58f1c22SToby Isaac   if (nroots >= 0) {
1293c58f1c22SToby Isaac     if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1294c58f1c22SToby Isaac     ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr);
1295c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1296c58f1c22SToby Isaac       ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr);
1297c58f1c22SToby Isaac     } else {
1298c58f1c22SToby Isaac       tmpOff = &(*gsection)->atlasDof[-pStart];
1299c58f1c22SToby Isaac     }
1300c58f1c22SToby Isaac   }
1301c58f1c22SToby Isaac   /* Mark ghost points with negative dof */
1302c58f1c22SToby Isaac   for (p = pStart; p < pEnd; ++p) {
1303c58f1c22SToby Isaac     PetscInt value;
1304c58f1c22SToby Isaac 
1305c58f1c22SToby Isaac     ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr);
1306c58f1c22SToby Isaac     if (value != labelValue) continue;
1307c58f1c22SToby Isaac     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1308c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr);
1309c58f1c22SToby Isaac     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1310c58f1c22SToby Isaac     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);}
1311c58f1c22SToby Isaac     if (neg) neg[p] = -(dof+1);
1312c58f1c22SToby Isaac   }
1313c58f1c22SToby Isaac   ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr);
1314c58f1c22SToby Isaac   if (nroots >= 0) {
1315c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1316c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1317c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1318c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1319c58f1c22SToby Isaac     }
1320c58f1c22SToby Isaac   }
1321c58f1c22SToby Isaac   /* Calculate new sizes, get proccess offset, and calculate point offsets */
1322c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1323c58f1c22SToby Isaac     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1324c58f1c22SToby Isaac     (*gsection)->atlasOff[p] = off;
1325c58f1c22SToby Isaac     off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1326c58f1c22SToby Isaac   }
1327c58f1c22SToby Isaac   ierr       = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1328c58f1c22SToby Isaac   globalOff -= off;
1329c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1330c58f1c22SToby Isaac     (*gsection)->atlasOff[p] += globalOff;
1331c58f1c22SToby Isaac     if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1332c58f1c22SToby Isaac   }
1333c58f1c22SToby Isaac   /* Put in negative offsets for ghost points */
1334c58f1c22SToby Isaac   if (nroots >= 0) {
1335c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1336c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1337c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1338c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1339c58f1c22SToby Isaac     }
1340c58f1c22SToby Isaac   }
1341c58f1c22SToby Isaac   if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);}
1342c58f1c22SToby Isaac   ierr = PetscFree(neg);CHKERRQ(ierr);
1343c58f1c22SToby Isaac   PetscFunctionReturn(0);
1344c58f1c22SToby Isaac }
1345c58f1c22SToby Isaac 
13465fdea053SToby Isaac typedef struct _n_PetscSectionSym_Label
13475fdea053SToby Isaac {
13485fdea053SToby Isaac   DMLabel           label;
13495fdea053SToby Isaac   PetscCopyMode     *modes;
13505fdea053SToby Isaac   PetscInt          *sizes;
13515fdea053SToby Isaac   const PetscInt    ***perms;
13525fdea053SToby Isaac   const PetscScalar ***rots;
13535fdea053SToby Isaac   PetscInt          (*minMaxOrients)[2];
13545fdea053SToby Isaac   PetscInt          numStrata; /* numStrata is only increasing, functions as a state */
13555fdea053SToby Isaac } PetscSectionSym_Label;
13565fdea053SToby Isaac 
13575fdea053SToby Isaac #undef __FUNCT__
13585fdea053SToby Isaac #define __FUNCT__ "PetscSectionSymLabelReset"
13595fdea053SToby Isaac static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
13605fdea053SToby Isaac {
13615fdea053SToby Isaac   PetscInt              i, j;
13625fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
13635fdea053SToby Isaac   PetscErrorCode        ierr;
13645fdea053SToby Isaac 
13655fdea053SToby Isaac   PetscFunctionBegin;
13665fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
13675fdea053SToby Isaac     if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
13685fdea053SToby Isaac       for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
13695fdea053SToby Isaac         if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);}
13705fdea053SToby Isaac         if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);}
13715fdea053SToby Isaac       }
13725fdea053SToby Isaac       if (sl->perms[i]) {
13735fdea053SToby Isaac         const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
13745fdea053SToby Isaac 
13755fdea053SToby Isaac         ierr = PetscFree(perms);CHKERRQ(ierr);
13765fdea053SToby Isaac       }
13775fdea053SToby Isaac       if (sl->rots[i]) {
13785fdea053SToby Isaac         const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
13795fdea053SToby Isaac 
13805fdea053SToby Isaac         ierr = PetscFree(rots);CHKERRQ(ierr);
13815fdea053SToby Isaac       }
13825fdea053SToby Isaac     }
13835fdea053SToby Isaac   }
13845fdea053SToby Isaac   ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr);
13855fdea053SToby Isaac   ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr);
13865fdea053SToby Isaac   sl->numStrata = 0;
13875fdea053SToby Isaac   PetscFunctionReturn(0);
13885fdea053SToby Isaac }
13895fdea053SToby Isaac 
13905fdea053SToby Isaac #undef __FUNCT__
13915fdea053SToby Isaac #define __FUNCT__ "PetscSectionSymDestroy_Label"
13925fdea053SToby Isaac static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
13935fdea053SToby Isaac {
13945fdea053SToby Isaac   PetscErrorCode ierr;
13955fdea053SToby Isaac 
13965fdea053SToby Isaac   PetscFunctionBegin;
13975fdea053SToby Isaac   ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);
13985fdea053SToby Isaac   ierr = PetscFree(sym->data);CHKERRQ(ierr);
13995fdea053SToby Isaac   PetscFunctionReturn(0);
14005fdea053SToby Isaac }
14015fdea053SToby Isaac 
14025fdea053SToby Isaac #undef __FUNCT__
14035fdea053SToby Isaac #define __FUNCT__ "PetscSectionSymView_Label"
14045fdea053SToby Isaac static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
14055fdea053SToby Isaac {
14065fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
14075fdea053SToby Isaac   PetscBool             isAscii;
14085fdea053SToby Isaac   DMLabel               label = sl->label;
14095fdea053SToby Isaac   PetscErrorCode        ierr;
14105fdea053SToby Isaac 
14115fdea053SToby Isaac   PetscFunctionBegin;
14125fdea053SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
14135fdea053SToby Isaac   if (isAscii) {
14145fdea053SToby Isaac     PetscInt          i, j, k;
14155fdea053SToby Isaac     PetscViewerFormat format;
14165fdea053SToby Isaac 
14175fdea053SToby Isaac     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
14185fdea053SToby Isaac     if (label) {
14195fdea053SToby Isaac       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
14205fdea053SToby Isaac       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
14215fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
14225fdea053SToby Isaac         ierr = DMLabelView(label, viewer);CHKERRQ(ierr);
14235fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
14245fdea053SToby Isaac       } else {
14255fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer,"  Label '%s'\n",sl->label->name);CHKERRQ(ierr);
14265fdea053SToby Isaac       }
14275fdea053SToby Isaac     } else {
14285fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr);
14295fdea053SToby Isaac     }
14305fdea053SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
14315fdea053SToby Isaac     for (i = 0; i <= sl->numStrata; i++) {
14325fdea053SToby Isaac       PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
14335fdea053SToby Isaac 
14345fdea053SToby Isaac       if (!(sl->perms[i] || sl->rots[i])) {
14355fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr);
14365fdea053SToby Isaac       } else {
14375fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr);
14385fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
14395fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr);
14405fdea053SToby Isaac         if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
14415fdea053SToby Isaac           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
14425fdea053SToby Isaac           for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
14435fdea053SToby Isaac             if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
14445fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr);
14455fdea053SToby Isaac             } else {
14465fdea053SToby Isaac               PetscInt tab;
14475fdea053SToby Isaac 
14485fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr);
14495fdea053SToby Isaac               ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
14505fdea053SToby Isaac               ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr);
14515fdea053SToby Isaac               if (sl->perms[i] && sl->perms[i][j]) {
14525fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr);
14535fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
14545fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);}
14555fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
14565fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
14575fdea053SToby Isaac               }
14585fdea053SToby Isaac               if (sl->rots[i] && sl->rots[i][j]) {
14595fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Rotations:  ");CHKERRQ(ierr);
14605fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
14615fdea053SToby Isaac #if defined(PETSC_USE_COMPLEX)
14625fdea053SToby 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);}
14635fdea053SToby Isaac #else
14645fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);}
14655fdea053SToby Isaac #endif
14665fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
14675fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
14685fdea053SToby Isaac               }
14695fdea053SToby Isaac               ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
14705fdea053SToby Isaac             }
14715fdea053SToby Isaac           }
14725fdea053SToby Isaac           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
14735fdea053SToby Isaac         }
14745fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
14755fdea053SToby Isaac       }
14765fdea053SToby Isaac     }
14775fdea053SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
14785fdea053SToby Isaac   }
14795fdea053SToby Isaac   PetscFunctionReturn(0);
14805fdea053SToby Isaac }
14815fdea053SToby Isaac 
14825fdea053SToby Isaac #undef __FUNCT__
14835fdea053SToby Isaac #define __FUNCT__ "PetscSectionSymLabelSetLabel"
14845fdea053SToby Isaac /*@
14855fdea053SToby Isaac   PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
14865fdea053SToby Isaac 
14875fdea053SToby Isaac   Logically collective on sym
14885fdea053SToby Isaac 
14895fdea053SToby Isaac   Input parameters:
14905fdea053SToby Isaac + sym - the section symmetries
14915fdea053SToby Isaac - label - the DMLabel describing the types of points
14925fdea053SToby Isaac 
14935fdea053SToby Isaac   Level: developer:
14945fdea053SToby Isaac 
14955fdea053SToby Isaac .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()
14965fdea053SToby Isaac @*/
14975fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
14985fdea053SToby Isaac {
14995fdea053SToby Isaac   PetscSectionSym_Label *sl;
15005fdea053SToby Isaac   PetscErrorCode ierr;
15015fdea053SToby Isaac 
15025fdea053SToby Isaac   PetscFunctionBegin;
15035fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
15045fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
15055fdea053SToby Isaac   if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);}
15065fdea053SToby Isaac   if (label) {
15075fdea053SToby Isaac     label->refct++;
15085fdea053SToby Isaac     sl->label = label;
15095fdea053SToby Isaac     ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr);
1510*1a834cf9SToby 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);
1511*1a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr);
1512*1a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr);
1513*1a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr);
1514*1a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr);
1515*1a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr);
15165fdea053SToby Isaac   }
15175fdea053SToby Isaac   PetscFunctionReturn(0);
15185fdea053SToby Isaac }
15195fdea053SToby Isaac 
15205fdea053SToby Isaac #undef __FUNCT__
15215fdea053SToby Isaac #define __FUNCT__ "PetscSectionSymLabelSetStratum"
15225fdea053SToby Isaac /*@C
15235fdea053SToby Isaac   PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
15245fdea053SToby Isaac 
15255fdea053SToby Isaac   Logically collective on PetscSectionSym
15265fdea053SToby Isaac 
15275fdea053SToby Isaac   InputParameters:
15285fdea053SToby Isaac + sys       - the section symmetries
15295fdea053SToby Isaac . stratum   - the stratum value in the label that we are assigning symmetries for
15305fdea053SToby Isaac . size      - the number of dofs for points in the stratum of the label
15315fdea053SToby Isaac . minOrient - the smallest orientation for a point in this stratum
15325fdea053SToby Isaac . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
15335fdea053SToby Isaac . mode      - how sym should copy the perms and rots arrays
15345fdea053SToby Isaac . perms     - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation.  A NULL permutation is the identity
15355fdea053SToby 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
15365fdea053SToby Isaac 
15375fdea053SToby Isaac   Level: developer
15385fdea053SToby Isaac 
15395fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()
15405fdea053SToby Isaac @*/
15415fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
15425fdea053SToby Isaac {
15435fdea053SToby Isaac   PetscInt       i, j, k;
15445fdea053SToby Isaac   PetscSectionSym_Label *sl;
15455fdea053SToby Isaac   PetscErrorCode ierr;
15465fdea053SToby Isaac 
15475fdea053SToby Isaac   PetscFunctionBegin;
15485fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
15495fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
15505fdea053SToby Isaac   if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet");
15515fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
15525fdea053SToby Isaac     PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
15535fdea053SToby Isaac 
15545fdea053SToby Isaac     if (stratum == value) break;
15555fdea053SToby Isaac   }
15565fdea053SToby 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);
15575fdea053SToby Isaac   sl->sizes[i] = size;
15585fdea053SToby Isaac   sl->modes[i] = mode;
15595fdea053SToby Isaac   sl->minMaxOrients[i][0] = minOrient;
15605fdea053SToby Isaac   sl->minMaxOrients[i][1] = maxOrient;
15615fdea053SToby Isaac   if (mode == PETSC_COPY_VALUES) {
15625fdea053SToby Isaac     if (perms) {
15635fdea053SToby Isaac       PetscInt    **ownPerms;
15645fdea053SToby Isaac 
15655fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr);
15665fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
15675fdea053SToby Isaac         if (perms[j]) {
15685fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr);
15695fdea053SToby Isaac           for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];}
15705fdea053SToby Isaac         }
15715fdea053SToby Isaac       }
15725fdea053SToby Isaac       sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient];
15735fdea053SToby Isaac     }
15745fdea053SToby Isaac     if (rots) {
15755fdea053SToby Isaac       PetscScalar **ownRots;
15765fdea053SToby Isaac 
15775fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr);
15785fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
15795fdea053SToby Isaac         if (rots[j]) {
15805fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr);
15815fdea053SToby Isaac           for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];}
15825fdea053SToby Isaac         }
15835fdea053SToby Isaac       }
15845fdea053SToby Isaac       sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient];
15855fdea053SToby Isaac     }
15865fdea053SToby Isaac   } else {
15875fdea053SToby Isaac     sl->perms[i] = perms ? &perms[-minOrient] : NULL;
15885fdea053SToby Isaac     sl->rots[i]  = rots ? &rots[-minOrient] : NULL;
15895fdea053SToby Isaac   }
15905fdea053SToby Isaac   PetscFunctionReturn(0);
15915fdea053SToby Isaac }
15925fdea053SToby Isaac 
15935fdea053SToby Isaac #undef __FUNCT__
15945fdea053SToby Isaac #define __FUNCT__ "PetscSectionSymGetPoints_Label"
15955fdea053SToby Isaac static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
15965fdea053SToby Isaac {
15975fdea053SToby Isaac   PetscInt              i, j, numStrata;
15985fdea053SToby Isaac   PetscSectionSym_Label *sl;
15995fdea053SToby Isaac   DMLabel               label;
16005fdea053SToby Isaac   PetscErrorCode        ierr;
16015fdea053SToby Isaac 
16025fdea053SToby Isaac   PetscFunctionBegin;
16035fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
16045fdea053SToby Isaac   numStrata = sl->numStrata;
16055fdea053SToby Isaac   label     = sl->label;
16065fdea053SToby Isaac   for (i = 0; i < numPoints; i++) {
16075fdea053SToby Isaac     PetscInt point = points[2*i];
16085fdea053SToby Isaac     PetscInt ornt  = points[2*i+1];
16095fdea053SToby Isaac 
16105fdea053SToby Isaac     for (j = 0; j < numStrata; j++) {
16115fdea053SToby Isaac       if (label->validIS[j]) {
16125fdea053SToby Isaac         PetscInt k;
16135fdea053SToby Isaac 
16145fdea053SToby Isaac         ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr);
16155fdea053SToby Isaac         if (k >= 0) break;
16165fdea053SToby Isaac       } else {
16175fdea053SToby Isaac         PetscBool has;
16185fdea053SToby Isaac 
16195fdea053SToby Isaac         PetscHashIHasKey(label->ht[j], point, has);
16205fdea053SToby Isaac         if (has) break;
16215fdea053SToby Isaac       }
16225fdea053SToby Isaac     }
16235fdea053SToby 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);
16245fdea053SToby Isaac     if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;}
16255fdea053SToby Isaac     if (rots) {rots[i]  = sl->rots[j] ? sl->rots[j][ornt] : NULL;}
16265fdea053SToby Isaac   }
16275fdea053SToby Isaac   PetscFunctionReturn(0);
16285fdea053SToby Isaac }
16295fdea053SToby Isaac 
16305fdea053SToby Isaac #undef __FUNCT__
16315fdea053SToby Isaac #define __FUNCT__ "PetscSectionSymCreate_Label"
16325fdea053SToby Isaac PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
16335fdea053SToby Isaac {
16345fdea053SToby Isaac   PetscSectionSym_Label *sl;
16355fdea053SToby Isaac   PetscErrorCode        ierr;
16365fdea053SToby Isaac 
16375fdea053SToby Isaac   PetscFunctionBegin;
16385fdea053SToby Isaac   ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr);
16395fdea053SToby Isaac   sym->ops->getpoints = PetscSectionSymGetPoints_Label;
16405fdea053SToby Isaac   sym->ops->view      = PetscSectionSymView_Label;
16415fdea053SToby Isaac   sym->ops->destroy   = PetscSectionSymDestroy_Label;
16425fdea053SToby Isaac   sym->data           = (void *) sl;
16435fdea053SToby Isaac   PetscFunctionReturn(0);
16445fdea053SToby Isaac }
16455fdea053SToby Isaac 
16465fdea053SToby Isaac #undef __FUNCT__
16475fdea053SToby Isaac #define __FUNCT__ "PetscSectionSymCreateLabel"
16485fdea053SToby Isaac /*@
16495fdea053SToby Isaac   PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
16505fdea053SToby Isaac 
16515fdea053SToby Isaac   Collective
16525fdea053SToby Isaac 
16535fdea053SToby Isaac   Input Parameters:
16545fdea053SToby Isaac + comm - the MPI communicator for the new symmetry
16555fdea053SToby Isaac - label - the label defining the strata
16565fdea053SToby Isaac 
16575fdea053SToby Isaac   Output Parameters:
16585fdea053SToby Isaac . sym - the section symmetries
16595fdea053SToby Isaac 
16605fdea053SToby Isaac   Level: developer
16615fdea053SToby Isaac 
16625fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()
16635fdea053SToby Isaac @*/
16645fdea053SToby Isaac PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
16655fdea053SToby Isaac {
16665fdea053SToby Isaac   PetscErrorCode        ierr;
16675fdea053SToby Isaac 
16685fdea053SToby Isaac   PetscFunctionBegin;
16695fdea053SToby Isaac   ierr = DMInitializePackage();CHKERRQ(ierr);
16705fdea053SToby Isaac   ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr);
16715fdea053SToby Isaac   ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr);
16725fdea053SToby Isaac   ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr);
16735fdea053SToby Isaac   PetscFunctionReturn(0);
16745fdea053SToby Isaac }
1675