1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 #include <ceed.h> 9 #include <ceed/backend.h> 10 #include <ceed/jit-tools.h> 11 #include <cuda.h> 12 #include <cuda_runtime.h> 13 14 #include "../cuda/ceed-cuda-common.h" 15 #include "../cuda/ceed-cuda-compile.h" 16 #include "ceed-cuda-ref.h" 17 18 //------------------------------------------------------------------------------ 19 // Basis apply - tensor 20 //------------------------------------------------------------------------------ 21 int CeedBasisApply_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 22 Ceed ceed; 23 CeedInt Q_1d, dim; 24 const CeedInt transpose = t_mode == CEED_TRANSPOSE; 25 const int max_block_size = 32; 26 const CeedScalar *d_u; 27 CeedScalar *d_v; 28 CeedBasis_Cuda *data; 29 30 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 31 CeedCallBackend(CeedBasisGetData(basis, &data)); 32 33 // Read vectors 34 if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 35 else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 36 CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 37 38 // Clear v for transpose operation 39 if (transpose) { 40 CeedSize length; 41 42 CeedCallBackend(CeedVectorGetLength(v, &length)); 43 CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar))); 44 } 45 CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 46 CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 47 48 // Basis action 49 switch (eval_mode) { 50 case CEED_EVAL_INTERP: { 51 void *interp_args[] = {(void *)&num_elem, (void *)&transpose, &data->d_interp_1d, &d_u, &d_v}; 52 const CeedInt block_size = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size); 53 54 CeedCallBackend(CeedRunKernel_Cuda(ceed, data->Interp, num_elem, block_size, interp_args)); 55 } break; 56 case CEED_EVAL_GRAD: { 57 void *grad_args[] = {(void *)&num_elem, (void *)&transpose, &data->d_interp_1d, &data->d_grad_1d, &d_u, &d_v}; 58 const CeedInt block_size = max_block_size; 59 60 CeedCallBackend(CeedRunKernel_Cuda(ceed, data->Grad, num_elem, block_size, grad_args)); 61 } break; 62 case CEED_EVAL_WEIGHT: { 63 void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v}; 64 const int block_size_x = Q_1d; 65 const int block_size_y = dim >= 2 ? Q_1d : 1; 66 67 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, num_elem, block_size_x, block_size_y, 1, weight_args)); 68 } break; 69 // LCOV_EXCL_START 70 case CEED_EVAL_DIV: 71 case CEED_EVAL_CURL: 72 return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 73 case CEED_EVAL_NONE: 74 return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_NONE does not make sense in this context"); 75 // LCOV_EXCL_STOP 76 } 77 78 // Restore vectors 79 if (eval_mode != CEED_EVAL_WEIGHT) { 80 CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 81 } 82 CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 83 return CEED_ERROR_SUCCESS; 84 } 85 86 //------------------------------------------------------------------------------ 87 // Basis apply - non-tensor 88 //------------------------------------------------------------------------------ 89 int CeedBasisApplyNonTensor_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, 90 CeedVector v) { 91 Ceed ceed; 92 CeedInt num_nodes, num_qpts; 93 const CeedInt transpose = t_mode == CEED_TRANSPOSE; 94 const int elems_per_block = 1; 95 const int grid = CeedDivUpInt(num_elem, elems_per_block); 96 const CeedScalar *d_u; 97 CeedScalar *d_v; 98 CeedBasisNonTensor_Cuda *data; 99 100 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 101 CeedCallBackend(CeedBasisGetData(basis, &data)); 102 CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 103 CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes)); 104 105 // Read vectors 106 if (eval_mode != CEED_EVAL_WEIGHT) { 107 CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 108 } 109 CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 110 111 // Clear v for transpose operation 112 if (transpose) { 113 CeedSize length; 114 115 CeedCallBackend(CeedVectorGetLength(v, &length)); 116 CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar))); 117 } 118 119 // Apply basis operation 120 switch (eval_mode) { 121 case CEED_EVAL_INTERP: { 122 void *interp_args[] = {(void *)&num_elem, &data->d_interp, &d_u, &d_v}; 123 const int block_size_x = transpose ? num_nodes : num_qpts; 124 125 if (transpose) { 126 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->InterpTranspose, grid, block_size_x, 1, elems_per_block, interp_args)); 127 } else { 128 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Interp, grid, block_size_x, 1, elems_per_block, interp_args)); 129 } 130 } break; 131 case CEED_EVAL_GRAD: { 132 void *grad_args[] = {(void *)&num_elem, &data->d_grad, &d_u, &d_v}; 133 const int block_size_x = transpose ? num_nodes : num_qpts; 134 135 if (transpose) { 136 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, grad_args)); 137 } else { 138 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, grad_args)); 139 } 140 } break; 141 case CEED_EVAL_DIV: { 142 void *div_args[] = {(void *)&num_elem, &data->d_div, &d_u, &d_v}; 143 const int block_size_x = transpose ? num_nodes : num_qpts; 144 145 if (transpose) { 146 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, div_args)); 147 } else { 148 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, div_args)); 149 } 150 } break; 151 case CEED_EVAL_CURL: { 152 void *curl_args[] = {(void *)&num_elem, &data->d_curl, &d_u, &d_v}; 153 const int block_size_x = transpose ? num_nodes : num_qpts; 154 155 if (transpose) { 156 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, curl_args)); 157 } else { 158 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, curl_args)); 159 } 160 } break; 161 case CEED_EVAL_WEIGHT: { 162 void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight, &d_v}; 163 164 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, num_qpts, 1, elems_per_block, weight_args)); 165 } break; 166 // LCOV_EXCL_START 167 case CEED_EVAL_NONE: 168 return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_NONE does not make sense in this context"); 169 // LCOV_EXCL_STOP 170 } 171 172 // Restore vectors 173 if (eval_mode != CEED_EVAL_WEIGHT) { 174 CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 175 } 176 CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 177 return CEED_ERROR_SUCCESS; 178 } 179 180 //------------------------------------------------------------------------------ 181 // Destroy tensor basis 182 //------------------------------------------------------------------------------ 183 static int CeedBasisDestroy_Cuda(CeedBasis basis) { 184 Ceed ceed; 185 CeedBasis_Cuda *data; 186 187 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 188 CeedCallBackend(CeedBasisGetData(basis, &data)); 189 CeedCallCuda(ceed, cuModuleUnload(data->module)); 190 CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d)); 191 CeedCallCuda(ceed, cudaFree(data->d_interp_1d)); 192 CeedCallCuda(ceed, cudaFree(data->d_grad_1d)); 193 CeedCallBackend(CeedFree(&data)); 194 return CEED_ERROR_SUCCESS; 195 } 196 197 //------------------------------------------------------------------------------ 198 // Destroy non-tensor basis 199 //------------------------------------------------------------------------------ 200 static int CeedBasisDestroyNonTensor_Cuda(CeedBasis basis) { 201 Ceed ceed; 202 CeedBasisNonTensor_Cuda *data; 203 204 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 205 CeedCallBackend(CeedBasisGetData(basis, &data)); 206 CeedCallCuda(ceed, cuModuleUnload(data->module)); 207 CeedCallCuda(ceed, cudaFree(data->d_q_weight)); 208 CeedCallCuda(ceed, cudaFree(data->d_interp)); 209 CeedCallCuda(ceed, cudaFree(data->d_grad)); 210 CeedCallCuda(ceed, cudaFree(data->d_div)); 211 CeedCallCuda(ceed, cudaFree(data->d_curl)); 212 CeedCallBackend(CeedFree(&data)); 213 return CEED_ERROR_SUCCESS; 214 } 215 216 //------------------------------------------------------------------------------ 217 // Create tensor 218 //------------------------------------------------------------------------------ 219 int CeedBasisCreateTensorH1_Cuda(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d, 220 const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) { 221 Ceed ceed; 222 char *basis_kernel_path, *basis_kernel_source; 223 CeedInt num_comp; 224 const CeedInt q_bytes = Q_1d * sizeof(CeedScalar); 225 const CeedInt interp_bytes = q_bytes * P_1d; 226 CeedBasis_Cuda *data; 227 228 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 229 CeedCallBackend(CeedCalloc(1, &data)); 230 231 // Copy data to GPU 232 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes)); 233 CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice)); 234 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes)); 235 CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 236 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes)); 237 CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice)); 238 239 // Compile basis kernels 240 CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 241 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-tensor.h", &basis_kernel_path)); 242 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 243 CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 244 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n"); 245 CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 7, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "BASIS_BUF_LEN", 246 num_comp * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, 247 "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim))); 248 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 249 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad)); 250 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 251 CeedCallBackend(CeedFree(&basis_kernel_path)); 252 CeedCallBackend(CeedFree(&basis_kernel_source)); 253 254 CeedCallBackend(CeedBasisSetData(basis, data)); 255 256 // Register backend functions 257 CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApply_Cuda)); 258 CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda)); 259 return CEED_ERROR_SUCCESS; 260 } 261 262 //------------------------------------------------------------------------------ 263 // Create non-tensor H^1 264 //------------------------------------------------------------------------------ 265 int CeedBasisCreateH1_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *grad, 266 const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 267 Ceed ceed; 268 char *basis_kernel_path, *basis_kernel_source; 269 CeedInt num_comp, q_comp_interp, q_comp_grad; 270 const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 271 CeedBasisNonTensor_Cuda *data; 272 273 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 274 CeedCallBackend(CeedCalloc(1, &data)); 275 276 // Copy basis data to GPU 277 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 278 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad)); 279 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes)); 280 CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice)); 281 if (interp) { 282 const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 283 284 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes)); 285 CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice)); 286 } 287 if (grad) { 288 const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad; 289 290 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad, grad_bytes)); 291 CeedCallCuda(ceed, cudaMemcpy(data->d_grad, grad, grad_bytes, cudaMemcpyHostToDevice)); 292 } 293 294 // Compile basis kernels 295 CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 296 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path)); 297 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 298 CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 299 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n"); 300 CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP", 301 q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_grad, "BASIS_NUM_COMP", num_comp)); 302 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 303 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 304 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv)); 305 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose)); 306 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 307 CeedCallBackend(CeedFree(&basis_kernel_path)); 308 CeedCallBackend(CeedFree(&basis_kernel_source)); 309 310 CeedCallBackend(CeedBasisSetData(basis, data)); 311 312 // Register backend functions 313 CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda)); 314 CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda)); 315 return CEED_ERROR_SUCCESS; 316 } 317 318 //------------------------------------------------------------------------------ 319 // Create non-tensor H(div) 320 //------------------------------------------------------------------------------ 321 int CeedBasisCreateHdiv_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *div, 322 const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 323 Ceed ceed; 324 char *basis_kernel_path, *basis_kernel_source; 325 CeedInt num_comp, q_comp_interp, q_comp_div; 326 const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 327 CeedBasisNonTensor_Cuda *data; 328 329 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 330 CeedCallBackend(CeedCalloc(1, &data)); 331 332 // Copy basis data to GPU 333 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 334 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div)); 335 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes)); 336 CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice)); 337 if (interp) { 338 const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 339 340 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes)); 341 CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice)); 342 } 343 if (div) { 344 const CeedInt div_bytes = q_bytes * num_nodes * q_comp_div; 345 346 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_div, div_bytes)); 347 CeedCallCuda(ceed, cudaMemcpy(data->d_div, div, div_bytes, cudaMemcpyHostToDevice)); 348 } 349 350 // Compile basis kernels 351 CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 352 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path)); 353 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 354 CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 355 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n"); 356 CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP", 357 q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_div, "BASIS_NUM_COMP", num_comp)); 358 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 359 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 360 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv)); 361 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose)); 362 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 363 CeedCallBackend(CeedFree(&basis_kernel_path)); 364 CeedCallBackend(CeedFree(&basis_kernel_source)); 365 366 CeedCallBackend(CeedBasisSetData(basis, data)); 367 368 // Register backend functions 369 CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda)); 370 CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda)); 371 return CEED_ERROR_SUCCESS; 372 } 373 374 //------------------------------------------------------------------------------ 375 // Create non-tensor H(curl) 376 //------------------------------------------------------------------------------ 377 int CeedBasisCreateHcurl_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 378 const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 379 Ceed ceed; 380 char *basis_kernel_path, *basis_kernel_source; 381 CeedInt num_comp, q_comp_interp, q_comp_curl; 382 const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 383 CeedBasisNonTensor_Cuda *data; 384 385 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 386 CeedCallBackend(CeedCalloc(1, &data)); 387 388 // Copy basis data to GPU 389 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 390 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl)); 391 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes)); 392 CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice)); 393 if (interp) { 394 const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 395 396 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes)); 397 CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice)); 398 } 399 if (curl) { 400 const CeedInt curl_bytes = q_bytes * num_nodes * q_comp_curl; 401 402 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_curl, curl_bytes)); 403 CeedCallCuda(ceed, cudaMemcpy(data->d_curl, curl, curl_bytes, cudaMemcpyHostToDevice)); 404 } 405 406 // Compile basis kernels 407 CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 408 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path)); 409 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 410 CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 411 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n"); 412 CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP", 413 q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_curl, "BASIS_NUM_COMP", num_comp)); 414 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 415 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 416 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv)); 417 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose)); 418 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 419 CeedCallBackend(CeedFree(&basis_kernel_path)); 420 CeedCallBackend(CeedFree(&basis_kernel_source)); 421 422 CeedCallBackend(CeedBasisSetData(basis, data)); 423 424 // Register backend functions 425 CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda)); 426 CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda)); 427 return CEED_ERROR_SUCCESS; 428 } 429 430 //------------------------------------------------------------------------------ 431