15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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; 23ca735530SJeremy L Thompson CeedInt Q_1d, dim; 247bbbfca3SJeremy L Thompson const CeedInt is_transpose = t_mode == CEED_TRANSPOSE; 25437930d1SJeremy L Thompson const int max_block_size = 32; 260d0321e0SJeremy L Thompson const CeedScalar *d_u; 270d0321e0SJeremy L Thompson CeedScalar *d_v; 28ca735530SJeremy L Thompson CeedBasis_Cuda *data; 29ca735530SJeremy L Thompson 30ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 31ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 32ca735530SJeremy L Thompson 339ea2cfd9SJeremy L Thompson // Get read/write access to u, v 346574a04fSJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 356574a04fSJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 362b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 370d0321e0SJeremy L Thompson 380d0321e0SJeremy L Thompson // Clear v for transpose operation 397bbbfca3SJeremy L Thompson if (is_transpose) { 401f9221feSJeremy L Thompson CeedSize length; 41ca735530SJeremy L Thompson 422b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorGetLength(v, &length)); 432b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar))); 440d0321e0SJeremy L Thompson } 452b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 462b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 470d0321e0SJeremy L Thompson 480d0321e0SJeremy L Thompson // Basis action 49437930d1SJeremy L Thompson switch (eval_mode) { 500d0321e0SJeremy L Thompson case CEED_EVAL_INTERP: { 517bbbfca3SJeremy L Thompson void *interp_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_interp_1d, &d_u, &d_v}; 52b2165e7aSSebastian Grimberg const CeedInt block_size = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size); 530d0321e0SJeremy L Thompson 54eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernel_Cuda(ceed, data->Interp, num_elem, block_size, interp_args)); 550d0321e0SJeremy L Thompson } break; 560d0321e0SJeremy L Thompson case CEED_EVAL_GRAD: { 577bbbfca3SJeremy L Thompson void *grad_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_interp_1d, &data->d_grad_1d, &d_u, &d_v}; 58b2165e7aSSebastian Grimberg const CeedInt block_size = max_block_size; 590d0321e0SJeremy L Thompson 60eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernel_Cuda(ceed, data->Grad, num_elem, block_size, grad_args)); 610d0321e0SJeremy L Thompson } break; 620d0321e0SJeremy L Thompson case CEED_EVAL_WEIGHT: { 63*097cc795SJames Wright CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights_1d not set", CeedEvalModes[eval_mode]); 64437930d1SJeremy L Thompson void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v}; 65b2165e7aSSebastian Grimberg const int block_size_x = Q_1d; 66b2165e7aSSebastian Grimberg const int block_size_y = dim >= 2 ? Q_1d : 1; 67b2165e7aSSebastian Grimberg 68b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, num_elem, block_size_x, block_size_y, 1, weight_args)); 690d0321e0SJeremy L Thompson } break; 709ea2cfd9SJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 719ea2cfd9SJeremy L Thompson break; 720d0321e0SJeremy L Thompson // LCOV_EXCL_START 730d0321e0SJeremy L Thompson case CEED_EVAL_DIV: 740d0321e0SJeremy L Thompson case CEED_EVAL_CURL: 75bcbe1c99SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 760d0321e0SJeremy L Thompson // LCOV_EXCL_STOP 770d0321e0SJeremy L Thompson } 780d0321e0SJeremy L Thompson 799ea2cfd9SJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 802b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 819ea2cfd9SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 829ea2cfd9SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 830d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 840d0321e0SJeremy L Thompson } 850d0321e0SJeremy L Thompson 860d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 870d0321e0SJeremy L Thompson // Basis apply - non-tensor 880d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 892b730f8bSJeremy L Thompson int CeedBasisApplyNonTensor_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, 902b730f8bSJeremy L Thompson CeedVector v) { 910d0321e0SJeremy L Thompson Ceed ceed; 92437930d1SJeremy L Thompson CeedInt num_nodes, num_qpts; 937bbbfca3SJeremy L Thompson const CeedInt is_transpose = t_mode == CEED_TRANSPOSE; 94d075f50bSSebastian Grimberg const int elems_per_block = 1; 95d075f50bSSebastian Grimberg const int grid = CeedDivUpInt(num_elem, elems_per_block); 960d0321e0SJeremy L Thompson const CeedScalar *d_u; 970d0321e0SJeremy L Thompson CeedScalar *d_v; 98ca735530SJeremy L Thompson CeedBasisNonTensor_Cuda *data; 99ca735530SJeremy L Thompson 100ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 101ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 102ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 103ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes)); 104ca735530SJeremy L Thompson 1059ea2cfd9SJeremy L Thompson // Get read/write access to u, v 1069ea2cfd9SJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 1079ea2cfd9SJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 1082b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 1090d0321e0SJeremy L Thompson 1100d0321e0SJeremy L Thompson // Clear v for transpose operation 1117bbbfca3SJeremy L Thompson if (is_transpose) { 1121f9221feSJeremy L Thompson CeedSize length; 113ca735530SJeremy L Thompson 1142b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorGetLength(v, &length)); 1152b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar))); 1160d0321e0SJeremy L Thompson } 1170d0321e0SJeremy L Thompson 1180d0321e0SJeremy L Thompson // Apply basis operation 119437930d1SJeremy L Thompson switch (eval_mode) { 1200d0321e0SJeremy L Thompson case CEED_EVAL_INTERP: { 121d075f50bSSebastian Grimberg void *interp_args[] = {(void *)&num_elem, &data->d_interp, &d_u, &d_v}; 1227bbbfca3SJeremy L Thompson const int block_size_x = is_transpose ? num_nodes : num_qpts; 123b2165e7aSSebastian Grimberg 1247bbbfca3SJeremy L Thompson if (is_transpose) { 125d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->InterpTranspose, grid, block_size_x, 1, elems_per_block, interp_args)); 126d075f50bSSebastian Grimberg } else { 127b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Interp, grid, block_size_x, 1, elems_per_block, interp_args)); 128d075f50bSSebastian Grimberg } 1290d0321e0SJeremy L Thompson } break; 1300d0321e0SJeremy L Thompson case CEED_EVAL_GRAD: { 131d075f50bSSebastian Grimberg void *grad_args[] = {(void *)&num_elem, &data->d_grad, &d_u, &d_v}; 1327bbbfca3SJeremy L Thompson const int block_size_x = is_transpose ? num_nodes : num_qpts; 133b2165e7aSSebastian Grimberg 1347bbbfca3SJeremy L Thompson if (is_transpose) { 135d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, grad_args)); 136d075f50bSSebastian Grimberg } else { 137d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, grad_args)); 138d075f50bSSebastian Grimberg } 139d075f50bSSebastian Grimberg } break; 140d075f50bSSebastian Grimberg case CEED_EVAL_DIV: { 141d075f50bSSebastian Grimberg void *div_args[] = {(void *)&num_elem, &data->d_div, &d_u, &d_v}; 1427bbbfca3SJeremy L Thompson const int block_size_x = is_transpose ? num_nodes : num_qpts; 143d075f50bSSebastian Grimberg 1447bbbfca3SJeremy L Thompson if (is_transpose) { 145d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, div_args)); 146d075f50bSSebastian Grimberg } else { 147d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, div_args)); 148d075f50bSSebastian Grimberg } 149d075f50bSSebastian Grimberg } break; 150d075f50bSSebastian Grimberg case CEED_EVAL_CURL: { 151d075f50bSSebastian Grimberg void *curl_args[] = {(void *)&num_elem, &data->d_curl, &d_u, &d_v}; 1527bbbfca3SJeremy L Thompson const int block_size_x = is_transpose ? num_nodes : num_qpts; 153d075f50bSSebastian Grimberg 1547bbbfca3SJeremy L Thompson if (is_transpose) { 155d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, curl_args)); 156d075f50bSSebastian Grimberg } else { 157d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, curl_args)); 158d075f50bSSebastian Grimberg } 1590d0321e0SJeremy L Thompson } break; 1600d0321e0SJeremy L Thompson case CEED_EVAL_WEIGHT: { 161*097cc795SJames Wright CeedCheck(data->d_q_weight, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights not set", CeedEvalModes[eval_mode]); 162437930d1SJeremy L Thompson void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight, &d_v}; 163b2165e7aSSebastian Grimberg 164eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, num_qpts, 1, elems_per_block, weight_args)); 1650d0321e0SJeremy L Thompson } break; 1669ea2cfd9SJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 1679ea2cfd9SJeremy L Thompson break; 1680d0321e0SJeremy L Thompson } 1690d0321e0SJeremy L Thompson 1709ea2cfd9SJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 1712b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 1729ea2cfd9SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 1739ea2cfd9SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 1740d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1750d0321e0SJeremy L Thompson } 1760d0321e0SJeremy L Thompson 1770d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1780d0321e0SJeremy L Thompson // Destroy tensor basis 1790d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1800d0321e0SJeremy L Thompson static int CeedBasisDestroy_Cuda(CeedBasis basis) { 1810d0321e0SJeremy L Thompson Ceed ceed; 1820d0321e0SJeremy L Thompson CeedBasis_Cuda *data; 183ca735530SJeremy L Thompson 184ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 1852b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 1862b730f8bSJeremy L Thompson CeedCallCuda(ceed, cuModuleUnload(data->module)); 187*097cc795SJames Wright if (data->d_q_weight_1d) CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d)); 1882b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_interp_1d)); 1892b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_grad_1d)); 1902b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&data)); 1910d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1920d0321e0SJeremy L Thompson } 1930d0321e0SJeremy L Thompson 1940d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1950d0321e0SJeremy L Thompson // Destroy non-tensor basis 1960d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1970d0321e0SJeremy L Thompson static int CeedBasisDestroyNonTensor_Cuda(CeedBasis basis) { 1980d0321e0SJeremy L Thompson Ceed ceed; 1990d0321e0SJeremy L Thompson CeedBasisNonTensor_Cuda *data; 200ca735530SJeremy L Thompson 201ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 2022b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 2032b730f8bSJeremy L Thompson CeedCallCuda(ceed, cuModuleUnload(data->module)); 204*097cc795SJames Wright if (data->d_q_weight) CeedCallCuda(ceed, cudaFree(data->d_q_weight)); 2052b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_interp)); 2062b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_grad)); 207d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaFree(data->d_div)); 208d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaFree(data->d_curl)); 2092b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&data)); 2100d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2110d0321e0SJeremy L Thompson } 2120d0321e0SJeremy L Thompson 2130d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 2140d0321e0SJeremy L Thompson // Create tensor 2150d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 2162b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d, 2172b730f8bSJeremy L Thompson const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) { 2180d0321e0SJeremy L Thompson Ceed ceed; 21922070f95SJeremy L Thompson char *basis_kernel_source; 22022070f95SJeremy L Thompson const char *basis_kernel_path; 221ca735530SJeremy L Thompson CeedInt num_comp; 222ca735530SJeremy L Thompson const CeedInt q_bytes = Q_1d * sizeof(CeedScalar); 223ca735530SJeremy L Thompson const CeedInt interp_bytes = q_bytes * P_1d; 2240d0321e0SJeremy L Thompson CeedBasis_Cuda *data; 225ca735530SJeremy L Thompson 226ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 2272b730f8bSJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 2280d0321e0SJeremy L Thompson 2290d0321e0SJeremy L Thompson // Copy data to GPU 230*097cc795SJames Wright if (q_weight_1d) { 2312b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes)); 2322b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice)); 233*097cc795SJames Wright } 2342b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes)); 2352b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 2362b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes)); 2372b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice)); 2380d0321e0SJeremy L Thompson 239ecc88aebSJeremy L Thompson // Compile basis kernels 2402b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 2412b730f8bSJeremy L Thompson CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-tensor.h", &basis_kernel_path)); 24223d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 2432b730f8bSJeremy L Thompson CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 24423d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n"); 245eb7e6cafSJeremy 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", 2462b730f8bSJeremy L Thompson num_comp * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, 2472b730f8bSJeremy L Thompson "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim))); 248eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 249eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad)); 250eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 2512b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&basis_kernel_path)); 2522b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&basis_kernel_source)); 253437930d1SJeremy L Thompson 2542b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 2550d0321e0SJeremy L Thompson 256d075f50bSSebastian Grimberg // Register backend functions 2572b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApply_Cuda)); 2582b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda)); 2590d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2600d0321e0SJeremy L Thompson } 2610d0321e0SJeremy L Thompson 2620d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 263d075f50bSSebastian Grimberg // Create non-tensor H^1 2640d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 2652b730f8bSJeremy L Thompson int CeedBasisCreateH1_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *grad, 26651475c7cSJeremy L Thompson const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 2670d0321e0SJeremy L Thompson Ceed ceed; 26822070f95SJeremy L Thompson char *basis_kernel_source; 26922070f95SJeremy L Thompson const char *basis_kernel_path; 270d075f50bSSebastian Grimberg CeedInt num_comp, q_comp_interp, q_comp_grad; 271ca735530SJeremy L Thompson const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 2720d0321e0SJeremy L Thompson CeedBasisNonTensor_Cuda *data; 273ca735530SJeremy L Thompson 274ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 2752b730f8bSJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 2760d0321e0SJeremy L Thompson 2770d0321e0SJeremy L Thompson // Copy basis data to GPU 278d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 279d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad)); 280*097cc795SJames Wright if (q_weight) { 2812b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes)); 2822b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice)); 283*097cc795SJames Wright } 284d075f50bSSebastian Grimberg if (interp) { 285d075f50bSSebastian Grimberg const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 286d075f50bSSebastian Grimberg 2872b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes)); 2882b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice)); 289d075f50bSSebastian Grimberg } 290d075f50bSSebastian Grimberg if (grad) { 291d075f50bSSebastian Grimberg const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad; 292d075f50bSSebastian Grimberg 2932b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad, grad_bytes)); 2942b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad, grad, grad_bytes, cudaMemcpyHostToDevice)); 295d075f50bSSebastian Grimberg } 2960d0321e0SJeremy L Thompson 2970d0321e0SJeremy L Thompson // Compile basis kernels 2982b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 2992b730f8bSJeremy L Thompson CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path)); 30023d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 3012b730f8bSJeremy L Thompson CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 30223d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n"); 303d075f50bSSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP", 304d075f50bSSebastian Grimberg q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_grad, "BASIS_NUM_COMP", num_comp)); 305d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 306d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 307d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv)); 308d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose)); 309d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 310d075f50bSSebastian Grimberg CeedCallBackend(CeedFree(&basis_kernel_path)); 311d075f50bSSebastian Grimberg CeedCallBackend(CeedFree(&basis_kernel_source)); 312d075f50bSSebastian Grimberg 313d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisSetData(basis, data)); 314d075f50bSSebastian Grimberg 315d075f50bSSebastian Grimberg // Register backend functions 316d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda)); 317d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda)); 318d075f50bSSebastian Grimberg return CEED_ERROR_SUCCESS; 319d075f50bSSebastian Grimberg } 320d075f50bSSebastian Grimberg 321d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 322d075f50bSSebastian Grimberg // Create non-tensor H(div) 323d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 324d075f50bSSebastian Grimberg int CeedBasisCreateHdiv_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *div, 325d075f50bSSebastian Grimberg const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 326d075f50bSSebastian Grimberg Ceed ceed; 32722070f95SJeremy L Thompson char *basis_kernel_source; 32822070f95SJeremy L Thompson const char *basis_kernel_path; 329d075f50bSSebastian Grimberg CeedInt num_comp, q_comp_interp, q_comp_div; 330d075f50bSSebastian Grimberg const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 331d075f50bSSebastian Grimberg CeedBasisNonTensor_Cuda *data; 332d075f50bSSebastian Grimberg 333d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 334d075f50bSSebastian Grimberg CeedCallBackend(CeedCalloc(1, &data)); 335d075f50bSSebastian Grimberg 336d075f50bSSebastian Grimberg // Copy basis data to GPU 337d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 338d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div)); 339*097cc795SJames Wright if (q_weight) { 340d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes)); 341d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice)); 342*097cc795SJames Wright } 343d075f50bSSebastian Grimberg if (interp) { 344d075f50bSSebastian Grimberg const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 345d075f50bSSebastian Grimberg 346d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes)); 347d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice)); 348d075f50bSSebastian Grimberg } 349d075f50bSSebastian Grimberg if (div) { 350d075f50bSSebastian Grimberg const CeedInt div_bytes = q_bytes * num_nodes * q_comp_div; 351d075f50bSSebastian Grimberg 352d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_div, div_bytes)); 353d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_div, div, div_bytes, cudaMemcpyHostToDevice)); 354d075f50bSSebastian Grimberg } 355d075f50bSSebastian Grimberg 356d075f50bSSebastian Grimberg // Compile basis kernels 357d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 358d075f50bSSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path)); 359d075f50bSSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 360d075f50bSSebastian Grimberg CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 361d075f50bSSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n"); 362d075f50bSSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP", 363d075f50bSSebastian Grimberg q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_div, "BASIS_NUM_COMP", num_comp)); 364d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 365d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 366d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv)); 367d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose)); 368d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 369d075f50bSSebastian Grimberg CeedCallBackend(CeedFree(&basis_kernel_path)); 370d075f50bSSebastian Grimberg CeedCallBackend(CeedFree(&basis_kernel_source)); 371d075f50bSSebastian Grimberg 372d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisSetData(basis, data)); 373d075f50bSSebastian Grimberg 374d075f50bSSebastian Grimberg // Register backend functions 375d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda)); 376d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda)); 377d075f50bSSebastian Grimberg return CEED_ERROR_SUCCESS; 378d075f50bSSebastian Grimberg } 379d075f50bSSebastian Grimberg 380d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 381d075f50bSSebastian Grimberg // Create non-tensor H(curl) 382d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 383d075f50bSSebastian Grimberg int CeedBasisCreateHcurl_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 384d075f50bSSebastian Grimberg const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 385d075f50bSSebastian Grimberg Ceed ceed; 38622070f95SJeremy L Thompson char *basis_kernel_source; 38722070f95SJeremy L Thompson const char *basis_kernel_path; 388d075f50bSSebastian Grimberg CeedInt num_comp, q_comp_interp, q_comp_curl; 389d075f50bSSebastian Grimberg const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 390d075f50bSSebastian Grimberg CeedBasisNonTensor_Cuda *data; 391d075f50bSSebastian Grimberg 392d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 393d075f50bSSebastian Grimberg CeedCallBackend(CeedCalloc(1, &data)); 394d075f50bSSebastian Grimberg 395d075f50bSSebastian Grimberg // Copy basis data to GPU 396d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 397d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl)); 398*097cc795SJames Wright if (q_weight) { 399d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes)); 400d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice)); 401*097cc795SJames Wright } 402d075f50bSSebastian Grimberg if (interp) { 403d075f50bSSebastian Grimberg const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 404d075f50bSSebastian Grimberg 405d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes)); 406d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice)); 407d075f50bSSebastian Grimberg } 408d075f50bSSebastian Grimberg if (curl) { 409d075f50bSSebastian Grimberg const CeedInt curl_bytes = q_bytes * num_nodes * q_comp_curl; 410d075f50bSSebastian Grimberg 411d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_curl, curl_bytes)); 412d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_curl, curl, curl_bytes, cudaMemcpyHostToDevice)); 413d075f50bSSebastian Grimberg } 414d075f50bSSebastian Grimberg 415d075f50bSSebastian Grimberg // Compile basis kernels 416d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 417d075f50bSSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path)); 418d075f50bSSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 419d075f50bSSebastian Grimberg CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 420d075f50bSSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n"); 421d075f50bSSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP", 422d075f50bSSebastian Grimberg q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_curl, "BASIS_NUM_COMP", num_comp)); 423d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 424d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 425d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv)); 426d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose)); 427d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 4282b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&basis_kernel_path)); 4292b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&basis_kernel_source)); 4300d0321e0SJeremy L Thompson 4312b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 4320d0321e0SJeremy L Thompson 4330d0321e0SJeremy L Thompson // Register backend functions 4342b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda)); 4352b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda)); 4360d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 4370d0321e0SJeremy L Thompson } 4382a86cc9dSSebastian Grimberg 4390d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 440