xref: /petsc/src/dm/interface/dmceed.c (revision 5962854d720b7b8d98c62edd758f00bbb980e600)
1f918ec44SMatthew G. Knepley #include <petsc/private/dmimpl.h> /*I      "petscdm.h"          I*/
2d2b2dc1eSMatthew G. Knepley #include <petscdmceed.h>
3f918ec44SMatthew G. Knepley 
4f918ec44SMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
5d2b2dc1eSMatthew G. Knepley   #include <petsc/private/dmpleximpl.h>
6d2b2dc1eSMatthew G. Knepley   #include <petscdmplexceed.h>
7d2b2dc1eSMatthew 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 
42d2b2dc1eSMatthew G. Knepley static CeedMemType PetscMemType2Ceed(PetscMemType mem_type)
43d2b2dc1eSMatthew 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));
58d2b2dc1eSMatthew 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));
66d2b2dc1eSMatthew 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));
79d2b2dc1eSMatthew 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));
87d2b2dc1eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
88d2b2dc1eSMatthew G. Knepley }
89d2b2dc1eSMatthew G. Knepley 
90d2b2dc1eSMatthew G. Knepley CEED_QFUNCTION(Geometry2D)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)
91d2b2dc1eSMatthew G. Knepley {
92d2b2dc1eSMatthew G. Knepley   const CeedScalar *x = in[0], *Jac = in[1], *w = in[2];
93d2b2dc1eSMatthew G. Knepley   CeedScalar       *qdata = out[0];
94d2b2dc1eSMatthew G. Knepley 
95d2b2dc1eSMatthew G. Knepley   CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i)
96d2b2dc1eSMatthew G. Knepley   {
97d2b2dc1eSMatthew G. Knepley     const CeedScalar J[2][2] = {
98d2b2dc1eSMatthew G. Knepley       {Jac[i + Q * 0], Jac[i + Q * 2]},
99d2b2dc1eSMatthew G. Knepley       {Jac[i + Q * 1], Jac[i + Q * 3]}
100d2b2dc1eSMatthew G. Knepley     };
101d2b2dc1eSMatthew G. Knepley     const CeedScalar det = J[0][0] * J[1][1] - J[0][1] * J[1][0];
102d2b2dc1eSMatthew G. Knepley 
103d2b2dc1eSMatthew G. Knepley     qdata[i + Q * 0] = det * w[i];
104d2b2dc1eSMatthew G. Knepley     qdata[i + Q * 1] = x[i + Q * 0];
105d2b2dc1eSMatthew G. Knepley     qdata[i + Q * 2] = x[i + Q * 1];
106d2b2dc1eSMatthew G. Knepley     qdata[i + Q * 3] = J[1][1] / det;
107d2b2dc1eSMatthew G. Knepley     qdata[i + Q * 4] = -J[1][0] / det;
108d2b2dc1eSMatthew G. Knepley     qdata[i + Q * 5] = -J[0][1] / det;
109d2b2dc1eSMatthew G. Knepley     qdata[i + Q * 6] = J[0][0] / det;
110d2b2dc1eSMatthew G. Knepley   }
111d2b2dc1eSMatthew G. Knepley   return CEED_ERROR_SUCCESS;
112d2b2dc1eSMatthew G. Knepley }
113d2b2dc1eSMatthew G. Knepley 
114d2b2dc1eSMatthew G. Knepley CEED_QFUNCTION(Geometry3D)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)
115d2b2dc1eSMatthew G. Knepley {
116d2b2dc1eSMatthew G. Knepley   const CeedScalar *Jac = in[1], *w = in[2];
117d2b2dc1eSMatthew G. Knepley   CeedScalar       *qdata = out[0];
118d2b2dc1eSMatthew G. Knepley 
119d2b2dc1eSMatthew G. Knepley   CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i)
120d2b2dc1eSMatthew G. Knepley   {
121d2b2dc1eSMatthew G. Knepley     const CeedScalar J[3][3] = {
122d2b2dc1eSMatthew G. Knepley       {Jac[i + Q * 0], Jac[i + Q * 3], Jac[i + Q * 6]},
123d2b2dc1eSMatthew G. Knepley       {Jac[i + Q * 1], Jac[i + Q * 4], Jac[i + Q * 7]},
124d2b2dc1eSMatthew G. Knepley       {Jac[i + Q * 2], Jac[i + Q * 5], Jac[i + Q * 8]}
125d2b2dc1eSMatthew G. Knepley     };
126d2b2dc1eSMatthew 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]);
127d2b2dc1eSMatthew G. Knepley 
128d2b2dc1eSMatthew G. Knepley     qdata[i + Q * 0] = det * w[i]; /* det J * weight */
129d2b2dc1eSMatthew G. Knepley   }
130d2b2dc1eSMatthew G. Knepley   return CEED_ERROR_SUCCESS;
131d2b2dc1eSMatthew G. Knepley }
132d2b2dc1eSMatthew G. Knepley 
133d2b2dc1eSMatthew G. Knepley static PetscErrorCode DMCeedCreateGeometry(DM dm, IS cellIS, PetscInt *Nqdata, CeedElemRestriction *erq, CeedVector *qd, DMCeed *soldata)
134d2b2dc1eSMatthew G. Knepley {
135d2b2dc1eSMatthew G. Knepley   Ceed              ceed;
136d2b2dc1eSMatthew G. Knepley   DMCeed            sd;
137d2b2dc1eSMatthew G. Knepley   PetscDS           ds;
138d2b2dc1eSMatthew G. Knepley   PetscFE           fe;
139d2b2dc1eSMatthew G. Knepley   CeedQFunctionUser geom     = NULL;
140d2b2dc1eSMatthew G. Knepley   const char       *geomName = NULL;
141d2b2dc1eSMatthew G. Knepley   const PetscInt   *cells;
142d2b2dc1eSMatthew G. Knepley   PetscInt          dim, cdim, cStart, cEnd, Ncell, Nq;
143d2b2dc1eSMatthew G. Knepley 
144d2b2dc1eSMatthew G. Knepley   PetscFunctionBegin;
145d2b2dc1eSMatthew G. Knepley   PetscCall(PetscCalloc1(1, &sd));
146d2b2dc1eSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
147d2b2dc1eSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
148d2b2dc1eSMatthew G. Knepley   PetscCall(DMGetCeed(dm, &ceed));
149d2b2dc1eSMatthew G. Knepley   PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
150d2b2dc1eSMatthew G. Knepley   Ncell = cEnd - cStart;
151d2b2dc1eSMatthew G. Knepley 
152d2b2dc1eSMatthew G. Knepley   PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &ds, NULL));
153d2b2dc1eSMatthew G. Knepley   PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe));
154d2b2dc1eSMatthew G. Knepley   PetscCall(PetscFEGetCeedBasis(fe, &sd->basis));
155d2b2dc1eSMatthew G. Knepley   PetscCall(CeedBasisGetNumQuadraturePoints(sd->basis, &Nq));
156d2b2dc1eSMatthew G. Knepley   PetscCall(DMPlexGetCeedRestriction(dm, NULL, 0, 0, 0, &sd->er));
157d2b2dc1eSMatthew G. Knepley 
158*5962854dSMatthew G. Knepley   *Nqdata = 1 + cdim + cdim * dim; // |J| * w_q, x, J^{-1}
159d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedElemRestrictionCreateStrided(ceed, Ncell, Nq, *Nqdata, Ncell * Nq * (*Nqdata), CEED_STRIDES_BACKEND, erq));
160d2b2dc1eSMatthew G. Knepley 
161d2b2dc1eSMatthew G. Knepley   switch (dim) {
162d2b2dc1eSMatthew G. Knepley   case 2:
163d2b2dc1eSMatthew G. Knepley     geom     = Geometry2D;
164d2b2dc1eSMatthew G. Knepley     geomName = Geometry2D_loc;
165d2b2dc1eSMatthew G. Knepley     break;
166d2b2dc1eSMatthew G. Knepley   case 3:
167d2b2dc1eSMatthew G. Knepley     geom     = Geometry3D;
168d2b2dc1eSMatthew G. Knepley     geomName = Geometry3D_loc;
169d2b2dc1eSMatthew G. Knepley     break;
170d2b2dc1eSMatthew G. Knepley   }
171d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedQFunctionCreateInterior(ceed, 1, geom, geomName, &sd->qf));
172d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedQFunctionAddInput(sd->qf, "x", cdim, CEED_EVAL_INTERP));
173d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedQFunctionAddInput(sd->qf, "dx", cdim * dim, CEED_EVAL_GRAD));
174d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedQFunctionAddInput(sd->qf, "weight", 1, CEED_EVAL_WEIGHT));
175d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedQFunctionAddOutput(sd->qf, "qdata", *Nqdata, CEED_EVAL_NONE));
176d2b2dc1eSMatthew G. Knepley 
177d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedOperatorCreate(ceed, sd->qf, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &sd->op));
178d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedOperatorSetField(sd->op, "x", sd->er, sd->basis, CEED_VECTOR_ACTIVE));
179d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedOperatorSetField(sd->op, "dx", sd->er, sd->basis, CEED_VECTOR_ACTIVE));
180d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedOperatorSetField(sd->op, "weight", CEED_ELEMRESTRICTION_NONE, sd->basis, CEED_VECTOR_NONE));
181d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedOperatorSetField(sd->op, "qdata", *erq, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
182d2b2dc1eSMatthew G. Knepley 
183d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedElemRestrictionCreateVector(*erq, qd, NULL));
184d2b2dc1eSMatthew G. Knepley   *soldata = sd;
185d2b2dc1eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
186d2b2dc1eSMatthew G. Knepley }
187d2b2dc1eSMatthew G. Knepley 
188d2b2dc1eSMatthew G. Knepley PetscErrorCode DMRefineHook_Ceed(DM coarse, DM fine, void *ctx)
189d2b2dc1eSMatthew G. Knepley {
190d2b2dc1eSMatthew G. Knepley   PetscFunctionBegin;
191d2b2dc1eSMatthew G. Knepley   if (coarse->dmceed) PetscCall(DMCeedCreate(fine, coarse->dmceed->geom ? PETSC_TRUE : PETSC_FALSE, coarse->dmceed->func, coarse->dmceed->funcSource));
192d2b2dc1eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
193d2b2dc1eSMatthew G. Knepley }
194d2b2dc1eSMatthew G. Knepley 
195d2b2dc1eSMatthew G. Knepley PetscErrorCode DMCeedCreate_Internal(DM dm, IS cellIS, PetscBool createGeometry, CeedQFunctionUser func, const char *func_source, DMCeed *soldata)
196d2b2dc1eSMatthew G. Knepley {
197d2b2dc1eSMatthew G. Knepley   PetscDS  ds;
198d2b2dc1eSMatthew G. Knepley   PetscFE  fe;
199d2b2dc1eSMatthew G. Knepley   DMCeed   sd;
200d2b2dc1eSMatthew G. Knepley   Ceed     ceed;
201d2b2dc1eSMatthew G. Knepley   PetscInt dim, Nc, Nq, Nqdata = 0;
202d2b2dc1eSMatthew G. Knepley 
203d2b2dc1eSMatthew G. Knepley   PetscFunctionBegin;
204d2b2dc1eSMatthew G. Knepley   PetscCall(PetscCalloc1(1, &sd));
205d2b2dc1eSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
206d2b2dc1eSMatthew G. Knepley   PetscCall(DMGetCeed(dm, &ceed));
207d2b2dc1eSMatthew G. Knepley   PetscCall(DMGetDS(dm, &ds));
208d2b2dc1eSMatthew G. Knepley   PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe));
209d2b2dc1eSMatthew G. Knepley   PetscCall(PetscFEGetCeedBasis(fe, &sd->basis));
210d2b2dc1eSMatthew G. Knepley   PetscCall(PetscFEGetNumComponents(fe, &Nc));
211d2b2dc1eSMatthew G. Knepley   PetscCall(CeedBasisGetNumQuadraturePoints(sd->basis, &Nq));
212d2b2dc1eSMatthew G. Knepley   PetscCall(DMPlexGetCeedRestriction(dm, NULL, 0, 0, 0, &sd->er));
213d2b2dc1eSMatthew G. Knepley 
214d2b2dc1eSMatthew G. Knepley   if (createGeometry) {
215d2b2dc1eSMatthew G. Knepley     DM cdm;
216d2b2dc1eSMatthew G. Knepley 
217d2b2dc1eSMatthew G. Knepley     PetscCall(DMGetCoordinateDM(dm, &cdm));
218d2b2dc1eSMatthew G. Knepley     PetscCall(DMCeedCreateGeometry(cdm, cellIS, &Nqdata, &sd->erq, &sd->qd, &sd->geom));
219d2b2dc1eSMatthew G. Knepley   }
220d2b2dc1eSMatthew G. Knepley 
221d2b2dc1eSMatthew G. Knepley   if (sd->geom) {
222d2b2dc1eSMatthew G. Knepley     PetscInt cdim, Nqx;
223d2b2dc1eSMatthew G. Knepley 
224d2b2dc1eSMatthew G. Knepley     PetscCallCEED(CeedBasisGetNumQuadraturePoints(sd->geom->basis, &Nqx));
225d2b2dc1eSMatthew 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);
226d2b2dc1eSMatthew G. Knepley     /* TODO Remove this limitation */
227d2b2dc1eSMatthew G. Knepley     PetscCall(DMGetCoordinateDim(dm, &cdim));
228d2b2dc1eSMatthew G. Knepley     PetscCheck(dim == cdim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_INCOMP, "Topological dimension %" PetscInt_FMT " != %" PetscInt_FMT " embedding dimension", dim, cdim);
229d2b2dc1eSMatthew G. Knepley   }
230d2b2dc1eSMatthew G. Knepley 
231d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedQFunctionCreateInterior(ceed, 1, func, func_source, &sd->qf));
232d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedQFunctionAddInput(sd->qf, "u", Nc, CEED_EVAL_INTERP));
233d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedQFunctionAddInput(sd->qf, "du", Nc * dim, CEED_EVAL_GRAD));
234d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedQFunctionAddInput(sd->qf, "qdata", Nqdata, CEED_EVAL_NONE));
235d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedQFunctionAddOutput(sd->qf, "v", Nc, CEED_EVAL_INTERP));
236d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedQFunctionAddOutput(sd->qf, "dv", Nc * dim, CEED_EVAL_GRAD));
237d2b2dc1eSMatthew G. Knepley 
238d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedOperatorCreate(ceed, sd->qf, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &sd->op));
239d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedOperatorSetField(sd->op, "u", sd->er, sd->basis, CEED_VECTOR_ACTIVE));
240d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedOperatorSetField(sd->op, "du", sd->er, sd->basis, CEED_VECTOR_ACTIVE));
241d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedOperatorSetField(sd->op, "qdata", sd->erq, CEED_BASIS_COLLOCATED, sd->qd));
242d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedOperatorSetField(sd->op, "v", sd->er, sd->basis, CEED_VECTOR_ACTIVE));
243d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedOperatorSetField(sd->op, "dv", sd->er, sd->basis, CEED_VECTOR_ACTIVE));
244d2b2dc1eSMatthew G. Knepley 
245d2b2dc1eSMatthew G. Knepley   // Handle refinement
246d2b2dc1eSMatthew G. Knepley   sd->func = func;
247d2b2dc1eSMatthew G. Knepley   PetscCall(PetscStrallocpy(func_source, &sd->funcSource));
248d2b2dc1eSMatthew G. Knepley   PetscCall(DMRefineHookAdd(dm, DMRefineHook_Ceed, NULL, NULL));
249d2b2dc1eSMatthew G. Knepley 
250d2b2dc1eSMatthew G. Knepley   *soldata = sd;
251d2b2dc1eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
252d2b2dc1eSMatthew G. Knepley }
253d2b2dc1eSMatthew G. Knepley 
254d2b2dc1eSMatthew G. Knepley PetscErrorCode DMCeedCreate(DM dm, PetscBool createGeometry, CeedQFunctionUser func, const char *func_source)
255d2b2dc1eSMatthew G. Knepley {
256d2b2dc1eSMatthew G. Knepley   DM plex;
257d2b2dc1eSMatthew G. Knepley   IS cellIS;
258d2b2dc1eSMatthew G. Knepley 
259d2b2dc1eSMatthew G. Knepley   PetscFunctionBegin;
260d2b2dc1eSMatthew G. Knepley   PetscCall(DMConvert(dm, DMPLEX, &plex));
261d2b2dc1eSMatthew G. Knepley   PetscCall(DMPlexGetAllCells_Internal(plex, &cellIS));
262d2b2dc1eSMatthew G. Knepley   #ifdef PETSC_HAVE_LIBCEED
263d2b2dc1eSMatthew G. Knepley   PetscCall(DMCeedCreate_Internal(dm, cellIS, createGeometry, func, func_source, &dm->dmceed));
264d2b2dc1eSMatthew G. Knepley   #endif
265d2b2dc1eSMatthew G. Knepley   PetscCall(ISDestroy(&cellIS));
266d2b2dc1eSMatthew G. Knepley   PetscCall(DMDestroy(&plex));
267d2b2dc1eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
2687ce467fcSMatthew G. Knepley }
2697ce467fcSMatthew G. Knepley 
270*5962854dSMatthew G. Knepley static PetscErrorCode DMCeedCreateGeometryFVM(DM dm, IS faceIS, PetscInt *Nqdata, CeedElemRestriction *erq, CeedVector *qd, DMCeed *soldata)
271*5962854dSMatthew G. Knepley {
272*5962854dSMatthew G. Knepley   Ceed            ceed;
273*5962854dSMatthew G. Knepley   DMCeed          sd;
274*5962854dSMatthew G. Knepley   const PetscInt *faces;
275*5962854dSMatthew G. Knepley   PetscInt        dim, cdim, fStart, fEnd, Nface, Nq = 1;
276*5962854dSMatthew G. Knepley 
277*5962854dSMatthew G. Knepley   PetscFunctionBegin;
278*5962854dSMatthew G. Knepley   PetscCall(PetscCalloc1(1, &sd));
279*5962854dSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
280*5962854dSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
281*5962854dSMatthew G. Knepley   PetscCall(DMGetCeed(dm, &ceed));
282*5962854dSMatthew G. Knepley   PetscCall(ISGetPointRange(faceIS, &fStart, &fEnd, &faces));
283*5962854dSMatthew G. Knepley   Nface = fEnd - fStart;
284*5962854dSMatthew G. Knepley 
285*5962854dSMatthew G. Knepley   *Nqdata = cdim + 2; // face normal and support cell volumes
286*5962854dSMatthew G. Knepley   PetscCallCEED(CeedElemRestrictionCreateStrided(ceed, Nface, Nq, *Nqdata, Nface * Nq * (*Nqdata), CEED_STRIDES_BACKEND, erq));
287*5962854dSMatthew G. Knepley   PetscCallCEED(CeedElemRestrictionCreateVector(*erq, qd, NULL));
288*5962854dSMatthew G. Knepley   *soldata = sd;
289*5962854dSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
290*5962854dSMatthew G. Knepley }
291*5962854dSMatthew G. Knepley 
292*5962854dSMatthew G. Knepley PetscErrorCode DMCeedCreateFVM_Internal(DM dm, IS faceIS, PetscBool createGeometry, CeedQFunctionUser func, const char *func_source, DMCeed *soldata, CeedQFunctionContext qfCtx)
293*5962854dSMatthew G. Knepley {
294*5962854dSMatthew G. Knepley   PetscDS  ds;
295*5962854dSMatthew G. Knepley   PetscFV  fv;
296*5962854dSMatthew G. Knepley   DMCeed   sd;
297*5962854dSMatthew G. Knepley   Ceed     ceed;
298*5962854dSMatthew G. Knepley   PetscInt dim, Nc, Nqdata = 0;
299*5962854dSMatthew G. Knepley 
300*5962854dSMatthew G. Knepley   PetscFunctionBegin;
301*5962854dSMatthew G. Knepley   PetscCall(PetscCalloc1(1, &sd));
302*5962854dSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
303*5962854dSMatthew G. Knepley   PetscCall(DMGetCeed(dm, &ceed));
304*5962854dSMatthew G. Knepley   PetscCall(DMGetDS(dm, &ds));
305*5962854dSMatthew G. Knepley   PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fv));
306*5962854dSMatthew G. Knepley   PetscCall(PetscFVGetNumComponents(fv, &Nc));
307*5962854dSMatthew G. Knepley   PetscCall(DMPlexCreateCeedRestrictionFVM(dm, &sd->erL, &sd->erR));
308*5962854dSMatthew G. Knepley 
309*5962854dSMatthew G. Knepley   if (createGeometry) {
310*5962854dSMatthew G. Knepley     DM cdm;
311*5962854dSMatthew G. Knepley 
312*5962854dSMatthew G. Knepley     PetscCall(DMGetCoordinateDM(dm, &cdm));
313*5962854dSMatthew G. Knepley     PetscCall(DMCeedCreateGeometryFVM(cdm, faceIS, &Nqdata, &sd->erq, &sd->qd, &sd->geom));
314*5962854dSMatthew G. Knepley   }
315*5962854dSMatthew G. Knepley 
316*5962854dSMatthew G. Knepley   PetscCallCEED(CeedQFunctionCreateInterior(ceed, 1, func, func_source, &sd->qf));
317*5962854dSMatthew G. Knepley   PetscCallCEED(CeedQFunctionAddInput(sd->qf, "uL", Nc, CEED_EVAL_NONE));
318*5962854dSMatthew G. Knepley   PetscCallCEED(CeedQFunctionAddInput(sd->qf, "uR", Nc, CEED_EVAL_NONE));
319*5962854dSMatthew G. Knepley   PetscCallCEED(CeedQFunctionAddInput(sd->qf, "geom", Nqdata, CEED_EVAL_NONE));
320*5962854dSMatthew G. Knepley   PetscCallCEED(CeedQFunctionAddOutput(sd->qf, "cL", Nc, CEED_EVAL_NONE));
321*5962854dSMatthew G. Knepley   PetscCallCEED(CeedQFunctionAddOutput(sd->qf, "cR", Nc, CEED_EVAL_NONE));
322*5962854dSMatthew G. Knepley 
323*5962854dSMatthew G. Knepley   PetscCallCEED(CeedQFunctionSetContext(sd->qf, qfCtx));
324*5962854dSMatthew G. Knepley 
325*5962854dSMatthew G. Knepley   PetscCallCEED(CeedOperatorCreate(ceed, sd->qf, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &sd->op));
326*5962854dSMatthew G. Knepley   PetscCallCEED(CeedOperatorSetField(sd->op, "uL", sd->erL, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
327*5962854dSMatthew G. Knepley   PetscCallCEED(CeedOperatorSetField(sd->op, "uR", sd->erR, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
328*5962854dSMatthew G. Knepley   PetscCallCEED(CeedOperatorSetField(sd->op, "geom", sd->erq, CEED_BASIS_COLLOCATED, sd->qd));
329*5962854dSMatthew G. Knepley   PetscCallCEED(CeedOperatorSetField(sd->op, "cL", sd->erL, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
330*5962854dSMatthew G. Knepley   PetscCallCEED(CeedOperatorSetField(sd->op, "cR", sd->erR, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
331*5962854dSMatthew G. Knepley 
332*5962854dSMatthew G. Knepley   // Handle refinement
333*5962854dSMatthew G. Knepley   sd->func = func;
334*5962854dSMatthew G. Knepley   PetscCall(PetscStrallocpy(func_source, &sd->funcSource));
335*5962854dSMatthew G. Knepley   PetscCall(DMRefineHookAdd(dm, DMRefineHook_Ceed, NULL, NULL));
336*5962854dSMatthew G. Knepley 
337*5962854dSMatthew G. Knepley   *soldata = sd;
338*5962854dSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
339*5962854dSMatthew G. Knepley }
340*5962854dSMatthew G. Knepley 
341*5962854dSMatthew G. Knepley PetscErrorCode DMCeedCreateFVM(DM dm, PetscBool createGeometry, CeedQFunctionUser func, const char *func_source, CeedQFunctionContext qfCtx)
342*5962854dSMatthew G. Knepley {
343*5962854dSMatthew G. Knepley   DM plex;
344*5962854dSMatthew G. Knepley   IS faceIS;
345*5962854dSMatthew G. Knepley 
346*5962854dSMatthew G. Knepley   PetscFunctionBegin;
347*5962854dSMatthew G. Knepley   PetscCall(DMConvert(dm, DMPLEX, &plex));
348*5962854dSMatthew G. Knepley   PetscCall(DMPlexGetAllFaces_Internal(plex, &faceIS));
349*5962854dSMatthew G. Knepley   #ifdef PETSC_HAVE_LIBCEED
350*5962854dSMatthew G. Knepley   PetscCall(DMCeedCreateFVM_Internal(dm, faceIS, createGeometry, func, func_source, &dm->dmceed, qfCtx));
351*5962854dSMatthew G. Knepley   #endif
352*5962854dSMatthew G. Knepley   PetscCall(ISDestroy(&faceIS));
353*5962854dSMatthew G. Knepley   PetscCall(DMDestroy(&plex));
354*5962854dSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
355*5962854dSMatthew G. Knepley }
356*5962854dSMatthew G. Knepley 
357f918ec44SMatthew G. Knepley #endif
358d2b2dc1eSMatthew G. Knepley 
359d2b2dc1eSMatthew G. Knepley PetscErrorCode DMCeedDestroy(DMCeed *pceed)
360d2b2dc1eSMatthew G. Knepley {
361d2b2dc1eSMatthew G. Knepley   DMCeed p = *pceed;
362d2b2dc1eSMatthew G. Knepley 
363d2b2dc1eSMatthew G. Knepley   PetscFunctionBegin;
364d2b2dc1eSMatthew G. Knepley   if (!p) PetscFunctionReturn(PETSC_SUCCESS);
365d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
366d2b2dc1eSMatthew G. Knepley   PetscCall(PetscFree(p->funcSource));
367d2b2dc1eSMatthew G. Knepley   if (p->qd) PetscCallCEED(CeedVectorDestroy(&p->qd));
368d2b2dc1eSMatthew G. Knepley   if (p->op) PetscCallCEED(CeedOperatorDestroy(&p->op));
369d2b2dc1eSMatthew G. Knepley   if (p->qf) PetscCallCEED(CeedQFunctionDestroy(&p->qf));
370*5962854dSMatthew G. Knepley   if (p->erL) PetscCallCEED(CeedElemRestrictionDestroy(&p->erL));
371*5962854dSMatthew G. Knepley   if (p->erR) PetscCallCEED(CeedElemRestrictionDestroy(&p->erR));
372d2b2dc1eSMatthew G. Knepley   if (p->erq) PetscCallCEED(CeedElemRestrictionDestroy(&p->erq));
373d2b2dc1eSMatthew G. Knepley   PetscCall(DMCeedDestroy(&p->geom));
374d2b2dc1eSMatthew G. Knepley #endif
375d2b2dc1eSMatthew G. Knepley   PetscCall(PetscFree(p));
376d2b2dc1eSMatthew G. Knepley   *pceed = NULL;
377d2b2dc1eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
378d2b2dc1eSMatthew G. Knepley }
379d2b2dc1eSMatthew G. Knepley 
380d2b2dc1eSMatthew G. Knepley PetscErrorCode DMCeedComputeGeometry(DM dm, DMCeed sd)
381d2b2dc1eSMatthew G. Knepley {
382d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
383d2b2dc1eSMatthew G. Knepley   Ceed       ceed;
384d2b2dc1eSMatthew G. Knepley   Vec        coords;
385d2b2dc1eSMatthew G. Knepley   CeedVector ccoords;
386d2b2dc1eSMatthew G. Knepley #endif
387d2b2dc1eSMatthew G. Knepley 
388d2b2dc1eSMatthew G. Knepley   PetscFunctionBegin;
389d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
390d2b2dc1eSMatthew G. Knepley   PetscCall(DMGetCeed(dm, &ceed));
391d2b2dc1eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocal(dm, &coords));
392d2b2dc1eSMatthew G. Knepley   PetscCall(VecGetCeedVectorRead(coords, ceed, &ccoords));
393*5962854dSMatthew G. Knepley   if (sd->geom->op) PetscCallCEED(CeedOperatorApply(sd->geom->op, ccoords, sd->qd, CEED_REQUEST_IMMEDIATE));
394*5962854dSMatthew G. Knepley   else PetscCall(DMPlexCeedComputeGeometryFVM(dm, sd->qd));
395*5962854dSMatthew G. Knepley   //PetscCallCEED(CeedVectorView(sd->qd, "%g", stdout));
396d2b2dc1eSMatthew G. Knepley   PetscCall(VecRestoreCeedVectorRead(coords, &ccoords));
397d2b2dc1eSMatthew G. Knepley #endif
398d2b2dc1eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
399d2b2dc1eSMatthew G. Knepley }
400