xref: /petsc/src/dm/label/dmlabel.c (revision fa8e8ae5f164b8f4ad3a50332e068fc0b5da93c6)
15fdea053SToby Isaac #include <petscdm.h>
2c58f1c22SToby Isaac #include <petsc/private/dmlabelimpl.h>   /*I      "petscdmlabel.h"   I*/
3c58f1c22SToby Isaac #include <petsc/private/isimpl.h>        /*I      "petscis.h"        I*/
4c58f1c22SToby Isaac #include <petscsf.h>
5c58f1c22SToby Isaac 
6c58f1c22SToby Isaac /*@C
7c58f1c22SToby Isaac   DMLabelCreate - Create a DMLabel object, which is a multimap
8c58f1c22SToby Isaac 
9c58f1c22SToby Isaac   Input parameter:
10c58f1c22SToby Isaac . name - The label name
11c58f1c22SToby Isaac 
12c58f1c22SToby Isaac   Output parameter:
13c58f1c22SToby Isaac . label - The DMLabel
14c58f1c22SToby Isaac 
15c58f1c22SToby Isaac   Level: beginner
16c58f1c22SToby Isaac 
17c58f1c22SToby Isaac .seealso: DMLabelDestroy()
18c58f1c22SToby Isaac @*/
19c58f1c22SToby Isaac PetscErrorCode DMLabelCreate(const char name[], DMLabel *label)
20c58f1c22SToby Isaac {
21c58f1c22SToby Isaac   PetscErrorCode ierr;
22c58f1c22SToby Isaac 
23c58f1c22SToby Isaac   PetscFunctionBegin;
24c58f1c22SToby Isaac   ierr = PetscNew(label);CHKERRQ(ierr);
25c58f1c22SToby Isaac   ierr = PetscStrallocpy(name, &(*label)->name);CHKERRQ(ierr);
26c58f1c22SToby Isaac 
27c58f1c22SToby Isaac   (*label)->refct          = 1;
28c58f1c22SToby Isaac   (*label)->state          = -1;
29c58f1c22SToby Isaac   (*label)->numStrata      = 0;
305aa44df4SToby Isaac   (*label)->defaultValue   = -1;
31c58f1c22SToby Isaac   (*label)->stratumValues  = NULL;
32ad8374ffSToby Isaac   (*label)->validIS        = NULL;
33c58f1c22SToby Isaac   (*label)->stratumSizes   = NULL;
34c58f1c22SToby Isaac   (*label)->points         = NULL;
35c58f1c22SToby Isaac   (*label)->ht             = NULL;
36c58f1c22SToby Isaac   (*label)->pStart         = -1;
37c58f1c22SToby Isaac   (*label)->pEnd           = -1;
38c58f1c22SToby Isaac   (*label)->bt             = NULL;
39c58f1c22SToby Isaac   PetscFunctionReturn(0);
40c58f1c22SToby Isaac }
41c58f1c22SToby Isaac 
42c58f1c22SToby Isaac /*
43c58f1c22SToby Isaac   DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format
44c58f1c22SToby Isaac 
45c58f1c22SToby Isaac   Input parameter:
46c58f1c22SToby Isaac + label - The DMLabel
47c58f1c22SToby Isaac - v - The stratum value
48c58f1c22SToby Isaac 
49c58f1c22SToby Isaac   Output parameter:
50c58f1c22SToby Isaac . label - The DMLabel with stratum in sorted list format
51c58f1c22SToby Isaac 
52c58f1c22SToby Isaac   Level: developer
53c58f1c22SToby Isaac 
54c58f1c22SToby Isaac .seealso: DMLabelCreate()
55c58f1c22SToby Isaac */
56c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v)
57c58f1c22SToby Isaac {
58c58f1c22SToby Isaac   PetscInt       off;
59ad8374ffSToby Isaac   PetscInt       *pointArray;
60c58f1c22SToby Isaac   PetscErrorCode ierr;
61c58f1c22SToby Isaac 
62ad8374ffSToby Isaac   if (label->validIS[v]) return 0;
63c58f1c22SToby Isaac   if (v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeValid_Private\n", v);
64c58f1c22SToby Isaac   PetscFunctionBegin;
65c58f1c22SToby Isaac   PetscHashISize(label->ht[v], label->stratumSizes[v]);
66c58f1c22SToby Isaac 
67ad8374ffSToby Isaac   ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr);
68c58f1c22SToby Isaac   off  = 0;
69ad8374ffSToby Isaac   ierr = PetscHashIGetKeys(label->ht[v], &off, pointArray);CHKERRQ(ierr);
70c58f1c22SToby Isaac   if (off != label->stratumSizes[v]) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of contributed points %D from value %D should be %D", off, label->stratumValues[v], label->stratumSizes[v]);
71c58f1c22SToby Isaac   PetscHashIClear(label->ht[v]);
72ad8374ffSToby Isaac   ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr);
73c58f1c22SToby Isaac   if (label->bt) {
74c58f1c22SToby Isaac     PetscInt p;
75c58f1c22SToby Isaac 
76c58f1c22SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
77ad8374ffSToby Isaac       const PetscInt point = pointArray[p];
78c58f1c22SToby Isaac 
79c58f1c22SToby Isaac       if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, label->pStart, label->pEnd);
80c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
81c58f1c22SToby Isaac     }
82c58f1c22SToby Isaac   }
83ad8374ffSToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF,label->stratumSizes[v],pointArray,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
84ad8374ffSToby Isaac   ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
85ad8374ffSToby Isaac   label->validIS[v] = PETSC_TRUE;
86c58f1c22SToby Isaac   ++label->state;
87c58f1c22SToby Isaac   PetscFunctionReturn(0);
88c58f1c22SToby Isaac }
89c58f1c22SToby Isaac 
90c58f1c22SToby Isaac /*
91c58f1c22SToby Isaac   DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format
92c58f1c22SToby Isaac 
93c58f1c22SToby Isaac   Input parameter:
94c58f1c22SToby Isaac . label - The DMLabel
95c58f1c22SToby Isaac 
96c58f1c22SToby Isaac   Output parameter:
97c58f1c22SToby Isaac . label - The DMLabel with all strata in sorted list format
98c58f1c22SToby Isaac 
99c58f1c22SToby Isaac   Level: developer
100c58f1c22SToby Isaac 
101c58f1c22SToby Isaac .seealso: DMLabelCreate()
102c58f1c22SToby Isaac */
103c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label)
104c58f1c22SToby Isaac {
105c58f1c22SToby Isaac   PetscInt       v;
106c58f1c22SToby Isaac   PetscErrorCode ierr;
107c58f1c22SToby Isaac 
108c58f1c22SToby Isaac   PetscFunctionBegin;
109c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; v++){
110c58f1c22SToby Isaac     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
111c58f1c22SToby Isaac   }
112c58f1c22SToby Isaac   PetscFunctionReturn(0);
113c58f1c22SToby Isaac }
114c58f1c22SToby Isaac 
115c58f1c22SToby Isaac /*
116c58f1c22SToby Isaac   DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format
117c58f1c22SToby Isaac 
118c58f1c22SToby Isaac   Input parameter:
119c58f1c22SToby Isaac + label - The DMLabel
120c58f1c22SToby Isaac - v - The stratum value
121c58f1c22SToby Isaac 
122c58f1c22SToby Isaac   Output parameter:
123c58f1c22SToby Isaac . label - The DMLabel with stratum in hash format
124c58f1c22SToby Isaac 
125c58f1c22SToby Isaac   Level: developer
126c58f1c22SToby Isaac 
127c58f1c22SToby Isaac .seealso: DMLabelCreate()
128c58f1c22SToby Isaac */
129c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v)
130c58f1c22SToby Isaac {
131c58f1c22SToby Isaac   PETSC_UNUSED PetscHashIIter ret, iter;
132c58f1c22SToby Isaac   PetscInt                    p;
133ad8374ffSToby Isaac   const PetscInt              *points;
134c58f1c22SToby Isaac   PetscErrorCode ierr;
135c58f1c22SToby Isaac 
136c58f1c22SToby Isaac   PetscFunctionBegin;
137ad8374ffSToby Isaac   if (!label->validIS[v]) PetscFunctionReturn(0);
138ad8374ffSToby Isaac   if (label->points[v]) {
139ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
140ad8374ffSToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) PetscHashIPut(label->ht[v], points[p], ret, iter);
141ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
142ad8374ffSToby Isaac     ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
143ad8374ffSToby Isaac   }
144ad8374ffSToby Isaac   label->validIS[v] = PETSC_FALSE;
145c58f1c22SToby Isaac   PetscFunctionReturn(0);
146c58f1c22SToby Isaac }
147c58f1c22SToby Isaac 
148c58f1c22SToby Isaac PetscErrorCode DMLabelGetState(DMLabel label, PetscObjectState *state)
149c58f1c22SToby Isaac {
150c58f1c22SToby Isaac   PetscFunctionBegin;
151c58f1c22SToby Isaac   PetscValidPointer(state, 2);
152c58f1c22SToby Isaac   *state = label->state;
153c58f1c22SToby Isaac   PetscFunctionReturn(0);
154c58f1c22SToby Isaac }
155c58f1c22SToby Isaac 
156c58f1c22SToby Isaac PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
157c58f1c22SToby Isaac {
158ad8374ffSToby Isaac   PetscInt    v, *tmpV, *tmpS;
159ad8374ffSToby Isaac   IS         *tmpP;
160c58f1c22SToby Isaac   PetscHashI *tmpH;
161c58f1c22SToby Isaac   PetscBool  *tmpB;
162c58f1c22SToby Isaac   PetscErrorCode ierr;
163c58f1c22SToby Isaac 
164c58f1c22SToby Isaac   PetscFunctionBegin;
165c58f1c22SToby Isaac 
166ad8374ffSToby Isaac   for (v = 0; v < label->numStrata; v++) {
167ad8374ffSToby Isaac     if (label->stratumValues[v] == value) PetscFunctionReturn(0);
168ad8374ffSToby Isaac   }
169c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpV);CHKERRQ(ierr);
170c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpS);CHKERRQ(ierr);
171c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpH);CHKERRQ(ierr);
172c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpP);CHKERRQ(ierr);
173c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpB);CHKERRQ(ierr);
174c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
175c58f1c22SToby Isaac     tmpV[v] = label->stratumValues[v];
176c58f1c22SToby Isaac     tmpS[v] = label->stratumSizes[v];
177c58f1c22SToby Isaac     tmpH[v] = label->ht[v];
178c58f1c22SToby Isaac     tmpP[v] = label->points[v];
179ad8374ffSToby Isaac     tmpB[v] = label->validIS[v];
180c58f1c22SToby Isaac   }
181c58f1c22SToby Isaac   tmpV[v] = value;
182c58f1c22SToby Isaac   tmpS[v] = 0;
183c58f1c22SToby Isaac   PetscHashICreate(tmpH[v]);
184ad8374ffSToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&tmpP[v]);CHKERRQ(ierr);
185c58f1c22SToby Isaac   tmpB[v] = PETSC_TRUE;
186c58f1c22SToby Isaac   ++label->numStrata;
187c58f1c22SToby Isaac   ierr = PetscFree(label->stratumValues);CHKERRQ(ierr);
188c58f1c22SToby Isaac   ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr);
189c58f1c22SToby Isaac   ierr = PetscFree(label->ht);CHKERRQ(ierr);
190c58f1c22SToby Isaac   ierr = PetscFree(label->points);CHKERRQ(ierr);
191ad8374ffSToby Isaac   ierr = PetscFree(label->validIS);CHKERRQ(ierr);
192c58f1c22SToby Isaac   label->stratumValues = tmpV;
193c58f1c22SToby Isaac   label->stratumSizes  = tmpS;
194c58f1c22SToby Isaac   label->ht            = tmpH;
195c58f1c22SToby Isaac   label->points        = tmpP;
196ad8374ffSToby Isaac   label->validIS       = tmpB;
197c58f1c22SToby Isaac 
198c58f1c22SToby Isaac   PetscFunctionReturn(0);
199c58f1c22SToby Isaac }
200c58f1c22SToby Isaac 
201c58f1c22SToby Isaac /*@C
202c58f1c22SToby Isaac   DMLabelGetName - Return the name of a DMLabel object
203c58f1c22SToby Isaac 
204c58f1c22SToby Isaac   Input parameter:
205c58f1c22SToby Isaac . label - The DMLabel
206c58f1c22SToby Isaac 
207c58f1c22SToby Isaac   Output parameter:
208c58f1c22SToby Isaac . name - The label name
209c58f1c22SToby Isaac 
210c58f1c22SToby Isaac   Level: beginner
211c58f1c22SToby Isaac 
212c58f1c22SToby Isaac .seealso: DMLabelCreate()
213c58f1c22SToby Isaac @*/
214c58f1c22SToby Isaac PetscErrorCode DMLabelGetName(DMLabel label, const char **name)
215c58f1c22SToby Isaac {
216c58f1c22SToby Isaac   PetscFunctionBegin;
217c58f1c22SToby Isaac   PetscValidPointer(name, 2);
218c58f1c22SToby Isaac   *name = label->name;
219c58f1c22SToby Isaac   PetscFunctionReturn(0);
220c58f1c22SToby Isaac }
221c58f1c22SToby Isaac 
222c58f1c22SToby Isaac static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
223c58f1c22SToby Isaac {
224c58f1c22SToby Isaac   PetscInt       v;
225c58f1c22SToby Isaac   PetscMPIInt    rank;
226c58f1c22SToby Isaac   PetscErrorCode ierr;
227c58f1c22SToby Isaac 
228c58f1c22SToby Isaac   PetscFunctionBegin;
229c58f1c22SToby Isaac   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr);
230c58f1c22SToby Isaac   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
231c58f1c22SToby Isaac   if (label) {
232c58f1c22SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", label->name);CHKERRQ(ierr);
233c58f1c22SToby Isaac     if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);}
234c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
235c58f1c22SToby Isaac       const PetscInt value = label->stratumValues[v];
236ad8374ffSToby Isaac       const PetscInt *points;
237c58f1c22SToby Isaac       PetscInt       p;
238c58f1c22SToby Isaac 
239ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
240c58f1c22SToby Isaac       for (p = 0; p < label->stratumSizes[v]; ++p) {
241ad8374ffSToby Isaac         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr);
242c58f1c22SToby Isaac       }
243ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
244c58f1c22SToby Isaac     }
245c58f1c22SToby Isaac   }
246c58f1c22SToby Isaac   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
247c58f1c22SToby Isaac   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
248c58f1c22SToby Isaac   PetscFunctionReturn(0);
249c58f1c22SToby Isaac }
250c58f1c22SToby Isaac 
251c58f1c22SToby Isaac /*@C
252c58f1c22SToby Isaac   DMLabelView - View the label
253c58f1c22SToby Isaac 
254c58f1c22SToby Isaac   Input Parameters:
255c58f1c22SToby Isaac + label - The DMLabel
256c58f1c22SToby Isaac - viewer - The PetscViewer
257c58f1c22SToby Isaac 
258c58f1c22SToby Isaac   Level: intermediate
259c58f1c22SToby Isaac 
260c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy()
261c58f1c22SToby Isaac @*/
262c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
263c58f1c22SToby Isaac {
264c58f1c22SToby Isaac   PetscBool      iascii;
265c58f1c22SToby Isaac   PetscErrorCode ierr;
266c58f1c22SToby Isaac 
267c58f1c22SToby Isaac   PetscFunctionBegin;
268c58f1c22SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
269c58f1c22SToby Isaac   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
270c58f1c22SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
271c58f1c22SToby Isaac   if (iascii) {
272c58f1c22SToby Isaac     ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr);
273c58f1c22SToby Isaac   }
274c58f1c22SToby Isaac   PetscFunctionReturn(0);
275c58f1c22SToby Isaac }
276c58f1c22SToby Isaac 
27784f0b6dfSMatthew G. Knepley /*@
27884f0b6dfSMatthew G. Knepley   DMLabelDestroy - Destroys a DMLabel
27984f0b6dfSMatthew G. Knepley 
28084f0b6dfSMatthew G. Knepley   Input Parameter:
28184f0b6dfSMatthew G. Knepley . label - The DMLabel
28284f0b6dfSMatthew G. Knepley 
28384f0b6dfSMatthew G. Knepley   Level: beginner
28484f0b6dfSMatthew G. Knepley 
28584f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate()
28684f0b6dfSMatthew G. Knepley @*/
287c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label)
288c58f1c22SToby Isaac {
289c58f1c22SToby Isaac   PetscInt       v;
290c58f1c22SToby Isaac   PetscErrorCode ierr;
291c58f1c22SToby Isaac 
292c58f1c22SToby Isaac   PetscFunctionBegin;
293c58f1c22SToby Isaac   if (!(*label)) PetscFunctionReturn(0);
294c58f1c22SToby Isaac   if (--(*label)->refct > 0) PetscFunctionReturn(0);
295c58f1c22SToby Isaac   ierr = PetscFree((*label)->name);CHKERRQ(ierr);
296c58f1c22SToby Isaac   ierr = PetscFree((*label)->stratumValues);CHKERRQ(ierr);
297c58f1c22SToby Isaac   ierr = PetscFree((*label)->stratumSizes);CHKERRQ(ierr);
298ad8374ffSToby Isaac   for (v = 0; v < (*label)->numStrata; ++v) {ierr = ISDestroy(&((*label)->points[v]));CHKERRQ(ierr);}
299c58f1c22SToby Isaac   ierr = PetscFree((*label)->points);CHKERRQ(ierr);
300ad8374ffSToby Isaac   ierr = PetscFree((*label)->validIS);CHKERRQ(ierr);
301c58f1c22SToby Isaac   if ((*label)->ht) {
302c58f1c22SToby Isaac     for (v = 0; v < (*label)->numStrata; ++v) {PetscHashIDestroy((*label)->ht[v]);}
303c58f1c22SToby Isaac     ierr = PetscFree((*label)->ht);CHKERRQ(ierr);
304c58f1c22SToby Isaac   }
305c58f1c22SToby Isaac   ierr = PetscBTDestroy(&(*label)->bt);CHKERRQ(ierr);
306c58f1c22SToby Isaac   ierr = PetscFree(*label);CHKERRQ(ierr);
307c58f1c22SToby Isaac   PetscFunctionReturn(0);
308c58f1c22SToby Isaac }
309c58f1c22SToby Isaac 
31084f0b6dfSMatthew G. Knepley /*@
31184f0b6dfSMatthew G. Knepley   DMLabelDuplicate - Duplicates a DMLabel
31284f0b6dfSMatthew G. Knepley 
31384f0b6dfSMatthew G. Knepley   Input Parameter:
31484f0b6dfSMatthew G. Knepley . label - The DMLabel
31584f0b6dfSMatthew G. Knepley 
31684f0b6dfSMatthew G. Knepley   Output Parameter:
31784f0b6dfSMatthew G. Knepley . labelnew - location to put new vector
31884f0b6dfSMatthew G. Knepley 
31984f0b6dfSMatthew G. Knepley   Level: intermediate
32084f0b6dfSMatthew G. Knepley 
32184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelDestroy()
32284f0b6dfSMatthew G. Knepley @*/
323c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
324c58f1c22SToby Isaac {
325ad8374ffSToby Isaac   PetscInt       v;
326c58f1c22SToby Isaac   PetscErrorCode ierr;
327c58f1c22SToby Isaac 
328c58f1c22SToby Isaac   PetscFunctionBegin;
329c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
330c58f1c22SToby Isaac   ierr = PetscNew(labelnew);CHKERRQ(ierr);
331c58f1c22SToby Isaac   ierr = PetscStrallocpy(label->name, &(*labelnew)->name);CHKERRQ(ierr);
332c58f1c22SToby Isaac 
333c58f1c22SToby Isaac   (*labelnew)->refct        = 1;
334c58f1c22SToby Isaac   (*labelnew)->numStrata    = label->numStrata;
3355aa44df4SToby Isaac   (*labelnew)->defaultValue = label->defaultValue;
336c58f1c22SToby Isaac   if (label->numStrata) {
337c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr);
338c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr);
339c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr);
340c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr);
341ad8374ffSToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr);
342c58f1c22SToby Isaac     /* Could eliminate unused space here */
343c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
344c58f1c22SToby Isaac       PetscHashICreate((*labelnew)->ht[v]);
345ad8374ffSToby Isaac       (*labelnew)->validIS[v]        = PETSC_TRUE;
346c58f1c22SToby Isaac       (*labelnew)->stratumValues[v]  = label->stratumValues[v];
347c58f1c22SToby Isaac       (*labelnew)->stratumSizes[v]   = label->stratumSizes[v];
348ad8374ffSToby Isaac       ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr);
349ad8374ffSToby Isaac       (*labelnew)->points[v]         = label->points[v];
350c58f1c22SToby Isaac     }
351c58f1c22SToby Isaac   }
352c58f1c22SToby Isaac   (*labelnew)->pStart = -1;
353c58f1c22SToby Isaac   (*labelnew)->pEnd   = -1;
354c58f1c22SToby Isaac   (*labelnew)->bt     = NULL;
355c58f1c22SToby Isaac   PetscFunctionReturn(0);
356c58f1c22SToby Isaac }
357c58f1c22SToby Isaac 
358c58f1c22SToby Isaac /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
359c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
360c58f1c22SToby Isaac {
361c58f1c22SToby Isaac   PetscInt       v;
362c58f1c22SToby Isaac   PetscErrorCode ierr;
363c58f1c22SToby Isaac 
364c58f1c22SToby Isaac   PetscFunctionBegin;
365c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
366c58f1c22SToby Isaac   if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);}
367c58f1c22SToby Isaac   label->pStart = pStart;
368c58f1c22SToby Isaac   label->pEnd   = pEnd;
369c58f1c22SToby Isaac   ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr);
370c58f1c22SToby Isaac   ierr = PetscBTMemzero(pEnd - pStart, label->bt);CHKERRQ(ierr);
371c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
372ad8374ffSToby Isaac     const PetscInt *points;
373c58f1c22SToby Isaac     PetscInt       i;
374c58f1c22SToby Isaac 
375ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
376c58f1c22SToby Isaac     for (i = 0; i < label->stratumSizes[v]; ++i) {
377ad8374ffSToby Isaac       const PetscInt point = points[i];
378c58f1c22SToby Isaac 
379c58f1c22SToby 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);
380c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr);
381c58f1c22SToby Isaac     }
382ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
383c58f1c22SToby Isaac   }
384c58f1c22SToby Isaac   PetscFunctionReturn(0);
385c58f1c22SToby Isaac }
386c58f1c22SToby Isaac 
387c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label)
388c58f1c22SToby Isaac {
389c58f1c22SToby Isaac   PetscErrorCode ierr;
390c58f1c22SToby Isaac 
391c58f1c22SToby Isaac   PetscFunctionBegin;
392c58f1c22SToby Isaac   label->pStart = -1;
393c58f1c22SToby Isaac   label->pEnd   = -1;
394c58f1c22SToby Isaac   if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);}
395c58f1c22SToby Isaac   PetscFunctionReturn(0);
396c58f1c22SToby Isaac }
397c58f1c22SToby Isaac 
398c58f1c22SToby Isaac /*@
399c58f1c22SToby Isaac   DMLabelHasValue - Determine whether a label assigns the value to any point
400c58f1c22SToby Isaac 
401c58f1c22SToby Isaac   Input Parameters:
402c58f1c22SToby Isaac + label - the DMLabel
403c58f1c22SToby Isaac - value - the value
404c58f1c22SToby Isaac 
405c58f1c22SToby Isaac   Output Parameter:
406c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point
407c58f1c22SToby Isaac 
408c58f1c22SToby Isaac   Level: developer
409c58f1c22SToby Isaac 
410c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
411c58f1c22SToby Isaac @*/
412c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
413c58f1c22SToby Isaac {
414c58f1c22SToby Isaac   PetscInt v;
415c58f1c22SToby Isaac 
416c58f1c22SToby Isaac   PetscFunctionBegin;
417c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
418c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
419c58f1c22SToby Isaac     if (value == label->stratumValues[v]) break;
420c58f1c22SToby Isaac   }
421c58f1c22SToby Isaac   *contains = (v < label->numStrata ? PETSC_TRUE : PETSC_FALSE);
422c58f1c22SToby Isaac   PetscFunctionReturn(0);
423c58f1c22SToby Isaac }
424c58f1c22SToby Isaac 
425c58f1c22SToby Isaac /*@
426c58f1c22SToby Isaac   DMLabelHasPoint - Determine whether a label assigns a value to a point
427c58f1c22SToby Isaac 
428c58f1c22SToby Isaac   Input Parameters:
429c58f1c22SToby Isaac + label - the DMLabel
430c58f1c22SToby Isaac - point - the point
431c58f1c22SToby Isaac 
432c58f1c22SToby Isaac   Output Parameter:
433c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value
434c58f1c22SToby Isaac 
435c58f1c22SToby Isaac   Note: The user must call DMLabelCreateIndex() before this function.
436c58f1c22SToby Isaac 
437c58f1c22SToby Isaac   Level: developer
438c58f1c22SToby Isaac 
439c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
440c58f1c22SToby Isaac @*/
441c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
442c58f1c22SToby Isaac {
443c58f1c22SToby Isaac   PetscErrorCode ierr;
444c58f1c22SToby Isaac 
445c58f1c22SToby Isaac   PetscFunctionBeginHot;
446c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
447c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
448c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG)
449c58f1c22SToby Isaac   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
450c58f1c22SToby 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);
451c58f1c22SToby Isaac #endif
452c58f1c22SToby Isaac   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
453c58f1c22SToby Isaac   PetscFunctionReturn(0);
454c58f1c22SToby Isaac }
455c58f1c22SToby Isaac 
456c58f1c22SToby Isaac /*@
457c58f1c22SToby Isaac   DMLabelStratumHasPoint - Return true if the stratum contains a point
458c58f1c22SToby Isaac 
459c58f1c22SToby Isaac   Input Parameters:
460c58f1c22SToby Isaac + label - the DMLabel
461c58f1c22SToby Isaac . value - the stratum value
462c58f1c22SToby Isaac - point - the point
463c58f1c22SToby Isaac 
464c58f1c22SToby Isaac   Output Parameter:
465c58f1c22SToby Isaac . contains - true if the stratum contains the point
466c58f1c22SToby Isaac 
467c58f1c22SToby Isaac   Level: intermediate
468c58f1c22SToby Isaac 
469c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
470c58f1c22SToby Isaac @*/
471c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
472c58f1c22SToby Isaac {
473c58f1c22SToby Isaac   PetscInt       v;
474c58f1c22SToby Isaac   PetscErrorCode ierr;
475c58f1c22SToby Isaac 
476c58f1c22SToby Isaac   PetscFunctionBegin;
477c58f1c22SToby Isaac   PetscValidPointer(contains, 4);
478c58f1c22SToby Isaac   *contains = PETSC_FALSE;
479c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
480c58f1c22SToby Isaac     if (label->stratumValues[v] == value) {
481ad8374ffSToby Isaac       if (label->validIS[v]) {
482c58f1c22SToby Isaac         PetscInt i;
483c58f1c22SToby Isaac 
484a2d74346SToby Isaac         ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr);
485c58f1c22SToby Isaac         if (i >= 0) {
486c58f1c22SToby Isaac           *contains = PETSC_TRUE;
487c58f1c22SToby Isaac           break;
488c58f1c22SToby Isaac         }
489c58f1c22SToby Isaac       } else {
490c58f1c22SToby Isaac         PetscBool has;
491c58f1c22SToby Isaac 
492c58f1c22SToby Isaac         PetscHashIHasKey(label->ht[v], point, has);
493c58f1c22SToby Isaac         if (has) {
494c58f1c22SToby Isaac           *contains = PETSC_TRUE;
495c58f1c22SToby Isaac           break;
496c58f1c22SToby Isaac         }
497c58f1c22SToby Isaac       }
498c58f1c22SToby Isaac     }
499c58f1c22SToby Isaac   }
500c58f1c22SToby Isaac   PetscFunctionReturn(0);
501c58f1c22SToby Isaac }
502c58f1c22SToby Isaac 
50384f0b6dfSMatthew G. Knepley /*@
5045aa44df4SToby Isaac   DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
5055aa44df4SToby Isaac   When a label is created, it is initialized to -1.
5065aa44df4SToby Isaac 
5075aa44df4SToby Isaac   Input parameter:
5085aa44df4SToby Isaac . label - a DMLabel object
5095aa44df4SToby Isaac 
5105aa44df4SToby Isaac   Output parameter:
5115aa44df4SToby Isaac . defaultValue - the default value
5125aa44df4SToby Isaac 
5135aa44df4SToby Isaac   Level: beginner
5145aa44df4SToby Isaac 
5155aa44df4SToby Isaac .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
51684f0b6dfSMatthew G. Knepley @*/
5175aa44df4SToby Isaac PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
5185aa44df4SToby Isaac {
5195aa44df4SToby Isaac   PetscFunctionBegin;
5205aa44df4SToby Isaac   *defaultValue = label->defaultValue;
5215aa44df4SToby Isaac   PetscFunctionReturn(0);
5225aa44df4SToby Isaac }
5235aa44df4SToby Isaac 
52484f0b6dfSMatthew G. Knepley /*@
5255aa44df4SToby Isaac   DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
5265aa44df4SToby Isaac   When a label is created, it is initialized to -1.
5275aa44df4SToby Isaac 
5285aa44df4SToby Isaac   Input parameter:
5295aa44df4SToby Isaac . label - a DMLabel object
5305aa44df4SToby Isaac 
5315aa44df4SToby Isaac   Output parameter:
5325aa44df4SToby Isaac . defaultValue - the default value
5335aa44df4SToby Isaac 
5345aa44df4SToby Isaac   Level: beginner
5355aa44df4SToby Isaac 
5365aa44df4SToby Isaac .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
53784f0b6dfSMatthew G. Knepley @*/
5385aa44df4SToby Isaac PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
5395aa44df4SToby Isaac {
5405aa44df4SToby Isaac   PetscFunctionBegin;
5415aa44df4SToby Isaac   label->defaultValue = defaultValue;
5425aa44df4SToby Isaac   PetscFunctionReturn(0);
5435aa44df4SToby Isaac }
5445aa44df4SToby Isaac 
545c58f1c22SToby Isaac /*@
5465aa44df4SToby 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())
547c58f1c22SToby Isaac 
548c58f1c22SToby Isaac   Input Parameters:
549c58f1c22SToby Isaac + label - the DMLabel
550c58f1c22SToby Isaac - point - the point
551c58f1c22SToby Isaac 
552c58f1c22SToby Isaac   Output Parameter:
553c58f1c22SToby Isaac . value - The point value, or -1
554c58f1c22SToby Isaac 
555c58f1c22SToby Isaac   Level: intermediate
556c58f1c22SToby Isaac 
5575aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
558c58f1c22SToby Isaac @*/
559c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
560c58f1c22SToby Isaac {
561c58f1c22SToby Isaac   PetscInt       v;
562c58f1c22SToby Isaac   PetscErrorCode ierr;
563c58f1c22SToby Isaac 
564c58f1c22SToby Isaac   PetscFunctionBegin;
565c58f1c22SToby Isaac   PetscValidPointer(value, 3);
5665aa44df4SToby Isaac   *value = label->defaultValue;
567c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
568ad8374ffSToby Isaac     if (label->validIS[v]) {
569c58f1c22SToby Isaac       PetscInt i;
570c58f1c22SToby Isaac 
571a2d74346SToby Isaac       ierr = ISLocate(label->points[v],point,&i);CHKERRQ(ierr);
572c58f1c22SToby Isaac       if (i >= 0) {
573c58f1c22SToby Isaac         *value = label->stratumValues[v];
574c58f1c22SToby Isaac         break;
575c58f1c22SToby Isaac       }
576c58f1c22SToby Isaac     } else {
577c58f1c22SToby Isaac       PetscBool has;
578c58f1c22SToby Isaac 
579c58f1c22SToby Isaac       PetscHashIHasKey(label->ht[v], point, has);
580c58f1c22SToby Isaac       if (has) {
581c58f1c22SToby Isaac         *value = label->stratumValues[v];
582c58f1c22SToby Isaac         break;
583c58f1c22SToby Isaac       }
584c58f1c22SToby Isaac     }
585c58f1c22SToby Isaac   }
586c58f1c22SToby Isaac   PetscFunctionReturn(0);
587c58f1c22SToby Isaac }
588c58f1c22SToby Isaac 
589c58f1c22SToby Isaac /*@
5905aa44df4SToby Isaac   DMLabelSetValue - Set the value a label assigns to a point.  If the value is the same as the label's default value (which is initially -1, and can be changed with DMLabelSetDefaultValue() to somethingg different), then this function will do nothing.
591c58f1c22SToby Isaac 
592c58f1c22SToby Isaac   Input Parameters:
593c58f1c22SToby Isaac + label - the DMLabel
594c58f1c22SToby Isaac . point - the point
595c58f1c22SToby Isaac - value - The point value
596c58f1c22SToby Isaac 
597c58f1c22SToby Isaac   Level: intermediate
598c58f1c22SToby Isaac 
5995aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
600c58f1c22SToby Isaac @*/
601c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
602c58f1c22SToby Isaac {
603c58f1c22SToby Isaac   PETSC_UNUSED PetscHashIIter iter, ret;
604c58f1c22SToby Isaac   PetscInt                    v;
605c58f1c22SToby Isaac   PetscErrorCode              ierr;
606c58f1c22SToby Isaac 
607c58f1c22SToby Isaac   PetscFunctionBegin;
608c58f1c22SToby Isaac   /* Find, or add, label value */
6095aa44df4SToby Isaac   if (value == label->defaultValue) PetscFunctionReturn(0);
610c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
611c58f1c22SToby Isaac     if (label->stratumValues[v] == value) break;
612c58f1c22SToby Isaac   }
613c58f1c22SToby Isaac   /* Create new table */
614c58f1c22SToby Isaac   if (v >= label->numStrata) {ierr = DMLabelAddStratum(label, value);CHKERRQ(ierr);}
615c58f1c22SToby Isaac   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
616c58f1c22SToby Isaac   /* Set key */
617c58f1c22SToby Isaac   PetscHashIPut(label->ht[v], point, ret, iter);
618c58f1c22SToby Isaac   PetscFunctionReturn(0);
619c58f1c22SToby Isaac }
620c58f1c22SToby Isaac 
621c58f1c22SToby Isaac /*@
622c58f1c22SToby Isaac   DMLabelClearValue - Clear the value a label assigns to a point
623c58f1c22SToby Isaac 
624c58f1c22SToby Isaac   Input Parameters:
625c58f1c22SToby Isaac + label - the DMLabel
626c58f1c22SToby Isaac . point - the point
627c58f1c22SToby Isaac - value - The point value
628c58f1c22SToby Isaac 
629c58f1c22SToby Isaac   Level: intermediate
630c58f1c22SToby Isaac 
631c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
632c58f1c22SToby Isaac @*/
633c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
634c58f1c22SToby Isaac {
635ad8374ffSToby Isaac   PetscInt       v;
636c58f1c22SToby Isaac   PetscErrorCode ierr;
637c58f1c22SToby Isaac 
638c58f1c22SToby Isaac   PetscFunctionBegin;
639c58f1c22SToby Isaac   /* Find label value */
640c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
641c58f1c22SToby Isaac     if (label->stratumValues[v] == value) break;
642c58f1c22SToby Isaac   }
643c58f1c22SToby Isaac   if (v >= label->numStrata) PetscFunctionReturn(0);
644eeed21e7SToby Isaac   if (label->validIS[v]) {
645eeed21e7SToby Isaac     ierr = DMLabelMakeInvalid_Private(label,v);CHKERRQ(ierr);
646eeed21e7SToby Isaac   }
647eeed21e7SToby Isaac   if (label->bt) {
648eeed21e7SToby 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);
649eeed21e7SToby Isaac     ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
650eeed21e7SToby Isaac   }
651c58f1c22SToby Isaac   ierr = PetscHashIDelKey(label->ht[v], point);CHKERRQ(ierr);
652c58f1c22SToby Isaac   PetscFunctionReturn(0);
653c58f1c22SToby Isaac }
654c58f1c22SToby Isaac 
655c58f1c22SToby Isaac /*@
656c58f1c22SToby Isaac   DMLabelInsertIS - Set all points in the IS to a value
657c58f1c22SToby Isaac 
658c58f1c22SToby Isaac   Input Parameters:
659c58f1c22SToby Isaac + label - the DMLabel
660c58f1c22SToby Isaac . is    - the point IS
661c58f1c22SToby Isaac - value - The point value
662c58f1c22SToby Isaac 
663c58f1c22SToby Isaac   Level: intermediate
664c58f1c22SToby Isaac 
665c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
666c58f1c22SToby Isaac @*/
667c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
668c58f1c22SToby Isaac {
669c58f1c22SToby Isaac   const PetscInt *points;
670c58f1c22SToby Isaac   PetscInt        n, p;
671c58f1c22SToby Isaac   PetscErrorCode  ierr;
672c58f1c22SToby Isaac 
673c58f1c22SToby Isaac   PetscFunctionBegin;
674c58f1c22SToby Isaac   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
675c58f1c22SToby Isaac   ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr);
676c58f1c22SToby Isaac   ierr = ISGetIndices(is, &points);CHKERRQ(ierr);
677c58f1c22SToby Isaac   for (p = 0; p < n; ++p) {ierr = DMLabelSetValue(label, points[p], value);CHKERRQ(ierr);}
678c58f1c22SToby Isaac   ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr);
679c58f1c22SToby Isaac   PetscFunctionReturn(0);
680c58f1c22SToby Isaac }
681c58f1c22SToby Isaac 
68284f0b6dfSMatthew G. Knepley /*@
68384f0b6dfSMatthew G. Knepley   DMLabelGetNumValues - Get the number of values that the DMLabel takes
68484f0b6dfSMatthew G. Knepley 
68584f0b6dfSMatthew G. Knepley   Input Parameter:
68684f0b6dfSMatthew G. Knepley . label - the DMLabel
68784f0b6dfSMatthew G. Knepley 
68884f0b6dfSMatthew G. Knepley   Output Paramater:
68984f0b6dfSMatthew G. Knepley . numValues - the number of values
69084f0b6dfSMatthew G. Knepley 
69184f0b6dfSMatthew G. Knepley   Level: intermediate
69284f0b6dfSMatthew G. Knepley 
69384f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
69484f0b6dfSMatthew G. Knepley @*/
695c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
696c58f1c22SToby Isaac {
697c58f1c22SToby Isaac   PetscFunctionBegin;
698c58f1c22SToby Isaac   PetscValidPointer(numValues, 2);
699c58f1c22SToby Isaac   *numValues = label->numStrata;
700c58f1c22SToby Isaac   PetscFunctionReturn(0);
701c58f1c22SToby Isaac }
702c58f1c22SToby Isaac 
70384f0b6dfSMatthew G. Knepley /*@
70484f0b6dfSMatthew G. Knepley   DMLabelGetValueIS - Get an IS of all values that the DMlabel takes
70584f0b6dfSMatthew G. Knepley 
70684f0b6dfSMatthew G. Knepley   Input Parameter:
70784f0b6dfSMatthew G. Knepley . label - the DMLabel
70884f0b6dfSMatthew G. Knepley 
70984f0b6dfSMatthew G. Knepley   Output Paramater:
71084f0b6dfSMatthew G. Knepley . is    - the value IS
71184f0b6dfSMatthew G. Knepley 
71284f0b6dfSMatthew G. Knepley   Level: intermediate
71384f0b6dfSMatthew G. Knepley 
71484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
71584f0b6dfSMatthew G. Knepley @*/
716c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
717c58f1c22SToby Isaac {
718c58f1c22SToby Isaac   PetscErrorCode ierr;
719c58f1c22SToby Isaac 
720c58f1c22SToby Isaac   PetscFunctionBegin;
721c58f1c22SToby Isaac   PetscValidPointer(values, 2);
722c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr);
723c58f1c22SToby Isaac   PetscFunctionReturn(0);
724c58f1c22SToby Isaac }
725c58f1c22SToby Isaac 
72684f0b6dfSMatthew G. Knepley /*@
72784f0b6dfSMatthew G. Knepley   DMLabelHasStratum - Determine whether points exist with the given value
72884f0b6dfSMatthew G. Knepley 
72984f0b6dfSMatthew G. Knepley   Input Parameters:
73084f0b6dfSMatthew G. Knepley + label - the DMLabel
73184f0b6dfSMatthew G. Knepley - value - the stratum value
73284f0b6dfSMatthew G. Knepley 
73384f0b6dfSMatthew G. Knepley   Output Paramater:
73484f0b6dfSMatthew G. Knepley . exists - Flag saying whether points exist
73584f0b6dfSMatthew G. Knepley 
73684f0b6dfSMatthew G. Knepley   Level: intermediate
73784f0b6dfSMatthew G. Knepley 
73884f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
73984f0b6dfSMatthew G. Knepley @*/
740fada774cSMatthew G. Knepley PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
741fada774cSMatthew G. Knepley {
742fada774cSMatthew G. Knepley   PetscInt v;
743fada774cSMatthew G. Knepley 
744fada774cSMatthew G. Knepley   PetscFunctionBegin;
745fada774cSMatthew G. Knepley   PetscValidPointer(exists, 3);
746fada774cSMatthew G. Knepley   *exists = PETSC_FALSE;
747fada774cSMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
748fada774cSMatthew G. Knepley     if (label->stratumValues[v] == value) {
749fada774cSMatthew G. Knepley       *exists = PETSC_TRUE;
750fada774cSMatthew G. Knepley       break;
751fada774cSMatthew G. Knepley     }
752fada774cSMatthew G. Knepley   }
753fada774cSMatthew G. Knepley   PetscFunctionReturn(0);
754fada774cSMatthew G. Knepley }
755fada774cSMatthew G. Knepley 
75684f0b6dfSMatthew G. Knepley /*@
75784f0b6dfSMatthew G. Knepley   DMLabelGetStratumSize - Get the size of a stratum
75884f0b6dfSMatthew G. Knepley 
75984f0b6dfSMatthew G. Knepley   Input Parameters:
76084f0b6dfSMatthew G. Knepley + label - the DMLabel
76184f0b6dfSMatthew G. Knepley - value - the stratum value
76284f0b6dfSMatthew G. Knepley 
76384f0b6dfSMatthew G. Knepley   Output Paramater:
76484f0b6dfSMatthew G. Knepley . size - The number of points in the stratum
76584f0b6dfSMatthew G. Knepley 
76684f0b6dfSMatthew G. Knepley   Level: intermediate
76784f0b6dfSMatthew G. Knepley 
76884f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
76984f0b6dfSMatthew G. Knepley @*/
770c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
771c58f1c22SToby Isaac {
772c58f1c22SToby Isaac   PetscInt       v;
773c58f1c22SToby Isaac   PetscErrorCode ierr;
774c58f1c22SToby Isaac 
775c58f1c22SToby Isaac   PetscFunctionBegin;
776c58f1c22SToby Isaac   PetscValidPointer(size, 3);
777c58f1c22SToby Isaac   *size = 0;
778c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
779c58f1c22SToby Isaac     if (label->stratumValues[v] == value) {
780c58f1c22SToby Isaac       ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
781c58f1c22SToby Isaac       *size = label->stratumSizes[v];
782c58f1c22SToby Isaac       break;
783c58f1c22SToby Isaac     }
784c58f1c22SToby Isaac   }
785c58f1c22SToby Isaac   PetscFunctionReturn(0);
786c58f1c22SToby Isaac }
787c58f1c22SToby Isaac 
78884f0b6dfSMatthew G. Knepley /*@
78984f0b6dfSMatthew G. Knepley   DMLabelGetStratumBounds - Get the largest and smallest point of a stratum
79084f0b6dfSMatthew G. Knepley 
79184f0b6dfSMatthew G. Knepley   Input Parameters:
79284f0b6dfSMatthew G. Knepley + label - the DMLabel
79384f0b6dfSMatthew G. Knepley - value - the stratum value
79484f0b6dfSMatthew G. Knepley 
79584f0b6dfSMatthew G. Knepley   Output Paramaters:
79684f0b6dfSMatthew G. Knepley + start - the smallest point in the stratum
79784f0b6dfSMatthew G. Knepley - end - the largest point in the stratum
79884f0b6dfSMatthew G. Knepley 
79984f0b6dfSMatthew G. Knepley   Level: intermediate
80084f0b6dfSMatthew G. Knepley 
80184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
80284f0b6dfSMatthew G. Knepley @*/
803c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
804c58f1c22SToby Isaac {
805c58f1c22SToby Isaac   PetscInt       v;
806c58f1c22SToby Isaac   PetscErrorCode ierr;
807c58f1c22SToby Isaac 
808c58f1c22SToby Isaac   PetscFunctionBegin;
809c58f1c22SToby Isaac   if (start) {PetscValidPointer(start, 3); *start = 0;}
810c58f1c22SToby Isaac   if (end)   {PetscValidPointer(end,   4); *end   = 0;}
811c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
812d6cb179aSToby Isaac     PetscInt min, max;
813ad8374ffSToby Isaac 
814c58f1c22SToby Isaac     if (label->stratumValues[v] != value) continue;
815c58f1c22SToby Isaac     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
816c58f1c22SToby Isaac     if (label->stratumSizes[v]  <= 0)     break;
817d6cb179aSToby Isaac     ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr);
818d6cb179aSToby Isaac     if (start) *start = min;
819d6cb179aSToby Isaac     if (end)   *end   = max+1;
820c58f1c22SToby Isaac     break;
821c58f1c22SToby Isaac   }
822c58f1c22SToby Isaac   PetscFunctionReturn(0);
823c58f1c22SToby Isaac }
824c58f1c22SToby Isaac 
82584f0b6dfSMatthew G. Knepley /*@
82684f0b6dfSMatthew G. Knepley   DMLabelGetStratumIS - Get an IS with the stratum points
82784f0b6dfSMatthew G. Knepley 
82884f0b6dfSMatthew G. Knepley   Input Parameters:
82984f0b6dfSMatthew G. Knepley + label - the DMLabel
83084f0b6dfSMatthew G. Knepley - value - the stratum value
83184f0b6dfSMatthew G. Knepley 
83284f0b6dfSMatthew G. Knepley   Output Paramater:
83384f0b6dfSMatthew G. Knepley . points - The stratum points
83484f0b6dfSMatthew G. Knepley 
83584f0b6dfSMatthew G. Knepley   Level: intermediate
83684f0b6dfSMatthew G. Knepley 
83784f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
83884f0b6dfSMatthew G. Knepley @*/
839c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
840c58f1c22SToby Isaac {
841c58f1c22SToby Isaac   PetscInt       v;
842c58f1c22SToby Isaac   PetscErrorCode ierr;
843c58f1c22SToby Isaac 
844c58f1c22SToby Isaac   PetscFunctionBegin;
845c58f1c22SToby Isaac   PetscValidPointer(points, 3);
846c58f1c22SToby Isaac   *points = NULL;
847c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
848c58f1c22SToby Isaac     if (label->stratumValues[v] == value) {
849c58f1c22SToby Isaac       ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
850ad8374ffSToby Isaac       if (label->validIS[v]) {
851ad8374ffSToby Isaac         ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr);
852ad8374ffSToby Isaac         *points = label->points[v];
853c58f1c22SToby Isaac       } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to implement this to speedup Stratify");
854c58f1c22SToby Isaac       break;
855c58f1c22SToby Isaac     }
856c58f1c22SToby Isaac   }
857c58f1c22SToby Isaac   PetscFunctionReturn(0);
858c58f1c22SToby Isaac }
859c58f1c22SToby Isaac 
86084f0b6dfSMatthew G. Knepley /*@
86184f0b6dfSMatthew G. Knepley   DMLabelSetStratumIS - Set the stratum points using an IS
86284f0b6dfSMatthew G. Knepley 
86384f0b6dfSMatthew G. Knepley   Input Parameters:
86484f0b6dfSMatthew G. Knepley + label - the DMLabel
86584f0b6dfSMatthew G. Knepley . value - the stratum value
86684f0b6dfSMatthew G. Knepley - points - The stratum points
86784f0b6dfSMatthew G. Knepley 
86884f0b6dfSMatthew G. Knepley   Level: intermediate
86984f0b6dfSMatthew G. Knepley 
87084f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
87184f0b6dfSMatthew G. Knepley @*/
8724de306b1SToby Isaac PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
8734de306b1SToby Isaac {
8744de306b1SToby Isaac   PetscInt       v, numStrata;
8754de306b1SToby Isaac   PetscErrorCode ierr;
8764de306b1SToby Isaac 
8774de306b1SToby Isaac   PetscFunctionBegin;
8784de306b1SToby Isaac   numStrata = label->numStrata;
8794de306b1SToby Isaac   for (v = 0; v < numStrata; v++) {
8804de306b1SToby Isaac     if (label->stratumValues[v] == value) break;
8814de306b1SToby Isaac   }
8824de306b1SToby Isaac   if (v >= numStrata) {ierr = DMLabelAddStratum(label,value);CHKERRQ(ierr);}
8834de306b1SToby Isaac   if (is == label->points[v]) PetscFunctionReturn(0);
8844de306b1SToby Isaac   ierr = DMLabelClearStratum(label,value);CHKERRQ(ierr);
8854de306b1SToby Isaac   ierr = ISGetLocalSize(is,&(label->stratumSizes[v]));CHKERRQ(ierr);
8864de306b1SToby Isaac   label->stratumValues[v] = value;
8874de306b1SToby Isaac   label->validIS[v] = PETSC_TRUE;
8884de306b1SToby Isaac   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
8894de306b1SToby Isaac   ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
8904de306b1SToby Isaac   if (label->bt) {
8914de306b1SToby Isaac     const PetscInt *points;
8924de306b1SToby Isaac     PetscInt p;
8934de306b1SToby Isaac 
8944de306b1SToby Isaac     ierr = ISGetIndices(is,&points);CHKERRQ(ierr);
8954de306b1SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
8964de306b1SToby Isaac       const PetscInt point = points[p];
8974de306b1SToby Isaac 
8984de306b1SToby 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);
8994de306b1SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
9004de306b1SToby Isaac     }
9014de306b1SToby Isaac   }
9024de306b1SToby Isaac   label->points[v] = is;
9034de306b1SToby Isaac   PetscFunctionReturn(0);
9044de306b1SToby Isaac }
9054de306b1SToby Isaac 
90684f0b6dfSMatthew G. Knepley /*@
90784f0b6dfSMatthew G. Knepley   DMLabelClearStratum - Remove a stratum
9084de306b1SToby Isaac 
90984f0b6dfSMatthew G. Knepley   Input Parameters:
91084f0b6dfSMatthew G. Knepley + label - the DMLabel
91184f0b6dfSMatthew G. Knepley - value - the stratum value
91284f0b6dfSMatthew G. Knepley 
91384f0b6dfSMatthew G. Knepley   Level: intermediate
91484f0b6dfSMatthew G. Knepley 
91584f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
91684f0b6dfSMatthew G. Knepley @*/
917c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
918c58f1c22SToby Isaac {
919c58f1c22SToby Isaac   PetscInt       v;
920c58f1c22SToby Isaac   PetscErrorCode ierr;
921c58f1c22SToby Isaac 
922c58f1c22SToby Isaac   PetscFunctionBegin;
923c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
924c58f1c22SToby Isaac     if (label->stratumValues[v] == value) break;
925c58f1c22SToby Isaac   }
926c58f1c22SToby Isaac   if (v >= label->numStrata) PetscFunctionReturn(0);
9274de306b1SToby Isaac   if (label->validIS[v]) {
9284de306b1SToby Isaac     if (label->bt) {
929c58f1c22SToby Isaac       PetscInt       i;
930ad8374ffSToby Isaac       const PetscInt *points;
931c58f1c22SToby Isaac 
932ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
933c58f1c22SToby Isaac       for (i = 0; i < label->stratumSizes[v]; ++i) {
934ad8374ffSToby Isaac         const PetscInt point = points[i];
935c58f1c22SToby Isaac 
936c58f1c22SToby 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);
937c58f1c22SToby Isaac         ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
938c58f1c22SToby Isaac       }
939ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
940c58f1c22SToby Isaac     }
941ad8374ffSToby Isaac     ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
942c58f1c22SToby Isaac     label->stratumSizes[v] = 0;
943ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
944ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
945c58f1c22SToby Isaac   } else {
946c58f1c22SToby Isaac     PetscHashIClear(label->ht[v]);
947c58f1c22SToby Isaac   }
948c58f1c22SToby Isaac   PetscFunctionReturn(0);
949c58f1c22SToby Isaac }
950c58f1c22SToby Isaac 
95184f0b6dfSMatthew G. Knepley /*@
95284f0b6dfSMatthew G. Knepley   DMLabelFilter - Remove all points outside of [start, end)
95384f0b6dfSMatthew G. Knepley 
95484f0b6dfSMatthew G. Knepley   Input Parameters:
95584f0b6dfSMatthew G. Knepley + label - the DMLabel
95684f0b6dfSMatthew G. Knepley . start - the first point
95784f0b6dfSMatthew G. Knepley - end - the last point
95884f0b6dfSMatthew G. Knepley 
95984f0b6dfSMatthew G. Knepley   Level: intermediate
96084f0b6dfSMatthew G. Knepley 
96184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
96284f0b6dfSMatthew G. Knepley @*/
963c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
964c58f1c22SToby Isaac {
965c58f1c22SToby Isaac   PetscInt       v;
966c58f1c22SToby Isaac   PetscErrorCode ierr;
967c58f1c22SToby Isaac 
968c58f1c22SToby Isaac   PetscFunctionBegin;
969c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
970c58f1c22SToby Isaac   label->pStart = start;
971c58f1c22SToby Isaac   label->pEnd   = end;
972c58f1c22SToby Isaac   if (label->bt) {ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);}
973c58f1c22SToby Isaac   /* Could squish offsets, but would only make sense if I reallocate the storage */
974c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
975c58f1c22SToby Isaac     PetscInt off, q;
976ad8374ffSToby Isaac     const PetscInt *points;
977ad8374ffSToby Isaac     PetscInt *pointsNew = NULL;
978c58f1c22SToby Isaac 
979ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
980c58f1c22SToby Isaac     for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) {
981ad8374ffSToby Isaac       const PetscInt point = points[q];
982c58f1c22SToby Isaac 
983ad8374ffSToby Isaac       if ((point < start) || (point >= end)) {
984ad8374ffSToby Isaac         if (!pointsNew) {
985ad8374ffSToby Isaac           ierr = PetscMalloc1(label->stratumSizes[v],&pointsNew);CHKERRQ(ierr);
986ad8374ffSToby Isaac           ierr = PetscMemcpy(pointsNew,points,(size_t) off * sizeof(PetscInt));CHKERRQ(ierr);
987ad8374ffSToby Isaac         }
988ad8374ffSToby Isaac         continue;
989ad8374ffSToby Isaac       }
990ad8374ffSToby Isaac       if (pointsNew) {
991ad8374ffSToby Isaac         pointsNew[off++] = point;
992ad8374ffSToby Isaac       }
993ad8374ffSToby Isaac     }
994ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
995ad8374ffSToby Isaac     if (pointsNew) {
996ad8374ffSToby Isaac       ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
997ad8374ffSToby Isaac       ierr = ISCreateGeneral(PETSC_COMM_SELF,off,pointsNew,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
998ad8374ffSToby Isaac       ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
999c58f1c22SToby Isaac     }
1000c58f1c22SToby Isaac     label->stratumSizes[v] = off;
1001c58f1c22SToby Isaac   }
1002c58f1c22SToby Isaac   ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr);
1003c58f1c22SToby Isaac   PetscFunctionReturn(0);
1004c58f1c22SToby Isaac }
1005c58f1c22SToby Isaac 
100684f0b6dfSMatthew G. Knepley /*@
100784f0b6dfSMatthew G. Knepley   DMLabelPermute - Create a new label with permuted points
100884f0b6dfSMatthew G. Knepley 
100984f0b6dfSMatthew G. Knepley   Input Parameters:
101084f0b6dfSMatthew G. Knepley + label - the DMLabel
101184f0b6dfSMatthew G. Knepley - permutation - the point permutation
101284f0b6dfSMatthew G. Knepley 
101384f0b6dfSMatthew G. Knepley   Output Parameter:
101484f0b6dfSMatthew G. Knepley . labelnew - the new label containing the permuted points
101584f0b6dfSMatthew G. Knepley 
101684f0b6dfSMatthew G. Knepley   Level: intermediate
101784f0b6dfSMatthew G. Knepley 
101884f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
101984f0b6dfSMatthew G. Knepley @*/
1020c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
1021c58f1c22SToby Isaac {
1022c58f1c22SToby Isaac   const PetscInt *perm;
1023c58f1c22SToby Isaac   PetscInt        numValues, numPoints, v, q;
1024c58f1c22SToby Isaac   PetscErrorCode  ierr;
1025c58f1c22SToby Isaac 
1026c58f1c22SToby Isaac   PetscFunctionBegin;
1027c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1028c58f1c22SToby Isaac   ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr);
1029c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr);
1030c58f1c22SToby Isaac   ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr);
1031c58f1c22SToby Isaac   ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr);
1032c58f1c22SToby Isaac   for (v = 0; v < numValues; ++v) {
1033c58f1c22SToby Isaac     const PetscInt size   = (*labelNew)->stratumSizes[v];
1034ad8374ffSToby Isaac     const PetscInt *points;
1035ad8374ffSToby Isaac     PetscInt *pointsNew;
1036c58f1c22SToby Isaac 
1037ad8374ffSToby Isaac     ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1038ad8374ffSToby Isaac     ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr);
1039c58f1c22SToby Isaac     for (q = 0; q < size; ++q) {
1040ad8374ffSToby Isaac       const PetscInt point = points[q];
1041c58f1c22SToby Isaac 
1042c58f1c22SToby 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);
1043ad8374ffSToby Isaac       pointsNew[q] = perm[point];
1044c58f1c22SToby Isaac     }
1045ad8374ffSToby Isaac     ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1046ad8374ffSToby Isaac     ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr);
1047ad8374ffSToby Isaac     ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr);
1048*fa8e8ae5SToby Isaac     if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) {
1049*fa8e8ae5SToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));CHKERRQ(ierr);
1050*fa8e8ae5SToby Isaac       ierr = PetscFree(pointsNew);CHKERRQ(ierr);
1051*fa8e8ae5SToby Isaac     } else {
1052ad8374ffSToby Isaac       ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr);
1053*fa8e8ae5SToby Isaac     }
1054ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr);
1055c58f1c22SToby Isaac   }
1056c58f1c22SToby Isaac   ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr);
1057c58f1c22SToby Isaac   if (label->bt) {
1058c58f1c22SToby Isaac     ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
1059c58f1c22SToby Isaac     ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
1060c58f1c22SToby Isaac   }
1061c58f1c22SToby Isaac   PetscFunctionReturn(0);
1062c58f1c22SToby Isaac }
1063c58f1c22SToby Isaac 
106426c55118SMichael Lange PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
106526c55118SMichael Lange {
106626c55118SMichael Lange   MPI_Comm       comm;
106726c55118SMichael Lange   PetscInt       s, l, nroots, nleaves, dof, offset, size;
106826c55118SMichael Lange   PetscInt      *remoteOffsets, *rootStrata, *rootIdx;
106926c55118SMichael Lange   PetscSection   rootSection;
107026c55118SMichael Lange   PetscSF        labelSF;
107126c55118SMichael Lange   PetscErrorCode ierr;
107226c55118SMichael Lange 
107326c55118SMichael Lange   PetscFunctionBegin;
107426c55118SMichael Lange   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
107526c55118SMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
107626c55118SMichael Lange   /* Build a section of stratum values per point, generate the according SF
107726c55118SMichael Lange      and distribute point-wise stratum values to leaves. */
107826c55118SMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
107926c55118SMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
108026c55118SMichael Lange   ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr);
108126c55118SMichael Lange   if (label) {
108226c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1083ad8374ffSToby Isaac       const PetscInt *points;
1084ad8374ffSToby Isaac 
1085ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
108626c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1087ad8374ffSToby Isaac         ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr);
1088ad8374ffSToby Isaac         ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr);
108926c55118SMichael Lange       }
1090ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
109126c55118SMichael Lange     }
109226c55118SMichael Lange   }
109326c55118SMichael Lange   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
109426c55118SMichael Lange   /* Create a point-wise array of stratum values */
109526c55118SMichael Lange   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
109626c55118SMichael Lange   ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr);
109726c55118SMichael Lange   ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr);
109826c55118SMichael Lange   if (label) {
109926c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1100ad8374ffSToby Isaac       const PetscInt *points;
1101ad8374ffSToby Isaac 
1102ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
110326c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1104ad8374ffSToby Isaac         const PetscInt p = points[l];
110526c55118SMichael Lange         ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr);
110626c55118SMichael Lange         rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
110726c55118SMichael Lange       }
1108ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
110926c55118SMichael Lange     }
111026c55118SMichael Lange   }
111126c55118SMichael Lange   /* Build SF that maps label points to remote processes */
111226c55118SMichael Lange   ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr);
111326c55118SMichael Lange   ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr);
111426c55118SMichael Lange   ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr);
111526c55118SMichael Lange   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
111626c55118SMichael Lange   /* Send the strata for each point over the derived SF */
111726c55118SMichael Lange   ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr);
111826c55118SMichael Lange   ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr);
111926c55118SMichael Lange   ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
112026c55118SMichael Lange   ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
112126c55118SMichael Lange   /* Clean up */
112226c55118SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
112326c55118SMichael Lange   ierr = PetscFree(rootIdx);CHKERRQ(ierr);
112426c55118SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
112526c55118SMichael Lange   ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr);
112626c55118SMichael Lange   PetscFunctionReturn(0);
112726c55118SMichael Lange }
112826c55118SMichael Lange 
112984f0b6dfSMatthew G. Knepley /*@
113084f0b6dfSMatthew G. Knepley   DMLabelDistribute - Create a new label pushed forward over the PetscSF
113184f0b6dfSMatthew G. Knepley 
113284f0b6dfSMatthew G. Knepley   Input Parameters:
113384f0b6dfSMatthew G. Knepley + label - the DMLabel
113484f0b6dfSMatthew G. Knepley - sf    - the map from old to new distribution
113584f0b6dfSMatthew G. Knepley 
113684f0b6dfSMatthew G. Knepley   Output Parameter:
113784f0b6dfSMatthew G. Knepley . labelnew - the new redistributed label
113884f0b6dfSMatthew G. Knepley 
113984f0b6dfSMatthew G. Knepley   Level: intermediate
114084f0b6dfSMatthew G. Knepley 
114184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
114284f0b6dfSMatthew G. Knepley @*/
1143c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1144c58f1c22SToby Isaac {
1145c58f1c22SToby Isaac   MPI_Comm       comm;
114626c55118SMichael Lange   PetscSection   leafSection;
114726c55118SMichael Lange   PetscInt       p, pStart, pEnd, s, size, dof, offset, stratum;
114826c55118SMichael Lange   PetscInt      *leafStrata, *strataIdx;
1149ad8374ffSToby Isaac   PetscInt     **points;
1150c58f1c22SToby Isaac   char          *name;
1151c58f1c22SToby Isaac   PetscInt       nameSize;
11525cbdf6fcSMichael Lange   PetscHashI     stratumHash;
11535cbdf6fcSMichael Lange   PETSC_UNUSED   PetscHashIIter ret, iter;
1154c58f1c22SToby Isaac   size_t         len = 0;
115526c55118SMichael Lange   PetscMPIInt    rank;
1156c58f1c22SToby Isaac   PetscErrorCode ierr;
1157c58f1c22SToby Isaac 
1158c58f1c22SToby Isaac   PetscFunctionBegin;
1159c58f1c22SToby Isaac   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
1160c58f1c22SToby Isaac   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1161c58f1c22SToby Isaac   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1162c58f1c22SToby Isaac   /* Bcast name */
1163c58f1c22SToby Isaac   if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);}
1164c58f1c22SToby Isaac   nameSize = len;
1165c58f1c22SToby Isaac   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1166c58f1c22SToby Isaac   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1167c58f1c22SToby Isaac   if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);}
1168c58f1c22SToby Isaac   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1169c58f1c22SToby Isaac   ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr);
1170c58f1c22SToby Isaac   ierr = PetscFree(name);CHKERRQ(ierr);
117177d236dfSMichael Lange   /* Bcast defaultValue */
117277d236dfSMichael Lange   if (!rank) (*labelNew)->defaultValue = label->defaultValue;
117377d236dfSMichael Lange   ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
117426c55118SMichael Lange   /* Distribute stratum values over the SF and get the point mapping on the receiver */
117526c55118SMichael Lange   ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr);
11765cbdf6fcSMichael Lange   /* Determine received stratum values and initialise new label*/
11775cbdf6fcSMichael Lange   PetscHashICreate(stratumHash);
117826c55118SMichael Lange   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
11795cbdf6fcSMichael Lange   for (p = 0; p < size; ++p) PetscHashIPut(stratumHash, leafStrata[p], ret, iter);
11805cbdf6fcSMichael Lange   PetscHashISize(stratumHash, (*labelNew)->numStrata);
1181ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr);
1182ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
11835cbdf6fcSMichael Lange   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr);
11845cbdf6fcSMichael Lange   /* Turn leafStrata into indices rather than stratum values */
11855cbdf6fcSMichael Lange   offset = 0;
11865cbdf6fcSMichael Lange   ierr = PetscHashIGetKeys(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr);
1187a205f1fdSToby Isaac   ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr);
11885cbdf6fcSMichael Lange   for (p = 0; p < size; ++p) {
1189231b9e6fSMatthew G. Knepley     for (s = 0; s < (*labelNew)->numStrata; ++s) {
1190231b9e6fSMatthew G. Knepley       if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
11915cbdf6fcSMichael Lange     }
11925cbdf6fcSMichael Lange   }
1193c58f1c22SToby Isaac   /* Rebuild the point strata on the receiver */
1194c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr);
1195c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1196c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1197c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1198c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1199c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1200c58f1c22SToby Isaac       (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1201c58f1c22SToby Isaac     }
1202c58f1c22SToby Isaac   }
1203c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr);
1204c58f1c22SToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr);
1205ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr);
1206c58f1c22SToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1207c58f1c22SToby Isaac     PetscHashICreate((*labelNew)->ht[s]);
1208ad8374ffSToby Isaac     ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr);
1209c58f1c22SToby Isaac   }
1210c58f1c22SToby Isaac   /* Insert points into new strata */
1211c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr);
1212c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1213c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1214c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1215c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1216c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1217c58f1c22SToby Isaac       stratum = leafStrata[offset+s];
1218ad8374ffSToby Isaac       points[stratum][strataIdx[stratum]++] = p;
1219c58f1c22SToby Isaac     }
1220c58f1c22SToby Isaac   }
1221ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; s++) {
1222ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr);
1223ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr);
1224ad8374ffSToby Isaac   }
1225ad8374ffSToby Isaac   ierr = PetscFree(points);CHKERRQ(ierr);
12265cbdf6fcSMichael Lange   PetscHashIDestroy(stratumHash);
1227c58f1c22SToby Isaac   ierr = PetscFree(leafStrata);CHKERRQ(ierr);
1228c58f1c22SToby Isaac   ierr = PetscFree(strataIdx);CHKERRQ(ierr);
1229c58f1c22SToby Isaac   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
1230c58f1c22SToby Isaac   PetscFunctionReturn(0);
1231c58f1c22SToby Isaac }
1232c58f1c22SToby Isaac 
12337937d9ceSMichael Lange /*@
12347937d9ceSMichael Lange   DMLabelGather - Gather all label values from leafs into roots
12357937d9ceSMichael Lange 
12367937d9ceSMichael Lange   Input Parameters:
12377937d9ceSMichael Lange + label - the DMLabel
123884f0b6dfSMatthew G. Knepley - sf - the Star Forest point communication map
12397937d9ceSMichael Lange 
124084f0b6dfSMatthew G. Knepley   Output Parameters:
124184f0b6dfSMatthew G. Knepley . labelNew - the new DMLabel with localised leaf values
12427937d9ceSMichael Lange 
12437937d9ceSMichael Lange   Level: developer
12447937d9ceSMichael Lange 
12457937d9ceSMichael Lange   Note: This is the inverse operation to DMLabelDistribute.
12467937d9ceSMichael Lange 
12477937d9ceSMichael Lange .seealso: DMLabelDistribute()
12487937d9ceSMichael Lange @*/
12497937d9ceSMichael Lange PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
12507937d9ceSMichael Lange {
12517937d9ceSMichael Lange   MPI_Comm       comm;
12527937d9ceSMichael Lange   PetscSection   rootSection;
12537937d9ceSMichael Lange   PetscSF        sfLabel;
12547937d9ceSMichael Lange   PetscSFNode   *rootPoints, *leafPoints;
12557937d9ceSMichael Lange   PetscInt       p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
12567937d9ceSMichael Lange   const PetscInt *rootDegree, *ilocal;
12577937d9ceSMichael Lange   PetscInt       *rootStrata;
12587937d9ceSMichael Lange   char          *name;
12597937d9ceSMichael Lange   PetscInt       nameSize;
12607937d9ceSMichael Lange   size_t         len = 0;
12619852e123SBarry Smith   PetscMPIInt    rank, size;
12627937d9ceSMichael Lange   PetscErrorCode ierr;
12637937d9ceSMichael Lange 
12647937d9ceSMichael Lange   PetscFunctionBegin;
12657937d9ceSMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
12667937d9ceSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
12679852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
12687937d9ceSMichael Lange   /* Bcast name */
12697937d9ceSMichael Lange   if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);}
12707937d9ceSMichael Lange   nameSize = len;
12717937d9ceSMichael Lange   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
12727937d9ceSMichael Lange   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
12737937d9ceSMichael Lange   if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);}
12747937d9ceSMichael Lange   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
12757937d9ceSMichael Lange   ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr);
12767937d9ceSMichael Lange   ierr = PetscFree(name);CHKERRQ(ierr);
12777937d9ceSMichael Lange   /* Gather rank/index pairs of leaves into local roots to build
12787937d9ceSMichael Lange      an inverse, multi-rooted SF. Note that this ignores local leaf
12797937d9ceSMichael Lange      indexing due to the use of the multiSF in PetscSFGather. */
12807937d9ceSMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr);
1281dc53bc9bSMatthew G. Knepley   ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr);
1282dc53bc9bSMatthew G. Knepley   for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
12837937d9ceSMichael Lange   for (p = 0; p < nleaves; p++) {
1284dc53bc9bSMatthew G. Knepley     leafPoints[ilocal[p]].index = ilocal[p];
1285dc53bc9bSMatthew G. Knepley     leafPoints[ilocal[p]].rank  = rank;
12867937d9ceSMichael Lange   }
12877937d9ceSMichael Lange   ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr);
12887937d9ceSMichael Lange   ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr);
12897937d9ceSMichael Lange   for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
12907937d9ceSMichael Lange   ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr);
12917937d9ceSMichael Lange   ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
12927937d9ceSMichael Lange   ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
12937937d9ceSMichael Lange   ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr);
12947937d9ceSMichael Lange   ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
12957937d9ceSMichael Lange   /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
12967937d9ceSMichael Lange   ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr);
12977937d9ceSMichael Lange   /* Rebuild the point strata on the receiver */
12987937d9ceSMichael Lange   for (p = 0, idx = 0; p < nroots; p++) {
12997937d9ceSMichael Lange     for (d = 0; d < rootDegree[p]; d++) {
13007937d9ceSMichael Lange       ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr);
13017937d9ceSMichael Lange       ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr);
13027937d9ceSMichael Lange       for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);}
13037937d9ceSMichael Lange     }
13047937d9ceSMichael Lange     idx += rootDegree[p];
13057937d9ceSMichael Lange   }
130677e0c0e7SMichael Lange   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
130777e0c0e7SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
130877e0c0e7SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
130977e0c0e7SMichael Lange   ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr);
13107937d9ceSMichael Lange   PetscFunctionReturn(0);
13117937d9ceSMichael Lange }
13127937d9ceSMichael Lange 
131384f0b6dfSMatthew G. Knepley /*@
131484f0b6dfSMatthew G. Knepley   DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label
131584f0b6dfSMatthew G. Knepley 
131684f0b6dfSMatthew G. Knepley   Input Parameter:
131784f0b6dfSMatthew G. Knepley . label - the DMLabel
131884f0b6dfSMatthew G. Knepley 
131984f0b6dfSMatthew G. Knepley   Output Parameters:
132084f0b6dfSMatthew G. Knepley + section - the section giving offsets for each stratum
132184f0b6dfSMatthew G. Knepley - is - An IS containing all the label points
132284f0b6dfSMatthew G. Knepley 
132384f0b6dfSMatthew G. Knepley   Level: developer
132484f0b6dfSMatthew G. Knepley 
132584f0b6dfSMatthew G. Knepley .seealso: DMLabelDistribute()
132684f0b6dfSMatthew G. Knepley @*/
1327c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1328c58f1c22SToby Isaac {
1329c58f1c22SToby Isaac   IS              vIS;
1330c58f1c22SToby Isaac   const PetscInt *values;
1331c58f1c22SToby Isaac   PetscInt       *points;
1332c58f1c22SToby Isaac   PetscInt        nV, vS = 0, vE = 0, v, N;
1333c58f1c22SToby Isaac   PetscErrorCode  ierr;
1334c58f1c22SToby Isaac 
1335c58f1c22SToby Isaac   PetscFunctionBegin;
1336c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr);
1337c58f1c22SToby Isaac   ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
1338c58f1c22SToby Isaac   ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr);
1339c58f1c22SToby Isaac   if (nV) {vS = values[0]; vE = values[0]+1;}
1340c58f1c22SToby Isaac   for (v = 1; v < nV; ++v) {
1341c58f1c22SToby Isaac     vS = PetscMin(vS, values[v]);
1342c58f1c22SToby Isaac     vE = PetscMax(vE, values[v]+1);
1343c58f1c22SToby Isaac   }
1344c58f1c22SToby Isaac   ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr);
1345c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr);
1346c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1347c58f1c22SToby Isaac     PetscInt n;
1348c58f1c22SToby Isaac 
1349c58f1c22SToby Isaac     ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr);
1350c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr);
1351c58f1c22SToby Isaac   }
1352c58f1c22SToby Isaac   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
1353c58f1c22SToby Isaac   ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr);
1354c58f1c22SToby Isaac   ierr = PetscMalloc1(N, &points);CHKERRQ(ierr);
1355c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1356c58f1c22SToby Isaac     IS              is;
1357c58f1c22SToby Isaac     const PetscInt *spoints;
1358c58f1c22SToby Isaac     PetscInt        dof, off, p;
1359c58f1c22SToby Isaac 
1360c58f1c22SToby Isaac     ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr);
1361c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr);
1362c58f1c22SToby Isaac     ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr);
1363c58f1c22SToby Isaac     ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr);
1364c58f1c22SToby Isaac     for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1365c58f1c22SToby Isaac     ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr);
1366c58f1c22SToby Isaac     ierr = ISDestroy(&is);CHKERRQ(ierr);
1367c58f1c22SToby Isaac   }
1368c58f1c22SToby Isaac   ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr);
1369c58f1c22SToby Isaac   ierr = ISDestroy(&vIS);CHKERRQ(ierr);
1370c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1371c58f1c22SToby Isaac   PetscFunctionReturn(0);
1372c58f1c22SToby Isaac }
1373c58f1c22SToby Isaac 
137484f0b6dfSMatthew G. Knepley /*@
1375c58f1c22SToby Isaac   PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1376c58f1c22SToby Isaac   the local section and an SF describing the section point overlap.
1377c58f1c22SToby Isaac 
1378c58f1c22SToby Isaac   Input Parameters:
1379c58f1c22SToby Isaac   + s - The PetscSection for the local field layout
1380c58f1c22SToby Isaac   . sf - The SF describing parallel layout of the section points
1381c58f1c22SToby Isaac   . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1382c58f1c22SToby Isaac   . label - The label specifying the points
1383c58f1c22SToby Isaac   - labelValue - The label stratum specifying the points
1384c58f1c22SToby Isaac 
1385c58f1c22SToby Isaac   Output Parameter:
1386c58f1c22SToby Isaac   . gsection - The PetscSection for the global field layout
1387c58f1c22SToby Isaac 
1388c58f1c22SToby Isaac   Note: This gives negative sizes and offsets to points not owned by this process
1389c58f1c22SToby Isaac 
1390c58f1c22SToby Isaac   Level: developer
1391c58f1c22SToby Isaac 
1392c58f1c22SToby Isaac .seealso: PetscSectionCreate()
1393c58f1c22SToby Isaac @*/
1394c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1395c58f1c22SToby Isaac {
1396c58f1c22SToby Isaac   PetscInt      *neg = NULL, *tmpOff = NULL;
1397c58f1c22SToby Isaac   PetscInt       pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1398c58f1c22SToby Isaac   PetscErrorCode ierr;
1399c58f1c22SToby Isaac 
1400c58f1c22SToby Isaac   PetscFunctionBegin;
1401c58f1c22SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr);
1402c58f1c22SToby Isaac   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1403c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr);
1404c58f1c22SToby Isaac   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1405c58f1c22SToby Isaac   if (nroots >= 0) {
1406c58f1c22SToby Isaac     if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1407c58f1c22SToby Isaac     ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr);
1408c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1409c58f1c22SToby Isaac       ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr);
1410c58f1c22SToby Isaac     } else {
1411c58f1c22SToby Isaac       tmpOff = &(*gsection)->atlasDof[-pStart];
1412c58f1c22SToby Isaac     }
1413c58f1c22SToby Isaac   }
1414c58f1c22SToby Isaac   /* Mark ghost points with negative dof */
1415c58f1c22SToby Isaac   for (p = pStart; p < pEnd; ++p) {
1416c58f1c22SToby Isaac     PetscInt value;
1417c58f1c22SToby Isaac 
1418c58f1c22SToby Isaac     ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr);
1419c58f1c22SToby Isaac     if (value != labelValue) continue;
1420c58f1c22SToby Isaac     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1421c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr);
1422c58f1c22SToby Isaac     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1423c58f1c22SToby Isaac     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);}
1424c58f1c22SToby Isaac     if (neg) neg[p] = -(dof+1);
1425c58f1c22SToby Isaac   }
1426c58f1c22SToby Isaac   ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr);
1427c58f1c22SToby Isaac   if (nroots >= 0) {
1428c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1429c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1430c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1431c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1432c58f1c22SToby Isaac     }
1433c58f1c22SToby Isaac   }
1434c58f1c22SToby Isaac   /* Calculate new sizes, get proccess offset, and calculate point offsets */
1435c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1436c58f1c22SToby Isaac     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1437c58f1c22SToby Isaac     (*gsection)->atlasOff[p] = off;
1438c58f1c22SToby Isaac     off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1439c58f1c22SToby Isaac   }
1440c58f1c22SToby Isaac   ierr       = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1441c58f1c22SToby Isaac   globalOff -= off;
1442c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1443c58f1c22SToby Isaac     (*gsection)->atlasOff[p] += globalOff;
1444c58f1c22SToby Isaac     if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1445c58f1c22SToby Isaac   }
1446c58f1c22SToby Isaac   /* Put in negative offsets for ghost points */
1447c58f1c22SToby Isaac   if (nroots >= 0) {
1448c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1449c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1450c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1451c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1452c58f1c22SToby Isaac     }
1453c58f1c22SToby Isaac   }
1454c58f1c22SToby Isaac   if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);}
1455c58f1c22SToby Isaac   ierr = PetscFree(neg);CHKERRQ(ierr);
1456c58f1c22SToby Isaac   PetscFunctionReturn(0);
1457c58f1c22SToby Isaac }
1458c58f1c22SToby Isaac 
14595fdea053SToby Isaac typedef struct _n_PetscSectionSym_Label
14605fdea053SToby Isaac {
14615fdea053SToby Isaac   DMLabel           label;
14625fdea053SToby Isaac   PetscCopyMode     *modes;
14635fdea053SToby Isaac   PetscInt          *sizes;
14645fdea053SToby Isaac   const PetscInt    ***perms;
14655fdea053SToby Isaac   const PetscScalar ***rots;
14665fdea053SToby Isaac   PetscInt          (*minMaxOrients)[2];
14675fdea053SToby Isaac   PetscInt          numStrata; /* numStrata is only increasing, functions as a state */
14685fdea053SToby Isaac } PetscSectionSym_Label;
14695fdea053SToby Isaac 
14705fdea053SToby Isaac static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
14715fdea053SToby Isaac {
14725fdea053SToby Isaac   PetscInt              i, j;
14735fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
14745fdea053SToby Isaac   PetscErrorCode        ierr;
14755fdea053SToby Isaac 
14765fdea053SToby Isaac   PetscFunctionBegin;
14775fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
14785fdea053SToby Isaac     if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
14795fdea053SToby Isaac       for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
14805fdea053SToby Isaac         if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);}
14815fdea053SToby Isaac         if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);}
14825fdea053SToby Isaac       }
14835fdea053SToby Isaac       if (sl->perms[i]) {
14845fdea053SToby Isaac         const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
14855fdea053SToby Isaac 
14865fdea053SToby Isaac         ierr = PetscFree(perms);CHKERRQ(ierr);
14875fdea053SToby Isaac       }
14885fdea053SToby Isaac       if (sl->rots[i]) {
14895fdea053SToby Isaac         const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
14905fdea053SToby Isaac 
14915fdea053SToby Isaac         ierr = PetscFree(rots);CHKERRQ(ierr);
14925fdea053SToby Isaac       }
14935fdea053SToby Isaac     }
14945fdea053SToby Isaac   }
14955fdea053SToby Isaac   ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr);
14965fdea053SToby Isaac   ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr);
14975fdea053SToby Isaac   sl->numStrata = 0;
14985fdea053SToby Isaac   PetscFunctionReturn(0);
14995fdea053SToby Isaac }
15005fdea053SToby Isaac 
15015fdea053SToby Isaac static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
15025fdea053SToby Isaac {
15035fdea053SToby Isaac   PetscErrorCode ierr;
15045fdea053SToby Isaac 
15055fdea053SToby Isaac   PetscFunctionBegin;
15065fdea053SToby Isaac   ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);
15075fdea053SToby Isaac   ierr = PetscFree(sym->data);CHKERRQ(ierr);
15085fdea053SToby Isaac   PetscFunctionReturn(0);
15095fdea053SToby Isaac }
15105fdea053SToby Isaac 
15115fdea053SToby Isaac static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
15125fdea053SToby Isaac {
15135fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
15145fdea053SToby Isaac   PetscBool             isAscii;
15155fdea053SToby Isaac   DMLabel               label = sl->label;
15165fdea053SToby Isaac   PetscErrorCode        ierr;
15175fdea053SToby Isaac 
15185fdea053SToby Isaac   PetscFunctionBegin;
15195fdea053SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
15205fdea053SToby Isaac   if (isAscii) {
15215fdea053SToby Isaac     PetscInt          i, j, k;
15225fdea053SToby Isaac     PetscViewerFormat format;
15235fdea053SToby Isaac 
15245fdea053SToby Isaac     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
15255fdea053SToby Isaac     if (label) {
15265fdea053SToby Isaac       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
15275fdea053SToby Isaac       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
15285fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15295fdea053SToby Isaac         ierr = DMLabelView(label, viewer);CHKERRQ(ierr);
15305fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
15315fdea053SToby Isaac       } else {
15325fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer,"  Label '%s'\n",sl->label->name);CHKERRQ(ierr);
15335fdea053SToby Isaac       }
15345fdea053SToby Isaac     } else {
15355fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr);
15365fdea053SToby Isaac     }
15375fdea053SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15385fdea053SToby Isaac     for (i = 0; i <= sl->numStrata; i++) {
15395fdea053SToby Isaac       PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
15405fdea053SToby Isaac 
15415fdea053SToby Isaac       if (!(sl->perms[i] || sl->rots[i])) {
15425fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr);
15435fdea053SToby Isaac       } else {
15445fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr);
15455fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15465fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr);
15475fdea053SToby Isaac         if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
15485fdea053SToby Isaac           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15495fdea053SToby Isaac           for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
15505fdea053SToby Isaac             if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
15515fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr);
15525fdea053SToby Isaac             } else {
15535fdea053SToby Isaac               PetscInt tab;
15545fdea053SToby Isaac 
15555fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr);
15565fdea053SToby Isaac               ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15575fdea053SToby Isaac               ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr);
15585fdea053SToby Isaac               if (sl->perms[i] && sl->perms[i][j]) {
15595fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr);
15605fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
15615fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);}
15625fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
15635fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
15645fdea053SToby Isaac               }
15655fdea053SToby Isaac               if (sl->rots[i] && sl->rots[i][j]) {
15665fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Rotations:  ");CHKERRQ(ierr);
15675fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
15685fdea053SToby Isaac #if defined(PETSC_USE_COMPLEX)
15695fdea053SToby 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);}
15705fdea053SToby Isaac #else
15715fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);}
15725fdea053SToby Isaac #endif
15735fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
15745fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
15755fdea053SToby Isaac               }
15765fdea053SToby Isaac               ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
15775fdea053SToby Isaac             }
15785fdea053SToby Isaac           }
15795fdea053SToby Isaac           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
15805fdea053SToby Isaac         }
15815fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
15825fdea053SToby Isaac       }
15835fdea053SToby Isaac     }
15845fdea053SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
15855fdea053SToby Isaac   }
15865fdea053SToby Isaac   PetscFunctionReturn(0);
15875fdea053SToby Isaac }
15885fdea053SToby Isaac 
15895fdea053SToby Isaac /*@
15905fdea053SToby Isaac   PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
15915fdea053SToby Isaac 
15925fdea053SToby Isaac   Logically collective on sym
15935fdea053SToby Isaac 
15945fdea053SToby Isaac   Input parameters:
15955fdea053SToby Isaac + sym - the section symmetries
15965fdea053SToby Isaac - label - the DMLabel describing the types of points
15975fdea053SToby Isaac 
15985fdea053SToby Isaac   Level: developer:
15995fdea053SToby Isaac 
16005fdea053SToby Isaac .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()
16015fdea053SToby Isaac @*/
16025fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
16035fdea053SToby Isaac {
16045fdea053SToby Isaac   PetscSectionSym_Label *sl;
16055fdea053SToby Isaac   PetscErrorCode ierr;
16065fdea053SToby Isaac 
16075fdea053SToby Isaac   PetscFunctionBegin;
16085fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
16095fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
16105fdea053SToby Isaac   if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);}
16115fdea053SToby Isaac   if (label) {
16125fdea053SToby Isaac     label->refct++;
16135fdea053SToby Isaac     sl->label = label;
16145fdea053SToby Isaac     ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr);
16151a834cf9SToby 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);
16161a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr);
16171a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr);
16181a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr);
16191a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr);
16201a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr);
16215fdea053SToby Isaac   }
16225fdea053SToby Isaac   PetscFunctionReturn(0);
16235fdea053SToby Isaac }
16245fdea053SToby Isaac 
16255fdea053SToby Isaac /*@C
16265fdea053SToby Isaac   PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
16275fdea053SToby Isaac 
16285fdea053SToby Isaac   Logically collective on PetscSectionSym
16295fdea053SToby Isaac 
16305fdea053SToby Isaac   InputParameters:
16315fdea053SToby Isaac + sys       - the section symmetries
16325fdea053SToby Isaac . stratum   - the stratum value in the label that we are assigning symmetries for
16335fdea053SToby Isaac . size      - the number of dofs for points in the stratum of the label
16345fdea053SToby Isaac . minOrient - the smallest orientation for a point in this stratum
16355fdea053SToby Isaac . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
16365fdea053SToby Isaac . mode      - how sym should copy the perms and rots arrays
16375fdea053SToby Isaac . perms     - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation.  A NULL permutation is the identity
16385fdea053SToby 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
16395fdea053SToby Isaac 
16405fdea053SToby Isaac   Level: developer
16415fdea053SToby Isaac 
16425fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()
16435fdea053SToby Isaac @*/
16445fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
16455fdea053SToby Isaac {
16465fdea053SToby Isaac   PetscInt       i, j, k;
16475fdea053SToby Isaac   PetscSectionSym_Label *sl;
16485fdea053SToby Isaac   PetscErrorCode ierr;
16495fdea053SToby Isaac 
16505fdea053SToby Isaac   PetscFunctionBegin;
16515fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
16525fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
16535fdea053SToby Isaac   if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet");
16545fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
16555fdea053SToby Isaac     PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
16565fdea053SToby Isaac 
16575fdea053SToby Isaac     if (stratum == value) break;
16585fdea053SToby Isaac   }
16595fdea053SToby 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);
16605fdea053SToby Isaac   sl->sizes[i] = size;
16615fdea053SToby Isaac   sl->modes[i] = mode;
16625fdea053SToby Isaac   sl->minMaxOrients[i][0] = minOrient;
16635fdea053SToby Isaac   sl->minMaxOrients[i][1] = maxOrient;
16645fdea053SToby Isaac   if (mode == PETSC_COPY_VALUES) {
16655fdea053SToby Isaac     if (perms) {
16665fdea053SToby Isaac       PetscInt    **ownPerms;
16675fdea053SToby Isaac 
16685fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr);
16695fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
16705fdea053SToby Isaac         if (perms[j]) {
16715fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr);
16725fdea053SToby Isaac           for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];}
16735fdea053SToby Isaac         }
16745fdea053SToby Isaac       }
16755fdea053SToby Isaac       sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient];
16765fdea053SToby Isaac     }
16775fdea053SToby Isaac     if (rots) {
16785fdea053SToby Isaac       PetscScalar **ownRots;
16795fdea053SToby Isaac 
16805fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr);
16815fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
16825fdea053SToby Isaac         if (rots[j]) {
16835fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr);
16845fdea053SToby Isaac           for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];}
16855fdea053SToby Isaac         }
16865fdea053SToby Isaac       }
16875fdea053SToby Isaac       sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient];
16885fdea053SToby Isaac     }
16895fdea053SToby Isaac   } else {
16905fdea053SToby Isaac     sl->perms[i] = perms ? &perms[-minOrient] : NULL;
16915fdea053SToby Isaac     sl->rots[i]  = rots ? &rots[-minOrient] : NULL;
16925fdea053SToby Isaac   }
16935fdea053SToby Isaac   PetscFunctionReturn(0);
16945fdea053SToby Isaac }
16955fdea053SToby Isaac 
16965fdea053SToby Isaac static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
16975fdea053SToby Isaac {
16985fdea053SToby Isaac   PetscInt              i, j, numStrata;
16995fdea053SToby Isaac   PetscSectionSym_Label *sl;
17005fdea053SToby Isaac   DMLabel               label;
17015fdea053SToby Isaac   PetscErrorCode        ierr;
17025fdea053SToby Isaac 
17035fdea053SToby Isaac   PetscFunctionBegin;
17045fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
17055fdea053SToby Isaac   numStrata = sl->numStrata;
17065fdea053SToby Isaac   label     = sl->label;
17075fdea053SToby Isaac   for (i = 0; i < numPoints; i++) {
17085fdea053SToby Isaac     PetscInt point = points[2*i];
17095fdea053SToby Isaac     PetscInt ornt  = points[2*i+1];
17105fdea053SToby Isaac 
17115fdea053SToby Isaac     for (j = 0; j < numStrata; j++) {
17125fdea053SToby Isaac       if (label->validIS[j]) {
17135fdea053SToby Isaac         PetscInt k;
17145fdea053SToby Isaac 
17155fdea053SToby Isaac         ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr);
17165fdea053SToby Isaac         if (k >= 0) break;
17175fdea053SToby Isaac       } else {
17185fdea053SToby Isaac         PetscBool has;
17195fdea053SToby Isaac 
17205fdea053SToby Isaac         PetscHashIHasKey(label->ht[j], point, has);
17215fdea053SToby Isaac         if (has) break;
17225fdea053SToby Isaac       }
17235fdea053SToby Isaac     }
17245fdea053SToby 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);
17255fdea053SToby Isaac     if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;}
17265fdea053SToby Isaac     if (rots) {rots[i]  = sl->rots[j] ? sl->rots[j][ornt] : NULL;}
17275fdea053SToby Isaac   }
17285fdea053SToby Isaac   PetscFunctionReturn(0);
17295fdea053SToby Isaac }
17305fdea053SToby Isaac 
17315fdea053SToby Isaac PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
17325fdea053SToby Isaac {
17335fdea053SToby Isaac   PetscSectionSym_Label *sl;
17345fdea053SToby Isaac   PetscErrorCode        ierr;
17355fdea053SToby Isaac 
17365fdea053SToby Isaac   PetscFunctionBegin;
17375fdea053SToby Isaac   ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr);
17385fdea053SToby Isaac   sym->ops->getpoints = PetscSectionSymGetPoints_Label;
17395fdea053SToby Isaac   sym->ops->view      = PetscSectionSymView_Label;
17405fdea053SToby Isaac   sym->ops->destroy   = PetscSectionSymDestroy_Label;
17415fdea053SToby Isaac   sym->data           = (void *) sl;
17425fdea053SToby Isaac   PetscFunctionReturn(0);
17435fdea053SToby Isaac }
17445fdea053SToby Isaac 
17455fdea053SToby Isaac /*@
17465fdea053SToby Isaac   PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
17475fdea053SToby Isaac 
17485fdea053SToby Isaac   Collective
17495fdea053SToby Isaac 
17505fdea053SToby Isaac   Input Parameters:
17515fdea053SToby Isaac + comm - the MPI communicator for the new symmetry
17525fdea053SToby Isaac - label - the label defining the strata
17535fdea053SToby Isaac 
17545fdea053SToby Isaac   Output Parameters:
17555fdea053SToby Isaac . sym - the section symmetries
17565fdea053SToby Isaac 
17575fdea053SToby Isaac   Level: developer
17585fdea053SToby Isaac 
17595fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()
17605fdea053SToby Isaac @*/
17615fdea053SToby Isaac PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
17625fdea053SToby Isaac {
17635fdea053SToby Isaac   PetscErrorCode        ierr;
17645fdea053SToby Isaac 
17655fdea053SToby Isaac   PetscFunctionBegin;
17665fdea053SToby Isaac   ierr = DMInitializePackage();CHKERRQ(ierr);
17675fdea053SToby Isaac   ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr);
17685fdea053SToby Isaac   ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr);
17695fdea053SToby Isaac   ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr);
17705fdea053SToby Isaac   PetscFunctionReturn(0);
17715fdea053SToby Isaac }
1772