xref: /petsc/src/dm/impls/plex/plexsection.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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 */
4baef929fSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionFields(DM dm, const PetscInt numComp[], PetscSection *section)
59cbb8f0dSMatthew G. Knepley {
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));
219cbb8f0dSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {isFE[f] = PETSC_TRUE;}
229cbb8f0dSMatthew G. Knepley     else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;}
239cbb8f0dSMatthew G. Knepley   }
24baef929fSMatthew G. Knepley 
259566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), section));
26baef929fSMatthew G. Knepley   if (Nf > 0) {
279566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetNumFields(*section, Nf));
289cbb8f0dSMatthew G. Knepley     if (numComp) {
29baef929fSMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
309566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetFieldComponents(*section, f, numComp[f]));
319cbb8f0dSMatthew G. Knepley         if (isFE[f]) {
329cbb8f0dSMatthew G. Knepley           PetscFE           fe;
339cbb8f0dSMatthew G. Knepley           PetscDualSpace    dspace;
349cbb8f0dSMatthew G. Knepley           const PetscInt    ***perms;
359cbb8f0dSMatthew G. Knepley           const PetscScalar ***flips;
369cbb8f0dSMatthew G. Knepley           const PetscInt    *numDof;
379cbb8f0dSMatthew G. Knepley 
389566063dSJacob Faibussowitsch           PetscCall(DMGetField(dm,f,NULL,(PetscObject *) &fe));
399566063dSJacob Faibussowitsch           PetscCall(PetscFEGetDualSpace(fe,&dspace));
409566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceGetSymmetries(dspace,&perms,&flips));
419566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceGetNumDof(dspace,&numDof));
429cbb8f0dSMatthew G. Knepley           if (perms || flips) {
439cbb8f0dSMatthew G. Knepley             DM              K;
444dff28b8SMatthew G. Knepley             PetscInt        sph, spdepth;
459cbb8f0dSMatthew G. Knepley             PetscSectionSym sym;
469cbb8f0dSMatthew G. Knepley 
479566063dSJacob Faibussowitsch             PetscCall(PetscDualSpaceGetDM(dspace,&K));
489566063dSJacob Faibussowitsch             PetscCall(DMPlexGetDepth(K, &spdepth));
499566063dSJacob Faibussowitsch             PetscCall(PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)*section),depthLabel,&sym));
504dff28b8SMatthew G. Knepley             for (sph = 0; sph <= spdepth; sph++) {
519cbb8f0dSMatthew G. Knepley               PetscDualSpace    hspace;
529cbb8f0dSMatthew G. Knepley               PetscInt          kStart, kEnd;
534dff28b8SMatthew G. Knepley               PetscInt          kConeSize, h = sph + (depth - spdepth);
549cbb8f0dSMatthew G. Knepley               const PetscInt    **perms0 = NULL;
559cbb8f0dSMatthew G. Knepley               const PetscScalar **flips0 = NULL;
569cbb8f0dSMatthew G. Knepley 
579566063dSJacob Faibussowitsch               PetscCall(PetscDualSpaceGetHeightSubspace(dspace, sph, &hspace));
589566063dSJacob Faibussowitsch               PetscCall(DMPlexGetHeightStratum(K, h, &kStart, &kEnd));
599cbb8f0dSMatthew G. Knepley               if (!hspace) continue;
609566063dSJacob Faibussowitsch               PetscCall(PetscDualSpaceGetSymmetries(hspace,&perms,&flips));
619cbb8f0dSMatthew G. Knepley               if (perms) perms0 = perms[0];
629cbb8f0dSMatthew G. Knepley               if (flips) flips0 = flips[0];
639cbb8f0dSMatthew G. Knepley               if (!(perms0 || flips0)) continue;
64b5a892a1SMatthew G. Knepley               {
65b5a892a1SMatthew G. Knepley                 DMPolytopeType ct;
66b5a892a1SMatthew G. Knepley                 /* The number of arrangements is no longer based on the number of faces */
679566063dSJacob Faibussowitsch                 PetscCall(DMPlexGetCellType(K, kStart, &ct));
68b5a892a1SMatthew G. Knepley                 kConeSize = DMPolytopeTypeGetNumArrangments(ct) / 2;
69b5a892a1SMatthew G. Knepley               }
709566063dSJacob Faibussowitsch               PetscCall(PetscSectionSymLabelSetStratum(sym,depth - h,numDof[depth - h],-kConeSize,kConeSize,PETSC_USE_POINTER,perms0 ? &perms0[-kConeSize] : NULL,flips0 ? &flips0[-kConeSize] : NULL));
719cbb8f0dSMatthew G. Knepley             }
729566063dSJacob Faibussowitsch             PetscCall(PetscSectionSetFieldSym(*section,f,sym));
739566063dSJacob Faibussowitsch             PetscCall(PetscSectionSymDestroy(&sym));
749cbb8f0dSMatthew G. Knepley           }
759cbb8f0dSMatthew G. Knepley         }
769cbb8f0dSMatthew G. Knepley       }
779cbb8f0dSMatthew G. Knepley     }
789cbb8f0dSMatthew G. Knepley   }
799566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
809566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(*section, pStart, pEnd));
819566063dSJacob Faibussowitsch   PetscCall(PetscFree(isFE));
82baef929fSMatthew G. Knepley   PetscFunctionReturn(0);
83baef929fSMatthew G. Knepley }
84baef929fSMatthew G. Knepley 
85baef929fSMatthew G. Knepley /* Set the number of dof on each point and separate by fields */
86baef929fSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionDof(DM dm, DMLabel label[],const PetscInt numDof[], PetscSection section)
87baef929fSMatthew G. Knepley {
88baef929fSMatthew G. Knepley   DMLabel        depthLabel;
89412e9a14SMatthew G. Knepley   DMPolytopeType ct;
90baef929fSMatthew G. Knepley   PetscInt       depth, cellHeight, pStart = 0, pEnd = 0;
91a55f9a55SMatthew G. Knepley   PetscInt       Nf, f, Nds, n, dim, d, dep, p;
925fedec97SMatthew G. Knepley   PetscBool     *isFE, hasCohesive = PETSC_FALSE;
93baef929fSMatthew G. Knepley 
94baef929fSMatthew G. Knepley   PetscFunctionBegin;
959566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
969566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
979566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(dm,&depthLabel));
989566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
999566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
100a55f9a55SMatthew G. Knepley   for (n = 0; n < Nds; ++n) {
101a55f9a55SMatthew G. Knepley     PetscDS   ds;
1025fedec97SMatthew G. Knepley     PetscBool isCohesive;
103a55f9a55SMatthew G. Knepley 
1049566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds));
1059566063dSJacob Faibussowitsch     PetscCall(PetscDSIsCohesive(ds, &isCohesive));
1065fedec97SMatthew G. Knepley     if (isCohesive) {hasCohesive = PETSC_TRUE; break;}
107a55f9a55SMatthew G. Knepley   }
1089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nf, &isFE));
109baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
110baef929fSMatthew G. Knepley     PetscObject  obj;
111baef929fSMatthew G. Knepley     PetscClassId id;
112baef929fSMatthew G. Knepley 
1139566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, f, NULL, &obj));
1149566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(obj, &id));
1151274e1a1SLawrence Mitchell     /* User is allowed to put a "placeholder" field in (c.f. DMCreateDS) */
1161274e1a1SLawrence Mitchell     isFE[f] = id == PETSCFE_CLASSID ? PETSC_TRUE : PETSC_FALSE;
117baef929fSMatthew G. Knepley   }
118baef929fSMatthew G. Knepley 
1199566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
120baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
121e0b68406SMatthew Knepley     PetscBool avoidTensor;
122e0b68406SMatthew Knepley 
1239566063dSJacob Faibussowitsch     PetscCall(DMGetFieldAvoidTensor(dm, f, &avoidTensor));
1245fedec97SMatthew G. Knepley     avoidTensor = (avoidTensor || hasCohesive) ? PETSC_TRUE : PETSC_FALSE;
125baef929fSMatthew G. Knepley     if (label && label[f]) {
126baef929fSMatthew G. Knepley       IS              pointIS;
127baef929fSMatthew G. Knepley       const PetscInt *points;
128baef929fSMatthew G. Knepley       PetscInt        n;
129baef929fSMatthew G. Knepley 
1309566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(label[f], 1, &pointIS));
131baef929fSMatthew G. Knepley       if (!pointIS) continue;
1329566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(pointIS, &n));
1339566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(pointIS, &points));
134baef929fSMatthew G. Knepley       for (p = 0; p < n; ++p) {
135baef929fSMatthew G. Knepley         const PetscInt point = points[p];
136baef929fSMatthew G. Knepley         PetscInt       dof, d;
137baef929fSMatthew G. Knepley 
1389566063dSJacob Faibussowitsch         PetscCall(DMPlexGetCellType(dm, point, &ct));
1399566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValue(depthLabel, point, &d));
140412e9a14SMatthew G. Knepley         /* If this is a tensor prism point, use dof for one dimension lower */
141412e9a14SMatthew G. Knepley         switch (ct) {
142412e9a14SMatthew G. Knepley           case DM_POLYTOPE_POINT_PRISM_TENSOR:
143412e9a14SMatthew G. Knepley           case DM_POLYTOPE_SEG_PRISM_TENSOR:
144412e9a14SMatthew G. Knepley           case DM_POLYTOPE_TRI_PRISM_TENSOR:
145412e9a14SMatthew G. Knepley           case DM_POLYTOPE_QUAD_PRISM_TENSOR:
1465fedec97SMatthew G. Knepley             if (hasCohesive) {--d;} break;
147412e9a14SMatthew G. Knepley           default: break;
148412e9a14SMatthew G. Knepley         }
149baef929fSMatthew G. Knepley         dof  = d < 0 ? 0 : numDof[f*(dim+1)+d];
1509566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetFieldDof(section, point, f, dof));
1519566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(section, point, dof));
152baef929fSMatthew G. Knepley       }
1539566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(pointIS, &points));
1549566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
155baef929fSMatthew G. Knepley     } else {
1569cbb8f0dSMatthew G. Knepley       for (dep = 0; dep <= depth - cellHeight; ++dep) {
1573fd2e703SMatthew G. Knepley         /* Cases: dim > depth (cell-vertex mesh), dim == depth (fully interpolated), dim < depth (interpolated submesh) */
1583fd2e703SMatthew G. Knepley         d    = dim <= depth ? dep : (!dep ? 0 : dim);
1599566063dSJacob Faibussowitsch         PetscCall(DMPlexGetDepthStratum(dm, dep, &pStart, &pEnd));
1609cbb8f0dSMatthew G. Knepley         for (p = pStart; p < pEnd; ++p) {
161baef929fSMatthew G. Knepley           const PetscInt dof = numDof[f*(dim+1)+d];
1629cbb8f0dSMatthew G. Knepley 
1639566063dSJacob Faibussowitsch           PetscCall(DMPlexGetCellType(dm, p, &ct));
164412e9a14SMatthew G. Knepley           switch (ct) {
165412e9a14SMatthew G. Knepley             case DM_POLYTOPE_POINT_PRISM_TENSOR:
166412e9a14SMatthew G. Knepley             case DM_POLYTOPE_SEG_PRISM_TENSOR:
167412e9a14SMatthew G. Knepley             case DM_POLYTOPE_TRI_PRISM_TENSOR:
168412e9a14SMatthew G. Knepley             case DM_POLYTOPE_QUAD_PRISM_TENSOR:
169e0b68406SMatthew Knepley               if (avoidTensor && isFE[f]) continue;
170412e9a14SMatthew G. Knepley             default: break;
171412e9a14SMatthew G. Knepley           }
1729566063dSJacob Faibussowitsch           PetscCall(PetscSectionSetFieldDof(section, p, f, dof));
1739566063dSJacob Faibussowitsch           PetscCall(PetscSectionAddDof(section, p, dof));
1749cbb8f0dSMatthew G. Knepley         }
175baef929fSMatthew G. Knepley       }
1769cbb8f0dSMatthew G. Knepley     }
1779cbb8f0dSMatthew G. Knepley   }
1789566063dSJacob Faibussowitsch   PetscCall(PetscFree(isFE));
1799cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
1809cbb8f0dSMatthew G. Knepley }
1819cbb8f0dSMatthew G. Knepley 
1829cbb8f0dSMatthew G. Knepley /* Set the number of dof on each point and separate by fields
1839cbb8f0dSMatthew G. Knepley    If bcComps is NULL or the IS is NULL, constrain every dof on the point
1849cbb8f0dSMatthew G. Knepley */
1859cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCDof(DM dm, PetscInt numBC, const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section)
1869cbb8f0dSMatthew G. Knepley {
187baef929fSMatthew G. Knepley   PetscInt       Nf;
1889cbb8f0dSMatthew G. Knepley   PetscInt       bc;
1899cbb8f0dSMatthew G. Knepley   PetscSection   aSec;
1909cbb8f0dSMatthew G. Knepley 
1919cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
1929566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &Nf));
1939cbb8f0dSMatthew G. Knepley   for (bc = 0; bc < numBC; ++bc) {
1949cbb8f0dSMatthew G. Knepley     PetscInt        field = 0;
1959cbb8f0dSMatthew G. Knepley     const PetscInt *comp;
1969cbb8f0dSMatthew G. Knepley     const PetscInt *idx;
197ce093827SMatthew G. Knepley     PetscInt        Nc = 0, cNc = -1, n, i;
1989cbb8f0dSMatthew G. Knepley 
199ce093827SMatthew G. Knepley     if (Nf) {
200ce093827SMatthew G. Knepley       field = bcField[bc];
2019566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldComponents(section, field, &Nc));
202ce093827SMatthew G. Knepley     }
2039566063dSJacob Faibussowitsch     if (bcComps && bcComps[bc]) PetscCall(ISGetLocalSize(bcComps[bc], &cNc));
2049566063dSJacob Faibussowitsch     if (bcComps && bcComps[bc]) PetscCall(ISGetIndices(bcComps[bc], &comp));
2059566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(bcPoints[bc], &n));
2069566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(bcPoints[bc], &idx));
2079cbb8f0dSMatthew G. Knepley     for (i = 0; i < n; ++i) {
2089cbb8f0dSMatthew G. Knepley       const PetscInt p = idx[i];
2099cbb8f0dSMatthew G. Knepley       PetscInt       numConst;
2109cbb8f0dSMatthew G. Knepley 
211baef929fSMatthew G. Knepley       if (Nf) {
2129566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetFieldDof(section, p, field, &numConst));
2139cbb8f0dSMatthew G. Knepley       } else {
2149566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, p, &numConst));
2159cbb8f0dSMatthew G. Knepley       }
216ce093827SMatthew G. Knepley       /* If Nc <= 0, constrain every dof on the point */
217ce093827SMatthew G. Knepley       if (cNc > 0) {
218ce093827SMatthew G. Knepley         /* We assume that a point may have multiple "nodes", which are collections of Nc dofs,
219ce093827SMatthew G. Knepley            and that those dofs are numbered n*Nc+c */
220ce093827SMatthew G. Knepley         if (Nf) {
221*08401ef6SPierre Jolivet           PetscCheck(numConst % Nc == 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %D has %D dof which is not divisible by %D field components", p, numConst, Nc);
222ce093827SMatthew G. Knepley           numConst = (numConst/Nc) * cNc;
223ce093827SMatthew G. Knepley         } else {
224ce093827SMatthew G. Knepley           numConst = PetscMin(numConst, cNc);
225ce093827SMatthew G. Knepley         }
226ce093827SMatthew G. Knepley       }
2279566063dSJacob Faibussowitsch       if (Nf) PetscCall(PetscSectionAddFieldConstraintDof(section, p, field, numConst));
2289566063dSJacob Faibussowitsch       PetscCall(PetscSectionAddConstraintDof(section, p, numConst));
2299cbb8f0dSMatthew G. Knepley     }
2309566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(bcPoints[bc], &idx));
2319566063dSJacob Faibussowitsch     if (bcComps && bcComps[bc]) PetscCall(ISRestoreIndices(bcComps[bc], &comp));
2329cbb8f0dSMatthew G. Knepley   }
2339566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAnchors(dm, &aSec, NULL));
2349cbb8f0dSMatthew G. Knepley   if (aSec) {
2359cbb8f0dSMatthew G. Knepley     PetscInt aStart, aEnd, a;
2369cbb8f0dSMatthew G. Knepley 
2379566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd));
2389cbb8f0dSMatthew G. Knepley     for (a = aStart; a < aEnd; a++) {
2399cbb8f0dSMatthew G. Knepley       PetscInt dof, f;
2409cbb8f0dSMatthew G. Knepley 
2419566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(aSec, a, &dof));
2429cbb8f0dSMatthew G. Knepley       if (dof) {
2439cbb8f0dSMatthew G. Knepley         /* if there are point-to-point constraints, then all dofs are constrained */
2449566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, a, &dof));
2459566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetConstraintDof(section, a, dof));
246baef929fSMatthew G. Knepley         for (f = 0; f < Nf; f++) {
2479566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldDof(section, a, f, &dof));
2489566063dSJacob Faibussowitsch           PetscCall(PetscSectionSetFieldConstraintDof(section, a, f, dof));
2499cbb8f0dSMatthew G. Knepley         }
2509cbb8f0dSMatthew G. Knepley       }
2519cbb8f0dSMatthew G. Knepley     }
2529cbb8f0dSMatthew G. Knepley   }
2539cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
2549cbb8f0dSMatthew G. Knepley }
2559cbb8f0dSMatthew G. Knepley 
2569cbb8f0dSMatthew G. Knepley /* Set the constrained field indices on each point
2579cbb8f0dSMatthew G. Knepley    If bcComps is NULL or the IS is NULL, constrain every dof on the point
2589cbb8f0dSMatthew G. Knepley */
2599cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCIndicesField(DM dm, PetscInt numBC,const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section)
2609cbb8f0dSMatthew G. Knepley {
2619cbb8f0dSMatthew G. Knepley   PetscSection   aSec;
2629cbb8f0dSMatthew G. Knepley   PetscInt      *indices;
263baef929fSMatthew G. Knepley   PetscInt       Nf, cdof, maxDof = 0, pStart, pEnd, p, bc, f, d;
2649cbb8f0dSMatthew G. Knepley 
2659cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
2669566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &Nf));
267baef929fSMatthew G. Knepley   if (!Nf) PetscFunctionReturn(0);
2689cbb8f0dSMatthew G. Knepley   /* Initialize all field indices to -1 */
2699566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
2709566063dSJacob Faibussowitsch   for (p = pStart; p < pEnd; ++p) {PetscCall(PetscSectionGetConstraintDof(section, p, &cdof)); maxDof = PetscMax(maxDof, cdof);}
2719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxDof, &indices));
2729cbb8f0dSMatthew G. Knepley   for (d = 0; d < maxDof; ++d) indices[d] = -1;
2739566063dSJacob Faibussowitsch   for (p = pStart; p < pEnd; ++p) for (f = 0; f < Nf; ++f) PetscCall(PetscSectionSetFieldConstraintIndices(section, p, f, indices));
2749cbb8f0dSMatthew G. Knepley   /* Handle BC constraints */
2759cbb8f0dSMatthew G. Knepley   for (bc = 0; bc < numBC; ++bc) {
2769cbb8f0dSMatthew G. Knepley     const PetscInt  field = bcField[bc];
2779cbb8f0dSMatthew G. Knepley     const PetscInt *comp, *idx;
278ce093827SMatthew G. Knepley     PetscInt        Nc, cNc = -1, n, i;
2799cbb8f0dSMatthew G. Knepley 
2809566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldComponents(section, field, &Nc));
2819566063dSJacob Faibussowitsch     if (bcComps && bcComps[bc]) PetscCall(ISGetLocalSize(bcComps[bc], &cNc));
2829566063dSJacob Faibussowitsch     if (bcComps && bcComps[bc]) PetscCall(ISGetIndices(bcComps[bc], &comp));
2839566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(bcPoints[bc], &n));
2849566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(bcPoints[bc], &idx));
2859cbb8f0dSMatthew G. Knepley     for (i = 0; i < n; ++i) {
2869cbb8f0dSMatthew G. Knepley       const PetscInt  p = idx[i];
2879cbb8f0dSMatthew G. Knepley       const PetscInt *find;
288ce093827SMatthew G. Knepley       PetscInt        fdof, fcdof, c, j;
2899cbb8f0dSMatthew G. Knepley 
2909566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldDof(section, p, field, &fdof));
2919cbb8f0dSMatthew G. Knepley       if (!fdof) continue;
292ce093827SMatthew G. Knepley       if (cNc < 0) {
2939cbb8f0dSMatthew G. Knepley         for (d = 0; d < fdof; ++d) indices[d] = d;
2949cbb8f0dSMatthew G. Knepley         fcdof = fdof;
2959cbb8f0dSMatthew G. Knepley       } else {
296ce093827SMatthew G. Knepley         /* We assume that a point may have multiple "nodes", which are collections of Nc dofs,
297ce093827SMatthew G. Knepley            and that those dofs are numbered n*Nc+c */
2989566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetFieldConstraintDof(section, p, field, &fcdof));
2999566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetFieldConstraintIndices(section, p, field, &find));
300ce093827SMatthew G. Knepley         /* Get indices constrained by previous bcs */
3019cbb8f0dSMatthew G. Knepley         for (d = 0; d < fcdof; ++d) {if (find[d] < 0) break; indices[d] = find[d];}
302ce093827SMatthew G. Knepley         for (j = 0; j < fdof/Nc; ++j) for (c = 0; c < cNc; ++c) indices[d++] = j*Nc + comp[c];
3039566063dSJacob Faibussowitsch         PetscCall(PetscSortRemoveDupsInt(&d, indices));
3049cbb8f0dSMatthew G. Knepley         for (c = d; c < fcdof; ++c) indices[c] = -1;
3059cbb8f0dSMatthew G. Knepley         fcdof = d;
3069cbb8f0dSMatthew G. Knepley       }
3079566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetFieldConstraintDof(section, p, field, fcdof));
3089566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetFieldConstraintIndices(section, p, field, indices));
3099cbb8f0dSMatthew G. Knepley     }
3109566063dSJacob Faibussowitsch     if (bcComps && bcComps[bc]) PetscCall(ISRestoreIndices(bcComps[bc], &comp));
3119566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(bcPoints[bc], &idx));
3129cbb8f0dSMatthew G. Knepley   }
3139cbb8f0dSMatthew G. Knepley   /* Handle anchors */
3149566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAnchors(dm, &aSec, NULL));
3159cbb8f0dSMatthew G. Knepley   if (aSec) {
3169cbb8f0dSMatthew G. Knepley     PetscInt aStart, aEnd, a;
3179cbb8f0dSMatthew G. Knepley 
3189cbb8f0dSMatthew G. Knepley     for (d = 0; d < maxDof; ++d) indices[d] = d;
3199566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd));
3209cbb8f0dSMatthew G. Knepley     for (a = aStart; a < aEnd; a++) {
3219cbb8f0dSMatthew G. Knepley       PetscInt dof, f;
3229cbb8f0dSMatthew G. Knepley 
3239566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(aSec, a, &dof));
3249cbb8f0dSMatthew G. Knepley       if (dof) {
3259cbb8f0dSMatthew G. Knepley         /* if there are point-to-point constraints, then all dofs are constrained */
326baef929fSMatthew G. Knepley         for (f = 0; f < Nf; f++) {
3279566063dSJacob Faibussowitsch           PetscCall(PetscSectionSetFieldConstraintIndices(section, a, f, indices));
3289cbb8f0dSMatthew G. Knepley         }
3299cbb8f0dSMatthew G. Knepley       }
3309cbb8f0dSMatthew G. Knepley     }
3319cbb8f0dSMatthew G. Knepley   }
3329566063dSJacob Faibussowitsch   PetscCall(PetscFree(indices));
3339cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
3349cbb8f0dSMatthew G. Knepley }
3359cbb8f0dSMatthew G. Knepley 
3369cbb8f0dSMatthew G. Knepley /* Set the constrained indices on each point */
3379cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCIndices(DM dm, PetscSection section)
3389cbb8f0dSMatthew G. Knepley {
3399cbb8f0dSMatthew G. Knepley   PetscInt      *indices;
340baef929fSMatthew G. Knepley   PetscInt       Nf, maxDof, pStart, pEnd, p, f, d;
3419cbb8f0dSMatthew G. Knepley 
3429cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
3439566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &Nf));
3449566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetMaxDof(section, &maxDof));
3459566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
3469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxDof, &indices));
3479cbb8f0dSMatthew G. Knepley   for (d = 0; d < maxDof; ++d) indices[d] = -1;
3489cbb8f0dSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
3499cbb8f0dSMatthew G. Knepley     PetscInt cdof, d;
3509cbb8f0dSMatthew G. Knepley 
3519566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
3529cbb8f0dSMatthew G. Knepley     if (cdof) {
353baef929fSMatthew G. Knepley       if (Nf) {
3549cbb8f0dSMatthew G. Knepley         PetscInt numConst = 0, foff = 0;
3559cbb8f0dSMatthew G. Knepley 
356baef929fSMatthew G. Knepley         for (f = 0; f < Nf; ++f) {
3579cbb8f0dSMatthew G. Knepley           const PetscInt *find;
3589cbb8f0dSMatthew G. Knepley           PetscInt        fcdof, fdof;
3599cbb8f0dSMatthew G. Knepley 
3609566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldDof(section, p, f, &fdof));
3619566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintDof(section, p, f, &fcdof));
3629cbb8f0dSMatthew G. Knepley           /* Change constraint numbering from field component to local dof number */
3639566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintIndices(section, p, f, &find));
3649cbb8f0dSMatthew G. Knepley           for (d = 0; d < fcdof; ++d) indices[numConst+d] = find[d] + foff;
3659cbb8f0dSMatthew G. Knepley           numConst += fcdof;
3669cbb8f0dSMatthew G. Knepley           foff     += fdof;
3679cbb8f0dSMatthew G. Knepley         }
3689566063dSJacob Faibussowitsch         if (cdof != numConst) PetscCall(PetscSectionSetConstraintDof(section, p, numConst));
3699cbb8f0dSMatthew G. Knepley       } else {
3709cbb8f0dSMatthew G. Knepley         for (d = 0; d < cdof; ++d) indices[d] = d;
3719cbb8f0dSMatthew G. Knepley       }
3729566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetConstraintIndices(section, p, indices));
3739cbb8f0dSMatthew G. Knepley     }
3749cbb8f0dSMatthew G. Knepley   }
3759566063dSJacob Faibussowitsch   PetscCall(PetscFree(indices));
3769cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
3779cbb8f0dSMatthew G. Knepley }
3789cbb8f0dSMatthew G. Knepley 
3799cbb8f0dSMatthew G. Knepley /*@C
3809cbb8f0dSMatthew G. Knepley   DMPlexCreateSection - Create a PetscSection based upon the dof layout specification provided.
3819cbb8f0dSMatthew G. Knepley 
3829cbb8f0dSMatthew G. Knepley   Not Collective
3839cbb8f0dSMatthew G. Knepley 
3849cbb8f0dSMatthew G. Knepley   Input Parameters:
3859cbb8f0dSMatthew G. Knepley + dm        - The DMPlex object
38692cfd99aSMartin Diehl . label     - The label indicating the mesh support of each field, or NULL for the whole mesh
3879cbb8f0dSMatthew G. Knepley . numComp   - An array of size numFields that holds the number of components for each field
3889cbb8f0dSMatthew G. Knepley . numDof    - An array of size numFields*(dim+1) which holds the number of dof for each field on a mesh piece of dimension d
3899cbb8f0dSMatthew G. Knepley . numBC     - The number of boundary conditions
3909cbb8f0dSMatthew G. Knepley . bcField   - An array of size numBC giving the field number for each boundry condition
3919cbb8f0dSMatthew G. Knepley . bcComps   - [Optional] An array of size numBC giving an IS holding the field components to which each boundary condition applies
3929cbb8f0dSMatthew G. Knepley . bcPoints  - An array of size numBC giving an IS holding the Plex points to which each boundary condition applies
3939cbb8f0dSMatthew G. Knepley - perm      - Optional permutation of the chart, or NULL
3949cbb8f0dSMatthew G. Knepley 
3959cbb8f0dSMatthew G. Knepley   Output Parameter:
3969cbb8f0dSMatthew G. Knepley . section - The PetscSection object
3979cbb8f0dSMatthew G. Knepley 
3989cbb8f0dSMatthew G. Knepley   Notes:
3999cbb8f0dSMatthew 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
4009cbb8f0dSMatthew G. Knepley   number of dof for field 0 on each edge.
4019cbb8f0dSMatthew G. Knepley 
4029cbb8f0dSMatthew G. Knepley   The chart permutation is the same one set using PetscSectionSetPermutation()
4039cbb8f0dSMatthew G. Knepley 
4049cbb8f0dSMatthew G. Knepley   Level: developer
4059cbb8f0dSMatthew G. Knepley 
4061bb6d2a8SBarry Smith   TODO: How is this related to DMCreateLocalSection()
4071bb6d2a8SBarry Smith 
4089cbb8f0dSMatthew G. Knepley .seealso: DMPlexCreate(), PetscSectionCreate(), PetscSectionSetPermutation()
4099cbb8f0dSMatthew G. Knepley @*/
410baef929fSMatthew G. Knepley 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)
4119cbb8f0dSMatthew G. Knepley {
4129cbb8f0dSMatthew G. Knepley   PetscSection   aSec;
4139cbb8f0dSMatthew G. Knepley 
4149cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
4159566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSectionFields(dm, numComp, section));
4169566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSectionDof(dm, label, numDof, *section));
4179566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSectionBCDof(dm, numBC, bcField, bcComps, bcPoints, *section));
4189566063dSJacob Faibussowitsch   if (perm) PetscCall(PetscSectionSetPermutation(*section, perm));
4199566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFromOptions(*section));
4209566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(*section));
4219566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAnchors(dm,&aSec,NULL));
4229cbb8f0dSMatthew G. Knepley   if (numBC || aSec) {
4239566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSectionBCIndicesField(dm, numBC, bcField, bcComps, bcPoints, *section));
4249566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSectionBCIndices(dm, *section));
4259cbb8f0dSMatthew G. Knepley   }
4269566063dSJacob Faibussowitsch   PetscCall(PetscSectionViewFromOptions(*section,NULL,"-section_view"));
4279cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
4289cbb8f0dSMatthew G. Knepley }
4299cbb8f0dSMatthew G. Knepley 
4301bb6d2a8SBarry Smith PetscErrorCode DMCreateLocalSection_Plex(DM dm)
4319cbb8f0dSMatthew G. Knepley {
4329cbb8f0dSMatthew G. Knepley   PetscSection   section;
433baef929fSMatthew G. Knepley   DMLabel       *labels;
4349cbb8f0dSMatthew G. Knepley   IS            *bcPoints, *bcComps;
4359cbb8f0dSMatthew G. Knepley   PetscBool     *isFE;
4369cbb8f0dSMatthew G. Knepley   PetscInt      *bcFields, *numComp, *numDof;
4379310035eSMatthew G. Knepley   PetscInt       depth, dim, numBC = 0, Nf, Nds, s, bc = 0, f;
4389cbb8f0dSMatthew G. Knepley   PetscInt       cStart, cEnd, cEndInterior;
4399cbb8f0dSMatthew G. Knepley 
4409cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
4419566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
4429566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
4439566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
4449cbb8f0dSMatthew G. Knepley   /* FE and FV boundary conditions are handled slightly differently */
4459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nf, &isFE));
446baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
4479cbb8f0dSMatthew G. Knepley     PetscObject  obj;
4489cbb8f0dSMatthew G. Knepley     PetscClassId id;
4499cbb8f0dSMatthew G. Knepley 
4509566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, f, NULL, &obj));
4519566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(obj, &id));
4529cbb8f0dSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {isFE[f] = PETSC_TRUE;}
4539cbb8f0dSMatthew G. Knepley     else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;}
45498921bdaSJacob Faibussowitsch     else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f);
4559cbb8f0dSMatthew G. Knepley   }
4569cbb8f0dSMatthew G. Knepley   /* Allocate boundary point storage for FEM boundaries */
4579566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
4589310035eSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
4599310035eSMatthew G. Knepley     PetscDS  dsBC;
4609310035eSMatthew G. Knepley     PetscInt numBd, bd;
4619310035eSMatthew G. Knepley 
4629566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC));
4639566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumBoundary(dsBC, &numBd));
464*08401ef6SPierre 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)");
4659cbb8f0dSMatthew G. Knepley     for (bd = 0; bd < numBd; ++bd) {
4669cbb8f0dSMatthew G. Knepley       PetscInt                field;
4679cbb8f0dSMatthew G. Knepley       DMBoundaryConditionType type;
4689cbb8f0dSMatthew G. Knepley       DMLabel                 label;
4699cbb8f0dSMatthew G. Knepley 
4709566063dSJacob Faibussowitsch       PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, &type, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL));
4719cbb8f0dSMatthew G. Knepley       if (label && isFE[field] && (type & DM_BC_ESSENTIAL)) ++numBC;
4729cbb8f0dSMatthew G. Knepley     }
4739310035eSMatthew G. Knepley   }
4749cbb8f0dSMatthew G. Knepley   /* Add ghost cell boundaries for FVM */
4759566063dSJacob Faibussowitsch   PetscCall(DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL));
476baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) if (!isFE[f] && cEndInterior >= 0) ++numBC;
4779566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(numBC, &bcFields, numBC, &bcPoints, numBC, &bcComps));
4789cbb8f0dSMatthew G. Knepley   /* Constrain ghost cells for FV */
479baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
4809cbb8f0dSMatthew G. Knepley     PetscInt *newidx, c;
4819cbb8f0dSMatthew G. Knepley 
4829cbb8f0dSMatthew G. Knepley     if (isFE[f] || cEndInterior < 0) continue;
4839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(cEnd-cEndInterior,&newidx));
4849cbb8f0dSMatthew G. Knepley     for (c = cEndInterior; c < cEnd; ++c) newidx[c-cEndInterior] = c;
4859cbb8f0dSMatthew G. Knepley     bcFields[bc] = f;
4869566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, cEnd-cEndInterior, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]));
4879cbb8f0dSMatthew G. Knepley   }
4889cbb8f0dSMatthew G. Knepley   /* Handle FEM Dirichlet boundaries */
4899566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
4909310035eSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
4919310035eSMatthew G. Knepley     PetscDS  dsBC;
4929310035eSMatthew G. Knepley     PetscInt numBd, bd;
4939310035eSMatthew G. Knepley 
4949566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC));
4959566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumBoundary(dsBC, &numBd));
4969cbb8f0dSMatthew G. Knepley     for (bd = 0; bd < numBd; ++bd) {
4979cbb8f0dSMatthew G. Knepley       DMLabel                 label;
4989cbb8f0dSMatthew G. Knepley       const PetscInt         *comps;
4999cbb8f0dSMatthew G. Knepley       const PetscInt         *values;
5009310035eSMatthew G. Knepley       PetscInt                bd2, field, numComps, numValues;
5019cbb8f0dSMatthew G. Knepley       DMBoundaryConditionType type;
5029310035eSMatthew G. Knepley       PetscBool               duplicate = PETSC_FALSE;
5039cbb8f0dSMatthew G. Knepley 
5049566063dSJacob Faibussowitsch       PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, &type, NULL, &label, &numValues, &values, &field, &numComps, &comps, NULL, NULL, NULL));
5059cbb8f0dSMatthew G. Knepley       if (!isFE[field] || !label) continue;
5069310035eSMatthew G. Knepley       /* Only want to modify label once */
5079310035eSMatthew G. Knepley       for (bd2 = 0; bd2 < bd; ++bd2) {
50845480ffeSMatthew G. Knepley         DMLabel l;
50945480ffeSMatthew G. Knepley 
5109566063dSJacob Faibussowitsch         PetscCall(PetscDSGetBoundary(dsBC, bd2, NULL, NULL, NULL, &l, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
51145480ffeSMatthew G. Knepley         duplicate = l == label ? PETSC_TRUE : PETSC_FALSE;
5129310035eSMatthew G. Knepley         if (duplicate) break;
5139310035eSMatthew G. Knepley       }
5149310035eSMatthew G. Knepley       if (!duplicate && (isFE[field])) {
5159310035eSMatthew G. Knepley         /* don't complete cells, which are just present to give orientation to the boundary */
5169566063dSJacob Faibussowitsch         PetscCall(DMPlexLabelComplete(dm, label));
5179310035eSMatthew G. Knepley       }
5189cbb8f0dSMatthew G. Knepley       /* Filter out cells, if you actually want to constrain cells you need to do things by hand right now */
5199cbb8f0dSMatthew G. Knepley       if (type & DM_BC_ESSENTIAL) {
5209cbb8f0dSMatthew G. Knepley         PetscInt       *newidx;
5219cbb8f0dSMatthew G. Knepley         PetscInt        n, newn = 0, p, v;
5229cbb8f0dSMatthew G. Knepley 
5239cbb8f0dSMatthew G. Knepley         bcFields[bc] = field;
5249566063dSJacob Faibussowitsch         if (numComps) PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject) dm), numComps, comps, PETSC_COPY_VALUES, &bcComps[bc]));
5259cbb8f0dSMatthew G. Knepley         for (v = 0; v < numValues; ++v) {
5269cbb8f0dSMatthew G. Knepley           IS              tmp;
5279cbb8f0dSMatthew G. Knepley           const PetscInt *idx;
5289cbb8f0dSMatthew G. Knepley 
5299566063dSJacob Faibussowitsch           PetscCall(DMLabelGetStratumIS(label, values[v], &tmp));
5309cbb8f0dSMatthew G. Knepley           if (!tmp) continue;
5319566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(tmp, &n));
5329566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(tmp, &idx));
5339cbb8f0dSMatthew G. Knepley           if (isFE[field]) {
5349cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) ++newn;
5359cbb8f0dSMatthew G. Knepley           } else {
5369cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) ++newn;
5379cbb8f0dSMatthew G. Knepley           }
5389566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(tmp, &idx));
5399566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&tmp));
5409cbb8f0dSMatthew G. Knepley         }
5419566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(newn, &newidx));
5429cbb8f0dSMatthew G. Knepley         newn = 0;
5439cbb8f0dSMatthew G. Knepley         for (v = 0; v < numValues; ++v) {
5449cbb8f0dSMatthew G. Knepley           IS              tmp;
5459cbb8f0dSMatthew G. Knepley           const PetscInt *idx;
5469cbb8f0dSMatthew G. Knepley 
5479566063dSJacob Faibussowitsch           PetscCall(DMLabelGetStratumIS(label, values[v], &tmp));
5489cbb8f0dSMatthew G. Knepley           if (!tmp) continue;
5499566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(tmp, &n));
5509566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(tmp, &idx));
5519cbb8f0dSMatthew G. Knepley           if (isFE[field]) {
5529cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) newidx[newn++] = idx[p];
5539cbb8f0dSMatthew G. Knepley           } else {
5549cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) newidx[newn++] = idx[p];
5559cbb8f0dSMatthew G. Knepley           }
5569566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(tmp, &idx));
5579566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&tmp));
5589cbb8f0dSMatthew G. Knepley         }
5599566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, newn, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]));
5609310035eSMatthew G. Knepley       }
5619cbb8f0dSMatthew G. Knepley     }
5629cbb8f0dSMatthew G. Knepley   }
5639cbb8f0dSMatthew G. Knepley   /* Handle discretization */
5649566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(Nf,&labels,Nf,&numComp,Nf*(dim+1),&numDof));
565baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
56644a7f3ddSMatthew G. Knepley     labels[f] = dm->fields[f].label;
5679cbb8f0dSMatthew G. Knepley     if (isFE[f]) {
56844a7f3ddSMatthew G. Knepley       PetscFE         fe = (PetscFE) dm->fields[f].disc;
5699cbb8f0dSMatthew G. Knepley       const PetscInt *numFieldDof;
570a3254110SMatthew Knepley       PetscInt        fedim, d;
5719cbb8f0dSMatthew G. Knepley 
5729566063dSJacob Faibussowitsch       PetscCall(PetscFEGetNumComponents(fe, &numComp[f]));
5739566063dSJacob Faibussowitsch       PetscCall(PetscFEGetNumDof(fe, &numFieldDof));
5749566063dSJacob Faibussowitsch       PetscCall(PetscFEGetSpatialDimension(fe, &fedim));
575a3254110SMatthew Knepley       for (d = 0; d < PetscMin(dim, fedim)+1; ++d) numDof[f*(dim+1)+d] = numFieldDof[d];
5769cbb8f0dSMatthew G. Knepley     } else {
57744a7f3ddSMatthew G. Knepley       PetscFV fv = (PetscFV) dm->fields[f].disc;
5789cbb8f0dSMatthew G. Knepley 
5799566063dSJacob Faibussowitsch       PetscCall(PetscFVGetNumComponents(fv, &numComp[f]));
5809cbb8f0dSMatthew G. Knepley       numDof[f*(dim+1)+dim] = numComp[f];
5819cbb8f0dSMatthew G. Knepley     }
5829cbb8f0dSMatthew G. Knepley   }
5839566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
584baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
5859cbb8f0dSMatthew G. Knepley     PetscInt d;
5869cbb8f0dSMatthew G. Knepley     for (d = 1; d < dim; ++d) {
5872c71b3e2SJacob Faibussowitsch       PetscCheckFalse((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.");
5889cbb8f0dSMatthew G. Knepley     }
5899cbb8f0dSMatthew G. Knepley   }
5909566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSection(dm, labels, numComp, numDof, numBC, bcFields, bcComps, bcPoints, NULL, &section));
591baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
5929cbb8f0dSMatthew G. Knepley     PetscFE     fe;
5939cbb8f0dSMatthew G. Knepley     const char *name;
5949cbb8f0dSMatthew G. Knepley 
5959566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, f, NULL, (PetscObject *) &fe));
596cdfe53dfSMatthew G. Knepley     if (!((PetscObject) fe)->name) continue;
5979566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject) fe, &name));
5989566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldName(section, f, name));
5999cbb8f0dSMatthew G. Knepley   }
6009566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, section));
6019566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&section));
6029566063dSJacob Faibussowitsch   for (bc = 0; bc < numBC; ++bc) {PetscCall(ISDestroy(&bcPoints[bc]));PetscCall(ISDestroy(&bcComps[bc]));}
6039566063dSJacob Faibussowitsch   PetscCall(PetscFree3(bcFields,bcPoints,bcComps));
6049566063dSJacob Faibussowitsch   PetscCall(PetscFree3(labels,numComp,numDof));
6059566063dSJacob Faibussowitsch   PetscCall(PetscFree(isFE));
6069cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
6079cbb8f0dSMatthew G. Knepley }
608