xref: /petsc/src/vec/is/section/interface/section.c (revision d5b43468fb8780a8feea140ccd6fa3e6a50411cc)
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 <petscsf.h>
7ea844a1aSMatthew Knepley 
8ea844a1aSMatthew Knepley PetscClassId PETSC_SECTION_CLASSID;
9ea844a1aSMatthew Knepley 
10ea844a1aSMatthew Knepley /*@
11cab54364SBarry Smith   PetscSectionCreate - Allocates `PetscSection` and sets the map contents to the default.
12ea844a1aSMatthew Knepley 
13ea844a1aSMatthew Knepley   Collective
14ea844a1aSMatthew Knepley 
15ea844a1aSMatthew Knepley   Input Parameters:
16ea844a1aSMatthew Knepley + comm - the MPI communicator
17ea844a1aSMatthew Knepley - s    - pointer to the section
18ea844a1aSMatthew Knepley 
19ea844a1aSMatthew Knepley   Level: beginner
20ea844a1aSMatthew Knepley 
21ea844a1aSMatthew Knepley   Notes:
22ea844a1aSMatthew Knepley   Typical calling sequence
23cab54364SBarry Smith .vb
24cab54364SBarry Smith        PetscSectionCreate(MPI_Comm,PetscSection *);!
25cab54364SBarry Smith        PetscSectionSetNumFields(PetscSection, numFields);
26cab54364SBarry Smith        PetscSectionSetChart(PetscSection,low,high);
27cab54364SBarry Smith        PetscSectionSetDof(PetscSection,point,numdof);
28cab54364SBarry Smith        PetscSectionSetUp(PetscSection);
29cab54364SBarry Smith        PetscSectionGetOffset(PetscSection,point,PetscInt *);
30cab54364SBarry Smith        PetscSectionDestroy(PetscSection);
31cab54364SBarry Smith .ve
32ea844a1aSMatthew Knepley 
33cab54364SBarry Smith   The `PetscSection` object and methods are intended to be used in the PETSc `Vec` and `Mat` implementations.
34ea844a1aSMatthew Knepley 
35cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionDestroy()`
36ea844a1aSMatthew Knepley @*/
37d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionCreate(MPI_Comm comm, PetscSection *s)
38d71ae5a4SJacob Faibussowitsch {
39ea844a1aSMatthew Knepley   PetscFunctionBegin;
40ea844a1aSMatthew Knepley   PetscValidPointer(s, 2);
419566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
42ea844a1aSMatthew Knepley 
439566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*s, PETSC_SECTION_CLASSID, "PetscSection", "Section", "IS", comm, PetscSectionDestroy, PetscSectionView));
44ea844a1aSMatthew Knepley 
45ea844a1aSMatthew Knepley   (*s)->pStart              = -1;
46ea844a1aSMatthew Knepley   (*s)->pEnd                = -1;
47ea844a1aSMatthew Knepley   (*s)->perm                = NULL;
48ea844a1aSMatthew Knepley   (*s)->pointMajor          = PETSC_TRUE;
4987e637c6Sksagiyam   (*s)->includesConstraints = PETSC_TRUE;
50ea844a1aSMatthew Knepley   (*s)->atlasDof            = NULL;
51ea844a1aSMatthew Knepley   (*s)->atlasOff            = NULL;
52ea844a1aSMatthew Knepley   (*s)->bc                  = NULL;
53ea844a1aSMatthew Knepley   (*s)->bcIndices           = NULL;
54ea844a1aSMatthew Knepley   (*s)->setup               = PETSC_FALSE;
55ea844a1aSMatthew Knepley   (*s)->numFields           = 0;
56ea844a1aSMatthew Knepley   (*s)->fieldNames          = NULL;
57ea844a1aSMatthew Knepley   (*s)->field               = NULL;
58ea844a1aSMatthew Knepley   (*s)->useFieldOff         = PETSC_FALSE;
59b778fa18SValeria Barra   (*s)->compNames           = NULL;
60ea844a1aSMatthew Knepley   (*s)->clObj               = NULL;
61c459fbc1SJed Brown   (*s)->clHash              = NULL;
62ea844a1aSMatthew Knepley   (*s)->clSection           = NULL;
63ea844a1aSMatthew Knepley   (*s)->clPoints            = NULL;
6469c11d05SVaclav Hapla   PetscCall(PetscSectionInvalidateMaxDof_Internal(*s));
65ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
66ea844a1aSMatthew Knepley }
67ea844a1aSMatthew Knepley 
68ea844a1aSMatthew Knepley /*@
69cab54364SBarry Smith   PetscSectionCopy - Creates a shallow (if possible) copy of the `PetscSection`
70ea844a1aSMatthew Knepley 
71ea844a1aSMatthew Knepley   Collective
72ea844a1aSMatthew Knepley 
73ea844a1aSMatthew Knepley   Input Parameter:
74cab54364SBarry Smith . section - the `PetscSection`
75ea844a1aSMatthew Knepley 
76ea844a1aSMatthew Knepley   Output Parameter:
77ea844a1aSMatthew Knepley . newSection - the copy
78ea844a1aSMatthew Knepley 
79ea844a1aSMatthew Knepley   Level: intermediate
80ea844a1aSMatthew Knepley 
81cab54364SBarry Smith   Developer Note:
82cab54364SBarry Smith   What exactly does shallow mean in this context?
83cab54364SBarry Smith 
84cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreate()`, `PetscSectionDestroy()`
85ea844a1aSMatthew Knepley @*/
86d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionCopy(PetscSection section, PetscSection newSection)
87d71ae5a4SJacob Faibussowitsch {
88ea844a1aSMatthew Knepley   PetscSectionSym sym;
89ea844a1aSMatthew Knepley   IS              perm;
90b778fa18SValeria Barra   PetscInt        numFields, f, c, pStart, pEnd, p;
91ea844a1aSMatthew Knepley 
92ea844a1aSMatthew Knepley   PetscFunctionBegin;
93ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1);
94ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(newSection, PETSC_SECTION_CLASSID, 2);
959566063dSJacob Faibussowitsch   PetscCall(PetscSectionReset(newSection));
969566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &numFields));
979566063dSJacob Faibussowitsch   if (numFields) PetscCall(PetscSectionSetNumFields(newSection, numFields));
98ea844a1aSMatthew Knepley   for (f = 0; f < numFields; ++f) {
99b778fa18SValeria Barra     const char *fieldName = NULL, *compName = NULL;
100ea844a1aSMatthew Knepley     PetscInt    numComp = 0;
101ea844a1aSMatthew Knepley 
1029566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldName(section, f, &fieldName));
1039566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldName(newSection, f, fieldName));
1049566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldComponents(section, f, &numComp));
1059566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldComponents(newSection, f, numComp));
106b778fa18SValeria Barra     for (c = 0; c < numComp; ++c) {
1079566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetComponentName(section, f, c, &compName));
1089566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetComponentName(newSection, f, c, compName));
109b778fa18SValeria Barra     }
1109566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldSym(section, f, &sym));
1119566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldSym(newSection, f, sym));
112ea844a1aSMatthew Knepley   }
1139566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
1149566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(newSection, pStart, pEnd));
1159566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetPermutation(section, &perm));
1169566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetPermutation(newSection, perm));
1179566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetSym(section, &sym));
1189566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetSym(newSection, sym));
119ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
120ea844a1aSMatthew Knepley     PetscInt dof, cdof, fcdof = 0;
121ea844a1aSMatthew Knepley 
1229566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
1239566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(newSection, p, dof));
1249566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
1259566063dSJacob Faibussowitsch     if (cdof) PetscCall(PetscSectionSetConstraintDof(newSection, p, cdof));
126ea844a1aSMatthew Knepley     for (f = 0; f < numFields; ++f) {
1279566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldDof(section, p, f, &dof));
1289566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetFieldDof(newSection, p, f, dof));
129ea844a1aSMatthew Knepley       if (cdof) {
1309566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetFieldConstraintDof(section, p, f, &fcdof));
1319566063dSJacob Faibussowitsch         if (fcdof) PetscCall(PetscSectionSetFieldConstraintDof(newSection, p, f, fcdof));
132ea844a1aSMatthew Knepley       }
133ea844a1aSMatthew Knepley     }
134ea844a1aSMatthew Knepley   }
1359566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(newSection));
136ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
137ea844a1aSMatthew Knepley     PetscInt        off, cdof, fcdof = 0;
138ea844a1aSMatthew Knepley     const PetscInt *cInd;
139ea844a1aSMatthew Knepley 
140ea844a1aSMatthew Knepley     /* Must set offsets in case they do not agree with the prefix sums */
1419566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
1429566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetOffset(newSection, p, off));
1439566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
144ea844a1aSMatthew Knepley     if (cdof) {
1459566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintIndices(section, p, &cInd));
1469566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetConstraintIndices(newSection, p, cInd));
147ea844a1aSMatthew Knepley       for (f = 0; f < numFields; ++f) {
1489566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetFieldConstraintDof(section, p, f, &fcdof));
149ea844a1aSMatthew Knepley         if (fcdof) {
1509566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintIndices(section, p, f, &cInd));
1519566063dSJacob Faibussowitsch           PetscCall(PetscSectionSetFieldConstraintIndices(newSection, p, f, cInd));
152ea844a1aSMatthew Knepley         }
153ea844a1aSMatthew Knepley       }
154ea844a1aSMatthew Knepley     }
155ea844a1aSMatthew Knepley   }
156ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
157ea844a1aSMatthew Knepley }
158ea844a1aSMatthew Knepley 
159ea844a1aSMatthew Knepley /*@
160cab54364SBarry Smith   PetscSectionClone - Creates a shallow (if possible) copy of the `PetscSection`
161ea844a1aSMatthew Knepley 
162ea844a1aSMatthew Knepley   Collective
163ea844a1aSMatthew Knepley 
164ea844a1aSMatthew Knepley   Input Parameter:
165cab54364SBarry Smith . section - the `PetscSection`
166ea844a1aSMatthew Knepley 
167ea844a1aSMatthew Knepley   Output Parameter:
168ea844a1aSMatthew Knepley . newSection - the copy
169ea844a1aSMatthew Knepley 
170ea844a1aSMatthew Knepley   Level: beginner
171ea844a1aSMatthew Knepley 
172cab54364SBarry Smith   Developer Note:
173cab54364SBarry Smith   With standard PETSc terminology this should be called `PetscSectionDuplicate()`
174cab54364SBarry Smith 
175cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreate()`, `PetscSectionDestroy()`, `PetscSectionCopy()`
176ea844a1aSMatthew Knepley @*/
177d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionClone(PetscSection section, PetscSection *newSection)
178d71ae5a4SJacob Faibussowitsch {
179ea844a1aSMatthew Knepley   PetscFunctionBegin;
180ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1);
181ea844a1aSMatthew Knepley   PetscValidPointer(newSection, 2);
1829566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)section), newSection));
1839566063dSJacob Faibussowitsch   PetscCall(PetscSectionCopy(section, *newSection));
184ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
185ea844a1aSMatthew Knepley }
186ea844a1aSMatthew Knepley 
187ea844a1aSMatthew Knepley /*@
188cab54364SBarry Smith   PetscSectionSetFromOptions - sets parameters in a `PetscSection` from the options database
189ea844a1aSMatthew Knepley 
19040750d41SVaclav Hapla   Collective
191ea844a1aSMatthew Knepley 
192ea844a1aSMatthew Knepley   Input Parameter:
193cab54364SBarry Smith . section - the `PetscSection`
194ea844a1aSMatthew Knepley 
195ea844a1aSMatthew Knepley   Options Database:
196cab54364SBarry Smith . -petscsection_point_major - `PETSC_TRUE` for point-major order
197ea844a1aSMatthew Knepley 
198ea844a1aSMatthew Knepley   Level: intermediate
199ea844a1aSMatthew Knepley 
200cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreate()`, `PetscSectionDestroy()`
201ea844a1aSMatthew Knepley @*/
202d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetFromOptions(PetscSection s)
203d71ae5a4SJacob Faibussowitsch {
204ea844a1aSMatthew Knepley   PetscFunctionBegin;
205ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
206d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)s);
2079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-petscsection_point_major", "The for ordering, either point major or field major", "PetscSectionSetPointMajor", s->pointMajor, &s->pointMajor, NULL));
208ea844a1aSMatthew Knepley   /* process any options handlers added with PetscObjectAddOptionsHandler() */
209dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)s, PetscOptionsObject));
210d0609cedSBarry Smith   PetscOptionsEnd();
2119566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)s, NULL, "-petscsection_view"));
212ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
213ea844a1aSMatthew Knepley }
214ea844a1aSMatthew Knepley 
215ea844a1aSMatthew Knepley /*@
216ea844a1aSMatthew Knepley   PetscSectionCompare - Compares two sections
217ea844a1aSMatthew Knepley 
21840750d41SVaclav Hapla   Collective
219ea844a1aSMatthew Knepley 
2207a7aea1fSJed Brown   Input Parameters:
221cab54364SBarry Smith + s1 - the first `PetscSection`
222cab54364SBarry Smith - s2 - the second `PetscSection`
223ea844a1aSMatthew Knepley 
224ea844a1aSMatthew Knepley   Output Parameter:
225cab54364SBarry Smith . congruent - `PETSC_TRUE` if the two sections are congruent, `PETSC_FALSE` otherwise
226ea844a1aSMatthew Knepley 
227ea844a1aSMatthew Knepley   Level: intermediate
228ea844a1aSMatthew Knepley 
229cab54364SBarry Smith   Note:
230ea844a1aSMatthew Knepley   Field names are disregarded.
231ea844a1aSMatthew Knepley 
232cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreate()`, `PetscSectionCopy()`, `PetscSectionClone()`
233ea844a1aSMatthew Knepley @*/
234d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionCompare(PetscSection s1, PetscSection s2, PetscBool *congruent)
235d71ae5a4SJacob Faibussowitsch {
236ea844a1aSMatthew Knepley   PetscInt        pStart, pEnd, nfields, ncdof, nfcdof, p, f, n1, n2;
237ea844a1aSMatthew Knepley   const PetscInt *idx1, *idx2;
238ea844a1aSMatthew Knepley   IS              perm1, perm2;
239ea844a1aSMatthew Knepley   PetscBool       flg;
240ea844a1aSMatthew Knepley   PetscMPIInt     mflg;
241ea844a1aSMatthew Knepley 
242ea844a1aSMatthew Knepley   PetscFunctionBegin;
243ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s1, PETSC_SECTION_CLASSID, 1);
244ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s2, PETSC_SECTION_CLASSID, 2);
245064a246eSJacob Faibussowitsch   PetscValidBoolPointer(congruent, 3);
246ea844a1aSMatthew Knepley   flg = PETSC_FALSE;
247ea844a1aSMatthew Knepley 
2489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)s1), PetscObjectComm((PetscObject)s2), &mflg));
249ea844a1aSMatthew Knepley   if (mflg != MPI_CONGRUENT && mflg != MPI_IDENT) {
250ea844a1aSMatthew Knepley     *congruent = PETSC_FALSE;
251ea844a1aSMatthew Knepley     PetscFunctionReturn(0);
252ea844a1aSMatthew Knepley   }
253ea844a1aSMatthew Knepley 
2549566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(s1, &pStart, &pEnd));
2559566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(s2, &n1, &n2));
256ea844a1aSMatthew Knepley   if (pStart != n1 || pEnd != n2) goto not_congruent;
257ea844a1aSMatthew Knepley 
2589566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetPermutation(s1, &perm1));
2599566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetPermutation(s2, &perm2));
260ea844a1aSMatthew Knepley   if (perm1 && perm2) {
2619566063dSJacob Faibussowitsch     PetscCall(ISEqual(perm1, perm2, congruent));
262ea844a1aSMatthew Knepley     if (!(*congruent)) goto not_congruent;
263ea844a1aSMatthew Knepley   } else if (perm1 != perm2) goto not_congruent;
264ea844a1aSMatthew Knepley 
265ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
2669566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(s1, p, &n1));
2679566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(s2, p, &n2));
268ea844a1aSMatthew Knepley     if (n1 != n2) goto not_congruent;
269ea844a1aSMatthew Knepley 
2709566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(s1, p, &n1));
2719566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(s2, p, &n2));
272ea844a1aSMatthew Knepley     if (n1 != n2) goto not_congruent;
273ea844a1aSMatthew Knepley 
2749566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(s1, p, &ncdof));
2759566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(s2, p, &n2));
276ea844a1aSMatthew Knepley     if (ncdof != n2) goto not_congruent;
277ea844a1aSMatthew Knepley 
2789566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintIndices(s1, p, &idx1));
2799566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintIndices(s2, p, &idx2));
2809566063dSJacob Faibussowitsch     PetscCall(PetscArraycmp(idx1, idx2, ncdof, congruent));
281ea844a1aSMatthew Knepley     if (!(*congruent)) goto not_congruent;
282ea844a1aSMatthew Knepley   }
283ea844a1aSMatthew Knepley 
2849566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(s1, &nfields));
2859566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(s2, &n2));
286ea844a1aSMatthew Knepley   if (nfields != n2) goto not_congruent;
287ea844a1aSMatthew Knepley 
288ea844a1aSMatthew Knepley   for (f = 0; f < nfields; ++f) {
2899566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldComponents(s1, f, &n1));
2909566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldComponents(s2, f, &n2));
291ea844a1aSMatthew Knepley     if (n1 != n2) goto not_congruent;
292ea844a1aSMatthew Knepley 
293ea844a1aSMatthew Knepley     for (p = pStart; p < pEnd; ++p) {
2949566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldOffset(s1, p, f, &n1));
2959566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldOffset(s2, p, f, &n2));
296ea844a1aSMatthew Knepley       if (n1 != n2) goto not_congruent;
297ea844a1aSMatthew Knepley 
2989566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldDof(s1, p, f, &n1));
2999566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldDof(s2, p, f, &n2));
300ea844a1aSMatthew Knepley       if (n1 != n2) goto not_congruent;
301ea844a1aSMatthew Knepley 
3029566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldConstraintDof(s1, p, f, &nfcdof));
3039566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldConstraintDof(s2, p, f, &n2));
304ea844a1aSMatthew Knepley       if (nfcdof != n2) goto not_congruent;
305ea844a1aSMatthew Knepley 
3069566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldConstraintIndices(s1, p, f, &idx1));
3079566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldConstraintIndices(s2, p, f, &idx2));
3089566063dSJacob Faibussowitsch       PetscCall(PetscArraycmp(idx1, idx2, nfcdof, congruent));
309ea844a1aSMatthew Knepley       if (!(*congruent)) goto not_congruent;
310ea844a1aSMatthew Knepley     }
311ea844a1aSMatthew Knepley   }
312ea844a1aSMatthew Knepley 
313ea844a1aSMatthew Knepley   flg = PETSC_TRUE;
314ea844a1aSMatthew Knepley not_congruent:
3151c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&flg, congruent, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)s1)));
316ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
317ea844a1aSMatthew Knepley }
318ea844a1aSMatthew Knepley 
319ea844a1aSMatthew Knepley /*@
320cab54364SBarry Smith   PetscSectionGetNumFields - Returns the number of fields in a `PetscSection`, or 0 if no fields were defined.
321ea844a1aSMatthew Knepley 
32240750d41SVaclav Hapla   Not Collective
323ea844a1aSMatthew Knepley 
324ea844a1aSMatthew Knepley   Input Parameter:
325cab54364SBarry Smith . s - the `PetscSection`
326ea844a1aSMatthew Knepley 
327ea844a1aSMatthew Knepley   Output Parameter:
328ea844a1aSMatthew Knepley . numFields - the number of fields defined, or 0 if none were defined
329ea844a1aSMatthew Knepley 
330ea844a1aSMatthew Knepley   Level: intermediate
331ea844a1aSMatthew Knepley 
332cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionSetNumFields()`
333ea844a1aSMatthew Knepley @*/
334d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetNumFields(PetscSection s, PetscInt *numFields)
335d71ae5a4SJacob Faibussowitsch {
336ea844a1aSMatthew Knepley   PetscFunctionBegin;
337ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
338dadcf809SJacob Faibussowitsch   PetscValidIntPointer(numFields, 2);
339ea844a1aSMatthew Knepley   *numFields = s->numFields;
340ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
341ea844a1aSMatthew Knepley }
342ea844a1aSMatthew Knepley 
343ea844a1aSMatthew Knepley /*@
344cab54364SBarry Smith   PetscSectionSetNumFields - Sets the number of fields in a `PetscSection`
345ea844a1aSMatthew Knepley 
34640750d41SVaclav Hapla   Not Collective
347ea844a1aSMatthew Knepley 
348ea844a1aSMatthew Knepley   Input Parameters:
349ea844a1aSMatthew Knepley + s - the PetscSection
350ea844a1aSMatthew Knepley - numFields - the number of fields
351ea844a1aSMatthew Knepley 
352ea844a1aSMatthew Knepley   Level: intermediate
353ea844a1aSMatthew Knepley 
354cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection()`, `PetscSectionGetNumFields()`
355ea844a1aSMatthew Knepley @*/
356d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetNumFields(PetscSection s, PetscInt numFields)
357d71ae5a4SJacob Faibussowitsch {
358ea844a1aSMatthew Knepley   PetscInt f;
359ea844a1aSMatthew Knepley 
360ea844a1aSMatthew Knepley   PetscFunctionBegin;
361ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
36208401ef6SPierre Jolivet   PetscCheck(numFields > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "The number of fields %" PetscInt_FMT " must be positive", numFields);
3639566063dSJacob Faibussowitsch   PetscCall(PetscSectionReset(s));
364ea844a1aSMatthew Knepley 
365ea844a1aSMatthew Knepley   s->numFields = numFields;
3669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s->numFields, &s->numFieldComponents));
3679566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s->numFields, &s->fieldNames));
3689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s->numFields, &s->compNames));
3699566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s->numFields, &s->field));
370ea844a1aSMatthew Knepley   for (f = 0; f < s->numFields; ++f) {
371ea844a1aSMatthew Knepley     char name[64];
372ea844a1aSMatthew Knepley 
373ea844a1aSMatthew Knepley     s->numFieldComponents[f] = 1;
374ea844a1aSMatthew Knepley 
3759566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)s), &s->field[f]));
3769566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, 64, "Field_%" PetscInt_FMT, f));
3779566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(name, (char **)&s->fieldNames[f]));
3789566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, 64, "Component_0"));
3799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(s->numFieldComponents[f], &s->compNames[f]));
3809566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(name, (char **)&s->compNames[f][0]));
381ea844a1aSMatthew Knepley   }
382ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
383ea844a1aSMatthew Knepley }
384ea844a1aSMatthew Knepley 
385ea844a1aSMatthew Knepley /*@C
386cab54364SBarry Smith   PetscSectionGetFieldName - Returns the name of a field in the `PetscSection`
387ea844a1aSMatthew Knepley 
388ea844a1aSMatthew Knepley   Not Collective
389ea844a1aSMatthew Knepley 
390ea844a1aSMatthew Knepley   Input Parameters:
391cab54364SBarry Smith + s     - the `PetscSection`
392ea844a1aSMatthew Knepley - field - the field number
393ea844a1aSMatthew Knepley 
394ea844a1aSMatthew Knepley   Output Parameter:
395ea844a1aSMatthew Knepley . fieldName - the field name
396ea844a1aSMatthew Knepley 
397ea844a1aSMatthew Knepley   Level: intermediate
398ea844a1aSMatthew Knepley 
399cab54364SBarry Smith   Note:
400cab54364SBarry Smith   Will error if the field number is out of range
401cab54364SBarry Smith 
402cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionSetFieldName()`, `PetscSectionSetNumFields()`, `PetscSectionGetNumFields()`
403ea844a1aSMatthew Knepley @*/
404d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetFieldName(PetscSection s, PetscInt field, const char *fieldName[])
405d71ae5a4SJacob Faibussowitsch {
406ea844a1aSMatthew Knepley   PetscFunctionBegin;
407ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
408ea844a1aSMatthew Knepley   PetscValidPointer(fieldName, 3);
4092abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
410ea844a1aSMatthew Knepley   *fieldName = s->fieldNames[field];
411ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
412ea844a1aSMatthew Knepley }
413ea844a1aSMatthew Knepley 
414ea844a1aSMatthew Knepley /*@C
415cab54364SBarry Smith   PetscSectionSetFieldName - Sets the name of a field in the `PetscSection`
416ea844a1aSMatthew Knepley 
417ea844a1aSMatthew Knepley   Not Collective
418ea844a1aSMatthew Knepley 
419ea844a1aSMatthew Knepley   Input Parameters:
420cab54364SBarry Smith + s     - the `PetscSection`
421ea844a1aSMatthew Knepley . field - the field number
422ea844a1aSMatthew Knepley - fieldName - the field name
423ea844a1aSMatthew Knepley 
424ea844a1aSMatthew Knepley   Level: intermediate
425ea844a1aSMatthew Knepley 
426cab54364SBarry Smith   Note:
427cab54364SBarry Smith   Will error if the field number is out of range
428cab54364SBarry Smith 
429cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionGetFieldName()`, `PetscSectionSetNumFields()`, `PetscSectionGetNumFields()`
430ea844a1aSMatthew Knepley @*/
431d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetFieldName(PetscSection s, PetscInt field, const char fieldName[])
432d71ae5a4SJacob Faibussowitsch {
433ea844a1aSMatthew Knepley   PetscFunctionBegin;
434ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
435ea844a1aSMatthew Knepley   if (fieldName) PetscValidCharPointer(fieldName, 3);
4362abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
4379566063dSJacob Faibussowitsch   PetscCall(PetscFree(s->fieldNames[field]));
4389566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(fieldName, (char **)&s->fieldNames[field]));
439ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
440ea844a1aSMatthew Knepley }
441ea844a1aSMatthew Knepley 
442b778fa18SValeria Barra /*@C
443cab54364SBarry Smith   PetscSectionGetComponentName - Gets the name of a field component in the `PetscSection`
444b778fa18SValeria Barra 
445b778fa18SValeria Barra   Not Collective
446b778fa18SValeria Barra 
447b778fa18SValeria Barra   Input Parameters:
448cab54364SBarry Smith + s     - the `PetscSection`
449b778fa18SValeria Barra . field - the field number
450cab54364SBarry Smith - comp  - the component number
451cab54364SBarry Smith 
452*d5b43468SJose E. Roman   Output Parameter:
453cab54364SBarry Smith . compName - the component name
454b778fa18SValeria Barra 
455b778fa18SValeria Barra   Level: intermediate
456b778fa18SValeria Barra 
457cab54364SBarry Smith   Note:
458cab54364SBarry Smith   Will error if the field or component number do not exist
459cab54364SBarry Smith 
460cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetFieldName()`, `PetscSectionSetNumFields()`, `PetscSectionGetNumFields()`,
461cab54364SBarry Smith           `PetscSectionSetComponentName()`, `PetscSectionSetFieldName()`, `PetscSectionGetFieldComponents()`, `PetscSectionSetFieldComponents()`
462b778fa18SValeria Barra @*/
463d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetComponentName(PetscSection s, PetscInt field, PetscInt comp, const char *compName[])
464d71ae5a4SJacob Faibussowitsch {
465b778fa18SValeria Barra   PetscFunctionBegin;
466b778fa18SValeria Barra   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
467064a246eSJacob Faibussowitsch   PetscValidPointer(compName, 4);
4682abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
4692abc8c78SJacob Faibussowitsch   PetscSectionCheckValidFieldComponent(comp, s->numFieldComponents[field]);
470b778fa18SValeria Barra   *compName = s->compNames[field][comp];
471b778fa18SValeria Barra   PetscFunctionReturn(0);
472b778fa18SValeria Barra }
473b778fa18SValeria Barra 
474b778fa18SValeria Barra /*@C
475cab54364SBarry Smith   PetscSectionSetComponentName - Sets the name of a field component in the `PetscSection`
476b778fa18SValeria Barra 
477b778fa18SValeria Barra   Not Collective
478b778fa18SValeria Barra 
479b778fa18SValeria Barra   Input Parameters:
480b778fa18SValeria Barra + s     - the PetscSection
481b778fa18SValeria Barra . field - the field number
482b778fa18SValeria Barra . comp  - the component number
483b778fa18SValeria Barra - compName - the component name
484b778fa18SValeria Barra 
485b778fa18SValeria Barra   Level: intermediate
486b778fa18SValeria Barra 
487cab54364SBarry Smith   Note:
488cab54364SBarry Smith   Will error if the field or component number do not exist
489cab54364SBarry Smith 
490cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetComponentName()`, `PetscSectionSetNumFields()`, `PetscSectionGetNumFields()`,
491cab54364SBarry Smith           `PetscSectionSetComponentName()`, `PetscSectionSetFieldName()`, `PetscSectionGetFieldComponents()`, `PetscSectionSetFieldComponents()`
492b778fa18SValeria Barra @*/
493d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetComponentName(PetscSection s, PetscInt field, PetscInt comp, const char compName[])
494d71ae5a4SJacob Faibussowitsch {
495b778fa18SValeria Barra   PetscFunctionBegin;
496b778fa18SValeria Barra   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
497064a246eSJacob Faibussowitsch   if (compName) PetscValidCharPointer(compName, 4);
4982abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
4992abc8c78SJacob Faibussowitsch   PetscSectionCheckValidFieldComponent(comp, s->numFieldComponents[field]);
5009566063dSJacob Faibussowitsch   PetscCall(PetscFree(s->compNames[field][comp]));
5019566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(compName, (char **)&s->compNames[field][comp]));
502b778fa18SValeria Barra   PetscFunctionReturn(0);
503b778fa18SValeria Barra }
504b778fa18SValeria Barra 
505ea844a1aSMatthew Knepley /*@
506ea844a1aSMatthew Knepley   PetscSectionGetFieldComponents - Returns the number of field components for the given field.
507ea844a1aSMatthew Knepley 
50840750d41SVaclav Hapla   Not Collective
509ea844a1aSMatthew Knepley 
510ea844a1aSMatthew Knepley   Input Parameters:
511cab54364SBarry Smith + s - the `PetscSection`
512ea844a1aSMatthew Knepley - field - the field number
513ea844a1aSMatthew Knepley 
514ea844a1aSMatthew Knepley   Output Parameter:
515ea844a1aSMatthew Knepley . numComp - the number of field components
516ea844a1aSMatthew Knepley 
517ea844a1aSMatthew Knepley   Level: intermediate
518ea844a1aSMatthew Knepley 
519cab54364SBarry Smith   Developer Note:
520cab54364SBarry Smith   This function is misnamed. There is a Num in `PetscSectionGetNumFields()` but not in this name
521cab54364SBarry Smith 
522cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionSetFieldComponents()`, `PetscSectionGetNumFields()`
523ea844a1aSMatthew Knepley @*/
524d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetFieldComponents(PetscSection s, PetscInt field, PetscInt *numComp)
525d71ae5a4SJacob Faibussowitsch {
526ea844a1aSMatthew Knepley   PetscFunctionBegin;
527ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
5282abc8c78SJacob Faibussowitsch   PetscValidIntPointer(numComp, 3);
5292abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
530ea844a1aSMatthew Knepley   *numComp = s->numFieldComponents[field];
531ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
532ea844a1aSMatthew Knepley }
533ea844a1aSMatthew Knepley 
534ea844a1aSMatthew Knepley /*@
535ea844a1aSMatthew Knepley   PetscSectionSetFieldComponents - Sets the number of field components for the given field.
536ea844a1aSMatthew Knepley 
53740750d41SVaclav Hapla   Not Collective
538ea844a1aSMatthew Knepley 
539ea844a1aSMatthew Knepley   Input Parameters:
540ea844a1aSMatthew Knepley + s - the PetscSection
541ea844a1aSMatthew Knepley . field - the field number
542ea844a1aSMatthew Knepley - numComp - the number of field components
543ea844a1aSMatthew Knepley 
544ea844a1aSMatthew Knepley   Level: intermediate
545ea844a1aSMatthew Knepley 
546cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetFieldComponents()`, `PetscSectionGetNumFields()`
547ea844a1aSMatthew Knepley @*/
548d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetFieldComponents(PetscSection s, PetscInt field, PetscInt numComp)
549d71ae5a4SJacob Faibussowitsch {
550b778fa18SValeria Barra   PetscInt c;
551b778fa18SValeria Barra 
552ea844a1aSMatthew Knepley   PetscFunctionBegin;
553ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
5542abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
555b778fa18SValeria Barra   if (s->compNames) {
55648a46eb9SPierre Jolivet     for (c = 0; c < s->numFieldComponents[field]; ++c) PetscCall(PetscFree(s->compNames[field][c]));
5579566063dSJacob Faibussowitsch     PetscCall(PetscFree(s->compNames[field]));
558b778fa18SValeria Barra   }
559b778fa18SValeria Barra 
560ea844a1aSMatthew Knepley   s->numFieldComponents[field] = numComp;
561b778fa18SValeria Barra   if (numComp) {
5629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numComp, (char ***)&s->compNames[field]));
563b778fa18SValeria Barra     for (c = 0; c < numComp; ++c) {
5642abc8c78SJacob Faibussowitsch       char name[64];
5652abc8c78SJacob Faibussowitsch 
5669566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(name, 64, "%" PetscInt_FMT, c));
5679566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(name, (char **)&s->compNames[field][c]));
568b778fa18SValeria Barra     }
569b778fa18SValeria Barra   }
570ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
571ea844a1aSMatthew Knepley }
572ea844a1aSMatthew Knepley 
573ea844a1aSMatthew Knepley /*@
574cab54364SBarry Smith   PetscSectionGetChart - Returns the range [pStart, pEnd) in which points (indices) lie for this `PetscSection`
575ea844a1aSMatthew Knepley 
57640750d41SVaclav Hapla   Not Collective
577ea844a1aSMatthew Knepley 
578ea844a1aSMatthew Knepley   Input Parameter:
579cab54364SBarry Smith . s - the `PetscSection`
580ea844a1aSMatthew Knepley 
581ea844a1aSMatthew Knepley   Output Parameters:
582ea844a1aSMatthew Knepley + pStart - the first point
583ea844a1aSMatthew Knepley - pEnd - one past the last point
584ea844a1aSMatthew Knepley 
585ea844a1aSMatthew Knepley   Level: intermediate
586ea844a1aSMatthew Knepley 
587cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionSetChart()`, `PetscSectionCreate()`
588ea844a1aSMatthew Knepley @*/
589d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetChart(PetscSection s, PetscInt *pStart, PetscInt *pEnd)
590d71ae5a4SJacob Faibussowitsch {
591ea844a1aSMatthew Knepley   PetscFunctionBegin;
592ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
593ea844a1aSMatthew Knepley   if (pStart) *pStart = s->pStart;
594ea844a1aSMatthew Knepley   if (pEnd) *pEnd = s->pEnd;
595ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
596ea844a1aSMatthew Knepley }
597ea844a1aSMatthew Knepley 
598ea844a1aSMatthew Knepley /*@
599cab54364SBarry Smith   PetscSectionSetChart - Sets the range [pStart, pEnd) in which points (indices) lie for this `PetscSection`
600ea844a1aSMatthew Knepley 
60140750d41SVaclav Hapla   Not Collective
602ea844a1aSMatthew Knepley 
603ea844a1aSMatthew Knepley   Input Parameters:
604ea844a1aSMatthew Knepley + s - the PetscSection
605ea844a1aSMatthew Knepley . pStart - the first point
606ea844a1aSMatthew Knepley - pEnd - one past the last point
607ea844a1aSMatthew Knepley 
608ea844a1aSMatthew Knepley   Level: intermediate
609ea844a1aSMatthew Knepley 
610cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetChart()`, `PetscSectionCreate()`
611ea844a1aSMatthew Knepley @*/
612d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetChart(PetscSection s, PetscInt pStart, PetscInt pEnd)
613d71ae5a4SJacob Faibussowitsch {
614ea844a1aSMatthew Knepley   PetscInt f;
615ea844a1aSMatthew Knepley 
616ea844a1aSMatthew Knepley   PetscFunctionBegin;
617ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
6189318fe57SMatthew G. Knepley   if (pStart == s->pStart && pEnd == s->pEnd) PetscFunctionReturn(0);
619ea844a1aSMatthew Knepley   /* Cannot Reset() because it destroys field information */
620ea844a1aSMatthew Knepley   s->setup = PETSC_FALSE;
6219566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&s->bc));
6229566063dSJacob Faibussowitsch   PetscCall(PetscFree(s->bcIndices));
6239566063dSJacob Faibussowitsch   PetscCall(PetscFree2(s->atlasDof, s->atlasOff));
624ea844a1aSMatthew Knepley 
625ea844a1aSMatthew Knepley   s->pStart = pStart;
626ea844a1aSMatthew Knepley   s->pEnd   = pEnd;
6279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2((pEnd - pStart), &s->atlasDof, (pEnd - pStart), &s->atlasOff));
6289566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(s->atlasDof, pEnd - pStart));
62948a46eb9SPierre Jolivet   for (f = 0; f < s->numFields; ++f) PetscCall(PetscSectionSetChart(s->field[f], pStart, pEnd));
630ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
631ea844a1aSMatthew Knepley }
632ea844a1aSMatthew Knepley 
633ea844a1aSMatthew Knepley /*@
634cab54364SBarry Smith   PetscSectionGetPermutation - Returns the permutation of [0, pEnd-pStart) or NULL that was set with `PetscSectionSetPermutation()`
635ea844a1aSMatthew Knepley 
63640750d41SVaclav Hapla   Not Collective
637ea844a1aSMatthew Knepley 
638ea844a1aSMatthew Knepley   Input Parameter:
639cab54364SBarry Smith . s - the `PetscSection`
640ea844a1aSMatthew Knepley 
641ea844a1aSMatthew Knepley   Output Parameters:
642cab54364SBarry Smith . perm - The permutation as an `IS`
643ea844a1aSMatthew Knepley 
644ea844a1aSMatthew Knepley   Level: intermediate
645ea844a1aSMatthew Knepley 
646cab54364SBarry Smith .seealso: [](sec_scatter), `IS`, `PetscSection`, `PetscSectionSetPermutation()`, `PetscSectionCreate()`
647ea844a1aSMatthew Knepley @*/
648d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetPermutation(PetscSection s, IS *perm)
649d71ae5a4SJacob Faibussowitsch {
650ea844a1aSMatthew Knepley   PetscFunctionBegin;
651ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
6529371c9d4SSatish Balay   if (perm) {
6539371c9d4SSatish Balay     PetscValidPointer(perm, 2);
6549371c9d4SSatish Balay     *perm = s->perm;
6559371c9d4SSatish Balay   }
656ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
657ea844a1aSMatthew Knepley }
658ea844a1aSMatthew Knepley 
659ea844a1aSMatthew Knepley /*@
660ea844a1aSMatthew Knepley   PetscSectionSetPermutation - Sets the permutation for [0, pEnd-pStart)
661ea844a1aSMatthew Knepley 
66240750d41SVaclav Hapla   Not Collective
663ea844a1aSMatthew Knepley 
664ea844a1aSMatthew Knepley   Input Parameters:
665ea844a1aSMatthew Knepley + s - the PetscSection
666ea844a1aSMatthew Knepley - perm - the permutation of points
667ea844a1aSMatthew Knepley 
668ea844a1aSMatthew Knepley   Level: intermediate
669ea844a1aSMatthew Knepley 
670cab54364SBarry Smith   Developer Note:
671cab54364SBarry Smith   What purpose does this permutation serve?
672cab54364SBarry Smith 
673cab54364SBarry Smith .seealso: [](sec_scatter), `IS`, `PetscSection`, `PetscSectionGetPermutation()`, `PetscSectionCreate()`
674ea844a1aSMatthew Knepley @*/
675d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetPermutation(PetscSection s, IS perm)
676d71ae5a4SJacob Faibussowitsch {
677ea844a1aSMatthew Knepley   PetscFunctionBegin;
678ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
679ea844a1aSMatthew Knepley   if (perm) PetscValidHeaderSpecific(perm, IS_CLASSID, 2);
68028b400f6SJacob Faibussowitsch   PetscCheck(!s->setup, PetscObjectComm((PetscObject)s), PETSC_ERR_ARG_WRONGSTATE, "Cannot set a permutation after the section is setup");
681ea844a1aSMatthew Knepley   if (s->perm != perm) {
6829566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&s->perm));
683ea844a1aSMatthew Knepley     if (perm) {
684ea844a1aSMatthew Knepley       s->perm = perm;
6859566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)s->perm));
686ea844a1aSMatthew Knepley     }
687ea844a1aSMatthew Knepley   }
688ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
689ea844a1aSMatthew Knepley }
690ea844a1aSMatthew Knepley 
691ea844a1aSMatthew Knepley /*@
692cab54364SBarry Smith   PetscSectionGetPointMajor - Returns the flag for dof ordering, `PETSC_TRUE` if it is point major, `PETSC_FALSE` if it is field major
693ea844a1aSMatthew Knepley 
69440750d41SVaclav Hapla   Not Collective
695ea844a1aSMatthew Knepley 
696ea844a1aSMatthew Knepley   Input Parameter:
697cab54364SBarry Smith . s - the `PetscSection`
698ea844a1aSMatthew Knepley 
699ea844a1aSMatthew Knepley   Output Parameter:
700ea844a1aSMatthew Knepley . pm - the flag for point major ordering
701ea844a1aSMatthew Knepley 
702ea844a1aSMatthew Knepley   Level: intermediate
703ea844a1aSMatthew Knepley 
704cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, PetscSectionSetPointMajor()`
705ea844a1aSMatthew Knepley @*/
706d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetPointMajor(PetscSection s, PetscBool *pm)
707d71ae5a4SJacob Faibussowitsch {
708ea844a1aSMatthew Knepley   PetscFunctionBegin;
709ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
710dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(pm, 2);
711ea844a1aSMatthew Knepley   *pm = s->pointMajor;
712ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
713ea844a1aSMatthew Knepley }
714ea844a1aSMatthew Knepley 
715ea844a1aSMatthew Knepley /*@
716cab54364SBarry Smith   PetscSectionSetPointMajor - Sets the flag for dof ordering, `PETSC_TRUE` for point major, otherwise it will be field major
717ea844a1aSMatthew Knepley 
71840750d41SVaclav Hapla   Not Collective
719ea844a1aSMatthew Knepley 
720ea844a1aSMatthew Knepley   Input Parameters:
721cab54364SBarry Smith + s  - the `PetscSection`
722ea844a1aSMatthew Knepley - pm - the flag for point major ordering
723ea844a1aSMatthew Knepley 
724ea844a1aSMatthew Knepley   Level: intermediate
725ea844a1aSMatthew Knepley 
726cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetPointMajor()`
727ea844a1aSMatthew Knepley @*/
728d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetPointMajor(PetscSection s, PetscBool pm)
729d71ae5a4SJacob Faibussowitsch {
730ea844a1aSMatthew Knepley   PetscFunctionBegin;
731ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
73228b400f6SJacob Faibussowitsch   PetscCheck(!s->setup, PetscObjectComm((PetscObject)s), PETSC_ERR_ARG_WRONGSTATE, "Cannot set the dof ordering after the section is setup");
733ea844a1aSMatthew Knepley   s->pointMajor = pm;
734ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
735ea844a1aSMatthew Knepley }
736ea844a1aSMatthew Knepley 
737ea844a1aSMatthew Knepley /*@
738cab54364SBarry Smith   PetscSectionGetIncludesConstraints - Returns the flag indicating if constrained dofs were included when computing offsets in the `PetscSection`.
739cab54364SBarry Smith   The value is set with `PetscSectionSetIncludesConstraints()`
74087e637c6Sksagiyam 
74140750d41SVaclav Hapla   Not Collective
74287e637c6Sksagiyam 
74387e637c6Sksagiyam   Input Parameter:
744cab54364SBarry Smith . s - the `PetscSection`
74587e637c6Sksagiyam 
74687e637c6Sksagiyam   Output Parameter:
74787e637c6Sksagiyam . includesConstraints - the flag indicating if constrained dofs were included when computing offsets
74887e637c6Sksagiyam 
74987e637c6Sksagiyam   Level: intermediate
75087e637c6Sksagiyam 
751cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionSetIncludesConstraints()`
75287e637c6Sksagiyam @*/
753d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetIncludesConstraints(PetscSection s, PetscBool *includesConstraints)
754d71ae5a4SJacob Faibussowitsch {
75587e637c6Sksagiyam   PetscFunctionBegin;
75687e637c6Sksagiyam   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
75787e637c6Sksagiyam   PetscValidBoolPointer(includesConstraints, 2);
75887e637c6Sksagiyam   *includesConstraints = s->includesConstraints;
75987e637c6Sksagiyam   PetscFunctionReturn(0);
76087e637c6Sksagiyam }
76187e637c6Sksagiyam 
76287e637c6Sksagiyam /*@
76387e637c6Sksagiyam   PetscSectionSetIncludesConstraints - Sets the flag indicating if constrained dofs are to be included when computing offsets
76487e637c6Sksagiyam 
76540750d41SVaclav Hapla   Not Collective
76687e637c6Sksagiyam 
76787e637c6Sksagiyam   Input Parameters:
76887e637c6Sksagiyam + s  - the PetscSection
76987e637c6Sksagiyam - includesConstraints - the flag indicating if constrained dofs are to be included when computing offsets
77087e637c6Sksagiyam 
77187e637c6Sksagiyam   Level: intermediate
77287e637c6Sksagiyam 
773cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetIncludesConstraints()`
77487e637c6Sksagiyam @*/
775d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetIncludesConstraints(PetscSection s, PetscBool includesConstraints)
776d71ae5a4SJacob Faibussowitsch {
77787e637c6Sksagiyam   PetscFunctionBegin;
77887e637c6Sksagiyam   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
77928b400f6SJacob Faibussowitsch   PetscCheck(!s->setup, PetscObjectComm((PetscObject)s), PETSC_ERR_ARG_WRONGSTATE, "Cannot set includesConstraints after the section is set up");
78087e637c6Sksagiyam   s->includesConstraints = includesConstraints;
78187e637c6Sksagiyam   PetscFunctionReturn(0);
78287e637c6Sksagiyam }
78387e637c6Sksagiyam 
78487e637c6Sksagiyam /*@
785ea844a1aSMatthew Knepley   PetscSectionGetDof - Return the number of degrees of freedom associated with a given point.
786ea844a1aSMatthew Knepley 
78740750d41SVaclav Hapla   Not Collective
788ea844a1aSMatthew Knepley 
789ea844a1aSMatthew Knepley   Input Parameters:
790cab54364SBarry Smith + s - the `PetscSection`
791ea844a1aSMatthew Knepley - point - the point
792ea844a1aSMatthew Knepley 
793ea844a1aSMatthew Knepley   Output Parameter:
794ea844a1aSMatthew Knepley . numDof - the number of dof
795ea844a1aSMatthew Knepley 
796ea844a1aSMatthew Knepley   Level: intermediate
797ea844a1aSMatthew Knepley 
798cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionSetDof()`, `PetscSectionCreate()`
799ea844a1aSMatthew Knepley @*/
800d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetDof(PetscSection s, PetscInt point, PetscInt *numDof)
801d71ae5a4SJacob Faibussowitsch {
80276bd3646SJed Brown   PetscFunctionBeginHot;
803ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
804dadcf809SJacob Faibussowitsch   PetscValidIntPointer(numDof, 3);
805c8bf3acbSVaclav Hapla   PetscAssert(point >= s->pStart && point < s->pEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %" PetscInt_FMT " should be in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, s->pStart, s->pEnd);
806ea844a1aSMatthew Knepley   *numDof = s->atlasDof[point - s->pStart];
807ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
808ea844a1aSMatthew Knepley }
809ea844a1aSMatthew Knepley 
810ea844a1aSMatthew Knepley /*@
811ea844a1aSMatthew Knepley   PetscSectionSetDof - Sets the number of degrees of freedom associated with a given point.
812ea844a1aSMatthew Knepley 
81340750d41SVaclav Hapla   Not Collective
814ea844a1aSMatthew Knepley 
815ea844a1aSMatthew Knepley   Input Parameters:
816cab54364SBarry Smith + s - the `PetscSection`
817ea844a1aSMatthew Knepley . point - the point
818ea844a1aSMatthew Knepley - numDof - the number of dof
819ea844a1aSMatthew Knepley 
820ea844a1aSMatthew Knepley   Level: intermediate
821ea844a1aSMatthew Knepley 
822cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetDof()`, `PetscSectionAddDof()`, `PetscSectionCreate()`
823ea844a1aSMatthew Knepley @*/
824d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetDof(PetscSection s, PetscInt point, PetscInt numDof)
825d71ae5a4SJacob Faibussowitsch {
826ea844a1aSMatthew Knepley   PetscFunctionBegin;
827ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
828c8bf3acbSVaclav Hapla   PetscAssert(point >= s->pStart && point < s->pEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %" PetscInt_FMT " should be in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, s->pStart, s->pEnd);
829ea844a1aSMatthew Knepley   s->atlasDof[point - s->pStart] = numDof;
83069c11d05SVaclav Hapla   PetscCall(PetscSectionInvalidateMaxDof_Internal(s));
831ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
832ea844a1aSMatthew Knepley }
833ea844a1aSMatthew Knepley 
834ea844a1aSMatthew Knepley /*@
835ea844a1aSMatthew Knepley   PetscSectionAddDof - Adds to the number of degrees of freedom associated with a given point.
836ea844a1aSMatthew Knepley 
83740750d41SVaclav Hapla   Not Collective
838ea844a1aSMatthew Knepley 
839ea844a1aSMatthew Knepley   Input Parameters:
840cab54364SBarry Smith + s - the `PetscSection`
841ea844a1aSMatthew Knepley . point - the point
842ea844a1aSMatthew Knepley - numDof - the number of additional dof
843ea844a1aSMatthew Knepley 
844ea844a1aSMatthew Knepley   Level: intermediate
845ea844a1aSMatthew Knepley 
846cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetDof()`, `PetscSectionSetDof()`, `PetscSectionCreate()`
847ea844a1aSMatthew Knepley @*/
848d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionAddDof(PetscSection s, PetscInt point, PetscInt numDof)
849d71ae5a4SJacob Faibussowitsch {
85076bd3646SJed Brown   PetscFunctionBeginHot;
851ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
852c8bf3acbSVaclav Hapla   PetscAssert(point >= s->pStart && point < s->pEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %" PetscInt_FMT " should be in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, s->pStart, s->pEnd);
853ea844a1aSMatthew Knepley   s->atlasDof[point - s->pStart] += numDof;
85469c11d05SVaclav Hapla   PetscCall(PetscSectionInvalidateMaxDof_Internal(s));
855ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
856ea844a1aSMatthew Knepley }
857ea844a1aSMatthew Knepley 
858ea844a1aSMatthew Knepley /*@
859ea844a1aSMatthew Knepley   PetscSectionGetFieldDof - Return the number of degrees of freedom associated with a field on a given point.
860ea844a1aSMatthew Knepley 
86140750d41SVaclav Hapla   Not Collective
862ea844a1aSMatthew Knepley 
863ea844a1aSMatthew Knepley   Input Parameters:
864cab54364SBarry Smith + s - the `PetscSection`
865ea844a1aSMatthew Knepley . point - the point
866ea844a1aSMatthew Knepley - field - the field
867ea844a1aSMatthew Knepley 
868ea844a1aSMatthew Knepley   Output Parameter:
869ea844a1aSMatthew Knepley . numDof - the number of dof
870ea844a1aSMatthew Knepley 
871ea844a1aSMatthew Knepley   Level: intermediate
872ea844a1aSMatthew Knepley 
873cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionSetFieldDof()`, `PetscSectionCreate()`
874ea844a1aSMatthew Knepley @*/
875d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt *numDof)
876d71ae5a4SJacob Faibussowitsch {
877ea844a1aSMatthew Knepley   PetscFunctionBegin;
878ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
8792abc8c78SJacob Faibussowitsch   PetscValidIntPointer(numDof, 4);
8802abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
8819566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetDof(s->field[field], point, numDof));
882ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
883ea844a1aSMatthew Knepley }
884ea844a1aSMatthew Knepley 
885ea844a1aSMatthew Knepley /*@
886ea844a1aSMatthew Knepley   PetscSectionSetFieldDof - Sets the number of degrees of freedom associated with a field on a given point.
887ea844a1aSMatthew Knepley 
88840750d41SVaclav Hapla   Not Collective
889ea844a1aSMatthew Knepley 
890ea844a1aSMatthew Knepley   Input Parameters:
891cab54364SBarry Smith + s - the `PetscSection`
892ea844a1aSMatthew Knepley . point - the point
893ea844a1aSMatthew Knepley . field - the field
894ea844a1aSMatthew Knepley - numDof - the number of dof
895ea844a1aSMatthew Knepley 
896ea844a1aSMatthew Knepley   Level: intermediate
897ea844a1aSMatthew Knepley 
898cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetFieldDof()`, `PetscSectionCreate()`
899ea844a1aSMatthew Knepley @*/
900d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof)
901d71ae5a4SJacob Faibussowitsch {
902ea844a1aSMatthew Knepley   PetscFunctionBegin;
903ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
9042abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
9059566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetDof(s->field[field], point, numDof));
906ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
907ea844a1aSMatthew Knepley }
908ea844a1aSMatthew Knepley 
909ea844a1aSMatthew Knepley /*@
910ea844a1aSMatthew Knepley   PetscSectionAddFieldDof - Adds a number of degrees of freedom associated with a field on a given point.
911ea844a1aSMatthew Knepley 
91240750d41SVaclav Hapla   Not Collective
913ea844a1aSMatthew Knepley 
914ea844a1aSMatthew Knepley   Input Parameters:
915cab54364SBarry Smith + s - the `PetscSection`
916ea844a1aSMatthew Knepley . point - the point
917ea844a1aSMatthew Knepley . field - the field
918ea844a1aSMatthew Knepley - numDof - the number of dof
919ea844a1aSMatthew Knepley 
920ea844a1aSMatthew Knepley   Level: intermediate
921ea844a1aSMatthew Knepley 
922cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionSetFieldDof()`, `PetscSectionGetFieldDof()`, `PetscSectionCreate()`
923ea844a1aSMatthew Knepley @*/
924d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionAddFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof)
925d71ae5a4SJacob Faibussowitsch {
926ea844a1aSMatthew Knepley   PetscFunctionBegin;
927ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
9282abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
9299566063dSJacob Faibussowitsch   PetscCall(PetscSectionAddDof(s->field[field], point, numDof));
930ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
931ea844a1aSMatthew Knepley }
932ea844a1aSMatthew Knepley 
933ea844a1aSMatthew Knepley /*@
934ea844a1aSMatthew Knepley   PetscSectionGetConstraintDof - Return the number of constrained degrees of freedom associated with a given point.
935ea844a1aSMatthew Knepley 
93640750d41SVaclav Hapla   Not Collective
937ea844a1aSMatthew Knepley 
938ea844a1aSMatthew Knepley   Input Parameters:
939cab54364SBarry Smith + s - the `PetscSection`
940ea844a1aSMatthew Knepley - point - the point
941ea844a1aSMatthew Knepley 
942ea844a1aSMatthew Knepley   Output Parameter:
943ea844a1aSMatthew Knepley . numDof - the number of dof which are fixed by constraints
944ea844a1aSMatthew Knepley 
945ea844a1aSMatthew Knepley   Level: intermediate
946ea844a1aSMatthew Knepley 
947cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetDof()`, `PetscSectionSetConstraintDof()`, `PetscSectionCreate()`
948ea844a1aSMatthew Knepley @*/
949d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetConstraintDof(PetscSection s, PetscInt point, PetscInt *numDof)
950d71ae5a4SJacob Faibussowitsch {
951ea844a1aSMatthew Knepley   PetscFunctionBegin;
952ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
953dadcf809SJacob Faibussowitsch   PetscValidIntPointer(numDof, 3);
954ea844a1aSMatthew Knepley   if (s->bc) {
9559566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(s->bc, point, numDof));
956ea844a1aSMatthew Knepley   } else *numDof = 0;
957ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
958ea844a1aSMatthew Knepley }
959ea844a1aSMatthew Knepley 
960ea844a1aSMatthew Knepley /*@
961ea844a1aSMatthew Knepley   PetscSectionSetConstraintDof - Set the number of constrained degrees of freedom associated with a given point.
962ea844a1aSMatthew Knepley 
96340750d41SVaclav Hapla   Not Collective
964ea844a1aSMatthew Knepley 
965ea844a1aSMatthew Knepley   Input Parameters:
966cab54364SBarry Smith + s - the `PetscSection`
967ea844a1aSMatthew Knepley . point - the point
968ea844a1aSMatthew Knepley - numDof - the number of dof which are fixed by constraints
969ea844a1aSMatthew Knepley 
970ea844a1aSMatthew Knepley   Level: intermediate
971ea844a1aSMatthew Knepley 
972cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionSetDof()`, `PetscSectionGetConstraintDof()`, `PetscSectionCreate()`
973ea844a1aSMatthew Knepley @*/
974d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetConstraintDof(PetscSection s, PetscInt point, PetscInt numDof)
975d71ae5a4SJacob Faibussowitsch {
976ea844a1aSMatthew Knepley   PetscFunctionBegin;
977ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
978ea844a1aSMatthew Knepley   if (numDof) {
9797c0883d5SVaclav Hapla     PetscCall(PetscSectionCheckConstraints_Private(s));
9809566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(s->bc, point, numDof));
981ea844a1aSMatthew Knepley   }
982ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
983ea844a1aSMatthew Knepley }
984ea844a1aSMatthew Knepley 
985ea844a1aSMatthew Knepley /*@
986ea844a1aSMatthew Knepley   PetscSectionAddConstraintDof - Increment the number of constrained degrees of freedom associated with a given point.
987ea844a1aSMatthew Knepley 
98840750d41SVaclav Hapla   Not Collective
989ea844a1aSMatthew Knepley 
990ea844a1aSMatthew Knepley   Input Parameters:
991cab54364SBarry Smith + s - the `PetscSection`
992ea844a1aSMatthew Knepley . point - the point
993ea844a1aSMatthew Knepley - numDof - the number of additional dof which are fixed by constraints
994ea844a1aSMatthew Knepley 
995ea844a1aSMatthew Knepley   Level: intermediate
996ea844a1aSMatthew Knepley 
997cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionAddDof()`, `PetscSectionGetConstraintDof()`, `PetscSectionCreate()`
998ea844a1aSMatthew Knepley @*/
999d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionAddConstraintDof(PetscSection s, PetscInt point, PetscInt numDof)
1000d71ae5a4SJacob Faibussowitsch {
1001ea844a1aSMatthew Knepley   PetscFunctionBegin;
1002ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1003ea844a1aSMatthew Knepley   if (numDof) {
10047c0883d5SVaclav Hapla     PetscCall(PetscSectionCheckConstraints_Private(s));
10059566063dSJacob Faibussowitsch     PetscCall(PetscSectionAddDof(s->bc, point, numDof));
1006ea844a1aSMatthew Knepley   }
1007ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1008ea844a1aSMatthew Knepley }
1009ea844a1aSMatthew Knepley 
1010ea844a1aSMatthew Knepley /*@
1011ea844a1aSMatthew Knepley   PetscSectionGetFieldConstraintDof - Return the number of constrained degrees of freedom associated with a given field on a point.
1012ea844a1aSMatthew Knepley 
101340750d41SVaclav Hapla   Not Collective
1014ea844a1aSMatthew Knepley 
1015ea844a1aSMatthew Knepley   Input Parameters:
1016cab54364SBarry Smith + s - the `PetscSection`
1017ea844a1aSMatthew Knepley . point - the point
1018ea844a1aSMatthew Knepley - field - the field
1019ea844a1aSMatthew Knepley 
1020ea844a1aSMatthew Knepley   Output Parameter:
1021ea844a1aSMatthew Knepley . numDof - the number of dof which are fixed by constraints
1022ea844a1aSMatthew Knepley 
1023ea844a1aSMatthew Knepley   Level: intermediate
1024ea844a1aSMatthew Knepley 
1025cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetDof()`, `PetscSectionSetFieldConstraintDof()`, `PetscSectionCreate()`
1026ea844a1aSMatthew Knepley @*/
1027d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetFieldConstraintDof(PetscSection s, PetscInt point, PetscInt field, PetscInt *numDof)
1028d71ae5a4SJacob Faibussowitsch {
1029ea844a1aSMatthew Knepley   PetscFunctionBegin;
1030ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
10312abc8c78SJacob Faibussowitsch   PetscValidIntPointer(numDof, 4);
10322abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
10339566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstraintDof(s->field[field], point, numDof));
1034ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1035ea844a1aSMatthew Knepley }
1036ea844a1aSMatthew Knepley 
1037ea844a1aSMatthew Knepley /*@
1038ea844a1aSMatthew Knepley   PetscSectionSetFieldConstraintDof - Set the number of constrained degrees of freedom associated with a given field on a point.
1039ea844a1aSMatthew Knepley 
104040750d41SVaclav Hapla   Not Collective
1041ea844a1aSMatthew Knepley 
1042ea844a1aSMatthew Knepley   Input Parameters:
1043cab54364SBarry Smith + s - the `PetscSection`
1044ea844a1aSMatthew Knepley . point - the point
1045ea844a1aSMatthew Knepley . field - the field
1046ea844a1aSMatthew Knepley - numDof - the number of dof which are fixed by constraints
1047ea844a1aSMatthew Knepley 
1048ea844a1aSMatthew Knepley   Level: intermediate
1049ea844a1aSMatthew Knepley 
1050cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionSetDof()`, `PetscSectionGetFieldConstraintDof()`, `PetscSectionCreate()`
1051ea844a1aSMatthew Knepley @*/
1052d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetFieldConstraintDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof)
1053d71ae5a4SJacob Faibussowitsch {
1054ea844a1aSMatthew Knepley   PetscFunctionBegin;
1055ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
10562abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
10579566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetConstraintDof(s->field[field], point, numDof));
1058ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1059ea844a1aSMatthew Knepley }
1060ea844a1aSMatthew Knepley 
1061ea844a1aSMatthew Knepley /*@
1062ea844a1aSMatthew Knepley   PetscSectionAddFieldConstraintDof - Increment the number of constrained degrees of freedom associated with a given field on a point.
1063ea844a1aSMatthew Knepley 
106440750d41SVaclav Hapla   Not Collective
1065ea844a1aSMatthew Knepley 
1066ea844a1aSMatthew Knepley   Input Parameters:
1067cab54364SBarry Smith + s - the `PetscSection`
1068ea844a1aSMatthew Knepley . point - the point
1069ea844a1aSMatthew Knepley . field - the field
1070ea844a1aSMatthew Knepley - numDof - the number of additional dof which are fixed by constraints
1071ea844a1aSMatthew Knepley 
1072ea844a1aSMatthew Knepley   Level: intermediate
1073ea844a1aSMatthew Knepley 
1074cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionAddDof()`, `PetscSectionGetFieldConstraintDof()`, `PetscSectionCreate()`
1075ea844a1aSMatthew Knepley @*/
1076d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionAddFieldConstraintDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof)
1077d71ae5a4SJacob Faibussowitsch {
1078ea844a1aSMatthew Knepley   PetscFunctionBegin;
1079ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
10802abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
10819566063dSJacob Faibussowitsch   PetscCall(PetscSectionAddConstraintDof(s->field[field], point, numDof));
1082ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1083ea844a1aSMatthew Knepley }
1084ea844a1aSMatthew Knepley 
1085ea844a1aSMatthew Knepley /*@
1086ea844a1aSMatthew Knepley   PetscSectionSetUpBC - Setup the subsections describing boundary conditions.
1087ea844a1aSMatthew Knepley 
108840750d41SVaclav Hapla   Not Collective
1089ea844a1aSMatthew Knepley 
1090ea844a1aSMatthew Knepley   Input Parameter:
1091cab54364SBarry Smith . s - the `PetscSection`
1092ea844a1aSMatthew Knepley 
1093ea844a1aSMatthew Knepley   Level: advanced
1094ea844a1aSMatthew Knepley 
1095cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionSetUp()`, `PetscSectionCreate()`
1096ea844a1aSMatthew Knepley @*/
1097d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetUpBC(PetscSection s)
1098d71ae5a4SJacob Faibussowitsch {
1099ea844a1aSMatthew Knepley   PetscFunctionBegin;
1100ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1101ea844a1aSMatthew Knepley   if (s->bc) {
1102ea844a1aSMatthew Knepley     const PetscInt last = (s->bc->pEnd - s->bc->pStart) - 1;
1103ea844a1aSMatthew Knepley 
11049566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(s->bc));
11058da24d32SBarry Smith     PetscCall(PetscMalloc1((last >= 0 ? s->bc->atlasOff[last] + s->bc->atlasDof[last] : 0), &s->bcIndices));
1106ea844a1aSMatthew Knepley   }
1107ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1108ea844a1aSMatthew Knepley }
1109ea844a1aSMatthew Knepley 
1110ea844a1aSMatthew Knepley /*@
1111ea844a1aSMatthew Knepley   PetscSectionSetUp - Calculate offsets based upon the number of degrees of freedom for each point.
1112ea844a1aSMatthew Knepley 
111340750d41SVaclav Hapla   Not Collective
1114ea844a1aSMatthew Knepley 
1115ea844a1aSMatthew Knepley   Input Parameter:
1116cab54364SBarry Smith . s - the `PetscSection`
1117ea844a1aSMatthew Knepley 
1118ea844a1aSMatthew Knepley   Level: intermediate
1119ea844a1aSMatthew Knepley 
1120cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreate()`
1121ea844a1aSMatthew Knepley @*/
1122d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetUp(PetscSection s)
1123d71ae5a4SJacob Faibussowitsch {
1124ea844a1aSMatthew Knepley   const PetscInt *pind   = NULL;
1125ea844a1aSMatthew Knepley   PetscInt        offset = 0, foff, p, f;
1126ea844a1aSMatthew Knepley 
1127ea844a1aSMatthew Knepley   PetscFunctionBegin;
1128ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1129ea844a1aSMatthew Knepley   if (s->setup) PetscFunctionReturn(0);
1130ea844a1aSMatthew Knepley   s->setup = PETSC_TRUE;
1131ea844a1aSMatthew Knepley   /* Set offsets and field offsets for all points */
1132ea844a1aSMatthew Knepley   /*   Assume that all fields have the same chart */
113328b400f6SJacob Faibussowitsch   PetscCheck(s->includesConstraints, PETSC_COMM_SELF, PETSC_ERR_SUP, "PetscSectionSetUp is currently unsupported for includesConstraints = PETSC_TRUE");
11349566063dSJacob Faibussowitsch   if (s->perm) PetscCall(ISGetIndices(s->perm, &pind));
1135ea844a1aSMatthew Knepley   if (s->pointMajor) {
1136ea844a1aSMatthew Knepley     for (p = 0; p < s->pEnd - s->pStart; ++p) {
1137ea844a1aSMatthew Knepley       const PetscInt q = pind ? pind[p] : p;
1138ea844a1aSMatthew Knepley 
1139ea844a1aSMatthew Knepley       /* Set point offset */
1140ea844a1aSMatthew Knepley       s->atlasOff[q] = offset;
1141ea844a1aSMatthew Knepley       offset += s->atlasDof[q];
1142ea844a1aSMatthew Knepley       /* Set field offset */
1143ea844a1aSMatthew Knepley       for (f = 0, foff = s->atlasOff[q]; f < s->numFields; ++f) {
1144ea844a1aSMatthew Knepley         PetscSection sf = s->field[f];
1145ea844a1aSMatthew Knepley 
1146ea844a1aSMatthew Knepley         sf->atlasOff[q] = foff;
1147ea844a1aSMatthew Knepley         foff += sf->atlasDof[q];
1148ea844a1aSMatthew Knepley       }
1149ea844a1aSMatthew Knepley     }
1150ea844a1aSMatthew Knepley   } else {
1151ea844a1aSMatthew Knepley     /* Set field offsets for all points */
1152ea844a1aSMatthew Knepley     for (f = 0; f < s->numFields; ++f) {
1153ea844a1aSMatthew Knepley       PetscSection sf = s->field[f];
1154ea844a1aSMatthew Knepley 
1155ea844a1aSMatthew Knepley       for (p = 0; p < s->pEnd - s->pStart; ++p) {
1156ea844a1aSMatthew Knepley         const PetscInt q = pind ? pind[p] : p;
1157ea844a1aSMatthew Knepley 
1158ea844a1aSMatthew Knepley         sf->atlasOff[q] = offset;
1159ea844a1aSMatthew Knepley         offset += sf->atlasDof[q];
1160ea844a1aSMatthew Knepley       }
1161ea844a1aSMatthew Knepley     }
1162ea844a1aSMatthew Knepley     /* Disable point offsets since these are unused */
1163ad540459SPierre Jolivet     for (p = 0; p < s->pEnd - s->pStart; ++p) s->atlasOff[p] = -1;
1164ea844a1aSMatthew Knepley   }
11659566063dSJacob Faibussowitsch   if (s->perm) PetscCall(ISRestoreIndices(s->perm, &pind));
1166ea844a1aSMatthew Knepley   /* Setup BC sections */
11679566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUpBC(s));
11689566063dSJacob Faibussowitsch   for (f = 0; f < s->numFields; ++f) PetscCall(PetscSectionSetUpBC(s->field[f]));
1169ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1170ea844a1aSMatthew Knepley }
1171ea844a1aSMatthew Knepley 
1172ea844a1aSMatthew Knepley /*@
1173ea844a1aSMatthew Knepley   PetscSectionGetMaxDof - Return the maximum number of degrees of freedom on any point in the chart
1174ea844a1aSMatthew Knepley 
117540750d41SVaclav Hapla   Not Collective
1176ea844a1aSMatthew Knepley 
1177ea844a1aSMatthew Knepley   Input Parameters:
1178cab54364SBarry Smith . s - the `PetscSection`
1179ea844a1aSMatthew Knepley 
1180ea844a1aSMatthew Knepley   Output Parameter:
1181ea844a1aSMatthew Knepley . maxDof - the maximum dof
1182ea844a1aSMatthew Knepley 
1183ea844a1aSMatthew Knepley   Level: intermediate
1184ea844a1aSMatthew Knepley 
1185cab54364SBarry Smith   Note:
1186cab54364SBarry Smith   The returned number is up-to-date without need for `PetscSectionSetUp()`.
118769c11d05SVaclav Hapla 
1188cab54364SBarry Smith   Developer Note:
118969c11d05SVaclav Hapla   The returned number is calculated lazily and stashed.
1190cab54364SBarry Smith 
1191cab54364SBarry Smith   A call to `PetscSectionInvalidateMaxDof_Internal()` invalidates the stashed value.
1192cab54364SBarry Smith 
1193cab54364SBarry Smith   `PetscSectionInvalidateMaxDof_Internal()` is called in `PetscSectionSetDof()`, `PetscSectionAddDof()` and `PetscSectionReset()`
1194cab54364SBarry Smith 
119569c11d05SVaclav Hapla   It should also be called every time atlasDof is modified directly.
119669c11d05SVaclav Hapla 
1197cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetDof()`, `PetscSectionSetDof()`, `PetscSectionAddDof()`, `PetscSectionCreate()`
1198ea844a1aSMatthew Knepley @*/
1199d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetMaxDof(PetscSection s, PetscInt *maxDof)
1200d71ae5a4SJacob Faibussowitsch {
120169c11d05SVaclav Hapla   PetscInt p;
120269c11d05SVaclav Hapla 
1203ea844a1aSMatthew Knepley   PetscFunctionBegin;
1204ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1205dadcf809SJacob Faibussowitsch   PetscValidIntPointer(maxDof, 2);
120669c11d05SVaclav Hapla   if (s->maxDof == PETSC_MIN_INT) {
120769c11d05SVaclav Hapla     s->maxDof = 0;
1208ad540459SPierre Jolivet     for (p = 0; p < s->pEnd - s->pStart; ++p) s->maxDof = PetscMax(s->maxDof, s->atlasDof[p]);
120969c11d05SVaclav Hapla   }
1210ea844a1aSMatthew Knepley   *maxDof = s->maxDof;
1211ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1212ea844a1aSMatthew Knepley }
1213ea844a1aSMatthew Knepley 
1214ea844a1aSMatthew Knepley /*@
1215ea844a1aSMatthew Knepley   PetscSectionGetStorageSize - Return the size of an array or local Vec capable of holding all the degrees of freedom.
1216ea844a1aSMatthew Knepley 
121740750d41SVaclav Hapla   Not Collective
1218ea844a1aSMatthew Knepley 
1219ea844a1aSMatthew Knepley   Input Parameter:
1220cab54364SBarry Smith . s - the `PetscSection`
1221ea844a1aSMatthew Knepley 
1222ea844a1aSMatthew Knepley   Output Parameter:
1223ea844a1aSMatthew Knepley . size - the size of an array which can hold all the dofs
1224ea844a1aSMatthew Knepley 
1225ea844a1aSMatthew Knepley   Level: intermediate
1226ea844a1aSMatthew Knepley 
1227cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetOffset()`, `PetscSectionGetConstrainedStorageSize()`, `PetscSectionCreate()`
1228ea844a1aSMatthew Knepley @*/
1229d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetStorageSize(PetscSection s, PetscInt *size)
1230d71ae5a4SJacob Faibussowitsch {
1231ea844a1aSMatthew Knepley   PetscInt p, n = 0;
1232ea844a1aSMatthew Knepley 
1233ea844a1aSMatthew Knepley   PetscFunctionBegin;
1234ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1235dadcf809SJacob Faibussowitsch   PetscValidIntPointer(size, 2);
1236ea844a1aSMatthew Knepley   for (p = 0; p < s->pEnd - s->pStart; ++p) n += s->atlasDof[p] > 0 ? s->atlasDof[p] : 0;
1237ea844a1aSMatthew Knepley   *size = n;
1238ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1239ea844a1aSMatthew Knepley }
1240ea844a1aSMatthew Knepley 
1241ea844a1aSMatthew Knepley /*@
1242cab54364SBarry Smith   PetscSectionGetConstrainedStorageSize - Return the size of an array or local `Vec` capable of holding all unconstrained degrees of freedom.
1243ea844a1aSMatthew Knepley 
124440750d41SVaclav Hapla   Not Collective
1245ea844a1aSMatthew Knepley 
1246064ec1f3Sprj-   Input Parameter:
1247cab54364SBarry Smith . s - the `PetscSection`
1248ea844a1aSMatthew Knepley 
1249ea844a1aSMatthew Knepley   Output Parameter:
1250ea844a1aSMatthew Knepley . size - the size of an array which can hold all unconstrained dofs
1251ea844a1aSMatthew Knepley 
1252ea844a1aSMatthew Knepley   Level: intermediate
1253ea844a1aSMatthew Knepley 
1254cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetStorageSize()`, `PetscSectionGetOffset()`, `PetscSectionCreate()`
1255ea844a1aSMatthew Knepley @*/
1256d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetConstrainedStorageSize(PetscSection s, PetscInt *size)
1257d71ae5a4SJacob Faibussowitsch {
1258ea844a1aSMatthew Knepley   PetscInt p, n = 0;
1259ea844a1aSMatthew Knepley 
1260ea844a1aSMatthew Knepley   PetscFunctionBegin;
1261ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1262dadcf809SJacob Faibussowitsch   PetscValidIntPointer(size, 2);
1263ea844a1aSMatthew Knepley   for (p = 0; p < s->pEnd - s->pStart; ++p) {
1264ea844a1aSMatthew Knepley     const PetscInt cdof = s->bc ? s->bc->atlasDof[p] : 0;
1265ea844a1aSMatthew Knepley     n += s->atlasDof[p] > 0 ? s->atlasDof[p] - cdof : 0;
1266ea844a1aSMatthew Knepley   }
1267ea844a1aSMatthew Knepley   *size = n;
1268ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1269ea844a1aSMatthew Knepley }
1270ea844a1aSMatthew Knepley 
1271ea844a1aSMatthew Knepley /*@
1272ea844a1aSMatthew Knepley   PetscSectionCreateGlobalSection - Create a section describing the global field layout using
1273cab54364SBarry Smith   the local section and a `PetscSF` describing the section point overlap.
1274ea844a1aSMatthew Knepley 
1275ea844a1aSMatthew Knepley   Input Parameters:
1276cab54364SBarry Smith + s - The `PetscSection` for the local field layout
1277cab54364SBarry Smith . sf - The `PetscSF` describing parallel layout of the section points (leaves are unowned local points)
1278cab54364SBarry Smith . includeConstraints - By default this is `PETSC_FALSE`, meaning that the global field vector will not possess constrained dofs
1279cab54364SBarry Smith - localOffsets - If `PETSC_TRUE`, use local rather than global offsets for the points
1280ea844a1aSMatthew Knepley 
1281ea844a1aSMatthew Knepley   Output Parameter:
1282cab54364SBarry Smith . gsection - The `PetscSection` for the global field layout
1283ea844a1aSMatthew Knepley 
1284ea844a1aSMatthew Knepley   Level: intermediate
1285ea844a1aSMatthew Knepley 
1286cab54364SBarry Smith   Note:
1287cab54364SBarry Smith   This gives negative sizes and offsets to points not owned by this process
1288cab54364SBarry Smith 
1289cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreate()`
1290ea844a1aSMatthew Knepley @*/
1291d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionCreateGlobalSection(PetscSection s, PetscSF sf, PetscBool includeConstraints, PetscBool localOffsets, PetscSection *gsection)
1292d71ae5a4SJacob Faibussowitsch {
1293ea844a1aSMatthew Knepley   PetscSection    gs;
1294ea844a1aSMatthew Knepley   const PetscInt *pind = NULL;
1295ea844a1aSMatthew Knepley   PetscInt       *recv = NULL, *neg = NULL;
1296046d2115Sksagiyam   PetscInt        pStart, pEnd, p, dof, cdof, off, foff, globalOff = 0, nroots, nlocal, maxleaf;
1297046d2115Sksagiyam   PetscInt        numFields, f, numComponents;
1298ea844a1aSMatthew Knepley 
1299ea844a1aSMatthew Knepley   PetscFunctionBegin;
1300ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1301ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1302ea844a1aSMatthew Knepley   PetscValidLogicalCollectiveBool(s, includeConstraints, 3);
1303ea844a1aSMatthew Knepley   PetscValidLogicalCollectiveBool(s, localOffsets, 4);
1304ea844a1aSMatthew Knepley   PetscValidPointer(gsection, 5);
130528b400f6SJacob Faibussowitsch   PetscCheck(s->pointMajor, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for field major ordering");
13069566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)s), &gs));
13079566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(s, &numFields));
13089566063dSJacob Faibussowitsch   if (numFields > 0) PetscCall(PetscSectionSetNumFields(gs, numFields));
13099566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
13109566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(gs, pStart, pEnd));
131187e637c6Sksagiyam   gs->includesConstraints = includeConstraints;
13129566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
1313ea844a1aSMatthew Knepley   nlocal = nroots; /* The local/leaf space matches global/root space */
1314ea844a1aSMatthew Knepley   /* Must allocate for all points visible to SF, which may be more than this section */
1315ea844a1aSMatthew Knepley   if (nroots >= 0) { /* nroots < 0 means that the graph has not been set, only happens in serial */
13169566063dSJacob Faibussowitsch     PetscCall(PetscSFGetLeafRange(sf, NULL, &maxleaf));
131708401ef6SPierre Jolivet     PetscCheck(nroots >= pEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "SF roots %" PetscInt_FMT " < pEnd %" PetscInt_FMT, nroots, pEnd);
131808401ef6SPierre Jolivet     PetscCheck(maxleaf < nroots, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Max local leaf %" PetscInt_FMT " >= nroots %" PetscInt_FMT, maxleaf, nroots);
13199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nroots, &neg, nlocal, &recv));
13209566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(neg, nroots));
1321ea844a1aSMatthew Knepley   }
1322ea844a1aSMatthew Knepley   /* Mark all local points with negative dof */
1323ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
13249566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(s, p, &dof));
13259566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(gs, p, dof));
13269566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
13279566063dSJacob Faibussowitsch     if (!includeConstraints && cdof > 0) PetscCall(PetscSectionSetConstraintDof(gs, p, cdof));
1328ea844a1aSMatthew Knepley     if (neg) neg[p] = -(dof + 1);
1329ea844a1aSMatthew Knepley   }
13309566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUpBC(gs));
13319566063dSJacob Faibussowitsch   if (gs->bcIndices) PetscCall(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]));
1332ea844a1aSMatthew Knepley   if (nroots >= 0) {
13339566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(recv, nlocal));
13349566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, recv, MPI_REPLACE));
13359566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, recv, MPI_REPLACE));
1336ea844a1aSMatthew Knepley     for (p = pStart; p < pEnd; ++p) {
1337ea844a1aSMatthew Knepley       if (recv[p] < 0) {
1338ea844a1aSMatthew Knepley         gs->atlasDof[p - pStart] = recv[p];
13399566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(s, p, &dof));
134008401ef6SPierre Jolivet         PetscCheck(-(recv[p] + 1) == dof, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global dof %" PetscInt_FMT " for point %" PetscInt_FMT " is not the unconstrained %" PetscInt_FMT, -(recv[p] + 1), p, dof);
1341ea844a1aSMatthew Knepley       }
1342ea844a1aSMatthew Knepley     }
1343ea844a1aSMatthew Knepley   }
1344ea844a1aSMatthew Knepley   /* Calculate new sizes, get process offset, and calculate point offsets */
13459566063dSJacob Faibussowitsch   if (s->perm) PetscCall(ISGetIndices(s->perm, &pind));
1346ea844a1aSMatthew Knepley   for (p = 0, off = 0; p < pEnd - pStart; ++p) {
1347ea844a1aSMatthew Knepley     const PetscInt q = pind ? pind[p] : p;
1348ea844a1aSMatthew Knepley 
1349ea844a1aSMatthew Knepley     cdof            = (!includeConstraints && s->bc) ? s->bc->atlasDof[q] : 0;
1350ea844a1aSMatthew Knepley     gs->atlasOff[q] = off;
1351ea844a1aSMatthew Knepley     off += gs->atlasDof[q] > 0 ? gs->atlasDof[q] - cdof : 0;
1352ea844a1aSMatthew Knepley   }
1353ea844a1aSMatthew Knepley   if (!localOffsets) {
13549566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)s)));
1355ea844a1aSMatthew Knepley     globalOff -= off;
1356ea844a1aSMatthew Knepley   }
1357ea844a1aSMatthew Knepley   for (p = pStart, off = 0; p < pEnd; ++p) {
1358ea844a1aSMatthew Knepley     gs->atlasOff[p - pStart] += globalOff;
1359ea844a1aSMatthew Knepley     if (neg) neg[p] = -(gs->atlasOff[p - pStart] + 1);
1360ea844a1aSMatthew Knepley   }
13619566063dSJacob Faibussowitsch   if (s->perm) PetscCall(ISRestoreIndices(s->perm, &pind));
1362ea844a1aSMatthew Knepley   /* Put in negative offsets for ghost points */
1363ea844a1aSMatthew Knepley   if (nroots >= 0) {
13649566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(recv, nlocal));
13659566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, recv, MPI_REPLACE));
13669566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, recv, MPI_REPLACE));
1367ea844a1aSMatthew Knepley     for (p = pStart; p < pEnd; ++p) {
1368ea844a1aSMatthew Knepley       if (recv[p] < 0) gs->atlasOff[p - pStart] = recv[p];
1369ea844a1aSMatthew Knepley     }
1370ea844a1aSMatthew Knepley   }
13719566063dSJacob Faibussowitsch   PetscCall(PetscFree2(neg, recv));
1372046d2115Sksagiyam   /* Set field dofs/offsets/constraints */
1373046d2115Sksagiyam   for (f = 0; f < numFields; ++f) {
1374046d2115Sksagiyam     gs->field[f]->includesConstraints = includeConstraints;
13759566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldComponents(s, f, &numComponents));
13769566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldComponents(gs, f, numComponents));
1377046d2115Sksagiyam   }
1378046d2115Sksagiyam   for (p = pStart; p < pEnd; ++p) {
13799566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(gs, p, &off));
1380046d2115Sksagiyam     for (f = 0, foff = off; f < numFields; ++f) {
13819566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof));
13829566063dSJacob Faibussowitsch       if (!includeConstraints && cdof > 0) PetscCall(PetscSectionSetFieldConstraintDof(gs, p, f, cdof));
13839566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldDof(s, p, f, &dof));
13849566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetFieldDof(gs, p, f, off < 0 ? -(dof + 1) : dof));
13859566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetFieldOffset(gs, p, f, foff));
13869566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldConstraintDof(gs, p, f, &cdof));
1387046d2115Sksagiyam       foff = off < 0 ? foff - (dof - cdof) : foff + (dof - cdof);
1388046d2115Sksagiyam     }
1389046d2115Sksagiyam   }
1390046d2115Sksagiyam   for (f = 0; f < numFields; ++f) {
1391046d2115Sksagiyam     PetscSection gfs = gs->field[f];
1392046d2115Sksagiyam 
13939566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUpBC(gfs));
13949566063dSJacob Faibussowitsch     if (gfs->bcIndices) PetscCall(PetscArraycpy(gfs->bcIndices, s->field[f]->bcIndices, gfs->bc->atlasOff[gfs->bc->pEnd - gfs->bc->pStart - 1] + gfs->bc->atlasDof[gfs->bc->pEnd - gfs->bc->pStart - 1]));
1395046d2115Sksagiyam   }
1396046d2115Sksagiyam   gs->setup = PETSC_TRUE;
13979566063dSJacob Faibussowitsch   PetscCall(PetscSectionViewFromOptions(gs, NULL, "-global_section_view"));
1398ea844a1aSMatthew Knepley   *gsection = gs;
1399ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1400ea844a1aSMatthew Knepley }
1401ea844a1aSMatthew Knepley 
1402ea844a1aSMatthew Knepley /*@
1403cab54364SBarry Smith   PetscSectionCreateGlobalSectionCensored - Create a `PetscSection` describing the global field layout using
1404cab54364SBarry Smith   the local section and an `PetscSF` describing the section point overlap.
1405ea844a1aSMatthew Knepley 
1406ea844a1aSMatthew Knepley   Input Parameters:
1407cab54364SBarry Smith + s - The `PetscSection` for the local field layout
1408cab54364SBarry Smith . sf - The `PetscSF` describing parallel layout of the section points
1409cab54364SBarry Smith . includeConstraints - By default this is `PETSC_FALSE`, meaning that the global field vector will not possess constrained dofs
1410ea844a1aSMatthew Knepley . numExcludes - The number of exclusion ranges
1411ea844a1aSMatthew Knepley - excludes - An array [start_0, end_0, start_1, end_1, ...] where there are numExcludes pairs
1412ea844a1aSMatthew Knepley 
1413ea844a1aSMatthew Knepley   Output Parameter:
1414cab54364SBarry Smith . gsection - The `PetscSection` for the global field layout
1415ea844a1aSMatthew Knepley 
1416ea844a1aSMatthew Knepley   Level: advanced
1417ea844a1aSMatthew Knepley 
1418cab54364SBarry Smith   Note:
1419cab54364SBarry Smith   This gives negative sizes and offsets to points not owned by this process
1420cab54364SBarry Smith 
1421cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreate()`
1422ea844a1aSMatthew Knepley @*/
1423d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionCreateGlobalSectionCensored(PetscSection s, PetscSF sf, PetscBool includeConstraints, PetscInt numExcludes, const PetscInt excludes[], PetscSection *gsection)
1424d71ae5a4SJacob Faibussowitsch {
1425ea844a1aSMatthew Knepley   const PetscInt *pind = NULL;
1426ea844a1aSMatthew Knepley   PetscInt       *neg = NULL, *tmpOff = NULL;
1427ea844a1aSMatthew Knepley   PetscInt        pStart, pEnd, p, e, dof, cdof, off, globalOff = 0, nroots;
1428ea844a1aSMatthew Knepley 
1429ea844a1aSMatthew Knepley   PetscFunctionBegin;
1430ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1431ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2);
1432ea844a1aSMatthew Knepley   PetscValidPointer(gsection, 6);
14339566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)s), gsection));
14349566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
14359566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(*gsection, pStart, pEnd));
14369566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
1437ea844a1aSMatthew Knepley   if (nroots >= 0) {
143808401ef6SPierre Jolivet     PetscCheck(nroots >= pEnd - pStart, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %" PetscInt_FMT " < %" PetscInt_FMT " section size", nroots, pEnd - pStart);
14399566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(nroots, &neg));
1440ea844a1aSMatthew Knepley     if (nroots > pEnd - pStart) {
14419566063dSJacob Faibussowitsch       PetscCall(PetscCalloc1(nroots, &tmpOff));
1442ea844a1aSMatthew Knepley     } else {
1443ea844a1aSMatthew Knepley       tmpOff = &(*gsection)->atlasDof[-pStart];
1444ea844a1aSMatthew Knepley     }
1445ea844a1aSMatthew Knepley   }
1446ea844a1aSMatthew Knepley   /* Mark ghost points with negative dof */
1447ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
1448ea844a1aSMatthew Knepley     for (e = 0; e < numExcludes; ++e) {
1449ea844a1aSMatthew Knepley       if ((p >= excludes[e * 2 + 0]) && (p < excludes[e * 2 + 1])) {
14509566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetDof(*gsection, p, 0));
1451ea844a1aSMatthew Knepley         break;
1452ea844a1aSMatthew Knepley       }
1453ea844a1aSMatthew Knepley     }
1454ea844a1aSMatthew Knepley     if (e < numExcludes) continue;
14559566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(s, p, &dof));
14569566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(*gsection, p, dof));
14579566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
14589566063dSJacob Faibussowitsch     if (!includeConstraints && cdof > 0) PetscCall(PetscSectionSetConstraintDof(*gsection, p, cdof));
1459ea844a1aSMatthew Knepley     if (neg) neg[p] = -(dof + 1);
1460ea844a1aSMatthew Knepley   }
14619566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUpBC(*gsection));
1462ea844a1aSMatthew Knepley   if (nroots >= 0) {
14639566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE));
14649566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE));
1465ea844a1aSMatthew Knepley     if (nroots > pEnd - pStart) {
14669371c9d4SSatish Balay       for (p = pStart; p < pEnd; ++p) {
14679371c9d4SSatish Balay         if (tmpOff[p] < 0) (*gsection)->atlasDof[p - pStart] = tmpOff[p];
14689371c9d4SSatish Balay       }
1469ea844a1aSMatthew Knepley     }
1470ea844a1aSMatthew Knepley   }
1471ea844a1aSMatthew Knepley   /* Calculate new sizes, get proccess offset, and calculate point offsets */
14729566063dSJacob Faibussowitsch   if (s->perm) PetscCall(ISGetIndices(s->perm, &pind));
1473ea844a1aSMatthew Knepley   for (p = 0, off = 0; p < pEnd - pStart; ++p) {
1474ea844a1aSMatthew Knepley     const PetscInt q = pind ? pind[p] : p;
1475ea844a1aSMatthew Knepley 
1476ea844a1aSMatthew Knepley     cdof                     = (!includeConstraints && s->bc) ? s->bc->atlasDof[q] : 0;
1477ea844a1aSMatthew Knepley     (*gsection)->atlasOff[q] = off;
1478ea844a1aSMatthew Knepley     off += (*gsection)->atlasDof[q] > 0 ? (*gsection)->atlasDof[q] - cdof : 0;
1479ea844a1aSMatthew Knepley   }
14809566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)s)));
1481ea844a1aSMatthew Knepley   globalOff -= off;
1482ea844a1aSMatthew Knepley   for (p = 0, off = 0; p < pEnd - pStart; ++p) {
1483ea844a1aSMatthew Knepley     (*gsection)->atlasOff[p] += globalOff;
1484ea844a1aSMatthew Knepley     if (neg) neg[p + pStart] = -((*gsection)->atlasOff[p] + 1);
1485ea844a1aSMatthew Knepley   }
14869566063dSJacob Faibussowitsch   if (s->perm) PetscCall(ISRestoreIndices(s->perm, &pind));
1487ea844a1aSMatthew Knepley   /* Put in negative offsets for ghost points */
1488ea844a1aSMatthew Knepley   if (nroots >= 0) {
1489ea844a1aSMatthew Knepley     if (nroots == pEnd - pStart) tmpOff = &(*gsection)->atlasOff[-pStart];
14909566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE));
14919566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff, MPI_REPLACE));
1492ea844a1aSMatthew Knepley     if (nroots > pEnd - pStart) {
14939371c9d4SSatish Balay       for (p = pStart; p < pEnd; ++p) {
14949371c9d4SSatish Balay         if (tmpOff[p] < 0) (*gsection)->atlasOff[p - pStart] = tmpOff[p];
14959371c9d4SSatish Balay       }
1496ea844a1aSMatthew Knepley     }
1497ea844a1aSMatthew Knepley   }
14989566063dSJacob Faibussowitsch   if (nroots >= 0 && nroots > pEnd - pStart) PetscCall(PetscFree(tmpOff));
14999566063dSJacob Faibussowitsch   PetscCall(PetscFree(neg));
1500ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1501ea844a1aSMatthew Knepley }
1502ea844a1aSMatthew Knepley 
1503ea844a1aSMatthew Knepley /*@
1504cab54364SBarry Smith   PetscSectionGetPointLayout - Get the `PetscLayout` associated with the section points
1505ea844a1aSMatthew Knepley 
1506cab54364SBarry Smith   Collective
1507ea844a1aSMatthew Knepley 
1508ea844a1aSMatthew Knepley   Input Parameters:
1509cab54364SBarry Smith + comm - The `MPI_Comm`
1510cab54364SBarry Smith - s    - The `PetscSection`
1511ea844a1aSMatthew Knepley 
1512ea844a1aSMatthew Knepley   Output Parameter:
1513ea844a1aSMatthew Knepley . layout - The point layout for the section
1514ea844a1aSMatthew Knepley 
1515ea844a1aSMatthew Knepley   Level: advanced
1516ea844a1aSMatthew Knepley 
1517cab54364SBarry Smith   Note:
1518cab54364SBarry Smith   This is usually called for the default global section.
1519cab54364SBarry Smith 
1520cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetValueLayout()`, `PetscSectionCreate()`
1521ea844a1aSMatthew Knepley @*/
1522d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetPointLayout(MPI_Comm comm, PetscSection s, PetscLayout *layout)
1523d71ae5a4SJacob Faibussowitsch {
1524ea844a1aSMatthew Knepley   PetscInt pStart, pEnd, p, localSize = 0;
1525ea844a1aSMatthew Knepley 
1526ea844a1aSMatthew Knepley   PetscFunctionBegin;
15279566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
1528ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
1529ea844a1aSMatthew Knepley     PetscInt dof;
1530ea844a1aSMatthew Knepley 
15319566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(s, p, &dof));
15329120ba30SVaclav Hapla     if (dof >= 0) ++localSize;
1533ea844a1aSMatthew Knepley   }
15349566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, layout));
15359566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(*layout, localSize));
15369566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(*layout, 1));
15379566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(*layout));
1538ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1539ea844a1aSMatthew Knepley }
1540ea844a1aSMatthew Knepley 
1541ea844a1aSMatthew Knepley /*@
1542cab54364SBarry Smith   PetscSectionGetValueLayout - Get the `PetscLayout` associated with the section dofs.
1543ea844a1aSMatthew Knepley 
1544cab54364SBarry Smith   Collective
1545ea844a1aSMatthew Knepley 
1546ea844a1aSMatthew Knepley   Input Parameters:
1547cab54364SBarry Smith + comm - The `MPI_Comm`
1548cab54364SBarry Smith - s    - The `PetscSection`
1549ea844a1aSMatthew Knepley 
1550ea844a1aSMatthew Knepley   Output Parameter:
1551ea844a1aSMatthew Knepley . layout - The dof layout for the section
1552ea844a1aSMatthew Knepley 
1553ea844a1aSMatthew Knepley   Level: advanced
1554ea844a1aSMatthew Knepley 
1555cab54364SBarry Smith   Note:
1556cab54364SBarry Smith   This is usually called for the default global section.
1557cab54364SBarry Smith 
1558cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetPointLayout()`, `PetscSectionCreate()`
1559ea844a1aSMatthew Knepley @*/
1560d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetValueLayout(MPI_Comm comm, PetscSection s, PetscLayout *layout)
1561d71ae5a4SJacob Faibussowitsch {
1562ea844a1aSMatthew Knepley   PetscInt pStart, pEnd, p, localSize = 0;
1563ea844a1aSMatthew Knepley 
1564ea844a1aSMatthew Knepley   PetscFunctionBegin;
1565ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 2);
1566ea844a1aSMatthew Knepley   PetscValidPointer(layout, 3);
15679566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
1568ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
1569ea844a1aSMatthew Knepley     PetscInt dof, cdof;
1570ea844a1aSMatthew Knepley 
15719566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(s, p, &dof));
15729566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
1573ea844a1aSMatthew Knepley     if (dof - cdof > 0) localSize += dof - cdof;
1574ea844a1aSMatthew Knepley   }
15759566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, layout));
15769566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(*layout, localSize));
15779566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(*layout, 1));
15789566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(*layout));
1579ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1580ea844a1aSMatthew Knepley }
1581ea844a1aSMatthew Knepley 
1582ea844a1aSMatthew Knepley /*@
1583cab54364SBarry Smith   PetscSectionGetOffset - Return the offset into an array or local `Vec` for the dof associated with the given point.
1584ea844a1aSMatthew Knepley 
158540750d41SVaclav Hapla   Not Collective
1586ea844a1aSMatthew Knepley 
1587ea844a1aSMatthew Knepley   Input Parameters:
1588cab54364SBarry Smith + s - the `PetscSection`
1589ea844a1aSMatthew Knepley - point - the point
1590ea844a1aSMatthew Knepley 
1591ea844a1aSMatthew Knepley   Output Parameter:
1592ea844a1aSMatthew Knepley . offset - the offset
1593ea844a1aSMatthew Knepley 
1594ea844a1aSMatthew Knepley   Level: intermediate
1595ea844a1aSMatthew Knepley 
1596cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetFieldOffset()`, `PetscSectionCreate()`
1597ea844a1aSMatthew Knepley @*/
1598d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetOffset(PetscSection s, PetscInt point, PetscInt *offset)
1599d71ae5a4SJacob Faibussowitsch {
1600ea844a1aSMatthew Knepley   PetscFunctionBegin;
1601ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1602dadcf809SJacob Faibussowitsch   PetscValidIntPointer(offset, 3);
1603ad540459SPierre Jolivet   if (PetscDefined(USE_DEBUG)) PetscCheck(!(point < s->pStart) && !(point >= s->pEnd), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %" PetscInt_FMT " should be in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, s->pStart, s->pEnd);
1604ea844a1aSMatthew Knepley   *offset = s->atlasOff[point - s->pStart];
1605ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1606ea844a1aSMatthew Knepley }
1607ea844a1aSMatthew Knepley 
1608ea844a1aSMatthew Knepley /*@
1609cab54364SBarry Smith   PetscSectionSetOffset - Set the offset into an array or local `Vec` for the dof associated with the given point.
1610ea844a1aSMatthew Knepley 
161140750d41SVaclav Hapla   Not Collective
1612ea844a1aSMatthew Knepley 
1613ea844a1aSMatthew Knepley   Input Parameters:
1614cab54364SBarry Smith + s - the `PetscSection`
1615ea844a1aSMatthew Knepley . point - the point
1616ea844a1aSMatthew Knepley - offset - the offset
1617ea844a1aSMatthew Knepley 
1618ea844a1aSMatthew Knepley   Level: intermediate
1619ea844a1aSMatthew Knepley 
1620cab54364SBarry Smith   Note:
1621cab54364SBarry Smith   The user usually does not call this function, but uses `PetscSectionSetUp()`
1622cab54364SBarry Smith 
1623cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetFieldOffset()`, `PetscSectionCreate()`, `PetscSectionSetUp()`
1624ea844a1aSMatthew Knepley @*/
1625d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetOffset(PetscSection s, PetscInt point, PetscInt offset)
1626d71ae5a4SJacob Faibussowitsch {
1627ea844a1aSMatthew Knepley   PetscFunctionBegin;
1628ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
162908401ef6SPierre Jolivet   PetscCheck(!(point < s->pStart) && !(point >= s->pEnd), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %" PetscInt_FMT " should be in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, s->pStart, s->pEnd);
1630ea844a1aSMatthew Knepley   s->atlasOff[point - s->pStart] = offset;
1631ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1632ea844a1aSMatthew Knepley }
1633ea844a1aSMatthew Knepley 
1634ea844a1aSMatthew Knepley /*@
1635cab54364SBarry Smith   PetscSectionGetFieldOffset - Return the offset into an array or local `Vec` for the dof associated with the given point.
1636ea844a1aSMatthew Knepley 
163740750d41SVaclav Hapla   Not Collective
1638ea844a1aSMatthew Knepley 
1639ea844a1aSMatthew Knepley   Input Parameters:
1640cab54364SBarry Smith + s - the `PetscSection`
1641ea844a1aSMatthew Knepley . point - the point
1642ea844a1aSMatthew Knepley - field - the field
1643ea844a1aSMatthew Knepley 
1644ea844a1aSMatthew Knepley   Output Parameter:
1645ea844a1aSMatthew Knepley . offset - the offset
1646ea844a1aSMatthew Knepley 
1647ea844a1aSMatthew Knepley   Level: intermediate
1648ea844a1aSMatthew Knepley 
1649cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetOffset()`, `PetscSectionCreate()`
1650ea844a1aSMatthew Knepley @*/
1651d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetFieldOffset(PetscSection s, PetscInt point, PetscInt field, PetscInt *offset)
1652d71ae5a4SJacob Faibussowitsch {
1653ea844a1aSMatthew Knepley   PetscFunctionBegin;
1654ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
16552abc8c78SJacob Faibussowitsch   PetscValidIntPointer(offset, 4);
16562abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
16579566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(s->field[field], point, offset));
1658ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1659ea844a1aSMatthew Knepley }
1660ea844a1aSMatthew Knepley 
1661ea844a1aSMatthew Knepley /*@
1662cab54364SBarry Smith   PetscSectionSetFieldOffset - Set the offset into an array or local `Vec` for the dof associated with the given field at a point.
1663ea844a1aSMatthew Knepley 
166440750d41SVaclav Hapla   Not Collective
1665ea844a1aSMatthew Knepley 
1666ea844a1aSMatthew Knepley   Input Parameters:
1667ea844a1aSMatthew Knepley + s - the PetscSection
1668ea844a1aSMatthew Knepley . point - the point
1669ea844a1aSMatthew Knepley . field - the field
1670ea844a1aSMatthew Knepley - offset - the offset
1671ea844a1aSMatthew Knepley 
1672cab54364SBarry Smith   Level: developer
1673ea844a1aSMatthew Knepley 
1674cab54364SBarry Smith   Note:
1675cab54364SBarry Smith   The user usually does not call this function, but uses `PetscSectionSetUp()`
1676ea844a1aSMatthew Knepley 
1677cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetFieldOffset()`, `PetscSectionSetOffset()`, `PetscSectionCreate()`, `PetscSectionSetUp()`
1678ea844a1aSMatthew Knepley @*/
1679d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetFieldOffset(PetscSection s, PetscInt point, PetscInt field, PetscInt offset)
1680d71ae5a4SJacob Faibussowitsch {
1681ea844a1aSMatthew Knepley   PetscFunctionBegin;
1682ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
16832abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
16849566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetOffset(s->field[field], point, offset));
1685ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1686ea844a1aSMatthew Knepley }
1687ea844a1aSMatthew Knepley 
1688ea844a1aSMatthew Knepley /*@
1689ea844a1aSMatthew Knepley   PetscSectionGetFieldPointOffset - Return the offset on the given point for the dof associated with the given point.
1690ea844a1aSMatthew Knepley 
169140750d41SVaclav Hapla   Not Collective
1692ea844a1aSMatthew Knepley 
1693ea844a1aSMatthew Knepley   Input Parameters:
1694cab54364SBarry Smith + s - the `PetscSection`
1695ea844a1aSMatthew Knepley . point - the point
1696ea844a1aSMatthew Knepley - field - the field
1697ea844a1aSMatthew Knepley 
1698ea844a1aSMatthew Knepley   Output Parameter:
1699ea844a1aSMatthew Knepley . offset - the offset
1700ea844a1aSMatthew Knepley 
1701ea844a1aSMatthew Knepley   Level: advanced
1702ea844a1aSMatthew Knepley 
1703cab54364SBarry Smith   Note:
1704cab54364SBarry Smith   This gives the offset on a point of the field, ignoring constraints, meaning starting at the first dof for
1705cab54364SBarry Smith   this point, what number is the first dof with this field.
1706cab54364SBarry Smith 
1707cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetOffset()`, `PetscSectionCreate()`
1708ea844a1aSMatthew Knepley @*/
1709d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetFieldPointOffset(PetscSection s, PetscInt point, PetscInt field, PetscInt *offset)
1710d71ae5a4SJacob Faibussowitsch {
1711ea844a1aSMatthew Knepley   PetscInt off, foff;
1712ea844a1aSMatthew Knepley 
1713ea844a1aSMatthew Knepley   PetscFunctionBegin;
1714ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
17152abc8c78SJacob Faibussowitsch   PetscValidIntPointer(offset, 4);
17162abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
17179566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(s, point, &off));
17189566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(s->field[field], point, &foff));
1719ea844a1aSMatthew Knepley   *offset = foff - off;
1720ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1721ea844a1aSMatthew Knepley }
1722ea844a1aSMatthew Knepley 
1723ea844a1aSMatthew Knepley /*@
1724ea844a1aSMatthew Knepley   PetscSectionGetOffsetRange - Return the full range of offsets [start, end)
1725ea844a1aSMatthew Knepley 
172640750d41SVaclav Hapla   Not Collective
1727ea844a1aSMatthew Knepley 
1728ea844a1aSMatthew Knepley   Input Parameter:
1729cab54364SBarry Smith . s - the `PetscSection`
1730ea844a1aSMatthew Knepley 
1731ea844a1aSMatthew Knepley   Output Parameters:
1732ea844a1aSMatthew Knepley + start - the minimum offset
1733ea844a1aSMatthew Knepley - end   - one more than the maximum offset
1734ea844a1aSMatthew Knepley 
1735ea844a1aSMatthew Knepley   Level: intermediate
1736ea844a1aSMatthew Knepley 
1737cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetOffset()`, `PetscSectionCreate()`
1738ea844a1aSMatthew Knepley @*/
1739d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetOffsetRange(PetscSection s, PetscInt *start, PetscInt *end)
1740d71ae5a4SJacob Faibussowitsch {
1741ea844a1aSMatthew Knepley   PetscInt os = 0, oe = 0, pStart, pEnd, p;
1742ea844a1aSMatthew Knepley 
1743ea844a1aSMatthew Knepley   PetscFunctionBegin;
1744ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
17459371c9d4SSatish Balay   if (s->atlasOff) {
17469371c9d4SSatish Balay     os = s->atlasOff[0];
17479371c9d4SSatish Balay     oe = s->atlasOff[0];
17489371c9d4SSatish Balay   }
17499566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
1750ea844a1aSMatthew Knepley   for (p = 0; p < pEnd - pStart; ++p) {
1751ea844a1aSMatthew Knepley     PetscInt dof = s->atlasDof[p], off = s->atlasOff[p];
1752ea844a1aSMatthew Knepley 
1753ea844a1aSMatthew Knepley     if (off >= 0) {
1754ea844a1aSMatthew Knepley       os = PetscMin(os, off);
1755ea844a1aSMatthew Knepley       oe = PetscMax(oe, off + dof);
1756ea844a1aSMatthew Knepley     }
1757ea844a1aSMatthew Knepley   }
1758ea844a1aSMatthew Knepley   if (start) *start = os;
1759ea844a1aSMatthew Knepley   if (end) *end = oe;
1760ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1761ea844a1aSMatthew Knepley }
1762ea844a1aSMatthew Knepley 
1763ea844a1aSMatthew Knepley /*@
1764ea844a1aSMatthew Knepley   PetscSectionCreateSubsection - Create a new, smaller section composed of only the selected fields
1765ea844a1aSMatthew Knepley 
176640750d41SVaclav Hapla   Collective
1767ea844a1aSMatthew Knepley 
1768d8d19677SJose E. Roman   Input Parameters:
1769cab54364SBarry Smith + s      - the `PetscSection`
1770ea844a1aSMatthew Knepley . len    - the number of subfields
1771ea844a1aSMatthew Knepley - fields - the subfield numbers
1772ea844a1aSMatthew Knepley 
1773ea844a1aSMatthew Knepley   Output Parameter:
1774ea844a1aSMatthew Knepley . subs   - the subsection
1775ea844a1aSMatthew Knepley 
1776ea844a1aSMatthew Knepley   Level: advanced
1777ea844a1aSMatthew Knepley 
1778cab54364SBarry Smith   Notes:
1779cab54364SBarry Smith   The section offsets now refer to a new, smaller vector.
1780cab54364SBarry Smith 
1781cab54364SBarry Smith   This will error if a fieldnumber is out of range
1782cab54364SBarry Smith 
1783cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreateSupersection()`, `PetscSectionCreate()`
1784ea844a1aSMatthew Knepley @*/
1785d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionCreateSubsection(PetscSection s, PetscInt len, const PetscInt fields[], PetscSection *subs)
1786d71ae5a4SJacob Faibussowitsch {
1787b778fa18SValeria Barra   PetscInt nF, f, c, pStart, pEnd, p, maxCdof = 0;
1788ea844a1aSMatthew Knepley 
1789ea844a1aSMatthew Knepley   PetscFunctionBegin;
1790ea844a1aSMatthew Knepley   if (!len) PetscFunctionReturn(0);
1791ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1792dadcf809SJacob Faibussowitsch   PetscValidIntPointer(fields, 3);
1793ea844a1aSMatthew Knepley   PetscValidPointer(subs, 4);
17949566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(s, &nF));
179508401ef6SPierre Jolivet   PetscCheck(len <= nF, PetscObjectComm((PetscObject)s), PETSC_ERR_ARG_WRONG, "Number of requested fields %" PetscInt_FMT " greater than number of fields %" PetscInt_FMT, len, nF);
17969566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)s), subs));
17979566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetNumFields(*subs, len));
1798ea844a1aSMatthew Knepley   for (f = 0; f < len; ++f) {
1799ea844a1aSMatthew Knepley     const char *name    = NULL;
1800ea844a1aSMatthew Knepley     PetscInt    numComp = 0;
1801ea844a1aSMatthew Knepley 
18029566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldName(s, fields[f], &name));
18039566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldName(*subs, f, name));
18049566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldComponents(s, fields[f], &numComp));
18059566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldComponents(*subs, f, numComp));
1806b778fa18SValeria Barra     for (c = 0; c < s->numFieldComponents[fields[f]]; ++c) {
18079566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetComponentName(s, fields[f], c, &name));
18089566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetComponentName(*subs, f, c, name));
1809b778fa18SValeria Barra     }
1810ea844a1aSMatthew Knepley   }
18119566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
18129566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(*subs, pStart, pEnd));
1813ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
1814ea844a1aSMatthew Knepley     PetscInt dof = 0, cdof = 0, fdof = 0, cfdof = 0;
1815ea844a1aSMatthew Knepley 
1816ea844a1aSMatthew Knepley     for (f = 0; f < len; ++f) {
18179566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof));
18189566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetFieldDof(*subs, p, f, fdof));
18199566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &cfdof));
18209566063dSJacob Faibussowitsch       if (cfdof) PetscCall(PetscSectionSetFieldConstraintDof(*subs, p, f, cfdof));
1821ea844a1aSMatthew Knepley       dof += fdof;
1822ea844a1aSMatthew Knepley       cdof += cfdof;
1823ea844a1aSMatthew Knepley     }
18249566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(*subs, p, dof));
18259566063dSJacob Faibussowitsch     if (cdof) PetscCall(PetscSectionSetConstraintDof(*subs, p, cdof));
1826ea844a1aSMatthew Knepley     maxCdof = PetscMax(cdof, maxCdof);
1827ea844a1aSMatthew Knepley   }
18289566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(*subs));
1829ea844a1aSMatthew Knepley   if (maxCdof) {
1830ea844a1aSMatthew Knepley     PetscInt *indices;
1831ea844a1aSMatthew Knepley 
18329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(maxCdof, &indices));
1833ea844a1aSMatthew Knepley     for (p = pStart; p < pEnd; ++p) {
1834ea844a1aSMatthew Knepley       PetscInt cdof;
1835ea844a1aSMatthew Knepley 
18369566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintDof(*subs, p, &cdof));
1837ea844a1aSMatthew Knepley       if (cdof) {
1838ea844a1aSMatthew Knepley         const PetscInt *oldIndices = NULL;
1839ea844a1aSMatthew Knepley         PetscInt        fdof = 0, cfdof = 0, fc, numConst = 0, fOff = 0;
1840ea844a1aSMatthew Knepley 
1841ea844a1aSMatthew Knepley         for (f = 0; f < len; ++f) {
18429566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof));
18439566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &cfdof));
18449566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintIndices(s, p, fields[f], &oldIndices));
18459566063dSJacob Faibussowitsch           PetscCall(PetscSectionSetFieldConstraintIndices(*subs, p, f, oldIndices));
18469574d0f0SMatthew G. Knepley           for (fc = 0; fc < cfdof; ++fc) indices[numConst + fc] = oldIndices[fc] + fOff;
1847ea844a1aSMatthew Knepley           numConst += cfdof;
1848ea844a1aSMatthew Knepley           fOff += fdof;
1849ea844a1aSMatthew Knepley         }
18509566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetConstraintIndices(*subs, p, indices));
1851ea844a1aSMatthew Knepley       }
1852ea844a1aSMatthew Knepley     }
18539566063dSJacob Faibussowitsch     PetscCall(PetscFree(indices));
1854ea844a1aSMatthew Knepley   }
1855ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1856ea844a1aSMatthew Knepley }
1857ea844a1aSMatthew Knepley 
1858ea844a1aSMatthew Knepley /*@
1859ea844a1aSMatthew Knepley   PetscSectionCreateSupersection - Create a new, larger section composed of the input sections
1860ea844a1aSMatthew Knepley 
186140750d41SVaclav Hapla   Collective
1862ea844a1aSMatthew Knepley 
1863ea844a1aSMatthew Knepley   Input Parameters:
1864ea844a1aSMatthew Knepley + s     - the input sections
1865ea844a1aSMatthew Knepley - len   - the number of input sections
1866ea844a1aSMatthew Knepley 
1867ea844a1aSMatthew Knepley   Output Parameter:
1868ea844a1aSMatthew Knepley . supers - the supersection
1869ea844a1aSMatthew Knepley 
1870ea844a1aSMatthew Knepley   Level: advanced
1871ea844a1aSMatthew Knepley 
1872cab54364SBarry Smith   Note:
1873cab54364SBarry Smith   The section offsets now refer to a new, larger vector.
1874cab54364SBarry Smith 
1875cab54364SBarry Smith   Developer Note:
1876cab54364SBarry Smith   Needs to explain how the sections are composed
1877cab54364SBarry Smith 
1878cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreateSubsection()`, `PetscSectionCreate()`
1879ea844a1aSMatthew Knepley @*/
1880d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionCreateSupersection(PetscSection s[], PetscInt len, PetscSection *supers)
1881d71ae5a4SJacob Faibussowitsch {
1882ea844a1aSMatthew Knepley   PetscInt Nf = 0, f, pStart = PETSC_MAX_INT, pEnd = 0, p, maxCdof = 0, i;
1883ea844a1aSMatthew Knepley 
1884ea844a1aSMatthew Knepley   PetscFunctionBegin;
1885ea844a1aSMatthew Knepley   if (!len) PetscFunctionReturn(0);
1886ea844a1aSMatthew Knepley   for (i = 0; i < len; ++i) {
1887ea844a1aSMatthew Knepley     PetscInt nf, pStarti, pEndi;
1888ea844a1aSMatthew Knepley 
18899566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetNumFields(s[i], &nf));
18909566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(s[i], &pStarti, &pEndi));
1891ea844a1aSMatthew Knepley     pStart = PetscMin(pStart, pStarti);
1892ea844a1aSMatthew Knepley     pEnd   = PetscMax(pEnd, pEndi);
1893ea844a1aSMatthew Knepley     Nf += nf;
1894ea844a1aSMatthew Knepley   }
18959566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)s[0]), supers));
18969566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetNumFields(*supers, Nf));
1897ea844a1aSMatthew Knepley   for (i = 0, f = 0; i < len; ++i) {
1898b778fa18SValeria Barra     PetscInt nf, fi, ci;
1899ea844a1aSMatthew Knepley 
19009566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetNumFields(s[i], &nf));
1901ea844a1aSMatthew Knepley     for (fi = 0; fi < nf; ++fi, ++f) {
1902ea844a1aSMatthew Knepley       const char *name    = NULL;
1903ea844a1aSMatthew Knepley       PetscInt    numComp = 0;
1904ea844a1aSMatthew Knepley 
19059566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldName(s[i], fi, &name));
19069566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetFieldName(*supers, f, name));
19079566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldComponents(s[i], fi, &numComp));
19089566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetFieldComponents(*supers, f, numComp));
1909b778fa18SValeria Barra       for (ci = 0; ci < s[i]->numFieldComponents[fi]; ++ci) {
19109566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetComponentName(s[i], fi, ci, &name));
19119566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetComponentName(*supers, f, ci, name));
1912b778fa18SValeria Barra       }
1913ea844a1aSMatthew Knepley     }
1914ea844a1aSMatthew Knepley   }
19159566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(*supers, pStart, pEnd));
1916ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
1917ea844a1aSMatthew Knepley     PetscInt dof = 0, cdof = 0;
1918ea844a1aSMatthew Knepley 
1919ea844a1aSMatthew Knepley     for (i = 0, f = 0; i < len; ++i) {
1920ea844a1aSMatthew Knepley       PetscInt nf, fi, pStarti, pEndi;
1921ea844a1aSMatthew Knepley       PetscInt fdof = 0, cfdof = 0;
1922ea844a1aSMatthew Knepley 
19239566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetNumFields(s[i], &nf));
19249566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetChart(s[i], &pStarti, &pEndi));
1925ea844a1aSMatthew Knepley       if ((p < pStarti) || (p >= pEndi)) continue;
1926ea844a1aSMatthew Knepley       for (fi = 0; fi < nf; ++fi, ++f) {
19279566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetFieldDof(s[i], p, fi, &fdof));
19289566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddFieldDof(*supers, p, f, fdof));
19299566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetFieldConstraintDof(s[i], p, fi, &cfdof));
19309566063dSJacob Faibussowitsch         if (cfdof) PetscCall(PetscSectionAddFieldConstraintDof(*supers, p, f, cfdof));
1931ea844a1aSMatthew Knepley         dof += fdof;
1932ea844a1aSMatthew Knepley         cdof += cfdof;
1933ea844a1aSMatthew Knepley       }
1934ea844a1aSMatthew Knepley     }
19359566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(*supers, p, dof));
19369566063dSJacob Faibussowitsch     if (cdof) PetscCall(PetscSectionSetConstraintDof(*supers, p, cdof));
1937ea844a1aSMatthew Knepley     maxCdof = PetscMax(cdof, maxCdof);
1938ea844a1aSMatthew Knepley   }
19399566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(*supers));
1940ea844a1aSMatthew Knepley   if (maxCdof) {
1941ea844a1aSMatthew Knepley     PetscInt *indices;
1942ea844a1aSMatthew Knepley 
19439566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(maxCdof, &indices));
1944ea844a1aSMatthew Knepley     for (p = pStart; p < pEnd; ++p) {
1945ea844a1aSMatthew Knepley       PetscInt cdof;
1946ea844a1aSMatthew Knepley 
19479566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintDof(*supers, p, &cdof));
1948ea844a1aSMatthew Knepley       if (cdof) {
1949ea844a1aSMatthew Knepley         PetscInt dof, numConst = 0, fOff = 0;
1950ea844a1aSMatthew Knepley 
1951ea844a1aSMatthew Knepley         for (i = 0, f = 0; i < len; ++i) {
1952ea844a1aSMatthew Knepley           const PetscInt *oldIndices = NULL;
1953ea844a1aSMatthew Knepley           PetscInt        nf, fi, pStarti, pEndi, fdof, cfdof, fc;
1954ea844a1aSMatthew Knepley 
19559566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetNumFields(s[i], &nf));
19569566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetChart(s[i], &pStarti, &pEndi));
1957ea844a1aSMatthew Knepley           if ((p < pStarti) || (p >= pEndi)) continue;
1958ea844a1aSMatthew Knepley           for (fi = 0; fi < nf; ++fi, ++f) {
19599566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(s[i], p, fi, &fdof));
19609566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldConstraintDof(s[i], p, fi, &cfdof));
19619566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldConstraintIndices(s[i], p, fi, &oldIndices));
1962fcd2503dSMatthew G. Knepley             for (fc = 0; fc < cfdof; ++fc) indices[numConst + fc] = oldIndices[fc];
19639566063dSJacob Faibussowitsch             PetscCall(PetscSectionSetFieldConstraintIndices(*supers, p, f, &indices[numConst]));
1964fcd2503dSMatthew G. Knepley             for (fc = 0; fc < cfdof; ++fc) indices[numConst + fc] += fOff;
1965ea844a1aSMatthew Knepley             numConst += cfdof;
1966ea844a1aSMatthew Knepley           }
19679566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetDof(s[i], p, &dof));
1968ea844a1aSMatthew Knepley           fOff += dof;
1969ea844a1aSMatthew Knepley         }
19709566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetConstraintIndices(*supers, p, indices));
1971ea844a1aSMatthew Knepley       }
1972ea844a1aSMatthew Knepley     }
19739566063dSJacob Faibussowitsch     PetscCall(PetscFree(indices));
1974ea844a1aSMatthew Knepley   }
1975ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
1976ea844a1aSMatthew Knepley }
1977ea844a1aSMatthew Knepley 
1978d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionCreateSubplexSection_Internal(PetscSection s, IS subpointMap, PetscBool renumberPoints, PetscSection *subs)
1979d71ae5a4SJacob Faibussowitsch {
1980ea844a1aSMatthew Knepley   const PetscInt *points = NULL, *indices = NULL;
198141f23ed0SMatthew G. Knepley   PetscInt        numFields, f, c, numSubpoints = 0, pStart, pEnd, p, spStart, spEnd, subp;
1982ea844a1aSMatthew Knepley 
1983ea844a1aSMatthew Knepley   PetscFunctionBegin;
1984ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
1985ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(subpointMap, IS_CLASSID, 2);
198641f23ed0SMatthew G. Knepley   PetscValidPointer(subs, 4);
19879566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(s, &numFields));
19889566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)s), subs));
19899566063dSJacob Faibussowitsch   if (numFields) PetscCall(PetscSectionSetNumFields(*subs, numFields));
1990ea844a1aSMatthew Knepley   for (f = 0; f < numFields; ++f) {
1991ea844a1aSMatthew Knepley     const char *name    = NULL;
1992ea844a1aSMatthew Knepley     PetscInt    numComp = 0;
1993ea844a1aSMatthew Knepley 
19949566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldName(s, f, &name));
19959566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldName(*subs, f, name));
19969566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldComponents(s, f, &numComp));
19979566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldComponents(*subs, f, numComp));
1998b778fa18SValeria Barra     for (c = 0; c < s->numFieldComponents[f]; ++c) {
19999566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetComponentName(s, f, c, &name));
20009566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetComponentName(*subs, f, c, name));
2001b778fa18SValeria Barra     }
2002ea844a1aSMatthew Knepley   }
2003ea844a1aSMatthew Knepley   /* For right now, we do not try to squeeze the subchart */
2004ea844a1aSMatthew Knepley   if (subpointMap) {
20059566063dSJacob Faibussowitsch     PetscCall(ISGetSize(subpointMap, &numSubpoints));
20069566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(subpointMap, &points));
2007ea844a1aSMatthew Knepley   }
20089566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
200941f23ed0SMatthew G. Knepley   if (renumberPoints) {
201041f23ed0SMatthew G. Knepley     spStart = 0;
201141f23ed0SMatthew G. Knepley     spEnd   = numSubpoints;
201241f23ed0SMatthew G. Knepley   } else {
201341f23ed0SMatthew G. Knepley     PetscCall(ISGetMinMax(subpointMap, &spStart, &spEnd));
201441f23ed0SMatthew G. Knepley     ++spEnd;
201541f23ed0SMatthew G. Knepley   }
201641f23ed0SMatthew G. Knepley   PetscCall(PetscSectionSetChart(*subs, spStart, spEnd));
2017ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
2018ea844a1aSMatthew Knepley     PetscInt dof, cdof, fdof = 0, cfdof = 0;
2019ea844a1aSMatthew Knepley 
20209566063dSJacob Faibussowitsch     PetscCall(PetscFindInt(p, numSubpoints, points, &subp));
2021ea844a1aSMatthew Knepley     if (subp < 0) continue;
202241f23ed0SMatthew G. Knepley     if (!renumberPoints) subp = p;
2023ea844a1aSMatthew Knepley     for (f = 0; f < numFields; ++f) {
20249566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof));
20259566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetFieldDof(*subs, subp, f, fdof));
20269566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cfdof));
20279566063dSJacob Faibussowitsch       if (cfdof) PetscCall(PetscSectionSetFieldConstraintDof(*subs, subp, f, cfdof));
2028ea844a1aSMatthew Knepley     }
20299566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(s, p, &dof));
20309566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(*subs, subp, dof));
20319566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
20329566063dSJacob Faibussowitsch     if (cdof) PetscCall(PetscSectionSetConstraintDof(*subs, subp, cdof));
2033ea844a1aSMatthew Knepley   }
20349566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(*subs));
2035ea844a1aSMatthew Knepley   /* Change offsets to original offsets */
2036ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
2037ea844a1aSMatthew Knepley     PetscInt off, foff = 0;
2038ea844a1aSMatthew Knepley 
20399566063dSJacob Faibussowitsch     PetscCall(PetscFindInt(p, numSubpoints, points, &subp));
2040ea844a1aSMatthew Knepley     if (subp < 0) continue;
204141f23ed0SMatthew G. Knepley     if (!renumberPoints) subp = p;
2042ea844a1aSMatthew Knepley     for (f = 0; f < numFields; ++f) {
20439566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldOffset(s, p, f, &foff));
20449566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetFieldOffset(*subs, subp, f, foff));
2045ea844a1aSMatthew Knepley     }
20469566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(s, p, &off));
20479566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetOffset(*subs, subp, off));
2048ea844a1aSMatthew Knepley   }
2049ea844a1aSMatthew Knepley   /* Copy constraint indices */
205041f23ed0SMatthew G. Knepley   for (subp = spStart; subp < spEnd; ++subp) {
2051ea844a1aSMatthew Knepley     PetscInt cdof;
2052ea844a1aSMatthew Knepley 
20539566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(*subs, subp, &cdof));
2054ea844a1aSMatthew Knepley     if (cdof) {
2055ea844a1aSMatthew Knepley       for (f = 0; f < numFields; ++f) {
205641f23ed0SMatthew G. Knepley         PetscCall(PetscSectionGetFieldConstraintIndices(s, points[subp - spStart], f, &indices));
20579566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetFieldConstraintIndices(*subs, subp, f, indices));
2058ea844a1aSMatthew Knepley       }
205941f23ed0SMatthew G. Knepley       PetscCall(PetscSectionGetConstraintIndices(s, points[subp - spStart], &indices));
20609566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetConstraintIndices(*subs, subp, indices));
2061ea844a1aSMatthew Knepley     }
2062ea844a1aSMatthew Knepley   }
20639566063dSJacob Faibussowitsch   if (subpointMap) PetscCall(ISRestoreIndices(subpointMap, &points));
2064ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2065ea844a1aSMatthew Knepley }
2066ea844a1aSMatthew Knepley 
206741f23ed0SMatthew G. Knepley /*@
206841f23ed0SMatthew G. Knepley   PetscSectionCreateSubmeshSection - Create a new, smaller section with support on the submesh
206941f23ed0SMatthew G. Knepley 
207041f23ed0SMatthew G. Knepley   Collective on s
207141f23ed0SMatthew G. Knepley 
207241f23ed0SMatthew G. Knepley   Input Parameters:
2073cab54364SBarry Smith + s           - the `PetscSection`
207441f23ed0SMatthew G. Knepley - subpointMap - a sorted list of points in the original mesh which are in the submesh
207541f23ed0SMatthew G. Knepley 
207641f23ed0SMatthew G. Knepley   Output Parameter:
207741f23ed0SMatthew G. Knepley . subs - the subsection
207841f23ed0SMatthew G. Knepley 
2079cab54364SBarry Smith   Level: advanced
2080cab54364SBarry Smith 
208141f23ed0SMatthew G. Knepley   Note:
208241f23ed0SMatthew G. Knepley   The points are renumbered from 0, and the section offsets now refer to a new, smaller vector.
208341f23ed0SMatthew G. Knepley 
2084cab54364SBarry Smith   Developer Note:
2085cab54364SBarry Smith   The use of the term Submesh is confusing and unneeded
208641f23ed0SMatthew G. Knepley 
2087cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreateSubdomainSection()`, `PetscSectionCreateSubsection()`, `DMPlexGetSubpointMap()`, `PetscSectionCreate()`
208841f23ed0SMatthew G. Knepley @*/
2089d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionCreateSubmeshSection(PetscSection s, IS subpointMap, PetscSection *subs)
2090d71ae5a4SJacob Faibussowitsch {
209141f23ed0SMatthew G. Knepley   PetscFunctionBegin;
209241f23ed0SMatthew G. Knepley   PetscCall(PetscSectionCreateSubplexSection_Internal(s, subpointMap, PETSC_TRUE, subs));
209341f23ed0SMatthew G. Knepley   PetscFunctionReturn(0);
209441f23ed0SMatthew G. Knepley }
209541f23ed0SMatthew G. Knepley 
209641f23ed0SMatthew G. Knepley /*@
209741f23ed0SMatthew G. Knepley   PetscSectionCreateSubdomainSection - Create a new, smaller section with support on a subdomain of the mesh
209841f23ed0SMatthew G. Knepley 
209941f23ed0SMatthew G. Knepley   Collective on s
210041f23ed0SMatthew G. Knepley 
210141f23ed0SMatthew G. Knepley   Input Parameters:
2102cab54364SBarry Smith + s           - the `PetscSection`
210341f23ed0SMatthew G. Knepley - subpointMap - a sorted list of points in the original mesh which are in the subdomain
210441f23ed0SMatthew G. Knepley 
210541f23ed0SMatthew G. Knepley   Output Parameter:
210641f23ed0SMatthew G. Knepley . subs - the subsection
210741f23ed0SMatthew G. Knepley 
2108cab54364SBarry Smith   Level: advanced
2109cab54364SBarry Smith 
211041f23ed0SMatthew G. Knepley   Note:
211141f23ed0SMatthew G. Knepley   The point numbers remain the same, but the section offsets now refer to a new, smaller vector.
211241f23ed0SMatthew G. Knepley 
2113cab54364SBarry Smith   Developer Notes:
2114cab54364SBarry Smith   It is unclear what the difference with `PetscSectionCreateSubmeshSection()` is.
211541f23ed0SMatthew G. Knepley 
2116cab54364SBarry Smith   The use of the term Subdomain is unneeded and confusing
2117cab54364SBarry Smith 
2118cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreateSubmeshSection()`, `PetscSectionCreateSubsection()`, `DMPlexGetSubpointMap()`, `PetscSectionCreate()`
211941f23ed0SMatthew G. Knepley @*/
2120d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionCreateSubdomainSection(PetscSection s, IS subpointMap, PetscSection *subs)
2121d71ae5a4SJacob Faibussowitsch {
212241f23ed0SMatthew G. Knepley   PetscFunctionBegin;
212341f23ed0SMatthew G. Knepley   PetscCall(PetscSectionCreateSubplexSection_Internal(s, subpointMap, PETSC_FALSE, subs));
212441f23ed0SMatthew G. Knepley   PetscFunctionReturn(0);
212541f23ed0SMatthew G. Knepley }
212641f23ed0SMatthew G. Knepley 
2127d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSectionView_ASCII(PetscSection s, PetscViewer viewer)
2128d71ae5a4SJacob Faibussowitsch {
2129ea844a1aSMatthew Knepley   PetscInt    p;
2130ea844a1aSMatthew Knepley   PetscMPIInt rank;
2131ea844a1aSMatthew Knepley 
2132ea844a1aSMatthew Knepley   PetscFunctionBegin;
21339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
21349566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
21359566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Process %d:\n", rank));
2136ea844a1aSMatthew Knepley   for (p = 0; p < s->pEnd - s->pStart; ++p) {
2137ea844a1aSMatthew Knepley     if ((s->bc) && (s->bc->atlasDof[p] > 0)) {
2138ea844a1aSMatthew Knepley       PetscInt b;
2139ea844a1aSMatthew Knepley 
21409566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  (%4" PetscInt_FMT ") dim %2" PetscInt_FMT " offset %3" PetscInt_FMT " constrained", p + s->pStart, s->atlasDof[p], s->atlasOff[p]));
2141083401c6SMatthew G. Knepley       if (s->bcIndices) {
214248a46eb9SPierre Jolivet         for (b = 0; b < s->bc->atlasDof[p]; ++b) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %" PetscInt_FMT, s->bcIndices[s->bc->atlasOff[p] + b]));
2143083401c6SMatthew G. Knepley       }
21449566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
2145ea844a1aSMatthew Knepley     } else {
21469566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  (%4" PetscInt_FMT ") dim %2" PetscInt_FMT " offset %3" PetscInt_FMT "\n", p + s->pStart, s->atlasDof[p], s->atlasOff[p]));
2147ea844a1aSMatthew Knepley     }
2148ea844a1aSMatthew Knepley   }
21499566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
21509566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2151ea844a1aSMatthew Knepley   if (s->sym) {
21529566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
21539566063dSJacob Faibussowitsch     PetscCall(PetscSectionSymView(s->sym, viewer));
21549566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
2155ea844a1aSMatthew Knepley   }
2156ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2157ea844a1aSMatthew Knepley }
2158ea844a1aSMatthew Knepley 
2159ea844a1aSMatthew Knepley /*@C
2160cab54364SBarry Smith    PetscSectionViewFromOptions - View the `PetscSection` based on values in the options database
2161fe2efc57SMark 
216240750d41SVaclav Hapla    Collective
2163fe2efc57SMark 
2164fe2efc57SMark    Input Parameters:
2165fe2efc57SMark +  A - the PetscSection object to view
2166cab54364SBarry Smith .  obj - Optional object that provides the options prefix used for the options
2167736c3998SJose E. Roman -  name - command line option
2168fe2efc57SMark 
2169fe2efc57SMark    Level: intermediate
2170cab54364SBarry Smith 
2171cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionView`, `PetscObjectViewFromOptions()`, `PetscSectionCreate()`, `PetscSectionView()`
2172fe2efc57SMark @*/
2173d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionViewFromOptions(PetscSection A, PetscObject obj, const char name[])
2174d71ae5a4SJacob Faibussowitsch {
2175fe2efc57SMark   PetscFunctionBegin;
2176fe2efc57SMark   PetscValidHeaderSpecific(A, PETSC_SECTION_CLASSID, 1);
21779566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
2178fe2efc57SMark   PetscFunctionReturn(0);
2179fe2efc57SMark }
2180fe2efc57SMark 
2181fe2efc57SMark /*@C
2182cab54364SBarry Smith   PetscSectionView - Views a `PetscSection`
2183ea844a1aSMatthew Knepley 
218440750d41SVaclav Hapla   Collective
2185ea844a1aSMatthew Knepley 
2186ea844a1aSMatthew Knepley   Input Parameters:
2187cab54364SBarry Smith + s - the `PetscSection` object to view
2188ea844a1aSMatthew Knepley - v - the viewer
2189ea844a1aSMatthew Knepley 
2190cab54364SBarry Smith   Level: beginner
2191cab54364SBarry Smith 
21921ddd528eSksagiyam   Note:
2193cab54364SBarry Smith   `PetscSectionView()`, when viewer is of type `PETSCVIEWERHDF5`, only saves
21941ddd528eSksagiyam   distribution independent data, such as dofs, offsets, constraint dofs,
21951ddd528eSksagiyam   and constraint indices. Points that have negative dofs, for instance,
21961ddd528eSksagiyam   are not saved as they represent points owned by other processes.
21971ddd528eSksagiyam   Point numbering and rank assignment is currently not stored.
2198cab54364SBarry Smith   The saved section can be loaded with `PetscSectionLoad()`.
21991ddd528eSksagiyam 
2200cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreate()`, `PetscSectionDestroy()`, `PetscSectionLoad()`, `PetscViewer`
2201ea844a1aSMatthew Knepley @*/
2202d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionView(PetscSection s, PetscViewer viewer)
2203d71ae5a4SJacob Faibussowitsch {
22041ddd528eSksagiyam   PetscBool isascii, ishdf5;
2205ea844a1aSMatthew Knepley   PetscInt  f;
2206ea844a1aSMatthew Knepley 
2207ea844a1aSMatthew Knepley   PetscFunctionBegin;
2208ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
22099566063dSJacob Faibussowitsch   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)s), &viewer));
2210ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
22119566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
22129566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2213ea844a1aSMatthew Knepley   if (isascii) {
22149566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)s, viewer));
2215ea844a1aSMatthew Knepley     if (s->numFields) {
22169566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " fields\n", s->numFields));
2217ea844a1aSMatthew Knepley       for (f = 0; f < s->numFields; ++f) {
22189566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  field %" PetscInt_FMT " with %" PetscInt_FMT " components\n", f, s->numFieldComponents[f]));
22199566063dSJacob Faibussowitsch         PetscCall(PetscSectionView_ASCII(s->field[f], viewer));
2220ea844a1aSMatthew Knepley       }
2221ea844a1aSMatthew Knepley     } else {
22229566063dSJacob Faibussowitsch       PetscCall(PetscSectionView_ASCII(s, viewer));
2223ea844a1aSMatthew Knepley     }
22241ddd528eSksagiyam   } else if (ishdf5) {
22251ddd528eSksagiyam #if PetscDefined(HAVE_HDF5)
22269566063dSJacob Faibussowitsch     PetscCall(PetscSectionView_HDF5_Internal(s, viewer));
22271ddd528eSksagiyam #else
22281ddd528eSksagiyam     SETERRQ(PetscObjectComm((PetscObject)s), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
22291ddd528eSksagiyam #endif
2230ea844a1aSMatthew Knepley   }
2231ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2232ea844a1aSMatthew Knepley }
2233ea844a1aSMatthew Knepley 
2234fde5e3acSksagiyam /*@C
2235cab54364SBarry Smith   PetscSectionLoad - Loads a `PetscSection`
2236fde5e3acSksagiyam 
223740750d41SVaclav Hapla   Collective
2238fde5e3acSksagiyam 
2239fde5e3acSksagiyam   Input Parameters:
2240cab54364SBarry Smith + s - the `PetscSection` object to load
2241fde5e3acSksagiyam - v - the viewer
2242fde5e3acSksagiyam 
2243cab54364SBarry Smith   Level: beginner
2244cab54364SBarry Smith 
2245fde5e3acSksagiyam   Note:
2246cab54364SBarry Smith   `PetscSectionLoad()`, when viewer is of type `PETSCVIEWERHDF5`, loads
2247cab54364SBarry Smith   a section saved with `PetscSectionView()`. The number of processes
2248fde5e3acSksagiyam   used here (N) does not need to be the same as that used when saving.
2249fde5e3acSksagiyam   After calling this function, the chart of s on rank i will be set
2250fde5e3acSksagiyam   to [0, E_i), where \sum_{i=0}^{N-1}E_i equals to the total number of
2251fde5e3acSksagiyam   saved section points.
2252fde5e3acSksagiyam 
2253cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreate()`, `PetscSectionDestroy()`, `PetscSectionView()`
2254fde5e3acSksagiyam @*/
2255d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionLoad(PetscSection s, PetscViewer viewer)
2256d71ae5a4SJacob Faibussowitsch {
2257fde5e3acSksagiyam   PetscBool ishdf5;
2258fde5e3acSksagiyam 
2259fde5e3acSksagiyam   PetscFunctionBegin;
2260fde5e3acSksagiyam   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
2261fde5e3acSksagiyam   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
22629566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2263fde5e3acSksagiyam   if (ishdf5) {
2264fde5e3acSksagiyam #if PetscDefined(HAVE_HDF5)
22659566063dSJacob Faibussowitsch     PetscCall(PetscSectionLoad_HDF5_Internal(s, viewer));
2266fde5e3acSksagiyam     PetscFunctionReturn(0);
2267fde5e3acSksagiyam #else
2268fde5e3acSksagiyam     SETERRQ(PetscObjectComm((PetscObject)s), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2269fde5e3acSksagiyam #endif
227098921bdaSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)s), PETSC_ERR_SUP, "Viewer type %s not yet supported for PetscSection loading", ((PetscObject)viewer)->type_name);
2271fde5e3acSksagiyam }
2272fde5e3acSksagiyam 
2273d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSectionResetClosurePermutation(PetscSection section)
2274d71ae5a4SJacob Faibussowitsch {
2275c459fbc1SJed Brown   PetscSectionClosurePermVal clVal;
2276c459fbc1SJed Brown 
2277c459fbc1SJed Brown   PetscFunctionBegin;
2278c459fbc1SJed Brown   if (!section->clHash) PetscFunctionReturn(0);
2279c459fbc1SJed Brown   kh_foreach_value(section->clHash, clVal, {
22809566063dSJacob Faibussowitsch     PetscCall(PetscFree(clVal.perm));
22819566063dSJacob Faibussowitsch     PetscCall(PetscFree(clVal.invPerm));
2282c459fbc1SJed Brown   });
2283c459fbc1SJed Brown   kh_destroy(ClPerm, section->clHash);
2284c459fbc1SJed Brown   section->clHash = NULL;
2285c459fbc1SJed Brown   PetscFunctionReturn(0);
2286c459fbc1SJed Brown }
2287c459fbc1SJed Brown 
2288ea844a1aSMatthew Knepley /*@
2289ea844a1aSMatthew Knepley   PetscSectionReset - Frees all section data.
2290ea844a1aSMatthew Knepley 
229140750d41SVaclav Hapla   Not Collective
2292ea844a1aSMatthew Knepley 
2293ea844a1aSMatthew Knepley   Input Parameters:
2294cab54364SBarry Smith . s - the `PetscSection`
2295ea844a1aSMatthew Knepley 
2296ea844a1aSMatthew Knepley   Level: beginner
2297ea844a1aSMatthew Knepley 
2298cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreate()`
2299ea844a1aSMatthew Knepley @*/
2300d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionReset(PetscSection s)
2301d71ae5a4SJacob Faibussowitsch {
2302b778fa18SValeria Barra   PetscInt f, c;
2303ea844a1aSMatthew Knepley 
2304ea844a1aSMatthew Knepley   PetscFunctionBegin;
2305ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
2306ea844a1aSMatthew Knepley   for (f = 0; f < s->numFields; ++f) {
23079566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&s->field[f]));
23089566063dSJacob Faibussowitsch     PetscCall(PetscFree(s->fieldNames[f]));
230948a46eb9SPierre Jolivet     for (c = 0; c < s->numFieldComponents[f]; ++c) PetscCall(PetscFree(s->compNames[f][c]));
23109566063dSJacob Faibussowitsch     PetscCall(PetscFree(s->compNames[f]));
2311ea844a1aSMatthew Knepley   }
23129566063dSJacob Faibussowitsch   PetscCall(PetscFree(s->numFieldComponents));
23139566063dSJacob Faibussowitsch   PetscCall(PetscFree(s->fieldNames));
23149566063dSJacob Faibussowitsch   PetscCall(PetscFree(s->compNames));
23159566063dSJacob Faibussowitsch   PetscCall(PetscFree(s->field));
23169566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&s->bc));
23179566063dSJacob Faibussowitsch   PetscCall(PetscFree(s->bcIndices));
23189566063dSJacob Faibussowitsch   PetscCall(PetscFree2(s->atlasDof, s->atlasOff));
23199566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&s->clSection));
23209566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&s->clPoints));
23219566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&s->perm));
23229566063dSJacob Faibussowitsch   PetscCall(PetscSectionResetClosurePermutation(s));
23239566063dSJacob Faibussowitsch   PetscCall(PetscSectionSymDestroy(&s->sym));
23249566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&s->clSection));
23259566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&s->clPoints));
232669c11d05SVaclav Hapla   PetscCall(PetscSectionInvalidateMaxDof_Internal(s));
2327ea844a1aSMatthew Knepley   s->pStart    = -1;
2328ea844a1aSMatthew Knepley   s->pEnd      = -1;
2329ea844a1aSMatthew Knepley   s->maxDof    = 0;
2330ea844a1aSMatthew Knepley   s->setup     = PETSC_FALSE;
2331ea844a1aSMatthew Knepley   s->numFields = 0;
2332ea844a1aSMatthew Knepley   s->clObj     = NULL;
2333ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2334ea844a1aSMatthew Knepley }
2335ea844a1aSMatthew Knepley 
2336ea844a1aSMatthew Knepley /*@
2337ea844a1aSMatthew Knepley   PetscSectionDestroy - Frees a section object and frees its range if that exists.
2338ea844a1aSMatthew Knepley 
233940750d41SVaclav Hapla   Not Collective
2340ea844a1aSMatthew Knepley 
2341ea844a1aSMatthew Knepley   Input Parameters:
2342cab54364SBarry Smith . s - the `PetscSection`
2343ea844a1aSMatthew Knepley 
2344ea844a1aSMatthew Knepley   Level: beginner
2345ea844a1aSMatthew Knepley 
2346cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionCreate()`, `PetscSectionReset()`
2347ea844a1aSMatthew Knepley @*/
2348d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionDestroy(PetscSection *s)
2349d71ae5a4SJacob Faibussowitsch {
2350ea844a1aSMatthew Knepley   PetscFunctionBegin;
2351ea844a1aSMatthew Knepley   if (!*s) PetscFunctionReturn(0);
235262f17a49SStefano Zampini   PetscValidHeaderSpecific(*s, PETSC_SECTION_CLASSID, 1);
2353ea844a1aSMatthew Knepley   if (--((PetscObject)(*s))->refct > 0) {
2354ea844a1aSMatthew Knepley     *s = NULL;
2355ea844a1aSMatthew Knepley     PetscFunctionReturn(0);
2356ea844a1aSMatthew Knepley   }
23579566063dSJacob Faibussowitsch   PetscCall(PetscSectionReset(*s));
23589566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(s));
2359ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2360ea844a1aSMatthew Knepley }
2361ea844a1aSMatthew Knepley 
2362d71ae5a4SJacob Faibussowitsch PetscErrorCode VecIntGetValuesSection(PetscInt *baseArray, PetscSection s, PetscInt point, const PetscInt **values)
2363d71ae5a4SJacob Faibussowitsch {
2364ea844a1aSMatthew Knepley   const PetscInt p = point - s->pStart;
2365ea844a1aSMatthew Knepley 
2366ea844a1aSMatthew Knepley   PetscFunctionBegin;
2367ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 2);
2368ea844a1aSMatthew Knepley   *values = &baseArray[s->atlasOff[p]];
2369ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2370ea844a1aSMatthew Knepley }
2371ea844a1aSMatthew Knepley 
2372d71ae5a4SJacob Faibussowitsch PetscErrorCode VecIntSetValuesSection(PetscInt *baseArray, PetscSection s, PetscInt point, const PetscInt values[], InsertMode mode)
2373d71ae5a4SJacob Faibussowitsch {
2374ea844a1aSMatthew Knepley   PetscInt      *array;
2375ea844a1aSMatthew Knepley   const PetscInt p           = point - s->pStart;
2376ea844a1aSMatthew Knepley   const PetscInt orientation = 0; /* Needs to be included for use in closure operations */
2377ea844a1aSMatthew Knepley   PetscInt       cDim        = 0;
2378ea844a1aSMatthew Knepley 
2379ea844a1aSMatthew Knepley   PetscFunctionBegin;
2380ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 2);
23819566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstraintDof(s, p, &cDim));
2382ea844a1aSMatthew Knepley   array = &baseArray[s->atlasOff[p]];
2383ea844a1aSMatthew Knepley   if (!cDim) {
2384ea844a1aSMatthew Knepley     if (orientation >= 0) {
2385ea844a1aSMatthew Knepley       const PetscInt dim = s->atlasDof[p];
2386ea844a1aSMatthew Knepley       PetscInt       i;
2387ea844a1aSMatthew Knepley 
2388ea844a1aSMatthew Knepley       if (mode == INSERT_VALUES) {
2389ea844a1aSMatthew Knepley         for (i = 0; i < dim; ++i) array[i] = values[i];
2390ea844a1aSMatthew Knepley       } else {
2391ea844a1aSMatthew Knepley         for (i = 0; i < dim; ++i) array[i] += values[i];
2392ea844a1aSMatthew Knepley       }
2393ea844a1aSMatthew Knepley     } else {
2394ea844a1aSMatthew Knepley       PetscInt offset = 0;
2395ea844a1aSMatthew Knepley       PetscInt j      = -1, field, i;
2396ea844a1aSMatthew Knepley 
2397ea844a1aSMatthew Knepley       for (field = 0; field < s->numFields; ++field) {
2398ea844a1aSMatthew Knepley         const PetscInt dim = s->field[field]->atlasDof[p];
2399ea844a1aSMatthew Knepley 
2400ea844a1aSMatthew Knepley         for (i = dim - 1; i >= 0; --i) array[++j] = values[i + offset];
2401ea844a1aSMatthew Knepley         offset += dim;
2402ea844a1aSMatthew Knepley       }
2403ea844a1aSMatthew Knepley     }
2404ea844a1aSMatthew Knepley   } else {
2405ea844a1aSMatthew Knepley     if (orientation >= 0) {
2406ea844a1aSMatthew Knepley       const PetscInt  dim  = s->atlasDof[p];
2407ea844a1aSMatthew Knepley       PetscInt        cInd = 0, i;
2408ea844a1aSMatthew Knepley       const PetscInt *cDof;
2409ea844a1aSMatthew Knepley 
24109566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintIndices(s, point, &cDof));
2411ea844a1aSMatthew Knepley       if (mode == INSERT_VALUES) {
2412ea844a1aSMatthew Knepley         for (i = 0; i < dim; ++i) {
24139371c9d4SSatish Balay           if ((cInd < cDim) && (i == cDof[cInd])) {
24149371c9d4SSatish Balay             ++cInd;
24159371c9d4SSatish Balay             continue;
24169371c9d4SSatish Balay           }
2417ea844a1aSMatthew Knepley           array[i] = values[i];
2418ea844a1aSMatthew Knepley         }
2419ea844a1aSMatthew Knepley       } else {
2420ea844a1aSMatthew Knepley         for (i = 0; i < dim; ++i) {
24219371c9d4SSatish Balay           if ((cInd < cDim) && (i == cDof[cInd])) {
24229371c9d4SSatish Balay             ++cInd;
24239371c9d4SSatish Balay             continue;
24249371c9d4SSatish Balay           }
2425ea844a1aSMatthew Knepley           array[i] += values[i];
2426ea844a1aSMatthew Knepley         }
2427ea844a1aSMatthew Knepley       }
2428ea844a1aSMatthew Knepley     } else {
2429ea844a1aSMatthew Knepley       const PetscInt *cDof;
2430ea844a1aSMatthew Knepley       PetscInt        offset  = 0;
2431ea844a1aSMatthew Knepley       PetscInt        cOffset = 0;
2432ea844a1aSMatthew Knepley       PetscInt        j       = 0, field;
2433ea844a1aSMatthew Knepley 
24349566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintIndices(s, point, &cDof));
2435ea844a1aSMatthew Knepley       for (field = 0; field < s->numFields; ++field) {
2436ea844a1aSMatthew Knepley         const PetscInt dim  = s->field[field]->atlasDof[p];     /* PetscSectionGetFieldDof() */
2437ea844a1aSMatthew Knepley         const PetscInt tDim = s->field[field]->bc->atlasDof[p]; /* PetscSectionGetFieldConstraintDof() */
2438ea844a1aSMatthew Knepley         const PetscInt sDim = dim - tDim;
2439ea844a1aSMatthew Knepley         PetscInt       cInd = 0, i, k;
2440ea844a1aSMatthew Knepley 
2441ea844a1aSMatthew Knepley         for (i = 0, k = dim + offset - 1; i < dim; ++i, ++j, --k) {
24429371c9d4SSatish Balay           if ((cInd < sDim) && (j == cDof[cInd + cOffset])) {
24439371c9d4SSatish Balay             ++cInd;
24449371c9d4SSatish Balay             continue;
24459371c9d4SSatish Balay           }
2446ea844a1aSMatthew Knepley           array[j] = values[k];
2447ea844a1aSMatthew Knepley         }
2448ea844a1aSMatthew Knepley         offset += dim;
2449ea844a1aSMatthew Knepley         cOffset += dim - tDim;
2450ea844a1aSMatthew Knepley       }
2451ea844a1aSMatthew Knepley     }
2452ea844a1aSMatthew Knepley   }
2453ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2454ea844a1aSMatthew Knepley }
2455ea844a1aSMatthew Knepley 
2456ea844a1aSMatthew Knepley /*@C
2457ea844a1aSMatthew Knepley   PetscSectionHasConstraints - Determine whether a section has constrained dofs
2458ea844a1aSMatthew Knepley 
245940750d41SVaclav Hapla   Not Collective
2460ea844a1aSMatthew Knepley 
2461ea844a1aSMatthew Knepley   Input Parameter:
2462cab54364SBarry Smith . s - The `PetscSection`
2463ea844a1aSMatthew Knepley 
2464ea844a1aSMatthew Knepley   Output Parameter:
2465ea844a1aSMatthew Knepley . hasConstraints - flag indicating that the section has constrained dofs
2466ea844a1aSMatthew Knepley 
2467ea844a1aSMatthew Knepley   Level: intermediate
2468ea844a1aSMatthew Knepley 
2469cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSetConstraintIndices()`, `PetscSectionGetConstraintDof()`, `PetscSection`
2470ea844a1aSMatthew Knepley @*/
2471d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionHasConstraints(PetscSection s, PetscBool *hasConstraints)
2472d71ae5a4SJacob Faibussowitsch {
2473ea844a1aSMatthew Knepley   PetscFunctionBegin;
2474ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
2475dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(hasConstraints, 2);
2476ea844a1aSMatthew Knepley   *hasConstraints = s->bc ? PETSC_TRUE : PETSC_FALSE;
2477ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2478ea844a1aSMatthew Knepley }
2479ea844a1aSMatthew Knepley 
2480ea844a1aSMatthew Knepley /*@C
2481ea844a1aSMatthew Knepley   PetscSectionGetConstraintIndices - Get the point dof numbers, in [0, dof), which are constrained
2482ea844a1aSMatthew Knepley 
248340750d41SVaclav Hapla   Not Collective
2484ea844a1aSMatthew Knepley 
2485ea844a1aSMatthew Knepley   Input Parameters:
2486cab54364SBarry Smith + s     - The `PetscSection`
2487ea844a1aSMatthew Knepley - point - The point
2488ea844a1aSMatthew Knepley 
2489ea844a1aSMatthew Knepley   Output Parameter:
2490ea844a1aSMatthew Knepley . indices - The constrained dofs
2491ea844a1aSMatthew Knepley 
2492ea844a1aSMatthew Knepley   Level: intermediate
2493ea844a1aSMatthew Knepley 
2494cab54364SBarry Smith   Fortran Note:
2495cab54364SBarry Smith   Call `PetscSectionGetConstraintIndicesF90()` and `PetscSectionRestoreConstraintIndicesF90()`
2496cab54364SBarry Smith 
2497cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSetConstraintIndices()`, `PetscSectionGetConstraintDof()`, `PetscSection`
2498ea844a1aSMatthew Knepley @*/
2499d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetConstraintIndices(PetscSection s, PetscInt point, const PetscInt **indices)
2500d71ae5a4SJacob Faibussowitsch {
2501ea844a1aSMatthew Knepley   PetscFunctionBegin;
2502ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
2503ea844a1aSMatthew Knepley   if (s->bc) {
25049566063dSJacob Faibussowitsch     PetscCall(VecIntGetValuesSection(s->bcIndices, s->bc, point, indices));
2505ea844a1aSMatthew Knepley   } else *indices = NULL;
2506ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2507ea844a1aSMatthew Knepley }
2508ea844a1aSMatthew Knepley 
2509ea844a1aSMatthew Knepley /*@C
2510ea844a1aSMatthew Knepley   PetscSectionSetConstraintIndices - Set the point dof numbers, in [0, dof), which are constrained
2511ea844a1aSMatthew Knepley 
251240750d41SVaclav Hapla   Not Collective
2513ea844a1aSMatthew Knepley 
2514ea844a1aSMatthew Knepley   Input Parameters:
2515cab54364SBarry Smith + s     - The `PetscSection`
2516ea844a1aSMatthew Knepley . point - The point
2517ea844a1aSMatthew Knepley - indices - The constrained dofs
2518ea844a1aSMatthew Knepley 
2519ea844a1aSMatthew Knepley   Level: intermediate
2520ea844a1aSMatthew Knepley 
2521cab54364SBarry Smith   Fortran Note:
2522cab54364SBarry Smith   Use `PetscSectionSetConstraintIndicesF90()`
2523cab54364SBarry Smith 
2524cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionGetConstraintIndices()`, `PetscSectionGetConstraintDof()`, `PetscSection`
2525ea844a1aSMatthew Knepley @*/
2526d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetConstraintIndices(PetscSection s, PetscInt point, const PetscInt indices[])
2527d71ae5a4SJacob Faibussowitsch {
2528ea844a1aSMatthew Knepley   PetscFunctionBegin;
2529ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
2530ea844a1aSMatthew Knepley   if (s->bc) {
2531d57bb9dbSMatthew G. Knepley     const PetscInt dof  = s->atlasDof[point];
2532d57bb9dbSMatthew G. Knepley     const PetscInt cdof = s->bc->atlasDof[point];
2533d57bb9dbSMatthew G. Knepley     PetscInt       d;
2534d57bb9dbSMatthew G. Knepley 
2535ad540459SPierre Jolivet     for (d = 0; d < cdof; ++d) PetscCheck(indices[d] < dof, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " dof %" PetscInt_FMT ", invalid constraint index[%" PetscInt_FMT "]: %" PetscInt_FMT, point, dof, d, indices[d]);
25369566063dSJacob Faibussowitsch     PetscCall(VecIntSetValuesSection(s->bcIndices, s->bc, point, indices, INSERT_VALUES));
2537ea844a1aSMatthew Knepley   }
2538ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2539ea844a1aSMatthew Knepley }
2540ea844a1aSMatthew Knepley 
2541ea844a1aSMatthew Knepley /*@C
2542ea844a1aSMatthew Knepley   PetscSectionGetFieldConstraintIndices - Get the field dof numbers, in [0, fdof), which are constrained
2543ea844a1aSMatthew Knepley 
254440750d41SVaclav Hapla   Not Collective
2545ea844a1aSMatthew Knepley 
2546ea844a1aSMatthew Knepley   Input Parameters:
2547cab54364SBarry Smith + s     - The `PetscSection`
2548ea844a1aSMatthew Knepley . field  - The field number
2549ea844a1aSMatthew Knepley - point - The point
2550ea844a1aSMatthew Knepley 
2551ea844a1aSMatthew Knepley   Output Parameter:
25529759eb26SJed Brown . indices - The constrained dofs sorted in ascending order
2553ea844a1aSMatthew Knepley 
2554ea844a1aSMatthew Knepley   Level: intermediate
2555ea844a1aSMatthew Knepley 
2556cab54364SBarry Smith   Note:
2557cab54364SBarry Smith   The indices array, which is provided by the caller, must have capacity to hold the number of constrained dofs, e.g., as returned by `PetscSectionGetConstraintDof()`.
2558cab54364SBarry Smith 
2559cab54364SBarry Smith   Fortran Note:
2560cab54364SBarry Smith   Call `PetscSectionGetFieldConstraintIndicesF90()` and `PetscSectionRestoreFieldConstraintIndicesF90()`
2561cab54364SBarry Smith 
2562cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSetFieldConstraintIndices()`, `PetscSectionGetConstraintIndices()`, `PetscSectionGetConstraintDof()`, `PetscSection`
2563ea844a1aSMatthew Knepley @*/
2564d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetFieldConstraintIndices(PetscSection s, PetscInt point, PetscInt field, const PetscInt **indices)
2565d71ae5a4SJacob Faibussowitsch {
2566ea844a1aSMatthew Knepley   PetscFunctionBegin;
2567ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
2568dadcf809SJacob Faibussowitsch   PetscValidPointer(indices, 4);
25692abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
25709566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstraintIndices(s->field[field], point, indices));
2571ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2572ea844a1aSMatthew Knepley }
2573ea844a1aSMatthew Knepley 
2574ea844a1aSMatthew Knepley /*@C
2575ea844a1aSMatthew Knepley   PetscSectionSetFieldConstraintIndices - Set the field dof numbers, in [0, fdof), which are constrained
2576ea844a1aSMatthew Knepley 
257740750d41SVaclav Hapla   Not Collective
2578ea844a1aSMatthew Knepley 
2579ea844a1aSMatthew Knepley   Input Parameters:
2580cab54364SBarry Smith + s       - The `PetscSection`
2581ea844a1aSMatthew Knepley . point   - The point
2582ea844a1aSMatthew Knepley . field   - The field number
2583ea844a1aSMatthew Knepley - indices - The constrained dofs
2584ea844a1aSMatthew Knepley 
2585ea844a1aSMatthew Knepley   Level: intermediate
2586ea844a1aSMatthew Knepley 
2587cab54364SBarry Smith   Fortran Note:
2588cab54364SBarry Smith   Use `PetscSectionSetFieldConstraintIndicesF90()`
2589cab54364SBarry Smith 
2590cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSetConstraintIndices()`, `PetscSectionGetFieldConstraintIndices()`, `PetscSectionGetConstraintDof()`, `PetscSection`
2591ea844a1aSMatthew Knepley @*/
2592d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetFieldConstraintIndices(PetscSection s, PetscInt point, PetscInt field, const PetscInt indices[])
2593d71ae5a4SJacob Faibussowitsch {
2594ea844a1aSMatthew Knepley   PetscFunctionBegin;
2595ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
25962abc8c78SJacob Faibussowitsch   if (PetscDefined(USE_DEBUG)) {
2597e2bfaee7SMatthew G. Knepley     PetscInt nfdof;
25982abc8c78SJacob Faibussowitsch 
25999566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldConstraintDof(s, point, field, &nfdof));
2600e2bfaee7SMatthew G. Knepley     if (nfdof) PetscValidIntPointer(indices, 4);
26012abc8c78SJacob Faibussowitsch   }
26022abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
26039566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetConstraintIndices(s->field[field], point, indices));
2604ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2605ea844a1aSMatthew Knepley }
2606ea844a1aSMatthew Knepley 
2607ea844a1aSMatthew Knepley /*@
2608ea844a1aSMatthew Knepley   PetscSectionPermute - Reorder the section according to the input point permutation
2609ea844a1aSMatthew Knepley 
261040750d41SVaclav Hapla   Collective
2611ea844a1aSMatthew Knepley 
2612d8d19677SJose E. Roman   Input Parameters:
2613cab54364SBarry Smith + section - The `PetscSection` object
2614ea844a1aSMatthew Knepley - perm - The point permutation, old point p becomes new point perm[p]
2615ea844a1aSMatthew Knepley 
2616ea844a1aSMatthew Knepley   Output Parameter:
2617cab54364SBarry Smith . sectionNew - The permuted `PetscSection`
2618ea844a1aSMatthew Knepley 
2619ea844a1aSMatthew Knepley   Level: intermediate
2620ea844a1aSMatthew Knepley 
2621cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `IS`, `PetscSection`, `MatPermute()`
2622ea844a1aSMatthew Knepley @*/
2623d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionPermute(PetscSection section, IS permutation, PetscSection *sectionNew)
2624d71ae5a4SJacob Faibussowitsch {
2625ea844a1aSMatthew Knepley   PetscSection    s = section, sNew;
2626ea844a1aSMatthew Knepley   const PetscInt *perm;
2627b778fa18SValeria Barra   PetscInt        numFields, f, c, numPoints, pStart, pEnd, p;
2628ea844a1aSMatthew Knepley 
2629ea844a1aSMatthew Knepley   PetscFunctionBegin;
2630ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1);
2631ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(permutation, IS_CLASSID, 2);
2632ea844a1aSMatthew Knepley   PetscValidPointer(sectionNew, 3);
26339566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)s), &sNew));
26349566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(s, &numFields));
26359566063dSJacob Faibussowitsch   if (numFields) PetscCall(PetscSectionSetNumFields(sNew, numFields));
2636ea844a1aSMatthew Knepley   for (f = 0; f < numFields; ++f) {
2637ea844a1aSMatthew Knepley     const char *name;
2638ea844a1aSMatthew Knepley     PetscInt    numComp;
2639ea844a1aSMatthew Knepley 
26409566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldName(s, f, &name));
26419566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldName(sNew, f, name));
26429566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldComponents(s, f, &numComp));
26439566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldComponents(sNew, f, numComp));
2644b778fa18SValeria Barra     for (c = 0; c < s->numFieldComponents[f]; ++c) {
26459566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetComponentName(s, f, c, &name));
26469566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetComponentName(sNew, f, c, name));
2647b778fa18SValeria Barra     }
2648ea844a1aSMatthew Knepley   }
26499566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(permutation, &numPoints));
26509566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(permutation, &perm));
26519566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
26529566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(sNew, pStart, pEnd));
265308401ef6SPierre Jolivet   PetscCheck(numPoints >= pEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Permutation size %" PetscInt_FMT " is less than largest Section point %" PetscInt_FMT, numPoints, pEnd);
2654ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
2655ea844a1aSMatthew Knepley     PetscInt dof, cdof;
2656ea844a1aSMatthew Knepley 
26579566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(s, p, &dof));
26589566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(sNew, perm[p], dof));
26599566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
26609566063dSJacob Faibussowitsch     if (cdof) PetscCall(PetscSectionSetConstraintDof(sNew, perm[p], cdof));
2661ea844a1aSMatthew Knepley     for (f = 0; f < numFields; ++f) {
26629566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldDof(s, p, f, &dof));
26639566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetFieldDof(sNew, perm[p], f, dof));
26649566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof));
26659566063dSJacob Faibussowitsch       if (cdof) PetscCall(PetscSectionSetFieldConstraintDof(sNew, perm[p], f, cdof));
2666ea844a1aSMatthew Knepley     }
2667ea844a1aSMatthew Knepley   }
26689566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(sNew));
2669ea844a1aSMatthew Knepley   for (p = pStart; p < pEnd; ++p) {
2670ea844a1aSMatthew Knepley     const PetscInt *cind;
2671ea844a1aSMatthew Knepley     PetscInt        cdof;
2672ea844a1aSMatthew Knepley 
26739566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
2674ea844a1aSMatthew Knepley     if (cdof) {
26759566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintIndices(s, p, &cind));
26769566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetConstraintIndices(sNew, perm[p], cind));
2677ea844a1aSMatthew Knepley     }
2678ea844a1aSMatthew Knepley     for (f = 0; f < numFields; ++f) {
26799566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof));
2680ea844a1aSMatthew Knepley       if (cdof) {
26819566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetFieldConstraintIndices(s, p, f, &cind));
26829566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetFieldConstraintIndices(sNew, perm[p], f, cind));
2683ea844a1aSMatthew Knepley       }
2684ea844a1aSMatthew Knepley     }
2685ea844a1aSMatthew Knepley   }
26869566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(permutation, &perm));
2687ea844a1aSMatthew Knepley   *sectionNew = sNew;
2688ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2689ea844a1aSMatthew Knepley }
2690ea844a1aSMatthew Knepley 
2691ea844a1aSMatthew Knepley /*@
2692ea844a1aSMatthew Knepley   PetscSectionSetClosureIndex - Set a cache of points in the closure of each point in the section
2693ea844a1aSMatthew Knepley 
269440750d41SVaclav Hapla   Collective
2695ea844a1aSMatthew Knepley 
2696ea844a1aSMatthew Knepley   Input Parameters:
2697cab54364SBarry Smith + section   - The `PetscSection`
2698cab54364SBarry Smith . obj       - A `PetscObject` which serves as the key for this index
2699cab54364SBarry Smith . clSection - `PetscSection` giving the size of the closure of each point
2700cab54364SBarry Smith - clPoints  - `IS` giving the points in each closure
2701ea844a1aSMatthew Knepley 
2702ea844a1aSMatthew Knepley   Level: advanced
2703ea844a1aSMatthew Knepley 
2704cab54364SBarry Smith   Note:
2705cab54364SBarry Smith   We compress out closure points with no dofs in this section
2706cab54364SBarry Smith 
2707cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionGetClosureIndex()`, `DMPlexCreateClosureIndex()`
2708ea844a1aSMatthew Knepley @*/
2709d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetClosureIndex(PetscSection section, PetscObject obj, PetscSection clSection, IS clPoints)
2710d71ae5a4SJacob Faibussowitsch {
2711ea844a1aSMatthew Knepley   PetscFunctionBegin;
27123e03ebd8SStefano Zampini   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1);
27133e03ebd8SStefano Zampini   PetscValidHeaderSpecific(clSection, PETSC_SECTION_CLASSID, 3);
27143e03ebd8SStefano Zampini   PetscValidHeaderSpecific(clPoints, IS_CLASSID, 4);
27159566063dSJacob Faibussowitsch   if (section->clObj != obj) PetscCall(PetscSectionResetClosurePermutation(section));
2716ea844a1aSMatthew Knepley   section->clObj = obj;
27179566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)clSection));
27189566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)clPoints));
27199566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&section->clSection));
27209566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&section->clPoints));
2721ea844a1aSMatthew Knepley   section->clSection = clSection;
2722ea844a1aSMatthew Knepley   section->clPoints  = clPoints;
2723ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2724ea844a1aSMatthew Knepley }
2725ea844a1aSMatthew Knepley 
2726ea844a1aSMatthew Knepley /*@
2727cab54364SBarry Smith   PetscSectionGetClosureIndex - Get the cache of points in the closure of each point in the section set with `PetscSectionSetClosureIndex()`
2728ea844a1aSMatthew Knepley 
272940750d41SVaclav Hapla   Collective
2730ea844a1aSMatthew Knepley 
2731ea844a1aSMatthew Knepley   Input Parameters:
2732cab54364SBarry Smith + section   - The `PetscSection`
2733cab54364SBarry Smith - obj       - A `PetscObject` which serves as the key for this index
2734ea844a1aSMatthew Knepley 
2735ea844a1aSMatthew Knepley   Output Parameters:
2736cab54364SBarry Smith + clSection - `PetscSection` giving the size of the closure of each point
2737cab54364SBarry Smith - clPoints  - `IS` giving the points in each closure
2738ea844a1aSMatthew Knepley 
2739ea844a1aSMatthew Knepley   Level: advanced
2740ea844a1aSMatthew Knepley 
2741cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSetClosureIndex()`, `DMPlexCreateClosureIndex()`
2742ea844a1aSMatthew Knepley @*/
2743d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetClosureIndex(PetscSection section, PetscObject obj, PetscSection *clSection, IS *clPoints)
2744d71ae5a4SJacob Faibussowitsch {
2745ea844a1aSMatthew Knepley   PetscFunctionBegin;
2746ea844a1aSMatthew Knepley   if (section->clObj == obj) {
2747ea844a1aSMatthew Knepley     if (clSection) *clSection = section->clSection;
2748ea844a1aSMatthew Knepley     if (clPoints) *clPoints = section->clPoints;
2749ea844a1aSMatthew Knepley   } else {
2750ea844a1aSMatthew Knepley     if (clSection) *clSection = NULL;
2751ea844a1aSMatthew Knepley     if (clPoints) *clPoints = NULL;
2752ea844a1aSMatthew Knepley   }
2753ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2754ea844a1aSMatthew Knepley }
2755ea844a1aSMatthew Knepley 
2756d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetClosurePermutation_Internal(PetscSection section, PetscObject obj, PetscInt depth, PetscInt clSize, PetscCopyMode mode, PetscInt *clPerm)
2757d71ae5a4SJacob Faibussowitsch {
2758ea844a1aSMatthew Knepley   PetscInt                    i;
2759c459fbc1SJed Brown   khiter_t                    iter;
2760c459fbc1SJed Brown   int                         new_entry;
2761c459fbc1SJed Brown   PetscSectionClosurePermKey  key = {depth, clSize};
2762c459fbc1SJed Brown   PetscSectionClosurePermVal *val;
2763ea844a1aSMatthew Knepley 
2764ea844a1aSMatthew Knepley   PetscFunctionBegin;
2765ea844a1aSMatthew Knepley   if (section->clObj != obj) {
27669566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&section->clSection));
27679566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&section->clPoints));
2768ea844a1aSMatthew Knepley   }
2769ea844a1aSMatthew Knepley   section->clObj = obj;
27709566063dSJacob Faibussowitsch   if (!section->clHash) PetscCall(PetscClPermCreate(&section->clHash));
2771c459fbc1SJed Brown   iter = kh_put(ClPerm, section->clHash, key, &new_entry);
2772c459fbc1SJed Brown   val  = &kh_val(section->clHash, iter);
2773c459fbc1SJed Brown   if (!new_entry) {
27749566063dSJacob Faibussowitsch     PetscCall(PetscFree(val->perm));
27759566063dSJacob Faibussowitsch     PetscCall(PetscFree(val->invPerm));
2776c459fbc1SJed Brown   }
2777ea844a1aSMatthew Knepley   if (mode == PETSC_COPY_VALUES) {
27789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(clSize, &val->perm));
27799566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(val->perm, clPerm, clSize));
2780ea844a1aSMatthew Knepley   } else if (mode == PETSC_OWN_POINTER) {
2781c459fbc1SJed Brown     val->perm = clPerm;
2782ea844a1aSMatthew Knepley   } else SETERRQ(PetscObjectComm(obj), PETSC_ERR_SUP, "Do not support borrowed arrays");
27839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(clSize, &val->invPerm));
2784c459fbc1SJed Brown   for (i = 0; i < clSize; ++i) val->invPerm[clPerm[i]] = i;
2785ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2786ea844a1aSMatthew Knepley }
2787ea844a1aSMatthew Knepley 
2788ea844a1aSMatthew Knepley /*@
2789c459fbc1SJed Brown   PetscSectionSetClosurePermutation - Set the dof permutation for the closure of each cell in the section, meaning clPerm[newIndex] = oldIndex.
2790ea844a1aSMatthew Knepley 
2791ea844a1aSMatthew Knepley   Not Collective
2792ea844a1aSMatthew Knepley 
2793ea844a1aSMatthew Knepley   Input Parameters:
2794cab54364SBarry Smith + section - The `PetscSection`
2795cab54364SBarry Smith . obj     - A `PetscObject` which serves as the key for this index (usually a `DM`)
2796c459fbc1SJed Brown . depth   - Depth of points on which to apply the given permutation
2797ea844a1aSMatthew Knepley - perm    - Permutation of the cell dof closure
2798ea844a1aSMatthew Knepley 
2799cab54364SBarry Smith   Level: intermediate
2800cab54364SBarry Smith 
2801cab54364SBarry Smith   Notes:
2802c459fbc1SJed Brown   The specified permutation will only be applied to points at depth whose closure size matches the length of perm.  In a
2803c459fbc1SJed Brown   mixed-topology or variable-degree finite element space, this function can be called multiple times at each depth for
2804c459fbc1SJed Brown   each topology and degree.
2805c459fbc1SJed Brown 
2806c459fbc1SJed Brown   This approach assumes that (depth, len(perm)) uniquely identifies the desired permutation; this might not be true for
2807c459fbc1SJed Brown   exotic/enriched spaces on mixed topology meshes.
2808ea844a1aSMatthew Knepley 
2809cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `IS`, `PetscSectionGetClosurePermutation()`, `PetscSectionGetClosureIndex()`, `DMPlexCreateClosureIndex()`, `PetscCopyMode`
2810ea844a1aSMatthew Knepley @*/
2811d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetClosurePermutation(PetscSection section, PetscObject obj, PetscInt depth, IS perm)
2812d71ae5a4SJacob Faibussowitsch {
2813ea844a1aSMatthew Knepley   const PetscInt *clPerm = NULL;
2814ea844a1aSMatthew Knepley   PetscInt        clSize = 0;
2815ea844a1aSMatthew Knepley 
2816ea844a1aSMatthew Knepley   PetscFunctionBegin;
2817ea844a1aSMatthew Knepley   if (perm) {
28189566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(perm, &clSize));
28199566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(perm, &clPerm));
2820ea844a1aSMatthew Knepley   }
28219566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetClosurePermutation_Internal(section, obj, depth, clSize, PETSC_COPY_VALUES, (PetscInt *)clPerm));
28229566063dSJacob Faibussowitsch   if (perm) PetscCall(ISRestoreIndices(perm, &clPerm));
2823ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2824ea844a1aSMatthew Knepley }
2825ea844a1aSMatthew Knepley 
2826d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetClosurePermutation_Internal(PetscSection section, PetscObject obj, PetscInt depth, PetscInt size, const PetscInt *perm[])
2827d71ae5a4SJacob Faibussowitsch {
2828ea844a1aSMatthew Knepley   PetscFunctionBegin;
2829ea844a1aSMatthew Knepley   if (section->clObj == obj) {
2830c459fbc1SJed Brown     PetscSectionClosurePermKey k = {depth, size};
2831c459fbc1SJed Brown     PetscSectionClosurePermVal v;
28329566063dSJacob Faibussowitsch     PetscCall(PetscClPermGet(section->clHash, k, &v));
2833c459fbc1SJed Brown     if (perm) *perm = v.perm;
2834ea844a1aSMatthew Knepley   } else {
2835ea844a1aSMatthew Knepley     if (perm) *perm = NULL;
2836ea844a1aSMatthew Knepley   }
2837ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2838ea844a1aSMatthew Knepley }
2839ea844a1aSMatthew Knepley 
2840ea844a1aSMatthew Knepley /*@
2841ea844a1aSMatthew Knepley   PetscSectionGetClosurePermutation - Get the dof permutation for the closure of each cell in the section, meaning clPerm[newIndex] = oldIndex.
2842ea844a1aSMatthew Knepley 
284340750d41SVaclav Hapla   Not Collective
2844ea844a1aSMatthew Knepley 
2845ea844a1aSMatthew Knepley   Input Parameters:
2846cab54364SBarry Smith + section   - The `PetscSection`
2847cab54364SBarry Smith . obj       - A `PetscObject` which serves as the key for this index (usually a DM)
2848c459fbc1SJed Brown . depth     - Depth stratum on which to obtain closure permutation
2849c459fbc1SJed Brown - clSize    - Closure size to be permuted (e.g., may vary with element topology and degree)
2850ea844a1aSMatthew Knepley 
2851ea844a1aSMatthew Knepley   Output Parameter:
2852ea844a1aSMatthew Knepley . perm - The dof closure permutation
2853ea844a1aSMatthew Knepley 
2854ea844a1aSMatthew Knepley   Level: intermediate
2855ea844a1aSMatthew Knepley 
2856cab54364SBarry Smith   Note:
2857cab54364SBarry Smith   The user must destroy the `IS` that is returned.
2858cab54364SBarry Smith 
2859cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `IS`, `PetscSectionSetClosurePermutation()`, `PetscSectionGetClosureInversePermutation()`, `PetscSectionGetClosureIndex()`, `PetscSectionSetClosureIndex()`, `DMPlexCreateClosureIndex()`
2860ea844a1aSMatthew Knepley @*/
2861d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetClosurePermutation(PetscSection section, PetscObject obj, PetscInt depth, PetscInt clSize, IS *perm)
2862d71ae5a4SJacob Faibussowitsch {
2863ea844a1aSMatthew Knepley   const PetscInt *clPerm;
2864ea844a1aSMatthew Knepley 
2865ea844a1aSMatthew Knepley   PetscFunctionBegin;
28669566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetClosurePermutation_Internal(section, obj, depth, clSize, &clPerm));
28679566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, clSize, clPerm, PETSC_USE_POINTER, perm));
2868ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2869ea844a1aSMatthew Knepley }
2870ea844a1aSMatthew Knepley 
2871d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetClosureInversePermutation_Internal(PetscSection section, PetscObject obj, PetscInt depth, PetscInt size, const PetscInt *perm[])
2872d71ae5a4SJacob Faibussowitsch {
2873ea844a1aSMatthew Knepley   PetscFunctionBegin;
2874c459fbc1SJed Brown   if (section->clObj == obj && section->clHash) {
2875c459fbc1SJed Brown     PetscSectionClosurePermKey k = {depth, size};
2876c459fbc1SJed Brown     PetscSectionClosurePermVal v;
28779566063dSJacob Faibussowitsch     PetscCall(PetscClPermGet(section->clHash, k, &v));
2878c459fbc1SJed Brown     if (perm) *perm = v.invPerm;
2879ea844a1aSMatthew Knepley   } else {
2880ea844a1aSMatthew Knepley     if (perm) *perm = NULL;
2881ea844a1aSMatthew Knepley   }
2882ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2883ea844a1aSMatthew Knepley }
2884ea844a1aSMatthew Knepley 
2885ea844a1aSMatthew Knepley /*@
2886ea844a1aSMatthew Knepley   PetscSectionGetClosureInversePermutation - Get the inverse dof permutation for the closure of each cell in the section, meaning clPerm[oldIndex] = newIndex.
2887ea844a1aSMatthew Knepley 
288840750d41SVaclav Hapla   Not Collective
2889ea844a1aSMatthew Knepley 
2890ea844a1aSMatthew Knepley   Input Parameters:
2891cab54364SBarry Smith + section   - The `PetscSection`
2892cab54364SBarry Smith . obj       - A `PetscObject` which serves as the key for this index (usually a `DM`)
2893c459fbc1SJed Brown . depth     - Depth stratum on which to obtain closure permutation
2894c459fbc1SJed Brown - clSize    - Closure size to be permuted (e.g., may vary with element topology and degree)
2895ea844a1aSMatthew Knepley 
2896ea844a1aSMatthew Knepley   Output Parameters:
2897c459fbc1SJed Brown . perm - The dof closure permutation
2898ea844a1aSMatthew Knepley 
2899ea844a1aSMatthew Knepley   Level: intermediate
2900ea844a1aSMatthew Knepley 
2901cab54364SBarry Smith   Note:
2902cab54364SBarry Smith   The user must destroy the `IS` that is returned.
2903cab54364SBarry Smith 
2904cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `IS`, `PetscSectionSetClosurePermutation()`, `PetscSectionGetClosureIndex()`, `PetscSectionSetClosureIndex()`, `DMPlexCreateClosureIndex()`
2905ea844a1aSMatthew Knepley @*/
2906d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetClosureInversePermutation(PetscSection section, PetscObject obj, PetscInt depth, PetscInt clSize, IS *perm)
2907d71ae5a4SJacob Faibussowitsch {
2908ea844a1aSMatthew Knepley   const PetscInt *clPerm;
2909ea844a1aSMatthew Knepley 
2910ea844a1aSMatthew Knepley   PetscFunctionBegin;
29119566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetClosureInversePermutation_Internal(section, obj, depth, clSize, &clPerm));
29129566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, clSize, clPerm, PETSC_USE_POINTER, perm));
2913ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2914ea844a1aSMatthew Knepley }
2915ea844a1aSMatthew Knepley 
2916ea844a1aSMatthew Knepley /*@
2917ea844a1aSMatthew Knepley   PetscSectionGetField - Get the subsection associated with a single field
2918ea844a1aSMatthew Knepley 
2919ea844a1aSMatthew Knepley   Input Parameters:
2920cab54364SBarry Smith + s     - The `PetscSection`
2921ea844a1aSMatthew Knepley - field - The field number
2922ea844a1aSMatthew Knepley 
2923ea844a1aSMatthew Knepley   Output Parameter:
2924ea844a1aSMatthew Knepley . subs  - The subsection for the given field
2925ea844a1aSMatthew Knepley 
2926ea844a1aSMatthew Knepley   Level: intermediate
2927ea844a1aSMatthew Knepley 
2928cab54364SBarry Smith   Note:
2929cab54364SBarry Smith   Does not increase the reference count of the selected sub-section. There is no matching `PetscSectionRestoreField()`
2930cab54364SBarry Smith 
2931cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `IS`, `PetscSectionSetNumFields()`
2932ea844a1aSMatthew Knepley @*/
2933d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetField(PetscSection s, PetscInt field, PetscSection *subs)
2934d71ae5a4SJacob Faibussowitsch {
2935ea844a1aSMatthew Knepley   PetscFunctionBegin;
2936ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
2937ea844a1aSMatthew Knepley   PetscValidPointer(subs, 3);
29382abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, s->numFields);
2939ea844a1aSMatthew Knepley   *subs = s->field[field];
2940ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2941ea844a1aSMatthew Knepley }
2942ea844a1aSMatthew Knepley 
2943ea844a1aSMatthew Knepley PetscClassId      PETSC_SECTION_SYM_CLASSID;
2944ea844a1aSMatthew Knepley PetscFunctionList PetscSectionSymList = NULL;
2945ea844a1aSMatthew Knepley 
2946ea844a1aSMatthew Knepley /*@
2947cab54364SBarry Smith   PetscSectionSymCreate - Creates an empty `PetscSectionSym` object.
2948ea844a1aSMatthew Knepley 
2949ea844a1aSMatthew Knepley   Collective
2950ea844a1aSMatthew Knepley 
2951ea844a1aSMatthew Knepley   Input Parameter:
2952ea844a1aSMatthew Knepley . comm - the MPI communicator
2953ea844a1aSMatthew Knepley 
2954ea844a1aSMatthew Knepley   Output Parameter:
2955ea844a1aSMatthew Knepley . sym - pointer to the new set of symmetries
2956ea844a1aSMatthew Knepley 
2957ea844a1aSMatthew Knepley   Level: developer
2958ea844a1aSMatthew Knepley 
2959cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSection`, `PetscSectionSym`, `PetscSectionSymDestroy()`
2960ea844a1aSMatthew Knepley @*/
2961d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSymCreate(MPI_Comm comm, PetscSectionSym *sym)
2962d71ae5a4SJacob Faibussowitsch {
2963ea844a1aSMatthew Knepley   PetscFunctionBegin;
2964ea844a1aSMatthew Knepley   PetscValidPointer(sym, 2);
29659566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
29669566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*sym, PETSC_SECTION_SYM_CLASSID, "PetscSectionSym", "Section Symmetry", "IS", comm, PetscSectionSymDestroy, PetscSectionSymView));
2967ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
2968ea844a1aSMatthew Knepley }
2969ea844a1aSMatthew Knepley 
2970ea844a1aSMatthew Knepley /*@C
2971cab54364SBarry Smith   PetscSectionSymSetType - Builds a `PetscSectionSym`, for a particular implementation.
2972ea844a1aSMatthew Knepley 
297340750d41SVaclav Hapla   Collective
2974ea844a1aSMatthew Knepley 
2975ea844a1aSMatthew Knepley   Input Parameters:
2976ea844a1aSMatthew Knepley + sym    - The section symmetry object
2977ea844a1aSMatthew Knepley - method - The name of the section symmetry type
2978ea844a1aSMatthew Knepley 
2979ea844a1aSMatthew Knepley   Level: developer
2980ea844a1aSMatthew Knepley 
2981cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSym`, `PetscSectionSymType`, `PetscSectionSymGetType()`, `PetscSectionSymCreate()`
2982ea844a1aSMatthew Knepley @*/
2983d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSymSetType(PetscSectionSym sym, PetscSectionSymType method)
2984d71ae5a4SJacob Faibussowitsch {
2985ea844a1aSMatthew Knepley   PetscErrorCode (*r)(PetscSectionSym);
2986ea844a1aSMatthew Knepley   PetscBool match;
2987ea844a1aSMatthew Knepley 
2988ea844a1aSMatthew Knepley   PetscFunctionBegin;
2989ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1);
29909566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)sym, method, &match));
2991ea844a1aSMatthew Knepley   if (match) PetscFunctionReturn(0);
2992ea844a1aSMatthew Knepley 
29939566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(PetscSectionSymList, method, &r));
299428b400f6SJacob Faibussowitsch   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscSectionSym type: %s", method);
2995dbbe0bcdSBarry Smith   PetscTryTypeMethod(sym, destroy);
2996ea844a1aSMatthew Knepley   sym->ops->destroy = NULL;
2997dbbe0bcdSBarry Smith 
29989566063dSJacob Faibussowitsch   PetscCall((*r)(sym));
29999566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)sym, method));
3000ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3001ea844a1aSMatthew Knepley }
3002ea844a1aSMatthew Knepley 
3003ea844a1aSMatthew Knepley /*@C
3004cab54364SBarry Smith   PetscSectionSymGetType - Gets the section symmetry type name (as a string) from the `PetscSectionSym`.
3005ea844a1aSMatthew Knepley 
3006ea844a1aSMatthew Knepley   Not Collective
3007ea844a1aSMatthew Knepley 
3008ea844a1aSMatthew Knepley   Input Parameter:
3009ea844a1aSMatthew Knepley . sym  - The section symmetry
3010ea844a1aSMatthew Knepley 
3011ea844a1aSMatthew Knepley   Output Parameter:
3012ea844a1aSMatthew Knepley . type - The index set type name
3013ea844a1aSMatthew Knepley 
3014ea844a1aSMatthew Knepley   Level: developer
3015ea844a1aSMatthew Knepley 
3016cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSym`, `PetscSectionSymType`, `PetscSectionSymSetType()`, `PetscSectionSymCreate()`
3017ea844a1aSMatthew Knepley @*/
3018d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSymGetType(PetscSectionSym sym, PetscSectionSymType *type)
3019d71ae5a4SJacob Faibussowitsch {
3020ea844a1aSMatthew Knepley   PetscFunctionBegin;
3021ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1);
3022dadcf809SJacob Faibussowitsch   PetscValidPointer(type, 2);
3023ea844a1aSMatthew Knepley   *type = ((PetscObject)sym)->type_name;
3024ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3025ea844a1aSMatthew Knepley }
3026ea844a1aSMatthew Knepley 
3027ea844a1aSMatthew Knepley /*@C
3028cab54364SBarry Smith   PetscSectionSymRegister - Registers a new section symmetry implementation
3029ea844a1aSMatthew Knepley 
3030ea844a1aSMatthew Knepley   Not Collective
3031ea844a1aSMatthew Knepley 
3032ea844a1aSMatthew Knepley   Input Parameters:
3033ea844a1aSMatthew Knepley + name        - The name of a new user-defined creation routine
3034ea844a1aSMatthew Knepley - create_func - The creation routine itself
3035ea844a1aSMatthew Knepley 
3036ea844a1aSMatthew Knepley   Level: developer
3037ea844a1aSMatthew Knepley 
3038cab54364SBarry Smith   Notes:
3039cab54364SBarry Smith   `PetscSectionSymRegister()` may be called multiple times to add several user-defined vectors
3040cab54364SBarry Smith 
3041cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSym`, `PetscSectionSymType`, `PetscSectionSymCreate()`, `PetscSectionSymSetType()`
3042ea844a1aSMatthew Knepley @*/
3043d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSymRegister(const char sname[], PetscErrorCode (*function)(PetscSectionSym))
3044d71ae5a4SJacob Faibussowitsch {
3045ea844a1aSMatthew Knepley   PetscFunctionBegin;
30469566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
30479566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PetscSectionSymList, sname, function));
3048ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3049ea844a1aSMatthew Knepley }
3050ea844a1aSMatthew Knepley 
3051ea844a1aSMatthew Knepley /*@
3052ea844a1aSMatthew Knepley    PetscSectionSymDestroy - Destroys a section symmetry.
3053ea844a1aSMatthew Knepley 
305440750d41SVaclav Hapla    Collective
3055ea844a1aSMatthew Knepley 
3056ea844a1aSMatthew Knepley    Input Parameters:
3057ea844a1aSMatthew Knepley .  sym - the section symmetry
3058ea844a1aSMatthew Knepley 
3059ea844a1aSMatthew Knepley    Level: developer
3060ea844a1aSMatthew Knepley 
3061cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection),  `PetscSectionSym`, `PetscSectionSymCreate()`, `PetscSectionSymDestroy()`
3062ea844a1aSMatthew Knepley @*/
3063d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSymDestroy(PetscSectionSym *sym)
3064d71ae5a4SJacob Faibussowitsch {
3065ea844a1aSMatthew Knepley   SymWorkLink link, next;
3066ea844a1aSMatthew Knepley 
3067ea844a1aSMatthew Knepley   PetscFunctionBegin;
3068ea844a1aSMatthew Knepley   if (!*sym) PetscFunctionReturn(0);
3069ea844a1aSMatthew Knepley   PetscValidHeaderSpecific((*sym), PETSC_SECTION_SYM_CLASSID, 1);
30709371c9d4SSatish Balay   if (--((PetscObject)(*sym))->refct > 0) {
30719371c9d4SSatish Balay     *sym = NULL;
30729371c9d4SSatish Balay     PetscFunctionReturn(0);
3073ea844a1aSMatthew Knepley   }
307448a46eb9SPierre Jolivet   if ((*sym)->ops->destroy) PetscCall((*(*sym)->ops->destroy)(*sym));
3075c9cc58a2SBarry Smith   PetscCheck(!(*sym)->workout, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Work array still checked out");
3076ea844a1aSMatthew Knepley   for (link = (*sym)->workin; link; link = next) {
30770c99d500SBarry Smith     PetscInt    **perms = (PetscInt **)link->perms;
30780c99d500SBarry Smith     PetscScalar **rots  = (PetscScalar **)link->rots;
30790c99d500SBarry Smith     PetscCall(PetscFree2(perms, rots));
3080ea844a1aSMatthew Knepley     next = link->next;
30819566063dSJacob Faibussowitsch     PetscCall(PetscFree(link));
3082ea844a1aSMatthew Knepley   }
3083ea844a1aSMatthew Knepley   (*sym)->workin = NULL;
30849566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(sym));
3085ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3086ea844a1aSMatthew Knepley }
3087ea844a1aSMatthew Knepley 
3088ea844a1aSMatthew Knepley /*@C
3089ea844a1aSMatthew Knepley    PetscSectionSymView - Displays a section symmetry
3090ea844a1aSMatthew Knepley 
309140750d41SVaclav Hapla    Collective
3092ea844a1aSMatthew Knepley 
3093ea844a1aSMatthew Knepley    Input Parameters:
3094ea844a1aSMatthew Knepley +  sym - the index set
3095cab54364SBarry Smith -  viewer - viewer used to display the set, for example `PETSC_VIEWER_STDOUT_SELF`.
3096ea844a1aSMatthew Knepley 
3097ea844a1aSMatthew Knepley    Level: developer
3098ea844a1aSMatthew Knepley 
3099cab54364SBarry Smith .seealso:  `PetscSectionSym`, `PetscViewer`, `PetscViewerASCIIOpen()`
3100ea844a1aSMatthew Knepley @*/
3101d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSymView(PetscSectionSym sym, PetscViewer viewer)
3102d71ae5a4SJacob Faibussowitsch {
3103ea844a1aSMatthew Knepley   PetscFunctionBegin;
3104ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1);
310548a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sym), &viewer));
3106ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
3107ea844a1aSMatthew Knepley   PetscCheckSameComm(sym, 1, viewer, 2);
31089566063dSJacob Faibussowitsch   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)sym, viewer));
3109dbbe0bcdSBarry Smith   PetscTryTypeMethod(sym, view, viewer);
3110ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3111ea844a1aSMatthew Knepley }
3112ea844a1aSMatthew Knepley 
3113ea844a1aSMatthew Knepley /*@
3114ea844a1aSMatthew Knepley   PetscSectionSetSym - Set the symmetries for the data referred to by the section
3115ea844a1aSMatthew Knepley 
311640750d41SVaclav Hapla   Collective
3117ea844a1aSMatthew Knepley 
3118ea844a1aSMatthew Knepley   Input Parameters:
3119ea844a1aSMatthew Knepley + section - the section describing data layout
3120ea844a1aSMatthew Knepley - sym - the symmetry describing the affect of orientation on the access of the data
3121ea844a1aSMatthew Knepley 
3122ea844a1aSMatthew Knepley   Level: developer
3123ea844a1aSMatthew Knepley 
3124cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection),  `PetscSectionSym`, `PetscSectionGetSym()`, `PetscSectionSymCreate()`
3125ea844a1aSMatthew Knepley @*/
3126d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetSym(PetscSection section, PetscSectionSym sym)
3127d71ae5a4SJacob Faibussowitsch {
3128ea844a1aSMatthew Knepley   PetscFunctionBegin;
3129ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1);
31309566063dSJacob Faibussowitsch   PetscCall(PetscSectionSymDestroy(&(section->sym)));
3131ea844a1aSMatthew Knepley   if (sym) {
3132ea844a1aSMatthew Knepley     PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 2);
3133ea844a1aSMatthew Knepley     PetscCheckSameComm(section, 1, sym, 2);
31349566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)sym));
3135ea844a1aSMatthew Knepley   }
3136ea844a1aSMatthew Knepley   section->sym = sym;
3137ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3138ea844a1aSMatthew Knepley }
3139ea844a1aSMatthew Knepley 
3140ea844a1aSMatthew Knepley /*@
3141ea844a1aSMatthew Knepley   PetscSectionGetSym - Get the symmetries for the data referred to by the section
3142ea844a1aSMatthew Knepley 
314340750d41SVaclav Hapla   Not Collective
3144ea844a1aSMatthew Knepley 
3145ea844a1aSMatthew Knepley   Input Parameters:
3146ea844a1aSMatthew Knepley . section - the section describing data layout
3147ea844a1aSMatthew Knepley 
3148ea844a1aSMatthew Knepley   Output Parameters:
3149ea844a1aSMatthew Knepley . sym - the symmetry describing the affect of orientation on the access of the data
3150ea844a1aSMatthew Knepley 
3151ea844a1aSMatthew Knepley   Level: developer
3152ea844a1aSMatthew Knepley 
3153cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSym`, `PetscSectionSetSym()`, `PetscSectionSymCreate()`
3154ea844a1aSMatthew Knepley @*/
3155d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetSym(PetscSection section, PetscSectionSym *sym)
3156d71ae5a4SJacob Faibussowitsch {
3157ea844a1aSMatthew Knepley   PetscFunctionBegin;
3158ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1);
3159ea844a1aSMatthew Knepley   *sym = section->sym;
3160ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3161ea844a1aSMatthew Knepley }
3162ea844a1aSMatthew Knepley 
3163ea844a1aSMatthew Knepley /*@
3164ea844a1aSMatthew Knepley   PetscSectionSetFieldSym - Set the symmetries for the data referred to by a field of the section
3165ea844a1aSMatthew Knepley 
316640750d41SVaclav Hapla   Collective
3167ea844a1aSMatthew Knepley 
3168ea844a1aSMatthew Knepley   Input Parameters:
3169ea844a1aSMatthew Knepley + section - the section describing data layout
3170ea844a1aSMatthew Knepley . field - the field number
3171ea844a1aSMatthew Knepley - sym - the symmetry describing the affect of orientation on the access of the data
3172ea844a1aSMatthew Knepley 
3173ea844a1aSMatthew Knepley   Level: developer
3174ea844a1aSMatthew Knepley 
3175cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSym`, `PetscSectionGetFieldSym()`, `PetscSectionSymCreate()`
3176ea844a1aSMatthew Knepley @*/
3177d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetFieldSym(PetscSection section, PetscInt field, PetscSectionSym sym)
3178d71ae5a4SJacob Faibussowitsch {
3179ea844a1aSMatthew Knepley   PetscFunctionBegin;
3180ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1);
31812abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, section->numFields);
31829566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetSym(section->field[field], sym));
3183ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3184ea844a1aSMatthew Knepley }
3185ea844a1aSMatthew Knepley 
3186ea844a1aSMatthew Knepley /*@
3187ea844a1aSMatthew Knepley   PetscSectionGetFieldSym - Get the symmetries for the data referred to by a field of the section
3188ea844a1aSMatthew Knepley 
318940750d41SVaclav Hapla   Collective
3190ea844a1aSMatthew Knepley 
3191ea844a1aSMatthew Knepley   Input Parameters:
3192ea844a1aSMatthew Knepley + section - the section describing data layout
3193ea844a1aSMatthew Knepley - field - the field number
3194ea844a1aSMatthew Knepley 
3195ea844a1aSMatthew Knepley   Output Parameters:
3196ea844a1aSMatthew Knepley . sym - the symmetry describing the affect of orientation on the access of the data
3197ea844a1aSMatthew Knepley 
3198ea844a1aSMatthew Knepley   Level: developer
3199ea844a1aSMatthew Knepley 
3200cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSym`, `PetscSectionSetFieldSym()`, `PetscSectionSymCreate()`
3201ea844a1aSMatthew Knepley @*/
3202d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetFieldSym(PetscSection section, PetscInt field, PetscSectionSym *sym)
3203d71ae5a4SJacob Faibussowitsch {
3204ea844a1aSMatthew Knepley   PetscFunctionBegin;
3205ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1);
32062abc8c78SJacob Faibussowitsch   PetscSectionCheckValidField(field, section->numFields);
3207ea844a1aSMatthew Knepley   *sym = section->field[field]->sym;
3208ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3209ea844a1aSMatthew Knepley }
3210ea844a1aSMatthew Knepley 
3211ea844a1aSMatthew Knepley /*@C
3212cab54364SBarry Smith   PetscSectionGetPointSyms - Get the symmetries for a set of points in a `PetscSection` under specific orientations.
3213ea844a1aSMatthew Knepley 
321440750d41SVaclav Hapla   Not Collective
3215ea844a1aSMatthew Knepley 
3216ea844a1aSMatthew Knepley   Input Parameters:
3217ea844a1aSMatthew Knepley + section - the section
3218ea844a1aSMatthew Knepley . numPoints - the number of points
3219ea844a1aSMatthew Knepley - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an
3220ea844a1aSMatthew Knepley     arbitrary integer: its interpretation is up to sym.  Orientations are used by DM: for their interpretation in that
3221ea844a1aSMatthew Knepley     context, see DMPlexGetConeOrientation()).
3222ea844a1aSMatthew Knepley 
3223d8d19677SJose E. Roman   Output Parameters:
3224ea844a1aSMatthew Knepley + perms - The permutations for the given orientations (or NULL if there is no symmetry or the permutation is the identity).
3225ea844a1aSMatthew Knepley - rots - The field rotations symmetries for the given orientations (or NULL if there is no symmetry or the rotations are all
3226ea844a1aSMatthew Knepley     identity).
3227ea844a1aSMatthew Knepley 
3228ea844a1aSMatthew Knepley   Example of usage, gathering dofs into a local array (lArray) from a section array (sArray):
3229ea844a1aSMatthew Knepley .vb
3230ea844a1aSMatthew Knepley      const PetscInt    **perms;
3231ea844a1aSMatthew Knepley      const PetscScalar **rots;
3232ea844a1aSMatthew Knepley      PetscInt            lOffset;
3233ea844a1aSMatthew Knepley 
3234ea844a1aSMatthew Knepley      PetscSectionGetPointSyms(section,numPoints,points,&perms,&rots);
3235ea844a1aSMatthew Knepley      for (i = 0, lOffset = 0; i < numPoints; i++) {
3236ea844a1aSMatthew Knepley        PetscInt           point = points[2*i], dof, sOffset;
3237ea844a1aSMatthew Knepley        const PetscInt    *perm  = perms ? perms[i] : NULL;
3238ea844a1aSMatthew Knepley        const PetscScalar *rot   = rots  ? rots[i]  : NULL;
3239ea844a1aSMatthew Knepley 
3240ea844a1aSMatthew Knepley        PetscSectionGetDof(section,point,&dof);
3241ea844a1aSMatthew Knepley        PetscSectionGetOffset(section,point,&sOffset);
3242ea844a1aSMatthew Knepley 
3243ea844a1aSMatthew Knepley        if (perm) {for (j = 0; j < dof; j++) {lArray[lOffset + perm[j]]  = sArray[sOffset + j];}}
3244ea844a1aSMatthew Knepley        else      {for (j = 0; j < dof; j++) {lArray[lOffset +      j ]  = sArray[sOffset + j];}}
3245ea844a1aSMatthew Knepley        if (rot)  {for (j = 0; j < dof; j++) {lArray[lOffset +      j ] *= rot[j];             }}
3246ea844a1aSMatthew Knepley        lOffset += dof;
3247ea844a1aSMatthew Knepley      }
3248ea844a1aSMatthew Knepley      PetscSectionRestorePointSyms(section,numPoints,points,&perms,&rots);
3249ea844a1aSMatthew Knepley .ve
3250ea844a1aSMatthew Knepley 
3251ea844a1aSMatthew Knepley   Example of usage, adding dofs into a section array (sArray) from a local array (lArray):
3252ea844a1aSMatthew Knepley .vb
3253ea844a1aSMatthew Knepley      const PetscInt    **perms;
3254ea844a1aSMatthew Knepley      const PetscScalar **rots;
3255ea844a1aSMatthew Knepley      PetscInt            lOffset;
3256ea844a1aSMatthew Knepley 
3257ea844a1aSMatthew Knepley      PetscSectionGetPointSyms(section,numPoints,points,&perms,&rots);
3258ea844a1aSMatthew Knepley      for (i = 0, lOffset = 0; i < numPoints; i++) {
3259ea844a1aSMatthew Knepley        PetscInt           point = points[2*i], dof, sOffset;
3260ea844a1aSMatthew Knepley        const PetscInt    *perm  = perms ? perms[i] : NULL;
3261ea844a1aSMatthew Knepley        const PetscScalar *rot   = rots  ? rots[i]  : NULL;
3262ea844a1aSMatthew Knepley 
3263ea844a1aSMatthew Knepley        PetscSectionGetDof(section,point,&dof);
3264ea844a1aSMatthew Knepley        PetscSectionGetOffset(section,point,&sOff);
3265ea844a1aSMatthew Knepley 
3266ea844a1aSMatthew Knepley        if (perm) {for (j = 0; j < dof; j++) {sArray[sOffset + j] += lArray[lOffset + perm[j]] * (rot ? PetscConj(rot[perm[j]]) : 1.);}}
3267ea844a1aSMatthew Knepley        else      {for (j = 0; j < dof; j++) {sArray[sOffset + j] += lArray[lOffset +      j ] * (rot ? PetscConj(rot[     j ]) : 1.);}}
3268ea844a1aSMatthew Knepley        offset += dof;
3269ea844a1aSMatthew Knepley      }
3270ea844a1aSMatthew Knepley      PetscSectionRestorePointSyms(section,numPoints,points,&perms,&rots);
3271ea844a1aSMatthew Knepley .ve
3272ea844a1aSMatthew Knepley 
3273ea844a1aSMatthew Knepley   Level: developer
3274ea844a1aSMatthew Knepley 
3275cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSym`, `PetscSectionRestorePointSyms()`, `PetscSectionSymCreate()`, `PetscSectionSetSym()`, `PetscSectionGetSym()`
3276ea844a1aSMatthew Knepley @*/
3277d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetPointSyms(PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots)
3278d71ae5a4SJacob Faibussowitsch {
3279ea844a1aSMatthew Knepley   PetscSectionSym sym;
3280ea844a1aSMatthew Knepley 
3281ea844a1aSMatthew Knepley   PetscFunctionBegin;
3282ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1);
3283ea844a1aSMatthew Knepley   if (numPoints) PetscValidIntPointer(points, 3);
3284ea844a1aSMatthew Knepley   if (perms) *perms = NULL;
3285ea844a1aSMatthew Knepley   if (rots) *rots = NULL;
3286ea844a1aSMatthew Knepley   sym = section->sym;
3287ea844a1aSMatthew Knepley   if (sym && (perms || rots)) {
3288ea844a1aSMatthew Knepley     SymWorkLink link;
3289ea844a1aSMatthew Knepley 
3290ea844a1aSMatthew Knepley     if (sym->workin) {
3291ea844a1aSMatthew Knepley       link        = sym->workin;
3292ea844a1aSMatthew Knepley       sym->workin = sym->workin->next;
3293ea844a1aSMatthew Knepley     } else {
32944dfa11a4SJacob Faibussowitsch       PetscCall(PetscNew(&link));
3295ea844a1aSMatthew Knepley     }
3296ea844a1aSMatthew Knepley     if (numPoints > link->numPoints) {
32970c99d500SBarry Smith       PetscInt    **perms = (PetscInt **)link->perms;
32980c99d500SBarry Smith       PetscScalar **rots  = (PetscScalar **)link->rots;
32990c99d500SBarry Smith       PetscCall(PetscFree2(perms, rots));
33000c99d500SBarry Smith       PetscCall(PetscMalloc2(numPoints, (PetscInt ***)&link->perms, numPoints, (PetscScalar ***)&link->rots));
3301ea844a1aSMatthew Knepley       link->numPoints = numPoints;
3302ea844a1aSMatthew Knepley     }
3303ea844a1aSMatthew Knepley     link->next   = sym->workout;
3304ea844a1aSMatthew Knepley     sym->workout = link;
33059566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero((PetscInt **)link->perms, numPoints));
33069566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero((PetscInt **)link->rots, numPoints));
33079566063dSJacob Faibussowitsch     PetscCall((*sym->ops->getpoints)(sym, section, numPoints, points, link->perms, link->rots));
3308ea844a1aSMatthew Knepley     if (perms) *perms = link->perms;
3309ea844a1aSMatthew Knepley     if (rots) *rots = link->rots;
3310ea844a1aSMatthew Knepley   }
3311ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3312ea844a1aSMatthew Knepley }
3313ea844a1aSMatthew Knepley 
3314ea844a1aSMatthew Knepley /*@C
3315cab54364SBarry Smith   PetscSectionRestorePointSyms - Restore the symmetries returned by `PetscSectionGetPointSyms()`
3316ea844a1aSMatthew Knepley 
331740750d41SVaclav Hapla   Not Collective
3318ea844a1aSMatthew Knepley 
3319ea844a1aSMatthew Knepley   Input Parameters:
3320ea844a1aSMatthew Knepley + section - the section
3321ea844a1aSMatthew Knepley . numPoints - the number of points
3322ea844a1aSMatthew Knepley - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an
3323cab54364SBarry Smith     arbitrary integer: its interpretation is up to sym.  Orientations are used by `DM`: for their interpretation in that
3324cab54364SBarry Smith     context, see `DMPlexGetConeOrientation()`).
3325ea844a1aSMatthew Knepley 
3326d8d19677SJose E. Roman   Output Parameters:
3327ea844a1aSMatthew Knepley + perms - The permutations for the given orientations: set to NULL at conclusion
3328ea844a1aSMatthew Knepley - rots - The field rotations symmetries for the given orientations: set to NULL at conclusion
3329ea844a1aSMatthew Knepley 
3330ea844a1aSMatthew Knepley   Level: developer
3331ea844a1aSMatthew Knepley 
3332cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSym`, `PetscSectionGetPointSyms()`, `PetscSectionSymCreate()`, `PetscSectionSetSym()`, `PetscSectionGetSym()`
3333ea844a1aSMatthew Knepley @*/
3334d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionRestorePointSyms(PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots)
3335d71ae5a4SJacob Faibussowitsch {
3336ea844a1aSMatthew Knepley   PetscSectionSym sym;
3337ea844a1aSMatthew Knepley 
3338ea844a1aSMatthew Knepley   PetscFunctionBegin;
3339ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1);
3340ea844a1aSMatthew Knepley   sym = section->sym;
3341ea844a1aSMatthew Knepley   if (sym && (perms || rots)) {
3342ea844a1aSMatthew Knepley     SymWorkLink *p, link;
3343ea844a1aSMatthew Knepley 
3344ea844a1aSMatthew Knepley     for (p = &sym->workout; (link = *p); p = &link->next) {
3345ea844a1aSMatthew Knepley       if ((perms && link->perms == *perms) || (rots && link->rots == *rots)) {
3346ea844a1aSMatthew Knepley         *p          = link->next;
3347ea844a1aSMatthew Knepley         link->next  = sym->workin;
3348ea844a1aSMatthew Knepley         sym->workin = link;
3349ea844a1aSMatthew Knepley         if (perms) *perms = NULL;
3350ea844a1aSMatthew Knepley         if (rots) *rots = NULL;
3351ea844a1aSMatthew Knepley         PetscFunctionReturn(0);
3352ea844a1aSMatthew Knepley       }
3353ea844a1aSMatthew Knepley     }
3354ea844a1aSMatthew Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Array was not checked out");
3355ea844a1aSMatthew Knepley   }
3356ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3357ea844a1aSMatthew Knepley }
3358ea844a1aSMatthew Knepley 
3359ea844a1aSMatthew Knepley /*@C
3360cab54364SBarry Smith   PetscSectionGetFieldPointSyms - Get the symmetries for a set of points in a field of a `PetscSection` under specific orientations.
3361ea844a1aSMatthew Knepley 
336240750d41SVaclav Hapla   Not Collective
3363ea844a1aSMatthew Knepley 
3364ea844a1aSMatthew Knepley   Input Parameters:
3365ea844a1aSMatthew Knepley + section - the section
3366ea844a1aSMatthew Knepley . field - the field of the section
3367ea844a1aSMatthew Knepley . numPoints - the number of points
3368ea844a1aSMatthew Knepley - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an
3369cab54364SBarry Smith     arbitrary integer: its interpretation is up to sym.  Orientations are used by `DM`: for their interpretation in that
3370cab54364SBarry Smith     context, see `DMPlexGetConeOrientation()`).
3371ea844a1aSMatthew Knepley 
3372d8d19677SJose E. Roman   Output Parameters:
3373ea844a1aSMatthew Knepley + perms - The permutations for the given orientations (or NULL if there is no symmetry or the permutation is the identity).
3374ea844a1aSMatthew Knepley - rots - The field rotations symmetries for the given orientations (or NULL if there is no symmetry or the rotations are all
3375ea844a1aSMatthew Knepley     identity).
3376ea844a1aSMatthew Knepley 
3377ea844a1aSMatthew Knepley   Level: developer
3378ea844a1aSMatthew Knepley 
3379cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSym`, `PetscSectionGetPointSyms()`, `PetscSectionRestoreFieldPointSyms()`
3380ea844a1aSMatthew Knepley @*/
3381d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetFieldPointSyms(PetscSection section, PetscInt field, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots)
3382d71ae5a4SJacob Faibussowitsch {
3383ea844a1aSMatthew Knepley   PetscFunctionBegin;
3384ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1);
338508401ef6SPierre Jolivet   PetscCheck(field <= section->numFields, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "field %" PetscInt_FMT " greater than number of fields (%" PetscInt_FMT ") in section", field, section->numFields);
33869566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetPointSyms(section->field[field], numPoints, points, perms, rots));
3387ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3388ea844a1aSMatthew Knepley }
3389ea844a1aSMatthew Knepley 
3390ea844a1aSMatthew Knepley /*@C
3391cab54364SBarry Smith   PetscSectionRestoreFieldPointSyms - Restore the symmetries returned by `PetscSectionGetFieldPointSyms()`
3392ea844a1aSMatthew Knepley 
339340750d41SVaclav Hapla   Not Collective
3394ea844a1aSMatthew Knepley 
3395ea844a1aSMatthew Knepley   Input Parameters:
3396ea844a1aSMatthew Knepley + section - the section
3397ea844a1aSMatthew Knepley . field - the field number
3398ea844a1aSMatthew Knepley . numPoints - the number of points
3399ea844a1aSMatthew Knepley - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an
3400cab54364SBarry Smith     arbitrary integer: its interpretation is up to sym.  Orientations are used by `DM`: for their interpretation in that
3401cab54364SBarry Smith     context, see `DMPlexGetConeOrientation()`).
3402ea844a1aSMatthew Knepley 
3403d8d19677SJose E. Roman   Output Parameters:
3404ea844a1aSMatthew Knepley + perms - The permutations for the given orientations: set to NULL at conclusion
3405ea844a1aSMatthew Knepley - rots - The field rotations symmetries for the given orientations: set to NULL at conclusion
3406ea844a1aSMatthew Knepley 
3407ea844a1aSMatthew Knepley   Level: developer
3408ea844a1aSMatthew Knepley 
3409cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSym`, `PetscSectionRestorePointSyms()`, `petscSectionGetFieldPointSyms()`, `PetscSectionSymCreate()`, `PetscSectionSetSym()`, `PetscSectionGetSym()`
3410ea844a1aSMatthew Knepley @*/
3411d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionRestoreFieldPointSyms(PetscSection section, PetscInt field, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots)
3412d71ae5a4SJacob Faibussowitsch {
3413ea844a1aSMatthew Knepley   PetscFunctionBegin;
3414ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1);
341508401ef6SPierre Jolivet   PetscCheck(field <= section->numFields, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "field %" PetscInt_FMT " greater than number of fields (%" PetscInt_FMT ") in section", field, section->numFields);
34169566063dSJacob Faibussowitsch   PetscCall(PetscSectionRestorePointSyms(section->field[field], numPoints, points, perms, rots));
3417ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3418ea844a1aSMatthew Knepley }
3419ea844a1aSMatthew Knepley 
3420ea844a1aSMatthew Knepley /*@
3421b004864fSMatthew G. Knepley   PetscSectionSymCopy - Copy the symmetries, assuming that the point structure is compatible
3422b004864fSMatthew G. Knepley 
342340750d41SVaclav Hapla   Not Collective
3424b004864fSMatthew G. Knepley 
3425b004864fSMatthew G. Knepley   Input Parameter:
3426cab54364SBarry Smith . sym - the `PetscSectionSym`
3427b004864fSMatthew G. Knepley 
3428b004864fSMatthew G. Knepley   Output Parameter:
3429b004864fSMatthew G. Knepley . nsym - the equivalent symmetries
3430b004864fSMatthew G. Knepley 
3431b004864fSMatthew G. Knepley   Level: developer
3432b004864fSMatthew G. Knepley 
3433cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSym`, `PetscSectionSymCreate()`, `PetscSectionSetSym()`, `PetscSectionGetSym()`, `PetscSectionSymLabelSetStratum()`, `PetscSectionGetPointSyms()`
3434b004864fSMatthew G. Knepley @*/
3435d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSymCopy(PetscSectionSym sym, PetscSectionSym nsym)
3436d71ae5a4SJacob Faibussowitsch {
3437b004864fSMatthew G. Knepley   PetscFunctionBegin;
3438b004864fSMatthew G. Knepley   PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1);
3439b004864fSMatthew G. Knepley   PetscValidHeaderSpecific(nsym, PETSC_SECTION_SYM_CLASSID, 2);
3440dbbe0bcdSBarry Smith   PetscTryTypeMethod(sym, copy, nsym);
3441b004864fSMatthew G. Knepley   PetscFunctionReturn(0);
3442b004864fSMatthew G. Knepley }
3443b004864fSMatthew G. Knepley 
3444b004864fSMatthew G. Knepley /*@
3445cab54364SBarry Smith   PetscSectionSymDistribute - Distribute the symmetries in accordance with the input `PetscSF`
3446b004864fSMatthew G. Knepley 
3447b004864fSMatthew G. Knepley   Collective
3448b004864fSMatthew G. Knepley 
3449b004864fSMatthew G. Knepley   Input Parameters:
3450cab54364SBarry Smith + sym - the `PetscSectionSym`
3451b004864fSMatthew G. Knepley - migrationSF - the distribution map from roots to leaves
3452b004864fSMatthew G. Knepley 
3453b004864fSMatthew G. Knepley   Output Parameters:
3454b004864fSMatthew G. Knepley . dsym - the redistributed symmetries
3455b004864fSMatthew G. Knepley 
3456b004864fSMatthew G. Knepley   Level: developer
3457b004864fSMatthew G. Knepley 
3458cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSym`, `PetscSectionSymCreate()`, `PetscSectionSetSym()`, `PetscSectionGetSym()`, `PetscSectionSymLabelSetStratum()`, `PetscSectionGetPointSyms()`
3459b004864fSMatthew G. Knepley @*/
3460d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSymDistribute(PetscSectionSym sym, PetscSF migrationSF, PetscSectionSym *dsym)
3461d71ae5a4SJacob Faibussowitsch {
3462b004864fSMatthew G. Knepley   PetscFunctionBegin;
3463b004864fSMatthew G. Knepley   PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1);
3464b004864fSMatthew G. Knepley   PetscValidHeaderSpecific(migrationSF, PETSCSF_CLASSID, 2);
3465b004864fSMatthew G. Knepley   PetscValidPointer(dsym, 3);
3466dbbe0bcdSBarry Smith   PetscTryTypeMethod(sym, distribute, migrationSF, dsym);
3467b004864fSMatthew G. Knepley   PetscFunctionReturn(0);
3468b004864fSMatthew G. Knepley }
3469b004864fSMatthew G. Knepley 
3470b004864fSMatthew G. Knepley /*@
3471ea844a1aSMatthew Knepley   PetscSectionGetUseFieldOffsets - Get the flag to use field offsets directly in a global section, rather than just the point offset
3472ea844a1aSMatthew Knepley 
347340750d41SVaclav Hapla   Not Collective
3474ea844a1aSMatthew Knepley 
3475ea844a1aSMatthew Knepley   Input Parameter:
3476cab54364SBarry Smith . s - the global `PetscSection`
3477ea844a1aSMatthew Knepley 
3478ea844a1aSMatthew Knepley   Output Parameters:
3479ea844a1aSMatthew Knepley . flg - the flag
3480ea844a1aSMatthew Knepley 
3481ea844a1aSMatthew Knepley   Level: developer
3482ea844a1aSMatthew Knepley 
3483cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSym`, `PetscSectionSetChart()`, `PetscSectionCreate()`
3484ea844a1aSMatthew Knepley @*/
3485d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionGetUseFieldOffsets(PetscSection s, PetscBool *flg)
3486d71ae5a4SJacob Faibussowitsch {
3487ea844a1aSMatthew Knepley   PetscFunctionBegin;
3488ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
3489ea844a1aSMatthew Knepley   *flg = s->useFieldOff;
3490ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3491ea844a1aSMatthew Knepley }
3492ea844a1aSMatthew Knepley 
3493ea844a1aSMatthew Knepley /*@
3494ea844a1aSMatthew Knepley   PetscSectionSetUseFieldOffsets - Set the flag to use field offsets directly in a global section, rather than just the point offset
3495ea844a1aSMatthew Knepley 
349640750d41SVaclav Hapla   Not Collective
3497ea844a1aSMatthew Knepley 
3498ea844a1aSMatthew Knepley   Input Parameters:
3499cab54364SBarry Smith + s   - the global `PetscSection`
3500ea844a1aSMatthew Knepley - flg - the flag
3501ea844a1aSMatthew Knepley 
3502ea844a1aSMatthew Knepley   Level: developer
3503ea844a1aSMatthew Knepley 
3504cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSym`, `PetscSectionGetUseFieldOffsets()`, `PetscSectionSetChart()`, `PetscSectionCreate()`
3505ea844a1aSMatthew Knepley @*/
3506d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionSetUseFieldOffsets(PetscSection s, PetscBool flg)
3507d71ae5a4SJacob Faibussowitsch {
3508ea844a1aSMatthew Knepley   PetscFunctionBegin;
3509ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1);
3510ea844a1aSMatthew Knepley   s->useFieldOff = flg;
3511ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3512ea844a1aSMatthew Knepley }
3513ea844a1aSMatthew Knepley 
3514ea844a1aSMatthew Knepley #define PetscSectionExpandPoints_Loop(TYPE) \
3515ea844a1aSMatthew Knepley   { \
3516ea844a1aSMatthew Knepley     PetscInt i, n, o0, o1, size; \
3517ea844a1aSMatthew Knepley     TYPE    *a0 = (TYPE *)origArray, *a1; \
35189566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetStorageSize(s, &size)); \
35199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &a1)); \
3520ea844a1aSMatthew Knepley     for (i = 0; i < npoints; i++) { \
35219566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(origSection, points_[i], &o0)); \
35229566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(s, i, &o1)); \
35239566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(s, i, &n)); \
35249566063dSJacob Faibussowitsch       PetscCall(PetscMemcpy(&a1[o1], &a0[o0], n *unitsize)); \
3525ea844a1aSMatthew Knepley     } \
3526ea844a1aSMatthew Knepley     *newArray = (void *)a1; \
3527ea844a1aSMatthew Knepley   }
3528ea844a1aSMatthew Knepley 
3529ea844a1aSMatthew Knepley /*@
3530ea844a1aSMatthew Knepley   PetscSectionExtractDofsFromArray - Extracts elements of an array corresponding to DOFs of specified points.
3531ea844a1aSMatthew Knepley 
353240750d41SVaclav Hapla   Not Collective
3533ea844a1aSMatthew Knepley 
3534ea844a1aSMatthew Knepley   Input Parameters:
3535cab54364SBarry Smith + origSection - the `PetscSection` describing the layout of the array
3536cab54364SBarry Smith . dataType - `MPI_Datatype` describing the data type of the array (currently only `MPIU_INT`, `MPIU_SCALAR`, `MPIU_REAL`)
3537ea844a1aSMatthew Knepley . origArray - the array; its size must be equal to the storage size of origSection
3538cab54364SBarry Smith - points - `IS` with points to extract; its indices must lie in the chart of origSection
3539ea844a1aSMatthew Knepley 
3540ea844a1aSMatthew Knepley   Output Parameters:
3541cab54364SBarry Smith + newSection - the new `PetscSection` desribing the layout of the new array (with points renumbered 0,1,... but preserving numbers of DOFs)
3542ea844a1aSMatthew Knepley - newArray - the array of the extracted DOFs; its size is the storage size of newSection
3543ea844a1aSMatthew Knepley 
3544ea844a1aSMatthew Knepley   Level: developer
3545ea844a1aSMatthew Knepley 
3546cab54364SBarry Smith .seealso: [PetscSection](sec_petscsection), `PetscSectionSym`, `PetscSectionGetChart()`, `PetscSectionGetDof()`, `PetscSectionGetStorageSize()`, `PetscSectionCreate()`
3547ea844a1aSMatthew Knepley @*/
3548d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSectionExtractDofsFromArray(PetscSection origSection, MPI_Datatype dataType, const void *origArray, IS points, PetscSection *newSection, void *newArray[])
3549d71ae5a4SJacob Faibussowitsch {
3550ea844a1aSMatthew Knepley   PetscSection    s;
3551ea844a1aSMatthew Knepley   const PetscInt *points_;
3552ea844a1aSMatthew Knepley   PetscInt        i, n, npoints, pStart, pEnd;
3553ea844a1aSMatthew Knepley   PetscMPIInt     unitsize;
3554ea844a1aSMatthew Knepley 
3555ea844a1aSMatthew Knepley   PetscFunctionBegin;
3556ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(origSection, PETSC_SECTION_CLASSID, 1);
3557ea844a1aSMatthew Knepley   PetscValidPointer(origArray, 3);
3558ea844a1aSMatthew Knepley   PetscValidHeaderSpecific(points, IS_CLASSID, 4);
3559ea844a1aSMatthew Knepley   if (newSection) PetscValidPointer(newSection, 5);
3560ea844a1aSMatthew Knepley   if (newArray) PetscValidPointer(newArray, 6);
35619566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_size(dataType, &unitsize));
35629566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(points, &npoints));
35639566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(points, &points_));
35649566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(origSection, &pStart, &pEnd));
35659566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &s));
35669566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(s, 0, npoints));
3567ea844a1aSMatthew Knepley   for (i = 0; i < npoints; i++) {
3568c9cc58a2SBarry Smith     PetscCheck(points_[i] >= pStart && points_[i] < pEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "point %" PetscInt_FMT " (index %" PetscInt_FMT ") in input IS out of input section's chart", points_[i], i);
35699566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(origSection, points_[i], &n));
35709566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(s, i, n));
3571ea844a1aSMatthew Knepley   }
35729566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(s));
3573ea844a1aSMatthew Knepley   if (newArray) {
35749371c9d4SSatish Balay     if (dataType == MPIU_INT) {
35759371c9d4SSatish Balay       PetscSectionExpandPoints_Loop(PetscInt);
35769371c9d4SSatish Balay     } else if (dataType == MPIU_SCALAR) {
35779371c9d4SSatish Balay       PetscSectionExpandPoints_Loop(PetscScalar);
35789371c9d4SSatish Balay     } else if (dataType == MPIU_REAL) {
35799371c9d4SSatish Balay       PetscSectionExpandPoints_Loop(PetscReal);
35809371c9d4SSatish Balay     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "not implemented for this MPI_Datatype");
3581ea844a1aSMatthew Knepley   }
3582ea844a1aSMatthew Knepley   if (newSection) {
3583ea844a1aSMatthew Knepley     *newSection = s;
3584ea844a1aSMatthew Knepley   } else {
35859566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&s));
3586ea844a1aSMatthew Knepley   }
35879566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(points, &points_));
3588ea844a1aSMatthew Knepley   PetscFunctionReturn(0);
3589ea844a1aSMatthew Knepley }
3590