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 50*4cbc44e0SJeremy 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 98*4cbc44e0SJeremy 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"); 300db2becc9SJeremy L Thompson if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 301db2becc9SJeremy L Thompson else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 3021dda9c1aSJeremy L Thompson 3031dda9c1aSJeremy L Thompson // Clear v for transpose operation 304db2becc9SJeremy L Thompson if (is_transpose && !apply_add) { 30519a04db8SJeremy L Thompson CeedInt num_comp, q_comp, num_nodes; 3061dda9c1aSJeremy L Thompson CeedSize length; 3071dda9c1aSJeremy L Thompson 30819a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 30919a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 31019a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes)); 31119a04db8SJeremy L Thompson length = 31219a04db8SJeremy L Thompson (CeedSize)num_elem * (CeedSize)num_comp * (t_mode == CEED_TRANSPOSE ? (CeedSize)num_nodes : ((CeedSize)max_num_points * (CeedSize)q_comp)); 3131dda9c1aSJeremy L Thompson CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar))); 3141dda9c1aSJeremy L Thompson } 3151dda9c1aSJeremy L Thompson 3161dda9c1aSJeremy L Thompson // Basis action 3171dda9c1aSJeremy L Thompson switch (eval_mode) { 3181dda9c1aSJeremy L Thompson case CEED_EVAL_INTERP: { 3199e1d4b82SJeremy L Thompson CeedInt P_1d, Q_1d; 3201dda9c1aSJeremy L Thompson 3219e1d4b82SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 3229e1d4b82SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 3239e1d4b82SJeremy L Thompson CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 3249e1d4b82SJeremy L Thompson 325aa4002adSJeremy L Thompson void *interp_args[] = {(void *)&num_elem, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v}; 3269e1d4b82SJeremy L Thompson 3279e1d4b82SJeremy L Thompson if (dim == 1) { 328dac82721SJeremy L Thompson // avoid >512 total threads 329dac82721SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1)); 330a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3319e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 3329e1d4b82SJeremy L Thompson 3339e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d, 1, 3349e1d4b82SJeremy L Thompson elems_per_block, shared_mem, interp_args)); 3359e1d4b82SJeremy L Thompson } else if (dim == 2) { 3369e1d4b82SJeremy L Thompson const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8}; 3379e1d4b82SJeremy L Thompson // elems_per_block must be at least 1 3389e1d4b82SJeremy L Thompson CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1); 339a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3409e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 3419e1d4b82SJeremy L Thompson 3429e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d, 3439e1d4b82SJeremy L Thompson thread_1d, elems_per_block, shared_mem, interp_args)); 3449e1d4b82SJeremy L Thompson } else if (dim == 3) { 3459e1d4b82SJeremy L Thompson CeedInt elems_per_block = 1; 346a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3479e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 3489e1d4b82SJeremy L Thompson 3499e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d, 3509e1d4b82SJeremy L Thompson thread_1d, elems_per_block, shared_mem, interp_args)); 3519e1d4b82SJeremy L Thompson } 3521dda9c1aSJeremy L Thompson } break; 3531dda9c1aSJeremy L Thompson case CEED_EVAL_GRAD: { 3549e1d4b82SJeremy L Thompson CeedInt P_1d, Q_1d; 3551dda9c1aSJeremy L Thompson 3569e1d4b82SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 3579e1d4b82SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 3589e1d4b82SJeremy L Thompson CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 3599e1d4b82SJeremy L Thompson 3609e1d4b82SJeremy L Thompson void *grad_args[] = {(void *)&num_elem, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v}; 3619e1d4b82SJeremy L Thompson 3629e1d4b82SJeremy L Thompson if (dim == 1) { 363dac82721SJeremy L Thompson // avoid >512 total threads 364dac82721SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1)); 365a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3669e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 3679e1d4b82SJeremy L Thompson 3689e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, 1, 3699e1d4b82SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 3709e1d4b82SJeremy L Thompson } else if (dim == 2) { 3719e1d4b82SJeremy L Thompson const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8}; 3729e1d4b82SJeremy L Thompson // elems_per_block must be at least 1 3739e1d4b82SJeremy L Thompson CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1); 374a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3759e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 3769e1d4b82SJeremy L Thompson 3779e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, thread_1d, 3789e1d4b82SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 3799e1d4b82SJeremy L Thompson } else if (dim == 3) { 3809e1d4b82SJeremy L Thompson CeedInt elems_per_block = 1; 381a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3829e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 3839e1d4b82SJeremy L Thompson 3849e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, thread_1d, 3859e1d4b82SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 3869e1d4b82SJeremy L Thompson } 3871dda9c1aSJeremy L Thompson } break; 3881dda9c1aSJeremy L Thompson case CEED_EVAL_WEIGHT: 3891dda9c1aSJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 3901dda9c1aSJeremy L Thompson break; 3911dda9c1aSJeremy L Thompson // LCOV_EXCL_START 3921dda9c1aSJeremy L Thompson case CEED_EVAL_DIV: 3931dda9c1aSJeremy L Thompson case CEED_EVAL_CURL: 3941dda9c1aSJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 3951dda9c1aSJeremy L Thompson // LCOV_EXCL_STOP 3961dda9c1aSJeremy L Thompson } 3971dda9c1aSJeremy L Thompson 3981dda9c1aSJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 3991dda9c1aSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArrayRead(x_ref, &d_x)); 4001dda9c1aSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 4011dda9c1aSJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 4021dda9c1aSJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 4039bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 4041dda9c1aSJeremy L Thompson return CEED_ERROR_SUCCESS; 4051dda9c1aSJeremy L Thompson } 4061dda9c1aSJeremy L Thompson 407db2becc9SJeremy L Thompson static int CeedBasisApplyAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, 408db2becc9SJeremy L Thompson CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) { 409db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda_shared(basis, false, num_elem, num_points, t_mode, eval_mode, x_ref, u, v)); 410db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 411db2becc9SJeremy L Thompson } 412db2becc9SJeremy L Thompson 413db2becc9SJeremy L Thompson static int CeedBasisApplyAddAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, 414db2becc9SJeremy L Thompson CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) { 415db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda_shared(basis, true, num_elem, num_points, t_mode, eval_mode, x_ref, u, v)); 416db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 417db2becc9SJeremy L Thompson } 418db2becc9SJeremy L Thompson 4191dda9c1aSJeremy L Thompson //------------------------------------------------------------------------------ 4209ff05d55SJeremy L Thompson // Apply non-tensor basis 4219ff05d55SJeremy L Thompson //------------------------------------------------------------------------------ 4229ff05d55SJeremy L Thompson static int CeedBasisApplyNonTensorCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode, 4239ff05d55SJeremy L Thompson CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 4249ff05d55SJeremy L Thompson Ceed ceed; 4259ff05d55SJeremy L Thompson Ceed_Cuda *ceed_Cuda; 4269ff05d55SJeremy L Thompson CeedInt dim; 4279ff05d55SJeremy L Thompson const CeedScalar *d_u; 4289ff05d55SJeremy L Thompson CeedScalar *d_v; 4299ff05d55SJeremy L Thompson CeedBasis_Cuda_shared *data; 4309ff05d55SJeremy L Thompson 4319ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 4329ff05d55SJeremy L Thompson CeedCallBackend(CeedGetData(ceed, &ceed_Cuda)); 4339ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 4349ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 4359ff05d55SJeremy L Thompson 4369ff05d55SJeremy L Thompson // Get read/write access to u, v 4379ff05d55SJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 4389ff05d55SJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 4399ff05d55SJeremy L Thompson if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 4409ff05d55SJeremy L Thompson else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 4419ff05d55SJeremy L Thompson 4429ff05d55SJeremy L Thompson // Apply basis operation 4439ff05d55SJeremy L Thompson switch (eval_mode) { 4449ff05d55SJeremy L Thompson case CEED_EVAL_INTERP: { 4459ff05d55SJeremy L Thompson CeedInt P, Q; 4469ff05d55SJeremy L Thompson 447*4cbc44e0SJeremy L Thompson CeedCheck(data->d_interp_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; interp not set", CeedEvalModes[eval_mode]); 4489ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes(basis, &P)); 4499ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q)); 4509ff05d55SJeremy L Thompson CeedInt thread = CeedIntMax(Q, P); 4519ff05d55SJeremy L Thompson 452aa4002adSJeremy L Thompson void *interp_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_u, &d_v}; 4539ff05d55SJeremy L Thompson 4549ff05d55SJeremy L Thompson { 4559ff05d55SJeremy L Thompson // avoid >512 total threads 4569ff05d55SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1)); 4579ff05d55SJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 4589ff05d55SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread * sizeof(CeedScalar); 4599ff05d55SJeremy L Thompson 4609ff05d55SJeremy L Thompson if (t_mode == CEED_TRANSPOSE) { 4619ff05d55SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread, 1, 4629ff05d55SJeremy L Thompson elems_per_block, shared_mem, interp_args)); 4639ff05d55SJeremy L Thompson } else { 4649ff05d55SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread, 1, elems_per_block, shared_mem, interp_args)); 4659ff05d55SJeremy L Thompson } 4669ff05d55SJeremy L Thompson } 4679ff05d55SJeremy L Thompson } break; 4689ff05d55SJeremy L Thompson case CEED_EVAL_GRAD: { 4699ff05d55SJeremy L Thompson CeedInt P, Q; 4709ff05d55SJeremy L Thompson 471*4cbc44e0SJeremy L Thompson CeedCheck(data->d_grad_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; grad not set", CeedEvalModes[eval_mode]); 4729ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes(basis, &P)); 4739ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q)); 4749ff05d55SJeremy L Thompson CeedInt thread = CeedIntMax(Q, P); 4759ff05d55SJeremy L Thompson 476aa4002adSJeremy L Thompson void *grad_args[] = {(void *)&num_elem, &data->d_grad_1d, &d_u, &d_v}; 4779ff05d55SJeremy L Thompson 4789ff05d55SJeremy L Thompson { 4799ff05d55SJeremy L Thompson // avoid >512 total threads 4809ff05d55SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1)); 4819ff05d55SJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 4829ff05d55SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread * sizeof(CeedScalar); 4839ff05d55SJeremy L Thompson 4849ff05d55SJeremy L Thompson if (t_mode == CEED_TRANSPOSE) { 4859ff05d55SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread, 1, 4869ff05d55SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 4879ff05d55SJeremy L Thompson } else { 4889ff05d55SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread, 1, elems_per_block, shared_mem, grad_args)); 4899ff05d55SJeremy L Thompson } 4909ff05d55SJeremy L Thompson } 4919ff05d55SJeremy L Thompson } break; 4929ff05d55SJeremy L Thompson case CEED_EVAL_WEIGHT: { 49397011eabSZach Atkins CeedInt P, Q; 4949ff05d55SJeremy L Thompson 495*4cbc44e0SJeremy L Thompson CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights not set", CeedEvalModes[eval_mode]); 49697011eabSZach Atkins CeedCallBackend(CeedBasisGetNumNodes(basis, &P)); 4979ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q)); 49897011eabSZach Atkins CeedInt thread = CeedIntMax(Q, P); 49997011eabSZach Atkins 5009ff05d55SJeremy L Thompson void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v}; 5019ff05d55SJeremy L Thompson 5029ff05d55SJeremy L Thompson { 5039ff05d55SJeremy L Thompson // avoid >512 total threads 50497011eabSZach Atkins CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1)); 5059ff05d55SJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 5069ff05d55SJeremy L Thompson 50797011eabSZach Atkins CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, thread, elems_per_block, 1, weight_args)); 5089ff05d55SJeremy L Thompson } 5099ff05d55SJeremy L Thompson } break; 5109ff05d55SJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 5119ff05d55SJeremy L Thompson break; 5129ff05d55SJeremy L Thompson // LCOV_EXCL_START 5139ff05d55SJeremy L Thompson case CEED_EVAL_DIV: 5149ff05d55SJeremy L Thompson case CEED_EVAL_CURL: 5159ff05d55SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 5169ff05d55SJeremy L Thompson // LCOV_EXCL_STOP 5179ff05d55SJeremy L Thompson } 5189ff05d55SJeremy L Thompson 5199ff05d55SJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 5209ff05d55SJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 5219ff05d55SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 5229ff05d55SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 5239ff05d55SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 5249ff05d55SJeremy L Thompson return CEED_ERROR_SUCCESS; 5259ff05d55SJeremy L Thompson } 5269ff05d55SJeremy L Thompson 5279ff05d55SJeremy L Thompson static int CeedBasisApplyNonTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, 5289ff05d55SJeremy L Thompson CeedVector u, CeedVector v) { 5299ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda_shared(basis, false, num_elem, t_mode, eval_mode, u, v)); 5309ff05d55SJeremy L Thompson return CEED_ERROR_SUCCESS; 5319ff05d55SJeremy L Thompson } 5329ff05d55SJeremy L Thompson 5339ff05d55SJeremy L Thompson static int CeedBasisApplyAddNonTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, 5349ff05d55SJeremy L Thompson CeedVector u, CeedVector v) { 5359ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda_shared(basis, true, num_elem, t_mode, eval_mode, u, v)); 5369ff05d55SJeremy L Thompson return CEED_ERROR_SUCCESS; 5379ff05d55SJeremy L Thompson } 5389ff05d55SJeremy L Thompson 5399ff05d55SJeremy L Thompson //------------------------------------------------------------------------------ 540ab213215SJeremy L Thompson // Destroy basis 541ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 542c532df63SYohann static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) { 543c532df63SYohann Ceed ceed; 544c532df63SYohann CeedBasis_Cuda_shared *data; 545ca735530SJeremy L Thompson 546ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 5472b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 5482b730f8bSJeremy L Thompson CeedCallCuda(ceed, cuModuleUnload(data->module)); 5491dda9c1aSJeremy L Thompson if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints)); 550097cc795SJames Wright if (data->d_q_weight_1d) CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d)); 551111870feSJeremy L Thompson CeedCallBackend(CeedFree(&data->h_points_per_elem)); 552111870feSJeremy L Thompson if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem)); 5532b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_interp_1d)); 5542b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_grad_1d)); 5552b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_collo_grad_1d)); 5561dda9c1aSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_chebyshev_interp_1d)); 5572b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&data)); 5589bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 559e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 560c532df63SYohann } 561c532df63SYohann 562ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 563ab213215SJeremy L Thompson // Create tensor basis 564ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 5652b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d, 5662b730f8bSJeremy L Thompson const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) { 567c532df63SYohann Ceed ceed; 568ca735530SJeremy L Thompson CeedInt num_comp; 569ca735530SJeremy L Thompson const CeedInt q_bytes = Q_1d * sizeof(CeedScalar); 570ca735530SJeremy L Thompson const CeedInt interp_bytes = q_bytes * P_1d; 571c532df63SYohann CeedBasis_Cuda_shared *data; 572ca735530SJeremy L Thompson 573ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 5742b730f8bSJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 575c532df63SYohann 576ab213215SJeremy L Thompson // Copy basis data to GPU 577097cc795SJames Wright if (q_weight_1d) { 5782b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes)); 5792b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice)); 580097cc795SJames Wright } 5812b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes)); 5822b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 5832b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes)); 5842b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice)); 585c532df63SYohann 586ab213215SJeremy L Thompson // Compute collocated gradient and copy to GPU 587437930d1SJeremy L Thompson data->d_collo_grad_1d = NULL; 5889e201c85SYohann bool has_collocated_grad = dim == 3 && Q_1d >= P_1d; 589ca735530SJeremy L Thompson 5909e201c85SYohann if (has_collocated_grad) { 591437930d1SJeremy L Thompson CeedScalar *collo_grad_1d; 592ca735530SJeremy L Thompson 5932b730f8bSJeremy L Thompson CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d)); 5942b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d)); 5952b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d)); 5962b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, cudaMemcpyHostToDevice)); 5972b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&collo_grad_1d)); 598ac421f39SYohann } 599ac421f39SYohann 600ab213215SJeremy L Thompson // Compile basis kernels 6019c25dd66SJeremy L Thompson const char basis_kernel_source[] = "// Tensor basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-tensor.h>\n"; 6029c25dd66SJeremy L Thompson 6032b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 604eb7e6cafSJeremy L Thompson CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D", 605eb7e6cafSJeremy L Thompson CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), 606eb7e6cafSJeremy L Thompson "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_HAS_COLLOCATED_GRAD", has_collocated_grad)); 607eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 608eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 609db2becc9SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd)); 610eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad)); 611eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose)); 612db2becc9SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd)); 613eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 614c532df63SYohann 6152b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 616ab213215SJeremy L Thompson 617ab213215SJeremy L Thompson // Register backend functions 6182b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared)); 619db2becc9SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddTensor_Cuda_shared)); 6201dda9c1aSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAtPoints", CeedBasisApplyAtPoints_Cuda_shared)); 621d474dafeSZach Atkins CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAddAtPoints", CeedBasisApplyAddAtPoints_Cuda_shared)); 6222b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared)); 6239bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 624e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 625c532df63SYohann } 6262a86cc9dSSebastian Grimberg 627ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 6289ff05d55SJeremy L Thompson // Create non-tensor basis 6299ff05d55SJeremy L Thompson //------------------------------------------------------------------------------ 6309ff05d55SJeremy L Thompson int CeedBasisCreateH1_Cuda_shared(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 6319ff05d55SJeremy L Thompson const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 6329ff05d55SJeremy L Thompson Ceed ceed; 6339ff05d55SJeremy L Thompson CeedInt num_comp, q_comp_interp, q_comp_grad; 6349ff05d55SJeremy L Thompson const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 6359ff05d55SJeremy L Thompson CeedBasis_Cuda_shared *data; 6369ff05d55SJeremy L Thompson 6379ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 6389ff05d55SJeremy L Thompson 639fda26546SJeremy L Thompson // Check shared memory size 640fda26546SJeremy L Thompson { 641fda26546SJeremy L Thompson Ceed_Cuda *cuda_data; 642fda26546SJeremy L Thompson 643fda26546SJeremy L Thompson CeedCallBackend(CeedGetData(ceed, &cuda_data)); 644fda26546SJeremy L Thompson if (((size_t)num_nodes * (size_t)num_qpts * (size_t)dim + (size_t)CeedIntMax(num_nodes, num_qpts)) * sizeof(CeedScalar) > 645fda26546SJeremy L Thompson cuda_data->device_prop.sharedMemPerBlock) { 646fda26546SJeremy L Thompson CeedCallBackend(CeedBasisCreateH1Fallback(ceed, topo, dim, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis)); 647fda26546SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 648fda26546SJeremy L Thompson return CEED_ERROR_SUCCESS; 649fda26546SJeremy L Thompson } 650fda26546SJeremy L Thompson } 651fda26546SJeremy L Thompson 652fda26546SJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 653aa4002adSJeremy L Thompson 6549ff05d55SJeremy L Thompson // Copy basis data to GPU 6559ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 6569ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad)); 6579ff05d55SJeremy L Thompson if (q_weight) { 6589ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes)); 6599ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight, q_bytes, cudaMemcpyHostToDevice)); 6609ff05d55SJeremy L Thompson } 6619ff05d55SJeremy L Thompson if (interp) { 6629ff05d55SJeremy L Thompson const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 6639ff05d55SJeremy L Thompson 6649ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes)); 6659ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp, interp_bytes, cudaMemcpyHostToDevice)); 6669ff05d55SJeremy L Thompson } 6679ff05d55SJeremy L Thompson if (grad) { 6689ff05d55SJeremy L Thompson const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad; 6699ff05d55SJeremy L Thompson 6709ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, grad_bytes)); 6719ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad, grad_bytes, cudaMemcpyHostToDevice)); 6729ff05d55SJeremy L Thompson } 6739ff05d55SJeremy L Thompson 6749ff05d55SJeremy L Thompson // Compile basis kernels 6759ff05d55SJeremy L Thompson const char basis_kernel_source[] = "// Non-tensor basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-nontensor.h>\n"; 6769ff05d55SJeremy L Thompson 6779ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 6789ff05d55SJeremy L Thompson CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "T_1D", 6799ff05d55SJeremy L Thompson CeedIntMax(num_qpts, num_nodes), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp)); 6809ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 6819ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 6829ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd)); 6839ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad)); 6849ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose)); 6859ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd)); 6869ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 6879ff05d55SJeremy L Thompson 6889ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 6899ff05d55SJeremy L Thompson 6909ff05d55SJeremy L Thompson // Register backend functions 6919ff05d55SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda_shared)); 6929ff05d55SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Cuda_shared)); 6939ff05d55SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared)); 6949ff05d55SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 6959ff05d55SJeremy L Thompson return CEED_ERROR_SUCCESS; 6969ff05d55SJeremy L Thompson } 6979ff05d55SJeremy L Thompson 6989ff05d55SJeremy L Thompson //------------------------------------------------------------------------------ 699