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> 13111870feSJeremy L Thompson #include <string.h> 142b730f8bSJeremy L Thompson 1549aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h" 160d0321e0SJeremy L Thompson #include "../cuda/ceed-cuda-compile.h" 172b730f8bSJeremy L Thompson #include "ceed-cuda-ref.h" 180d0321e0SJeremy L Thompson 190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 200d0321e0SJeremy L Thompson // Basis apply - tensor 210d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 22db2becc9SJeremy L Thompson static int CeedBasisApplyCore_Cuda(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, 23db2becc9SJeremy L Thompson CeedVector u, CeedVector v) { 240d0321e0SJeremy L Thompson Ceed ceed; 25ca735530SJeremy L Thompson CeedInt Q_1d, dim; 267bbbfca3SJeremy L Thompson const CeedInt is_transpose = t_mode == CEED_TRANSPOSE; 27437930d1SJeremy L Thompson const int max_block_size = 32; 280d0321e0SJeremy L Thompson const CeedScalar *d_u; 290d0321e0SJeremy L Thompson CeedScalar *d_v; 30ca735530SJeremy L Thompson CeedBasis_Cuda *data; 31ca735530SJeremy L Thompson 32ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 33ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 34ca735530SJeremy L Thompson 359ea2cfd9SJeremy L Thompson // Get read/write access to u, v 366574a04fSJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 376574a04fSJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 38db2becc9SJeremy L Thompson if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 39db2becc9SJeremy L Thompson else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 400d0321e0SJeremy L Thompson 410d0321e0SJeremy L Thompson // Clear v for transpose operation 42db2becc9SJeremy L Thompson if (is_transpose && !apply_add) { 4319a04db8SJeremy L Thompson CeedInt num_comp, q_comp, num_nodes, num_qpts; 441f9221feSJeremy L Thompson CeedSize length; 45ca735530SJeremy L Thompson 4619a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 4719a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 4819a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes)); 4919a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 5019a04db8SJeremy L Thompson length = (CeedSize)num_elem * (CeedSize)num_comp * (t_mode == CEED_TRANSPOSE ? (CeedSize)num_nodes : ((CeedSize)num_qpts * (CeedSize)q_comp)); 512b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar))); 520d0321e0SJeremy L Thompson } 532b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 542b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 550d0321e0SJeremy L Thompson 560d0321e0SJeremy L Thompson // Basis action 57437930d1SJeremy L Thompson switch (eval_mode) { 580d0321e0SJeremy L Thompson case CEED_EVAL_INTERP: { 597bbbfca3SJeremy L Thompson void *interp_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_interp_1d, &d_u, &d_v}; 60b2165e7aSSebastian Grimberg const CeedInt block_size = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size); 610d0321e0SJeremy L Thompson 62eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernel_Cuda(ceed, data->Interp, num_elem, block_size, interp_args)); 630d0321e0SJeremy L Thompson } break; 640d0321e0SJeremy L Thompson case CEED_EVAL_GRAD: { 657bbbfca3SJeremy L Thompson void *grad_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_interp_1d, &data->d_grad_1d, &d_u, &d_v}; 66b2165e7aSSebastian Grimberg const CeedInt block_size = max_block_size; 670d0321e0SJeremy L Thompson 68eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernel_Cuda(ceed, data->Grad, num_elem, block_size, grad_args)); 690d0321e0SJeremy L Thompson } break; 700d0321e0SJeremy L Thompson case CEED_EVAL_WEIGHT: { 71097cc795SJames Wright CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights_1d not set", CeedEvalModes[eval_mode]); 72437930d1SJeremy L Thompson void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v}; 73b2165e7aSSebastian Grimberg const int block_size_x = Q_1d; 74b2165e7aSSebastian Grimberg const int block_size_y = dim >= 2 ? Q_1d : 1; 75b2165e7aSSebastian Grimberg 76b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, num_elem, block_size_x, block_size_y, 1, weight_args)); 770d0321e0SJeremy L Thompson } break; 789ea2cfd9SJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 799ea2cfd9SJeremy L Thompson break; 800d0321e0SJeremy L Thompson // LCOV_EXCL_START 810d0321e0SJeremy L Thompson case CEED_EVAL_DIV: 820d0321e0SJeremy L Thompson case CEED_EVAL_CURL: 83bcbe1c99SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 840d0321e0SJeremy L Thompson // LCOV_EXCL_STOP 850d0321e0SJeremy L Thompson } 860d0321e0SJeremy L Thompson 879ea2cfd9SJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 882b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 899ea2cfd9SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 909ea2cfd9SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 919bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 920d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 930d0321e0SJeremy L Thompson } 940d0321e0SJeremy L Thompson 95db2becc9SJeremy L Thompson static int CeedBasisApply_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, 96db2becc9SJeremy L Thompson CeedVector v) { 97db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyCore_Cuda(basis, false, num_elem, t_mode, eval_mode, u, v)); 98db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 99db2becc9SJeremy L Thompson } 100db2becc9SJeremy L Thompson 101db2becc9SJeremy L Thompson static int CeedBasisApplyAdd_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, 102db2becc9SJeremy L Thompson CeedVector v) { 103db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyCore_Cuda(basis, true, num_elem, t_mode, eval_mode, u, v)); 104db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 105db2becc9SJeremy L Thompson } 106db2becc9SJeremy L Thompson 1070d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 10834d14614SJeremy L Thompson // Basis apply - tensor AtPoints 10934d14614SJeremy L Thompson //------------------------------------------------------------------------------ 110db2becc9SJeremy L Thompson static int CeedBasisApplyAtPointsCore_Cuda(CeedBasis basis, bool apply_add, const CeedInt num_elem, const CeedInt *num_points, 111db2becc9SJeremy L Thompson CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) { 11234d14614SJeremy L Thompson Ceed ceed; 11334d14614SJeremy L Thompson CeedInt Q_1d, dim, max_num_points = num_points[0]; 11434d14614SJeremy L Thompson const CeedInt is_transpose = t_mode == CEED_TRANSPOSE; 11534d14614SJeremy L Thompson const int max_block_size = 32; 11634d14614SJeremy L Thompson const CeedScalar *d_x, *d_u; 11734d14614SJeremy L Thompson CeedScalar *d_v; 11834d14614SJeremy L Thompson CeedBasis_Cuda *data; 11934d14614SJeremy L Thompson 12034d14614SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 12134d14614SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 12234d14614SJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 12334d14614SJeremy L Thompson 12434d14614SJeremy L Thompson // Weight handled separately 12534d14614SJeremy L Thompson if (eval_mode == CEED_EVAL_WEIGHT) { 1265a5594ffSJeremy L Thompson CeedCallBackend(CeedVectorSetValue(v, 1.0)); 12734d14614SJeremy L Thompson return CEED_ERROR_SUCCESS; 12834d14614SJeremy L Thompson } 12934d14614SJeremy L Thompson 1309bc66399SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 1319bc66399SJeremy L Thompson 132111870feSJeremy L Thompson // Check padded to uniform number of points per elem 133111870feSJeremy L Thompson for (CeedInt i = 1; i < num_elem; i++) max_num_points = CeedIntMax(max_num_points, num_points[i]); 134111870feSJeremy L Thompson { 135111870feSJeremy L Thompson CeedInt num_comp, q_comp; 136111870feSJeremy L Thompson CeedSize len, len_required; 137111870feSJeremy L Thompson 138111870feSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 139111870feSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 140111870feSJeremy L Thompson CeedCallBackend(CeedVectorGetLength(is_transpose ? u : v, &len)); 141111870feSJeremy L Thompson len_required = (CeedSize)num_comp * (CeedSize)q_comp * (CeedSize)num_elem * (CeedSize)max_num_points; 142111870feSJeremy L Thompson CeedCheck(len >= len_required, ceed, CEED_ERROR_BACKEND, 143111870feSJeremy L Thompson "Vector at points must be padded to the same number of points in each element for BasisApplyAtPoints on GPU backends." 144111870feSJeremy L Thompson " Found %" CeedSize_FMT ", Required %" CeedSize_FMT, 145111870feSJeremy L Thompson len, len_required); 146111870feSJeremy L Thompson } 147111870feSJeremy L Thompson 148111870feSJeremy L Thompson // Move num_points array to device 149111870feSJeremy L Thompson if (is_transpose) { 150111870feSJeremy L Thompson const CeedInt num_bytes = num_elem * sizeof(CeedInt); 151111870feSJeremy L Thompson 152111870feSJeremy L Thompson if (num_elem != data->num_elem_at_points) { 153111870feSJeremy L Thompson data->num_elem_at_points = num_elem; 154111870feSJeremy L Thompson 155111870feSJeremy L Thompson if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem)); 156111870feSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_points_per_elem, num_bytes)); 157111870feSJeremy L Thompson CeedCallBackend(CeedFree(&data->h_points_per_elem)); 158111870feSJeremy L Thompson CeedCallBackend(CeedCalloc(num_elem, &data->h_points_per_elem)); 159111870feSJeremy L Thompson } 1609e511c80SJeremy L Thompson if (memcmp(data->h_points_per_elem, num_points, num_bytes)) { 161111870feSJeremy L Thompson memcpy(data->h_points_per_elem, num_points, num_bytes); 162111870feSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_points_per_elem, num_points, num_bytes, cudaMemcpyHostToDevice)); 163111870feSJeremy L Thompson } 164111870feSJeremy L Thompson } 165111870feSJeremy L Thompson 16634d14614SJeremy L Thompson // Build kernels if needed 16734d14614SJeremy L Thompson if (data->num_points != max_num_points) { 16834d14614SJeremy L Thompson CeedInt P_1d; 16934d14614SJeremy L Thompson 17034d14614SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 17134d14614SJeremy L Thompson data->num_points = max_num_points; 17234d14614SJeremy L Thompson 17334d14614SJeremy L Thompson // -- Create interp matrix to Chebyshev coefficients 17434d14614SJeremy L Thompson if (!data->d_chebyshev_interp_1d) { 17534d14614SJeremy L Thompson CeedSize interp_bytes; 17634d14614SJeremy L Thompson CeedScalar *chebyshev_interp_1d; 17734d14614SJeremy L Thompson 17834d14614SJeremy L Thompson interp_bytes = P_1d * Q_1d * sizeof(CeedScalar); 17934d14614SJeremy L Thompson CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d)); 1805a5594ffSJeremy L Thompson CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d)); 18134d14614SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_chebyshev_interp_1d, interp_bytes)); 18234d14614SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 18334d14614SJeremy L Thompson CeedCallBackend(CeedFree(&chebyshev_interp_1d)); 18434d14614SJeremy L Thompson } 18534d14614SJeremy L Thompson 18634d14614SJeremy L Thompson // -- Compile kernels 1879c25dd66SJeremy L Thompson const char basis_kernel_source[] = "// AtPoints basis source\n#include <ceed/jit-source/cuda/cuda-ref-basis-tensor-at-points.h>\n"; 18834d14614SJeremy L Thompson CeedInt num_comp; 18934d14614SJeremy L Thompson 19034d14614SJeremy L Thompson if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints)); 19134d14614SJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 19234d14614SJeremy L Thompson CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->moduleAtPoints, 9, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "BASIS_BUF_LEN", 193f7c9815fSJeremy L Thompson Q_1d * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim - 1), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, 19434d14614SJeremy L Thompson "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_NUM_PTS", 195f7c9815fSJeremy L Thompson max_num_points, "POINTS_BUFF_LEN", CeedIntPow(Q_1d, dim - 1))); 19634d14614SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpAtPoints", &data->InterpAtPoints)); 197*81ae6159SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpTransposeAtPoints", &data->InterpTransposeAtPoints)); 19834d14614SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradAtPoints", &data->GradAtPoints)); 199*81ae6159SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradTransposeAtPoints", &data->GradTransposeAtPoints)); 20034d14614SJeremy L Thompson } 20134d14614SJeremy L Thompson 20234d14614SJeremy L Thompson // Get read/write access to u, v 20334d14614SJeremy L Thompson CeedCallBackend(CeedVectorGetArrayRead(x_ref, CEED_MEM_DEVICE, &d_x)); 20434d14614SJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 20534d14614SJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 206db2becc9SJeremy L Thompson if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 207db2becc9SJeremy L Thompson else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 20834d14614SJeremy L Thompson 20934d14614SJeremy L Thompson // Clear v for transpose operation 210db2becc9SJeremy L Thompson if (is_transpose && !apply_add) { 21119a04db8SJeremy L Thompson CeedInt num_comp, q_comp, num_nodes; 21234d14614SJeremy L Thompson CeedSize length; 21334d14614SJeremy L Thompson 21419a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 21519a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 21619a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes)); 21719a04db8SJeremy L Thompson length = 21819a04db8SJeremy L Thompson (CeedSize)num_elem * (CeedSize)num_comp * (t_mode == CEED_TRANSPOSE ? (CeedSize)num_nodes : ((CeedSize)max_num_points * (CeedSize)q_comp)); 21934d14614SJeremy L Thompson CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar))); 22034d14614SJeremy L Thompson } 22134d14614SJeremy L Thompson 22234d14614SJeremy L Thompson // Basis action 22334d14614SJeremy L Thompson switch (eval_mode) { 22434d14614SJeremy L Thompson case CEED_EVAL_INTERP: { 225*81ae6159SJeremy L Thompson void *interp_args[] = {(void *)&num_elem, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v}; 22634d14614SJeremy L Thompson const CeedInt block_size = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size); 22734d14614SJeremy L Thompson 228*81ae6159SJeremy L Thompson CeedCallBackend( 229*81ae6159SJeremy L Thompson CeedRunKernel_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, num_elem, block_size, interp_args)); 23034d14614SJeremy L Thompson } break; 23134d14614SJeremy L Thompson case CEED_EVAL_GRAD: { 232*81ae6159SJeremy L Thompson void *grad_args[] = {(void *)&num_elem, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v}; 2332d10e82cSJeremy L Thompson const CeedInt block_size = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size); 23434d14614SJeremy L Thompson 235*81ae6159SJeremy L Thompson CeedCallBackend(CeedRunKernel_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, num_elem, block_size, grad_args)); 23634d14614SJeremy L Thompson } break; 23734d14614SJeremy L Thompson case CEED_EVAL_WEIGHT: 23834d14614SJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 23934d14614SJeremy L Thompson break; 24034d14614SJeremy L Thompson // LCOV_EXCL_START 24134d14614SJeremy L Thompson case CEED_EVAL_DIV: 24234d14614SJeremy L Thompson case CEED_EVAL_CURL: 24334d14614SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 24434d14614SJeremy L Thompson // LCOV_EXCL_STOP 24534d14614SJeremy L Thompson } 24634d14614SJeremy L Thompson 24734d14614SJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 24834d14614SJeremy L Thompson CeedCallBackend(CeedVectorRestoreArrayRead(x_ref, &d_x)); 24934d14614SJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 25034d14614SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 25134d14614SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 2529bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 25334d14614SJeremy L Thompson return CEED_ERROR_SUCCESS; 25434d14614SJeremy L Thompson } 25534d14614SJeremy L Thompson 256db2becc9SJeremy L Thompson static int CeedBasisApplyAtPoints_Cuda(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, 257db2becc9SJeremy L Thompson CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) { 258db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda(basis, false, num_elem, num_points, t_mode, eval_mode, x_ref, u, v)); 259db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 260db2becc9SJeremy L Thompson } 261db2becc9SJeremy L Thompson 262db2becc9SJeremy L Thompson static int CeedBasisApplyAddAtPoints_Cuda(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, 263db2becc9SJeremy L Thompson CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) { 264db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda(basis, true, num_elem, num_points, t_mode, eval_mode, x_ref, u, v)); 265db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 266db2becc9SJeremy L Thompson } 267db2becc9SJeremy L Thompson 26834d14614SJeremy L Thompson //------------------------------------------------------------------------------ 2690d0321e0SJeremy L Thompson // Basis apply - non-tensor 2700d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 271db2becc9SJeremy L Thompson static int CeedBasisApplyNonTensorCore_Cuda(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, 272db2becc9SJeremy L Thompson CeedVector u, CeedVector v) { 2730d0321e0SJeremy L Thompson Ceed ceed; 274437930d1SJeremy L Thompson CeedInt num_nodes, num_qpts; 2757bbbfca3SJeremy L Thompson const CeedInt is_transpose = t_mode == CEED_TRANSPOSE; 276d075f50bSSebastian Grimberg const int elems_per_block = 1; 277d075f50bSSebastian Grimberg const int grid = CeedDivUpInt(num_elem, elems_per_block); 2780d0321e0SJeremy L Thompson const CeedScalar *d_u; 2790d0321e0SJeremy L Thompson CeedScalar *d_v; 280ca735530SJeremy L Thompson CeedBasisNonTensor_Cuda *data; 281ca735530SJeremy L Thompson 282ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 283ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 284ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 285ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes)); 286ca735530SJeremy L Thompson 2879ea2cfd9SJeremy L Thompson // Get read/write access to u, v 2889ea2cfd9SJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 2899ea2cfd9SJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 290db2becc9SJeremy L Thompson if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 291db2becc9SJeremy L Thompson else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 2920d0321e0SJeremy L Thompson 2930d0321e0SJeremy L Thompson // Clear v for transpose operation 294db2becc9SJeremy L Thompson if (is_transpose && !apply_add) { 29519a04db8SJeremy L Thompson CeedInt num_comp, q_comp; 2961f9221feSJeremy L Thompson CeedSize length; 297ca735530SJeremy L Thompson 29819a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 29919a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 30019a04db8SJeremy L Thompson length = (CeedSize)num_elem * (CeedSize)num_comp * (t_mode == CEED_TRANSPOSE ? (CeedSize)num_nodes : ((CeedSize)num_qpts * (CeedSize)q_comp)); 3012b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar))); 3020d0321e0SJeremy L Thompson } 3030d0321e0SJeremy L Thompson 3040d0321e0SJeremy L Thompson // Apply basis operation 305437930d1SJeremy L Thompson switch (eval_mode) { 3060d0321e0SJeremy L Thompson case CEED_EVAL_INTERP: { 307d075f50bSSebastian Grimberg void *interp_args[] = {(void *)&num_elem, &data->d_interp, &d_u, &d_v}; 3087bbbfca3SJeremy L Thompson const int block_size_x = is_transpose ? num_nodes : num_qpts; 309b2165e7aSSebastian Grimberg 3107bbbfca3SJeremy L Thompson if (is_transpose) { 311d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->InterpTranspose, grid, block_size_x, 1, elems_per_block, interp_args)); 312d075f50bSSebastian Grimberg } else { 313b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Interp, grid, block_size_x, 1, elems_per_block, interp_args)); 314d075f50bSSebastian Grimberg } 3150d0321e0SJeremy L Thompson } break; 3160d0321e0SJeremy L Thompson case CEED_EVAL_GRAD: { 317d075f50bSSebastian Grimberg void *grad_args[] = {(void *)&num_elem, &data->d_grad, &d_u, &d_v}; 3187bbbfca3SJeremy L Thompson const int block_size_x = is_transpose ? num_nodes : num_qpts; 319b2165e7aSSebastian Grimberg 3207bbbfca3SJeremy L Thompson if (is_transpose) { 321d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, grad_args)); 322d075f50bSSebastian Grimberg } else { 323d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, grad_args)); 324d075f50bSSebastian Grimberg } 325d075f50bSSebastian Grimberg } break; 326d075f50bSSebastian Grimberg case CEED_EVAL_DIV: { 327d075f50bSSebastian Grimberg void *div_args[] = {(void *)&num_elem, &data->d_div, &d_u, &d_v}; 3287bbbfca3SJeremy L Thompson const int block_size_x = is_transpose ? num_nodes : num_qpts; 329d075f50bSSebastian Grimberg 3307bbbfca3SJeremy L Thompson if (is_transpose) { 331d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, div_args)); 332d075f50bSSebastian Grimberg } else { 333d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, div_args)); 334d075f50bSSebastian Grimberg } 335d075f50bSSebastian Grimberg } break; 336d075f50bSSebastian Grimberg case CEED_EVAL_CURL: { 337d075f50bSSebastian Grimberg void *curl_args[] = {(void *)&num_elem, &data->d_curl, &d_u, &d_v}; 3387bbbfca3SJeremy L Thompson const int block_size_x = is_transpose ? num_nodes : num_qpts; 339d075f50bSSebastian Grimberg 3407bbbfca3SJeremy L Thompson if (is_transpose) { 341d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, curl_args)); 342d075f50bSSebastian Grimberg } else { 343d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, curl_args)); 344d075f50bSSebastian Grimberg } 3450d0321e0SJeremy L Thompson } break; 3460d0321e0SJeremy L Thompson case CEED_EVAL_WEIGHT: { 347097cc795SJames Wright CeedCheck(data->d_q_weight, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights not set", CeedEvalModes[eval_mode]); 348437930d1SJeremy L Thompson void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight, &d_v}; 349b2165e7aSSebastian Grimberg 350eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, num_qpts, 1, elems_per_block, weight_args)); 3510d0321e0SJeremy L Thompson } break; 3529ea2cfd9SJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 3539ea2cfd9SJeremy L Thompson break; 3540d0321e0SJeremy L Thompson } 3550d0321e0SJeremy L Thompson 3569ea2cfd9SJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 3572b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 3589ea2cfd9SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 3599ea2cfd9SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 3609bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 3610d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3620d0321e0SJeremy L Thompson } 3630d0321e0SJeremy L Thompson 364db2becc9SJeremy L Thompson static int CeedBasisApplyNonTensor_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, 365db2becc9SJeremy L Thompson CeedVector v) { 366db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda(basis, false, num_elem, t_mode, eval_mode, u, v)); 367db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 368db2becc9SJeremy L Thompson } 369db2becc9SJeremy L Thompson 370db2becc9SJeremy L Thompson static int CeedBasisApplyAddNonTensor_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, 371db2becc9SJeremy L Thompson CeedVector v) { 372db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda(basis, true, num_elem, t_mode, eval_mode, u, v)); 373db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 374db2becc9SJeremy L Thompson } 375db2becc9SJeremy L Thompson 3760d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 3770d0321e0SJeremy L Thompson // Destroy tensor basis 3780d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 3790d0321e0SJeremy L Thompson static int CeedBasisDestroy_Cuda(CeedBasis basis) { 3800d0321e0SJeremy L Thompson Ceed ceed; 3810d0321e0SJeremy L Thompson CeedBasis_Cuda *data; 382ca735530SJeremy L Thompson 383ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 3842b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 3852b730f8bSJeremy L Thompson CeedCallCuda(ceed, cuModuleUnload(data->module)); 38634d14614SJeremy L Thompson if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints)); 387097cc795SJames Wright if (data->d_q_weight_1d) CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d)); 388111870feSJeremy L Thompson CeedCallBackend(CeedFree(&data->h_points_per_elem)); 389111870feSJeremy L Thompson if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem)); 3902b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_interp_1d)); 3912b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_grad_1d)); 39234d14614SJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_chebyshev_interp_1d)); 3932b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&data)); 3949bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 3950d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3960d0321e0SJeremy L Thompson } 3970d0321e0SJeremy L Thompson 3980d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 3990d0321e0SJeremy L Thompson // Destroy non-tensor basis 4000d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 4010d0321e0SJeremy L Thompson static int CeedBasisDestroyNonTensor_Cuda(CeedBasis basis) { 4020d0321e0SJeremy L Thompson Ceed ceed; 4030d0321e0SJeremy L Thompson CeedBasisNonTensor_Cuda *data; 404ca735530SJeremy L Thompson 405ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 4062b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 4072b730f8bSJeremy L Thompson CeedCallCuda(ceed, cuModuleUnload(data->module)); 408097cc795SJames Wright if (data->d_q_weight) CeedCallCuda(ceed, cudaFree(data->d_q_weight)); 4092b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_interp)); 4102b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_grad)); 411d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaFree(data->d_div)); 412d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaFree(data->d_curl)); 4132b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&data)); 4149bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 4150d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 4160d0321e0SJeremy L Thompson } 4170d0321e0SJeremy L Thompson 4180d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 4190d0321e0SJeremy L Thompson // Create tensor 4200d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 4212b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d, 4222b730f8bSJeremy L Thompson const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) { 4230d0321e0SJeremy L Thompson Ceed ceed; 424ca735530SJeremy L Thompson CeedInt num_comp; 425ca735530SJeremy L Thompson const CeedInt q_bytes = Q_1d * sizeof(CeedScalar); 426ca735530SJeremy L Thompson const CeedInt interp_bytes = q_bytes * P_1d; 4270d0321e0SJeremy L Thompson CeedBasis_Cuda *data; 428ca735530SJeremy L Thompson 429ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 4302b730f8bSJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 4310d0321e0SJeremy L Thompson 4320d0321e0SJeremy L Thompson // Copy data to GPU 433097cc795SJames Wright if (q_weight_1d) { 4342b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes)); 4352b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice)); 436097cc795SJames Wright } 4372b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes)); 4382b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 4392b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes)); 4402b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice)); 4410d0321e0SJeremy L Thompson 442ecc88aebSJeremy L Thompson // Compile basis kernels 4439c25dd66SJeremy L Thompson const char basis_kernel_source[] = "// Tensor basis source\n#include <ceed/jit-source/cuda/cuda-ref-basis-tensor.h>\n"; 4449c25dd66SJeremy L Thompson 4452b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 446eb7e6cafSJeremy 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", 447f7c9815fSJeremy L Thompson Q_1d * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim - 1), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, 4482b730f8bSJeremy L Thompson "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim))); 449eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 450eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad)); 451eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 452437930d1SJeremy L Thompson 4532b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 4540d0321e0SJeremy L Thompson 455d075f50bSSebastian Grimberg // Register backend functions 4562b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApply_Cuda)); 457db2becc9SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAdd_Cuda)); 45834d14614SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAtPoints", CeedBasisApplyAtPoints_Cuda)); 459db2becc9SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAddAtPoints", CeedBasisApplyAddAtPoints_Cuda)); 4602b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda)); 4619bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 4620d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 4630d0321e0SJeremy L Thompson } 4640d0321e0SJeremy L Thompson 4650d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 466d075f50bSSebastian Grimberg // Create non-tensor H^1 4670d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 4682b730f8bSJeremy L Thompson int CeedBasisCreateH1_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *grad, 46951475c7cSJeremy L Thompson const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 4700d0321e0SJeremy L Thompson Ceed ceed; 471d075f50bSSebastian Grimberg CeedInt num_comp, q_comp_interp, q_comp_grad; 472ca735530SJeremy L Thompson const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 4730d0321e0SJeremy L Thompson CeedBasisNonTensor_Cuda *data; 474ca735530SJeremy L Thompson 475ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 4762b730f8bSJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 4770d0321e0SJeremy L Thompson 4780d0321e0SJeremy L Thompson // Copy basis data to GPU 479d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 480d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad)); 481097cc795SJames Wright if (q_weight) { 4822b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes)); 4832b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice)); 484097cc795SJames Wright } 485d075f50bSSebastian Grimberg if (interp) { 486d075f50bSSebastian Grimberg const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 487d075f50bSSebastian Grimberg 4882b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes)); 4892b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice)); 490d075f50bSSebastian Grimberg } 491d075f50bSSebastian Grimberg if (grad) { 492d075f50bSSebastian Grimberg const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad; 493d075f50bSSebastian Grimberg 4942b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad, grad_bytes)); 4952b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad, grad, grad_bytes, cudaMemcpyHostToDevice)); 496d075f50bSSebastian Grimberg } 4970d0321e0SJeremy L Thompson 4980d0321e0SJeremy L Thompson // Compile basis kernels 4999c25dd66SJeremy L Thompson const char basis_kernel_source[] = "// Nontensor basis source\n#include <ceed/jit-source/cuda/cuda-ref-basis-nontensor.h>\n"; 5009c25dd66SJeremy L Thompson 5012b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 502d075f50bSSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP", 503d075f50bSSebastian Grimberg q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_grad, "BASIS_NUM_COMP", num_comp)); 504d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 505d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 506d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv)); 507d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose)); 508d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 509d075f50bSSebastian Grimberg 510d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisSetData(basis, data)); 511d075f50bSSebastian Grimberg 512d075f50bSSebastian Grimberg // Register backend functions 513d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda)); 514db2becc9SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Cuda)); 515d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda)); 5169bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 517d075f50bSSebastian Grimberg return CEED_ERROR_SUCCESS; 518d075f50bSSebastian Grimberg } 519d075f50bSSebastian Grimberg 520d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 521d075f50bSSebastian Grimberg // Create non-tensor H(div) 522d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 523d075f50bSSebastian Grimberg int CeedBasisCreateHdiv_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *div, 524d075f50bSSebastian Grimberg const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 525d075f50bSSebastian Grimberg Ceed ceed; 526d075f50bSSebastian Grimberg CeedInt num_comp, q_comp_interp, q_comp_div; 527d075f50bSSebastian Grimberg const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 528d075f50bSSebastian Grimberg CeedBasisNonTensor_Cuda *data; 529d075f50bSSebastian Grimberg 530d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 531d075f50bSSebastian Grimberg CeedCallBackend(CeedCalloc(1, &data)); 532d075f50bSSebastian Grimberg 533d075f50bSSebastian Grimberg // Copy basis data to GPU 534d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 535d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div)); 536097cc795SJames Wright if (q_weight) { 537d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes)); 538d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice)); 539097cc795SJames Wright } 540d075f50bSSebastian Grimberg if (interp) { 541d075f50bSSebastian Grimberg const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 542d075f50bSSebastian Grimberg 543d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes)); 544d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice)); 545d075f50bSSebastian Grimberg } 546d075f50bSSebastian Grimberg if (div) { 547d075f50bSSebastian Grimberg const CeedInt div_bytes = q_bytes * num_nodes * q_comp_div; 548d075f50bSSebastian Grimberg 549d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_div, div_bytes)); 550d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_div, div, div_bytes, cudaMemcpyHostToDevice)); 551d075f50bSSebastian Grimberg } 552d075f50bSSebastian Grimberg 553d075f50bSSebastian Grimberg // Compile basis kernels 5549c25dd66SJeremy L Thompson const char basis_kernel_source[] = "// Nontensor basis source\n#include <ceed/jit-source/cuda/cuda-ref-basis-nontensor.h>\n"; 5559c25dd66SJeremy L Thompson 556d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 557d075f50bSSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP", 558d075f50bSSebastian Grimberg q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_div, "BASIS_NUM_COMP", num_comp)); 559d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 560d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 561d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv)); 562d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose)); 563d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 564d075f50bSSebastian Grimberg 565d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisSetData(basis, data)); 566d075f50bSSebastian Grimberg 567d075f50bSSebastian Grimberg // Register backend functions 568d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda)); 569db2becc9SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Cuda)); 570d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda)); 5719bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 572d075f50bSSebastian Grimberg return CEED_ERROR_SUCCESS; 573d075f50bSSebastian Grimberg } 574d075f50bSSebastian Grimberg 575d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 576d075f50bSSebastian Grimberg // Create non-tensor H(curl) 577d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 578d075f50bSSebastian Grimberg int CeedBasisCreateHcurl_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 579d075f50bSSebastian Grimberg const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 580d075f50bSSebastian Grimberg Ceed ceed; 581d075f50bSSebastian Grimberg CeedInt num_comp, q_comp_interp, q_comp_curl; 582d075f50bSSebastian Grimberg const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 583d075f50bSSebastian Grimberg CeedBasisNonTensor_Cuda *data; 584d075f50bSSebastian Grimberg 585d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 586d075f50bSSebastian Grimberg CeedCallBackend(CeedCalloc(1, &data)); 587d075f50bSSebastian Grimberg 588d075f50bSSebastian Grimberg // Copy basis data to GPU 589d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 590d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl)); 591097cc795SJames Wright if (q_weight) { 592d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes)); 593d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice)); 594097cc795SJames Wright } 595d075f50bSSebastian Grimberg if (interp) { 596d075f50bSSebastian Grimberg const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 597d075f50bSSebastian Grimberg 598d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes)); 599d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice)); 600d075f50bSSebastian Grimberg } 601d075f50bSSebastian Grimberg if (curl) { 602d075f50bSSebastian Grimberg const CeedInt curl_bytes = q_bytes * num_nodes * q_comp_curl; 603d075f50bSSebastian Grimberg 604d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_curl, curl_bytes)); 605d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_curl, curl, curl_bytes, cudaMemcpyHostToDevice)); 606d075f50bSSebastian Grimberg } 607d075f50bSSebastian Grimberg 608d075f50bSSebastian Grimberg // Compile basis kernels 6099c25dd66SJeremy L Thompson const char basis_kernel_source[] = "// Nontensor basis source\n#include <ceed/jit-source/cuda/cuda-ref-basis-nontensor.h>\n"; 6109c25dd66SJeremy L Thompson 611d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 612d075f50bSSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP", 613d075f50bSSebastian Grimberg q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_curl, "BASIS_NUM_COMP", num_comp)); 614d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 615d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 616d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv)); 617d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose)); 618d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 6190d0321e0SJeremy L Thompson 6202b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 6210d0321e0SJeremy L Thompson 6220d0321e0SJeremy L Thompson // Register backend functions 6232b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda)); 624db2becc9SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Cuda)); 6252b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda)); 6269bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 6270d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 6280d0321e0SJeremy L Thompson } 6292a86cc9dSSebastian Grimberg 6300d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 631