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. 3c532df63SYohann // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5c532df63SYohann // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7c532df63SYohann 849aac155SJeremy L Thompson #include <ceed.h> 9ec3da8bcSJed Brown #include <ceed/backend.h> 10437930d1SJeremy L Thompson #include <ceed/jit-tools.h> 113d576824SJeremy L Thompson #include <cuda.h> 123d576824SJeremy L Thompson #include <cuda_runtime.h> 1349aac155SJeremy L Thompson #include <stdbool.h> 143d576824SJeremy L Thompson #include <stddef.h> 15111870feSJeremy L Thompson #include <string.h> 16c532df63SYohann 1749aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h" 182b730f8bSJeremy L Thompson #include "../cuda/ceed-cuda-compile.h" 192b730f8bSJeremy L Thompson #include "ceed-cuda-shared.h" 20c532df63SYohann 21ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 229ff05d55SJeremy L Thompson // Apply tensor basis 23ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 24db2becc9SJeremy L Thompson static int CeedBasisApplyTensorCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode, 25db2becc9SJeremy L Thompson CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 26c532df63SYohann Ceed ceed; 276dbfb411Snbeams Ceed_Cuda *ceed_Cuda; 28437930d1SJeremy L Thompson CeedInt dim, num_comp; 29ca735530SJeremy L Thompson const CeedScalar *d_u; 30ca735530SJeremy L Thompson CeedScalar *d_v; 31ca735530SJeremy L Thompson CeedBasis_Cuda_shared *data; 32ca735530SJeremy L Thompson 33ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 34ca735530SJeremy L Thompson CeedCallBackend(CeedGetData(ceed, &ceed_Cuda)); 35ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 362b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 372b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 38c532df63SYohann 399ea2cfd9SJeremy L Thompson // Get read/write access to u, v 406574a04fSJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 416574a04fSJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 42db2becc9SJeremy L Thompson if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 43db2becc9SJeremy L Thompson else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 44c532df63SYohann 45ab213215SJeremy L Thompson // Apply basis operation 46437930d1SJeremy L Thompson switch (eval_mode) { 47ab213215SJeremy L Thompson case CEED_EVAL_INTERP: { 48437930d1SJeremy L Thompson CeedInt P_1d, Q_1d; 49ca735530SJeremy L Thompson 504cbc44e0SJeremy L Thompson CeedCheck(data->d_interp_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; interp_1d not set", CeedEvalModes[eval_mode]); 512b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 522b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 53437930d1SJeremy L Thompson CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 54ca735530SJeremy L Thompson 55aa4002adSJeremy L Thompson void *interp_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_u, &d_v}; 56ca735530SJeremy L Thompson 574d537eeaSYohann if (dim == 1) { 58dac82721SJeremy L Thompson // avoid >512 total threads 59dac82721SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1)); 60a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 61437930d1SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 62b2165e7aSSebastian Grimberg 639e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 64db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, 1, 65db2becc9SJeremy L Thompson elems_per_block, shared_mem, interp_args)); 669e201c85SYohann } else { 67eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args)); 689e201c85SYohann } 69074be161SYohann Dudouit } else if (dim == 2) { 70437930d1SJeremy L Thompson const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8}; 71437930d1SJeremy L Thompson // elems_per_block must be at least 1 722b730f8bSJeremy L Thompson CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1); 73a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 742b730f8bSJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 75b2165e7aSSebastian Grimberg 769e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 77db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, thread_1d, 78db2becc9SJeremy L Thompson elems_per_block, shared_mem, interp_args)); 799e201c85SYohann } else { 80eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 819e201c85SYohann } 82074be161SYohann Dudouit } else if (dim == 3) { 83437930d1SJeremy L Thompson CeedInt elems_per_block = 1; 84a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 852b730f8bSJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 86b2165e7aSSebastian Grimberg 879e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 88db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, thread_1d, 89db2becc9SJeremy L Thompson elems_per_block, shared_mem, interp_args)); 909e201c85SYohann } else { 91eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 92074be161SYohann Dudouit } 939e201c85SYohann } 94ab213215SJeremy L Thompson } break; 95ab213215SJeremy L Thompson case CEED_EVAL_GRAD: { 96437930d1SJeremy L Thompson CeedInt P_1d, Q_1d; 97ca735530SJeremy L Thompson 984cbc44e0SJeremy L Thompson CeedCheck(data->d_grad_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; grad_1d not set", CeedEvalModes[eval_mode]); 992b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 1002b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 101437930d1SJeremy L Thompson CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 102aa4002adSJeremy L Thompson CeedScalar *d_grad_1d = data->d_grad_1d; 103ca735530SJeremy L Thompson 1049e201c85SYohann if (data->d_collo_grad_1d) { 105aa4002adSJeremy L Thompson d_grad_1d = data->d_collo_grad_1d; 1069e201c85SYohann } 107aa4002adSJeremy L Thompson void *grad_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_grad_1d, &d_u, &d_v}; 108aa4002adSJeremy L Thompson 1094d537eeaSYohann if (dim == 1) { 110dac82721SJeremy L Thompson // avoid >512 total threads 111dac82721SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1)); 112a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 113437930d1SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 114b2165e7aSSebastian Grimberg 1159e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 116db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, 1, 117db2becc9SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 1189e201c85SYohann } else { 119eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args)); 1209e201c85SYohann } 121074be161SYohann Dudouit } else if (dim == 2) { 122437930d1SJeremy L Thompson const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8}; 123437930d1SJeremy L Thompson // elems_per_block must be at least 1 1242b730f8bSJeremy L Thompson CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1); 125a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 1262b730f8bSJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 127b2165e7aSSebastian Grimberg 1289e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 129db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, thread_1d, 130db2becc9SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 1319e201c85SYohann } else { 132eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 1339e201c85SYohann } 134074be161SYohann Dudouit } else if (dim == 3) { 135437930d1SJeremy L Thompson CeedInt elems_per_block = 1; 136a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 1372b730f8bSJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 138b2165e7aSSebastian Grimberg 1399e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 140db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, thread_1d, 141db2becc9SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 1429e201c85SYohann } else { 143eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 1449e201c85SYohann } 145074be161SYohann Dudouit } 146ab213215SJeremy L Thompson } break; 147ab213215SJeremy L Thompson case CEED_EVAL_WEIGHT: { 148437930d1SJeremy L Thompson CeedInt Q_1d; 149397164e9SSebastian Grimberg CeedInt block_size = 32; 150ca735530SJeremy L Thompson 151097cc795SJames Wright CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights_1d not set", CeedEvalModes[eval_mode]); 1522b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 153437930d1SJeremy L Thompson void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v}; 154074be161SYohann Dudouit if (dim == 1) { 155397164e9SSebastian Grimberg const CeedInt elems_per_block = block_size / Q_1d; 156a8d440fbSJeremy L Thompson const CeedInt grid_size = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 157b2165e7aSSebastian Grimberg 158b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, elems_per_block, 1, weight_args)); 159074be161SYohann Dudouit } else if (dim == 2) { 160397164e9SSebastian Grimberg const CeedInt opt_elems = block_size / (Q_1d * Q_1d); 161437930d1SJeremy L Thompson const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1; 162a8d440fbSJeremy L Thompson const CeedInt grid_size = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 163b2165e7aSSebastian Grimberg 164b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args)); 165074be161SYohann Dudouit } else if (dim == 3) { 166397164e9SSebastian Grimberg const CeedInt opt_elems = block_size / (Q_1d * Q_1d); 1679e201c85SYohann const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1; 168a8d440fbSJeremy L Thompson const CeedInt grid_size = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 169b2165e7aSSebastian Grimberg 170b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args)); 171074be161SYohann Dudouit } 172ab213215SJeremy L Thompson } break; 1739ea2cfd9SJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 1749ea2cfd9SJeremy L Thompson break; 175ab213215SJeremy L Thompson // LCOV_EXCL_START 176ab213215SJeremy L Thompson case CEED_EVAL_DIV: 177ab213215SJeremy L Thompson case CEED_EVAL_CURL: 178bcbe1c99SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 179ab213215SJeremy L Thompson // LCOV_EXCL_STOP 180c532df63SYohann } 181c532df63SYohann 1829ea2cfd9SJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 1832b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 1849ea2cfd9SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 1859ea2cfd9SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 1869bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 187e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 188c532df63SYohann } 189c532df63SYohann 190db2becc9SJeremy L Thompson static int CeedBasisApplyTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, 191db2becc9SJeremy L Thompson CeedVector v) { 192db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyTensorCore_Cuda_shared(basis, false, num_elem, t_mode, eval_mode, u, v)); 193db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 194db2becc9SJeremy L Thompson } 195db2becc9SJeremy L Thompson 196db2becc9SJeremy L Thompson static int CeedBasisApplyAddTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, 197db2becc9SJeremy L Thompson CeedVector u, CeedVector v) { 198db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyTensorCore_Cuda_shared(basis, true, num_elem, t_mode, eval_mode, u, v)); 199db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 200db2becc9SJeremy L Thompson } 201db2becc9SJeremy L Thompson 202ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 2031dda9c1aSJeremy L Thompson // Basis apply - tensor AtPoints 2041dda9c1aSJeremy L Thompson //------------------------------------------------------------------------------ 205db2becc9SJeremy L Thompson static int CeedBasisApplyAtPointsCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, const CeedInt *num_points, 206db2becc9SJeremy L Thompson CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) { 2071dda9c1aSJeremy L Thompson Ceed ceed; 2089e1d4b82SJeremy L Thompson Ceed_Cuda *ceed_Cuda; 2099e1d4b82SJeremy L Thompson CeedInt Q_1d, dim, num_comp, max_num_points = num_points[0]; 2101dda9c1aSJeremy L Thompson const CeedInt is_transpose = t_mode == CEED_TRANSPOSE; 2111dda9c1aSJeremy L Thompson const CeedScalar *d_x, *d_u; 2121dda9c1aSJeremy L Thompson CeedScalar *d_v; 2131dda9c1aSJeremy L Thompson CeedBasis_Cuda_shared *data; 2141dda9c1aSJeremy L Thompson 2151dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 2161dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 2171dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 2189e1d4b82SJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 2191dda9c1aSJeremy L Thompson 2201dda9c1aSJeremy L Thompson // Weight handled separately 2211dda9c1aSJeremy L Thompson if (eval_mode == CEED_EVAL_WEIGHT) { 2225a5594ffSJeremy L Thompson CeedCallBackend(CeedVectorSetValue(v, 1.0)); 2231dda9c1aSJeremy L Thompson return CEED_ERROR_SUCCESS; 2241dda9c1aSJeremy L Thompson } 2251dda9c1aSJeremy L Thompson 2269bc66399SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 2279e1d4b82SJeremy L Thompson CeedCallBackend(CeedGetData(ceed, &ceed_Cuda)); 2289bc66399SJeremy L Thompson 229111870feSJeremy L Thompson // Check padded to uniform number of points per elem 230111870feSJeremy L Thompson for (CeedInt i = 1; i < num_elem; i++) max_num_points = CeedIntMax(max_num_points, num_points[i]); 231111870feSJeremy L Thompson { 2329e1d4b82SJeremy L Thompson CeedInt q_comp; 233111870feSJeremy L Thompson CeedSize len, len_required; 234111870feSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 235111870feSJeremy L Thompson CeedCallBackend(CeedVectorGetLength(is_transpose ? u : v, &len)); 236111870feSJeremy L Thompson len_required = (CeedSize)num_comp * (CeedSize)q_comp * (CeedSize)num_elem * (CeedSize)max_num_points; 237111870feSJeremy L Thompson CeedCheck(len >= len_required, ceed, CEED_ERROR_BACKEND, 238111870feSJeremy L Thompson "Vector at points must be padded to the same number of points in each element for BasisApplyAtPoints on GPU backends." 239111870feSJeremy L Thompson " Found %" CeedSize_FMT ", Required %" CeedSize_FMT, 240111870feSJeremy L Thompson len, len_required); 241111870feSJeremy L Thompson } 242111870feSJeremy L Thompson 243111870feSJeremy L Thompson // Move num_points array to device 244111870feSJeremy L Thompson if (is_transpose) { 245111870feSJeremy L Thompson const CeedInt num_bytes = num_elem * sizeof(CeedInt); 246111870feSJeremy L Thompson 247111870feSJeremy L Thompson if (num_elem != data->num_elem_at_points) { 248111870feSJeremy L Thompson data->num_elem_at_points = num_elem; 249111870feSJeremy L Thompson 250111870feSJeremy L Thompson if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem)); 251111870feSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_points_per_elem, num_bytes)); 252111870feSJeremy L Thompson CeedCallBackend(CeedFree(&data->h_points_per_elem)); 253111870feSJeremy L Thompson CeedCallBackend(CeedCalloc(num_elem, &data->h_points_per_elem)); 254111870feSJeremy L Thompson } 2559e511c80SJeremy L Thompson if (memcmp(data->h_points_per_elem, num_points, num_bytes)) { 256111870feSJeremy L Thompson memcpy(data->h_points_per_elem, num_points, num_bytes); 257111870feSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_points_per_elem, num_points, num_bytes, cudaMemcpyHostToDevice)); 258111870feSJeremy L Thompson } 259111870feSJeremy L Thompson } 260111870feSJeremy L Thompson 2611dda9c1aSJeremy L Thompson // Build kernels if needed 2621dda9c1aSJeremy L Thompson if (data->num_points != max_num_points) { 2631dda9c1aSJeremy L Thompson CeedInt P_1d; 2641dda9c1aSJeremy L Thompson 2651dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 2661dda9c1aSJeremy L Thompson data->num_points = max_num_points; 2671dda9c1aSJeremy L Thompson 2681dda9c1aSJeremy L Thompson // -- Create interp matrix to Chebyshev coefficients 2691dda9c1aSJeremy L Thompson if (!data->d_chebyshev_interp_1d) { 2701dda9c1aSJeremy L Thompson CeedSize interp_bytes; 2711dda9c1aSJeremy L Thompson CeedScalar *chebyshev_interp_1d; 2721dda9c1aSJeremy L Thompson 2731dda9c1aSJeremy L Thompson interp_bytes = P_1d * Q_1d * sizeof(CeedScalar); 2741dda9c1aSJeremy L Thompson CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d)); 2755a5594ffSJeremy L Thompson CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d)); 2761dda9c1aSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_chebyshev_interp_1d, interp_bytes)); 2771dda9c1aSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 2781dda9c1aSJeremy L Thompson CeedCallBackend(CeedFree(&chebyshev_interp_1d)); 2791dda9c1aSJeremy L Thompson } 2801dda9c1aSJeremy L Thompson 2811dda9c1aSJeremy L Thompson // -- Compile kernels 2829e1d4b82SJeremy L Thompson const char basis_kernel_source[] = "// AtPoints basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-at-points.h>\n"; 2831dda9c1aSJeremy L Thompson CeedInt num_comp; 2841dda9c1aSJeremy L Thompson 2851dda9c1aSJeremy L Thompson if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints)); 2861dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 2879e1d4b82SJeremy L Thompson CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->moduleAtPoints, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D", 2889e1d4b82SJeremy L Thompson CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), 2899e1d4b82SJeremy L Thompson "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_NUM_PTS", max_num_points)); 2901dda9c1aSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpAtPoints", &data->InterpAtPoints)); 29181ae6159SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpTransposeAtPoints", &data->InterpTransposeAtPoints)); 2921dda9c1aSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradAtPoints", &data->GradAtPoints)); 29381ae6159SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradTransposeAtPoints", &data->GradTransposeAtPoints)); 2941dda9c1aSJeremy L Thompson } 2951dda9c1aSJeremy L Thompson 2961dda9c1aSJeremy L Thompson // Get read/write access to u, v 2971dda9c1aSJeremy L Thompson CeedCallBackend(CeedVectorGetArrayRead(x_ref, CEED_MEM_DEVICE, &d_x)); 2981dda9c1aSJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 2991dda9c1aSJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 300*9dafd6dfSZach Atkins if (apply_add) { 301*9dafd6dfSZach Atkins CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 302*9dafd6dfSZach Atkins } else { 3031dda9c1aSJeremy L Thompson // Clear v for transpose operation 304*9dafd6dfSZach Atkins if (is_transpose) CeedCallBackend(CeedVectorSetValue(v, 0.0)); 305*9dafd6dfSZach Atkins CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 3061dda9c1aSJeremy L Thompson } 3071dda9c1aSJeremy L Thompson 3081dda9c1aSJeremy L Thompson // Basis action 3091dda9c1aSJeremy L Thompson switch (eval_mode) { 3101dda9c1aSJeremy L Thompson case CEED_EVAL_INTERP: { 3119e1d4b82SJeremy L Thompson CeedInt P_1d, Q_1d; 3121dda9c1aSJeremy L Thompson 3139e1d4b82SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 3149e1d4b82SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 3159e1d4b82SJeremy L Thompson CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 3169e1d4b82SJeremy L Thompson 317aa4002adSJeremy L Thompson void *interp_args[] = {(void *)&num_elem, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v}; 3189e1d4b82SJeremy L Thompson 3199e1d4b82SJeremy L Thompson if (dim == 1) { 320dac82721SJeremy L Thompson // avoid >512 total threads 321dac82721SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1)); 322a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3239e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 3249e1d4b82SJeremy L Thompson 3259e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d, 1, 3269e1d4b82SJeremy L Thompson elems_per_block, shared_mem, interp_args)); 3279e1d4b82SJeremy L Thompson } else if (dim == 2) { 3289e1d4b82SJeremy L Thompson const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8}; 3299e1d4b82SJeremy L Thompson // elems_per_block must be at least 1 3309e1d4b82SJeremy L Thompson CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1); 331a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3329e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 3339e1d4b82SJeremy L Thompson 3349e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d, 3359e1d4b82SJeremy L Thompson thread_1d, elems_per_block, shared_mem, interp_args)); 3369e1d4b82SJeremy L Thompson } else if (dim == 3) { 3379e1d4b82SJeremy L Thompson CeedInt elems_per_block = 1; 338a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3399e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 3409e1d4b82SJeremy L Thompson 3419e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d, 3429e1d4b82SJeremy L Thompson thread_1d, elems_per_block, shared_mem, interp_args)); 3439e1d4b82SJeremy L Thompson } 3441dda9c1aSJeremy L Thompson } break; 3451dda9c1aSJeremy L Thompson case CEED_EVAL_GRAD: { 3469e1d4b82SJeremy L Thompson CeedInt P_1d, Q_1d; 3471dda9c1aSJeremy L Thompson 3489e1d4b82SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 3499e1d4b82SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 3509e1d4b82SJeremy L Thompson CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 3519e1d4b82SJeremy L Thompson 3529e1d4b82SJeremy L Thompson void *grad_args[] = {(void *)&num_elem, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v}; 3539e1d4b82SJeremy L Thompson 3549e1d4b82SJeremy L Thompson if (dim == 1) { 355dac82721SJeremy L Thompson // avoid >512 total threads 356dac82721SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1)); 357a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3589e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 3599e1d4b82SJeremy L Thompson 3609e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, 1, 3619e1d4b82SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 3629e1d4b82SJeremy L Thompson } else if (dim == 2) { 3639e1d4b82SJeremy L Thompson const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8}; 3649e1d4b82SJeremy L Thompson // elems_per_block must be at least 1 3659e1d4b82SJeremy L Thompson CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1); 366a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3679e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 3689e1d4b82SJeremy L Thompson 3699e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, thread_1d, 3709e1d4b82SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 3719e1d4b82SJeremy L Thompson } else if (dim == 3) { 3729e1d4b82SJeremy L Thompson CeedInt elems_per_block = 1; 373a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3749e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 3759e1d4b82SJeremy L Thompson 3769e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, thread_1d, 3779e1d4b82SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 3789e1d4b82SJeremy L Thompson } 3791dda9c1aSJeremy L Thompson } break; 3801dda9c1aSJeremy L Thompson case CEED_EVAL_WEIGHT: 3811dda9c1aSJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 3821dda9c1aSJeremy L Thompson break; 3831dda9c1aSJeremy L Thompson // LCOV_EXCL_START 3841dda9c1aSJeremy L Thompson case CEED_EVAL_DIV: 3851dda9c1aSJeremy L Thompson case CEED_EVAL_CURL: 3861dda9c1aSJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 3871dda9c1aSJeremy L Thompson // LCOV_EXCL_STOP 3881dda9c1aSJeremy L Thompson } 3891dda9c1aSJeremy L Thompson 3901dda9c1aSJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 3911dda9c1aSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArrayRead(x_ref, &d_x)); 3921dda9c1aSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 3931dda9c1aSJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 3941dda9c1aSJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 3959bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 3961dda9c1aSJeremy L Thompson return CEED_ERROR_SUCCESS; 3971dda9c1aSJeremy L Thompson } 3981dda9c1aSJeremy L Thompson 399db2becc9SJeremy L Thompson static int CeedBasisApplyAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, 400db2becc9SJeremy L Thompson CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) { 401db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda_shared(basis, false, num_elem, num_points, t_mode, eval_mode, x_ref, u, v)); 402db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 403db2becc9SJeremy L Thompson } 404db2becc9SJeremy L Thompson 405db2becc9SJeremy L Thompson static int CeedBasisApplyAddAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, 406db2becc9SJeremy L Thompson CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) { 407db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda_shared(basis, true, num_elem, num_points, t_mode, eval_mode, x_ref, u, v)); 408db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 409db2becc9SJeremy L Thompson } 410db2becc9SJeremy L Thompson 4111dda9c1aSJeremy L Thompson //------------------------------------------------------------------------------ 4129ff05d55SJeremy L Thompson // Apply non-tensor basis 4139ff05d55SJeremy L Thompson //------------------------------------------------------------------------------ 4149ff05d55SJeremy L Thompson static int CeedBasisApplyNonTensorCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode, 4159ff05d55SJeremy L Thompson CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 4169ff05d55SJeremy L Thompson Ceed ceed; 4179ff05d55SJeremy L Thompson Ceed_Cuda *ceed_Cuda; 4189ff05d55SJeremy L Thompson CeedInt dim; 4199ff05d55SJeremy L Thompson const CeedScalar *d_u; 4209ff05d55SJeremy L Thompson CeedScalar *d_v; 4219ff05d55SJeremy L Thompson CeedBasis_Cuda_shared *data; 4229ff05d55SJeremy L Thompson 4239ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 4249ff05d55SJeremy L Thompson CeedCallBackend(CeedGetData(ceed, &ceed_Cuda)); 4259ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 4269ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 4279ff05d55SJeremy L Thompson 4289ff05d55SJeremy L Thompson // Get read/write access to u, v 4299ff05d55SJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 4309ff05d55SJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 4319ff05d55SJeremy L Thompson if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 4329ff05d55SJeremy L Thompson else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 4339ff05d55SJeremy L Thompson 4349ff05d55SJeremy L Thompson // Apply basis operation 4359ff05d55SJeremy L Thompson switch (eval_mode) { 4369ff05d55SJeremy L Thompson case CEED_EVAL_INTERP: { 4379ff05d55SJeremy L Thompson CeedInt P, Q; 4389ff05d55SJeremy L Thompson 4394cbc44e0SJeremy L Thompson CeedCheck(data->d_interp_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; interp not set", CeedEvalModes[eval_mode]); 4409ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes(basis, &P)); 4419ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q)); 4429ff05d55SJeremy L Thompson CeedInt thread = CeedIntMax(Q, P); 4439ff05d55SJeremy L Thompson 444aa4002adSJeremy L Thompson void *interp_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_u, &d_v}; 4459ff05d55SJeremy L Thompson 4469ff05d55SJeremy L Thompson { 4479ff05d55SJeremy L Thompson // avoid >512 total threads 4489ff05d55SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1)); 4499ff05d55SJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 4509ff05d55SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread * sizeof(CeedScalar); 4519ff05d55SJeremy L Thompson 4529ff05d55SJeremy L Thompson if (t_mode == CEED_TRANSPOSE) { 4539ff05d55SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread, 1, 4549ff05d55SJeremy L Thompson elems_per_block, shared_mem, interp_args)); 4559ff05d55SJeremy L Thompson } else { 4569ff05d55SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread, 1, elems_per_block, shared_mem, interp_args)); 4579ff05d55SJeremy L Thompson } 4589ff05d55SJeremy L Thompson } 4599ff05d55SJeremy L Thompson } break; 4609ff05d55SJeremy L Thompson case CEED_EVAL_GRAD: { 4619ff05d55SJeremy L Thompson CeedInt P, Q; 4629ff05d55SJeremy L Thompson 4634cbc44e0SJeremy L Thompson CeedCheck(data->d_grad_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; grad not set", CeedEvalModes[eval_mode]); 4649ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes(basis, &P)); 4659ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q)); 4669ff05d55SJeremy L Thompson CeedInt thread = CeedIntMax(Q, P); 4679ff05d55SJeremy L Thompson 468aa4002adSJeremy L Thompson void *grad_args[] = {(void *)&num_elem, &data->d_grad_1d, &d_u, &d_v}; 4699ff05d55SJeremy L Thompson 4709ff05d55SJeremy L Thompson { 4719ff05d55SJeremy L Thompson // avoid >512 total threads 4729ff05d55SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1)); 4739ff05d55SJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 4749ff05d55SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread * sizeof(CeedScalar); 4759ff05d55SJeremy L Thompson 4769ff05d55SJeremy L Thompson if (t_mode == CEED_TRANSPOSE) { 4779ff05d55SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread, 1, 4789ff05d55SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 4799ff05d55SJeremy L Thompson } else { 4809ff05d55SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread, 1, elems_per_block, shared_mem, grad_args)); 4819ff05d55SJeremy L Thompson } 4829ff05d55SJeremy L Thompson } 4839ff05d55SJeremy L Thompson } break; 4849ff05d55SJeremy L Thompson case CEED_EVAL_WEIGHT: { 48597011eabSZach Atkins CeedInt P, Q; 4869ff05d55SJeremy L Thompson 4874cbc44e0SJeremy L Thompson CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights not set", CeedEvalModes[eval_mode]); 48897011eabSZach Atkins CeedCallBackend(CeedBasisGetNumNodes(basis, &P)); 4899ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q)); 49097011eabSZach Atkins CeedInt thread = CeedIntMax(Q, P); 49197011eabSZach Atkins 4929ff05d55SJeremy L Thompson void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v}; 4939ff05d55SJeremy L Thompson 4949ff05d55SJeremy L Thompson { 4959ff05d55SJeremy L Thompson // avoid >512 total threads 49697011eabSZach Atkins CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1)); 4979ff05d55SJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 4989ff05d55SJeremy L Thompson 49997011eabSZach Atkins CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, thread, elems_per_block, 1, weight_args)); 5009ff05d55SJeremy L Thompson } 5019ff05d55SJeremy L Thompson } break; 5029ff05d55SJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 5039ff05d55SJeremy L Thompson break; 5049ff05d55SJeremy L Thompson // LCOV_EXCL_START 5059ff05d55SJeremy L Thompson case CEED_EVAL_DIV: 5069ff05d55SJeremy L Thompson case CEED_EVAL_CURL: 5079ff05d55SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 5089ff05d55SJeremy L Thompson // LCOV_EXCL_STOP 5099ff05d55SJeremy L Thompson } 5109ff05d55SJeremy L Thompson 5119ff05d55SJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 5129ff05d55SJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 5139ff05d55SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 5149ff05d55SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 5159ff05d55SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 5169ff05d55SJeremy L Thompson return CEED_ERROR_SUCCESS; 5179ff05d55SJeremy L Thompson } 5189ff05d55SJeremy L Thompson 5199ff05d55SJeremy L Thompson static int CeedBasisApplyNonTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, 5209ff05d55SJeremy L Thompson CeedVector u, CeedVector v) { 5219ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda_shared(basis, false, num_elem, t_mode, eval_mode, u, v)); 5229ff05d55SJeremy L Thompson return CEED_ERROR_SUCCESS; 5239ff05d55SJeremy L Thompson } 5249ff05d55SJeremy L Thompson 5259ff05d55SJeremy L Thompson static int CeedBasisApplyAddNonTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, 5269ff05d55SJeremy L Thompson CeedVector u, CeedVector v) { 5279ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda_shared(basis, true, num_elem, t_mode, eval_mode, u, v)); 5289ff05d55SJeremy L Thompson return CEED_ERROR_SUCCESS; 5299ff05d55SJeremy L Thompson } 5309ff05d55SJeremy L Thompson 5319ff05d55SJeremy L Thompson //------------------------------------------------------------------------------ 532ab213215SJeremy L Thompson // Destroy basis 533ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 534c532df63SYohann static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) { 535c532df63SYohann Ceed ceed; 536c532df63SYohann CeedBasis_Cuda_shared *data; 537ca735530SJeremy L Thompson 538ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 5392b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 5402b730f8bSJeremy L Thompson CeedCallCuda(ceed, cuModuleUnload(data->module)); 5411dda9c1aSJeremy L Thompson if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints)); 542097cc795SJames Wright if (data->d_q_weight_1d) CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d)); 543111870feSJeremy L Thompson CeedCallBackend(CeedFree(&data->h_points_per_elem)); 544111870feSJeremy L Thompson if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem)); 5452b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_interp_1d)); 5462b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_grad_1d)); 5472b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_collo_grad_1d)); 5481dda9c1aSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_chebyshev_interp_1d)); 5492b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&data)); 5509bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 551e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 552c532df63SYohann } 553c532df63SYohann 554ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 555ab213215SJeremy L Thompson // Create tensor basis 556ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 5572b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d, 5582b730f8bSJeremy L Thompson const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) { 559c532df63SYohann Ceed ceed; 560ca735530SJeremy L Thompson CeedInt num_comp; 561ca735530SJeremy L Thompson const CeedInt q_bytes = Q_1d * sizeof(CeedScalar); 562ca735530SJeremy L Thompson const CeedInt interp_bytes = q_bytes * P_1d; 563c532df63SYohann CeedBasis_Cuda_shared *data; 564ca735530SJeremy L Thompson 565ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 5662b730f8bSJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 567c532df63SYohann 568ab213215SJeremy L Thompson // Copy basis data to GPU 569097cc795SJames Wright if (q_weight_1d) { 5702b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes)); 5712b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice)); 572097cc795SJames Wright } 5732b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes)); 5742b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 5752b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes)); 5762b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice)); 577c532df63SYohann 578ab213215SJeremy L Thompson // Compute collocated gradient and copy to GPU 579437930d1SJeremy L Thompson data->d_collo_grad_1d = NULL; 5809e201c85SYohann bool has_collocated_grad = dim == 3 && Q_1d >= P_1d; 581ca735530SJeremy L Thompson 5829e201c85SYohann if (has_collocated_grad) { 583437930d1SJeremy L Thompson CeedScalar *collo_grad_1d; 584ca735530SJeremy L Thompson 5852b730f8bSJeremy L Thompson CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d)); 5862b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d)); 5872b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d)); 5882b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, cudaMemcpyHostToDevice)); 5892b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&collo_grad_1d)); 590ac421f39SYohann } 591ac421f39SYohann 592ab213215SJeremy L Thompson // Compile basis kernels 5939c25dd66SJeremy L Thompson const char basis_kernel_source[] = "// Tensor basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-tensor.h>\n"; 5949c25dd66SJeremy L Thompson 5952b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 596eb7e6cafSJeremy L Thompson CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D", 597eb7e6cafSJeremy L Thompson CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), 598eb7e6cafSJeremy L Thompson "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_HAS_COLLOCATED_GRAD", has_collocated_grad)); 599eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 600eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 601db2becc9SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd)); 602eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad)); 603eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose)); 604db2becc9SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd)); 605eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 606c532df63SYohann 6072b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 608ab213215SJeremy L Thompson 609ab213215SJeremy L Thompson // Register backend functions 6102b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared)); 611db2becc9SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddTensor_Cuda_shared)); 6121dda9c1aSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAtPoints", CeedBasisApplyAtPoints_Cuda_shared)); 613d474dafeSZach Atkins CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAddAtPoints", CeedBasisApplyAddAtPoints_Cuda_shared)); 6142b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared)); 6159bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 616e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 617c532df63SYohann } 6182a86cc9dSSebastian Grimberg 619ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 6209ff05d55SJeremy L Thompson // Create non-tensor basis 6219ff05d55SJeremy L Thompson //------------------------------------------------------------------------------ 6229ff05d55SJeremy L Thompson int CeedBasisCreateH1_Cuda_shared(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 6239ff05d55SJeremy L Thompson const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 6249ff05d55SJeremy L Thompson Ceed ceed; 6259ff05d55SJeremy L Thompson CeedInt num_comp, q_comp_interp, q_comp_grad; 6269ff05d55SJeremy L Thompson const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 6279ff05d55SJeremy L Thompson CeedBasis_Cuda_shared *data; 6289ff05d55SJeremy L Thompson 6299ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 6309ff05d55SJeremy L Thompson 631fda26546SJeremy L Thompson // Check shared memory size 632fda26546SJeremy L Thompson { 633fda26546SJeremy L Thompson Ceed_Cuda *cuda_data; 634fda26546SJeremy L Thompson 635fda26546SJeremy L Thompson CeedCallBackend(CeedGetData(ceed, &cuda_data)); 636fda26546SJeremy L Thompson if (((size_t)num_nodes * (size_t)num_qpts * (size_t)dim + (size_t)CeedIntMax(num_nodes, num_qpts)) * sizeof(CeedScalar) > 637fda26546SJeremy L Thompson cuda_data->device_prop.sharedMemPerBlock) { 638fda26546SJeremy L Thompson CeedCallBackend(CeedBasisCreateH1Fallback(ceed, topo, dim, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis)); 639fda26546SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 640fda26546SJeremy L Thompson return CEED_ERROR_SUCCESS; 641fda26546SJeremy L Thompson } 642fda26546SJeremy L Thompson } 643fda26546SJeremy L Thompson 644fda26546SJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 645aa4002adSJeremy L Thompson 6469ff05d55SJeremy L Thompson // Copy basis data to GPU 6479ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 6489ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad)); 6499ff05d55SJeremy L Thompson if (q_weight) { 6509ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes)); 6519ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight, q_bytes, cudaMemcpyHostToDevice)); 6529ff05d55SJeremy L Thompson } 6539ff05d55SJeremy L Thompson if (interp) { 6549ff05d55SJeremy L Thompson const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 6559ff05d55SJeremy L Thompson 6569ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes)); 6579ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp, interp_bytes, cudaMemcpyHostToDevice)); 6589ff05d55SJeremy L Thompson } 6599ff05d55SJeremy L Thompson if (grad) { 6609ff05d55SJeremy L Thompson const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad; 6619ff05d55SJeremy L Thompson 6629ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, grad_bytes)); 6639ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad, grad_bytes, cudaMemcpyHostToDevice)); 6649ff05d55SJeremy L Thompson } 6659ff05d55SJeremy L Thompson 6669ff05d55SJeremy L Thompson // Compile basis kernels 6679ff05d55SJeremy L Thompson const char basis_kernel_source[] = "// Non-tensor basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-nontensor.h>\n"; 6689ff05d55SJeremy L Thompson 6699ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 6709ff05d55SJeremy L Thompson CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "T_1D", 6719ff05d55SJeremy L Thompson CeedIntMax(num_qpts, num_nodes), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp)); 6729ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 6739ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 6749ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd)); 6759ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad)); 6769ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose)); 6779ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd)); 6789ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 6799ff05d55SJeremy L Thompson 6809ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 6819ff05d55SJeremy L Thompson 6829ff05d55SJeremy L Thompson // Register backend functions 6839ff05d55SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda_shared)); 6849ff05d55SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Cuda_shared)); 6859ff05d55SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared)); 6869ff05d55SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 6879ff05d55SJeremy L Thompson return CEED_ERROR_SUCCESS; 6889ff05d55SJeremy L Thompson } 6899ff05d55SJeremy L Thompson 6909ff05d55SJeremy L Thompson //------------------------------------------------------------------------------ 691