xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-shared/ceed-cuda-shared-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.
3c532df63SYohann //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5c532df63SYohann //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7c532df63SYohann 
8ec3da8bcSJed Brown #include <ceed/backend.h>
9*2b730f8bSJeremy L Thompson #include <ceed/ceed.h>
10437930d1SJeremy L Thompson #include <ceed/jit-tools.h>
113d576824SJeremy L Thompson #include <cuda.h>
123d576824SJeremy L Thompson #include <cuda_runtime.h>
133d576824SJeremy L Thompson #include <stddef.h>
14c532df63SYohann 
15*2b730f8bSJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
16*2b730f8bSJeremy L Thompson #include "ceed-cuda-shared.h"
17c532df63SYohann 
18ab213215SJeremy L Thompson //------------------------------------------------------------------------------
19ab213215SJeremy L Thompson // Device initalization
20ab213215SJeremy L Thompson //------------------------------------------------------------------------------
21*2b730f8bSJeremy L Thompson int CeedCudaInitInterp(CeedScalar *d_B, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B);
22*2b730f8bSJeremy L Thompson int CeedCudaInitGrad(CeedScalar *d_B, CeedScalar *d_G, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B_ptr, CeedScalar **c_G_ptr);
23*2b730f8bSJeremy L Thompson int CeedCudaInitCollocatedGrad(CeedScalar *d_B, CeedScalar *d_G, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B_ptr, CeedScalar **c_G_ptr);
24c532df63SYohann 
25ab213215SJeremy L Thompson //------------------------------------------------------------------------------
26ab213215SJeremy L Thompson // Apply basis
27ab213215SJeremy L Thompson //------------------------------------------------------------------------------
28*2b730f8bSJeremy L Thompson int CeedBasisApplyTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
297f823360Sjeremylt                                      CeedVector v) {
30c532df63SYohann   Ceed ceed;
31*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
326dbfb411Snbeams   Ceed_Cuda *ceed_Cuda;
33*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
34c532df63SYohann   CeedBasis_Cuda_shared *data;
35*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
36437930d1SJeremy L Thompson   CeedInt dim, num_comp;
37*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
38*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
39c532df63SYohann 
40ab213215SJeremy L Thompson   // Read vectors
41c532df63SYohann   const CeedScalar *d_u;
42c532df63SYohann   CeedScalar       *d_v;
43437930d1SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
44*2b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
45c532df63SYohann   }
46*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
47c532df63SYohann 
48ab213215SJeremy L Thompson   // Apply basis operation
49437930d1SJeremy L Thompson   switch (eval_mode) {
50ab213215SJeremy L Thompson     case CEED_EVAL_INTERP: {
51437930d1SJeremy L Thompson       CeedInt P_1d, Q_1d;
52*2b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
53*2b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
54437930d1SJeremy L Thompson       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
55*2b730f8bSJeremy L Thompson       CeedCallBackend(CeedCudaInitInterp(data->d_interp_1d, P_1d, Q_1d, &data->c_B));
56*2b730f8bSJeremy L Thompson       void *interp_args[] = {(void *)&num_elem, &data->c_B, &d_u, &d_v};
574d537eeaSYohann       if (dim == 1) {
58*2b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d,
59e6f67ff7SJed Brown                                                                                                  1));  // avoid >512 total threads
60*2b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
61437930d1SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
629e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
63*2b730f8bSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->InterpTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
649e201c85SYohann         } else {
65*2b730f8bSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Interp, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
669e201c85SYohann         }
67074be161SYohann Dudouit       } else if (dim == 2) {
68437930d1SJeremy L Thompson         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
69437930d1SJeremy L Thompson         // elems_per_block must be at least 1
70*2b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
71*2b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
72*2b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
739e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
74*2b730f8bSJeremy L Thompson           CeedCallBackend(
75*2b730f8bSJeremy L Thompson               CeedRunKernelDimSharedCuda(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
769e201c85SYohann         } else {
77*2b730f8bSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
789e201c85SYohann         }
79074be161SYohann Dudouit       } else if (dim == 3) {
80437930d1SJeremy L Thompson         CeedInt elems_per_block = 1;
81*2b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
82*2b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
839e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
84*2b730f8bSJeremy L Thompson           CeedCallBackend(
85*2b730f8bSJeremy L Thompson               CeedRunKernelDimSharedCuda(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
869e201c85SYohann         } else {
87*2b730f8bSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
88074be161SYohann Dudouit         }
899e201c85SYohann       }
90ab213215SJeremy L Thompson     } break;
91ab213215SJeremy L Thompson     case CEED_EVAL_GRAD: {
92437930d1SJeremy L Thompson       CeedInt P_1d, Q_1d;
93*2b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
94*2b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
95437930d1SJeremy L Thompson       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
969e201c85SYohann       if (data->d_collo_grad_1d) {
97*2b730f8bSJeremy L Thompson         CeedCallBackend(CeedCudaInitCollocatedGrad(data->d_interp_1d, data->d_collo_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G));
989e201c85SYohann       } else {
99*2b730f8bSJeremy L Thompson         CeedCallBackend(CeedCudaInitGrad(data->d_interp_1d, data->d_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G));
1009e201c85SYohann       }
101*2b730f8bSJeremy L Thompson       void *grad_args[] = {(void *)&num_elem, &data->c_B, &data->c_G, &d_u, &d_v};
1024d537eeaSYohann       if (dim == 1) {
103*2b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d,
104e6f67ff7SJed Brown                                                                                                  1));  // avoid >512 total threads
105*2b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
106437930d1SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
1079e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
108*2b730f8bSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->GradTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
1099e201c85SYohann         } else {
110*2b730f8bSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
1119e201c85SYohann         }
112074be161SYohann Dudouit       } else if (dim == 2) {
113437930d1SJeremy L Thompson         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
114437930d1SJeremy L Thompson         // elems_per_block must be at least 1
115*2b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
116*2b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
117*2b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
1189e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
119*2b730f8bSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
1209e201c85SYohann         } else {
121*2b730f8bSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
1229e201c85SYohann         }
123074be161SYohann Dudouit       } else if (dim == 3) {
124437930d1SJeremy L Thompson         CeedInt elems_per_block = 1;
125*2b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
126*2b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
1279e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
128*2b730f8bSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
1299e201c85SYohann         } else {
130*2b730f8bSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
1319e201c85SYohann         }
132074be161SYohann Dudouit       }
133ab213215SJeremy L Thompson     } break;
134ab213215SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
135437930d1SJeremy L Thompson       CeedInt Q_1d;
136*2b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
137437930d1SJeremy L Thompson       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
138074be161SYohann Dudouit       if (dim == 1) {
139437930d1SJeremy L Thompson         const CeedInt elems_per_block = 32 / Q_1d;
140*2b730f8bSJeremy L Thompson         const CeedInt gridsize        = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
141*2b730f8bSJeremy L Thompson         CeedCallBackend(CeedRunKernelDimCuda(ceed, data->Weight, gridsize, Q_1d, elems_per_block, 1, weight_args));
142074be161SYohann Dudouit       } else if (dim == 2) {
143437930d1SJeremy L Thompson         const CeedInt opt_elems       = 32 / (Q_1d * Q_1d);
144437930d1SJeremy L Thompson         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
145*2b730f8bSJeremy L Thompson         const CeedInt gridsize        = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
146*2b730f8bSJeremy L Thompson         CeedCallBackend(CeedRunKernelDimCuda(ceed, data->Weight, gridsize, Q_1d, Q_1d, elems_per_block, weight_args));
147074be161SYohann Dudouit       } else if (dim == 3) {
1489e201c85SYohann         const CeedInt opt_elems       = 32 / (Q_1d * Q_1d);
1499e201c85SYohann         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
150*2b730f8bSJeremy L Thompson         const CeedInt gridsize        = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
151*2b730f8bSJeremy L Thompson         CeedCallBackend(CeedRunKernelDimCuda(ceed, data->Weight, gridsize, Q_1d, Q_1d, elems_per_block, weight_args));
152074be161SYohann Dudouit       }
153ab213215SJeremy L Thompson     } break;
154ab213215SJeremy L Thompson     // LCOV_EXCL_START
155ab213215SJeremy L Thompson     // Evaluate the divergence to/from the quadrature points
156ab213215SJeremy L Thompson     case CEED_EVAL_DIV:
157e15f9bd0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_DIV not supported");
158ab213215SJeremy L Thompson     // Evaluate the curl to/from the quadrature points
159ab213215SJeremy L Thompson     case CEED_EVAL_CURL:
160e15f9bd0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_CURL not supported");
161ab213215SJeremy L Thompson     // Take no action, BasisApply should not have been called
162ab213215SJeremy L Thompson     case CEED_EVAL_NONE:
163*2b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_NONE does not make sense in this context");
164ab213215SJeremy L Thompson       // LCOV_EXCL_STOP
165c532df63SYohann   }
166c532df63SYohann 
167ab213215SJeremy L Thompson   // Restore vectors
168437930d1SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
169*2b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
170c532df63SYohann   }
171*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
172e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
173c532df63SYohann }
174c532df63SYohann 
175ab213215SJeremy L Thompson //------------------------------------------------------------------------------
176ab213215SJeremy L Thompson // Destroy basis
177ab213215SJeremy L Thompson //------------------------------------------------------------------------------
178c532df63SYohann static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) {
179c532df63SYohann   Ceed ceed;
180*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
181c532df63SYohann 
182c532df63SYohann   CeedBasis_Cuda_shared *data;
183*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
184c532df63SYohann 
185*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
186c532df63SYohann 
187*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d));
188*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp_1d));
189*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad_1d));
190*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_collo_grad_1d));
191c532df63SYohann 
192*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
193c532df63SYohann 
194e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
195c532df63SYohann }
196c532df63SYohann 
197ab213215SJeremy L Thompson //------------------------------------------------------------------------------
198ab213215SJeremy L Thompson // Create tensor basis
199ab213215SJeremy L Thompson //------------------------------------------------------------------------------
200*2b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
201*2b730f8bSJeremy L Thompson                                         const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
202c532df63SYohann   Ceed ceed;
203*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
204c532df63SYohann   CeedBasis_Cuda_shared *data;
205*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
206c532df63SYohann 
207ab213215SJeremy L Thompson   // Copy basis data to GPU
208437930d1SJeremy L Thompson   const CeedInt q_bytes = Q_1d * sizeof(CeedScalar);
209*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
210*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice));
211c532df63SYohann 
212437930d1SJeremy L Thompson   const CeedInt interp_bytes = q_bytes * P_1d;
213*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
214*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice));
215c532df63SYohann 
216*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes));
217*2b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice));
218c532df63SYohann 
219ab213215SJeremy L Thompson   // Compute collocated gradient and copy to GPU
220437930d1SJeremy L Thompson   data->d_collo_grad_1d    = NULL;
2219e201c85SYohann   bool has_collocated_grad = dim == 3 && Q_1d >= P_1d;
2229e201c85SYohann   if (has_collocated_grad) {
223437930d1SJeremy L Thompson     CeedScalar *collo_grad_1d;
224*2b730f8bSJeremy L Thompson     CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d));
225*2b730f8bSJeremy L Thompson     CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d));
226*2b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d));
227*2b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, cudaMemcpyHostToDevice));
228*2b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&collo_grad_1d));
229ac421f39SYohann   }
230ac421f39SYohann 
231ab213215SJeremy L Thompson   // Compile basis kernels
232437930d1SJeremy L Thompson   CeedInt num_comp;
233*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
234437930d1SJeremy L Thompson   char *basis_kernel_path, *basis_kernel_source;
235*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-shared-basis-tensor.h", &basis_kernel_path));
23646dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading Basis Kernel Source -----\n");
237*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
23846dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading Basis Kernel Source Complete -----\n");
239*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedCompileCuda(ceed, basis_kernel_source, &data->module, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D", CeedIntMax(Q_1d, P_1d),
240*2b730f8bSJeremy L Thompson                                   "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS",
241*2b730f8bSJeremy L Thompson                                   CeedIntPow(Q_1d, dim), "BASIS_HAS_COLLOCATED_GRAD", has_collocated_grad));
242*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "Interp", &data->Interp));
243*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
244*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "Grad", &data->Grad));
245*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "GradTranspose", &data->GradTranspose));
246*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "Weight", &data->Weight));
247*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_path));
248*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_source));
249c532df63SYohann 
250*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
251ab213215SJeremy L Thompson 
252ab213215SJeremy L Thompson   // Register backend functions
253*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared));
254*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared));
255e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
256c532df63SYohann }
257ab213215SJeremy L Thompson //------------------------------------------------------------------------------
258