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 */ 49371c9d4SSatish Balay static PetscErrorCode DMPlexCreateSectionFields(DM dm, const PetscInt numComp[], PetscSection *section) { 5baef929fSMatthew G. Knepley DMLabel depthLabel; 6baef929fSMatthew G. Knepley PetscInt depth, Nf, f, pStart, pEnd; 79cbb8f0dSMatthew G. Knepley PetscBool *isFE; 89cbb8f0dSMatthew G. Knepley 99cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 109566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 119566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(dm, &depthLabel)); 129566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 139566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nf, &isFE)); 14baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 159cbb8f0dSMatthew G. Knepley PetscObject obj; 169cbb8f0dSMatthew G. Knepley PetscClassId id; 179cbb8f0dSMatthew G. Knepley 189566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, f, NULL, &obj)); 199566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(obj, &id)); 209371c9d4SSatish Balay if (id == PETSCFE_CLASSID) { 219371c9d4SSatish Balay isFE[f] = PETSC_TRUE; 229371c9d4SSatish Balay } else if (id == PETSCFV_CLASSID) { 239371c9d4SSatish Balay isFE[f] = PETSC_FALSE; 249371c9d4SSatish Balay } 259cbb8f0dSMatthew G. Knepley } 26baef929fSMatthew G. Knepley 279566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), section)); 28baef929fSMatthew G. Knepley if (Nf > 0) { 299566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(*section, Nf)); 309cbb8f0dSMatthew G. Knepley if (numComp) { 31baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 329566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(*section, f, numComp[f])); 339cbb8f0dSMatthew G. Knepley if (isFE[f]) { 349cbb8f0dSMatthew G. Knepley PetscFE fe; 359cbb8f0dSMatthew G. Knepley PetscDualSpace dspace; 369cbb8f0dSMatthew G. Knepley const PetscInt ***perms; 379cbb8f0dSMatthew G. Knepley const PetscScalar ***flips; 389cbb8f0dSMatthew G. Knepley const PetscInt *numDof; 399cbb8f0dSMatthew G. Knepley 409566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, f, NULL, (PetscObject *)&fe)); 419566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(fe, &dspace)); 429566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetSymmetries(dspace, &perms, &flips)); 439566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetNumDof(dspace, &numDof)); 449cbb8f0dSMatthew G. Knepley if (perms || flips) { 459cbb8f0dSMatthew G. Knepley DM K; 464dff28b8SMatthew G. Knepley PetscInt sph, spdepth; 479cbb8f0dSMatthew G. Knepley PetscSectionSym sym; 489cbb8f0dSMatthew G. Knepley 499566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(dspace, &K)); 509566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(K, &spdepth)); 519566063dSJacob Faibussowitsch PetscCall(PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)*section), depthLabel, &sym)); 524dff28b8SMatthew G. Knepley for (sph = 0; sph <= spdepth; sph++) { 539cbb8f0dSMatthew G. Knepley PetscDualSpace hspace; 549cbb8f0dSMatthew G. Knepley PetscInt kStart, kEnd; 554dff28b8SMatthew G. Knepley PetscInt kConeSize, h = sph + (depth - spdepth); 569cbb8f0dSMatthew G. Knepley const PetscInt **perms0 = NULL; 579cbb8f0dSMatthew G. Knepley const PetscScalar **flips0 = NULL; 589cbb8f0dSMatthew G. Knepley 599566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetHeightSubspace(dspace, sph, &hspace)); 609566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(K, h, &kStart, &kEnd)); 619cbb8f0dSMatthew G. Knepley if (!hspace) continue; 629566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetSymmetries(hspace, &perms, &flips)); 639cbb8f0dSMatthew G. Knepley if (perms) perms0 = perms[0]; 649cbb8f0dSMatthew G. Knepley if (flips) flips0 = flips[0]; 659cbb8f0dSMatthew G. Knepley if (!(perms0 || flips0)) continue; 66b5a892a1SMatthew G. Knepley { 67b5a892a1SMatthew G. Knepley DMPolytopeType ct; 68b5a892a1SMatthew G. Knepley /* The number of arrangements is no longer based on the number of faces */ 699566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(K, kStart, &ct)); 70b5a892a1SMatthew G. Knepley kConeSize = DMPolytopeTypeGetNumArrangments(ct) / 2; 71b5a892a1SMatthew G. Knepley } 729566063dSJacob Faibussowitsch PetscCall(PetscSectionSymLabelSetStratum(sym, depth - h, numDof[depth - h], -kConeSize, kConeSize, PETSC_USE_POINTER, perms0 ? &perms0[-kConeSize] : NULL, flips0 ? &flips0[-kConeSize] : NULL)); 739cbb8f0dSMatthew G. Knepley } 749566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldSym(*section, f, sym)); 759566063dSJacob Faibussowitsch PetscCall(PetscSectionSymDestroy(&sym)); 769cbb8f0dSMatthew G. Knepley } 779cbb8f0dSMatthew G. Knepley } 789cbb8f0dSMatthew G. Knepley } 799cbb8f0dSMatthew G. Knepley } 809cbb8f0dSMatthew G. Knepley } 819566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 829566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(*section, pStart, pEnd)); 839566063dSJacob Faibussowitsch PetscCall(PetscFree(isFE)); 84baef929fSMatthew G. Knepley PetscFunctionReturn(0); 85baef929fSMatthew G. Knepley } 86baef929fSMatthew G. Knepley 87baef929fSMatthew G. Knepley /* Set the number of dof on each point and separate by fields */ 889371c9d4SSatish Balay static PetscErrorCode DMPlexCreateSectionDof(DM dm, DMLabel label[], const PetscInt numDof[], PetscSection section) { 89baef929fSMatthew G. Knepley DMLabel depthLabel; 90412e9a14SMatthew G. Knepley DMPolytopeType ct; 91baef929fSMatthew G. Knepley PetscInt depth, cellHeight, pStart = 0, pEnd = 0; 92a55f9a55SMatthew G. Knepley PetscInt Nf, f, Nds, n, dim, d, dep, p; 935fedec97SMatthew G. Knepley PetscBool *isFE, hasCohesive = PETSC_FALSE; 94baef929fSMatthew G. Knepley 95baef929fSMatthew G. Knepley PetscFunctionBegin; 969566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 979566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 989566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(dm, &depthLabel)); 999566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 1009566063dSJacob Faibussowitsch PetscCall(DMGetNumDS(dm, &Nds)); 101a55f9a55SMatthew G. Knepley for (n = 0; n < Nds; ++n) { 102a55f9a55SMatthew G. Knepley PetscDS ds; 1035fedec97SMatthew G. Knepley PetscBool isCohesive; 104a55f9a55SMatthew G. Knepley 1059566063dSJacob Faibussowitsch PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds)); 1069566063dSJacob Faibussowitsch PetscCall(PetscDSIsCohesive(ds, &isCohesive)); 1079371c9d4SSatish Balay if (isCohesive) { 1089371c9d4SSatish Balay hasCohesive = PETSC_TRUE; 1099371c9d4SSatish Balay break; 1109371c9d4SSatish Balay } 111a55f9a55SMatthew G. Knepley } 1129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nf, &isFE)); 113baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 114baef929fSMatthew G. Knepley PetscObject obj; 115baef929fSMatthew G. Knepley PetscClassId id; 116baef929fSMatthew G. Knepley 1179566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, f, NULL, &obj)); 1189566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(obj, &id)); 1191274e1a1SLawrence Mitchell /* User is allowed to put a "placeholder" field in (c.f. DMCreateDS) */ 1201274e1a1SLawrence Mitchell isFE[f] = id == PETSCFE_CLASSID ? PETSC_TRUE : PETSC_FALSE; 121baef929fSMatthew G. Knepley } 122baef929fSMatthew G. Knepley 1239566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 124baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 125e0b68406SMatthew Knepley PetscBool avoidTensor; 126e0b68406SMatthew Knepley 1279566063dSJacob Faibussowitsch PetscCall(DMGetFieldAvoidTensor(dm, f, &avoidTensor)); 1285fedec97SMatthew G. Knepley avoidTensor = (avoidTensor || hasCohesive) ? PETSC_TRUE : PETSC_FALSE; 129baef929fSMatthew G. Knepley if (label && label[f]) { 130baef929fSMatthew G. Knepley IS pointIS; 131baef929fSMatthew G. Knepley const PetscInt *points; 132baef929fSMatthew G. Knepley PetscInt n; 133baef929fSMatthew G. Knepley 1349566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label[f], 1, &pointIS)); 135baef929fSMatthew G. Knepley if (!pointIS) continue; 1369566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointIS, &n)); 1379566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pointIS, &points)); 138baef929fSMatthew G. Knepley for (p = 0; p < n; ++p) { 139baef929fSMatthew G. Knepley const PetscInt point = points[p]; 140baef929fSMatthew G. Knepley PetscInt dof, d; 141baef929fSMatthew G. Knepley 1429566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, point, &ct)); 1439566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(depthLabel, point, &d)); 144412e9a14SMatthew G. Knepley /* If this is a tensor prism point, use dof for one dimension lower */ 145412e9a14SMatthew G. Knepley switch (ct) { 146412e9a14SMatthew G. Knepley case DM_POLYTOPE_POINT_PRISM_TENSOR: 147412e9a14SMatthew G. Knepley case DM_POLYTOPE_SEG_PRISM_TENSOR: 148412e9a14SMatthew G. Knepley case DM_POLYTOPE_TRI_PRISM_TENSOR: 149412e9a14SMatthew G. Knepley case DM_POLYTOPE_QUAD_PRISM_TENSOR: 1509371c9d4SSatish Balay if (hasCohesive) { --d; } 1519371c9d4SSatish Balay break; 152412e9a14SMatthew G. Knepley default: break; 153412e9a14SMatthew G. Knepley } 154baef929fSMatthew G. Knepley dof = d < 0 ? 0 : numDof[f * (dim + 1) + d]; 1559566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(section, point, f, dof)); 1569566063dSJacob Faibussowitsch PetscCall(PetscSectionAddDof(section, point, dof)); 157baef929fSMatthew G. Knepley } 1589566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pointIS, &points)); 1599566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointIS)); 160baef929fSMatthew G. Knepley } else { 1619cbb8f0dSMatthew G. Knepley for (dep = 0; dep <= depth - cellHeight; ++dep) { 1623fd2e703SMatthew G. Knepley /* Cases: dim > depth (cell-vertex mesh), dim == depth (fully interpolated), dim < depth (interpolated submesh) */ 1633fd2e703SMatthew G. Knepley d = dim <= depth ? dep : (!dep ? 0 : dim); 1649566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, dep, &pStart, &pEnd)); 1659cbb8f0dSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 166baef929fSMatthew G. Knepley const PetscInt dof = numDof[f * (dim + 1) + d]; 1679cbb8f0dSMatthew G. Knepley 1689566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, p, &ct)); 169412e9a14SMatthew G. Knepley switch (ct) { 170412e9a14SMatthew G. Knepley case DM_POLYTOPE_POINT_PRISM_TENSOR: 171412e9a14SMatthew G. Knepley case DM_POLYTOPE_SEG_PRISM_TENSOR: 172412e9a14SMatthew G. Knepley case DM_POLYTOPE_TRI_PRISM_TENSOR: 173412e9a14SMatthew G. Knepley case DM_POLYTOPE_QUAD_PRISM_TENSOR: 174e0b68406SMatthew Knepley if (avoidTensor && isFE[f]) continue; 175412e9a14SMatthew G. Knepley default: break; 176412e9a14SMatthew G. Knepley } 1779566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(section, p, f, dof)); 1789566063dSJacob Faibussowitsch PetscCall(PetscSectionAddDof(section, p, dof)); 1799cbb8f0dSMatthew G. Knepley } 180baef929fSMatthew G. Knepley } 1819cbb8f0dSMatthew G. Knepley } 1829cbb8f0dSMatthew G. Knepley } 1839566063dSJacob Faibussowitsch PetscCall(PetscFree(isFE)); 1849cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 1859cbb8f0dSMatthew G. Knepley } 1869cbb8f0dSMatthew G. Knepley 1879cbb8f0dSMatthew G. Knepley /* Set the number of dof on each point and separate by fields 1889cbb8f0dSMatthew G. Knepley If bcComps is NULL or the IS is NULL, constrain every dof on the point 1899cbb8f0dSMatthew G. Knepley */ 1909371c9d4SSatish Balay static PetscErrorCode DMPlexCreateSectionBCDof(DM dm, PetscInt numBC, const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section) { 191baef929fSMatthew G. Knepley PetscInt Nf; 1929cbb8f0dSMatthew G. Knepley PetscInt bc; 1939cbb8f0dSMatthew G. Knepley PetscSection aSec; 1949cbb8f0dSMatthew G. Knepley 1959cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 1969566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &Nf)); 1979cbb8f0dSMatthew G. Knepley for (bc = 0; bc < numBC; ++bc) { 1989cbb8f0dSMatthew G. Knepley PetscInt field = 0; 1999cbb8f0dSMatthew G. Knepley const PetscInt *comp; 2009cbb8f0dSMatthew G. Knepley const PetscInt *idx; 201ce093827SMatthew G. Knepley PetscInt Nc = 0, cNc = -1, n, i; 2029cbb8f0dSMatthew G. Knepley 203ce093827SMatthew G. Knepley if (Nf) { 204ce093827SMatthew G. Knepley field = bcField[bc]; 2059566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldComponents(section, field, &Nc)); 206ce093827SMatthew G. Knepley } 2079566063dSJacob Faibussowitsch if (bcComps && bcComps[bc]) PetscCall(ISGetLocalSize(bcComps[bc], &cNc)); 2089566063dSJacob Faibussowitsch if (bcComps && bcComps[bc]) PetscCall(ISGetIndices(bcComps[bc], &comp)); 20924cfc5f6SToby Isaac if (bcPoints[bc]) { 2109566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(bcPoints[bc], &n)); 2119566063dSJacob Faibussowitsch PetscCall(ISGetIndices(bcPoints[bc], &idx)); 2129cbb8f0dSMatthew G. Knepley for (i = 0; i < n; ++i) { 2139cbb8f0dSMatthew G. Knepley const PetscInt p = idx[i]; 2149cbb8f0dSMatthew G. Knepley PetscInt numConst; 2159cbb8f0dSMatthew G. Knepley 216baef929fSMatthew G. Knepley if (Nf) { 2179566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, p, field, &numConst)); 2189cbb8f0dSMatthew G. Knepley } else { 2199566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &numConst)); 2209cbb8f0dSMatthew G. Knepley } 221ce093827SMatthew G. Knepley /* If Nc <= 0, constrain every dof on the point */ 222ce093827SMatthew G. Knepley if (cNc > 0) { 223ce093827SMatthew G. Knepley /* We assume that a point may have multiple "nodes", which are collections of Nc dofs, 224ce093827SMatthew G. Knepley and that those dofs are numbered n*Nc+c */ 225ce093827SMatthew G. Knepley if (Nf) { 22663a3b9bcSJacob 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); 227ce093827SMatthew G. Knepley numConst = (numConst / Nc) * cNc; 228ce093827SMatthew G. Knepley } else { 229ce093827SMatthew G. Knepley numConst = PetscMin(numConst, cNc); 230ce093827SMatthew G. Knepley } 231ce093827SMatthew G. Knepley } 2329566063dSJacob Faibussowitsch if (Nf) PetscCall(PetscSectionAddFieldConstraintDof(section, p, field, numConst)); 2339566063dSJacob Faibussowitsch PetscCall(PetscSectionAddConstraintDof(section, p, numConst)); 2349cbb8f0dSMatthew G. Knepley } 2359566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(bcPoints[bc], &idx)); 23624cfc5f6SToby Isaac } 2379566063dSJacob Faibussowitsch if (bcComps && bcComps[bc]) PetscCall(ISRestoreIndices(bcComps[bc], &comp)); 2389cbb8f0dSMatthew G. Knepley } 2399566063dSJacob Faibussowitsch PetscCall(DMPlexGetAnchors(dm, &aSec, NULL)); 2409cbb8f0dSMatthew G. Knepley if (aSec) { 2419cbb8f0dSMatthew G. Knepley PetscInt aStart, aEnd, a; 2429cbb8f0dSMatthew G. Knepley 2439566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd)); 2449cbb8f0dSMatthew G. Knepley for (a = aStart; a < aEnd; a++) { 2459cbb8f0dSMatthew G. Knepley PetscInt dof, f; 2469cbb8f0dSMatthew G. Knepley 2479566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(aSec, a, &dof)); 2489cbb8f0dSMatthew G. Knepley if (dof) { 2499cbb8f0dSMatthew G. Knepley /* if there are point-to-point constraints, then all dofs are constrained */ 2509566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, a, &dof)); 2519566063dSJacob Faibussowitsch PetscCall(PetscSectionSetConstraintDof(section, a, dof)); 252baef929fSMatthew G. Knepley for (f = 0; f < Nf; f++) { 2539566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, a, f, &dof)); 2549566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldConstraintDof(section, a, f, dof)); 2559cbb8f0dSMatthew G. Knepley } 2569cbb8f0dSMatthew G. Knepley } 2579cbb8f0dSMatthew G. Knepley } 2589cbb8f0dSMatthew G. Knepley } 2599cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 2609cbb8f0dSMatthew G. Knepley } 2619cbb8f0dSMatthew G. Knepley 2629cbb8f0dSMatthew G. Knepley /* Set the constrained field indices on each point 2639cbb8f0dSMatthew G. Knepley If bcComps is NULL or the IS is NULL, constrain every dof on the point 2649cbb8f0dSMatthew G. Knepley */ 2659371c9d4SSatish Balay static PetscErrorCode DMPlexCreateSectionBCIndicesField(DM dm, PetscInt numBC, const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section) { 2669cbb8f0dSMatthew G. Knepley PetscSection aSec; 2679cbb8f0dSMatthew G. Knepley PetscInt *indices; 268baef929fSMatthew G. Knepley PetscInt Nf, cdof, maxDof = 0, pStart, pEnd, p, bc, f, d; 2699cbb8f0dSMatthew G. Knepley 2709cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 2719566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &Nf)); 272baef929fSMatthew G. Knepley if (!Nf) PetscFunctionReturn(0); 2739cbb8f0dSMatthew G. Knepley /* Initialize all field indices to -1 */ 2749566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 2759371c9d4SSatish Balay for (p = pStart; p < pEnd; ++p) { 2769371c9d4SSatish Balay PetscCall(PetscSectionGetConstraintDof(section, p, &cdof)); 2779371c9d4SSatish Balay maxDof = PetscMax(maxDof, cdof); 2789371c9d4SSatish Balay } 2799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxDof, &indices)); 2809cbb8f0dSMatthew G. Knepley for (d = 0; d < maxDof; ++d) indices[d] = -1; 2819371c9d4SSatish Balay for (p = pStart; p < pEnd; ++p) 2829371c9d4SSatish Balay for (f = 0; f < Nf; ++f) PetscCall(PetscSectionSetFieldConstraintIndices(section, p, f, indices)); 2839cbb8f0dSMatthew G. Knepley /* Handle BC constraints */ 2849cbb8f0dSMatthew G. Knepley for (bc = 0; bc < numBC; ++bc) { 2859cbb8f0dSMatthew G. Knepley const PetscInt field = bcField[bc]; 2869cbb8f0dSMatthew G. Knepley const PetscInt *comp, *idx; 287ce093827SMatthew G. Knepley PetscInt Nc, cNc = -1, n, i; 2889cbb8f0dSMatthew G. Knepley 2899566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldComponents(section, field, &Nc)); 2909566063dSJacob Faibussowitsch if (bcComps && bcComps[bc]) PetscCall(ISGetLocalSize(bcComps[bc], &cNc)); 2919566063dSJacob Faibussowitsch if (bcComps && bcComps[bc]) PetscCall(ISGetIndices(bcComps[bc], &comp)); 29224cfc5f6SToby Isaac if (bcPoints[bc]) { 2939566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(bcPoints[bc], &n)); 2949566063dSJacob Faibussowitsch PetscCall(ISGetIndices(bcPoints[bc], &idx)); 2959cbb8f0dSMatthew G. Knepley for (i = 0; i < n; ++i) { 2969cbb8f0dSMatthew G. Knepley const PetscInt p = idx[i]; 2979cbb8f0dSMatthew G. Knepley const PetscInt *find; 298ce093827SMatthew G. Knepley PetscInt fdof, fcdof, c, j; 2999cbb8f0dSMatthew G. Knepley 3009566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, p, field, &fdof)); 3019cbb8f0dSMatthew G. Knepley if (!fdof) continue; 302ce093827SMatthew G. Knepley if (cNc < 0) { 3039cbb8f0dSMatthew G. Knepley for (d = 0; d < fdof; ++d) indices[d] = d; 3049cbb8f0dSMatthew G. Knepley fcdof = fdof; 3059cbb8f0dSMatthew G. Knepley } else { 306ce093827SMatthew G. Knepley /* We assume that a point may have multiple "nodes", which are collections of Nc dofs, 307ce093827SMatthew G. Knepley and that those dofs are numbered n*Nc+c */ 3089566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(section, p, field, &fcdof)); 3099566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintIndices(section, p, field, &find)); 310ce093827SMatthew G. Knepley /* Get indices constrained by previous bcs */ 3119371c9d4SSatish Balay for (d = 0; d < fcdof; ++d) { 3129371c9d4SSatish Balay if (find[d] < 0) break; 3139371c9d4SSatish Balay indices[d] = find[d]; 3149371c9d4SSatish Balay } 3159371c9d4SSatish Balay for (j = 0; j < fdof / Nc; ++j) 3169371c9d4SSatish Balay for (c = 0; c < cNc; ++c) indices[d++] = j * Nc + comp[c]; 3179566063dSJacob Faibussowitsch PetscCall(PetscSortRemoveDupsInt(&d, indices)); 3189cbb8f0dSMatthew G. Knepley for (c = d; c < fcdof; ++c) indices[c] = -1; 3199cbb8f0dSMatthew G. Knepley fcdof = d; 3209cbb8f0dSMatthew G. Knepley } 3219566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldConstraintDof(section, p, field, fcdof)); 3229566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldConstraintIndices(section, p, field, indices)); 3239cbb8f0dSMatthew G. Knepley } 3249566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(bcPoints[bc], &idx)); 3259cbb8f0dSMatthew G. Knepley } 32624cfc5f6SToby Isaac if (bcComps && bcComps[bc]) PetscCall(ISRestoreIndices(bcComps[bc], &comp)); 32724cfc5f6SToby Isaac } 3289cbb8f0dSMatthew G. Knepley /* Handle anchors */ 3299566063dSJacob Faibussowitsch PetscCall(DMPlexGetAnchors(dm, &aSec, NULL)); 3309cbb8f0dSMatthew G. Knepley if (aSec) { 3319cbb8f0dSMatthew G. Knepley PetscInt aStart, aEnd, a; 3329cbb8f0dSMatthew G. Knepley 3339cbb8f0dSMatthew G. Knepley for (d = 0; d < maxDof; ++d) indices[d] = d; 3349566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd)); 3359cbb8f0dSMatthew G. Knepley for (a = aStart; a < aEnd; a++) { 3369cbb8f0dSMatthew G. Knepley PetscInt dof, f; 3379cbb8f0dSMatthew G. Knepley 3389566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(aSec, a, &dof)); 3399cbb8f0dSMatthew G. Knepley if (dof) { 3409cbb8f0dSMatthew G. Knepley /* if there are point-to-point constraints, then all dofs are constrained */ 341*48a46eb9SPierre Jolivet for (f = 0; f < Nf; f++) PetscCall(PetscSectionSetFieldConstraintIndices(section, a, f, indices)); 3429cbb8f0dSMatthew G. Knepley } 3439cbb8f0dSMatthew G. Knepley } 3449cbb8f0dSMatthew G. Knepley } 3459566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 3469cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 3479cbb8f0dSMatthew G. Knepley } 3489cbb8f0dSMatthew G. Knepley 3499cbb8f0dSMatthew G. Knepley /* Set the constrained indices on each point */ 3509371c9d4SSatish Balay static PetscErrorCode DMPlexCreateSectionBCIndices(DM dm, PetscSection section) { 3519cbb8f0dSMatthew G. Knepley PetscInt *indices; 352baef929fSMatthew G. Knepley PetscInt Nf, maxDof, pStart, pEnd, p, f, d; 3539cbb8f0dSMatthew G. Knepley 3549cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 3559566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &Nf)); 3569566063dSJacob Faibussowitsch PetscCall(PetscSectionGetMaxDof(section, &maxDof)); 3579566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 3589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxDof, &indices)); 3599cbb8f0dSMatthew G. Knepley for (d = 0; d < maxDof; ++d) indices[d] = -1; 3609cbb8f0dSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 3619cbb8f0dSMatthew G. Knepley PetscInt cdof, d; 3629cbb8f0dSMatthew G. Knepley 3639566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, p, &cdof)); 3649cbb8f0dSMatthew G. Knepley if (cdof) { 365baef929fSMatthew G. Knepley if (Nf) { 3669cbb8f0dSMatthew G. Knepley PetscInt numConst = 0, foff = 0; 3679cbb8f0dSMatthew G. Knepley 368baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 3699cbb8f0dSMatthew G. Knepley const PetscInt *find; 3709cbb8f0dSMatthew G. Knepley PetscInt fcdof, fdof; 3719cbb8f0dSMatthew G. Knepley 3729566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, p, f, &fdof)); 3739566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(section, p, f, &fcdof)); 3749cbb8f0dSMatthew G. Knepley /* Change constraint numbering from field component to local dof number */ 3759566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintIndices(section, p, f, &find)); 3769cbb8f0dSMatthew G. Knepley for (d = 0; d < fcdof; ++d) indices[numConst + d] = find[d] + foff; 3779cbb8f0dSMatthew G. Knepley numConst += fcdof; 3789cbb8f0dSMatthew G. Knepley foff += fdof; 3799cbb8f0dSMatthew G. Knepley } 3809566063dSJacob Faibussowitsch if (cdof != numConst) PetscCall(PetscSectionSetConstraintDof(section, p, numConst)); 3819cbb8f0dSMatthew G. Knepley } else { 3829cbb8f0dSMatthew G. Knepley for (d = 0; d < cdof; ++d) indices[d] = d; 3839cbb8f0dSMatthew G. Knepley } 3849566063dSJacob Faibussowitsch PetscCall(PetscSectionSetConstraintIndices(section, p, indices)); 3859cbb8f0dSMatthew G. Knepley } 3869cbb8f0dSMatthew G. Knepley } 3879566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 3889cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 3899cbb8f0dSMatthew G. Knepley } 3909cbb8f0dSMatthew G. Knepley 3919cbb8f0dSMatthew G. Knepley /*@C 3929cbb8f0dSMatthew G. Knepley DMPlexCreateSection - Create a PetscSection based upon the dof layout specification provided. 3939cbb8f0dSMatthew G. Knepley 3949cbb8f0dSMatthew G. Knepley Not Collective 3959cbb8f0dSMatthew G. Knepley 3969cbb8f0dSMatthew G. Knepley Input Parameters: 3979cbb8f0dSMatthew G. Knepley + dm - The DMPlex object 39892cfd99aSMartin Diehl . label - The label indicating the mesh support of each field, or NULL for the whole mesh 3999cbb8f0dSMatthew G. Knepley . numComp - An array of size numFields that holds the number of components for each field 4009cbb8f0dSMatthew 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 4019cbb8f0dSMatthew G. Knepley . numBC - The number of boundary conditions 4029cbb8f0dSMatthew G. Knepley . bcField - An array of size numBC giving the field number for each boundry condition 4039cbb8f0dSMatthew G. Knepley . bcComps - [Optional] An array of size numBC giving an IS holding the field components to which each boundary condition applies 4049cbb8f0dSMatthew G. Knepley . bcPoints - An array of size numBC giving an IS holding the Plex points to which each boundary condition applies 4059cbb8f0dSMatthew G. Knepley - perm - Optional permutation of the chart, or NULL 4069cbb8f0dSMatthew G. Knepley 4079cbb8f0dSMatthew G. Knepley Output Parameter: 4089cbb8f0dSMatthew G. Knepley . section - The PetscSection object 4099cbb8f0dSMatthew G. Knepley 4109cbb8f0dSMatthew G. Knepley Notes: 4119cbb8f0dSMatthew 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 4129cbb8f0dSMatthew G. Knepley number of dof for field 0 on each edge. 4139cbb8f0dSMatthew G. Knepley 4149cbb8f0dSMatthew G. Knepley The chart permutation is the same one set using PetscSectionSetPermutation() 4159cbb8f0dSMatthew G. Knepley 4169cbb8f0dSMatthew G. Knepley Level: developer 4179cbb8f0dSMatthew G. Knepley 4181bb6d2a8SBarry Smith TODO: How is this related to DMCreateLocalSection() 4191bb6d2a8SBarry Smith 420db781477SPatrick Sanan .seealso: `DMPlexCreate()`, `PetscSectionCreate()`, `PetscSectionSetPermutation()` 4219cbb8f0dSMatthew G. Knepley @*/ 4229371c9d4SSatish Balay 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) { 4239cbb8f0dSMatthew G. Knepley PetscSection aSec; 4249cbb8f0dSMatthew G. Knepley 4259cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 4269566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSectionFields(dm, numComp, section)); 4279566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSectionDof(dm, label, numDof, *section)); 4289566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSectionBCDof(dm, numBC, bcField, bcComps, bcPoints, *section)); 4299566063dSJacob Faibussowitsch if (perm) PetscCall(PetscSectionSetPermutation(*section, perm)); 4309566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFromOptions(*section)); 4319566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(*section)); 4329566063dSJacob Faibussowitsch PetscCall(DMPlexGetAnchors(dm, &aSec, NULL)); 4339cbb8f0dSMatthew G. Knepley if (numBC || aSec) { 4349566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSectionBCIndicesField(dm, numBC, bcField, bcComps, bcPoints, *section)); 4359566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSectionBCIndices(dm, *section)); 4369cbb8f0dSMatthew G. Knepley } 4379566063dSJacob Faibussowitsch PetscCall(PetscSectionViewFromOptions(*section, NULL, "-section_view")); 4389cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 4399cbb8f0dSMatthew G. Knepley } 4409cbb8f0dSMatthew G. Knepley 4419371c9d4SSatish Balay PetscErrorCode DMCreateLocalSection_Plex(DM dm) { 4429cbb8f0dSMatthew G. Knepley PetscSection section; 443baef929fSMatthew G. Knepley DMLabel *labels; 4449cbb8f0dSMatthew G. Knepley IS *bcPoints, *bcComps; 4459cbb8f0dSMatthew G. Knepley PetscBool *isFE; 4469cbb8f0dSMatthew G. Knepley PetscInt *bcFields, *numComp, *numDof; 4479310035eSMatthew G. Knepley PetscInt depth, dim, numBC = 0, Nf, Nds, s, bc = 0, f; 4489cbb8f0dSMatthew G. Knepley PetscInt cStart, cEnd, cEndInterior; 4499cbb8f0dSMatthew G. Knepley 4509cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 4519566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 4529566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4539566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 4549cbb8f0dSMatthew G. Knepley /* FE and FV boundary conditions are handled slightly differently */ 4559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nf, &isFE)); 456baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 4579cbb8f0dSMatthew G. Knepley PetscObject obj; 4589cbb8f0dSMatthew G. Knepley PetscClassId id; 4599cbb8f0dSMatthew G. Knepley 4609566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, f, NULL, &obj)); 4619566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(obj, &id)); 4629371c9d4SSatish Balay if (id == PETSCFE_CLASSID) { 4639371c9d4SSatish Balay isFE[f] = PETSC_TRUE; 4649371c9d4SSatish Balay } else if (id == PETSCFV_CLASSID) { 4659371c9d4SSatish Balay isFE[f] = PETSC_FALSE; 4669371c9d4SSatish Balay } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f); 4679cbb8f0dSMatthew G. Knepley } 4689cbb8f0dSMatthew G. Knepley /* Allocate boundary point storage for FEM boundaries */ 4699566063dSJacob Faibussowitsch PetscCall(DMGetNumDS(dm, &Nds)); 4709310035eSMatthew G. Knepley for (s = 0; s < Nds; ++s) { 4719310035eSMatthew G. Knepley PetscDS dsBC; 4729310035eSMatthew G. Knepley PetscInt numBd, bd; 4739310035eSMatthew G. Knepley 4749566063dSJacob Faibussowitsch PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC)); 4759566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumBoundary(dsBC, &numBd)); 47608401ef6SPierre 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)"); 4779cbb8f0dSMatthew G. Knepley for (bd = 0; bd < numBd; ++bd) { 4789cbb8f0dSMatthew G. Knepley PetscInt field; 4799cbb8f0dSMatthew G. Knepley DMBoundaryConditionType type; 4809cbb8f0dSMatthew G. Knepley DMLabel label; 4819cbb8f0dSMatthew G. Knepley 4829566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, &type, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL)); 4839cbb8f0dSMatthew G. Knepley if (label && isFE[field] && (type & DM_BC_ESSENTIAL)) ++numBC; 4849cbb8f0dSMatthew G. Knepley } 4859310035eSMatthew G. Knepley } 4869cbb8f0dSMatthew G. Knepley /* Add ghost cell boundaries for FVM */ 4879566063dSJacob Faibussowitsch PetscCall(DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL)); 4889371c9d4SSatish Balay for (f = 0; f < Nf; ++f) 4899371c9d4SSatish Balay if (!isFE[f] && cEndInterior >= 0) ++numBC; 4909566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(numBC, &bcFields, numBC, &bcPoints, numBC, &bcComps)); 4919cbb8f0dSMatthew G. Knepley /* Constrain ghost cells for FV */ 492baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 4939cbb8f0dSMatthew G. Knepley PetscInt *newidx, c; 4949cbb8f0dSMatthew G. Knepley 4959cbb8f0dSMatthew G. Knepley if (isFE[f] || cEndInterior < 0) continue; 4969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cEnd - cEndInterior, &newidx)); 4979cbb8f0dSMatthew G. Knepley for (c = cEndInterior; c < cEnd; ++c) newidx[c - cEndInterior] = c; 4989cbb8f0dSMatthew G. Knepley bcFields[bc] = f; 4999566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, cEnd - cEndInterior, newidx, PETSC_OWN_POINTER, &bcPoints[bc++])); 5009cbb8f0dSMatthew G. Knepley } 501f985c9a8SMatthew G. Knepley /* Complete labels for boundary conditions */ 502799db056SMatthew G. Knepley PetscCall(DMCompleteBCLabels_Internal(dm)); 5039cbb8f0dSMatthew G. Knepley /* Handle FEM Dirichlet boundaries */ 5049566063dSJacob Faibussowitsch PetscCall(DMGetNumDS(dm, &Nds)); 5059310035eSMatthew G. Knepley for (s = 0; s < Nds; ++s) { 5069310035eSMatthew G. Knepley PetscDS dsBC; 5079310035eSMatthew G. Knepley PetscInt numBd, bd; 5089310035eSMatthew G. Knepley 5099566063dSJacob Faibussowitsch PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC)); 5109566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumBoundary(dsBC, &numBd)); 5119cbb8f0dSMatthew G. Knepley for (bd = 0; bd < numBd; ++bd) { 5129cbb8f0dSMatthew G. Knepley DMLabel label; 5139cbb8f0dSMatthew G. Knepley const PetscInt *comps; 5149cbb8f0dSMatthew G. Knepley const PetscInt *values; 5159310035eSMatthew G. Knepley PetscInt bd2, field, numComps, numValues; 5169cbb8f0dSMatthew G. Knepley DMBoundaryConditionType type; 5179310035eSMatthew G. Knepley PetscBool duplicate = PETSC_FALSE; 5189cbb8f0dSMatthew G. Knepley 5199566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, &type, NULL, &label, &numValues, &values, &field, &numComps, &comps, NULL, NULL, NULL)); 5209cbb8f0dSMatthew G. Knepley if (!isFE[field] || !label) continue; 5219310035eSMatthew G. Knepley /* Only want to modify label once */ 5229310035eSMatthew G. Knepley for (bd2 = 0; bd2 < bd; ++bd2) { 52345480ffeSMatthew G. Knepley DMLabel l; 52445480ffeSMatthew G. Knepley 5259566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(dsBC, bd2, NULL, NULL, NULL, &l, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 52645480ffeSMatthew G. Knepley duplicate = l == label ? PETSC_TRUE : PETSC_FALSE; 5279310035eSMatthew G. Knepley if (duplicate) break; 5289310035eSMatthew G. Knepley } 5299cbb8f0dSMatthew G. Knepley /* Filter out cells, if you actually want to constrain cells you need to do things by hand right now */ 5309cbb8f0dSMatthew G. Knepley if (type & DM_BC_ESSENTIAL) { 5319cbb8f0dSMatthew G. Knepley PetscInt *newidx; 5329cbb8f0dSMatthew G. Knepley PetscInt n, newn = 0, p, v; 5339cbb8f0dSMatthew G. Knepley 5349cbb8f0dSMatthew G. Knepley bcFields[bc] = field; 53540cbb1a0SMatthew G. Knepley if (numComps) PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numComps, comps, PETSC_COPY_VALUES, &bcComps[bc])); 5369cbb8f0dSMatthew G. Knepley for (v = 0; v < numValues; ++v) { 5379cbb8f0dSMatthew G. Knepley IS tmp; 5389cbb8f0dSMatthew G. Knepley const PetscInt *idx; 5399cbb8f0dSMatthew G. Knepley 5409566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, values[v], &tmp)); 5419cbb8f0dSMatthew G. Knepley if (!tmp) continue; 5429566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(tmp, &n)); 5439566063dSJacob Faibussowitsch PetscCall(ISGetIndices(tmp, &idx)); 5449cbb8f0dSMatthew G. Knepley if (isFE[field]) { 5459371c9d4SSatish Balay for (p = 0; p < n; ++p) 5469371c9d4SSatish Balay if ((idx[p] < cStart) || (idx[p] >= cEnd)) ++newn; 5479cbb8f0dSMatthew G. Knepley } else { 5489371c9d4SSatish Balay for (p = 0; p < n; ++p) 5499371c9d4SSatish Balay if ((idx[p] >= cStart) || (idx[p] < cEnd)) ++newn; 5509cbb8f0dSMatthew G. Knepley } 5519566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(tmp, &idx)); 5529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tmp)); 5539cbb8f0dSMatthew G. Knepley } 5549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(newn, &newidx)); 5559cbb8f0dSMatthew G. Knepley newn = 0; 5569cbb8f0dSMatthew G. Knepley for (v = 0; v < numValues; ++v) { 5579cbb8f0dSMatthew G. Knepley IS tmp; 5589cbb8f0dSMatthew G. Knepley const PetscInt *idx; 5599cbb8f0dSMatthew G. Knepley 5609566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, values[v], &tmp)); 5619cbb8f0dSMatthew G. Knepley if (!tmp) continue; 5629566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(tmp, &n)); 5639566063dSJacob Faibussowitsch PetscCall(ISGetIndices(tmp, &idx)); 5649cbb8f0dSMatthew G. Knepley if (isFE[field]) { 5659371c9d4SSatish Balay for (p = 0; p < n; ++p) 5669371c9d4SSatish Balay if ((idx[p] < cStart) || (idx[p] >= cEnd)) newidx[newn++] = idx[p]; 5679cbb8f0dSMatthew G. Knepley } else { 5689371c9d4SSatish Balay for (p = 0; p < n; ++p) 5699371c9d4SSatish Balay if ((idx[p] >= cStart) || (idx[p] < cEnd)) newidx[newn++] = idx[p]; 5709cbb8f0dSMatthew G. Knepley } 5719566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(tmp, &idx)); 5729566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tmp)); 5739cbb8f0dSMatthew G. Knepley } 5749566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, newn, newidx, PETSC_OWN_POINTER, &bcPoints[bc++])); 5759310035eSMatthew G. Knepley } 5769cbb8f0dSMatthew G. Knepley } 5779cbb8f0dSMatthew G. Knepley } 5789cbb8f0dSMatthew G. Knepley /* Handle discretization */ 5799566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(Nf, &labels, Nf, &numComp, Nf * (dim + 1), &numDof)); 580baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 58144a7f3ddSMatthew G. Knepley labels[f] = dm->fields[f].label; 5829cbb8f0dSMatthew G. Knepley if (isFE[f]) { 58344a7f3ddSMatthew G. Knepley PetscFE fe = (PetscFE)dm->fields[f].disc; 5849cbb8f0dSMatthew G. Knepley const PetscInt *numFieldDof; 585a3254110SMatthew Knepley PetscInt fedim, d; 5869cbb8f0dSMatthew G. Knepley 5879566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumComponents(fe, &numComp[f])); 5889566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumDof(fe, &numFieldDof)); 5899566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &fedim)); 590a3254110SMatthew Knepley for (d = 0; d < PetscMin(dim, fedim) + 1; ++d) numDof[f * (dim + 1) + d] = numFieldDof[d]; 5919cbb8f0dSMatthew G. Knepley } else { 59244a7f3ddSMatthew G. Knepley PetscFV fv = (PetscFV)dm->fields[f].disc; 5939cbb8f0dSMatthew G. Knepley 5949566063dSJacob Faibussowitsch PetscCall(PetscFVGetNumComponents(fv, &numComp[f])); 5959cbb8f0dSMatthew G. Knepley numDof[f * (dim + 1) + dim] = numComp[f]; 5969cbb8f0dSMatthew G. Knepley } 5979cbb8f0dSMatthew G. Knepley } 5989566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 599baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 6009cbb8f0dSMatthew G. Knepley PetscInt d; 6019371c9d4SSatish Balay for (d = 1; d < dim; ++d) { PetscCheck(numDof[f * (dim + 1) + d] <= 0 || depth >= dim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated when unknowns are specified on edges or faces."); } 6029cbb8f0dSMatthew G. Knepley } 6039566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSection(dm, labels, numComp, numDof, numBC, bcFields, bcComps, bcPoints, NULL, §ion)); 604baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 6059cbb8f0dSMatthew G. Knepley PetscFE fe; 6069cbb8f0dSMatthew G. Knepley const char *name; 6079cbb8f0dSMatthew G. Knepley 6089566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, f, NULL, (PetscObject *)&fe)); 609cdfe53dfSMatthew G. Knepley if (!((PetscObject)fe)->name) continue; 6109566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)fe, &name)); 6119566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(section, f, name)); 6129cbb8f0dSMatthew G. Knepley } 6139566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(dm, section)); 6149566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ion)); 6159371c9d4SSatish Balay for (bc = 0; bc < numBC; ++bc) { 6169371c9d4SSatish Balay PetscCall(ISDestroy(&bcPoints[bc])); 6179371c9d4SSatish Balay PetscCall(ISDestroy(&bcComps[bc])); 6189371c9d4SSatish Balay } 6199566063dSJacob Faibussowitsch PetscCall(PetscFree3(bcFields, bcPoints, bcComps)); 6209566063dSJacob Faibussowitsch PetscCall(PetscFree3(labels, numComp, numDof)); 6219566063dSJacob Faibussowitsch PetscCall(PetscFree(isFE)); 6229cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 6239cbb8f0dSMatthew G. Knepley } 624