xref: /petsc/src/dm/impls/plex/plexsection.c (revision 5fedec97950c19de564efaecd0f125b1a6cb2b20)
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   PetscErrorCode ierr;
109cbb8f0dSMatthew G. Knepley 
119cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
12baef929fSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
13baef929fSMatthew G. Knepley   ierr = DMPlexGetDepthLabel(dm,&depthLabel);CHKERRQ(ierr);
149cbb8f0dSMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
15baef929fSMatthew G. Knepley   ierr = PetscCalloc1(Nf, &isFE);CHKERRQ(ierr);
16baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
179cbb8f0dSMatthew G. Knepley     PetscObject  obj;
189cbb8f0dSMatthew G. Knepley     PetscClassId id;
199cbb8f0dSMatthew G. Knepley 
2044a7f3ddSMatthew G. Knepley     ierr = DMGetField(dm, f, NULL, &obj);CHKERRQ(ierr);
219cbb8f0dSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
229cbb8f0dSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {isFE[f] = PETSC_TRUE;}
239cbb8f0dSMatthew G. Knepley     else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;}
249cbb8f0dSMatthew G. Knepley   }
25baef929fSMatthew G. Knepley 
269cbb8f0dSMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), section);CHKERRQ(ierr);
27baef929fSMatthew G. Knepley   if (Nf > 0) {
28baef929fSMatthew G. Knepley     ierr = PetscSectionSetNumFields(*section, Nf);CHKERRQ(ierr);
299cbb8f0dSMatthew G. Knepley     if (numComp) {
30baef929fSMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
319cbb8f0dSMatthew G. Knepley         ierr = PetscSectionSetFieldComponents(*section, f, numComp[f]);CHKERRQ(ierr);
329cbb8f0dSMatthew G. Knepley         if (isFE[f]) {
339cbb8f0dSMatthew G. Knepley           PetscFE           fe;
349cbb8f0dSMatthew G. Knepley           PetscDualSpace    dspace;
359cbb8f0dSMatthew G. Knepley           const PetscInt    ***perms;
369cbb8f0dSMatthew G. Knepley           const PetscScalar ***flips;
379cbb8f0dSMatthew G. Knepley           const PetscInt    *numDof;
389cbb8f0dSMatthew G. Knepley 
3944a7f3ddSMatthew G. Knepley           ierr = DMGetField(dm,f,NULL,(PetscObject *) &fe);CHKERRQ(ierr);
409cbb8f0dSMatthew G. Knepley           ierr = PetscFEGetDualSpace(fe,&dspace);CHKERRQ(ierr);
419cbb8f0dSMatthew G. Knepley           ierr = PetscDualSpaceGetSymmetries(dspace,&perms,&flips);CHKERRQ(ierr);
429cbb8f0dSMatthew G. Knepley           ierr = PetscDualSpaceGetNumDof(dspace,&numDof);CHKERRQ(ierr);
439cbb8f0dSMatthew G. Knepley           if (perms || flips) {
449cbb8f0dSMatthew G. Knepley             DM              K;
454dff28b8SMatthew G. Knepley             PetscInt        sph, spdepth;
469cbb8f0dSMatthew G. Knepley             PetscSectionSym sym;
479cbb8f0dSMatthew G. Knepley 
489cbb8f0dSMatthew G. Knepley             ierr = PetscDualSpaceGetDM(dspace,&K);CHKERRQ(ierr);
494dff28b8SMatthew G. Knepley             ierr = DMPlexGetDepth(K, &spdepth);CHKERRQ(ierr);
509cbb8f0dSMatthew G. Knepley             ierr = PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)*section),depthLabel,&sym);CHKERRQ(ierr);
514dff28b8SMatthew G. Knepley             for (sph = 0; sph <= spdepth; sph++) {
529cbb8f0dSMatthew G. Knepley               PetscDualSpace    hspace;
539cbb8f0dSMatthew G. Knepley               PetscInt          kStart, kEnd;
544dff28b8SMatthew G. Knepley               PetscInt          kConeSize, h = sph + (depth - spdepth);
559cbb8f0dSMatthew G. Knepley               const PetscInt    **perms0 = NULL;
569cbb8f0dSMatthew G. Knepley               const PetscScalar **flips0 = NULL;
579cbb8f0dSMatthew G. Knepley 
584dff28b8SMatthew G. Knepley               ierr = PetscDualSpaceGetHeightSubspace(dspace, sph, &hspace);CHKERRQ(ierr);
599cbb8f0dSMatthew G. Knepley               ierr = DMPlexGetHeightStratum(K, h, &kStart, &kEnd);CHKERRQ(ierr);
609cbb8f0dSMatthew G. Knepley               if (!hspace) continue;
619cbb8f0dSMatthew G. Knepley               ierr = PetscDualSpaceGetSymmetries(hspace,&perms,&flips);CHKERRQ(ierr);
629cbb8f0dSMatthew G. Knepley               if (perms) perms0 = perms[0];
639cbb8f0dSMatthew G. Knepley               if (flips) flips0 = flips[0];
649cbb8f0dSMatthew G. Knepley               if (!(perms0 || flips0)) continue;
65b5a892a1SMatthew G. Knepley               {
66b5a892a1SMatthew G. Knepley                 DMPolytopeType ct;
67b5a892a1SMatthew G. Knepley                 /* The number of arrangements is no longer based on the number of faces */
68b5a892a1SMatthew G. Knepley                 ierr = DMPlexGetCellType(K, kStart, &ct);CHKERRQ(ierr);
69b5a892a1SMatthew G. Knepley                 kConeSize = DMPolytopeTypeGetNumArrangments(ct) / 2;
70b5a892a1SMatthew G. Knepley               }
719cbb8f0dSMatthew G. Knepley               ierr = PetscSectionSymLabelSetStratum(sym,depth - h,numDof[depth - h],-kConeSize,kConeSize,PETSC_USE_POINTER,perms0 ? &perms0[-kConeSize] : NULL,flips0 ? &flips0[-kConeSize] : NULL);CHKERRQ(ierr);
729cbb8f0dSMatthew G. Knepley             }
739cbb8f0dSMatthew G. Knepley             ierr = PetscSectionSetFieldSym(*section,f,sym);CHKERRQ(ierr);
749cbb8f0dSMatthew G. Knepley             ierr = PetscSectionSymDestroy(&sym);CHKERRQ(ierr);
759cbb8f0dSMatthew G. Knepley           }
769cbb8f0dSMatthew G. Knepley         }
779cbb8f0dSMatthew G. Knepley       }
789cbb8f0dSMatthew G. Knepley     }
799cbb8f0dSMatthew G. Knepley   }
809cbb8f0dSMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
819cbb8f0dSMatthew G. Knepley   ierr = PetscSectionSetChart(*section, pStart, pEnd);CHKERRQ(ierr);
82baef929fSMatthew G. Knepley   ierr = PetscFree(isFE);CHKERRQ(ierr);
83baef929fSMatthew G. Knepley   PetscFunctionReturn(0);
84baef929fSMatthew G. Knepley }
85baef929fSMatthew G. Knepley 
86baef929fSMatthew G. Knepley /* Set the number of dof on each point and separate by fields */
87baef929fSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionDof(DM dm, DMLabel label[],const PetscInt numDof[], PetscSection section)
88baef929fSMatthew G. Knepley {
89baef929fSMatthew G. Knepley   DMLabel        depthLabel;
90412e9a14SMatthew G. Knepley   DMPolytopeType ct;
91baef929fSMatthew G. Knepley   PetscInt       depth, cellHeight, pStart = 0, pEnd = 0;
92a55f9a55SMatthew G. Knepley   PetscInt       Nf, f, Nds, n, dim, d, dep, p;
93*5fedec97SMatthew G. Knepley   PetscBool     *isFE, hasCohesive = PETSC_FALSE;
94baef929fSMatthew G. Knepley   PetscErrorCode ierr;
95baef929fSMatthew G. Knepley 
96baef929fSMatthew G. Knepley   PetscFunctionBegin;
97baef929fSMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
989cbb8f0dSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
99baef929fSMatthew G. Knepley   ierr = DMPlexGetDepthLabel(dm,&depthLabel);CHKERRQ(ierr);
100baef929fSMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
101a55f9a55SMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
102a55f9a55SMatthew G. Knepley   for (n = 0; n < Nds; ++n) {
103a55f9a55SMatthew G. Knepley     PetscDS   ds;
104*5fedec97SMatthew G. Knepley     PetscBool isCohesive;
105a55f9a55SMatthew G. Knepley 
106a55f9a55SMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, n, NULL, NULL, &ds);CHKERRQ(ierr);
107*5fedec97SMatthew G. Knepley     ierr = PetscDSIsCohesive(ds, &isCohesive);CHKERRQ(ierr);
108*5fedec97SMatthew G. Knepley     if (isCohesive) {hasCohesive = PETSC_TRUE; break;}
109a55f9a55SMatthew G. Knepley   }
110baef929fSMatthew G. Knepley   ierr = PetscMalloc1(Nf, &isFE);CHKERRQ(ierr);
111baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
112baef929fSMatthew G. Knepley     PetscObject  obj;
113baef929fSMatthew G. Knepley     PetscClassId id;
114baef929fSMatthew G. Knepley 
11544a7f3ddSMatthew G. Knepley     ierr = DMGetField(dm, f, NULL, &obj);CHKERRQ(ierr);
116baef929fSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1171274e1a1SLawrence Mitchell     /* User is allowed to put a "placeholder" field in (c.f. DMCreateDS) */
1181274e1a1SLawrence Mitchell     isFE[f] = id == PETSCFE_CLASSID ? PETSC_TRUE : PETSC_FALSE;
119baef929fSMatthew G. Knepley   }
120baef929fSMatthew G. Knepley 
1219cbb8f0dSMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
122baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
123e0b68406SMatthew Knepley     PetscBool avoidTensor;
124e0b68406SMatthew Knepley 
125e0b68406SMatthew Knepley     ierr = DMGetFieldAvoidTensor(dm, f, &avoidTensor);CHKERRQ(ierr);
126*5fedec97SMatthew G. Knepley     avoidTensor = (avoidTensor || hasCohesive) ? PETSC_TRUE : PETSC_FALSE;
127baef929fSMatthew G. Knepley     if (label && label[f]) {
128baef929fSMatthew G. Knepley       IS              pointIS;
129baef929fSMatthew G. Knepley       const PetscInt *points;
130baef929fSMatthew G. Knepley       PetscInt        n;
131baef929fSMatthew G. Knepley 
132baef929fSMatthew G. Knepley       ierr = DMLabelGetStratumIS(label[f], 1, &pointIS);CHKERRQ(ierr);
133baef929fSMatthew G. Knepley       if (!pointIS) continue;
134baef929fSMatthew G. Knepley       ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
135baef929fSMatthew G. Knepley       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
136baef929fSMatthew G. Knepley       for (p = 0; p < n; ++p) {
137baef929fSMatthew G. Knepley         const PetscInt point = points[p];
138baef929fSMatthew G. Knepley         PetscInt       dof, d;
139baef929fSMatthew G. Knepley 
140412e9a14SMatthew G. Knepley         ierr = DMPlexGetCellType(dm, point, &ct);CHKERRQ(ierr);
141baef929fSMatthew G. Knepley         ierr = DMLabelGetValue(depthLabel, point, &d);CHKERRQ(ierr);
142412e9a14SMatthew G. Knepley         /* If this is a tensor prism point, use dof for one dimension lower */
143412e9a14SMatthew G. Knepley         switch (ct) {
144412e9a14SMatthew G. Knepley           case DM_POLYTOPE_POINT_PRISM_TENSOR:
145412e9a14SMatthew G. Knepley           case DM_POLYTOPE_SEG_PRISM_TENSOR:
146412e9a14SMatthew G. Knepley           case DM_POLYTOPE_TRI_PRISM_TENSOR:
147412e9a14SMatthew G. Knepley           case DM_POLYTOPE_QUAD_PRISM_TENSOR:
148*5fedec97SMatthew G. Knepley             if (hasCohesive) {--d;} break;
149412e9a14SMatthew G. Knepley           default: break;
150412e9a14SMatthew G. Knepley         }
151baef929fSMatthew G. Knepley         dof  = d < 0 ? 0 : numDof[f*(dim+1)+d];
152baef929fSMatthew G. Knepley         ierr = PetscSectionSetFieldDof(section, point, f, dof);CHKERRQ(ierr);
153baef929fSMatthew G. Knepley         ierr = PetscSectionAddDof(section, point, dof);CHKERRQ(ierr);
154baef929fSMatthew G. Knepley       }
155baef929fSMatthew G. Knepley       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
156baef929fSMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
157baef929fSMatthew G. Knepley     } else {
1589cbb8f0dSMatthew G. Knepley       for (dep = 0; dep <= depth - cellHeight; ++dep) {
1593fd2e703SMatthew G. Knepley         /* Cases: dim > depth (cell-vertex mesh), dim == depth (fully interpolated), dim < depth (interpolated submesh) */
1603fd2e703SMatthew G. Knepley         d    = dim <= depth ? dep : (!dep ? 0 : dim);
1619cbb8f0dSMatthew G. Knepley         ierr = DMPlexGetDepthStratum(dm, dep, &pStart, &pEnd);CHKERRQ(ierr);
1629cbb8f0dSMatthew G. Knepley         for (p = pStart; p < pEnd; ++p) {
163baef929fSMatthew G. Knepley           const PetscInt dof = numDof[f*(dim+1)+d];
1649cbb8f0dSMatthew G. Knepley 
165412e9a14SMatthew G. Knepley           ierr = DMPlexGetCellType(dm, p, &ct);CHKERRQ(ierr);
166412e9a14SMatthew G. Knepley           switch (ct) {
167412e9a14SMatthew G. Knepley             case DM_POLYTOPE_POINT_PRISM_TENSOR:
168412e9a14SMatthew G. Knepley             case DM_POLYTOPE_SEG_PRISM_TENSOR:
169412e9a14SMatthew G. Knepley             case DM_POLYTOPE_TRI_PRISM_TENSOR:
170412e9a14SMatthew G. Knepley             case DM_POLYTOPE_QUAD_PRISM_TENSOR:
171e0b68406SMatthew Knepley               if (avoidTensor && isFE[f]) continue;
172412e9a14SMatthew G. Knepley             default: break;
173412e9a14SMatthew G. Knepley           }
174baef929fSMatthew G. Knepley           ierr = PetscSectionSetFieldDof(section, p, f, dof);CHKERRQ(ierr);
175baef929fSMatthew G. Knepley           ierr = PetscSectionAddDof(section, p, dof);CHKERRQ(ierr);
1769cbb8f0dSMatthew G. Knepley         }
177baef929fSMatthew G. Knepley       }
1789cbb8f0dSMatthew G. Knepley     }
1799cbb8f0dSMatthew G. Knepley   }
1809cbb8f0dSMatthew G. Knepley   ierr = PetscFree(isFE);CHKERRQ(ierr);
1819cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
1829cbb8f0dSMatthew G. Knepley }
1839cbb8f0dSMatthew G. Knepley 
1849cbb8f0dSMatthew G. Knepley /* Set the number of dof on each point and separate by fields
1859cbb8f0dSMatthew G. Knepley    If bcComps is NULL or the IS is NULL, constrain every dof on the point
1869cbb8f0dSMatthew G. Knepley */
1879cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCDof(DM dm, PetscInt numBC, const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section)
1889cbb8f0dSMatthew G. Knepley {
189baef929fSMatthew G. Knepley   PetscInt       Nf;
1909cbb8f0dSMatthew G. Knepley   PetscInt       bc;
1919cbb8f0dSMatthew G. Knepley   PetscSection   aSec;
1929cbb8f0dSMatthew G. Knepley   PetscErrorCode ierr;
1939cbb8f0dSMatthew G. Knepley 
1949cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
195baef929fSMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1969cbb8f0dSMatthew G. Knepley   for (bc = 0; bc < numBC; ++bc) {
1979cbb8f0dSMatthew G. Knepley     PetscInt        field = 0;
1989cbb8f0dSMatthew G. Knepley     const PetscInt *comp;
1999cbb8f0dSMatthew G. Knepley     const PetscInt *idx;
200ce093827SMatthew G. Knepley     PetscInt        Nc = 0, cNc = -1, n, i;
2019cbb8f0dSMatthew G. Knepley 
202ce093827SMatthew G. Knepley     if (Nf) {
203ce093827SMatthew G. Knepley       field = bcField[bc];
204ce093827SMatthew G. Knepley       ierr = PetscSectionGetFieldComponents(section, field, &Nc);CHKERRQ(ierr);
205ce093827SMatthew G. Knepley     }
206ce093827SMatthew G. Knepley     if (bcComps && bcComps[bc]) {ierr = ISGetLocalSize(bcComps[bc], &cNc);CHKERRQ(ierr);}
2079cbb8f0dSMatthew G. Knepley     if (bcComps && bcComps[bc]) {ierr = ISGetIndices(bcComps[bc], &comp);CHKERRQ(ierr);}
2089cbb8f0dSMatthew G. Knepley     ierr = ISGetLocalSize(bcPoints[bc], &n);CHKERRQ(ierr);
2099cbb8f0dSMatthew G. Knepley     ierr = ISGetIndices(bcPoints[bc], &idx);CHKERRQ(ierr);
2109cbb8f0dSMatthew G. Knepley     for (i = 0; i < n; ++i) {
2119cbb8f0dSMatthew G. Knepley       const PetscInt p = idx[i];
2129cbb8f0dSMatthew G. Knepley       PetscInt       numConst;
2139cbb8f0dSMatthew G. Knepley 
214baef929fSMatthew G. Knepley       if (Nf) {
2159cbb8f0dSMatthew G. Knepley         ierr = PetscSectionGetFieldDof(section, p, field, &numConst);CHKERRQ(ierr);
2169cbb8f0dSMatthew G. Knepley       } else {
2179cbb8f0dSMatthew G. Knepley         ierr = PetscSectionGetDof(section, p, &numConst);CHKERRQ(ierr);
2189cbb8f0dSMatthew G. Knepley       }
219ce093827SMatthew G. Knepley       /* If Nc <= 0, constrain every dof on the point */
220ce093827SMatthew G. Knepley       if (cNc > 0) {
221ce093827SMatthew G. Knepley         /* We assume that a point may have multiple "nodes", which are collections of Nc dofs,
222ce093827SMatthew G. Knepley            and that those dofs are numbered n*Nc+c */
223ce093827SMatthew G. Knepley         if (Nf) {
224ce093827SMatthew G. Knepley           if (numConst % Nc) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %D has %D dof which is not divisible by %D field components", p, numConst, Nc);
225ce093827SMatthew G. Knepley           numConst = (numConst/Nc) * cNc;
226ce093827SMatthew G. Knepley         } else {
227ce093827SMatthew G. Knepley           numConst = PetscMin(numConst, cNc);
228ce093827SMatthew G. Knepley         }
229ce093827SMatthew G. Knepley       }
230baef929fSMatthew G. Knepley       if (Nf) {ierr = PetscSectionAddFieldConstraintDof(section, p, field, numConst);CHKERRQ(ierr);}
2319cbb8f0dSMatthew G. Knepley       ierr = PetscSectionAddConstraintDof(section, p, numConst);CHKERRQ(ierr);
2329cbb8f0dSMatthew G. Knepley     }
2339cbb8f0dSMatthew G. Knepley     ierr = ISRestoreIndices(bcPoints[bc], &idx);CHKERRQ(ierr);
2349cbb8f0dSMatthew G. Knepley     if (bcComps && bcComps[bc]) {ierr = ISRestoreIndices(bcComps[bc], &comp);CHKERRQ(ierr);}
2359cbb8f0dSMatthew G. Knepley   }
2369cbb8f0dSMatthew G. Knepley   ierr = DMPlexGetAnchors(dm, &aSec, NULL);CHKERRQ(ierr);
2379cbb8f0dSMatthew G. Knepley   if (aSec) {
2389cbb8f0dSMatthew G. Knepley     PetscInt aStart, aEnd, a;
2399cbb8f0dSMatthew G. Knepley 
2409cbb8f0dSMatthew G. Knepley     ierr = PetscSectionGetChart(aSec, &aStart, &aEnd);CHKERRQ(ierr);
2419cbb8f0dSMatthew G. Knepley     for (a = aStart; a < aEnd; a++) {
2429cbb8f0dSMatthew G. Knepley       PetscInt dof, f;
2439cbb8f0dSMatthew G. Knepley 
2449cbb8f0dSMatthew G. Knepley       ierr = PetscSectionGetDof(aSec, a, &dof);CHKERRQ(ierr);
2459cbb8f0dSMatthew G. Knepley       if (dof) {
2469cbb8f0dSMatthew G. Knepley         /* if there are point-to-point constraints, then all dofs are constrained */
2479cbb8f0dSMatthew G. Knepley         ierr = PetscSectionGetDof(section, a, &dof);CHKERRQ(ierr);
2489cbb8f0dSMatthew G. Knepley         ierr = PetscSectionSetConstraintDof(section, a, dof);CHKERRQ(ierr);
249baef929fSMatthew G. Knepley         for (f = 0; f < Nf; f++) {
2509cbb8f0dSMatthew G. Knepley           ierr = PetscSectionGetFieldDof(section, a, f, &dof);CHKERRQ(ierr);
2519cbb8f0dSMatthew G. Knepley           ierr = PetscSectionSetFieldConstraintDof(section, a, f, dof);CHKERRQ(ierr);
2529cbb8f0dSMatthew G. Knepley         }
2539cbb8f0dSMatthew G. Knepley       }
2549cbb8f0dSMatthew G. Knepley     }
2559cbb8f0dSMatthew G. Knepley   }
2569cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
2579cbb8f0dSMatthew G. Knepley }
2589cbb8f0dSMatthew G. Knepley 
2599cbb8f0dSMatthew G. Knepley /* Set the constrained field indices on each point
2609cbb8f0dSMatthew G. Knepley    If bcComps is NULL or the IS is NULL, constrain every dof on the point
2619cbb8f0dSMatthew G. Knepley */
2629cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCIndicesField(DM dm, PetscInt numBC,const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section)
2639cbb8f0dSMatthew G. Knepley {
2649cbb8f0dSMatthew G. Knepley   PetscSection   aSec;
2659cbb8f0dSMatthew G. Knepley   PetscInt      *indices;
266baef929fSMatthew G. Knepley   PetscInt       Nf, cdof, maxDof = 0, pStart, pEnd, p, bc, f, d;
2679cbb8f0dSMatthew G. Knepley   PetscErrorCode ierr;
2689cbb8f0dSMatthew G. Knepley 
2699cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
270baef929fSMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
271baef929fSMatthew G. Knepley   if (!Nf) PetscFunctionReturn(0);
2729cbb8f0dSMatthew G. Knepley   /* Initialize all field indices to -1 */
2739cbb8f0dSMatthew G. Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
2749cbb8f0dSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); maxDof = PetscMax(maxDof, cdof);}
2759cbb8f0dSMatthew G. Knepley   ierr = PetscMalloc1(maxDof, &indices);CHKERRQ(ierr);
2769cbb8f0dSMatthew G. Knepley   for (d = 0; d < maxDof; ++d) indices[d] = -1;
277baef929fSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) for (f = 0; f < Nf; ++f) {ierr = PetscSectionSetFieldConstraintIndices(section, p, f, indices);CHKERRQ(ierr);}
2789cbb8f0dSMatthew G. Knepley   /* Handle BC constraints */
2799cbb8f0dSMatthew G. Knepley   for (bc = 0; bc < numBC; ++bc) {
2809cbb8f0dSMatthew G. Knepley     const PetscInt  field = bcField[bc];
2819cbb8f0dSMatthew G. Knepley     const PetscInt *comp, *idx;
282ce093827SMatthew G. Knepley     PetscInt        Nc, cNc = -1, n, i;
2839cbb8f0dSMatthew G. Knepley 
284ce093827SMatthew G. Knepley     ierr = PetscSectionGetFieldComponents(section, field, &Nc);CHKERRQ(ierr);
285ce093827SMatthew G. Knepley     if (bcComps && bcComps[bc]) {ierr = ISGetLocalSize(bcComps[bc], &cNc);CHKERRQ(ierr);}
2869cbb8f0dSMatthew G. Knepley     if (bcComps && bcComps[bc]) {ierr = ISGetIndices(bcComps[bc], &comp);CHKERRQ(ierr);}
2879cbb8f0dSMatthew G. Knepley     ierr = ISGetLocalSize(bcPoints[bc], &n);CHKERRQ(ierr);
2889cbb8f0dSMatthew G. Knepley     ierr = ISGetIndices(bcPoints[bc], &idx);CHKERRQ(ierr);
2899cbb8f0dSMatthew G. Knepley     for (i = 0; i < n; ++i) {
2909cbb8f0dSMatthew G. Knepley       const PetscInt  p = idx[i];
2919cbb8f0dSMatthew G. Knepley       const PetscInt *find;
292ce093827SMatthew G. Knepley       PetscInt        fdof, fcdof, c, j;
2939cbb8f0dSMatthew G. Knepley 
2949cbb8f0dSMatthew G. Knepley       ierr = PetscSectionGetFieldDof(section, p, field, &fdof);CHKERRQ(ierr);
2959cbb8f0dSMatthew G. Knepley       if (!fdof) continue;
296ce093827SMatthew G. Knepley       if (cNc < 0) {
2979cbb8f0dSMatthew G. Knepley         for (d = 0; d < fdof; ++d) indices[d] = d;
2989cbb8f0dSMatthew G. Knepley         fcdof = fdof;
2999cbb8f0dSMatthew G. Knepley       } else {
300ce093827SMatthew G. Knepley         /* We assume that a point may have multiple "nodes", which are collections of Nc dofs,
301ce093827SMatthew G. Knepley            and that those dofs are numbered n*Nc+c */
3029cbb8f0dSMatthew G. Knepley         ierr = PetscSectionGetFieldConstraintDof(section, p, field, &fcdof);CHKERRQ(ierr);
3039cbb8f0dSMatthew G. Knepley         ierr = PetscSectionGetFieldConstraintIndices(section, p, field, &find);CHKERRQ(ierr);
304ce093827SMatthew G. Knepley         /* Get indices constrained by previous bcs */
3059cbb8f0dSMatthew G. Knepley         for (d = 0; d < fcdof; ++d) {if (find[d] < 0) break; indices[d] = find[d];}
306ce093827SMatthew G. Knepley         for (j = 0; j < fdof/Nc; ++j) for (c = 0; c < cNc; ++c) indices[d++] = j*Nc + comp[c];
3079cbb8f0dSMatthew G. Knepley         ierr = PetscSortRemoveDupsInt(&d, indices);CHKERRQ(ierr);
3089cbb8f0dSMatthew G. Knepley         for (c = d; c < fcdof; ++c) indices[c] = -1;
3099cbb8f0dSMatthew G. Knepley         fcdof = d;
3109cbb8f0dSMatthew G. Knepley       }
3119cbb8f0dSMatthew G. Knepley       ierr = PetscSectionSetFieldConstraintDof(section, p, field, fcdof);CHKERRQ(ierr);
3129cbb8f0dSMatthew G. Knepley       ierr = PetscSectionSetFieldConstraintIndices(section, p, field, indices);CHKERRQ(ierr);
3139cbb8f0dSMatthew G. Knepley     }
3149cbb8f0dSMatthew G. Knepley     if (bcComps && bcComps[bc]) {ierr = ISRestoreIndices(bcComps[bc], &comp);CHKERRQ(ierr);}
3159cbb8f0dSMatthew G. Knepley     ierr = ISRestoreIndices(bcPoints[bc], &idx);CHKERRQ(ierr);
3169cbb8f0dSMatthew G. Knepley   }
3179cbb8f0dSMatthew G. Knepley   /* Handle anchors */
3189cbb8f0dSMatthew G. Knepley   ierr = DMPlexGetAnchors(dm, &aSec, NULL);CHKERRQ(ierr);
3199cbb8f0dSMatthew G. Knepley   if (aSec) {
3209cbb8f0dSMatthew G. Knepley     PetscInt aStart, aEnd, a;
3219cbb8f0dSMatthew G. Knepley 
3229cbb8f0dSMatthew G. Knepley     for (d = 0; d < maxDof; ++d) indices[d] = d;
3239cbb8f0dSMatthew G. Knepley     ierr = PetscSectionGetChart(aSec, &aStart, &aEnd);CHKERRQ(ierr);
3249cbb8f0dSMatthew G. Knepley     for (a = aStart; a < aEnd; a++) {
3259cbb8f0dSMatthew G. Knepley       PetscInt dof, f;
3269cbb8f0dSMatthew G. Knepley 
3279cbb8f0dSMatthew G. Knepley       ierr = PetscSectionGetDof(aSec, a, &dof);CHKERRQ(ierr);
3289cbb8f0dSMatthew G. Knepley       if (dof) {
3299cbb8f0dSMatthew G. Knepley         /* if there are point-to-point constraints, then all dofs are constrained */
330baef929fSMatthew G. Knepley         for (f = 0; f < Nf; f++) {
3319cbb8f0dSMatthew G. Knepley           ierr = PetscSectionSetFieldConstraintIndices(section, a, f, indices);CHKERRQ(ierr);
3329cbb8f0dSMatthew G. Knepley         }
3339cbb8f0dSMatthew G. Knepley       }
3349cbb8f0dSMatthew G. Knepley     }
3359cbb8f0dSMatthew G. Knepley   }
3369cbb8f0dSMatthew G. Knepley   ierr = PetscFree(indices);CHKERRQ(ierr);
3379cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
3389cbb8f0dSMatthew G. Knepley }
3399cbb8f0dSMatthew G. Knepley 
3409cbb8f0dSMatthew G. Knepley /* Set the constrained indices on each point */
3419cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCIndices(DM dm, PetscSection section)
3429cbb8f0dSMatthew G. Knepley {
3439cbb8f0dSMatthew G. Knepley   PetscInt      *indices;
344baef929fSMatthew G. Knepley   PetscInt       Nf, maxDof, pStart, pEnd, p, f, d;
3459cbb8f0dSMatthew G. Knepley   PetscErrorCode ierr;
3469cbb8f0dSMatthew G. Knepley 
3479cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
348baef929fSMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
3499cbb8f0dSMatthew G. Knepley   ierr = PetscSectionGetMaxDof(section, &maxDof);CHKERRQ(ierr);
3509cbb8f0dSMatthew G. Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
3519cbb8f0dSMatthew G. Knepley   ierr = PetscMalloc1(maxDof, &indices);CHKERRQ(ierr);
3529cbb8f0dSMatthew G. Knepley   for (d = 0; d < maxDof; ++d) indices[d] = -1;
3539cbb8f0dSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
3549cbb8f0dSMatthew G. Knepley     PetscInt cdof, d;
3559cbb8f0dSMatthew G. Knepley 
3569cbb8f0dSMatthew G. Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
3579cbb8f0dSMatthew G. Knepley     if (cdof) {
358baef929fSMatthew G. Knepley       if (Nf) {
3599cbb8f0dSMatthew G. Knepley         PetscInt numConst = 0, foff = 0;
3609cbb8f0dSMatthew G. Knepley 
361baef929fSMatthew G. Knepley         for (f = 0; f < Nf; ++f) {
3629cbb8f0dSMatthew G. Knepley           const PetscInt *find;
3639cbb8f0dSMatthew G. Knepley           PetscInt        fcdof, fdof;
3649cbb8f0dSMatthew G. Knepley 
3659cbb8f0dSMatthew G. Knepley           ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr);
3669cbb8f0dSMatthew G. Knepley           ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr);
3679cbb8f0dSMatthew G. Knepley           /* Change constraint numbering from field component to local dof number */
3689cbb8f0dSMatthew G. Knepley           ierr = PetscSectionGetFieldConstraintIndices(section, p, f, &find);CHKERRQ(ierr);
3699cbb8f0dSMatthew G. Knepley           for (d = 0; d < fcdof; ++d) indices[numConst+d] = find[d] + foff;
3709cbb8f0dSMatthew G. Knepley           numConst += fcdof;
3719cbb8f0dSMatthew G. Knepley           foff     += fdof;
3729cbb8f0dSMatthew G. Knepley         }
3739cbb8f0dSMatthew G. Knepley         if (cdof != numConst) {ierr = PetscSectionSetConstraintDof(section, p, numConst);CHKERRQ(ierr);}
3749cbb8f0dSMatthew G. Knepley       } else {
3759cbb8f0dSMatthew G. Knepley         for (d = 0; d < cdof; ++d) indices[d] = d;
3769cbb8f0dSMatthew G. Knepley       }
3779cbb8f0dSMatthew G. Knepley       ierr = PetscSectionSetConstraintIndices(section, p, indices);CHKERRQ(ierr);
3789cbb8f0dSMatthew G. Knepley     }
3799cbb8f0dSMatthew G. Knepley   }
3809cbb8f0dSMatthew G. Knepley   ierr = PetscFree(indices);CHKERRQ(ierr);
3819cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
3829cbb8f0dSMatthew G. Knepley }
3839cbb8f0dSMatthew G. Knepley 
3849cbb8f0dSMatthew G. Knepley /*@C
3859cbb8f0dSMatthew G. Knepley   DMPlexCreateSection - Create a PetscSection based upon the dof layout specification provided.
3869cbb8f0dSMatthew G. Knepley 
3879cbb8f0dSMatthew G. Knepley   Not Collective
3889cbb8f0dSMatthew G. Knepley 
3899cbb8f0dSMatthew G. Knepley   Input Parameters:
3909cbb8f0dSMatthew G. Knepley + dm        - The DMPlex object
39192cfd99aSMartin Diehl . label     - The label indicating the mesh support of each field, or NULL for the whole mesh
3929cbb8f0dSMatthew G. Knepley . numComp   - An array of size numFields that holds the number of components for each field
3939cbb8f0dSMatthew 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
3949cbb8f0dSMatthew G. Knepley . numBC     - The number of boundary conditions
3959cbb8f0dSMatthew G. Knepley . bcField   - An array of size numBC giving the field number for each boundry condition
3969cbb8f0dSMatthew G. Knepley . bcComps   - [Optional] An array of size numBC giving an IS holding the field components to which each boundary condition applies
3979cbb8f0dSMatthew G. Knepley . bcPoints  - An array of size numBC giving an IS holding the Plex points to which each boundary condition applies
3989cbb8f0dSMatthew G. Knepley - perm      - Optional permutation of the chart, or NULL
3999cbb8f0dSMatthew G. Knepley 
4009cbb8f0dSMatthew G. Knepley   Output Parameter:
4019cbb8f0dSMatthew G. Knepley . section - The PetscSection object
4029cbb8f0dSMatthew G. Knepley 
4039cbb8f0dSMatthew G. Knepley   Notes:
4049cbb8f0dSMatthew 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
4059cbb8f0dSMatthew G. Knepley   number of dof for field 0 on each edge.
4069cbb8f0dSMatthew G. Knepley 
4079cbb8f0dSMatthew G. Knepley   The chart permutation is the same one set using PetscSectionSetPermutation()
4089cbb8f0dSMatthew G. Knepley 
4099cbb8f0dSMatthew G. Knepley   Level: developer
4109cbb8f0dSMatthew G. Knepley 
4111bb6d2a8SBarry Smith   TODO: How is this related to DMCreateLocalSection()
4121bb6d2a8SBarry Smith 
4139cbb8f0dSMatthew G. Knepley .seealso: DMPlexCreate(), PetscSectionCreate(), PetscSectionSetPermutation()
4149cbb8f0dSMatthew G. Knepley @*/
415baef929fSMatthew 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)
4169cbb8f0dSMatthew G. Knepley {
4179cbb8f0dSMatthew G. Knepley   PetscSection   aSec;
4189cbb8f0dSMatthew G. Knepley   PetscErrorCode ierr;
4199cbb8f0dSMatthew G. Knepley 
4209cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
421baef929fSMatthew G. Knepley   ierr = DMPlexCreateSectionFields(dm, numComp, section);CHKERRQ(ierr);
422baef929fSMatthew G. Knepley   ierr = DMPlexCreateSectionDof(dm, label, numDof, *section);CHKERRQ(ierr);
4239cbb8f0dSMatthew G. Knepley   ierr = DMPlexCreateSectionBCDof(dm, numBC, bcField, bcComps, bcPoints, *section);CHKERRQ(ierr);
4249cbb8f0dSMatthew G. Knepley   if (perm) {ierr = PetscSectionSetPermutation(*section, perm);CHKERRQ(ierr);}
4259cbb8f0dSMatthew G. Knepley   ierr = PetscSectionSetFromOptions(*section);CHKERRQ(ierr);
4269cbb8f0dSMatthew G. Knepley   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
4279cbb8f0dSMatthew G. Knepley   ierr = DMPlexGetAnchors(dm,&aSec,NULL);CHKERRQ(ierr);
4289cbb8f0dSMatthew G. Knepley   if (numBC || aSec) {
4299cbb8f0dSMatthew G. Knepley     ierr = DMPlexCreateSectionBCIndicesField(dm, numBC, bcField, bcComps, bcPoints, *section);CHKERRQ(ierr);
4309cbb8f0dSMatthew G. Knepley     ierr = DMPlexCreateSectionBCIndices(dm, *section);CHKERRQ(ierr);
4319cbb8f0dSMatthew G. Knepley   }
4329cbb8f0dSMatthew G. Knepley   ierr = PetscSectionViewFromOptions(*section,NULL,"-section_view");CHKERRQ(ierr);
4339cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
4349cbb8f0dSMatthew G. Knepley }
4359cbb8f0dSMatthew G. Knepley 
4361bb6d2a8SBarry Smith PetscErrorCode DMCreateLocalSection_Plex(DM dm)
4379cbb8f0dSMatthew G. Knepley {
4389cbb8f0dSMatthew G. Knepley   PetscSection   section;
439baef929fSMatthew G. Knepley   DMLabel       *labels;
4409cbb8f0dSMatthew G. Knepley   IS            *bcPoints, *bcComps;
4419cbb8f0dSMatthew G. Knepley   PetscBool     *isFE;
4429cbb8f0dSMatthew G. Knepley   PetscInt      *bcFields, *numComp, *numDof;
4439310035eSMatthew G. Knepley   PetscInt       depth, dim, numBC = 0, Nf, Nds, s, bc = 0, f;
4449cbb8f0dSMatthew G. Knepley   PetscInt       cStart, cEnd, cEndInterior;
4459cbb8f0dSMatthew G. Knepley   PetscErrorCode ierr;
4469cbb8f0dSMatthew G. Knepley 
4479cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
448baef929fSMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
4499310035eSMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
4509310035eSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
4519cbb8f0dSMatthew G. Knepley   /* FE and FV boundary conditions are handled slightly differently */
452baef929fSMatthew G. Knepley   ierr = PetscMalloc1(Nf, &isFE);CHKERRQ(ierr);
453baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
4549cbb8f0dSMatthew G. Knepley     PetscObject  obj;
4559cbb8f0dSMatthew G. Knepley     PetscClassId id;
4569cbb8f0dSMatthew G. Knepley 
45744a7f3ddSMatthew G. Knepley     ierr = DMGetField(dm, f, NULL, &obj);CHKERRQ(ierr);
4589cbb8f0dSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
4599cbb8f0dSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {isFE[f] = PETSC_TRUE;}
4609cbb8f0dSMatthew G. Knepley     else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;}
4619cbb8f0dSMatthew G. Knepley     else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f);
4629cbb8f0dSMatthew G. Knepley   }
4639cbb8f0dSMatthew G. Knepley   /* Allocate boundary point storage for FEM boundaries */
4649310035eSMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
4659310035eSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
4669310035eSMatthew G. Knepley     PetscDS  dsBC;
4679310035eSMatthew G. Knepley     PetscInt numBd, bd;
4689310035eSMatthew G. Knepley 
4699310035eSMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC);CHKERRQ(ierr);
4709310035eSMatthew G. Knepley     ierr = PetscDSGetNumBoundary(dsBC, &numBd);CHKERRQ(ierr);
471baef929fSMatthew G. Knepley     if (!Nf && numBd) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "number of fields is zero and number of boundary conditions is nonzero (this should never happen)");
4729cbb8f0dSMatthew G. Knepley     for (bd = 0; bd < numBd; ++bd) {
4739cbb8f0dSMatthew G. Knepley       PetscInt                field;
4749cbb8f0dSMatthew G. Knepley       DMBoundaryConditionType type;
4759cbb8f0dSMatthew G. Knepley       DMLabel                 label;
4769cbb8f0dSMatthew G. Knepley 
47745480ffeSMatthew G. Knepley       ierr = PetscDSGetBoundary(dsBC, bd, NULL, &type, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
4789cbb8f0dSMatthew G. Knepley       if (label && isFE[field] && (type & DM_BC_ESSENTIAL)) ++numBC;
4799cbb8f0dSMatthew G. Knepley     }
4809310035eSMatthew G. Knepley   }
4819cbb8f0dSMatthew G. Knepley   /* Add ghost cell boundaries for FVM */
4829310035eSMatthew G. Knepley   ierr = DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);CHKERRQ(ierr);
483baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) if (!isFE[f] && cEndInterior >= 0) ++numBC;
4849cbb8f0dSMatthew G. Knepley   ierr = PetscCalloc3(numBC, &bcFields, numBC, &bcPoints, numBC, &bcComps);CHKERRQ(ierr);
4859cbb8f0dSMatthew G. Knepley   /* Constrain ghost cells for FV */
486baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
4879cbb8f0dSMatthew G. Knepley     PetscInt *newidx, c;
4889cbb8f0dSMatthew G. Knepley 
4899cbb8f0dSMatthew G. Knepley     if (isFE[f] || cEndInterior < 0) continue;
4909cbb8f0dSMatthew G. Knepley     ierr = PetscMalloc1(cEnd-cEndInterior,&newidx);CHKERRQ(ierr);
4919cbb8f0dSMatthew G. Knepley     for (c = cEndInterior; c < cEnd; ++c) newidx[c-cEndInterior] = c;
4929cbb8f0dSMatthew G. Knepley     bcFields[bc] = f;
49369293d4bSJed Brown     ierr = ISCreateGeneral(PETSC_COMM_SELF, cEnd-cEndInterior, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]);CHKERRQ(ierr);
4949cbb8f0dSMatthew G. Knepley   }
4959cbb8f0dSMatthew G. Knepley   /* Handle FEM Dirichlet boundaries */
4969310035eSMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
4979310035eSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
4989310035eSMatthew G. Knepley     PetscDS  dsBC;
4999310035eSMatthew G. Knepley     PetscInt numBd, bd;
5009310035eSMatthew G. Knepley 
5019310035eSMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC);CHKERRQ(ierr);
5029310035eSMatthew G. Knepley     ierr = PetscDSGetNumBoundary(dsBC, &numBd);CHKERRQ(ierr);
5039cbb8f0dSMatthew G. Knepley     for (bd = 0; bd < numBd; ++bd) {
5049cbb8f0dSMatthew G. Knepley       DMLabel                 label;
5059cbb8f0dSMatthew G. Knepley       const PetscInt         *comps;
5069cbb8f0dSMatthew G. Knepley       const PetscInt         *values;
5079310035eSMatthew G. Knepley       PetscInt                bd2, field, numComps, numValues;
5089cbb8f0dSMatthew G. Knepley       DMBoundaryConditionType type;
5099310035eSMatthew G. Knepley       PetscBool               duplicate = PETSC_FALSE;
5109cbb8f0dSMatthew G. Knepley 
51145480ffeSMatthew G. Knepley       ierr = PetscDSGetBoundary(dsBC, bd, NULL, &type, NULL, &label, &numValues, &values, &field, &numComps, &comps, NULL, NULL, NULL);CHKERRQ(ierr);
5129cbb8f0dSMatthew G. Knepley       if (!isFE[field] || !label) continue;
5139310035eSMatthew G. Knepley       /* Only want to modify label once */
5149310035eSMatthew G. Knepley       for (bd2 = 0; bd2 < bd; ++bd2) {
51545480ffeSMatthew G. Knepley         DMLabel l;
51645480ffeSMatthew G. Knepley 
51745480ffeSMatthew G. Knepley         ierr = PetscDSGetBoundary(dsBC, bd2, NULL, NULL, NULL, &l, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
51845480ffeSMatthew G. Knepley         duplicate = l == label ? PETSC_TRUE : PETSC_FALSE;
5199310035eSMatthew G. Knepley         if (duplicate) break;
5209310035eSMatthew G. Knepley       }
5219310035eSMatthew G. Knepley       if (!duplicate && (isFE[field])) {
5229310035eSMatthew G. Knepley         /* don't complete cells, which are just present to give orientation to the boundary */
5239310035eSMatthew G. Knepley         ierr = DMPlexLabelComplete(dm, label);CHKERRQ(ierr);
5249310035eSMatthew G. Knepley       }
5259cbb8f0dSMatthew G. Knepley       /* Filter out cells, if you actually want to constrain cells you need to do things by hand right now */
5269cbb8f0dSMatthew G. Knepley       if (type & DM_BC_ESSENTIAL) {
5279cbb8f0dSMatthew G. Knepley         PetscInt       *newidx;
5289cbb8f0dSMatthew G. Knepley         PetscInt        n, newn = 0, p, v;
5299cbb8f0dSMatthew G. Knepley 
5309cbb8f0dSMatthew G. Knepley         bcFields[bc] = field;
5319310035eSMatthew G. Knepley         if (numComps) {ierr = ISCreateGeneral(PetscObjectComm((PetscObject) dm), numComps, comps, PETSC_COPY_VALUES, &bcComps[bc]);CHKERRQ(ierr);}
5329cbb8f0dSMatthew G. Knepley         for (v = 0; v < numValues; ++v) {
5339cbb8f0dSMatthew G. Knepley           IS              tmp;
5349cbb8f0dSMatthew G. Knepley           const PetscInt *idx;
5359cbb8f0dSMatthew G. Knepley 
53645480ffeSMatthew G. Knepley           ierr = DMLabelGetStratumIS(label, values[v], &tmp);CHKERRQ(ierr);
5379cbb8f0dSMatthew G. Knepley           if (!tmp) continue;
5389cbb8f0dSMatthew G. Knepley           ierr = ISGetLocalSize(tmp, &n);CHKERRQ(ierr);
5399cbb8f0dSMatthew G. Knepley           ierr = ISGetIndices(tmp, &idx);CHKERRQ(ierr);
5409cbb8f0dSMatthew G. Knepley           if (isFE[field]) {
5419cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) ++newn;
5429cbb8f0dSMatthew G. Knepley           } else {
5439cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) ++newn;
5449cbb8f0dSMatthew G. Knepley           }
5459cbb8f0dSMatthew G. Knepley           ierr = ISRestoreIndices(tmp, &idx);CHKERRQ(ierr);
5469cbb8f0dSMatthew G. Knepley           ierr = ISDestroy(&tmp);CHKERRQ(ierr);
5479cbb8f0dSMatthew G. Knepley         }
5489cbb8f0dSMatthew G. Knepley         ierr = PetscMalloc1(newn, &newidx);CHKERRQ(ierr);
5499cbb8f0dSMatthew G. Knepley         newn = 0;
5509cbb8f0dSMatthew G. Knepley         for (v = 0; v < numValues; ++v) {
5519cbb8f0dSMatthew G. Knepley           IS              tmp;
5529cbb8f0dSMatthew G. Knepley           const PetscInt *idx;
5539cbb8f0dSMatthew G. Knepley 
55445480ffeSMatthew G. Knepley           ierr = DMLabelGetStratumIS(label, values[v], &tmp);CHKERRQ(ierr);
5559cbb8f0dSMatthew G. Knepley           if (!tmp) continue;
5569cbb8f0dSMatthew G. Knepley           ierr = ISGetLocalSize(tmp, &n);CHKERRQ(ierr);
5579cbb8f0dSMatthew G. Knepley           ierr = ISGetIndices(tmp, &idx);CHKERRQ(ierr);
5589cbb8f0dSMatthew G. Knepley           if (isFE[field]) {
5599cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) newidx[newn++] = idx[p];
5609cbb8f0dSMatthew G. Knepley           } else {
5619cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) newidx[newn++] = idx[p];
5629cbb8f0dSMatthew G. Knepley           }
5639cbb8f0dSMatthew G. Knepley           ierr = ISRestoreIndices(tmp, &idx);CHKERRQ(ierr);
5649cbb8f0dSMatthew G. Knepley           ierr = ISDestroy(&tmp);CHKERRQ(ierr);
5659cbb8f0dSMatthew G. Knepley         }
566665f567fSMatthew G. Knepley         ierr = ISCreateGeneral(PETSC_COMM_SELF, newn, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]);CHKERRQ(ierr);
5679310035eSMatthew G. Knepley       }
5689cbb8f0dSMatthew G. Knepley     }
5699cbb8f0dSMatthew G. Knepley   }
5709cbb8f0dSMatthew G. Knepley   /* Handle discretization */
571baef929fSMatthew G. Knepley   ierr = PetscCalloc3(Nf,&labels,Nf,&numComp,Nf*(dim+1),&numDof);CHKERRQ(ierr);
572baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
57344a7f3ddSMatthew G. Knepley     labels[f] = dm->fields[f].label;
5749cbb8f0dSMatthew G. Knepley     if (isFE[f]) {
57544a7f3ddSMatthew G. Knepley       PetscFE         fe = (PetscFE) dm->fields[f].disc;
5769cbb8f0dSMatthew G. Knepley       const PetscInt *numFieldDof;
577a3254110SMatthew Knepley       PetscInt        fedim, d;
5789cbb8f0dSMatthew G. Knepley 
5799cbb8f0dSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr);
5809cbb8f0dSMatthew G. Knepley       ierr = PetscFEGetNumDof(fe, &numFieldDof);CHKERRQ(ierr);
581a3254110SMatthew Knepley       ierr = PetscFEGetSpatialDimension(fe, &fedim);CHKERRQ(ierr);
582a3254110SMatthew Knepley       for (d = 0; d < PetscMin(dim, fedim)+1; ++d) numDof[f*(dim+1)+d] = numFieldDof[d];
5839cbb8f0dSMatthew G. Knepley     } else {
58444a7f3ddSMatthew G. Knepley       PetscFV fv = (PetscFV) dm->fields[f].disc;
5859cbb8f0dSMatthew G. Knepley 
5869cbb8f0dSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr);
5879cbb8f0dSMatthew G. Knepley       numDof[f*(dim+1)+dim] = numComp[f];
5889cbb8f0dSMatthew G. Knepley     }
5899cbb8f0dSMatthew G. Knepley   }
5909310035eSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
591baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
5929cbb8f0dSMatthew G. Knepley     PetscInt d;
5939cbb8f0dSMatthew G. Knepley     for (d = 1; d < dim; ++d) {
5949cbb8f0dSMatthew G. Knepley       if ((numDof[f*(dim+1)+d] > 0) && (depth < dim)) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated when unknowns are specified on edges or faces.");
5959cbb8f0dSMatthew G. Knepley     }
5969cbb8f0dSMatthew G. Knepley   }
597baef929fSMatthew G. Knepley   ierr = DMPlexCreateSection(dm, labels, numComp, numDof, numBC, bcFields, bcComps, bcPoints, NULL, &section);CHKERRQ(ierr);
598baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
5999cbb8f0dSMatthew G. Knepley     PetscFE     fe;
6009cbb8f0dSMatthew G. Knepley     const char *name;
6019cbb8f0dSMatthew G. Knepley 
60244a7f3ddSMatthew G. Knepley     ierr = DMGetField(dm, f, NULL, (PetscObject *) &fe);CHKERRQ(ierr);
603cdfe53dfSMatthew G. Knepley     if (!((PetscObject) fe)->name) continue;
6049cbb8f0dSMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) fe, &name);CHKERRQ(ierr);
6059cbb8f0dSMatthew G. Knepley     ierr = PetscSectionSetFieldName(section, f, name);CHKERRQ(ierr);
6069cbb8f0dSMatthew G. Knepley   }
60792fd8e1eSJed Brown   ierr = DMSetLocalSection(dm, section);CHKERRQ(ierr);
6089cbb8f0dSMatthew G. Knepley   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
6099cbb8f0dSMatthew G. Knepley   for (bc = 0; bc < numBC; ++bc) {ierr = ISDestroy(&bcPoints[bc]);CHKERRQ(ierr);ierr = ISDestroy(&bcComps[bc]);CHKERRQ(ierr);}
6109cbb8f0dSMatthew G. Knepley   ierr = PetscFree3(bcFields,bcPoints,bcComps);CHKERRQ(ierr);
611baef929fSMatthew G. Knepley   ierr = PetscFree3(labels,numComp,numDof);CHKERRQ(ierr);
6129cbb8f0dSMatthew G. Knepley   ierr = PetscFree(isFE);CHKERRQ(ierr);
6139cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
6149cbb8f0dSMatthew G. Knepley }
615