1f918ec44SMatthew G. Knepley #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 2*d2b2dc1eSMatthew G. Knepley #include <petscdmceed.h> 3f918ec44SMatthew G. Knepley 4f918ec44SMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED 5*d2b2dc1eSMatthew G. Knepley #include <petsc/private/dmpleximpl.h> 6*d2b2dc1eSMatthew G. Knepley #include <petscdmplexceed.h> 7*d2b2dc1eSMatthew G. Knepley #include <petscfeceed.h> 8f918ec44SMatthew G. Knepley 9f918ec44SMatthew G. Knepley /*@C 1020f4b53cSBarry Smith DMGetCeed - Get the LibCEED context associated with this `DM` 11f918ec44SMatthew G. Knepley 1220f4b53cSBarry Smith Not Collective 13f918ec44SMatthew G. Knepley 14f918ec44SMatthew G. Knepley Input Parameter: 1520f4b53cSBarry Smith . DM - The `DM` 16f918ec44SMatthew G. Knepley 17f918ec44SMatthew G. Knepley Output Parameter: 18f918ec44SMatthew G. Knepley . ceed - The LibCEED context 19f918ec44SMatthew G. Knepley 20f918ec44SMatthew G. Knepley Level: intermediate 21f918ec44SMatthew G. Knepley 2220f4b53cSBarry Smith .seealso: `DM`, `DMCreate()` 23f918ec44SMatthew G. Knepley @*/ 24d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCeed(DM dm, Ceed *ceed) 25d71ae5a4SJacob Faibussowitsch { 26f918ec44SMatthew G. Knepley PetscFunctionBegin; 27f918ec44SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 284f572ea9SToby Isaac PetscAssertPointer(ceed, 2); 29f918ec44SMatthew G. Knepley if (!dm->ceed) { 30f918ec44SMatthew G. Knepley char ceedresource[PETSC_MAX_PATH_LEN]; /* libCEED resource specifier */ 31f918ec44SMatthew G. Knepley const char *prefix; 32f918ec44SMatthew G. Knepley 33c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(ceedresource, "/cpu/self", sizeof(ceedresource))); 349566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, prefix, "-dm_ceed", ceedresource, sizeof(ceedresource), NULL)); 369566063dSJacob Faibussowitsch PetscCallCEED(CeedInit(ceedresource, &dm->ceed)); 37f918ec44SMatthew G. Knepley } 38f918ec44SMatthew G. Knepley *ceed = dm->ceed; 393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40f918ec44SMatthew G. Knepley } 41f918ec44SMatthew G. Knepley 42*d2b2dc1eSMatthew G. Knepley static CeedMemType PetscMemType2Ceed(PetscMemType mem_type) 43*d2b2dc1eSMatthew G. Knepley { 447ce467fcSMatthew G. Knepley return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; 457ce467fcSMatthew G. Knepley } 467ce467fcSMatthew G. Knepley 477ce467fcSMatthew G. Knepley PetscErrorCode VecGetCeedVector(Vec X, Ceed ceed, CeedVector *cx) 487ce467fcSMatthew G. Knepley { 497ce467fcSMatthew G. Knepley PetscMemType memtype; 507ce467fcSMatthew G. Knepley PetscScalar *x; 517ce467fcSMatthew G. Knepley PetscInt n; 527ce467fcSMatthew G. Knepley 537ce467fcSMatthew G. Knepley PetscFunctionBegin; 547ce467fcSMatthew G. Knepley PetscCall(VecGetLocalSize(X, &n)); 557ce467fcSMatthew G. Knepley PetscCall(VecGetArrayAndMemType(X, &x, &memtype)); 567ce467fcSMatthew G. Knepley PetscCallCEED(CeedVectorCreate(ceed, n, cx)); 577ce467fcSMatthew G. Knepley PetscCallCEED(CeedVectorSetArray(*cx, PetscMemType2Ceed(memtype), CEED_USE_POINTER, x)); 58*d2b2dc1eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 597ce467fcSMatthew G. Knepley } 607ce467fcSMatthew G. Knepley 617ce467fcSMatthew G. Knepley PetscErrorCode VecRestoreCeedVector(Vec X, CeedVector *cx) 627ce467fcSMatthew G. Knepley { 637ce467fcSMatthew G. Knepley PetscFunctionBegin; 647ce467fcSMatthew G. Knepley PetscCall(VecRestoreArrayAndMemType(X, NULL)); 657ce467fcSMatthew G. Knepley PetscCallCEED(CeedVectorDestroy(cx)); 66*d2b2dc1eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 677ce467fcSMatthew G. Knepley } 687ce467fcSMatthew G. Knepley 697ce467fcSMatthew G. Knepley PetscErrorCode VecGetCeedVectorRead(Vec X, Ceed ceed, CeedVector *cx) 707ce467fcSMatthew G. Knepley { 717ce467fcSMatthew G. Knepley PetscMemType memtype; 727ce467fcSMatthew G. Knepley const PetscScalar *x; 737ce467fcSMatthew G. Knepley PetscInt n; 747ce467fcSMatthew G. Knepley PetscFunctionBegin; 757ce467fcSMatthew G. Knepley PetscCall(VecGetLocalSize(X, &n)); 767ce467fcSMatthew G. Knepley PetscCall(VecGetArrayReadAndMemType(X, &x, &memtype)); 777ce467fcSMatthew G. Knepley PetscCallCEED(CeedVectorCreate(ceed, n, cx)); 787ce467fcSMatthew G. Knepley PetscCallCEED(CeedVectorSetArray(*cx, PetscMemType2Ceed(memtype), CEED_USE_POINTER, (PetscScalar *)x)); 79*d2b2dc1eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 807ce467fcSMatthew G. Knepley } 817ce467fcSMatthew G. Knepley 827ce467fcSMatthew G. Knepley PetscErrorCode VecRestoreCeedVectorRead(Vec X, CeedVector *cx) 837ce467fcSMatthew G. Knepley { 847ce467fcSMatthew G. Knepley PetscFunctionBegin; 857ce467fcSMatthew G. Knepley PetscCall(VecRestoreArrayReadAndMemType(X, NULL)); 867ce467fcSMatthew G. Knepley PetscCallCEED(CeedVectorDestroy(cx)); 87*d2b2dc1eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 88*d2b2dc1eSMatthew G. Knepley } 89*d2b2dc1eSMatthew G. Knepley 90*d2b2dc1eSMatthew G. Knepley CEED_QFUNCTION(Geometry2D)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) 91*d2b2dc1eSMatthew G. Knepley { 92*d2b2dc1eSMatthew G. Knepley const CeedScalar *x = in[0], *Jac = in[1], *w = in[2]; 93*d2b2dc1eSMatthew G. Knepley CeedScalar *qdata = out[0]; 94*d2b2dc1eSMatthew G. Knepley 95*d2b2dc1eSMatthew G. Knepley CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) 96*d2b2dc1eSMatthew G. Knepley { 97*d2b2dc1eSMatthew G. Knepley const CeedScalar J[2][2] = { 98*d2b2dc1eSMatthew G. Knepley {Jac[i + Q * 0], Jac[i + Q * 2]}, 99*d2b2dc1eSMatthew G. Knepley {Jac[i + Q * 1], Jac[i + Q * 3]} 100*d2b2dc1eSMatthew G. Knepley }; 101*d2b2dc1eSMatthew G. Knepley const CeedScalar det = J[0][0] * J[1][1] - J[0][1] * J[1][0]; 102*d2b2dc1eSMatthew G. Knepley 103*d2b2dc1eSMatthew G. Knepley qdata[i + Q * 0] = det * w[i]; 104*d2b2dc1eSMatthew G. Knepley qdata[i + Q * 1] = x[i + Q * 0]; 105*d2b2dc1eSMatthew G. Knepley qdata[i + Q * 2] = x[i + Q * 1]; 106*d2b2dc1eSMatthew G. Knepley qdata[i + Q * 3] = J[1][1] / det; 107*d2b2dc1eSMatthew G. Knepley qdata[i + Q * 4] = -J[1][0] / det; 108*d2b2dc1eSMatthew G. Knepley qdata[i + Q * 5] = -J[0][1] / det; 109*d2b2dc1eSMatthew G. Knepley qdata[i + Q * 6] = J[0][0] / det; 110*d2b2dc1eSMatthew G. Knepley } 111*d2b2dc1eSMatthew G. Knepley return CEED_ERROR_SUCCESS; 112*d2b2dc1eSMatthew G. Knepley } 113*d2b2dc1eSMatthew G. Knepley 114*d2b2dc1eSMatthew G. Knepley CEED_QFUNCTION(Geometry3D)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) 115*d2b2dc1eSMatthew G. Knepley { 116*d2b2dc1eSMatthew G. Knepley const CeedScalar *Jac = in[1], *w = in[2]; 117*d2b2dc1eSMatthew G. Knepley CeedScalar *qdata = out[0]; 118*d2b2dc1eSMatthew G. Knepley 119*d2b2dc1eSMatthew G. Knepley CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) 120*d2b2dc1eSMatthew G. Knepley { 121*d2b2dc1eSMatthew G. Knepley const CeedScalar J[3][3] = { 122*d2b2dc1eSMatthew G. Knepley {Jac[i + Q * 0], Jac[i + Q * 3], Jac[i + Q * 6]}, 123*d2b2dc1eSMatthew G. Knepley {Jac[i + Q * 1], Jac[i + Q * 4], Jac[i + Q * 7]}, 124*d2b2dc1eSMatthew G. Knepley {Jac[i + Q * 2], Jac[i + Q * 5], Jac[i + Q * 8]} 125*d2b2dc1eSMatthew G. Knepley }; 126*d2b2dc1eSMatthew G. Knepley const CeedScalar det = J[0][0] * (J[1][1] * J[2][2] - J[1][2] * J[2][1]) + J[0][1] * (J[1][2] * J[2][0] - J[1][0] * J[2][2]) + J[0][2] * (J[1][0] * J[2][1] - J[1][1] * J[2][0]); 127*d2b2dc1eSMatthew G. Knepley 128*d2b2dc1eSMatthew G. Knepley qdata[i + Q * 0] = det * w[i]; /* det J * weight */ 129*d2b2dc1eSMatthew G. Knepley } 130*d2b2dc1eSMatthew G. Knepley return CEED_ERROR_SUCCESS; 131*d2b2dc1eSMatthew G. Knepley } 132*d2b2dc1eSMatthew G. Knepley 133*d2b2dc1eSMatthew G. Knepley static PetscErrorCode DMCeedCreateGeometry(DM dm, IS cellIS, PetscInt *Nqdata, CeedElemRestriction *erq, CeedVector *qd, DMCeed *soldata) 134*d2b2dc1eSMatthew G. Knepley { 135*d2b2dc1eSMatthew G. Knepley Ceed ceed; 136*d2b2dc1eSMatthew G. Knepley DMCeed sd; 137*d2b2dc1eSMatthew G. Knepley PetscDS ds; 138*d2b2dc1eSMatthew G. Knepley PetscFE fe; 139*d2b2dc1eSMatthew G. Knepley CeedQFunctionUser geom = NULL; 140*d2b2dc1eSMatthew G. Knepley const char *geomName = NULL; 141*d2b2dc1eSMatthew G. Knepley const PetscInt *cells; 142*d2b2dc1eSMatthew G. Knepley PetscInt dim, cdim, cStart, cEnd, Ncell, Nq; 143*d2b2dc1eSMatthew G. Knepley 144*d2b2dc1eSMatthew G. Knepley PetscFunctionBegin; 145*d2b2dc1eSMatthew G. Knepley PetscCall(PetscCalloc1(1, &sd)); 146*d2b2dc1eSMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 147*d2b2dc1eSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &cdim)); 148*d2b2dc1eSMatthew G. Knepley PetscCall(DMGetCeed(dm, &ceed)); 149*d2b2dc1eSMatthew G. Knepley PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells)); 150*d2b2dc1eSMatthew G. Knepley Ncell = cEnd - cStart; 151*d2b2dc1eSMatthew G. Knepley 152*d2b2dc1eSMatthew G. Knepley PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &ds, NULL)); 153*d2b2dc1eSMatthew G. Knepley PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe)); 154*d2b2dc1eSMatthew G. Knepley PetscCall(PetscFEGetCeedBasis(fe, &sd->basis)); 155*d2b2dc1eSMatthew G. Knepley PetscCall(CeedBasisGetNumQuadraturePoints(sd->basis, &Nq)); 156*d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexGetCeedRestriction(dm, NULL, 0, 0, 0, &sd->er)); 157*d2b2dc1eSMatthew G. Knepley 158*d2b2dc1eSMatthew G. Knepley *Nqdata = 1 + cdim + cdim * dim; 159*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedElemRestrictionCreateStrided(ceed, Ncell, Nq, *Nqdata, Ncell * Nq * (*Nqdata), CEED_STRIDES_BACKEND, erq)); 160*d2b2dc1eSMatthew G. Knepley 161*d2b2dc1eSMatthew G. Knepley switch (dim) { 162*d2b2dc1eSMatthew G. Knepley case 2: 163*d2b2dc1eSMatthew G. Knepley geom = Geometry2D; 164*d2b2dc1eSMatthew G. Knepley geomName = Geometry2D_loc; 165*d2b2dc1eSMatthew G. Knepley break; 166*d2b2dc1eSMatthew G. Knepley case 3: 167*d2b2dc1eSMatthew G. Knepley geom = Geometry3D; 168*d2b2dc1eSMatthew G. Knepley geomName = Geometry3D_loc; 169*d2b2dc1eSMatthew G. Knepley break; 170*d2b2dc1eSMatthew G. Knepley } 171*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedQFunctionCreateInterior(ceed, 1, geom, geomName, &sd->qf)); 172*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedQFunctionAddInput(sd->qf, "x", cdim, CEED_EVAL_INTERP)); 173*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedQFunctionAddInput(sd->qf, "dx", cdim * dim, CEED_EVAL_GRAD)); 174*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedQFunctionAddInput(sd->qf, "weight", 1, CEED_EVAL_WEIGHT)); 175*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedQFunctionAddOutput(sd->qf, "qdata", *Nqdata, CEED_EVAL_NONE)); 176*d2b2dc1eSMatthew G. Knepley 177*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedOperatorCreate(ceed, sd->qf, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &sd->op)); 178*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedOperatorSetField(sd->op, "x", sd->er, sd->basis, CEED_VECTOR_ACTIVE)); 179*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedOperatorSetField(sd->op, "dx", sd->er, sd->basis, CEED_VECTOR_ACTIVE)); 180*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedOperatorSetField(sd->op, "weight", CEED_ELEMRESTRICTION_NONE, sd->basis, CEED_VECTOR_NONE)); 181*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedOperatorSetField(sd->op, "qdata", *erq, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE)); 182*d2b2dc1eSMatthew G. Knepley 183*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedElemRestrictionCreateVector(*erq, qd, NULL)); 184*d2b2dc1eSMatthew G. Knepley *soldata = sd; 185*d2b2dc1eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 186*d2b2dc1eSMatthew G. Knepley } 187*d2b2dc1eSMatthew G. Knepley 188*d2b2dc1eSMatthew G. Knepley PetscErrorCode DMRefineHook_Ceed(DM coarse, DM fine, void *ctx) 189*d2b2dc1eSMatthew G. Knepley { 190*d2b2dc1eSMatthew G. Knepley PetscFunctionBegin; 191*d2b2dc1eSMatthew G. Knepley if (coarse->dmceed) PetscCall(DMCeedCreate(fine, coarse->dmceed->geom ? PETSC_TRUE : PETSC_FALSE, coarse->dmceed->func, coarse->dmceed->funcSource)); 192*d2b2dc1eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 193*d2b2dc1eSMatthew G. Knepley } 194*d2b2dc1eSMatthew G. Knepley 195*d2b2dc1eSMatthew G. Knepley PetscErrorCode DMCeedCreate_Internal(DM dm, IS cellIS, PetscBool createGeometry, CeedQFunctionUser func, const char *func_source, DMCeed *soldata) 196*d2b2dc1eSMatthew G. Knepley { 197*d2b2dc1eSMatthew G. Knepley PetscDS ds; 198*d2b2dc1eSMatthew G. Knepley PetscFE fe; 199*d2b2dc1eSMatthew G. Knepley DMCeed sd; 200*d2b2dc1eSMatthew G. Knepley Ceed ceed; 201*d2b2dc1eSMatthew G. Knepley PetscInt dim, Nc, Nq, Nqdata = 0; 202*d2b2dc1eSMatthew G. Knepley 203*d2b2dc1eSMatthew G. Knepley PetscFunctionBegin; 204*d2b2dc1eSMatthew G. Knepley PetscCall(PetscCalloc1(1, &sd)); 205*d2b2dc1eSMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 206*d2b2dc1eSMatthew G. Knepley PetscCall(DMGetCeed(dm, &ceed)); 207*d2b2dc1eSMatthew G. Knepley PetscCall(DMGetDS(dm, &ds)); 208*d2b2dc1eSMatthew G. Knepley PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe)); 209*d2b2dc1eSMatthew G. Knepley PetscCall(PetscFEGetCeedBasis(fe, &sd->basis)); 210*d2b2dc1eSMatthew G. Knepley PetscCall(PetscFEGetNumComponents(fe, &Nc)); 211*d2b2dc1eSMatthew G. Knepley PetscCall(CeedBasisGetNumQuadraturePoints(sd->basis, &Nq)); 212*d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexGetCeedRestriction(dm, NULL, 0, 0, 0, &sd->er)); 213*d2b2dc1eSMatthew G. Knepley 214*d2b2dc1eSMatthew G. Knepley if (createGeometry) { 215*d2b2dc1eSMatthew G. Knepley DM cdm; 216*d2b2dc1eSMatthew G. Knepley 217*d2b2dc1eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); 218*d2b2dc1eSMatthew G. Knepley PetscCall(DMCeedCreateGeometry(cdm, cellIS, &Nqdata, &sd->erq, &sd->qd, &sd->geom)); 219*d2b2dc1eSMatthew G. Knepley } 220*d2b2dc1eSMatthew G. Knepley 221*d2b2dc1eSMatthew G. Knepley if (sd->geom) { 222*d2b2dc1eSMatthew G. Knepley PetscInt cdim, Nqx; 223*d2b2dc1eSMatthew G. Knepley 224*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedBasisGetNumQuadraturePoints(sd->geom->basis, &Nqx)); 225*d2b2dc1eSMatthew G. Knepley PetscCheck(Nqx == Nq, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_INCOMP, "Number of qpoints for solution %" PetscInt_FMT " != %" PetscInt_FMT " Number of qpoints for coordinates", Nq, Nqx); 226*d2b2dc1eSMatthew G. Knepley /* TODO Remove this limitation */ 227*d2b2dc1eSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &cdim)); 228*d2b2dc1eSMatthew G. Knepley PetscCheck(dim == cdim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_INCOMP, "Topological dimension %" PetscInt_FMT " != %" PetscInt_FMT " embedding dimension", dim, cdim); 229*d2b2dc1eSMatthew G. Knepley } 230*d2b2dc1eSMatthew G. Knepley 231*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedQFunctionCreateInterior(ceed, 1, func, func_source, &sd->qf)); 232*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedQFunctionAddInput(sd->qf, "u", Nc, CEED_EVAL_INTERP)); 233*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedQFunctionAddInput(sd->qf, "du", Nc * dim, CEED_EVAL_GRAD)); 234*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedQFunctionAddInput(sd->qf, "qdata", Nqdata, CEED_EVAL_NONE)); 235*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedQFunctionAddOutput(sd->qf, "v", Nc, CEED_EVAL_INTERP)); 236*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedQFunctionAddOutput(sd->qf, "dv", Nc * dim, CEED_EVAL_GRAD)); 237*d2b2dc1eSMatthew G. Knepley 238*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedOperatorCreate(ceed, sd->qf, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &sd->op)); 239*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedOperatorSetField(sd->op, "u", sd->er, sd->basis, CEED_VECTOR_ACTIVE)); 240*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedOperatorSetField(sd->op, "du", sd->er, sd->basis, CEED_VECTOR_ACTIVE)); 241*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedOperatorSetField(sd->op, "qdata", sd->erq, CEED_BASIS_COLLOCATED, sd->qd)); 242*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedOperatorSetField(sd->op, "v", sd->er, sd->basis, CEED_VECTOR_ACTIVE)); 243*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedOperatorSetField(sd->op, "dv", sd->er, sd->basis, CEED_VECTOR_ACTIVE)); 244*d2b2dc1eSMatthew G. Knepley 245*d2b2dc1eSMatthew G. Knepley // Handle refinement 246*d2b2dc1eSMatthew G. Knepley sd->func = func; 247*d2b2dc1eSMatthew G. Knepley PetscCall(PetscStrallocpy(func_source, &sd->funcSource)); 248*d2b2dc1eSMatthew G. Knepley PetscCall(DMRefineHookAdd(dm, DMRefineHook_Ceed, NULL, NULL)); 249*d2b2dc1eSMatthew G. Knepley 250*d2b2dc1eSMatthew G. Knepley *soldata = sd; 251*d2b2dc1eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 252*d2b2dc1eSMatthew G. Knepley } 253*d2b2dc1eSMatthew G. Knepley 254*d2b2dc1eSMatthew G. Knepley PetscErrorCode DMCeedCreate(DM dm, PetscBool createGeometry, CeedQFunctionUser func, const char *func_source) 255*d2b2dc1eSMatthew G. Knepley { 256*d2b2dc1eSMatthew G. Knepley DM plex; 257*d2b2dc1eSMatthew G. Knepley IS cellIS; 258*d2b2dc1eSMatthew G. Knepley 259*d2b2dc1eSMatthew G. Knepley PetscFunctionBegin; 260*d2b2dc1eSMatthew G. Knepley PetscCall(DMConvert(dm, DMPLEX, &plex)); 261*d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexGetAllCells_Internal(plex, &cellIS)); 262*d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED 263*d2b2dc1eSMatthew G. Knepley PetscCall(DMCeedCreate_Internal(dm, cellIS, createGeometry, func, func_source, &dm->dmceed)); 264*d2b2dc1eSMatthew G. Knepley #endif 265*d2b2dc1eSMatthew G. Knepley PetscCall(ISDestroy(&cellIS)); 266*d2b2dc1eSMatthew G. Knepley PetscCall(DMDestroy(&plex)); 267*d2b2dc1eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 2687ce467fcSMatthew G. Knepley } 2697ce467fcSMatthew G. Knepley 270f918ec44SMatthew G. Knepley #endif 271*d2b2dc1eSMatthew G. Knepley 272*d2b2dc1eSMatthew G. Knepley PetscErrorCode DMCeedDestroy(DMCeed *pceed) 273*d2b2dc1eSMatthew G. Knepley { 274*d2b2dc1eSMatthew G. Knepley DMCeed p = *pceed; 275*d2b2dc1eSMatthew G. Knepley 276*d2b2dc1eSMatthew G. Knepley PetscFunctionBegin; 277*d2b2dc1eSMatthew G. Knepley if (!p) PetscFunctionReturn(PETSC_SUCCESS); 278*d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED 279*d2b2dc1eSMatthew G. Knepley PetscCall(PetscFree(p->funcSource)); 280*d2b2dc1eSMatthew G. Knepley if (p->qd) PetscCallCEED(CeedVectorDestroy(&p->qd)); 281*d2b2dc1eSMatthew G. Knepley if (p->op) PetscCallCEED(CeedOperatorDestroy(&p->op)); 282*d2b2dc1eSMatthew G. Knepley if (p->qf) PetscCallCEED(CeedQFunctionDestroy(&p->qf)); 283*d2b2dc1eSMatthew G. Knepley if (p->erq) PetscCallCEED(CeedElemRestrictionDestroy(&p->erq)); 284*d2b2dc1eSMatthew G. Knepley PetscCall(DMCeedDestroy(&p->geom)); 285*d2b2dc1eSMatthew G. Knepley #endif 286*d2b2dc1eSMatthew G. Knepley PetscCall(PetscFree(p)); 287*d2b2dc1eSMatthew G. Knepley *pceed = NULL; 288*d2b2dc1eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 289*d2b2dc1eSMatthew G. Knepley } 290*d2b2dc1eSMatthew G. Knepley 291*d2b2dc1eSMatthew G. Knepley PetscErrorCode DMCeedComputeGeometry(DM dm, DMCeed sd) 292*d2b2dc1eSMatthew G. Knepley { 293*d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED 294*d2b2dc1eSMatthew G. Knepley Ceed ceed; 295*d2b2dc1eSMatthew G. Knepley Vec coords; 296*d2b2dc1eSMatthew G. Knepley CeedVector ccoords; 297*d2b2dc1eSMatthew G. Knepley #endif 298*d2b2dc1eSMatthew G. Knepley 299*d2b2dc1eSMatthew G. Knepley PetscFunctionBegin; 300*d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED 301*d2b2dc1eSMatthew G. Knepley PetscCall(DMGetCeed(dm, &ceed)); 302*d2b2dc1eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocal(dm, &coords)); 303*d2b2dc1eSMatthew G. Knepley PetscCall(VecGetCeedVectorRead(coords, ceed, &ccoords)); 304*d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedOperatorApply(sd->geom->op, ccoords, sd->qd, CEED_REQUEST_IMMEDIATE)); 305*d2b2dc1eSMatthew G. Knepley PetscCall(VecRestoreCeedVectorRead(coords, &ccoords)); 306*d2b2dc1eSMatthew G. Knepley #endif 307*d2b2dc1eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 308*d2b2dc1eSMatthew G. Knepley } 309