xref: /petsc/src/vec/is/section/interface/section.c (revision fe2efc57b7777594bce7568e90822861480cbdc8)
1ea844a1aSMatthew Knepley /*
2ea844a1aSMatthew Knepley    This file contains routines for basic section object implementation.
3ea844a1aSMatthew Knepley */
4ea844a1aSMatthew Knepley 
5ea844a1aSMatthew Knepley #include <petsc/private/sectionimpl.h>   /*I  "petscsection.h"   I*/
6ea844a1aSMatthew Knepley #include <petscsection.h>
7ea844a1aSMatthew Knepley #include <petscsf.h>
8ea844a1aSMatthew Knepley #include <petscviewer.h>
9ea844a1aSMatthew Knepley 
10ea844a1aSMatthew Knepley PetscClassId PETSC_SECTION_CLASSID;
11ea844a1aSMatthew Knepley 
12ea844a1aSMatthew Knepley /*@
13ea844a1aSMatthew Knepley   PetscSectionCreate - Allocates PetscSection space and sets the map contents to the default.
14ea844a1aSMatthew Knepley 
15ea844a1aSMatthew Knepley   Collective
16ea844a1aSMatthew Knepley 
17ea844a1aSMatthew Knepley   Input Parameters:
18ea844a1aSMatthew Knepley + comm - the MPI communicator
19ea844a1aSMatthew Knepley - s    - pointer to the section
20ea844a1aSMatthew Knepley 
21ea844a1aSMatthew Knepley   Level: beginner
22ea844a1aSMatthew Knepley 
23ea844a1aSMatthew Knepley   Notes:
24ea844a1aSMatthew Knepley   Typical calling sequence
25ea844a1aSMatthew Knepley $       PetscSectionCreate(MPI_Comm,PetscSection *);
26ea844a1aSMatthew Knepley $       PetscSectionSetNumFields(PetscSection, numFields);
27ea844a1aSMatthew Knepley $       PetscSectionSetChart(PetscSection,low,high);
28ea844a1aSMatthew Knepley $       PetscSectionSetDof(PetscSection,point,numdof);
29ea844a1aSMatthew Knepley $       PetscSectionSetUp(PetscSection);
30ea844a1aSMatthew Knepley $       PetscSectionGetOffset(PetscSection,point,PetscInt *);
31ea844a1aSMatthew Knepley $       PetscSectionDestroy(PetscSection);
32ea844a1aSMatthew Knepley 
33ea844a1aSMatthew Knepley   The PetscSection object and methods are intended to be used in the PETSc Vec and Mat implementions; it is
34ea844a1aSMatthew Knepley   recommended they not be used in user codes unless you really gain something in their use.
35ea844a1aSMatthew Knepley 
36ea844a1aSMatthew Knepley .seealso: PetscSection, PetscSectionDestroy()
37ea844a1aSMatthew Knepley @*/
38ea844a1aSMatthew Knepley PetscErrorCode PetscSectionCreate(MPI_Comm comm, PetscSection *s)
39ea844a1aSMatthew Knepley {
40ea844a1aSMatthew Knepley   PetscErrorCode ierr;
41ea844a1aSMatthew Knepley 
42ea844a1aSMatthew Knepley   PetscFunctionBegin;
43ea844a1aSMatthew Knepley   PetscValidPointer(s,2);
44ea844a1aSMatthew Knepley   ierr = ISInitializePackage();CHKERRQ(ierr);
45ea844a1aSMatthew Knepley 
46ea844a1aSMatthew Knepley   ierr = PetscHeaderCreate(*s,PETSC_SECTION_CLASSID,"PetscSection","Section","IS",comm,PetscSectionDestroy,PetscSectionView);CHKERRQ(ierr);
47ea844a1aSMatthew Knepley 
48ea844a1aSMatthew Knepley   (*s)->pStart             = -1;
49ea844a1aSMatthew Knepley   (*s)->pEnd               = -1;
50ea844a1aSMatthew Knepley   (*s)->perm               = NULL;
51ea844a1aSMatthew Knepley   (*s)->pointMajor         = PETSC_TRUE;
52ea844a1aSMatthew Knepley   (*s)->maxDof             = 0;
53ea844a1aSMatthew Knepley   (*s)->atlasDof           = NULL;
54ea844a1aSMatthew Knepley   (*s)->atlasOff           = NULL;
55ea844a1aSMatthew Knepley   (*s)->bc                 = NULL;
56ea844a1aSMatthew Knepley   (*s)->bcIndices          = NULL;
57ea844a1aSMatthew Knepley   (*s)->setup              = PETSC_FALSE;
58ea844a1aSMatthew Knepley   (*s)->numFields          = 0;
59ea844a1aSMatthew Knepley   (*s)->fieldNames         = NULL;
60ea844a1aSMatthew Knepley   (*s)->field              = NULL;
61ea844a1aSMatthew Knepley   (*s)->useFieldOff        = PETSC_FALSE;
62ea844a1aSMatthew Knepley   (*s)->clObj              = NULL;
63ea844a1aSMatthew Knepley   (*s)->clSection          = NULL;
64ea844a1aSMatthew Knepley   (*s)->clPoints           = NULL;
65ea844a1aSMatthew Knepley   (*s)->clSize             = 0;
66ea844a1aSMatthew Knepley   (*s)->clPerm             = NULL;
67ea844a1aSMatthew Knepley   (*s)->clInvPerm          = NULL;
68ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
69ea844a1aSMatthew Knepley }
70ea844a1aSMatthew Knepley 
71ea844a1aSMatthew Knepley /*@
72ea844a1aSMatthew Knepley   PetscSectionCopy - Creates a shallow (if possible) copy of the PetscSection
73ea844a1aSMatthew Knepley 
74ea844a1aSMatthew Knepley   Collective
75ea844a1aSMatthew Knepley 
76ea844a1aSMatthew Knepley   Input Parameter:
77ea844a1aSMatthew Knepley . section - the PetscSection
78ea844a1aSMatthew Knepley 
79ea844a1aSMatthew Knepley   Output Parameter:
80ea844a1aSMatthew Knepley . newSection - the copy
81ea844a1aSMatthew Knepley 
82ea844a1aSMatthew Knepley   Level: intermediate
83ea844a1aSMatthew Knepley 
84ea844a1aSMatthew Knepley .seealso: PetscSection, PetscSectionCreate(), PetscSectionDestroy()
85ea844a1aSMatthew Knepley @*/
86ea844a1aSMatthew Knepley PetscErrorCode PetscSectionCopy(PetscSection section, PetscSection newSection)
87ea844a1aSMatthew Knepley {
88ea844a1aSMatthew Knepley   PetscSectionSym sym;
89ea844a1aSMatthew Knepley   IS              perm;
90ea844a1aSMatthew Knepley   PetscInt        numFields, f, pStart, pEnd, p;
91ea844a1aSMatthew Knepley   PetscErrorCode  ierr;
92ea844a1aSMatthew Knepley 
93ea844a1aSMatthew Knepley   PetscFunctionBegin;
94ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1);
95ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(newSection, PETSC_SECTION_CLASSID, 2);
96ea844a1aSMatthew Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
97ea844a1aSMatthew Knepley   if (numFields) {ierr = PetscSectionSetNumFields(newSection, numFields);CHKERRQ(ierr);}
98ea844a1aSMatthew Knepley   for (f = 0; f < numFields; ++f) {
99ea844a1aSMatthew Knepley     const char *name   = NULL;
100ea844a1aSMatthew Knepley     PetscInt   numComp = 0;
101ea844a1aSMatthew Knepley 
102ea844a1aSMatthew Knepley     ierr = PetscSectionGetFieldName(section, f, &name);CHKERRQ(ierr);
103ea844a1aSMatthew Knepley     ierr = PetscSectionSetFieldName(newSection, f, name);CHKERRQ(ierr);
104ea844a1aSMatthew Knepley     ierr = PetscSectionGetFieldComponents(section, f, &numComp);CHKERRQ(ierr);
105ea844a1aSMatthew Knepley     ierr = PetscSectionSetFieldComponents(newSection, f, numComp);CHKERRQ(ierr);
106ea844a1aSMatthew Knepley     ierr = PetscSectionGetFieldSym(section, f, &sym);CHKERRQ(ierr);
107ea844a1aSMatthew Knepley     ierr = PetscSectionSetFieldSym(newSection, f, sym);CHKERRQ(ierr);
108ea844a1aSMatthew Knepley   }
109ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
110ea844a1aSMatthew Knepley   ierr = PetscSectionSetChart(newSection, pStart, pEnd);CHKERRQ(ierr);
111ea844a1aSMatthew Knepley   ierr = PetscSectionGetPermutation(section, &perm);CHKERRQ(ierr);
112ea844a1aSMatthew Knepley   ierr = PetscSectionSetPermutation(newSection, perm);CHKERRQ(ierr);
113ea844a1aSMatthew Knepley   ierr = PetscSectionGetSym(section, &sym);CHKERRQ(ierr);
114ea844a1aSMatthew Knepley   ierr = PetscSectionSetSym(newSection, sym);CHKERRQ(ierr);
115ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
116ea844a1aSMatthew Knepley     PetscInt dof, cdof, fcdof = 0;
117ea844a1aSMatthew Knepley 
118ea844a1aSMatthew Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
119ea844a1aSMatthew Knepley     ierr = PetscSectionSetDof(newSection, p, dof);CHKERRQ(ierr);
120ea844a1aSMatthew Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
121ea844a1aSMatthew Knepley     if (cdof) {ierr = PetscSectionSetConstraintDof(newSection, p, cdof);CHKERRQ(ierr);}
122ea844a1aSMatthew Knepley     for (f = 0; f < numFields; ++f) {
123ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldDof(section, p, f, &dof);CHKERRQ(ierr);
124ea844a1aSMatthew Knepley       ierr = PetscSectionSetFieldDof(newSection, p, f, dof);CHKERRQ(ierr);
125ea844a1aSMatthew Knepley       if (cdof) {
126ea844a1aSMatthew Knepley         ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr);
127ea844a1aSMatthew Knepley         if (fcdof) {ierr = PetscSectionSetFieldConstraintDof(newSection, p, f, fcdof);CHKERRQ(ierr);}
128ea844a1aSMatthew Knepley       }
129ea844a1aSMatthew Knepley     }
130ea844a1aSMatthew Knepley   }
131ea844a1aSMatthew Knepley   ierr = PetscSectionSetUp(newSection);CHKERRQ(ierr);
132ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
133ea844a1aSMatthew Knepley     PetscInt       off, cdof, fcdof = 0;
134ea844a1aSMatthew Knepley     const PetscInt *cInd;
135ea844a1aSMatthew Knepley 
136ea844a1aSMatthew Knepley     /* Must set offsets in case they do not agree with the prefix sums */
137ea844a1aSMatthew Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
138ea844a1aSMatthew Knepley     ierr = PetscSectionSetOffset(newSection, p, off);CHKERRQ(ierr);
139ea844a1aSMatthew Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
140ea844a1aSMatthew Knepley     if (cdof) {
141ea844a1aSMatthew Knepley       ierr = PetscSectionGetConstraintIndices(section, p, &cInd);CHKERRQ(ierr);
142ea844a1aSMatthew Knepley       ierr = PetscSectionSetConstraintIndices(newSection, p, cInd);CHKERRQ(ierr);
143ea844a1aSMatthew Knepley       for (f = 0; f < numFields; ++f) {
144ea844a1aSMatthew Knepley         ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr);
145ea844a1aSMatthew Knepley         if (fcdof) {
146ea844a1aSMatthew Knepley           ierr = PetscSectionGetFieldConstraintIndices(section, p, f, &cInd);CHKERRQ(ierr);
147ea844a1aSMatthew Knepley           ierr = PetscSectionSetFieldConstraintIndices(newSection, p, f, cInd);CHKERRQ(ierr);
148ea844a1aSMatthew Knepley         }
149ea844a1aSMatthew Knepley       }
150ea844a1aSMatthew Knepley     }
151ea844a1aSMatthew Knepley   }
152ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
153ea844a1aSMatthew Knepley }
154ea844a1aSMatthew Knepley 
155ea844a1aSMatthew Knepley /*@
156ea844a1aSMatthew Knepley   PetscSectionClone - Creates a shallow (if possible) copy of the PetscSection
157ea844a1aSMatthew Knepley 
158ea844a1aSMatthew Knepley   Collective
159ea844a1aSMatthew Knepley 
160ea844a1aSMatthew Knepley   Input Parameter:
161ea844a1aSMatthew Knepley . section - the PetscSection
162ea844a1aSMatthew Knepley 
163ea844a1aSMatthew Knepley   Output Parameter:
164ea844a1aSMatthew Knepley . newSection - the copy
165ea844a1aSMatthew Knepley 
166ea844a1aSMatthew Knepley   Level: beginner
167ea844a1aSMatthew Knepley 
168ea844a1aSMatthew Knepley .seealso: PetscSection, PetscSectionCreate(), PetscSectionDestroy()
169ea844a1aSMatthew Knepley @*/
170ea844a1aSMatthew Knepley PetscErrorCode PetscSectionClone(PetscSection section, PetscSection *newSection)
171ea844a1aSMatthew Knepley {
172ea844a1aSMatthew Knepley   PetscErrorCode ierr;
173ea844a1aSMatthew Knepley 
174ea844a1aSMatthew Knepley   PetscFunctionBegin;
175ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1);
176ea844a1aSMatthew Knepley   PetscValidPointer(newSection, 2);
177ea844a1aSMatthew Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) section), newSection);CHKERRQ(ierr);
178ea844a1aSMatthew Knepley   ierr = PetscSectionCopy(section, *newSection);CHKERRQ(ierr);
179ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
180ea844a1aSMatthew Knepley }
181ea844a1aSMatthew Knepley 
182ea844a1aSMatthew Knepley /*@
183ea844a1aSMatthew Knepley   PetscSectionSetFromOptions - sets parameters in a PetscSection from the options database
184ea844a1aSMatthew Knepley 
185ea844a1aSMatthew Knepley   Collective on PetscSection
186ea844a1aSMatthew Knepley 
187ea844a1aSMatthew Knepley   Input Parameter:
188ea844a1aSMatthew Knepley . section - the PetscSection
189ea844a1aSMatthew Knepley 
190ea844a1aSMatthew Knepley   Options Database:
191ea844a1aSMatthew Knepley . -petscsection_point_major the dof order
192ea844a1aSMatthew Knepley 
193ea844a1aSMatthew Knepley   Level: intermediate
194ea844a1aSMatthew Knepley 
195ea844a1aSMatthew Knepley .seealso: PetscSection, PetscSectionCreate(), PetscSectionDestroy()
196ea844a1aSMatthew Knepley @*/
197ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetFromOptions(PetscSection s)
198ea844a1aSMatthew Knepley {
199ea844a1aSMatthew Knepley   PetscErrorCode ierr;
200ea844a1aSMatthew Knepley 
201ea844a1aSMatthew Knepley   PetscFunctionBegin;
202ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
203ea844a1aSMatthew Knepley   ierr = PetscObjectOptionsBegin((PetscObject) s);CHKERRQ(ierr);
204ea844a1aSMatthew Knepley   ierr = PetscOptionsBool("-petscsection_point_major", "The for ordering, either point major or field major", "PetscSectionSetPointMajor", s->pointMajor, &s->pointMajor, NULL);CHKERRQ(ierr);
205ea844a1aSMatthew Knepley   /* process any options handlers added with PetscObjectAddOptionsHandler() */
206ea844a1aSMatthew Knepley   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) s);CHKERRQ(ierr);
207ea844a1aSMatthew Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
208ea844a1aSMatthew Knepley   ierr = PetscObjectViewFromOptions((PetscObject) s, NULL, "-petscsection_view");CHKERRQ(ierr);
209ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
210ea844a1aSMatthew Knepley }
211ea844a1aSMatthew Knepley 
212ea844a1aSMatthew Knepley /*@
213ea844a1aSMatthew Knepley   PetscSectionCompare - Compares two sections
214ea844a1aSMatthew Knepley 
215ea844a1aSMatthew Knepley   Collective on PetscSection
216ea844a1aSMatthew Knepley 
217ea844a1aSMatthew Knepley   Input Parameterss:
218ea844a1aSMatthew Knepley + s1 - the first PetscSection
219ea844a1aSMatthew Knepley - s2 - the second PetscSection
220ea844a1aSMatthew Knepley 
221ea844a1aSMatthew Knepley   Output Parameter:
222ea844a1aSMatthew Knepley . congruent - PETSC_TRUE if the two sections are congruent, PETSC_FALSE otherwise
223ea844a1aSMatthew Knepley 
224ea844a1aSMatthew Knepley   Level: intermediate
225ea844a1aSMatthew Knepley 
226ea844a1aSMatthew Knepley   Notes:
227ea844a1aSMatthew Knepley   Field names are disregarded.
228ea844a1aSMatthew Knepley 
229ea844a1aSMatthew Knepley .seealso: PetscSection, PetscSectionCreate(), PetscSectionCopy(), PetscSectionClone()
230ea844a1aSMatthew Knepley @*/
231ea844a1aSMatthew Knepley PetscErrorCode PetscSectionCompare(PetscSection s1, PetscSection s2, PetscBool *congruent)
232ea844a1aSMatthew Knepley {
233ea844a1aSMatthew Knepley   PetscInt pStart, pEnd, nfields, ncdof, nfcdof, p, f, n1, n2;
234ea844a1aSMatthew Knepley   const PetscInt *idx1, *idx2;
235ea844a1aSMatthew Knepley   IS perm1, perm2;
236ea844a1aSMatthew Knepley   PetscBool flg;
237ea844a1aSMatthew Knepley   PetscMPIInt mflg;
238ea844a1aSMatthew Knepley   PetscErrorCode ierr;
239ea844a1aSMatthew Knepley 
240ea844a1aSMatthew Knepley   PetscFunctionBegin;
241ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s1, PETSC_SECTION_CLASSID, 1);
242ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s2, PETSC_SECTION_CLASSID, 2);
243ea844a1aSMatthew Knepley   PetscValidIntPointer(congruent,3);
244ea844a1aSMatthew Knepley   flg = PETSC_FALSE;
245ea844a1aSMatthew Knepley 
246ea844a1aSMatthew Knepley   ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)s1),PetscObjectComm((PetscObject)s2),&mflg);CHKERRQ(ierr);
247ea844a1aSMatthew Knepley   if (mflg != MPI_CONGRUENT && mflg != MPI_IDENT) {
248ea844a1aSMatthew Knepley     *congruent = PETSC_FALSE;
249ea844a1aSMatthew Knepley     PetscFunctionReturn(0);
250ea844a1aSMatthew Knepley   }
251ea844a1aSMatthew Knepley 
252ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(s1, &pStart, &pEnd);CHKERRQ(ierr);
253ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(s2, &n1, &n2);CHKERRQ(ierr);
254ea844a1aSMatthew Knepley   if (pStart != n1 || pEnd != n2) goto not_congruent;
255ea844a1aSMatthew Knepley 
256ea844a1aSMatthew Knepley   ierr = PetscSectionGetPermutation(s1, &perm1);CHKERRQ(ierr);
257ea844a1aSMatthew Knepley   ierr = PetscSectionGetPermutation(s2, &perm2);CHKERRQ(ierr);
258ea844a1aSMatthew Knepley   if (perm1 && perm2) {
259ea844a1aSMatthew Knepley     ierr = ISEqual(perm1, perm2, congruent);CHKERRQ(ierr);
260ea844a1aSMatthew Knepley     if (!(*congruent)) goto not_congruent;
261ea844a1aSMatthew Knepley   } else if (perm1 != perm2) goto not_congruent;
262ea844a1aSMatthew Knepley 
263ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
264ea844a1aSMatthew Knepley     ierr = PetscSectionGetOffset(s1, p, &n1);CHKERRQ(ierr);
265ea844a1aSMatthew Knepley     ierr = PetscSectionGetOffset(s2, p, &n2);CHKERRQ(ierr);
266ea844a1aSMatthew Knepley     if (n1 != n2) goto not_congruent;
267ea844a1aSMatthew Knepley 
268ea844a1aSMatthew Knepley     ierr = PetscSectionGetDof(s1, p, &n1);CHKERRQ(ierr);
269ea844a1aSMatthew Knepley     ierr = PetscSectionGetDof(s2, p, &n2);CHKERRQ(ierr);
270ea844a1aSMatthew Knepley     if (n1 != n2) goto not_congruent;
271ea844a1aSMatthew Knepley 
272ea844a1aSMatthew Knepley     ierr = PetscSectionGetConstraintDof(s1, p, &ncdof);CHKERRQ(ierr);
273ea844a1aSMatthew Knepley     ierr = PetscSectionGetConstraintDof(s2, p, &n2);CHKERRQ(ierr);
274ea844a1aSMatthew Knepley     if (ncdof != n2) goto not_congruent;
275ea844a1aSMatthew Knepley 
276ea844a1aSMatthew Knepley     ierr = PetscSectionGetConstraintIndices(s1, p, &idx1);CHKERRQ(ierr);
277ea844a1aSMatthew Knepley     ierr = PetscSectionGetConstraintIndices(s2, p, &idx2);CHKERRQ(ierr);
278ea844a1aSMatthew Knepley     ierr = PetscArraycmp(idx1, idx2, ncdof, congruent);CHKERRQ(ierr);
279ea844a1aSMatthew Knepley     if (!(*congruent)) goto not_congruent;
280ea844a1aSMatthew Knepley   }
281ea844a1aSMatthew Knepley 
282ea844a1aSMatthew Knepley   ierr = PetscSectionGetNumFields(s1, &nfields);CHKERRQ(ierr);
283ea844a1aSMatthew Knepley   ierr = PetscSectionGetNumFields(s2, &n2);CHKERRQ(ierr);
284ea844a1aSMatthew Knepley   if (nfields != n2) goto not_congruent;
285ea844a1aSMatthew Knepley 
286ea844a1aSMatthew Knepley   for (f = 0; f < nfields; ++f) {
287ea844a1aSMatthew Knepley     ierr = PetscSectionGetFieldComponents(s1, f, &n1);CHKERRQ(ierr);
288ea844a1aSMatthew Knepley     ierr = PetscSectionGetFieldComponents(s2, f, &n2);CHKERRQ(ierr);
289ea844a1aSMatthew Knepley     if (n1 != n2) goto not_congruent;
290ea844a1aSMatthew Knepley 
291ea844a1aSMatthew Knepley     for (p = pStart; p < pEnd; ++p) {
292ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldOffset(s1, p, f, &n1);CHKERRQ(ierr);
293ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldOffset(s2, p, f, &n2);CHKERRQ(ierr);
294ea844a1aSMatthew Knepley       if (n1 != n2) goto not_congruent;
295ea844a1aSMatthew Knepley 
296ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldDof(s1, p, f, &n1);CHKERRQ(ierr);
297ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldDof(s2, p, f, &n2);CHKERRQ(ierr);
298ea844a1aSMatthew Knepley       if (n1 != n2) goto not_congruent;
299ea844a1aSMatthew Knepley 
300ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldConstraintDof(s1, p, f, &nfcdof);CHKERRQ(ierr);
301ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldConstraintDof(s2, p, f, &n2);CHKERRQ(ierr);
302ea844a1aSMatthew Knepley       if (nfcdof != n2) goto not_congruent;
303ea844a1aSMatthew Knepley 
304ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldConstraintIndices(s1, p, f, &idx1);CHKERRQ(ierr);
305ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldConstraintIndices(s2, p, f, &idx2);CHKERRQ(ierr);
306ea844a1aSMatthew Knepley       ierr = PetscArraycmp(idx1, idx2, nfcdof, congruent);CHKERRQ(ierr);
307ea844a1aSMatthew Knepley       if (!(*congruent)) goto not_congruent;
308ea844a1aSMatthew Knepley     }
309ea844a1aSMatthew Knepley   }
310ea844a1aSMatthew Knepley 
311ea844a1aSMatthew Knepley   flg = PETSC_TRUE;
312ea844a1aSMatthew Knepley not_congruent:
313ea844a1aSMatthew Knepley   ierr = MPIU_Allreduce(&flg,congruent,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)s1));CHKERRQ(ierr);
314ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
315ea844a1aSMatthew Knepley }
316ea844a1aSMatthew Knepley 
317ea844a1aSMatthew Knepley /*@
318ea844a1aSMatthew Knepley   PetscSectionGetNumFields - Returns the number of fields, or 0 if no fields were defined.
319ea844a1aSMatthew Knepley 
320ea844a1aSMatthew Knepley   Not collective
321ea844a1aSMatthew Knepley 
322ea844a1aSMatthew Knepley   Input Parameter:
323ea844a1aSMatthew Knepley . s - the PetscSection
324ea844a1aSMatthew Knepley 
325ea844a1aSMatthew Knepley   Output Parameter:
326ea844a1aSMatthew Knepley . numFields - the number of fields defined, or 0 if none were defined
327ea844a1aSMatthew Knepley 
328ea844a1aSMatthew Knepley   Level: intermediate
329ea844a1aSMatthew Knepley 
330ea844a1aSMatthew Knepley .seealso: PetscSectionSetNumFields()
331ea844a1aSMatthew Knepley @*/
332ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetNumFields(PetscSection s, PetscInt *numFields)
333ea844a1aSMatthew Knepley {
334ea844a1aSMatthew Knepley   PetscFunctionBegin;
335ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
336ea844a1aSMatthew Knepley   PetscValidPointer(numFields,2);
337ea844a1aSMatthew Knepley   *numFields = s->numFields;
338ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
339ea844a1aSMatthew Knepley }
340ea844a1aSMatthew Knepley 
341ea844a1aSMatthew Knepley /*@
342ea844a1aSMatthew Knepley   PetscSectionSetNumFields - Sets the number of fields.
343ea844a1aSMatthew Knepley 
344ea844a1aSMatthew Knepley   Not collective
345ea844a1aSMatthew Knepley 
346ea844a1aSMatthew Knepley   Input Parameters:
347ea844a1aSMatthew Knepley + s - the PetscSection
348ea844a1aSMatthew Knepley - numFields - the number of fields
349ea844a1aSMatthew Knepley 
350ea844a1aSMatthew Knepley   Level: intermediate
351ea844a1aSMatthew Knepley 
352ea844a1aSMatthew Knepley .seealso: PetscSectionGetNumFields()
353ea844a1aSMatthew Knepley @*/
354ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetNumFields(PetscSection s, PetscInt numFields)
355ea844a1aSMatthew Knepley {
356ea844a1aSMatthew Knepley   PetscInt       f;
357ea844a1aSMatthew Knepley   PetscErrorCode ierr;
358ea844a1aSMatthew Knepley 
359ea844a1aSMatthew Knepley   PetscFunctionBegin;
360ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
361ea844a1aSMatthew Knepley   if (numFields <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "The number of fields %D must be positive", numFields);
362ea844a1aSMatthew Knepley   ierr = PetscSectionReset(s);CHKERRQ(ierr);
363ea844a1aSMatthew Knepley 
364ea844a1aSMatthew Knepley   s->numFields = numFields;
365ea844a1aSMatthew Knepley   ierr = PetscMalloc1(s->numFields, &s->numFieldComponents);CHKERRQ(ierr);
366ea844a1aSMatthew Knepley   ierr = PetscMalloc1(s->numFields, &s->fieldNames);CHKERRQ(ierr);
367ea844a1aSMatthew Knepley   ierr = PetscMalloc1(s->numFields, &s->field);CHKERRQ(ierr);
368ea844a1aSMatthew Knepley   for (f = 0; f < s->numFields; ++f) {
369ea844a1aSMatthew Knepley     char name[64];
370ea844a1aSMatthew Knepley 
371ea844a1aSMatthew Knepley     s->numFieldComponents[f] = 1;
372ea844a1aSMatthew Knepley 
373ea844a1aSMatthew Knepley     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), &s->field[f]);CHKERRQ(ierr);
374ea844a1aSMatthew Knepley     ierr = PetscSNPrintf(name, 64, "Field_%D", f);CHKERRQ(ierr);
375ea844a1aSMatthew Knepley     ierr = PetscStrallocpy(name, (char **) &s->fieldNames[f]);CHKERRQ(ierr);
376ea844a1aSMatthew Knepley   }
377ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
378ea844a1aSMatthew Knepley }
379ea844a1aSMatthew Knepley 
380ea844a1aSMatthew Knepley /*@C
381ea844a1aSMatthew Knepley   PetscSectionGetFieldName - Returns the name of a field in the PetscSection
382ea844a1aSMatthew Knepley 
383ea844a1aSMatthew Knepley   Not Collective
384ea844a1aSMatthew Knepley 
385ea844a1aSMatthew Knepley   Input Parameters:
386ea844a1aSMatthew Knepley + s     - the PetscSection
387ea844a1aSMatthew Knepley - field - the field number
388ea844a1aSMatthew Knepley 
389ea844a1aSMatthew Knepley   Output Parameter:
390ea844a1aSMatthew Knepley . fieldName - the field name
391ea844a1aSMatthew Knepley 
392ea844a1aSMatthew Knepley   Level: intermediate
393ea844a1aSMatthew Knepley 
394ea844a1aSMatthew Knepley .seealso: PetscSectionSetFieldName()
395ea844a1aSMatthew Knepley @*/
396ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetFieldName(PetscSection s, PetscInt field, const char *fieldName[])
397ea844a1aSMatthew Knepley {
398ea844a1aSMatthew Knepley   PetscFunctionBegin;
399ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
400ea844a1aSMatthew Knepley   PetscValidPointer(fieldName,3);
401ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
402ea844a1aSMatthew Knepley   *fieldName = s->fieldNames[field];
403ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
404ea844a1aSMatthew Knepley }
405ea844a1aSMatthew Knepley 
406ea844a1aSMatthew Knepley /*@C
407ea844a1aSMatthew Knepley   PetscSectionSetFieldName - Sets the name of a field in the PetscSection
408ea844a1aSMatthew Knepley 
409ea844a1aSMatthew Knepley   Not Collective
410ea844a1aSMatthew Knepley 
411ea844a1aSMatthew Knepley   Input Parameters:
412ea844a1aSMatthew Knepley + s     - the PetscSection
413ea844a1aSMatthew Knepley . field - the field number
414ea844a1aSMatthew Knepley - fieldName - the field name
415ea844a1aSMatthew Knepley 
416ea844a1aSMatthew Knepley   Level: intermediate
417ea844a1aSMatthew Knepley 
418ea844a1aSMatthew Knepley .seealso: PetscSectionGetFieldName()
419ea844a1aSMatthew Knepley @*/
420ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetFieldName(PetscSection s, PetscInt field, const char fieldName[])
421ea844a1aSMatthew Knepley {
422ea844a1aSMatthew Knepley   PetscErrorCode ierr;
423ea844a1aSMatthew Knepley 
424ea844a1aSMatthew Knepley   PetscFunctionBegin;
425ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
426ea844a1aSMatthew Knepley   if (fieldName) PetscValidCharPointer(fieldName,3);
427ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
428ea844a1aSMatthew Knepley   ierr = PetscFree(s->fieldNames[field]);CHKERRQ(ierr);
429ea844a1aSMatthew Knepley   ierr = PetscStrallocpy(fieldName, (char**) &s->fieldNames[field]);CHKERRQ(ierr);
430ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
431ea844a1aSMatthew Knepley }
432ea844a1aSMatthew Knepley 
433ea844a1aSMatthew Knepley /*@
434ea844a1aSMatthew Knepley   PetscSectionGetFieldComponents - Returns the number of field components for the given field.
435ea844a1aSMatthew Knepley 
436ea844a1aSMatthew Knepley   Not collective
437ea844a1aSMatthew Knepley 
438ea844a1aSMatthew Knepley   Input Parameters:
439ea844a1aSMatthew Knepley + s - the PetscSection
440ea844a1aSMatthew Knepley - field - the field number
441ea844a1aSMatthew Knepley 
442ea844a1aSMatthew Knepley   Output Parameter:
443ea844a1aSMatthew Knepley . numComp - the number of field components
444ea844a1aSMatthew Knepley 
445ea844a1aSMatthew Knepley   Level: intermediate
446ea844a1aSMatthew Knepley 
447ea844a1aSMatthew Knepley .seealso: PetscSectionSetFieldComponents(), PetscSectionGetNumFields()
448ea844a1aSMatthew Knepley @*/
449ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetFieldComponents(PetscSection s, PetscInt field, PetscInt *numComp)
450ea844a1aSMatthew Knepley {
451ea844a1aSMatthew Knepley   PetscFunctionBegin;
452ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
453ea844a1aSMatthew Knepley   PetscValidPointer(numComp, 3);
454ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
455ea844a1aSMatthew Knepley   *numComp = s->numFieldComponents[field];
456ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
457ea844a1aSMatthew Knepley }
458ea844a1aSMatthew Knepley 
459ea844a1aSMatthew Knepley /*@
460ea844a1aSMatthew Knepley   PetscSectionSetFieldComponents - Sets the number of field components for the given field.
461ea844a1aSMatthew Knepley 
462ea844a1aSMatthew Knepley   Not collective
463ea844a1aSMatthew Knepley 
464ea844a1aSMatthew Knepley   Input Parameters:
465ea844a1aSMatthew Knepley + s - the PetscSection
466ea844a1aSMatthew Knepley . field - the field number
467ea844a1aSMatthew Knepley - numComp - the number of field components
468ea844a1aSMatthew Knepley 
469ea844a1aSMatthew Knepley   Level: intermediate
470ea844a1aSMatthew Knepley 
471ea844a1aSMatthew Knepley .seealso: PetscSectionGetNumFieldComponents(), PetscSectionGetNumFields()
472ea844a1aSMatthew Knepley @*/
473ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetFieldComponents(PetscSection s, PetscInt field, PetscInt numComp)
474ea844a1aSMatthew Knepley {
475ea844a1aSMatthew Knepley   PetscFunctionBegin;
476ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
477ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
478ea844a1aSMatthew Knepley   s->numFieldComponents[field] = numComp;
479ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
480ea844a1aSMatthew Knepley }
481ea844a1aSMatthew Knepley 
482ea844a1aSMatthew Knepley static PetscErrorCode PetscSectionCheckConstraints_Static(PetscSection s)
483ea844a1aSMatthew Knepley {
484ea844a1aSMatthew Knepley   PetscErrorCode ierr;
485ea844a1aSMatthew Knepley 
486ea844a1aSMatthew Knepley   PetscFunctionBegin;
487ea844a1aSMatthew Knepley   if (!s->bc) {
488ea844a1aSMatthew Knepley     ierr = PetscSectionCreate(PETSC_COMM_SELF, &s->bc);CHKERRQ(ierr);
489ea844a1aSMatthew Knepley     ierr = PetscSectionSetChart(s->bc, s->pStart, s->pEnd);CHKERRQ(ierr);
490ea844a1aSMatthew Knepley   }
491ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
492ea844a1aSMatthew Knepley }
493ea844a1aSMatthew Knepley 
494ea844a1aSMatthew Knepley /*@
495ea844a1aSMatthew Knepley   PetscSectionGetChart - Returns the range [pStart, pEnd) in which points lie.
496ea844a1aSMatthew Knepley 
497ea844a1aSMatthew Knepley   Not collective
498ea844a1aSMatthew Knepley 
499ea844a1aSMatthew Knepley   Input Parameter:
500ea844a1aSMatthew Knepley . s - the PetscSection
501ea844a1aSMatthew Knepley 
502ea844a1aSMatthew Knepley   Output Parameters:
503ea844a1aSMatthew Knepley + pStart - the first point
504ea844a1aSMatthew Knepley - pEnd - one past the last point
505ea844a1aSMatthew Knepley 
506ea844a1aSMatthew Knepley   Level: intermediate
507ea844a1aSMatthew Knepley 
508ea844a1aSMatthew Knepley .seealso: PetscSectionSetChart(), PetscSectionCreate()
509ea844a1aSMatthew Knepley @*/
510ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetChart(PetscSection s, PetscInt *pStart, PetscInt *pEnd)
511ea844a1aSMatthew Knepley {
512ea844a1aSMatthew Knepley   PetscFunctionBegin;
513ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
514ea844a1aSMatthew Knepley   if (pStart) *pStart = s->pStart;
515ea844a1aSMatthew Knepley   if (pEnd)   *pEnd   = s->pEnd;
516ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
517ea844a1aSMatthew Knepley }
518ea844a1aSMatthew Knepley 
519ea844a1aSMatthew Knepley /*@
520ea844a1aSMatthew Knepley   PetscSectionSetChart - Sets the range [pStart, pEnd) in which points lie.
521ea844a1aSMatthew Knepley 
522ea844a1aSMatthew Knepley   Not collective
523ea844a1aSMatthew Knepley 
524ea844a1aSMatthew Knepley   Input Parameters:
525ea844a1aSMatthew Knepley + s - the PetscSection
526ea844a1aSMatthew Knepley . pStart - the first point
527ea844a1aSMatthew Knepley - pEnd - one past the last point
528ea844a1aSMatthew Knepley 
529ea844a1aSMatthew Knepley   Level: intermediate
530ea844a1aSMatthew Knepley 
531ea844a1aSMatthew Knepley .seealso: PetscSectionGetChart(), PetscSectionCreate()
532ea844a1aSMatthew Knepley @*/
533ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetChart(PetscSection s, PetscInt pStart, PetscInt pEnd)
534ea844a1aSMatthew Knepley {
535ea844a1aSMatthew Knepley   PetscInt       f;
536ea844a1aSMatthew Knepley   PetscErrorCode ierr;
537ea844a1aSMatthew Knepley 
538ea844a1aSMatthew Knepley   PetscFunctionBegin;
539ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
540ea844a1aSMatthew Knepley   /* Cannot Reset() because it destroys field information */
541ea844a1aSMatthew Knepley   s->setup = PETSC_FALSE;
542ea844a1aSMatthew Knepley   ierr = PetscSectionDestroy(&s->bc);CHKERRQ(ierr);
543ea844a1aSMatthew Knepley   ierr = PetscFree(s->bcIndices);CHKERRQ(ierr);
544ea844a1aSMatthew Knepley   ierr = PetscFree2(s->atlasDof, s->atlasOff);CHKERRQ(ierr);
545ea844a1aSMatthew Knepley 
546ea844a1aSMatthew Knepley   s->pStart = pStart;
547ea844a1aSMatthew Knepley   s->pEnd   = pEnd;
548ea844a1aSMatthew Knepley   ierr = PetscMalloc2((pEnd - pStart), &s->atlasDof, (pEnd - pStart), &s->atlasOff);CHKERRQ(ierr);
549ea844a1aSMatthew Knepley   ierr = PetscArrayzero(s->atlasDof, pEnd - pStart);CHKERRQ(ierr);
550ea844a1aSMatthew Knepley   for (f = 0; f < s->numFields; ++f) {
551ea844a1aSMatthew Knepley     ierr = PetscSectionSetChart(s->field[f], pStart, pEnd);CHKERRQ(ierr);
552ea844a1aSMatthew Knepley   }
553ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
554ea844a1aSMatthew Knepley }
555ea844a1aSMatthew Knepley 
556ea844a1aSMatthew Knepley /*@
557ea844a1aSMatthew Knepley   PetscSectionGetPermutation - Returns the permutation of [0, pEnd-pStart) or NULL
558ea844a1aSMatthew Knepley 
559ea844a1aSMatthew Knepley   Not collective
560ea844a1aSMatthew Knepley 
561ea844a1aSMatthew Knepley   Input Parameter:
562ea844a1aSMatthew Knepley . s - the PetscSection
563ea844a1aSMatthew Knepley 
564ea844a1aSMatthew Knepley   Output Parameters:
565ea844a1aSMatthew Knepley . perm - The permutation as an IS
566ea844a1aSMatthew Knepley 
567ea844a1aSMatthew Knepley   Level: intermediate
568ea844a1aSMatthew Knepley 
569ea844a1aSMatthew Knepley .seealso: PetscSectionSetPermutation(), PetscSectionCreate()
570ea844a1aSMatthew Knepley @*/
571ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetPermutation(PetscSection s, IS *perm)
572ea844a1aSMatthew Knepley {
573ea844a1aSMatthew Knepley   PetscFunctionBegin;
574ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
575ea844a1aSMatthew Knepley   if (perm) {PetscValidPointer(perm, 2); *perm = s->perm;}
576ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
577ea844a1aSMatthew Knepley }
578ea844a1aSMatthew Knepley 
579ea844a1aSMatthew Knepley /*@
580ea844a1aSMatthew Knepley   PetscSectionSetPermutation - Sets the permutation for [0, pEnd-pStart)
581ea844a1aSMatthew Knepley 
582ea844a1aSMatthew Knepley   Not collective
583ea844a1aSMatthew Knepley 
584ea844a1aSMatthew Knepley   Input Parameters:
585ea844a1aSMatthew Knepley + s - the PetscSection
586ea844a1aSMatthew Knepley - perm - the permutation of points
587ea844a1aSMatthew Knepley 
588ea844a1aSMatthew Knepley   Level: intermediate
589ea844a1aSMatthew Knepley 
590ea844a1aSMatthew Knepley .seealso: PetscSectionGetPermutation(), PetscSectionCreate()
591ea844a1aSMatthew Knepley @*/
592ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetPermutation(PetscSection s, IS perm)
593ea844a1aSMatthew Knepley {
594ea844a1aSMatthew Knepley   PetscErrorCode ierr;
595ea844a1aSMatthew Knepley 
596ea844a1aSMatthew Knepley   PetscFunctionBegin;
597ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
598ea844a1aSMatthew Knepley   if (perm) PetscValidHeaderSpecific(perm, IS_CLASSID, 2);
599ea844a1aSMatthew Knepley   if (s->setup) SETERRQ(PetscObjectComm((PetscObject) s), PETSC_ERR_ARG_WRONGSTATE, "Cannot set a permutation after the section is setup");
600ea844a1aSMatthew Knepley   if (s->perm != perm) {
601ea844a1aSMatthew Knepley     ierr = ISDestroy(&s->perm);CHKERRQ(ierr);
602ea844a1aSMatthew Knepley     if (perm) {
603ea844a1aSMatthew Knepley       s->perm = perm;
604ea844a1aSMatthew Knepley       ierr = PetscObjectReference((PetscObject) s->perm);CHKERRQ(ierr);
605ea844a1aSMatthew Knepley     }
606ea844a1aSMatthew Knepley   }
607ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
608ea844a1aSMatthew Knepley }
609ea844a1aSMatthew Knepley 
610ea844a1aSMatthew Knepley /*@
611ea844a1aSMatthew Knepley   PetscSectionGetPointMajor - Returns the flag for dof ordering, true if it is point major, otherwise field major
612ea844a1aSMatthew Knepley 
613ea844a1aSMatthew Knepley   Not collective
614ea844a1aSMatthew Knepley 
615ea844a1aSMatthew Knepley   Input Parameter:
616ea844a1aSMatthew Knepley . s - the PetscSection
617ea844a1aSMatthew Knepley 
618ea844a1aSMatthew Knepley   Output Parameter:
619ea844a1aSMatthew Knepley . pm - the flag for point major ordering
620ea844a1aSMatthew Knepley 
621ea844a1aSMatthew Knepley   Level: intermediate
622ea844a1aSMatthew Knepley 
623ea844a1aSMatthew Knepley .seealso: PetscSectionSetPointMajor()
624ea844a1aSMatthew Knepley @*/
625ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetPointMajor(PetscSection s, PetscBool *pm)
626ea844a1aSMatthew Knepley {
627ea844a1aSMatthew Knepley   PetscFunctionBegin;
628ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
629ea844a1aSMatthew Knepley   PetscValidPointer(pm,2);
630ea844a1aSMatthew Knepley   *pm = s->pointMajor;
631ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
632ea844a1aSMatthew Knepley }
633ea844a1aSMatthew Knepley 
634ea844a1aSMatthew Knepley /*@
635ea844a1aSMatthew Knepley   PetscSectionSetPointMajor - Sets the flag for dof ordering, true if it is point major, otherwise field major
636ea844a1aSMatthew Knepley 
637ea844a1aSMatthew Knepley   Not collective
638ea844a1aSMatthew Knepley 
639ea844a1aSMatthew Knepley   Input Parameters:
640ea844a1aSMatthew Knepley + s  - the PetscSection
641ea844a1aSMatthew Knepley - pm - the flag for point major ordering
642ea844a1aSMatthew Knepley 
643ea844a1aSMatthew Knepley   Not collective
644ea844a1aSMatthew Knepley 
645ea844a1aSMatthew Knepley   Level: intermediate
646ea844a1aSMatthew Knepley 
647ea844a1aSMatthew Knepley .seealso: PetscSectionGetPointMajor()
648ea844a1aSMatthew Knepley @*/
649ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetPointMajor(PetscSection s, PetscBool pm)
650ea844a1aSMatthew Knepley {
651ea844a1aSMatthew Knepley   PetscFunctionBegin;
652ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
653ea844a1aSMatthew Knepley   if (s->setup) SETERRQ(PetscObjectComm((PetscObject) s), PETSC_ERR_ARG_WRONGSTATE, "Cannot set the dof ordering after the section is setup");
654ea844a1aSMatthew Knepley   s->pointMajor = pm;
655ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
656ea844a1aSMatthew Knepley }
657ea844a1aSMatthew Knepley 
658ea844a1aSMatthew Knepley /*@
659ea844a1aSMatthew Knepley   PetscSectionGetDof - Return the number of degrees of freedom associated with a given point.
660ea844a1aSMatthew Knepley 
661ea844a1aSMatthew Knepley   Not collective
662ea844a1aSMatthew Knepley 
663ea844a1aSMatthew Knepley   Input Parameters:
664ea844a1aSMatthew Knepley + s - the PetscSection
665ea844a1aSMatthew Knepley - point - the point
666ea844a1aSMatthew Knepley 
667ea844a1aSMatthew Knepley   Output Parameter:
668ea844a1aSMatthew Knepley . numDof - the number of dof
669ea844a1aSMatthew Knepley 
670ea844a1aSMatthew Knepley   Level: intermediate
671ea844a1aSMatthew Knepley 
672ea844a1aSMatthew Knepley .seealso: PetscSectionSetDof(), PetscSectionCreate()
673ea844a1aSMatthew Knepley @*/
674ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetDof(PetscSection s, PetscInt point, PetscInt *numDof)
675ea844a1aSMatthew Knepley {
676ea844a1aSMatthew Knepley   PetscFunctionBegin;
677ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
678ea844a1aSMatthew Knepley   PetscValidPointer(numDof, 3);
679ea844a1aSMatthew Knepley #if defined(PETSC_USE_DEBUG)
680ea844a1aSMatthew Knepley   if ((point < s->pStart) || (point >= s->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %D should be in [%D, %D)", point, s->pStart, s->pEnd);
681ea844a1aSMatthew Knepley #endif
682ea844a1aSMatthew Knepley   *numDof = s->atlasDof[point - s->pStart];
683ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
684ea844a1aSMatthew Knepley }
685ea844a1aSMatthew Knepley 
686ea844a1aSMatthew Knepley /*@
687ea844a1aSMatthew Knepley   PetscSectionSetDof - Sets the number of degrees of freedom associated with a given point.
688ea844a1aSMatthew Knepley 
689ea844a1aSMatthew Knepley   Not collective
690ea844a1aSMatthew Knepley 
691ea844a1aSMatthew Knepley   Input Parameters:
692ea844a1aSMatthew Knepley + s - the PetscSection
693ea844a1aSMatthew Knepley . point - the point
694ea844a1aSMatthew Knepley - numDof - the number of dof
695ea844a1aSMatthew Knepley 
696ea844a1aSMatthew Knepley   Level: intermediate
697ea844a1aSMatthew Knepley 
698ea844a1aSMatthew Knepley .seealso: PetscSectionGetDof(), PetscSectionAddDof(), PetscSectionCreate()
699ea844a1aSMatthew Knepley @*/
700ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetDof(PetscSection s, PetscInt point, PetscInt numDof)
701ea844a1aSMatthew Knepley {
702ea844a1aSMatthew Knepley   PetscFunctionBegin;
703ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
704ea844a1aSMatthew Knepley   if ((point < s->pStart) || (point >= s->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %D should be in [%D, %D)", point, s->pStart, s->pEnd);
705ea844a1aSMatthew Knepley   s->atlasDof[point - s->pStart] = numDof;
706ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
707ea844a1aSMatthew Knepley }
708ea844a1aSMatthew Knepley 
709ea844a1aSMatthew Knepley /*@
710ea844a1aSMatthew Knepley   PetscSectionAddDof - Adds to the number of degrees of freedom associated with a given point.
711ea844a1aSMatthew Knepley 
712ea844a1aSMatthew Knepley   Not collective
713ea844a1aSMatthew Knepley 
714ea844a1aSMatthew Knepley   Input Parameters:
715ea844a1aSMatthew Knepley + s - the PetscSection
716ea844a1aSMatthew Knepley . point - the point
717ea844a1aSMatthew Knepley - numDof - the number of additional dof
718ea844a1aSMatthew Knepley 
719ea844a1aSMatthew Knepley   Level: intermediate
720ea844a1aSMatthew Knepley 
721ea844a1aSMatthew Knepley .seealso: PetscSectionGetDof(), PetscSectionSetDof(), PetscSectionCreate()
722ea844a1aSMatthew Knepley @*/
723ea844a1aSMatthew Knepley PetscErrorCode PetscSectionAddDof(PetscSection s, PetscInt point, PetscInt numDof)
724ea844a1aSMatthew Knepley {
725ea844a1aSMatthew Knepley   PetscFunctionBegin;
726ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
727ea844a1aSMatthew Knepley #if defined(PETSC_USE_DEBUG)
728ea844a1aSMatthew Knepley   if ((point < s->pStart) || (point >= s->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %D should be in [%D, %D)", point, s->pStart, s->pEnd);
729ea844a1aSMatthew Knepley #endif
730ea844a1aSMatthew Knepley   s->atlasDof[point - s->pStart] += numDof;
731ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
732ea844a1aSMatthew Knepley }
733ea844a1aSMatthew Knepley 
734ea844a1aSMatthew Knepley /*@
735ea844a1aSMatthew Knepley   PetscSectionGetFieldDof - Return the number of degrees of freedom associated with a field on a given point.
736ea844a1aSMatthew Knepley 
737ea844a1aSMatthew Knepley   Not collective
738ea844a1aSMatthew Knepley 
739ea844a1aSMatthew Knepley   Input Parameters:
740ea844a1aSMatthew Knepley + s - the PetscSection
741ea844a1aSMatthew Knepley . point - the point
742ea844a1aSMatthew Knepley - field - the field
743ea844a1aSMatthew Knepley 
744ea844a1aSMatthew Knepley   Output Parameter:
745ea844a1aSMatthew Knepley . numDof - the number of dof
746ea844a1aSMatthew Knepley 
747ea844a1aSMatthew Knepley   Level: intermediate
748ea844a1aSMatthew Knepley 
749ea844a1aSMatthew Knepley .seealso: PetscSectionSetFieldDof(), PetscSectionCreate()
750ea844a1aSMatthew Knepley @*/
751ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt *numDof)
752ea844a1aSMatthew Knepley {
753ea844a1aSMatthew Knepley   PetscErrorCode ierr;
754ea844a1aSMatthew Knepley 
755ea844a1aSMatthew Knepley   PetscFunctionBegin;
756ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
757ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
758ea844a1aSMatthew Knepley   ierr = PetscSectionGetDof(s->field[field], point, numDof);CHKERRQ(ierr);
759ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
760ea844a1aSMatthew Knepley }
761ea844a1aSMatthew Knepley 
762ea844a1aSMatthew Knepley /*@
763ea844a1aSMatthew Knepley   PetscSectionSetFieldDof - Sets the number of degrees of freedom associated with a field on a given point.
764ea844a1aSMatthew Knepley 
765ea844a1aSMatthew Knepley   Not collective
766ea844a1aSMatthew Knepley 
767ea844a1aSMatthew Knepley   Input Parameters:
768ea844a1aSMatthew Knepley + s - the PetscSection
769ea844a1aSMatthew Knepley . point - the point
770ea844a1aSMatthew Knepley . field - the field
771ea844a1aSMatthew Knepley - numDof - the number of dof
772ea844a1aSMatthew Knepley 
773ea844a1aSMatthew Knepley   Level: intermediate
774ea844a1aSMatthew Knepley 
775ea844a1aSMatthew Knepley .seealso: PetscSectionGetFieldDof(), PetscSectionCreate()
776ea844a1aSMatthew Knepley @*/
777ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof)
778ea844a1aSMatthew Knepley {
779ea844a1aSMatthew Knepley   PetscErrorCode ierr;
780ea844a1aSMatthew Knepley 
781ea844a1aSMatthew Knepley   PetscFunctionBegin;
782ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
783ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
784ea844a1aSMatthew Knepley   ierr = PetscSectionSetDof(s->field[field], point, numDof);CHKERRQ(ierr);
785ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
786ea844a1aSMatthew Knepley }
787ea844a1aSMatthew Knepley 
788ea844a1aSMatthew Knepley /*@
789ea844a1aSMatthew Knepley   PetscSectionAddFieldDof - Adds a number of degrees of freedom associated with a field on a given point.
790ea844a1aSMatthew Knepley 
791ea844a1aSMatthew Knepley   Not collective
792ea844a1aSMatthew Knepley 
793ea844a1aSMatthew Knepley   Input Parameters:
794ea844a1aSMatthew Knepley + s - the PetscSection
795ea844a1aSMatthew Knepley . point - the point
796ea844a1aSMatthew Knepley . field - the field
797ea844a1aSMatthew Knepley - numDof - the number of dof
798ea844a1aSMatthew Knepley 
799ea844a1aSMatthew Knepley   Level: intermediate
800ea844a1aSMatthew Knepley 
801ea844a1aSMatthew Knepley .seealso: PetscSectionSetFieldDof(), PetscSectionGetFieldDof(), PetscSectionCreate()
802ea844a1aSMatthew Knepley @*/
803ea844a1aSMatthew Knepley PetscErrorCode PetscSectionAddFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof)
804ea844a1aSMatthew Knepley {
805ea844a1aSMatthew Knepley   PetscErrorCode ierr;
806ea844a1aSMatthew Knepley 
807ea844a1aSMatthew Knepley   PetscFunctionBegin;
808ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
809ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
810ea844a1aSMatthew Knepley   ierr = PetscSectionAddDof(s->field[field], point, numDof);CHKERRQ(ierr);
811ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
812ea844a1aSMatthew Knepley }
813ea844a1aSMatthew Knepley 
814ea844a1aSMatthew Knepley /*@
815ea844a1aSMatthew Knepley   PetscSectionGetConstraintDof - Return the number of constrained degrees of freedom associated with a given point.
816ea844a1aSMatthew Knepley 
817ea844a1aSMatthew Knepley   Not collective
818ea844a1aSMatthew Knepley 
819ea844a1aSMatthew Knepley   Input Parameters:
820ea844a1aSMatthew Knepley + s - the PetscSection
821ea844a1aSMatthew Knepley - point - the point
822ea844a1aSMatthew Knepley 
823ea844a1aSMatthew Knepley   Output Parameter:
824ea844a1aSMatthew Knepley . numDof - the number of dof which are fixed by constraints
825ea844a1aSMatthew Knepley 
826ea844a1aSMatthew Knepley   Level: intermediate
827ea844a1aSMatthew Knepley 
828ea844a1aSMatthew Knepley .seealso: PetscSectionGetDof(), PetscSectionSetConstraintDof(), PetscSectionCreate()
829ea844a1aSMatthew Knepley @*/
830ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetConstraintDof(PetscSection s, PetscInt point, PetscInt *numDof)
831ea844a1aSMatthew Knepley {
832ea844a1aSMatthew Knepley   PetscErrorCode ierr;
833ea844a1aSMatthew Knepley 
834ea844a1aSMatthew Knepley   PetscFunctionBegin;
835ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
836ea844a1aSMatthew Knepley   PetscValidPointer(numDof, 3);
837ea844a1aSMatthew Knepley   if (s->bc) {
838ea844a1aSMatthew Knepley     ierr = PetscSectionGetDof(s->bc, point, numDof);CHKERRQ(ierr);
839ea844a1aSMatthew Knepley   } else *numDof = 0;
840ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
841ea844a1aSMatthew Knepley }
842ea844a1aSMatthew Knepley 
843ea844a1aSMatthew Knepley /*@
844ea844a1aSMatthew Knepley   PetscSectionSetConstraintDof - Set the number of constrained degrees of freedom associated with a given point.
845ea844a1aSMatthew Knepley 
846ea844a1aSMatthew Knepley   Not collective
847ea844a1aSMatthew Knepley 
848ea844a1aSMatthew Knepley   Input Parameters:
849ea844a1aSMatthew Knepley + s - the PetscSection
850ea844a1aSMatthew Knepley . point - the point
851ea844a1aSMatthew Knepley - numDof - the number of dof which are fixed by constraints
852ea844a1aSMatthew Knepley 
853ea844a1aSMatthew Knepley   Level: intermediate
854ea844a1aSMatthew Knepley 
855ea844a1aSMatthew Knepley .seealso: PetscSectionSetDof(), PetscSectionGetConstraintDof(), PetscSectionCreate()
856ea844a1aSMatthew Knepley @*/
857ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetConstraintDof(PetscSection s, PetscInt point, PetscInt numDof)
858ea844a1aSMatthew Knepley {
859ea844a1aSMatthew Knepley   PetscErrorCode ierr;
860ea844a1aSMatthew Knepley 
861ea844a1aSMatthew Knepley   PetscFunctionBegin;
862ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
863ea844a1aSMatthew Knepley   if (numDof) {
864ea844a1aSMatthew Knepley     ierr = PetscSectionCheckConstraints_Static(s);CHKERRQ(ierr);
865ea844a1aSMatthew Knepley     ierr = PetscSectionSetDof(s->bc, point, numDof);CHKERRQ(ierr);
866ea844a1aSMatthew Knepley   }
867ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
868ea844a1aSMatthew Knepley }
869ea844a1aSMatthew Knepley 
870ea844a1aSMatthew Knepley /*@
871ea844a1aSMatthew Knepley   PetscSectionAddConstraintDof - Increment the number of constrained degrees of freedom associated with a given point.
872ea844a1aSMatthew Knepley 
873ea844a1aSMatthew Knepley   Not collective
874ea844a1aSMatthew Knepley 
875ea844a1aSMatthew Knepley   Input Parameters:
876ea844a1aSMatthew Knepley + s - the PetscSection
877ea844a1aSMatthew Knepley . point - the point
878ea844a1aSMatthew Knepley - numDof - the number of additional dof which are fixed by constraints
879ea844a1aSMatthew Knepley 
880ea844a1aSMatthew Knepley   Level: intermediate
881ea844a1aSMatthew Knepley 
882ea844a1aSMatthew Knepley .seealso: PetscSectionAddDof(), PetscSectionGetConstraintDof(), PetscSectionCreate()
883ea844a1aSMatthew Knepley @*/
884ea844a1aSMatthew Knepley PetscErrorCode PetscSectionAddConstraintDof(PetscSection s, PetscInt point, PetscInt numDof)
885ea844a1aSMatthew Knepley {
886ea844a1aSMatthew Knepley   PetscErrorCode ierr;
887ea844a1aSMatthew Knepley 
888ea844a1aSMatthew Knepley   PetscFunctionBegin;
889ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
890ea844a1aSMatthew Knepley   if (numDof) {
891ea844a1aSMatthew Knepley     ierr = PetscSectionCheckConstraints_Static(s);CHKERRQ(ierr);
892ea844a1aSMatthew Knepley     ierr = PetscSectionAddDof(s->bc, point, numDof);CHKERRQ(ierr);
893ea844a1aSMatthew Knepley   }
894ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
895ea844a1aSMatthew Knepley }
896ea844a1aSMatthew Knepley 
897ea844a1aSMatthew Knepley /*@
898ea844a1aSMatthew Knepley   PetscSectionGetFieldConstraintDof - Return the number of constrained degrees of freedom associated with a given field on a point.
899ea844a1aSMatthew Knepley 
900ea844a1aSMatthew Knepley   Not collective
901ea844a1aSMatthew Knepley 
902ea844a1aSMatthew Knepley   Input Parameters:
903ea844a1aSMatthew Knepley + s - the PetscSection
904ea844a1aSMatthew Knepley . point - the point
905ea844a1aSMatthew Knepley - field - the field
906ea844a1aSMatthew Knepley 
907ea844a1aSMatthew Knepley   Output Parameter:
908ea844a1aSMatthew Knepley . numDof - the number of dof which are fixed by constraints
909ea844a1aSMatthew Knepley 
910ea844a1aSMatthew Knepley   Level: intermediate
911ea844a1aSMatthew Knepley 
912ea844a1aSMatthew Knepley .seealso: PetscSectionGetDof(), PetscSectionSetFieldConstraintDof(), PetscSectionCreate()
913ea844a1aSMatthew Knepley @*/
914ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetFieldConstraintDof(PetscSection s, PetscInt point, PetscInt field, PetscInt *numDof)
915ea844a1aSMatthew Knepley {
916ea844a1aSMatthew Knepley   PetscErrorCode ierr;
917ea844a1aSMatthew Knepley 
918ea844a1aSMatthew Knepley   PetscFunctionBegin;
919ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
920ea844a1aSMatthew Knepley   PetscValidPointer(numDof, 4);
921ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
922ea844a1aSMatthew Knepley   ierr = PetscSectionGetConstraintDof(s->field[field], point, numDof);CHKERRQ(ierr);
923ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
924ea844a1aSMatthew Knepley }
925ea844a1aSMatthew Knepley 
926ea844a1aSMatthew Knepley /*@
927ea844a1aSMatthew Knepley   PetscSectionSetFieldConstraintDof - Set the number of constrained degrees of freedom associated with a given field on a point.
928ea844a1aSMatthew Knepley 
929ea844a1aSMatthew Knepley   Not collective
930ea844a1aSMatthew Knepley 
931ea844a1aSMatthew Knepley   Input Parameters:
932ea844a1aSMatthew Knepley + s - the PetscSection
933ea844a1aSMatthew Knepley . point - the point
934ea844a1aSMatthew Knepley . field - the field
935ea844a1aSMatthew Knepley - numDof - the number of dof which are fixed by constraints
936ea844a1aSMatthew Knepley 
937ea844a1aSMatthew Knepley   Level: intermediate
938ea844a1aSMatthew Knepley 
939ea844a1aSMatthew Knepley .seealso: PetscSectionSetDof(), PetscSectionGetFieldConstraintDof(), PetscSectionCreate()
940ea844a1aSMatthew Knepley @*/
941ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetFieldConstraintDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof)
942ea844a1aSMatthew Knepley {
943ea844a1aSMatthew Knepley   PetscErrorCode ierr;
944ea844a1aSMatthew Knepley 
945ea844a1aSMatthew Knepley   PetscFunctionBegin;
946ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
947ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
948ea844a1aSMatthew Knepley   ierr = PetscSectionSetConstraintDof(s->field[field], point, numDof);CHKERRQ(ierr);
949ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
950ea844a1aSMatthew Knepley }
951ea844a1aSMatthew Knepley 
952ea844a1aSMatthew Knepley /*@
953ea844a1aSMatthew Knepley   PetscSectionAddFieldConstraintDof - Increment the number of constrained degrees of freedom associated with a given field on a point.
954ea844a1aSMatthew Knepley 
955ea844a1aSMatthew Knepley   Not collective
956ea844a1aSMatthew Knepley 
957ea844a1aSMatthew Knepley   Input Parameters:
958ea844a1aSMatthew Knepley + s - the PetscSection
959ea844a1aSMatthew Knepley . point - the point
960ea844a1aSMatthew Knepley . field - the field
961ea844a1aSMatthew Knepley - numDof - the number of additional dof which are fixed by constraints
962ea844a1aSMatthew Knepley 
963ea844a1aSMatthew Knepley   Level: intermediate
964ea844a1aSMatthew Knepley 
965ea844a1aSMatthew Knepley .seealso: PetscSectionAddDof(), PetscSectionGetFieldConstraintDof(), PetscSectionCreate()
966ea844a1aSMatthew Knepley @*/
967ea844a1aSMatthew Knepley PetscErrorCode PetscSectionAddFieldConstraintDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof)
968ea844a1aSMatthew Knepley {
969ea844a1aSMatthew Knepley   PetscErrorCode ierr;
970ea844a1aSMatthew Knepley 
971ea844a1aSMatthew Knepley   PetscFunctionBegin;
972ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
973ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
974ea844a1aSMatthew Knepley   ierr = PetscSectionAddConstraintDof(s->field[field], point, numDof);CHKERRQ(ierr);
975ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
976ea844a1aSMatthew Knepley }
977ea844a1aSMatthew Knepley 
978ea844a1aSMatthew Knepley /*@
979ea844a1aSMatthew Knepley   PetscSectionSetUpBC - Setup the subsections describing boundary conditions.
980ea844a1aSMatthew Knepley 
981ea844a1aSMatthew Knepley   Not collective
982ea844a1aSMatthew Knepley 
983ea844a1aSMatthew Knepley   Input Parameter:
984ea844a1aSMatthew Knepley . s - the PetscSection
985ea844a1aSMatthew Knepley 
986ea844a1aSMatthew Knepley   Level: advanced
987ea844a1aSMatthew Knepley 
988ea844a1aSMatthew Knepley .seealso: PetscSectionSetUp(), PetscSectionCreate()
989ea844a1aSMatthew Knepley @*/
990ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetUpBC(PetscSection s)
991ea844a1aSMatthew Knepley {
992ea844a1aSMatthew Knepley   PetscErrorCode ierr;
993ea844a1aSMatthew Knepley 
994ea844a1aSMatthew Knepley   PetscFunctionBegin;
995ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
996ea844a1aSMatthew Knepley   if (s->bc) {
997ea844a1aSMatthew Knepley     const PetscInt last = (s->bc->pEnd-s->bc->pStart) - 1;
998ea844a1aSMatthew Knepley 
999ea844a1aSMatthew Knepley     ierr = PetscSectionSetUp(s->bc);CHKERRQ(ierr);
100028d58a37SPierre Jolivet     ierr = PetscMalloc1(last >= 0 ? s->bc->atlasOff[last] + s->bc->atlasDof[last] : 0, &s->bcIndices);CHKERRQ(ierr);
1001ea844a1aSMatthew Knepley   }
1002ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1003ea844a1aSMatthew Knepley }
1004ea844a1aSMatthew Knepley 
1005ea844a1aSMatthew Knepley /*@
1006ea844a1aSMatthew Knepley   PetscSectionSetUp - Calculate offsets based upon the number of degrees of freedom for each point.
1007ea844a1aSMatthew Knepley 
1008ea844a1aSMatthew Knepley   Not collective
1009ea844a1aSMatthew Knepley 
1010ea844a1aSMatthew Knepley   Input Parameter:
1011ea844a1aSMatthew Knepley . s - the PetscSection
1012ea844a1aSMatthew Knepley 
1013ea844a1aSMatthew Knepley   Level: intermediate
1014ea844a1aSMatthew Knepley 
1015ea844a1aSMatthew Knepley .seealso: PetscSectionCreate()
1016ea844a1aSMatthew Knepley @*/
1017ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetUp(PetscSection s)
1018ea844a1aSMatthew Knepley {
1019ea844a1aSMatthew Knepley   const PetscInt *pind   = NULL;
1020ea844a1aSMatthew Knepley   PetscInt        offset = 0, foff, p, f;
1021ea844a1aSMatthew Knepley   PetscErrorCode  ierr;
1022ea844a1aSMatthew Knepley 
1023ea844a1aSMatthew Knepley   PetscFunctionBegin;
1024ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1025ea844a1aSMatthew Knepley   if (s->setup) PetscFunctionReturn(0);
1026ea844a1aSMatthew Knepley   s->setup = PETSC_TRUE;
1027ea844a1aSMatthew Knepley   /* Set offsets and field offsets for all points */
1028ea844a1aSMatthew Knepley   /*   Assume that all fields have the same chart */
1029ea844a1aSMatthew Knepley   if (s->perm) {ierr = ISGetIndices(s->perm, &pind);CHKERRQ(ierr);}
1030ea844a1aSMatthew Knepley   if (s->pointMajor) {
1031ea844a1aSMatthew Knepley     for (p = 0; p < s->pEnd - s->pStart; ++p) {
1032ea844a1aSMatthew Knepley       const PetscInt q = pind ? pind[p] : p;
1033ea844a1aSMatthew Knepley 
1034ea844a1aSMatthew Knepley       /* Set point offset */
1035ea844a1aSMatthew Knepley       s->atlasOff[q] = offset;
1036ea844a1aSMatthew Knepley       offset        += s->atlasDof[q];
1037ea844a1aSMatthew Knepley       s->maxDof      = PetscMax(s->maxDof, s->atlasDof[q]);
1038ea844a1aSMatthew Knepley       /* Set field offset */
1039ea844a1aSMatthew Knepley       for (f = 0, foff = s->atlasOff[q]; f < s->numFields; ++f) {
1040ea844a1aSMatthew Knepley         PetscSection sf = s->field[f];
1041ea844a1aSMatthew Knepley 
1042ea844a1aSMatthew Knepley         sf->atlasOff[q] = foff;
1043ea844a1aSMatthew Knepley         foff += sf->atlasDof[q];
1044ea844a1aSMatthew Knepley       }
1045ea844a1aSMatthew Knepley     }
1046ea844a1aSMatthew Knepley   } else {
1047ea844a1aSMatthew Knepley     /* Set field offsets for all points */
1048ea844a1aSMatthew Knepley     for (f = 0; f < s->numFields; ++f) {
1049ea844a1aSMatthew Knepley       PetscSection sf = s->field[f];
1050ea844a1aSMatthew Knepley 
1051ea844a1aSMatthew Knepley       for (p = 0; p < s->pEnd - s->pStart; ++p) {
1052ea844a1aSMatthew Knepley         const PetscInt q = pind ? pind[p] : p;
1053ea844a1aSMatthew Knepley 
1054ea844a1aSMatthew Knepley         sf->atlasOff[q] = offset;
1055ea844a1aSMatthew Knepley         offset += sf->atlasDof[q];
1056ea844a1aSMatthew Knepley       }
1057ea844a1aSMatthew Knepley     }
1058ea844a1aSMatthew Knepley     /* Disable point offsets since these are unused */
1059ea844a1aSMatthew Knepley     for (p = 0; p < s->pEnd - s->pStart; ++p) {
1060ea844a1aSMatthew Knepley       s->atlasOff[p] = -1;
1061ea844a1aSMatthew Knepley       s->maxDof      = PetscMax(s->maxDof, s->atlasDof[p]);
1062ea844a1aSMatthew Knepley     }
1063ea844a1aSMatthew Knepley   }
1064ea844a1aSMatthew Knepley   if (s->perm) {ierr = ISRestoreIndices(s->perm, &pind);CHKERRQ(ierr);}
1065ea844a1aSMatthew Knepley   /* Setup BC sections */
1066ea844a1aSMatthew Knepley   ierr = PetscSectionSetUpBC(s);CHKERRQ(ierr);
1067ea844a1aSMatthew Knepley   for (f = 0; f < s->numFields; ++f) {ierr = PetscSectionSetUpBC(s->field[f]);CHKERRQ(ierr);}
1068ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1069ea844a1aSMatthew Knepley }
1070ea844a1aSMatthew Knepley 
1071ea844a1aSMatthew Knepley /*@
1072ea844a1aSMatthew Knepley   PetscSectionGetMaxDof - Return the maximum number of degrees of freedom on any point in the chart
1073ea844a1aSMatthew Knepley 
1074ea844a1aSMatthew Knepley   Not collective
1075ea844a1aSMatthew Knepley 
1076ea844a1aSMatthew Knepley   Input Parameters:
1077ea844a1aSMatthew Knepley . s - the PetscSection
1078ea844a1aSMatthew Knepley 
1079ea844a1aSMatthew Knepley   Output Parameter:
1080ea844a1aSMatthew Knepley . maxDof - the maximum dof
1081ea844a1aSMatthew Knepley 
1082ea844a1aSMatthew Knepley   Level: intermediate
1083ea844a1aSMatthew Knepley 
1084ea844a1aSMatthew Knepley .seealso: PetscSectionGetDof(), PetscSectionSetDof(), PetscSectionCreate()
1085ea844a1aSMatthew Knepley @*/
1086ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetMaxDof(PetscSection s, PetscInt *maxDof)
1087ea844a1aSMatthew Knepley {
1088ea844a1aSMatthew Knepley   PetscFunctionBegin;
1089ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1090ea844a1aSMatthew Knepley   PetscValidPointer(maxDof, 2);
1091ea844a1aSMatthew Knepley   *maxDof = s->maxDof;
1092ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1093ea844a1aSMatthew Knepley }
1094ea844a1aSMatthew Knepley 
1095ea844a1aSMatthew Knepley /*@
1096ea844a1aSMatthew Knepley   PetscSectionGetStorageSize - Return the size of an array or local Vec capable of holding all the degrees of freedom.
1097ea844a1aSMatthew Knepley 
1098ea844a1aSMatthew Knepley   Not collective
1099ea844a1aSMatthew Knepley 
1100ea844a1aSMatthew Knepley   Input Parameter:
1101ea844a1aSMatthew Knepley . s - the PetscSection
1102ea844a1aSMatthew Knepley 
1103ea844a1aSMatthew Knepley   Output Parameter:
1104ea844a1aSMatthew Knepley . size - the size of an array which can hold all the dofs
1105ea844a1aSMatthew Knepley 
1106ea844a1aSMatthew Knepley   Level: intermediate
1107ea844a1aSMatthew Knepley 
1108ea844a1aSMatthew Knepley .seealso: PetscSectionGetOffset(), PetscSectionGetConstrainedStorageSize(), PetscSectionCreate()
1109ea844a1aSMatthew Knepley @*/
1110ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetStorageSize(PetscSection s, PetscInt *size)
1111ea844a1aSMatthew Knepley {
1112ea844a1aSMatthew Knepley   PetscInt p, n = 0;
1113ea844a1aSMatthew Knepley 
1114ea844a1aSMatthew Knepley   PetscFunctionBegin;
1115ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1116ea844a1aSMatthew Knepley   PetscValidPointer(size, 2);
1117ea844a1aSMatthew Knepley   for (p = 0; p < s->pEnd - s->pStart; ++p) n += s->atlasDof[p] > 0 ? s->atlasDof[p] : 0;
1118ea844a1aSMatthew Knepley   *size = n;
1119ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1120ea844a1aSMatthew Knepley }
1121ea844a1aSMatthew Knepley 
1122ea844a1aSMatthew Knepley /*@
1123ea844a1aSMatthew Knepley   PetscSectionGetConstrainedStorageSize - Return the size of an array or local Vec capable of holding all unconstrained degrees of freedom.
1124ea844a1aSMatthew Knepley 
1125ea844a1aSMatthew Knepley   Not collective
1126ea844a1aSMatthew Knepley 
1127ea844a1aSMatthew Knepley   Input Parameters:
1128ea844a1aSMatthew Knepley + s - the PetscSection
1129ea844a1aSMatthew Knepley - point - the point
1130ea844a1aSMatthew Knepley 
1131ea844a1aSMatthew Knepley   Output Parameter:
1132ea844a1aSMatthew Knepley . size - the size of an array which can hold all unconstrained dofs
1133ea844a1aSMatthew Knepley 
1134ea844a1aSMatthew Knepley   Level: intermediate
1135ea844a1aSMatthew Knepley 
1136ea844a1aSMatthew Knepley .seealso: PetscSectionGetStorageSize(), PetscSectionGetOffset(), PetscSectionCreate()
1137ea844a1aSMatthew Knepley @*/
1138ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetConstrainedStorageSize(PetscSection s, PetscInt *size)
1139ea844a1aSMatthew Knepley {
1140ea844a1aSMatthew Knepley   PetscInt p, n = 0;
1141ea844a1aSMatthew Knepley 
1142ea844a1aSMatthew Knepley   PetscFunctionBegin;
1143ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1144ea844a1aSMatthew Knepley   PetscValidPointer(size, 2);
1145ea844a1aSMatthew Knepley   for (p = 0; p < s->pEnd - s->pStart; ++p) {
1146ea844a1aSMatthew Knepley     const PetscInt cdof = s->bc ? s->bc->atlasDof[p] : 0;
1147ea844a1aSMatthew Knepley     n += s->atlasDof[p] > 0 ? s->atlasDof[p] - cdof : 0;
1148ea844a1aSMatthew Knepley   }
1149ea844a1aSMatthew Knepley   *size = n;
1150ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1151ea844a1aSMatthew Knepley }
1152ea844a1aSMatthew Knepley 
1153ea844a1aSMatthew Knepley /*@
1154ea844a1aSMatthew Knepley   PetscSectionCreateGlobalSection - Create a section describing the global field layout using
1155ea844a1aSMatthew Knepley   the local section and an SF describing the section point overlap.
1156ea844a1aSMatthew Knepley 
1157ea844a1aSMatthew Knepley   Input Parameters:
1158ea844a1aSMatthew Knepley + s - The PetscSection for the local field layout
1159ea844a1aSMatthew Knepley . sf - The SF describing parallel layout of the section points (leaves are unowned local points)
1160ea844a1aSMatthew Knepley . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1161ea844a1aSMatthew Knepley - localOffsets - If PETSC_TRUE, use local rather than global offsets for the points
1162ea844a1aSMatthew Knepley 
1163ea844a1aSMatthew Knepley   Output Parameter:
1164ea844a1aSMatthew Knepley . gsection - The PetscSection for the global field layout
1165ea844a1aSMatthew Knepley 
1166ea844a1aSMatthew Knepley   Note: This gives negative sizes and offsets to points not owned by this process
1167ea844a1aSMatthew Knepley 
1168ea844a1aSMatthew Knepley   Level: intermediate
1169ea844a1aSMatthew Knepley 
1170ea844a1aSMatthew Knepley .seealso: PetscSectionCreate()
1171ea844a1aSMatthew Knepley @*/
1172ea844a1aSMatthew Knepley PetscErrorCode PetscSectionCreateGlobalSection(PetscSection s, PetscSF sf, PetscBool includeConstraints, PetscBool localOffsets, PetscSection *gsection)
1173ea844a1aSMatthew Knepley {
1174ea844a1aSMatthew Knepley   PetscSection    gs;
1175ea844a1aSMatthew Knepley   const PetscInt *pind = NULL;
1176ea844a1aSMatthew Knepley   PetscInt       *recv = NULL, *neg = NULL;
1177ea844a1aSMatthew Knepley   PetscInt        pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots, nlocal, maxleaf;
1178ea844a1aSMatthew Knepley   PetscErrorCode  ierr;
1179ea844a1aSMatthew Knepley 
1180ea844a1aSMatthew Knepley   PetscFunctionBegin;
1181ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1182ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1183ea844a1aSMatthew Knepley   PetscValidLogicalCollectiveBool(s, includeConstraints, 3);
1184ea844a1aSMatthew Knepley   PetscValidLogicalCollectiveBool(s, localOffsets, 4);
1185ea844a1aSMatthew Knepley   PetscValidPointer(gsection, 5);
1186ea844a1aSMatthew Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), &gs);CHKERRQ(ierr);
1187ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1188ea844a1aSMatthew Knepley   ierr = PetscSectionSetChart(gs, pStart, pEnd);CHKERRQ(ierr);
1189ea844a1aSMatthew Knepley   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1190ea844a1aSMatthew Knepley   nlocal = nroots;              /* The local/leaf space matches global/root space */
1191ea844a1aSMatthew Knepley   /* Must allocate for all points visible to SF, which may be more than this section */
1192ea844a1aSMatthew Knepley   if (nroots >= 0) {             /* nroots < 0 means that the graph has not been set, only happens in serial */
1193ea844a1aSMatthew Knepley     ierr = PetscSFGetLeafRange(sf, NULL, &maxleaf);CHKERRQ(ierr);
1194ea844a1aSMatthew Knepley     if (nroots < pEnd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "SF roots %D < pEnd %D", nroots, pEnd);
1195ea844a1aSMatthew Knepley     if (maxleaf >= nroots) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Max local leaf %D >= nroots %D", maxleaf, nroots);
1196ea844a1aSMatthew Knepley     ierr = PetscMalloc2(nroots,&neg,nlocal,&recv);CHKERRQ(ierr);
1197ea844a1aSMatthew Knepley     ierr = PetscArrayzero(neg,nroots);CHKERRQ(ierr);
1198ea844a1aSMatthew Knepley   }
1199ea844a1aSMatthew Knepley   /* Mark all local points with negative dof */
1200ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
1201ea844a1aSMatthew Knepley     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1202ea844a1aSMatthew Knepley     ierr = PetscSectionSetDof(gs, p, dof);CHKERRQ(ierr);
1203ea844a1aSMatthew Knepley     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1204ea844a1aSMatthew Knepley     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(gs, p, cdof);CHKERRQ(ierr);}
1205ea844a1aSMatthew Knepley     if (neg) neg[p] = -(dof+1);
1206ea844a1aSMatthew Knepley   }
1207ea844a1aSMatthew Knepley   ierr = PetscSectionSetUpBC(gs);CHKERRQ(ierr);
1208ea844a1aSMatthew Knepley   if (gs->bcIndices) {ierr = PetscArraycpy(gs->bcIndices, s->bcIndices,gs->bc->atlasOff[gs->bc->pEnd-gs->bc->pStart-1] + gs->bc->atlasDof[gs->bc->pEnd-gs->bc->pStart-1]);CHKERRQ(ierr);}
1209ea844a1aSMatthew Knepley   if (nroots >= 0) {
1210ea844a1aSMatthew Knepley     ierr = PetscArrayzero(recv,nlocal);CHKERRQ(ierr);
1211ea844a1aSMatthew Knepley     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, recv);CHKERRQ(ierr);
1212ea844a1aSMatthew Knepley     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, recv);CHKERRQ(ierr);
1213ea844a1aSMatthew Knepley     for (p = pStart; p < pEnd; ++p) {
1214ea844a1aSMatthew Knepley       if (recv[p] < 0) {
1215ea844a1aSMatthew Knepley         gs->atlasDof[p-pStart] = recv[p];
1216ea844a1aSMatthew Knepley         ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1217ea844a1aSMatthew Knepley         if (-(recv[p]+1) != dof) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global dof %D for point %D is not the unconstrained %D", -(recv[p]+1), p, dof);
1218ea844a1aSMatthew Knepley       }
1219ea844a1aSMatthew Knepley     }
1220ea844a1aSMatthew Knepley   }
1221ea844a1aSMatthew Knepley   /* Calculate new sizes, get process offset, and calculate point offsets */
1222ea844a1aSMatthew Knepley   if (s->perm) {ierr = ISGetIndices(s->perm, &pind);CHKERRQ(ierr);}
1223ea844a1aSMatthew Knepley   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1224ea844a1aSMatthew Knepley     const PetscInt q = pind ? pind[p] : p;
1225ea844a1aSMatthew Knepley 
1226ea844a1aSMatthew Knepley     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[q] : 0;
1227ea844a1aSMatthew Knepley     gs->atlasOff[q] = off;
1228ea844a1aSMatthew Knepley     off += gs->atlasDof[q] > 0 ? gs->atlasDof[q]-cdof : 0;
1229ea844a1aSMatthew Knepley   }
1230ea844a1aSMatthew Knepley   if (!localOffsets) {
1231ea844a1aSMatthew Knepley     ierr = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1232ea844a1aSMatthew Knepley     globalOff -= off;
1233ea844a1aSMatthew Knepley   }
1234ea844a1aSMatthew Knepley   for (p = pStart, off = 0; p < pEnd; ++p) {
1235ea844a1aSMatthew Knepley     gs->atlasOff[p-pStart] += globalOff;
1236ea844a1aSMatthew Knepley     if (neg) neg[p] = -(gs->atlasOff[p-pStart]+1);
1237ea844a1aSMatthew Knepley   }
1238ea844a1aSMatthew Knepley   if (s->perm) {ierr = ISRestoreIndices(s->perm, &pind);CHKERRQ(ierr);}
1239ea844a1aSMatthew Knepley   /* Put in negative offsets for ghost points */
1240ea844a1aSMatthew Knepley   if (nroots >= 0) {
1241ea844a1aSMatthew Knepley     ierr = PetscArrayzero(recv,nlocal);CHKERRQ(ierr);
1242ea844a1aSMatthew Knepley     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, recv);CHKERRQ(ierr);
1243ea844a1aSMatthew Knepley     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, recv);CHKERRQ(ierr);
1244ea844a1aSMatthew Knepley     for (p = pStart; p < pEnd; ++p) {
1245ea844a1aSMatthew Knepley       if (recv[p] < 0) gs->atlasOff[p-pStart] = recv[p];
1246ea844a1aSMatthew Knepley     }
1247ea844a1aSMatthew Knepley   }
1248ea844a1aSMatthew Knepley   ierr = PetscFree2(neg,recv);CHKERRQ(ierr);
1249ea844a1aSMatthew Knepley   ierr = PetscSectionViewFromOptions(gs, NULL, "-global_section_view");CHKERRQ(ierr);
1250ea844a1aSMatthew Knepley   *gsection = gs;
1251ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1252ea844a1aSMatthew Knepley }
1253ea844a1aSMatthew Knepley 
1254ea844a1aSMatthew Knepley /*@
1255ea844a1aSMatthew Knepley   PetscSectionCreateGlobalSectionCensored - Create a section describing the global field layout using
1256ea844a1aSMatthew Knepley   the local section and an SF describing the section point overlap.
1257ea844a1aSMatthew Knepley 
1258ea844a1aSMatthew Knepley   Input Parameters:
1259ea844a1aSMatthew Knepley   + s - The PetscSection for the local field layout
1260ea844a1aSMatthew Knepley   . sf - The SF describing parallel layout of the section points
1261ea844a1aSMatthew Knepley   . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
1262ea844a1aSMatthew Knepley   . numExcludes - The number of exclusion ranges
1263ea844a1aSMatthew Knepley   - excludes - An array [start_0, end_0, start_1, end_1, ...] where there are numExcludes pairs
1264ea844a1aSMatthew Knepley 
1265ea844a1aSMatthew Knepley   Output Parameter:
1266ea844a1aSMatthew Knepley   . gsection - The PetscSection for the global field layout
1267ea844a1aSMatthew Knepley 
1268ea844a1aSMatthew Knepley   Note: This gives negative sizes and offsets to points not owned by this process
1269ea844a1aSMatthew Knepley 
1270ea844a1aSMatthew Knepley   Level: advanced
1271ea844a1aSMatthew Knepley 
1272ea844a1aSMatthew Knepley .seealso: PetscSectionCreate()
1273ea844a1aSMatthew Knepley @*/
1274ea844a1aSMatthew Knepley PetscErrorCode PetscSectionCreateGlobalSectionCensored(PetscSection s, PetscSF sf, PetscBool includeConstraints, PetscInt numExcludes, const PetscInt excludes[], PetscSection *gsection)
1275ea844a1aSMatthew Knepley {
1276ea844a1aSMatthew Knepley   const PetscInt *pind = NULL;
1277ea844a1aSMatthew Knepley   PetscInt       *neg  = NULL, *tmpOff = NULL;
1278ea844a1aSMatthew Knepley   PetscInt        pStart, pEnd, p, e, dof, cdof, off, globalOff = 0, nroots;
1279ea844a1aSMatthew Knepley   PetscErrorCode  ierr;
1280ea844a1aSMatthew Knepley 
1281ea844a1aSMatthew Knepley   PetscFunctionBegin;
1282ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1283ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1284ea844a1aSMatthew Knepley   PetscValidPointer(gsection, 6);
1285ea844a1aSMatthew Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection);CHKERRQ(ierr);
1286ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1287ea844a1aSMatthew Knepley   ierr = PetscSectionSetChart(*gsection, pStart, pEnd);CHKERRQ(ierr);
1288ea844a1aSMatthew Knepley   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1289ea844a1aSMatthew Knepley   if (nroots >= 0) {
1290ea844a1aSMatthew Knepley     if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %D < %D section size", nroots, pEnd-pStart);
1291ea844a1aSMatthew Knepley     ierr = PetscCalloc1(nroots, &neg);CHKERRQ(ierr);
1292ea844a1aSMatthew Knepley     if (nroots > pEnd-pStart) {
1293ea844a1aSMatthew Knepley       ierr = PetscCalloc1(nroots, &tmpOff);CHKERRQ(ierr);
1294ea844a1aSMatthew Knepley     } else {
1295ea844a1aSMatthew Knepley       tmpOff = &(*gsection)->atlasDof[-pStart];
1296ea844a1aSMatthew Knepley     }
1297ea844a1aSMatthew Knepley   }
1298ea844a1aSMatthew Knepley   /* Mark ghost points with negative dof */
1299ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
1300ea844a1aSMatthew Knepley     for (e = 0; e < numExcludes; ++e) {
1301ea844a1aSMatthew Knepley       if ((p >= excludes[e*2+0]) && (p < excludes[e*2+1])) {
1302ea844a1aSMatthew Knepley         ierr = PetscSectionSetDof(*gsection, p, 0);CHKERRQ(ierr);
1303ea844a1aSMatthew Knepley         break;
1304ea844a1aSMatthew Knepley       }
1305ea844a1aSMatthew Knepley     }
1306ea844a1aSMatthew Knepley     if (e < numExcludes) continue;
1307ea844a1aSMatthew Knepley     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1308ea844a1aSMatthew Knepley     ierr = PetscSectionSetDof(*gsection, p, dof);CHKERRQ(ierr);
1309ea844a1aSMatthew Knepley     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1310ea844a1aSMatthew Knepley     if (!includeConstraints && cdof > 0) {ierr = PetscSectionSetConstraintDof(*gsection, p, cdof);CHKERRQ(ierr);}
1311ea844a1aSMatthew Knepley     if (neg) neg[p] = -(dof+1);
1312ea844a1aSMatthew Knepley   }
1313ea844a1aSMatthew Knepley   ierr = PetscSectionSetUpBC(*gsection);CHKERRQ(ierr);
1314ea844a1aSMatthew Knepley   if (nroots >= 0) {
1315ea844a1aSMatthew Knepley     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1316ea844a1aSMatthew Knepley     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1317ea844a1aSMatthew Knepley     if (nroots > pEnd - pStart) {
1318ea844a1aSMatthew Knepley       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
1319ea844a1aSMatthew Knepley     }
1320ea844a1aSMatthew Knepley   }
1321ea844a1aSMatthew Knepley   /* Calculate new sizes, get proccess offset, and calculate point offsets */
1322ea844a1aSMatthew Knepley   if (s->perm) {ierr = ISGetIndices(s->perm, &pind);CHKERRQ(ierr);}
1323ea844a1aSMatthew Knepley   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1324ea844a1aSMatthew Knepley     const PetscInt q = pind ? pind[p] : p;
1325ea844a1aSMatthew Knepley 
1326ea844a1aSMatthew Knepley     cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[q] : 0;
1327ea844a1aSMatthew Knepley     (*gsection)->atlasOff[q] = off;
1328ea844a1aSMatthew Knepley     off += (*gsection)->atlasDof[q] > 0 ? (*gsection)->atlasDof[q]-cdof : 0;
1329ea844a1aSMatthew Knepley   }
1330ea844a1aSMatthew Knepley   ierr = MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s));CHKERRQ(ierr);
1331ea844a1aSMatthew Knepley   globalOff -= off;
1332ea844a1aSMatthew Knepley   for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1333ea844a1aSMatthew Knepley     (*gsection)->atlasOff[p] += globalOff;
1334ea844a1aSMatthew Knepley     if (neg) neg[p+pStart] = -((*gsection)->atlasOff[p]+1);
1335ea844a1aSMatthew Knepley   }
1336ea844a1aSMatthew Knepley   if (s->perm) {ierr = ISRestoreIndices(s->perm, &pind);CHKERRQ(ierr);}
1337ea844a1aSMatthew Knepley   /* Put in negative offsets for ghost points */
1338ea844a1aSMatthew Knepley   if (nroots >= 0) {
1339ea844a1aSMatthew Knepley     if (nroots == pEnd-pStart) tmpOff = &(*gsection)->atlasOff[-pStart];
1340ea844a1aSMatthew Knepley     ierr = PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1341ea844a1aSMatthew Knepley     ierr = PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);CHKERRQ(ierr);
1342ea844a1aSMatthew Knepley     if (nroots > pEnd - pStart) {
1343ea844a1aSMatthew Knepley       for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1344ea844a1aSMatthew Knepley     }
1345ea844a1aSMatthew Knepley   }
1346ea844a1aSMatthew Knepley   if (nroots >= 0 && nroots > pEnd-pStart) {ierr = PetscFree(tmpOff);CHKERRQ(ierr);}
1347ea844a1aSMatthew Knepley   ierr = PetscFree(neg);CHKERRQ(ierr);
1348ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1349ea844a1aSMatthew Knepley }
1350ea844a1aSMatthew Knepley 
1351ea844a1aSMatthew Knepley /*@
1352ea844a1aSMatthew Knepley   PetscSectionGetPointLayout - Get the PetscLayout associated with the section points
1353ea844a1aSMatthew Knepley 
1354ea844a1aSMatthew Knepley   Collective on comm
1355ea844a1aSMatthew Knepley 
1356ea844a1aSMatthew Knepley   Input Parameters:
1357ea844a1aSMatthew Knepley + comm - The MPI_Comm
1358ea844a1aSMatthew Knepley - s    - The PetscSection
1359ea844a1aSMatthew Knepley 
1360ea844a1aSMatthew Knepley   Output Parameter:
1361ea844a1aSMatthew Knepley . layout - The point layout for the section
1362ea844a1aSMatthew Knepley 
1363ea844a1aSMatthew Knepley   Note: This is usually caleld for the default global section.
1364ea844a1aSMatthew Knepley 
1365ea844a1aSMatthew Knepley   Level: advanced
1366ea844a1aSMatthew Knepley 
1367ea844a1aSMatthew Knepley .seealso: PetscSectionGetValueLayout(), PetscSectionCreate()
1368ea844a1aSMatthew Knepley @*/
1369ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetPointLayout(MPI_Comm comm, PetscSection s, PetscLayout *layout)
1370ea844a1aSMatthew Knepley {
1371ea844a1aSMatthew Knepley   PetscInt       pStart, pEnd, p, localSize = 0;
1372ea844a1aSMatthew Knepley   PetscErrorCode ierr;
1373ea844a1aSMatthew Knepley 
1374ea844a1aSMatthew Knepley   PetscFunctionBegin;
1375ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1376ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
1377ea844a1aSMatthew Knepley     PetscInt dof;
1378ea844a1aSMatthew Knepley 
1379ea844a1aSMatthew Knepley     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1380ea844a1aSMatthew Knepley     if (dof > 0) ++localSize;
1381ea844a1aSMatthew Knepley   }
1382ea844a1aSMatthew Knepley   ierr = PetscLayoutCreate(comm, layout);CHKERRQ(ierr);
1383ea844a1aSMatthew Knepley   ierr = PetscLayoutSetLocalSize(*layout, localSize);CHKERRQ(ierr);
1384ea844a1aSMatthew Knepley   ierr = PetscLayoutSetBlockSize(*layout, 1);CHKERRQ(ierr);
1385ea844a1aSMatthew Knepley   ierr = PetscLayoutSetUp(*layout);CHKERRQ(ierr);
1386ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1387ea844a1aSMatthew Knepley }
1388ea844a1aSMatthew Knepley 
1389ea844a1aSMatthew Knepley /*@
1390ea844a1aSMatthew Knepley   PetscSectionGetValueLayout - Get the PetscLayout associated with the section dofs.
1391ea844a1aSMatthew Knepley 
1392ea844a1aSMatthew Knepley   Collective on comm
1393ea844a1aSMatthew Knepley 
1394ea844a1aSMatthew Knepley   Input Parameters:
1395ea844a1aSMatthew Knepley + comm - The MPI_Comm
1396ea844a1aSMatthew Knepley - s    - The PetscSection
1397ea844a1aSMatthew Knepley 
1398ea844a1aSMatthew Knepley   Output Parameter:
1399ea844a1aSMatthew Knepley . layout - The dof layout for the section
1400ea844a1aSMatthew Knepley 
1401ea844a1aSMatthew Knepley   Note: This is usually called for the default global section.
1402ea844a1aSMatthew Knepley 
1403ea844a1aSMatthew Knepley   Level: advanced
1404ea844a1aSMatthew Knepley 
1405ea844a1aSMatthew Knepley .seealso: PetscSectionGetPointLayout(), PetscSectionCreate()
1406ea844a1aSMatthew Knepley @*/
1407ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetValueLayout(MPI_Comm comm, PetscSection s, PetscLayout *layout)
1408ea844a1aSMatthew Knepley {
1409ea844a1aSMatthew Knepley   PetscInt       pStart, pEnd, p, localSize = 0;
1410ea844a1aSMatthew Knepley   PetscErrorCode ierr;
1411ea844a1aSMatthew Knepley 
1412ea844a1aSMatthew Knepley   PetscFunctionBegin;
1413ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 2);
1414ea844a1aSMatthew Knepley   PetscValidPointer(layout, 3);
1415ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1416ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
1417ea844a1aSMatthew Knepley     PetscInt dof,cdof;
1418ea844a1aSMatthew Knepley 
1419ea844a1aSMatthew Knepley     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1420ea844a1aSMatthew Knepley     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1421ea844a1aSMatthew Knepley     if (dof-cdof > 0) localSize += dof-cdof;
1422ea844a1aSMatthew Knepley   }
1423ea844a1aSMatthew Knepley   ierr = PetscLayoutCreate(comm, layout);CHKERRQ(ierr);
1424ea844a1aSMatthew Knepley   ierr = PetscLayoutSetLocalSize(*layout, localSize);CHKERRQ(ierr);
1425ea844a1aSMatthew Knepley   ierr = PetscLayoutSetBlockSize(*layout, 1);CHKERRQ(ierr);
1426ea844a1aSMatthew Knepley   ierr = PetscLayoutSetUp(*layout);CHKERRQ(ierr);
1427ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1428ea844a1aSMatthew Knepley }
1429ea844a1aSMatthew Knepley 
1430ea844a1aSMatthew Knepley /*@
1431ea844a1aSMatthew Knepley   PetscSectionGetOffset - Return the offset into an array or local Vec for the dof associated with the given point.
1432ea844a1aSMatthew Knepley 
1433ea844a1aSMatthew Knepley   Not collective
1434ea844a1aSMatthew Knepley 
1435ea844a1aSMatthew Knepley   Input Parameters:
1436ea844a1aSMatthew Knepley + s - the PetscSection
1437ea844a1aSMatthew Knepley - point - the point
1438ea844a1aSMatthew Knepley 
1439ea844a1aSMatthew Knepley   Output Parameter:
1440ea844a1aSMatthew Knepley . offset - the offset
1441ea844a1aSMatthew Knepley 
1442ea844a1aSMatthew Knepley   Level: intermediate
1443ea844a1aSMatthew Knepley 
1444ea844a1aSMatthew Knepley .seealso: PetscSectionGetFieldOffset(), PetscSectionCreate()
1445ea844a1aSMatthew Knepley @*/
1446ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetOffset(PetscSection s, PetscInt point, PetscInt *offset)
1447ea844a1aSMatthew Knepley {
1448ea844a1aSMatthew Knepley   PetscFunctionBegin;
1449ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1450ea844a1aSMatthew Knepley   PetscValidPointer(offset, 3);
1451ea844a1aSMatthew Knepley #if defined(PETSC_USE_DEBUG)
1452ea844a1aSMatthew Knepley   if ((point < s->pStart) || (point >= s->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %D should be in [%D, %D)", point, s->pStart, s->pEnd);
1453ea844a1aSMatthew Knepley #endif
1454ea844a1aSMatthew Knepley   *offset = s->atlasOff[point - s->pStart];
1455ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1456ea844a1aSMatthew Knepley }
1457ea844a1aSMatthew Knepley 
1458ea844a1aSMatthew Knepley /*@
1459ea844a1aSMatthew Knepley   PetscSectionSetOffset - Set the offset into an array or local Vec for the dof associated with the given point.
1460ea844a1aSMatthew Knepley 
1461ea844a1aSMatthew Knepley   Not collective
1462ea844a1aSMatthew Knepley 
1463ea844a1aSMatthew Knepley   Input Parameters:
1464ea844a1aSMatthew Knepley + s - the PetscSection
1465ea844a1aSMatthew Knepley . point - the point
1466ea844a1aSMatthew Knepley - offset - the offset
1467ea844a1aSMatthew Knepley 
1468ea844a1aSMatthew Knepley   Note: The user usually does not call this function, but uses PetscSectionSetUp()
1469ea844a1aSMatthew Knepley 
1470ea844a1aSMatthew Knepley   Level: intermediate
1471ea844a1aSMatthew Knepley 
1472ea844a1aSMatthew Knepley .seealso: PetscSectionGetFieldOffset(), PetscSectionCreate(), PetscSectionSetUp()
1473ea844a1aSMatthew Knepley @*/
1474ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetOffset(PetscSection s, PetscInt point, PetscInt offset)
1475ea844a1aSMatthew Knepley {
1476ea844a1aSMatthew Knepley   PetscFunctionBegin;
1477ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1478ea844a1aSMatthew Knepley   if ((point < s->pStart) || (point >= s->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %D should be in [%D, %D)", point, s->pStart, s->pEnd);
1479ea844a1aSMatthew Knepley   s->atlasOff[point - s->pStart] = offset;
1480ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1481ea844a1aSMatthew Knepley }
1482ea844a1aSMatthew Knepley 
1483ea844a1aSMatthew Knepley /*@
1484ea844a1aSMatthew Knepley   PetscSectionGetFieldOffset - Return the offset into an array or local Vec for the dof associated with the given point.
1485ea844a1aSMatthew Knepley 
1486ea844a1aSMatthew Knepley   Not collective
1487ea844a1aSMatthew Knepley 
1488ea844a1aSMatthew Knepley   Input Parameters:
1489ea844a1aSMatthew Knepley + s - the PetscSection
1490ea844a1aSMatthew Knepley . point - the point
1491ea844a1aSMatthew Knepley - field - the field
1492ea844a1aSMatthew Knepley 
1493ea844a1aSMatthew Knepley   Output Parameter:
1494ea844a1aSMatthew Knepley . offset - the offset
1495ea844a1aSMatthew Knepley 
1496ea844a1aSMatthew Knepley   Level: intermediate
1497ea844a1aSMatthew Knepley 
1498ea844a1aSMatthew Knepley .seealso: PetscSectionGetOffset(), PetscSectionCreate()
1499ea844a1aSMatthew Knepley @*/
1500ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetFieldOffset(PetscSection s, PetscInt point, PetscInt field, PetscInt *offset)
1501ea844a1aSMatthew Knepley {
1502ea844a1aSMatthew Knepley   PetscErrorCode ierr;
1503ea844a1aSMatthew Knepley 
1504ea844a1aSMatthew Knepley   PetscFunctionBegin;
1505ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1506ea844a1aSMatthew Knepley   PetscValidPointer(offset, 4);
1507ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
1508ea844a1aSMatthew Knepley   ierr = PetscSectionGetOffset(s->field[field], point, offset);CHKERRQ(ierr);
1509ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1510ea844a1aSMatthew Knepley }
1511ea844a1aSMatthew Knepley 
1512ea844a1aSMatthew Knepley /*@
1513ea844a1aSMatthew Knepley   PetscSectionSetFieldOffset - Set the offset into an array or local Vec for the dof associated with the given point.
1514ea844a1aSMatthew Knepley 
1515ea844a1aSMatthew Knepley   Not collective
1516ea844a1aSMatthew Knepley 
1517ea844a1aSMatthew Knepley   Input Parameters:
1518ea844a1aSMatthew Knepley + s - the PetscSection
1519ea844a1aSMatthew Knepley . point - the point
1520ea844a1aSMatthew Knepley . field - the field
1521ea844a1aSMatthew Knepley - offset - the offset
1522ea844a1aSMatthew Knepley 
1523ea844a1aSMatthew Knepley   Note: The user usually does not call this function, but uses PetscSectionSetUp()
1524ea844a1aSMatthew Knepley 
1525ea844a1aSMatthew Knepley   Level: intermediate
1526ea844a1aSMatthew Knepley 
1527ea844a1aSMatthew Knepley .seealso: PetscSectionGetOffset(), PetscSectionCreate(), PetscSectionSetUp()
1528ea844a1aSMatthew Knepley @*/
1529ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetFieldOffset(PetscSection s, PetscInt point, PetscInt field, PetscInt offset)
1530ea844a1aSMatthew Knepley {
1531ea844a1aSMatthew Knepley   PetscErrorCode ierr;
1532ea844a1aSMatthew Knepley 
1533ea844a1aSMatthew Knepley   PetscFunctionBegin;
1534ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1535ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
1536ea844a1aSMatthew Knepley   ierr = PetscSectionSetOffset(s->field[field], point, offset);CHKERRQ(ierr);
1537ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1538ea844a1aSMatthew Knepley }
1539ea844a1aSMatthew Knepley 
1540ea844a1aSMatthew Knepley /*@
1541ea844a1aSMatthew Knepley   PetscSectionGetFieldPointOffset - Return the offset on the given point for the dof associated with the given point.
1542ea844a1aSMatthew Knepley 
1543ea844a1aSMatthew Knepley   Not collective
1544ea844a1aSMatthew Knepley 
1545ea844a1aSMatthew Knepley   Input Parameters:
1546ea844a1aSMatthew Knepley + s - the PetscSection
1547ea844a1aSMatthew Knepley . point - the point
1548ea844a1aSMatthew Knepley - field - the field
1549ea844a1aSMatthew Knepley 
1550ea844a1aSMatthew Knepley   Output Parameter:
1551ea844a1aSMatthew Knepley . offset - the offset
1552ea844a1aSMatthew Knepley 
1553ea844a1aSMatthew Knepley   Note: This gives the offset on a point of the field, ignoring constraints, meaning starting at the first dof for
1554ea844a1aSMatthew Knepley         this point, what number is the first dof with this field.
1555ea844a1aSMatthew Knepley 
1556ea844a1aSMatthew Knepley   Level: advanced
1557ea844a1aSMatthew Knepley 
1558ea844a1aSMatthew Knepley .seealso: PetscSectionGetOffset(), PetscSectionCreate()
1559ea844a1aSMatthew Knepley @*/
1560ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetFieldPointOffset(PetscSection s, PetscInt point, PetscInt field, PetscInt *offset)
1561ea844a1aSMatthew Knepley {
1562ea844a1aSMatthew Knepley   PetscInt       off, foff;
1563ea844a1aSMatthew Knepley   PetscErrorCode ierr;
1564ea844a1aSMatthew Knepley 
1565ea844a1aSMatthew Knepley   PetscFunctionBegin;
1566ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1567ea844a1aSMatthew Knepley   PetscValidPointer(offset, 4);
1568ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
1569ea844a1aSMatthew Knepley   ierr = PetscSectionGetOffset(s, point, &off);CHKERRQ(ierr);
1570ea844a1aSMatthew Knepley   ierr = PetscSectionGetOffset(s->field[field], point, &foff);CHKERRQ(ierr);
1571ea844a1aSMatthew Knepley   *offset = foff - off;
1572ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1573ea844a1aSMatthew Knepley }
1574ea844a1aSMatthew Knepley 
1575ea844a1aSMatthew Knepley /*@
1576ea844a1aSMatthew Knepley   PetscSectionGetOffsetRange - Return the full range of offsets [start, end)
1577ea844a1aSMatthew Knepley 
1578ea844a1aSMatthew Knepley   Not collective
1579ea844a1aSMatthew Knepley 
1580ea844a1aSMatthew Knepley   Input Parameter:
1581ea844a1aSMatthew Knepley . s - the PetscSection
1582ea844a1aSMatthew Knepley 
1583ea844a1aSMatthew Knepley   Output Parameters:
1584ea844a1aSMatthew Knepley + start - the minimum offset
1585ea844a1aSMatthew Knepley - end   - one more than the maximum offset
1586ea844a1aSMatthew Knepley 
1587ea844a1aSMatthew Knepley   Level: intermediate
1588ea844a1aSMatthew Knepley 
1589ea844a1aSMatthew Knepley .seealso: PetscSectionGetOffset(), PetscSectionCreate()
1590ea844a1aSMatthew Knepley @*/
1591ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetOffsetRange(PetscSection s, PetscInt *start, PetscInt *end)
1592ea844a1aSMatthew Knepley {
1593ea844a1aSMatthew Knepley   PetscInt       os = 0, oe = 0, pStart, pEnd, p;
1594ea844a1aSMatthew Knepley   PetscErrorCode ierr;
1595ea844a1aSMatthew Knepley 
1596ea844a1aSMatthew Knepley   PetscFunctionBegin;
1597ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1598ea844a1aSMatthew Knepley   if (s->atlasOff) {os = s->atlasOff[0]; oe = s->atlasOff[0];}
1599ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1600ea844a1aSMatthew Knepley   for (p = 0; p < pEnd-pStart; ++p) {
1601ea844a1aSMatthew Knepley     PetscInt dof = s->atlasDof[p], off = s->atlasOff[p];
1602ea844a1aSMatthew Knepley 
1603ea844a1aSMatthew Knepley     if (off >= 0) {
1604ea844a1aSMatthew Knepley       os = PetscMin(os, off);
1605ea844a1aSMatthew Knepley       oe = PetscMax(oe, off+dof);
1606ea844a1aSMatthew Knepley     }
1607ea844a1aSMatthew Knepley   }
1608ea844a1aSMatthew Knepley   if (start) *start = os;
1609ea844a1aSMatthew Knepley   if (end)   *end   = oe;
1610ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1611ea844a1aSMatthew Knepley }
1612ea844a1aSMatthew Knepley 
1613ea844a1aSMatthew Knepley /*@
1614ea844a1aSMatthew Knepley   PetscSectionCreateSubsection - Create a new, smaller section composed of only the selected fields
1615ea844a1aSMatthew Knepley 
1616ea844a1aSMatthew Knepley   Collective on s
1617ea844a1aSMatthew Knepley 
1618ea844a1aSMatthew Knepley   Input Parameter:
1619ea844a1aSMatthew Knepley + s      - the PetscSection
1620ea844a1aSMatthew Knepley . len    - the number of subfields
1621ea844a1aSMatthew Knepley - fields - the subfield numbers
1622ea844a1aSMatthew Knepley 
1623ea844a1aSMatthew Knepley   Output Parameter:
1624ea844a1aSMatthew Knepley . subs   - the subsection
1625ea844a1aSMatthew Knepley 
1626ea844a1aSMatthew Knepley   Note: The section offsets now refer to a new, smaller vector.
1627ea844a1aSMatthew Knepley 
1628ea844a1aSMatthew Knepley   Level: advanced
1629ea844a1aSMatthew Knepley 
1630ea844a1aSMatthew Knepley .seealso: PetscSectionCreateSupersection(), PetscSectionCreate()
1631ea844a1aSMatthew Knepley @*/
1632ea844a1aSMatthew Knepley PetscErrorCode PetscSectionCreateSubsection(PetscSection s, PetscInt len, const PetscInt fields[], PetscSection *subs)
1633ea844a1aSMatthew Knepley {
1634ea844a1aSMatthew Knepley   PetscInt       nF, f, pStart, pEnd, p, maxCdof = 0;
1635ea844a1aSMatthew Knepley   PetscErrorCode ierr;
1636ea844a1aSMatthew Knepley 
1637ea844a1aSMatthew Knepley   PetscFunctionBegin;
1638ea844a1aSMatthew Knepley   if (!len) PetscFunctionReturn(0);
1639ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1640ea844a1aSMatthew Knepley   PetscValidPointer(fields, 3);
1641ea844a1aSMatthew Knepley   PetscValidPointer(subs, 4);
1642ea844a1aSMatthew Knepley   ierr = PetscSectionGetNumFields(s, &nF);CHKERRQ(ierr);
1643ea844a1aSMatthew Knepley   if (len > nF) SETERRQ2(PetscObjectComm((PetscObject) s), PETSC_ERR_ARG_WRONG, "Number of requested fields %D greater than number of fields %D", len, nF);
1644ea844a1aSMatthew Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), subs);CHKERRQ(ierr);
1645ea844a1aSMatthew Knepley   ierr = PetscSectionSetNumFields(*subs, len);CHKERRQ(ierr);
1646ea844a1aSMatthew Knepley   for (f = 0; f < len; ++f) {
1647ea844a1aSMatthew Knepley     const char *name   = NULL;
1648ea844a1aSMatthew Knepley     PetscInt   numComp = 0;
1649ea844a1aSMatthew Knepley 
1650ea844a1aSMatthew Knepley     ierr = PetscSectionGetFieldName(s, fields[f], &name);CHKERRQ(ierr);
1651ea844a1aSMatthew Knepley     ierr = PetscSectionSetFieldName(*subs, f, name);CHKERRQ(ierr);
1652ea844a1aSMatthew Knepley     ierr = PetscSectionGetFieldComponents(s, fields[f], &numComp);CHKERRQ(ierr);
1653ea844a1aSMatthew Knepley     ierr = PetscSectionSetFieldComponents(*subs, f, numComp);CHKERRQ(ierr);
1654ea844a1aSMatthew Knepley   }
1655ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1656ea844a1aSMatthew Knepley   ierr = PetscSectionSetChart(*subs, pStart, pEnd);CHKERRQ(ierr);
1657ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
1658ea844a1aSMatthew Knepley     PetscInt dof = 0, cdof = 0, fdof = 0, cfdof = 0;
1659ea844a1aSMatthew Knepley 
1660ea844a1aSMatthew Knepley     for (f = 0; f < len; ++f) {
1661ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldDof(s, p, fields[f], &fdof);CHKERRQ(ierr);
1662ea844a1aSMatthew Knepley       ierr = PetscSectionSetFieldDof(*subs, p, f, fdof);CHKERRQ(ierr);
1663ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldConstraintDof(s, p, fields[f], &cfdof);CHKERRQ(ierr);
1664ea844a1aSMatthew Knepley       if (cfdof) {ierr = PetscSectionSetFieldConstraintDof(*subs, p, f, cfdof);CHKERRQ(ierr);}
1665ea844a1aSMatthew Knepley       dof  += fdof;
1666ea844a1aSMatthew Knepley       cdof += cfdof;
1667ea844a1aSMatthew Knepley     }
1668ea844a1aSMatthew Knepley     ierr = PetscSectionSetDof(*subs, p, dof);CHKERRQ(ierr);
1669ea844a1aSMatthew Knepley     if (cdof) {ierr = PetscSectionSetConstraintDof(*subs, p, cdof);CHKERRQ(ierr);}
1670ea844a1aSMatthew Knepley     maxCdof = PetscMax(cdof, maxCdof);
1671ea844a1aSMatthew Knepley   }
1672ea844a1aSMatthew Knepley   ierr = PetscSectionSetUp(*subs);CHKERRQ(ierr);
1673ea844a1aSMatthew Knepley   if (maxCdof) {
1674ea844a1aSMatthew Knepley     PetscInt *indices;
1675ea844a1aSMatthew Knepley 
1676ea844a1aSMatthew Knepley     ierr = PetscMalloc1(maxCdof, &indices);CHKERRQ(ierr);
1677ea844a1aSMatthew Knepley     for (p = pStart; p < pEnd; ++p) {
1678ea844a1aSMatthew Knepley       PetscInt cdof;
1679ea844a1aSMatthew Knepley 
1680ea844a1aSMatthew Knepley       ierr = PetscSectionGetConstraintDof(*subs, p, &cdof);CHKERRQ(ierr);
1681ea844a1aSMatthew Knepley       if (cdof) {
1682ea844a1aSMatthew Knepley         const PetscInt *oldIndices = NULL;
1683ea844a1aSMatthew Knepley         PetscInt       fdof = 0, cfdof = 0, fc, numConst = 0, fOff = 0;
1684ea844a1aSMatthew Knepley 
1685ea844a1aSMatthew Knepley         for (f = 0; f < len; ++f) {
1686ea844a1aSMatthew Knepley           PetscInt oldFoff = 0, oldf;
1687ea844a1aSMatthew Knepley 
1688ea844a1aSMatthew Knepley           ierr = PetscSectionGetFieldDof(s, p, fields[f], &fdof);CHKERRQ(ierr);
1689ea844a1aSMatthew Knepley           ierr = PetscSectionGetFieldConstraintDof(s, p, fields[f], &cfdof);CHKERRQ(ierr);
1690ea844a1aSMatthew Knepley           ierr = PetscSectionGetFieldConstraintIndices(s, p, fields[f], &oldIndices);CHKERRQ(ierr);
1691ea844a1aSMatthew Knepley           /* This can be sped up if we assume sorted fields */
1692ea844a1aSMatthew Knepley           for (oldf = 0; oldf < fields[f]; ++oldf) {
1693ea844a1aSMatthew Knepley             PetscInt oldfdof = 0;
1694ea844a1aSMatthew Knepley             ierr = PetscSectionGetFieldDof(s, p, oldf, &oldfdof);CHKERRQ(ierr);
1695ea844a1aSMatthew Knepley             oldFoff += oldfdof;
1696ea844a1aSMatthew Knepley           }
1697ea844a1aSMatthew Knepley           for (fc = 0; fc < cfdof; ++fc) indices[numConst+fc] = oldIndices[fc] + (fOff - oldFoff);
1698ea844a1aSMatthew Knepley           ierr = PetscSectionSetFieldConstraintIndices(*subs, p, f, &indices[numConst]);CHKERRQ(ierr);
1699ea844a1aSMatthew Knepley           numConst += cfdof;
1700ea844a1aSMatthew Knepley           fOff     += fdof;
1701ea844a1aSMatthew Knepley         }
1702ea844a1aSMatthew Knepley         ierr = PetscSectionSetConstraintIndices(*subs, p, indices);CHKERRQ(ierr);
1703ea844a1aSMatthew Knepley       }
1704ea844a1aSMatthew Knepley     }
1705ea844a1aSMatthew Knepley     ierr = PetscFree(indices);CHKERRQ(ierr);
1706ea844a1aSMatthew Knepley   }
1707ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1708ea844a1aSMatthew Knepley }
1709ea844a1aSMatthew Knepley 
1710ea844a1aSMatthew Knepley /*@
1711ea844a1aSMatthew Knepley   PetscSectionCreateSupersection - Create a new, larger section composed of the input sections
1712ea844a1aSMatthew Knepley 
1713ea844a1aSMatthew Knepley   Collective on s
1714ea844a1aSMatthew Knepley 
1715ea844a1aSMatthew Knepley   Input Parameters:
1716ea844a1aSMatthew Knepley + s     - the input sections
1717ea844a1aSMatthew Knepley - len   - the number of input sections
1718ea844a1aSMatthew Knepley 
1719ea844a1aSMatthew Knepley   Output Parameter:
1720ea844a1aSMatthew Knepley . supers - the supersection
1721ea844a1aSMatthew Knepley 
1722ea844a1aSMatthew Knepley   Note: The section offsets now refer to a new, larger vector.
1723ea844a1aSMatthew Knepley 
1724ea844a1aSMatthew Knepley   Level: advanced
1725ea844a1aSMatthew Knepley 
1726ea844a1aSMatthew Knepley .seealso: PetscSectionCreateSubsection(), PetscSectionCreate()
1727ea844a1aSMatthew Knepley @*/
1728ea844a1aSMatthew Knepley PetscErrorCode PetscSectionCreateSupersection(PetscSection s[], PetscInt len, PetscSection *supers)
1729ea844a1aSMatthew Knepley {
1730ea844a1aSMatthew Knepley   PetscInt       Nf = 0, f, pStart = PETSC_MAX_INT, pEnd = 0, p, maxCdof = 0, i;
1731ea844a1aSMatthew Knepley   PetscErrorCode ierr;
1732ea844a1aSMatthew Knepley 
1733ea844a1aSMatthew Knepley   PetscFunctionBegin;
1734ea844a1aSMatthew Knepley   if (!len) PetscFunctionReturn(0);
1735ea844a1aSMatthew Knepley   for (i = 0; i < len; ++i) {
1736ea844a1aSMatthew Knepley     PetscInt nf, pStarti, pEndi;
1737ea844a1aSMatthew Knepley 
1738ea844a1aSMatthew Knepley     ierr = PetscSectionGetNumFields(s[i], &nf);CHKERRQ(ierr);
1739ea844a1aSMatthew Knepley     ierr = PetscSectionGetChart(s[i], &pStarti, &pEndi);CHKERRQ(ierr);
1740ea844a1aSMatthew Knepley     pStart = PetscMin(pStart, pStarti);
1741ea844a1aSMatthew Knepley     pEnd   = PetscMax(pEnd,   pEndi);
1742ea844a1aSMatthew Knepley     Nf += nf;
1743ea844a1aSMatthew Knepley   }
1744ea844a1aSMatthew Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s[0]), supers);CHKERRQ(ierr);
1745ea844a1aSMatthew Knepley   ierr = PetscSectionSetNumFields(*supers, Nf);CHKERRQ(ierr);
1746ea844a1aSMatthew Knepley   for (i = 0, f = 0; i < len; ++i) {
1747ea844a1aSMatthew Knepley     PetscInt nf, fi;
1748ea844a1aSMatthew Knepley 
1749ea844a1aSMatthew Knepley     ierr = PetscSectionGetNumFields(s[i], &nf);CHKERRQ(ierr);
1750ea844a1aSMatthew Knepley     for (fi = 0; fi < nf; ++fi, ++f) {
1751ea844a1aSMatthew Knepley       const char *name   = NULL;
1752ea844a1aSMatthew Knepley       PetscInt   numComp = 0;
1753ea844a1aSMatthew Knepley 
1754ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldName(s[i], fi, &name);CHKERRQ(ierr);
1755ea844a1aSMatthew Knepley       ierr = PetscSectionSetFieldName(*supers, f, name);CHKERRQ(ierr);
1756ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldComponents(s[i], fi, &numComp);CHKERRQ(ierr);
1757ea844a1aSMatthew Knepley       ierr = PetscSectionSetFieldComponents(*supers, f, numComp);CHKERRQ(ierr);
1758ea844a1aSMatthew Knepley     }
1759ea844a1aSMatthew Knepley   }
1760ea844a1aSMatthew Knepley   ierr = PetscSectionSetChart(*supers, pStart, pEnd);CHKERRQ(ierr);
1761ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
1762ea844a1aSMatthew Knepley     PetscInt dof = 0, cdof = 0;
1763ea844a1aSMatthew Knepley 
1764ea844a1aSMatthew Knepley     for (i = 0, f = 0; i < len; ++i) {
1765ea844a1aSMatthew Knepley       PetscInt nf, fi, pStarti, pEndi;
1766ea844a1aSMatthew Knepley       PetscInt fdof = 0, cfdof = 0;
1767ea844a1aSMatthew Knepley 
1768ea844a1aSMatthew Knepley       ierr = PetscSectionGetNumFields(s[i], &nf);CHKERRQ(ierr);
1769ea844a1aSMatthew Knepley       ierr = PetscSectionGetChart(s[i], &pStarti, &pEndi);CHKERRQ(ierr);
1770ea844a1aSMatthew Knepley       if ((p < pStarti) || (p >= pEndi)) continue;
1771ea844a1aSMatthew Knepley       for (fi = 0; fi < nf; ++fi, ++f) {
1772ea844a1aSMatthew Knepley         ierr = PetscSectionGetFieldDof(s[i], p, fi, &fdof);CHKERRQ(ierr);
1773ea844a1aSMatthew Knepley         ierr = PetscSectionAddFieldDof(*supers, p, f, fdof);CHKERRQ(ierr);
1774ea844a1aSMatthew Knepley         ierr = PetscSectionGetFieldConstraintDof(s[i], p, fi, &cfdof);CHKERRQ(ierr);
1775ea844a1aSMatthew Knepley         if (cfdof) {ierr = PetscSectionAddFieldConstraintDof(*supers, p, f, cfdof);CHKERRQ(ierr);}
1776ea844a1aSMatthew Knepley         dof  += fdof;
1777ea844a1aSMatthew Knepley         cdof += cfdof;
1778ea844a1aSMatthew Knepley       }
1779ea844a1aSMatthew Knepley     }
1780ea844a1aSMatthew Knepley     ierr = PetscSectionSetDof(*supers, p, dof);CHKERRQ(ierr);
1781ea844a1aSMatthew Knepley     if (cdof) {ierr = PetscSectionSetConstraintDof(*supers, p, cdof);CHKERRQ(ierr);}
1782ea844a1aSMatthew Knepley     maxCdof = PetscMax(cdof, maxCdof);
1783ea844a1aSMatthew Knepley   }
1784ea844a1aSMatthew Knepley   ierr = PetscSectionSetUp(*supers);CHKERRQ(ierr);
1785ea844a1aSMatthew Knepley   if (maxCdof) {
1786ea844a1aSMatthew Knepley     PetscInt *indices;
1787ea844a1aSMatthew Knepley 
1788ea844a1aSMatthew Knepley     ierr = PetscMalloc1(maxCdof, &indices);CHKERRQ(ierr);
1789ea844a1aSMatthew Knepley     for (p = pStart; p < pEnd; ++p) {
1790ea844a1aSMatthew Knepley       PetscInt cdof;
1791ea844a1aSMatthew Knepley 
1792ea844a1aSMatthew Knepley       ierr = PetscSectionGetConstraintDof(*supers, p, &cdof);CHKERRQ(ierr);
1793ea844a1aSMatthew Knepley       if (cdof) {
1794ea844a1aSMatthew Knepley         PetscInt dof, numConst = 0, fOff = 0;
1795ea844a1aSMatthew Knepley 
1796ea844a1aSMatthew Knepley         for (i = 0, f = 0; i < len; ++i) {
1797ea844a1aSMatthew Knepley           const PetscInt *oldIndices = NULL;
1798ea844a1aSMatthew Knepley           PetscInt        nf, fi, pStarti, pEndi, fdof, cfdof, fc;
1799ea844a1aSMatthew Knepley 
1800ea844a1aSMatthew Knepley           ierr = PetscSectionGetNumFields(s[i], &nf);CHKERRQ(ierr);
1801ea844a1aSMatthew Knepley           ierr = PetscSectionGetChart(s[i], &pStarti, &pEndi);CHKERRQ(ierr);
1802ea844a1aSMatthew Knepley           if ((p < pStarti) || (p >= pEndi)) continue;
1803ea844a1aSMatthew Knepley           for (fi = 0; fi < nf; ++fi, ++f) {
1804ea844a1aSMatthew Knepley             ierr = PetscSectionGetFieldDof(s[i], p, fi, &fdof);CHKERRQ(ierr);
1805ea844a1aSMatthew Knepley             ierr = PetscSectionGetFieldConstraintDof(s[i], p, fi, &cfdof);CHKERRQ(ierr);
1806ea844a1aSMatthew Knepley             ierr = PetscSectionGetFieldConstraintIndices(s[i], p, fi, &oldIndices);CHKERRQ(ierr);
1807ea844a1aSMatthew Knepley             for (fc = 0; fc < cfdof; ++fc) indices[numConst+fc] = oldIndices[fc] + fOff;
1808ea844a1aSMatthew Knepley             ierr = PetscSectionSetFieldConstraintIndices(*supers, p, f, &indices[numConst]);CHKERRQ(ierr);
1809ea844a1aSMatthew Knepley             numConst += cfdof;
1810ea844a1aSMatthew Knepley           }
1811ea844a1aSMatthew Knepley           ierr = PetscSectionGetDof(s[i], p, &dof);CHKERRQ(ierr);
1812ea844a1aSMatthew Knepley           fOff += dof;
1813ea844a1aSMatthew Knepley         }
1814ea844a1aSMatthew Knepley         ierr = PetscSectionSetConstraintIndices(*supers, p, indices);CHKERRQ(ierr);
1815ea844a1aSMatthew Knepley       }
1816ea844a1aSMatthew Knepley     }
1817ea844a1aSMatthew Knepley     ierr = PetscFree(indices);CHKERRQ(ierr);
1818ea844a1aSMatthew Knepley   }
1819ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1820ea844a1aSMatthew Knepley }
1821ea844a1aSMatthew Knepley 
1822ea844a1aSMatthew Knepley /*@
1823ea844a1aSMatthew Knepley   PetscSectionCreateSubmeshSection - Create a new, smaller section with support on the submesh
1824ea844a1aSMatthew Knepley 
1825ea844a1aSMatthew Knepley   Collective on s
1826ea844a1aSMatthew Knepley 
1827ea844a1aSMatthew Knepley   Input Parameters:
1828ea844a1aSMatthew Knepley + s           - the PetscSection
1829ea844a1aSMatthew Knepley - subpointMap - a sorted list of points in the original mesh which are in the submesh
1830ea844a1aSMatthew Knepley 
1831ea844a1aSMatthew Knepley   Output Parameter:
1832ea844a1aSMatthew Knepley . subs - the subsection
1833ea844a1aSMatthew Knepley 
1834ea844a1aSMatthew Knepley   Note: The section offsets now refer to a new, smaller vector.
1835ea844a1aSMatthew Knepley 
1836ea844a1aSMatthew Knepley   Level: advanced
1837ea844a1aSMatthew Knepley 
1838ea844a1aSMatthew Knepley .seealso: PetscSectionCreateSubsection(), DMPlexGetSubpointMap(), PetscSectionCreate()
1839ea844a1aSMatthew Knepley @*/
1840ea844a1aSMatthew Knepley PetscErrorCode PetscSectionCreateSubmeshSection(PetscSection s, IS subpointMap, PetscSection *subs)
1841ea844a1aSMatthew Knepley {
1842ea844a1aSMatthew Knepley   const PetscInt *points = NULL, *indices = NULL;
1843ea844a1aSMatthew Knepley   PetscInt       numFields, f, numSubpoints = 0, pStart, pEnd, p, subp;
1844ea844a1aSMatthew Knepley   PetscErrorCode ierr;
1845ea844a1aSMatthew Knepley 
1846ea844a1aSMatthew Knepley   PetscFunctionBegin;
1847ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1848ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(subpointMap, IS_CLASSID, 2);
1849ea844a1aSMatthew Knepley   PetscValidPointer(subs, 3);
1850ea844a1aSMatthew Knepley   ierr = PetscSectionGetNumFields(s, &numFields);CHKERRQ(ierr);
1851ea844a1aSMatthew Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), subs);CHKERRQ(ierr);
1852ea844a1aSMatthew Knepley   if (numFields) {ierr = PetscSectionSetNumFields(*subs, numFields);CHKERRQ(ierr);}
1853ea844a1aSMatthew Knepley   for (f = 0; f < numFields; ++f) {
1854ea844a1aSMatthew Knepley     const char *name   = NULL;
1855ea844a1aSMatthew Knepley     PetscInt   numComp = 0;
1856ea844a1aSMatthew Knepley 
1857ea844a1aSMatthew Knepley     ierr = PetscSectionGetFieldName(s, f, &name);CHKERRQ(ierr);
1858ea844a1aSMatthew Knepley     ierr = PetscSectionSetFieldName(*subs, f, name);CHKERRQ(ierr);
1859ea844a1aSMatthew Knepley     ierr = PetscSectionGetFieldComponents(s, f, &numComp);CHKERRQ(ierr);
1860ea844a1aSMatthew Knepley     ierr = PetscSectionSetFieldComponents(*subs, f, numComp);CHKERRQ(ierr);
1861ea844a1aSMatthew Knepley   }
1862ea844a1aSMatthew Knepley   /* For right now, we do not try to squeeze the subchart */
1863ea844a1aSMatthew Knepley   if (subpointMap) {
1864ea844a1aSMatthew Knepley     ierr = ISGetSize(subpointMap, &numSubpoints);CHKERRQ(ierr);
1865ea844a1aSMatthew Knepley     ierr = ISGetIndices(subpointMap, &points);CHKERRQ(ierr);
1866ea844a1aSMatthew Knepley   }
1867ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1868ea844a1aSMatthew Knepley   ierr = PetscSectionSetChart(*subs, 0, numSubpoints);CHKERRQ(ierr);
1869ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
1870ea844a1aSMatthew Knepley     PetscInt dof, cdof, fdof = 0, cfdof = 0;
1871ea844a1aSMatthew Knepley 
1872ea844a1aSMatthew Knepley     ierr = PetscFindInt(p, numSubpoints, points, &subp);CHKERRQ(ierr);
1873ea844a1aSMatthew Knepley     if (subp < 0) continue;
1874ea844a1aSMatthew Knepley     for (f = 0; f < numFields; ++f) {
1875ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
1876ea844a1aSMatthew Knepley       ierr = PetscSectionSetFieldDof(*subs, subp, f, fdof);CHKERRQ(ierr);
1877ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldConstraintDof(s, p, f, &cfdof);CHKERRQ(ierr);
1878ea844a1aSMatthew Knepley       if (cfdof) {ierr = PetscSectionSetFieldConstraintDof(*subs, subp, f, cfdof);CHKERRQ(ierr);}
1879ea844a1aSMatthew Knepley     }
1880ea844a1aSMatthew Knepley     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
1881ea844a1aSMatthew Knepley     ierr = PetscSectionSetDof(*subs, subp, dof);CHKERRQ(ierr);
1882ea844a1aSMatthew Knepley     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1883ea844a1aSMatthew Knepley     if (cdof) {ierr = PetscSectionSetConstraintDof(*subs, subp, cdof);CHKERRQ(ierr);}
1884ea844a1aSMatthew Knepley   }
1885ea844a1aSMatthew Knepley   ierr = PetscSectionSetUp(*subs);CHKERRQ(ierr);
1886ea844a1aSMatthew Knepley   /* Change offsets to original offsets */
1887ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
1888ea844a1aSMatthew Knepley     PetscInt off, foff = 0;
1889ea844a1aSMatthew Knepley 
1890ea844a1aSMatthew Knepley     ierr = PetscFindInt(p, numSubpoints, points, &subp);CHKERRQ(ierr);
1891ea844a1aSMatthew Knepley     if (subp < 0) continue;
1892ea844a1aSMatthew Knepley     for (f = 0; f < numFields; ++f) {
1893ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr);
1894ea844a1aSMatthew Knepley       ierr = PetscSectionSetFieldOffset(*subs, subp, f, foff);CHKERRQ(ierr);
1895ea844a1aSMatthew Knepley     }
1896ea844a1aSMatthew Knepley     ierr = PetscSectionGetOffset(s, p, &off);CHKERRQ(ierr);
1897ea844a1aSMatthew Knepley     ierr = PetscSectionSetOffset(*subs, subp, off);CHKERRQ(ierr);
1898ea844a1aSMatthew Knepley   }
1899ea844a1aSMatthew Knepley   /* Copy constraint indices */
1900ea844a1aSMatthew Knepley   for (subp = 0; subp < numSubpoints; ++subp) {
1901ea844a1aSMatthew Knepley     PetscInt cdof;
1902ea844a1aSMatthew Knepley 
1903ea844a1aSMatthew Knepley     ierr = PetscSectionGetConstraintDof(*subs, subp, &cdof);CHKERRQ(ierr);
1904ea844a1aSMatthew Knepley     if (cdof) {
1905ea844a1aSMatthew Knepley       for (f = 0; f < numFields; ++f) {
1906ea844a1aSMatthew Knepley         ierr = PetscSectionGetFieldConstraintIndices(s, points[subp], f, &indices);CHKERRQ(ierr);
1907ea844a1aSMatthew Knepley         ierr = PetscSectionSetFieldConstraintIndices(*subs, subp, f, indices);CHKERRQ(ierr);
1908ea844a1aSMatthew Knepley       }
1909ea844a1aSMatthew Knepley       ierr = PetscSectionGetConstraintIndices(s, points[subp], &indices);CHKERRQ(ierr);
1910ea844a1aSMatthew Knepley       ierr = PetscSectionSetConstraintIndices(*subs, subp, indices);CHKERRQ(ierr);
1911ea844a1aSMatthew Knepley     }
1912ea844a1aSMatthew Knepley   }
1913ea844a1aSMatthew Knepley   if (subpointMap) {ierr = ISRestoreIndices(subpointMap, &points);CHKERRQ(ierr);}
1914ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1915ea844a1aSMatthew Knepley }
1916ea844a1aSMatthew Knepley 
1917ea844a1aSMatthew Knepley static PetscErrorCode PetscSectionView_ASCII(PetscSection s, PetscViewer viewer)
1918ea844a1aSMatthew Knepley {
1919ea844a1aSMatthew Knepley   PetscInt       p;
1920ea844a1aSMatthew Knepley   PetscMPIInt    rank;
1921ea844a1aSMatthew Knepley   PetscErrorCode ierr;
1922ea844a1aSMatthew Knepley 
1923ea844a1aSMatthew Knepley   PetscFunctionBegin;
1924ea844a1aSMatthew Knepley   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);CHKERRQ(ierr);
1925ea844a1aSMatthew Knepley   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
1926ea844a1aSMatthew Knepley   ierr = PetscViewerASCIISynchronizedPrintf(viewer, "Process %d:\n", rank);CHKERRQ(ierr);
1927ea844a1aSMatthew Knepley   for (p = 0; p < s->pEnd - s->pStart; ++p) {
1928ea844a1aSMatthew Knepley     if ((s->bc) && (s->bc->atlasDof[p] > 0)) {
1929ea844a1aSMatthew Knepley       PetscInt b;
1930ea844a1aSMatthew Knepley 
1931ea844a1aSMatthew Knepley       ierr = PetscViewerASCIISynchronizedPrintf(viewer, "  (%4D) dim %2D offset %3D constrained", p+s->pStart, s->atlasDof[p], s->atlasOff[p]);CHKERRQ(ierr);
1932ea844a1aSMatthew Knepley       for (b = 0; b < s->bc->atlasDof[p]; ++b) {
1933ea844a1aSMatthew Knepley         ierr = PetscViewerASCIISynchronizedPrintf(viewer, " %D", s->bcIndices[s->bc->atlasOff[p]+b]);CHKERRQ(ierr);
1934ea844a1aSMatthew Knepley       }
1935ea844a1aSMatthew Knepley       ierr = PetscViewerASCIISynchronizedPrintf(viewer, "\n");CHKERRQ(ierr);
1936ea844a1aSMatthew Knepley     } else {
1937ea844a1aSMatthew Knepley       ierr = PetscViewerASCIISynchronizedPrintf(viewer, "  (%4D) dim %2D offset %3D\n", p+s->pStart, s->atlasDof[p], s->atlasOff[p]);CHKERRQ(ierr);
1938ea844a1aSMatthew Knepley     }
1939ea844a1aSMatthew Knepley   }
1940ea844a1aSMatthew Knepley   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1941ea844a1aSMatthew Knepley   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
1942ea844a1aSMatthew Knepley   if (s->sym) {
1943ea844a1aSMatthew Knepley     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1944ea844a1aSMatthew Knepley     ierr = PetscSectionSymView(s->sym,viewer);CHKERRQ(ierr);
1945ea844a1aSMatthew Knepley     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1946ea844a1aSMatthew Knepley   }
1947ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1948ea844a1aSMatthew Knepley }
1949ea844a1aSMatthew Knepley 
1950ea844a1aSMatthew Knepley /*@C
1951*fe2efc57SMark    PetscSectionViewFromOptions - View from Options
1952*fe2efc57SMark 
1953*fe2efc57SMark    Collective on PetscSection
1954*fe2efc57SMark 
1955*fe2efc57SMark    Input Parameters:
1956*fe2efc57SMark +  A - the PetscSection object to view
1957*fe2efc57SMark -  obj - Optional object
1958*fe2efc57SMark .  name - command line option
1959*fe2efc57SMark 
1960*fe2efc57SMark    Level: intermediate
1961*fe2efc57SMark .seealso:  PetscSection, PetscSectionView, PetscObjectViewFromOptions(), PetscSectionCreate()
1962*fe2efc57SMark @*/
1963*fe2efc57SMark PetscErrorCode  PetscSectionViewFromOptions(PetscSection A,PetscObject obj,const char name[])
1964*fe2efc57SMark {
1965*fe2efc57SMark   PetscErrorCode ierr;
1966*fe2efc57SMark 
1967*fe2efc57SMark   PetscFunctionBegin;
1968*fe2efc57SMark   PetscValidHeaderSpecific(A,PETSC_SECTION_CLASSID,1);
1969*fe2efc57SMark   ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr);
1970*fe2efc57SMark   PetscFunctionReturn(0);
1971*fe2efc57SMark }
1972*fe2efc57SMark 
1973*fe2efc57SMark /*@C
1974ea844a1aSMatthew Knepley   PetscSectionView - Views a PetscSection
1975ea844a1aSMatthew Knepley 
1976ea844a1aSMatthew Knepley   Collective on PetscSection
1977ea844a1aSMatthew Knepley 
1978ea844a1aSMatthew Knepley   Input Parameters:
1979ea844a1aSMatthew Knepley + s - the PetscSection object to view
1980ea844a1aSMatthew Knepley - v - the viewer
1981ea844a1aSMatthew Knepley 
1982ea844a1aSMatthew Knepley   Level: beginner
1983ea844a1aSMatthew Knepley 
1984ea844a1aSMatthew Knepley .seealso PetscSectionCreate(), PetscSectionDestroy()
1985ea844a1aSMatthew Knepley @*/
1986ea844a1aSMatthew Knepley PetscErrorCode PetscSectionView(PetscSection s, PetscViewer viewer)
1987ea844a1aSMatthew Knepley {
1988ea844a1aSMatthew Knepley   PetscBool      isascii;
1989ea844a1aSMatthew Knepley   PetscInt       f;
1990ea844a1aSMatthew Knepley   PetscErrorCode ierr;
1991ea844a1aSMatthew Knepley 
1992ea844a1aSMatthew Knepley   PetscFunctionBegin;
1993ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1994ea844a1aSMatthew Knepley   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)s), &viewer);CHKERRQ(ierr);}
1995ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1996ea844a1aSMatthew Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
1997ea844a1aSMatthew Knepley   if (isascii) {
1998ea844a1aSMatthew Knepley     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)s,viewer);CHKERRQ(ierr);
1999ea844a1aSMatthew Knepley     if (s->numFields) {
2000ea844a1aSMatthew Knepley       ierr = PetscViewerASCIIPrintf(viewer, "%D fields\n", s->numFields);CHKERRQ(ierr);
2001ea844a1aSMatthew Knepley       for (f = 0; f < s->numFields; ++f) {
2002ea844a1aSMatthew Knepley         ierr = PetscViewerASCIIPrintf(viewer, "  field %D with %D components\n", f, s->numFieldComponents[f]);CHKERRQ(ierr);
2003ea844a1aSMatthew Knepley         ierr = PetscSectionView_ASCII(s->field[f], viewer);CHKERRQ(ierr);
2004ea844a1aSMatthew Knepley       }
2005ea844a1aSMatthew Knepley     } else {
2006ea844a1aSMatthew Knepley       ierr = PetscSectionView_ASCII(s, viewer);CHKERRQ(ierr);
2007ea844a1aSMatthew Knepley     }
2008ea844a1aSMatthew Knepley   }
2009ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2010ea844a1aSMatthew Knepley }
2011ea844a1aSMatthew Knepley 
2012ea844a1aSMatthew Knepley /*@
2013ea844a1aSMatthew Knepley   PetscSectionReset - Frees all section data.
2014ea844a1aSMatthew Knepley 
2015ea844a1aSMatthew Knepley   Not collective
2016ea844a1aSMatthew Knepley 
2017ea844a1aSMatthew Knepley   Input Parameters:
2018ea844a1aSMatthew Knepley . s - the PetscSection
2019ea844a1aSMatthew Knepley 
2020ea844a1aSMatthew Knepley   Level: beginner
2021ea844a1aSMatthew Knepley 
2022ea844a1aSMatthew Knepley .seealso: PetscSection, PetscSectionCreate()
2023ea844a1aSMatthew Knepley @*/
2024ea844a1aSMatthew Knepley PetscErrorCode PetscSectionReset(PetscSection s)
2025ea844a1aSMatthew Knepley {
2026ea844a1aSMatthew Knepley   PetscInt       f;
2027ea844a1aSMatthew Knepley   PetscErrorCode ierr;
2028ea844a1aSMatthew Knepley 
2029ea844a1aSMatthew Knepley   PetscFunctionBegin;
2030ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
2031ea844a1aSMatthew Knepley   ierr = PetscFree(s->numFieldComponents);CHKERRQ(ierr);
2032ea844a1aSMatthew Knepley   for (f = 0; f < s->numFields; ++f) {
2033ea844a1aSMatthew Knepley     ierr = PetscSectionDestroy(&s->field[f]);CHKERRQ(ierr);
2034ea844a1aSMatthew Knepley     ierr = PetscFree(s->fieldNames[f]);CHKERRQ(ierr);
2035ea844a1aSMatthew Knepley   }
2036ea844a1aSMatthew Knepley   ierr = PetscFree(s->fieldNames);CHKERRQ(ierr);
2037ea844a1aSMatthew Knepley   ierr = PetscFree(s->field);CHKERRQ(ierr);
2038ea844a1aSMatthew Knepley   ierr = PetscSectionDestroy(&s->bc);CHKERRQ(ierr);
2039ea844a1aSMatthew Knepley   ierr = PetscFree(s->bcIndices);CHKERRQ(ierr);
2040ea844a1aSMatthew Knepley   ierr = PetscFree2(s->atlasDof, s->atlasOff);CHKERRQ(ierr);
2041ea844a1aSMatthew Knepley   ierr = PetscSectionDestroy(&s->clSection);CHKERRQ(ierr);
2042ea844a1aSMatthew Knepley   ierr = ISDestroy(&s->clPoints);CHKERRQ(ierr);
2043ea844a1aSMatthew Knepley   ierr = ISDestroy(&s->perm);CHKERRQ(ierr);
2044ea844a1aSMatthew Knepley   ierr = PetscFree(s->clPerm);CHKERRQ(ierr);
2045ea844a1aSMatthew Knepley   ierr = PetscFree(s->clInvPerm);CHKERRQ(ierr);
2046ea844a1aSMatthew Knepley   ierr = PetscSectionSymDestroy(&s->sym);CHKERRQ(ierr);
2047ea844a1aSMatthew Knepley 
2048ea844a1aSMatthew Knepley   s->pStart    = -1;
2049ea844a1aSMatthew Knepley   s->pEnd      = -1;
2050ea844a1aSMatthew Knepley   s->maxDof    = 0;
2051ea844a1aSMatthew Knepley   s->setup     = PETSC_FALSE;
2052ea844a1aSMatthew Knepley   s->numFields = 0;
2053ea844a1aSMatthew Knepley   s->clObj     = NULL;
2054ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2055ea844a1aSMatthew Knepley }
2056ea844a1aSMatthew Knepley 
2057ea844a1aSMatthew Knepley /*@
2058ea844a1aSMatthew Knepley   PetscSectionDestroy - Frees a section object and frees its range if that exists.
2059ea844a1aSMatthew Knepley 
2060ea844a1aSMatthew Knepley   Not collective
2061ea844a1aSMatthew Knepley 
2062ea844a1aSMatthew Knepley   Input Parameters:
2063ea844a1aSMatthew Knepley . s - the PetscSection
2064ea844a1aSMatthew Knepley 
2065ea844a1aSMatthew Knepley   Level: beginner
2066ea844a1aSMatthew Knepley 
2067ea844a1aSMatthew Knepley .seealso: PetscSection, PetscSectionCreate()
2068ea844a1aSMatthew Knepley @*/
2069ea844a1aSMatthew Knepley PetscErrorCode PetscSectionDestroy(PetscSection *s)
2070ea844a1aSMatthew Knepley {
2071ea844a1aSMatthew Knepley   PetscErrorCode ierr;
2072ea844a1aSMatthew Knepley 
2073ea844a1aSMatthew Knepley   PetscFunctionBegin;
2074ea844a1aSMatthew Knepley   if (!*s) PetscFunctionReturn(0);
207562f17a49SStefano Zampini   PetscValidHeaderSpecific(*s,PETSC_SECTION_CLASSID,1);
2076ea844a1aSMatthew Knepley   if (--((PetscObject)(*s))->refct > 0) {
2077ea844a1aSMatthew Knepley     *s = NULL;
2078ea844a1aSMatthew Knepley     PetscFunctionReturn(0);
2079ea844a1aSMatthew Knepley   }
2080ea844a1aSMatthew Knepley   ierr = PetscSectionReset(*s);CHKERRQ(ierr);
2081ea844a1aSMatthew Knepley   ierr = PetscHeaderDestroy(s);CHKERRQ(ierr);
2082ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2083ea844a1aSMatthew Knepley }
2084ea844a1aSMatthew Knepley 
2085ea844a1aSMatthew Knepley PetscErrorCode VecIntGetValuesSection(PetscInt *baseArray, PetscSection s, PetscInt point, const PetscInt **values)
2086ea844a1aSMatthew Knepley {
2087ea844a1aSMatthew Knepley   const PetscInt p = point - s->pStart;
2088ea844a1aSMatthew Knepley 
2089ea844a1aSMatthew Knepley   PetscFunctionBegin;
2090ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 2);
2091ea844a1aSMatthew Knepley   *values = &baseArray[s->atlasOff[p]];
2092ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2093ea844a1aSMatthew Knepley }
2094ea844a1aSMatthew Knepley 
2095ea844a1aSMatthew Knepley PetscErrorCode VecIntSetValuesSection(PetscInt *baseArray, PetscSection s, PetscInt point, const PetscInt values[], InsertMode mode)
2096ea844a1aSMatthew Knepley {
2097ea844a1aSMatthew Knepley   PetscInt       *array;
2098ea844a1aSMatthew Knepley   const PetscInt p           = point - s->pStart;
2099ea844a1aSMatthew Knepley   const PetscInt orientation = 0; /* Needs to be included for use in closure operations */
2100ea844a1aSMatthew Knepley   PetscInt       cDim        = 0;
2101ea844a1aSMatthew Knepley   PetscErrorCode ierr;
2102ea844a1aSMatthew Knepley 
2103ea844a1aSMatthew Knepley   PetscFunctionBegin;
2104ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 2);
2105ea844a1aSMatthew Knepley   ierr  = PetscSectionGetConstraintDof(s, p, &cDim);CHKERRQ(ierr);
2106ea844a1aSMatthew Knepley   array = &baseArray[s->atlasOff[p]];
2107ea844a1aSMatthew Knepley   if (!cDim) {
2108ea844a1aSMatthew Knepley     if (orientation >= 0) {
2109ea844a1aSMatthew Knepley       const PetscInt dim = s->atlasDof[p];
2110ea844a1aSMatthew Knepley       PetscInt       i;
2111ea844a1aSMatthew Knepley 
2112ea844a1aSMatthew Knepley       if (mode == INSERT_VALUES) {
2113ea844a1aSMatthew Knepley         for (i = 0; i < dim; ++i) array[i] = values[i];
2114ea844a1aSMatthew Knepley       } else {
2115ea844a1aSMatthew Knepley         for (i = 0; i < dim; ++i) array[i] += values[i];
2116ea844a1aSMatthew Knepley       }
2117ea844a1aSMatthew Knepley     } else {
2118ea844a1aSMatthew Knepley       PetscInt offset = 0;
2119ea844a1aSMatthew Knepley       PetscInt j      = -1, field, i;
2120ea844a1aSMatthew Knepley 
2121ea844a1aSMatthew Knepley       for (field = 0; field < s->numFields; ++field) {
2122ea844a1aSMatthew Knepley         const PetscInt dim = s->field[field]->atlasDof[p];
2123ea844a1aSMatthew Knepley 
2124ea844a1aSMatthew Knepley         for (i = dim-1; i >= 0; --i) array[++j] = values[i+offset];
2125ea844a1aSMatthew Knepley         offset += dim;
2126ea844a1aSMatthew Knepley       }
2127ea844a1aSMatthew Knepley     }
2128ea844a1aSMatthew Knepley   } else {
2129ea844a1aSMatthew Knepley     if (orientation >= 0) {
2130ea844a1aSMatthew Knepley       const PetscInt dim  = s->atlasDof[p];
2131ea844a1aSMatthew Knepley       PetscInt       cInd = 0, i;
2132ea844a1aSMatthew Knepley       const PetscInt *cDof;
2133ea844a1aSMatthew Knepley 
2134ea844a1aSMatthew Knepley       ierr = PetscSectionGetConstraintIndices(s, point, &cDof);CHKERRQ(ierr);
2135ea844a1aSMatthew Knepley       if (mode == INSERT_VALUES) {
2136ea844a1aSMatthew Knepley         for (i = 0; i < dim; ++i) {
2137ea844a1aSMatthew Knepley           if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;}
2138ea844a1aSMatthew Knepley           array[i] = values[i];
2139ea844a1aSMatthew Knepley         }
2140ea844a1aSMatthew Knepley       } else {
2141ea844a1aSMatthew Knepley         for (i = 0; i < dim; ++i) {
2142ea844a1aSMatthew Knepley           if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;}
2143ea844a1aSMatthew Knepley           array[i] += values[i];
2144ea844a1aSMatthew Knepley         }
2145ea844a1aSMatthew Knepley       }
2146ea844a1aSMatthew Knepley     } else {
2147ea844a1aSMatthew Knepley       const PetscInt *cDof;
2148ea844a1aSMatthew Knepley       PetscInt       offset  = 0;
2149ea844a1aSMatthew Knepley       PetscInt       cOffset = 0;
2150ea844a1aSMatthew Knepley       PetscInt       j       = 0, field;
2151ea844a1aSMatthew Knepley 
2152ea844a1aSMatthew Knepley       ierr = PetscSectionGetConstraintIndices(s, point, &cDof);CHKERRQ(ierr);
2153ea844a1aSMatthew Knepley       for (field = 0; field < s->numFields; ++field) {
2154ea844a1aSMatthew Knepley         const PetscInt dim  = s->field[field]->atlasDof[p];     /* PetscSectionGetFieldDof() */
2155ea844a1aSMatthew Knepley         const PetscInt tDim = s->field[field]->bc->atlasDof[p]; /* PetscSectionGetFieldConstraintDof() */
2156ea844a1aSMatthew Knepley         const PetscInt sDim = dim - tDim;
2157ea844a1aSMatthew Knepley         PetscInt       cInd = 0, i,k;
2158ea844a1aSMatthew Knepley 
2159ea844a1aSMatthew Knepley         for (i = 0, k = dim+offset-1; i < dim; ++i, ++j, --k) {
2160ea844a1aSMatthew Knepley           if ((cInd < sDim) && (j == cDof[cInd+cOffset])) {++cInd; continue;}
2161ea844a1aSMatthew Knepley           array[j] = values[k];
2162ea844a1aSMatthew Knepley         }
2163ea844a1aSMatthew Knepley         offset  += dim;
2164ea844a1aSMatthew Knepley         cOffset += dim - tDim;
2165ea844a1aSMatthew Knepley       }
2166ea844a1aSMatthew Knepley     }
2167ea844a1aSMatthew Knepley   }
2168ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2169ea844a1aSMatthew Knepley }
2170ea844a1aSMatthew Knepley 
2171ea844a1aSMatthew Knepley /*@C
2172ea844a1aSMatthew Knepley   PetscSectionHasConstraints - Determine whether a section has constrained dofs
2173ea844a1aSMatthew Knepley 
2174ea844a1aSMatthew Knepley   Not collective
2175ea844a1aSMatthew Knepley 
2176ea844a1aSMatthew Knepley   Input Parameter:
2177ea844a1aSMatthew Knepley . s - The PetscSection
2178ea844a1aSMatthew Knepley 
2179ea844a1aSMatthew Knepley   Output Parameter:
2180ea844a1aSMatthew Knepley . hasConstraints - flag indicating that the section has constrained dofs
2181ea844a1aSMatthew Knepley 
2182ea844a1aSMatthew Knepley   Level: intermediate
2183ea844a1aSMatthew Knepley 
2184ea844a1aSMatthew Knepley .seealso: PetscSectionSetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection
2185ea844a1aSMatthew Knepley @*/
2186ea844a1aSMatthew Knepley PetscErrorCode PetscSectionHasConstraints(PetscSection s, PetscBool *hasConstraints)
2187ea844a1aSMatthew Knepley {
2188ea844a1aSMatthew Knepley   PetscFunctionBegin;
2189ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
2190ea844a1aSMatthew Knepley   PetscValidPointer(hasConstraints, 2);
2191ea844a1aSMatthew Knepley   *hasConstraints = s->bc ? PETSC_TRUE : PETSC_FALSE;
2192ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2193ea844a1aSMatthew Knepley }
2194ea844a1aSMatthew Knepley 
2195ea844a1aSMatthew Knepley /*@C
2196ea844a1aSMatthew Knepley   PetscSectionGetConstraintIndices - Get the point dof numbers, in [0, dof), which are constrained
2197ea844a1aSMatthew Knepley 
2198ea844a1aSMatthew Knepley   Not collective
2199ea844a1aSMatthew Knepley 
2200ea844a1aSMatthew Knepley   Input Parameters:
2201ea844a1aSMatthew Knepley + s     - The PetscSection
2202ea844a1aSMatthew Knepley - point - The point
2203ea844a1aSMatthew Knepley 
2204ea844a1aSMatthew Knepley   Output Parameter:
2205ea844a1aSMatthew Knepley . indices - The constrained dofs
2206ea844a1aSMatthew Knepley 
2207ea844a1aSMatthew Knepley   Note: In Fortran, you call PetscSectionGetConstraintIndicesF90() and PetscSectionRestoreConstraintIndicesF90()
2208ea844a1aSMatthew Knepley 
2209ea844a1aSMatthew Knepley   Level: intermediate
2210ea844a1aSMatthew Knepley 
2211ea844a1aSMatthew Knepley .seealso: PetscSectionSetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection
2212ea844a1aSMatthew Knepley @*/
2213ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetConstraintIndices(PetscSection s, PetscInt point, const PetscInt **indices)
2214ea844a1aSMatthew Knepley {
2215ea844a1aSMatthew Knepley   PetscErrorCode ierr;
2216ea844a1aSMatthew Knepley 
2217ea844a1aSMatthew Knepley   PetscFunctionBegin;
2218ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
2219ea844a1aSMatthew Knepley   if (s->bc) {
2220ea844a1aSMatthew Knepley     ierr = VecIntGetValuesSection(s->bcIndices, s->bc, point, indices);CHKERRQ(ierr);
2221ea844a1aSMatthew Knepley   } else *indices = NULL;
2222ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2223ea844a1aSMatthew Knepley }
2224ea844a1aSMatthew Knepley 
2225ea844a1aSMatthew Knepley /*@C
2226ea844a1aSMatthew Knepley   PetscSectionSetConstraintIndices - Set the point dof numbers, in [0, dof), which are constrained
2227ea844a1aSMatthew Knepley 
2228ea844a1aSMatthew Knepley   Not collective
2229ea844a1aSMatthew Knepley 
2230ea844a1aSMatthew Knepley   Input Parameters:
2231ea844a1aSMatthew Knepley + s     - The PetscSection
2232ea844a1aSMatthew Knepley . point - The point
2233ea844a1aSMatthew Knepley - indices - The constrained dofs
2234ea844a1aSMatthew Knepley 
2235ea844a1aSMatthew Knepley   Note: The Fortran is PetscSectionSetConstraintIndicesF90()
2236ea844a1aSMatthew Knepley 
2237ea844a1aSMatthew Knepley   Level: intermediate
2238ea844a1aSMatthew Knepley 
2239ea844a1aSMatthew Knepley .seealso: PetscSectionGetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection
2240ea844a1aSMatthew Knepley @*/
2241ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetConstraintIndices(PetscSection s, PetscInt point, const PetscInt indices[])
2242ea844a1aSMatthew Knepley {
2243ea844a1aSMatthew Knepley   PetscErrorCode ierr;
2244ea844a1aSMatthew Knepley 
2245ea844a1aSMatthew Knepley   PetscFunctionBegin;
2246ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
2247ea844a1aSMatthew Knepley   if (s->bc) {
2248ea844a1aSMatthew Knepley     ierr = VecIntSetValuesSection(s->bcIndices, s->bc, point, indices, INSERT_VALUES);CHKERRQ(ierr);
2249ea844a1aSMatthew Knepley   }
2250ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2251ea844a1aSMatthew Knepley }
2252ea844a1aSMatthew Knepley 
2253ea844a1aSMatthew Knepley /*@C
2254ea844a1aSMatthew Knepley   PetscSectionGetFieldConstraintIndices - Get the field dof numbers, in [0, fdof), which are constrained
2255ea844a1aSMatthew Knepley 
2256ea844a1aSMatthew Knepley   Not collective
2257ea844a1aSMatthew Knepley 
2258ea844a1aSMatthew Knepley   Input Parameters:
2259ea844a1aSMatthew Knepley + s     - The PetscSection
2260ea844a1aSMatthew Knepley . field  - The field number
2261ea844a1aSMatthew Knepley - point - The point
2262ea844a1aSMatthew Knepley 
2263ea844a1aSMatthew Knepley   Output Parameter:
2264ea844a1aSMatthew Knepley . indices - The constrained dofs
2265ea844a1aSMatthew Knepley 
2266ea844a1aSMatthew Knepley   Note: In Fortran, you call PetscSectionGetFieldConstraintIndicesF90() and PetscSectionRestoreFieldConstraintIndicesF90()
2267ea844a1aSMatthew Knepley 
2268ea844a1aSMatthew Knepley   Level: intermediate
2269ea844a1aSMatthew Knepley 
2270ea844a1aSMatthew Knepley .seealso: PetscSectionSetFieldConstraintIndices(), PetscSectionGetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection
2271ea844a1aSMatthew Knepley @*/
2272ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetFieldConstraintIndices(PetscSection s, PetscInt point, PetscInt field, const PetscInt **indices)
2273ea844a1aSMatthew Knepley {
2274ea844a1aSMatthew Knepley   PetscErrorCode ierr;
2275ea844a1aSMatthew Knepley 
2276ea844a1aSMatthew Knepley   PetscFunctionBegin;
2277ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
2278ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
2279ea844a1aSMatthew Knepley   ierr = PetscSectionGetConstraintIndices(s->field[field], point, indices);CHKERRQ(ierr);
2280ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2281ea844a1aSMatthew Knepley }
2282ea844a1aSMatthew Knepley 
2283ea844a1aSMatthew Knepley /*@C
2284ea844a1aSMatthew Knepley   PetscSectionSetFieldConstraintIndices - Set the field dof numbers, in [0, fdof), which are constrained
2285ea844a1aSMatthew Knepley 
2286ea844a1aSMatthew Knepley   Not collective
2287ea844a1aSMatthew Knepley 
2288ea844a1aSMatthew Knepley   Input Parameters:
2289ea844a1aSMatthew Knepley + s       - The PetscSection
2290ea844a1aSMatthew Knepley . point   - The point
2291ea844a1aSMatthew Knepley . field   - The field number
2292ea844a1aSMatthew Knepley - indices - The constrained dofs
2293ea844a1aSMatthew Knepley 
2294ea844a1aSMatthew Knepley   Note: The Fortran is PetscSectionSetFieldConstraintIndicesF90()
2295ea844a1aSMatthew Knepley 
2296ea844a1aSMatthew Knepley   Level: intermediate
2297ea844a1aSMatthew Knepley 
2298ea844a1aSMatthew Knepley .seealso: PetscSectionSetConstraintIndices(), PetscSectionGetFieldConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection
2299ea844a1aSMatthew Knepley @*/
2300ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetFieldConstraintIndices(PetscSection s, PetscInt point, PetscInt field, const PetscInt indices[])
2301ea844a1aSMatthew Knepley {
2302ea844a1aSMatthew Knepley   PetscErrorCode ierr;
2303ea844a1aSMatthew Knepley 
2304ea844a1aSMatthew Knepley   PetscFunctionBegin;
2305ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
2306ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
2307ea844a1aSMatthew Knepley   ierr = PetscSectionSetConstraintIndices(s->field[field], point, indices);CHKERRQ(ierr);
2308ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2309ea844a1aSMatthew Knepley }
2310ea844a1aSMatthew Knepley 
2311ea844a1aSMatthew Knepley /*@
2312ea844a1aSMatthew Knepley   PetscSectionPermute - Reorder the section according to the input point permutation
2313ea844a1aSMatthew Knepley 
2314ea844a1aSMatthew Knepley   Collective on PetscSection
2315ea844a1aSMatthew Knepley 
2316ea844a1aSMatthew Knepley   Input Parameter:
2317ea844a1aSMatthew Knepley + section - The PetscSection object
2318ea844a1aSMatthew Knepley - perm - The point permutation, old point p becomes new point perm[p]
2319ea844a1aSMatthew Knepley 
2320ea844a1aSMatthew Knepley   Output Parameter:
2321ea844a1aSMatthew Knepley . sectionNew - The permuted PetscSection
2322ea844a1aSMatthew Knepley 
2323ea844a1aSMatthew Knepley   Level: intermediate
2324ea844a1aSMatthew Knepley 
2325ea844a1aSMatthew Knepley .seealso: MatPermute()
2326ea844a1aSMatthew Knepley @*/
2327ea844a1aSMatthew Knepley PetscErrorCode PetscSectionPermute(PetscSection section, IS permutation, PetscSection *sectionNew)
2328ea844a1aSMatthew Knepley {
2329ea844a1aSMatthew Knepley   PetscSection    s = section, sNew;
2330ea844a1aSMatthew Knepley   const PetscInt *perm;
2331ea844a1aSMatthew Knepley   PetscInt        numFields, f, numPoints, pStart, pEnd, p;
2332ea844a1aSMatthew Knepley   PetscErrorCode  ierr;
2333ea844a1aSMatthew Knepley 
2334ea844a1aSMatthew Knepley   PetscFunctionBegin;
2335ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1);
2336ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(permutation, IS_CLASSID, 2);
2337ea844a1aSMatthew Knepley   PetscValidPointer(sectionNew, 3);
2338ea844a1aSMatthew Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), &sNew);CHKERRQ(ierr);
2339ea844a1aSMatthew Knepley   ierr = PetscSectionGetNumFields(s, &numFields);CHKERRQ(ierr);
2340ea844a1aSMatthew Knepley   if (numFields) {ierr = PetscSectionSetNumFields(sNew, numFields);CHKERRQ(ierr);}
2341ea844a1aSMatthew Knepley   for (f = 0; f < numFields; ++f) {
2342ea844a1aSMatthew Knepley     const char *name;
2343ea844a1aSMatthew Knepley     PetscInt    numComp;
2344ea844a1aSMatthew Knepley 
2345ea844a1aSMatthew Knepley     ierr = PetscSectionGetFieldName(s, f, &name);CHKERRQ(ierr);
2346ea844a1aSMatthew Knepley     ierr = PetscSectionSetFieldName(sNew, f, name);CHKERRQ(ierr);
2347ea844a1aSMatthew Knepley     ierr = PetscSectionGetFieldComponents(s, f, &numComp);CHKERRQ(ierr);
2348ea844a1aSMatthew Knepley     ierr = PetscSectionSetFieldComponents(sNew, f, numComp);CHKERRQ(ierr);
2349ea844a1aSMatthew Knepley   }
2350ea844a1aSMatthew Knepley   ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr);
2351ea844a1aSMatthew Knepley   ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr);
2352ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
2353ea844a1aSMatthew Knepley   ierr = PetscSectionSetChart(sNew, pStart, pEnd);CHKERRQ(ierr);
2354ea844a1aSMatthew Knepley   if (numPoints < pEnd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Permutation size %D is less than largest Section point %D", numPoints, pEnd);
2355ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
2356ea844a1aSMatthew Knepley     PetscInt dof, cdof;
2357ea844a1aSMatthew Knepley 
2358ea844a1aSMatthew Knepley     ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr);
2359ea844a1aSMatthew Knepley     ierr = PetscSectionSetDof(sNew, perm[p], dof);CHKERRQ(ierr);
2360ea844a1aSMatthew Knepley     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
2361ea844a1aSMatthew Knepley     if (cdof) {ierr = PetscSectionSetConstraintDof(sNew, perm[p], cdof);CHKERRQ(ierr);}
2362ea844a1aSMatthew Knepley     for (f = 0; f < numFields; ++f) {
2363ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldDof(s, p, f, &dof);CHKERRQ(ierr);
2364ea844a1aSMatthew Knepley       ierr = PetscSectionSetFieldDof(sNew, perm[p], f, dof);CHKERRQ(ierr);
2365ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldConstraintDof(s, p, f, &cdof);CHKERRQ(ierr);
2366ea844a1aSMatthew Knepley       if (cdof) {ierr = PetscSectionSetFieldConstraintDof(sNew, perm[p], f, cdof);CHKERRQ(ierr);}
2367ea844a1aSMatthew Knepley     }
2368ea844a1aSMatthew Knepley   }
2369ea844a1aSMatthew Knepley   ierr = PetscSectionSetUp(sNew);CHKERRQ(ierr);
2370ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
2371ea844a1aSMatthew Knepley     const PetscInt *cind;
2372ea844a1aSMatthew Knepley     PetscInt        cdof;
2373ea844a1aSMatthew Knepley 
2374ea844a1aSMatthew Knepley     ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
2375ea844a1aSMatthew Knepley     if (cdof) {
2376ea844a1aSMatthew Knepley       ierr = PetscSectionGetConstraintIndices(s, p, &cind);CHKERRQ(ierr);
2377ea844a1aSMatthew Knepley       ierr = PetscSectionSetConstraintIndices(sNew, perm[p], cind);CHKERRQ(ierr);
2378ea844a1aSMatthew Knepley     }
2379ea844a1aSMatthew Knepley     for (f = 0; f < numFields; ++f) {
2380ea844a1aSMatthew Knepley       ierr = PetscSectionGetFieldConstraintDof(s, p, f, &cdof);CHKERRQ(ierr);
2381ea844a1aSMatthew Knepley       if (cdof) {
2382ea844a1aSMatthew Knepley         ierr = PetscSectionGetFieldConstraintIndices(s, p, f, &cind);CHKERRQ(ierr);
2383ea844a1aSMatthew Knepley         ierr = PetscSectionSetFieldConstraintIndices(sNew, perm[p], f, cind);CHKERRQ(ierr);
2384ea844a1aSMatthew Knepley       }
2385ea844a1aSMatthew Knepley     }
2386ea844a1aSMatthew Knepley   }
2387ea844a1aSMatthew Knepley   ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr);
2388ea844a1aSMatthew Knepley   *sectionNew = sNew;
2389ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2390ea844a1aSMatthew Knepley }
2391ea844a1aSMatthew Knepley 
2392ea844a1aSMatthew Knepley /* TODO: the next three functions should be moved to sf/utils */
2393ea844a1aSMatthew Knepley #include <petsc/private/sfimpl.h>
2394ea844a1aSMatthew Knepley 
2395ea844a1aSMatthew Knepley /*@C
2396ea844a1aSMatthew Knepley   PetscSFDistributeSection - Create a new PetscSection reorganized, moving from the root to the leaves of the SF
2397ea844a1aSMatthew Knepley 
2398ea844a1aSMatthew Knepley   Collective on sf
2399ea844a1aSMatthew Knepley 
2400ea844a1aSMatthew Knepley   Input Parameters:
2401ea844a1aSMatthew Knepley + sf - The SF
2402ea844a1aSMatthew Knepley - rootSection - Section defined on root space
2403ea844a1aSMatthew Knepley 
2404ea844a1aSMatthew Knepley   Output Parameters:
2405ea844a1aSMatthew Knepley + remoteOffsets - root offsets in leaf storage, or NULL
2406ea844a1aSMatthew Knepley - leafSection - Section defined on the leaf space
2407ea844a1aSMatthew Knepley 
2408ea844a1aSMatthew Knepley   Level: advanced
2409ea844a1aSMatthew Knepley 
2410ea844a1aSMatthew Knepley .seealso: PetscSFCreate()
2411ea844a1aSMatthew Knepley @*/
2412ea844a1aSMatthew Knepley PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection)
2413ea844a1aSMatthew Knepley {
2414ea844a1aSMatthew Knepley   PetscSF        embedSF;
2415ea844a1aSMatthew Knepley   const PetscInt *indices;
2416ea844a1aSMatthew Knepley   IS             selected;
2417ea844a1aSMatthew Knepley   PetscInt       numFields, nroots, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, f;
241828d58a37SPierre Jolivet   PetscBool      *sub, hasc;
2419ea844a1aSMatthew Knepley   PetscErrorCode ierr;
2420ea844a1aSMatthew Knepley 
2421ea844a1aSMatthew Knepley   PetscFunctionBegin;
2422ea844a1aSMatthew Knepley   ierr = PetscLogEventBegin(PETSCSF_DistSect,sf,0,0,0);CHKERRQ(ierr);
2423ea844a1aSMatthew Knepley   ierr = PetscSectionGetNumFields(rootSection, &numFields);CHKERRQ(ierr);
2424ea844a1aSMatthew Knepley   if (numFields) {ierr = PetscSectionSetNumFields(leafSection, numFields);CHKERRQ(ierr);}
242528d58a37SPierre Jolivet   ierr = PetscMalloc1(numFields+2, &sub);CHKERRQ(ierr);
242628d58a37SPierre Jolivet   sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE;
2427ea844a1aSMatthew Knepley   for (f = 0; f < numFields; ++f) {
242828d58a37SPierre Jolivet     PetscSectionSym sym;
242928d58a37SPierre Jolivet     const char      *name   = NULL;
2430ea844a1aSMatthew Knepley     PetscInt        numComp = 0;
243128d58a37SPierre Jolivet 
243228d58a37SPierre Jolivet     sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE;
2433ea844a1aSMatthew Knepley     ierr = PetscSectionGetFieldComponents(rootSection, f, &numComp);CHKERRQ(ierr);
243428d58a37SPierre Jolivet     ierr = PetscSectionGetFieldName(rootSection, f, &name);CHKERRQ(ierr);
243528d58a37SPierre Jolivet     ierr = PetscSectionGetFieldSym(rootSection, f, &sym);CHKERRQ(ierr);
2436ea844a1aSMatthew Knepley     ierr = PetscSectionSetFieldComponents(leafSection, f, numComp);CHKERRQ(ierr);
243728d58a37SPierre Jolivet     ierr = PetscSectionSetFieldName(leafSection, f, name);CHKERRQ(ierr);
243828d58a37SPierre Jolivet     ierr = PetscSectionSetFieldSym(leafSection, f, sym);CHKERRQ(ierr);
2439ea844a1aSMatthew Knepley   }
2440ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(rootSection, &rpStart, &rpEnd);CHKERRQ(ierr);
2441ea844a1aSMatthew Knepley   ierr = PetscSFGetGraph(sf,&nroots,NULL,NULL,NULL);CHKERRQ(ierr);
2442ea844a1aSMatthew Knepley   rpEnd = PetscMin(rpEnd,nroots);
2443ea844a1aSMatthew Knepley   rpEnd = PetscMax(rpStart,rpEnd);
2444ea844a1aSMatthew Knepley   /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */
244528d58a37SPierre Jolivet   sub[0] = (PetscBool)(nroots != rpEnd - rpStart);
244628d58a37SPierre Jolivet   ierr = MPIU_Allreduce(MPI_IN_PLACE, sub, 2+numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf));CHKERRQ(ierr);
244728d58a37SPierre Jolivet   if (sub[0]) {
2448ea844a1aSMatthew Knepley     ierr = ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected);CHKERRQ(ierr);
2449ea844a1aSMatthew Knepley     ierr = ISGetIndices(selected, &indices);CHKERRQ(ierr);
2450ea844a1aSMatthew Knepley     ierr = PetscSFCreateEmbeddedSF(sf, rpEnd - rpStart, indices, &embedSF);CHKERRQ(ierr);
2451ea844a1aSMatthew Knepley     ierr = ISRestoreIndices(selected, &indices);CHKERRQ(ierr);
2452ea844a1aSMatthew Knepley     ierr = ISDestroy(&selected);CHKERRQ(ierr);
2453ea844a1aSMatthew Knepley   } else {
2454ea844a1aSMatthew Knepley     ierr = PetscObjectReference((PetscObject)sf);CHKERRQ(ierr);
2455ea844a1aSMatthew Knepley     embedSF = sf;
2456ea844a1aSMatthew Knepley   }
2457ea844a1aSMatthew Knepley   ierr = PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd);CHKERRQ(ierr);
2458ea844a1aSMatthew Knepley   lpEnd++;
2459ea844a1aSMatthew Knepley 
2460ea844a1aSMatthew Knepley   ierr = PetscSectionSetChart(leafSection, lpStart, lpEnd);CHKERRQ(ierr);
246128d58a37SPierre Jolivet 
246228d58a37SPierre Jolivet   /* Constrained dof section */
246328d58a37SPierre Jolivet   hasc = sub[1];
246428d58a37SPierre Jolivet   for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2+f]);
246528d58a37SPierre Jolivet 
246628d58a37SPierre Jolivet   /* Could fuse these at the cost of copies and extra allocation */
2467ea844a1aSMatthew Knepley   ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart]);CHKERRQ(ierr);
2468ea844a1aSMatthew Knepley   ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart]);CHKERRQ(ierr);
246928d58a37SPierre Jolivet   if (sub[1]) {
247028d58a37SPierre Jolivet     ierr = PetscSectionCheckConstraints_Static(rootSection);CHKERRQ(ierr);
247128d58a37SPierre Jolivet     ierr = PetscSectionCheckConstraints_Static(leafSection);CHKERRQ(ierr);
247228d58a37SPierre Jolivet     ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart]);CHKERRQ(ierr);
247328d58a37SPierre Jolivet     ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart]);CHKERRQ(ierr);
247428d58a37SPierre Jolivet   }
2475ea844a1aSMatthew Knepley   for (f = 0; f < numFields; ++f) {
2476ea844a1aSMatthew Knepley     ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart]);CHKERRQ(ierr);
2477ea844a1aSMatthew Knepley     ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart]);CHKERRQ(ierr);
247828d58a37SPierre Jolivet     if (sub[2+f]) {
247928d58a37SPierre Jolivet       ierr = PetscSectionCheckConstraints_Static(rootSection->field[f]);CHKERRQ(ierr);
248028d58a37SPierre Jolivet       ierr = PetscSectionCheckConstraints_Static(leafSection->field[f]);CHKERRQ(ierr);
248128d58a37SPierre Jolivet       ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart]);CHKERRQ(ierr);
248228d58a37SPierre Jolivet       ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart]);CHKERRQ(ierr);
248328d58a37SPierre Jolivet     }
2484ea844a1aSMatthew Knepley   }
2485ea844a1aSMatthew Knepley   if (remoteOffsets) {
2486ea844a1aSMatthew Knepley     ierr = PetscMalloc1(lpEnd - lpStart, remoteOffsets);CHKERRQ(ierr);
2487ea844a1aSMatthew Knepley     ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart]);CHKERRQ(ierr);
2488ea844a1aSMatthew Knepley     ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart]);CHKERRQ(ierr);
2489ea844a1aSMatthew Knepley   }
2490ea844a1aSMatthew Knepley   ierr = PetscSectionSetUp(leafSection);CHKERRQ(ierr);
249128d58a37SPierre Jolivet   if (hasc) { /* need to communicate bcIndices */
249228d58a37SPierre Jolivet     PetscSF  bcSF;
249328d58a37SPierre Jolivet     PetscInt *rOffBc;
249428d58a37SPierre Jolivet 
249528d58a37SPierre Jolivet     ierr = PetscMalloc1(lpEnd - lpStart, &rOffBc);CHKERRQ(ierr);
249628d58a37SPierre Jolivet     if (sub[1]) {
249728d58a37SPierre Jolivet       ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart]);CHKERRQ(ierr);
249828d58a37SPierre Jolivet       ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart]);CHKERRQ(ierr);
249928d58a37SPierre Jolivet       ierr = PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF);CHKERRQ(ierr);
250028d58a37SPierre Jolivet       ierr = PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices);CHKERRQ(ierr);
250128d58a37SPierre Jolivet       ierr = PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices);CHKERRQ(ierr);
250228d58a37SPierre Jolivet       ierr = PetscSFDestroy(&bcSF);CHKERRQ(ierr);
250328d58a37SPierre Jolivet     }
250428d58a37SPierre Jolivet     for (f = 0; f < numFields; ++f) {
250528d58a37SPierre Jolivet       if (sub[2+f]) {
250628d58a37SPierre Jolivet         ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart]);CHKERRQ(ierr);
250728d58a37SPierre Jolivet         ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart]);CHKERRQ(ierr);
250828d58a37SPierre Jolivet         ierr = PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF);CHKERRQ(ierr);
250928d58a37SPierre Jolivet         ierr = PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices);CHKERRQ(ierr);
251028d58a37SPierre Jolivet         ierr = PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices);CHKERRQ(ierr);
251128d58a37SPierre Jolivet         ierr = PetscSFDestroy(&bcSF);CHKERRQ(ierr);
251228d58a37SPierre Jolivet       }
251328d58a37SPierre Jolivet     }
251428d58a37SPierre Jolivet     ierr = PetscFree(rOffBc);CHKERRQ(ierr);
251528d58a37SPierre Jolivet   }
251628d58a37SPierre Jolivet   ierr = PetscSFDestroy(&embedSF);CHKERRQ(ierr);
251728d58a37SPierre Jolivet   ierr = PetscFree(sub);CHKERRQ(ierr);
2518ea844a1aSMatthew Knepley   ierr = PetscLogEventEnd(PETSCSF_DistSect,sf,0,0,0);CHKERRQ(ierr);
2519ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2520ea844a1aSMatthew Knepley }
2521ea844a1aSMatthew Knepley 
2522ea844a1aSMatthew Knepley /*@C
2523ea844a1aSMatthew Knepley   PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes
2524ea844a1aSMatthew Knepley 
2525ea844a1aSMatthew Knepley   Collective on sf
2526ea844a1aSMatthew Knepley 
2527ea844a1aSMatthew Knepley   Input Parameters:
2528ea844a1aSMatthew Knepley + sf          - The SF
2529ea844a1aSMatthew Knepley . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots)
2530ea844a1aSMatthew Knepley - leafSection - Data layout of local points for incoming data  (this is layout for SF leaves)
2531ea844a1aSMatthew Knepley 
2532ea844a1aSMatthew Knepley   Output Parameter:
2533ea844a1aSMatthew Knepley . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL
2534ea844a1aSMatthew Knepley 
2535ea844a1aSMatthew Knepley   Level: developer
2536ea844a1aSMatthew Knepley 
2537ea844a1aSMatthew Knepley .seealso: PetscSFCreate()
2538ea844a1aSMatthew Knepley @*/
2539ea844a1aSMatthew Knepley PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets)
2540ea844a1aSMatthew Knepley {
2541ea844a1aSMatthew Knepley   PetscSF         embedSF;
2542ea844a1aSMatthew Knepley   const PetscInt *indices;
2543ea844a1aSMatthew Knepley   IS              selected;
2544ea844a1aSMatthew Knepley   PetscInt        numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0;
2545ea844a1aSMatthew Knepley   PetscErrorCode  ierr;
2546ea844a1aSMatthew Knepley 
2547ea844a1aSMatthew Knepley   PetscFunctionBegin;
2548ea844a1aSMatthew Knepley   *remoteOffsets = NULL;
2549ea844a1aSMatthew Knepley   ierr = PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL);CHKERRQ(ierr);
2550ea844a1aSMatthew Knepley   if (numRoots < 0) PetscFunctionReturn(0);
2551ea844a1aSMatthew Knepley   ierr = PetscLogEventBegin(PETSCSF_RemoteOff,sf,0,0,0);CHKERRQ(ierr);
2552ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(rootSection, &rpStart, &rpEnd);CHKERRQ(ierr);
2553ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(leafSection, &lpStart, &lpEnd);CHKERRQ(ierr);
2554ea844a1aSMatthew Knepley   ierr = ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected);CHKERRQ(ierr);
2555ea844a1aSMatthew Knepley   ierr = ISGetIndices(selected, &indices);CHKERRQ(ierr);
2556ea844a1aSMatthew Knepley   ierr = PetscSFCreateEmbeddedSF(sf, rpEnd - rpStart, indices, &embedSF);CHKERRQ(ierr);
2557ea844a1aSMatthew Knepley   ierr = ISRestoreIndices(selected, &indices);CHKERRQ(ierr);
2558ea844a1aSMatthew Knepley   ierr = ISDestroy(&selected);CHKERRQ(ierr);
2559ea844a1aSMatthew Knepley   ierr = PetscCalloc1(lpEnd - lpStart, remoteOffsets);CHKERRQ(ierr);
2560ea844a1aSMatthew Knepley   ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart]);CHKERRQ(ierr);
2561ea844a1aSMatthew Knepley   ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart]);CHKERRQ(ierr);
2562ea844a1aSMatthew Knepley   ierr = PetscSFDestroy(&embedSF);CHKERRQ(ierr);
2563ea844a1aSMatthew Knepley   ierr = PetscLogEventEnd(PETSCSF_RemoteOff,sf,0,0,0);CHKERRQ(ierr);
2564ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2565ea844a1aSMatthew Knepley }
2566ea844a1aSMatthew Knepley 
2567ea844a1aSMatthew Knepley /*@C
2568ea844a1aSMatthew Knepley   PetscSFCreateSectionSF - Create an expanded SF of dofs, assuming the input SF relates points
2569ea844a1aSMatthew Knepley 
2570ea844a1aSMatthew Knepley   Collective on sf
2571ea844a1aSMatthew Knepley 
2572ea844a1aSMatthew Knepley   Input Parameters:
2573ea844a1aSMatthew Knepley + sf - The SF
2574ea844a1aSMatthew Knepley . rootSection - Data layout of remote points for outgoing data (this is usually the serial section)
2575ea844a1aSMatthew Knepley . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL
2576ea844a1aSMatthew Knepley - leafSection - Data layout of local points for incoming data  (this is the distributed section)
2577ea844a1aSMatthew Knepley 
2578ea844a1aSMatthew Knepley   Output Parameters:
2579ea844a1aSMatthew Knepley - sectionSF - The new SF
2580ea844a1aSMatthew Knepley 
2581ea844a1aSMatthew Knepley   Note: Either rootSection or remoteOffsets can be specified
2582ea844a1aSMatthew Knepley 
2583ea844a1aSMatthew Knepley   Level: advanced
2584ea844a1aSMatthew Knepley 
2585ea844a1aSMatthew Knepley .seealso: PetscSFCreate()
2586ea844a1aSMatthew Knepley @*/
2587ea844a1aSMatthew Knepley PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF)
2588ea844a1aSMatthew Knepley {
2589ea844a1aSMatthew Knepley   MPI_Comm          comm;
2590ea844a1aSMatthew Knepley   const PetscInt    *localPoints;
2591ea844a1aSMatthew Knepley   const PetscSFNode *remotePoints;
2592ea844a1aSMatthew Knepley   PetscInt          lpStart, lpEnd;
2593ea844a1aSMatthew Knepley   PetscInt          numRoots, numSectionRoots, numPoints, numIndices = 0;
2594ea844a1aSMatthew Knepley   PetscInt          *localIndices;
2595ea844a1aSMatthew Knepley   PetscSFNode       *remoteIndices;
2596ea844a1aSMatthew Knepley   PetscInt          i, ind;
2597ea844a1aSMatthew Knepley   PetscErrorCode    ierr;
2598ea844a1aSMatthew Knepley 
2599ea844a1aSMatthew Knepley   PetscFunctionBegin;
2600ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
2601ea844a1aSMatthew Knepley   PetscValidPointer(rootSection,2);
2602ea844a1aSMatthew Knepley   /* Cannot check PetscValidIntPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */
2603ea844a1aSMatthew Knepley   PetscValidPointer(leafSection,4);
2604ea844a1aSMatthew Knepley   PetscValidPointer(sectionSF,5);
2605ea844a1aSMatthew Knepley   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
2606ea844a1aSMatthew Knepley   ierr = PetscSFCreate(comm, sectionSF);CHKERRQ(ierr);
2607ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(leafSection, &lpStart, &lpEnd);CHKERRQ(ierr);
2608ea844a1aSMatthew Knepley   ierr = PetscSectionGetStorageSize(rootSection, &numSectionRoots);CHKERRQ(ierr);
2609ea844a1aSMatthew Knepley   ierr = PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints);CHKERRQ(ierr);
2610ea844a1aSMatthew Knepley   if (numRoots < 0) PetscFunctionReturn(0);
2611ea844a1aSMatthew Knepley   ierr = PetscLogEventBegin(PETSCSF_SectSF,sf,0,0,0);CHKERRQ(ierr);
2612ea844a1aSMatthew Knepley   for (i = 0; i < numPoints; ++i) {
2613ea844a1aSMatthew Knepley     PetscInt localPoint = localPoints ? localPoints[i] : i;
2614ea844a1aSMatthew Knepley     PetscInt dof;
2615ea844a1aSMatthew Knepley 
2616ea844a1aSMatthew Knepley     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
2617ea844a1aSMatthew Knepley       ierr = PetscSectionGetDof(leafSection, localPoint, &dof);CHKERRQ(ierr);
2618ea844a1aSMatthew Knepley       numIndices += dof;
2619ea844a1aSMatthew Knepley     }
2620ea844a1aSMatthew Knepley   }
2621ea844a1aSMatthew Knepley   ierr = PetscMalloc1(numIndices, &localIndices);CHKERRQ(ierr);
2622ea844a1aSMatthew Knepley   ierr = PetscMalloc1(numIndices, &remoteIndices);CHKERRQ(ierr);
2623ea844a1aSMatthew Knepley   /* Create new index graph */
2624ea844a1aSMatthew Knepley   for (i = 0, ind = 0; i < numPoints; ++i) {
2625ea844a1aSMatthew Knepley     PetscInt localPoint = localPoints ? localPoints[i] : i;
2626ea844a1aSMatthew Knepley     PetscInt rank       = remotePoints[i].rank;
2627ea844a1aSMatthew Knepley 
2628ea844a1aSMatthew Knepley     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
2629ea844a1aSMatthew Knepley       PetscInt remoteOffset = remoteOffsets[localPoint-lpStart];
2630ea844a1aSMatthew Knepley       PetscInt loff, dof, d;
2631ea844a1aSMatthew Knepley 
2632ea844a1aSMatthew Knepley       ierr = PetscSectionGetOffset(leafSection, localPoint, &loff);CHKERRQ(ierr);
2633ea844a1aSMatthew Knepley       ierr = PetscSectionGetDof(leafSection, localPoint, &dof);CHKERRQ(ierr);
2634ea844a1aSMatthew Knepley       for (d = 0; d < dof; ++d, ++ind) {
2635ea844a1aSMatthew Knepley         localIndices[ind]        = loff+d;
2636ea844a1aSMatthew Knepley         remoteIndices[ind].rank  = rank;
2637ea844a1aSMatthew Knepley         remoteIndices[ind].index = remoteOffset+d;
2638ea844a1aSMatthew Knepley       }
2639ea844a1aSMatthew Knepley     }
2640ea844a1aSMatthew Knepley   }
2641ea844a1aSMatthew Knepley   if (numIndices != ind) SETERRQ2(comm, PETSC_ERR_PLIB, "Inconsistency in indices, %D should be %D", ind, numIndices);
2642ea844a1aSMatthew Knepley   ierr = PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER);CHKERRQ(ierr);
2643ea844a1aSMatthew Knepley   ierr = PetscSFSetUp(*sectionSF);CHKERRQ(ierr);
2644ea844a1aSMatthew Knepley   ierr = PetscLogEventEnd(PETSCSF_SectSF,sf,0,0,0);CHKERRQ(ierr);
2645ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2646ea844a1aSMatthew Knepley }
2647ea844a1aSMatthew Knepley 
2648ea844a1aSMatthew Knepley /*@
2649ea844a1aSMatthew Knepley   PetscSectionSetClosureIndex - Set a cache of points in the closure of each point in the section
2650ea844a1aSMatthew Knepley 
2651ea844a1aSMatthew Knepley   Collective on section
2652ea844a1aSMatthew Knepley 
2653ea844a1aSMatthew Knepley   Input Parameters:
2654ea844a1aSMatthew Knepley + section   - The PetscSection
2655ea844a1aSMatthew Knepley . obj       - A PetscObject which serves as the key for this index
2656ea844a1aSMatthew Knepley . clSection - Section giving the size of the closure of each point
2657ea844a1aSMatthew Knepley - clPoints  - IS giving the points in each closure
2658ea844a1aSMatthew Knepley 
2659ea844a1aSMatthew Knepley   Note: We compress out closure points with no dofs in this section
2660ea844a1aSMatthew Knepley 
2661ea844a1aSMatthew Knepley   Level: advanced
2662ea844a1aSMatthew Knepley 
2663ea844a1aSMatthew Knepley .seealso: PetscSectionGetClosureIndex(), DMPlexCreateClosureIndex()
2664ea844a1aSMatthew Knepley @*/
2665ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetClosureIndex(PetscSection section, PetscObject obj, PetscSection clSection, IS clPoints)
2666ea844a1aSMatthew Knepley {
2667ea844a1aSMatthew Knepley   PetscErrorCode ierr;
2668ea844a1aSMatthew Knepley 
2669ea844a1aSMatthew Knepley   PetscFunctionBegin;
2670ea844a1aSMatthew Knepley   if (section->clObj != obj) {ierr = PetscFree(section->clPerm);CHKERRQ(ierr);ierr = PetscFree(section->clInvPerm);CHKERRQ(ierr);}
2671ea844a1aSMatthew Knepley   section->clObj     = obj;
2672ea844a1aSMatthew Knepley   ierr = PetscSectionDestroy(&section->clSection);CHKERRQ(ierr);
2673ea844a1aSMatthew Knepley   ierr = ISDestroy(&section->clPoints);CHKERRQ(ierr);
2674ea844a1aSMatthew Knepley   section->clSection = clSection;
2675ea844a1aSMatthew Knepley   section->clPoints  = clPoints;
2676ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2677ea844a1aSMatthew Knepley }
2678ea844a1aSMatthew Knepley 
2679ea844a1aSMatthew Knepley /*@
2680ea844a1aSMatthew Knepley   PetscSectionGetClosureIndex - Get the cache of points in the closure of each point in the section
2681ea844a1aSMatthew Knepley 
2682ea844a1aSMatthew Knepley   Collective on section
2683ea844a1aSMatthew Knepley 
2684ea844a1aSMatthew Knepley   Input Parameters:
2685ea844a1aSMatthew Knepley + section   - The PetscSection
2686ea844a1aSMatthew Knepley - obj       - A PetscObject which serves as the key for this index
2687ea844a1aSMatthew Knepley 
2688ea844a1aSMatthew Knepley   Output Parameters:
2689ea844a1aSMatthew Knepley + clSection - Section giving the size of the closure of each point
2690ea844a1aSMatthew Knepley - clPoints  - IS giving the points in each closure
2691ea844a1aSMatthew Knepley 
2692ea844a1aSMatthew Knepley   Note: We compress out closure points with no dofs in this section
2693ea844a1aSMatthew Knepley 
2694ea844a1aSMatthew Knepley   Level: advanced
2695ea844a1aSMatthew Knepley 
2696ea844a1aSMatthew Knepley .seealso: PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex()
2697ea844a1aSMatthew Knepley @*/
2698ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetClosureIndex(PetscSection section, PetscObject obj, PetscSection *clSection, IS *clPoints)
2699ea844a1aSMatthew Knepley {
2700ea844a1aSMatthew Knepley   PetscFunctionBegin;
2701ea844a1aSMatthew Knepley   if (section->clObj == obj) {
2702ea844a1aSMatthew Knepley     if (clSection) *clSection = section->clSection;
2703ea844a1aSMatthew Knepley     if (clPoints)  *clPoints  = section->clPoints;
2704ea844a1aSMatthew Knepley   } else {
2705ea844a1aSMatthew Knepley     if (clSection) *clSection = NULL;
2706ea844a1aSMatthew Knepley     if (clPoints)  *clPoints  = NULL;
2707ea844a1aSMatthew Knepley   }
2708ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2709ea844a1aSMatthew Knepley }
2710ea844a1aSMatthew Knepley 
2711ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetClosurePermutation_Internal(PetscSection section, PetscObject obj, PetscInt clSize, PetscCopyMode mode, PetscInt *clPerm)
2712ea844a1aSMatthew Knepley {
2713ea844a1aSMatthew Knepley   PetscInt       i;
2714ea844a1aSMatthew Knepley   PetscErrorCode ierr;
2715ea844a1aSMatthew Knepley 
2716ea844a1aSMatthew Knepley   PetscFunctionBegin;
2717ea844a1aSMatthew Knepley   if (section->clObj != obj) {
2718ea844a1aSMatthew Knepley     ierr = PetscSectionDestroy(&section->clSection);CHKERRQ(ierr);
2719ea844a1aSMatthew Knepley     ierr = ISDestroy(&section->clPoints);CHKERRQ(ierr);
2720ea844a1aSMatthew Knepley   }
2721ea844a1aSMatthew Knepley   section->clObj  = obj;
2722ea844a1aSMatthew Knepley   ierr = PetscFree(section->clPerm);CHKERRQ(ierr);
2723ea844a1aSMatthew Knepley   ierr = PetscFree(section->clInvPerm);CHKERRQ(ierr);
2724ea844a1aSMatthew Knepley   section->clSize = clSize;
2725ea844a1aSMatthew Knepley   if (mode == PETSC_COPY_VALUES) {
2726ea844a1aSMatthew Knepley     ierr = PetscMalloc1(clSize, &section->clPerm);CHKERRQ(ierr);
2727ea844a1aSMatthew Knepley     ierr = PetscLogObjectMemory((PetscObject) obj, clSize*sizeof(PetscInt));CHKERRQ(ierr);
2728ea844a1aSMatthew Knepley     ierr = PetscArraycpy(section->clPerm, clPerm, clSize);CHKERRQ(ierr);
2729ea844a1aSMatthew Knepley   } else if (mode == PETSC_OWN_POINTER) {
2730ea844a1aSMatthew Knepley     section->clPerm = clPerm;
2731ea844a1aSMatthew Knepley   } else SETERRQ(PetscObjectComm(obj), PETSC_ERR_SUP, "Do not support borrowed arrays");
2732ea844a1aSMatthew Knepley   ierr = PetscMalloc1(clSize, &section->clInvPerm);CHKERRQ(ierr);
2733ea844a1aSMatthew Knepley   for (i = 0; i < clSize; ++i) section->clInvPerm[section->clPerm[i]] = i;
2734ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2735ea844a1aSMatthew Knepley }
2736ea844a1aSMatthew Knepley 
2737ea844a1aSMatthew Knepley /*@
2738ea844a1aSMatthew Knepley   PetscSectionSetClosurePermutation - Get the dof permutation for the closure of each cell in the section, meaning clPerm[newIndex] = oldIndex.
2739ea844a1aSMatthew Knepley 
2740ea844a1aSMatthew Knepley   Not Collective
2741ea844a1aSMatthew Knepley 
2742ea844a1aSMatthew Knepley   Input Parameters:
2743ea844a1aSMatthew Knepley + section - The PetscSection
2744ea844a1aSMatthew Knepley . obj     - A PetscObject which serves as the key for this index
2745ea844a1aSMatthew Knepley - perm    - Permutation of the cell dof closure
2746ea844a1aSMatthew Knepley 
2747ea844a1aSMatthew Knepley   Note: This strategy only works when all cells have the same size dof closure, and no closures are retrieved for
2748ea844a1aSMatthew Knepley   other points (like faces).
2749ea844a1aSMatthew Knepley 
2750ea844a1aSMatthew Knepley   Level: intermediate
2751ea844a1aSMatthew Knepley 
2752ea844a1aSMatthew Knepley .seealso: PetscSectionGetClosurePermutation(), PetscSectionGetClosureIndex(), DMPlexCreateClosureIndex(), PetscCopyMode
2753ea844a1aSMatthew Knepley @*/
2754ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetClosurePermutation(PetscSection section, PetscObject obj, IS perm)
2755ea844a1aSMatthew Knepley {
2756ea844a1aSMatthew Knepley   const PetscInt *clPerm = NULL;
2757ea844a1aSMatthew Knepley   PetscInt        clSize = 0;
2758ea844a1aSMatthew Knepley   PetscErrorCode  ierr;
2759ea844a1aSMatthew Knepley 
2760ea844a1aSMatthew Knepley   PetscFunctionBegin;
2761ea844a1aSMatthew Knepley   if (perm) {
2762ea844a1aSMatthew Knepley     ierr = ISGetLocalSize(perm, &clSize);CHKERRQ(ierr);
2763ea844a1aSMatthew Knepley     ierr = ISGetIndices(perm, &clPerm);CHKERRQ(ierr);
2764ea844a1aSMatthew Knepley   }
2765ea844a1aSMatthew Knepley   ierr = PetscSectionSetClosurePermutation_Internal(section, obj, clSize, PETSC_COPY_VALUES, (PetscInt *) clPerm);CHKERRQ(ierr);
2766ea844a1aSMatthew Knepley   if (perm) {ierr = ISRestoreIndices(perm, &clPerm);CHKERRQ(ierr);}
2767ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2768ea844a1aSMatthew Knepley }
2769ea844a1aSMatthew Knepley 
2770ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetClosurePermutation_Internal(PetscSection section, PetscObject obj, PetscInt *size, const PetscInt *perm[])
2771ea844a1aSMatthew Knepley {
2772ea844a1aSMatthew Knepley   PetscFunctionBegin;
2773ea844a1aSMatthew Knepley   if (section->clObj == obj) {
2774ea844a1aSMatthew Knepley     if (size) *size = section->clSize;
2775ea844a1aSMatthew Knepley     if (perm) *perm = section->clPerm;
2776ea844a1aSMatthew Knepley   } else {
2777ea844a1aSMatthew Knepley     if (size) *size = 0;
2778ea844a1aSMatthew Knepley     if (perm) *perm = NULL;
2779ea844a1aSMatthew Knepley   }
2780ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2781ea844a1aSMatthew Knepley }
2782ea844a1aSMatthew Knepley 
2783ea844a1aSMatthew Knepley /*@
2784ea844a1aSMatthew Knepley   PetscSectionGetClosurePermutation - Get the dof permutation for the closure of each cell in the section, meaning clPerm[newIndex] = oldIndex.
2785ea844a1aSMatthew Knepley 
2786ea844a1aSMatthew Knepley   Not collective
2787ea844a1aSMatthew Knepley 
2788ea844a1aSMatthew Knepley   Input Parameters:
2789ea844a1aSMatthew Knepley + section   - The PetscSection
2790ea844a1aSMatthew Knepley - obj       - A PetscObject which serves as the key for this index
2791ea844a1aSMatthew Knepley 
2792ea844a1aSMatthew Knepley   Output Parameter:
2793ea844a1aSMatthew Knepley . perm - The dof closure permutation
2794ea844a1aSMatthew Knepley 
2795ea844a1aSMatthew Knepley   Note: This strategy only works when all cells have the same size dof closure, and no closures are retrieved for
2796ea844a1aSMatthew Knepley   other points (like faces).
2797ea844a1aSMatthew Knepley 
2798ea844a1aSMatthew Knepley   The user must destroy the IS that is returned.
2799ea844a1aSMatthew Knepley 
2800ea844a1aSMatthew Knepley   Level: intermediate
2801ea844a1aSMatthew Knepley 
2802ea844a1aSMatthew Knepley .seealso: PetscSectionSetClosurePermutation(), PetscSectionGetClosureInversePermutation(), PetscSectionGetClosureIndex(), PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex()
2803ea844a1aSMatthew Knepley @*/
2804ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetClosurePermutation(PetscSection section, PetscObject obj, IS *perm)
2805ea844a1aSMatthew Knepley {
2806ea844a1aSMatthew Knepley   const PetscInt *clPerm;
2807ea844a1aSMatthew Knepley   PetscInt        clSize;
2808ea844a1aSMatthew Knepley   PetscErrorCode  ierr;
2809ea844a1aSMatthew Knepley 
2810ea844a1aSMatthew Knepley   PetscFunctionBegin;
2811ea844a1aSMatthew Knepley   ierr = PetscSectionGetClosurePermutation_Internal(section, obj, &clSize, &clPerm);CHKERRQ(ierr);
2812ea844a1aSMatthew Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, clSize, clPerm, PETSC_USE_POINTER, perm);CHKERRQ(ierr);
2813ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2814ea844a1aSMatthew Knepley }
2815ea844a1aSMatthew Knepley 
2816ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetClosureInversePermutation_Internal(PetscSection section, PetscObject obj, PetscInt *size, const PetscInt *perm[])
2817ea844a1aSMatthew Knepley {
2818ea844a1aSMatthew Knepley   PetscFunctionBegin;
2819ea844a1aSMatthew Knepley   if (section->clObj == obj) {
2820ea844a1aSMatthew Knepley     if (size) *size = section->clSize;
2821ea844a1aSMatthew Knepley     if (perm) *perm = section->clInvPerm;
2822ea844a1aSMatthew Knepley   } else {
2823ea844a1aSMatthew Knepley     if (size) *size = 0;
2824ea844a1aSMatthew Knepley     if (perm) *perm = NULL;
2825ea844a1aSMatthew Knepley   }
2826ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2827ea844a1aSMatthew Knepley }
2828ea844a1aSMatthew Knepley 
2829ea844a1aSMatthew Knepley /*@
2830ea844a1aSMatthew Knepley   PetscSectionGetClosureInversePermutation - Get the inverse dof permutation for the closure of each cell in the section, meaning clPerm[oldIndex] = newIndex.
2831ea844a1aSMatthew Knepley 
2832ea844a1aSMatthew Knepley   Not collective
2833ea844a1aSMatthew Knepley 
2834ea844a1aSMatthew Knepley   Input Parameters:
2835ea844a1aSMatthew Knepley + section   - The PetscSection
2836ea844a1aSMatthew Knepley - obj       - A PetscObject which serves as the key for this index
2837ea844a1aSMatthew Knepley 
2838ea844a1aSMatthew Knepley   Output Parameters:
2839ea844a1aSMatthew Knepley + size - The dof closure size
2840ea844a1aSMatthew Knepley - perm - The dof closure permutation
2841ea844a1aSMatthew Knepley 
2842ea844a1aSMatthew Knepley   Note: This strategy only works when all cells have the same size dof closure, and no closures are retrieved for
2843ea844a1aSMatthew Knepley   other points (like faces).
2844ea844a1aSMatthew Knepley 
2845ea844a1aSMatthew Knepley   The user must destroy the IS that is returned.
2846ea844a1aSMatthew Knepley 
2847ea844a1aSMatthew Knepley   Level: intermediate
2848ea844a1aSMatthew Knepley 
2849ea844a1aSMatthew Knepley .seealso: PetscSectionSetClosurePermutation(), PetscSectionGetClosureIndex(), PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex()
2850ea844a1aSMatthew Knepley @*/
2851ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetClosureInversePermutation(PetscSection section, PetscObject obj, IS *perm)
2852ea844a1aSMatthew Knepley {
2853ea844a1aSMatthew Knepley   const PetscInt *clPerm;
2854ea844a1aSMatthew Knepley   PetscInt        clSize;
2855ea844a1aSMatthew Knepley   PetscErrorCode  ierr;
2856ea844a1aSMatthew Knepley 
2857ea844a1aSMatthew Knepley   PetscFunctionBegin;
2858ea844a1aSMatthew Knepley   ierr = PetscSectionGetClosureInversePermutation_Internal(section, obj, &clSize, &clPerm);CHKERRQ(ierr);
2859ea844a1aSMatthew Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, clSize, clPerm, PETSC_USE_POINTER, perm);CHKERRQ(ierr);
2860ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2861ea844a1aSMatthew Knepley }
2862ea844a1aSMatthew Knepley 
2863ea844a1aSMatthew Knepley /*@
2864ea844a1aSMatthew Knepley   PetscSectionGetField - Get the subsection associated with a single field
2865ea844a1aSMatthew Knepley 
2866ea844a1aSMatthew Knepley   Input Parameters:
2867ea844a1aSMatthew Knepley + s     - The PetscSection
2868ea844a1aSMatthew Knepley - field - The field number
2869ea844a1aSMatthew Knepley 
2870ea844a1aSMatthew Knepley   Output Parameter:
2871ea844a1aSMatthew Knepley . subs  - The subsection for the given field
2872ea844a1aSMatthew Knepley 
2873ea844a1aSMatthew Knepley   Level: intermediate
2874ea844a1aSMatthew Knepley 
2875ea844a1aSMatthew Knepley .seealso: PetscSectionSetNumFields()
2876ea844a1aSMatthew Knepley @*/
2877ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetField(PetscSection s, PetscInt field, PetscSection *subs)
2878ea844a1aSMatthew Knepley {
2879ea844a1aSMatthew Knepley   PetscFunctionBegin;
2880ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s,PETSC_SECTION_CLASSID,1);
2881ea844a1aSMatthew Knepley   PetscValidPointer(subs,3);
2882ea844a1aSMatthew Knepley   if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields);
2883ea844a1aSMatthew Knepley   *subs = s->field[field];
2884ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2885ea844a1aSMatthew Knepley }
2886ea844a1aSMatthew Knepley 
2887ea844a1aSMatthew Knepley PetscClassId      PETSC_SECTION_SYM_CLASSID;
2888ea844a1aSMatthew Knepley PetscFunctionList PetscSectionSymList = NULL;
2889ea844a1aSMatthew Knepley 
2890ea844a1aSMatthew Knepley /*@
2891ea844a1aSMatthew Knepley   PetscSectionSymCreate - Creates an empty PetscSectionSym object.
2892ea844a1aSMatthew Knepley 
2893ea844a1aSMatthew Knepley   Collective
2894ea844a1aSMatthew Knepley 
2895ea844a1aSMatthew Knepley   Input Parameter:
2896ea844a1aSMatthew Knepley . comm - the MPI communicator
2897ea844a1aSMatthew Knepley 
2898ea844a1aSMatthew Knepley   Output Parameter:
2899ea844a1aSMatthew Knepley . sym - pointer to the new set of symmetries
2900ea844a1aSMatthew Knepley 
2901ea844a1aSMatthew Knepley   Level: developer
2902ea844a1aSMatthew Knepley 
2903ea844a1aSMatthew Knepley .seealso: PetscSectionSym, PetscSectionSymDestroy()
2904ea844a1aSMatthew Knepley @*/
2905ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSymCreate(MPI_Comm comm, PetscSectionSym *sym)
2906ea844a1aSMatthew Knepley {
2907ea844a1aSMatthew Knepley   PetscErrorCode ierr;
2908ea844a1aSMatthew Knepley 
2909ea844a1aSMatthew Knepley   PetscFunctionBegin;
2910ea844a1aSMatthew Knepley   PetscValidPointer(sym,2);
2911ea844a1aSMatthew Knepley   ierr = ISInitializePackage();CHKERRQ(ierr);
2912ea844a1aSMatthew Knepley   ierr = PetscHeaderCreate(*sym,PETSC_SECTION_SYM_CLASSID,"PetscSectionSym","Section Symmetry","IS",comm,PetscSectionSymDestroy,PetscSectionSymView);CHKERRQ(ierr);
2913ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2914ea844a1aSMatthew Knepley }
2915ea844a1aSMatthew Knepley 
2916ea844a1aSMatthew Knepley /*@C
2917ea844a1aSMatthew Knepley   PetscSectionSymSetType - Builds a PetscSection symmetry, for a particular implementation.
2918ea844a1aSMatthew Knepley 
2919ea844a1aSMatthew Knepley   Collective on PetscSectionSym
2920ea844a1aSMatthew Knepley 
2921ea844a1aSMatthew Knepley   Input Parameters:
2922ea844a1aSMatthew Knepley + sym    - The section symmetry object
2923ea844a1aSMatthew Knepley - method - The name of the section symmetry type
2924ea844a1aSMatthew Knepley 
2925ea844a1aSMatthew Knepley   Level: developer
2926ea844a1aSMatthew Knepley 
2927ea844a1aSMatthew Knepley .seealso: PetscSectionSymGetType(), PetscSectionSymCreate()
2928ea844a1aSMatthew Knepley @*/
2929ea844a1aSMatthew Knepley PetscErrorCode  PetscSectionSymSetType(PetscSectionSym sym, PetscSectionSymType method)
2930ea844a1aSMatthew Knepley {
2931ea844a1aSMatthew Knepley   PetscErrorCode (*r)(PetscSectionSym);
2932ea844a1aSMatthew Knepley   PetscBool      match;
2933ea844a1aSMatthew Knepley   PetscErrorCode ierr;
2934ea844a1aSMatthew Knepley 
2935ea844a1aSMatthew Knepley   PetscFunctionBegin;
2936ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID,1);
2937ea844a1aSMatthew Knepley   ierr = PetscObjectTypeCompare((PetscObject) sym, method, &match);CHKERRQ(ierr);
2938ea844a1aSMatthew Knepley   if (match) PetscFunctionReturn(0);
2939ea844a1aSMatthew Knepley 
2940ea844a1aSMatthew Knepley   ierr = PetscFunctionListFind(PetscSectionSymList,method,&r);CHKERRQ(ierr);
2941ea844a1aSMatthew Knepley   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscSectionSym type: %s", method);
2942ea844a1aSMatthew Knepley   if (sym->ops->destroy) {
2943ea844a1aSMatthew Knepley     ierr = (*sym->ops->destroy)(sym);CHKERRQ(ierr);
2944ea844a1aSMatthew Knepley     sym->ops->destroy = NULL;
2945ea844a1aSMatthew Knepley   }
2946ea844a1aSMatthew Knepley   ierr = (*r)(sym);CHKERRQ(ierr);
2947ea844a1aSMatthew Knepley   ierr = PetscObjectChangeTypeName((PetscObject)sym,method);CHKERRQ(ierr);
2948ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2949ea844a1aSMatthew Knepley }
2950ea844a1aSMatthew Knepley 
2951ea844a1aSMatthew Knepley 
2952ea844a1aSMatthew Knepley /*@C
2953ea844a1aSMatthew Knepley   PetscSectionSymGetType - Gets the section symmetry type name (as a string) from the PetscSectionSym.
2954ea844a1aSMatthew Knepley 
2955ea844a1aSMatthew Knepley   Not Collective
2956ea844a1aSMatthew Knepley 
2957ea844a1aSMatthew Knepley   Input Parameter:
2958ea844a1aSMatthew Knepley . sym  - The section symmetry
2959ea844a1aSMatthew Knepley 
2960ea844a1aSMatthew Knepley   Output Parameter:
2961ea844a1aSMatthew Knepley . type - The index set type name
2962ea844a1aSMatthew Knepley 
2963ea844a1aSMatthew Knepley   Level: developer
2964ea844a1aSMatthew Knepley 
2965ea844a1aSMatthew Knepley .seealso: PetscSectionSymSetType(), PetscSectionSymCreate()
2966ea844a1aSMatthew Knepley @*/
2967ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSymGetType(PetscSectionSym sym, PetscSectionSymType *type)
2968ea844a1aSMatthew Knepley {
2969ea844a1aSMatthew Knepley   PetscFunctionBegin;
2970ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID,1);
2971ea844a1aSMatthew Knepley   PetscValidCharPointer(type,2);
2972ea844a1aSMatthew Knepley   *type = ((PetscObject)sym)->type_name;
2973ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2974ea844a1aSMatthew Knepley }
2975ea844a1aSMatthew Knepley 
2976ea844a1aSMatthew Knepley /*@C
2977ea844a1aSMatthew Knepley   PetscSectionSymRegister - Adds a new section symmetry implementation
2978ea844a1aSMatthew Knepley 
2979ea844a1aSMatthew Knepley   Not Collective
2980ea844a1aSMatthew Knepley 
2981ea844a1aSMatthew Knepley   Input Parameters:
2982ea844a1aSMatthew Knepley + name        - The name of a new user-defined creation routine
2983ea844a1aSMatthew Knepley - create_func - The creation routine itself
2984ea844a1aSMatthew Knepley 
2985ea844a1aSMatthew Knepley   Notes:
2986ea844a1aSMatthew Knepley   PetscSectionSymRegister() may be called multiple times to add several user-defined vectors
2987ea844a1aSMatthew Knepley 
2988ea844a1aSMatthew Knepley   Level: developer
2989ea844a1aSMatthew Knepley 
2990ea844a1aSMatthew Knepley .seealso: PetscSectionSymCreate(), PetscSectionSymSetType()
2991ea844a1aSMatthew Knepley @*/
2992ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSymRegister(const char sname[], PetscErrorCode (*function)(PetscSectionSym))
2993ea844a1aSMatthew Knepley {
2994ea844a1aSMatthew Knepley   PetscErrorCode ierr;
2995ea844a1aSMatthew Knepley 
2996ea844a1aSMatthew Knepley   PetscFunctionBegin;
2997ea844a1aSMatthew Knepley   ierr = ISInitializePackage();CHKERRQ(ierr);
2998ea844a1aSMatthew Knepley   ierr = PetscFunctionListAdd(&PetscSectionSymList,sname,function);CHKERRQ(ierr);
2999ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3000ea844a1aSMatthew Knepley }
3001ea844a1aSMatthew Knepley 
3002ea844a1aSMatthew Knepley /*@
3003ea844a1aSMatthew Knepley    PetscSectionSymDestroy - Destroys a section symmetry.
3004ea844a1aSMatthew Knepley 
3005ea844a1aSMatthew Knepley    Collective on PetscSectionSym
3006ea844a1aSMatthew Knepley 
3007ea844a1aSMatthew Knepley    Input Parameters:
3008ea844a1aSMatthew Knepley .  sym - the section symmetry
3009ea844a1aSMatthew Knepley 
3010ea844a1aSMatthew Knepley    Level: developer
3011ea844a1aSMatthew Knepley 
3012ea844a1aSMatthew Knepley .seealso: PetscSectionSymCreate(), PetscSectionSymDestroy()
3013ea844a1aSMatthew Knepley @*/
3014ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSymDestroy(PetscSectionSym *sym)
3015ea844a1aSMatthew Knepley {
3016ea844a1aSMatthew Knepley   SymWorkLink    link,next;
3017ea844a1aSMatthew Knepley   PetscErrorCode ierr;
3018ea844a1aSMatthew Knepley 
3019ea844a1aSMatthew Knepley   PetscFunctionBegin;
3020ea844a1aSMatthew Knepley   if (!*sym) PetscFunctionReturn(0);
3021ea844a1aSMatthew Knepley   PetscValidHeaderSpecific((*sym),PETSC_SECTION_SYM_CLASSID,1);
3022ea844a1aSMatthew Knepley   if (--((PetscObject)(*sym))->refct > 0) {*sym = 0; PetscFunctionReturn(0);}
3023ea844a1aSMatthew Knepley   if ((*sym)->ops->destroy) {
3024ea844a1aSMatthew Knepley     ierr = (*(*sym)->ops->destroy)(*sym);CHKERRQ(ierr);
3025ea844a1aSMatthew Knepley   }
3026ea844a1aSMatthew Knepley   if ((*sym)->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Work array still checked out");
3027ea844a1aSMatthew Knepley   for (link=(*sym)->workin; link; link=next) {
3028ea844a1aSMatthew Knepley     next = link->next;
3029ea844a1aSMatthew Knepley     ierr = PetscFree2(*(PetscInt***)&link->perms,*(PetscScalar***)&link->rots);CHKERRQ(ierr);
3030ea844a1aSMatthew Knepley     ierr = PetscFree(link);CHKERRQ(ierr);
3031ea844a1aSMatthew Knepley   }
3032ea844a1aSMatthew Knepley   (*sym)->workin = NULL;
3033ea844a1aSMatthew Knepley   ierr = PetscHeaderDestroy(sym);CHKERRQ(ierr);
3034ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3035ea844a1aSMatthew Knepley }
3036ea844a1aSMatthew Knepley 
3037ea844a1aSMatthew Knepley /*@C
3038ea844a1aSMatthew Knepley    PetscSectionSymView - Displays a section symmetry
3039ea844a1aSMatthew Knepley 
3040ea844a1aSMatthew Knepley    Collective on PetscSectionSym
3041ea844a1aSMatthew Knepley 
3042ea844a1aSMatthew Knepley    Input Parameters:
3043ea844a1aSMatthew Knepley +  sym - the index set
3044ea844a1aSMatthew Knepley -  viewer - viewer used to display the set, for example PETSC_VIEWER_STDOUT_SELF.
3045ea844a1aSMatthew Knepley 
3046ea844a1aSMatthew Knepley    Level: developer
3047ea844a1aSMatthew Knepley 
3048ea844a1aSMatthew Knepley .seealso: PetscViewerASCIIOpen()
3049ea844a1aSMatthew Knepley @*/
3050ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSymView(PetscSectionSym sym,PetscViewer viewer)
3051ea844a1aSMatthew Knepley {
3052ea844a1aSMatthew Knepley   PetscErrorCode ierr;
3053ea844a1aSMatthew Knepley 
3054ea844a1aSMatthew Knepley   PetscFunctionBegin;
3055ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1);
3056ea844a1aSMatthew Knepley   if (!viewer) {
3057ea844a1aSMatthew Knepley     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sym),&viewer);CHKERRQ(ierr);
3058ea844a1aSMatthew Knepley   }
3059ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
3060ea844a1aSMatthew Knepley   PetscCheckSameComm(sym,1,viewer,2);
3061ea844a1aSMatthew Knepley   ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sym,viewer);CHKERRQ(ierr);
3062ea844a1aSMatthew Knepley   if (sym->ops->view) {
3063ea844a1aSMatthew Knepley     ierr = (*sym->ops->view)(sym,viewer);CHKERRQ(ierr);
3064ea844a1aSMatthew Knepley   }
3065ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3066ea844a1aSMatthew Knepley }
3067ea844a1aSMatthew Knepley 
3068ea844a1aSMatthew Knepley /*@
3069ea844a1aSMatthew Knepley   PetscSectionSetSym - Set the symmetries for the data referred to by the section
3070ea844a1aSMatthew Knepley 
3071ea844a1aSMatthew Knepley   Collective on PetscSection
3072ea844a1aSMatthew Knepley 
3073ea844a1aSMatthew Knepley   Input Parameters:
3074ea844a1aSMatthew Knepley + section - the section describing data layout
3075ea844a1aSMatthew Knepley - sym - the symmetry describing the affect of orientation on the access of the data
3076ea844a1aSMatthew Knepley 
3077ea844a1aSMatthew Knepley   Level: developer
3078ea844a1aSMatthew Knepley 
3079ea844a1aSMatthew Knepley .seealso: PetscSectionGetSym(), PetscSectionSymCreate()
3080ea844a1aSMatthew Knepley @*/
3081ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetSym(PetscSection section, PetscSectionSym sym)
3082ea844a1aSMatthew Knepley {
3083ea844a1aSMatthew Knepley   PetscErrorCode ierr;
3084ea844a1aSMatthew Knepley 
3085ea844a1aSMatthew Knepley   PetscFunctionBegin;
3086ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1);
3087ea844a1aSMatthew Knepley   ierr = PetscSectionSymDestroy(&(section->sym));CHKERRQ(ierr);
3088ea844a1aSMatthew Knepley   if (sym) {
3089ea844a1aSMatthew Knepley     PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,2);
3090ea844a1aSMatthew Knepley     PetscCheckSameComm(section,1,sym,2);
3091ea844a1aSMatthew Knepley     ierr = PetscObjectReference((PetscObject) sym);CHKERRQ(ierr);
3092ea844a1aSMatthew Knepley   }
3093ea844a1aSMatthew Knepley   section->sym = sym;
3094ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3095ea844a1aSMatthew Knepley }
3096ea844a1aSMatthew Knepley 
3097ea844a1aSMatthew Knepley /*@
3098ea844a1aSMatthew Knepley   PetscSectionGetSym - Get the symmetries for the data referred to by the section
3099ea844a1aSMatthew Knepley 
3100ea844a1aSMatthew Knepley   Not collective
3101ea844a1aSMatthew Knepley 
3102ea844a1aSMatthew Knepley   Input Parameters:
3103ea844a1aSMatthew Knepley . section - the section describing data layout
3104ea844a1aSMatthew Knepley 
3105ea844a1aSMatthew Knepley   Output Parameters:
3106ea844a1aSMatthew Knepley . sym - the symmetry describing the affect of orientation on the access of the data
3107ea844a1aSMatthew Knepley 
3108ea844a1aSMatthew Knepley   Level: developer
3109ea844a1aSMatthew Knepley 
3110ea844a1aSMatthew Knepley .seealso: PetscSectionSetSym(), PetscSectionSymCreate()
3111ea844a1aSMatthew Knepley @*/
3112ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetSym(PetscSection section, PetscSectionSym *sym)
3113ea844a1aSMatthew Knepley {
3114ea844a1aSMatthew Knepley   PetscFunctionBegin;
3115ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1);
3116ea844a1aSMatthew Knepley   *sym = section->sym;
3117ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3118ea844a1aSMatthew Knepley }
3119ea844a1aSMatthew Knepley 
3120ea844a1aSMatthew Knepley /*@
3121ea844a1aSMatthew Knepley   PetscSectionSetFieldSym - Set the symmetries for the data referred to by a field of the section
3122ea844a1aSMatthew Knepley 
3123ea844a1aSMatthew Knepley   Collective on PetscSection
3124ea844a1aSMatthew Knepley 
3125ea844a1aSMatthew Knepley   Input Parameters:
3126ea844a1aSMatthew Knepley + section - the section describing data layout
3127ea844a1aSMatthew Knepley . field - the field number
3128ea844a1aSMatthew Knepley - sym - the symmetry describing the affect of orientation on the access of the data
3129ea844a1aSMatthew Knepley 
3130ea844a1aSMatthew Knepley   Level: developer
3131ea844a1aSMatthew Knepley 
3132ea844a1aSMatthew Knepley .seealso: PetscSectionGetFieldSym(), PetscSectionSymCreate()
3133ea844a1aSMatthew Knepley @*/
3134ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetFieldSym(PetscSection section, PetscInt field, PetscSectionSym sym)
3135ea844a1aSMatthew Knepley {
3136ea844a1aSMatthew Knepley   PetscErrorCode ierr;
3137ea844a1aSMatthew Knepley 
3138ea844a1aSMatthew Knepley   PetscFunctionBegin;
3139ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1);
3140ea844a1aSMatthew Knepley   if (field < 0 || field >= section->numFields) SETERRQ2(PetscObjectComm((PetscObject)section),PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number %D (not in [0,%D)", field, section->numFields);
3141ea844a1aSMatthew Knepley   ierr = PetscSectionSetSym(section->field[field],sym);CHKERRQ(ierr);
3142ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3143ea844a1aSMatthew Knepley }
3144ea844a1aSMatthew Knepley 
3145ea844a1aSMatthew Knepley /*@
3146ea844a1aSMatthew Knepley   PetscSectionGetFieldSym - Get the symmetries for the data referred to by a field of the section
3147ea844a1aSMatthew Knepley 
3148ea844a1aSMatthew Knepley   Collective on PetscSection
3149ea844a1aSMatthew Knepley 
3150ea844a1aSMatthew Knepley   Input Parameters:
3151ea844a1aSMatthew Knepley + section - the section describing data layout
3152ea844a1aSMatthew Knepley - field - the field number
3153ea844a1aSMatthew Knepley 
3154ea844a1aSMatthew Knepley   Output Parameters:
3155ea844a1aSMatthew Knepley . sym - the symmetry describing the affect of orientation on the access of the data
3156ea844a1aSMatthew Knepley 
3157ea844a1aSMatthew Knepley   Level: developer
3158ea844a1aSMatthew Knepley 
3159ea844a1aSMatthew Knepley .seealso: PetscSectionSetFieldSym(), PetscSectionSymCreate()
3160ea844a1aSMatthew Knepley @*/
3161ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetFieldSym(PetscSection section, PetscInt field, PetscSectionSym *sym)
3162ea844a1aSMatthew Knepley {
3163ea844a1aSMatthew Knepley   PetscFunctionBegin;
3164ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1);
3165ea844a1aSMatthew Knepley   if (field < 0 || field >= section->numFields) SETERRQ2(PetscObjectComm((PetscObject)section),PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number %D (not in [0,%D)", field, section->numFields);
3166ea844a1aSMatthew Knepley   *sym = section->field[field]->sym;
3167ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3168ea844a1aSMatthew Knepley }
3169ea844a1aSMatthew Knepley 
3170ea844a1aSMatthew Knepley /*@C
3171ea844a1aSMatthew Knepley   PetscSectionGetPointSyms - Get the symmetries for a set of points in a PetscSection under specific orientations.
3172ea844a1aSMatthew Knepley 
3173ea844a1aSMatthew Knepley   Not collective
3174ea844a1aSMatthew Knepley 
3175ea844a1aSMatthew Knepley   Input Parameters:
3176ea844a1aSMatthew Knepley + section - the section
3177ea844a1aSMatthew Knepley . numPoints - the number of points
3178ea844a1aSMatthew Knepley - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an
3179ea844a1aSMatthew Knepley     arbitrary integer: its interpretation is up to sym.  Orientations are used by DM: for their interpretation in that
3180ea844a1aSMatthew Knepley     context, see DMPlexGetConeOrientation()).
3181ea844a1aSMatthew Knepley 
3182ea844a1aSMatthew Knepley   Output Parameter:
3183ea844a1aSMatthew Knepley + perms - The permutations for the given orientations (or NULL if there is no symmetry or the permutation is the identity).
3184ea844a1aSMatthew Knepley - rots - The field rotations symmetries for the given orientations (or NULL if there is no symmetry or the rotations are all
3185ea844a1aSMatthew Knepley     identity).
3186ea844a1aSMatthew Knepley 
3187ea844a1aSMatthew Knepley   Example of usage, gathering dofs into a local array (lArray) from a section array (sArray):
3188ea844a1aSMatthew Knepley .vb
3189ea844a1aSMatthew Knepley      const PetscInt    **perms;
3190ea844a1aSMatthew Knepley      const PetscScalar **rots;
3191ea844a1aSMatthew Knepley      PetscInt            lOffset;
3192ea844a1aSMatthew Knepley 
3193ea844a1aSMatthew Knepley      PetscSectionGetPointSyms(section,numPoints,points,&perms,&rots);
3194ea844a1aSMatthew Knepley      for (i = 0, lOffset = 0; i < numPoints; i++) {
3195ea844a1aSMatthew Knepley        PetscInt           point = points[2*i], dof, sOffset;
3196ea844a1aSMatthew Knepley        const PetscInt    *perm  = perms ? perms[i] : NULL;
3197ea844a1aSMatthew Knepley        const PetscScalar *rot   = rots  ? rots[i]  : NULL;
3198ea844a1aSMatthew Knepley 
3199ea844a1aSMatthew Knepley        PetscSectionGetDof(section,point,&dof);
3200ea844a1aSMatthew Knepley        PetscSectionGetOffset(section,point,&sOffset);
3201ea844a1aSMatthew Knepley 
3202ea844a1aSMatthew Knepley        if (perm) {for (j = 0; j < dof; j++) {lArray[lOffset + perm[j]]  = sArray[sOffset + j];}}
3203ea844a1aSMatthew Knepley        else      {for (j = 0; j < dof; j++) {lArray[lOffset +      j ]  = sArray[sOffset + j];}}
3204ea844a1aSMatthew Knepley        if (rot)  {for (j = 0; j < dof; j++) {lArray[lOffset +      j ] *= rot[j];             }}
3205ea844a1aSMatthew Knepley        lOffset += dof;
3206ea844a1aSMatthew Knepley      }
3207ea844a1aSMatthew Knepley      PetscSectionRestorePointSyms(section,numPoints,points,&perms,&rots);
3208ea844a1aSMatthew Knepley .ve
3209ea844a1aSMatthew Knepley 
3210ea844a1aSMatthew Knepley   Example of usage, adding dofs into a section array (sArray) from a local array (lArray):
3211ea844a1aSMatthew Knepley .vb
3212ea844a1aSMatthew Knepley      const PetscInt    **perms;
3213ea844a1aSMatthew Knepley      const PetscScalar **rots;
3214ea844a1aSMatthew Knepley      PetscInt            lOffset;
3215ea844a1aSMatthew Knepley 
3216ea844a1aSMatthew Knepley      PetscSectionGetPointSyms(section,numPoints,points,&perms,&rots);
3217ea844a1aSMatthew Knepley      for (i = 0, lOffset = 0; i < numPoints; i++) {
3218ea844a1aSMatthew Knepley        PetscInt           point = points[2*i], dof, sOffset;
3219ea844a1aSMatthew Knepley        const PetscInt    *perm  = perms ? perms[i] : NULL;
3220ea844a1aSMatthew Knepley        const PetscScalar *rot   = rots  ? rots[i]  : NULL;
3221ea844a1aSMatthew Knepley 
3222ea844a1aSMatthew Knepley        PetscSectionGetDof(section,point,&dof);
3223ea844a1aSMatthew Knepley        PetscSectionGetOffset(section,point,&sOff);
3224ea844a1aSMatthew Knepley 
3225ea844a1aSMatthew Knepley        if (perm) {for (j = 0; j < dof; j++) {sArray[sOffset + j] += lArray[lOffset + perm[j]] * (rot ? PetscConj(rot[perm[j]]) : 1.);}}
3226ea844a1aSMatthew Knepley        else      {for (j = 0; j < dof; j++) {sArray[sOffset + j] += lArray[lOffset +      j ] * (rot ? PetscConj(rot[     j ]) : 1.);}}
3227ea844a1aSMatthew Knepley        offset += dof;
3228ea844a1aSMatthew Knepley      }
3229ea844a1aSMatthew Knepley      PetscSectionRestorePointSyms(section,numPoints,points,&perms,&rots);
3230ea844a1aSMatthew Knepley .ve
3231ea844a1aSMatthew Knepley 
3232ea844a1aSMatthew Knepley   Level: developer
3233ea844a1aSMatthew Knepley 
3234ea844a1aSMatthew Knepley .seealso: PetscSectionRestorePointSyms(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym()
3235ea844a1aSMatthew Knepley @*/
3236ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetPointSyms(PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots)
3237ea844a1aSMatthew Knepley {
3238ea844a1aSMatthew Knepley   PetscSectionSym sym;
3239ea844a1aSMatthew Knepley   PetscErrorCode  ierr;
3240ea844a1aSMatthew Knepley 
3241ea844a1aSMatthew Knepley   PetscFunctionBegin;
3242ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1);
3243ea844a1aSMatthew Knepley   if (numPoints) PetscValidIntPointer(points,3);
3244ea844a1aSMatthew Knepley   if (perms) *perms = NULL;
3245ea844a1aSMatthew Knepley   if (rots)  *rots  = NULL;
3246ea844a1aSMatthew Knepley   sym = section->sym;
3247ea844a1aSMatthew Knepley   if (sym && (perms || rots)) {
3248ea844a1aSMatthew Knepley     SymWorkLink link;
3249ea844a1aSMatthew Knepley 
3250ea844a1aSMatthew Knepley     if (sym->workin) {
3251ea844a1aSMatthew Knepley       link        = sym->workin;
3252ea844a1aSMatthew Knepley       sym->workin = sym->workin->next;
3253ea844a1aSMatthew Knepley     } else {
3254ea844a1aSMatthew Knepley       ierr = PetscNewLog(sym,&link);CHKERRQ(ierr);
3255ea844a1aSMatthew Knepley     }
3256ea844a1aSMatthew Knepley     if (numPoints > link->numPoints) {
3257ea844a1aSMatthew Knepley       ierr = PetscFree2(*(PetscInt***)&link->perms,*(PetscInt***)&link->rots);CHKERRQ(ierr);
3258ea844a1aSMatthew Knepley       ierr = PetscMalloc2(numPoints,(PetscInt***)&link->perms,numPoints,(PetscScalar***)&link->rots);CHKERRQ(ierr);
3259ea844a1aSMatthew Knepley       link->numPoints = numPoints;
3260ea844a1aSMatthew Knepley     }
3261ea844a1aSMatthew Knepley     link->next   = sym->workout;
3262ea844a1aSMatthew Knepley     sym->workout = link;
3263ea844a1aSMatthew Knepley     ierr = PetscArrayzero((PetscInt**)link->perms,numPoints);CHKERRQ(ierr);
3264ea844a1aSMatthew Knepley     ierr = PetscArrayzero((PetscInt**)link->rots,numPoints);CHKERRQ(ierr);
3265ea844a1aSMatthew Knepley     ierr = (*sym->ops->getpoints) (sym, section, numPoints, points, link->perms, link->rots);CHKERRQ(ierr);
3266ea844a1aSMatthew Knepley     if (perms) *perms = link->perms;
3267ea844a1aSMatthew Knepley     if (rots)  *rots  = link->rots;
3268ea844a1aSMatthew Knepley   }
3269ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3270ea844a1aSMatthew Knepley }
3271ea844a1aSMatthew Knepley 
3272ea844a1aSMatthew Knepley /*@C
3273ea844a1aSMatthew Knepley   PetscSectionRestorePointSyms - Restore the symmetries returned by PetscSectionGetPointSyms()
3274ea844a1aSMatthew Knepley 
3275ea844a1aSMatthew Knepley   Not collective
3276ea844a1aSMatthew Knepley 
3277ea844a1aSMatthew Knepley   Input Parameters:
3278ea844a1aSMatthew Knepley + section - the section
3279ea844a1aSMatthew Knepley . numPoints - the number of points
3280ea844a1aSMatthew Knepley - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an
3281ea844a1aSMatthew Knepley     arbitrary integer: its interpretation is up to sym.  Orientations are used by DM: for their interpretation in that
3282ea844a1aSMatthew Knepley     context, see DMPlexGetConeOrientation()).
3283ea844a1aSMatthew Knepley 
3284ea844a1aSMatthew Knepley   Output Parameter:
3285ea844a1aSMatthew Knepley + perms - The permutations for the given orientations: set to NULL at conclusion
3286ea844a1aSMatthew Knepley - rots - The field rotations symmetries for the given orientations: set to NULL at conclusion
3287ea844a1aSMatthew Knepley 
3288ea844a1aSMatthew Knepley   Level: developer
3289ea844a1aSMatthew Knepley 
3290ea844a1aSMatthew Knepley .seealso: PetscSectionGetPointSyms(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym()
3291ea844a1aSMatthew Knepley @*/
3292ea844a1aSMatthew Knepley PetscErrorCode PetscSectionRestorePointSyms(PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots)
3293ea844a1aSMatthew Knepley {
3294ea844a1aSMatthew Knepley   PetscSectionSym sym;
3295ea844a1aSMatthew Knepley 
3296ea844a1aSMatthew Knepley   PetscFunctionBegin;
3297ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1);
3298ea844a1aSMatthew Knepley   sym = section->sym;
3299ea844a1aSMatthew Knepley   if (sym && (perms || rots)) {
3300ea844a1aSMatthew Knepley     SymWorkLink *p,link;
3301ea844a1aSMatthew Knepley 
3302ea844a1aSMatthew Knepley     for (p=&sym->workout; (link=*p); p=&link->next) {
3303ea844a1aSMatthew Knepley       if ((perms && link->perms == *perms) || (rots && link->rots == *rots)) {
3304ea844a1aSMatthew Knepley         *p          = link->next;
3305ea844a1aSMatthew Knepley         link->next  = sym->workin;
3306ea844a1aSMatthew Knepley         sym->workin = link;
3307ea844a1aSMatthew Knepley         if (perms) *perms = NULL;
3308ea844a1aSMatthew Knepley         if (rots)  *rots  = NULL;
3309ea844a1aSMatthew Knepley         PetscFunctionReturn(0);
3310ea844a1aSMatthew Knepley       }
3311ea844a1aSMatthew Knepley     }
3312ea844a1aSMatthew Knepley     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Array was not checked out");
3313ea844a1aSMatthew Knepley   }
3314ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3315ea844a1aSMatthew Knepley }
3316ea844a1aSMatthew Knepley 
3317ea844a1aSMatthew Knepley /*@C
3318ea844a1aSMatthew Knepley   PetscSectionGetFieldPointSyms - Get the symmetries for a set of points in a field of a PetscSection under specific orientations.
3319ea844a1aSMatthew Knepley 
3320ea844a1aSMatthew Knepley   Not collective
3321ea844a1aSMatthew Knepley 
3322ea844a1aSMatthew Knepley   Input Parameters:
3323ea844a1aSMatthew Knepley + section - the section
3324ea844a1aSMatthew Knepley . field - the field of the section
3325ea844a1aSMatthew Knepley . numPoints - the number of points
3326ea844a1aSMatthew Knepley - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an
3327ea844a1aSMatthew Knepley     arbitrary integer: its interpretation is up to sym.  Orientations are used by DM: for their interpretation in that
3328ea844a1aSMatthew Knepley     context, see DMPlexGetConeOrientation()).
3329ea844a1aSMatthew Knepley 
3330ea844a1aSMatthew Knepley   Output Parameter:
3331ea844a1aSMatthew Knepley + perms - The permutations for the given orientations (or NULL if there is no symmetry or the permutation is the identity).
3332ea844a1aSMatthew Knepley - rots - The field rotations symmetries for the given orientations (or NULL if there is no symmetry or the rotations are all
3333ea844a1aSMatthew Knepley     identity).
3334ea844a1aSMatthew Knepley 
3335ea844a1aSMatthew Knepley   Level: developer
3336ea844a1aSMatthew Knepley 
3337ea844a1aSMatthew Knepley .seealso: PetscSectionGetPointSyms(), PetscSectionRestoreFieldPointSyms()
3338ea844a1aSMatthew Knepley @*/
3339ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetFieldPointSyms(PetscSection section, PetscInt field, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots)
3340ea844a1aSMatthew Knepley {
3341ea844a1aSMatthew Knepley   PetscErrorCode ierr;
3342ea844a1aSMatthew Knepley 
3343ea844a1aSMatthew Knepley   PetscFunctionBegin;
3344ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1);
3345ea844a1aSMatthew Knepley   if (field > section->numFields) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"field %D greater than number of fields (%D) in section",field,section->numFields);
3346ea844a1aSMatthew Knepley   ierr = PetscSectionGetPointSyms(section->field[field],numPoints,points,perms,rots);CHKERRQ(ierr);
3347ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3348ea844a1aSMatthew Knepley }
3349ea844a1aSMatthew Knepley 
3350ea844a1aSMatthew Knepley /*@C
3351ea844a1aSMatthew Knepley   PetscSectionRestoreFieldPointSyms - Restore the symmetries returned by PetscSectionGetFieldPointSyms()
3352ea844a1aSMatthew Knepley 
3353ea844a1aSMatthew Knepley   Not collective
3354ea844a1aSMatthew Knepley 
3355ea844a1aSMatthew Knepley   Input Parameters:
3356ea844a1aSMatthew Knepley + section - the section
3357ea844a1aSMatthew Knepley . field - the field number
3358ea844a1aSMatthew Knepley . numPoints - the number of points
3359ea844a1aSMatthew Knepley - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an
3360ea844a1aSMatthew Knepley     arbitrary integer: its interpretation is up to sym.  Orientations are used by DM: for their interpretation in that
3361ea844a1aSMatthew Knepley     context, see DMPlexGetConeOrientation()).
3362ea844a1aSMatthew Knepley 
3363ea844a1aSMatthew Knepley   Output Parameter:
3364ea844a1aSMatthew Knepley + perms - The permutations for the given orientations: set to NULL at conclusion
3365ea844a1aSMatthew Knepley - rots - The field rotations symmetries for the given orientations: set to NULL at conclusion
3366ea844a1aSMatthew Knepley 
3367ea844a1aSMatthew Knepley   Level: developer
3368ea844a1aSMatthew Knepley 
3369ea844a1aSMatthew Knepley .seealso: PetscSectionRestorePointSyms(), petscSectionGetFieldPointSyms(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym()
3370ea844a1aSMatthew Knepley @*/
3371ea844a1aSMatthew Knepley PetscErrorCode PetscSectionRestoreFieldPointSyms(PetscSection section, PetscInt field, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots)
3372ea844a1aSMatthew Knepley {
3373ea844a1aSMatthew Knepley   PetscErrorCode ierr;
3374ea844a1aSMatthew Knepley 
3375ea844a1aSMatthew Knepley   PetscFunctionBegin;
3376ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1);
3377ea844a1aSMatthew Knepley   if (field > section->numFields) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"field %D greater than number of fields (%D) in section",field,section->numFields);
3378ea844a1aSMatthew Knepley   ierr = PetscSectionRestorePointSyms(section->field[field],numPoints,points,perms,rots);CHKERRQ(ierr);
3379ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3380ea844a1aSMatthew Knepley }
3381ea844a1aSMatthew Knepley 
3382ea844a1aSMatthew Knepley /*@
3383ea844a1aSMatthew Knepley   PetscSectionGetUseFieldOffsets - Get the flag to use field offsets directly in a global section, rather than just the point offset
3384ea844a1aSMatthew Knepley 
3385ea844a1aSMatthew Knepley   Not collective
3386ea844a1aSMatthew Knepley 
3387ea844a1aSMatthew Knepley   Input Parameter:
3388ea844a1aSMatthew Knepley . s - the global PetscSection
3389ea844a1aSMatthew Knepley 
3390ea844a1aSMatthew Knepley   Output Parameters:
3391ea844a1aSMatthew Knepley . flg - the flag
3392ea844a1aSMatthew Knepley 
3393ea844a1aSMatthew Knepley   Level: developer
3394ea844a1aSMatthew Knepley 
3395ea844a1aSMatthew Knepley .seealso: PetscSectionSetChart(), PetscSectionCreate()
3396ea844a1aSMatthew Knepley @*/
3397ea844a1aSMatthew Knepley PetscErrorCode PetscSectionGetUseFieldOffsets(PetscSection s, PetscBool *flg)
3398ea844a1aSMatthew Knepley {
3399ea844a1aSMatthew Knepley   PetscFunctionBegin;
3400ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
3401ea844a1aSMatthew Knepley   *flg = s->useFieldOff;
3402ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3403ea844a1aSMatthew Knepley }
3404ea844a1aSMatthew Knepley 
3405ea844a1aSMatthew Knepley /*@
3406ea844a1aSMatthew Knepley   PetscSectionSetUseFieldOffsets - Set the flag to use field offsets directly in a global section, rather than just the point offset
3407ea844a1aSMatthew Knepley 
3408ea844a1aSMatthew Knepley   Not collective
3409ea844a1aSMatthew Knepley 
3410ea844a1aSMatthew Knepley   Input Parameters:
3411ea844a1aSMatthew Knepley + s   - the global PetscSection
3412ea844a1aSMatthew Knepley - flg - the flag
3413ea844a1aSMatthew Knepley 
3414ea844a1aSMatthew Knepley   Level: developer
3415ea844a1aSMatthew Knepley 
3416ea844a1aSMatthew Knepley .seealso: PetscSectionGetUseFieldOffsets(), PetscSectionSetChart(), PetscSectionCreate()
3417ea844a1aSMatthew Knepley @*/
3418ea844a1aSMatthew Knepley PetscErrorCode PetscSectionSetUseFieldOffsets(PetscSection s, PetscBool flg)
3419ea844a1aSMatthew Knepley {
3420ea844a1aSMatthew Knepley   PetscFunctionBegin;
3421ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
3422ea844a1aSMatthew Knepley   s->useFieldOff = flg;
3423ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3424ea844a1aSMatthew Knepley }
3425ea844a1aSMatthew Knepley 
3426ea844a1aSMatthew Knepley #define PetscSectionExpandPoints_Loop(TYPE) \
3427ea844a1aSMatthew Knepley { \
3428ea844a1aSMatthew Knepley   PetscInt i, n, o0, o1, size; \
3429ea844a1aSMatthew Knepley   TYPE *a0 = (TYPE*)origArray, *a1; \
3430ea844a1aSMatthew Knepley   ierr = PetscSectionGetStorageSize(s, &size);CHKERRQ(ierr); \
3431ea844a1aSMatthew Knepley   ierr = PetscMalloc1(size, &a1);CHKERRQ(ierr); \
3432ea844a1aSMatthew Knepley   for (i=0; i<npoints; i++) { \
3433ea844a1aSMatthew Knepley     ierr = PetscSectionGetOffset(origSection, points_[i], &o0);CHKERRQ(ierr); \
3434ea844a1aSMatthew Knepley     ierr = PetscSectionGetOffset(s, i, &o1);CHKERRQ(ierr); \
3435ea844a1aSMatthew Knepley     ierr = PetscSectionGetDof(s, i, &n);CHKERRQ(ierr); \
3436ea844a1aSMatthew Knepley     ierr = PetscMemcpy(&a1[o1], &a0[o0], n*unitsize);CHKERRQ(ierr); \
3437ea844a1aSMatthew Knepley   } \
3438ea844a1aSMatthew Knepley   *newArray = (void*)a1; \
3439ea844a1aSMatthew Knepley }
3440ea844a1aSMatthew Knepley 
3441ea844a1aSMatthew Knepley /*@
3442ea844a1aSMatthew Knepley   PetscSectionExtractDofsFromArray - Extracts elements of an array corresponding to DOFs of specified points.
3443ea844a1aSMatthew Knepley 
3444ea844a1aSMatthew Knepley   Not collective
3445ea844a1aSMatthew Knepley 
3446ea844a1aSMatthew Knepley   Input Parameters:
3447ea844a1aSMatthew Knepley + origSection - the PetscSection describing the layout of the array
3448ea844a1aSMatthew Knepley . dataType - MPI_Datatype describing the data type of the array (currently only MPIU_INT, MPIU_SCALAR, MPIU_REAL)
3449ea844a1aSMatthew Knepley . origArray - the array; its size must be equal to the storage size of origSection
3450ea844a1aSMatthew Knepley - points - IS with points to extract; its indices must lie in the chart of origSection
3451ea844a1aSMatthew Knepley 
3452ea844a1aSMatthew Knepley   Output Parameters:
3453ea844a1aSMatthew Knepley + newSection - the new PetscSection desribing the layout of the new array (with points renumbered 0,1,... but preserving numbers of DOFs)
3454ea844a1aSMatthew Knepley - newArray - the array of the extracted DOFs; its size is the storage size of newSection
3455ea844a1aSMatthew Knepley 
3456ea844a1aSMatthew Knepley   Level: developer
3457ea844a1aSMatthew Knepley 
3458ea844a1aSMatthew Knepley .seealso: PetscSectionGetChart(), PetscSectionGetDof(), PetscSectionGetStorageSize(), PetscSectionCreate()
3459ea844a1aSMatthew Knepley @*/
3460ea844a1aSMatthew Knepley PetscErrorCode PetscSectionExtractDofsFromArray(PetscSection origSection, MPI_Datatype dataType, const void *origArray, IS points, PetscSection *newSection, void *newArray[])
3461ea844a1aSMatthew Knepley {
3462ea844a1aSMatthew Knepley   PetscSection        s;
3463ea844a1aSMatthew Knepley   const PetscInt      *points_;
3464ea844a1aSMatthew Knepley   PetscInt            i, n, npoints, pStart, pEnd;
3465ea844a1aSMatthew Knepley   PetscMPIInt         unitsize;
3466ea844a1aSMatthew Knepley   PetscErrorCode      ierr;
3467ea844a1aSMatthew Knepley 
3468ea844a1aSMatthew Knepley   PetscFunctionBegin;
3469ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(origSection, PETSC_SECTION_CLASSID, 1);
3470ea844a1aSMatthew Knepley   PetscValidPointer(origArray, 3);
3471ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(points, IS_CLASSID, 4);
3472ea844a1aSMatthew Knepley   if (newSection) PetscValidPointer(newSection, 5);
3473ea844a1aSMatthew Knepley   if (newArray) PetscValidPointer(newArray, 6);
3474ea844a1aSMatthew Knepley   ierr = MPI_Type_size(dataType, &unitsize);CHKERRQ(ierr);
3475ea844a1aSMatthew Knepley   ierr = ISGetLocalSize(points, &npoints);CHKERRQ(ierr);
3476ea844a1aSMatthew Knepley   ierr = ISGetIndices(points, &points_);CHKERRQ(ierr);
3477ea844a1aSMatthew Knepley   ierr = PetscSectionGetChart(origSection, &pStart, &pEnd);CHKERRQ(ierr);
3478ea844a1aSMatthew Knepley   ierr = PetscSectionCreate(PETSC_COMM_SELF, &s);CHKERRQ(ierr);
3479ea844a1aSMatthew Knepley   ierr = PetscSectionSetChart(s, 0, npoints);CHKERRQ(ierr);
3480ea844a1aSMatthew Knepley   for (i=0; i<npoints; i++) {
3481ea844a1aSMatthew Knepley     if (PetscUnlikely(points_[i] < pStart || points_[i] >= pEnd)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "point %d (index %d) in input IS out of input section's chart", points_[i], i);
3482ea844a1aSMatthew Knepley     ierr = PetscSectionGetDof(origSection, points_[i], &n);CHKERRQ(ierr);
3483ea844a1aSMatthew Knepley     ierr = PetscSectionSetDof(s, i, n);CHKERRQ(ierr);
3484ea844a1aSMatthew Knepley   }
3485ea844a1aSMatthew Knepley   ierr = PetscSectionSetUp(s);CHKERRQ(ierr);
3486ea844a1aSMatthew Knepley   if (newArray) {
3487ea844a1aSMatthew Knepley     if (dataType == MPIU_INT)           {PetscSectionExpandPoints_Loop(PetscInt);}
3488ea844a1aSMatthew Knepley     else if (dataType == MPIU_SCALAR)   {PetscSectionExpandPoints_Loop(PetscScalar);}
3489ea844a1aSMatthew Knepley     else if (dataType == MPIU_REAL)     {PetscSectionExpandPoints_Loop(PetscReal);}
3490ea844a1aSMatthew Knepley     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "not implemented for this MPI_Datatype");
3491ea844a1aSMatthew Knepley   }
3492ea844a1aSMatthew Knepley   if (newSection) {
3493ea844a1aSMatthew Knepley     *newSection = s;
3494ea844a1aSMatthew Knepley   } else {
3495ea844a1aSMatthew Knepley     ierr = PetscSectionDestroy(&s);CHKERRQ(ierr);
3496ea844a1aSMatthew Knepley   }
3497ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3498ea844a1aSMatthew Knepley }
3499