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)); 910d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 920d0321e0SJeremy L Thompson } 930d0321e0SJeremy L Thompson 94db2becc9SJeremy L Thompson static int CeedBasisApply_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, 95db2becc9SJeremy L Thompson CeedVector v) { 96db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyCore_Cuda(basis, false, num_elem, t_mode, eval_mode, u, v)); 97db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 98db2becc9SJeremy L Thompson } 99db2becc9SJeremy L Thompson 100db2becc9SJeremy L Thompson static int CeedBasisApplyAdd_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, 101db2becc9SJeremy L Thompson CeedVector v) { 102db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyCore_Cuda(basis, true, num_elem, t_mode, eval_mode, u, v)); 103db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 104db2becc9SJeremy L Thompson } 105db2becc9SJeremy L Thompson 1060d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 10734d14614SJeremy L Thompson // Basis apply - tensor AtPoints 10834d14614SJeremy L Thompson //------------------------------------------------------------------------------ 109db2becc9SJeremy L Thompson static int CeedBasisApplyAtPointsCore_Cuda(CeedBasis basis, bool apply_add, const CeedInt num_elem, const CeedInt *num_points, 110db2becc9SJeremy L Thompson CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) { 11134d14614SJeremy L Thompson Ceed ceed; 11234d14614SJeremy L Thompson CeedInt Q_1d, dim, max_num_points = num_points[0]; 11334d14614SJeremy L Thompson const CeedInt is_transpose = t_mode == CEED_TRANSPOSE; 11434d14614SJeremy L Thompson const int max_block_size = 32; 11534d14614SJeremy L Thompson const CeedScalar *d_x, *d_u; 11634d14614SJeremy L Thompson CeedScalar *d_v; 11734d14614SJeremy L Thompson CeedBasis_Cuda *data; 11834d14614SJeremy L Thompson 11934d14614SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 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 130111870feSJeremy L Thompson // Check padded to uniform number of points per elem 131111870feSJeremy L Thompson for (CeedInt i = 1; i < num_elem; i++) max_num_points = CeedIntMax(max_num_points, num_points[i]); 132111870feSJeremy L Thompson { 133111870feSJeremy L Thompson CeedInt num_comp, q_comp; 134111870feSJeremy L Thompson CeedSize len, len_required; 135111870feSJeremy L Thompson 136111870feSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 137111870feSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 138111870feSJeremy L Thompson CeedCallBackend(CeedVectorGetLength(is_transpose ? u : v, &len)); 139111870feSJeremy L Thompson len_required = (CeedSize)num_comp * (CeedSize)q_comp * (CeedSize)num_elem * (CeedSize)max_num_points; 140111870feSJeremy L Thompson CeedCheck(len >= len_required, ceed, CEED_ERROR_BACKEND, 141111870feSJeremy L Thompson "Vector at points must be padded to the same number of points in each element for BasisApplyAtPoints on GPU backends." 142111870feSJeremy L Thompson " Found %" CeedSize_FMT ", Required %" CeedSize_FMT, 143111870feSJeremy L Thompson len, len_required); 144111870feSJeremy L Thompson } 145111870feSJeremy L Thompson 146111870feSJeremy L Thompson // Move num_points array to device 147111870feSJeremy L Thompson if (is_transpose) { 148111870feSJeremy L Thompson const CeedInt num_bytes = num_elem * sizeof(CeedInt); 149111870feSJeremy L Thompson 150111870feSJeremy L Thompson if (num_elem != data->num_elem_at_points) { 151111870feSJeremy L Thompson data->num_elem_at_points = num_elem; 152111870feSJeremy L Thompson 153111870feSJeremy L Thompson if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem)); 154111870feSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_points_per_elem, num_bytes)); 155111870feSJeremy L Thompson CeedCallBackend(CeedFree(&data->h_points_per_elem)); 156111870feSJeremy L Thompson CeedCallBackend(CeedCalloc(num_elem, &data->h_points_per_elem)); 157111870feSJeremy L Thompson } 1589e511c80SJeremy L Thompson if (memcmp(data->h_points_per_elem, num_points, num_bytes)) { 159111870feSJeremy L Thompson memcpy(data->h_points_per_elem, num_points, num_bytes); 160111870feSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_points_per_elem, num_points, num_bytes, cudaMemcpyHostToDevice)); 161111870feSJeremy L Thompson } 162111870feSJeremy L Thompson } 163111870feSJeremy L Thompson 16434d14614SJeremy L Thompson // Build kernels if needed 16534d14614SJeremy L Thompson if (data->num_points != max_num_points) { 16634d14614SJeremy L Thompson CeedInt P_1d; 16734d14614SJeremy L Thompson 16834d14614SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 16934d14614SJeremy L Thompson data->num_points = max_num_points; 17034d14614SJeremy L Thompson 17134d14614SJeremy L Thompson // -- Create interp matrix to Chebyshev coefficients 17234d14614SJeremy L Thompson if (!data->d_chebyshev_interp_1d) { 17334d14614SJeremy L Thompson CeedSize interp_bytes; 17434d14614SJeremy L Thompson CeedScalar *chebyshev_interp_1d; 17534d14614SJeremy L Thompson 17634d14614SJeremy L Thompson interp_bytes = P_1d * Q_1d * sizeof(CeedScalar); 17734d14614SJeremy L Thompson CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d)); 1785a5594ffSJeremy L Thompson CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d)); 17934d14614SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_chebyshev_interp_1d, interp_bytes)); 18034d14614SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 18134d14614SJeremy L Thompson CeedCallBackend(CeedFree(&chebyshev_interp_1d)); 18234d14614SJeremy L Thompson } 18334d14614SJeremy L Thompson 18434d14614SJeremy L Thompson // -- Compile kernels 185*9c25dd66SJeremy L Thompson const char basis_kernel_source[] = "// AtPoints basis source\n#include <ceed/jit-source/cuda/cuda-ref-basis-tensor-at-points.h>\n"; 18634d14614SJeremy L Thompson CeedInt num_comp; 18734d14614SJeremy L Thompson 18834d14614SJeremy L Thompson if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints)); 18934d14614SJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 19034d14614SJeremy 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", 191f7c9815fSJeremy L Thompson Q_1d * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim - 1), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, 19234d14614SJeremy L Thompson "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_NUM_PTS", 193f7c9815fSJeremy L Thompson max_num_points, "POINTS_BUFF_LEN", CeedIntPow(Q_1d, dim - 1))); 19434d14614SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpAtPoints", &data->InterpAtPoints)); 19534d14614SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradAtPoints", &data->GradAtPoints)); 19634d14614SJeremy L Thompson } 19734d14614SJeremy L Thompson 19834d14614SJeremy L Thompson // Get read/write access to u, v 19934d14614SJeremy L Thompson CeedCallBackend(CeedVectorGetArrayRead(x_ref, CEED_MEM_DEVICE, &d_x)); 20034d14614SJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 20134d14614SJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 202db2becc9SJeremy L Thompson if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 203db2becc9SJeremy L Thompson else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 20434d14614SJeremy L Thompson 20534d14614SJeremy L Thompson // Clear v for transpose operation 206db2becc9SJeremy L Thompson if (is_transpose && !apply_add) { 20719a04db8SJeremy L Thompson CeedInt num_comp, q_comp, num_nodes; 20834d14614SJeremy L Thompson CeedSize length; 20934d14614SJeremy L Thompson 21019a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 21119a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 21219a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes)); 21319a04db8SJeremy L Thompson length = 21419a04db8SJeremy L Thompson (CeedSize)num_elem * (CeedSize)num_comp * (t_mode == CEED_TRANSPOSE ? (CeedSize)num_nodes : ((CeedSize)max_num_points * (CeedSize)q_comp)); 21534d14614SJeremy L Thompson CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar))); 21634d14614SJeremy L Thompson } 21734d14614SJeremy L Thompson 21834d14614SJeremy L Thompson // Basis action 21934d14614SJeremy L Thompson switch (eval_mode) { 22034d14614SJeremy L Thompson case CEED_EVAL_INTERP: { 221111870feSJeremy L Thompson void *interp_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v}; 22234d14614SJeremy L Thompson const CeedInt block_size = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size); 22334d14614SJeremy L Thompson 22434d14614SJeremy L Thompson CeedCallBackend(CeedRunKernel_Cuda(ceed, data->InterpAtPoints, num_elem, block_size, interp_args)); 22534d14614SJeremy L Thompson } break; 22634d14614SJeremy L Thompson case CEED_EVAL_GRAD: { 227111870feSJeremy L Thompson void *grad_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v}; 2282d10e82cSJeremy L Thompson const CeedInt block_size = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size); 22934d14614SJeremy L Thompson 23034d14614SJeremy L Thompson CeedCallBackend(CeedRunKernel_Cuda(ceed, data->GradAtPoints, num_elem, block_size, grad_args)); 23134d14614SJeremy L Thompson } break; 23234d14614SJeremy L Thompson case CEED_EVAL_WEIGHT: 23334d14614SJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 23434d14614SJeremy L Thompson break; 23534d14614SJeremy L Thompson // LCOV_EXCL_START 23634d14614SJeremy L Thompson case CEED_EVAL_DIV: 23734d14614SJeremy L Thompson case CEED_EVAL_CURL: 23834d14614SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 23934d14614SJeremy L Thompson // LCOV_EXCL_STOP 24034d14614SJeremy L Thompson } 24134d14614SJeremy L Thompson 24234d14614SJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 24334d14614SJeremy L Thompson CeedCallBackend(CeedVectorRestoreArrayRead(x_ref, &d_x)); 24434d14614SJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 24534d14614SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 24634d14614SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 24734d14614SJeremy L Thompson return CEED_ERROR_SUCCESS; 24834d14614SJeremy L Thompson } 24934d14614SJeremy L Thompson 250db2becc9SJeremy L Thompson static int CeedBasisApplyAtPoints_Cuda(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, 251db2becc9SJeremy L Thompson CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) { 252db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda(basis, false, num_elem, num_points, t_mode, eval_mode, x_ref, u, v)); 253db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 254db2becc9SJeremy L Thompson } 255db2becc9SJeremy L Thompson 256db2becc9SJeremy L Thompson static int CeedBasisApplyAddAtPoints_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, true, 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 26234d14614SJeremy L Thompson //------------------------------------------------------------------------------ 2630d0321e0SJeremy L Thompson // Basis apply - non-tensor 2640d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 265db2becc9SJeremy L Thompson static int CeedBasisApplyNonTensorCore_Cuda(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, 266db2becc9SJeremy L Thompson CeedVector u, CeedVector v) { 2670d0321e0SJeremy L Thompson Ceed ceed; 268437930d1SJeremy L Thompson CeedInt num_nodes, num_qpts; 2697bbbfca3SJeremy L Thompson const CeedInt is_transpose = t_mode == CEED_TRANSPOSE; 270d075f50bSSebastian Grimberg const int elems_per_block = 1; 271d075f50bSSebastian Grimberg const int grid = CeedDivUpInt(num_elem, elems_per_block); 2720d0321e0SJeremy L Thompson const CeedScalar *d_u; 2730d0321e0SJeremy L Thompson CeedScalar *d_v; 274ca735530SJeremy L Thompson CeedBasisNonTensor_Cuda *data; 275ca735530SJeremy L Thompson 276ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 277ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 278ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 279ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes)); 280ca735530SJeremy L Thompson 2819ea2cfd9SJeremy L Thompson // Get read/write access to u, v 2829ea2cfd9SJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 2839ea2cfd9SJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 284db2becc9SJeremy L Thompson if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 285db2becc9SJeremy L Thompson else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 2860d0321e0SJeremy L Thompson 2870d0321e0SJeremy L Thompson // Clear v for transpose operation 288db2becc9SJeremy L Thompson if (is_transpose && !apply_add) { 28919a04db8SJeremy L Thompson CeedInt num_comp, q_comp; 2901f9221feSJeremy L Thompson CeedSize length; 291ca735530SJeremy L Thompson 29219a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 29319a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 29419a04db8SJeremy L Thompson length = (CeedSize)num_elem * (CeedSize)num_comp * (t_mode == CEED_TRANSPOSE ? (CeedSize)num_nodes : ((CeedSize)num_qpts * (CeedSize)q_comp)); 2952b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar))); 2960d0321e0SJeremy L Thompson } 2970d0321e0SJeremy L Thompson 2980d0321e0SJeremy L Thompson // Apply basis operation 299437930d1SJeremy L Thompson switch (eval_mode) { 3000d0321e0SJeremy L Thompson case CEED_EVAL_INTERP: { 301d075f50bSSebastian Grimberg void *interp_args[] = {(void *)&num_elem, &data->d_interp, &d_u, &d_v}; 3027bbbfca3SJeremy L Thompson const int block_size_x = is_transpose ? num_nodes : num_qpts; 303b2165e7aSSebastian Grimberg 3047bbbfca3SJeremy L Thompson if (is_transpose) { 305d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->InterpTranspose, grid, block_size_x, 1, elems_per_block, interp_args)); 306d075f50bSSebastian Grimberg } else { 307b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Interp, grid, block_size_x, 1, elems_per_block, interp_args)); 308d075f50bSSebastian Grimberg } 3090d0321e0SJeremy L Thompson } break; 3100d0321e0SJeremy L Thompson case CEED_EVAL_GRAD: { 311d075f50bSSebastian Grimberg void *grad_args[] = {(void *)&num_elem, &data->d_grad, &d_u, &d_v}; 3127bbbfca3SJeremy L Thompson const int block_size_x = is_transpose ? num_nodes : num_qpts; 313b2165e7aSSebastian Grimberg 3147bbbfca3SJeremy L Thompson if (is_transpose) { 315d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, grad_args)); 316d075f50bSSebastian Grimberg } else { 317d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, grad_args)); 318d075f50bSSebastian Grimberg } 319d075f50bSSebastian Grimberg } break; 320d075f50bSSebastian Grimberg case CEED_EVAL_DIV: { 321d075f50bSSebastian Grimberg void *div_args[] = {(void *)&num_elem, &data->d_div, &d_u, &d_v}; 3227bbbfca3SJeremy L Thompson const int block_size_x = is_transpose ? num_nodes : num_qpts; 323d075f50bSSebastian Grimberg 3247bbbfca3SJeremy L Thompson if (is_transpose) { 325d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, div_args)); 326d075f50bSSebastian Grimberg } else { 327d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, div_args)); 328d075f50bSSebastian Grimberg } 329d075f50bSSebastian Grimberg } break; 330d075f50bSSebastian Grimberg case CEED_EVAL_CURL: { 331d075f50bSSebastian Grimberg void *curl_args[] = {(void *)&num_elem, &data->d_curl, &d_u, &d_v}; 3327bbbfca3SJeremy L Thompson const int block_size_x = is_transpose ? num_nodes : num_qpts; 333d075f50bSSebastian Grimberg 3347bbbfca3SJeremy L Thompson if (is_transpose) { 335d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, curl_args)); 336d075f50bSSebastian Grimberg } else { 337d075f50bSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, curl_args)); 338d075f50bSSebastian Grimberg } 3390d0321e0SJeremy L Thompson } break; 3400d0321e0SJeremy L Thompson case CEED_EVAL_WEIGHT: { 341097cc795SJames Wright CeedCheck(data->d_q_weight, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights not set", CeedEvalModes[eval_mode]); 342437930d1SJeremy L Thompson void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight, &d_v}; 343b2165e7aSSebastian Grimberg 344eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, num_qpts, 1, elems_per_block, weight_args)); 3450d0321e0SJeremy L Thompson } break; 3469ea2cfd9SJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 3479ea2cfd9SJeremy L Thompson break; 3480d0321e0SJeremy L Thompson } 3490d0321e0SJeremy L Thompson 3509ea2cfd9SJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 3512b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 3529ea2cfd9SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 3539ea2cfd9SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 3540d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3550d0321e0SJeremy L Thompson } 3560d0321e0SJeremy L Thompson 357db2becc9SJeremy L Thompson static int CeedBasisApplyNonTensor_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, 358db2becc9SJeremy L Thompson CeedVector v) { 359db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda(basis, false, num_elem, t_mode, eval_mode, u, v)); 360db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 361db2becc9SJeremy L Thompson } 362db2becc9SJeremy L Thompson 363db2becc9SJeremy L Thompson static int CeedBasisApplyAddNonTensor_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, 364db2becc9SJeremy L Thompson CeedVector v) { 365db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda(basis, true, num_elem, t_mode, eval_mode, u, v)); 366db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 367db2becc9SJeremy L Thompson } 368db2becc9SJeremy L Thompson 3690d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 3700d0321e0SJeremy L Thompson // Destroy tensor basis 3710d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 3720d0321e0SJeremy L Thompson static int CeedBasisDestroy_Cuda(CeedBasis basis) { 3730d0321e0SJeremy L Thompson Ceed ceed; 3740d0321e0SJeremy L Thompson CeedBasis_Cuda *data; 375ca735530SJeremy L Thompson 376ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 3772b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 3782b730f8bSJeremy L Thompson CeedCallCuda(ceed, cuModuleUnload(data->module)); 37934d14614SJeremy L Thompson if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints)); 380097cc795SJames Wright if (data->d_q_weight_1d) CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d)); 381111870feSJeremy L Thompson CeedCallBackend(CeedFree(&data->h_points_per_elem)); 382111870feSJeremy L Thompson if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem)); 3832b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_interp_1d)); 3842b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_grad_1d)); 38534d14614SJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_chebyshev_interp_1d)); 3862b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&data)); 3870d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3880d0321e0SJeremy L Thompson } 3890d0321e0SJeremy L Thompson 3900d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 3910d0321e0SJeremy L Thompson // Destroy non-tensor basis 3920d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 3930d0321e0SJeremy L Thompson static int CeedBasisDestroyNonTensor_Cuda(CeedBasis basis) { 3940d0321e0SJeremy L Thompson Ceed ceed; 3950d0321e0SJeremy L Thompson CeedBasisNonTensor_Cuda *data; 396ca735530SJeremy L Thompson 397ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 3982b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 3992b730f8bSJeremy L Thompson CeedCallCuda(ceed, cuModuleUnload(data->module)); 400097cc795SJames Wright if (data->d_q_weight) CeedCallCuda(ceed, cudaFree(data->d_q_weight)); 4012b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_interp)); 4022b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_grad)); 403d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaFree(data->d_div)); 404d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaFree(data->d_curl)); 4052b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&data)); 4060d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 4070d0321e0SJeremy L Thompson } 4080d0321e0SJeremy L Thompson 4090d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 4100d0321e0SJeremy L Thompson // Create tensor 4110d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 4122b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d, 4132b730f8bSJeremy L Thompson const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) { 4140d0321e0SJeremy L Thompson Ceed ceed; 415ca735530SJeremy L Thompson CeedInt num_comp; 416ca735530SJeremy L Thompson const CeedInt q_bytes = Q_1d * sizeof(CeedScalar); 417ca735530SJeremy L Thompson const CeedInt interp_bytes = q_bytes * P_1d; 4180d0321e0SJeremy L Thompson CeedBasis_Cuda *data; 419ca735530SJeremy L Thompson 420ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 4212b730f8bSJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 4220d0321e0SJeremy L Thompson 4230d0321e0SJeremy L Thompson // Copy data to GPU 424097cc795SJames Wright if (q_weight_1d) { 4252b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes)); 4262b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice)); 427097cc795SJames Wright } 4282b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes)); 4292b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 4302b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes)); 4312b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice)); 4320d0321e0SJeremy L Thompson 433ecc88aebSJeremy L Thompson // Compile basis kernels 434*9c25dd66SJeremy L Thompson const char basis_kernel_source[] = "// Tensor basis source\n#include <ceed/jit-source/cuda/cuda-ref-basis-tensor.h>\n"; 435*9c25dd66SJeremy L Thompson 4362b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 437eb7e6cafSJeremy 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", 438f7c9815fSJeremy L Thompson Q_1d * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim - 1), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, 4392b730f8bSJeremy L Thompson "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim))); 440eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 441eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad)); 442eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 443437930d1SJeremy L Thompson 4442b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 4450d0321e0SJeremy L Thompson 446d075f50bSSebastian Grimberg // Register backend functions 4472b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApply_Cuda)); 448db2becc9SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAdd_Cuda)); 44934d14614SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAtPoints", CeedBasisApplyAtPoints_Cuda)); 450db2becc9SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAddAtPoints", CeedBasisApplyAddAtPoints_Cuda)); 4512b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda)); 4520d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 4530d0321e0SJeremy L Thompson } 4540d0321e0SJeremy L Thompson 4550d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 456d075f50bSSebastian Grimberg // Create non-tensor H^1 4570d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 4582b730f8bSJeremy L Thompson int CeedBasisCreateH1_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *grad, 45951475c7cSJeremy L Thompson const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 4600d0321e0SJeremy L Thompson Ceed ceed; 461d075f50bSSebastian Grimberg CeedInt num_comp, q_comp_interp, q_comp_grad; 462ca735530SJeremy L Thompson const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 4630d0321e0SJeremy L Thompson CeedBasisNonTensor_Cuda *data; 464ca735530SJeremy L Thompson 465ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 4662b730f8bSJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 4670d0321e0SJeremy L Thompson 4680d0321e0SJeremy L Thompson // Copy basis data to GPU 469d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 470d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad)); 471097cc795SJames Wright if (q_weight) { 4722b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes)); 4732b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice)); 474097cc795SJames Wright } 475d075f50bSSebastian Grimberg if (interp) { 476d075f50bSSebastian Grimberg const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 477d075f50bSSebastian Grimberg 4782b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes)); 4792b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice)); 480d075f50bSSebastian Grimberg } 481d075f50bSSebastian Grimberg if (grad) { 482d075f50bSSebastian Grimberg const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad; 483d075f50bSSebastian Grimberg 4842b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad, grad_bytes)); 4852b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad, grad, grad_bytes, cudaMemcpyHostToDevice)); 486d075f50bSSebastian Grimberg } 4870d0321e0SJeremy L Thompson 4880d0321e0SJeremy L Thompson // Compile basis kernels 489*9c25dd66SJeremy L Thompson const char basis_kernel_source[] = "// Nontensor basis source\n#include <ceed/jit-source/cuda/cuda-ref-basis-nontensor.h>\n"; 490*9c25dd66SJeremy L Thompson 4912b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 492d075f50bSSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP", 493d075f50bSSebastian Grimberg q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_grad, "BASIS_NUM_COMP", num_comp)); 494d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 495d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 496d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv)); 497d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose)); 498d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 499d075f50bSSebastian Grimberg 500d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisSetData(basis, data)); 501d075f50bSSebastian Grimberg 502d075f50bSSebastian Grimberg // Register backend functions 503d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda)); 504db2becc9SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Cuda)); 505d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda)); 506d075f50bSSebastian Grimberg return CEED_ERROR_SUCCESS; 507d075f50bSSebastian Grimberg } 508d075f50bSSebastian Grimberg 509d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 510d075f50bSSebastian Grimberg // Create non-tensor H(div) 511d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 512d075f50bSSebastian Grimberg int CeedBasisCreateHdiv_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *div, 513d075f50bSSebastian Grimberg const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 514d075f50bSSebastian Grimberg Ceed ceed; 515d075f50bSSebastian Grimberg CeedInt num_comp, q_comp_interp, q_comp_div; 516d075f50bSSebastian Grimberg const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 517d075f50bSSebastian Grimberg CeedBasisNonTensor_Cuda *data; 518d075f50bSSebastian Grimberg 519d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 520d075f50bSSebastian Grimberg CeedCallBackend(CeedCalloc(1, &data)); 521d075f50bSSebastian Grimberg 522d075f50bSSebastian Grimberg // Copy basis data to GPU 523d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 524d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div)); 525097cc795SJames Wright if (q_weight) { 526d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes)); 527d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice)); 528097cc795SJames Wright } 529d075f50bSSebastian Grimberg if (interp) { 530d075f50bSSebastian Grimberg const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 531d075f50bSSebastian Grimberg 532d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes)); 533d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice)); 534d075f50bSSebastian Grimberg } 535d075f50bSSebastian Grimberg if (div) { 536d075f50bSSebastian Grimberg const CeedInt div_bytes = q_bytes * num_nodes * q_comp_div; 537d075f50bSSebastian Grimberg 538d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_div, div_bytes)); 539d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_div, div, div_bytes, cudaMemcpyHostToDevice)); 540d075f50bSSebastian Grimberg } 541d075f50bSSebastian Grimberg 542d075f50bSSebastian Grimberg // Compile basis kernels 543*9c25dd66SJeremy L Thompson const char basis_kernel_source[] = "// Nontensor basis source\n#include <ceed/jit-source/cuda/cuda-ref-basis-nontensor.h>\n"; 544*9c25dd66SJeremy L Thompson 545d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 546d075f50bSSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP", 547d075f50bSSebastian Grimberg q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_div, "BASIS_NUM_COMP", num_comp)); 548d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 549d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 550d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv)); 551d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose)); 552d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 553d075f50bSSebastian Grimberg 554d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisSetData(basis, data)); 555d075f50bSSebastian Grimberg 556d075f50bSSebastian Grimberg // Register backend functions 557d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda)); 558db2becc9SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Cuda)); 559d075f50bSSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda)); 560d075f50bSSebastian Grimberg return CEED_ERROR_SUCCESS; 561d075f50bSSebastian Grimberg } 562d075f50bSSebastian Grimberg 563d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 564d075f50bSSebastian Grimberg // Create non-tensor H(curl) 565d075f50bSSebastian Grimberg //------------------------------------------------------------------------------ 566d075f50bSSebastian Grimberg int CeedBasisCreateHcurl_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 567d075f50bSSebastian Grimberg const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 568d075f50bSSebastian Grimberg Ceed ceed; 569d075f50bSSebastian Grimberg CeedInt num_comp, q_comp_interp, q_comp_curl; 570d075f50bSSebastian Grimberg const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 571d075f50bSSebastian Grimberg CeedBasisNonTensor_Cuda *data; 572d075f50bSSebastian Grimberg 573d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 574d075f50bSSebastian Grimberg CeedCallBackend(CeedCalloc(1, &data)); 575d075f50bSSebastian Grimberg 576d075f50bSSebastian Grimberg // Copy basis data to GPU 577d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 578d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl)); 579097cc795SJames Wright if (q_weight) { 580d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes)); 581d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice)); 582097cc795SJames Wright } 583d075f50bSSebastian Grimberg if (interp) { 584d075f50bSSebastian Grimberg const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 585d075f50bSSebastian Grimberg 586d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes)); 587d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice)); 588d075f50bSSebastian Grimberg } 589d075f50bSSebastian Grimberg if (curl) { 590d075f50bSSebastian Grimberg const CeedInt curl_bytes = q_bytes * num_nodes * q_comp_curl; 591d075f50bSSebastian Grimberg 592d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&data->d_curl, curl_bytes)); 593d075f50bSSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(data->d_curl, curl, curl_bytes, cudaMemcpyHostToDevice)); 594d075f50bSSebastian Grimberg } 595d075f50bSSebastian Grimberg 596d075f50bSSebastian Grimberg // Compile basis kernels 597*9c25dd66SJeremy L Thompson const char basis_kernel_source[] = "// Nontensor basis source\n#include <ceed/jit-source/cuda/cuda-ref-basis-nontensor.h>\n"; 598*9c25dd66SJeremy L Thompson 599d075f50bSSebastian Grimberg CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 600d075f50bSSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP", 601d075f50bSSebastian Grimberg q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_curl, "BASIS_NUM_COMP", num_comp)); 602d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 603d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 604d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv)); 605d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose)); 606d075f50bSSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 6070d0321e0SJeremy L Thompson 6082b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 6090d0321e0SJeremy L Thompson 6100d0321e0SJeremy L Thompson // Register backend functions 6112b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda)); 612db2becc9SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Cuda)); 6132b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda)); 6140d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 6150d0321e0SJeremy L Thompson } 6162a86cc9dSSebastian Grimberg 6170d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 618