xref: /petsc/src/dm/impls/plex/plexsection.c (revision 09b68a49ed2854d1e4985cc2aa6af33c7c4e69b3)
19cbb8f0dSMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
29cbb8f0dSMatthew G. Knepley 
39cbb8f0dSMatthew G. Knepley /* Set the number of dof on each point and separate by fields */
DMPlexCreateSectionFields(DM dm,const PetscInt numComp[],PetscSection * section)4d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateSectionFields(DM dm, const PetscInt numComp[], PetscSection *section)
5d71ae5a4SJacob Faibussowitsch {
6baef929fSMatthew G. Knepley   DMLabel    depthLabel;
7baef929fSMatthew G. Knepley   PetscInt   depth, Nf, f, pStart, pEnd;
89cbb8f0dSMatthew G. Knepley   PetscBool *isFE;
99cbb8f0dSMatthew G. Knepley 
109cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
119566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
129566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
139566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
149566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(Nf, &isFE));
15baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
169cbb8f0dSMatthew G. Knepley     PetscObject  obj;
179cbb8f0dSMatthew G. Knepley     PetscClassId id;
189cbb8f0dSMatthew G. Knepley 
199566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, f, NULL, &obj));
209566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(obj, &id));
219371c9d4SSatish Balay     if (id == PETSCFE_CLASSID) {
229371c9d4SSatish Balay       isFE[f] = PETSC_TRUE;
239371c9d4SSatish Balay     } else if (id == PETSCFV_CLASSID) {
249371c9d4SSatish Balay       isFE[f] = PETSC_FALSE;
259371c9d4SSatish Balay     }
269cbb8f0dSMatthew G. Knepley   }
27baef929fSMatthew G. Knepley 
289566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), section));
29baef929fSMatthew G. Knepley   if (Nf > 0) {
309566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetNumFields(*section, Nf));
319cbb8f0dSMatthew G. Knepley     if (numComp) {
32baef929fSMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
339566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetFieldComponents(*section, f, numComp[f]));
349cbb8f0dSMatthew G. Knepley         if (isFE[f]) {
359cbb8f0dSMatthew G. Knepley           PetscFE              fe;
369cbb8f0dSMatthew G. Knepley           PetscDualSpace       dspace;
379cbb8f0dSMatthew G. Knepley           const PetscInt    ***perms;
389cbb8f0dSMatthew G. Knepley           const PetscScalar ***flips;
399cbb8f0dSMatthew G. Knepley           const PetscInt      *numDof;
409cbb8f0dSMatthew G. Knepley 
419566063dSJacob Faibussowitsch           PetscCall(DMGetField(dm, f, NULL, (PetscObject *)&fe));
429566063dSJacob Faibussowitsch           PetscCall(PetscFEGetDualSpace(fe, &dspace));
439566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceGetSymmetries(dspace, &perms, &flips));
449566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceGetNumDof(dspace, &numDof));
459cbb8f0dSMatthew G. Knepley           if (perms || flips) {
469cbb8f0dSMatthew G. Knepley             DM              K;
474dff28b8SMatthew G. Knepley             PetscInt        sph, spdepth;
489cbb8f0dSMatthew G. Knepley             PetscSectionSym sym;
499cbb8f0dSMatthew G. Knepley 
509566063dSJacob Faibussowitsch             PetscCall(PetscDualSpaceGetDM(dspace, &K));
519566063dSJacob Faibussowitsch             PetscCall(DMPlexGetDepth(K, &spdepth));
529566063dSJacob Faibussowitsch             PetscCall(PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)*section), depthLabel, &sym));
534dff28b8SMatthew G. Knepley             for (sph = 0; sph <= spdepth; sph++) {
549cbb8f0dSMatthew G. Knepley               PetscDualSpace      hspace;
559cbb8f0dSMatthew G. Knepley               PetscInt            kStart, kEnd;
564dff28b8SMatthew G. Knepley               PetscInt            kConeSize, h = sph + (depth - spdepth);
579cbb8f0dSMatthew G. Knepley               const PetscInt    **perms0 = NULL;
589cbb8f0dSMatthew G. Knepley               const PetscScalar **flips0 = NULL;
599cbb8f0dSMatthew G. Knepley 
609566063dSJacob Faibussowitsch               PetscCall(PetscDualSpaceGetHeightSubspace(dspace, sph, &hspace));
619566063dSJacob Faibussowitsch               PetscCall(DMPlexGetHeightStratum(K, h, &kStart, &kEnd));
629cbb8f0dSMatthew G. Knepley               if (!hspace) continue;
639566063dSJacob Faibussowitsch               PetscCall(PetscDualSpaceGetSymmetries(hspace, &perms, &flips));
649cbb8f0dSMatthew G. Knepley               if (perms) perms0 = perms[0];
659cbb8f0dSMatthew G. Knepley               if (flips) flips0 = flips[0];
669cbb8f0dSMatthew G. Knepley               if (!(perms0 || flips0)) continue;
67b5a892a1SMatthew G. Knepley               {
68b5a892a1SMatthew G. Knepley                 DMPolytopeType ct;
69b5a892a1SMatthew G. Knepley                 /* The number of arrangements is no longer based on the number of faces */
709566063dSJacob Faibussowitsch                 PetscCall(DMPlexGetCellType(K, kStart, &ct));
7185036b15SMatthew G. Knepley                 kConeSize = DMPolytopeTypeGetNumArrangements(ct) / 2;
72b5a892a1SMatthew G. Knepley               }
738e3a54c0SPierre Jolivet               PetscCall(PetscSectionSymLabelSetStratum(sym, depth - h, numDof[depth - h], -kConeSize, kConeSize, PETSC_USE_POINTER, PetscSafePointerPlusOffset(perms0, -kConeSize), PetscSafePointerPlusOffset(flips0, -kConeSize)));
749cbb8f0dSMatthew G. Knepley             }
759566063dSJacob Faibussowitsch             PetscCall(PetscSectionSetFieldSym(*section, f, sym));
769566063dSJacob Faibussowitsch             PetscCall(PetscSectionSymDestroy(&sym));
779cbb8f0dSMatthew G. Knepley           }
789cbb8f0dSMatthew G. Knepley         }
799cbb8f0dSMatthew G. Knepley       }
809cbb8f0dSMatthew G. Knepley     }
819cbb8f0dSMatthew G. Knepley   }
829566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
839566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(*section, pStart, pEnd));
849566063dSJacob Faibussowitsch   PetscCall(PetscFree(isFE));
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86baef929fSMatthew G. Knepley }
87baef929fSMatthew G. Knepley 
88baef929fSMatthew G. Knepley /* Set the number of dof on each point and separate by fields */
DMPlexCreateSectionDof(DM dm,DMLabel label[],const PetscInt numDof[],PetscSection section)89d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateSectionDof(DM dm, DMLabel label[], const PetscInt numDof[], PetscSection section)
90d71ae5a4SJacob Faibussowitsch {
91baef929fSMatthew G. Knepley   DMLabel        depthLabel;
92412e9a14SMatthew G. Knepley   DMPolytopeType ct;
93baef929fSMatthew G. Knepley   PetscInt       depth, cellHeight, pStart = 0, pEnd = 0;
94a55f9a55SMatthew G. Knepley   PetscInt       Nf, f, Nds, n, dim, d, dep, p;
955fedec97SMatthew G. Knepley   PetscBool     *isFE, hasCohesive = PETSC_FALSE;
96baef929fSMatthew G. Knepley 
97baef929fSMatthew G. Knepley   PetscFunctionBegin;
989566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
999566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
1009566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
1019566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
1029566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
103a55f9a55SMatthew G. Knepley   for (n = 0; n < Nds; ++n) {
104a55f9a55SMatthew G. Knepley     PetscDS   ds;
1055fedec97SMatthew G. Knepley     PetscBool isCohesive;
106a55f9a55SMatthew G. Knepley 
10707218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds, NULL));
1089566063dSJacob Faibussowitsch     PetscCall(PetscDSIsCohesive(ds, &isCohesive));
1099371c9d4SSatish Balay     if (isCohesive) {
1109371c9d4SSatish Balay       hasCohesive = PETSC_TRUE;
1119371c9d4SSatish Balay       break;
1129371c9d4SSatish Balay     }
113a55f9a55SMatthew G. Knepley   }
1149566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nf, &isFE));
115baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
116baef929fSMatthew G. Knepley     PetscObject  obj;
117baef929fSMatthew G. Knepley     PetscClassId id;
118baef929fSMatthew G. Knepley 
1199566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, f, NULL, &obj));
1209566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(obj, &id));
1211274e1a1SLawrence Mitchell     /* User is allowed to put a "placeholder" field in (c.f. DMCreateDS) */
1221274e1a1SLawrence Mitchell     isFE[f] = id == PETSCFE_CLASSID ? PETSC_TRUE : PETSC_FALSE;
123baef929fSMatthew G. Knepley   }
124baef929fSMatthew G. Knepley 
1259566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
126baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
127e0b68406SMatthew Knepley     PetscBool avoidTensor;
128e0b68406SMatthew Knepley 
1299566063dSJacob Faibussowitsch     PetscCall(DMGetFieldAvoidTensor(dm, f, &avoidTensor));
1305fedec97SMatthew G. Knepley     avoidTensor = (avoidTensor || hasCohesive) ? PETSC_TRUE : PETSC_FALSE;
131baef929fSMatthew G. Knepley     if (label && label[f]) {
132baef929fSMatthew G. Knepley       IS              pointIS;
133baef929fSMatthew G. Knepley       const PetscInt *points;
134baef929fSMatthew G. Knepley       PetscInt        n;
135baef929fSMatthew G. Knepley 
1369566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(label[f], 1, &pointIS));
137baef929fSMatthew G. Knepley       if (!pointIS) continue;
1389566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(pointIS, &n));
1399566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(pointIS, &points));
140baef929fSMatthew G. Knepley       for (p = 0; p < n; ++p) {
141baef929fSMatthew G. Knepley         const PetscInt point = points[p];
142baef929fSMatthew G. Knepley         PetscInt       dof, d;
143baef929fSMatthew G. Knepley 
1449566063dSJacob Faibussowitsch         PetscCall(DMPlexGetCellType(dm, point, &ct));
1459566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValue(depthLabel, point, &d));
146412e9a14SMatthew G. Knepley         /* If this is a tensor prism point, use dof for one dimension lower */
147412e9a14SMatthew G. Knepley         switch (ct) {
148412e9a14SMatthew G. Knepley         case DM_POLYTOPE_POINT_PRISM_TENSOR:
149412e9a14SMatthew G. Knepley         case DM_POLYTOPE_SEG_PRISM_TENSOR:
150412e9a14SMatthew G. Knepley         case DM_POLYTOPE_TRI_PRISM_TENSOR:
151412e9a14SMatthew G. Knepley         case DM_POLYTOPE_QUAD_PRISM_TENSOR:
152ad540459SPierre Jolivet           if (hasCohesive) --d;
1539371c9d4SSatish Balay           break;
154d71ae5a4SJacob Faibussowitsch         default:
155d71ae5a4SJacob Faibussowitsch           break;
156412e9a14SMatthew G. Knepley         }
157baef929fSMatthew G. Knepley         dof = d < 0 ? 0 : numDof[f * (dim + 1) + d];
1589566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetFieldDof(section, point, f, dof));
1599566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(section, point, dof));
160baef929fSMatthew G. Knepley       }
1619566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(pointIS, &points));
1629566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
163baef929fSMatthew G. Knepley     } else {
1649cbb8f0dSMatthew G. Knepley       for (dep = 0; dep <= depth - cellHeight; ++dep) {
1653fd2e703SMatthew G. Knepley         /* Cases: dim > depth (cell-vertex mesh), dim == depth (fully interpolated), dim < depth (interpolated submesh) */
1663fd2e703SMatthew G. Knepley         d = dim <= depth ? dep : (!dep ? 0 : dim);
1679566063dSJacob Faibussowitsch         PetscCall(DMPlexGetDepthStratum(dm, dep, &pStart, &pEnd));
1689cbb8f0dSMatthew G. Knepley         for (p = pStart; p < pEnd; ++p) {
169baef929fSMatthew G. Knepley           const PetscInt dof = numDof[f * (dim + 1) + d];
1709cbb8f0dSMatthew G. Knepley 
1719566063dSJacob Faibussowitsch           PetscCall(DMPlexGetCellType(dm, p, &ct));
172412e9a14SMatthew G. Knepley           switch (ct) {
173412e9a14SMatthew G. Knepley           case DM_POLYTOPE_POINT_PRISM_TENSOR:
174412e9a14SMatthew G. Knepley           case DM_POLYTOPE_SEG_PRISM_TENSOR:
175412e9a14SMatthew G. Knepley           case DM_POLYTOPE_TRI_PRISM_TENSOR:
176412e9a14SMatthew G. Knepley           case DM_POLYTOPE_QUAD_PRISM_TENSOR:
177e0b68406SMatthew Knepley             if (avoidTensor && isFE[f]) continue;
178d71ae5a4SJacob Faibussowitsch           default:
179d71ae5a4SJacob Faibussowitsch             break;
180412e9a14SMatthew G. Knepley           }
1819566063dSJacob Faibussowitsch           PetscCall(PetscSectionSetFieldDof(section, p, f, dof));
1829566063dSJacob Faibussowitsch           PetscCall(PetscSectionAddDof(section, p, dof));
1839cbb8f0dSMatthew G. Knepley         }
184baef929fSMatthew G. Knepley       }
1859cbb8f0dSMatthew G. Knepley     }
1869cbb8f0dSMatthew G. Knepley   }
1879566063dSJacob Faibussowitsch   PetscCall(PetscFree(isFE));
1883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1899cbb8f0dSMatthew G. Knepley }
1909cbb8f0dSMatthew G. Knepley 
1919cbb8f0dSMatthew G. Knepley /* Set the number of dof on each point and separate by fields
1929cbb8f0dSMatthew G. Knepley    If bcComps is NULL or the IS is NULL, constrain every dof on the point
1939cbb8f0dSMatthew G. Knepley */
DMPlexCreateSectionBCDof(DM dm,PetscInt numBC,const PetscInt bcField[],const IS bcComps[],const IS bcPoints[],PetscSection section)194d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateSectionBCDof(DM dm, PetscInt numBC, const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section)
195d71ae5a4SJacob Faibussowitsch {
196baef929fSMatthew G. Knepley   PetscInt     Nf;
1979cbb8f0dSMatthew G. Knepley   PetscInt     bc;
1989cbb8f0dSMatthew G. Knepley   PetscSection aSec;
1999cbb8f0dSMatthew G. Knepley 
2009cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
2019566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &Nf));
2029cbb8f0dSMatthew G. Knepley   for (bc = 0; bc < numBC; ++bc) {
2039cbb8f0dSMatthew G. Knepley     PetscInt        field = 0;
2049cbb8f0dSMatthew G. Knepley     const PetscInt *comp;
2059cbb8f0dSMatthew G. Knepley     const PetscInt *idx;
206ce093827SMatthew G. Knepley     PetscInt        Nc = 0, cNc = -1, n, i;
2079cbb8f0dSMatthew G. Knepley 
208ce093827SMatthew G. Knepley     if (Nf) {
209ce093827SMatthew G. Knepley       field = bcField[bc];
2109566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldComponents(section, field, &Nc));
211ce093827SMatthew G. Knepley     }
2129566063dSJacob Faibussowitsch     if (bcComps && bcComps[bc]) PetscCall(ISGetLocalSize(bcComps[bc], &cNc));
2139566063dSJacob Faibussowitsch     if (bcComps && bcComps[bc]) PetscCall(ISGetIndices(bcComps[bc], &comp));
21424cfc5f6SToby Isaac     if (bcPoints[bc]) {
2159566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(bcPoints[bc], &n));
2169566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(bcPoints[bc], &idx));
2179cbb8f0dSMatthew G. Knepley       for (i = 0; i < n; ++i) {
2189cbb8f0dSMatthew G. Knepley         const PetscInt p = idx[i];
2199cbb8f0dSMatthew G. Knepley         PetscInt       numConst;
2209cbb8f0dSMatthew G. Knepley 
221baef929fSMatthew G. Knepley         if (Nf) {
2229566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldDof(section, p, field, &numConst));
2239cbb8f0dSMatthew G. Knepley         } else {
2249566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetDof(section, p, &numConst));
2259cbb8f0dSMatthew G. Knepley         }
226ce093827SMatthew G. Knepley         /* If Nc <= 0, constrain every dof on the point */
227ce093827SMatthew G. Knepley         if (cNc > 0) {
228ce093827SMatthew G. Knepley           /* We assume that a point may have multiple "nodes", which are collections of Nc dofs,
229ce093827SMatthew G. Knepley              and that those dofs are numbered n*Nc+c */
230ce093827SMatthew G. Knepley           if (Nf) {
23163a3b9bcSJacob Faibussowitsch             PetscCheck(numConst % Nc == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT " has %" PetscInt_FMT " dof which is not divisible by %" PetscInt_FMT " field components", p, numConst, Nc);
232ce093827SMatthew G. Knepley             numConst = (numConst / Nc) * cNc;
233ce093827SMatthew G. Knepley           } else {
234ce093827SMatthew G. Knepley             numConst = PetscMin(numConst, cNc);
235ce093827SMatthew G. Knepley           }
236ce093827SMatthew G. Knepley         }
2379566063dSJacob Faibussowitsch         if (Nf) PetscCall(PetscSectionAddFieldConstraintDof(section, p, field, numConst));
2389566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddConstraintDof(section, p, numConst));
2399cbb8f0dSMatthew G. Knepley       }
2409566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(bcPoints[bc], &idx));
24124cfc5f6SToby Isaac     }
2429566063dSJacob Faibussowitsch     if (bcComps && bcComps[bc]) PetscCall(ISRestoreIndices(bcComps[bc], &comp));
2439cbb8f0dSMatthew G. Knepley   }
2449566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAnchors(dm, &aSec, NULL));
2459cbb8f0dSMatthew G. Knepley   if (aSec) {
2469cbb8f0dSMatthew G. Knepley     PetscInt aStart, aEnd, a;
2479cbb8f0dSMatthew G. Knepley 
2489566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd));
2499cbb8f0dSMatthew G. Knepley     for (a = aStart; a < aEnd; a++) {
2509cbb8f0dSMatthew G. Knepley       PetscInt dof, f;
2519cbb8f0dSMatthew G. Knepley 
2529566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(aSec, a, &dof));
2539cbb8f0dSMatthew G. Knepley       if (dof) {
2549cbb8f0dSMatthew G. Knepley         /* if there are point-to-point constraints, then all dofs are constrained */
2559566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, a, &dof));
2569566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetConstraintDof(section, a, dof));
257baef929fSMatthew G. Knepley         for (f = 0; f < Nf; f++) {
2589566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldDof(section, a, f, &dof));
2599566063dSJacob Faibussowitsch           PetscCall(PetscSectionSetFieldConstraintDof(section, a, f, dof));
2609cbb8f0dSMatthew G. Knepley         }
2619cbb8f0dSMatthew G. Knepley       }
2629cbb8f0dSMatthew G. Knepley     }
2639cbb8f0dSMatthew G. Knepley   }
2643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2659cbb8f0dSMatthew G. Knepley }
2669cbb8f0dSMatthew G. Knepley 
2679cbb8f0dSMatthew G. Knepley /* Set the constrained field indices on each point
2689cbb8f0dSMatthew G. Knepley    If bcComps is NULL or the IS is NULL, constrain every dof on the point
2699cbb8f0dSMatthew G. Knepley */
DMPlexCreateSectionBCIndicesField(DM dm,PetscInt numBC,const PetscInt bcField[],const IS bcComps[],const IS bcPoints[],PetscSection section)270d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateSectionBCIndicesField(DM dm, PetscInt numBC, const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section)
271d71ae5a4SJacob Faibussowitsch {
2729cbb8f0dSMatthew G. Knepley   PetscSection aSec;
2739cbb8f0dSMatthew G. Knepley   PetscInt    *indices;
274baef929fSMatthew G. Knepley   PetscInt     Nf, cdof, maxDof = 0, pStart, pEnd, p, bc, f, d;
2759cbb8f0dSMatthew G. Knepley 
2769cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
2779566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &Nf));
2783ba16761SJacob Faibussowitsch   if (!Nf) PetscFunctionReturn(PETSC_SUCCESS);
2799cbb8f0dSMatthew G. Knepley   /* Initialize all field indices to -1 */
2809566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
2819371c9d4SSatish Balay   for (p = pStart; p < pEnd; ++p) {
2829371c9d4SSatish Balay     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
2839371c9d4SSatish Balay     maxDof = PetscMax(maxDof, cdof);
2849371c9d4SSatish Balay   }
2859566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxDof, &indices));
2869cbb8f0dSMatthew G. Knepley   for (d = 0; d < maxDof; ++d) indices[d] = -1;
2879371c9d4SSatish Balay   for (p = pStart; p < pEnd; ++p)
2889371c9d4SSatish Balay     for (f = 0; f < Nf; ++f) PetscCall(PetscSectionSetFieldConstraintIndices(section, p, f, indices));
2899cbb8f0dSMatthew G. Knepley   /* Handle BC constraints */
2909cbb8f0dSMatthew G. Knepley   for (bc = 0; bc < numBC; ++bc) {
2919cbb8f0dSMatthew G. Knepley     const PetscInt  field = bcField[bc];
2929cbb8f0dSMatthew G. Knepley     const PetscInt *comp, *idx;
293ce093827SMatthew G. Knepley     PetscInt        Nc, cNc = -1, n, i;
2949cbb8f0dSMatthew G. Knepley 
2959566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldComponents(section, field, &Nc));
2969566063dSJacob Faibussowitsch     if (bcComps && bcComps[bc]) PetscCall(ISGetLocalSize(bcComps[bc], &cNc));
2979566063dSJacob Faibussowitsch     if (bcComps && bcComps[bc]) PetscCall(ISGetIndices(bcComps[bc], &comp));
29824cfc5f6SToby Isaac     if (bcPoints[bc]) {
2999566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(bcPoints[bc], &n));
3009566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(bcPoints[bc], &idx));
3019cbb8f0dSMatthew G. Knepley       for (i = 0; i < n; ++i) {
3029cbb8f0dSMatthew G. Knepley         const PetscInt  p = idx[i];
3039cbb8f0dSMatthew G. Knepley         const PetscInt *find;
304ce093827SMatthew G. Knepley         PetscInt        fdof, fcdof, c, j;
3059cbb8f0dSMatthew G. Knepley 
3069566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetFieldDof(section, p, field, &fdof));
3079cbb8f0dSMatthew G. Knepley         if (!fdof) continue;
308ce093827SMatthew G. Knepley         if (cNc < 0) {
3099cbb8f0dSMatthew G. Knepley           for (d = 0; d < fdof; ++d) indices[d] = d;
3109cbb8f0dSMatthew G. Knepley           fcdof = fdof;
3119cbb8f0dSMatthew G. Knepley         } else {
312ce093827SMatthew G. Knepley           /* We assume that a point may have multiple "nodes", which are collections of Nc dofs,
313ce093827SMatthew G. Knepley              and that those dofs are numbered n*Nc+c */
3149566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintDof(section, p, field, &fcdof));
3159566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintIndices(section, p, field, &find));
316ce093827SMatthew G. Knepley           /* Get indices constrained by previous bcs */
3179371c9d4SSatish Balay           for (d = 0; d < fcdof; ++d) {
3189371c9d4SSatish Balay             if (find[d] < 0) break;
3199371c9d4SSatish Balay             indices[d] = find[d];
3209371c9d4SSatish Balay           }
3219371c9d4SSatish Balay           for (j = 0; j < fdof / Nc; ++j)
3229371c9d4SSatish Balay             for (c = 0; c < cNc; ++c) indices[d++] = j * Nc + comp[c];
3239566063dSJacob Faibussowitsch           PetscCall(PetscSortRemoveDupsInt(&d, indices));
3249cbb8f0dSMatthew G. Knepley           for (c = d; c < fcdof; ++c) indices[c] = -1;
3259cbb8f0dSMatthew G. Knepley           fcdof = d;
3269cbb8f0dSMatthew G. Knepley         }
3279566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetFieldConstraintDof(section, p, field, fcdof));
3289566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetFieldConstraintIndices(section, p, field, indices));
3299cbb8f0dSMatthew G. Knepley       }
3309566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(bcPoints[bc], &idx));
3319cbb8f0dSMatthew G. Knepley     }
33224cfc5f6SToby Isaac     if (bcComps && bcComps[bc]) PetscCall(ISRestoreIndices(bcComps[bc], &comp));
33324cfc5f6SToby Isaac   }
3349cbb8f0dSMatthew G. Knepley   /* Handle anchors */
3359566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAnchors(dm, &aSec, NULL));
3369cbb8f0dSMatthew G. Knepley   if (aSec) {
3379cbb8f0dSMatthew G. Knepley     PetscInt aStart, aEnd, a;
3389cbb8f0dSMatthew G. Knepley 
3399cbb8f0dSMatthew G. Knepley     for (d = 0; d < maxDof; ++d) indices[d] = d;
3409566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd));
3419cbb8f0dSMatthew G. Knepley     for (a = aStart; a < aEnd; a++) {
3429cbb8f0dSMatthew G. Knepley       PetscInt dof, f;
3439cbb8f0dSMatthew G. Knepley 
3449566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(aSec, a, &dof));
3459cbb8f0dSMatthew G. Knepley       if (dof) {
3469cbb8f0dSMatthew G. Knepley         /* if there are point-to-point constraints, then all dofs are constrained */
34748a46eb9SPierre Jolivet         for (f = 0; f < Nf; f++) PetscCall(PetscSectionSetFieldConstraintIndices(section, a, f, indices));
3489cbb8f0dSMatthew G. Knepley       }
3499cbb8f0dSMatthew G. Knepley     }
3509cbb8f0dSMatthew G. Knepley   }
3519566063dSJacob Faibussowitsch   PetscCall(PetscFree(indices));
3523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3539cbb8f0dSMatthew G. Knepley }
3549cbb8f0dSMatthew G. Knepley 
3559cbb8f0dSMatthew G. Knepley /* Set the constrained indices on each point */
DMPlexCreateSectionBCIndices(DM dm,PetscSection section)356d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateSectionBCIndices(DM dm, PetscSection section)
357d71ae5a4SJacob Faibussowitsch {
3589cbb8f0dSMatthew G. Knepley   PetscInt *indices;
359baef929fSMatthew G. Knepley   PetscInt  Nf, maxDof, pStart, pEnd, p, f, d;
3609cbb8f0dSMatthew G. Knepley 
3619cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
3629566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &Nf));
3639566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetMaxDof(section, &maxDof));
3649566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
3659566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxDof, &indices));
3669cbb8f0dSMatthew G. Knepley   for (d = 0; d < maxDof; ++d) indices[d] = -1;
3679cbb8f0dSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
3689cbb8f0dSMatthew G. Knepley     PetscInt cdof, d;
3699cbb8f0dSMatthew G. Knepley 
3709566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
3719cbb8f0dSMatthew G. Knepley     if (cdof) {
372baef929fSMatthew G. Knepley       if (Nf) {
3739cbb8f0dSMatthew G. Knepley         PetscInt numConst = 0, foff = 0;
3749cbb8f0dSMatthew G. Knepley 
375baef929fSMatthew G. Knepley         for (f = 0; f < Nf; ++f) {
3769cbb8f0dSMatthew G. Knepley           const PetscInt *find;
3779cbb8f0dSMatthew G. Knepley           PetscInt        fcdof, fdof;
3789cbb8f0dSMatthew G. Knepley 
3799566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldDof(section, p, f, &fdof));
3809566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintDof(section, p, f, &fcdof));
3819cbb8f0dSMatthew G. Knepley           /* Change constraint numbering from field component to local dof number */
3829566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintIndices(section, p, f, &find));
3839cbb8f0dSMatthew G. Knepley           for (d = 0; d < fcdof; ++d) indices[numConst + d] = find[d] + foff;
3849cbb8f0dSMatthew G. Knepley           numConst += fcdof;
3859cbb8f0dSMatthew G. Knepley           foff += fdof;
3869cbb8f0dSMatthew G. Knepley         }
3879566063dSJacob Faibussowitsch         if (cdof != numConst) PetscCall(PetscSectionSetConstraintDof(section, p, numConst));
3889cbb8f0dSMatthew G. Knepley       } else {
3899cbb8f0dSMatthew G. Knepley         for (d = 0; d < cdof; ++d) indices[d] = d;
3909cbb8f0dSMatthew G. Knepley       }
3919566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetConstraintIndices(section, p, indices));
3929cbb8f0dSMatthew G. Knepley     }
3939cbb8f0dSMatthew G. Knepley   }
3949566063dSJacob Faibussowitsch   PetscCall(PetscFree(indices));
3953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3969cbb8f0dSMatthew G. Knepley }
3979cbb8f0dSMatthew G. Knepley 
3989cbb8f0dSMatthew G. Knepley /*@C
399a1cb98faSBarry Smith   DMPlexCreateSection - Create a `PetscSection` based upon the dof layout specification provided.
4009cbb8f0dSMatthew G. Knepley 
4019cbb8f0dSMatthew G. Knepley   Not Collective
4029cbb8f0dSMatthew G. Knepley 
4039cbb8f0dSMatthew G. Knepley   Input Parameters:
404a1cb98faSBarry Smith + dm       - The `DMPLEX` object
405*bfe80ac4SPierre Jolivet . label    - An array of `DMLabel` of length `numFields` indicating the mesh support of each field, or `NULL` for the whole mesh
406ce78bad3SBarry Smith . numComp  - An array of size `numFields` that holds the number of components for each field
407ce78bad3SBarry Smith . numDof   - An array of size $ numFields \time (dim+1)$ which holds the number of dof for each field on a mesh piece of dimension d
4089cbb8f0dSMatthew G. Knepley . numBC    - The number of boundary conditions
409ce78bad3SBarry Smith . bcField  - An array of size `numBC` giving the field number for each boundary condition
410ce78bad3SBarry Smith . bcComps  - [Optional] An array of size `numBC` of `IS` holding the field components to which each boundary condition applies
411ce78bad3SBarry Smith . bcPoints - An array of size `numBC` of `IS` holding the `DMPLEX` points to which each boundary condition applies
412ce78bad3SBarry Smith - perm     - Optional permutation of the chart, or `NULL`
4139cbb8f0dSMatthew G. Knepley 
4149cbb8f0dSMatthew G. Knepley   Output Parameter:
415a1cb98faSBarry Smith . section - The `PetscSection` object
416a1cb98faSBarry Smith 
417a1cb98faSBarry Smith   Level: developer
4189cbb8f0dSMatthew G. Knepley 
4199cbb8f0dSMatthew G. Knepley   Notes:
4209cbb8f0dSMatthew G. Knepley   numDof[f*(dim+1)+d] gives the number of dof for field f on points of dimension d. For instance, numDof[1] is the
4219cbb8f0dSMatthew G. Knepley   number of dof for field 0 on each edge.
4229cbb8f0dSMatthew G. Knepley 
423a1cb98faSBarry Smith   The chart permutation is the same one set using `PetscSectionSetPermutation()`
4249cbb8f0dSMatthew G. Knepley 
4251cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `PetscSectionCreate()`, `PetscSectionSetPermutation()`
4269cbb8f0dSMatthew G. Knepley @*/
DMPlexCreateSection(DM dm,DMLabel label[],const PetscInt numComp[],const PetscInt numDof[],PetscInt numBC,const PetscInt bcField[],const IS bcComps[],const IS bcPoints[],IS perm,PetscSection * section)427d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateSection(DM dm, DMLabel label[], const PetscInt numComp[], const PetscInt numDof[], PetscInt numBC, const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], IS perm, PetscSection *section)
428d71ae5a4SJacob Faibussowitsch {
4299cbb8f0dSMatthew G. Knepley   PetscSection aSec;
4309cbb8f0dSMatthew G. Knepley 
4319cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
4329566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSectionFields(dm, numComp, section));
4339566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSectionDof(dm, label, numDof, *section));
4349566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSectionBCDof(dm, numBC, bcField, bcComps, bcPoints, *section));
4359566063dSJacob Faibussowitsch   if (perm) PetscCall(PetscSectionSetPermutation(*section, perm));
4369566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFromOptions(*section));
4379566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(*section));
4389566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAnchors(dm, &aSec, NULL));
4399cbb8f0dSMatthew G. Knepley   if (numBC || aSec) {
4409566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSectionBCIndicesField(dm, numBC, bcField, bcComps, bcPoints, *section));
4419566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSectionBCIndices(dm, *section));
4429cbb8f0dSMatthew G. Knepley   }
4439566063dSJacob Faibussowitsch   PetscCall(PetscSectionViewFromOptions(*section, NULL, "-section_view"));
4443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4459cbb8f0dSMatthew G. Knepley }
4469cbb8f0dSMatthew G. Knepley 
DMCreateLocalSection_Plex(DM dm)447d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateLocalSection_Plex(DM dm)
448d71ae5a4SJacob Faibussowitsch {
4499cbb8f0dSMatthew G. Knepley   PetscSection section;
450baef929fSMatthew G. Knepley   DMLabel     *labels;
451d02c7345SMatthew G. Knepley   IS          *bcPoints, *bcComps, permIS;
452d02c7345SMatthew G. Knepley   PetscBT      blockStarts;
4539cbb8f0dSMatthew G. Knepley   PetscBool   *isFE;
4549cbb8f0dSMatthew G. Knepley   PetscInt    *bcFields, *numComp, *numDof;
4559310035eSMatthew G. Knepley   PetscInt     depth, dim, numBC = 0, Nf, Nds, s, bc = 0, f;
4569cbb8f0dSMatthew G. Knepley   PetscInt     cStart, cEnd, cEndInterior;
4579cbb8f0dSMatthew G. Knepley 
4589cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
4599566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
4609566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
4619566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
4629cbb8f0dSMatthew G. Knepley   /* FE and FV boundary conditions are handled slightly differently */
4639566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nf, &isFE));
464baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
4659cbb8f0dSMatthew G. Knepley     PetscObject  obj;
4669cbb8f0dSMatthew G. Knepley     PetscClassId id;
4679cbb8f0dSMatthew G. Knepley 
4689566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, f, NULL, &obj));
4699566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(obj, &id));
4709371c9d4SSatish Balay     if (id == PETSCFE_CLASSID) {
4719371c9d4SSatish Balay       isFE[f] = PETSC_TRUE;
4729371c9d4SSatish Balay     } else if (id == PETSCFV_CLASSID) {
4739371c9d4SSatish Balay       isFE[f] = PETSC_FALSE;
4749371c9d4SSatish Balay     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
4759cbb8f0dSMatthew G. Knepley   }
4769cbb8f0dSMatthew G. Knepley   /* Allocate boundary point storage for FEM boundaries */
4779566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
4789310035eSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
4799310035eSMatthew G. Knepley     PetscDS  dsBC;
4809310035eSMatthew G. Knepley     PetscInt numBd, bd;
4819310035eSMatthew G. Knepley 
48207218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC, NULL));
4839566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumBoundary(dsBC, &numBd));
48408401ef6SPierre Jolivet     PetscCheck(Nf || !numBd, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "number of fields is zero and number of boundary conditions is nonzero (this should never happen)");
4859cbb8f0dSMatthew G. Knepley     for (bd = 0; bd < numBd; ++bd) {
4869cbb8f0dSMatthew G. Knepley       PetscInt                field;
4879cbb8f0dSMatthew G. Knepley       DMBoundaryConditionType type;
4889cbb8f0dSMatthew G. Knepley       DMLabel                 label;
4899cbb8f0dSMatthew G. Knepley 
4909566063dSJacob Faibussowitsch       PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, &type, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL));
4919cbb8f0dSMatthew G. Knepley       if (label && isFE[field] && (type & DM_BC_ESSENTIAL)) ++numBC;
4929cbb8f0dSMatthew G. Knepley     }
4939310035eSMatthew G. Knepley   }
4949cbb8f0dSMatthew G. Knepley   /* Add ghost cell boundaries for FVM */
4952827ebadSStefano Zampini   PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &cEndInterior, NULL));
4969371c9d4SSatish Balay   for (f = 0; f < Nf; ++f)
4979371c9d4SSatish Balay     if (!isFE[f] && cEndInterior >= 0) ++numBC;
4989566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(numBC, &bcFields, numBC, &bcPoints, numBC, &bcComps));
4999cbb8f0dSMatthew G. Knepley   /* Constrain ghost cells for FV */
500baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
5019cbb8f0dSMatthew G. Knepley     PetscInt *newidx, c;
5029cbb8f0dSMatthew G. Knepley 
5039cbb8f0dSMatthew G. Knepley     if (isFE[f] || cEndInterior < 0) continue;
5049566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(cEnd - cEndInterior, &newidx));
5059cbb8f0dSMatthew G. Knepley     for (c = cEndInterior; c < cEnd; ++c) newidx[c - cEndInterior] = c;
5069cbb8f0dSMatthew G. Knepley     bcFields[bc] = f;
5079566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, cEnd - cEndInterior, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]));
5089cbb8f0dSMatthew G. Knepley   }
509f985c9a8SMatthew G. Knepley   /* Complete labels for boundary conditions */
510799db056SMatthew G. Knepley   PetscCall(DMCompleteBCLabels_Internal(dm));
5119cbb8f0dSMatthew G. Knepley   /* Handle FEM Dirichlet boundaries */
5129566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
5139310035eSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
5149310035eSMatthew G. Knepley     PetscDS  dsBC;
5159310035eSMatthew G. Knepley     PetscInt numBd, bd;
5169310035eSMatthew G. Knepley 
51707218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC, NULL));
5189566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumBoundary(dsBC, &numBd));
5199cbb8f0dSMatthew G. Knepley     for (bd = 0; bd < numBd; ++bd) {
5209cbb8f0dSMatthew G. Knepley       DMLabel                 label;
5219cbb8f0dSMatthew G. Knepley       const PetscInt         *comps;
5229cbb8f0dSMatthew G. Knepley       const PetscInt         *values;
5239310035eSMatthew G. Knepley       PetscInt                bd2, field, numComps, numValues;
5249cbb8f0dSMatthew G. Knepley       DMBoundaryConditionType type;
5259310035eSMatthew G. Knepley       PetscBool               duplicate = PETSC_FALSE;
5269cbb8f0dSMatthew G. Knepley 
5279566063dSJacob Faibussowitsch       PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, &type, NULL, &label, &numValues, &values, &field, &numComps, &comps, NULL, NULL, NULL));
5289cbb8f0dSMatthew G. Knepley       if (!isFE[field] || !label) continue;
5299310035eSMatthew G. Knepley       /* Only want to modify label once */
5309310035eSMatthew G. Knepley       for (bd2 = 0; bd2 < bd; ++bd2) {
53145480ffeSMatthew G. Knepley         DMLabel l;
53245480ffeSMatthew G. Knepley 
5339566063dSJacob Faibussowitsch         PetscCall(PetscDSGetBoundary(dsBC, bd2, NULL, NULL, NULL, &l, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
53445480ffeSMatthew G. Knepley         duplicate = l == label ? PETSC_TRUE : PETSC_FALSE;
5359310035eSMatthew G. Knepley         if (duplicate) break;
5369310035eSMatthew G. Knepley       }
5379cbb8f0dSMatthew G. Knepley       /* Filter out cells, if you actually want to constrain cells you need to do things by hand right now */
5389cbb8f0dSMatthew G. Knepley       if (type & DM_BC_ESSENTIAL) {
5399cbb8f0dSMatthew G. Knepley         PetscInt *newidx;
5409cbb8f0dSMatthew G. Knepley         PetscInt  n, newn = 0, p, v;
5419cbb8f0dSMatthew G. Knepley 
5429cbb8f0dSMatthew G. Knepley         bcFields[bc] = field;
54340cbb1a0SMatthew G. Knepley         if (numComps) PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numComps, comps, PETSC_COPY_VALUES, &bcComps[bc]));
5449cbb8f0dSMatthew G. Knepley         for (v = 0; v < numValues; ++v) {
5459cbb8f0dSMatthew G. Knepley           IS              tmp;
5469cbb8f0dSMatthew G. Knepley           const PetscInt *idx;
5479cbb8f0dSMatthew G. Knepley 
5489566063dSJacob Faibussowitsch           PetscCall(DMLabelGetStratumIS(label, values[v], &tmp));
5499cbb8f0dSMatthew G. Knepley           if (!tmp) continue;
5509566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(tmp, &n));
5519566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(tmp, &idx));
5529cbb8f0dSMatthew G. Knepley           if (isFE[field]) {
5539371c9d4SSatish Balay             for (p = 0; p < n; ++p)
5549371c9d4SSatish Balay               if ((idx[p] < cStart) || (idx[p] >= cEnd)) ++newn;
5559cbb8f0dSMatthew G. Knepley           } else {
5569371c9d4SSatish Balay             for (p = 0; p < n; ++p)
5579371c9d4SSatish Balay               if ((idx[p] >= cStart) || (idx[p] < cEnd)) ++newn;
5589cbb8f0dSMatthew G. Knepley           }
5599566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(tmp, &idx));
5609566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&tmp));
5619cbb8f0dSMatthew G. Knepley         }
5629566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(newn, &newidx));
5639cbb8f0dSMatthew G. Knepley         newn = 0;
5649cbb8f0dSMatthew G. Knepley         for (v = 0; v < numValues; ++v) {
5659cbb8f0dSMatthew G. Knepley           IS              tmp;
5669cbb8f0dSMatthew G. Knepley           const PetscInt *idx;
5679cbb8f0dSMatthew G. Knepley 
5689566063dSJacob Faibussowitsch           PetscCall(DMLabelGetStratumIS(label, values[v], &tmp));
5699cbb8f0dSMatthew G. Knepley           if (!tmp) continue;
5709566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(tmp, &n));
5719566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(tmp, &idx));
5729cbb8f0dSMatthew G. Knepley           if (isFE[field]) {
5739371c9d4SSatish Balay             for (p = 0; p < n; ++p)
5749371c9d4SSatish Balay               if ((idx[p] < cStart) || (idx[p] >= cEnd)) newidx[newn++] = idx[p];
5759cbb8f0dSMatthew G. Knepley           } else {
5769371c9d4SSatish Balay             for (p = 0; p < n; ++p)
5779371c9d4SSatish Balay               if ((idx[p] >= cStart) || (idx[p] < cEnd)) newidx[newn++] = idx[p];
5789cbb8f0dSMatthew G. Knepley           }
5799566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(tmp, &idx));
5809566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&tmp));
5819cbb8f0dSMatthew G. Knepley         }
5829566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, newn, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]));
5839310035eSMatthew G. Knepley       }
5849cbb8f0dSMatthew G. Knepley     }
5859cbb8f0dSMatthew G. Knepley   }
5869cbb8f0dSMatthew G. Knepley   /* Handle discretization */
5879566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(Nf, &labels, Nf, &numComp, Nf * (dim + 1), &numDof));
588baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
58944a7f3ddSMatthew G. Knepley     labels[f] = dm->fields[f].label;
5909cbb8f0dSMatthew G. Knepley     if (isFE[f]) {
59144a7f3ddSMatthew G. Knepley       PetscFE         fe = (PetscFE)dm->fields[f].disc;
5929cbb8f0dSMatthew G. Knepley       const PetscInt *numFieldDof;
593a3254110SMatthew Knepley       PetscInt        fedim, d;
5949cbb8f0dSMatthew G. Knepley 
5959566063dSJacob Faibussowitsch       PetscCall(PetscFEGetNumComponents(fe, &numComp[f]));
5969566063dSJacob Faibussowitsch       PetscCall(PetscFEGetNumDof(fe, &numFieldDof));
5979566063dSJacob Faibussowitsch       PetscCall(PetscFEGetSpatialDimension(fe, &fedim));
598a3254110SMatthew Knepley       for (d = 0; d < PetscMin(dim, fedim) + 1; ++d) numDof[f * (dim + 1) + d] = numFieldDof[d];
5999cbb8f0dSMatthew G. Knepley     } else {
60044a7f3ddSMatthew G. Knepley       PetscFV fv = (PetscFV)dm->fields[f].disc;
6019cbb8f0dSMatthew G. Knepley 
6029566063dSJacob Faibussowitsch       PetscCall(PetscFVGetNumComponents(fv, &numComp[f]));
6039cbb8f0dSMatthew G. Knepley       numDof[f * (dim + 1) + dim] = numComp[f];
6049cbb8f0dSMatthew G. Knepley     }
6059cbb8f0dSMatthew G. Knepley   }
6069566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
607baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
6089cbb8f0dSMatthew G. Knepley     PetscInt d;
609ad540459SPierre Jolivet     for (d = 1; d < dim; ++d) PetscCheck(numDof[f * (dim + 1) + d] <= 0 || depth >= dim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated when unknowns are specified on edges or faces.");
6109cbb8f0dSMatthew G. Knepley   }
611adc21957SMatthew G. Knepley   PetscCall(DMCreateSectionPermutation(dm, &permIS, &blockStarts));
612d02c7345SMatthew G. Knepley   PetscCall(DMPlexCreateSection(dm, labels, numComp, numDof, numBC, bcFields, bcComps, bcPoints, permIS, &section));
613d02c7345SMatthew G. Knepley   section->blockStarts = blockStarts;
614d02c7345SMatthew G. Knepley   PetscCall(ISDestroy(&permIS));
615baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
6169cbb8f0dSMatthew G. Knepley     PetscFE     fe;
6179cbb8f0dSMatthew G. Knepley     const char *name;
6189cbb8f0dSMatthew G. Knepley 
6199566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, f, NULL, (PetscObject *)&fe));
620cdfe53dfSMatthew G. Knepley     if (!((PetscObject)fe)->name) continue;
6219566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)fe, &name));
6229566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldName(section, f, name));
6239cbb8f0dSMatthew G. Knepley   }
6249566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, section));
6259566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&section));
6269371c9d4SSatish Balay   for (bc = 0; bc < numBC; ++bc) {
6279371c9d4SSatish Balay     PetscCall(ISDestroy(&bcPoints[bc]));
6289371c9d4SSatish Balay     PetscCall(ISDestroy(&bcComps[bc]));
6299371c9d4SSatish Balay   }
6309566063dSJacob Faibussowitsch   PetscCall(PetscFree3(bcFields, bcPoints, bcComps));
6319566063dSJacob Faibussowitsch   PetscCall(PetscFree3(labels, numComp, numDof));
6329566063dSJacob Faibussowitsch   PetscCall(PetscFree(isFE));
633d2b2dc1eSMatthew G. Knepley   /* Checking for CEED usage */
634d2b2dc1eSMatthew G. Knepley   {
635d2b2dc1eSMatthew G. Knepley     PetscBool useCeed;
636d2b2dc1eSMatthew G. Knepley 
637d2b2dc1eSMatthew G. Knepley     PetscCall(DMPlexGetUseCeed(dm, &useCeed));
638d2b2dc1eSMatthew G. Knepley     if (useCeed) {
639d2b2dc1eSMatthew G. Knepley       PetscCall(DMPlexSetUseMatClosurePermutation(dm, PETSC_FALSE));
640d2b2dc1eSMatthew G. Knepley       PetscCall(DMUseTensorOrder(dm, PETSC_TRUE));
641d2b2dc1eSMatthew G. Knepley     }
642d2b2dc1eSMatthew G. Knepley   }
6433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6449cbb8f0dSMatthew G. Knepley }
645