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 //------------------------------------------------------------------------------ 232b730f8bSJeremy L Thompson int CeedCudaInitInterp(CeedScalar *d_B, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B); 242b730f8bSJeremy L Thompson int CeedCudaInitGrad(CeedScalar *d_B, CeedScalar *d_G, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B_ptr, CeedScalar **c_G_ptr); 252b730f8bSJeremy L Thompson int CeedCudaInitCollocatedGrad(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; 45437930d1SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) { 462b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 47c532df63SYohann } 482b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 49c532df63SYohann 50ab213215SJeremy L Thompson // Apply basis operation 51437930d1SJeremy L Thompson switch (eval_mode) { 52ab213215SJeremy L Thompson case CEED_EVAL_INTERP: { 53437930d1SJeremy L Thompson CeedInt P_1d, Q_1d; 542b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 552b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 56437930d1SJeremy L Thompson CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 572b730f8bSJeremy L Thompson CeedCallBackend(CeedCudaInitInterp(data->d_interp_1d, P_1d, Q_1d, &data->c_B)); 582b730f8bSJeremy L Thompson void *interp_args[] = {(void *)&num_elem, &data->c_B, &d_u, &d_v}; 594d537eeaSYohann if (dim == 1) { 602b730f8bSJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 61e6f67ff7SJed Brown 1)); // avoid >512 total threads 622b730f8bSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 63437930d1SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 649e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 652b730f8bSJeremy L Thompson CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->InterpTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args)); 669e201c85SYohann } else { 672b730f8bSJeremy L Thompson CeedCallBackend(CeedRunKernelDimSharedCuda(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); 759e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 762b730f8bSJeremy L Thompson CeedCallBackend( 772b730f8bSJeremy L Thompson CeedRunKernelDimSharedCuda(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 789e201c85SYohann } else { 792b730f8bSJeremy L Thompson CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 809e201c85SYohann } 81074be161SYohann Dudouit } else if (dim == 3) { 82437930d1SJeremy L Thompson CeedInt elems_per_block = 1; 832b730f8bSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 842b730f8bSJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 859e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 862b730f8bSJeremy L Thompson CeedCallBackend( 872b730f8bSJeremy L Thompson CeedRunKernelDimSharedCuda(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 889e201c85SYohann } else { 892b730f8bSJeremy L Thompson CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 90074be161SYohann Dudouit } 919e201c85SYohann } 92ab213215SJeremy L Thompson } break; 93ab213215SJeremy L Thompson case CEED_EVAL_GRAD: { 94437930d1SJeremy L Thompson CeedInt P_1d, Q_1d; 952b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 962b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 97437930d1SJeremy L Thompson CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 989e201c85SYohann if (data->d_collo_grad_1d) { 992b730f8bSJeremy L Thompson CeedCallBackend(CeedCudaInitCollocatedGrad(data->d_interp_1d, data->d_collo_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G)); 1009e201c85SYohann } else { 1012b730f8bSJeremy L Thompson CeedCallBackend(CeedCudaInitGrad(data->d_interp_1d, data->d_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G)); 1029e201c85SYohann } 1032b730f8bSJeremy L Thompson void *grad_args[] = {(void *)&num_elem, &data->c_B, &data->c_G, &d_u, &d_v}; 1044d537eeaSYohann if (dim == 1) { 1052b730f8bSJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 106e6f67ff7SJed Brown 1)); // avoid >512 total threads 1072b730f8bSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 108437930d1SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 1099e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 1102b730f8bSJeremy L Thompson CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->GradTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args)); 1119e201c85SYohann } else { 1122b730f8bSJeremy L Thompson CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args)); 1139e201c85SYohann } 114074be161SYohann Dudouit } else if (dim == 2) { 115437930d1SJeremy L Thompson const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8}; 116437930d1SJeremy L Thompson // elems_per_block must be at least 1 1172b730f8bSJeremy L Thompson CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1); 1182b730f8bSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 1192b730f8bSJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 1209e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 1212b730f8bSJeremy L Thompson CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 1229e201c85SYohann } else { 1232b730f8bSJeremy L Thompson CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 1249e201c85SYohann } 125074be161SYohann Dudouit } else if (dim == 3) { 126437930d1SJeremy L Thompson CeedInt elems_per_block = 1; 1272b730f8bSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 1282b730f8bSJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 1299e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 1302b730f8bSJeremy L Thompson CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 1319e201c85SYohann } else { 1322b730f8bSJeremy L Thompson CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 1339e201c85SYohann } 134074be161SYohann Dudouit } 135ab213215SJeremy L Thompson } break; 136ab213215SJeremy L Thompson case CEED_EVAL_WEIGHT: { 137437930d1SJeremy L Thompson CeedInt Q_1d; 1382b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 139437930d1SJeremy L Thompson void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v}; 140074be161SYohann Dudouit if (dim == 1) { 141437930d1SJeremy L Thompson const CeedInt elems_per_block = 32 / Q_1d; 1422b730f8bSJeremy L Thompson const CeedInt gridsize = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 1432b730f8bSJeremy L Thompson CeedCallBackend(CeedRunKernelDimCuda(ceed, data->Weight, gridsize, Q_1d, elems_per_block, 1, weight_args)); 144074be161SYohann Dudouit } else if (dim == 2) { 145437930d1SJeremy L Thompson const CeedInt opt_elems = 32 / (Q_1d * Q_1d); 146437930d1SJeremy L Thompson const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1; 1472b730f8bSJeremy L Thompson const CeedInt gridsize = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 1482b730f8bSJeremy L Thompson CeedCallBackend(CeedRunKernelDimCuda(ceed, data->Weight, gridsize, Q_1d, Q_1d, elems_per_block, weight_args)); 149074be161SYohann Dudouit } else if (dim == 3) { 1509e201c85SYohann const CeedInt opt_elems = 32 / (Q_1d * Q_1d); 1519e201c85SYohann const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1; 1522b730f8bSJeremy L Thompson const CeedInt gridsize = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 1532b730f8bSJeremy L Thompson CeedCallBackend(CeedRunKernelDimCuda(ceed, data->Weight, gridsize, Q_1d, Q_1d, elems_per_block, weight_args)); 154074be161SYohann Dudouit } 155ab213215SJeremy L Thompson } break; 156ab213215SJeremy L Thompson // LCOV_EXCL_START 157ab213215SJeremy L Thompson // Evaluate the divergence to/from the quadrature points 158ab213215SJeremy L Thompson case CEED_EVAL_DIV: 159e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_DIV not supported"); 160ab213215SJeremy L Thompson // Evaluate the curl to/from the quadrature points 161ab213215SJeremy L Thompson case CEED_EVAL_CURL: 162e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_CURL not supported"); 163ab213215SJeremy L Thompson // Take no action, BasisApply should not have been called 164ab213215SJeremy L Thompson case CEED_EVAL_NONE: 1652b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_NONE does not make sense in this context"); 166ab213215SJeremy L Thompson // LCOV_EXCL_STOP 167c532df63SYohann } 168c532df63SYohann 169ab213215SJeremy L Thompson // Restore vectors 170437930d1SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) { 1712b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 172c532df63SYohann } 1732b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 174e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 175c532df63SYohann } 176c532df63SYohann 177ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 178ab213215SJeremy L Thompson // Destroy basis 179ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 180c532df63SYohann static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) { 181c532df63SYohann Ceed ceed; 1822b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 183c532df63SYohann 184c532df63SYohann CeedBasis_Cuda_shared *data; 1852b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 186c532df63SYohann 1872b730f8bSJeremy L Thompson CeedCallCuda(ceed, cuModuleUnload(data->module)); 188c532df63SYohann 1892b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d)); 1902b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_interp_1d)); 1912b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_grad_1d)); 1922b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_collo_grad_1d)); 193c532df63SYohann 1942b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&data)); 195c532df63SYohann 196e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 197c532df63SYohann } 198c532df63SYohann 199ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 200ab213215SJeremy L Thompson // Create tensor basis 201ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 2022b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d, 2032b730f8bSJeremy L Thompson const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) { 204c532df63SYohann Ceed ceed; 2052b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 206c532df63SYohann CeedBasis_Cuda_shared *data; 2072b730f8bSJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 208c532df63SYohann 209ab213215SJeremy L Thompson // Copy basis data to GPU 210437930d1SJeremy L Thompson const CeedInt q_bytes = Q_1d * sizeof(CeedScalar); 2112b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes)); 2122b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice)); 213c532df63SYohann 214437930d1SJeremy L Thompson const CeedInt interp_bytes = q_bytes * P_1d; 2152b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes)); 2162b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 217c532df63SYohann 2182b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes)); 2192b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice)); 220c532df63SYohann 221ab213215SJeremy L Thompson // Compute collocated gradient and copy to GPU 222437930d1SJeremy L Thompson data->d_collo_grad_1d = NULL; 2239e201c85SYohann bool has_collocated_grad = dim == 3 && Q_1d >= P_1d; 2249e201c85SYohann if (has_collocated_grad) { 225437930d1SJeremy L Thompson CeedScalar *collo_grad_1d; 2262b730f8bSJeremy L Thompson CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d)); 2272b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d)); 2282b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d)); 2292b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, cudaMemcpyHostToDevice)); 2302b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&collo_grad_1d)); 231ac421f39SYohann } 232ac421f39SYohann 233ab213215SJeremy L Thompson // Compile basis kernels 234437930d1SJeremy L Thompson CeedInt num_comp; 2352b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 236437930d1SJeremy L Thompson char *basis_kernel_path, *basis_kernel_source; 2372b730f8bSJeremy L Thompson CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-shared-basis-tensor.h", &basis_kernel_path)); 23846dc0734SJeremy L Thompson CeedDebug256(ceed, 2, "----- Loading Basis Kernel Source -----\n"); 2392b730f8bSJeremy L Thompson CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 24046dc0734SJeremy L Thompson CeedDebug256(ceed, 2, "----- Loading Basis Kernel Source Complete -----\n"); 2412b730f8bSJeremy 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), 2422b730f8bSJeremy L Thompson "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", 2432b730f8bSJeremy L Thompson CeedIntPow(Q_1d, dim), "BASIS_HAS_COLLOCATED_GRAD", has_collocated_grad)); 2442b730f8bSJeremy L Thompson CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "Interp", &data->Interp)); 2452b730f8bSJeremy L Thompson CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 2462b730f8bSJeremy L Thompson CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "Grad", &data->Grad)); 2472b730f8bSJeremy L Thompson CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "GradTranspose", &data->GradTranspose)); 2482b730f8bSJeremy L Thompson CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "Weight", &data->Weight)); 2492b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&basis_kernel_path)); 2502b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&basis_kernel_source)); 251c532df63SYohann 2522b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 253ab213215SJeremy L Thompson 254ab213215SJeremy L Thompson // Register backend functions 2552b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared)); 2562b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared)); 257e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 258c532df63SYohann } 259*2a86cc9dSSebastian Grimberg 260ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 261