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; 45baef929fSMatthew G. Knepley PetscInt h; 469cbb8f0dSMatthew G. Knepley PetscSectionSym sym; 479cbb8f0dSMatthew G. Knepley 489cbb8f0dSMatthew G. Knepley ierr = PetscDualSpaceGetDM(dspace,&K);CHKERRQ(ierr); 499cbb8f0dSMatthew G. Knepley ierr = PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)*section),depthLabel,&sym);CHKERRQ(ierr); 509cbb8f0dSMatthew G. Knepley for (h = 0; h <= depth; h++) { 519cbb8f0dSMatthew G. Knepley PetscDualSpace hspace; 529cbb8f0dSMatthew G. Knepley PetscInt kStart, kEnd; 539cbb8f0dSMatthew G. Knepley PetscInt kConeSize; 549cbb8f0dSMatthew G. Knepley const PetscInt **perms0 = NULL; 559cbb8f0dSMatthew G. Knepley const PetscScalar **flips0 = NULL; 569cbb8f0dSMatthew G. Knepley 579cbb8f0dSMatthew G. Knepley ierr = PetscDualSpaceGetHeightSubspace(dspace,h,&hspace);CHKERRQ(ierr); 589cbb8f0dSMatthew G. Knepley ierr = DMPlexGetHeightStratum(K,h,&kStart,&kEnd);CHKERRQ(ierr); 599cbb8f0dSMatthew G. Knepley if (!hspace) continue; 609cbb8f0dSMatthew G. Knepley ierr = PetscDualSpaceGetSymmetries(hspace,&perms,&flips);CHKERRQ(ierr); 619cbb8f0dSMatthew G. Knepley if (perms) perms0 = perms[0]; 629cbb8f0dSMatthew G. Knepley if (flips) flips0 = flips[0]; 639cbb8f0dSMatthew G. Knepley if (!(perms0 || flips0)) continue; 649cbb8f0dSMatthew G. Knepley ierr = DMPlexGetConeSize(K,kStart,&kConeSize);CHKERRQ(ierr); 659cbb8f0dSMatthew 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); 669cbb8f0dSMatthew G. Knepley } 679cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFieldSym(*section,f,sym);CHKERRQ(ierr); 689cbb8f0dSMatthew G. Knepley ierr = PetscSectionSymDestroy(&sym);CHKERRQ(ierr); 699cbb8f0dSMatthew G. Knepley } 709cbb8f0dSMatthew G. Knepley } 719cbb8f0dSMatthew G. Knepley } 729cbb8f0dSMatthew G. Knepley } 739cbb8f0dSMatthew G. Knepley } 749cbb8f0dSMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 759cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetChart(*section, pStart, pEnd);CHKERRQ(ierr); 76baef929fSMatthew G. Knepley ierr = PetscFree(isFE);CHKERRQ(ierr); 77baef929fSMatthew G. Knepley PetscFunctionReturn(0); 78baef929fSMatthew G. Knepley } 79baef929fSMatthew G. Knepley 80baef929fSMatthew G. Knepley /* Set the number of dof on each point and separate by fields */ 81baef929fSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionDof(DM dm, DMLabel label[],const PetscInt numDof[], PetscSection section) 82baef929fSMatthew G. Knepley { 83baef929fSMatthew G. Knepley DMLabel depthLabel; 84*412e9a14SMatthew G. Knepley DMPolytopeType ct; 85baef929fSMatthew G. Knepley PetscInt depth, cellHeight, pStart = 0, pEnd = 0; 86baef929fSMatthew G. Knepley PetscInt Nf, f, dim, d, dep, p; 87baef929fSMatthew G. Knepley PetscBool *isFE; 88baef929fSMatthew G. Knepley PetscErrorCode ierr; 89baef929fSMatthew G. Knepley 90baef929fSMatthew G. Knepley PetscFunctionBegin; 91baef929fSMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 929cbb8f0dSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 93baef929fSMatthew G. Knepley ierr = DMPlexGetDepthLabel(dm,&depthLabel);CHKERRQ(ierr); 94baef929fSMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 95baef929fSMatthew G. Knepley ierr = PetscMalloc1(Nf, &isFE);CHKERRQ(ierr); 96baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 97baef929fSMatthew G. Knepley PetscObject obj; 98baef929fSMatthew G. Knepley PetscClassId id; 99baef929fSMatthew G. Knepley 10044a7f3ddSMatthew G. Knepley ierr = DMGetField(dm, f, NULL, &obj);CHKERRQ(ierr); 101baef929fSMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1021274e1a1SLawrence Mitchell /* User is allowed to put a "placeholder" field in (c.f. DMCreateDS) */ 1031274e1a1SLawrence Mitchell isFE[f] = id == PETSCFE_CLASSID ? PETSC_TRUE : PETSC_FALSE; 104baef929fSMatthew G. Knepley } 105baef929fSMatthew G. Knepley 1069cbb8f0dSMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 107baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 108baef929fSMatthew G. Knepley if (label && label[f]) { 109baef929fSMatthew G. Knepley IS pointIS; 110baef929fSMatthew G. Knepley const PetscInt *points; 111baef929fSMatthew G. Knepley PetscInt n; 112baef929fSMatthew G. Knepley 113baef929fSMatthew G. Knepley ierr = DMLabelGetStratumIS(label[f], 1, &pointIS);CHKERRQ(ierr); 114baef929fSMatthew G. Knepley if (!pointIS) continue; 115baef929fSMatthew G. Knepley ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr); 116baef929fSMatthew G. Knepley ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 117baef929fSMatthew G. Knepley for (p = 0; p < n; ++p) { 118baef929fSMatthew G. Knepley const PetscInt point = points[p]; 119baef929fSMatthew G. Knepley PetscInt dof, d; 120baef929fSMatthew G. Knepley 121*412e9a14SMatthew G. Knepley ierr = DMPlexGetCellType(dm, point, &ct);CHKERRQ(ierr); 122baef929fSMatthew G. Knepley ierr = DMLabelGetValue(depthLabel, point, &d);CHKERRQ(ierr); 123*412e9a14SMatthew G. Knepley /* If this is a tensor prism point, use dof for one dimension lower */ 124*412e9a14SMatthew G. Knepley switch (ct) { 125*412e9a14SMatthew G. Knepley case DM_POLYTOPE_POINT_PRISM_TENSOR: 126*412e9a14SMatthew G. Knepley case DM_POLYTOPE_SEG_PRISM_TENSOR: 127*412e9a14SMatthew G. Knepley case DM_POLYTOPE_TRI_PRISM_TENSOR: 128*412e9a14SMatthew G. Knepley case DM_POLYTOPE_QUAD_PRISM_TENSOR: 129*412e9a14SMatthew G. Knepley --d;break; 130*412e9a14SMatthew G. Knepley default: break; 131*412e9a14SMatthew G. Knepley } 132baef929fSMatthew G. Knepley dof = d < 0 ? 0 : numDof[f*(dim+1)+d]; 133baef929fSMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, point, f, dof);CHKERRQ(ierr); 134baef929fSMatthew G. Knepley ierr = PetscSectionAddDof(section, point, dof);CHKERRQ(ierr); 135baef929fSMatthew G. Knepley } 136baef929fSMatthew G. Knepley ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 137baef929fSMatthew G. Knepley ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 138baef929fSMatthew G. Knepley } else { 1399cbb8f0dSMatthew G. Knepley for (dep = 0; dep <= depth - cellHeight; ++dep) { 1403fd2e703SMatthew G. Knepley /* Cases: dim > depth (cell-vertex mesh), dim == depth (fully interpolated), dim < depth (interpolated submesh) */ 1413fd2e703SMatthew G. Knepley d = dim <= depth ? dep : (!dep ? 0 : dim); 1429cbb8f0dSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, dep, &pStart, &pEnd);CHKERRQ(ierr); 1439cbb8f0dSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 144baef929fSMatthew G. Knepley const PetscInt dof = numDof[f*(dim+1)+d]; 1459cbb8f0dSMatthew G. Knepley 146*412e9a14SMatthew G. Knepley ierr = DMPlexGetCellType(dm, p, &ct);CHKERRQ(ierr); 147*412e9a14SMatthew G. Knepley switch (ct) { 148*412e9a14SMatthew G. Knepley case DM_POLYTOPE_POINT_PRISM_TENSOR: 149*412e9a14SMatthew G. Knepley case DM_POLYTOPE_SEG_PRISM_TENSOR: 150*412e9a14SMatthew G. Knepley case DM_POLYTOPE_TRI_PRISM_TENSOR: 151*412e9a14SMatthew G. Knepley case DM_POLYTOPE_QUAD_PRISM_TENSOR: 152*412e9a14SMatthew G. Knepley if (isFE[f]) continue; 153*412e9a14SMatthew G. Knepley default: break; 154*412e9a14SMatthew G. Knepley } 155baef929fSMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, p, f, dof);CHKERRQ(ierr); 156baef929fSMatthew G. Knepley ierr = PetscSectionAddDof(section, p, dof);CHKERRQ(ierr); 1579cbb8f0dSMatthew G. Knepley } 158baef929fSMatthew G. Knepley } 1599cbb8f0dSMatthew G. Knepley } 1609cbb8f0dSMatthew G. Knepley } 1619cbb8f0dSMatthew G. Knepley ierr = PetscFree(isFE);CHKERRQ(ierr); 1629cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 1639cbb8f0dSMatthew G. Knepley } 1649cbb8f0dSMatthew G. Knepley 1659cbb8f0dSMatthew G. Knepley /* Set the number of dof on each point and separate by fields 1669cbb8f0dSMatthew G. Knepley If bcComps is NULL or the IS is NULL, constrain every dof on the point 1679cbb8f0dSMatthew G. Knepley */ 1689cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCDof(DM dm, PetscInt numBC, const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section) 1699cbb8f0dSMatthew G. Knepley { 170baef929fSMatthew G. Knepley PetscInt Nf; 1719cbb8f0dSMatthew G. Knepley PetscInt bc; 1729cbb8f0dSMatthew G. Knepley PetscSection aSec; 1739cbb8f0dSMatthew G. Knepley PetscErrorCode ierr; 1749cbb8f0dSMatthew G. Knepley 1759cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 176baef929fSMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 1779cbb8f0dSMatthew G. Knepley for (bc = 0; bc < numBC; ++bc) { 1789cbb8f0dSMatthew G. Knepley PetscInt field = 0; 1799cbb8f0dSMatthew G. Knepley const PetscInt *comp; 1809cbb8f0dSMatthew G. Knepley const PetscInt *idx; 181ce093827SMatthew G. Knepley PetscInt Nc = 0, cNc = -1, n, i; 1829cbb8f0dSMatthew G. Knepley 183ce093827SMatthew G. Knepley if (Nf) { 184ce093827SMatthew G. Knepley field = bcField[bc]; 185ce093827SMatthew G. Knepley ierr = PetscSectionGetFieldComponents(section, field, &Nc);CHKERRQ(ierr); 186ce093827SMatthew G. Knepley } 187ce093827SMatthew G. Knepley if (bcComps && bcComps[bc]) {ierr = ISGetLocalSize(bcComps[bc], &cNc);CHKERRQ(ierr);} 1889cbb8f0dSMatthew G. Knepley if (bcComps && bcComps[bc]) {ierr = ISGetIndices(bcComps[bc], &comp);CHKERRQ(ierr);} 1899cbb8f0dSMatthew G. Knepley ierr = ISGetLocalSize(bcPoints[bc], &n);CHKERRQ(ierr); 1909cbb8f0dSMatthew G. Knepley ierr = ISGetIndices(bcPoints[bc], &idx);CHKERRQ(ierr); 1919cbb8f0dSMatthew G. Knepley for (i = 0; i < n; ++i) { 1929cbb8f0dSMatthew G. Knepley const PetscInt p = idx[i]; 1939cbb8f0dSMatthew G. Knepley PetscInt numConst; 1949cbb8f0dSMatthew G. Knepley 195baef929fSMatthew G. Knepley if (Nf) { 1969cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldDof(section, p, field, &numConst);CHKERRQ(ierr); 1979cbb8f0dSMatthew G. Knepley } else { 1989cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetDof(section, p, &numConst);CHKERRQ(ierr); 1999cbb8f0dSMatthew G. Knepley } 200ce093827SMatthew G. Knepley /* If Nc <= 0, constrain every dof on the point */ 201ce093827SMatthew G. Knepley if (cNc > 0) { 202ce093827SMatthew G. Knepley /* We assume that a point may have multiple "nodes", which are collections of Nc dofs, 203ce093827SMatthew G. Knepley and that those dofs are numbered n*Nc+c */ 204ce093827SMatthew G. Knepley if (Nf) { 205ce093827SMatthew 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); 206ce093827SMatthew G. Knepley numConst = (numConst/Nc) * cNc; 207ce093827SMatthew G. Knepley } else { 208ce093827SMatthew G. Knepley numConst = PetscMin(numConst, cNc); 209ce093827SMatthew G. Knepley } 210ce093827SMatthew G. Knepley } 211baef929fSMatthew G. Knepley if (Nf) {ierr = PetscSectionAddFieldConstraintDof(section, p, field, numConst);CHKERRQ(ierr);} 2129cbb8f0dSMatthew G. Knepley ierr = PetscSectionAddConstraintDof(section, p, numConst);CHKERRQ(ierr); 2139cbb8f0dSMatthew G. Knepley } 2149cbb8f0dSMatthew G. Knepley ierr = ISRestoreIndices(bcPoints[bc], &idx);CHKERRQ(ierr); 2159cbb8f0dSMatthew G. Knepley if (bcComps && bcComps[bc]) {ierr = ISRestoreIndices(bcComps[bc], &comp);CHKERRQ(ierr);} 2169cbb8f0dSMatthew G. Knepley } 2179cbb8f0dSMatthew G. Knepley ierr = DMPlexGetAnchors(dm, &aSec, NULL);CHKERRQ(ierr); 2189cbb8f0dSMatthew G. Knepley if (aSec) { 2199cbb8f0dSMatthew G. Knepley PetscInt aStart, aEnd, a; 2209cbb8f0dSMatthew G. Knepley 2219cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetChart(aSec, &aStart, &aEnd);CHKERRQ(ierr); 2229cbb8f0dSMatthew G. Knepley for (a = aStart; a < aEnd; a++) { 2239cbb8f0dSMatthew G. Knepley PetscInt dof, f; 2249cbb8f0dSMatthew G. Knepley 2259cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetDof(aSec, a, &dof);CHKERRQ(ierr); 2269cbb8f0dSMatthew G. Knepley if (dof) { 2279cbb8f0dSMatthew G. Knepley /* if there are point-to-point constraints, then all dofs are constrained */ 2289cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetDof(section, a, &dof);CHKERRQ(ierr); 2299cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetConstraintDof(section, a, dof);CHKERRQ(ierr); 230baef929fSMatthew G. Knepley for (f = 0; f < Nf; f++) { 2319cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldDof(section, a, f, &dof);CHKERRQ(ierr); 2329cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFieldConstraintDof(section, a, f, dof);CHKERRQ(ierr); 2339cbb8f0dSMatthew G. Knepley } 2349cbb8f0dSMatthew G. Knepley } 2359cbb8f0dSMatthew G. Knepley } 2369cbb8f0dSMatthew G. Knepley } 2379cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 2389cbb8f0dSMatthew G. Knepley } 2399cbb8f0dSMatthew G. Knepley 2409cbb8f0dSMatthew G. Knepley /* Set the constrained field indices on each point 2419cbb8f0dSMatthew G. Knepley If bcComps is NULL or the IS is NULL, constrain every dof on the point 2429cbb8f0dSMatthew G. Knepley */ 2439cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCIndicesField(DM dm, PetscInt numBC,const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section) 2449cbb8f0dSMatthew G. Knepley { 2459cbb8f0dSMatthew G. Knepley PetscSection aSec; 2469cbb8f0dSMatthew G. Knepley PetscInt *indices; 247baef929fSMatthew G. Knepley PetscInt Nf, cdof, maxDof = 0, pStart, pEnd, p, bc, f, d; 2489cbb8f0dSMatthew G. Knepley PetscErrorCode ierr; 2499cbb8f0dSMatthew G. Knepley 2509cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 251baef929fSMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 252baef929fSMatthew G. Knepley if (!Nf) PetscFunctionReturn(0); 2539cbb8f0dSMatthew G. Knepley /* Initialize all field indices to -1 */ 2549cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 2559cbb8f0dSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) {ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); maxDof = PetscMax(maxDof, cdof);} 2569cbb8f0dSMatthew G. Knepley ierr = PetscMalloc1(maxDof, &indices);CHKERRQ(ierr); 2579cbb8f0dSMatthew G. Knepley for (d = 0; d < maxDof; ++d) indices[d] = -1; 258baef929fSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) for (f = 0; f < Nf; ++f) {ierr = PetscSectionSetFieldConstraintIndices(section, p, f, indices);CHKERRQ(ierr);} 2599cbb8f0dSMatthew G. Knepley /* Handle BC constraints */ 2609cbb8f0dSMatthew G. Knepley for (bc = 0; bc < numBC; ++bc) { 2619cbb8f0dSMatthew G. Knepley const PetscInt field = bcField[bc]; 2629cbb8f0dSMatthew G. Knepley const PetscInt *comp, *idx; 263ce093827SMatthew G. Knepley PetscInt Nc, cNc = -1, n, i; 2649cbb8f0dSMatthew G. Knepley 265ce093827SMatthew G. Knepley ierr = PetscSectionGetFieldComponents(section, field, &Nc);CHKERRQ(ierr); 266ce093827SMatthew G. Knepley if (bcComps && bcComps[bc]) {ierr = ISGetLocalSize(bcComps[bc], &cNc);CHKERRQ(ierr);} 2679cbb8f0dSMatthew G. Knepley if (bcComps && bcComps[bc]) {ierr = ISGetIndices(bcComps[bc], &comp);CHKERRQ(ierr);} 2689cbb8f0dSMatthew G. Knepley ierr = ISGetLocalSize(bcPoints[bc], &n);CHKERRQ(ierr); 2699cbb8f0dSMatthew G. Knepley ierr = ISGetIndices(bcPoints[bc], &idx);CHKERRQ(ierr); 2709cbb8f0dSMatthew G. Knepley for (i = 0; i < n; ++i) { 2719cbb8f0dSMatthew G. Knepley const PetscInt p = idx[i]; 2729cbb8f0dSMatthew G. Knepley const PetscInt *find; 273ce093827SMatthew G. Knepley PetscInt fdof, fcdof, c, j; 2749cbb8f0dSMatthew G. Knepley 2759cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldDof(section, p, field, &fdof);CHKERRQ(ierr); 2769cbb8f0dSMatthew G. Knepley if (!fdof) continue; 277ce093827SMatthew G. Knepley if (cNc < 0) { 2789cbb8f0dSMatthew G. Knepley for (d = 0; d < fdof; ++d) indices[d] = d; 2799cbb8f0dSMatthew G. Knepley fcdof = fdof; 2809cbb8f0dSMatthew G. Knepley } else { 281ce093827SMatthew G. Knepley /* We assume that a point may have multiple "nodes", which are collections of Nc dofs, 282ce093827SMatthew G. Knepley and that those dofs are numbered n*Nc+c */ 2839cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintDof(section, p, field, &fcdof);CHKERRQ(ierr); 2849cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintIndices(section, p, field, &find);CHKERRQ(ierr); 285ce093827SMatthew G. Knepley /* Get indices constrained by previous bcs */ 2869cbb8f0dSMatthew G. Knepley for (d = 0; d < fcdof; ++d) {if (find[d] < 0) break; indices[d] = find[d];} 287ce093827SMatthew G. Knepley for (j = 0; j < fdof/Nc; ++j) for (c = 0; c < cNc; ++c) indices[d++] = j*Nc + comp[c]; 2889cbb8f0dSMatthew G. Knepley ierr = PetscSortRemoveDupsInt(&d, indices);CHKERRQ(ierr); 2899cbb8f0dSMatthew G. Knepley for (c = d; c < fcdof; ++c) indices[c] = -1; 2909cbb8f0dSMatthew G. Knepley fcdof = d; 2919cbb8f0dSMatthew G. Knepley } 2929cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFieldConstraintDof(section, p, field, fcdof);CHKERRQ(ierr); 2939cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFieldConstraintIndices(section, p, field, indices);CHKERRQ(ierr); 2949cbb8f0dSMatthew G. Knepley } 2959cbb8f0dSMatthew G. Knepley if (bcComps && bcComps[bc]) {ierr = ISRestoreIndices(bcComps[bc], &comp);CHKERRQ(ierr);} 2969cbb8f0dSMatthew G. Knepley ierr = ISRestoreIndices(bcPoints[bc], &idx);CHKERRQ(ierr); 2979cbb8f0dSMatthew G. Knepley } 2989cbb8f0dSMatthew G. Knepley /* Handle anchors */ 2999cbb8f0dSMatthew G. Knepley ierr = DMPlexGetAnchors(dm, &aSec, NULL);CHKERRQ(ierr); 3009cbb8f0dSMatthew G. Knepley if (aSec) { 3019cbb8f0dSMatthew G. Knepley PetscInt aStart, aEnd, a; 3029cbb8f0dSMatthew G. Knepley 3039cbb8f0dSMatthew G. Knepley for (d = 0; d < maxDof; ++d) indices[d] = d; 3049cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetChart(aSec, &aStart, &aEnd);CHKERRQ(ierr); 3059cbb8f0dSMatthew G. Knepley for (a = aStart; a < aEnd; a++) { 3069cbb8f0dSMatthew G. Knepley PetscInt dof, f; 3079cbb8f0dSMatthew G. Knepley 3089cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetDof(aSec, a, &dof);CHKERRQ(ierr); 3099cbb8f0dSMatthew G. Knepley if (dof) { 3109cbb8f0dSMatthew G. Knepley /* if there are point-to-point constraints, then all dofs are constrained */ 311baef929fSMatthew G. Knepley for (f = 0; f < Nf; f++) { 3129cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFieldConstraintIndices(section, a, f, indices);CHKERRQ(ierr); 3139cbb8f0dSMatthew G. Knepley } 3149cbb8f0dSMatthew G. Knepley } 3159cbb8f0dSMatthew G. Knepley } 3169cbb8f0dSMatthew G. Knepley } 3179cbb8f0dSMatthew G. Knepley ierr = PetscFree(indices);CHKERRQ(ierr); 3189cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 3199cbb8f0dSMatthew G. Knepley } 3209cbb8f0dSMatthew G. Knepley 3219cbb8f0dSMatthew G. Knepley /* Set the constrained indices on each point */ 3229cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCIndices(DM dm, PetscSection section) 3239cbb8f0dSMatthew G. Knepley { 3249cbb8f0dSMatthew G. Knepley PetscInt *indices; 325baef929fSMatthew G. Knepley PetscInt Nf, maxDof, pStart, pEnd, p, f, d; 3269cbb8f0dSMatthew G. Knepley PetscErrorCode ierr; 3279cbb8f0dSMatthew G. Knepley 3289cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 329baef929fSMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 3309cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetMaxDof(section, &maxDof);CHKERRQ(ierr); 3319cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 3329cbb8f0dSMatthew G. Knepley ierr = PetscMalloc1(maxDof, &indices);CHKERRQ(ierr); 3339cbb8f0dSMatthew G. Knepley for (d = 0; d < maxDof; ++d) indices[d] = -1; 3349cbb8f0dSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 3359cbb8f0dSMatthew G. Knepley PetscInt cdof, d; 3369cbb8f0dSMatthew G. Knepley 3379cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 3389cbb8f0dSMatthew G. Knepley if (cdof) { 339baef929fSMatthew G. Knepley if (Nf) { 3409cbb8f0dSMatthew G. Knepley PetscInt numConst = 0, foff = 0; 3419cbb8f0dSMatthew G. Knepley 342baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 3439cbb8f0dSMatthew G. Knepley const PetscInt *find; 3449cbb8f0dSMatthew G. Knepley PetscInt fcdof, fdof; 3459cbb8f0dSMatthew G. Knepley 3469cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr); 3479cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr); 3489cbb8f0dSMatthew G. Knepley /* Change constraint numbering from field component to local dof number */ 3499cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintIndices(section, p, f, &find);CHKERRQ(ierr); 3509cbb8f0dSMatthew G. Knepley for (d = 0; d < fcdof; ++d) indices[numConst+d] = find[d] + foff; 3519cbb8f0dSMatthew G. Knepley numConst += fcdof; 3529cbb8f0dSMatthew G. Knepley foff += fdof; 3539cbb8f0dSMatthew G. Knepley } 3549cbb8f0dSMatthew G. Knepley if (cdof != numConst) {ierr = PetscSectionSetConstraintDof(section, p, numConst);CHKERRQ(ierr);} 3559cbb8f0dSMatthew G. Knepley } else { 3569cbb8f0dSMatthew G. Knepley for (d = 0; d < cdof; ++d) indices[d] = d; 3579cbb8f0dSMatthew G. Knepley } 3589cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetConstraintIndices(section, p, indices);CHKERRQ(ierr); 3599cbb8f0dSMatthew G. Knepley } 3609cbb8f0dSMatthew G. Knepley } 3619cbb8f0dSMatthew G. Knepley ierr = PetscFree(indices);CHKERRQ(ierr); 3629cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 3639cbb8f0dSMatthew G. Knepley } 3649cbb8f0dSMatthew G. Knepley 3659cbb8f0dSMatthew G. Knepley /*@C 3669cbb8f0dSMatthew G. Knepley DMPlexCreateSection - Create a PetscSection based upon the dof layout specification provided. 3679cbb8f0dSMatthew G. Knepley 3689cbb8f0dSMatthew G. Knepley Not Collective 3699cbb8f0dSMatthew G. Knepley 3709cbb8f0dSMatthew G. Knepley Input Parameters: 3719cbb8f0dSMatthew G. Knepley + dm - The DMPlex object 37292cfd99aSMartin Diehl . label - The label indicating the mesh support of each field, or NULL for the whole mesh 3739cbb8f0dSMatthew G. Knepley . numComp - An array of size numFields that holds the number of components for each field 3749cbb8f0dSMatthew 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 3759cbb8f0dSMatthew G. Knepley . numBC - The number of boundary conditions 3769cbb8f0dSMatthew G. Knepley . bcField - An array of size numBC giving the field number for each boundry condition 3779cbb8f0dSMatthew G. Knepley . bcComps - [Optional] An array of size numBC giving an IS holding the field components to which each boundary condition applies 3789cbb8f0dSMatthew G. Knepley . bcPoints - An array of size numBC giving an IS holding the Plex points to which each boundary condition applies 3799cbb8f0dSMatthew G. Knepley - perm - Optional permutation of the chart, or NULL 3809cbb8f0dSMatthew G. Knepley 3819cbb8f0dSMatthew G. Knepley Output Parameter: 3829cbb8f0dSMatthew G. Knepley . section - The PetscSection object 3839cbb8f0dSMatthew G. Knepley 3849cbb8f0dSMatthew G. Knepley Notes: 3859cbb8f0dSMatthew 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 3869cbb8f0dSMatthew G. Knepley number of dof for field 0 on each edge. 3879cbb8f0dSMatthew G. Knepley 3889cbb8f0dSMatthew G. Knepley The chart permutation is the same one set using PetscSectionSetPermutation() 3899cbb8f0dSMatthew G. Knepley 3909cbb8f0dSMatthew G. Knepley Level: developer 3919cbb8f0dSMatthew G. Knepley 3921bb6d2a8SBarry Smith TODO: How is this related to DMCreateLocalSection() 3931bb6d2a8SBarry Smith 3949cbb8f0dSMatthew G. Knepley .seealso: DMPlexCreate(), PetscSectionCreate(), PetscSectionSetPermutation() 3959cbb8f0dSMatthew G. Knepley @*/ 396baef929fSMatthew 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) 3979cbb8f0dSMatthew G. Knepley { 3989cbb8f0dSMatthew G. Knepley PetscSection aSec; 3999cbb8f0dSMatthew G. Knepley PetscErrorCode ierr; 4009cbb8f0dSMatthew G. Knepley 4019cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 402baef929fSMatthew G. Knepley ierr = DMPlexCreateSectionFields(dm, numComp, section);CHKERRQ(ierr); 403baef929fSMatthew G. Knepley ierr = DMPlexCreateSectionDof(dm, label, numDof, *section);CHKERRQ(ierr); 4049cbb8f0dSMatthew G. Knepley ierr = DMPlexCreateSectionBCDof(dm, numBC, bcField, bcComps, bcPoints, *section);CHKERRQ(ierr); 4059cbb8f0dSMatthew G. Knepley if (perm) {ierr = PetscSectionSetPermutation(*section, perm);CHKERRQ(ierr);} 4069cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFromOptions(*section);CHKERRQ(ierr); 4079cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 4089cbb8f0dSMatthew G. Knepley ierr = DMPlexGetAnchors(dm,&aSec,NULL);CHKERRQ(ierr); 4099cbb8f0dSMatthew G. Knepley if (numBC || aSec) { 4109cbb8f0dSMatthew G. Knepley ierr = DMPlexCreateSectionBCIndicesField(dm, numBC, bcField, bcComps, bcPoints, *section);CHKERRQ(ierr); 4119cbb8f0dSMatthew G. Knepley ierr = DMPlexCreateSectionBCIndices(dm, *section);CHKERRQ(ierr); 4129cbb8f0dSMatthew G. Knepley } 4139cbb8f0dSMatthew G. Knepley ierr = PetscSectionViewFromOptions(*section,NULL,"-section_view");CHKERRQ(ierr); 4149cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 4159cbb8f0dSMatthew G. Knepley } 4169cbb8f0dSMatthew G. Knepley 4171bb6d2a8SBarry Smith PetscErrorCode DMCreateLocalSection_Plex(DM dm) 4189cbb8f0dSMatthew G. Knepley { 419e5e52638SMatthew G. Knepley PetscDS probBC; 4209cbb8f0dSMatthew G. Knepley PetscSection section; 421baef929fSMatthew G. Knepley DMLabel *labels; 4229cbb8f0dSMatthew G. Knepley IS *bcPoints, *bcComps; 4239cbb8f0dSMatthew G. Knepley PetscBool *isFE; 4249cbb8f0dSMatthew G. Knepley PetscInt *bcFields, *numComp, *numDof; 425baef929fSMatthew G. Knepley PetscInt depth, dim, numBd, numBC = 0, Nf, bd, bc = 0, f; 4269cbb8f0dSMatthew G. Knepley PetscInt cStart, cEnd, cEndInterior; 4279cbb8f0dSMatthew G. Knepley PetscErrorCode ierr; 4289cbb8f0dSMatthew G. Knepley 4299cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 430baef929fSMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 4319cbb8f0dSMatthew G. Knepley /* FE and FV boundary conditions are handled slightly differently */ 432baef929fSMatthew G. Knepley ierr = PetscMalloc1(Nf, &isFE);CHKERRQ(ierr); 433baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 4349cbb8f0dSMatthew G. Knepley PetscObject obj; 4359cbb8f0dSMatthew G. Knepley PetscClassId id; 4369cbb8f0dSMatthew G. Knepley 43744a7f3ddSMatthew G. Knepley ierr = DMGetField(dm, f, NULL, &obj);CHKERRQ(ierr); 4389cbb8f0dSMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 4399cbb8f0dSMatthew G. Knepley if (id == PETSCFE_CLASSID) {isFE[f] = PETSC_TRUE;} 4409cbb8f0dSMatthew G. Knepley else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;} 4419cbb8f0dSMatthew G. Knepley else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f); 4429cbb8f0dSMatthew G. Knepley } 4439cbb8f0dSMatthew G. Knepley /* Allocate boundary point storage for FEM boundaries */ 4449cbb8f0dSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 4459cbb8f0dSMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 4469cbb8f0dSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 447485ad865SMatthew G. Knepley ierr = DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);CHKERRQ(ierr); 448e5e52638SMatthew G. Knepley ierr = DMGetDS(dm, &probBC);CHKERRQ(ierr); 449e5e52638SMatthew G. Knepley ierr = PetscDSGetNumBoundary(probBC, &numBd);CHKERRQ(ierr); 450baef929fSMatthew 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)"); 4519cbb8f0dSMatthew G. Knepley for (bd = 0; bd < numBd; ++bd) { 4529cbb8f0dSMatthew G. Knepley PetscInt field; 4539cbb8f0dSMatthew G. Knepley DMBoundaryConditionType type; 4549cbb8f0dSMatthew G. Knepley const char *labelName; 4559cbb8f0dSMatthew G. Knepley DMLabel label; 4569cbb8f0dSMatthew G. Knepley 457e5e52638SMatthew G. Knepley ierr = PetscDSGetBoundary(probBC, bd, &type, NULL, &labelName, &field, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 4589cbb8f0dSMatthew G. Knepley ierr = DMGetLabel(dm,labelName,&label);CHKERRQ(ierr); 4599cbb8f0dSMatthew G. Knepley if (label && isFE[field] && (type & DM_BC_ESSENTIAL)) ++numBC; 4609cbb8f0dSMatthew G. Knepley } 4619cbb8f0dSMatthew G. Knepley /* Add ghost cell boundaries for FVM */ 462baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) if (!isFE[f] && cEndInterior >= 0) ++numBC; 4639cbb8f0dSMatthew G. Knepley ierr = PetscCalloc3(numBC,&bcFields,numBC,&bcPoints,numBC,&bcComps);CHKERRQ(ierr); 4649cbb8f0dSMatthew G. Knepley /* Constrain ghost cells for FV */ 465baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 4669cbb8f0dSMatthew G. Knepley PetscInt *newidx, c; 4679cbb8f0dSMatthew G. Knepley 4689cbb8f0dSMatthew G. Knepley if (isFE[f] || cEndInterior < 0) continue; 4699cbb8f0dSMatthew G. Knepley ierr = PetscMalloc1(cEnd-cEndInterior,&newidx);CHKERRQ(ierr); 4709cbb8f0dSMatthew G. Knepley for (c = cEndInterior; c < cEnd; ++c) newidx[c-cEndInterior] = c; 4719cbb8f0dSMatthew G. Knepley bcFields[bc] = f; 47269293d4bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF, cEnd-cEndInterior, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]);CHKERRQ(ierr); 4739cbb8f0dSMatthew G. Knepley } 4749cbb8f0dSMatthew G. Knepley /* Handle FEM Dirichlet boundaries */ 4759cbb8f0dSMatthew G. Knepley for (bd = 0; bd < numBd; ++bd) { 4769cbb8f0dSMatthew G. Knepley const char *bdLabel; 4779cbb8f0dSMatthew G. Knepley DMLabel label; 4789cbb8f0dSMatthew G. Knepley const PetscInt *comps; 4799cbb8f0dSMatthew G. Knepley const PetscInt *values; 4809cbb8f0dSMatthew G. Knepley PetscInt bd2, field, numComps, numValues; 4819cbb8f0dSMatthew G. Knepley DMBoundaryConditionType type; 4829cbb8f0dSMatthew G. Knepley PetscBool duplicate = PETSC_FALSE; 4839cbb8f0dSMatthew G. Knepley 484e5e52638SMatthew G. Knepley ierr = PetscDSGetBoundary(probBC, bd, &type, NULL, &bdLabel, &field, &numComps, &comps, NULL, &numValues, &values, NULL);CHKERRQ(ierr); 4859cbb8f0dSMatthew G. Knepley ierr = DMGetLabel(dm, bdLabel, &label);CHKERRQ(ierr); 4869cbb8f0dSMatthew G. Knepley if (!isFE[field] || !label) continue; 4879cbb8f0dSMatthew G. Knepley /* Only want to modify label once */ 4889cbb8f0dSMatthew G. Knepley for (bd2 = 0; bd2 < bd; ++bd2) { 4899cbb8f0dSMatthew G. Knepley const char *bdname; 490e5e52638SMatthew G. Knepley ierr = PetscDSGetBoundary(probBC, bd2, NULL, NULL, &bdname, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 4919cbb8f0dSMatthew G. Knepley ierr = PetscStrcmp(bdname, bdLabel, &duplicate);CHKERRQ(ierr); 4929cbb8f0dSMatthew G. Knepley if (duplicate) break; 4939cbb8f0dSMatthew G. Knepley } 4949cbb8f0dSMatthew G. Knepley if (!duplicate && (isFE[field])) { 4959cbb8f0dSMatthew G. Knepley /* don't complete cells, which are just present to give orientation to the boundary */ 4969cbb8f0dSMatthew G. Knepley ierr = DMPlexLabelComplete(dm, label);CHKERRQ(ierr); 4979cbb8f0dSMatthew G. Knepley } 4989cbb8f0dSMatthew G. Knepley /* Filter out cells, if you actually want to constrain cells you need to do things by hand right now */ 4999cbb8f0dSMatthew G. Knepley if (type & DM_BC_ESSENTIAL) { 5009cbb8f0dSMatthew G. Knepley PetscInt *newidx; 5019cbb8f0dSMatthew G. Knepley PetscInt n, newn = 0, p, v; 5029cbb8f0dSMatthew G. Knepley 5039cbb8f0dSMatthew G. Knepley bcFields[bc] = field; 50469293d4bSJed Brown if (numComps) {ierr = ISCreateGeneral(PETSC_COMM_SELF, numComps, comps, PETSC_COPY_VALUES, &bcComps[bc]);CHKERRQ(ierr);} 5059cbb8f0dSMatthew G. Knepley for (v = 0; v < numValues; ++v) { 5069cbb8f0dSMatthew G. Knepley IS tmp; 5079cbb8f0dSMatthew G. Knepley const PetscInt *idx; 5089cbb8f0dSMatthew G. Knepley 5099cbb8f0dSMatthew G. Knepley ierr = DMGetStratumIS(dm, bdLabel, values[v], &tmp);CHKERRQ(ierr); 5109cbb8f0dSMatthew G. Knepley if (!tmp) continue; 5119cbb8f0dSMatthew G. Knepley ierr = ISGetLocalSize(tmp, &n);CHKERRQ(ierr); 5129cbb8f0dSMatthew G. Knepley ierr = ISGetIndices(tmp, &idx);CHKERRQ(ierr); 5139cbb8f0dSMatthew G. Knepley if (isFE[field]) { 5149cbb8f0dSMatthew G. Knepley for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) ++newn; 5159cbb8f0dSMatthew G. Knepley } else { 5169cbb8f0dSMatthew G. Knepley for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) ++newn; 5179cbb8f0dSMatthew G. Knepley } 5189cbb8f0dSMatthew G. Knepley ierr = ISRestoreIndices(tmp, &idx);CHKERRQ(ierr); 5199cbb8f0dSMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 5209cbb8f0dSMatthew G. Knepley } 5219cbb8f0dSMatthew G. Knepley ierr = PetscMalloc1(newn,&newidx);CHKERRQ(ierr); 5229cbb8f0dSMatthew G. Knepley newn = 0; 5239cbb8f0dSMatthew G. Knepley for (v = 0; v < numValues; ++v) { 5249cbb8f0dSMatthew G. Knepley IS tmp; 5259cbb8f0dSMatthew G. Knepley const PetscInt *idx; 5269cbb8f0dSMatthew G. Knepley 5279cbb8f0dSMatthew G. Knepley ierr = DMGetStratumIS(dm, bdLabel, values[v], &tmp);CHKERRQ(ierr); 5289cbb8f0dSMatthew G. Knepley if (!tmp) continue; 5299cbb8f0dSMatthew G. Knepley ierr = ISGetLocalSize(tmp, &n);CHKERRQ(ierr); 5309cbb8f0dSMatthew G. Knepley ierr = ISGetIndices(tmp, &idx);CHKERRQ(ierr); 5319cbb8f0dSMatthew G. Knepley if (isFE[field]) { 5329cbb8f0dSMatthew G. Knepley for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) newidx[newn++] = idx[p]; 5339cbb8f0dSMatthew G. Knepley } else { 5349cbb8f0dSMatthew G. Knepley for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) newidx[newn++] = idx[p]; 5359cbb8f0dSMatthew G. Knepley } 5369cbb8f0dSMatthew G. Knepley ierr = ISRestoreIndices(tmp, &idx);CHKERRQ(ierr); 5379cbb8f0dSMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 5389cbb8f0dSMatthew G. Knepley } 53969293d4bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF, newn, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]);CHKERRQ(ierr); 5409cbb8f0dSMatthew G. Knepley } 5419cbb8f0dSMatthew G. Knepley } 5429cbb8f0dSMatthew G. Knepley /* Handle discretization */ 543baef929fSMatthew G. Knepley ierr = PetscCalloc3(Nf,&labels,Nf,&numComp,Nf*(dim+1),&numDof);CHKERRQ(ierr); 544baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 54544a7f3ddSMatthew G. Knepley labels[f] = dm->fields[f].label; 5469cbb8f0dSMatthew G. Knepley if (isFE[f]) { 54744a7f3ddSMatthew G. Knepley PetscFE fe = (PetscFE) dm->fields[f].disc; 5489cbb8f0dSMatthew G. Knepley const PetscInt *numFieldDof; 549a3254110SMatthew Knepley PetscInt fedim, d; 5509cbb8f0dSMatthew G. Knepley 5519cbb8f0dSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr); 5529cbb8f0dSMatthew G. Knepley ierr = PetscFEGetNumDof(fe, &numFieldDof);CHKERRQ(ierr); 553a3254110SMatthew Knepley ierr = PetscFEGetSpatialDimension(fe, &fedim);CHKERRQ(ierr); 554a3254110SMatthew Knepley for (d = 0; d < PetscMin(dim, fedim)+1; ++d) numDof[f*(dim+1)+d] = numFieldDof[d]; 5559cbb8f0dSMatthew G. Knepley } else { 55644a7f3ddSMatthew G. Knepley PetscFV fv = (PetscFV) dm->fields[f].disc; 5579cbb8f0dSMatthew G. Knepley 5589cbb8f0dSMatthew G. Knepley ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr); 5599cbb8f0dSMatthew G. Knepley numDof[f*(dim+1)+dim] = numComp[f]; 5609cbb8f0dSMatthew G. Knepley } 5619cbb8f0dSMatthew G. Knepley } 562baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 5639cbb8f0dSMatthew G. Knepley PetscInt d; 5649cbb8f0dSMatthew G. Knepley for (d = 1; d < dim; ++d) { 5659cbb8f0dSMatthew 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."); 5669cbb8f0dSMatthew G. Knepley } 5679cbb8f0dSMatthew G. Knepley } 568baef929fSMatthew G. Knepley ierr = DMPlexCreateSection(dm, labels, numComp, numDof, numBC, bcFields, bcComps, bcPoints, NULL, §ion);CHKERRQ(ierr); 569baef929fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 5709cbb8f0dSMatthew G. Knepley PetscFE fe; 5719cbb8f0dSMatthew G. Knepley const char *name; 5729cbb8f0dSMatthew G. Knepley 57344a7f3ddSMatthew G. Knepley if (isFE[f]) { 57444a7f3ddSMatthew G. Knepley ierr = DMGetField(dm, f, NULL, (PetscObject *) &fe);CHKERRQ(ierr); 5759cbb8f0dSMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) fe, &name);CHKERRQ(ierr); 5769cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFieldName(section, f, name);CHKERRQ(ierr); 5779cbb8f0dSMatthew G. Knepley } 57844a7f3ddSMatthew G. Knepley } 57992fd8e1eSJed Brown ierr = DMSetLocalSection(dm, section);CHKERRQ(ierr); 5809cbb8f0dSMatthew G. Knepley ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 5819cbb8f0dSMatthew G. Knepley for (bc = 0; bc < numBC; ++bc) {ierr = ISDestroy(&bcPoints[bc]);CHKERRQ(ierr);ierr = ISDestroy(&bcComps[bc]);CHKERRQ(ierr);} 5829cbb8f0dSMatthew G. Knepley ierr = PetscFree3(bcFields,bcPoints,bcComps);CHKERRQ(ierr); 583baef929fSMatthew G. Knepley ierr = PetscFree3(labels,numComp,numDof);CHKERRQ(ierr); 5849cbb8f0dSMatthew G. Knepley ierr = PetscFree(isFE);CHKERRQ(ierr); 5859cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 5869cbb8f0dSMatthew G. Knepley } 587