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. 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; 24437930d1SJeremy L Thompson const CeedInt 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 33*9ea2cfd9SJeremy 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 39d075f50bSSebastian Grimberg if (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: { 512b730f8bSJeremy L Thompson void *interp_args[] = {(void *)&num_elem, (void *)&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: { 572b730f8bSJeremy L Thompson void *grad_args[] = {(void *)&num_elem, (void *)&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: { 63437930d1SJeremy L Thompson void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v}; 64b2165e7aSSebastian Grimberg const int block_size_x = Q_1d; 65b2165e7aSSebastian Grimberg const int block_size_y = dim >= 2 ? Q_1d : 1; 66b2165e7aSSebastian Grimberg 67b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, num_elem, block_size_x, block_size_y, 1, weight_args)); 680d0321e0SJeremy L Thompson } break; 69*9ea2cfd9SJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 70*9ea2cfd9SJeremy L Thompson break; 710d0321e0SJeremy L Thompson // LCOV_EXCL_START 720d0321e0SJeremy L Thompson case CEED_EVAL_DIV: 730d0321e0SJeremy L Thompson case CEED_EVAL_CURL: 74bcbe1c99SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 750d0321e0SJeremy L Thompson // LCOV_EXCL_STOP 760d0321e0SJeremy L Thompson } 770d0321e0SJeremy L Thompson 78*9ea2cfd9SJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 792b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 80*9ea2cfd9SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 81*9ea2cfd9SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 820d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 830d0321e0SJeremy L Thompson } 840d0321e0SJeremy L Thompson 850d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 860d0321e0SJeremy L Thompson // Basis apply - non-tensor 870d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 882b730f8bSJeremy L Thompson int CeedBasisApplyNonTensor_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, 892b730f8bSJeremy L Thompson CeedVector v) { 900d0321e0SJeremy L Thompson Ceed ceed; 91437930d1SJeremy L Thompson CeedInt num_nodes, num_qpts; 92437930d1SJeremy L Thompson const CeedInt transpose = t_mode == CEED_TRANSPOSE; 93d075f50bSSebastian Grimberg const int elems_per_block = 1; 94d075f50bSSebastian Grimberg const int grid = CeedDivUpInt(num_elem, elems_per_block); 950d0321e0SJeremy L Thompson const CeedScalar *d_u; 960d0321e0SJeremy L Thompson CeedScalar *d_v; 97ca735530SJeremy L Thompson CeedBasisNonTensor_Cuda *data; 98ca735530SJeremy L Thompson 99ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 100ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 101ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 102ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes)); 103ca735530SJeremy L Thompson 104*9ea2cfd9SJeremy L Thompson // Get read/write access to u, v 105*9ea2cfd9SJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 106*9ea2cfd9SJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 1072b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 1080d0321e0SJeremy L Thompson 1090d0321e0SJeremy L Thompson // Clear v for transpose operation 110d075f50bSSebastian Grimberg if (transpose) { 1111f9221feSJeremy L Thompson CeedSize length; 112ca735530SJeremy L Thompson 1132b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorGetLength(v, &length)); 1142b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar))); 1150d0321e0SJeremy L Thompson } 1160d0321e0SJeremy L Thompson 1170d0321e0SJeremy L Thompson // Apply basis operation 118437930d1SJeremy L Thompson switch (eval_mode) { 1190d0321e0SJeremy L Thompson case CEED_EVAL_INTERP: { 120d075f50bSSebastian Grimberg void *interp_args[] = {(void *)&num_elem, &data->d_interp, &d_u, &d_v}; 121b2165e7aSSebastian Grimberg const int block_size_x = transpose ? num_nodes : num_qpts; 122b2165e7aSSebastian Grimberg 123d075f50bSSebastian Grimberg if (transpose) { 124d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->InterpTranspose, grid, block_size_x, 1, elems_per_block, interp_args)); 125d075f50bSSebastian Grimberg } else { 126b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Interp, grid, block_size_x, 1, elems_per_block, interp_args)); 127d075f50bSSebastian Grimberg } 1280d0321e0SJeremy L Thompson } break; 1290d0321e0SJeremy L Thompson case CEED_EVAL_GRAD: { 130d075f50bSSebastian Grimberg void *grad_args[] = {(void *)&num_elem, &data->d_grad, &d_u, &d_v}; 131b2165e7aSSebastian Grimberg const int block_size_x = transpose ? num_nodes : num_qpts; 132b2165e7aSSebastian Grimberg 133d075f50bSSebastian Grimberg if (transpose) { 134d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, grad_args)); 135d075f50bSSebastian Grimberg } else { 136d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, grad_args)); 137d075f50bSSebastian Grimberg } 138d075f50bSSebastian Grimberg } break; 139d075f50bSSebastian Grimberg case CEED_EVAL_DIV: { 140d075f50bSSebastian Grimberg void *div_args[] = {(void *)&num_elem, &data->d_div, &d_u, &d_v}; 141d075f50bSSebastian Grimberg const int block_size_x = transpose ? num_nodes : num_qpts; 142d075f50bSSebastian Grimberg 143d075f50bSSebastian Grimberg if (transpose) { 144d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, div_args)); 145d075f50bSSebastian Grimberg } else { 146d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, div_args)); 147d075f50bSSebastian Grimberg } 148d075f50bSSebastian Grimberg } break; 149d075f50bSSebastian Grimberg case CEED_EVAL_CURL: { 150d075f50bSSebastian Grimberg void *curl_args[] = {(void *)&num_elem, &data->d_curl, &d_u, &d_v}; 151d075f50bSSebastian Grimberg const int block_size_x = transpose ? num_nodes : num_qpts; 152d075f50bSSebastian Grimberg 153d075f50bSSebastian Grimberg if (transpose) { 154d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, curl_args)); 155d075f50bSSebastian Grimberg } else { 156d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, curl_args)); 157d075f50bSSebastian Grimberg } 1580d0321e0SJeremy L Thompson } break; 1590d0321e0SJeremy L Thompson case CEED_EVAL_WEIGHT: { 160437930d1SJeremy L Thompson void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight, &d_v}; 161b2165e7aSSebastian Grimberg 162eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, num_qpts, 1, elems_per_block, weight_args)); 1630d0321e0SJeremy L Thompson } break; 164*9ea2cfd9SJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 165*9ea2cfd9SJeremy L Thompson break; 1660d0321e0SJeremy L Thompson } 1670d0321e0SJeremy L Thompson 168*9ea2cfd9SJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 1692b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 170*9ea2cfd9SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 171*9ea2cfd9SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 1720d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1730d0321e0SJeremy L Thompson } 1740d0321e0SJeremy L Thompson 1750d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1760d0321e0SJeremy L Thompson // Destroy tensor basis 1770d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1780d0321e0SJeremy L Thompson static int CeedBasisDestroy_Cuda(CeedBasis basis) { 1790d0321e0SJeremy L Thompson Ceed ceed; 1800d0321e0SJeremy L Thompson CeedBasis_Cuda *data; 181ca735530SJeremy L Thompson 182ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 1832b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 1842b730f8bSJeremy L Thompson CeedCallCuda(ceed, cuModuleUnload(data->module)); 1852b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d)); 1862b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_interp_1d)); 1872b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_grad_1d)); 1882b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&data)); 1890d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1900d0321e0SJeremy L Thompson } 1910d0321e0SJeremy L Thompson 1920d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1930d0321e0SJeremy L Thompson // Destroy non-tensor basis 1940d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1950d0321e0SJeremy L Thompson static int CeedBasisDestroyNonTensor_Cuda(CeedBasis basis) { 1960d0321e0SJeremy L Thompson Ceed ceed; 1970d0321e0SJeremy L Thompson CeedBasisNonTensor_Cuda *data; 198ca735530SJeremy L Thompson 199ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 2002b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 2012b730f8bSJeremy L Thompson CeedCallCuda(ceed, cuModuleUnload(data->module)); 2022b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_q_weight)); 2032b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_interp)); 2042b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_grad)); 205d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaFree(data->d_div)); 206d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaFree(data->d_curl)); 2072b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&data)); 2080d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2090d0321e0SJeremy L Thompson } 2100d0321e0SJeremy L Thompson 2110d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 2120d0321e0SJeremy L Thompson // Create tensor 2130d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 2142b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d, 2152b730f8bSJeremy L Thompson const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) { 2160d0321e0SJeremy L Thompson Ceed ceed; 217ca735530SJeremy L Thompson char *basis_kernel_path, *basis_kernel_source; 218ca735530SJeremy L Thompson CeedInt num_comp; 219ca735530SJeremy L Thompson const CeedInt q_bytes = Q_1d * sizeof(CeedScalar); 220ca735530SJeremy L Thompson const CeedInt interp_bytes = q_bytes * P_1d; 2210d0321e0SJeremy L Thompson CeedBasis_Cuda *data; 222ca735530SJeremy L Thompson 223ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 2242b730f8bSJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 2250d0321e0SJeremy L Thompson 2260d0321e0SJeremy L Thompson // Copy data to GPU 2272b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes)); 2282b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice)); 2292b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes)); 2302b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 2312b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes)); 2322b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice)); 2330d0321e0SJeremy L Thompson 234ecc88aebSJeremy L Thompson // Compile basis kernels 2352b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 2362b730f8bSJeremy L Thompson CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-tensor.h", &basis_kernel_path)); 23723d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 2382b730f8bSJeremy L Thompson CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 23923d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n"); 240eb7e6cafSJeremy 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", 2412b730f8bSJeremy L Thompson num_comp * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, 2422b730f8bSJeremy L Thompson "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim))); 243eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 244eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad)); 245eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 2462b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&basis_kernel_path)); 2472b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&basis_kernel_source)); 248437930d1SJeremy L Thompson 2492b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 2500d0321e0SJeremy L Thompson 251d075f50bSSebastian Grimberg // Register backend functions 2522b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApply_Cuda)); 2532b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda)); 2540d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2550d0321e0SJeremy L Thompson } 2560d0321e0SJeremy L Thompson 2570d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 258d075f50bSSebastian Grimberg // Create non-tensor H^1 2590d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 2602b730f8bSJeremy L Thompson int CeedBasisCreateH1_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *grad, 26151475c7cSJeremy L Thompson const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 2620d0321e0SJeremy L Thompson Ceed ceed; 263ca735530SJeremy L Thompson char *basis_kernel_path, *basis_kernel_source; 264d075f50bSSebastian Grimberg CeedInt num_comp, q_comp_interp, q_comp_grad; 265ca735530SJeremy L Thompson const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 2660d0321e0SJeremy L Thompson CeedBasisNonTensor_Cuda *data; 267ca735530SJeremy L Thompson 268ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 2692b730f8bSJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 2700d0321e0SJeremy L Thompson 2710d0321e0SJeremy L Thompson // Copy basis data to GPU 272d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 273d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad)); 2742b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes)); 2752b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice)); 276d075f50bSSebastian Grimberg if (interp) { 277d075f50bSSebastian Grimberg const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 278d075f50bSSebastian Grimberg 2792b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes)); 2802b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice)); 281d075f50bSSebastian Grimberg } 282d075f50bSSebastian Grimberg if (grad) { 283d075f50bSSebastian Grimberg const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad; 284d075f50bSSebastian Grimberg 2852b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad, grad_bytes)); 2862b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad, grad, grad_bytes, cudaMemcpyHostToDevice)); 287d075f50bSSebastian Grimberg } 2880d0321e0SJeremy L Thompson 2890d0321e0SJeremy L Thompson // Compile basis kernels 2902b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 2912b730f8bSJeremy L Thompson CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path)); 29223d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 2932b730f8bSJeremy L Thompson CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 29423d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n"); 295d075f50bSSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP", 296d075f50bSSebastian Grimberg q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_grad, "BASIS_NUM_COMP", num_comp)); 297d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 298d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 299d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv)); 300d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose)); 301d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 302d075f50bSSebastian Grimberg CeedCallBackend(CeedFree(&basis_kernel_path)); 303d075f50bSSebastian Grimberg CeedCallBackend(CeedFree(&basis_kernel_source)); 304d075f50bSSebastian Grimberg 305d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisSetData(basis, data)); 306d075f50bSSebastian Grimberg 307d075f50bSSebastian Grimberg // Register backend functions 308d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda)); 309d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda)); 310d075f50bSSebastian Grimberg return CEED_ERROR_SUCCESS; 311d075f50bSSebastian Grimberg } 312d075f50bSSebastian Grimberg 313d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 314d075f50bSSebastian Grimberg // Create non-tensor H(div) 315d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 316d075f50bSSebastian Grimberg int CeedBasisCreateHdiv_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *div, 317d075f50bSSebastian Grimberg const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 318d075f50bSSebastian Grimberg Ceed ceed; 319d075f50bSSebastian Grimberg char *basis_kernel_path, *basis_kernel_source; 320d075f50bSSebastian Grimberg CeedInt num_comp, q_comp_interp, q_comp_div; 321d075f50bSSebastian Grimberg const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 322d075f50bSSebastian Grimberg CeedBasisNonTensor_Cuda *data; 323d075f50bSSebastian Grimberg 324d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 325d075f50bSSebastian Grimberg CeedCallBackend(CeedCalloc(1, &data)); 326d075f50bSSebastian Grimberg 327d075f50bSSebastian Grimberg // Copy basis data to GPU 328d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 329d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div)); 330d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes)); 331d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice)); 332d075f50bSSebastian Grimberg if (interp) { 333d075f50bSSebastian Grimberg const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 334d075f50bSSebastian Grimberg 335d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes)); 336d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice)); 337d075f50bSSebastian Grimberg } 338d075f50bSSebastian Grimberg if (div) { 339d075f50bSSebastian Grimberg const CeedInt div_bytes = q_bytes * num_nodes * q_comp_div; 340d075f50bSSebastian Grimberg 341d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_div, div_bytes)); 342d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_div, div, div_bytes, cudaMemcpyHostToDevice)); 343d075f50bSSebastian Grimberg } 344d075f50bSSebastian Grimberg 345d075f50bSSebastian Grimberg // Compile basis kernels 346d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 347d075f50bSSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path)); 348d075f50bSSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 349d075f50bSSebastian Grimberg CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 350d075f50bSSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n"); 351d075f50bSSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP", 352d075f50bSSebastian Grimberg q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_div, "BASIS_NUM_COMP", num_comp)); 353d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 354d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 355d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv)); 356d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose)); 357d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 358d075f50bSSebastian Grimberg CeedCallBackend(CeedFree(&basis_kernel_path)); 359d075f50bSSebastian Grimberg CeedCallBackend(CeedFree(&basis_kernel_source)); 360d075f50bSSebastian Grimberg 361d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisSetData(basis, data)); 362d075f50bSSebastian Grimberg 363d075f50bSSebastian Grimberg // Register backend functions 364d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda)); 365d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda)); 366d075f50bSSebastian Grimberg return CEED_ERROR_SUCCESS; 367d075f50bSSebastian Grimberg } 368d075f50bSSebastian Grimberg 369d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 370d075f50bSSebastian Grimberg // Create non-tensor H(curl) 371d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 372d075f50bSSebastian Grimberg int CeedBasisCreateHcurl_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 373d075f50bSSebastian Grimberg const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 374d075f50bSSebastian Grimberg Ceed ceed; 375d075f50bSSebastian Grimberg char *basis_kernel_path, *basis_kernel_source; 376d075f50bSSebastian Grimberg CeedInt num_comp, q_comp_interp, q_comp_curl; 377d075f50bSSebastian Grimberg const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 378d075f50bSSebastian Grimberg CeedBasisNonTensor_Cuda *data; 379d075f50bSSebastian Grimberg 380d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 381d075f50bSSebastian Grimberg CeedCallBackend(CeedCalloc(1, &data)); 382d075f50bSSebastian Grimberg 383d075f50bSSebastian Grimberg // Copy basis data to GPU 384d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 385d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl)); 386d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes)); 387d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice)); 388d075f50bSSebastian Grimberg if (interp) { 389d075f50bSSebastian Grimberg const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 390d075f50bSSebastian Grimberg 391d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes)); 392d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice)); 393d075f50bSSebastian Grimberg } 394d075f50bSSebastian Grimberg if (curl) { 395d075f50bSSebastian Grimberg const CeedInt curl_bytes = q_bytes * num_nodes * q_comp_curl; 396d075f50bSSebastian Grimberg 397d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_curl, curl_bytes)); 398d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_curl, curl, curl_bytes, cudaMemcpyHostToDevice)); 399d075f50bSSebastian Grimberg } 400d075f50bSSebastian Grimberg 401d075f50bSSebastian Grimberg // Compile basis kernels 402d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 403d075f50bSSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path)); 404d075f50bSSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 405d075f50bSSebastian Grimberg CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 406d075f50bSSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n"); 407d075f50bSSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP", 408d075f50bSSebastian Grimberg q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_curl, "BASIS_NUM_COMP", num_comp)); 409d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 410d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 411d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv)); 412d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose)); 413d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 4142b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&basis_kernel_path)); 4152b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&basis_kernel_source)); 4160d0321e0SJeremy L Thompson 4172b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 4180d0321e0SJeremy L Thompson 4190d0321e0SJeremy L Thompson // Register backend functions 4202b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda)); 4212b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda)); 4220d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 4230d0321e0SJeremy L Thompson } 4242a86cc9dSSebastian Grimberg 4250d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 426