xref: /libCEED/backends/cuda-shared/ceed-cuda-shared-basis.c (revision b2165e7a2e371018feedcb47974a027ed47e0487)
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 
849aac155SJeremy L Thompson #include <ceed.h>
9ec3da8bcSJed Brown #include <ceed/backend.h>
10437930d1SJeremy L Thompson #include <ceed/jit-tools.h>
113d576824SJeremy L Thompson #include <cuda.h>
123d576824SJeremy L Thompson #include <cuda_runtime.h>
1349aac155SJeremy L Thompson #include <stdbool.h>
143d576824SJeremy L Thompson #include <stddef.h>
15c532df63SYohann 
1649aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
172b730f8bSJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
182b730f8bSJeremy L Thompson #include "ceed-cuda-shared.h"
19c532df63SYohann 
20ab213215SJeremy L Thompson //------------------------------------------------------------------------------
21ab213215SJeremy L Thompson // Device initalization
22ab213215SJeremy L Thompson //------------------------------------------------------------------------------
23eb7e6cafSJeremy L Thompson int CeedInit_CudaInterp(CeedScalar *d_B, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B);
24eb7e6cafSJeremy L Thompson int CeedInit_CudaGrad(CeedScalar *d_B, CeedScalar *d_G, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B_ptr, CeedScalar **c_G_ptr);
25eb7e6cafSJeremy L Thompson int CeedInit_CudaCollocatedGrad(CeedScalar *d_B, CeedScalar *d_G, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B_ptr, CeedScalar **c_G_ptr);
26c532df63SYohann 
27ab213215SJeremy L Thompson //------------------------------------------------------------------------------
28ab213215SJeremy L Thompson // Apply basis
29ab213215SJeremy L Thompson //------------------------------------------------------------------------------
302b730f8bSJeremy L Thompson int CeedBasisApplyTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
317f823360Sjeremylt                                      CeedVector v) {
32c532df63SYohann   Ceed ceed;
332b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
346dbfb411Snbeams   Ceed_Cuda *ceed_Cuda;
352b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
36c532df63SYohann   CeedBasis_Cuda_shared *data;
372b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
38437930d1SJeremy L Thompson   CeedInt dim, num_comp;
392b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
402b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
41c532df63SYohann 
42ab213215SJeremy L Thompson   // Read vectors
43c532df63SYohann   const CeedScalar *d_u;
44c532df63SYohann   CeedScalar       *d_v;
456574a04fSJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
466574a04fSJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
472b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
48c532df63SYohann 
49ab213215SJeremy L Thompson   // Apply basis operation
50437930d1SJeremy L Thompson   switch (eval_mode) {
51ab213215SJeremy L Thompson     case CEED_EVAL_INTERP: {
52437930d1SJeremy L Thompson       CeedInt P_1d, Q_1d;
532b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
542b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
55437930d1SJeremy L Thompson       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
56eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedInit_CudaInterp(data->d_interp_1d, P_1d, Q_1d, &data->c_B));
572b730f8bSJeremy L Thompson       void *interp_args[] = {(void *)&num_elem, &data->c_B, &d_u, &d_v};
584d537eeaSYohann       if (dim == 1) {
592b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d,
60e6f67ff7SJed Brown                                                                                                  1));  // avoid >512 total threads
612b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
62437930d1SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
63*b2165e7aSSebastian Grimberg 
649e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
65eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->InterpTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
669e201c85SYohann         } else {
67eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
689e201c85SYohann         }
69074be161SYohann Dudouit       } else if (dim == 2) {
70437930d1SJeremy L Thompson         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
71437930d1SJeremy L Thompson         // elems_per_block must be at least 1
722b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
732b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
742b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
75*b2165e7aSSebastian Grimberg 
769e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
772b730f8bSJeremy L Thompson           CeedCallBackend(
78eb7e6cafSJeremy L Thompson               CeedRunKernelDimShared_Cuda(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
799e201c85SYohann         } else {
80eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
819e201c85SYohann         }
82074be161SYohann Dudouit       } else if (dim == 3) {
83437930d1SJeremy L Thompson         CeedInt elems_per_block = 1;
842b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
852b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
86*b2165e7aSSebastian Grimberg 
879e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
882b730f8bSJeremy L Thompson           CeedCallBackend(
89eb7e6cafSJeremy L Thompson               CeedRunKernelDimShared_Cuda(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
909e201c85SYohann         } else {
91eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
92074be161SYohann Dudouit         }
939e201c85SYohann       }
94ab213215SJeremy L Thompson     } break;
95ab213215SJeremy L Thompson     case CEED_EVAL_GRAD: {
96437930d1SJeremy L Thompson       CeedInt P_1d, Q_1d;
972b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
982b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
99437930d1SJeremy L Thompson       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
1009e201c85SYohann       if (data->d_collo_grad_1d) {
101eb7e6cafSJeremy L Thompson         CeedCallBackend(CeedInit_CudaCollocatedGrad(data->d_interp_1d, data->d_collo_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G));
1029e201c85SYohann       } else {
103eb7e6cafSJeremy L Thompson         CeedCallBackend(CeedInit_CudaGrad(data->d_interp_1d, data->d_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G));
1049e201c85SYohann       }
1052b730f8bSJeremy L Thompson       void *grad_args[] = {(void *)&num_elem, &data->c_B, &data->c_G, &d_u, &d_v};
1064d537eeaSYohann       if (dim == 1) {
1072b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d,
108e6f67ff7SJed Brown                                                                                                  1));  // avoid >512 total threads
1092b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
110437930d1SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
111*b2165e7aSSebastian Grimberg 
1129e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
113eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->GradTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
1149e201c85SYohann         } else {
115eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
1169e201c85SYohann         }
117074be161SYohann Dudouit       } else if (dim == 2) {
118437930d1SJeremy L Thompson         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
119437930d1SJeremy L Thompson         // elems_per_block must be at least 1
1202b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
1212b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
1222b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
123*b2165e7aSSebastian Grimberg 
1249e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
125eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
1269e201c85SYohann         } else {
127eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
1289e201c85SYohann         }
129074be161SYohann Dudouit       } else if (dim == 3) {
130437930d1SJeremy L Thompson         CeedInt elems_per_block = 1;
1312b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
1322b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
133*b2165e7aSSebastian Grimberg 
1349e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
135eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
1369e201c85SYohann         } else {
137eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
1389e201c85SYohann         }
139074be161SYohann Dudouit       }
140ab213215SJeremy L Thompson     } break;
141ab213215SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
142437930d1SJeremy L Thompson       CeedInt Q_1d;
1432b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
144437930d1SJeremy L Thompson       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
145074be161SYohann Dudouit       if (dim == 1) {
146437930d1SJeremy L Thompson         const CeedInt elems_per_block = 32 / Q_1d;
147*b2165e7aSSebastian Grimberg         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
148*b2165e7aSSebastian Grimberg 
149*b2165e7aSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, elems_per_block, 1, weight_args));
150074be161SYohann Dudouit       } else if (dim == 2) {
151437930d1SJeremy L Thompson         const CeedInt opt_elems       = 32 / (Q_1d * Q_1d);
152437930d1SJeremy L Thompson         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
153*b2165e7aSSebastian Grimberg         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
154*b2165e7aSSebastian Grimberg 
155*b2165e7aSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
156074be161SYohann Dudouit       } else if (dim == 3) {
1579e201c85SYohann         const CeedInt opt_elems       = 32 / (Q_1d * Q_1d);
1589e201c85SYohann         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
159*b2165e7aSSebastian Grimberg         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
160*b2165e7aSSebastian Grimberg 
161*b2165e7aSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
162074be161SYohann Dudouit       }
163ab213215SJeremy L Thompson     } break;
164ab213215SJeremy L Thompson     // LCOV_EXCL_START
165ab213215SJeremy L Thompson     // Evaluate the divergence to/from the quadrature points
166ab213215SJeremy L Thompson     case CEED_EVAL_DIV:
167e15f9bd0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_DIV not supported");
168ab213215SJeremy L Thompson     // Evaluate the curl to/from the quadrature points
169ab213215SJeremy L Thompson     case CEED_EVAL_CURL:
170e15f9bd0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_CURL not supported");
171ab213215SJeremy L Thompson     // Take no action, BasisApply should not have been called
172ab213215SJeremy L Thompson     case CEED_EVAL_NONE:
1732b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_NONE does not make sense in this context");
174ab213215SJeremy L Thompson       // LCOV_EXCL_STOP
175c532df63SYohann   }
176c532df63SYohann 
177ab213215SJeremy L Thompson   // Restore vectors
178437930d1SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
1792b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
180c532df63SYohann   }
1812b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
182e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
183c532df63SYohann }
184c532df63SYohann 
185ab213215SJeremy L Thompson //------------------------------------------------------------------------------
186ab213215SJeremy L Thompson // Destroy basis
187ab213215SJeremy L Thompson //------------------------------------------------------------------------------
188c532df63SYohann static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) {
189c532df63SYohann   Ceed ceed;
1902b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
191c532df63SYohann 
192c532df63SYohann   CeedBasis_Cuda_shared *data;
1932b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
194c532df63SYohann 
1952b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
196c532df63SYohann 
1972b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d));
1982b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp_1d));
1992b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad_1d));
2002b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_collo_grad_1d));
201c532df63SYohann 
2022b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
203c532df63SYohann 
204e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
205c532df63SYohann }
206c532df63SYohann 
207ab213215SJeremy L Thompson //------------------------------------------------------------------------------
208ab213215SJeremy L Thompson // Create tensor basis
209ab213215SJeremy L Thompson //------------------------------------------------------------------------------
2102b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
2112b730f8bSJeremy L Thompson                                         const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
212c532df63SYohann   Ceed ceed;
2132b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
214c532df63SYohann   CeedBasis_Cuda_shared *data;
2152b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
216c532df63SYohann 
217ab213215SJeremy L Thompson   // Copy basis data to GPU
218437930d1SJeremy L Thompson   const CeedInt q_bytes = Q_1d * sizeof(CeedScalar);
2192b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
2202b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice));
221c532df63SYohann 
222437930d1SJeremy L Thompson   const CeedInt interp_bytes = q_bytes * P_1d;
2232b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
2242b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice));
225c532df63SYohann 
2262b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes));
2272b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice));
228c532df63SYohann 
229ab213215SJeremy L Thompson   // Compute collocated gradient and copy to GPU
230437930d1SJeremy L Thompson   data->d_collo_grad_1d    = NULL;
2319e201c85SYohann   bool has_collocated_grad = dim == 3 && Q_1d >= P_1d;
2329e201c85SYohann   if (has_collocated_grad) {
233437930d1SJeremy L Thompson     CeedScalar *collo_grad_1d;
2342b730f8bSJeremy L Thompson     CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d));
2352b730f8bSJeremy L Thompson     CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d));
2362b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d));
2372b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, cudaMemcpyHostToDevice));
2382b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&collo_grad_1d));
239ac421f39SYohann   }
240ac421f39SYohann 
241ab213215SJeremy L Thompson   // Compile basis kernels
242437930d1SJeremy L Thompson   CeedInt num_comp;
2432b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
244437930d1SJeremy L Thompson   char *basis_kernel_path, *basis_kernel_source;
2452b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-shared-basis-tensor.h", &basis_kernel_path));
24623d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
2472b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
24823d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete -----\n");
249eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D",
250eb7e6cafSJeremy L Thompson                                    CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim),
251eb7e6cafSJeremy L Thompson                                    "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_HAS_COLLOCATED_GRAD", has_collocated_grad));
252eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
253eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
254eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
255eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose));
256eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
2572b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_path));
2582b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_source));
259c532df63SYohann 
2602b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
261ab213215SJeremy L Thompson 
262ab213215SJeremy L Thompson   // Register backend functions
2632b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared));
2642b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared));
265e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
266c532df63SYohann }
2672a86cc9dSSebastian Grimberg 
268ab213215SJeremy L Thompson //------------------------------------------------------------------------------
269