1*9cbb8f0dSMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2*9cbb8f0dSMatthew G. Knepley 3*9cbb8f0dSMatthew G. Knepley /* Set the number of dof on each point and separate by fields */ 4*9cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionInitial(DM dm, PetscInt dim, PetscInt numFields,const PetscInt numComp[],const PetscInt numDof[], PetscSection *section) 5*9cbb8f0dSMatthew G. Knepley { 6*9cbb8f0dSMatthew G. Knepley PetscInt *pMax; 7*9cbb8f0dSMatthew G. Knepley PetscInt depth, cellHeight, pStart = 0, pEnd = 0; 8*9cbb8f0dSMatthew G. Knepley PetscInt Nf, p, d, dep, f; 9*9cbb8f0dSMatthew G. Knepley PetscBool *isFE; 10*9cbb8f0dSMatthew G. Knepley PetscErrorCode ierr; 11*9cbb8f0dSMatthew G. Knepley 12*9cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 13*9cbb8f0dSMatthew G. Knepley ierr = PetscMalloc1(numFields, &isFE);CHKERRQ(ierr); 14*9cbb8f0dSMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 15*9cbb8f0dSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 16*9cbb8f0dSMatthew G. Knepley PetscObject obj; 17*9cbb8f0dSMatthew G. Knepley PetscClassId id; 18*9cbb8f0dSMatthew G. Knepley 19*9cbb8f0dSMatthew G. Knepley isFE[f] = PETSC_FALSE; 20*9cbb8f0dSMatthew G. Knepley if (f >= Nf) continue; 21*9cbb8f0dSMatthew G. Knepley ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr); 22*9cbb8f0dSMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 23*9cbb8f0dSMatthew G. Knepley if (id == PETSCFE_CLASSID) {isFE[f] = PETSC_TRUE;} 24*9cbb8f0dSMatthew G. Knepley else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;} 25*9cbb8f0dSMatthew G. Knepley } 26*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), section);CHKERRQ(ierr); 27*9cbb8f0dSMatthew G. Knepley if (numFields > 0) { 28*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetNumFields(*section, numFields);CHKERRQ(ierr); 29*9cbb8f0dSMatthew G. Knepley if (numComp) { 30*9cbb8f0dSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 31*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFieldComponents(*section, f, numComp[f]);CHKERRQ(ierr); 32*9cbb8f0dSMatthew G. Knepley if (isFE[f]) { 33*9cbb8f0dSMatthew G. Knepley PetscFE fe; 34*9cbb8f0dSMatthew G. Knepley PetscDualSpace dspace; 35*9cbb8f0dSMatthew G. Knepley const PetscInt ***perms; 36*9cbb8f0dSMatthew G. Knepley const PetscScalar ***flips; 37*9cbb8f0dSMatthew G. Knepley const PetscInt *numDof; 38*9cbb8f0dSMatthew G. Knepley 39*9cbb8f0dSMatthew G. Knepley ierr = DMGetField(dm,f,(PetscObject *) &fe);CHKERRQ(ierr); 40*9cbb8f0dSMatthew G. Knepley ierr = PetscFEGetDualSpace(fe,&dspace);CHKERRQ(ierr); 41*9cbb8f0dSMatthew G. Knepley ierr = PetscDualSpaceGetSymmetries(dspace,&perms,&flips);CHKERRQ(ierr); 42*9cbb8f0dSMatthew G. Knepley ierr = PetscDualSpaceGetNumDof(dspace,&numDof);CHKERRQ(ierr); 43*9cbb8f0dSMatthew G. Knepley if (perms || flips) { 44*9cbb8f0dSMatthew G. Knepley DM K; 45*9cbb8f0dSMatthew G. Knepley DMLabel depthLabel; 46*9cbb8f0dSMatthew G. Knepley PetscInt depth, h; 47*9cbb8f0dSMatthew G. Knepley PetscSectionSym sym; 48*9cbb8f0dSMatthew G. Knepley 49*9cbb8f0dSMatthew G. Knepley ierr = PetscDualSpaceGetDM(dspace,&K);CHKERRQ(ierr); 50*9cbb8f0dSMatthew G. Knepley ierr = DMPlexGetDepthLabel(dm,&depthLabel);CHKERRQ(ierr); 51*9cbb8f0dSMatthew G. Knepley ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 52*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)*section),depthLabel,&sym);CHKERRQ(ierr); 53*9cbb8f0dSMatthew G. Knepley for (h = 0; h <= depth; h++) { 54*9cbb8f0dSMatthew G. Knepley PetscDualSpace hspace; 55*9cbb8f0dSMatthew G. Knepley PetscInt kStart, kEnd; 56*9cbb8f0dSMatthew G. Knepley PetscInt kConeSize; 57*9cbb8f0dSMatthew G. Knepley const PetscInt **perms0 = NULL; 58*9cbb8f0dSMatthew G. Knepley const PetscScalar **flips0 = NULL; 59*9cbb8f0dSMatthew G. Knepley 60*9cbb8f0dSMatthew G. Knepley ierr = PetscDualSpaceGetHeightSubspace(dspace,h,&hspace);CHKERRQ(ierr); 61*9cbb8f0dSMatthew G. Knepley ierr = DMPlexGetHeightStratum(K,h,&kStart,&kEnd);CHKERRQ(ierr); 62*9cbb8f0dSMatthew G. Knepley if (!hspace) continue; 63*9cbb8f0dSMatthew G. Knepley ierr = PetscDualSpaceGetSymmetries(hspace,&perms,&flips);CHKERRQ(ierr); 64*9cbb8f0dSMatthew G. Knepley if (perms) perms0 = perms[0]; 65*9cbb8f0dSMatthew G. Knepley if (flips) flips0 = flips[0]; 66*9cbb8f0dSMatthew G. Knepley if (!(perms0 || flips0)) continue; 67*9cbb8f0dSMatthew G. Knepley ierr = DMPlexGetConeSize(K,kStart,&kConeSize);CHKERRQ(ierr); 68*9cbb8f0dSMatthew 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); 69*9cbb8f0dSMatthew G. Knepley } 70*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFieldSym(*section,f,sym);CHKERRQ(ierr); 71*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSymDestroy(&sym);CHKERRQ(ierr); 72*9cbb8f0dSMatthew G. Knepley } 73*9cbb8f0dSMatthew G. Knepley } 74*9cbb8f0dSMatthew G. Knepley } 75*9cbb8f0dSMatthew G. Knepley } 76*9cbb8f0dSMatthew G. Knepley } 77*9cbb8f0dSMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 78*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetChart(*section, pStart, pEnd);CHKERRQ(ierr); 79*9cbb8f0dSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 80*9cbb8f0dSMatthew G. Knepley ierr = PetscMalloc1(depth+1,&pMax);CHKERRQ(ierr); 81*9cbb8f0dSMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, depth >= 0 ? &pMax[depth] : NULL, depth>1 ? &pMax[depth-1] : NULL, depth>2 ? &pMax[1] : NULL, &pMax[0]);CHKERRQ(ierr); 82*9cbb8f0dSMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 83*9cbb8f0dSMatthew G. Knepley for (dep = 0; dep <= depth - cellHeight; ++dep) { 84*9cbb8f0dSMatthew G. Knepley d = dim == depth ? dep : (!dep ? 0 : dim); 85*9cbb8f0dSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, dep, &pStart, &pEnd);CHKERRQ(ierr); 86*9cbb8f0dSMatthew G. Knepley pMax[dep] = pMax[dep] < 0 ? pEnd : pMax[dep]; 87*9cbb8f0dSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 88*9cbb8f0dSMatthew G. Knepley PetscInt tot = 0; 89*9cbb8f0dSMatthew G. Knepley 90*9cbb8f0dSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 91*9cbb8f0dSMatthew G. Knepley if (isFE[f] && p >= pMax[dep]) continue; 92*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFieldDof(*section, p, f, numDof[f*(dim+1)+d]);CHKERRQ(ierr); 93*9cbb8f0dSMatthew G. Knepley tot += numDof[f*(dim+1)+d]; 94*9cbb8f0dSMatthew G. Knepley } 95*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetDof(*section, p, tot);CHKERRQ(ierr); 96*9cbb8f0dSMatthew G. Knepley } 97*9cbb8f0dSMatthew G. Knepley } 98*9cbb8f0dSMatthew G. Knepley ierr = PetscFree(pMax);CHKERRQ(ierr); 99*9cbb8f0dSMatthew G. Knepley ierr = PetscFree(isFE);CHKERRQ(ierr); 100*9cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 101*9cbb8f0dSMatthew G. Knepley } 102*9cbb8f0dSMatthew G. Knepley 103*9cbb8f0dSMatthew G. Knepley /* Set the number of dof on each point and separate by fields 104*9cbb8f0dSMatthew G. Knepley If bcComps is NULL or the IS is NULL, constrain every dof on the point 105*9cbb8f0dSMatthew G. Knepley */ 106*9cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCDof(DM dm, PetscInt numBC, const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section) 107*9cbb8f0dSMatthew G. Knepley { 108*9cbb8f0dSMatthew G. Knepley PetscInt numFields; 109*9cbb8f0dSMatthew G. Knepley PetscInt bc; 110*9cbb8f0dSMatthew G. Knepley PetscSection aSec; 111*9cbb8f0dSMatthew G. Knepley PetscErrorCode ierr; 112*9cbb8f0dSMatthew G. Knepley 113*9cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 114*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 115*9cbb8f0dSMatthew G. Knepley for (bc = 0; bc < numBC; ++bc) { 116*9cbb8f0dSMatthew G. Knepley PetscInt field = 0; 117*9cbb8f0dSMatthew G. Knepley const PetscInt *comp; 118*9cbb8f0dSMatthew G. Knepley const PetscInt *idx; 119*9cbb8f0dSMatthew G. Knepley PetscInt Nc = -1, n, i; 120*9cbb8f0dSMatthew G. Knepley 121*9cbb8f0dSMatthew G. Knepley if (numFields) field = bcField[bc]; 122*9cbb8f0dSMatthew G. Knepley if (bcComps && bcComps[bc]) {ierr = ISGetLocalSize(bcComps[bc], &Nc);CHKERRQ(ierr);} 123*9cbb8f0dSMatthew G. Knepley if (bcComps && bcComps[bc]) {ierr = ISGetIndices(bcComps[bc], &comp);CHKERRQ(ierr);} 124*9cbb8f0dSMatthew G. Knepley ierr = ISGetLocalSize(bcPoints[bc], &n);CHKERRQ(ierr); 125*9cbb8f0dSMatthew G. Knepley ierr = ISGetIndices(bcPoints[bc], &idx);CHKERRQ(ierr); 126*9cbb8f0dSMatthew G. Knepley for (i = 0; i < n; ++i) { 127*9cbb8f0dSMatthew G. Knepley const PetscInt p = idx[i]; 128*9cbb8f0dSMatthew G. Knepley PetscInt numConst; 129*9cbb8f0dSMatthew G. Knepley 130*9cbb8f0dSMatthew G. Knepley if (numFields) { 131*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldDof(section, p, field, &numConst);CHKERRQ(ierr); 132*9cbb8f0dSMatthew G. Knepley } else { 133*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetDof(section, p, &numConst);CHKERRQ(ierr); 134*9cbb8f0dSMatthew G. Knepley } 135*9cbb8f0dSMatthew G. Knepley /* If Nc < 0, constrain every dof on the point */ 136*9cbb8f0dSMatthew G. Knepley /* TODO: Matt, this only works if there is one node on the point. We need to handle numDofs > NumComponents */ 137*9cbb8f0dSMatthew G. Knepley if (Nc > 0) numConst = PetscMin(numConst, Nc); 138*9cbb8f0dSMatthew G. Knepley if (numFields) {ierr = PetscSectionAddFieldConstraintDof(section, p, field, numConst);CHKERRQ(ierr);} 139*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionAddConstraintDof(section, p, numConst);CHKERRQ(ierr); 140*9cbb8f0dSMatthew G. Knepley } 141*9cbb8f0dSMatthew G. Knepley ierr = ISRestoreIndices(bcPoints[bc], &idx);CHKERRQ(ierr); 142*9cbb8f0dSMatthew G. Knepley if (bcComps && bcComps[bc]) {ierr = ISRestoreIndices(bcComps[bc], &comp);CHKERRQ(ierr);} 143*9cbb8f0dSMatthew G. Knepley } 144*9cbb8f0dSMatthew G. Knepley ierr = DMPlexGetAnchors(dm, &aSec, NULL);CHKERRQ(ierr); 145*9cbb8f0dSMatthew G. Knepley if (aSec) { 146*9cbb8f0dSMatthew G. Knepley PetscInt aStart, aEnd, a; 147*9cbb8f0dSMatthew G. Knepley 148*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetChart(aSec, &aStart, &aEnd);CHKERRQ(ierr); 149*9cbb8f0dSMatthew G. Knepley for (a = aStart; a < aEnd; a++) { 150*9cbb8f0dSMatthew G. Knepley PetscInt dof, f; 151*9cbb8f0dSMatthew G. Knepley 152*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetDof(aSec, a, &dof);CHKERRQ(ierr); 153*9cbb8f0dSMatthew G. Knepley if (dof) { 154*9cbb8f0dSMatthew G. Knepley /* if there are point-to-point constraints, then all dofs are constrained */ 155*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetDof(section, a, &dof);CHKERRQ(ierr); 156*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetConstraintDof(section, a, dof);CHKERRQ(ierr); 157*9cbb8f0dSMatthew G. Knepley for (f = 0; f < numFields; f++) { 158*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldDof(section, a, f, &dof);CHKERRQ(ierr); 159*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFieldConstraintDof(section, a, f, dof);CHKERRQ(ierr); 160*9cbb8f0dSMatthew G. Knepley } 161*9cbb8f0dSMatthew G. Knepley } 162*9cbb8f0dSMatthew G. Knepley } 163*9cbb8f0dSMatthew G. Knepley } 164*9cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 165*9cbb8f0dSMatthew G. Knepley } 166*9cbb8f0dSMatthew G. Knepley 167*9cbb8f0dSMatthew G. Knepley /* Set the constrained field indices on each point 168*9cbb8f0dSMatthew G. Knepley If bcComps is NULL or the IS is NULL, constrain every dof on the point 169*9cbb8f0dSMatthew G. Knepley */ 170*9cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCIndicesField(DM dm, PetscInt numBC,const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section) 171*9cbb8f0dSMatthew G. Knepley { 172*9cbb8f0dSMatthew G. Knepley PetscSection aSec; 173*9cbb8f0dSMatthew G. Knepley PetscInt *indices; 174*9cbb8f0dSMatthew G. Knepley PetscInt numFields, cdof, maxDof = 0, pStart, pEnd, p, bc, f, d; 175*9cbb8f0dSMatthew G. Knepley PetscErrorCode ierr; 176*9cbb8f0dSMatthew G. Knepley 177*9cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 178*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 179*9cbb8f0dSMatthew G. Knepley if (!numFields) PetscFunctionReturn(0); 180*9cbb8f0dSMatthew G. Knepley /* Initialize all field indices to -1 */ 181*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 182*9cbb8f0dSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) {ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); maxDof = PetscMax(maxDof, cdof);} 183*9cbb8f0dSMatthew G. Knepley ierr = PetscMalloc1(maxDof, &indices);CHKERRQ(ierr); 184*9cbb8f0dSMatthew G. Knepley for (d = 0; d < maxDof; ++d) indices[d] = -1; 185*9cbb8f0dSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) for (f = 0; f < numFields; ++f) {ierr = PetscSectionSetFieldConstraintIndices(section, p, f, indices);CHKERRQ(ierr);} 186*9cbb8f0dSMatthew G. Knepley /* Handle BC constraints */ 187*9cbb8f0dSMatthew G. Knepley for (bc = 0; bc < numBC; ++bc) { 188*9cbb8f0dSMatthew G. Knepley const PetscInt field = bcField[bc]; 189*9cbb8f0dSMatthew G. Knepley const PetscInt *comp, *idx; 190*9cbb8f0dSMatthew G. Knepley PetscInt Nc = -1, n, i; 191*9cbb8f0dSMatthew G. Knepley 192*9cbb8f0dSMatthew G. Knepley if (bcComps && bcComps[bc]) {ierr = ISGetLocalSize(bcComps[bc], &Nc);CHKERRQ(ierr);} 193*9cbb8f0dSMatthew G. Knepley if (bcComps && bcComps[bc]) {ierr = ISGetIndices(bcComps[bc], &comp);CHKERRQ(ierr);} 194*9cbb8f0dSMatthew G. Knepley ierr = ISGetLocalSize(bcPoints[bc], &n);CHKERRQ(ierr); 195*9cbb8f0dSMatthew G. Knepley ierr = ISGetIndices(bcPoints[bc], &idx);CHKERRQ(ierr); 196*9cbb8f0dSMatthew G. Knepley for (i = 0; i < n; ++i) { 197*9cbb8f0dSMatthew G. Knepley const PetscInt p = idx[i]; 198*9cbb8f0dSMatthew G. Knepley const PetscInt *find; 199*9cbb8f0dSMatthew G. Knepley PetscInt fdof, fcdof, c; 200*9cbb8f0dSMatthew G. Knepley 201*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldDof(section, p, field, &fdof);CHKERRQ(ierr); 202*9cbb8f0dSMatthew G. Knepley if (!fdof) continue; 203*9cbb8f0dSMatthew G. Knepley if (Nc < 0) { 204*9cbb8f0dSMatthew G. Knepley for (d = 0; d < fdof; ++d) indices[d] = d; 205*9cbb8f0dSMatthew G. Knepley fcdof = fdof; 206*9cbb8f0dSMatthew G. Knepley } else { 207*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintDof(section, p, field, &fcdof);CHKERRQ(ierr); 208*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintIndices(section, p, field, &find);CHKERRQ(ierr); 209*9cbb8f0dSMatthew G. Knepley for (d = 0; d < fcdof; ++d) {if (find[d] < 0) break; indices[d] = find[d];} 210*9cbb8f0dSMatthew G. Knepley for (c = 0; c < Nc; ++c) indices[d++] = comp[c]; 211*9cbb8f0dSMatthew G. Knepley ierr = PetscSortRemoveDupsInt(&d, indices);CHKERRQ(ierr); 212*9cbb8f0dSMatthew G. Knepley for (c = d; c < fcdof; ++c) indices[c] = -1; 213*9cbb8f0dSMatthew G. Knepley fcdof = d; 214*9cbb8f0dSMatthew G. Knepley } 215*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFieldConstraintDof(section, p, field, fcdof);CHKERRQ(ierr); 216*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFieldConstraintIndices(section, p, field, indices);CHKERRQ(ierr); 217*9cbb8f0dSMatthew G. Knepley } 218*9cbb8f0dSMatthew G. Knepley if (bcComps && bcComps[bc]) {ierr = ISRestoreIndices(bcComps[bc], &comp);CHKERRQ(ierr);} 219*9cbb8f0dSMatthew G. Knepley ierr = ISRestoreIndices(bcPoints[bc], &idx);CHKERRQ(ierr); 220*9cbb8f0dSMatthew G. Knepley } 221*9cbb8f0dSMatthew G. Knepley /* Handle anchors */ 222*9cbb8f0dSMatthew G. Knepley ierr = DMPlexGetAnchors(dm, &aSec, NULL);CHKERRQ(ierr); 223*9cbb8f0dSMatthew G. Knepley if (aSec) { 224*9cbb8f0dSMatthew G. Knepley PetscInt aStart, aEnd, a; 225*9cbb8f0dSMatthew G. Knepley 226*9cbb8f0dSMatthew G. Knepley for (d = 0; d < maxDof; ++d) indices[d] = d; 227*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetChart(aSec, &aStart, &aEnd);CHKERRQ(ierr); 228*9cbb8f0dSMatthew G. Knepley for (a = aStart; a < aEnd; a++) { 229*9cbb8f0dSMatthew G. Knepley PetscInt dof, f; 230*9cbb8f0dSMatthew G. Knepley 231*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetDof(aSec, a, &dof);CHKERRQ(ierr); 232*9cbb8f0dSMatthew G. Knepley if (dof) { 233*9cbb8f0dSMatthew G. Knepley /* if there are point-to-point constraints, then all dofs are constrained */ 234*9cbb8f0dSMatthew G. Knepley for (f = 0; f < numFields; f++) { 235*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFieldConstraintIndices(section, a, f, indices);CHKERRQ(ierr); 236*9cbb8f0dSMatthew G. Knepley } 237*9cbb8f0dSMatthew G. Knepley } 238*9cbb8f0dSMatthew G. Knepley } 239*9cbb8f0dSMatthew G. Knepley } 240*9cbb8f0dSMatthew G. Knepley ierr = PetscFree(indices);CHKERRQ(ierr); 241*9cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 242*9cbb8f0dSMatthew G. Knepley } 243*9cbb8f0dSMatthew G. Knepley 244*9cbb8f0dSMatthew G. Knepley /* Set the constrained indices on each point */ 245*9cbb8f0dSMatthew G. Knepley static PetscErrorCode DMPlexCreateSectionBCIndices(DM dm, PetscSection section) 246*9cbb8f0dSMatthew G. Knepley { 247*9cbb8f0dSMatthew G. Knepley PetscInt *indices; 248*9cbb8f0dSMatthew G. Knepley PetscInt numFields, maxDof, pStart, pEnd, p, f, d; 249*9cbb8f0dSMatthew G. Knepley PetscErrorCode ierr; 250*9cbb8f0dSMatthew G. Knepley 251*9cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 252*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 253*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetMaxDof(section, &maxDof);CHKERRQ(ierr); 254*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 255*9cbb8f0dSMatthew G. Knepley ierr = PetscMalloc1(maxDof, &indices);CHKERRQ(ierr); 256*9cbb8f0dSMatthew G. Knepley for (d = 0; d < maxDof; ++d) indices[d] = -1; 257*9cbb8f0dSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 258*9cbb8f0dSMatthew G. Knepley PetscInt cdof, d; 259*9cbb8f0dSMatthew G. Knepley 260*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 261*9cbb8f0dSMatthew G. Knepley if (cdof) { 262*9cbb8f0dSMatthew G. Knepley if (numFields) { 263*9cbb8f0dSMatthew G. Knepley PetscInt numConst = 0, foff = 0; 264*9cbb8f0dSMatthew G. Knepley 265*9cbb8f0dSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 266*9cbb8f0dSMatthew G. Knepley const PetscInt *find; 267*9cbb8f0dSMatthew G. Knepley PetscInt fcdof, fdof; 268*9cbb8f0dSMatthew G. Knepley 269*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr); 270*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr); 271*9cbb8f0dSMatthew G. Knepley /* Change constraint numbering from field component to local dof number */ 272*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintIndices(section, p, f, &find);CHKERRQ(ierr); 273*9cbb8f0dSMatthew G. Knepley for (d = 0; d < fcdof; ++d) indices[numConst+d] = find[d] + foff; 274*9cbb8f0dSMatthew G. Knepley numConst += fcdof; 275*9cbb8f0dSMatthew G. Knepley foff += fdof; 276*9cbb8f0dSMatthew G. Knepley } 277*9cbb8f0dSMatthew G. Knepley if (cdof != numConst) {ierr = PetscSectionSetConstraintDof(section, p, numConst);CHKERRQ(ierr);} 278*9cbb8f0dSMatthew G. Knepley } else { 279*9cbb8f0dSMatthew G. Knepley for (d = 0; d < cdof; ++d) indices[d] = d; 280*9cbb8f0dSMatthew G. Knepley } 281*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetConstraintIndices(section, p, indices);CHKERRQ(ierr); 282*9cbb8f0dSMatthew G. Knepley } 283*9cbb8f0dSMatthew G. Knepley } 284*9cbb8f0dSMatthew G. Knepley ierr = PetscFree(indices);CHKERRQ(ierr); 285*9cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 286*9cbb8f0dSMatthew G. Knepley } 287*9cbb8f0dSMatthew G. Knepley 288*9cbb8f0dSMatthew G. Knepley /*@C 289*9cbb8f0dSMatthew G. Knepley DMPlexCreateSection - Create a PetscSection based upon the dof layout specification provided. 290*9cbb8f0dSMatthew G. Knepley 291*9cbb8f0dSMatthew G. Knepley Not Collective 292*9cbb8f0dSMatthew G. Knepley 293*9cbb8f0dSMatthew G. Knepley Input Parameters: 294*9cbb8f0dSMatthew G. Knepley + dm - The DMPlex object 295*9cbb8f0dSMatthew G. Knepley . dim - The spatial dimension of the problem 296*9cbb8f0dSMatthew G. Knepley . numFields - The number of fields in the problem 297*9cbb8f0dSMatthew G. Knepley . numComp - An array of size numFields that holds the number of components for each field 298*9cbb8f0dSMatthew 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 299*9cbb8f0dSMatthew G. Knepley . numBC - The number of boundary conditions 300*9cbb8f0dSMatthew G. Knepley . bcField - An array of size numBC giving the field number for each boundry condition 301*9cbb8f0dSMatthew G. Knepley . bcComps - [Optional] An array of size numBC giving an IS holding the field components to which each boundary condition applies 302*9cbb8f0dSMatthew G. Knepley . bcPoints - An array of size numBC giving an IS holding the Plex points to which each boundary condition applies 303*9cbb8f0dSMatthew G. Knepley - perm - Optional permutation of the chart, or NULL 304*9cbb8f0dSMatthew G. Knepley 305*9cbb8f0dSMatthew G. Knepley Output Parameter: 306*9cbb8f0dSMatthew G. Knepley . section - The PetscSection object 307*9cbb8f0dSMatthew G. Knepley 308*9cbb8f0dSMatthew G. Knepley Notes: 309*9cbb8f0dSMatthew 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 310*9cbb8f0dSMatthew G. Knepley number of dof for field 0 on each edge. 311*9cbb8f0dSMatthew G. Knepley 312*9cbb8f0dSMatthew G. Knepley The chart permutation is the same one set using PetscSectionSetPermutation() 313*9cbb8f0dSMatthew G. Knepley 314*9cbb8f0dSMatthew G. Knepley Level: developer 315*9cbb8f0dSMatthew G. Knepley 316*9cbb8f0dSMatthew G. Knepley Fortran Notes: 317*9cbb8f0dSMatthew G. Knepley A Fortran 90 version is available as DMPlexCreateSectionF90() 318*9cbb8f0dSMatthew G. Knepley 319*9cbb8f0dSMatthew G. Knepley .keywords: mesh, elements 320*9cbb8f0dSMatthew G. Knepley .seealso: DMPlexCreate(), PetscSectionCreate(), PetscSectionSetPermutation() 321*9cbb8f0dSMatthew G. Knepley @*/ 322*9cbb8f0dSMatthew G. Knepley PetscErrorCode DMPlexCreateSection(DM dm, PetscInt dim, PetscInt numFields,const PetscInt numComp[],const PetscInt numDof[], PetscInt numBC,const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], IS perm, PetscSection *section) 323*9cbb8f0dSMatthew G. Knepley { 324*9cbb8f0dSMatthew G. Knepley PetscSection aSec; 325*9cbb8f0dSMatthew G. Knepley PetscErrorCode ierr; 326*9cbb8f0dSMatthew G. Knepley 327*9cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 328*9cbb8f0dSMatthew G. Knepley ierr = DMPlexCreateSectionInitial(dm, dim, numFields, numComp, numDof, section);CHKERRQ(ierr); 329*9cbb8f0dSMatthew G. Knepley ierr = DMPlexCreateSectionBCDof(dm, numBC, bcField, bcComps, bcPoints, *section);CHKERRQ(ierr); 330*9cbb8f0dSMatthew G. Knepley if (perm) {ierr = PetscSectionSetPermutation(*section, perm);CHKERRQ(ierr);} 331*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFromOptions(*section);CHKERRQ(ierr); 332*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 333*9cbb8f0dSMatthew G. Knepley ierr = DMPlexGetAnchors(dm,&aSec,NULL);CHKERRQ(ierr); 334*9cbb8f0dSMatthew G. Knepley if (numBC || aSec) { 335*9cbb8f0dSMatthew G. Knepley ierr = DMPlexCreateSectionBCIndicesField(dm, numBC, bcField, bcComps, bcPoints, *section);CHKERRQ(ierr); 336*9cbb8f0dSMatthew G. Knepley ierr = DMPlexCreateSectionBCIndices(dm, *section);CHKERRQ(ierr); 337*9cbb8f0dSMatthew G. Knepley } 338*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionViewFromOptions(*section,NULL,"-section_view");CHKERRQ(ierr); 339*9cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 340*9cbb8f0dSMatthew G. Knepley } 341*9cbb8f0dSMatthew G. Knepley 342*9cbb8f0dSMatthew G. Knepley PetscErrorCode DMCreateDefaultSection_Plex(DM dm) 343*9cbb8f0dSMatthew G. Knepley { 344*9cbb8f0dSMatthew G. Knepley PetscSection section; 345*9cbb8f0dSMatthew G. Knepley IS *bcPoints, *bcComps; 346*9cbb8f0dSMatthew G. Knepley PetscBool *isFE; 347*9cbb8f0dSMatthew G. Knepley PetscInt *bcFields, *numComp, *numDof; 348*9cbb8f0dSMatthew G. Knepley PetscInt depth, dim, numBd, numBC = 0, numFields, bd, bc = 0, f; 349*9cbb8f0dSMatthew G. Knepley PetscInt cStart, cEnd, cEndInterior; 350*9cbb8f0dSMatthew G. Knepley PetscErrorCode ierr; 351*9cbb8f0dSMatthew G. Knepley 352*9cbb8f0dSMatthew G. Knepley PetscFunctionBegin; 353*9cbb8f0dSMatthew G. Knepley ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr); 354*9cbb8f0dSMatthew G. Knepley /* FE and FV boundary conditions are handled slightly differently */ 355*9cbb8f0dSMatthew G. Knepley ierr = PetscMalloc1(numFields, &isFE);CHKERRQ(ierr); 356*9cbb8f0dSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 357*9cbb8f0dSMatthew G. Knepley PetscObject obj; 358*9cbb8f0dSMatthew G. Knepley PetscClassId id; 359*9cbb8f0dSMatthew G. Knepley 360*9cbb8f0dSMatthew G. Knepley ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr); 361*9cbb8f0dSMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 362*9cbb8f0dSMatthew G. Knepley if (id == PETSCFE_CLASSID) {isFE[f] = PETSC_TRUE;} 363*9cbb8f0dSMatthew G. Knepley else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;} 364*9cbb8f0dSMatthew G. Knepley else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f); 365*9cbb8f0dSMatthew G. Knepley } 366*9cbb8f0dSMatthew G. Knepley /* Allocate boundary point storage for FEM boundaries */ 367*9cbb8f0dSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 368*9cbb8f0dSMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 369*9cbb8f0dSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 370*9cbb8f0dSMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 371*9cbb8f0dSMatthew G. Knepley ierr = PetscDSGetNumBoundary(dm->prob, &numBd);CHKERRQ(ierr); 372*9cbb8f0dSMatthew G. Knepley if (!numFields && numBd) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "number of fields is zero and number of boundary conditions is nonzero (this should never happen)"); 373*9cbb8f0dSMatthew G. Knepley for (bd = 0; bd < numBd; ++bd) { 374*9cbb8f0dSMatthew G. Knepley PetscInt field; 375*9cbb8f0dSMatthew G. Knepley DMBoundaryConditionType type; 376*9cbb8f0dSMatthew G. Knepley const char *labelName; 377*9cbb8f0dSMatthew G. Knepley DMLabel label; 378*9cbb8f0dSMatthew G. Knepley 379*9cbb8f0dSMatthew G. Knepley ierr = PetscDSGetBoundary(dm->prob, bd, &type, NULL, &labelName, &field, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 380*9cbb8f0dSMatthew G. Knepley ierr = DMGetLabel(dm,labelName,&label);CHKERRQ(ierr); 381*9cbb8f0dSMatthew G. Knepley if (label && isFE[field] && (type & DM_BC_ESSENTIAL)) ++numBC; 382*9cbb8f0dSMatthew G. Knepley } 383*9cbb8f0dSMatthew G. Knepley /* Add ghost cell boundaries for FVM */ 384*9cbb8f0dSMatthew G. Knepley for (f = 0; f < numFields; ++f) if (!isFE[f] && cEndInterior >= 0) ++numBC; 385*9cbb8f0dSMatthew G. Knepley ierr = PetscCalloc3(numBC,&bcFields,numBC,&bcPoints,numBC,&bcComps);CHKERRQ(ierr); 386*9cbb8f0dSMatthew G. Knepley /* Constrain ghost cells for FV */ 387*9cbb8f0dSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 388*9cbb8f0dSMatthew G. Knepley PetscInt *newidx, c; 389*9cbb8f0dSMatthew G. Knepley 390*9cbb8f0dSMatthew G. Knepley if (isFE[f] || cEndInterior < 0) continue; 391*9cbb8f0dSMatthew G. Knepley ierr = PetscMalloc1(cEnd-cEndInterior,&newidx);CHKERRQ(ierr); 392*9cbb8f0dSMatthew G. Knepley for (c = cEndInterior; c < cEnd; ++c) newidx[c-cEndInterior] = c; 393*9cbb8f0dSMatthew G. Knepley bcFields[bc] = f; 394*9cbb8f0dSMatthew G. Knepley ierr = ISCreateGeneral(PetscObjectComm((PetscObject) dm), cEnd-cEndInterior, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]);CHKERRQ(ierr); 395*9cbb8f0dSMatthew G. Knepley } 396*9cbb8f0dSMatthew G. Knepley /* Handle FEM Dirichlet boundaries */ 397*9cbb8f0dSMatthew G. Knepley for (bd = 0; bd < numBd; ++bd) { 398*9cbb8f0dSMatthew G. Knepley const char *bdLabel; 399*9cbb8f0dSMatthew G. Knepley DMLabel label; 400*9cbb8f0dSMatthew G. Knepley const PetscInt *comps; 401*9cbb8f0dSMatthew G. Knepley const PetscInt *values; 402*9cbb8f0dSMatthew G. Knepley PetscInt bd2, field, numComps, numValues; 403*9cbb8f0dSMatthew G. Knepley DMBoundaryConditionType type; 404*9cbb8f0dSMatthew G. Knepley PetscBool duplicate = PETSC_FALSE; 405*9cbb8f0dSMatthew G. Knepley 406*9cbb8f0dSMatthew G. Knepley ierr = PetscDSGetBoundary(dm->prob, bd, &type, NULL, &bdLabel, &field, &numComps, &comps, NULL, &numValues, &values, NULL);CHKERRQ(ierr); 407*9cbb8f0dSMatthew G. Knepley ierr = DMGetLabel(dm, bdLabel, &label);CHKERRQ(ierr); 408*9cbb8f0dSMatthew G. Knepley if (!isFE[field] || !label) continue; 409*9cbb8f0dSMatthew G. Knepley /* Only want to modify label once */ 410*9cbb8f0dSMatthew G. Knepley for (bd2 = 0; bd2 < bd; ++bd2) { 411*9cbb8f0dSMatthew G. Knepley const char *bdname; 412*9cbb8f0dSMatthew G. Knepley ierr = PetscDSGetBoundary(dm->prob, bd2, NULL, NULL, &bdname, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 413*9cbb8f0dSMatthew G. Knepley ierr = PetscStrcmp(bdname, bdLabel, &duplicate);CHKERRQ(ierr); 414*9cbb8f0dSMatthew G. Knepley if (duplicate) break; 415*9cbb8f0dSMatthew G. Knepley } 416*9cbb8f0dSMatthew G. Knepley if (!duplicate && (isFE[field])) { 417*9cbb8f0dSMatthew G. Knepley /* don't complete cells, which are just present to give orientation to the boundary */ 418*9cbb8f0dSMatthew G. Knepley ierr = DMPlexLabelComplete(dm, label);CHKERRQ(ierr); 419*9cbb8f0dSMatthew G. Knepley } 420*9cbb8f0dSMatthew G. Knepley /* Filter out cells, if you actually want to constrain cells you need to do things by hand right now */ 421*9cbb8f0dSMatthew G. Knepley if (type & DM_BC_ESSENTIAL) { 422*9cbb8f0dSMatthew G. Knepley PetscInt *newidx; 423*9cbb8f0dSMatthew G. Knepley PetscInt n, newn = 0, p, v; 424*9cbb8f0dSMatthew G. Knepley 425*9cbb8f0dSMatthew G. Knepley bcFields[bc] = field; 426*9cbb8f0dSMatthew G. Knepley if (numComps) {ierr = ISCreateGeneral(PetscObjectComm((PetscObject) dm), numComps, comps, PETSC_COPY_VALUES, &bcComps[bc]);CHKERRQ(ierr);} 427*9cbb8f0dSMatthew G. Knepley for (v = 0; v < numValues; ++v) { 428*9cbb8f0dSMatthew G. Knepley IS tmp; 429*9cbb8f0dSMatthew G. Knepley const PetscInt *idx; 430*9cbb8f0dSMatthew G. Knepley 431*9cbb8f0dSMatthew G. Knepley ierr = DMGetStratumIS(dm, bdLabel, values[v], &tmp);CHKERRQ(ierr); 432*9cbb8f0dSMatthew G. Knepley if (!tmp) continue; 433*9cbb8f0dSMatthew G. Knepley ierr = ISGetLocalSize(tmp, &n);CHKERRQ(ierr); 434*9cbb8f0dSMatthew G. Knepley ierr = ISGetIndices(tmp, &idx);CHKERRQ(ierr); 435*9cbb8f0dSMatthew G. Knepley if (isFE[field]) { 436*9cbb8f0dSMatthew G. Knepley for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) ++newn; 437*9cbb8f0dSMatthew G. Knepley } else { 438*9cbb8f0dSMatthew G. Knepley for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) ++newn; 439*9cbb8f0dSMatthew G. Knepley } 440*9cbb8f0dSMatthew G. Knepley ierr = ISRestoreIndices(tmp, &idx);CHKERRQ(ierr); 441*9cbb8f0dSMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 442*9cbb8f0dSMatthew G. Knepley } 443*9cbb8f0dSMatthew G. Knepley ierr = PetscMalloc1(newn,&newidx);CHKERRQ(ierr); 444*9cbb8f0dSMatthew G. Knepley newn = 0; 445*9cbb8f0dSMatthew G. Knepley for (v = 0; v < numValues; ++v) { 446*9cbb8f0dSMatthew G. Knepley IS tmp; 447*9cbb8f0dSMatthew G. Knepley const PetscInt *idx; 448*9cbb8f0dSMatthew G. Knepley 449*9cbb8f0dSMatthew G. Knepley ierr = DMGetStratumIS(dm, bdLabel, values[v], &tmp);CHKERRQ(ierr); 450*9cbb8f0dSMatthew G. Knepley if (!tmp) continue; 451*9cbb8f0dSMatthew G. Knepley ierr = ISGetLocalSize(tmp, &n);CHKERRQ(ierr); 452*9cbb8f0dSMatthew G. Knepley ierr = ISGetIndices(tmp, &idx);CHKERRQ(ierr); 453*9cbb8f0dSMatthew G. Knepley if (isFE[field]) { 454*9cbb8f0dSMatthew G. Knepley for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) newidx[newn++] = idx[p]; 455*9cbb8f0dSMatthew G. Knepley } else { 456*9cbb8f0dSMatthew G. Knepley for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) newidx[newn++] = idx[p]; 457*9cbb8f0dSMatthew G. Knepley } 458*9cbb8f0dSMatthew G. Knepley ierr = ISRestoreIndices(tmp, &idx);CHKERRQ(ierr); 459*9cbb8f0dSMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 460*9cbb8f0dSMatthew G. Knepley } 461*9cbb8f0dSMatthew G. Knepley ierr = ISCreateGeneral(PetscObjectComm((PetscObject) dm), newn, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]);CHKERRQ(ierr); 462*9cbb8f0dSMatthew G. Knepley } 463*9cbb8f0dSMatthew G. Knepley } 464*9cbb8f0dSMatthew G. Knepley /* Handle discretization */ 465*9cbb8f0dSMatthew G. Knepley ierr = PetscCalloc2(numFields,&numComp,numFields*(dim+1),&numDof);CHKERRQ(ierr); 466*9cbb8f0dSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 467*9cbb8f0dSMatthew G. Knepley PetscObject obj; 468*9cbb8f0dSMatthew G. Knepley 469*9cbb8f0dSMatthew G. Knepley ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr); 470*9cbb8f0dSMatthew G. Knepley if (isFE[f]) { 471*9cbb8f0dSMatthew G. Knepley PetscFE fe = (PetscFE) obj; 472*9cbb8f0dSMatthew G. Knepley const PetscInt *numFieldDof; 473*9cbb8f0dSMatthew G. Knepley PetscInt d; 474*9cbb8f0dSMatthew G. Knepley 475*9cbb8f0dSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr); 476*9cbb8f0dSMatthew G. Knepley ierr = PetscFEGetNumDof(fe, &numFieldDof);CHKERRQ(ierr); 477*9cbb8f0dSMatthew G. Knepley for (d = 0; d < dim+1; ++d) numDof[f*(dim+1)+d] = numFieldDof[d]; 478*9cbb8f0dSMatthew G. Knepley } else { 479*9cbb8f0dSMatthew G. Knepley PetscFV fv = (PetscFV) obj; 480*9cbb8f0dSMatthew G. Knepley 481*9cbb8f0dSMatthew G. Knepley ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr); 482*9cbb8f0dSMatthew G. Knepley numDof[f*(dim+1)+dim] = numComp[f]; 483*9cbb8f0dSMatthew G. Knepley } 484*9cbb8f0dSMatthew G. Knepley } 485*9cbb8f0dSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 486*9cbb8f0dSMatthew G. Knepley PetscInt d; 487*9cbb8f0dSMatthew G. Knepley for (d = 1; d < dim; ++d) { 488*9cbb8f0dSMatthew 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."); 489*9cbb8f0dSMatthew G. Knepley } 490*9cbb8f0dSMatthew G. Knepley } 491*9cbb8f0dSMatthew G. Knepley ierr = DMPlexCreateSection(dm, dim, numFields, numComp, numDof, numBC, bcFields, bcComps, bcPoints, NULL, §ion);CHKERRQ(ierr); 492*9cbb8f0dSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 493*9cbb8f0dSMatthew G. Knepley PetscFE fe; 494*9cbb8f0dSMatthew G. Knepley const char *name; 495*9cbb8f0dSMatthew G. Knepley 496*9cbb8f0dSMatthew G. Knepley ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr); 497*9cbb8f0dSMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) fe, &name);CHKERRQ(ierr); 498*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionSetFieldName(section, f, name);CHKERRQ(ierr); 499*9cbb8f0dSMatthew G. Knepley } 500*9cbb8f0dSMatthew G. Knepley ierr = DMSetSection(dm, section);CHKERRQ(ierr); 501*9cbb8f0dSMatthew G. Knepley ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 502*9cbb8f0dSMatthew G. Knepley for (bc = 0; bc < numBC; ++bc) {ierr = ISDestroy(&bcPoints[bc]);CHKERRQ(ierr);ierr = ISDestroy(&bcComps[bc]);CHKERRQ(ierr);} 503*9cbb8f0dSMatthew G. Knepley ierr = PetscFree3(bcFields,bcPoints,bcComps);CHKERRQ(ierr); 504*9cbb8f0dSMatthew G. Knepley ierr = PetscFree2(numComp,numDof);CHKERRQ(ierr); 505*9cbb8f0dSMatthew G. Knepley ierr = PetscFree(isFE);CHKERRQ(ierr); 506*9cbb8f0dSMatthew G. Knepley PetscFunctionReturn(0); 507*9cbb8f0dSMatthew G. Knepley } 508