xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-ref/ceed-cuda-ref-basis.c (revision 2b730f8b5a9c809740a0b3b302db43a719c636b1)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
30d0321e0SJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
50d0321e0SJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
70d0321e0SJeremy L Thompson 
80d0321e0SJeremy L Thompson #include <ceed/backend.h>
9*2b730f8bSJeremy L Thompson #include <ceed/ceed.h>
10437930d1SJeremy L Thompson #include <ceed/jit-tools.h>
110d0321e0SJeremy L Thompson #include <cuda.h>
120d0321e0SJeremy L Thompson #include <cuda_runtime.h>
13*2b730f8bSJeremy L Thompson 
140d0321e0SJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
15*2b730f8bSJeremy L Thompson #include "ceed-cuda-ref.h"
160d0321e0SJeremy L Thompson 
170d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
180d0321e0SJeremy L Thompson // Basis apply - tensor
190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
20*2b730f8bSJeremy L Thompson int CeedBasisApply_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
210d0321e0SJeremy L Thompson   Ceed ceed;
22*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
230d0321e0SJeremy L Thompson   Ceed_Cuda *ceed_Cuda;
24*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
250d0321e0SJeremy L Thompson   CeedBasis_Cuda *data;
26*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
27437930d1SJeremy L Thompson   const CeedInt transpose      = t_mode == CEED_TRANSPOSE;
28437930d1SJeremy L Thompson   const int     max_block_size = 32;
290d0321e0SJeremy L Thompson 
300d0321e0SJeremy L Thompson   // Read vectors
310d0321e0SJeremy L Thompson   const CeedScalar *d_u;
320d0321e0SJeremy L Thompson   CeedScalar       *d_v;
33437930d1SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
34*2b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
350d0321e0SJeremy L Thompson   }
36*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
370d0321e0SJeremy L Thompson 
380d0321e0SJeremy L Thompson   // Clear v for transpose operation
39437930d1SJeremy L Thompson   if (t_mode == CEED_TRANSPOSE) {
401f9221feSJeremy L Thompson     CeedSize length;
41*2b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(v, &length));
42*2b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
430d0321e0SJeremy L Thompson   }
44437930d1SJeremy L Thompson   CeedInt Q_1d, dim;
45*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
46*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
470d0321e0SJeremy L Thompson 
480d0321e0SJeremy L Thompson   // Basis action
49437930d1SJeremy L Thompson   switch (eval_mode) {
500d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP: {
51*2b730f8bSJeremy L Thompson       void   *interp_args[] = {(void *)&num_elem, (void *)&transpose, &data->d_interp_1d, &d_u, &d_v};
52437930d1SJeremy L Thompson       CeedInt block_size    = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
530d0321e0SJeremy L Thompson 
54*2b730f8bSJeremy L Thompson       CeedCallBackend(CeedRunKernelCuda(ceed, data->Interp, num_elem, block_size, interp_args));
550d0321e0SJeremy L Thompson     } break;
560d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD: {
57*2b730f8bSJeremy L Thompson       void   *grad_args[] = {(void *)&num_elem, (void *)&transpose, &data->d_interp_1d, &data->d_grad_1d, &d_u, &d_v};
58437930d1SJeremy L Thompson       CeedInt block_size  = max_block_size;
590d0321e0SJeremy L Thompson 
60*2b730f8bSJeremy L Thompson       CeedCallBackend(CeedRunKernelCuda(ceed, data->Grad, num_elem, block_size, grad_args));
610d0321e0SJeremy L Thompson     } break;
620d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
63437930d1SJeremy L Thompson       void     *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
64437930d1SJeremy L Thompson       const int grid_size     = num_elem;
65*2b730f8bSJeremy L Thompson       CeedCallBackend(CeedRunKernelDimCuda(ceed, data->Weight, grid_size, Q_1d, dim >= 2 ? Q_1d : 1, 1, weight_args));
660d0321e0SJeremy L Thompson     } break;
670d0321e0SJeremy L Thompson     // LCOV_EXCL_START
680d0321e0SJeremy L Thompson     // Evaluate the divergence to/from the quadrature points
690d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
700d0321e0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_DIV not supported");
710d0321e0SJeremy L Thompson     // Evaluate the curl to/from the quadrature points
720d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
730d0321e0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_CURL not supported");
740d0321e0SJeremy L Thompson     // Take no action, BasisApply should not have been called
750d0321e0SJeremy L Thompson     case CEED_EVAL_NONE:
76*2b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_NONE does not make sense in this context");
770d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
780d0321e0SJeremy L Thompson   }
790d0321e0SJeremy L Thompson 
800d0321e0SJeremy L Thompson   // Restore vectors
81437930d1SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
82*2b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
830d0321e0SJeremy L Thompson   }
84*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
850d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
860d0321e0SJeremy L Thompson }
870d0321e0SJeremy L Thompson 
880d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
890d0321e0SJeremy L Thompson // Basis apply - non-tensor
900d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
91*2b730f8bSJeremy L Thompson int CeedBasisApplyNonTensor_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
92*2b730f8bSJeremy L Thompson                                  CeedVector v) {
930d0321e0SJeremy L Thompson   Ceed ceed;
94*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
950d0321e0SJeremy L Thompson   Ceed_Cuda *ceed_Cuda;
96*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
970d0321e0SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
98*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
99437930d1SJeremy L Thompson   CeedInt num_nodes, num_qpts;
100*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
101*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes));
102437930d1SJeremy L Thompson   const CeedInt transpose       = t_mode == CEED_TRANSPOSE;
103437930d1SJeremy L Thompson   int           elems_per_block = 1;
104*2b730f8bSJeremy L Thompson   int           grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
1050d0321e0SJeremy L Thompson 
1060d0321e0SJeremy L Thompson   // Read vectors
1070d0321e0SJeremy L Thompson   const CeedScalar *d_u;
1080d0321e0SJeremy L Thompson   CeedScalar       *d_v;
109437930d1SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
110*2b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
1110d0321e0SJeremy L Thompson   }
112*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
1130d0321e0SJeremy L Thompson 
1140d0321e0SJeremy L Thompson   // Clear v for transpose operation
115437930d1SJeremy L Thompson   if (t_mode == CEED_TRANSPOSE) {
1161f9221feSJeremy L Thompson     CeedSize length;
117*2b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(v, &length));
118*2b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
1190d0321e0SJeremy L Thompson   }
1200d0321e0SJeremy L Thompson 
1210d0321e0SJeremy L Thompson   // Apply basis operation
122437930d1SJeremy L Thompson   switch (eval_mode) {
1230d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP: {
124*2b730f8bSJeremy L Thompson       void *interp_args[] = {(void *)&num_elem, (void *)&transpose, &data->d_interp, &d_u, &d_v};
125437930d1SJeremy L Thompson       if (transpose) {
126*2b730f8bSJeremy L Thompson         CeedCallBackend(CeedRunKernelDimCuda(ceed, data->Interp, grid, num_nodes, 1, elems_per_block, interp_args));
1270d0321e0SJeremy L Thompson       } else {
128*2b730f8bSJeremy L Thompson         CeedCallBackend(CeedRunKernelDimCuda(ceed, data->Interp, grid, num_qpts, 1, elems_per_block, interp_args));
1290d0321e0SJeremy L Thompson       }
1300d0321e0SJeremy L Thompson     } break;
1310d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD: {
132*2b730f8bSJeremy L Thompson       void *grad_args[] = {(void *)&num_elem, (void *)&transpose, &data->d_grad, &d_u, &d_v};
133437930d1SJeremy L Thompson       if (transpose) {
134*2b730f8bSJeremy L Thompson         CeedCallBackend(CeedRunKernelDimCuda(ceed, data->Grad, grid, num_nodes, 1, elems_per_block, grad_args));
1350d0321e0SJeremy L Thompson       } else {
136*2b730f8bSJeremy L Thompson         CeedCallBackend(CeedRunKernelDimCuda(ceed, data->Grad, grid, num_qpts, 1, elems_per_block, grad_args));
1370d0321e0SJeremy L Thompson       }
1380d0321e0SJeremy L Thompson     } break;
1390d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
140437930d1SJeremy L Thompson       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight, &d_v};
141*2b730f8bSJeremy L Thompson       CeedCallBackend(CeedRunKernelDimCuda(ceed, data->Weight, grid, num_qpts, 1, elems_per_block, weight_args));
1420d0321e0SJeremy L Thompson     } break;
1430d0321e0SJeremy L Thompson     // LCOV_EXCL_START
1440d0321e0SJeremy L Thompson     // Evaluate the divergence to/from the quadrature points
1450d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
1460d0321e0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_DIV not supported");
1470d0321e0SJeremy L Thompson     // Evaluate the curl to/from the quadrature points
1480d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
1490d0321e0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_CURL not supported");
1500d0321e0SJeremy L Thompson     // Take no action, BasisApply should not have been called
1510d0321e0SJeremy L Thompson     case CEED_EVAL_NONE:
152*2b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_NONE does not make sense in this context");
1530d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
1540d0321e0SJeremy L Thompson   }
1550d0321e0SJeremy L Thompson 
1560d0321e0SJeremy L Thompson   // Restore vectors
157437930d1SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
158*2b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
1590d0321e0SJeremy L Thompson   }
160*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
1610d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1620d0321e0SJeremy L Thompson }
1630d0321e0SJeremy L Thompson 
1640d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1650d0321e0SJeremy L Thompson // Destroy tensor basis
1660d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1670d0321e0SJeremy L Thompson static int CeedBasisDestroy_Cuda(CeedBasis basis) {
1680d0321e0SJeremy L Thompson   Ceed ceed;
169*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
1700d0321e0SJeremy L Thompson 
1710d0321e0SJeremy L Thompson   CeedBasis_Cuda *data;
172*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
1730d0321e0SJeremy L Thompson 
174*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
1750d0321e0SJeremy L Thompson 
176*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d));
177*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp_1d));
178*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad_1d));
179*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
180437930d1SJeremy L Thompson 
1810d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1820d0321e0SJeremy L Thompson }
1830d0321e0SJeremy L Thompson 
1840d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1850d0321e0SJeremy L Thompson // Destroy non-tensor basis
1860d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1870d0321e0SJeremy L Thompson static int CeedBasisDestroyNonTensor_Cuda(CeedBasis basis) {
1880d0321e0SJeremy L Thompson   Ceed ceed;
189*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
1900d0321e0SJeremy L Thompson 
1910d0321e0SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
192*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
1930d0321e0SJeremy L Thompson 
194*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
1950d0321e0SJeremy L Thompson 
196*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_q_weight));
197*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp));
198*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad));
199*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
200437930d1SJeremy L Thompson 
2010d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2020d0321e0SJeremy L Thompson }
2030d0321e0SJeremy L Thompson 
2040d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2050d0321e0SJeremy L Thompson // Create tensor
2060d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
207*2b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
208*2b730f8bSJeremy L Thompson                                  const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
2090d0321e0SJeremy L Thompson   Ceed ceed;
210*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
2110d0321e0SJeremy L Thompson   CeedBasis_Cuda *data;
212*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
2130d0321e0SJeremy L Thompson 
2140d0321e0SJeremy L Thompson   // Copy data to GPU
215437930d1SJeremy L Thompson   const CeedInt q_bytes = Q_1d * sizeof(CeedScalar);
216*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
217*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice));
2180d0321e0SJeremy L Thompson 
219437930d1SJeremy L Thompson   const CeedInt interp_bytes = q_bytes * P_1d;
220*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
221*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice));
2220d0321e0SJeremy L Thompson 
223*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes));
224*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice));
2250d0321e0SJeremy L Thompson 
2260d0321e0SJeremy L Thompson   // Complie basis kernels
227437930d1SJeremy L Thompson   CeedInt num_comp;
228*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
229437930d1SJeremy L Thompson   char *basis_kernel_path, *basis_kernel_source;
230*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-tensor.h", &basis_kernel_path));
23146dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading Basis Kernel Source -----\n");
232*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
23346dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading Basis Kernel Source Complete! -----\n");
234*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedCompileCuda(ceed, basis_kernel_source, &data->module, 7, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "BASIS_BUF_LEN",
235*2b730f8bSJeremy L Thompson                                   num_comp * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp,
236*2b730f8bSJeremy L Thompson                                   "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim)));
237*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "Interp", &data->Interp));
238*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "Grad", &data->Grad));
239*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "Weight", &data->Weight));
240*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_path));
241*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_source));
242437930d1SJeremy L Thompson 
243*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
2440d0321e0SJeremy L Thompson 
245*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApply_Cuda));
246*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda));
2470d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2480d0321e0SJeremy L Thompson }
2490d0321e0SJeremy L Thompson 
2500d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2510d0321e0SJeremy L Thompson // Create non-tensor
2520d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
253*2b730f8bSJeremy L Thompson int CeedBasisCreateH1_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *grad,
254*2b730f8bSJeremy L Thompson                            const CeedScalar *qref, const CeedScalar *q_weight, CeedBasis basis) {
2550d0321e0SJeremy L Thompson   Ceed ceed;
256*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
2570d0321e0SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
258*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
2590d0321e0SJeremy L Thompson 
2600d0321e0SJeremy L Thompson   // Copy basis data to GPU
261437930d1SJeremy L Thompson   const CeedInt q_bytes = num_qpts * sizeof(CeedScalar);
262*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes));
263*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice));
2640d0321e0SJeremy L Thompson 
265437930d1SJeremy L Thompson   const CeedInt interp_bytes = q_bytes * num_nodes;
266*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes));
267*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice));
2680d0321e0SJeremy L Thompson 
269437930d1SJeremy L Thompson   const CeedInt grad_bytes = q_bytes * num_nodes * dim;
270*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad, grad_bytes));
271*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_grad, grad, grad_bytes, cudaMemcpyHostToDevice));
2720d0321e0SJeremy L Thompson 
2730d0321e0SJeremy L Thompson   // Compile basis kernels
274437930d1SJeremy L Thompson   CeedInt num_comp;
275*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
276437930d1SJeremy L Thompson   char *basis_kernel_path, *basis_kernel_source;
277*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path));
27846dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading Basis Kernel Source -----\n");
279*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
28046dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading Basis Kernel Source Complete! -----\n");
281*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, CeedCompileCuda(ceed, basis_kernel_source, &data->module, 4, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_DIM", dim,
282*2b730f8bSJeremy L Thompson                                      "BASIS_NUM_COMP", num_comp));
283*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, CeedGetKernelCuda(ceed, data->module, "Interp", &data->Interp));
284*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, CeedGetKernelCuda(ceed, data->module, "Grad", &data->Grad));
285*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, CeedGetKernelCuda(ceed, data->module, "Weight", &data->Weight));
286*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_path));
287*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_source));
2880d0321e0SJeremy L Thompson 
289*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
2900d0321e0SJeremy L Thompson 
2910d0321e0SJeremy L Thompson   // Register backend functions
292*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda));
293*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda));
2940d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2950d0321e0SJeremy L Thompson }
2960d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
297