xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-ref/ceed-cuda-ref-basis.c (revision ca7355308a14805110a7d98dabf1f658d5ec16d5)
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 
849aac155SJeremy L Thompson #include <ceed.h>
90d0321e0SJeremy L Thompson #include <ceed/backend.h>
10437930d1SJeremy L Thompson #include <ceed/jit-tools.h>
110d0321e0SJeremy L Thompson #include <cuda.h>
120d0321e0SJeremy L Thompson #include <cuda_runtime.h>
132b730f8bSJeremy L Thompson 
1449aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
150d0321e0SJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
162b730f8bSJeremy L Thompson #include "ceed-cuda-ref.h"
170d0321e0SJeremy L Thompson 
180d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
190d0321e0SJeremy L Thompson // Basis apply - tensor
200d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
212b730f8bSJeremy L Thompson int CeedBasisApply_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
220d0321e0SJeremy L Thompson   Ceed              ceed;
230d0321e0SJeremy L Thompson   Ceed_Cuda        *ceed_Cuda;
24*ca735530SJeremy L Thompson   CeedInt           Q_1d, dim;
25437930d1SJeremy L Thompson   const CeedInt     transpose      = t_mode == CEED_TRANSPOSE;
26437930d1SJeremy L Thompson   const int         max_block_size = 32;
270d0321e0SJeremy L Thompson   const CeedScalar *d_u;
280d0321e0SJeremy L Thompson   CeedScalar       *d_v;
29*ca735530SJeremy L Thompson   CeedBasis_Cuda   *data;
30*ca735530SJeremy L Thompson 
31*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
32*ca735530SJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
33*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
34*ca735530SJeremy L Thompson 
35*ca735530SJeremy L Thompson   // Read vectors
366574a04fSJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
376574a04fSJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
382b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
390d0321e0SJeremy L Thompson 
400d0321e0SJeremy L Thompson   // Clear v for transpose operation
41437930d1SJeremy L Thompson   if (t_mode == CEED_TRANSPOSE) {
421f9221feSJeremy L Thompson     CeedSize length;
43*ca735530SJeremy L Thompson 
442b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(v, &length));
452b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
460d0321e0SJeremy L Thompson   }
472b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
482b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
490d0321e0SJeremy L Thompson 
500d0321e0SJeremy L Thompson   // Basis action
51437930d1SJeremy L Thompson   switch (eval_mode) {
520d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP: {
532b730f8bSJeremy L Thompson       void         *interp_args[] = {(void *)&num_elem, (void *)&transpose, &data->d_interp_1d, &d_u, &d_v};
54b2165e7aSSebastian Grimberg       const CeedInt block_size    = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
550d0321e0SJeremy L Thompson 
56eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->Interp, num_elem, block_size, interp_args));
570d0321e0SJeremy L Thompson     } break;
580d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD: {
592b730f8bSJeremy L Thompson       void         *grad_args[] = {(void *)&num_elem, (void *)&transpose, &data->d_interp_1d, &data->d_grad_1d, &d_u, &d_v};
60b2165e7aSSebastian Grimberg       const CeedInt block_size  = max_block_size;
610d0321e0SJeremy L Thompson 
62eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->Grad, num_elem, block_size, grad_args));
630d0321e0SJeremy L Thompson     } break;
640d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
65437930d1SJeremy L Thompson       void     *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
66b2165e7aSSebastian Grimberg       const int block_size_x  = Q_1d;
67b2165e7aSSebastian Grimberg       const int block_size_y  = dim >= 2 ? Q_1d : 1;
68b2165e7aSSebastian Grimberg 
69b2165e7aSSebastian Grimberg       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, num_elem, block_size_x, block_size_y, 1, weight_args));
700d0321e0SJeremy L Thompson     } break;
710d0321e0SJeremy L Thompson     // LCOV_EXCL_START
720d0321e0SJeremy L Thompson     // Evaluate the divergence to/from the quadrature points
730d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
740d0321e0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_DIV not supported");
750d0321e0SJeremy L Thompson     // Evaluate the curl to/from the quadrature points
760d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
770d0321e0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_CURL not supported");
780d0321e0SJeremy L Thompson     // Take no action, BasisApply should not have been called
790d0321e0SJeremy L Thompson     case CEED_EVAL_NONE:
802b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_NONE does not make sense in this context");
810d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
820d0321e0SJeremy L Thompson   }
830d0321e0SJeremy L Thompson 
840d0321e0SJeremy L Thompson   // Restore vectors
85437930d1SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
862b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
870d0321e0SJeremy L Thompson   }
882b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
890d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
900d0321e0SJeremy L Thompson }
910d0321e0SJeremy L Thompson 
920d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
930d0321e0SJeremy L Thompson // Basis apply - non-tensor
940d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
952b730f8bSJeremy L Thompson int CeedBasisApplyNonTensor_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
962b730f8bSJeremy L Thompson                                  CeedVector v) {
970d0321e0SJeremy L Thompson   Ceed                     ceed;
980d0321e0SJeremy L Thompson   Ceed_Cuda               *ceed_Cuda;
99437930d1SJeremy L Thompson   CeedInt                  num_nodes, num_qpts;
100437930d1SJeremy L Thompson   const CeedInt            transpose       = t_mode == CEED_TRANSPOSE;
101437930d1SJeremy L Thompson   int                      elems_per_block = 1;
1022b730f8bSJeremy L Thompson   int                      grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
1030d0321e0SJeremy L Thompson   const CeedScalar        *d_u;
1040d0321e0SJeremy L Thompson   CeedScalar              *d_v;
105*ca735530SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
106*ca735530SJeremy L Thompson 
107*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
108*ca735530SJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
109*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
110*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
111*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes));
112*ca735530SJeremy L Thompson 
113*ca735530SJeremy L Thompson   // Read vectors
114437930d1SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
1152b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
1160d0321e0SJeremy L Thompson   }
1172b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
1180d0321e0SJeremy L Thompson 
1190d0321e0SJeremy L Thompson   // Clear v for transpose operation
120437930d1SJeremy L Thompson   if (t_mode == CEED_TRANSPOSE) {
1211f9221feSJeremy L Thompson     CeedSize length;
122*ca735530SJeremy L Thompson 
1232b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(v, &length));
1242b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
1250d0321e0SJeremy L Thompson   }
1260d0321e0SJeremy L Thompson 
1270d0321e0SJeremy L Thompson   // Apply basis operation
128437930d1SJeremy L Thompson   switch (eval_mode) {
1290d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP: {
1302b730f8bSJeremy L Thompson       void     *interp_args[] = {(void *)&num_elem, (void *)&transpose, &data->d_interp, &d_u, &d_v};
131b2165e7aSSebastian Grimberg       const int block_size_x  = transpose ? num_nodes : num_qpts;
132b2165e7aSSebastian Grimberg 
133b2165e7aSSebastian Grimberg       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Interp, grid, block_size_x, 1, elems_per_block, interp_args));
1340d0321e0SJeremy L Thompson     } break;
1350d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD: {
1362b730f8bSJeremy L Thompson       void     *grad_args[]  = {(void *)&num_elem, (void *)&transpose, &data->d_grad, &d_u, &d_v};
137b2165e7aSSebastian Grimberg       const int block_size_x = transpose ? num_nodes : num_qpts;
138b2165e7aSSebastian Grimberg 
139b2165e7aSSebastian Grimberg       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Grad, grid, block_size_x, 1, elems_per_block, grad_args));
1400d0321e0SJeremy L Thompson     } break;
1410d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
142437930d1SJeremy L Thompson       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight, &d_v};
143b2165e7aSSebastian Grimberg 
144eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, num_qpts, 1, elems_per_block, weight_args));
1450d0321e0SJeremy L Thompson     } break;
1460d0321e0SJeremy L Thompson     // LCOV_EXCL_START
1470d0321e0SJeremy L Thompson     // Evaluate the divergence to/from the quadrature points
1480d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
1490d0321e0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_DIV not supported");
1500d0321e0SJeremy L Thompson     // Evaluate the curl to/from the quadrature points
1510d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
1520d0321e0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_CURL not supported");
1530d0321e0SJeremy L Thompson     // Take no action, BasisApply should not have been called
1540d0321e0SJeremy L Thompson     case CEED_EVAL_NONE:
1552b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_NONE does not make sense in this context");
1560d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
1570d0321e0SJeremy L Thompson   }
1580d0321e0SJeremy L Thompson 
1590d0321e0SJeremy L Thompson   // Restore vectors
160437930d1SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
1612b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
1620d0321e0SJeremy L Thompson   }
1632b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
1640d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1650d0321e0SJeremy L Thompson }
1660d0321e0SJeremy L Thompson 
1670d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1680d0321e0SJeremy L Thompson // Destroy tensor basis
1690d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1700d0321e0SJeremy L Thompson static int CeedBasisDestroy_Cuda(CeedBasis basis) {
1710d0321e0SJeremy L Thompson   Ceed            ceed;
1720d0321e0SJeremy L Thompson   CeedBasis_Cuda *data;
173*ca735530SJeremy L Thompson 
174*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
1752b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
1762b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
1772b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d));
1782b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp_1d));
1792b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad_1d));
1802b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
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;
1890d0321e0SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
190*ca735530SJeremy L Thompson 
191*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
1922b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
1932b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
1942b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_q_weight));
1952b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp));
1962b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad));
1972b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
1980d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1990d0321e0SJeremy L Thompson }
2000d0321e0SJeremy L Thompson 
2010d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2020d0321e0SJeremy L Thompson // Create tensor
2030d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2042b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
2052b730f8bSJeremy L Thompson                                  const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
2060d0321e0SJeremy L Thompson   Ceed            ceed;
207*ca735530SJeremy L Thompson   char           *basis_kernel_path, *basis_kernel_source;
208*ca735530SJeremy L Thompson   CeedInt         num_comp;
209*ca735530SJeremy L Thompson   const CeedInt   q_bytes      = Q_1d * sizeof(CeedScalar);
210*ca735530SJeremy L Thompson   const CeedInt   interp_bytes = q_bytes * P_1d;
2110d0321e0SJeremy L Thompson   CeedBasis_Cuda *data;
212*ca735530SJeremy L Thompson 
213*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
2142b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
2150d0321e0SJeremy L Thompson 
2160d0321e0SJeremy L Thompson   // Copy data to GPU
2172b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
2182b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice));
2192b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
2202b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice));
2212b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes));
2222b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice));
2230d0321e0SJeremy L Thompson 
224ecc88aebSJeremy L Thompson   // Compile basis kernels
2252b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
2262b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-tensor.h", &basis_kernel_path));
22723d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
2282b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
22923d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
230eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 7, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "BASIS_BUF_LEN",
2312b730f8bSJeremy L Thompson                                    num_comp * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp,
2322b730f8bSJeremy L Thompson                                    "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim)));
233eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
234eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
235eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
2362b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_path));
2372b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_source));
238437930d1SJeremy L Thompson 
2392b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
2400d0321e0SJeremy L Thompson 
2412b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApply_Cuda));
2422b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda));
2430d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2440d0321e0SJeremy L Thompson }
2450d0321e0SJeremy L Thompson 
2460d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2470d0321e0SJeremy L Thompson // Create non-tensor
2480d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2492b730f8bSJeremy L Thompson int CeedBasisCreateH1_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *grad,
25051475c7cSJeremy L Thompson                            const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
2510d0321e0SJeremy L Thompson   Ceed                     ceed;
252*ca735530SJeremy L Thompson   char                    *basis_kernel_path, *basis_kernel_source;
253*ca735530SJeremy L Thompson   CeedInt                  num_comp;
254*ca735530SJeremy L Thompson   const CeedInt            q_bytes      = num_qpts * sizeof(CeedScalar);
255*ca735530SJeremy L Thompson   const CeedInt            interp_bytes = q_bytes * num_nodes;
256*ca735530SJeremy L Thompson   const CeedInt            grad_bytes   = q_bytes * num_nodes * dim;
2570d0321e0SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
258*ca735530SJeremy L Thompson 
259*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
2602b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
2610d0321e0SJeremy L Thompson 
2620d0321e0SJeremy L Thompson   // Copy basis data to GPU
2632b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes));
2642b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice));
2652b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes));
2662b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice));
2672b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad, grad_bytes));
2682b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_grad, grad, grad_bytes, cudaMemcpyHostToDevice));
2690d0321e0SJeremy L Thompson 
2700d0321e0SJeremy L Thompson   // Compile basis kernels
2712b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
2722b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path));
27323d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
2742b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
27523d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
276eb7e6cafSJeremy L Thompson   CeedCallCuda(ceed, CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 4, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_DIM", dim,
2772b730f8bSJeremy L Thompson                                       "BASIS_NUM_COMP", num_comp));
278eb7e6cafSJeremy L Thompson   CeedCallCuda(ceed, CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
279eb7e6cafSJeremy L Thompson   CeedCallCuda(ceed, CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
280eb7e6cafSJeremy L Thompson   CeedCallCuda(ceed, CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
2812b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_path));
2822b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_source));
2830d0321e0SJeremy L Thompson 
2842b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
2850d0321e0SJeremy L Thompson 
2860d0321e0SJeremy L Thompson   // Register backend functions
2872b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda));
2882b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda));
2890d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2900d0321e0SJeremy L Thompson }
2912a86cc9dSSebastian Grimberg 
2920d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
293