xref: /petsc/src/dm/impls/plex/plexsection.c (revision e0b68406e04bbc41e6819f8c0ff59eea3b2c814e)
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;
659cbb8f0dSMatthew G. Knepley               ierr = DMPlexGetConeSize(K,kStart,&kConeSize);CHKERRQ(ierr);
669cbb8f0dSMatthew 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);
679cbb8f0dSMatthew G. Knepley             }
689cbb8f0dSMatthew G. Knepley             ierr = PetscSectionSetFieldSym(*section,f,sym);CHKERRQ(ierr);
699cbb8f0dSMatthew G. Knepley             ierr = PetscSectionSymDestroy(&sym);CHKERRQ(ierr);
709cbb8f0dSMatthew G. Knepley           }
719cbb8f0dSMatthew G. Knepley         }
729cbb8f0dSMatthew G. Knepley       }
739cbb8f0dSMatthew G. Knepley     }
749cbb8f0dSMatthew G. Knepley   }
759cbb8f0dSMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
769cbb8f0dSMatthew G. Knepley   ierr = PetscSectionSetChart(*section, pStart, pEnd);CHKERRQ(ierr);
77baef929fSMatthew G. Knepley   ierr = PetscFree(isFE);CHKERRQ(ierr);
78baef929fSMatthew G. Knepley   PetscFunctionReturn(0);
79baef929fSMatthew G. Knepley }
80baef929fSMatthew G. Knepley 
81baef929fSMatthew G. Knepley /* Set the number of dof on each point and separate by fields */
82baef929fSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionDof(DM dm, DMLabel label[],const PetscInt numDof[], PetscSection section)
83baef929fSMatthew G. Knepley {
84baef929fSMatthew G. Knepley   DMLabel        depthLabel;
85412e9a14SMatthew G. Knepley   DMPolytopeType ct;
86baef929fSMatthew G. Knepley   PetscInt       depth, cellHeight, pStart = 0, pEnd = 0;
87a55f9a55SMatthew G. Knepley   PetscInt       Nf, f, Nds, n, dim, d, dep, p;
88a55f9a55SMatthew G. Knepley   PetscBool     *isFE, hasHybrid = PETSC_FALSE;
89baef929fSMatthew G. Knepley   PetscErrorCode ierr;
90baef929fSMatthew G. Knepley 
91baef929fSMatthew G. Knepley   PetscFunctionBegin;
92baef929fSMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
939cbb8f0dSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
94baef929fSMatthew G. Knepley   ierr = DMPlexGetDepthLabel(dm,&depthLabel);CHKERRQ(ierr);
95baef929fSMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
96a55f9a55SMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
97a55f9a55SMatthew G. Knepley   for (n = 0; n < Nds; ++n) {
98a55f9a55SMatthew G. Knepley     PetscDS   ds;
99a55f9a55SMatthew G. Knepley     PetscBool isHybrid;
100a55f9a55SMatthew G. Knepley 
101a55f9a55SMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, n, NULL, NULL, &ds);CHKERRQ(ierr);
102a55f9a55SMatthew G. Knepley     ierr = PetscDSGetHybrid(ds, &isHybrid);CHKERRQ(ierr);
103a55f9a55SMatthew G. Knepley     if (isHybrid) {hasHybrid = PETSC_TRUE; break;}
104a55f9a55SMatthew G. Knepley   }
105baef929fSMatthew G. Knepley   ierr = PetscMalloc1(Nf, &isFE);CHKERRQ(ierr);
106baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
107baef929fSMatthew G. Knepley     PetscObject  obj;
108baef929fSMatthew G. Knepley     PetscClassId id;
109baef929fSMatthew G. Knepley 
11044a7f3ddSMatthew G. Knepley     ierr = DMGetField(dm, f, NULL, &obj);CHKERRQ(ierr);
111baef929fSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1121274e1a1SLawrence Mitchell     /* User is allowed to put a "placeholder" field in (c.f. DMCreateDS) */
1131274e1a1SLawrence Mitchell     isFE[f] = id == PETSCFE_CLASSID ? PETSC_TRUE : PETSC_FALSE;
114baef929fSMatthew G. Knepley   }
115baef929fSMatthew G. Knepley 
1169cbb8f0dSMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
117baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
118*e0b68406SMatthew Knepley     PetscBool avoidTensor;
119*e0b68406SMatthew Knepley 
120*e0b68406SMatthew Knepley     ierr = DMGetFieldAvoidTensor(dm, f, &avoidTensor);CHKERRQ(ierr);
121*e0b68406SMatthew Knepley     avoidTensor = (avoidTensor || hasHybrid) ? PETSC_TRUE : PETSC_FALSE;
122baef929fSMatthew G. Knepley     if (label && label[f]) {
123baef929fSMatthew G. Knepley       IS              pointIS;
124baef929fSMatthew G. Knepley       const PetscInt *points;
125baef929fSMatthew G. Knepley       PetscInt        n;
126baef929fSMatthew G. Knepley 
127baef929fSMatthew G. Knepley       ierr = DMLabelGetStratumIS(label[f], 1, &pointIS);CHKERRQ(ierr);
128baef929fSMatthew G. Knepley       if (!pointIS) continue;
129baef929fSMatthew G. Knepley       ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
130baef929fSMatthew G. Knepley       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
131baef929fSMatthew G. Knepley       for (p = 0; p < n; ++p) {
132baef929fSMatthew G. Knepley         const PetscInt point = points[p];
133baef929fSMatthew G. Knepley         PetscInt       dof, d;
134baef929fSMatthew G. Knepley 
135412e9a14SMatthew G. Knepley         ierr = DMPlexGetCellType(dm, point, &ct);CHKERRQ(ierr);
136baef929fSMatthew G. Knepley         ierr = DMLabelGetValue(depthLabel, point, &d);CHKERRQ(ierr);
137412e9a14SMatthew G. Knepley         /* If this is a tensor prism point, use dof for one dimension lower */
138412e9a14SMatthew G. Knepley         switch (ct) {
139412e9a14SMatthew G. Knepley           case DM_POLYTOPE_POINT_PRISM_TENSOR:
140412e9a14SMatthew G. Knepley           case DM_POLYTOPE_SEG_PRISM_TENSOR:
141412e9a14SMatthew G. Knepley           case DM_POLYTOPE_TRI_PRISM_TENSOR:
142412e9a14SMatthew G. Knepley           case DM_POLYTOPE_QUAD_PRISM_TENSOR:
143a55f9a55SMatthew G. Knepley             if (hasHybrid) {--d;} break;
144412e9a14SMatthew G. Knepley           default: break;
145412e9a14SMatthew G. Knepley         }
146baef929fSMatthew G. Knepley         dof  = d < 0 ? 0 : numDof[f*(dim+1)+d];
147baef929fSMatthew G. Knepley         ierr = PetscSectionSetFieldDof(section, point, f, dof);CHKERRQ(ierr);
148baef929fSMatthew G. Knepley         ierr = PetscSectionAddDof(section, point, dof);CHKERRQ(ierr);
149baef929fSMatthew G. Knepley       }
150baef929fSMatthew G. Knepley       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
151baef929fSMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
152baef929fSMatthew G. Knepley     } else {
1539cbb8f0dSMatthew G. Knepley       for (dep = 0; dep <= depth - cellHeight; ++dep) {
1543fd2e703SMatthew G. Knepley         /* Cases: dim > depth (cell-vertex mesh), dim == depth (fully interpolated), dim < depth (interpolated submesh) */
1553fd2e703SMatthew G. Knepley         d    = dim <= depth ? dep : (!dep ? 0 : dim);
1569cbb8f0dSMatthew G. Knepley         ierr = DMPlexGetDepthStratum(dm, dep, &pStart, &pEnd);CHKERRQ(ierr);
1579cbb8f0dSMatthew G. Knepley         for (p = pStart; p < pEnd; ++p) {
158baef929fSMatthew G. Knepley           const PetscInt dof = numDof[f*(dim+1)+d];
1599cbb8f0dSMatthew G. Knepley 
160412e9a14SMatthew G. Knepley           ierr = DMPlexGetCellType(dm, p, &ct);CHKERRQ(ierr);
161412e9a14SMatthew G. Knepley           switch (ct) {
162412e9a14SMatthew G. Knepley             case DM_POLYTOPE_POINT_PRISM_TENSOR:
163412e9a14SMatthew G. Knepley             case DM_POLYTOPE_SEG_PRISM_TENSOR:
164412e9a14SMatthew G. Knepley             case DM_POLYTOPE_TRI_PRISM_TENSOR:
165412e9a14SMatthew G. Knepley             case DM_POLYTOPE_QUAD_PRISM_TENSOR:
166*e0b68406SMatthew Knepley               if (avoidTensor && isFE[f]) continue;
167412e9a14SMatthew G. Knepley             default: break;
168412e9a14SMatthew G. Knepley           }
169baef929fSMatthew G. Knepley           ierr = PetscSectionSetFieldDof(section, p, f, dof);CHKERRQ(ierr);
170baef929fSMatthew G. Knepley           ierr = PetscSectionAddDof(section, p, dof);CHKERRQ(ierr);
1719cbb8f0dSMatthew G. Knepley         }
172baef929fSMatthew G. Knepley       }
1739cbb8f0dSMatthew G. Knepley     }
1749cbb8f0dSMatthew G. Knepley   }
1759cbb8f0dSMatthew G. Knepley   ierr = PetscFree(isFE);CHKERRQ(ierr);
1769cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
1779cbb8f0dSMatthew G. Knepley }
1789cbb8f0dSMatthew G. Knepley 
1799cbb8f0dSMatthew G. Knepley /* Set the number of dof on each point and separate by fields
1809cbb8f0dSMatthew G. Knepley    If bcComps is NULL or the IS is NULL, constrain every dof on the point
1819cbb8f0dSMatthew G. Knepley */
1829cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCDof(DM dm, PetscInt numBC, const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section)
1839cbb8f0dSMatthew G. Knepley {
184baef929fSMatthew G. Knepley   PetscInt       Nf;
1859cbb8f0dSMatthew G. Knepley   PetscInt       bc;
1869cbb8f0dSMatthew G. Knepley   PetscSection   aSec;
1879cbb8f0dSMatthew G. Knepley   PetscErrorCode ierr;
1889cbb8f0dSMatthew G. Knepley 
1899cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
190baef929fSMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1919cbb8f0dSMatthew G. Knepley   for (bc = 0; bc < numBC; ++bc) {
1929cbb8f0dSMatthew G. Knepley     PetscInt        field = 0;
1939cbb8f0dSMatthew G. Knepley     const PetscInt *comp;
1949cbb8f0dSMatthew G. Knepley     const PetscInt *idx;
195ce093827SMatthew G. Knepley     PetscInt        Nc = 0, cNc = -1, n, i;
1969cbb8f0dSMatthew G. Knepley 
197ce093827SMatthew G. Knepley     if (Nf) {
198ce093827SMatthew G. Knepley       field = bcField[bc];
199ce093827SMatthew G. Knepley       ierr = PetscSectionGetFieldComponents(section, field, &Nc);CHKERRQ(ierr);
200ce093827SMatthew G. Knepley     }
201ce093827SMatthew G. Knepley     if (bcComps && bcComps[bc]) {ierr = ISGetLocalSize(bcComps[bc], &cNc);CHKERRQ(ierr);}
2029cbb8f0dSMatthew G. Knepley     if (bcComps && bcComps[bc]) {ierr = ISGetIndices(bcComps[bc], &comp);CHKERRQ(ierr);}
2039cbb8f0dSMatthew G. Knepley     ierr = ISGetLocalSize(bcPoints[bc], &n);CHKERRQ(ierr);
2049cbb8f0dSMatthew G. Knepley     ierr = ISGetIndices(bcPoints[bc], &idx);CHKERRQ(ierr);
2059cbb8f0dSMatthew G. Knepley     for (i = 0; i < n; ++i) {
2069cbb8f0dSMatthew G. Knepley       const PetscInt p = idx[i];
2079cbb8f0dSMatthew G. Knepley       PetscInt       numConst;
2089cbb8f0dSMatthew G. Knepley 
209baef929fSMatthew G. Knepley       if (Nf) {
2109cbb8f0dSMatthew G. Knepley         ierr = PetscSectionGetFieldDof(section, p, field, &numConst);CHKERRQ(ierr);
2119cbb8f0dSMatthew G. Knepley       } else {
2129cbb8f0dSMatthew G. Knepley         ierr = PetscSectionGetDof(section, p, &numConst);CHKERRQ(ierr);
2139cbb8f0dSMatthew G. Knepley       }
214ce093827SMatthew G. Knepley       /* If Nc <= 0, constrain every dof on the point */
215ce093827SMatthew G. Knepley       if (cNc > 0) {
216ce093827SMatthew G. Knepley         /* We assume that a point may have multiple "nodes", which are collections of Nc dofs,
217ce093827SMatthew G. Knepley            and that those dofs are numbered n*Nc+c */
218ce093827SMatthew G. Knepley         if (Nf) {
219ce093827SMatthew 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);
220ce093827SMatthew G. Knepley           numConst = (numConst/Nc) * cNc;
221ce093827SMatthew G. Knepley         } else {
222ce093827SMatthew G. Knepley           numConst = PetscMin(numConst, cNc);
223ce093827SMatthew G. Knepley         }
224ce093827SMatthew G. Knepley       }
225baef929fSMatthew G. Knepley       if (Nf) {ierr = PetscSectionAddFieldConstraintDof(section, p, field, numConst);CHKERRQ(ierr);}
2269cbb8f0dSMatthew G. Knepley       ierr = PetscSectionAddConstraintDof(section, p, numConst);CHKERRQ(ierr);
2279cbb8f0dSMatthew G. Knepley     }
2289cbb8f0dSMatthew G. Knepley     ierr = ISRestoreIndices(bcPoints[bc], &idx);CHKERRQ(ierr);
2299cbb8f0dSMatthew G. Knepley     if (bcComps && bcComps[bc]) {ierr = ISRestoreIndices(bcComps[bc], &comp);CHKERRQ(ierr);}
2309cbb8f0dSMatthew G. Knepley   }
2319cbb8f0dSMatthew G. Knepley   ierr = DMPlexGetAnchors(dm, &aSec, NULL);CHKERRQ(ierr);
2329cbb8f0dSMatthew G. Knepley   if (aSec) {
2339cbb8f0dSMatthew G. Knepley     PetscInt aStart, aEnd, a;
2349cbb8f0dSMatthew G. Knepley 
2359cbb8f0dSMatthew G. Knepley     ierr = PetscSectionGetChart(aSec, &aStart, &aEnd);CHKERRQ(ierr);
2369cbb8f0dSMatthew G. Knepley     for (a = aStart; a < aEnd; a++) {
2379cbb8f0dSMatthew G. Knepley       PetscInt dof, f;
2389cbb8f0dSMatthew G. Knepley 
2399cbb8f0dSMatthew G. Knepley       ierr = PetscSectionGetDof(aSec, a, &dof);CHKERRQ(ierr);
2409cbb8f0dSMatthew G. Knepley       if (dof) {
2419cbb8f0dSMatthew G. Knepley         /* if there are point-to-point constraints, then all dofs are constrained */
2429cbb8f0dSMatthew G. Knepley         ierr = PetscSectionGetDof(section, a, &dof);CHKERRQ(ierr);
2439cbb8f0dSMatthew G. Knepley         ierr = PetscSectionSetConstraintDof(section, a, dof);CHKERRQ(ierr);
244baef929fSMatthew G. Knepley         for (f = 0; f < Nf; f++) {
2459cbb8f0dSMatthew G. Knepley           ierr = PetscSectionGetFieldDof(section, a, f, &dof);CHKERRQ(ierr);
2469cbb8f0dSMatthew G. Knepley           ierr = PetscSectionSetFieldConstraintDof(section, a, f, dof);CHKERRQ(ierr);
2479cbb8f0dSMatthew G. Knepley         }
2489cbb8f0dSMatthew G. Knepley       }
2499cbb8f0dSMatthew G. Knepley     }
2509cbb8f0dSMatthew G. Knepley   }
2519cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
2529cbb8f0dSMatthew G. Knepley }
2539cbb8f0dSMatthew G. Knepley 
2549cbb8f0dSMatthew G. Knepley /* Set the constrained field indices on each point
2559cbb8f0dSMatthew G. Knepley    If bcComps is NULL or the IS is NULL, constrain every dof on the point
2569cbb8f0dSMatthew G. Knepley */
2579cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCIndicesField(DM dm, PetscInt numBC,const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section)
2589cbb8f0dSMatthew G. Knepley {
2599cbb8f0dSMatthew G. Knepley   PetscSection   aSec;
2609cbb8f0dSMatthew G. Knepley   PetscInt      *indices;
261baef929fSMatthew G. Knepley   PetscInt       Nf, cdof, maxDof = 0, pStart, pEnd, p, bc, f, d;
2629cbb8f0dSMatthew G. Knepley   PetscErrorCode ierr;
2639cbb8f0dSMatthew G. Knepley 
2649cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
265baef929fSMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
266baef929fSMatthew G. Knepley   if (!Nf) PetscFunctionReturn(0);
2679cbb8f0dSMatthew G. Knepley   /* Initialize all field indices to -1 */
2689cbb8f0dSMatthew G. Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
2699cbb8f0dSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); maxDof = PetscMax(maxDof, cdof);}
2709cbb8f0dSMatthew G. Knepley   ierr = PetscMalloc1(maxDof, &indices);CHKERRQ(ierr);
2719cbb8f0dSMatthew G. Knepley   for (d = 0; d < maxDof; ++d) indices[d] = -1;
272baef929fSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) for (f = 0; f < Nf; ++f) {ierr = PetscSectionSetFieldConstraintIndices(section, p, f, indices);CHKERRQ(ierr);}
2739cbb8f0dSMatthew G. Knepley   /* Handle BC constraints */
2749cbb8f0dSMatthew G. Knepley   for (bc = 0; bc < numBC; ++bc) {
2759cbb8f0dSMatthew G. Knepley     const PetscInt  field = bcField[bc];
2769cbb8f0dSMatthew G. Knepley     const PetscInt *comp, *idx;
277ce093827SMatthew G. Knepley     PetscInt        Nc, cNc = -1, n, i;
2789cbb8f0dSMatthew G. Knepley 
279ce093827SMatthew G. Knepley     ierr = PetscSectionGetFieldComponents(section, field, &Nc);CHKERRQ(ierr);
280ce093827SMatthew G. Knepley     if (bcComps && bcComps[bc]) {ierr = ISGetLocalSize(bcComps[bc], &cNc);CHKERRQ(ierr);}
2819cbb8f0dSMatthew G. Knepley     if (bcComps && bcComps[bc]) {ierr = ISGetIndices(bcComps[bc], &comp);CHKERRQ(ierr);}
2829cbb8f0dSMatthew G. Knepley     ierr = ISGetLocalSize(bcPoints[bc], &n);CHKERRQ(ierr);
2839cbb8f0dSMatthew G. Knepley     ierr = ISGetIndices(bcPoints[bc], &idx);CHKERRQ(ierr);
2849cbb8f0dSMatthew G. Knepley     for (i = 0; i < n; ++i) {
2859cbb8f0dSMatthew G. Knepley       const PetscInt  p = idx[i];
2869cbb8f0dSMatthew G. Knepley       const PetscInt *find;
287ce093827SMatthew G. Knepley       PetscInt        fdof, fcdof, c, j;
2889cbb8f0dSMatthew G. Knepley 
2899cbb8f0dSMatthew G. Knepley       ierr = PetscSectionGetFieldDof(section, p, field, &fdof);CHKERRQ(ierr);
2909cbb8f0dSMatthew G. Knepley       if (!fdof) continue;
291ce093827SMatthew G. Knepley       if (cNc < 0) {
2929cbb8f0dSMatthew G. Knepley         for (d = 0; d < fdof; ++d) indices[d] = d;
2939cbb8f0dSMatthew G. Knepley         fcdof = fdof;
2949cbb8f0dSMatthew G. Knepley       } else {
295ce093827SMatthew G. Knepley         /* We assume that a point may have multiple "nodes", which are collections of Nc dofs,
296ce093827SMatthew G. Knepley            and that those dofs are numbered n*Nc+c */
2979cbb8f0dSMatthew G. Knepley         ierr = PetscSectionGetFieldConstraintDof(section, p, field, &fcdof);CHKERRQ(ierr);
2989cbb8f0dSMatthew G. Knepley         ierr = PetscSectionGetFieldConstraintIndices(section, p, field, &find);CHKERRQ(ierr);
299ce093827SMatthew G. Knepley         /* Get indices constrained by previous bcs */
3009cbb8f0dSMatthew G. Knepley         for (d = 0; d < fcdof; ++d) {if (find[d] < 0) break; indices[d] = find[d];}
301ce093827SMatthew G. Knepley         for (j = 0; j < fdof/Nc; ++j) for (c = 0; c < cNc; ++c) indices[d++] = j*Nc + comp[c];
3029cbb8f0dSMatthew G. Knepley         ierr = PetscSortRemoveDupsInt(&d, indices);CHKERRQ(ierr);
3039cbb8f0dSMatthew G. Knepley         for (c = d; c < fcdof; ++c) indices[c] = -1;
3049cbb8f0dSMatthew G. Knepley         fcdof = d;
3059cbb8f0dSMatthew G. Knepley       }
3069cbb8f0dSMatthew G. Knepley       ierr = PetscSectionSetFieldConstraintDof(section, p, field, fcdof);CHKERRQ(ierr);
3079cbb8f0dSMatthew G. Knepley       ierr = PetscSectionSetFieldConstraintIndices(section, p, field, indices);CHKERRQ(ierr);
3089cbb8f0dSMatthew G. Knepley     }
3099cbb8f0dSMatthew G. Knepley     if (bcComps && bcComps[bc]) {ierr = ISRestoreIndices(bcComps[bc], &comp);CHKERRQ(ierr);}
3109cbb8f0dSMatthew G. Knepley     ierr = ISRestoreIndices(bcPoints[bc], &idx);CHKERRQ(ierr);
3119cbb8f0dSMatthew G. Knepley   }
3129cbb8f0dSMatthew G. Knepley   /* Handle anchors */
3139cbb8f0dSMatthew G. Knepley   ierr = DMPlexGetAnchors(dm, &aSec, NULL);CHKERRQ(ierr);
3149cbb8f0dSMatthew G. Knepley   if (aSec) {
3159cbb8f0dSMatthew G. Knepley     PetscInt aStart, aEnd, a;
3169cbb8f0dSMatthew G. Knepley 
3179cbb8f0dSMatthew G. Knepley     for (d = 0; d < maxDof; ++d) indices[d] = d;
3189cbb8f0dSMatthew G. Knepley     ierr = PetscSectionGetChart(aSec, &aStart, &aEnd);CHKERRQ(ierr);
3199cbb8f0dSMatthew G. Knepley     for (a = aStart; a < aEnd; a++) {
3209cbb8f0dSMatthew G. Knepley       PetscInt dof, f;
3219cbb8f0dSMatthew G. Knepley 
3229cbb8f0dSMatthew G. Knepley       ierr = PetscSectionGetDof(aSec, a, &dof);CHKERRQ(ierr);
3239cbb8f0dSMatthew G. Knepley       if (dof) {
3249cbb8f0dSMatthew G. Knepley         /* if there are point-to-point constraints, then all dofs are constrained */
325baef929fSMatthew G. Knepley         for (f = 0; f < Nf; f++) {
3269cbb8f0dSMatthew G. Knepley           ierr = PetscSectionSetFieldConstraintIndices(section, a, f, indices);CHKERRQ(ierr);
3279cbb8f0dSMatthew G. Knepley         }
3289cbb8f0dSMatthew G. Knepley       }
3299cbb8f0dSMatthew G. Knepley     }
3309cbb8f0dSMatthew G. Knepley   }
3319cbb8f0dSMatthew G. Knepley   ierr = PetscFree(indices);CHKERRQ(ierr);
3329cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
3339cbb8f0dSMatthew G. Knepley }
3349cbb8f0dSMatthew G. Knepley 
3359cbb8f0dSMatthew G. Knepley /* Set the constrained indices on each point */
3369cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCIndices(DM dm, PetscSection section)
3379cbb8f0dSMatthew G. Knepley {
3389cbb8f0dSMatthew G. Knepley   PetscInt      *indices;
339baef929fSMatthew G. Knepley   PetscInt       Nf, maxDof, pStart, pEnd, p, f, d;
3409cbb8f0dSMatthew G. Knepley   PetscErrorCode ierr;
3419cbb8f0dSMatthew G. Knepley 
3429cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
343baef929fSMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
3449cbb8f0dSMatthew G. Knepley   ierr = PetscSectionGetMaxDof(section, &maxDof);CHKERRQ(ierr);
3459cbb8f0dSMatthew G. Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
3469cbb8f0dSMatthew G. Knepley   ierr = PetscMalloc1(maxDof, &indices);CHKERRQ(ierr);
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 
3519cbb8f0dSMatthew G. Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
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 
3609cbb8f0dSMatthew G. Knepley           ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr);
3619cbb8f0dSMatthew G. Knepley           ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr);
3629cbb8f0dSMatthew G. Knepley           /* Change constraint numbering from field component to local dof number */
3639cbb8f0dSMatthew G. Knepley           ierr = PetscSectionGetFieldConstraintIndices(section, p, f, &find);CHKERRQ(ierr);
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         }
3689cbb8f0dSMatthew G. Knepley         if (cdof != numConst) {ierr = PetscSectionSetConstraintDof(section, p, numConst);CHKERRQ(ierr);}
3699cbb8f0dSMatthew G. Knepley       } else {
3709cbb8f0dSMatthew G. Knepley         for (d = 0; d < cdof; ++d) indices[d] = d;
3719cbb8f0dSMatthew G. Knepley       }
3729cbb8f0dSMatthew G. Knepley       ierr = PetscSectionSetConstraintIndices(section, p, indices);CHKERRQ(ierr);
3739cbb8f0dSMatthew G. Knepley     }
3749cbb8f0dSMatthew G. Knepley   }
3759cbb8f0dSMatthew G. Knepley   ierr = PetscFree(indices);CHKERRQ(ierr);
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   PetscErrorCode ierr;
4149cbb8f0dSMatthew G. Knepley 
4159cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
416baef929fSMatthew G. Knepley   ierr = DMPlexCreateSectionFields(dm, numComp, section);CHKERRQ(ierr);
417baef929fSMatthew G. Knepley   ierr = DMPlexCreateSectionDof(dm, label, numDof, *section);CHKERRQ(ierr);
4189cbb8f0dSMatthew G. Knepley   ierr = DMPlexCreateSectionBCDof(dm, numBC, bcField, bcComps, bcPoints, *section);CHKERRQ(ierr);
4199cbb8f0dSMatthew G. Knepley   if (perm) {ierr = PetscSectionSetPermutation(*section, perm);CHKERRQ(ierr);}
4209cbb8f0dSMatthew G. Knepley   ierr = PetscSectionSetFromOptions(*section);CHKERRQ(ierr);
4219cbb8f0dSMatthew G. Knepley   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
4229cbb8f0dSMatthew G. Knepley   ierr = DMPlexGetAnchors(dm,&aSec,NULL);CHKERRQ(ierr);
4239cbb8f0dSMatthew G. Knepley   if (numBC || aSec) {
4249cbb8f0dSMatthew G. Knepley     ierr = DMPlexCreateSectionBCIndicesField(dm, numBC, bcField, bcComps, bcPoints, *section);CHKERRQ(ierr);
4259cbb8f0dSMatthew G. Knepley     ierr = DMPlexCreateSectionBCIndices(dm, *section);CHKERRQ(ierr);
4269cbb8f0dSMatthew G. Knepley   }
4279cbb8f0dSMatthew G. Knepley   ierr = PetscSectionViewFromOptions(*section,NULL,"-section_view");CHKERRQ(ierr);
4289cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
4299cbb8f0dSMatthew G. Knepley }
4309cbb8f0dSMatthew G. Knepley 
4311bb6d2a8SBarry Smith PetscErrorCode DMCreateLocalSection_Plex(DM dm)
4329cbb8f0dSMatthew G. Knepley {
4339cbb8f0dSMatthew G. Knepley   PetscSection   section;
434baef929fSMatthew G. Knepley   DMLabel       *labels;
4359cbb8f0dSMatthew G. Knepley   IS            *bcPoints, *bcComps;
4369cbb8f0dSMatthew G. Knepley   PetscBool     *isFE;
4379cbb8f0dSMatthew G. Knepley   PetscInt      *bcFields, *numComp, *numDof;
4389310035eSMatthew G. Knepley   PetscInt       depth, dim, numBC = 0, Nf, Nds, s, bc = 0, f;
4399cbb8f0dSMatthew G. Knepley   PetscInt       cStart, cEnd, cEndInterior;
4409cbb8f0dSMatthew G. Knepley   PetscErrorCode ierr;
4419cbb8f0dSMatthew G. Knepley 
4429cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
443baef929fSMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
4449310035eSMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
4459310035eSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
4469cbb8f0dSMatthew G. Knepley   /* FE and FV boundary conditions are handled slightly differently */
447baef929fSMatthew G. Knepley   ierr = PetscMalloc1(Nf, &isFE);CHKERRQ(ierr);
448baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
4499cbb8f0dSMatthew G. Knepley     PetscObject  obj;
4509cbb8f0dSMatthew G. Knepley     PetscClassId id;
4519cbb8f0dSMatthew G. Knepley 
45244a7f3ddSMatthew G. Knepley     ierr = DMGetField(dm, f, NULL, &obj);CHKERRQ(ierr);
4539cbb8f0dSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
4549cbb8f0dSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {isFE[f] = PETSC_TRUE;}
4559cbb8f0dSMatthew G. Knepley     else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;}
4569cbb8f0dSMatthew G. Knepley     else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f);
4579cbb8f0dSMatthew G. Knepley   }
4589cbb8f0dSMatthew G. Knepley   /* Allocate boundary point storage for FEM boundaries */
4599310035eSMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
4609310035eSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
4619310035eSMatthew G. Knepley     PetscDS  dsBC;
4629310035eSMatthew G. Knepley     PetscInt numBd, bd;
4639310035eSMatthew G. Knepley 
4649310035eSMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC);CHKERRQ(ierr);
4659310035eSMatthew G. Knepley     ierr = PetscDSGetNumBoundary(dsBC, &numBd);CHKERRQ(ierr);
466baef929fSMatthew 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)");
4679cbb8f0dSMatthew G. Knepley     for (bd = 0; bd < numBd; ++bd) {
4689cbb8f0dSMatthew G. Knepley       PetscInt                field;
4699cbb8f0dSMatthew G. Knepley       DMBoundaryConditionType type;
4709cbb8f0dSMatthew G. Knepley       const char             *labelName;
4719cbb8f0dSMatthew G. Knepley       DMLabel                 label;
4729cbb8f0dSMatthew G. Knepley 
47356cf3b9cSMatthew G. Knepley       ierr = PetscDSGetBoundary(dsBC, bd, &type, NULL, &labelName, &field, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
4749cbb8f0dSMatthew G. Knepley       ierr = DMGetLabel(dm, labelName, &label);CHKERRQ(ierr);
4759cbb8f0dSMatthew G. Knepley       if (label && isFE[field] && (type & DM_BC_ESSENTIAL)) ++numBC;
4769cbb8f0dSMatthew G. Knepley     }
4779310035eSMatthew G. Knepley   }
4789cbb8f0dSMatthew G. Knepley   /* Add ghost cell boundaries for FVM */
4799310035eSMatthew G. Knepley   ierr = DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);CHKERRQ(ierr);
480baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) if (!isFE[f] && cEndInterior >= 0) ++numBC;
4819cbb8f0dSMatthew G. Knepley   ierr = PetscCalloc3(numBC, &bcFields, numBC, &bcPoints, numBC, &bcComps);CHKERRQ(ierr);
4829cbb8f0dSMatthew G. Knepley   /* Constrain ghost cells for FV */
483baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
4849cbb8f0dSMatthew G. Knepley     PetscInt *newidx, c;
4859cbb8f0dSMatthew G. Knepley 
4869cbb8f0dSMatthew G. Knepley     if (isFE[f] || cEndInterior < 0) continue;
4879cbb8f0dSMatthew G. Knepley     ierr = PetscMalloc1(cEnd-cEndInterior,&newidx);CHKERRQ(ierr);
4889cbb8f0dSMatthew G. Knepley     for (c = cEndInterior; c < cEnd; ++c) newidx[c-cEndInterior] = c;
4899cbb8f0dSMatthew G. Knepley     bcFields[bc] = f;
49069293d4bSJed Brown     ierr = ISCreateGeneral(PETSC_COMM_SELF, cEnd-cEndInterior, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]);CHKERRQ(ierr);
4919cbb8f0dSMatthew G. Knepley   }
4929cbb8f0dSMatthew G. Knepley   /* Handle FEM Dirichlet boundaries */
4939310035eSMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
4949310035eSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
4959310035eSMatthew G. Knepley     PetscDS  dsBC;
4969310035eSMatthew G. Knepley     PetscInt numBd, bd;
4979310035eSMatthew G. Knepley 
4989310035eSMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC);CHKERRQ(ierr);
4999310035eSMatthew G. Knepley     ierr = PetscDSGetNumBoundary(dsBC, &numBd);CHKERRQ(ierr);
5009cbb8f0dSMatthew G. Knepley     for (bd = 0; bd < numBd; ++bd) {
5019cbb8f0dSMatthew G. Knepley       const char             *bdLabel;
5029cbb8f0dSMatthew G. Knepley       DMLabel                 label;
5039cbb8f0dSMatthew G. Knepley       const PetscInt         *comps;
5049cbb8f0dSMatthew G. Knepley       const PetscInt         *values;
5059310035eSMatthew G. Knepley       PetscInt                bd2, field, numComps, numValues;
5069cbb8f0dSMatthew G. Knepley       DMBoundaryConditionType type;
5079310035eSMatthew G. Knepley       PetscBool               duplicate = PETSC_FALSE;
5089cbb8f0dSMatthew G. Knepley 
50956cf3b9cSMatthew G. Knepley       ierr = PetscDSGetBoundary(dsBC, bd, &type, NULL, &bdLabel, &field, &numComps, &comps, NULL, NULL, &numValues, &values, NULL);CHKERRQ(ierr);
5109cbb8f0dSMatthew G. Knepley       ierr = DMGetLabel(dm, bdLabel, &label);CHKERRQ(ierr);
5119cbb8f0dSMatthew G. Knepley       if (!isFE[field] || !label) continue;
5129310035eSMatthew G. Knepley       /* Only want to modify label once */
5139310035eSMatthew G. Knepley       for (bd2 = 0; bd2 < bd; ++bd2) {
5149310035eSMatthew G. Knepley         const char *bdname;
51556cf3b9cSMatthew G. Knepley         ierr = PetscDSGetBoundary(dsBC, bd2, NULL, NULL, &bdname, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
5169310035eSMatthew G. Knepley         ierr = PetscStrcmp(bdname, bdLabel, &duplicate);CHKERRQ(ierr);
5179310035eSMatthew G. Knepley         if (duplicate) break;
5189310035eSMatthew G. Knepley       }
5199310035eSMatthew G. Knepley       if (!duplicate && (isFE[field])) {
5209310035eSMatthew G. Knepley         /* don't complete cells, which are just present to give orientation to the boundary */
5219310035eSMatthew G. Knepley         ierr = DMPlexLabelComplete(dm, label);CHKERRQ(ierr);
5229310035eSMatthew G. Knepley       }
5239cbb8f0dSMatthew G. Knepley       /* Filter out cells, if you actually want to constrain cells you need to do things by hand right now */
5249cbb8f0dSMatthew G. Knepley       if (type & DM_BC_ESSENTIAL) {
5259cbb8f0dSMatthew G. Knepley         PetscInt       *newidx;
5269cbb8f0dSMatthew G. Knepley         PetscInt        n, newn = 0, p, v;
5279cbb8f0dSMatthew G. Knepley 
5289cbb8f0dSMatthew G. Knepley         bcFields[bc] = field;
5299310035eSMatthew G. Knepley         if (numComps) {ierr = ISCreateGeneral(PetscObjectComm((PetscObject) dm), numComps, comps, PETSC_COPY_VALUES, &bcComps[bc]);CHKERRQ(ierr);}
5309cbb8f0dSMatthew G. Knepley         for (v = 0; v < numValues; ++v) {
5319cbb8f0dSMatthew G. Knepley           IS              tmp;
5329cbb8f0dSMatthew G. Knepley           const PetscInt *idx;
5339cbb8f0dSMatthew G. Knepley 
5349cbb8f0dSMatthew G. Knepley           ierr = DMGetStratumIS(dm, bdLabel, values[v], &tmp);CHKERRQ(ierr);
5359cbb8f0dSMatthew G. Knepley           if (!tmp) continue;
5369cbb8f0dSMatthew G. Knepley           ierr = ISGetLocalSize(tmp, &n);CHKERRQ(ierr);
5379cbb8f0dSMatthew G. Knepley           ierr = ISGetIndices(tmp, &idx);CHKERRQ(ierr);
5389cbb8f0dSMatthew G. Knepley           if (isFE[field]) {
5399cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) ++newn;
5409cbb8f0dSMatthew G. Knepley           } else {
5419cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) ++newn;
5429cbb8f0dSMatthew G. Knepley           }
5439cbb8f0dSMatthew G. Knepley           ierr = ISRestoreIndices(tmp, &idx);CHKERRQ(ierr);
5449cbb8f0dSMatthew G. Knepley           ierr = ISDestroy(&tmp);CHKERRQ(ierr);
5459cbb8f0dSMatthew G. Knepley         }
5469cbb8f0dSMatthew G. Knepley         ierr = PetscMalloc1(newn, &newidx);CHKERRQ(ierr);
5479cbb8f0dSMatthew G. Knepley         newn = 0;
5489cbb8f0dSMatthew G. Knepley         for (v = 0; v < numValues; ++v) {
5499cbb8f0dSMatthew G. Knepley           IS              tmp;
5509cbb8f0dSMatthew G. Knepley           const PetscInt *idx;
5519cbb8f0dSMatthew G. Knepley 
5529cbb8f0dSMatthew G. Knepley           ierr = DMGetStratumIS(dm, bdLabel, values[v], &tmp);CHKERRQ(ierr);
5539cbb8f0dSMatthew G. Knepley           if (!tmp) continue;
5549cbb8f0dSMatthew G. Knepley           ierr = ISGetLocalSize(tmp, &n);CHKERRQ(ierr);
5559cbb8f0dSMatthew G. Knepley           ierr = ISGetIndices(tmp, &idx);CHKERRQ(ierr);
5569cbb8f0dSMatthew G. Knepley           if (isFE[field]) {
5579cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) newidx[newn++] = idx[p];
5589cbb8f0dSMatthew G. Knepley           } else {
5599cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) newidx[newn++] = idx[p];
5609cbb8f0dSMatthew G. Knepley           }
5619cbb8f0dSMatthew G. Knepley           ierr = ISRestoreIndices(tmp, &idx);CHKERRQ(ierr);
5629cbb8f0dSMatthew G. Knepley           ierr = ISDestroy(&tmp);CHKERRQ(ierr);
5639cbb8f0dSMatthew G. Knepley         }
564665f567fSMatthew G. Knepley         ierr = ISCreateGeneral(PETSC_COMM_SELF, newn, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]);CHKERRQ(ierr);
5659310035eSMatthew G. Knepley       }
5669cbb8f0dSMatthew G. Knepley     }
5679cbb8f0dSMatthew G. Knepley   }
5689cbb8f0dSMatthew G. Knepley   /* Handle discretization */
569baef929fSMatthew G. Knepley   ierr = PetscCalloc3(Nf,&labels,Nf,&numComp,Nf*(dim+1),&numDof);CHKERRQ(ierr);
570baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
57144a7f3ddSMatthew G. Knepley     labels[f] = dm->fields[f].label;
5729cbb8f0dSMatthew G. Knepley     if (isFE[f]) {
57344a7f3ddSMatthew G. Knepley       PetscFE         fe = (PetscFE) dm->fields[f].disc;
5749cbb8f0dSMatthew G. Knepley       const PetscInt *numFieldDof;
575a3254110SMatthew Knepley       PetscInt        fedim, d;
5769cbb8f0dSMatthew G. Knepley 
5779cbb8f0dSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr);
5789cbb8f0dSMatthew G. Knepley       ierr = PetscFEGetNumDof(fe, &numFieldDof);CHKERRQ(ierr);
579a3254110SMatthew Knepley       ierr = PetscFEGetSpatialDimension(fe, &fedim);CHKERRQ(ierr);
580a3254110SMatthew Knepley       for (d = 0; d < PetscMin(dim, fedim)+1; ++d) numDof[f*(dim+1)+d] = numFieldDof[d];
5819cbb8f0dSMatthew G. Knepley     } else {
58244a7f3ddSMatthew G. Knepley       PetscFV fv = (PetscFV) dm->fields[f].disc;
5839cbb8f0dSMatthew G. Knepley 
5849cbb8f0dSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr);
5859cbb8f0dSMatthew G. Knepley       numDof[f*(dim+1)+dim] = numComp[f];
5869cbb8f0dSMatthew G. Knepley     }
5879cbb8f0dSMatthew G. Knepley   }
5889310035eSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
589baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
5909cbb8f0dSMatthew G. Knepley     PetscInt d;
5919cbb8f0dSMatthew G. Knepley     for (d = 1; d < dim; ++d) {
5929cbb8f0dSMatthew 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.");
5939cbb8f0dSMatthew G. Knepley     }
5949cbb8f0dSMatthew G. Knepley   }
595baef929fSMatthew G. Knepley   ierr = DMPlexCreateSection(dm, labels, numComp, numDof, numBC, bcFields, bcComps, bcPoints, NULL, &section);CHKERRQ(ierr);
596baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
5979cbb8f0dSMatthew G. Knepley     PetscFE     fe;
5989cbb8f0dSMatthew G. Knepley     const char *name;
5999cbb8f0dSMatthew G. Knepley 
60044a7f3ddSMatthew G. Knepley     ierr = DMGetField(dm, f, NULL, (PetscObject *) &fe);CHKERRQ(ierr);
601cdfe53dfSMatthew G. Knepley     if (!((PetscObject) fe)->name) continue;
6029cbb8f0dSMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) fe, &name);CHKERRQ(ierr);
6039cbb8f0dSMatthew G. Knepley     ierr = PetscSectionSetFieldName(section, f, name);CHKERRQ(ierr);
6049cbb8f0dSMatthew G. Knepley   }
60592fd8e1eSJed Brown   ierr = DMSetLocalSection(dm, section);CHKERRQ(ierr);
6069cbb8f0dSMatthew G. Knepley   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
6079cbb8f0dSMatthew G. Knepley   for (bc = 0; bc < numBC; ++bc) {ierr = ISDestroy(&bcPoints[bc]);CHKERRQ(ierr);ierr = ISDestroy(&bcComps[bc]);CHKERRQ(ierr);}
6089cbb8f0dSMatthew G. Knepley   ierr = PetscFree3(bcFields,bcPoints,bcComps);CHKERRQ(ierr);
609baef929fSMatthew G. Knepley   ierr = PetscFree3(labels,numComp,numDof);CHKERRQ(ierr);
6109cbb8f0dSMatthew G. Knepley   ierr = PetscFree(isFE);CHKERRQ(ierr);
6119cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
6129cbb8f0dSMatthew G. Knepley }
613