14411cf47Sjeremylt /// @file 24411cf47Sjeremylt /// Test creation, evaluation, and destruction for qfunction 34411cf47Sjeremylt /// \test Test creation, evaluation, and destruction for qfunction 457c64913Sjeremylt #include <ceed.h> 557c64913Sjeremylt 657c64913Sjeremylt static int setup(void *ctx, CeedInt Q, const CeedScalar *const *in, 757c64913Sjeremylt CeedScalar *const *out) { 857c64913Sjeremylt const CeedScalar *w = in[0]; 957c64913Sjeremylt CeedScalar *qdata = out[0]; 1057c64913Sjeremylt for (CeedInt i=0; i<Q; i++) { 1157c64913Sjeremylt qdata[i] = w[i]; 1257c64913Sjeremylt } 1357c64913Sjeremylt return 0; 1457c64913Sjeremylt } 1557c64913Sjeremylt 1657c64913Sjeremylt static int mass(void *ctx, CeedInt Q, const CeedScalar *const *in, 1757c64913Sjeremylt CeedScalar *const *out) { 1857c64913Sjeremylt const CeedScalar *qdata = in[0], *u = in[1]; 1957c64913Sjeremylt CeedScalar *v = out[0]; 2057c64913Sjeremylt for (CeedInt i=0; i<Q; i++) { 2157c64913Sjeremylt v[i] = qdata[i] * u[i]; 2257c64913Sjeremylt } 2357c64913Sjeremylt return 0; 2457c64913Sjeremylt } 2557c64913Sjeremylt 2657c64913Sjeremylt int main(int argc, char **argv) { 2757c64913Sjeremylt Ceed ceed; 28aedaa0e5Sjeremylt CeedVector in[16], out[16]; 29aedaa0e5Sjeremylt CeedVector Qdata, W, U, V; 3057c64913Sjeremylt CeedQFunction qf_setup, qf_mass; 3157c64913Sjeremylt CeedInt Q = 8; 32aedaa0e5Sjeremylt const CeedScalar *vv; 33aedaa0e5Sjeremylt CeedScalar w[Q], u[Q], v[Q]; 3457c64913Sjeremylt 3557c64913Sjeremylt 3657c64913Sjeremylt CeedInit(argv[1], &ceed); 37aedaa0e5Sjeremylt 3857c64913Sjeremylt CeedQFunctionCreateInterior(ceed, 1, setup, __FILE__ ":setup", &qf_setup); 3957c64913Sjeremylt CeedQFunctionAddInput(qf_setup, "w", 1, CEED_EVAL_INTERP); 4057c64913Sjeremylt CeedQFunctionAddOutput(qf_setup, "qdata", 1, CEED_EVAL_INTERP); 4157c64913Sjeremylt 4257c64913Sjeremylt CeedQFunctionCreateInterior(ceed, 1, mass, __FILE__ ":mass", &qf_mass); 4357c64913Sjeremylt CeedQFunctionAddInput(qf_mass, "qdata", 1, CEED_EVAL_INTERP); 4457c64913Sjeremylt CeedQFunctionAddInput(qf_mass, "u", 1, CEED_EVAL_INTERP); 4557c64913Sjeremylt CeedQFunctionAddOutput(qf_mass, "v", 1, CEED_EVAL_INTERP); 4657c64913Sjeremylt 4757c64913Sjeremylt for (CeedInt i=0; i<Q; i++) { 4857c64913Sjeremylt CeedScalar x = 2.*i/(Q-1) - 1; 4957c64913Sjeremylt w[i] = 1 - x*x; 5057c64913Sjeremylt u[i] = 2 + 3*x + 5*x*x; 5157c64913Sjeremylt v[i] = w[i] * u[i]; 5257c64913Sjeremylt } 53aedaa0e5Sjeremylt 54aedaa0e5Sjeremylt CeedVectorCreate(ceed, Q, &W); 55aedaa0e5Sjeremylt CeedVectorSetArray(W, CEED_MEM_HOST, CEED_USE_POINTER, (CeedScalar *)&w); 56aedaa0e5Sjeremylt CeedVectorCreate(ceed, Q, &U); 57aedaa0e5Sjeremylt CeedVectorSetArray(U, CEED_MEM_HOST, CEED_USE_POINTER, (CeedScalar *)&u); 58aedaa0e5Sjeremylt CeedVectorCreate(ceed, Q, &V); 59aedaa0e5Sjeremylt CeedVectorSetValue(V, 0); 60aedaa0e5Sjeremylt CeedVectorCreate(ceed, Q, &Qdata); 61aedaa0e5Sjeremylt CeedVectorSetValue(Qdata, 0); 62aedaa0e5Sjeremylt 6357c64913Sjeremylt { 64aedaa0e5Sjeremylt in[0] = W; 65aedaa0e5Sjeremylt out[0] = Qdata; 6657c64913Sjeremylt CeedQFunctionApply(qf_setup, Q, in, out); 6757c64913Sjeremylt } 6857c64913Sjeremylt { 69aedaa0e5Sjeremylt in[0] = W; 70aedaa0e5Sjeremylt in[1] = U; 71aedaa0e5Sjeremylt out[0] = V; 7257c64913Sjeremylt CeedQFunctionApply(qf_mass, Q, in, out); 7357c64913Sjeremylt } 74aedaa0e5Sjeremylt 75aedaa0e5Sjeremylt CeedVectorGetArrayRead(V, CEED_MEM_HOST, &vv); 76a2546046Sjeremylt for (CeedInt i=0; i<Q; i++) 77a2546046Sjeremylt if (v[i] != vv[i]) 78a2546046Sjeremylt // LCOV_EXCL_START 79a2546046Sjeremylt printf("[%d] v %f != vv %f\n",i, v[i], vv[i]); 80*de996c55Sjeremylt // LCOV_EXCL_STOP 81aedaa0e5Sjeremylt CeedVectorRestoreArrayRead(V, &vv); 82aedaa0e5Sjeremylt 83aedaa0e5Sjeremylt CeedVectorDestroy(&W); 84aedaa0e5Sjeremylt CeedVectorDestroy(&U); 85aedaa0e5Sjeremylt CeedVectorDestroy(&V); 86aedaa0e5Sjeremylt CeedVectorDestroy(&Qdata); 8757c64913Sjeremylt CeedQFunctionDestroy(&qf_setup); 8857c64913Sjeremylt CeedQFunctionDestroy(&qf_mass); 8957c64913Sjeremylt CeedDestroy(&ceed); 9057c64913Sjeremylt return 0; 9157c64913Sjeremylt } 92