xref: /petsc/src/dm/impls/plex/plexsection.c (revision f985c9a8db4e1c917cfec8e832913fbbb999d435)
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));
20524cfc5f6SToby Isaac     if (bcPoints[bc]) {
2069566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(bcPoints[bc], &n));
2079566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(bcPoints[bc], &idx));
2089cbb8f0dSMatthew G. Knepley       for (i = 0; i < n; ++i) {
2099cbb8f0dSMatthew G. Knepley         const PetscInt p = idx[i];
2109cbb8f0dSMatthew G. Knepley         PetscInt       numConst;
2119cbb8f0dSMatthew G. Knepley 
212baef929fSMatthew G. Knepley         if (Nf) {
2139566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldDof(section, p, field, &numConst));
2149cbb8f0dSMatthew G. Knepley         } else {
2159566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetDof(section, p, &numConst));
2169cbb8f0dSMatthew G. Knepley         }
217ce093827SMatthew G. Knepley         /* If Nc <= 0, constrain every dof on the point */
218ce093827SMatthew G. Knepley         if (cNc > 0) {
219ce093827SMatthew G. Knepley           /* We assume that a point may have multiple "nodes", which are collections of Nc dofs,
220ce093827SMatthew G. Knepley              and that those dofs are numbered n*Nc+c */
221ce093827SMatthew G. Knepley           if (Nf) {
22263a3b9bcSJacob Faibussowitsch             PetscCheck(numConst % Nc == 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT " has %" PetscInt_FMT " dof which is not divisible by %" PetscInt_FMT " field components", p, numConst, Nc);
223ce093827SMatthew G. Knepley             numConst = (numConst/Nc) * cNc;
224ce093827SMatthew G. Knepley           } else {
225ce093827SMatthew G. Knepley             numConst = PetscMin(numConst, cNc);
226ce093827SMatthew G. Knepley           }
227ce093827SMatthew G. Knepley         }
2289566063dSJacob Faibussowitsch         if (Nf) PetscCall(PetscSectionAddFieldConstraintDof(section, p, field, numConst));
2299566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddConstraintDof(section, p, numConst));
2309cbb8f0dSMatthew G. Knepley       }
2319566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(bcPoints[bc], &idx));
23224cfc5f6SToby Isaac     }
2339566063dSJacob Faibussowitsch     if (bcComps && bcComps[bc]) PetscCall(ISRestoreIndices(bcComps[bc], &comp));
2349cbb8f0dSMatthew G. Knepley   }
2359566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAnchors(dm, &aSec, NULL));
2369cbb8f0dSMatthew G. Knepley   if (aSec) {
2379cbb8f0dSMatthew G. Knepley     PetscInt aStart, aEnd, a;
2389cbb8f0dSMatthew G. Knepley 
2399566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd));
2409cbb8f0dSMatthew G. Knepley     for (a = aStart; a < aEnd; a++) {
2419cbb8f0dSMatthew G. Knepley       PetscInt dof, f;
2429cbb8f0dSMatthew G. Knepley 
2439566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(aSec, a, &dof));
2449cbb8f0dSMatthew G. Knepley       if (dof) {
2459cbb8f0dSMatthew G. Knepley         /* if there are point-to-point constraints, then all dofs are constrained */
2469566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, a, &dof));
2479566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetConstraintDof(section, a, dof));
248baef929fSMatthew G. Knepley         for (f = 0; f < Nf; f++) {
2499566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldDof(section, a, f, &dof));
2509566063dSJacob Faibussowitsch           PetscCall(PetscSectionSetFieldConstraintDof(section, a, f, dof));
2519cbb8f0dSMatthew G. Knepley         }
2529cbb8f0dSMatthew G. Knepley       }
2539cbb8f0dSMatthew G. Knepley     }
2549cbb8f0dSMatthew G. Knepley   }
2559cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
2569cbb8f0dSMatthew G. Knepley }
2579cbb8f0dSMatthew G. Knepley 
2589cbb8f0dSMatthew G. Knepley /* Set the constrained field indices on each point
2599cbb8f0dSMatthew G. Knepley    If bcComps is NULL or the IS is NULL, constrain every dof on the point
2609cbb8f0dSMatthew G. Knepley */
2619cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCIndicesField(DM dm, PetscInt numBC,const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section)
2629cbb8f0dSMatthew G. Knepley {
2639cbb8f0dSMatthew G. Knepley   PetscSection   aSec;
2649cbb8f0dSMatthew G. Knepley   PetscInt      *indices;
265baef929fSMatthew G. Knepley   PetscInt       Nf, cdof, maxDof = 0, pStart, pEnd, p, bc, f, d;
2669cbb8f0dSMatthew G. Knepley 
2679cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
2689566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &Nf));
269baef929fSMatthew G. Knepley   if (!Nf) PetscFunctionReturn(0);
2709cbb8f0dSMatthew G. Knepley   /* Initialize all field indices to -1 */
2719566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
2729566063dSJacob Faibussowitsch   for (p = pStart; p < pEnd; ++p) {PetscCall(PetscSectionGetConstraintDof(section, p, &cdof)); maxDof = PetscMax(maxDof, cdof);}
2739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxDof, &indices));
2749cbb8f0dSMatthew G. Knepley   for (d = 0; d < maxDof; ++d) indices[d] = -1;
2759566063dSJacob Faibussowitsch   for (p = pStart; p < pEnd; ++p) for (f = 0; f < Nf; ++f) PetscCall(PetscSectionSetFieldConstraintIndices(section, p, f, indices));
2769cbb8f0dSMatthew G. Knepley   /* Handle BC constraints */
2779cbb8f0dSMatthew G. Knepley   for (bc = 0; bc < numBC; ++bc) {
2789cbb8f0dSMatthew G. Knepley     const PetscInt  field = bcField[bc];
2799cbb8f0dSMatthew G. Knepley     const PetscInt *comp, *idx;
280ce093827SMatthew G. Knepley     PetscInt        Nc, cNc = -1, n, i;
2819cbb8f0dSMatthew G. Knepley 
2829566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldComponents(section, field, &Nc));
2839566063dSJacob Faibussowitsch     if (bcComps && bcComps[bc]) PetscCall(ISGetLocalSize(bcComps[bc], &cNc));
2849566063dSJacob Faibussowitsch     if (bcComps && bcComps[bc]) PetscCall(ISGetIndices(bcComps[bc], &comp));
28524cfc5f6SToby Isaac     if (bcPoints[bc]) {
2869566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(bcPoints[bc], &n));
2879566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(bcPoints[bc], &idx));
2889cbb8f0dSMatthew G. Knepley       for (i = 0; i < n; ++i) {
2899cbb8f0dSMatthew G. Knepley         const PetscInt  p = idx[i];
2909cbb8f0dSMatthew G. Knepley         const PetscInt *find;
291ce093827SMatthew G. Knepley         PetscInt        fdof, fcdof, c, j;
2929cbb8f0dSMatthew G. Knepley 
2939566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetFieldDof(section, p, field, &fdof));
2949cbb8f0dSMatthew G. Knepley         if (!fdof) continue;
295ce093827SMatthew G. Knepley         if (cNc < 0) {
2969cbb8f0dSMatthew G. Knepley           for (d = 0; d < fdof; ++d) indices[d] = d;
2979cbb8f0dSMatthew G. Knepley           fcdof = fdof;
2989cbb8f0dSMatthew G. Knepley         } else {
299ce093827SMatthew G. Knepley           /* We assume that a point may have multiple "nodes", which are collections of Nc dofs,
300ce093827SMatthew G. Knepley              and that those dofs are numbered n*Nc+c */
3019566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintDof(section, p, field, &fcdof));
3029566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintIndices(section, p, field, &find));
303ce093827SMatthew G. Knepley           /* Get indices constrained by previous bcs */
3049cbb8f0dSMatthew G. Knepley           for (d = 0; d < fcdof; ++d) {if (find[d] < 0) break; indices[d] = find[d];}
305ce093827SMatthew G. Knepley           for (j = 0; j < fdof/Nc; ++j) for (c = 0; c < cNc; ++c) indices[d++] = j*Nc + comp[c];
3069566063dSJacob Faibussowitsch           PetscCall(PetscSortRemoveDupsInt(&d, indices));
3079cbb8f0dSMatthew G. Knepley           for (c = d; c < fcdof; ++c) indices[c] = -1;
3089cbb8f0dSMatthew G. Knepley           fcdof = d;
3099cbb8f0dSMatthew G. Knepley         }
3109566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetFieldConstraintDof(section, p, field, fcdof));
3119566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetFieldConstraintIndices(section, p, field, indices));
3129cbb8f0dSMatthew G. Knepley       }
3139566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(bcPoints[bc], &idx));
3149cbb8f0dSMatthew G. Knepley     }
31524cfc5f6SToby Isaac     if (bcComps && bcComps[bc]) PetscCall(ISRestoreIndices(bcComps[bc], &comp));
31624cfc5f6SToby Isaac   }
3179cbb8f0dSMatthew G. Knepley   /* Handle anchors */
3189566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAnchors(dm, &aSec, NULL));
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;
3239566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd));
3249cbb8f0dSMatthew G. Knepley     for (a = aStart; a < aEnd; a++) {
3259cbb8f0dSMatthew G. Knepley       PetscInt dof, f;
3269cbb8f0dSMatthew G. Knepley 
3279566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(aSec, a, &dof));
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++) {
3319566063dSJacob Faibussowitsch           PetscCall(PetscSectionSetFieldConstraintIndices(section, a, f, indices));
3329cbb8f0dSMatthew G. Knepley         }
3339cbb8f0dSMatthew G. Knepley       }
3349cbb8f0dSMatthew G. Knepley     }
3359cbb8f0dSMatthew G. Knepley   }
3369566063dSJacob Faibussowitsch   PetscCall(PetscFree(indices));
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 
3469cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
3479566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &Nf));
3489566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetMaxDof(section, &maxDof));
3499566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
3509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxDof, &indices));
3519cbb8f0dSMatthew G. Knepley   for (d = 0; d < maxDof; ++d) indices[d] = -1;
3529cbb8f0dSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
3539cbb8f0dSMatthew G. Knepley     PetscInt cdof, d;
3549cbb8f0dSMatthew G. Knepley 
3559566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
3569cbb8f0dSMatthew G. Knepley     if (cdof) {
357baef929fSMatthew G. Knepley       if (Nf) {
3589cbb8f0dSMatthew G. Knepley         PetscInt numConst = 0, foff = 0;
3599cbb8f0dSMatthew G. Knepley 
360baef929fSMatthew G. Knepley         for (f = 0; f < Nf; ++f) {
3619cbb8f0dSMatthew G. Knepley           const PetscInt *find;
3629cbb8f0dSMatthew G. Knepley           PetscInt        fcdof, fdof;
3639cbb8f0dSMatthew G. Knepley 
3649566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldDof(section, p, f, &fdof));
3659566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintDof(section, p, f, &fcdof));
3669cbb8f0dSMatthew G. Knepley           /* Change constraint numbering from field component to local dof number */
3679566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintIndices(section, p, f, &find));
3689cbb8f0dSMatthew G. Knepley           for (d = 0; d < fcdof; ++d) indices[numConst+d] = find[d] + foff;
3699cbb8f0dSMatthew G. Knepley           numConst += fcdof;
3709cbb8f0dSMatthew G. Knepley           foff     += fdof;
3719cbb8f0dSMatthew G. Knepley         }
3729566063dSJacob Faibussowitsch         if (cdof != numConst) PetscCall(PetscSectionSetConstraintDof(section, p, numConst));
3739cbb8f0dSMatthew G. Knepley       } else {
3749cbb8f0dSMatthew G. Knepley         for (d = 0; d < cdof; ++d) indices[d] = d;
3759cbb8f0dSMatthew G. Knepley       }
3769566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetConstraintIndices(section, p, indices));
3779cbb8f0dSMatthew G. Knepley     }
3789cbb8f0dSMatthew G. Knepley   }
3799566063dSJacob Faibussowitsch   PetscCall(PetscFree(indices));
3809cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
3819cbb8f0dSMatthew G. Knepley }
3829cbb8f0dSMatthew G. Knepley 
3839cbb8f0dSMatthew G. Knepley /*@C
3849cbb8f0dSMatthew G. Knepley   DMPlexCreateSection - Create a PetscSection based upon the dof layout specification provided.
3859cbb8f0dSMatthew G. Knepley 
3869cbb8f0dSMatthew G. Knepley   Not Collective
3879cbb8f0dSMatthew G. Knepley 
3889cbb8f0dSMatthew G. Knepley   Input Parameters:
3899cbb8f0dSMatthew G. Knepley + dm        - The DMPlex object
39092cfd99aSMartin Diehl . label     - The label indicating the mesh support of each field, or NULL for the whole mesh
3919cbb8f0dSMatthew G. Knepley . numComp   - An array of size numFields that holds the number of components for each field
3929cbb8f0dSMatthew 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
3939cbb8f0dSMatthew G. Knepley . numBC     - The number of boundary conditions
3949cbb8f0dSMatthew G. Knepley . bcField   - An array of size numBC giving the field number for each boundry condition
3959cbb8f0dSMatthew G. Knepley . bcComps   - [Optional] An array of size numBC giving an IS holding the field components to which each boundary condition applies
3969cbb8f0dSMatthew G. Knepley . bcPoints  - An array of size numBC giving an IS holding the Plex points to which each boundary condition applies
3979cbb8f0dSMatthew G. Knepley - perm      - Optional permutation of the chart, or NULL
3989cbb8f0dSMatthew G. Knepley 
3999cbb8f0dSMatthew G. Knepley   Output Parameter:
4009cbb8f0dSMatthew G. Knepley . section - The PetscSection object
4019cbb8f0dSMatthew G. Knepley 
4029cbb8f0dSMatthew G. Knepley   Notes:
4039cbb8f0dSMatthew 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
4049cbb8f0dSMatthew G. Knepley   number of dof for field 0 on each edge.
4059cbb8f0dSMatthew G. Knepley 
4069cbb8f0dSMatthew G. Knepley   The chart permutation is the same one set using PetscSectionSetPermutation()
4079cbb8f0dSMatthew G. Knepley 
4089cbb8f0dSMatthew G. Knepley   Level: developer
4099cbb8f0dSMatthew G. Knepley 
4101bb6d2a8SBarry Smith   TODO: How is this related to DMCreateLocalSection()
4111bb6d2a8SBarry Smith 
412db781477SPatrick Sanan .seealso: `DMPlexCreate()`, `PetscSectionCreate()`, `PetscSectionSetPermutation()`
4139cbb8f0dSMatthew G. Knepley @*/
414baef929fSMatthew 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)
4159cbb8f0dSMatthew G. Knepley {
4169cbb8f0dSMatthew G. Knepley   PetscSection   aSec;
4179cbb8f0dSMatthew G. Knepley 
4189cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
4199566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSectionFields(dm, numComp, section));
4209566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSectionDof(dm, label, numDof, *section));
4219566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSectionBCDof(dm, numBC, bcField, bcComps, bcPoints, *section));
4229566063dSJacob Faibussowitsch   if (perm) PetscCall(PetscSectionSetPermutation(*section, perm));
4239566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFromOptions(*section));
4249566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(*section));
4259566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAnchors(dm,&aSec,NULL));
4269cbb8f0dSMatthew G. Knepley   if (numBC || aSec) {
4279566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSectionBCIndicesField(dm, numBC, bcField, bcComps, bcPoints, *section));
4289566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSectionBCIndices(dm, *section));
4299cbb8f0dSMatthew G. Knepley   }
4309566063dSJacob Faibussowitsch   PetscCall(PetscSectionViewFromOptions(*section,NULL,"-section_view"));
4319cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
4329cbb8f0dSMatthew G. Knepley }
4339cbb8f0dSMatthew G. Knepley 
434*f985c9a8SMatthew G. Knepley static PetscErrorCode DMCompleteBCLabels_Private(DM dm, const PetscBool isFE[])
435*f985c9a8SMatthew G. Knepley {
436*f985c9a8SMatthew G. Knepley   DMLabel     *labels, *glabels;
437*f985c9a8SMatthew G. Knepley   const char **names;
438*f985c9a8SMatthew G. Knepley   char        *sendNames, *recvNames;
439*f985c9a8SMatthew G. Knepley   PetscInt     Nds, s, maxLabels = 0, maxLen = 0, gmaxLen, Nl = 0, gNl, l, gl, m;
440*f985c9a8SMatthew G. Knepley   size_t       len;
441*f985c9a8SMatthew G. Knepley   MPI_Comm     comm;
442*f985c9a8SMatthew G. Knepley   PetscMPIInt  rank, size, p, *counts, *displs;
443*f985c9a8SMatthew G. Knepley 
444*f985c9a8SMatthew G. Knepley   PetscFunctionBegin;
445*f985c9a8SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
446*f985c9a8SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(comm, &size));
447*f985c9a8SMatthew G. Knepley   PetscCallMPI(MPI_Comm_rank(comm, &rank));
448*f985c9a8SMatthew G. Knepley   PetscCall(DMGetNumDS(dm, &Nds));
449*f985c9a8SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
450*f985c9a8SMatthew G. Knepley     PetscDS  dsBC;
451*f985c9a8SMatthew G. Knepley     PetscInt numBd;
452*f985c9a8SMatthew G. Knepley 
453*f985c9a8SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC));
454*f985c9a8SMatthew G. Knepley     PetscCall(PetscDSGetNumBoundary(dsBC, &numBd));
455*f985c9a8SMatthew G. Knepley     maxLabels += numBd;
456*f985c9a8SMatthew G. Knepley   }
457*f985c9a8SMatthew G. Knepley   PetscCall(PetscCalloc1(maxLabels, &labels));
458*f985c9a8SMatthew G. Knepley   /* Get list of labels to be completed */
459*f985c9a8SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
460*f985c9a8SMatthew G. Knepley     PetscDS  dsBC;
461*f985c9a8SMatthew G. Knepley     PetscInt numBd, bd;
462*f985c9a8SMatthew G. Knepley 
463*f985c9a8SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC));
464*f985c9a8SMatthew G. Knepley     PetscCall(PetscDSGetNumBoundary(dsBC, &numBd));
465*f985c9a8SMatthew G. Knepley     for (bd = 0; bd < numBd; ++bd) {
466*f985c9a8SMatthew G. Knepley       DMLabel  label;
467*f985c9a8SMatthew G. Knepley       PetscInt field;
468*f985c9a8SMatthew G. Knepley 
469*f985c9a8SMatthew G. Knepley       PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, NULL, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL));
470*f985c9a8SMatthew G. Knepley       if (!isFE[field] || !label) continue;
471*f985c9a8SMatthew G. Knepley       for (l = 0; l < Nl; ++l) if (labels[l] == label) break;
472*f985c9a8SMatthew G. Knepley       if (l == Nl) labels[Nl++] = label;
473*f985c9a8SMatthew G. Knepley     }
474*f985c9a8SMatthew G. Knepley   }
475*f985c9a8SMatthew G. Knepley   /* Get label names */
476*f985c9a8SMatthew G. Knepley   PetscCall(PetscMalloc1(Nl, &names));
477*f985c9a8SMatthew G. Knepley   for (l = 0; l < Nl; ++l) PetscCall(PetscObjectGetName((PetscObject) labels[l], &names[l]));
478*f985c9a8SMatthew G. Knepley   for (l = 0; l < Nl; ++l) {PetscCall(PetscStrlen(names[l], &len)); maxLen = PetscMax(maxLen, len+2);}
479*f985c9a8SMatthew G. Knepley   PetscCall(PetscFree(labels));
480*f985c9a8SMatthew G. Knepley   PetscCallMPI(MPI_Allreduce(&maxLen, &gmaxLen, 1, MPIU_INT, MPI_MAX, comm));
481*f985c9a8SMatthew G. Knepley   PetscCall(PetscMalloc1(Nl * gmaxLen, &sendNames));
482*f985c9a8SMatthew G. Knepley   for (l = 0; l < Nl; ++l) PetscCall(PetscStrcpy(&sendNames[gmaxLen*l], names[l]));
483*f985c9a8SMatthew G. Knepley   PetscCall(PetscFree(names));
484*f985c9a8SMatthew G. Knepley   /* Put all names on all processes */
485*f985c9a8SMatthew G. Knepley   PetscCall(PetscCalloc2(size, &counts, size+1, &displs));
486*f985c9a8SMatthew G. Knepley   PetscCallMPI(MPI_Allgather(&Nl, 1, MPIU_INT, counts, 1, MPIU_INT, comm));
487*f985c9a8SMatthew G. Knepley   for (p = 0; p < size; ++p) displs[p+1] = displs[p] + counts[p];
488*f985c9a8SMatthew G. Knepley   gNl = displs[size];
489*f985c9a8SMatthew G. Knepley   for (p = 0; p < size; ++p) {counts[p] *= gmaxLen; displs[p] *= gmaxLen;}
490*f985c9a8SMatthew G. Knepley   PetscCall(PetscMalloc2(gNl * gmaxLen, &recvNames, gNl, &glabels));
491*f985c9a8SMatthew G. Knepley   PetscCallMPI(MPI_Allgatherv(sendNames, counts[rank], MPI_CHAR, recvNames, counts, displs, MPI_CHAR, comm));
492*f985c9a8SMatthew G. Knepley   PetscCall(PetscFree2(counts, displs));
493*f985c9a8SMatthew G. Knepley   PetscCall(PetscFree(sendNames));
494*f985c9a8SMatthew G. Knepley   for (l = 0, gl = 0; l < gNl; ++l) {
495*f985c9a8SMatthew G. Knepley     PetscCall(DMGetLabel(dm, &recvNames[l*gmaxLen], &glabels[gl]));
496*f985c9a8SMatthew G. Knepley     PetscCheck(glabels[gl], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Label %s missing on rank %d", &recvNames[l*gmaxLen], rank);
497*f985c9a8SMatthew G. Knepley     for (m = 0; m < gl; ++m) if (glabels[m] == glabels[gl]) continue;
498*f985c9a8SMatthew G. Knepley     PetscCall(DMPlexLabelComplete(dm, glabels[gl]));
499*f985c9a8SMatthew G. Knepley     ++gl;
500*f985c9a8SMatthew G. Knepley   }
501*f985c9a8SMatthew G. Knepley   PetscCall(PetscFree2(recvNames, glabels));
502*f985c9a8SMatthew G. Knepley   PetscFunctionReturn(0);
503*f985c9a8SMatthew G. Knepley }
504*f985c9a8SMatthew G. Knepley 
5051bb6d2a8SBarry Smith PetscErrorCode DMCreateLocalSection_Plex(DM dm)
5069cbb8f0dSMatthew G. Knepley {
5079cbb8f0dSMatthew G. Knepley   PetscSection   section;
508baef929fSMatthew G. Knepley   DMLabel       *labels;
5099cbb8f0dSMatthew G. Knepley   IS            *bcPoints, *bcComps;
5109cbb8f0dSMatthew G. Knepley   PetscBool     *isFE;
5119cbb8f0dSMatthew G. Knepley   PetscInt      *bcFields, *numComp, *numDof;
5129310035eSMatthew G. Knepley   PetscInt       depth, dim, numBC = 0, Nf, Nds, s, bc = 0, f;
5139cbb8f0dSMatthew G. Knepley   PetscInt       cStart, cEnd, cEndInterior;
5149cbb8f0dSMatthew G. Knepley 
5159cbb8f0dSMatthew G. Knepley   PetscFunctionBegin;
5169566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
5179566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
5189566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
5199cbb8f0dSMatthew G. Knepley   /* FE and FV boundary conditions are handled slightly differently */
5209566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nf, &isFE));
521baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
5229cbb8f0dSMatthew G. Knepley     PetscObject  obj;
5239cbb8f0dSMatthew G. Knepley     PetscClassId id;
5249cbb8f0dSMatthew G. Knepley 
5259566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, f, NULL, &obj));
5269566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(obj, &id));
5279cbb8f0dSMatthew G. Knepley     if (id == PETSCFE_CLASSID)      {isFE[f] = PETSC_TRUE;}
5289cbb8f0dSMatthew G. Knepley     else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;}
52963a3b9bcSJacob Faibussowitsch     else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
5309cbb8f0dSMatthew G. Knepley   }
5319cbb8f0dSMatthew G. Knepley   /* Allocate boundary point storage for FEM boundaries */
5329566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
5339310035eSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
5349310035eSMatthew G. Knepley     PetscDS  dsBC;
5359310035eSMatthew G. Knepley     PetscInt numBd, bd;
5369310035eSMatthew G. Knepley 
5379566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC));
5389566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumBoundary(dsBC, &numBd));
53908401ef6SPierre 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)");
5409cbb8f0dSMatthew G. Knepley     for (bd = 0; bd < numBd; ++bd) {
5419cbb8f0dSMatthew G. Knepley       PetscInt                field;
5429cbb8f0dSMatthew G. Knepley       DMBoundaryConditionType type;
5439cbb8f0dSMatthew G. Knepley       DMLabel                 label;
5449cbb8f0dSMatthew G. Knepley 
5459566063dSJacob Faibussowitsch       PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, &type, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL));
5469cbb8f0dSMatthew G. Knepley       if (label && isFE[field] && (type & DM_BC_ESSENTIAL)) ++numBC;
5479cbb8f0dSMatthew G. Knepley     }
5489310035eSMatthew G. Knepley   }
5499cbb8f0dSMatthew G. Knepley   /* Add ghost cell boundaries for FVM */
5509566063dSJacob Faibussowitsch   PetscCall(DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL));
551baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) if (!isFE[f] && cEndInterior >= 0) ++numBC;
5529566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(numBC, &bcFields, numBC, &bcPoints, numBC, &bcComps));
5539cbb8f0dSMatthew G. Knepley   /* Constrain ghost cells for FV */
554baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
5559cbb8f0dSMatthew G. Knepley     PetscInt *newidx, c;
5569cbb8f0dSMatthew G. Knepley 
5579cbb8f0dSMatthew G. Knepley     if (isFE[f] || cEndInterior < 0) continue;
5589566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(cEnd-cEndInterior,&newidx));
5599cbb8f0dSMatthew G. Knepley     for (c = cEndInterior; c < cEnd; ++c) newidx[c-cEndInterior] = c;
5609cbb8f0dSMatthew G. Knepley     bcFields[bc] = f;
5619566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, cEnd-cEndInterior, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]));
5629cbb8f0dSMatthew G. Knepley   }
563*f985c9a8SMatthew G. Knepley   /* Complete labels for boundary conditions */
564*f985c9a8SMatthew G. Knepley   PetscCall(DMCompleteBCLabels_Private(dm, isFE));
5659cbb8f0dSMatthew G. Knepley   /* Handle FEM Dirichlet boundaries */
5669566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
5679310035eSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
5689310035eSMatthew G. Knepley     PetscDS  dsBC;
5699310035eSMatthew G. Knepley     PetscInt numBd, bd;
5709310035eSMatthew G. Knepley 
5719566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC));
5729566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumBoundary(dsBC, &numBd));
5739cbb8f0dSMatthew G. Knepley     for (bd = 0; bd < numBd; ++bd) {
5749cbb8f0dSMatthew G. Knepley       DMLabel                 label;
5759cbb8f0dSMatthew G. Knepley       const PetscInt         *comps;
5769cbb8f0dSMatthew G. Knepley       const PetscInt         *values;
5779310035eSMatthew G. Knepley       PetscInt                bd2, field, numComps, numValues;
5789cbb8f0dSMatthew G. Knepley       DMBoundaryConditionType type;
5799310035eSMatthew G. Knepley       PetscBool               duplicate = PETSC_FALSE;
5809cbb8f0dSMatthew G. Knepley 
5819566063dSJacob Faibussowitsch       PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, &type, NULL, &label, &numValues, &values, &field, &numComps, &comps, NULL, NULL, NULL));
5829cbb8f0dSMatthew G. Knepley       if (!isFE[field] || !label) continue;
5839310035eSMatthew G. Knepley       /* Only want to modify label once */
5849310035eSMatthew G. Knepley       for (bd2 = 0; bd2 < bd; ++bd2) {
58545480ffeSMatthew G. Knepley         DMLabel l;
58645480ffeSMatthew G. Knepley 
5879566063dSJacob Faibussowitsch         PetscCall(PetscDSGetBoundary(dsBC, bd2, NULL, NULL, NULL, &l, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
58845480ffeSMatthew G. Knepley         duplicate = l == label ? PETSC_TRUE : PETSC_FALSE;
5899310035eSMatthew G. Knepley         if (duplicate) break;
5909310035eSMatthew G. Knepley       }
5919cbb8f0dSMatthew G. Knepley       /* Filter out cells, if you actually want to constrain cells you need to do things by hand right now */
5929cbb8f0dSMatthew G. Knepley       if (type & DM_BC_ESSENTIAL) {
5939cbb8f0dSMatthew G. Knepley         PetscInt       *newidx;
5949cbb8f0dSMatthew G. Knepley         PetscInt        n, newn = 0, p, v;
5959cbb8f0dSMatthew G. Knepley 
5969cbb8f0dSMatthew G. Knepley         bcFields[bc] = field;
5979566063dSJacob Faibussowitsch         if (numComps) PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject) dm), numComps, comps, PETSC_COPY_VALUES, &bcComps[bc]));
5989cbb8f0dSMatthew G. Knepley         for (v = 0; v < numValues; ++v) {
5999cbb8f0dSMatthew G. Knepley           IS              tmp;
6009cbb8f0dSMatthew G. Knepley           const PetscInt *idx;
6019cbb8f0dSMatthew G. Knepley 
6029566063dSJacob Faibussowitsch           PetscCall(DMLabelGetStratumIS(label, values[v], &tmp));
6039cbb8f0dSMatthew G. Knepley           if (!tmp) continue;
6049566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(tmp, &n));
6059566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(tmp, &idx));
6069cbb8f0dSMatthew G. Knepley           if (isFE[field]) {
6079cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) ++newn;
6089cbb8f0dSMatthew G. Knepley           } else {
6099cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) ++newn;
6109cbb8f0dSMatthew G. Knepley           }
6119566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(tmp, &idx));
6129566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&tmp));
6139cbb8f0dSMatthew G. Knepley         }
6149566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(newn, &newidx));
6159cbb8f0dSMatthew G. Knepley         newn = 0;
6169cbb8f0dSMatthew G. Knepley         for (v = 0; v < numValues; ++v) {
6179cbb8f0dSMatthew G. Knepley           IS              tmp;
6189cbb8f0dSMatthew G. Knepley           const PetscInt *idx;
6199cbb8f0dSMatthew G. Knepley 
6209566063dSJacob Faibussowitsch           PetscCall(DMLabelGetStratumIS(label, values[v], &tmp));
6219cbb8f0dSMatthew G. Knepley           if (!tmp) continue;
6229566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(tmp, &n));
6239566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(tmp, &idx));
6249cbb8f0dSMatthew G. Knepley           if (isFE[field]) {
6259cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) newidx[newn++] = idx[p];
6269cbb8f0dSMatthew G. Knepley           } else {
6279cbb8f0dSMatthew G. Knepley             for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) newidx[newn++] = idx[p];
6289cbb8f0dSMatthew G. Knepley           }
6299566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(tmp, &idx));
6309566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&tmp));
6319cbb8f0dSMatthew G. Knepley         }
6329566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, newn, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]));
6339310035eSMatthew G. Knepley       }
6349cbb8f0dSMatthew G. Knepley     }
6359cbb8f0dSMatthew G. Knepley   }
6369cbb8f0dSMatthew G. Knepley   /* Handle discretization */
6379566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(Nf,&labels,Nf,&numComp,Nf*(dim+1),&numDof));
638baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
63944a7f3ddSMatthew G. Knepley     labels[f] = dm->fields[f].label;
6409cbb8f0dSMatthew G. Knepley     if (isFE[f]) {
64144a7f3ddSMatthew G. Knepley       PetscFE         fe = (PetscFE) dm->fields[f].disc;
6429cbb8f0dSMatthew G. Knepley       const PetscInt *numFieldDof;
643a3254110SMatthew Knepley       PetscInt        fedim, d;
6449cbb8f0dSMatthew G. Knepley 
6459566063dSJacob Faibussowitsch       PetscCall(PetscFEGetNumComponents(fe, &numComp[f]));
6469566063dSJacob Faibussowitsch       PetscCall(PetscFEGetNumDof(fe, &numFieldDof));
6479566063dSJacob Faibussowitsch       PetscCall(PetscFEGetSpatialDimension(fe, &fedim));
648a3254110SMatthew Knepley       for (d = 0; d < PetscMin(dim, fedim)+1; ++d) numDof[f*(dim+1)+d] = numFieldDof[d];
6499cbb8f0dSMatthew G. Knepley     } else {
65044a7f3ddSMatthew G. Knepley       PetscFV fv = (PetscFV) dm->fields[f].disc;
6519cbb8f0dSMatthew G. Knepley 
6529566063dSJacob Faibussowitsch       PetscCall(PetscFVGetNumComponents(fv, &numComp[f]));
6539cbb8f0dSMatthew G. Knepley       numDof[f*(dim+1)+dim] = numComp[f];
6549cbb8f0dSMatthew G. Knepley     }
6559cbb8f0dSMatthew G. Knepley   }
6569566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
657baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
6589cbb8f0dSMatthew G. Knepley     PetscInt d;
6599cbb8f0dSMatthew G. Knepley     for (d = 1; d < dim; ++d) {
6601dca8a05SBarry Smith       PetscCheck(numDof[f*(dim+1)+d] <= 0 || depth >= dim,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated when unknowns are specified on edges or faces.");
6619cbb8f0dSMatthew G. Knepley     }
6629cbb8f0dSMatthew G. Knepley   }
6639566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSection(dm, labels, numComp, numDof, numBC, bcFields, bcComps, bcPoints, NULL, &section));
664baef929fSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
6659cbb8f0dSMatthew G. Knepley     PetscFE     fe;
6669cbb8f0dSMatthew G. Knepley     const char *name;
6679cbb8f0dSMatthew G. Knepley 
6689566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, f, NULL, (PetscObject *) &fe));
669cdfe53dfSMatthew G. Knepley     if (!((PetscObject) fe)->name) continue;
6709566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject) fe, &name));
6719566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldName(section, f, name));
6729cbb8f0dSMatthew G. Knepley   }
6739566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, section));
6749566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&section));
6759566063dSJacob Faibussowitsch   for (bc = 0; bc < numBC; ++bc) {PetscCall(ISDestroy(&bcPoints[bc]));PetscCall(ISDestroy(&bcComps[bc]));}
6769566063dSJacob Faibussowitsch   PetscCall(PetscFree3(bcFields,bcPoints,bcComps));
6779566063dSJacob Faibussowitsch   PetscCall(PetscFree3(labels,numComp,numDof));
6789566063dSJacob Faibussowitsch   PetscCall(PetscFree(isFE));
6799cbb8f0dSMatthew G. Knepley   PetscFunctionReturn(0);
6809cbb8f0dSMatthew G. Knepley }
681