1*e83e87a5Sjeremylt #include "../include/libceedsetup.h" 2*e83e87a5Sjeremylt #include "../include/petscutils.h" 3*e83e87a5Sjeremylt #include <stdio.h> 4*e83e87a5Sjeremylt 5*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 6*e83e87a5Sjeremylt // Destroy libCEED operator objects 7*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 8*e83e87a5Sjeremylt PetscErrorCode CeedDataDestroy(CeedInt i, CeedData data) { 9*e83e87a5Sjeremylt int ierr; 10*e83e87a5Sjeremylt 11*e83e87a5Sjeremylt CeedVectorDestroy(&data->qdata); 12*e83e87a5Sjeremylt CeedVectorDestroy(&data->Xceed); 13*e83e87a5Sjeremylt CeedVectorDestroy(&data->Yceed); 14*e83e87a5Sjeremylt CeedBasisDestroy(&data->basisx); 15*e83e87a5Sjeremylt CeedBasisDestroy(&data->basisu); 16*e83e87a5Sjeremylt CeedElemRestrictionDestroy(&data->Erestrictu); 17*e83e87a5Sjeremylt CeedElemRestrictionDestroy(&data->Erestrictx); 18*e83e87a5Sjeremylt CeedElemRestrictionDestroy(&data->Erestrictui); 19*e83e87a5Sjeremylt CeedElemRestrictionDestroy(&data->Erestrictqdi); 20*e83e87a5Sjeremylt CeedQFunctionDestroy(&data->qfApply); 21*e83e87a5Sjeremylt CeedOperatorDestroy(&data->opApply); 22*e83e87a5Sjeremylt if (i > 0) { 23*e83e87a5Sjeremylt CeedOperatorDestroy(&data->opProlong); 24*e83e87a5Sjeremylt CeedBasisDestroy(&data->basisctof); 25*e83e87a5Sjeremylt CeedOperatorDestroy(&data->opRestrict); 26*e83e87a5Sjeremylt } 27*e83e87a5Sjeremylt ierr = PetscFree(data); CHKERRQ(ierr); 28*e83e87a5Sjeremylt 29*e83e87a5Sjeremylt PetscFunctionReturn(0); 30*e83e87a5Sjeremylt }; 31*e83e87a5Sjeremylt 32*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 33*e83e87a5Sjeremylt // Set up libCEED for a given degree 34*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 35*e83e87a5Sjeremylt PetscErrorCode SetupLibceedByDegree(DM dm, Ceed ceed, CeedInt degree, 36*e83e87a5Sjeremylt CeedInt topodim, CeedInt qextra, 37*e83e87a5Sjeremylt PetscInt ncompx, PetscInt ncompu, 38*e83e87a5Sjeremylt PetscInt gsize, PetscInt xlsize, 39*e83e87a5Sjeremylt bpData bpData, CeedData data, 40*e83e87a5Sjeremylt PetscBool setup_rhs, CeedVector rhsceed, 41*e83e87a5Sjeremylt CeedVector *target) { 42*e83e87a5Sjeremylt int ierr; 43*e83e87a5Sjeremylt DM dmcoord; 44*e83e87a5Sjeremylt Vec coords; 45*e83e87a5Sjeremylt const PetscScalar *coordArray; 46*e83e87a5Sjeremylt CeedBasis basisx, basisu; 47*e83e87a5Sjeremylt CeedElemRestriction Erestrictx, Erestrictu, Erestrictui, Erestrictqdi; 48*e83e87a5Sjeremylt CeedQFunction qfSetupGeo, qfApply; 49*e83e87a5Sjeremylt CeedOperator opSetupGeo, opApply; 50*e83e87a5Sjeremylt CeedVector xcoord, qdata, Xceed, Yceed; 51*e83e87a5Sjeremylt CeedInt P, Q, nqpts, cStart, cEnd, nelem, qdatasize = bpData.qdatasize; 52*e83e87a5Sjeremylt CeedScalar R = 1, // radius of the sphere 53*e83e87a5Sjeremylt l = 1.0/PetscSqrtReal(3.0); // half edge of the inscribed cube 54*e83e87a5Sjeremylt 55*e83e87a5Sjeremylt // CEED bases 56*e83e87a5Sjeremylt P = degree + 1; 57*e83e87a5Sjeremylt Q = P + qextra; 58*e83e87a5Sjeremylt CeedBasisCreateTensorH1Lagrange(ceed, topodim, ncompu, P, Q, bpData.qmode, 59*e83e87a5Sjeremylt &basisu); 60*e83e87a5Sjeremylt CeedBasisCreateTensorH1Lagrange(ceed, topodim, ncompx, 2, Q, bpData.qmode, 61*e83e87a5Sjeremylt &basisx); 62*e83e87a5Sjeremylt CeedBasisGetNumQuadraturePoints(basisu, &nqpts); 63*e83e87a5Sjeremylt 64*e83e87a5Sjeremylt // CEED restrictions 65*e83e87a5Sjeremylt ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr); 66*e83e87a5Sjeremylt ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL); 67*e83e87a5Sjeremylt CHKERRQ(ierr); 68*e83e87a5Sjeremylt ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, topodim, 0, 0, 0, 69*e83e87a5Sjeremylt &Erestrictx); 70*e83e87a5Sjeremylt CHKERRQ(ierr); 71*e83e87a5Sjeremylt ierr = CreateRestrictionFromPlex(ceed, dm, P, topodim, 0, 0, 0, &Erestrictu); 72*e83e87a5Sjeremylt CHKERRQ(ierr); 73*e83e87a5Sjeremylt 74*e83e87a5Sjeremylt ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr); 75*e83e87a5Sjeremylt nelem = cEnd - cStart; 76*e83e87a5Sjeremylt 77*e83e87a5Sjeremylt CeedElemRestrictionCreateStrided(ceed, nelem, nqpts, ncompu, ncompu*nelem*nqpts, 78*e83e87a5Sjeremylt CEED_STRIDES_BACKEND, &Erestrictui); 79*e83e87a5Sjeremylt CeedElemRestrictionCreateStrided(ceed, nelem, nqpts, qdatasize, 80*e83e87a5Sjeremylt qdatasize*nelem*nqpts, 81*e83e87a5Sjeremylt CEED_STRIDES_BACKEND, &Erestrictqdi); 82*e83e87a5Sjeremylt 83*e83e87a5Sjeremylt // Element coordinates 84*e83e87a5Sjeremylt ierr = DMGetCoordinatesLocal(dm, &coords); CHKERRQ(ierr); 85*e83e87a5Sjeremylt ierr = VecGetArrayRead(coords, &coordArray); CHKERRQ(ierr); 86*e83e87a5Sjeremylt 87*e83e87a5Sjeremylt CeedElemRestrictionCreateVector(Erestrictx, &xcoord, NULL); 88*e83e87a5Sjeremylt CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_COPY_VALUES, 89*e83e87a5Sjeremylt (PetscScalar *)coordArray); 90*e83e87a5Sjeremylt ierr = VecRestoreArrayRead(coords, &coordArray); 91*e83e87a5Sjeremylt 92*e83e87a5Sjeremylt // Create the persistent vectors that will be needed in setup and apply 93*e83e87a5Sjeremylt CeedVectorCreate(ceed, qdatasize*nelem*nqpts, &qdata); 94*e83e87a5Sjeremylt CeedVectorCreate(ceed, xlsize, &Xceed); 95*e83e87a5Sjeremylt CeedVectorCreate(ceed, xlsize, &Yceed); 96*e83e87a5Sjeremylt 97*e83e87a5Sjeremylt // Create the QFunction that builds the context data 98*e83e87a5Sjeremylt CeedQFunctionCreateInterior(ceed, 1, bpData.setupgeo, bpData.setupgeofname, 99*e83e87a5Sjeremylt &qfSetupGeo); 100*e83e87a5Sjeremylt CeedQFunctionAddInput(qfSetupGeo, "x", ncompx, CEED_EVAL_INTERP); 101*e83e87a5Sjeremylt CeedQFunctionAddInput(qfSetupGeo, "dx", ncompx*topodim, CEED_EVAL_GRAD); 102*e83e87a5Sjeremylt CeedQFunctionAddInput(qfSetupGeo, "weight", 1, CEED_EVAL_WEIGHT); 103*e83e87a5Sjeremylt CeedQFunctionAddOutput(qfSetupGeo, "qdata", qdatasize, CEED_EVAL_NONE); 104*e83e87a5Sjeremylt 105*e83e87a5Sjeremylt // Create the operator that builds the quadrature data 106*e83e87a5Sjeremylt CeedOperatorCreate(ceed, qfSetupGeo, NULL, NULL, &opSetupGeo); 107*e83e87a5Sjeremylt CeedOperatorSetField(opSetupGeo, "x", Erestrictx, basisx, CEED_VECTOR_ACTIVE); 108*e83e87a5Sjeremylt CeedOperatorSetField(opSetupGeo, "dx", Erestrictx, basisx, CEED_VECTOR_ACTIVE); 109*e83e87a5Sjeremylt CeedOperatorSetField(opSetupGeo, "weight", CEED_ELEMRESTRICTION_NONE, basisx, 110*e83e87a5Sjeremylt CEED_VECTOR_NONE); 111*e83e87a5Sjeremylt CeedOperatorSetField(opSetupGeo, "qdata", Erestrictqdi, 112*e83e87a5Sjeremylt CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 113*e83e87a5Sjeremylt 114*e83e87a5Sjeremylt // Setup qdata 115*e83e87a5Sjeremylt CeedOperatorApply(opSetupGeo, xcoord, qdata, CEED_REQUEST_IMMEDIATE); 116*e83e87a5Sjeremylt 117*e83e87a5Sjeremylt // Set up PDE operator 118*e83e87a5Sjeremylt CeedInt inscale = bpData.inmode == CEED_EVAL_GRAD ? topodim : 1; 119*e83e87a5Sjeremylt CeedInt outscale = bpData.outmode == CEED_EVAL_GRAD ? topodim : 1; 120*e83e87a5Sjeremylt CeedQFunctionCreateInterior(ceed, 1, bpData.apply, bpData.applyfname, &qfApply); 121*e83e87a5Sjeremylt CeedQFunctionAddInput(qfApply, "u", ncompu*inscale, bpData.inmode); 122*e83e87a5Sjeremylt CeedQFunctionAddInput(qfApply, "qdata", qdatasize, CEED_EVAL_NONE); 123*e83e87a5Sjeremylt CeedQFunctionAddOutput(qfApply, "v", ncompu*outscale, bpData.outmode); 124*e83e87a5Sjeremylt 125*e83e87a5Sjeremylt // Create the mass or diff operator 126*e83e87a5Sjeremylt CeedOperatorCreate(ceed, qfApply, NULL, NULL, &opApply); 127*e83e87a5Sjeremylt CeedOperatorSetField(opApply, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE); 128*e83e87a5Sjeremylt CeedOperatorSetField(opApply, "qdata", Erestrictqdi, CEED_BASIS_COLLOCATED, 129*e83e87a5Sjeremylt qdata); 130*e83e87a5Sjeremylt CeedOperatorSetField(opApply, "v", Erestrictu, basisu, CEED_VECTOR_ACTIVE); 131*e83e87a5Sjeremylt 132*e83e87a5Sjeremylt // Set up RHS if needed 133*e83e87a5Sjeremylt if (setup_rhs) { 134*e83e87a5Sjeremylt CeedQFunction qfSetupRHS; 135*e83e87a5Sjeremylt CeedOperator opSetupRHS; 136*e83e87a5Sjeremylt CeedVectorCreate(ceed, nelem*nqpts*ncompu, target); 137*e83e87a5Sjeremylt 138*e83e87a5Sjeremylt // Create the q-function that sets up the RHS and true solution 139*e83e87a5Sjeremylt CeedQFunctionCreateInterior(ceed, 1, bpData.setuprhs, bpData.setuprhsfname, 140*e83e87a5Sjeremylt &qfSetupRHS); 141*e83e87a5Sjeremylt CeedQFunctionAddInput(qfSetupRHS, "x", ncompx, CEED_EVAL_INTERP); 142*e83e87a5Sjeremylt CeedQFunctionAddInput(qfSetupRHS, "qdata", qdatasize, CEED_EVAL_NONE); 143*e83e87a5Sjeremylt CeedQFunctionAddOutput(qfSetupRHS, "true_soln", ncompu, CEED_EVAL_NONE); 144*e83e87a5Sjeremylt CeedQFunctionAddOutput(qfSetupRHS, "rhs", ncompu, CEED_EVAL_INTERP); 145*e83e87a5Sjeremylt 146*e83e87a5Sjeremylt // Create the operator that builds the RHS and true solution 147*e83e87a5Sjeremylt CeedOperatorCreate(ceed, qfSetupRHS, NULL, NULL, &opSetupRHS); 148*e83e87a5Sjeremylt CeedOperatorSetField(opSetupRHS, "x", Erestrictx, basisx, CEED_VECTOR_ACTIVE); 149*e83e87a5Sjeremylt CeedOperatorSetField(opSetupRHS, "qdata", Erestrictqdi, CEED_BASIS_COLLOCATED, 150*e83e87a5Sjeremylt qdata); 151*e83e87a5Sjeremylt CeedOperatorSetField(opSetupRHS, "true_soln", Erestrictui, 152*e83e87a5Sjeremylt CEED_BASIS_COLLOCATED, *target); 153*e83e87a5Sjeremylt CeedOperatorSetField(opSetupRHS, "rhs", Erestrictu, basisu, 154*e83e87a5Sjeremylt CEED_VECTOR_ACTIVE); 155*e83e87a5Sjeremylt 156*e83e87a5Sjeremylt // Set up the libCEED context 157*e83e87a5Sjeremylt CeedQFunctionContext rhsSetupCtx; 158*e83e87a5Sjeremylt CeedQFunctionContextCreate(ceed, &rhsSetupCtx); 159*e83e87a5Sjeremylt CeedScalar rhsSetupData[2] = {R, l}; 160*e83e87a5Sjeremylt CeedQFunctionContextSetData(rhsSetupCtx, CEED_MEM_HOST, CEED_COPY_VALUES, 161*e83e87a5Sjeremylt sizeof rhsSetupData, &rhsSetupData); 162*e83e87a5Sjeremylt CeedQFunctionSetContext(qfSetupRHS, rhsSetupCtx); 163*e83e87a5Sjeremylt CeedQFunctionContextDestroy(&rhsSetupCtx); 164*e83e87a5Sjeremylt 165*e83e87a5Sjeremylt // Setup RHS and target 166*e83e87a5Sjeremylt CeedOperatorApply(opSetupRHS, xcoord, rhsceed, CEED_REQUEST_IMMEDIATE); 167*e83e87a5Sjeremylt 168*e83e87a5Sjeremylt // Cleanup 169*e83e87a5Sjeremylt CeedQFunctionDestroy(&qfSetupRHS); 170*e83e87a5Sjeremylt CeedOperatorDestroy(&opSetupRHS); 171*e83e87a5Sjeremylt } 172*e83e87a5Sjeremylt 173*e83e87a5Sjeremylt // Cleanup 174*e83e87a5Sjeremylt CeedQFunctionDestroy(&qfSetupGeo); 175*e83e87a5Sjeremylt CeedOperatorDestroy(&opSetupGeo); 176*e83e87a5Sjeremylt CeedVectorDestroy(&xcoord); 177*e83e87a5Sjeremylt 178*e83e87a5Sjeremylt // Save libCEED data required for level 179*e83e87a5Sjeremylt data->basisx = basisx; data->basisu = basisu; 180*e83e87a5Sjeremylt data->Erestrictx = Erestrictx; 181*e83e87a5Sjeremylt data->Erestrictu = Erestrictu; 182*e83e87a5Sjeremylt data->Erestrictui = Erestrictui; 183*e83e87a5Sjeremylt data->Erestrictqdi = Erestrictqdi; 184*e83e87a5Sjeremylt data->qfApply = qfApply; 185*e83e87a5Sjeremylt data->opApply = opApply; 186*e83e87a5Sjeremylt data->qdata = qdata; 187*e83e87a5Sjeremylt data->Xceed = Xceed; 188*e83e87a5Sjeremylt data->Yceed = Yceed; 189*e83e87a5Sjeremylt 190*e83e87a5Sjeremylt PetscFunctionReturn(0); 191*e83e87a5Sjeremylt }; 192*e83e87a5Sjeremylt 193*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 194*e83e87a5Sjeremylt // Setup libCEED level transfer operator objects 195*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 196*e83e87a5Sjeremylt PetscErrorCode CeedLevelTransferSetup(Ceed ceed, CeedInt numlevels, 197*e83e87a5Sjeremylt CeedInt ncompu, 198*e83e87a5Sjeremylt CeedData *data, CeedInt *leveldegrees, 199*e83e87a5Sjeremylt CeedQFunction qfrestrict, CeedQFunction qfprolong) { 200*e83e87a5Sjeremylt // Return early if numlevels=1 201*e83e87a5Sjeremylt if (numlevels==1) 202*e83e87a5Sjeremylt PetscFunctionReturn(0); 203*e83e87a5Sjeremylt 204*e83e87a5Sjeremylt // Set up each level 205*e83e87a5Sjeremylt for (CeedInt i=1; i<numlevels; i++) { 206*e83e87a5Sjeremylt // P coarse and P fine 207*e83e87a5Sjeremylt CeedInt Pc = leveldegrees[i-1] + 1; 208*e83e87a5Sjeremylt CeedInt Pf = leveldegrees[i] + 1; 209*e83e87a5Sjeremylt 210*e83e87a5Sjeremylt // Restriction - Fine to corse 211*e83e87a5Sjeremylt CeedBasis basisctof; 212*e83e87a5Sjeremylt CeedOperator opRestrict; 213*e83e87a5Sjeremylt 214*e83e87a5Sjeremylt // Basis 215*e83e87a5Sjeremylt CeedBasisCreateTensorH1Lagrange(ceed, 3, ncompu, Pc, Pf, 216*e83e87a5Sjeremylt CEED_GAUSS_LOBATTO, &basisctof); 217*e83e87a5Sjeremylt 218*e83e87a5Sjeremylt // Create the restriction operator 219*e83e87a5Sjeremylt CeedOperatorCreate(ceed, qfrestrict, CEED_QFUNCTION_NONE, 220*e83e87a5Sjeremylt CEED_QFUNCTION_NONE, &opRestrict); 221*e83e87a5Sjeremylt CeedOperatorSetField(opRestrict, "input", data[i]->Erestrictu, 222*e83e87a5Sjeremylt CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 223*e83e87a5Sjeremylt CeedOperatorSetField(opRestrict, "output", data[i-1]->Erestrictu, 224*e83e87a5Sjeremylt basisctof, CEED_VECTOR_ACTIVE); 225*e83e87a5Sjeremylt 226*e83e87a5Sjeremylt // Save libCEED data required for level 227*e83e87a5Sjeremylt data[i]->basisctof = basisctof; 228*e83e87a5Sjeremylt data[i]->opRestrict = opRestrict; 229*e83e87a5Sjeremylt 230*e83e87a5Sjeremylt // Interpolation - Corse to fine 231*e83e87a5Sjeremylt CeedOperator opProlong; 232*e83e87a5Sjeremylt 233*e83e87a5Sjeremylt // Create the prolongation operator 234*e83e87a5Sjeremylt CeedOperatorCreate(ceed, qfprolong, CEED_QFUNCTION_NONE, 235*e83e87a5Sjeremylt CEED_QFUNCTION_NONE, &opProlong); 236*e83e87a5Sjeremylt CeedOperatorSetField(opProlong, "input", data[i-1]->Erestrictu, 237*e83e87a5Sjeremylt basisctof, CEED_VECTOR_ACTIVE); 238*e83e87a5Sjeremylt CeedOperatorSetField(opProlong, "output", data[i]->Erestrictu, 239*e83e87a5Sjeremylt CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 240*e83e87a5Sjeremylt 241*e83e87a5Sjeremylt // Save libCEED data required for level 242*e83e87a5Sjeremylt data[i]->opProlong = opProlong; 243*e83e87a5Sjeremylt } 244*e83e87a5Sjeremylt 245*e83e87a5Sjeremylt PetscFunctionReturn(0); 246*e83e87a5Sjeremylt }; 247*e83e87a5Sjeremylt 248*e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 249