xref: /petsc/src/dm/label/dmlabel.c (revision c6a43d28382f502ae1ec57d51dd44f0da9601430)
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 
9d67d17b1SMatthew G. Knepley   Input parameters:
10d67d17b1SMatthew G. Knepley + comm - The communicator, usually PETSC_COMM_SELF
11d67d17b1SMatthew G. Knepley - name - The label name
12c58f1c22SToby Isaac 
13c58f1c22SToby Isaac   Output parameter:
14c58f1c22SToby Isaac . label - The DMLabel
15c58f1c22SToby Isaac 
16c58f1c22SToby Isaac   Level: beginner
17c58f1c22SToby Isaac 
18c58f1c22SToby Isaac .seealso: DMLabelDestroy()
19c58f1c22SToby Isaac @*/
20d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelCreate(MPI_Comm comm, const char name[], DMLabel *label)
21c58f1c22SToby Isaac {
22c58f1c22SToby Isaac   PetscErrorCode ierr;
23c58f1c22SToby Isaac 
24c58f1c22SToby Isaac   PetscFunctionBegin;
25d67d17b1SMatthew G. Knepley   PetscValidPointer(label,2);
26d67d17b1SMatthew G. Knepley   ierr = DMInitializePackage();CHKERRQ(ierr);
27c58f1c22SToby Isaac 
28d67d17b1SMatthew G. Knepley   ierr = PetscHeaderCreate(*label,DMLABEL_CLASSID,"DMLabel","DMLabel","DM",comm,DMLabelDestroy,DMLabelView);CHKERRQ(ierr);
29d67d17b1SMatthew G. Knepley 
30c58f1c22SToby Isaac   (*label)->numStrata      = 0;
315aa44df4SToby Isaac   (*label)->defaultValue   = -1;
32c58f1c22SToby Isaac   (*label)->stratumValues  = NULL;
33ad8374ffSToby Isaac   (*label)->validIS        = NULL;
34c58f1c22SToby Isaac   (*label)->stratumSizes   = NULL;
35c58f1c22SToby Isaac   (*label)->points         = NULL;
36c58f1c22SToby Isaac   (*label)->ht             = NULL;
37c58f1c22SToby Isaac   (*label)->pStart         = -1;
38c58f1c22SToby Isaac   (*label)->pEnd           = -1;
39c58f1c22SToby Isaac   (*label)->bt             = NULL;
40d67d17b1SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *label, name);CHKERRQ(ierr);
41c58f1c22SToby Isaac   PetscFunctionReturn(0);
42c58f1c22SToby Isaac }
43c58f1c22SToby Isaac 
44c58f1c22SToby Isaac /*
45c58f1c22SToby Isaac   DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format
46c58f1c22SToby Isaac 
47c58f1c22SToby Isaac   Input parameter:
48c58f1c22SToby Isaac + label - The DMLabel
49c58f1c22SToby Isaac - v - The stratum value
50c58f1c22SToby Isaac 
51c58f1c22SToby Isaac   Output parameter:
52c58f1c22SToby Isaac . label - The DMLabel with stratum in sorted list format
53c58f1c22SToby Isaac 
54c58f1c22SToby Isaac   Level: developer
55c58f1c22SToby Isaac 
56c58f1c22SToby Isaac .seealso: DMLabelCreate()
57c58f1c22SToby Isaac */
58c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v)
59c58f1c22SToby Isaac {
600c3c4a36SLisandro Dalcin   PetscInt       off = 0;
61ad8374ffSToby Isaac   PetscInt       *pointArray;
62c58f1c22SToby Isaac   PetscErrorCode ierr;
63c58f1c22SToby Isaac 
64ad8374ffSToby Isaac   if (label->validIS[v]) return 0;
65c58f1c22SToby Isaac   PetscFunctionBegin;
660c3c4a36SLisandro 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);
67e8f14785SLisandro Dalcin   ierr = PetscHSetIGetSize(label->ht[v], &label->stratumSizes[v]);CHKERRQ(ierr);
68c58f1c22SToby Isaac 
69ad8374ffSToby Isaac   ierr = PetscMalloc1(label->stratumSizes[v], &pointArray);CHKERRQ(ierr);
70e8f14785SLisandro Dalcin   ierr = PetscHSetIGetElems(label->ht[v], &off, pointArray);CHKERRQ(ierr);
71c58f1c22SToby 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]);
72e8f14785SLisandro Dalcin   PetscHSetIClear(label->ht[v]);
73ad8374ffSToby Isaac   ierr = PetscSortInt(label->stratumSizes[v], pointArray);CHKERRQ(ierr);
74c58f1c22SToby Isaac   if (label->bt) {
75c58f1c22SToby Isaac     PetscInt p;
76c58f1c22SToby Isaac 
77c58f1c22SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
78ad8374ffSToby Isaac       const PetscInt point = pointArray[p];
79c58f1c22SToby Isaac 
80c58f1c22SToby 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);
81c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
82c58f1c22SToby Isaac     }
83c58f1c22SToby Isaac   }
84ad8374ffSToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF,label->stratumSizes[v],pointArray,PETSC_OWN_POINTER,&(label->points[v]));CHKERRQ(ierr);
85ad8374ffSToby Isaac   ierr = PetscObjectSetName((PetscObject) (label->points[v]), "indices");CHKERRQ(ierr);
86ad8374ffSToby Isaac   label->validIS[v] = PETSC_TRUE;
87d67d17b1SMatthew G. Knepley   ierr = PetscObjectStateIncrease((PetscObject) label);CHKERRQ(ierr);
88c58f1c22SToby Isaac   PetscFunctionReturn(0);
89c58f1c22SToby Isaac }
90c58f1c22SToby Isaac 
91c58f1c22SToby Isaac /*
92c58f1c22SToby Isaac   DMLabelMakeAllValid_Private - Transfer all strata from the hash format to the sorted list format
93c58f1c22SToby Isaac 
94c58f1c22SToby Isaac   Input parameter:
95c58f1c22SToby Isaac . label - The DMLabel
96c58f1c22SToby Isaac 
97c58f1c22SToby Isaac   Output parameter:
98c58f1c22SToby Isaac . label - The DMLabel with all strata in sorted list format
99c58f1c22SToby Isaac 
100c58f1c22SToby Isaac   Level: developer
101c58f1c22SToby Isaac 
102c58f1c22SToby Isaac .seealso: DMLabelCreate()
103c58f1c22SToby Isaac */
104c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeAllValid_Private(DMLabel label)
105c58f1c22SToby Isaac {
106c58f1c22SToby Isaac   PetscInt       v;
107c58f1c22SToby Isaac   PetscErrorCode ierr;
108c58f1c22SToby Isaac 
109c58f1c22SToby Isaac   PetscFunctionBegin;
110c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; v++){
111c58f1c22SToby Isaac     ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
112c58f1c22SToby Isaac   }
113c58f1c22SToby Isaac   PetscFunctionReturn(0);
114c58f1c22SToby Isaac }
115c58f1c22SToby Isaac 
116c58f1c22SToby Isaac /*
117c58f1c22SToby Isaac   DMLabelMakeInvalid_Private - Transfer stratum data from the sorted list format to the hash format
118c58f1c22SToby Isaac 
119c58f1c22SToby Isaac   Input parameter:
120c58f1c22SToby Isaac + label - The DMLabel
121c58f1c22SToby Isaac - v - The stratum value
122c58f1c22SToby Isaac 
123c58f1c22SToby Isaac   Output parameter:
124c58f1c22SToby Isaac . label - The DMLabel with stratum in hash format
125c58f1c22SToby Isaac 
126c58f1c22SToby Isaac   Level: developer
127c58f1c22SToby Isaac 
128c58f1c22SToby Isaac .seealso: DMLabelCreate()
129c58f1c22SToby Isaac */
130c58f1c22SToby Isaac static PetscErrorCode DMLabelMakeInvalid_Private(DMLabel label, PetscInt v)
131c58f1c22SToby Isaac {
132c58f1c22SToby Isaac   PetscInt       p;
133ad8374ffSToby Isaac   const PetscInt *points;
134c58f1c22SToby Isaac   PetscErrorCode ierr;
135c58f1c22SToby Isaac 
1360c3c4a36SLisandro Dalcin   if (!label->validIS[v]) return 0;
137c58f1c22SToby Isaac   PetscFunctionBegin;
1380c3c4a36SLisandro 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);
139ad8374ffSToby Isaac   if (label->points[v]) {
140ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
141e8f14785SLisandro Dalcin     for (p = 0; p < label->stratumSizes[v]; ++p) {
142e8f14785SLisandro Dalcin       ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);
143e8f14785SLisandro Dalcin     }
144ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
145ad8374ffSToby Isaac     ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
146ad8374ffSToby Isaac   }
147ad8374ffSToby Isaac   label->validIS[v] = PETSC_FALSE;
148c58f1c22SToby Isaac   PetscFunctionReturn(0);
149c58f1c22SToby Isaac }
150c58f1c22SToby Isaac 
1510c3c4a36SLisandro Dalcin PETSC_STATIC_INLINE PetscErrorCode DMLabelLookupStratum(DMLabel label, PetscInt value, PetscInt *index)
1520c3c4a36SLisandro Dalcin {
1530c3c4a36SLisandro Dalcin   PetscInt v;
1540e79e033SBarry Smith 
1550c3c4a36SLisandro Dalcin   PetscFunctionBegin;
1560e79e033SBarry Smith   *index = -1;
157d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1580e79e033SBarry Smith   for (v = 0; v < label->numStrata; ++v) {
1590e79e033SBarry Smith     if (label->stratumValues[v] == value) {
1600e79e033SBarry Smith       *index = v; break;
1610e79e033SBarry Smith     }
1620e79e033SBarry Smith   }
1630c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
1640c3c4a36SLisandro Dalcin }
1650c3c4a36SLisandro Dalcin 
1660c3c4a36SLisandro Dalcin static PetscErrorCode DMLabelNewStratum(DMLabel label, PetscInt value, PetscInt *index)
167c58f1c22SToby Isaac {
168ad8374ffSToby Isaac   PetscInt       v, *tmpV, *tmpS;
169ad8374ffSToby Isaac   IS             *tmpP;
170e8f14785SLisandro Dalcin   PetscHSetI     *tmpH;
171c58f1c22SToby Isaac   PetscBool      *tmpB;
172c58f1c22SToby Isaac   PetscErrorCode ierr;
173c58f1c22SToby Isaac 
174c58f1c22SToby Isaac   PetscFunctionBegin;
175d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
176c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpV);CHKERRQ(ierr);
177c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpS);CHKERRQ(ierr);
178c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpH);CHKERRQ(ierr);
179c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpP);CHKERRQ(ierr);
180c58f1c22SToby Isaac   ierr = PetscMalloc1((label->numStrata+1), &tmpB);CHKERRQ(ierr);
181c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
182c58f1c22SToby Isaac     tmpV[v] = label->stratumValues[v];
183c58f1c22SToby Isaac     tmpS[v] = label->stratumSizes[v];
184c58f1c22SToby Isaac     tmpH[v] = label->ht[v];
185c58f1c22SToby Isaac     tmpP[v] = label->points[v];
186ad8374ffSToby Isaac     tmpB[v] = label->validIS[v];
187c58f1c22SToby Isaac   }
188c58f1c22SToby Isaac   tmpV[v] = value;
189c58f1c22SToby Isaac   tmpS[v] = 0;
190e8f14785SLisandro Dalcin   ierr = PetscHSetICreate(&tmpH[v]);CHKERRQ(ierr);
191ecb9a371SToby Isaac   ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&tmpP[v]);CHKERRQ(ierr);
192c58f1c22SToby Isaac   tmpB[v] = PETSC_TRUE;
1930c3c4a36SLisandro Dalcin 
194c58f1c22SToby Isaac   ++label->numStrata;
195c58f1c22SToby Isaac   ierr = PetscFree(label->stratumValues);CHKERRQ(ierr);
196c58f1c22SToby Isaac   ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr);
197c58f1c22SToby Isaac   ierr = PetscFree(label->ht);CHKERRQ(ierr);
198c58f1c22SToby Isaac   ierr = PetscFree(label->points);CHKERRQ(ierr);
199ad8374ffSToby Isaac   ierr = PetscFree(label->validIS);CHKERRQ(ierr);
200c58f1c22SToby Isaac   label->stratumValues = tmpV;
201c58f1c22SToby Isaac   label->stratumSizes  = tmpS;
202c58f1c22SToby Isaac   label->ht            = tmpH;
203c58f1c22SToby Isaac   label->points        = tmpP;
204ad8374ffSToby Isaac   label->validIS       = tmpB;
205c58f1c22SToby Isaac 
2060c3c4a36SLisandro Dalcin   *index = v;
2070c3c4a36SLisandro Dalcin   PetscFunctionReturn(0);
2080c3c4a36SLisandro Dalcin }
2090c3c4a36SLisandro Dalcin 
2100c3c4a36SLisandro Dalcin PetscErrorCode DMLabelAddStratum(DMLabel label, PetscInt value)
2110c3c4a36SLisandro Dalcin {
2120c3c4a36SLisandro Dalcin   PetscInt       v;
2130c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
2140c3c4a36SLisandro Dalcin 
2150c3c4a36SLisandro Dalcin   PetscFunctionBegin;
216d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
2170c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
2180c3c4a36SLisandro Dalcin   if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);}
219c58f1c22SToby Isaac   PetscFunctionReturn(0);
220c58f1c22SToby Isaac }
221c58f1c22SToby Isaac 
222c58f1c22SToby Isaac static PetscErrorCode DMLabelView_Ascii(DMLabel label, PetscViewer viewer)
223c58f1c22SToby Isaac {
224c58f1c22SToby Isaac   PetscInt       v;
225c58f1c22SToby Isaac   PetscMPIInt    rank;
226c58f1c22SToby Isaac   PetscErrorCode ierr;
227c58f1c22SToby Isaac 
228c58f1c22SToby Isaac   PetscFunctionBegin;
229c58f1c22SToby Isaac   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr);
230c58f1c22SToby Isaac   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
231c58f1c22SToby Isaac   if (label) {
232d67d17b1SMatthew G. Knepley     const char *name;
233d67d17b1SMatthew G. Knepley 
234d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
235d67d17b1SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "Label '%s':\n", name);CHKERRQ(ierr);
236c58f1c22SToby Isaac     if (label->bt) {ierr = PetscViewerASCIIPrintf(viewer, "  Index has been calculated in [%D, %D)\n", label->pStart, label->pEnd);CHKERRQ(ierr);}
237c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
238c58f1c22SToby Isaac       const PetscInt value = label->stratumValues[v];
239ad8374ffSToby Isaac       const PetscInt *points;
240c58f1c22SToby Isaac       PetscInt       p;
241c58f1c22SToby Isaac 
242ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v],&points);CHKERRQ(ierr);
243c58f1c22SToby Isaac       for (p = 0; p < label->stratumSizes[v]; ++p) {
244ad8374ffSToby Isaac         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D (%D)\n", rank, points[p], value);CHKERRQ(ierr);
245c58f1c22SToby Isaac       }
246ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
247c58f1c22SToby Isaac     }
248c58f1c22SToby Isaac   }
249c58f1c22SToby Isaac   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
250c58f1c22SToby Isaac   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
251c58f1c22SToby Isaac   PetscFunctionReturn(0);
252c58f1c22SToby Isaac }
253c58f1c22SToby Isaac 
254c58f1c22SToby Isaac /*@C
255c58f1c22SToby Isaac   DMLabelView - View the label
256c58f1c22SToby Isaac 
257c58f1c22SToby Isaac   Input Parameters:
258c58f1c22SToby Isaac + label - The DMLabel
259c58f1c22SToby Isaac - viewer - The PetscViewer
260c58f1c22SToby Isaac 
261c58f1c22SToby Isaac   Level: intermediate
262c58f1c22SToby Isaac 
263c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelDestroy()
264c58f1c22SToby Isaac @*/
265c58f1c22SToby Isaac PetscErrorCode DMLabelView(DMLabel label, PetscViewer viewer)
266c58f1c22SToby Isaac {
267c58f1c22SToby Isaac   PetscBool      iascii;
268c58f1c22SToby Isaac   PetscErrorCode ierr;
269c58f1c22SToby Isaac 
270c58f1c22SToby Isaac   PetscFunctionBegin;
271d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
272c58f1c22SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
273c58f1c22SToby Isaac   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
274c58f1c22SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
275c58f1c22SToby Isaac   if (iascii) {
276c58f1c22SToby Isaac     ierr = DMLabelView_Ascii(label, viewer);CHKERRQ(ierr);
277c58f1c22SToby Isaac   }
278c58f1c22SToby Isaac   PetscFunctionReturn(0);
279c58f1c22SToby Isaac }
280c58f1c22SToby Isaac 
28184f0b6dfSMatthew G. Knepley /*@
282d67d17b1SMatthew G. Knepley   DMLabelReset - Destroys internal data structures in a DMLabel
283d67d17b1SMatthew G. Knepley 
284d67d17b1SMatthew G. Knepley   Input Parameter:
285d67d17b1SMatthew G. Knepley . label - The DMLabel
286d67d17b1SMatthew G. Knepley 
287d67d17b1SMatthew G. Knepley   Level: beginner
288d67d17b1SMatthew G. Knepley 
289d67d17b1SMatthew G. Knepley .seealso: DMLabelDestroy(), DMLabelCreate()
290d67d17b1SMatthew G. Knepley @*/
291d67d17b1SMatthew G. Knepley PetscErrorCode DMLabelReset(DMLabel label)
292d67d17b1SMatthew G. Knepley {
293d67d17b1SMatthew G. Knepley   PetscInt       v;
294d67d17b1SMatthew G. Knepley   PetscErrorCode ierr;
295d67d17b1SMatthew G. Knepley 
296d67d17b1SMatthew G. Knepley   PetscFunctionBegin;
297d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
298d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->stratumValues);CHKERRQ(ierr);
299d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->stratumSizes);CHKERRQ(ierr);
300d67d17b1SMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
301d67d17b1SMatthew G. Knepley     ierr = PetscHSetIDestroy(&label->ht[v]);CHKERRQ(ierr);
302d67d17b1SMatthew G. Knepley     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
303d67d17b1SMatthew G. Knepley   }
304d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->ht);CHKERRQ(ierr);
305d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->points);CHKERRQ(ierr);
306d67d17b1SMatthew G. Knepley   ierr = PetscFree(label->validIS);CHKERRQ(ierr);
307d67d17b1SMatthew G. Knepley   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
308d67d17b1SMatthew G. Knepley   PetscFunctionReturn(0);
309d67d17b1SMatthew G. Knepley }
310d67d17b1SMatthew G. Knepley 
311d67d17b1SMatthew G. Knepley /*@
31284f0b6dfSMatthew G. Knepley   DMLabelDestroy - Destroys a DMLabel
31384f0b6dfSMatthew G. Knepley 
31484f0b6dfSMatthew G. Knepley   Input Parameter:
31584f0b6dfSMatthew G. Knepley . label - The DMLabel
31684f0b6dfSMatthew G. Knepley 
31784f0b6dfSMatthew G. Knepley   Level: beginner
31884f0b6dfSMatthew G. Knepley 
319d67d17b1SMatthew G. Knepley .seealso: DMLabelReset(), DMLabelCreate()
32084f0b6dfSMatthew G. Knepley @*/
321c58f1c22SToby Isaac PetscErrorCode DMLabelDestroy(DMLabel *label)
322c58f1c22SToby Isaac {
323c58f1c22SToby Isaac   PetscErrorCode ierr;
324c58f1c22SToby Isaac 
325c58f1c22SToby Isaac   PetscFunctionBegin;
326d67d17b1SMatthew G. Knepley   if (!*label) PetscFunctionReturn(0);
327d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific((*label),DMLABEL_CLASSID,1);
328d67d17b1SMatthew G. Knepley   if (--((PetscObject)(*label))->refct > 0) {
329d67d17b1SMatthew G. Knepley     *label = NULL;
330d67d17b1SMatthew G. Knepley     PetscFunctionReturn(0);
3310c3c4a36SLisandro Dalcin   }
332d67d17b1SMatthew G. Knepley   ierr = DMLabelReset(*label);CHKERRQ(ierr);
333d67d17b1SMatthew G. Knepley   ierr = PetscHeaderDestroy(label);CHKERRQ(ierr);
334c58f1c22SToby Isaac   PetscFunctionReturn(0);
335c58f1c22SToby Isaac }
336c58f1c22SToby Isaac 
33784f0b6dfSMatthew G. Knepley /*@
33884f0b6dfSMatthew G. Knepley   DMLabelDuplicate - Duplicates a DMLabel
33984f0b6dfSMatthew G. Knepley 
34084f0b6dfSMatthew G. Knepley   Input Parameter:
34184f0b6dfSMatthew G. Knepley . label - The DMLabel
34284f0b6dfSMatthew G. Knepley 
34384f0b6dfSMatthew G. Knepley   Output Parameter:
34484f0b6dfSMatthew G. Knepley . labelnew - location to put new vector
34584f0b6dfSMatthew G. Knepley 
34684f0b6dfSMatthew G. Knepley   Level: intermediate
34784f0b6dfSMatthew G. Knepley 
34884f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelDestroy()
34984f0b6dfSMatthew G. Knepley @*/
350c58f1c22SToby Isaac PetscErrorCode DMLabelDuplicate(DMLabel label, DMLabel *labelnew)
351c58f1c22SToby Isaac {
352d67d17b1SMatthew G. Knepley   const char    *name;
353ad8374ffSToby Isaac   PetscInt       v;
354c58f1c22SToby Isaac   PetscErrorCode ierr;
355c58f1c22SToby Isaac 
356c58f1c22SToby Isaac   PetscFunctionBegin;
357d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
358c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
359d67d17b1SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) label, &name);CHKERRQ(ierr);
360d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PetscObjectComm((PetscObject) label), name, labelnew);CHKERRQ(ierr);
361c58f1c22SToby Isaac 
362c58f1c22SToby Isaac   (*labelnew)->numStrata    = label->numStrata;
3635aa44df4SToby Isaac   (*labelnew)->defaultValue = label->defaultValue;
364c58f1c22SToby Isaac   if (label->numStrata) {
365c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumValues);CHKERRQ(ierr);
366c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->stratumSizes);CHKERRQ(ierr);
367c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->ht);CHKERRQ(ierr);
368c58f1c22SToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->points);CHKERRQ(ierr);
369ad8374ffSToby Isaac     ierr = PetscMalloc1(label->numStrata, &(*labelnew)->validIS);CHKERRQ(ierr);
370c58f1c22SToby Isaac     /* Could eliminate unused space here */
371c58f1c22SToby Isaac     for (v = 0; v < label->numStrata; ++v) {
372e8f14785SLisandro Dalcin       ierr = PetscHSetICreate(&(*labelnew)->ht[v]);CHKERRQ(ierr);
373ad8374ffSToby Isaac       (*labelnew)->validIS[v]        = PETSC_TRUE;
374c58f1c22SToby Isaac       (*labelnew)->stratumValues[v]  = label->stratumValues[v];
375c58f1c22SToby Isaac       (*labelnew)->stratumSizes[v]   = label->stratumSizes[v];
376ad8374ffSToby Isaac       ierr = PetscObjectReference((PetscObject) (label->points[v]));CHKERRQ(ierr);
377ad8374ffSToby Isaac       (*labelnew)->points[v]         = label->points[v];
378c58f1c22SToby Isaac     }
379c58f1c22SToby Isaac   }
380c58f1c22SToby Isaac   (*labelnew)->pStart = -1;
381c58f1c22SToby Isaac   (*labelnew)->pEnd   = -1;
382c58f1c22SToby Isaac   (*labelnew)->bt     = NULL;
383c58f1c22SToby Isaac   PetscFunctionReturn(0);
384c58f1c22SToby Isaac }
385c58f1c22SToby Isaac 
386*c6a43d28SMatthew G. Knepley /*@
387*c6a43d28SMatthew G. Knepley   DMLabelComputeIndex - Create an index structure for membership determination, automatically determining the bounds
388*c6a43d28SMatthew G. Knepley 
389*c6a43d28SMatthew G. Knepley   Input Parameter:
390*c6a43d28SMatthew G. Knepley . label  - The DMLabel
391*c6a43d28SMatthew G. Knepley 
392*c6a43d28SMatthew G. Knepley   Level: intermediate
393*c6a43d28SMatthew G. Knepley 
394*c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
395*c6a43d28SMatthew G. Knepley @*/
396*c6a43d28SMatthew G. Knepley PetscErrorCode DMLabelComputeIndex(DMLabel label)
397*c6a43d28SMatthew G. Knepley {
398*c6a43d28SMatthew G. Knepley   PetscInt       pStart = PETSC_MAX_INT, pEnd = -1, v;
399*c6a43d28SMatthew G. Knepley   PetscErrorCode ierr;
400*c6a43d28SMatthew G. Knepley 
401*c6a43d28SMatthew G. Knepley   PetscFunctionBegin;
402*c6a43d28SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
403*c6a43d28SMatthew G. Knepley   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
404*c6a43d28SMatthew G. Knepley   for (v = 0; v < label->numStrata; ++v) {
405*c6a43d28SMatthew G. Knepley     const PetscInt *points;
406*c6a43d28SMatthew G. Knepley     PetscInt       i;
407*c6a43d28SMatthew G. Knepley 
408*c6a43d28SMatthew G. Knepley     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
409*c6a43d28SMatthew G. Knepley     for (i = 0; i < label->stratumSizes[v]; ++i) {
410*c6a43d28SMatthew G. Knepley       const PetscInt point = points[i];
411*c6a43d28SMatthew G. Knepley 
412*c6a43d28SMatthew G. Knepley       pStart = PetscMin(point,   pStart);
413*c6a43d28SMatthew G. Knepley       pEnd   = PetscMax(point+1, pEnd);
414*c6a43d28SMatthew G. Knepley     }
415*c6a43d28SMatthew G. Knepley     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
416*c6a43d28SMatthew G. Knepley   }
417*c6a43d28SMatthew G. Knepley   label->pStart = pStart == PETSC_MAX_INT ? -1 : pStart;
418*c6a43d28SMatthew G. Knepley   label->pEnd   = pEnd;
419*c6a43d28SMatthew G. Knepley   ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
420*c6a43d28SMatthew G. Knepley   PetscFunctionReturn(0);
421*c6a43d28SMatthew G. Knepley }
422*c6a43d28SMatthew G. Knepley 
423*c6a43d28SMatthew G. Knepley /*@
424*c6a43d28SMatthew G. Knepley   DMLabelCreateIndex - Create an index structure for membership determination
425*c6a43d28SMatthew G. Knepley 
426*c6a43d28SMatthew G. Knepley   Input Parameters:
427*c6a43d28SMatthew G. Knepley + label  - The DMLabel
428*c6a43d28SMatthew G. Knepley . pStart - The smallest point
429*c6a43d28SMatthew G. Knepley - pEnd   - The largest point + 1
430*c6a43d28SMatthew G. Knepley 
431*c6a43d28SMatthew G. Knepley   Level: intermediate
432*c6a43d28SMatthew G. Knepley 
433*c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelComputeIndex(), DMLabelDestroyIndex(), DMLabelGetValue(), DMLabelSetValue()
434*c6a43d28SMatthew G. Knepley @*/
435c58f1c22SToby Isaac PetscErrorCode DMLabelCreateIndex(DMLabel label, PetscInt pStart, PetscInt pEnd)
436c58f1c22SToby Isaac {
437c58f1c22SToby Isaac   PetscInt       v;
438c58f1c22SToby Isaac   PetscErrorCode ierr;
439c58f1c22SToby Isaac 
440c58f1c22SToby Isaac   PetscFunctionBegin;
441d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
4420c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
443c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
444c58f1c22SToby Isaac   label->pStart = pStart;
445c58f1c22SToby Isaac   label->pEnd   = pEnd;
446*c6a43d28SMatthew G. Knepley   /* This can be hooked into SetValue(),  ClearValue(), etc. for updating */
447c58f1c22SToby Isaac   ierr = PetscBTCreate(pEnd - pStart, &label->bt);CHKERRQ(ierr);
448c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
449ad8374ffSToby Isaac     const PetscInt *points;
450c58f1c22SToby Isaac     PetscInt       i;
451c58f1c22SToby Isaac 
452ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
453c58f1c22SToby Isaac     for (i = 0; i < label->stratumSizes[v]; ++i) {
454ad8374ffSToby Isaac       const PetscInt point = points[i];
455c58f1c22SToby Isaac 
456c58f1c22SToby 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);
457c58f1c22SToby Isaac       ierr = PetscBTSet(label->bt, point - pStart);CHKERRQ(ierr);
458c58f1c22SToby Isaac     }
459ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
460c58f1c22SToby Isaac   }
461c58f1c22SToby Isaac   PetscFunctionReturn(0);
462c58f1c22SToby Isaac }
463c58f1c22SToby Isaac 
464*c6a43d28SMatthew G. Knepley /*@
465*c6a43d28SMatthew G. Knepley   DMLabelDestroyIndex - Destroy the index structure
466*c6a43d28SMatthew G. Knepley 
467*c6a43d28SMatthew G. Knepley   Input Parameter:
468*c6a43d28SMatthew G. Knepley . label - the DMLabel
469*c6a43d28SMatthew G. Knepley 
470*c6a43d28SMatthew G. Knepley   Level: intermediate
471*c6a43d28SMatthew G. Knepley 
472*c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
473*c6a43d28SMatthew G. Knepley @*/
474c58f1c22SToby Isaac PetscErrorCode DMLabelDestroyIndex(DMLabel label)
475c58f1c22SToby Isaac {
476c58f1c22SToby Isaac   PetscErrorCode ierr;
477c58f1c22SToby Isaac 
478c58f1c22SToby Isaac   PetscFunctionBegin;
479d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
480c58f1c22SToby Isaac   label->pStart = -1;
481c58f1c22SToby Isaac   label->pEnd   = -1;
4820c3c4a36SLisandro Dalcin   ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
483c58f1c22SToby Isaac   PetscFunctionReturn(0);
484c58f1c22SToby Isaac }
485c58f1c22SToby Isaac 
486c58f1c22SToby Isaac /*@
487*c6a43d28SMatthew G. Knepley   DMLabelGetBounds - Return the smallest and largest point in the label
488*c6a43d28SMatthew G. Knepley 
489*c6a43d28SMatthew G. Knepley   Input Parameter:
490*c6a43d28SMatthew G. Knepley . label - the DMLabel
491*c6a43d28SMatthew G. Knepley 
492*c6a43d28SMatthew G. Knepley   Output Parameters:
493*c6a43d28SMatthew G. Knepley + pStart - The smallest point
494*c6a43d28SMatthew G. Knepley - pEnd   - The largest point + 1
495*c6a43d28SMatthew G. Knepley 
496*c6a43d28SMatthew G. Knepley   Note: This will compute an index for the label if one does not exist.
497*c6a43d28SMatthew G. Knepley 
498*c6a43d28SMatthew G. Knepley   Level: intermediate
499*c6a43d28SMatthew G. Knepley 
500*c6a43d28SMatthew G. Knepley .seealso: DMLabelHasPoint(), DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
501*c6a43d28SMatthew G. Knepley @*/
502*c6a43d28SMatthew G. Knepley PetscErrorCode DMLabelGetBounds(DMLabel label, PetscInt *pStart, PetscInt *pEnd)
503*c6a43d28SMatthew G. Knepley {
504*c6a43d28SMatthew G. Knepley   PetscErrorCode ierr;
505*c6a43d28SMatthew G. Knepley 
506*c6a43d28SMatthew G. Knepley   PetscFunctionBegin;
507*c6a43d28SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
508*c6a43d28SMatthew G. Knepley   if ((label->pStart == -1) && (label->pEnd == -1)) {ierr = DMLabelComputeIndex(label);CHKERRQ(ierr);}
509*c6a43d28SMatthew G. Knepley   if (pStart) {
510*c6a43d28SMatthew G. Knepley     PetscValidPointer(pStart, 2);
511*c6a43d28SMatthew G. Knepley     *pStart = label->pStart;
512*c6a43d28SMatthew G. Knepley   }
513*c6a43d28SMatthew G. Knepley   if (pEnd) {
514*c6a43d28SMatthew G. Knepley     PetscValidPointer(pEnd, 3);
515*c6a43d28SMatthew G. Knepley     *pEnd = label->pEnd;
516*c6a43d28SMatthew G. Knepley   }
517*c6a43d28SMatthew G. Knepley   PetscFunctionReturn(0);
518*c6a43d28SMatthew G. Knepley }
519*c6a43d28SMatthew G. Knepley 
520*c6a43d28SMatthew G. Knepley /*@
521c58f1c22SToby Isaac   DMLabelHasValue - Determine whether a label assigns the value to any point
522c58f1c22SToby Isaac 
523c58f1c22SToby Isaac   Input Parameters:
524c58f1c22SToby Isaac + label - the DMLabel
525c58f1c22SToby Isaac - value - the value
526c58f1c22SToby Isaac 
527c58f1c22SToby Isaac   Output Parameter:
528c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this value to any point
529c58f1c22SToby Isaac 
530c58f1c22SToby Isaac   Level: developer
531c58f1c22SToby Isaac 
532c58f1c22SToby Isaac .seealso: DMLabelHasPoint(), DMLabelGetValue(), DMLabelSetValue()
533c58f1c22SToby Isaac @*/
534c58f1c22SToby Isaac PetscErrorCode DMLabelHasValue(DMLabel label, PetscInt value, PetscBool *contains)
535c58f1c22SToby Isaac {
536c58f1c22SToby Isaac   PetscInt v;
5370c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
538c58f1c22SToby Isaac 
539c58f1c22SToby Isaac   PetscFunctionBegin;
540d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
541c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
5420c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
5430c3c4a36SLisandro Dalcin   *contains = v < 0 ? PETSC_FALSE : PETSC_TRUE;
544c58f1c22SToby Isaac   PetscFunctionReturn(0);
545c58f1c22SToby Isaac }
546c58f1c22SToby Isaac 
547c58f1c22SToby Isaac /*@
548c58f1c22SToby Isaac   DMLabelHasPoint - Determine whether a label assigns a value to a point
549c58f1c22SToby Isaac 
550c58f1c22SToby Isaac   Input Parameters:
551c58f1c22SToby Isaac + label - the DMLabel
552c58f1c22SToby Isaac - point - the point
553c58f1c22SToby Isaac 
554c58f1c22SToby Isaac   Output Parameter:
555c58f1c22SToby Isaac . contains - Flag indicating whether the label maps this point to a value
556c58f1c22SToby Isaac 
557c58f1c22SToby Isaac   Note: The user must call DMLabelCreateIndex() before this function.
558c58f1c22SToby Isaac 
559c58f1c22SToby Isaac   Level: developer
560c58f1c22SToby Isaac 
561c58f1c22SToby Isaac .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue()
562c58f1c22SToby Isaac @*/
563c58f1c22SToby Isaac PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains)
564c58f1c22SToby Isaac {
565c58f1c22SToby Isaac   PetscErrorCode ierr;
566c58f1c22SToby Isaac 
567c58f1c22SToby Isaac   PetscFunctionBeginHot;
568d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
569c58f1c22SToby Isaac   PetscValidPointer(contains, 3);
570c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
571c58f1c22SToby Isaac #if defined(PETSC_USE_DEBUG)
572c58f1c22SToby Isaac   if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()");
573c58f1c22SToby 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);
574c58f1c22SToby Isaac #endif
575c58f1c22SToby Isaac   *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE;
576c58f1c22SToby Isaac   PetscFunctionReturn(0);
577c58f1c22SToby Isaac }
578c58f1c22SToby Isaac 
579c58f1c22SToby Isaac /*@
580c58f1c22SToby Isaac   DMLabelStratumHasPoint - Return true if the stratum contains a point
581c58f1c22SToby Isaac 
582c58f1c22SToby Isaac   Input Parameters:
583c58f1c22SToby Isaac + label - the DMLabel
584c58f1c22SToby Isaac . value - the stratum value
585c58f1c22SToby Isaac - point - the point
586c58f1c22SToby Isaac 
587c58f1c22SToby Isaac   Output Parameter:
588c58f1c22SToby Isaac . contains - true if the stratum contains the point
589c58f1c22SToby Isaac 
590c58f1c22SToby Isaac   Level: intermediate
591c58f1c22SToby Isaac 
592c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue()
593c58f1c22SToby Isaac @*/
594c58f1c22SToby Isaac PetscErrorCode DMLabelStratumHasPoint(DMLabel label, PetscInt value, PetscInt point, PetscBool *contains)
595c58f1c22SToby Isaac {
596c58f1c22SToby Isaac   PetscInt       v;
597c58f1c22SToby Isaac   PetscErrorCode ierr;
598c58f1c22SToby Isaac 
5990c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
600d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
601c58f1c22SToby Isaac   PetscValidPointer(contains, 4);
602c58f1c22SToby Isaac   *contains = PETSC_FALSE;
6030c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
6040c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
6050c3c4a36SLisandro Dalcin 
606ad8374ffSToby Isaac   if (label->validIS[v]) {
607c58f1c22SToby Isaac     PetscInt i;
608c58f1c22SToby Isaac 
609a2d74346SToby Isaac     ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
6100c3c4a36SLisandro Dalcin     if (i >= 0) *contains = PETSC_TRUE;
611c58f1c22SToby Isaac   } else {
612c58f1c22SToby Isaac     PetscBool has;
613c58f1c22SToby Isaac 
614e8f14785SLisandro Dalcin     PetscHSetIHas(label->ht[v], point, &has);
6150c3c4a36SLisandro Dalcin     if (has) *contains = PETSC_TRUE;
616c58f1c22SToby Isaac   }
617c58f1c22SToby Isaac   PetscFunctionReturn(0);
618c58f1c22SToby Isaac }
619c58f1c22SToby Isaac 
62084f0b6dfSMatthew G. Knepley /*@
6215aa44df4SToby Isaac   DMLabelGetDefaultValue - Get the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
6225aa44df4SToby Isaac   When a label is created, it is initialized to -1.
6235aa44df4SToby Isaac 
6245aa44df4SToby Isaac   Input parameter:
6255aa44df4SToby Isaac . label - a DMLabel object
6265aa44df4SToby Isaac 
6275aa44df4SToby Isaac   Output parameter:
6285aa44df4SToby Isaac . defaultValue - the default value
6295aa44df4SToby Isaac 
6305aa44df4SToby Isaac   Level: beginner
6315aa44df4SToby Isaac 
6325aa44df4SToby Isaac .seealso: DMLabelSetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
63384f0b6dfSMatthew G. Knepley @*/
6345aa44df4SToby Isaac PetscErrorCode DMLabelGetDefaultValue(DMLabel label, PetscInt *defaultValue)
6355aa44df4SToby Isaac {
6365aa44df4SToby Isaac   PetscFunctionBegin;
637d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
6385aa44df4SToby Isaac   *defaultValue = label->defaultValue;
6395aa44df4SToby Isaac   PetscFunctionReturn(0);
6405aa44df4SToby Isaac }
6415aa44df4SToby Isaac 
64284f0b6dfSMatthew G. Knepley /*@
6435aa44df4SToby Isaac   DMLabelSetDefaultValue - Set the default value returned by DMLabelGetValue() if a point has not been explicitly given a value.
6445aa44df4SToby Isaac   When a label is created, it is initialized to -1.
6455aa44df4SToby Isaac 
6465aa44df4SToby Isaac   Input parameter:
6475aa44df4SToby Isaac . label - a DMLabel object
6485aa44df4SToby Isaac 
6495aa44df4SToby Isaac   Output parameter:
6505aa44df4SToby Isaac . defaultValue - the default value
6515aa44df4SToby Isaac 
6525aa44df4SToby Isaac   Level: beginner
6535aa44df4SToby Isaac 
6545aa44df4SToby Isaac .seealso: DMLabelGetDefaultValue(), DMLabelGetValue(), DMLabelSetValue()
65584f0b6dfSMatthew G. Knepley @*/
6565aa44df4SToby Isaac PetscErrorCode DMLabelSetDefaultValue(DMLabel label, PetscInt defaultValue)
6575aa44df4SToby Isaac {
6585aa44df4SToby Isaac   PetscFunctionBegin;
659d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
6605aa44df4SToby Isaac   label->defaultValue = defaultValue;
6615aa44df4SToby Isaac   PetscFunctionReturn(0);
6625aa44df4SToby Isaac }
6635aa44df4SToby Isaac 
664c58f1c22SToby Isaac /*@
6655aa44df4SToby 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())
666c58f1c22SToby Isaac 
667c58f1c22SToby Isaac   Input Parameters:
668c58f1c22SToby Isaac + label - the DMLabel
669c58f1c22SToby Isaac - point - the point
670c58f1c22SToby Isaac 
671c58f1c22SToby Isaac   Output Parameter:
6728e68afcfSMatthew G. Knepley . value - The point value, or the default value (-1 by default)
673c58f1c22SToby Isaac 
674c58f1c22SToby Isaac   Level: intermediate
675c58f1c22SToby Isaac 
6765aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelSetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
677c58f1c22SToby Isaac @*/
678c58f1c22SToby Isaac PetscErrorCode DMLabelGetValue(DMLabel label, PetscInt point, PetscInt *value)
679c58f1c22SToby Isaac {
680c58f1c22SToby Isaac   PetscInt       v;
681c58f1c22SToby Isaac   PetscErrorCode ierr;
682c58f1c22SToby Isaac 
6830c3c4a36SLisandro Dalcin   PetscFunctionBeginHot;
684d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
685c58f1c22SToby Isaac   PetscValidPointer(value, 3);
6865aa44df4SToby Isaac   *value = label->defaultValue;
687c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
688ad8374ffSToby Isaac     if (label->validIS[v]) {
689c58f1c22SToby Isaac       PetscInt i;
690c58f1c22SToby Isaac 
691a2d74346SToby Isaac       ierr = ISLocate(label->points[v], point, &i);CHKERRQ(ierr);
692c58f1c22SToby Isaac       if (i >= 0) {
693c58f1c22SToby Isaac         *value = label->stratumValues[v];
694c58f1c22SToby Isaac         break;
695c58f1c22SToby Isaac       }
696c58f1c22SToby Isaac     } else {
697c58f1c22SToby Isaac       PetscBool has;
698c58f1c22SToby Isaac 
699e8f14785SLisandro Dalcin       PetscHSetIHas(label->ht[v], point, &has);
700c58f1c22SToby Isaac       if (has) {
701c58f1c22SToby Isaac         *value = label->stratumValues[v];
702c58f1c22SToby Isaac         break;
703c58f1c22SToby Isaac       }
704c58f1c22SToby Isaac     }
705c58f1c22SToby Isaac   }
706c58f1c22SToby Isaac   PetscFunctionReturn(0);
707c58f1c22SToby Isaac }
708c58f1c22SToby Isaac 
709c58f1c22SToby Isaac /*@
710367003a6SStefano Zampini   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 something different), then this function will do nothing.
711c58f1c22SToby Isaac 
712c58f1c22SToby Isaac   Input Parameters:
713c58f1c22SToby Isaac + label - the DMLabel
714c58f1c22SToby Isaac . point - the point
715c58f1c22SToby Isaac - value - The point value
716c58f1c22SToby Isaac 
717c58f1c22SToby Isaac   Level: intermediate
718c58f1c22SToby Isaac 
7195aa44df4SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelClearValue(), DMLabelGetDefaultValue(), DMLabelSetDefaultValue()
720c58f1c22SToby Isaac @*/
721c58f1c22SToby Isaac PetscErrorCode DMLabelSetValue(DMLabel label, PetscInt point, PetscInt value)
722c58f1c22SToby Isaac {
723c58f1c22SToby Isaac   PetscInt       v;
724c58f1c22SToby Isaac   PetscErrorCode ierr;
725c58f1c22SToby Isaac 
726c58f1c22SToby Isaac   PetscFunctionBegin;
727d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
7280c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
7295aa44df4SToby Isaac   if (value == label->defaultValue) PetscFunctionReturn(0);
7300c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
7310c3c4a36SLisandro Dalcin   if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);}
732c58f1c22SToby Isaac   /* Set key */
7330c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
734e8f14785SLisandro Dalcin   ierr = PetscHSetIAdd(label->ht[v], point);CHKERRQ(ierr);
735c58f1c22SToby Isaac   PetscFunctionReturn(0);
736c58f1c22SToby Isaac }
737c58f1c22SToby Isaac 
738c58f1c22SToby Isaac /*@
739c58f1c22SToby Isaac   DMLabelClearValue - Clear the value a label assigns to a point
740c58f1c22SToby Isaac 
741c58f1c22SToby Isaac   Input Parameters:
742c58f1c22SToby Isaac + label - the DMLabel
743c58f1c22SToby Isaac . point - the point
744c58f1c22SToby Isaac - value - The point value
745c58f1c22SToby Isaac 
746c58f1c22SToby Isaac   Level: intermediate
747c58f1c22SToby Isaac 
748c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue()
749c58f1c22SToby Isaac @*/
750c58f1c22SToby Isaac PetscErrorCode DMLabelClearValue(DMLabel label, PetscInt point, PetscInt value)
751c58f1c22SToby Isaac {
752ad8374ffSToby Isaac   PetscInt       v;
753c58f1c22SToby Isaac   PetscErrorCode ierr;
754c58f1c22SToby Isaac 
755c58f1c22SToby Isaac   PetscFunctionBegin;
756d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
757c58f1c22SToby Isaac   /* Find label value */
7580c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
7590c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
7600c3c4a36SLisandro Dalcin 
761eeed21e7SToby Isaac   if (label->bt) {
762eeed21e7SToby 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);
763eeed21e7SToby Isaac     ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
764eeed21e7SToby Isaac   }
7650c3c4a36SLisandro Dalcin 
7660c3c4a36SLisandro Dalcin   /* Delete key */
7670c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
768e8f14785SLisandro Dalcin   ierr = PetscHSetIDel(label->ht[v], point);CHKERRQ(ierr);
769c58f1c22SToby Isaac   PetscFunctionReturn(0);
770c58f1c22SToby Isaac }
771c58f1c22SToby Isaac 
772c58f1c22SToby Isaac /*@
773c58f1c22SToby Isaac   DMLabelInsertIS - Set all points in the IS to a value
774c58f1c22SToby Isaac 
775c58f1c22SToby Isaac   Input Parameters:
776c58f1c22SToby Isaac + label - the DMLabel
777c58f1c22SToby Isaac . is    - the point IS
778c58f1c22SToby Isaac - value - The point value
779c58f1c22SToby Isaac 
780c58f1c22SToby Isaac   Level: intermediate
781c58f1c22SToby Isaac 
782c58f1c22SToby Isaac .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
783c58f1c22SToby Isaac @*/
784c58f1c22SToby Isaac PetscErrorCode DMLabelInsertIS(DMLabel label, IS is, PetscInt value)
785c58f1c22SToby Isaac {
7860c3c4a36SLisandro Dalcin   PetscInt        v, n, p;
787c58f1c22SToby Isaac   const PetscInt *points;
788c58f1c22SToby Isaac   PetscErrorCode  ierr;
789c58f1c22SToby Isaac 
790c58f1c22SToby Isaac   PetscFunctionBegin;
791d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
792c58f1c22SToby Isaac   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
7930c3c4a36SLisandro Dalcin   /* Find label value, add new entry if needed */
7940c3c4a36SLisandro Dalcin   if (value == label->defaultValue) PetscFunctionReturn(0);
7950c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
7960c3c4a36SLisandro Dalcin   if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);}
7970c3c4a36SLisandro Dalcin   /* Set keys */
7980c3c4a36SLisandro Dalcin   ierr = DMLabelMakeInvalid_Private(label, v);CHKERRQ(ierr);
799c58f1c22SToby Isaac   ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr);
800c58f1c22SToby Isaac   ierr = ISGetIndices(is, &points);CHKERRQ(ierr);
801e8f14785SLisandro Dalcin   for (p = 0; p < n; ++p) {ierr = PetscHSetIAdd(label->ht[v], points[p]);CHKERRQ(ierr);}
802c58f1c22SToby Isaac   ierr = ISRestoreIndices(is, &points);CHKERRQ(ierr);
803c58f1c22SToby Isaac   PetscFunctionReturn(0);
804c58f1c22SToby Isaac }
805c58f1c22SToby Isaac 
80684f0b6dfSMatthew G. Knepley /*@
80784f0b6dfSMatthew G. Knepley   DMLabelGetNumValues - Get the number of values that the DMLabel takes
80884f0b6dfSMatthew G. Knepley 
80984f0b6dfSMatthew G. Knepley   Input Parameter:
81084f0b6dfSMatthew G. Knepley . label - the DMLabel
81184f0b6dfSMatthew G. Knepley 
81284f0b6dfSMatthew G. Knepley   Output Paramater:
81384f0b6dfSMatthew G. Knepley . numValues - the number of values
81484f0b6dfSMatthew G. Knepley 
81584f0b6dfSMatthew G. Knepley   Level: intermediate
81684f0b6dfSMatthew G. Knepley 
81784f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
81884f0b6dfSMatthew G. Knepley @*/
819c58f1c22SToby Isaac PetscErrorCode DMLabelGetNumValues(DMLabel label, PetscInt *numValues)
820c58f1c22SToby Isaac {
821c58f1c22SToby Isaac   PetscFunctionBegin;
822d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
823c58f1c22SToby Isaac   PetscValidPointer(numValues, 2);
824c58f1c22SToby Isaac   *numValues = label->numStrata;
825c58f1c22SToby Isaac   PetscFunctionReturn(0);
826c58f1c22SToby Isaac }
827c58f1c22SToby Isaac 
82884f0b6dfSMatthew G. Knepley /*@
82984f0b6dfSMatthew G. Knepley   DMLabelGetValueIS - Get an IS of all values that the DMlabel takes
83084f0b6dfSMatthew G. Knepley 
83184f0b6dfSMatthew G. Knepley   Input Parameter:
83284f0b6dfSMatthew G. Knepley . label - the DMLabel
83384f0b6dfSMatthew G. Knepley 
83484f0b6dfSMatthew G. Knepley   Output Paramater:
83584f0b6dfSMatthew G. Knepley . is    - the value IS
83684f0b6dfSMatthew G. Knepley 
83784f0b6dfSMatthew G. Knepley   Level: intermediate
83884f0b6dfSMatthew G. Knepley 
83984f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
84084f0b6dfSMatthew G. Knepley @*/
841c58f1c22SToby Isaac PetscErrorCode DMLabelGetValueIS(DMLabel label, IS *values)
842c58f1c22SToby Isaac {
843c58f1c22SToby Isaac   PetscErrorCode ierr;
844c58f1c22SToby Isaac 
845c58f1c22SToby Isaac   PetscFunctionBegin;
846d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
847c58f1c22SToby Isaac   PetscValidPointer(values, 2);
848c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, label->numStrata, label->stratumValues, PETSC_USE_POINTER, values);CHKERRQ(ierr);
849c58f1c22SToby Isaac   PetscFunctionReturn(0);
850c58f1c22SToby Isaac }
851c58f1c22SToby Isaac 
85284f0b6dfSMatthew G. Knepley /*@
85384f0b6dfSMatthew G. Knepley   DMLabelHasStratum - Determine whether points exist with the given value
85484f0b6dfSMatthew G. Knepley 
85584f0b6dfSMatthew G. Knepley   Input Parameters:
85684f0b6dfSMatthew G. Knepley + label - the DMLabel
85784f0b6dfSMatthew G. Knepley - value - the stratum value
85884f0b6dfSMatthew G. Knepley 
85984f0b6dfSMatthew G. Knepley   Output Paramater:
86084f0b6dfSMatthew G. Knepley . exists - Flag saying whether points exist
86184f0b6dfSMatthew G. Knepley 
86284f0b6dfSMatthew G. Knepley   Level: intermediate
86384f0b6dfSMatthew G. Knepley 
86484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
86584f0b6dfSMatthew G. Knepley @*/
866fada774cSMatthew G. Knepley PetscErrorCode DMLabelHasStratum(DMLabel label, PetscInt value, PetscBool *exists)
867fada774cSMatthew G. Knepley {
868fada774cSMatthew G. Knepley   PetscInt       v;
8690c3c4a36SLisandro Dalcin   PetscErrorCode ierr;
870fada774cSMatthew G. Knepley 
871fada774cSMatthew G. Knepley   PetscFunctionBegin;
872d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
873fada774cSMatthew G. Knepley   PetscValidPointer(exists, 3);
8740c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
8750c3c4a36SLisandro Dalcin   *exists = v < 0 ? PETSC_FALSE : PETSC_TRUE;
876fada774cSMatthew G. Knepley   PetscFunctionReturn(0);
877fada774cSMatthew G. Knepley }
878fada774cSMatthew G. Knepley 
87984f0b6dfSMatthew G. Knepley /*@
88084f0b6dfSMatthew G. Knepley   DMLabelGetStratumSize - Get the size of a stratum
88184f0b6dfSMatthew G. Knepley 
88284f0b6dfSMatthew G. Knepley   Input Parameters:
88384f0b6dfSMatthew G. Knepley + label - the DMLabel
88484f0b6dfSMatthew G. Knepley - value - the stratum value
88584f0b6dfSMatthew G. Knepley 
88684f0b6dfSMatthew G. Knepley   Output Paramater:
88784f0b6dfSMatthew G. Knepley . size - The number of points in the stratum
88884f0b6dfSMatthew G. Knepley 
88984f0b6dfSMatthew G. Knepley   Level: intermediate
89084f0b6dfSMatthew G. Knepley 
89184f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
89284f0b6dfSMatthew G. Knepley @*/
893c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumSize(DMLabel label, PetscInt value, PetscInt *size)
894c58f1c22SToby Isaac {
895c58f1c22SToby Isaac   PetscInt       v;
896c58f1c22SToby Isaac   PetscErrorCode ierr;
897c58f1c22SToby Isaac 
898c58f1c22SToby Isaac   PetscFunctionBegin;
899d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
900c58f1c22SToby Isaac   PetscValidPointer(size, 3);
901c58f1c22SToby Isaac   *size = 0;
9020c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
9030c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
904c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
905c58f1c22SToby Isaac   *size = label->stratumSizes[v];
906c58f1c22SToby Isaac   PetscFunctionReturn(0);
907c58f1c22SToby Isaac }
908c58f1c22SToby Isaac 
90984f0b6dfSMatthew G. Knepley /*@
91084f0b6dfSMatthew G. Knepley   DMLabelGetStratumBounds - Get the largest and smallest point of a stratum
91184f0b6dfSMatthew G. Knepley 
91284f0b6dfSMatthew G. Knepley   Input Parameters:
91384f0b6dfSMatthew G. Knepley + label - the DMLabel
91484f0b6dfSMatthew G. Knepley - value - the stratum value
91584f0b6dfSMatthew G. Knepley 
91684f0b6dfSMatthew G. Knepley   Output Paramaters:
91784f0b6dfSMatthew G. Knepley + start - the smallest point in the stratum
91884f0b6dfSMatthew G. Knepley - end - the largest point in the stratum
91984f0b6dfSMatthew G. Knepley 
92084f0b6dfSMatthew G. Knepley   Level: intermediate
92184f0b6dfSMatthew G. Knepley 
92284f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
92384f0b6dfSMatthew G. Knepley @*/
924c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumBounds(DMLabel label, PetscInt value, PetscInt *start, PetscInt *end)
925c58f1c22SToby Isaac {
9260c3c4a36SLisandro Dalcin   PetscInt       v, min, max;
927c58f1c22SToby Isaac   PetscErrorCode ierr;
928c58f1c22SToby Isaac 
929c58f1c22SToby Isaac   PetscFunctionBegin;
930d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
931c58f1c22SToby Isaac   if (start) {PetscValidPointer(start, 3); *start = 0;}
932c58f1c22SToby Isaac   if (end)   {PetscValidPointer(end,   4); *end   = 0;}
9330c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
9340c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
935c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
9360c3c4a36SLisandro Dalcin   if (label->stratumSizes[v] <= 0) PetscFunctionReturn(0);
937d6cb179aSToby Isaac   ierr = ISGetMinMax(label->points[v], &min, &max);CHKERRQ(ierr);
938d6cb179aSToby Isaac   if (start) *start = min;
939d6cb179aSToby Isaac   if (end)   *end   = max+1;
940c58f1c22SToby Isaac   PetscFunctionReturn(0);
941c58f1c22SToby Isaac }
942c58f1c22SToby Isaac 
94384f0b6dfSMatthew G. Knepley /*@
94484f0b6dfSMatthew G. Knepley   DMLabelGetStratumIS - Get an IS with the stratum points
94584f0b6dfSMatthew G. Knepley 
94684f0b6dfSMatthew G. Knepley   Input Parameters:
94784f0b6dfSMatthew G. Knepley + label - the DMLabel
94884f0b6dfSMatthew G. Knepley - value - the stratum value
94984f0b6dfSMatthew G. Knepley 
95084f0b6dfSMatthew G. Knepley   Output Paramater:
95184f0b6dfSMatthew G. Knepley . points - The stratum points
95284f0b6dfSMatthew G. Knepley 
95384f0b6dfSMatthew G. Knepley   Level: intermediate
95484f0b6dfSMatthew G. Knepley 
95584f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
95684f0b6dfSMatthew G. Knepley @*/
957c58f1c22SToby Isaac PetscErrorCode DMLabelGetStratumIS(DMLabel label, PetscInt value, IS *points)
958c58f1c22SToby Isaac {
959c58f1c22SToby Isaac   PetscInt       v;
960c58f1c22SToby Isaac   PetscErrorCode ierr;
961c58f1c22SToby Isaac 
962c58f1c22SToby Isaac   PetscFunctionBegin;
963d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
964c58f1c22SToby Isaac   PetscValidPointer(points, 3);
965c58f1c22SToby Isaac   *points = NULL;
9660c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
9670c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
968c58f1c22SToby Isaac   ierr = DMLabelMakeValid_Private(label, v);CHKERRQ(ierr);
969ad8374ffSToby Isaac   ierr = PetscObjectReference((PetscObject) label->points[v]);CHKERRQ(ierr);
970ad8374ffSToby Isaac   *points = label->points[v];
971c58f1c22SToby Isaac   PetscFunctionReturn(0);
972c58f1c22SToby Isaac }
973c58f1c22SToby Isaac 
97484f0b6dfSMatthew G. Knepley /*@
97584f0b6dfSMatthew G. Knepley   DMLabelSetStratumIS - Set the stratum points using an IS
97684f0b6dfSMatthew G. Knepley 
97784f0b6dfSMatthew G. Knepley   Input Parameters:
97884f0b6dfSMatthew G. Knepley + label - the DMLabel
97984f0b6dfSMatthew G. Knepley . value - the stratum value
98084f0b6dfSMatthew G. Knepley - points - The stratum points
98184f0b6dfSMatthew G. Knepley 
98284f0b6dfSMatthew G. Knepley   Level: intermediate
98384f0b6dfSMatthew G. Knepley 
98484f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
98584f0b6dfSMatthew G. Knepley @*/
9864de306b1SToby Isaac PetscErrorCode DMLabelSetStratumIS(DMLabel label, PetscInt value, IS is)
9874de306b1SToby Isaac {
9880c3c4a36SLisandro Dalcin   PetscInt       v;
9894de306b1SToby Isaac   PetscErrorCode ierr;
9904de306b1SToby Isaac 
9914de306b1SToby Isaac   PetscFunctionBegin;
992d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
993d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
9940c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
9950c3c4a36SLisandro Dalcin   if (v < 0) {ierr = DMLabelNewStratum(label, value, &v);CHKERRQ(ierr);}
9964de306b1SToby Isaac   if (is == label->points[v]) PetscFunctionReturn(0);
9974de306b1SToby Isaac   ierr = DMLabelClearStratum(label, value);CHKERRQ(ierr);
9984de306b1SToby Isaac   ierr = ISGetLocalSize(is, &(label->stratumSizes[v]));CHKERRQ(ierr);
9994de306b1SToby Isaac   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
10004de306b1SToby Isaac   ierr = ISDestroy(&(label->points[v]));CHKERRQ(ierr);
10010c3c4a36SLisandro Dalcin   label->points[v] = is;
10020c3c4a36SLisandro Dalcin   label->validIS[v] = PETSC_TRUE;
10034de306b1SToby Isaac   if (label->bt) {
10044de306b1SToby Isaac     const PetscInt *points;
10054de306b1SToby Isaac     PetscInt p;
10064de306b1SToby Isaac 
10074de306b1SToby Isaac     ierr = ISGetIndices(is,&points);CHKERRQ(ierr);
10084de306b1SToby Isaac     for (p = 0; p < label->stratumSizes[v]; ++p) {
10094de306b1SToby Isaac       const PetscInt point = points[p];
10104de306b1SToby Isaac 
10114de306b1SToby 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);
10124de306b1SToby Isaac       ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr);
10134de306b1SToby Isaac     }
10144de306b1SToby Isaac   }
10154de306b1SToby Isaac   PetscFunctionReturn(0);
10164de306b1SToby Isaac }
10174de306b1SToby Isaac 
101884f0b6dfSMatthew G. Knepley /*@
101984f0b6dfSMatthew G. Knepley   DMLabelClearStratum - Remove a stratum
10204de306b1SToby Isaac 
102184f0b6dfSMatthew G. Knepley   Input Parameters:
102284f0b6dfSMatthew G. Knepley + label - the DMLabel
102384f0b6dfSMatthew G. Knepley - value - the stratum value
102484f0b6dfSMatthew G. Knepley 
102584f0b6dfSMatthew G. Knepley   Level: intermediate
102684f0b6dfSMatthew G. Knepley 
102784f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
102884f0b6dfSMatthew G. Knepley @*/
1029c58f1c22SToby Isaac PetscErrorCode DMLabelClearStratum(DMLabel label, PetscInt value)
1030c58f1c22SToby Isaac {
1031c58f1c22SToby Isaac   PetscInt       v;
1032c58f1c22SToby Isaac   PetscErrorCode ierr;
1033c58f1c22SToby Isaac 
1034c58f1c22SToby Isaac   PetscFunctionBegin;
1035d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
10360c3c4a36SLisandro Dalcin   ierr = DMLabelLookupStratum(label, value, &v);CHKERRQ(ierr);
10370c3c4a36SLisandro Dalcin   if (v < 0) PetscFunctionReturn(0);
10384de306b1SToby Isaac   if (label->validIS[v]) {
10394de306b1SToby Isaac     if (label->bt) {
1040c58f1c22SToby Isaac       PetscInt       i;
1041ad8374ffSToby Isaac       const PetscInt *points;
1042c58f1c22SToby Isaac 
1043ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
1044c58f1c22SToby Isaac       for (i = 0; i < label->stratumSizes[v]; ++i) {
1045ad8374ffSToby Isaac         const PetscInt point = points[i];
1046c58f1c22SToby Isaac 
1047c58f1c22SToby 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);
1048c58f1c22SToby Isaac         ierr = PetscBTClear(label->bt, point - label->pStart);CHKERRQ(ierr);
1049c58f1c22SToby Isaac       }
1050ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[v], &points);CHKERRQ(ierr);
1051c58f1c22SToby Isaac     }
1052c58f1c22SToby Isaac     label->stratumSizes[v] = 0;
10530c3c4a36SLisandro Dalcin     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
10540c3c4a36SLisandro Dalcin     ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr);
10550c3c4a36SLisandro Dalcin     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
1056c58f1c22SToby Isaac   } else {
1057e8f14785SLisandro Dalcin     PetscHSetIClear(label->ht[v]);
1058c58f1c22SToby Isaac   }
1059c58f1c22SToby Isaac   PetscFunctionReturn(0);
1060c58f1c22SToby Isaac }
1061c58f1c22SToby Isaac 
106284f0b6dfSMatthew G. Knepley /*@
106384f0b6dfSMatthew G. Knepley   DMLabelFilter - Remove all points outside of [start, end)
106484f0b6dfSMatthew G. Knepley 
106584f0b6dfSMatthew G. Knepley   Input Parameters:
106684f0b6dfSMatthew G. Knepley + label - the DMLabel
106784f0b6dfSMatthew G. Knepley . start - the first point
106884f0b6dfSMatthew G. Knepley - end - the last point
106984f0b6dfSMatthew G. Knepley 
107084f0b6dfSMatthew G. Knepley   Level: intermediate
107184f0b6dfSMatthew G. Knepley 
107284f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
107384f0b6dfSMatthew G. Knepley @*/
1074c58f1c22SToby Isaac PetscErrorCode DMLabelFilter(DMLabel label, PetscInt start, PetscInt end)
1075c58f1c22SToby Isaac {
1076c58f1c22SToby Isaac   PetscInt       v;
1077c58f1c22SToby Isaac   PetscErrorCode ierr;
1078c58f1c22SToby Isaac 
1079c58f1c22SToby Isaac   PetscFunctionBegin;
1080d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
10810c3c4a36SLisandro Dalcin   ierr = DMLabelDestroyIndex(label);CHKERRQ(ierr);
1082c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1083c58f1c22SToby Isaac   for (v = 0; v < label->numStrata; ++v) {
1084c58f1c22SToby Isaac     PetscInt off, q;
1085ad8374ffSToby Isaac     const PetscInt *points;
1086033405d5SLisandro Dalcin     PetscInt numPointsNew = 0, *pointsNew = NULL;
1087c58f1c22SToby Isaac 
1088ad8374ffSToby Isaac     ierr = ISGetIndices(label->points[v], &points);CHKERRQ(ierr);
1089033405d5SLisandro Dalcin     for (q = 0; q < label->stratumSizes[v]; ++q)
1090033405d5SLisandro Dalcin       if (points[q] >= start && points[q] < end)
1091033405d5SLisandro Dalcin         numPointsNew++;
1092033405d5SLisandro Dalcin     ierr = PetscMalloc1(numPointsNew, &pointsNew);CHKERRQ(ierr);
1093c58f1c22SToby Isaac     for (off = 0, q = 0; q < label->stratumSizes[v]; ++q) {
1094033405d5SLisandro Dalcin       if (points[q] >= start && points[q] < end)
1095033405d5SLisandro Dalcin         pointsNew[off++] = points[q];
1096ad8374ffSToby Isaac     }
1097ad8374ffSToby Isaac     ierr = ISRestoreIndices(label->points[v],&points);CHKERRQ(ierr);
1098033405d5SLisandro Dalcin 
1099033405d5SLisandro Dalcin     label->stratumSizes[v] = numPointsNew;
1100033405d5SLisandro Dalcin     ierr = ISDestroy(&label->points[v]);CHKERRQ(ierr);
1101033405d5SLisandro Dalcin     ierr = ISCreateGeneral(PETSC_COMM_SELF,numPointsNew, pointsNew, PETSC_OWN_POINTER, &label->points[v]);CHKERRQ(ierr);
1102033405d5SLisandro Dalcin     ierr = PetscObjectSetName((PetscObject) label->points[v], "indices");CHKERRQ(ierr);
1103c58f1c22SToby Isaac   }
1104c58f1c22SToby Isaac   ierr = DMLabelCreateIndex(label, start, end);CHKERRQ(ierr);
1105c58f1c22SToby Isaac   PetscFunctionReturn(0);
1106c58f1c22SToby Isaac }
1107c58f1c22SToby Isaac 
110884f0b6dfSMatthew G. Knepley /*@
110984f0b6dfSMatthew G. Knepley   DMLabelPermute - Create a new label with permuted points
111084f0b6dfSMatthew G. Knepley 
111184f0b6dfSMatthew G. Knepley   Input Parameters:
111284f0b6dfSMatthew G. Knepley + label - the DMLabel
111384f0b6dfSMatthew G. Knepley - permutation - the point permutation
111484f0b6dfSMatthew G. Knepley 
111584f0b6dfSMatthew G. Knepley   Output Parameter:
111684f0b6dfSMatthew G. Knepley . labelnew - the new label containing the permuted points
111784f0b6dfSMatthew G. Knepley 
111884f0b6dfSMatthew G. Knepley   Level: intermediate
111984f0b6dfSMatthew G. Knepley 
112084f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
112184f0b6dfSMatthew G. Knepley @*/
1122c58f1c22SToby Isaac PetscErrorCode DMLabelPermute(DMLabel label, IS permutation, DMLabel *labelNew)
1123c58f1c22SToby Isaac {
1124c58f1c22SToby Isaac   const PetscInt *perm;
1125c58f1c22SToby Isaac   PetscInt        numValues, numPoints, v, q;
1126c58f1c22SToby Isaac   PetscErrorCode  ierr;
1127c58f1c22SToby Isaac 
1128c58f1c22SToby Isaac   PetscFunctionBegin;
1129d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1130d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(permutation, IS_CLASSID, 2);
1131c58f1c22SToby Isaac   ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1132c58f1c22SToby Isaac   ierr = DMLabelDuplicate(label, labelNew);CHKERRQ(ierr);
1133c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(*labelNew, &numValues);CHKERRQ(ierr);
1134c58f1c22SToby Isaac   ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr);
1135c58f1c22SToby Isaac   ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr);
1136c58f1c22SToby Isaac   for (v = 0; v < numValues; ++v) {
1137c58f1c22SToby Isaac     const PetscInt size   = (*labelNew)->stratumSizes[v];
1138ad8374ffSToby Isaac     const PetscInt *points;
1139ad8374ffSToby Isaac     PetscInt *pointsNew;
1140c58f1c22SToby Isaac 
1141ad8374ffSToby Isaac     ierr = ISGetIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1142ad8374ffSToby Isaac     ierr = PetscMalloc1(size,&pointsNew);CHKERRQ(ierr);
1143c58f1c22SToby Isaac     for (q = 0; q < size; ++q) {
1144ad8374ffSToby Isaac       const PetscInt point = points[q];
1145c58f1c22SToby Isaac 
1146c58f1c22SToby 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);
1147ad8374ffSToby Isaac       pointsNew[q] = perm[point];
1148c58f1c22SToby Isaac     }
1149ad8374ffSToby Isaac     ierr = ISRestoreIndices((*labelNew)->points[v],&points);CHKERRQ(ierr);
1150ad8374ffSToby Isaac     ierr = PetscSortInt(size, pointsNew);CHKERRQ(ierr);
1151ad8374ffSToby Isaac     ierr = ISDestroy(&((*labelNew)->points[v]));CHKERRQ(ierr);
1152fa8e8ae5SToby Isaac     if (size > 0 && pointsNew[size - 1] == pointsNew[0] + size - 1) {
1153fa8e8ae5SToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,size,pointsNew[0],1,&((*labelNew)->points[v]));CHKERRQ(ierr);
1154fa8e8ae5SToby Isaac       ierr = PetscFree(pointsNew);CHKERRQ(ierr);
1155fa8e8ae5SToby Isaac     } else {
1156ad8374ffSToby Isaac       ierr = ISCreateGeneral(PETSC_COMM_SELF,size,pointsNew,PETSC_OWN_POINTER,&((*labelNew)->points[v]));CHKERRQ(ierr);
1157fa8e8ae5SToby Isaac     }
1158ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject) ((*labelNew)->points[v]), "indices");CHKERRQ(ierr);
1159c58f1c22SToby Isaac   }
1160c58f1c22SToby Isaac   ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr);
1161c58f1c22SToby Isaac   if (label->bt) {
1162c58f1c22SToby Isaac     ierr = PetscBTDestroy(&label->bt);CHKERRQ(ierr);
1163c58f1c22SToby Isaac     ierr = DMLabelCreateIndex(label, label->pStart, label->pEnd);CHKERRQ(ierr);
1164c58f1c22SToby Isaac   }
1165c58f1c22SToby Isaac   PetscFunctionReturn(0);
1166c58f1c22SToby Isaac }
1167c58f1c22SToby Isaac 
116826c55118SMichael Lange PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata)
116926c55118SMichael Lange {
117026c55118SMichael Lange   MPI_Comm       comm;
117126c55118SMichael Lange   PetscInt       s, l, nroots, nleaves, dof, offset, size;
117226c55118SMichael Lange   PetscInt      *remoteOffsets, *rootStrata, *rootIdx;
117326c55118SMichael Lange   PetscSection   rootSection;
117426c55118SMichael Lange   PetscSF        labelSF;
117526c55118SMichael Lange   PetscErrorCode ierr;
117626c55118SMichael Lange 
117726c55118SMichael Lange   PetscFunctionBegin;
117826c55118SMichael Lange   if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);}
117926c55118SMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
118026c55118SMichael Lange   /* Build a section of stratum values per point, generate the according SF
118126c55118SMichael Lange      and distribute point-wise stratum values to leaves. */
118226c55118SMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
118326c55118SMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
118426c55118SMichael Lange   ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr);
118526c55118SMichael Lange   if (label) {
118626c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1187ad8374ffSToby Isaac       const PetscInt *points;
1188ad8374ffSToby Isaac 
1189ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
119026c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1191ad8374ffSToby Isaac         ierr = PetscSectionGetDof(rootSection, points[l], &dof);CHKERRQ(ierr);
1192ad8374ffSToby Isaac         ierr = PetscSectionSetDof(rootSection, points[l], dof+1);CHKERRQ(ierr);
119326c55118SMichael Lange       }
1194ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
119526c55118SMichael Lange     }
119626c55118SMichael Lange   }
119726c55118SMichael Lange   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
119826c55118SMichael Lange   /* Create a point-wise array of stratum values */
119926c55118SMichael Lange   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
120026c55118SMichael Lange   ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr);
120126c55118SMichael Lange   ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr);
120226c55118SMichael Lange   if (label) {
120326c55118SMichael Lange     for (s = 0; s < label->numStrata; ++s) {
1204ad8374ffSToby Isaac       const PetscInt *points;
1205ad8374ffSToby Isaac 
1206ad8374ffSToby Isaac       ierr = ISGetIndices(label->points[s], &points);CHKERRQ(ierr);
120726c55118SMichael Lange       for (l = 0; l < label->stratumSizes[s]; l++) {
1208ad8374ffSToby Isaac         const PetscInt p = points[l];
120926c55118SMichael Lange         ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr);
121026c55118SMichael Lange         rootStrata[offset+rootIdx[p]++] = label->stratumValues[s];
121126c55118SMichael Lange       }
1212ad8374ffSToby Isaac       ierr = ISRestoreIndices(label->points[s], &points);CHKERRQ(ierr);
121326c55118SMichael Lange     }
121426c55118SMichael Lange   }
121526c55118SMichael Lange   /* Build SF that maps label points to remote processes */
121626c55118SMichael Lange   ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr);
121726c55118SMichael Lange   ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr);
121826c55118SMichael Lange   ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr);
121926c55118SMichael Lange   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
122026c55118SMichael Lange   /* Send the strata for each point over the derived SF */
122126c55118SMichael Lange   ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr);
122226c55118SMichael Lange   ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr);
122326c55118SMichael Lange   ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
122426c55118SMichael Lange   ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr);
122526c55118SMichael Lange   /* Clean up */
122626c55118SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
122726c55118SMichael Lange   ierr = PetscFree(rootIdx);CHKERRQ(ierr);
122826c55118SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
122926c55118SMichael Lange   ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr);
123026c55118SMichael Lange   PetscFunctionReturn(0);
123126c55118SMichael Lange }
123226c55118SMichael Lange 
123384f0b6dfSMatthew G. Knepley /*@
123484f0b6dfSMatthew G. Knepley   DMLabelDistribute - Create a new label pushed forward over the PetscSF
123584f0b6dfSMatthew G. Knepley 
123684f0b6dfSMatthew G. Knepley   Input Parameters:
123784f0b6dfSMatthew G. Knepley + label - the DMLabel
123884f0b6dfSMatthew G. Knepley - sf    - the map from old to new distribution
123984f0b6dfSMatthew G. Knepley 
124084f0b6dfSMatthew G. Knepley   Output Parameter:
124184f0b6dfSMatthew G. Knepley . labelnew - the new redistributed label
124284f0b6dfSMatthew G. Knepley 
124384f0b6dfSMatthew G. Knepley   Level: intermediate
124484f0b6dfSMatthew G. Knepley 
124584f0b6dfSMatthew G. Knepley .seealso: DMLabelCreate(), DMLabelGetValue(), DMLabelSetValue(), DMLabelClearValue()
124684f0b6dfSMatthew G. Knepley @*/
1247c58f1c22SToby Isaac PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew)
1248c58f1c22SToby Isaac {
1249c58f1c22SToby Isaac   MPI_Comm       comm;
125026c55118SMichael Lange   PetscSection   leafSection;
125126c55118SMichael Lange   PetscInt       p, pStart, pEnd, s, size, dof, offset, stratum;
125226c55118SMichael Lange   PetscInt      *leafStrata, *strataIdx;
1253ad8374ffSToby Isaac   PetscInt     **points;
1254d67d17b1SMatthew G. Knepley   const char    *lname = NULL;
1255c58f1c22SToby Isaac   char          *name;
1256c58f1c22SToby Isaac   PetscInt       nameSize;
1257e8f14785SLisandro Dalcin   PetscHSetI     stratumHash;
1258c58f1c22SToby Isaac   size_t         len = 0;
125926c55118SMichael Lange   PetscMPIInt    rank;
1260c58f1c22SToby Isaac   PetscErrorCode ierr;
1261c58f1c22SToby Isaac 
1262c58f1c22SToby Isaac   PetscFunctionBegin;
1263d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1264f018e600SMatthew G. Knepley   if (label) {
1265f018e600SMatthew G. Knepley     PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1266f018e600SMatthew G. Knepley     ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);
1267f018e600SMatthew G. Knepley   }
1268c58f1c22SToby Isaac   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
1269c58f1c22SToby Isaac   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1270c58f1c22SToby Isaac   /* Bcast name */
1271d67d17b1SMatthew G. Knepley   if (!rank) {
1272d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1273d67d17b1SMatthew G. Knepley     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1274d67d17b1SMatthew G. Knepley   }
1275c58f1c22SToby Isaac   nameSize = len;
1276c58f1c22SToby Isaac   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1277c58f1c22SToby Isaac   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1278d67d17b1SMatthew G. Knepley   if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);}
1279c58f1c22SToby Isaac   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1280d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
1281c58f1c22SToby Isaac   ierr = PetscFree(name);CHKERRQ(ierr);
128277d236dfSMichael Lange   /* Bcast defaultValue */
128377d236dfSMichael Lange   if (!rank) (*labelNew)->defaultValue = label->defaultValue;
128477d236dfSMichael Lange   ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
128526c55118SMichael Lange   /* Distribute stratum values over the SF and get the point mapping on the receiver */
128626c55118SMichael Lange   ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr);
12875cbdf6fcSMichael Lange   /* Determine received stratum values and initialise new label*/
1288e8f14785SLisandro Dalcin   ierr = PetscHSetICreate(&stratumHash);CHKERRQ(ierr);
128926c55118SMichael Lange   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
1290e8f14785SLisandro Dalcin   for (p = 0; p < size; ++p) {ierr = PetscHSetIAdd(stratumHash, leafStrata[p]);CHKERRQ(ierr);}
1291e8f14785SLisandro Dalcin   ierr = PetscHSetIGetSize(stratumHash, &(*labelNew)->numStrata);CHKERRQ(ierr);
1292ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->validIS);CHKERRQ(ierr);
1293ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->validIS[s] = PETSC_TRUE;
12945cbdf6fcSMichael Lange   ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr);
12955cbdf6fcSMichael Lange   /* Turn leafStrata into indices rather than stratum values */
12965cbdf6fcSMichael Lange   offset = 0;
1297e8f14785SLisandro Dalcin   ierr = PetscHSetIGetElems(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr);
1298a205f1fdSToby Isaac   ierr = PetscSortInt((*labelNew)->numStrata,(*labelNew)->stratumValues);CHKERRQ(ierr);
12995cbdf6fcSMichael Lange   for (p = 0; p < size; ++p) {
1300231b9e6fSMatthew G. Knepley     for (s = 0; s < (*labelNew)->numStrata; ++s) {
1301231b9e6fSMatthew G. Knepley       if (leafStrata[p] == (*labelNew)->stratumValues[s]) {leafStrata[p] = s; break;}
13025cbdf6fcSMichael Lange     }
13035cbdf6fcSMichael Lange   }
1304c58f1c22SToby Isaac   /* Rebuild the point strata on the receiver */
1305c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr);
1306c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1307c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1308c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1309c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1310c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1311c58f1c22SToby Isaac       (*labelNew)->stratumSizes[leafStrata[offset+s]]++;
1312c58f1c22SToby Isaac     }
1313c58f1c22SToby Isaac   }
1314c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr);
1315c58f1c22SToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr);
1316ad8374ffSToby Isaac   ierr = PetscMalloc1((*labelNew)->numStrata,&points);CHKERRQ(ierr);
1317c58f1c22SToby Isaac   for (s = 0; s < (*labelNew)->numStrata; ++s) {
1318e8f14785SLisandro Dalcin     ierr = PetscHSetICreate(&(*labelNew)->ht[s]);CHKERRQ(ierr);
1319ad8374ffSToby Isaac     ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(points[s]));CHKERRQ(ierr);
1320c58f1c22SToby Isaac   }
1321c58f1c22SToby Isaac   /* Insert points into new strata */
1322c58f1c22SToby Isaac   ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr);
1323c58f1c22SToby Isaac   ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr);
1324c58f1c22SToby Isaac   for (p=pStart; p<pEnd; p++) {
1325c58f1c22SToby Isaac     ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr);
1326c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr);
1327c58f1c22SToby Isaac     for (s=0; s<dof; s++) {
1328c58f1c22SToby Isaac       stratum = leafStrata[offset+s];
1329ad8374ffSToby Isaac       points[stratum][strataIdx[stratum]++] = p;
1330c58f1c22SToby Isaac     }
1331c58f1c22SToby Isaac   }
1332ad8374ffSToby Isaac   for (s = 0; s < (*labelNew)->numStrata; s++) {
1333ad8374ffSToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF,(*labelNew)->stratumSizes[s],&(points[s][0]),PETSC_OWN_POINTER,&((*labelNew)->points[s]));CHKERRQ(ierr);
1334ad8374ffSToby Isaac     ierr = PetscObjectSetName((PetscObject)((*labelNew)->points[s]),"indices");CHKERRQ(ierr);
1335ad8374ffSToby Isaac   }
1336ad8374ffSToby Isaac   ierr = PetscFree(points);CHKERRQ(ierr);
1337e8f14785SLisandro Dalcin   ierr = PetscHSetIDestroy(&stratumHash);CHKERRQ(ierr);
1338c58f1c22SToby Isaac   ierr = PetscFree(leafStrata);CHKERRQ(ierr);
1339c58f1c22SToby Isaac   ierr = PetscFree(strataIdx);CHKERRQ(ierr);
1340c58f1c22SToby Isaac   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
1341c58f1c22SToby Isaac   PetscFunctionReturn(0);
1342c58f1c22SToby Isaac }
1343c58f1c22SToby Isaac 
13447937d9ceSMichael Lange /*@
13457937d9ceSMichael Lange   DMLabelGather - Gather all label values from leafs into roots
13467937d9ceSMichael Lange 
13477937d9ceSMichael Lange   Input Parameters:
13487937d9ceSMichael Lange + label - the DMLabel
134984f0b6dfSMatthew G. Knepley - sf - the Star Forest point communication map
13507937d9ceSMichael Lange 
135184f0b6dfSMatthew G. Knepley   Output Parameters:
135284f0b6dfSMatthew G. Knepley . labelNew - the new DMLabel with localised leaf values
13537937d9ceSMichael Lange 
13547937d9ceSMichael Lange   Level: developer
13557937d9ceSMichael Lange 
13567937d9ceSMichael Lange   Note: This is the inverse operation to DMLabelDistribute.
13577937d9ceSMichael Lange 
13587937d9ceSMichael Lange .seealso: DMLabelDistribute()
13597937d9ceSMichael Lange @*/
13607937d9ceSMichael Lange PetscErrorCode DMLabelGather(DMLabel label, PetscSF sf, DMLabel *labelNew)
13617937d9ceSMichael Lange {
13627937d9ceSMichael Lange   MPI_Comm       comm;
13637937d9ceSMichael Lange   PetscSection   rootSection;
13647937d9ceSMichael Lange   PetscSF        sfLabel;
13657937d9ceSMichael Lange   PetscSFNode   *rootPoints, *leafPoints;
13667937d9ceSMichael Lange   PetscInt       p, s, d, nroots, nleaves, nmultiroots, idx, dof, offset;
13677937d9ceSMichael Lange   const PetscInt *rootDegree, *ilocal;
13687937d9ceSMichael Lange   PetscInt       *rootStrata;
1369d67d17b1SMatthew G. Knepley   const char    *lname;
13707937d9ceSMichael Lange   char          *name;
13717937d9ceSMichael Lange   PetscInt       nameSize;
13727937d9ceSMichael Lange   size_t         len = 0;
13739852e123SBarry Smith   PetscMPIInt    rank, size;
13747937d9ceSMichael Lange   PetscErrorCode ierr;
13757937d9ceSMichael Lange 
13767937d9ceSMichael Lange   PetscFunctionBegin;
1377d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1378d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
13797937d9ceSMichael Lange   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
13807937d9ceSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
13819852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
13827937d9ceSMichael Lange   /* Bcast name */
1383d67d17b1SMatthew G. Knepley   if (!rank) {
1384d67d17b1SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) label, &lname);CHKERRQ(ierr);
1385d67d17b1SMatthew G. Knepley     ierr = PetscStrlen(lname, &len);CHKERRQ(ierr);
1386d67d17b1SMatthew G. Knepley   }
13877937d9ceSMichael Lange   nameSize = len;
13887937d9ceSMichael Lange   ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
13897937d9ceSMichael Lange   ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr);
1390d67d17b1SMatthew G. Knepley   if (!rank) {ierr = PetscMemcpy(name, lname, nameSize+1);CHKERRQ(ierr);}
13917937d9ceSMichael Lange   ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1392d67d17b1SMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, name, labelNew);CHKERRQ(ierr);
13937937d9ceSMichael Lange   ierr = PetscFree(name);CHKERRQ(ierr);
13947937d9ceSMichael Lange   /* Gather rank/index pairs of leaves into local roots to build
13957937d9ceSMichael Lange      an inverse, multi-rooted SF. Note that this ignores local leaf
13967937d9ceSMichael Lange      indexing due to the use of the multiSF in PetscSFGather. */
13977937d9ceSMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL);CHKERRQ(ierr);
1398dc53bc9bSMatthew G. Knepley   ierr = PetscMalloc1(nroots, &leafPoints);CHKERRQ(ierr);
1399dc53bc9bSMatthew G. Knepley   for (p = 0; p < nroots; ++p) leafPoints[p].rank = leafPoints[p].index = -1;
14007937d9ceSMichael Lange   for (p = 0; p < nleaves; p++) {
14018212dd46SStefano Zampini     PetscInt ilp = ilocal ? ilocal[p] : p;
14028212dd46SStefano Zampini 
14038212dd46SStefano Zampini     leafPoints[ilp].index = ilp;
14048212dd46SStefano Zampini     leafPoints[ilp].rank  = rank;
14057937d9ceSMichael Lange   }
14067937d9ceSMichael Lange   ierr = PetscSFComputeDegreeBegin(sf, &rootDegree);CHKERRQ(ierr);
14077937d9ceSMichael Lange   ierr = PetscSFComputeDegreeEnd(sf, &rootDegree);CHKERRQ(ierr);
14087937d9ceSMichael Lange   for (p = 0, nmultiroots = 0; p < nroots; ++p) nmultiroots += rootDegree[p];
14097937d9ceSMichael Lange   ierr = PetscMalloc1(nmultiroots, &rootPoints);CHKERRQ(ierr);
14107937d9ceSMichael Lange   ierr = PetscSFGatherBegin(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
14117937d9ceSMichael Lange   ierr = PetscSFGatherEnd(sf, MPIU_2INT, leafPoints, rootPoints);CHKERRQ(ierr);
14127937d9ceSMichael Lange   ierr = PetscSFCreate(comm,& sfLabel);CHKERRQ(ierr);
14137937d9ceSMichael Lange   ierr = PetscSFSetGraph(sfLabel, nroots, nmultiroots, NULL, PETSC_OWN_POINTER, rootPoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
14147937d9ceSMichael Lange   /* Migrate label over inverted SF to pull stratum values at leaves into roots. */
14157937d9ceSMichael Lange   ierr = DMLabelDistribute_Internal(label, sfLabel, &rootSection, &rootStrata);CHKERRQ(ierr);
14167937d9ceSMichael Lange   /* Rebuild the point strata on the receiver */
14177937d9ceSMichael Lange   for (p = 0, idx = 0; p < nroots; p++) {
14187937d9ceSMichael Lange     for (d = 0; d < rootDegree[p]; d++) {
14197937d9ceSMichael Lange       ierr = PetscSectionGetDof(rootSection, idx+d, &dof);CHKERRQ(ierr);
14207937d9ceSMichael Lange       ierr = PetscSectionGetOffset(rootSection, idx+d, &offset);CHKERRQ(ierr);
14217937d9ceSMichael Lange       for (s = 0; s < dof; s++) {ierr = DMLabelSetValue(*labelNew, p, rootStrata[offset+s]);CHKERRQ(ierr);}
14227937d9ceSMichael Lange     }
14237937d9ceSMichael Lange     idx += rootDegree[p];
14247937d9ceSMichael Lange   }
142577e0c0e7SMichael Lange   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
142677e0c0e7SMichael Lange   ierr = PetscFree(rootStrata);CHKERRQ(ierr);
142777e0c0e7SMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
142877e0c0e7SMichael Lange   ierr = PetscSFDestroy(&sfLabel);CHKERRQ(ierr);
14297937d9ceSMichael Lange   PetscFunctionReturn(0);
14307937d9ceSMichael Lange }
14317937d9ceSMichael Lange 
143284f0b6dfSMatthew G. Knepley /*@
143384f0b6dfSMatthew G. Knepley   DMLabelConvertToSection - Make a PetscSection/IS pair that encodes the label
143484f0b6dfSMatthew G. Knepley 
143584f0b6dfSMatthew G. Knepley   Input Parameter:
143684f0b6dfSMatthew G. Knepley . label - the DMLabel
143784f0b6dfSMatthew G. Knepley 
143884f0b6dfSMatthew G. Knepley   Output Parameters:
143984f0b6dfSMatthew G. Knepley + section - the section giving offsets for each stratum
144084f0b6dfSMatthew G. Knepley - is - An IS containing all the label points
144184f0b6dfSMatthew G. Knepley 
144284f0b6dfSMatthew G. Knepley   Level: developer
144384f0b6dfSMatthew G. Knepley 
144484f0b6dfSMatthew G. Knepley .seealso: DMLabelDistribute()
144584f0b6dfSMatthew G. Knepley @*/
1446c58f1c22SToby Isaac PetscErrorCode DMLabelConvertToSection(DMLabel label, PetscSection *section, IS *is)
1447c58f1c22SToby Isaac {
1448c58f1c22SToby Isaac   IS              vIS;
1449c58f1c22SToby Isaac   const PetscInt *values;
1450c58f1c22SToby Isaac   PetscInt       *points;
1451c58f1c22SToby Isaac   PetscInt        nV, vS = 0, vE = 0, v, N;
1452c58f1c22SToby Isaac   PetscErrorCode  ierr;
1453c58f1c22SToby Isaac 
1454c58f1c22SToby Isaac   PetscFunctionBegin;
1455d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 1);
1456c58f1c22SToby Isaac   ierr = DMLabelGetNumValues(label, &nV);CHKERRQ(ierr);
1457c58f1c22SToby Isaac   ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
1458c58f1c22SToby Isaac   ierr = ISGetIndices(vIS, &values);CHKERRQ(ierr);
1459c58f1c22SToby Isaac   if (nV) {vS = values[0]; vE = values[0]+1;}
1460c58f1c22SToby Isaac   for (v = 1; v < nV; ++v) {
1461c58f1c22SToby Isaac     vS = PetscMin(vS, values[v]);
1462c58f1c22SToby Isaac     vE = PetscMax(vE, values[v]+1);
1463c58f1c22SToby Isaac   }
1464c58f1c22SToby Isaac   ierr = PetscSectionCreate(PETSC_COMM_SELF, section);CHKERRQ(ierr);
1465c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*section, vS, vE);CHKERRQ(ierr);
1466c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1467c58f1c22SToby Isaac     PetscInt n;
1468c58f1c22SToby Isaac 
1469c58f1c22SToby Isaac     ierr = DMLabelGetStratumSize(label, values[v], &n);CHKERRQ(ierr);
1470c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*section, values[v], n);CHKERRQ(ierr);
1471c58f1c22SToby Isaac   }
1472c58f1c22SToby Isaac   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
1473c58f1c22SToby Isaac   ierr = PetscSectionGetStorageSize(*section, &N);CHKERRQ(ierr);
1474c58f1c22SToby Isaac   ierr = PetscMalloc1(N, &points);CHKERRQ(ierr);
1475c58f1c22SToby Isaac   for (v = 0; v < nV; ++v) {
1476c58f1c22SToby Isaac     IS              is;
1477c58f1c22SToby Isaac     const PetscInt *spoints;
1478c58f1c22SToby Isaac     PetscInt        dof, off, p;
1479c58f1c22SToby Isaac 
1480c58f1c22SToby Isaac     ierr = PetscSectionGetDof(*section, values[v], &dof);CHKERRQ(ierr);
1481c58f1c22SToby Isaac     ierr = PetscSectionGetOffset(*section, values[v], &off);CHKERRQ(ierr);
1482c58f1c22SToby Isaac     ierr = DMLabelGetStratumIS(label, values[v], &is);CHKERRQ(ierr);
1483c58f1c22SToby Isaac     ierr = ISGetIndices(is, &spoints);CHKERRQ(ierr);
1484c58f1c22SToby Isaac     for (p = 0; p < dof; ++p) points[off+p] = spoints[p];
1485c58f1c22SToby Isaac     ierr = ISRestoreIndices(is, &spoints);CHKERRQ(ierr);
1486c58f1c22SToby Isaac     ierr = ISDestroy(&is);CHKERRQ(ierr);
1487c58f1c22SToby Isaac   }
1488c58f1c22SToby Isaac   ierr = ISRestoreIndices(vIS, &values);CHKERRQ(ierr);
1489c58f1c22SToby Isaac   ierr = ISDestroy(&vIS);CHKERRQ(ierr);
1490c58f1c22SToby Isaac   ierr = ISCreateGeneral(PETSC_COMM_SELF, N, points, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1491c58f1c22SToby Isaac   PetscFunctionReturn(0);
1492c58f1c22SToby Isaac }
1493c58f1c22SToby Isaac 
149484f0b6dfSMatthew G. Knepley /*@
1495c58f1c22SToby Isaac   PetscSectionCreateGlobalSectionLabel - Create a section describing the global field layout using
1496c58f1c22SToby Isaac   the local section and an SF describing the section point overlap.
1497c58f1c22SToby Isaac 
1498c58f1c22SToby Isaac   Input Parameters:
1499c58f1c22SToby Isaac   + s - The PetscSection for the local field layout
1500c58f1c22SToby Isaac   . sf - The SF describing parallel layout of the section points
1501c58f1c22SToby Isaac   . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1502c58f1c22SToby Isaac   . label - The label specifying the points
1503c58f1c22SToby Isaac   - labelValue - The label stratum specifying the points
1504c58f1c22SToby Isaac 
1505c58f1c22SToby Isaac   Output Parameter:
1506c58f1c22SToby Isaac   . gsection - The PetscSection for the global field layout
1507c58f1c22SToby Isaac 
1508c58f1c22SToby Isaac   Note: This gives negative sizes and offsets to points not owned by this process
1509c58f1c22SToby Isaac 
1510c58f1c22SToby Isaac   Level: developer
1511c58f1c22SToby Isaac 
1512c58f1c22SToby Isaac .seealso: PetscSectionCreate()
1513c58f1c22SToby Isaac @*/
1514c58f1c22SToby Isaac PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection s, PetscSF sf, PetscBool includeConstraints, DMLabel label, PetscInt labelValue, PetscSection *gsection)
1515c58f1c22SToby Isaac {
1516c58f1c22SToby Isaac   PetscInt      *neg = NULL, *tmpOff = NULL;
1517c58f1c22SToby Isaac   PetscInt       pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots;
1518c58f1c22SToby Isaac   PetscErrorCode ierr;
1519c58f1c22SToby Isaac 
1520c58f1c22SToby Isaac   PetscFunctionBegin;
1521d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1522d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1523d67d17b1SMatthew G. Knepley   PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 4);
1524c58f1c22SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr);
1525c58f1c22SToby Isaac   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1526c58f1c22SToby Isaac   ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr);
1527c58f1c22SToby Isaac   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1528c58f1c22SToby Isaac   if (nroots >= 0) {
1529c58f1c22SToby Isaac     if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
1530c58f1c22SToby Isaac     ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr);
1531c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1532c58f1c22SToby Isaac       ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr);
1533c58f1c22SToby Isaac     } else {
1534c58f1c22SToby Isaac       tmpOff = &(*gsection)->atlasDof[-pStart];
1535c58f1c22SToby Isaac     }
1536c58f1c22SToby Isaac   }
1537c58f1c22SToby Isaac   /* Mark ghost points with negative dof */
1538c58f1c22SToby Isaac   for (p = pStart; p < pEnd; ++p) {
1539c58f1c22SToby Isaac     PetscInt value;
1540c58f1c22SToby Isaac 
1541c58f1c22SToby Isaac     ierr = DMLabelGetValue(label, p, &value);CHKERRQ(ierr);
1542c58f1c22SToby Isaac     if (value != labelValue) continue;
1543c58f1c22SToby Isaac     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1544c58f1c22SToby Isaac     ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr);
1545c58f1c22SToby Isaac     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1546c58f1c22SToby Isaac     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);}
1547c58f1c22SToby Isaac     if (neg) neg[p] = -(dof+1);
1548c58f1c22SToby Isaac   }
1549c58f1c22SToby Isaac   ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr);
1550c58f1c22SToby Isaac   if (nroots >= 0) {
1551c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1552c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1553c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1554c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1555c58f1c22SToby Isaac     }
1556c58f1c22SToby Isaac   }
1557c58f1c22SToby Isaac   /* Calculate new sizes, get proccess offset, and calculate point offsets */
1558c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1559c58f1c22SToby Isaac     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
1560c58f1c22SToby Isaac     (*gsection)->atlasOff[p] = off;
1561c58f1c22SToby Isaac     off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
1562c58f1c22SToby Isaac   }
1563c58f1c22SToby Isaac   ierr       = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1564c58f1c22SToby Isaac   globalOff -= off;
1565c58f1c22SToby Isaac   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1566c58f1c22SToby Isaac     (*gsection)->atlasOff[p] += globalOff;
1567c58f1c22SToby Isaac     if (neg) neg[p] = -((*gsection)->atlasOff[p]+1);
1568c58f1c22SToby Isaac   }
1569c58f1c22SToby Isaac   /* Put in negative offsets for ghost points */
1570c58f1c22SToby Isaac   if (nroots >= 0) {
1571c58f1c22SToby Isaac     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1572c58f1c22SToby Isaac     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1573c58f1c22SToby Isaac     if (nroots > pEnd-pStart) {
1574c58f1c22SToby Isaac       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1575c58f1c22SToby Isaac     }
1576c58f1c22SToby Isaac   }
1577c58f1c22SToby Isaac   if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);}
1578c58f1c22SToby Isaac   ierr = PetscFree(neg);CHKERRQ(ierr);
1579c58f1c22SToby Isaac   PetscFunctionReturn(0);
1580c58f1c22SToby Isaac }
1581c58f1c22SToby Isaac 
15825fdea053SToby Isaac typedef struct _n_PetscSectionSym_Label
15835fdea053SToby Isaac {
15845fdea053SToby Isaac   DMLabel           label;
15855fdea053SToby Isaac   PetscCopyMode     *modes;
15865fdea053SToby Isaac   PetscInt          *sizes;
15875fdea053SToby Isaac   const PetscInt    ***perms;
15885fdea053SToby Isaac   const PetscScalar ***rots;
15895fdea053SToby Isaac   PetscInt          (*minMaxOrients)[2];
15905fdea053SToby Isaac   PetscInt          numStrata; /* numStrata is only increasing, functions as a state */
15915fdea053SToby Isaac } PetscSectionSym_Label;
15925fdea053SToby Isaac 
15935fdea053SToby Isaac static PetscErrorCode PetscSectionSymLabelReset(PetscSectionSym sym)
15945fdea053SToby Isaac {
15955fdea053SToby Isaac   PetscInt              i, j;
15965fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
15975fdea053SToby Isaac   PetscErrorCode        ierr;
15985fdea053SToby Isaac 
15995fdea053SToby Isaac   PetscFunctionBegin;
16005fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
16015fdea053SToby Isaac     if (sl->modes[i] == PETSC_OWN_POINTER || sl->modes[i] == PETSC_COPY_VALUES) {
16025fdea053SToby Isaac       for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
16035fdea053SToby Isaac         if (sl->perms[i]) {ierr = PetscFree(sl->perms[i][j]);CHKERRQ(ierr);}
16045fdea053SToby Isaac         if (sl->rots[i]) {ierr = PetscFree(sl->rots[i][j]);CHKERRQ(ierr);}
16055fdea053SToby Isaac       }
16065fdea053SToby Isaac       if (sl->perms[i]) {
16075fdea053SToby Isaac         const PetscInt **perms = &sl->perms[i][sl->minMaxOrients[i][0]];
16085fdea053SToby Isaac 
16095fdea053SToby Isaac         ierr = PetscFree(perms);CHKERRQ(ierr);
16105fdea053SToby Isaac       }
16115fdea053SToby Isaac       if (sl->rots[i]) {
16125fdea053SToby Isaac         const PetscScalar **rots = &sl->rots[i][sl->minMaxOrients[i][0]];
16135fdea053SToby Isaac 
16145fdea053SToby Isaac         ierr = PetscFree(rots);CHKERRQ(ierr);
16155fdea053SToby Isaac       }
16165fdea053SToby Isaac     }
16175fdea053SToby Isaac   }
16185fdea053SToby Isaac   ierr = PetscFree5(sl->modes,sl->sizes,sl->perms,sl->rots,sl->minMaxOrients);CHKERRQ(ierr);
16195fdea053SToby Isaac   ierr = DMLabelDestroy(&sl->label);CHKERRQ(ierr);
16205fdea053SToby Isaac   sl->numStrata = 0;
16215fdea053SToby Isaac   PetscFunctionReturn(0);
16225fdea053SToby Isaac }
16235fdea053SToby Isaac 
16245fdea053SToby Isaac static PetscErrorCode PetscSectionSymDestroy_Label(PetscSectionSym sym)
16255fdea053SToby Isaac {
16265fdea053SToby Isaac   PetscErrorCode ierr;
16275fdea053SToby Isaac 
16285fdea053SToby Isaac   PetscFunctionBegin;
16295fdea053SToby Isaac   ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);
16305fdea053SToby Isaac   ierr = PetscFree(sym->data);CHKERRQ(ierr);
16315fdea053SToby Isaac   PetscFunctionReturn(0);
16325fdea053SToby Isaac }
16335fdea053SToby Isaac 
16345fdea053SToby Isaac static PetscErrorCode PetscSectionSymView_Label(PetscSectionSym sym, PetscViewer viewer)
16355fdea053SToby Isaac {
16365fdea053SToby Isaac   PetscSectionSym_Label *sl = (PetscSectionSym_Label *) sym->data;
16375fdea053SToby Isaac   PetscBool             isAscii;
16385fdea053SToby Isaac   DMLabel               label = sl->label;
1639d67d17b1SMatthew G. Knepley   const char           *name;
16405fdea053SToby Isaac   PetscErrorCode        ierr;
16415fdea053SToby Isaac 
16425fdea053SToby Isaac   PetscFunctionBegin;
16435fdea053SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
16445fdea053SToby Isaac   if (isAscii) {
16455fdea053SToby Isaac     PetscInt          i, j, k;
16465fdea053SToby Isaac     PetscViewerFormat format;
16475fdea053SToby Isaac 
16485fdea053SToby Isaac     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
16495fdea053SToby Isaac     if (label) {
16505fdea053SToby Isaac       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
16515fdea053SToby Isaac       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
16525fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
16535fdea053SToby Isaac         ierr = DMLabelView(label, viewer);CHKERRQ(ierr);
16545fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
16555fdea053SToby Isaac       } else {
1656d67d17b1SMatthew G. Knepley         ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1657d67d17b1SMatthew G. Knepley         ierr = PetscViewerASCIIPrintf(viewer,"  Label '%s'\n",name);CHKERRQ(ierr);
16585fdea053SToby Isaac       }
16595fdea053SToby Isaac     } else {
16605fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "No label given\n");CHKERRQ(ierr);
16615fdea053SToby Isaac     }
16625fdea053SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
16635fdea053SToby Isaac     for (i = 0; i <= sl->numStrata; i++) {
16645fdea053SToby Isaac       PetscInt value = i < sl->numStrata ? label->stratumValues[i] : label->defaultValue;
16655fdea053SToby Isaac 
16665fdea053SToby Isaac       if (!(sl->perms[i] || sl->rots[i])) {
16675fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point): no symmetries\n", value, sl->sizes[i]);CHKERRQ(ierr);
16685fdea053SToby Isaac       } else {
16695fdea053SToby Isaac       ierr = PetscViewerASCIIPrintf(viewer, "Symmetry for stratum value %D (%D dofs per point):\n", value, sl->sizes[i]);CHKERRQ(ierr);
16705fdea053SToby Isaac         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
16715fdea053SToby Isaac         ierr = PetscViewerASCIIPrintf(viewer, "Orientation range: [%D, %D)\n", sl->minMaxOrients[i][0], sl->minMaxOrients[i][1]);CHKERRQ(ierr);
16725fdea053SToby Isaac         if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
16735fdea053SToby Isaac           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
16745fdea053SToby Isaac           for (j = sl->minMaxOrients[i][0]; j < sl->minMaxOrients[i][1]; j++) {
16755fdea053SToby Isaac             if (!((sl->perms[i] && sl->perms[i][j]) || (sl->rots[i] && sl->rots[i][j]))) {
16765fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D: identity\n",j);CHKERRQ(ierr);
16775fdea053SToby Isaac             } else {
16785fdea053SToby Isaac               PetscInt tab;
16795fdea053SToby Isaac 
16805fdea053SToby Isaac               ierr = PetscViewerASCIIPrintf(viewer, "Orientation %D:\n",j);CHKERRQ(ierr);
16815fdea053SToby Isaac               ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
16825fdea053SToby Isaac               ierr = PetscViewerASCIIGetTab(viewer,&tab);CHKERRQ(ierr);
16835fdea053SToby Isaac               if (sl->perms[i] && sl->perms[i][j]) {
16845fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Permutation:");CHKERRQ(ierr);
16855fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
16865fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %D",sl->perms[i][j][k]);CHKERRQ(ierr);}
16875fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
16885fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
16895fdea053SToby Isaac               }
16905fdea053SToby Isaac               if (sl->rots[i] && sl->rots[i][j]) {
16915fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"Rotations:  ");CHKERRQ(ierr);
16925fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,0);CHKERRQ(ierr);
16935fdea053SToby Isaac #if defined(PETSC_USE_COMPLEX)
16945fdea053SToby 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);}
16955fdea053SToby Isaac #else
16965fdea053SToby Isaac                 for (k = 0; k < sl->sizes[i]; k++) {ierr = PetscViewerASCIIPrintf(viewer," %+f",sl->rots[i][j][k]);CHKERRQ(ierr);}
16975fdea053SToby Isaac #endif
16985fdea053SToby Isaac                 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
16995fdea053SToby Isaac                 ierr = PetscViewerASCIISetTab(viewer,tab);CHKERRQ(ierr);
17005fdea053SToby Isaac               }
17015fdea053SToby Isaac               ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
17025fdea053SToby Isaac             }
17035fdea053SToby Isaac           }
17045fdea053SToby Isaac           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
17055fdea053SToby Isaac         }
17065fdea053SToby Isaac         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
17075fdea053SToby Isaac       }
17085fdea053SToby Isaac     }
17095fdea053SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
17105fdea053SToby Isaac   }
17115fdea053SToby Isaac   PetscFunctionReturn(0);
17125fdea053SToby Isaac }
17135fdea053SToby Isaac 
17145fdea053SToby Isaac /*@
17155fdea053SToby Isaac   PetscSectionSymLabelSetLabel - set the label whose strata will define the points that receive symmetries
17165fdea053SToby Isaac 
17175fdea053SToby Isaac   Logically collective on sym
17185fdea053SToby Isaac 
17195fdea053SToby Isaac   Input parameters:
17205fdea053SToby Isaac + sym - the section symmetries
17215fdea053SToby Isaac - label - the DMLabel describing the types of points
17225fdea053SToby Isaac 
17235fdea053SToby Isaac   Level: developer:
17245fdea053SToby Isaac 
17255fdea053SToby Isaac .seealso: PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()
17265fdea053SToby Isaac @*/
17275fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetLabel(PetscSectionSym sym, DMLabel label)
17285fdea053SToby Isaac {
17295fdea053SToby Isaac   PetscSectionSym_Label *sl;
17305fdea053SToby Isaac   PetscErrorCode        ierr;
17315fdea053SToby Isaac 
17325fdea053SToby Isaac   PetscFunctionBegin;
17335fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
17345fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
17355fdea053SToby Isaac   if (sl->label && sl->label != label) {ierr = PetscSectionSymLabelReset(sym);CHKERRQ(ierr);}
17365fdea053SToby Isaac   if (label) {
17375fdea053SToby Isaac     sl->label = label;
1738d67d17b1SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) label);CHKERRQ(ierr);
17395fdea053SToby Isaac     ierr = DMLabelGetNumValues(label,&sl->numStrata);CHKERRQ(ierr);
17401a834cf9SToby 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);
17411a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->modes,(sl->numStrata+1)*sizeof(PetscCopyMode));CHKERRQ(ierr);
17421a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->sizes,(sl->numStrata+1)*sizeof(PetscInt));CHKERRQ(ierr);
17431a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->perms,(sl->numStrata+1)*sizeof(const PetscInt **));CHKERRQ(ierr);
17441a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->rots,(sl->numStrata+1)*sizeof(const PetscScalar **));CHKERRQ(ierr);
17451a834cf9SToby Isaac     ierr = PetscMemzero((void *) sl->minMaxOrients,(sl->numStrata+1)*sizeof(PetscInt[2]));CHKERRQ(ierr);
17465fdea053SToby Isaac   }
17475fdea053SToby Isaac   PetscFunctionReturn(0);
17485fdea053SToby Isaac }
17495fdea053SToby Isaac 
17505fdea053SToby Isaac /*@C
17515fdea053SToby Isaac   PetscSectionSymLabelSetStratum - set the symmetries for the orientations of a stratum
17525fdea053SToby Isaac 
17535fdea053SToby Isaac   Logically collective on PetscSectionSym
17545fdea053SToby Isaac 
17555fdea053SToby Isaac   InputParameters:
17565fdea053SToby Isaac + sys       - the section symmetries
17575fdea053SToby Isaac . stratum   - the stratum value in the label that we are assigning symmetries for
17585fdea053SToby Isaac . size      - the number of dofs for points in the stratum of the label
17595fdea053SToby Isaac . minOrient - the smallest orientation for a point in this stratum
17605fdea053SToby Isaac . maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
17615fdea053SToby Isaac . mode      - how sym should copy the perms and rots arrays
17625fdea053SToby Isaac . perms     - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation.  A NULL permutation is the identity
17635fdea053SToby 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
17645fdea053SToby Isaac 
17655fdea053SToby Isaac   Level: developer
17665fdea053SToby Isaac 
17675fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()
17685fdea053SToby Isaac @*/
17695fdea053SToby Isaac PetscErrorCode PetscSectionSymLabelSetStratum(PetscSectionSym sym, PetscInt stratum, PetscInt size, PetscInt minOrient, PetscInt maxOrient, PetscCopyMode mode, const PetscInt **perms, const PetscScalar **rots)
17705fdea053SToby Isaac {
17715fdea053SToby Isaac   PetscSectionSym_Label *sl;
1772d67d17b1SMatthew G. Knepley   const char            *name;
1773d67d17b1SMatthew G. Knepley   PetscInt               i, j, k;
17745fdea053SToby Isaac   PetscErrorCode         ierr;
17755fdea053SToby Isaac 
17765fdea053SToby Isaac   PetscFunctionBegin;
17775fdea053SToby Isaac   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
17785fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
17795fdea053SToby Isaac   if (!sl->label) SETERRQ(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_WRONGSTATE,"No label set yet");
17805fdea053SToby Isaac   for (i = 0; i <= sl->numStrata; i++) {
17815fdea053SToby Isaac     PetscInt value = (i < sl->numStrata) ? sl->label->stratumValues[i] : sl->label->defaultValue;
17825fdea053SToby Isaac 
17835fdea053SToby Isaac     if (stratum == value) break;
17845fdea053SToby Isaac   }
1785d67d17b1SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) sl->label, &name);CHKERRQ(ierr);
1786d67d17b1SMatthew G. Knepley   if (i > sl->numStrata) SETERRQ2(PetscObjectComm((PetscObject)sym),PETSC_ERR_ARG_OUTOFRANGE,"Stratum %D not found in label %s\n",stratum,name);
17875fdea053SToby Isaac   sl->sizes[i] = size;
17885fdea053SToby Isaac   sl->modes[i] = mode;
17895fdea053SToby Isaac   sl->minMaxOrients[i][0] = minOrient;
17905fdea053SToby Isaac   sl->minMaxOrients[i][1] = maxOrient;
17915fdea053SToby Isaac   if (mode == PETSC_COPY_VALUES) {
17925fdea053SToby Isaac     if (perms) {
17935fdea053SToby Isaac       PetscInt    **ownPerms;
17945fdea053SToby Isaac 
17955fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownPerms);CHKERRQ(ierr);
17965fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
17975fdea053SToby Isaac         if (perms[j]) {
17985fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownPerms[j]);CHKERRQ(ierr);
17995fdea053SToby Isaac           for (k = 0; k < size; k++) {ownPerms[j][k] = perms[j][k];}
18005fdea053SToby Isaac         }
18015fdea053SToby Isaac       }
18025fdea053SToby Isaac       sl->perms[i] = (const PetscInt **) &ownPerms[-minOrient];
18035fdea053SToby Isaac     }
18045fdea053SToby Isaac     if (rots) {
18055fdea053SToby Isaac       PetscScalar **ownRots;
18065fdea053SToby Isaac 
18075fdea053SToby Isaac       ierr = PetscCalloc1(maxOrient - minOrient,&ownRots);CHKERRQ(ierr);
18085fdea053SToby Isaac       for (j = 0; j < maxOrient-minOrient; j++) {
18095fdea053SToby Isaac         if (rots[j]) {
18105fdea053SToby Isaac           ierr = PetscMalloc1(size,&ownRots[j]);CHKERRQ(ierr);
18115fdea053SToby Isaac           for (k = 0; k < size; k++) {ownRots[j][k] = rots[j][k];}
18125fdea053SToby Isaac         }
18135fdea053SToby Isaac       }
18145fdea053SToby Isaac       sl->rots[i] = (const PetscScalar **) &ownRots[-minOrient];
18155fdea053SToby Isaac     }
18165fdea053SToby Isaac   } else {
18175fdea053SToby Isaac     sl->perms[i] = perms ? &perms[-minOrient] : NULL;
18185fdea053SToby Isaac     sl->rots[i]  = rots ? &rots[-minOrient] : NULL;
18195fdea053SToby Isaac   }
18205fdea053SToby Isaac   PetscFunctionReturn(0);
18215fdea053SToby Isaac }
18225fdea053SToby Isaac 
18235fdea053SToby Isaac static PetscErrorCode PetscSectionSymGetPoints_Label(PetscSectionSym sym, PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt **perms, const PetscScalar **rots)
18245fdea053SToby Isaac {
18255fdea053SToby Isaac   PetscInt              i, j, numStrata;
18265fdea053SToby Isaac   PetscSectionSym_Label *sl;
18275fdea053SToby Isaac   DMLabel               label;
18285fdea053SToby Isaac   PetscErrorCode        ierr;
18295fdea053SToby Isaac 
18305fdea053SToby Isaac   PetscFunctionBegin;
18315fdea053SToby Isaac   sl = (PetscSectionSym_Label *) sym->data;
18325fdea053SToby Isaac   numStrata = sl->numStrata;
18335fdea053SToby Isaac   label     = sl->label;
18345fdea053SToby Isaac   for (i = 0; i < numPoints; i++) {
18355fdea053SToby Isaac     PetscInt point = points[2*i];
18365fdea053SToby Isaac     PetscInt ornt  = points[2*i+1];
18375fdea053SToby Isaac 
18385fdea053SToby Isaac     for (j = 0; j < numStrata; j++) {
18395fdea053SToby Isaac       if (label->validIS[j]) {
18405fdea053SToby Isaac         PetscInt k;
18415fdea053SToby Isaac 
18425fdea053SToby Isaac         ierr = ISLocate(label->points[j],point,&k);CHKERRQ(ierr);
18435fdea053SToby Isaac         if (k >= 0) break;
18445fdea053SToby Isaac       } else {
18455fdea053SToby Isaac         PetscBool has;
18465fdea053SToby Isaac 
1847e8f14785SLisandro Dalcin         PetscHSetIHas(label->ht[j], point, &has);
18485fdea053SToby Isaac         if (has) break;
18495fdea053SToby Isaac       }
18505fdea053SToby Isaac     }
18515fdea053SToby 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);
18525fdea053SToby Isaac     if (perms) {perms[i] = sl->perms[j] ? sl->perms[j][ornt] : NULL;}
18535fdea053SToby Isaac     if (rots) {rots[i]  = sl->rots[j] ? sl->rots[j][ornt] : NULL;}
18545fdea053SToby Isaac   }
18555fdea053SToby Isaac   PetscFunctionReturn(0);
18565fdea053SToby Isaac }
18575fdea053SToby Isaac 
18585fdea053SToby Isaac PetscErrorCode PetscSectionSymCreate_Label(PetscSectionSym sym)
18595fdea053SToby Isaac {
18605fdea053SToby Isaac   PetscSectionSym_Label *sl;
18615fdea053SToby Isaac   PetscErrorCode        ierr;
18625fdea053SToby Isaac 
18635fdea053SToby Isaac   PetscFunctionBegin;
18645fdea053SToby Isaac   ierr = PetscNewLog(sym,&sl);CHKERRQ(ierr);
18655fdea053SToby Isaac   sym->ops->getpoints = PetscSectionSymGetPoints_Label;
18665fdea053SToby Isaac   sym->ops->view      = PetscSectionSymView_Label;
18675fdea053SToby Isaac   sym->ops->destroy   = PetscSectionSymDestroy_Label;
18685fdea053SToby Isaac   sym->data           = (void *) sl;
18695fdea053SToby Isaac   PetscFunctionReturn(0);
18705fdea053SToby Isaac }
18715fdea053SToby Isaac 
18725fdea053SToby Isaac /*@
18735fdea053SToby Isaac   PetscSectionSymCreateLabel - Create a section symmetry that assigns one symmetry to each stratum of a label
18745fdea053SToby Isaac 
18755fdea053SToby Isaac   Collective
18765fdea053SToby Isaac 
18775fdea053SToby Isaac   Input Parameters:
18785fdea053SToby Isaac + comm - the MPI communicator for the new symmetry
18795fdea053SToby Isaac - label - the label defining the strata
18805fdea053SToby Isaac 
18815fdea053SToby Isaac   Output Parameters:
18825fdea053SToby Isaac . sym - the section symmetries
18835fdea053SToby Isaac 
18845fdea053SToby Isaac   Level: developer
18855fdea053SToby Isaac 
18865fdea053SToby Isaac .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()
18875fdea053SToby Isaac @*/
18885fdea053SToby Isaac PetscErrorCode PetscSectionSymCreateLabel(MPI_Comm comm, DMLabel label, PetscSectionSym *sym)
18895fdea053SToby Isaac {
18905fdea053SToby Isaac   PetscErrorCode        ierr;
18915fdea053SToby Isaac 
18925fdea053SToby Isaac   PetscFunctionBegin;
18935fdea053SToby Isaac   ierr = DMInitializePackage();CHKERRQ(ierr);
18945fdea053SToby Isaac   ierr = PetscSectionSymCreate(comm,sym);CHKERRQ(ierr);
18955fdea053SToby Isaac   ierr = PetscSectionSymSetType(*sym,PETSCSECTIONSYMLABEL);CHKERRQ(ierr);
18965fdea053SToby Isaac   ierr = PetscSectionSymLabelSetLabel(*sym,label);CHKERRQ(ierr);
18975fdea053SToby Isaac   PetscFunctionReturn(0);
18985fdea053SToby Isaac }
1899