xref: /petsc/src/dm/label/dmlabel.c (revision 033405d5d18a18b6458511dfad85b23321f5ae26)
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 {
580c3c4a36SLisandro Dalcin   PetscInt       off = 0;
59ad8374ffSToby Isaac   PetscInt       *pointArray;
60c58f1c22SToby Isaac   PetscErrorCode ierr;
61c58f1c22SToby Isaac 
62ad8374ffSToby Isaac   if (label->validIS[v]) return 0;
63c58f1c22SToby Isaac   PetscFunctionBegin;
640c3c4a36SLisandro Dalcin   if (v < 0 || v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeValid_Private\n", v);
65c58f1c22SToby Isaac   PetscHashISize(label->ht[v], label->stratumSizes[v]);
66c58f1c22SToby Isaac 
67ad8374ffSToby Isaac   ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr);
68ad8374ffSToby Isaac   ierr = PetscHashIGetKeys(label->ht[v], &off, pointArray);CHKERRQ(ierr);
69c58f1c22SToby 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]);
70c58f1c22SToby Isaac   PetscHashIClear(label->ht[v]);
71ad8374ffSToby Isaac   ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr);
72c58f1c22SToby Isaac   if (label->bt) {
73c58f1c22SToby Isaac     PetscInt p;
74c58f1c22SToby Isaac 
75c58f1c22SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
76ad8374ffSToby Isaac       const PetscInt point = pointArray[p];
77c58f1c22SToby Isaac 
78c58f1c22SToby 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);
79c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
80c58f1c22SToby Isaac     }
81c58f1c22SToby Isaac   }
82ad8374ffSToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF,label->stratumSizes[v],pointArray,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
83ad8374ffSToby Isaac   ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
84ad8374ffSToby Isaac   label->validIS[v] = PETSC_TRUE;
85c58f1c22SToby Isaac   ++label->state;
86c58f1c22SToby Isaac   PetscFunctionReturn(0);
87c58f1c22SToby Isaac }
88c58f1c22SToby Isaac 
89c58f1c22SToby Isaac /*
90c58f1c22SToby Isaac   DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format
91c58f1c22SToby Isaac 
92c58f1c22SToby Isaac   Input parameter:
93c58f1c22SToby Isaac . label - The DMLabel
94c58f1c22SToby Isaac 
95c58f1c22SToby Isaac   Output parameter:
96c58f1c22SToby Isaac . label - The DMLabel with all strata in sorted list format
97c58f1c22SToby Isaac 
98c58f1c22SToby Isaac   Level: developer
99c58f1c22SToby Isaac 
100c58f1c22SToby Isaac .seealso: DMLabelCreate()
101c58f1c22SToby Isaac */
102c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label)
103c58f1c22SToby Isaac {
104c58f1c22SToby Isaac   PetscInt       v;
105c58f1c22SToby Isaac   PetscErrorCode ierr;
106c58f1c22SToby Isaac 
107c58f1c22SToby Isaac   PetscFunctionBegin;
108c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; v++){
109c58f1c22SToby Isaac     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
110c58f1c22SToby Isaac   }
111c58f1c22SToby Isaac   PetscFunctionReturn(0);
112c58f1c22SToby Isaac }
113c58f1c22SToby Isaac 
114c58f1c22SToby Isaac /*
115c58f1c22SToby Isaac   DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format
116c58f1c22SToby Isaac 
117c58f1c22SToby Isaac   Input parameter:
118c58f1c22SToby Isaac + label - The DMLabel
119c58f1c22SToby Isaac - v - The stratum value
120c58f1c22SToby Isaac 
121c58f1c22SToby Isaac   Output parameter:
122c58f1c22SToby Isaac . label - The DMLabel with stratum in hash format
123c58f1c22SToby Isaac 
124c58f1c22SToby Isaac   Level: developer
125c58f1c22SToby Isaac 
126c58f1c22SToby Isaac .seealso: DMLabelCreate()
127c58f1c22SToby Isaac */
128c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v)
129c58f1c22SToby Isaac {
130c58f1c22SToby Isaac   PETSC_UNUSED PetscHashIIter ret, iter;
131c58f1c22SToby Isaac   PetscInt                    p;
132ad8374ffSToby Isaac   const PetscInt              *points;
133c58f1c22SToby Isaac   PetscErrorCode ierr;
134c58f1c22SToby Isaac 
1350c3c4a36SLisandro Dalcin   if (!label->validIS[v]) return 0;
136c58f1c22SToby Isaac   PetscFunctionBegin;
1370c3c4a36SLisandro Dalcin   if (v < 0 || v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeInvalid_Private\n", v);
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 
1560c3c4a36SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index)
1570c3c4a36SLisandro Dalcin {
1580c3c4a36SLisandro Dalcin   PetscInt v;
1590c3c4a36SLisandro Dalcin   PetscFunctionBegin;
1600c3c4a36SLisandro Dalcin   for (*index = -1, v = 0; v < label->numStrata; ++v)
1610c3c4a36SLisandro Dalcin     if (label->stratumValues[v] == value) { *index = v; break; }
1620c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
1630c3c4a36SLisandro Dalcin }
1640c3c4a36SLisandro Dalcin 
1650c3c4a36SLisandro Dalcin static PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index)
166c58f1c22SToby Isaac {
167ad8374ffSToby Isaac   PetscInt   v, *tmpV, *tmpS;
168ad8374ffSToby Isaac   IS         *tmpP;
169c58f1c22SToby Isaac   PetscHashI *tmpH;
170c58f1c22SToby Isaac   PetscBool  *tmpB;
171c58f1c22SToby Isaac   PetscErrorCode ierr;
172c58f1c22SToby Isaac 
173c58f1c22SToby Isaac   PetscFunctionBegin;
174c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpV);CHKERRQ(ierr);
175c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpS);CHKERRQ(ierr);
176c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpH);CHKERRQ(ierr);
177c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpP);CHKERRQ(ierr);
178c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpB);CHKERRQ(ierr);
179c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
180c58f1c22SToby Isaac     tmpV[v] = label->stratumValues[v];
181c58f1c22SToby Isaac     tmpS[v] = label->stratumSizes[v];
182c58f1c22SToby Isaac     tmpH[v] = label->ht[v];
183c58f1c22SToby Isaac     tmpP[v] = label->points[v];
184ad8374ffSToby Isaac     tmpB[v] = label->validIS[v];
185c58f1c22SToby Isaac   }
186c58f1c22SToby Isaac   tmpV[v] = value;
187c58f1c22SToby Isaac   tmpS[v] = 0;
188c58f1c22SToby Isaac   PetscHashICreate(tmpH[v]);
189ad8374ffSToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF,0,NULL,PETSC_OWN_POINTER,&tmpP[v]);CHKERRQ(ierr);
190c58f1c22SToby Isaac   tmpB[v] = PETSC_TRUE;
1910c3c4a36SLisandro Dalcin 
192c58f1c22SToby Isaac   ++label->numStrata;
193c58f1c22SToby Isaac   ierr = PetscFree(label->stratumValues);CHKERRQ(ierr);
194c58f1c22SToby Isaac   ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr);
195c58f1c22SToby Isaac   ierr = PetscFree(label->ht);CHKERRQ(ierr);
196c58f1c22SToby Isaac   ierr = PetscFree(label->points);CHKERRQ(ierr);
197ad8374ffSToby Isaac   ierr = PetscFree(label->validIS);CHKERRQ(ierr);
198c58f1c22SToby Isaac   label->stratumValues = tmpV;
199c58f1c22SToby Isaac   label->stratumSizes  = tmpS;
200c58f1c22SToby Isaac   label->ht            = tmpH;
201c58f1c22SToby Isaac   label->points        = tmpP;
202ad8374ffSToby Isaac   label->validIS       = tmpB;
203c58f1c22SToby Isaac 
2040c3c4a36SLisandro Dalcin   *index = v;
2050c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
2060c3c4a36SLisandro Dalcin }
2070c3c4a36SLisandro Dalcin 
2080c3c4a36SLisandro Dalcin PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
2090c3c4a36SLisandro Dalcin {
2100c3c4a36SLisandro Dalcin   PetscInt       v;
2110c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
2120c3c4a36SLisandro Dalcin 
2130c3c4a36SLisandro Dalcin   PetscFunctionBegin;
2140c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
2150c3c4a36SLisandro Dalcin   if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);}
216c58f1c22SToby Isaac   PetscFunctionReturn(0);
217c58f1c22SToby Isaac }
218c58f1c22SToby Isaac 
219c58f1c22SToby Isaac /*@C
220c58f1c22SToby Isaac   DMLabelGetName - Return the name of a DMLabel object
221c58f1c22SToby Isaac 
222c58f1c22SToby Isaac   Input parameter:
223c58f1c22SToby Isaac . label - The DMLabel
224c58f1c22SToby Isaac 
225c58f1c22SToby Isaac   Output parameter:
226c58f1c22SToby Isaac . name - The label name
227c58f1c22SToby Isaac 
228c58f1c22SToby Isaac   Level: beginner
229c58f1c22SToby Isaac 
230c58f1c22SToby Isaac .seealso: DMLabelCreate()
231c58f1c22SToby Isaac @*/
232c58f1c22SToby Isaac PetscErrorCode DMLabelGetName(DMLabel label, const char **name)
233c58f1c22SToby Isaac {
234c58f1c22SToby Isaac   PetscFunctionBegin;
235c58f1c22SToby Isaac   PetscValidPointer(name, 2);
236c58f1c22SToby Isaac   *name = label->name;
237c58f1c22SToby Isaac   PetscFunctionReturn(0);
238c58f1c22SToby Isaac }
239c58f1c22SToby Isaac 
240c58f1c22SToby Isaac static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
241c58f1c22SToby Isaac {
242c58f1c22SToby Isaac   PetscInt       v;
243c58f1c22SToby Isaac   PetscMPIInt    rank;
244c58f1c22SToby Isaac   PetscErrorCode ierr;
245c58f1c22SToby Isaac 
246c58f1c22SToby Isaac   PetscFunctionBegin;
247c58f1c22SToby Isaac   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr);
248c58f1c22SToby Isaac   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
249c58f1c22SToby Isaac   if (label) {
250c58f1c22SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", label->name);CHKERRQ(ierr);
251c58f1c22SToby Isaac     if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);}
252c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
253c58f1c22SToby Isaac       const PetscInt value = label->stratumValues[v];
254ad8374ffSToby Isaac       const PetscInt *points;
255c58f1c22SToby Isaac       PetscInt       p;
256c58f1c22SToby Isaac 
257ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
258c58f1c22SToby Isaac       for (p = 0; p < label->stratumSizes[v]; ++p) {
259ad8374ffSToby Isaac         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr);
260c58f1c22SToby Isaac       }
261ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
262c58f1c22SToby Isaac     }
263c58f1c22SToby Isaac   }
264c58f1c22SToby Isaac   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
265c58f1c22SToby Isaac   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
266c58f1c22SToby Isaac   PetscFunctionReturn(0);
267c58f1c22SToby Isaac }
268c58f1c22SToby Isaac 
269c58f1c22SToby Isaac /*@C
270c58f1c22SToby Isaac   DMLabelView - View the label
271c58f1c22SToby Isaac 
272c58f1c22SToby Isaac   Input Parameters:
273c58f1c22SToby Isaac + label - The DMLabel
274c58f1c22SToby Isaac - viewer - The PetscViewer
275c58f1c22SToby Isaac 
276c58f1c22SToby Isaac   Level: intermediate
277c58f1c22SToby Isaac 
278c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy()
279c58f1c22SToby Isaac @*/
280c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
281c58f1c22SToby Isaac {
282c58f1c22SToby Isaac   PetscBool      iascii;
283c58f1c22SToby Isaac   PetscErrorCode ierr;
284c58f1c22SToby Isaac 
285c58f1c22SToby Isaac   PetscFunctionBegin;
286c58f1c22SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
287c58f1c22SToby Isaac   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
288c58f1c22SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
289c58f1c22SToby Isaac   if (iascii) {
290c58f1c22SToby Isaac     ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr);
291c58f1c22SToby Isaac   }
292c58f1c22SToby Isaac   PetscFunctionReturn(0);
293c58f1c22SToby Isaac }
294c58f1c22SToby Isaac 
29584f0b6dfSMatthew G. Knepley /*@
29684f0b6dfSMatthew G. Knepley   DMLabelDestroy - Destroys a DMLabel
29784f0b6dfSMatthew G. Knepley 
29884f0b6dfSMatthew G. Knepley   Input Parameter:
29984f0b6dfSMatthew G. Knepley . label - The DMLabel
30084f0b6dfSMatthew G. Knepley 
30184f0b6dfSMatthew G. Knepley   Level: beginner
30284f0b6dfSMatthew G. Knepley 
30384f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate()
30484f0b6dfSMatthew G. Knepley @*/
305c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label)
306c58f1c22SToby Isaac {
307c58f1c22SToby Isaac   PetscInt       v;
308c58f1c22SToby Isaac   PetscErrorCode ierr;
309c58f1c22SToby Isaac 
310c58f1c22SToby Isaac   PetscFunctionBegin;
311c58f1c22SToby Isaac   if (!(*label)) PetscFunctionReturn(0);
312c58f1c22SToby Isaac   if (--(*label)->refct > 0) PetscFunctionReturn(0);
313c58f1c22SToby Isaac   ierr = PetscFree((*label)->name);CHKERRQ(ierr);
314c58f1c22SToby Isaac   ierr = PetscFree((*label)->stratumValues);CHKERRQ(ierr);
315c58f1c22SToby Isaac   ierr = PetscFree((*label)->stratumSizes);CHKERRQ(ierr);
3160c3c4a36SLisandro Dalcin   for (v = 0; v < (*label)->numStrata; ++v) {
3170c3c4a36SLisandro Dalcin     PetscHashIDestroy((*label)->ht[v]);
3180c3c4a36SLisandro Dalcin     ierr = ISDestroy(&((*label)->points[v]));CHKERRQ(ierr);
3190c3c4a36SLisandro Dalcin   }
3200c3c4a36SLisandro Dalcin   ierr = PetscFree((*label)->ht);CHKERRQ(ierr);
321c58f1c22SToby Isaac   ierr = PetscFree((*label)->points);CHKERRQ(ierr);
322ad8374ffSToby Isaac   ierr = PetscFree((*label)->validIS);CHKERRQ(ierr);
323c58f1c22SToby Isaac   ierr = PetscBTDestroy(&(*label)->bt);CHKERRQ(ierr);
324c58f1c22SToby Isaac   ierr = PetscFree(*label);CHKERRQ(ierr);
325c58f1c22SToby Isaac   PetscFunctionReturn(0);
326c58f1c22SToby Isaac }
327c58f1c22SToby Isaac 
32884f0b6dfSMatthew G. Knepley /*@
32984f0b6dfSMatthew G. Knepley   DMLabelDuplicate - Duplicates a DMLabel
33084f0b6dfSMatthew G. Knepley 
33184f0b6dfSMatthew G. Knepley   Input Parameter:
33284f0b6dfSMatthew G. Knepley . label - The DMLabel
33384f0b6dfSMatthew G. Knepley 
33484f0b6dfSMatthew G. Knepley   Output Parameter:
33584f0b6dfSMatthew G. Knepley . labelnew - location to put new vector
33684f0b6dfSMatthew G. Knepley 
33784f0b6dfSMatthew G. Knepley   Level: intermediate
33884f0b6dfSMatthew G. Knepley 
33984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelDestroy()
34084f0b6dfSMatthew G. Knepley @*/
341c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
342c58f1c22SToby Isaac {
343ad8374ffSToby Isaac   PetscInt       v;
344c58f1c22SToby Isaac   PetscErrorCode ierr;
345c58f1c22SToby Isaac 
346c58f1c22SToby Isaac   PetscFunctionBegin;
347c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
348c58f1c22SToby Isaac   ierr = PetscNew(labelnew);CHKERRQ(ierr);
349c58f1c22SToby Isaac   ierr = PetscStrallocpy(label->name, &(*labelnew)->name);CHKERRQ(ierr);
350c58f1c22SToby Isaac 
351c58f1c22SToby Isaac   (*labelnew)->refct        = 1;
352c58f1c22SToby Isaac   (*labelnew)->numStrata    = label->numStrata;
3535aa44df4SToby Isaac   (*labelnew)->defaultValue = label->defaultValue;
354c58f1c22SToby Isaac   if (label->numStrata) {
355c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr);
356c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr);
357c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr);
358c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr);
359ad8374ffSToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr);
360c58f1c22SToby Isaac     /* Could eliminate unused space here */
361c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
362c58f1c22SToby Isaac       PetscHashICreate((*labelnew)->ht[v]);
363ad8374ffSToby Isaac       (*labelnew)->validIS[v]        = PETSC_TRUE;
364c58f1c22SToby Isaac       (*labelnew)->stratumValues[v]  = label->stratumValues[v];
365c58f1c22SToby Isaac       (*labelnew)->stratumSizes[v]   = label->stratumSizes[v];
366ad8374ffSToby Isaac       ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr);
367ad8374ffSToby Isaac       (*labelnew)->points[v]         = label->points[v];
368c58f1c22SToby Isaac     }
369c58f1c22SToby Isaac   }
370c58f1c22SToby Isaac   (*labelnew)->pStart = -1;
371c58f1c22SToby Isaac   (*labelnew)->pEnd   = -1;
372c58f1c22SToby Isaac   (*labelnew)->bt     = NULL;
373c58f1c22SToby Isaac   PetscFunctionReturn(0);
374c58f1c22SToby Isaac }
375c58f1c22SToby Isaac 
376c58f1c22SToby Isaac /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
377c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
378c58f1c22SToby Isaac {
379c58f1c22SToby Isaac   PetscInt       v;
380c58f1c22SToby Isaac   PetscErrorCode ierr;
381c58f1c22SToby Isaac 
382c58f1c22SToby Isaac   PetscFunctionBegin;
3830c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
384c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
385c58f1c22SToby Isaac   label->pStart = pStart;
386c58f1c22SToby Isaac   label->pEnd   = pEnd;
387c58f1c22SToby Isaac   ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr);
388c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
389ad8374ffSToby Isaac     const PetscInt *points;
390c58f1c22SToby Isaac     PetscInt       i;
391c58f1c22SToby Isaac 
392ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
393c58f1c22SToby Isaac     for (i = 0; i < label->stratumSizes[v]; ++i) {
394ad8374ffSToby Isaac       const PetscInt point = points[i];
395c58f1c22SToby Isaac 
396c58f1c22SToby 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);
397c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr);
398c58f1c22SToby Isaac     }
399ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
400c58f1c22SToby Isaac   }
401c58f1c22SToby Isaac   PetscFunctionReturn(0);
402c58f1c22SToby Isaac }
403c58f1c22SToby Isaac 
404c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label)
405c58f1c22SToby Isaac {
406c58f1c22SToby Isaac   PetscErrorCode ierr;
407c58f1c22SToby Isaac 
408c58f1c22SToby Isaac   PetscFunctionBegin;
409c58f1c22SToby Isaac   label->pStart = -1;
410c58f1c22SToby Isaac   label->pEnd   = -1;
4110c3c4a36SLisandro Dalcin   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
412c58f1c22SToby Isaac   PetscFunctionReturn(0);
413c58f1c22SToby Isaac }
414c58f1c22SToby Isaac 
415c58f1c22SToby Isaac /*@
416c58f1c22SToby Isaac   DMLabelHasValue - Determine whether a label assigns the value to any point
417c58f1c22SToby Isaac 
418c58f1c22SToby Isaac   Input Parameters:
419c58f1c22SToby Isaac + label - the DMLabel
420c58f1c22SToby Isaac - value - the value
421c58f1c22SToby Isaac 
422c58f1c22SToby Isaac   Output Parameter:
423c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point
424c58f1c22SToby Isaac 
425c58f1c22SToby Isaac   Level: developer
426c58f1c22SToby Isaac 
427c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
428c58f1c22SToby Isaac @*/
429c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
430c58f1c22SToby Isaac {
431c58f1c22SToby Isaac   PetscInt v;
4320c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
433c58f1c22SToby Isaac 
434c58f1c22SToby Isaac   PetscFunctionBegin;
435c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
4360c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
4370c3c4a36SLisandro Dalcin   *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE;
438c58f1c22SToby Isaac   PetscFunctionReturn(0);
439c58f1c22SToby Isaac }
440c58f1c22SToby Isaac 
441c58f1c22SToby Isaac /*@
442c58f1c22SToby Isaac   DMLabelHasPoint - Determine whether a label assigns a value to a point
443c58f1c22SToby Isaac 
444c58f1c22SToby Isaac   Input Parameters:
445c58f1c22SToby Isaac + label - the DMLabel
446c58f1c22SToby Isaac - point - the point
447c58f1c22SToby Isaac 
448c58f1c22SToby Isaac   Output Parameter:
449c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value
450c58f1c22SToby Isaac 
451c58f1c22SToby Isaac   Note: The user must call DMLabelCreateIndex() before this function.
452c58f1c22SToby Isaac 
453c58f1c22SToby Isaac   Level: developer
454c58f1c22SToby Isaac 
455c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
456c58f1c22SToby Isaac @*/
457c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
458c58f1c22SToby Isaac {
459c58f1c22SToby Isaac   PetscErrorCode ierr;
460c58f1c22SToby Isaac 
461c58f1c22SToby Isaac   PetscFunctionBeginHot;
462c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
463c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
464c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG)
465c58f1c22SToby Isaac   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
466c58f1c22SToby 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);
467c58f1c22SToby Isaac #endif
468c58f1c22SToby Isaac   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
469c58f1c22SToby Isaac   PetscFunctionReturn(0);
470c58f1c22SToby Isaac }
471c58f1c22SToby Isaac 
472c58f1c22SToby Isaac /*@
473c58f1c22SToby Isaac   DMLabelStratumHasPoint - Return true if the stratum contains a point
474c58f1c22SToby Isaac 
475c58f1c22SToby Isaac   Input Parameters:
476c58f1c22SToby Isaac + label - the DMLabel
477c58f1c22SToby Isaac . value - the stratum value
478c58f1c22SToby Isaac - point - the point
479c58f1c22SToby Isaac 
480c58f1c22SToby Isaac   Output Parameter:
481c58f1c22SToby Isaac . contains - true if the stratum contains the point
482c58f1c22SToby Isaac 
483c58f1c22SToby Isaac   Level: intermediate
484c58f1c22SToby Isaac 
485c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
486c58f1c22SToby Isaac @*/
487c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
488c58f1c22SToby Isaac {
489c58f1c22SToby Isaac   PetscInt       v;
490c58f1c22SToby Isaac   PetscErrorCode ierr;
491c58f1c22SToby Isaac 
4920c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
493c58f1c22SToby Isaac   PetscValidPointer(contains, 4);
494c58f1c22SToby Isaac   *contains = PETSC_FALSE;
4950c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
4960c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
4970c3c4a36SLisandro Dalcin 
498ad8374ffSToby Isaac   if (label->validIS[v]) {
499c58f1c22SToby Isaac     PetscInt i;
500c58f1c22SToby Isaac 
501a2d74346SToby Isaac     ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
5020c3c4a36SLisandro Dalcin     if (i >= 0) *contains = PETSC_TRUE;
503c58f1c22SToby Isaac   } else {
504c58f1c22SToby Isaac     PetscBool has;
505c58f1c22SToby Isaac 
506c58f1c22SToby Isaac     PetscHashIHasKey(label->ht[v], point, has);
5070c3c4a36SLisandro Dalcin     if (has) *contains = PETSC_TRUE;
508c58f1c22SToby Isaac   }
509c58f1c22SToby Isaac   PetscFunctionReturn(0);
510c58f1c22SToby Isaac }
511c58f1c22SToby Isaac 
51284f0b6dfSMatthew G. Knepley /*@
5135aa44df4SToby Isaac   DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
5145aa44df4SToby Isaac   When a label is created, it is initialized to -1.
5155aa44df4SToby Isaac 
5165aa44df4SToby Isaac   Input parameter:
5175aa44df4SToby Isaac . label - a DMLabel object
5185aa44df4SToby Isaac 
5195aa44df4SToby Isaac   Output parameter:
5205aa44df4SToby Isaac . defaultValue - the default value
5215aa44df4SToby Isaac 
5225aa44df4SToby Isaac   Level: beginner
5235aa44df4SToby Isaac 
5245aa44df4SToby Isaac .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
52584f0b6dfSMatthew G. Knepley @*/
5265aa44df4SToby Isaac PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
5275aa44df4SToby Isaac {
5285aa44df4SToby Isaac   PetscFunctionBegin;
5295aa44df4SToby Isaac   *defaultValue = label->defaultValue;
5305aa44df4SToby Isaac   PetscFunctionReturn(0);
5315aa44df4SToby Isaac }
5325aa44df4SToby Isaac 
53384f0b6dfSMatthew G. Knepley /*@
5345aa44df4SToby Isaac   DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
5355aa44df4SToby Isaac   When a label is created, it is initialized to -1.
5365aa44df4SToby Isaac 
5375aa44df4SToby Isaac   Input parameter:
5385aa44df4SToby Isaac . label - a DMLabel object
5395aa44df4SToby Isaac 
5405aa44df4SToby Isaac   Output parameter:
5415aa44df4SToby Isaac . defaultValue - the default value
5425aa44df4SToby Isaac 
5435aa44df4SToby Isaac   Level: beginner
5445aa44df4SToby Isaac 
5455aa44df4SToby Isaac .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
54684f0b6dfSMatthew G. Knepley @*/
5475aa44df4SToby Isaac PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
5485aa44df4SToby Isaac {
5495aa44df4SToby Isaac   PetscFunctionBegin;
5505aa44df4SToby Isaac   label->defaultValue = defaultValue;
5515aa44df4SToby Isaac   PetscFunctionReturn(0);
5525aa44df4SToby Isaac }
5535aa44df4SToby Isaac 
554c58f1c22SToby Isaac /*@
5555aa44df4SToby 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())
556c58f1c22SToby Isaac 
557c58f1c22SToby Isaac   Input Parameters:
558c58f1c22SToby Isaac + label - the DMLabel
559c58f1c22SToby Isaac - point - the point
560c58f1c22SToby Isaac 
561c58f1c22SToby Isaac   Output Parameter:
562c58f1c22SToby Isaac . value - The point value, or -1
563c58f1c22SToby Isaac 
564c58f1c22SToby Isaac   Level: intermediate
565c58f1c22SToby Isaac 
5665aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
567c58f1c22SToby Isaac @*/
568c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
569c58f1c22SToby Isaac {
570c58f1c22SToby Isaac   PetscInt       v;
571c58f1c22SToby Isaac   PetscErrorCode ierr;
572c58f1c22SToby Isaac 
5730c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
574c58f1c22SToby Isaac   PetscValidPointer(value, 3);
5755aa44df4SToby Isaac   *value = label->defaultValue;
576c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
577ad8374ffSToby Isaac     if (label->validIS[v]) {
578c58f1c22SToby Isaac       PetscInt i;
579c58f1c22SToby Isaac 
580a2d74346SToby Isaac       ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
581c58f1c22SToby Isaac       if (i >= 0) {
582c58f1c22SToby Isaac         *value = label->stratumValues[v];
583c58f1c22SToby Isaac         break;
584c58f1c22SToby Isaac       }
585c58f1c22SToby Isaac     } else {
586c58f1c22SToby Isaac       PetscBool has;
587c58f1c22SToby Isaac 
588c58f1c22SToby Isaac       PetscHashIHasKey(label->ht[v], point, has);
589c58f1c22SToby Isaac       if (has) {
590c58f1c22SToby Isaac         *value = label->stratumValues[v];
591c58f1c22SToby Isaac         break;
592c58f1c22SToby Isaac       }
593c58f1c22SToby Isaac     }
594c58f1c22SToby Isaac   }
595c58f1c22SToby Isaac   PetscFunctionReturn(0);
596c58f1c22SToby Isaac }
597c58f1c22SToby Isaac 
598c58f1c22SToby Isaac /*@
5995aa44df4SToby 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.
600c58f1c22SToby Isaac 
601c58f1c22SToby Isaac   Input Parameters:
602c58f1c22SToby Isaac + label - the DMLabel
603c58f1c22SToby Isaac . point - the point
604c58f1c22SToby Isaac - value - The point value
605c58f1c22SToby Isaac 
606c58f1c22SToby Isaac   Level: intermediate
607c58f1c22SToby Isaac 
6085aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
609c58f1c22SToby Isaac @*/
610c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
611c58f1c22SToby Isaac {
612c58f1c22SToby Isaac   PETSC_UNUSED PetscHashIIter iter, ret;
613c58f1c22SToby Isaac   PetscInt                    v;
614c58f1c22SToby Isaac   PetscErrorCode              ierr;
615c58f1c22SToby Isaac 
616c58f1c22SToby Isaac   PetscFunctionBegin;
6170c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
6185aa44df4SToby Isaac   if (value == label->defaultValue) PetscFunctionReturn(0);
6190c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
6200c3c4a36SLisandro Dalcin   if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);}
621c58f1c22SToby Isaac   /* Set key */
6220c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
623c58f1c22SToby Isaac   PetscHashIPut(label->ht[v], point, ret, iter);
624c58f1c22SToby Isaac   PetscFunctionReturn(0);
625c58f1c22SToby Isaac }
626c58f1c22SToby Isaac 
627c58f1c22SToby Isaac /*@
628c58f1c22SToby Isaac   DMLabelClearValue - Clear the value a label assigns to a point
629c58f1c22SToby Isaac 
630c58f1c22SToby Isaac   Input Parameters:
631c58f1c22SToby Isaac + label - the DMLabel
632c58f1c22SToby Isaac . point - the point
633c58f1c22SToby Isaac - value - The point value
634c58f1c22SToby Isaac 
635c58f1c22SToby Isaac   Level: intermediate
636c58f1c22SToby Isaac 
637c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
638c58f1c22SToby Isaac @*/
639c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
640c58f1c22SToby Isaac {
641ad8374ffSToby Isaac   PetscInt       v;
642c58f1c22SToby Isaac   PetscErrorCode ierr;
643c58f1c22SToby Isaac 
644c58f1c22SToby Isaac   PetscFunctionBegin;
645c58f1c22SToby Isaac   /* Find label value */
6460c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
6470c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
6480c3c4a36SLisandro Dalcin 
649eeed21e7SToby Isaac   if (label->bt) {
650eeed21e7SToby 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);
651eeed21e7SToby Isaac     ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
652eeed21e7SToby Isaac   }
6530c3c4a36SLisandro Dalcin 
6540c3c4a36SLisandro Dalcin   /* Delete key */
6550c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
656c58f1c22SToby Isaac   ierr = PetscHashIDelKey(label->ht[v], point);CHKERRQ(ierr);
657c58f1c22SToby Isaac   PetscFunctionReturn(0);
658c58f1c22SToby Isaac }
659c58f1c22SToby Isaac 
660c58f1c22SToby Isaac /*@
661c58f1c22SToby Isaac   DMLabelInsertIS - Set all points in the IS to a value
662c58f1c22SToby Isaac 
663c58f1c22SToby Isaac   Input Parameters:
664c58f1c22SToby Isaac + label - the DMLabel
665c58f1c22SToby Isaac . is    - the point IS
666c58f1c22SToby Isaac - value - The point value
667c58f1c22SToby Isaac 
668c58f1c22SToby Isaac   Level: intermediate
669c58f1c22SToby Isaac 
670c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
671c58f1c22SToby Isaac @*/
672c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
673c58f1c22SToby Isaac {
6740c3c4a36SLisandro Dalcin   PETSC_UNUSED PetscHashIIter iter, ret;
6750c3c4a36SLisandro Dalcin   PetscInt        v, n, p;
676c58f1c22SToby Isaac   const PetscInt *points;
677c58f1c22SToby Isaac   PetscErrorCode  ierr;
678c58f1c22SToby Isaac 
679c58f1c22SToby Isaac   PetscFunctionBegin;
680c58f1c22SToby Isaac   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
6810c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
6820c3c4a36SLisandro Dalcin   if (value == label->defaultValue) PetscFunctionReturn(0);
6830c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
6840c3c4a36SLisandro Dalcin   if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);}
6850c3c4a36SLisandro Dalcin   /* Set keys */
6860c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
687c58f1c22SToby Isaac   ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr);
688c58f1c22SToby Isaac   ierr = ISGetIndices(is, &points);CHKERRQ(ierr);
6890c3c4a36SLisandro Dalcin   for (p = 0; p < n; ++p) PetscHashIPut(label->ht[v], points[p], ret, iter);
690c58f1c22SToby Isaac   ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr);
691c58f1c22SToby Isaac   PetscFunctionReturn(0);
692c58f1c22SToby Isaac }
693c58f1c22SToby Isaac 
69484f0b6dfSMatthew G. Knepley /*@
69584f0b6dfSMatthew G. Knepley   DMLabelGetNumValues - Get the number of values that the DMLabel takes
69684f0b6dfSMatthew G. Knepley 
69784f0b6dfSMatthew G. Knepley   Input Parameter:
69884f0b6dfSMatthew G. Knepley . label - the DMLabel
69984f0b6dfSMatthew G. Knepley 
70084f0b6dfSMatthew G. Knepley   Output Paramater:
70184f0b6dfSMatthew G. Knepley . numValues - the number of values
70284f0b6dfSMatthew G. Knepley 
70384f0b6dfSMatthew G. Knepley   Level: intermediate
70484f0b6dfSMatthew G. Knepley 
70584f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
70684f0b6dfSMatthew G. Knepley @*/
707c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
708c58f1c22SToby Isaac {
709c58f1c22SToby Isaac   PetscFunctionBegin;
710c58f1c22SToby Isaac   PetscValidPointer(numValues, 2);
711c58f1c22SToby Isaac   *numValues = label->numStrata;
712c58f1c22SToby Isaac   PetscFunctionReturn(0);
713c58f1c22SToby Isaac }
714c58f1c22SToby Isaac 
71584f0b6dfSMatthew G. Knepley /*@
71684f0b6dfSMatthew G. Knepley   DMLabelGetValueIS - Get an IS of all values that the DMlabel takes
71784f0b6dfSMatthew G. Knepley 
71884f0b6dfSMatthew G. Knepley   Input Parameter:
71984f0b6dfSMatthew G. Knepley . label - the DMLabel
72084f0b6dfSMatthew G. Knepley 
72184f0b6dfSMatthew G. Knepley   Output Paramater:
72284f0b6dfSMatthew G. Knepley . is    - the value IS
72384f0b6dfSMatthew G. Knepley 
72484f0b6dfSMatthew G. Knepley   Level: intermediate
72584f0b6dfSMatthew G. Knepley 
72684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
72784f0b6dfSMatthew G. Knepley @*/
728c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
729c58f1c22SToby Isaac {
730c58f1c22SToby Isaac   PetscErrorCode ierr;
731c58f1c22SToby Isaac 
732c58f1c22SToby Isaac   PetscFunctionBegin;
733c58f1c22SToby Isaac   PetscValidPointer(values, 2);
734c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr);
735c58f1c22SToby Isaac   PetscFunctionReturn(0);
736c58f1c22SToby Isaac }
737c58f1c22SToby Isaac 
73884f0b6dfSMatthew G. Knepley /*@
73984f0b6dfSMatthew G. Knepley   DMLabelHasStratum - Determine whether points exist with the given value
74084f0b6dfSMatthew G. Knepley 
74184f0b6dfSMatthew G. Knepley   Input Parameters:
74284f0b6dfSMatthew G. Knepley + label - the DMLabel
74384f0b6dfSMatthew G. Knepley - value - the stratum value
74484f0b6dfSMatthew G. Knepley 
74584f0b6dfSMatthew G. Knepley   Output Paramater:
74684f0b6dfSMatthew G. Knepley . exists - Flag saying whether points exist
74784f0b6dfSMatthew G. Knepley 
74884f0b6dfSMatthew G. Knepley   Level: intermediate
74984f0b6dfSMatthew G. Knepley 
75084f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
75184f0b6dfSMatthew G. Knepley @*/
752fada774cSMatthew G. Knepley PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
753fada774cSMatthew G. Knepley {
754fada774cSMatthew G. Knepley   PetscInt v;
7550c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
756fada774cSMatthew G. Knepley 
757fada774cSMatthew G. Knepley   PetscFunctionBegin;
758fada774cSMatthew G. Knepley   PetscValidPointer(exists, 3);
7590c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
7600c3c4a36SLisandro Dalcin   *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE;
761fada774cSMatthew G. Knepley   PetscFunctionReturn(0);
762fada774cSMatthew G. Knepley }
763fada774cSMatthew G. Knepley 
76484f0b6dfSMatthew G. Knepley /*@
76584f0b6dfSMatthew G. Knepley   DMLabelGetStratumSize - Get the size of a stratum
76684f0b6dfSMatthew G. Knepley 
76784f0b6dfSMatthew G. Knepley   Input Parameters:
76884f0b6dfSMatthew G. Knepley + label - the DMLabel
76984f0b6dfSMatthew G. Knepley - value - the stratum value
77084f0b6dfSMatthew G. Knepley 
77184f0b6dfSMatthew G. Knepley   Output Paramater:
77284f0b6dfSMatthew G. Knepley . size - The number of points in the stratum
77384f0b6dfSMatthew G. Knepley 
77484f0b6dfSMatthew G. Knepley   Level: intermediate
77584f0b6dfSMatthew G. Knepley 
77684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
77784f0b6dfSMatthew G. Knepley @*/
778c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
779c58f1c22SToby Isaac {
780c58f1c22SToby Isaac   PetscInt       v;
781c58f1c22SToby Isaac   PetscErrorCode ierr;
782c58f1c22SToby Isaac 
783c58f1c22SToby Isaac   PetscFunctionBegin;
784c58f1c22SToby Isaac   PetscValidPointer(size, 3);
785c58f1c22SToby Isaac   *size = 0;
7860c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
7870c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
788c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
789c58f1c22SToby Isaac   *size = label->stratumSizes[v];
790c58f1c22SToby Isaac   PetscFunctionReturn(0);
791c58f1c22SToby Isaac }
792c58f1c22SToby Isaac 
79384f0b6dfSMatthew G. Knepley /*@
79484f0b6dfSMatthew G. Knepley   DMLabelGetStratumBounds - Get the largest and smallest point of a stratum
79584f0b6dfSMatthew G. Knepley 
79684f0b6dfSMatthew G. Knepley   Input Parameters:
79784f0b6dfSMatthew G. Knepley + label - the DMLabel
79884f0b6dfSMatthew G. Knepley - value - the stratum value
79984f0b6dfSMatthew G. Knepley 
80084f0b6dfSMatthew G. Knepley   Output Paramaters:
80184f0b6dfSMatthew G. Knepley + start - the smallest point in the stratum
80284f0b6dfSMatthew G. Knepley - end - the largest point in the stratum
80384f0b6dfSMatthew G. Knepley 
80484f0b6dfSMatthew G. Knepley   Level: intermediate
80584f0b6dfSMatthew G. Knepley 
80684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
80784f0b6dfSMatthew G. Knepley @*/
808c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
809c58f1c22SToby Isaac {
8100c3c4a36SLisandro Dalcin   PetscInt       v, min, max;
811c58f1c22SToby Isaac   PetscErrorCode ierr;
812c58f1c22SToby Isaac 
813c58f1c22SToby Isaac   PetscFunctionBegin;
814c58f1c22SToby Isaac   if (start) {PetscValidPointer(start, 3); *start = 0;}
815c58f1c22SToby Isaac   if (end)   {PetscValidPointer(end,   4); *end   = 0;}
8160c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
8170c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
818c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
8190c3c4a36SLisandro Dalcin   if (label->stratumSizes[v] <= 0) PetscFunctionReturn(0);
820d6cb179aSToby Isaac   ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr);
821d6cb179aSToby Isaac   if (start) *start = min;
822d6cb179aSToby Isaac   if (end)   *end   = max+1;
823c58f1c22SToby Isaac   PetscFunctionReturn(0);
824c58f1c22SToby Isaac }
825c58f1c22SToby Isaac 
82684f0b6dfSMatthew G. Knepley /*@
82784f0b6dfSMatthew G. Knepley   DMLabelGetStratumIS - Get an IS with the stratum points
82884f0b6dfSMatthew G. Knepley 
82984f0b6dfSMatthew G. Knepley   Input Parameters:
83084f0b6dfSMatthew G. Knepley + label - the DMLabel
83184f0b6dfSMatthew G. Knepley - value - the stratum value
83284f0b6dfSMatthew G. Knepley 
83384f0b6dfSMatthew G. Knepley   Output Paramater:
83484f0b6dfSMatthew G. Knepley . points - The stratum points
83584f0b6dfSMatthew G. Knepley 
83684f0b6dfSMatthew G. Knepley   Level: intermediate
83784f0b6dfSMatthew G. Knepley 
83884f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
83984f0b6dfSMatthew G. Knepley @*/
840c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
841c58f1c22SToby Isaac {
842c58f1c22SToby Isaac   PetscInt       v;
843c58f1c22SToby Isaac   PetscErrorCode ierr;
844c58f1c22SToby Isaac 
845c58f1c22SToby Isaac   PetscFunctionBegin;
846c58f1c22SToby Isaac   PetscValidPointer(points, 3);
847c58f1c22SToby Isaac   *points = NULL;
8480c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
8490c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
850c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
851ad8374ffSToby Isaac   ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr);
852ad8374ffSToby Isaac   *points = label->points[v];
853c58f1c22SToby Isaac   PetscFunctionReturn(0);
854c58f1c22SToby Isaac }
855c58f1c22SToby Isaac 
85684f0b6dfSMatthew G. Knepley /*@
85784f0b6dfSMatthew G. Knepley   DMLabelSetStratumIS - Set the stratum points using an IS
85884f0b6dfSMatthew G. Knepley 
85984f0b6dfSMatthew G. Knepley   Input Parameters:
86084f0b6dfSMatthew G. Knepley + label - the DMLabel
86184f0b6dfSMatthew G. Knepley . value - the stratum value
86284f0b6dfSMatthew G. Knepley - points - The stratum points
86384f0b6dfSMatthew G. Knepley 
86484f0b6dfSMatthew G. Knepley   Level: intermediate
86584f0b6dfSMatthew G. Knepley 
86684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
86784f0b6dfSMatthew G. Knepley @*/
8684de306b1SToby Isaac PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
8694de306b1SToby Isaac {
8700c3c4a36SLisandro Dalcin   PetscInt       v;
8714de306b1SToby Isaac   PetscErrorCode ierr;
8724de306b1SToby Isaac 
8734de306b1SToby Isaac   PetscFunctionBegin;
8740c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
8750c3c4a36SLisandro Dalcin   if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);}
8764de306b1SToby Isaac   if (is == label->points[v]) PetscFunctionReturn(0);
8774de306b1SToby Isaac   ierr = DMLabelClearStratum(label, value);CHKERRQ(ierr);
8784de306b1SToby Isaac   ierr = ISGetLocalSize(is, &(label->stratumSizes[v]));CHKERRQ(ierr);
8794de306b1SToby Isaac   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
8804de306b1SToby Isaac   ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
8810c3c4a36SLisandro Dalcin   label->points[v] = is;
8820c3c4a36SLisandro Dalcin   label->validIS[v] = PETSC_TRUE;
8834de306b1SToby Isaac   if (label->bt) {
8844de306b1SToby Isaac     const PetscInt *points;
8854de306b1SToby Isaac     PetscInt p;
8864de306b1SToby Isaac 
8874de306b1SToby Isaac     ierr = ISGetIndices(is,&points);CHKERRQ(ierr);
8884de306b1SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
8894de306b1SToby Isaac       const PetscInt point = points[p];
8904de306b1SToby Isaac 
8914de306b1SToby 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);
8924de306b1SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
8934de306b1SToby Isaac     }
8944de306b1SToby Isaac   }
8954de306b1SToby Isaac   PetscFunctionReturn(0);
8964de306b1SToby Isaac }
8974de306b1SToby Isaac 
89884f0b6dfSMatthew G. Knepley /*@
89984f0b6dfSMatthew G. Knepley   DMLabelClearStratum - Remove a stratum
9004de306b1SToby Isaac 
90184f0b6dfSMatthew G. Knepley   Input Parameters:
90284f0b6dfSMatthew G. Knepley + label - the DMLabel
90384f0b6dfSMatthew G. Knepley - value - the stratum value
90484f0b6dfSMatthew G. Knepley 
90584f0b6dfSMatthew G. Knepley   Level: intermediate
90684f0b6dfSMatthew G. Knepley 
90784f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
90884f0b6dfSMatthew G. Knepley @*/
909c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
910c58f1c22SToby Isaac {
911c58f1c22SToby Isaac   PetscInt       v;
912c58f1c22SToby Isaac   PetscErrorCode ierr;
913c58f1c22SToby Isaac 
914c58f1c22SToby Isaac   PetscFunctionBegin;
9150c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
9160c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
9174de306b1SToby Isaac   if (label->validIS[v]) {
9184de306b1SToby Isaac     if (label->bt) {
919c58f1c22SToby Isaac       PetscInt       i;
920ad8374ffSToby Isaac       const PetscInt *points;
921c58f1c22SToby Isaac 
922ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
923c58f1c22SToby Isaac       for (i = 0; i < label->stratumSizes[v]; ++i) {
924ad8374ffSToby Isaac         const PetscInt point = points[i];
925c58f1c22SToby Isaac 
926c58f1c22SToby 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);
927c58f1c22SToby Isaac         ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
928c58f1c22SToby Isaac       }
929ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
930c58f1c22SToby Isaac     }
931c58f1c22SToby Isaac     label->stratumSizes[v] = 0;
9320c3c4a36SLisandro Dalcin     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
9330c3c4a36SLisandro Dalcin     ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr);
9340c3c4a36SLisandro Dalcin     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
935c58f1c22SToby Isaac   } else {
936c58f1c22SToby Isaac     PetscHashIClear(label->ht[v]);
937c58f1c22SToby Isaac   }
938c58f1c22SToby Isaac   PetscFunctionReturn(0);
939c58f1c22SToby Isaac }
940c58f1c22SToby Isaac 
94184f0b6dfSMatthew G. Knepley /*@
94284f0b6dfSMatthew G. Knepley   DMLabelFilter - Remove all points outside of [start, end)
94384f0b6dfSMatthew G. Knepley 
94484f0b6dfSMatthew G. Knepley   Input Parameters:
94584f0b6dfSMatthew G. Knepley + label - the DMLabel
94684f0b6dfSMatthew G. Knepley . start - the first point
94784f0b6dfSMatthew G. Knepley - end - the last point
94884f0b6dfSMatthew G. Knepley 
94984f0b6dfSMatthew G. Knepley   Level: intermediate
95084f0b6dfSMatthew G. Knepley 
95184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
95284f0b6dfSMatthew G. Knepley @*/
953c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
954c58f1c22SToby Isaac {
955c58f1c22SToby Isaac   PetscInt       v;
956c58f1c22SToby Isaac   PetscErrorCode ierr;
957c58f1c22SToby Isaac 
958c58f1c22SToby Isaac   PetscFunctionBegin;
9590c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
960c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
961c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
962c58f1c22SToby Isaac     PetscInt off, q;
963ad8374ffSToby Isaac     const PetscInt *points;
964*033405d5SLisandro Dalcin     PetscInt numPointsNew = 0, *pointsNew = NULL;
965c58f1c22SToby Isaac 
966ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
967*033405d5SLisandro Dalcin     for (q = 0; q < label->stratumSizes[v]; ++q)
968*033405d5SLisandro Dalcin       if (points[q] >= start && points[q] < end)
969*033405d5SLisandro Dalcin         numPointsNew++;
970*033405d5SLisandro Dalcin     ierr = PetscMalloc1(numPointsNew, &pointsNew);CHKERRQ(ierr);
971c58f1c22SToby Isaac     for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) {
972*033405d5SLisandro Dalcin       if (points[q] >= start && points[q] < end)
973*033405d5SLisandro Dalcin         pointsNew[off++] = points[q];
974ad8374ffSToby Isaac     }
975ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
976*033405d5SLisandro Dalcin 
977*033405d5SLisandro Dalcin     label->stratumSizes[v] = numPointsNew;
978*033405d5SLisandro Dalcin     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
979*033405d5SLisandro Dalcin     ierr = ISCreateGeneral(PETSC_COMM_SELF,numPointsNew, pointsNew, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr);
980*033405d5SLisandro Dalcin     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
981c58f1c22SToby Isaac   }
982c58f1c22SToby Isaac   ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr);
983c58f1c22SToby Isaac   PetscFunctionReturn(0);
984c58f1c22SToby Isaac }
985c58f1c22SToby Isaac 
98684f0b6dfSMatthew G. Knepley /*@
98784f0b6dfSMatthew G. Knepley   DMLabelPermute - Create a new label with permuted points
98884f0b6dfSMatthew G. Knepley 
98984f0b6dfSMatthew G. Knepley   Input Parameters:
99084f0b6dfSMatthew G. Knepley + label - the DMLabel
99184f0b6dfSMatthew G. Knepley - permutation - the point permutation
99284f0b6dfSMatthew G. Knepley 
99384f0b6dfSMatthew G. Knepley   Output Parameter:
99484f0b6dfSMatthew G. Knepley . labelnew - the new label containing the permuted points
99584f0b6dfSMatthew G. Knepley 
99684f0b6dfSMatthew G. Knepley   Level: intermediate
99784f0b6dfSMatthew G. Knepley 
99884f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
99984f0b6dfSMatthew G. Knepley @*/
1000c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
1001c58f1c22SToby Isaac {
1002c58f1c22SToby Isaac   const PetscInt *perm;
1003c58f1c22SToby Isaac   PetscInt        numValues, numPoints, v, q;
1004c58f1c22SToby Isaac   PetscErrorCode  ierr;
1005c58f1c22SToby Isaac 
1006c58f1c22SToby Isaac   PetscFunctionBegin;
1007c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1008c58f1c22SToby Isaac   ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr);
1009c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr);
1010c58f1c22SToby Isaac   ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr);
1011c58f1c22SToby Isaac   ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr);
1012c58f1c22SToby Isaac   for (v = 0; v < numValues; ++v) {
1013c58f1c22SToby Isaac     const PetscInt size   = (*labelNew)->stratumSizes[v];
1014ad8374ffSToby Isaac     const PetscInt *points;
1015ad8374ffSToby Isaac     PetscInt *pointsNew;
1016c58f1c22SToby Isaac 
1017ad8374ffSToby Isaac     ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1018ad8374ffSToby Isaac     ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr);
1019c58f1c22SToby Isaac     for (q = 0; q < size; ++q) {
1020ad8374ffSToby Isaac       const PetscInt point = points[q];
1021c58f1c22SToby Isaac 
1022c58f1c22SToby 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);
1023ad8374ffSToby Isaac       pointsNew[q] = perm[point];
1024c58f1c22SToby Isaac     }
1025ad8374ffSToby Isaac     ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1026ad8374ffSToby Isaac     ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr);
1027ad8374ffSToby Isaac     ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr);
1028ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr);
1029ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr);
1030c58f1c22SToby Isaac   }
1031c58f1c22SToby Isaac   ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr);
1032c58f1c22SToby Isaac   if (label->bt) {
1033c58f1c22SToby Isaac     ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
1034c58f1c22SToby Isaac     ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
1035c58f1c22SToby Isaac   }
1036c58f1c22SToby Isaac   PetscFunctionReturn(0);
1037c58f1c22SToby Isaac }
1038c58f1c22SToby Isaac 
103926c55118SMichael Lange PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
104026c55118SMichael Lange {
104126c55118SMichael Lange   MPI_Comm       comm;
104226c55118SMichael Lange   PetscInt       s, l, nroots, nleaves, dof, offset, size;
104326c55118SMichael Lange   PetscInt      *remoteOffsets, *rootStrata, *rootIdx;
104426c55118SMichael Lange   PetscSection   rootSection;
104526c55118SMichael Lange   PetscSF        labelSF;
104626c55118SMichael Lange   PetscErrorCode ierr;
104726c55118SMichael Lange 
104826c55118SMichael Lange   PetscFunctionBegin;
104926c55118SMichael Lange   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
105026c55118SMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
105126c55118SMichael Lange   /* Build a section of stratum values per point, generate the according SF
105226c55118SMichael Lange      and distribute point-wise stratum values to leaves. */
105326c55118SMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
105426c55118SMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
105526c55118SMichael Lange   ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr);
105626c55118SMichael Lange   if (label) {
105726c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1058ad8374ffSToby Isaac       const PetscInt *points;
1059ad8374ffSToby Isaac 
1060ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
106126c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1062ad8374ffSToby Isaac         ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr);
1063ad8374ffSToby Isaac         ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr);
106426c55118SMichael Lange       }
1065ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
106626c55118SMichael Lange     }
106726c55118SMichael Lange   }
106826c55118SMichael Lange   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
106926c55118SMichael Lange   /* Create a point-wise array of stratum values */
107026c55118SMichael Lange   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
107126c55118SMichael Lange   ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr);
107226c55118SMichael Lange   ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr);
107326c55118SMichael Lange   if (label) {
107426c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1075ad8374ffSToby Isaac       const PetscInt *points;
1076ad8374ffSToby Isaac 
1077ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
107826c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1079ad8374ffSToby Isaac         const PetscInt p = points[l];
108026c55118SMichael Lange         ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr);
108126c55118SMichael Lange         rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
108226c55118SMichael Lange       }
1083ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
108426c55118SMichael Lange     }
108526c55118SMichael Lange   }
108626c55118SMichael Lange   /* Build SF that maps label points to remote processes */
108726c55118SMichael Lange   ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr);
108826c55118SMichael Lange   ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr);
108926c55118SMichael Lange   ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr);
109026c55118SMichael Lange   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
109126c55118SMichael Lange   /* Send the strata for each point over the derived SF */
109226c55118SMichael Lange   ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr);
109326c55118SMichael Lange   ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr);
109426c55118SMichael Lange   ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
109526c55118SMichael Lange   ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
109626c55118SMichael Lange   /* Clean up */
109726c55118SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
109826c55118SMichael Lange   ierr = PetscFree(rootIdx);CHKERRQ(ierr);
109926c55118SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
110026c55118SMichael Lange   ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr);
110126c55118SMichael Lange   PetscFunctionReturn(0);
110226c55118SMichael Lange }
110326c55118SMichael Lange 
110484f0b6dfSMatthew G. Knepley /*@
110584f0b6dfSMatthew G. Knepley   DMLabelDistribute - Create a new label pushed forward over the PetscSF
110684f0b6dfSMatthew G. Knepley 
110784f0b6dfSMatthew G. Knepley   Input Parameters:
110884f0b6dfSMatthew G. Knepley + label - the DMLabel
110984f0b6dfSMatthew G. Knepley - sf    - the map from old to new distribution
111084f0b6dfSMatthew G. Knepley 
111184f0b6dfSMatthew G. Knepley   Output Parameter:
111284f0b6dfSMatthew G. Knepley . labelnew - the new redistributed label
111384f0b6dfSMatthew G. Knepley 
111484f0b6dfSMatthew G. Knepley   Level: intermediate
111584f0b6dfSMatthew G. Knepley 
111684f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
111784f0b6dfSMatthew G. Knepley @*/
1118c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1119c58f1c22SToby Isaac {
1120c58f1c22SToby Isaac   MPI_Comm       comm;
112126c55118SMichael Lange   PetscSection   leafSection;
112226c55118SMichael Lange   PetscInt       p, pStart, pEnd, s, size, dof, offset, stratum;
112326c55118SMichael Lange   PetscInt      *leafStrata, *strataIdx;
1124ad8374ffSToby Isaac   PetscInt     **points;
1125c58f1c22SToby Isaac   char          *name;
1126c58f1c22SToby Isaac   PetscInt       nameSize;
11275cbdf6fcSMichael Lange   PetscHashI     stratumHash;
11285cbdf6fcSMichael Lange   PETSC_UNUSED   PetscHashIIter ret, iter;
1129c58f1c22SToby Isaac   size_t         len = 0;
113026c55118SMichael Lange   PetscMPIInt    rank;
1131c58f1c22SToby Isaac   PetscErrorCode ierr;
1132c58f1c22SToby Isaac 
1133c58f1c22SToby Isaac   PetscFunctionBegin;
1134c58f1c22SToby Isaac   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
1135c58f1c22SToby Isaac   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1136c58f1c22SToby Isaac   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1137c58f1c22SToby Isaac   /* Bcast name */
1138c58f1c22SToby Isaac   if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);}
1139c58f1c22SToby Isaac   nameSize = len;
1140c58f1c22SToby Isaac   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1141c58f1c22SToby Isaac   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1142c58f1c22SToby Isaac   if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);}
1143c58f1c22SToby Isaac   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1144c58f1c22SToby Isaac   ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr);
1145c58f1c22SToby Isaac   ierr = PetscFree(name);CHKERRQ(ierr);
114677d236dfSMichael Lange   /* Bcast defaultValue */
114777d236dfSMichael Lange   if (!rank) (*labelNew)->defaultValue = label->defaultValue;
114877d236dfSMichael Lange   ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
114926c55118SMichael Lange   /* Distribute stratum values over the SF and get the point mapping on the receiver */
115026c55118SMichael Lange   ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr);
11515cbdf6fcSMichael Lange   /* Determine received stratum values and initialise new label*/
11525cbdf6fcSMichael Lange   PetscHashICreate(stratumHash);
115326c55118SMichael Lange   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
11545cbdf6fcSMichael Lange   for (p = 0; p < size; ++p) PetscHashIPut(stratumHash, leafStrata[p], ret, iter);
11555cbdf6fcSMichael Lange   PetscHashISize(stratumHash, (*labelNew)->numStrata);
1156ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr);
1157ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
11585cbdf6fcSMichael Lange   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr);
11595cbdf6fcSMichael Lange   /* Turn leafStrata into indices rather than stratum values */
11605cbdf6fcSMichael Lange   offset = 0;
11615cbdf6fcSMichael Lange   ierr = PetscHashIGetKeys(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr);
1162a205f1fdSToby Isaac   ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr);
11635cbdf6fcSMichael Lange   for (p = 0; p < size; ++p) {
1164231b9e6fSMatthew G. Knepley     for (s = 0; s < (*labelNew)->numStrata; ++s) {
1165231b9e6fSMatthew G. Knepley       if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
11665cbdf6fcSMichael Lange     }
11675cbdf6fcSMichael Lange   }
1168c58f1c22SToby Isaac   /* Rebuild the point strata on the receiver */
1169c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr);
1170c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1171c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1172c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1173c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1174c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1175c58f1c22SToby Isaac       (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1176c58f1c22SToby Isaac     }
1177c58f1c22SToby Isaac   }
1178c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr);
1179c58f1c22SToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr);
1180ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr);
1181c58f1c22SToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1182c58f1c22SToby Isaac     PetscHashICreate((*labelNew)->ht[s]);
1183ad8374ffSToby Isaac     ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr);
1184c58f1c22SToby Isaac   }
1185c58f1c22SToby Isaac   /* Insert points into new strata */
1186c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr);
1187c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1188c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1189c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1190c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1191c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1192c58f1c22SToby Isaac       stratum = leafStrata[offset+s];
1193ad8374ffSToby Isaac       points[stratum][strataIdx[stratum]++] = p;
1194c58f1c22SToby Isaac     }
1195c58f1c22SToby Isaac   }
1196ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; s++) {
1197ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr);
1198ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr);
1199ad8374ffSToby Isaac   }
1200ad8374ffSToby Isaac   ierr = PetscFree(points);CHKERRQ(ierr);
12015cbdf6fcSMichael Lange   PetscHashIDestroy(stratumHash);
1202c58f1c22SToby Isaac   ierr = PetscFree(leafStrata);CHKERRQ(ierr);
1203c58f1c22SToby Isaac   ierr = PetscFree(strataIdx);CHKERRQ(ierr);
1204c58f1c22SToby Isaac   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
1205c58f1c22SToby Isaac   PetscFunctionReturn(0);
1206c58f1c22SToby Isaac }
1207c58f1c22SToby Isaac 
12087937d9ceSMichael Lange /*@
12097937d9ceSMichael Lange   DMLabelGather - Gather all label values from leafs into roots
12107937d9ceSMichael Lange 
12117937d9ceSMichael Lange   Input Parameters:
12127937d9ceSMichael Lange + label - the DMLabel
121384f0b6dfSMatthew G. Knepley - sf - the Star Forest point communication map
12147937d9ceSMichael Lange 
121584f0b6dfSMatthew G. Knepley   Output Parameters:
121684f0b6dfSMatthew G. Knepley . labelNew - the new DMLabel with localised leaf values
12177937d9ceSMichael Lange 
12187937d9ceSMichael Lange   Level: developer
12197937d9ceSMichael Lange 
12207937d9ceSMichael Lange   Note: This is the inverse operation to DMLabelDistribute.
12217937d9ceSMichael Lange 
12227937d9ceSMichael Lange .seealso: DMLabelDistribute()
12237937d9ceSMichael Lange @*/
12247937d9ceSMichael Lange PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
12257937d9ceSMichael Lange {
12267937d9ceSMichael Lange   MPI_Comm       comm;
12277937d9ceSMichael Lange   PetscSection   rootSection;
12287937d9ceSMichael Lange   PetscSF        sfLabel;
12297937d9ceSMichael Lange   PetscSFNode   *rootPoints, *leafPoints;
12307937d9ceSMichael Lange   PetscInt       p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
12317937d9ceSMichael Lange   const PetscInt *rootDegree, *ilocal;
12327937d9ceSMichael Lange   PetscInt       *rootStrata;
12337937d9ceSMichael Lange   char          *name;
12347937d9ceSMichael Lange   PetscInt       nameSize;
12357937d9ceSMichael Lange   size_t         len = 0;
12369852e123SBarry Smith   PetscMPIInt    rank, size;
12377937d9ceSMichael Lange   PetscErrorCode ierr;
12387937d9ceSMichael Lange 
12397937d9ceSMichael Lange   PetscFunctionBegin;
12407937d9ceSMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
12417937d9ceSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
12429852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
12437937d9ceSMichael Lange   /* Bcast name */
12447937d9ceSMichael Lange   if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);}
12457937d9ceSMichael Lange   nameSize = len;
12467937d9ceSMichael Lange   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
12477937d9ceSMichael Lange   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
12487937d9ceSMichael Lange   if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);}
12497937d9ceSMichael Lange   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
12507937d9ceSMichael Lange   ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr);
12517937d9ceSMichael Lange   ierr = PetscFree(name);CHKERRQ(ierr);
12527937d9ceSMichael Lange   /* Gather rank/index pairs of leaves into local roots to build
12537937d9ceSMichael Lange      an inverse, multi-rooted SF. Note that this ignores local leaf
12547937d9ceSMichael Lange      indexing due to the use of the multiSF in PetscSFGather. */
12557937d9ceSMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr);
1256dc53bc9bSMatthew G. Knepley   ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr);
1257dc53bc9bSMatthew G. Knepley   for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
12587937d9ceSMichael Lange   for (p = 0; p < nleaves; p++) {
1259dc53bc9bSMatthew G. Knepley     leafPoints[ilocal[p]].index = ilocal[p];
1260dc53bc9bSMatthew G. Knepley     leafPoints[ilocal[p]].rank  = rank;
12617937d9ceSMichael Lange   }
12627937d9ceSMichael Lange   ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr);
12637937d9ceSMichael Lange   ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr);
12647937d9ceSMichael Lange   for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
12657937d9ceSMichael Lange   ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr);
12667937d9ceSMichael Lange   ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
12677937d9ceSMichael Lange   ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
12687937d9ceSMichael Lange   ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr);
12697937d9ceSMichael Lange   ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
12707937d9ceSMichael Lange   /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
12717937d9ceSMichael Lange   ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr);
12727937d9ceSMichael Lange   /* Rebuild the point strata on the receiver */
12737937d9ceSMichael Lange   for (p = 0, idx = 0; p < nroots; p++) {
12747937d9ceSMichael Lange     for (d = 0; d < rootDegree[p]; d++) {
12757937d9ceSMichael Lange       ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr);
12767937d9ceSMichael Lange       ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr);
12777937d9ceSMichael Lange       for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);}
12787937d9ceSMichael Lange     }
12797937d9ceSMichael Lange     idx += rootDegree[p];
12807937d9ceSMichael Lange   }
128177e0c0e7SMichael Lange   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
128277e0c0e7SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
128377e0c0e7SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
128477e0c0e7SMichael Lange   ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr);
12857937d9ceSMichael Lange   PetscFunctionReturn(0);
12867937d9ceSMichael Lange }
12877937d9ceSMichael Lange 
128884f0b6dfSMatthew G. Knepley /*@
128984f0b6dfSMatthew G. Knepley   DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label
129084f0b6dfSMatthew G. Knepley 
129184f0b6dfSMatthew G. Knepley   Input Parameter:
129284f0b6dfSMatthew G. Knepley . label - the DMLabel
129384f0b6dfSMatthew G. Knepley 
129484f0b6dfSMatthew G. Knepley   Output Parameters:
129584f0b6dfSMatthew G. Knepley + section - the section giving offsets for each stratum
129684f0b6dfSMatthew G. Knepley - is - An IS containing all the label points
129784f0b6dfSMatthew G. Knepley 
129884f0b6dfSMatthew G. Knepley   Level: developer
129984f0b6dfSMatthew G. Knepley 
130084f0b6dfSMatthew G. Knepley .seealso: DMLabelDistribute()
130184f0b6dfSMatthew G. Knepley @*/
1302c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1303c58f1c22SToby Isaac {
1304c58f1c22SToby Isaac   IS              vIS;
1305c58f1c22SToby Isaac   const PetscInt *values;
1306c58f1c22SToby Isaac   PetscInt       *points;
1307c58f1c22SToby Isaac   PetscInt        nV, vS = 0, vE = 0, v, N;
1308c58f1c22SToby Isaac   PetscErrorCode  ierr;
1309c58f1c22SToby Isaac 
1310c58f1c22SToby Isaac   PetscFunctionBegin;
1311c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr);
1312c58f1c22SToby Isaac   ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
1313c58f1c22SToby Isaac   ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr);
1314c58f1c22SToby Isaac   if (nV) {vS = values[0]; vE = values[0]+1;}
1315c58f1c22SToby Isaac   for (v = 1; v < nV; ++v) {
1316c58f1c22SToby Isaac     vS = PetscMin(vS, values[v]);
1317c58f1c22SToby Isaac     vE = PetscMax(vE, values[v]+1);
1318c58f1c22SToby Isaac   }
1319c58f1c22SToby Isaac   ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr);
1320c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr);
1321c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1322c58f1c22SToby Isaac     PetscInt n;
1323c58f1c22SToby Isaac 
1324c58f1c22SToby Isaac     ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr);
1325c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr);
1326c58f1c22SToby Isaac   }
1327c58f1c22SToby Isaac   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
1328c58f1c22SToby Isaac   ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr);
1329c58f1c22SToby Isaac   ierr = PetscMalloc1(N, &points);CHKERRQ(ierr);
1330c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1331c58f1c22SToby Isaac     IS              is;
1332c58f1c22SToby Isaac     const PetscInt *spoints;
1333c58f1c22SToby Isaac     PetscInt        dof, off, p;
1334c58f1c22SToby Isaac 
1335c58f1c22SToby Isaac     ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr);
1336c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr);
1337c58f1c22SToby Isaac     ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr);
1338c58f1c22SToby Isaac     ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr);
1339c58f1c22SToby Isaac     for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1340c58f1c22SToby Isaac     ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr);
1341c58f1c22SToby Isaac     ierr = ISDestroy(&is);CHKERRQ(ierr);
1342c58f1c22SToby Isaac   }
1343c58f1c22SToby Isaac   ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr);
1344c58f1c22SToby Isaac   ierr = ISDestroy(&vIS);CHKERRQ(ierr);
1345c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1346c58f1c22SToby Isaac   PetscFunctionReturn(0);
1347c58f1c22SToby Isaac }
1348c58f1c22SToby Isaac 
134984f0b6dfSMatthew G. Knepley /*@
1350c58f1c22SToby Isaac   PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1351c58f1c22SToby Isaac   the local section and an SF describing the section point overlap.
1352c58f1c22SToby Isaac 
1353c58f1c22SToby Isaac   Input Parameters:
1354c58f1c22SToby Isaac   + s - The PetscSection for the local field layout
1355c58f1c22SToby Isaac   . sf - The SF describing parallel layout of the section points
1356c58f1c22SToby Isaac   . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1357c58f1c22SToby Isaac   . label - The label specifying the points
1358c58f1c22SToby Isaac   - labelValue - The label stratum specifying the points
1359c58f1c22SToby Isaac 
1360c58f1c22SToby Isaac   Output Parameter:
1361c58f1c22SToby Isaac   . gsection - The PetscSection for the global field layout
1362c58f1c22SToby Isaac 
1363c58f1c22SToby Isaac   Note: This gives negative sizes and offsets to points not owned by this process
1364c58f1c22SToby Isaac 
1365c58f1c22SToby Isaac   Level: developer
1366c58f1c22SToby Isaac 
1367c58f1c22SToby Isaac .seealso: PetscSectionCreate()
1368c58f1c22SToby Isaac @*/
1369c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1370c58f1c22SToby Isaac {
1371c58f1c22SToby Isaac   PetscInt      *neg = NULL, *tmpOff = NULL;
1372c58f1c22SToby Isaac   PetscInt       pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1373c58f1c22SToby Isaac   PetscErrorCode ierr;
1374c58f1c22SToby Isaac 
1375c58f1c22SToby Isaac   PetscFunctionBegin;
1376c58f1c22SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr);
1377c58f1c22SToby Isaac   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1378c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr);
1379c58f1c22SToby Isaac   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1380c58f1c22SToby Isaac   if (nroots >= 0) {
1381c58f1c22SToby Isaac     if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1382c58f1c22SToby Isaac     ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr);
1383c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1384c58f1c22SToby Isaac       ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr);
1385c58f1c22SToby Isaac     } else {
1386c58f1c22SToby Isaac       tmpOff = &(*gsection)->atlasDof[-pStart];
1387c58f1c22SToby Isaac     }
1388c58f1c22SToby Isaac   }
1389c58f1c22SToby Isaac   /* Mark ghost points with negative dof */
1390c58f1c22SToby Isaac   for (p = pStart; p < pEnd; ++p) {
1391c58f1c22SToby Isaac     PetscInt value;
1392c58f1c22SToby Isaac 
1393c58f1c22SToby Isaac     ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr);
1394c58f1c22SToby Isaac     if (value != labelValue) continue;
1395c58f1c22SToby Isaac     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1396c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr);
1397c58f1c22SToby Isaac     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1398c58f1c22SToby Isaac     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);}
1399c58f1c22SToby Isaac     if (neg) neg[p] = -(dof+1);
1400c58f1c22SToby Isaac   }
1401c58f1c22SToby Isaac   ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr);
1402c58f1c22SToby Isaac   if (nroots >= 0) {
1403c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1404c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1405c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1406c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1407c58f1c22SToby Isaac     }
1408c58f1c22SToby Isaac   }
1409c58f1c22SToby Isaac   /* Calculate new sizes, get proccess offset, and calculate point offsets */
1410c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1411c58f1c22SToby Isaac     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1412c58f1c22SToby Isaac     (*gsection)->atlasOff[p] = off;
1413c58f1c22SToby Isaac     off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1414c58f1c22SToby Isaac   }
1415c58f1c22SToby Isaac   ierr       = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1416c58f1c22SToby Isaac   globalOff -= off;
1417c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1418c58f1c22SToby Isaac     (*gsection)->atlasOff[p] += globalOff;
1419c58f1c22SToby Isaac     if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1420c58f1c22SToby Isaac   }
1421c58f1c22SToby Isaac   /* Put in negative offsets for ghost points */
1422c58f1c22SToby Isaac   if (nroots >= 0) {
1423c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1424c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1425c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1426c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1427c58f1c22SToby Isaac     }
1428c58f1c22SToby Isaac   }
1429c58f1c22SToby Isaac   if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);}
1430c58f1c22SToby Isaac   ierr = PetscFree(neg);CHKERRQ(ierr);
1431c58f1c22SToby Isaac   PetscFunctionReturn(0);
1432c58f1c22SToby Isaac }
1433c58f1c22SToby Isaac 
14345fdea053SToby Isaac typedef struct _n_PetscSectionSym_Label
14355fdea053SToby Isaac {
14365fdea053SToby Isaac   DMLabel           label;
14375fdea053SToby Isaac   PetscCopyMode     *modes;
14385fdea053SToby Isaac   PetscInt          *sizes;
14395fdea053SToby Isaac   const PetscInt    ***perms;
14405fdea053SToby Isaac   const PetscScalar ***rots;
14415fdea053SToby Isaac   PetscInt          (*minMaxOrients)[2];
14425fdea053SToby Isaac   PetscInt          numStrata; /* numStrata is only increasing, functions as a state */
14435fdea053SToby Isaac } PetscSectionSym_Label;
14445fdea053SToby Isaac 
14455fdea053SToby Isaac static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
14465fdea053SToby Isaac {
14475fdea053SToby Isaac   PetscInt              i, j;
14485fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
14495fdea053SToby Isaac   PetscErrorCode        ierr;
14505fdea053SToby Isaac 
14515fdea053SToby Isaac   PetscFunctionBegin;
14525fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
14535fdea053SToby Isaac     if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
14545fdea053SToby Isaac       for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
14555fdea053SToby Isaac         if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);}
14565fdea053SToby Isaac         if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);}
14575fdea053SToby Isaac       }
14585fdea053SToby Isaac       if (sl->perms[i]) {
14595fdea053SToby Isaac         const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
14605fdea053SToby Isaac 
14615fdea053SToby Isaac         ierr = PetscFree(perms);CHKERRQ(ierr);
14625fdea053SToby Isaac       }
14635fdea053SToby Isaac       if (sl->rots[i]) {
14645fdea053SToby Isaac         const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
14655fdea053SToby Isaac 
14665fdea053SToby Isaac         ierr = PetscFree(rots);CHKERRQ(ierr);
14675fdea053SToby Isaac       }
14685fdea053SToby Isaac     }
14695fdea053SToby Isaac   }
14705fdea053SToby Isaac   ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr);
14715fdea053SToby Isaac   ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr);
14725fdea053SToby Isaac   sl->numStrata = 0;
14735fdea053SToby Isaac   PetscFunctionReturn(0);
14745fdea053SToby Isaac }
14755fdea053SToby Isaac 
14765fdea053SToby Isaac static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
14775fdea053SToby Isaac {
14785fdea053SToby Isaac   PetscErrorCode ierr;
14795fdea053SToby Isaac 
14805fdea053SToby Isaac   PetscFunctionBegin;
14815fdea053SToby Isaac   ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);
14825fdea053SToby Isaac   ierr = PetscFree(sym->data);CHKERRQ(ierr);
14835fdea053SToby Isaac   PetscFunctionReturn(0);
14845fdea053SToby Isaac }
14855fdea053SToby Isaac 
14865fdea053SToby Isaac static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
14875fdea053SToby Isaac {
14885fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
14895fdea053SToby Isaac   PetscBool             isAscii;
14905fdea053SToby Isaac   DMLabel               label = sl->label;
14915fdea053SToby Isaac   PetscErrorCode        ierr;
14925fdea053SToby Isaac 
14935fdea053SToby Isaac   PetscFunctionBegin;
14945fdea053SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
14955fdea053SToby Isaac   if (isAscii) {
14965fdea053SToby Isaac     PetscInt          i, j, k;
14975fdea053SToby Isaac     PetscViewerFormat format;
14985fdea053SToby Isaac 
14995fdea053SToby Isaac     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
15005fdea053SToby Isaac     if (label) {
15015fdea053SToby Isaac       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
15025fdea053SToby Isaac       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
15035fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15045fdea053SToby Isaac         ierr = DMLabelView(label, viewer);CHKERRQ(ierr);
15055fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
15065fdea053SToby Isaac       } else {
15075fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer,"  Label '%s'\n",sl->label->name);CHKERRQ(ierr);
15085fdea053SToby Isaac       }
15095fdea053SToby Isaac     } else {
15105fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr);
15115fdea053SToby Isaac     }
15125fdea053SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15135fdea053SToby Isaac     for (i = 0; i <= sl->numStrata; i++) {
15145fdea053SToby Isaac       PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
15155fdea053SToby Isaac 
15165fdea053SToby Isaac       if (!(sl->perms[i] || sl->rots[i])) {
15175fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr);
15185fdea053SToby Isaac       } else {
15195fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr);
15205fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15215fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr);
15225fdea053SToby Isaac         if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
15235fdea053SToby Isaac           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15245fdea053SToby Isaac           for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
15255fdea053SToby Isaac             if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
15265fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr);
15275fdea053SToby Isaac             } else {
15285fdea053SToby Isaac               PetscInt tab;
15295fdea053SToby Isaac 
15305fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr);
15315fdea053SToby Isaac               ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15325fdea053SToby Isaac               ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr);
15335fdea053SToby Isaac               if (sl->perms[i] && sl->perms[i][j]) {
15345fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr);
15355fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
15365fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);}
15375fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
15385fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
15395fdea053SToby Isaac               }
15405fdea053SToby Isaac               if (sl->rots[i] && sl->rots[i][j]) {
15415fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Rotations:  ");CHKERRQ(ierr);
15425fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
15435fdea053SToby Isaac #if defined(PETSC_USE_COMPLEX)
15445fdea053SToby 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);}
15455fdea053SToby Isaac #else
15465fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);}
15475fdea053SToby Isaac #endif
15485fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
15495fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
15505fdea053SToby Isaac               }
15515fdea053SToby Isaac               ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
15525fdea053SToby Isaac             }
15535fdea053SToby Isaac           }
15545fdea053SToby Isaac           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
15555fdea053SToby Isaac         }
15565fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
15575fdea053SToby Isaac       }
15585fdea053SToby Isaac     }
15595fdea053SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
15605fdea053SToby Isaac   }
15615fdea053SToby Isaac   PetscFunctionReturn(0);
15625fdea053SToby Isaac }
15635fdea053SToby Isaac 
15645fdea053SToby Isaac /*@
15655fdea053SToby Isaac   PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
15665fdea053SToby Isaac 
15675fdea053SToby Isaac   Logically collective on sym
15685fdea053SToby Isaac 
15695fdea053SToby Isaac   Input parameters:
15705fdea053SToby Isaac + sym - the section symmetries
15715fdea053SToby Isaac - label - the DMLabel describing the types of points
15725fdea053SToby Isaac 
15735fdea053SToby Isaac   Level: developer:
15745fdea053SToby Isaac 
15755fdea053SToby Isaac .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()
15765fdea053SToby Isaac @*/
15775fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
15785fdea053SToby Isaac {
15795fdea053SToby Isaac   PetscSectionSym_Label *sl;
15805fdea053SToby Isaac   PetscErrorCode ierr;
15815fdea053SToby Isaac 
15825fdea053SToby Isaac   PetscFunctionBegin;
15835fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
15845fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
15855fdea053SToby Isaac   if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);}
15865fdea053SToby Isaac   if (label) {
15875fdea053SToby Isaac     label->refct++;
15885fdea053SToby Isaac     sl->label = label;
15895fdea053SToby Isaac     ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr);
15901a834cf9SToby 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);
15911a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr);
15921a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr);
15931a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr);
15941a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr);
15951a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr);
15965fdea053SToby Isaac   }
15975fdea053SToby Isaac   PetscFunctionReturn(0);
15985fdea053SToby Isaac }
15995fdea053SToby Isaac 
16005fdea053SToby Isaac /*@C
16015fdea053SToby Isaac   PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
16025fdea053SToby Isaac 
16035fdea053SToby Isaac   Logically collective on PetscSectionSym
16045fdea053SToby Isaac 
16055fdea053SToby Isaac   InputParameters:
16065fdea053SToby Isaac + sys       - the section symmetries
16075fdea053SToby Isaac . stratum   - the stratum value in the label that we are assigning symmetries for
16085fdea053SToby Isaac . size      - the number of dofs for points in the stratum of the label
16095fdea053SToby Isaac . minOrient - the smallest orientation for a point in this stratum
16105fdea053SToby Isaac . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
16115fdea053SToby Isaac . mode      - how sym should copy the perms and rots arrays
16125fdea053SToby Isaac . perms     - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation.  A NULL permutation is the identity
16135fdea053SToby 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
16145fdea053SToby Isaac 
16155fdea053SToby Isaac   Level: developer
16165fdea053SToby Isaac 
16175fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()
16185fdea053SToby Isaac @*/
16195fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
16205fdea053SToby Isaac {
16215fdea053SToby Isaac   PetscInt       i, j, k;
16225fdea053SToby Isaac   PetscSectionSym_Label *sl;
16235fdea053SToby Isaac   PetscErrorCode ierr;
16245fdea053SToby Isaac 
16255fdea053SToby Isaac   PetscFunctionBegin;
16265fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
16275fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
16285fdea053SToby Isaac   if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet");
16295fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
16305fdea053SToby Isaac     PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
16315fdea053SToby Isaac 
16325fdea053SToby Isaac     if (stratum == value) break;
16335fdea053SToby Isaac   }
16345fdea053SToby 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);
16355fdea053SToby Isaac   sl->sizes[i] = size;
16365fdea053SToby Isaac   sl->modes[i] = mode;
16375fdea053SToby Isaac   sl->minMaxOrients[i][0] = minOrient;
16385fdea053SToby Isaac   sl->minMaxOrients[i][1] = maxOrient;
16395fdea053SToby Isaac   if (mode == PETSC_COPY_VALUES) {
16405fdea053SToby Isaac     if (perms) {
16415fdea053SToby Isaac       PetscInt    **ownPerms;
16425fdea053SToby Isaac 
16435fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr);
16445fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
16455fdea053SToby Isaac         if (perms[j]) {
16465fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr);
16475fdea053SToby Isaac           for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];}
16485fdea053SToby Isaac         }
16495fdea053SToby Isaac       }
16505fdea053SToby Isaac       sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient];
16515fdea053SToby Isaac     }
16525fdea053SToby Isaac     if (rots) {
16535fdea053SToby Isaac       PetscScalar **ownRots;
16545fdea053SToby Isaac 
16555fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr);
16565fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
16575fdea053SToby Isaac         if (rots[j]) {
16585fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr);
16595fdea053SToby Isaac           for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];}
16605fdea053SToby Isaac         }
16615fdea053SToby Isaac       }
16625fdea053SToby Isaac       sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient];
16635fdea053SToby Isaac     }
16645fdea053SToby Isaac   } else {
16655fdea053SToby Isaac     sl->perms[i] = perms ? &perms[-minOrient] : NULL;
16665fdea053SToby Isaac     sl->rots[i]  = rots ? &rots[-minOrient] : NULL;
16675fdea053SToby Isaac   }
16685fdea053SToby Isaac   PetscFunctionReturn(0);
16695fdea053SToby Isaac }
16705fdea053SToby Isaac 
16715fdea053SToby Isaac static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
16725fdea053SToby Isaac {
16735fdea053SToby Isaac   PetscInt              i, j, numStrata;
16745fdea053SToby Isaac   PetscSectionSym_Label *sl;
16755fdea053SToby Isaac   DMLabel               label;
16765fdea053SToby Isaac   PetscErrorCode        ierr;
16775fdea053SToby Isaac 
16785fdea053SToby Isaac   PetscFunctionBegin;
16795fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
16805fdea053SToby Isaac   numStrata = sl->numStrata;
16815fdea053SToby Isaac   label     = sl->label;
16825fdea053SToby Isaac   for (i = 0; i < numPoints; i++) {
16835fdea053SToby Isaac     PetscInt point = points[2*i];
16845fdea053SToby Isaac     PetscInt ornt  = points[2*i+1];
16855fdea053SToby Isaac 
16865fdea053SToby Isaac     for (j = 0; j < numStrata; j++) {
16875fdea053SToby Isaac       if (label->validIS[j]) {
16885fdea053SToby Isaac         PetscInt k;
16895fdea053SToby Isaac 
16905fdea053SToby Isaac         ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr);
16915fdea053SToby Isaac         if (k >= 0) break;
16925fdea053SToby Isaac       } else {
16935fdea053SToby Isaac         PetscBool has;
16945fdea053SToby Isaac 
16955fdea053SToby Isaac         PetscHashIHasKey(label->ht[j], point, has);
16965fdea053SToby Isaac         if (has) break;
16975fdea053SToby Isaac       }
16985fdea053SToby Isaac     }
16995fdea053SToby 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);
17005fdea053SToby Isaac     if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;}
17015fdea053SToby Isaac     if (rots) {rots[i]  = sl->rots[j] ? sl->rots[j][ornt] : NULL;}
17025fdea053SToby Isaac   }
17035fdea053SToby Isaac   PetscFunctionReturn(0);
17045fdea053SToby Isaac }
17055fdea053SToby Isaac 
17065fdea053SToby Isaac PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
17075fdea053SToby Isaac {
17085fdea053SToby Isaac   PetscSectionSym_Label *sl;
17095fdea053SToby Isaac   PetscErrorCode        ierr;
17105fdea053SToby Isaac 
17115fdea053SToby Isaac   PetscFunctionBegin;
17125fdea053SToby Isaac   ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr);
17135fdea053SToby Isaac   sym->ops->getpoints = PetscSectionSymGetPoints_Label;
17145fdea053SToby Isaac   sym->ops->view      = PetscSectionSymView_Label;
17155fdea053SToby Isaac   sym->ops->destroy   = PetscSectionSymDestroy_Label;
17165fdea053SToby Isaac   sym->data           = (void *) sl;
17175fdea053SToby Isaac   PetscFunctionReturn(0);
17185fdea053SToby Isaac }
17195fdea053SToby Isaac 
17205fdea053SToby Isaac /*@
17215fdea053SToby Isaac   PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
17225fdea053SToby Isaac 
17235fdea053SToby Isaac   Collective
17245fdea053SToby Isaac 
17255fdea053SToby Isaac   Input Parameters:
17265fdea053SToby Isaac + comm - the MPI communicator for the new symmetry
17275fdea053SToby Isaac - label - the label defining the strata
17285fdea053SToby Isaac 
17295fdea053SToby Isaac   Output Parameters:
17305fdea053SToby Isaac . sym - the section symmetries
17315fdea053SToby Isaac 
17325fdea053SToby Isaac   Level: developer
17335fdea053SToby Isaac 
17345fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()
17355fdea053SToby Isaac @*/
17365fdea053SToby Isaac PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
17375fdea053SToby Isaac {
17385fdea053SToby Isaac   PetscErrorCode        ierr;
17395fdea053SToby Isaac 
17405fdea053SToby Isaac   PetscFunctionBegin;
17415fdea053SToby Isaac   ierr = DMInitializePackage();CHKERRQ(ierr);
17425fdea053SToby Isaac   ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr);
17435fdea053SToby Isaac   ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr);
17445fdea053SToby Isaac   ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr);
17455fdea053SToby Isaac   PetscFunctionReturn(0);
17465fdea053SToby Isaac }
1747