1*e83e87a5Sjeremylt #include "../include/petscutils.h" 2*e83e87a5Sjeremylt 3*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 4*e83e87a5Sjeremylt // Convert PETSc MemType to libCEED MemType 5*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 6*e83e87a5Sjeremylt CeedMemType MemTypeP2C(PetscMemType mtype) { 7*e83e87a5Sjeremylt return PetscMemTypeDevice(mtype) ? CEED_MEM_DEVICE : CEED_MEM_HOST; 8*e83e87a5Sjeremylt } 9*e83e87a5Sjeremylt 10*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 11*e83e87a5Sjeremylt // Utility function taken from petsc/src/dm/impls/plex/examples/tutorials/ex7.c 12*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 13*e83e87a5Sjeremylt PetscErrorCode ProjectToUnitSphere(DM dm) { 14*e83e87a5Sjeremylt Vec coordinates; 15*e83e87a5Sjeremylt PetscScalar *coords; 16*e83e87a5Sjeremylt PetscInt Nv, v, dim, d; 17*e83e87a5Sjeremylt PetscErrorCode ierr; 18*e83e87a5Sjeremylt 19*e83e87a5Sjeremylt PetscFunctionBeginUser; 20*e83e87a5Sjeremylt ierr = DMGetCoordinatesLocal(dm, &coordinates); CHKERRQ(ierr); 21*e83e87a5Sjeremylt ierr = VecGetLocalSize(coordinates, &Nv); CHKERRQ(ierr); 22*e83e87a5Sjeremylt ierr = VecGetBlockSize(coordinates, &dim); CHKERRQ(ierr); 23*e83e87a5Sjeremylt Nv /= dim; 24*e83e87a5Sjeremylt ierr = VecGetArray(coordinates, &coords); CHKERRQ(ierr); 25*e83e87a5Sjeremylt for (v = 0; v < Nv; ++v) { 26*e83e87a5Sjeremylt PetscReal r = 0.0; 27*e83e87a5Sjeremylt 28*e83e87a5Sjeremylt for (d = 0; d < dim; ++d) r += PetscSqr(PetscRealPart(coords[v*dim+d])); 29*e83e87a5Sjeremylt r = PetscSqrtReal(r); 30*e83e87a5Sjeremylt for (d = 0; d < dim; ++d) coords[v*dim+d] /= r; 31*e83e87a5Sjeremylt } 32*e83e87a5Sjeremylt ierr = VecRestoreArray(coordinates, &coords); CHKERRQ(ierr); 33*e83e87a5Sjeremylt PetscFunctionReturn(0); 34*e83e87a5Sjeremylt }; 35*e83e87a5Sjeremylt 36*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 37*e83e87a5Sjeremylt // PETSc FE Boilerplate 38*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 39*e83e87a5Sjeremylt PetscErrorCode PetscFECreateByDegree(DM dm, PetscInt dim, PetscInt Nc, 40*e83e87a5Sjeremylt PetscBool isSimplex, const char prefix[], 41*e83e87a5Sjeremylt PetscInt order, PetscFE *fem) { 42*e83e87a5Sjeremylt PetscQuadrature q, fq; 43*e83e87a5Sjeremylt DM K; 44*e83e87a5Sjeremylt PetscSpace P; 45*e83e87a5Sjeremylt PetscDualSpace Q; 46*e83e87a5Sjeremylt PetscInt quadPointsPerEdge; 47*e83e87a5Sjeremylt PetscBool tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 48*e83e87a5Sjeremylt PetscErrorCode ierr; 49*e83e87a5Sjeremylt 50*e83e87a5Sjeremylt PetscFunctionBeginUser; 51*e83e87a5Sjeremylt /* Create space */ 52*e83e87a5Sjeremylt ierr = PetscSpaceCreate(PetscObjectComm((PetscObject) dm), &P); CHKERRQ(ierr); 53*e83e87a5Sjeremylt ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix); CHKERRQ(ierr); 54*e83e87a5Sjeremylt ierr = PetscSpacePolynomialSetTensor(P, tensor); CHKERRQ(ierr); 55*e83e87a5Sjeremylt ierr = PetscSpaceSetFromOptions(P); CHKERRQ(ierr); 56*e83e87a5Sjeremylt ierr = PetscSpaceSetNumComponents(P, Nc); CHKERRQ(ierr); 57*e83e87a5Sjeremylt ierr = PetscSpaceSetNumVariables(P, dim); CHKERRQ(ierr); 58*e83e87a5Sjeremylt ierr = PetscSpaceSetDegree(P, order, order); CHKERRQ(ierr); 59*e83e87a5Sjeremylt ierr = PetscSpaceSetUp(P); CHKERRQ(ierr); 60*e83e87a5Sjeremylt ierr = PetscSpacePolynomialGetTensor(P, &tensor); CHKERRQ(ierr); 61*e83e87a5Sjeremylt /* Create dual space */ 62*e83e87a5Sjeremylt ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) dm), &Q); 63*e83e87a5Sjeremylt CHKERRQ(ierr); 64*e83e87a5Sjeremylt ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE); CHKERRQ(ierr); 65*e83e87a5Sjeremylt ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix); CHKERRQ(ierr); 66*e83e87a5Sjeremylt ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K); CHKERRQ(ierr); 67*e83e87a5Sjeremylt ierr = PetscDualSpaceSetDM(Q, K); CHKERRQ(ierr); 68*e83e87a5Sjeremylt ierr = DMDestroy(&K); CHKERRQ(ierr); 69*e83e87a5Sjeremylt ierr = PetscDualSpaceSetNumComponents(Q, Nc); CHKERRQ(ierr); 70*e83e87a5Sjeremylt ierr = PetscDualSpaceSetOrder(Q, order); CHKERRQ(ierr); 71*e83e87a5Sjeremylt ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor); CHKERRQ(ierr); 72*e83e87a5Sjeremylt ierr = PetscDualSpaceSetFromOptions(Q); CHKERRQ(ierr); 73*e83e87a5Sjeremylt ierr = PetscDualSpaceSetUp(Q); CHKERRQ(ierr); 74*e83e87a5Sjeremylt /* Create element */ 75*e83e87a5Sjeremylt ierr = PetscFECreate(PetscObjectComm((PetscObject) dm), fem); CHKERRQ(ierr); 76*e83e87a5Sjeremylt ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix); CHKERRQ(ierr); 77*e83e87a5Sjeremylt ierr = PetscFESetFromOptions(*fem); CHKERRQ(ierr); 78*e83e87a5Sjeremylt ierr = PetscFESetBasisSpace(*fem, P); CHKERRQ(ierr); 79*e83e87a5Sjeremylt ierr = PetscFESetDualSpace(*fem, Q); CHKERRQ(ierr); 80*e83e87a5Sjeremylt ierr = PetscFESetNumComponents(*fem, Nc); CHKERRQ(ierr); 81*e83e87a5Sjeremylt ierr = PetscFESetUp(*fem); CHKERRQ(ierr); 82*e83e87a5Sjeremylt ierr = PetscSpaceDestroy(&P); CHKERRQ(ierr); 83*e83e87a5Sjeremylt ierr = PetscDualSpaceDestroy(&Q); CHKERRQ(ierr); 84*e83e87a5Sjeremylt /* Create quadrature */ 85*e83e87a5Sjeremylt quadPointsPerEdge = PetscMax(order + 1,1); 86*e83e87a5Sjeremylt if (isSimplex) { 87*e83e87a5Sjeremylt ierr = PetscDTStroudConicalQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, 88*e83e87a5Sjeremylt &q); CHKERRQ(ierr); 89*e83e87a5Sjeremylt ierr = PetscDTStroudConicalQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, 90*e83e87a5Sjeremylt &fq); CHKERRQ(ierr); 91*e83e87a5Sjeremylt } else { 92*e83e87a5Sjeremylt ierr = PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, 93*e83e87a5Sjeremylt &q); CHKERRQ(ierr); 94*e83e87a5Sjeremylt ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, 95*e83e87a5Sjeremylt &fq); CHKERRQ(ierr); 96*e83e87a5Sjeremylt } 97*e83e87a5Sjeremylt ierr = PetscFESetQuadrature(*fem, q); CHKERRQ(ierr); 98*e83e87a5Sjeremylt ierr = PetscFESetFaceQuadrature(*fem, fq); CHKERRQ(ierr); 99*e83e87a5Sjeremylt ierr = PetscQuadratureDestroy(&q); CHKERRQ(ierr); 100*e83e87a5Sjeremylt ierr = PetscQuadratureDestroy(&fq); CHKERRQ(ierr); 101*e83e87a5Sjeremylt 102*e83e87a5Sjeremylt PetscFunctionReturn(0); 103*e83e87a5Sjeremylt }; 104*e83e87a5Sjeremylt 105*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 106*e83e87a5Sjeremylt // Create BC label 107*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 108*e83e87a5Sjeremylt static PetscErrorCode CreateBCLabel(DM dm, const char name[]) { 109*e83e87a5Sjeremylt int ierr; 110*e83e87a5Sjeremylt DMLabel label; 111*e83e87a5Sjeremylt 112*e83e87a5Sjeremylt PetscFunctionBeginUser; 113*e83e87a5Sjeremylt 114*e83e87a5Sjeremylt ierr = DMCreateLabel(dm, name); CHKERRQ(ierr); 115*e83e87a5Sjeremylt ierr = DMGetLabel(dm, name, &label); CHKERRQ(ierr); 116*e83e87a5Sjeremylt ierr = DMPlexMarkBoundaryFaces(dm, 1, label); CHKERRQ(ierr); 117*e83e87a5Sjeremylt ierr = DMPlexLabelComplete(dm, label); CHKERRQ(ierr); 118*e83e87a5Sjeremylt 119*e83e87a5Sjeremylt PetscFunctionReturn(0); 120*e83e87a5Sjeremylt }; 121*e83e87a5Sjeremylt 122*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 123*e83e87a5Sjeremylt // This function sets up a DM for a given degree 124*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 125*e83e87a5Sjeremylt PetscErrorCode SetupDMByDegree(DM dm, PetscInt degree, PetscInt ncompu, 126*e83e87a5Sjeremylt PetscInt dim, bool enforcebc, BCFunction bcsfunc) { 127*e83e87a5Sjeremylt PetscInt ierr, marker_ids[1] = {1}; 128*e83e87a5Sjeremylt PetscFE fe; 129*e83e87a5Sjeremylt 130*e83e87a5Sjeremylt PetscFunctionBeginUser; 131*e83e87a5Sjeremylt 132*e83e87a5Sjeremylt // Setup FE 133*e83e87a5Sjeremylt ierr = PetscFECreateByDegree(dm, dim, ncompu, PETSC_FALSE, NULL, degree, &fe); 134*e83e87a5Sjeremylt CHKERRQ(ierr); 135*e83e87a5Sjeremylt ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 136*e83e87a5Sjeremylt ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr); 137*e83e87a5Sjeremylt 138*e83e87a5Sjeremylt // Setup DM 139*e83e87a5Sjeremylt ierr = DMCreateDS(dm); CHKERRQ(ierr); 140*e83e87a5Sjeremylt if (enforcebc) { 141*e83e87a5Sjeremylt PetscBool hasLabel; 142*e83e87a5Sjeremylt DMHasLabel(dm, "marker", &hasLabel); 143*e83e87a5Sjeremylt if (!hasLabel) {CreateBCLabel(dm, "marker");} 144*e83e87a5Sjeremylt ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, 145*e83e87a5Sjeremylt (void(*)(void))bcsfunc, NULL, 1, marker_ids, NULL); 146*e83e87a5Sjeremylt CHKERRQ(ierr); 147*e83e87a5Sjeremylt } 148*e83e87a5Sjeremylt ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 149*e83e87a5Sjeremylt CHKERRQ(ierr); 150*e83e87a5Sjeremylt ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 151*e83e87a5Sjeremylt 152*e83e87a5Sjeremylt PetscFunctionReturn(0); 153*e83e87a5Sjeremylt }; 154*e83e87a5Sjeremylt 155*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 156*e83e87a5Sjeremylt // Utility function - essential BC dofs are encoded in closure indices as -(i+1) 157*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 158*e83e87a5Sjeremylt PetscInt Involute(PetscInt i) { 159*e83e87a5Sjeremylt return i >= 0 ? i : -(i + 1); 160*e83e87a5Sjeremylt }; 161*e83e87a5Sjeremylt 162*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 163*e83e87a5Sjeremylt // Get CEED restriction data from DMPlex 164*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 165*e83e87a5Sjeremylt PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P, 166*e83e87a5Sjeremylt CeedInt topodim, CeedInt height, DMLabel domainLabel, CeedInt value, 167*e83e87a5Sjeremylt CeedElemRestriction *Erestrict) { 168*e83e87a5Sjeremylt PetscSection section; 169*e83e87a5Sjeremylt PetscInt p, Nelem, Ndof, *erestrict, eoffset, nfields, dim, depth; 170*e83e87a5Sjeremylt DMLabel depthLabel; 171*e83e87a5Sjeremylt IS depthIS, iterIS; 172*e83e87a5Sjeremylt Vec Uloc; 173*e83e87a5Sjeremylt const PetscInt *iterIndices; 174*e83e87a5Sjeremylt PetscErrorCode ierr; 175*e83e87a5Sjeremylt 176*e83e87a5Sjeremylt PetscFunctionBeginUser; 177*e83e87a5Sjeremylt 178*e83e87a5Sjeremylt ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 179*e83e87a5Sjeremylt dim -= height; 180*e83e87a5Sjeremylt ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 181*e83e87a5Sjeremylt ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr); 182*e83e87a5Sjeremylt PetscInt ncomp[nfields], fieldoff[nfields+1]; 183*e83e87a5Sjeremylt fieldoff[0] = 0; 184*e83e87a5Sjeremylt for (PetscInt f = 0; f < nfields; f++) { 185*e83e87a5Sjeremylt ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr); 186*e83e87a5Sjeremylt fieldoff[f+1] = fieldoff[f] + ncomp[f]; 187*e83e87a5Sjeremylt } 188*e83e87a5Sjeremylt 189*e83e87a5Sjeremylt ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr); 190*e83e87a5Sjeremylt ierr = DMPlexGetDepthLabel(dm, &depthLabel); CHKERRQ(ierr); 191*e83e87a5Sjeremylt ierr = DMLabelGetStratumIS(depthLabel, depth - height, &depthIS); CHKERRQ(ierr); 192*e83e87a5Sjeremylt if (domainLabel) { 193*e83e87a5Sjeremylt IS domainIS; 194*e83e87a5Sjeremylt ierr = DMLabelGetStratumIS(domainLabel, value, &domainIS); CHKERRQ(ierr); 195*e83e87a5Sjeremylt if (domainIS) { // domainIS is non-empty 196*e83e87a5Sjeremylt ierr = ISIntersect(depthIS, domainIS, &iterIS); CHKERRQ(ierr); 197*e83e87a5Sjeremylt ierr = ISDestroy(&domainIS); CHKERRQ(ierr); 198*e83e87a5Sjeremylt } else { // domainIS is NULL (empty) 199*e83e87a5Sjeremylt iterIS = NULL; 200*e83e87a5Sjeremylt } 201*e83e87a5Sjeremylt ierr = ISDestroy(&depthIS); CHKERRQ(ierr); 202*e83e87a5Sjeremylt } else { 203*e83e87a5Sjeremylt iterIS = depthIS; 204*e83e87a5Sjeremylt } 205*e83e87a5Sjeremylt if (iterIS) { 206*e83e87a5Sjeremylt ierr = ISGetLocalSize(iterIS, &Nelem); CHKERRQ(ierr); 207*e83e87a5Sjeremylt ierr = ISGetIndices(iterIS, &iterIndices); CHKERRQ(ierr); 208*e83e87a5Sjeremylt } else { 209*e83e87a5Sjeremylt Nelem = 0; 210*e83e87a5Sjeremylt iterIndices = NULL; 211*e83e87a5Sjeremylt } 212*e83e87a5Sjeremylt ierr = PetscMalloc1(Nelem*PetscPowInt(P, topodim), &erestrict); CHKERRQ(ierr); 213*e83e87a5Sjeremylt for (p = 0, eoffset = 0; p < Nelem; p++) { 214*e83e87a5Sjeremylt PetscInt c = iterIndices[p]; 215*e83e87a5Sjeremylt PetscInt numindices, *indices, nnodes; 216*e83e87a5Sjeremylt ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, 217*e83e87a5Sjeremylt &numindices, &indices, NULL, NULL); 218*e83e87a5Sjeremylt CHKERRQ(ierr); 219*e83e87a5Sjeremylt bool flip = false; 220*e83e87a5Sjeremylt if (height > 0) { 221*e83e87a5Sjeremylt PetscInt numCells, numFaces, start = -1; 222*e83e87a5Sjeremylt const PetscInt *orients, *faces, *cells; 223*e83e87a5Sjeremylt ierr = DMPlexGetSupport(dm, c, &cells); CHKERRQ(ierr); 224*e83e87a5Sjeremylt ierr = DMPlexGetSupportSize(dm, c, &numCells); CHKERRQ(ierr); 225*e83e87a5Sjeremylt if (numCells != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 226*e83e87a5Sjeremylt "Expected one cell in support of exterior face, but got %D cells", 227*e83e87a5Sjeremylt numCells); 228*e83e87a5Sjeremylt ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr); 229*e83e87a5Sjeremylt ierr = DMPlexGetConeSize(dm, cells[0], &numFaces); CHKERRQ(ierr); 230*e83e87a5Sjeremylt for (PetscInt i=0; i<numFaces; i++) {if (faces[i] == c) start = i;} 231*e83e87a5Sjeremylt if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, 232*e83e87a5Sjeremylt "Could not find face %D in cone of its support", 233*e83e87a5Sjeremylt c); 234*e83e87a5Sjeremylt ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr); 235*e83e87a5Sjeremylt if (orients[start] < 0) flip = true; 236*e83e87a5Sjeremylt } 237*e83e87a5Sjeremylt if (numindices % fieldoff[nfields]) SETERRQ1(PETSC_COMM_SELF, 238*e83e87a5Sjeremylt PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D", 239*e83e87a5Sjeremylt c); 240*e83e87a5Sjeremylt nnodes = numindices / fieldoff[nfields]; 241*e83e87a5Sjeremylt for (PetscInt i = 0; i < nnodes; i++) { 242*e83e87a5Sjeremylt PetscInt ii = i; 243*e83e87a5Sjeremylt if (flip) { 244*e83e87a5Sjeremylt if (P == nnodes) ii = nnodes - 1 - i; 245*e83e87a5Sjeremylt else if (P*P == nnodes) { 246*e83e87a5Sjeremylt PetscInt row = i / P, col = i % P; 247*e83e87a5Sjeremylt ii = row + col * P; 248*e83e87a5Sjeremylt } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, 249*e83e87a5Sjeremylt "No support for flipping point with %D nodes != P (%D) or P^2", 250*e83e87a5Sjeremylt nnodes, P); 251*e83e87a5Sjeremylt } 252*e83e87a5Sjeremylt // Check that indices are blocked by node and thus can be coalesced as a single field with 253*e83e87a5Sjeremylt // fieldoff[nfields] = sum(ncomp) components. 254*e83e87a5Sjeremylt for (PetscInt f = 0; f < nfields; f++) { 255*e83e87a5Sjeremylt for (PetscInt j = 0; j < ncomp[f]; j++) { 256*e83e87a5Sjeremylt if (Involute(indices[fieldoff[f]*nnodes + ii*ncomp[f] + j]) 257*e83e87a5Sjeremylt != Involute(indices[ii*ncomp[0]]) + fieldoff[f] + j) 258*e83e87a5Sjeremylt SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 259*e83e87a5Sjeremylt "Cell %D closure indices not interlaced for node %D field %D component %D", 260*e83e87a5Sjeremylt c, ii, f, j); 261*e83e87a5Sjeremylt } 262*e83e87a5Sjeremylt } 263*e83e87a5Sjeremylt // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode. 264*e83e87a5Sjeremylt PetscInt loc = Involute(indices[ii*ncomp[0]]); 265*e83e87a5Sjeremylt erestrict[eoffset++] = loc; 266*e83e87a5Sjeremylt } 267*e83e87a5Sjeremylt ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, 268*e83e87a5Sjeremylt &numindices, &indices, NULL, NULL); 269*e83e87a5Sjeremylt CHKERRQ(ierr); 270*e83e87a5Sjeremylt } 271*e83e87a5Sjeremylt if (eoffset != Nelem*PetscPowInt(P, topodim)) 272*e83e87a5Sjeremylt SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB, 273*e83e87a5Sjeremylt "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem, 274*e83e87a5Sjeremylt PetscPowInt(P, topodim),eoffset); 275*e83e87a5Sjeremylt if (iterIS) { 276*e83e87a5Sjeremylt ierr = ISRestoreIndices(iterIS, &iterIndices); CHKERRQ(ierr); 277*e83e87a5Sjeremylt } 278*e83e87a5Sjeremylt ierr = ISDestroy(&iterIS); CHKERRQ(ierr); 279*e83e87a5Sjeremylt 280*e83e87a5Sjeremylt ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr); 281*e83e87a5Sjeremylt ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr); 282*e83e87a5Sjeremylt ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr); 283*e83e87a5Sjeremylt CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, topodim), 284*e83e87a5Sjeremylt fieldoff[nfields], 285*e83e87a5Sjeremylt 1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict, 286*e83e87a5Sjeremylt Erestrict); 287*e83e87a5Sjeremylt ierr = PetscFree(erestrict); CHKERRQ(ierr); 288*e83e87a5Sjeremylt PetscFunctionReturn(0); 289*e83e87a5Sjeremylt }; 290*e83e87a5Sjeremylt 291*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 292