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 502b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 512b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 52437930d1SJeremy L Thompson CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 53ca735530SJeremy L Thompson 54*aa4002adSJeremy L Thompson void *interp_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_u, &d_v}; 55ca735530SJeremy L Thompson 564d537eeaSYohann if (dim == 1) { 57dac82721SJeremy L Thompson // avoid >512 total threads 58dac82721SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1)); 59a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 60437930d1SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 61b2165e7aSSebastian Grimberg 629e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 63db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, 1, 64db2becc9SJeremy L Thompson elems_per_block, shared_mem, interp_args)); 659e201c85SYohann } else { 66eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args)); 679e201c85SYohann } 68074be161SYohann Dudouit } else if (dim == 2) { 69437930d1SJeremy L Thompson const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8}; 70437930d1SJeremy L Thompson // elems_per_block must be at least 1 712b730f8bSJeremy L Thompson CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1); 72a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 732b730f8bSJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 74b2165e7aSSebastian Grimberg 759e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 76db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, thread_1d, 77db2becc9SJeremy L Thompson elems_per_block, shared_mem, interp_args)); 789e201c85SYohann } else { 79eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 809e201c85SYohann } 81074be161SYohann Dudouit } else if (dim == 3) { 82437930d1SJeremy L Thompson CeedInt elems_per_block = 1; 83a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 842b730f8bSJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 85b2165e7aSSebastian Grimberg 869e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 87db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, thread_1d, 88db2becc9SJeremy L Thompson elems_per_block, shared_mem, interp_args)); 899e201c85SYohann } else { 90eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 91074be161SYohann Dudouit } 929e201c85SYohann } 93ab213215SJeremy L Thompson } break; 94ab213215SJeremy L Thompson case CEED_EVAL_GRAD: { 95437930d1SJeremy L Thompson CeedInt P_1d, Q_1d; 96ca735530SJeremy L Thompson 972b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 982b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 99437930d1SJeremy L Thompson CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 100*aa4002adSJeremy L Thompson CeedScalar *d_grad_1d = data->d_grad_1d; 101ca735530SJeremy L Thompson 1029e201c85SYohann if (data->d_collo_grad_1d) { 103*aa4002adSJeremy L Thompson d_grad_1d = data->d_collo_grad_1d; 1049e201c85SYohann } 105*aa4002adSJeremy L Thompson void *grad_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_grad_1d, &d_u, &d_v}; 106*aa4002adSJeremy L Thompson 1074d537eeaSYohann if (dim == 1) { 108dac82721SJeremy L Thompson // avoid >512 total threads 109dac82721SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1)); 110a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 111437930d1SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 112b2165e7aSSebastian Grimberg 1139e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 114db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, 1, 115db2becc9SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 1169e201c85SYohann } else { 117eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args)); 1189e201c85SYohann } 119074be161SYohann Dudouit } else if (dim == 2) { 120437930d1SJeremy L Thompson const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8}; 121437930d1SJeremy L Thompson // elems_per_block must be at least 1 1222b730f8bSJeremy L Thompson CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1); 123a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 1242b730f8bSJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 125b2165e7aSSebastian Grimberg 1269e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 127db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, thread_1d, 128db2becc9SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 1299e201c85SYohann } else { 130eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 1319e201c85SYohann } 132074be161SYohann Dudouit } else if (dim == 3) { 133437930d1SJeremy L Thompson CeedInt elems_per_block = 1; 134a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 1352b730f8bSJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 136b2165e7aSSebastian Grimberg 1379e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 138db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, thread_1d, 139db2becc9SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 1409e201c85SYohann } else { 141eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 1429e201c85SYohann } 143074be161SYohann Dudouit } 144ab213215SJeremy L Thompson } break; 145ab213215SJeremy L Thompson case CEED_EVAL_WEIGHT: { 146437930d1SJeremy L Thompson CeedInt Q_1d; 147397164e9SSebastian Grimberg CeedInt block_size = 32; 148ca735530SJeremy L Thompson 149097cc795SJames Wright CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights_1d not set", CeedEvalModes[eval_mode]); 1502b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 151437930d1SJeremy L Thompson void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v}; 152074be161SYohann Dudouit if (dim == 1) { 153397164e9SSebastian Grimberg const CeedInt elems_per_block = block_size / Q_1d; 154a8d440fbSJeremy L Thompson const CeedInt grid_size = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 155b2165e7aSSebastian Grimberg 156b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, elems_per_block, 1, weight_args)); 157074be161SYohann Dudouit } else if (dim == 2) { 158397164e9SSebastian Grimberg const CeedInt opt_elems = block_size / (Q_1d * Q_1d); 159437930d1SJeremy L Thompson const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1; 160a8d440fbSJeremy L Thompson const CeedInt grid_size = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 161b2165e7aSSebastian Grimberg 162b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args)); 163074be161SYohann Dudouit } else if (dim == 3) { 164397164e9SSebastian Grimberg const CeedInt opt_elems = block_size / (Q_1d * Q_1d); 1659e201c85SYohann const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1; 166a8d440fbSJeremy L Thompson const CeedInt grid_size = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 167b2165e7aSSebastian Grimberg 168b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args)); 169074be161SYohann Dudouit } 170ab213215SJeremy L Thompson } break; 1719ea2cfd9SJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 1729ea2cfd9SJeremy L Thompson break; 173ab213215SJeremy L Thompson // LCOV_EXCL_START 174ab213215SJeremy L Thompson case CEED_EVAL_DIV: 175ab213215SJeremy L Thompson case CEED_EVAL_CURL: 176bcbe1c99SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 177ab213215SJeremy L Thompson // LCOV_EXCL_STOP 178c532df63SYohann } 179c532df63SYohann 1809ea2cfd9SJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 1812b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 1829ea2cfd9SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 1839ea2cfd9SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 1849bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 185e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 186c532df63SYohann } 187c532df63SYohann 188db2becc9SJeremy L Thompson static int CeedBasisApplyTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, 189db2becc9SJeremy L Thompson CeedVector v) { 190db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyTensorCore_Cuda_shared(basis, false, num_elem, t_mode, eval_mode, u, v)); 191db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 192db2becc9SJeremy L Thompson } 193db2becc9SJeremy L Thompson 194db2becc9SJeremy L Thompson static int CeedBasisApplyAddTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, 195db2becc9SJeremy L Thompson CeedVector u, CeedVector v) { 196db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyTensorCore_Cuda_shared(basis, true, num_elem, t_mode, eval_mode, u, v)); 197db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 198db2becc9SJeremy L Thompson } 199db2becc9SJeremy L Thompson 200ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 2011dda9c1aSJeremy L Thompson // Basis apply - tensor AtPoints 2021dda9c1aSJeremy L Thompson //------------------------------------------------------------------------------ 203db2becc9SJeremy L Thompson static int CeedBasisApplyAtPointsCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, const CeedInt *num_points, 204db2becc9SJeremy L Thompson CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) { 2051dda9c1aSJeremy L Thompson Ceed ceed; 2069e1d4b82SJeremy L Thompson Ceed_Cuda *ceed_Cuda; 2079e1d4b82SJeremy L Thompson CeedInt Q_1d, dim, num_comp, max_num_points = num_points[0]; 2081dda9c1aSJeremy L Thompson const CeedInt is_transpose = t_mode == CEED_TRANSPOSE; 2091dda9c1aSJeremy L Thompson const CeedScalar *d_x, *d_u; 2101dda9c1aSJeremy L Thompson CeedScalar *d_v; 2111dda9c1aSJeremy L Thompson CeedBasis_Cuda_shared *data; 2121dda9c1aSJeremy L Thompson 2131dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 2141dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 2151dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 2169e1d4b82SJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 2171dda9c1aSJeremy L Thompson 2181dda9c1aSJeremy L Thompson // Weight handled separately 2191dda9c1aSJeremy L Thompson if (eval_mode == CEED_EVAL_WEIGHT) { 2205a5594ffSJeremy L Thompson CeedCallBackend(CeedVectorSetValue(v, 1.0)); 2211dda9c1aSJeremy L Thompson return CEED_ERROR_SUCCESS; 2221dda9c1aSJeremy L Thompson } 2231dda9c1aSJeremy L Thompson 2249bc66399SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 2259e1d4b82SJeremy L Thompson CeedCallBackend(CeedGetData(ceed, &ceed_Cuda)); 2269bc66399SJeremy L Thompson 227111870feSJeremy L Thompson // Check padded to uniform number of points per elem 228111870feSJeremy L Thompson for (CeedInt i = 1; i < num_elem; i++) max_num_points = CeedIntMax(max_num_points, num_points[i]); 229111870feSJeremy L Thompson { 2309e1d4b82SJeremy L Thompson CeedInt q_comp; 231111870feSJeremy L Thompson CeedSize len, len_required; 232111870feSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 233111870feSJeremy L Thompson CeedCallBackend(CeedVectorGetLength(is_transpose ? u : v, &len)); 234111870feSJeremy L Thompson len_required = (CeedSize)num_comp * (CeedSize)q_comp * (CeedSize)num_elem * (CeedSize)max_num_points; 235111870feSJeremy L Thompson CeedCheck(len >= len_required, ceed, CEED_ERROR_BACKEND, 236111870feSJeremy L Thompson "Vector at points must be padded to the same number of points in each element for BasisApplyAtPoints on GPU backends." 237111870feSJeremy L Thompson " Found %" CeedSize_FMT ", Required %" CeedSize_FMT, 238111870feSJeremy L Thompson len, len_required); 239111870feSJeremy L Thompson } 240111870feSJeremy L Thompson 241111870feSJeremy L Thompson // Move num_points array to device 242111870feSJeremy L Thompson if (is_transpose) { 243111870feSJeremy L Thompson const CeedInt num_bytes = num_elem * sizeof(CeedInt); 244111870feSJeremy L Thompson 245111870feSJeremy L Thompson if (num_elem != data->num_elem_at_points) { 246111870feSJeremy L Thompson data->num_elem_at_points = num_elem; 247111870feSJeremy L Thompson 248111870feSJeremy L Thompson if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem)); 249111870feSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_points_per_elem, num_bytes)); 250111870feSJeremy L Thompson CeedCallBackend(CeedFree(&data->h_points_per_elem)); 251111870feSJeremy L Thompson CeedCallBackend(CeedCalloc(num_elem, &data->h_points_per_elem)); 252111870feSJeremy L Thompson } 2539e511c80SJeremy L Thompson if (memcmp(data->h_points_per_elem, num_points, num_bytes)) { 254111870feSJeremy L Thompson memcpy(data->h_points_per_elem, num_points, num_bytes); 255111870feSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_points_per_elem, num_points, num_bytes, cudaMemcpyHostToDevice)); 256111870feSJeremy L Thompson } 257111870feSJeremy L Thompson } 258111870feSJeremy L Thompson 2591dda9c1aSJeremy L Thompson // Build kernels if needed 2601dda9c1aSJeremy L Thompson if (data->num_points != max_num_points) { 2611dda9c1aSJeremy L Thompson CeedInt P_1d; 2621dda9c1aSJeremy L Thompson 2631dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 2641dda9c1aSJeremy L Thompson data->num_points = max_num_points; 2651dda9c1aSJeremy L Thompson 2661dda9c1aSJeremy L Thompson // -- Create interp matrix to Chebyshev coefficients 2671dda9c1aSJeremy L Thompson if (!data->d_chebyshev_interp_1d) { 2681dda9c1aSJeremy L Thompson CeedSize interp_bytes; 2691dda9c1aSJeremy L Thompson CeedScalar *chebyshev_interp_1d; 2701dda9c1aSJeremy L Thompson 2711dda9c1aSJeremy L Thompson interp_bytes = P_1d * Q_1d * sizeof(CeedScalar); 2721dda9c1aSJeremy L Thompson CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d)); 2735a5594ffSJeremy L Thompson CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d)); 2741dda9c1aSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_chebyshev_interp_1d, interp_bytes)); 2751dda9c1aSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 2761dda9c1aSJeremy L Thompson CeedCallBackend(CeedFree(&chebyshev_interp_1d)); 2771dda9c1aSJeremy L Thompson } 2781dda9c1aSJeremy L Thompson 2791dda9c1aSJeremy L Thompson // -- Compile kernels 2809e1d4b82SJeremy L Thompson const char basis_kernel_source[] = "// AtPoints basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-at-points.h>\n"; 2811dda9c1aSJeremy L Thompson CeedInt num_comp; 2821dda9c1aSJeremy L Thompson 2831dda9c1aSJeremy L Thompson if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints)); 2841dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 2859e1d4b82SJeremy L Thompson CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->moduleAtPoints, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D", 2869e1d4b82SJeremy L Thompson CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), 2879e1d4b82SJeremy L Thompson "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_NUM_PTS", max_num_points)); 2881dda9c1aSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpAtPoints", &data->InterpAtPoints)); 28981ae6159SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpTransposeAtPoints", &data->InterpTransposeAtPoints)); 2901dda9c1aSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradAtPoints", &data->GradAtPoints)); 29181ae6159SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradTransposeAtPoints", &data->GradTransposeAtPoints)); 2921dda9c1aSJeremy L Thompson } 2931dda9c1aSJeremy L Thompson 2941dda9c1aSJeremy L Thompson // Get read/write access to u, v 2951dda9c1aSJeremy L Thompson CeedCallBackend(CeedVectorGetArrayRead(x_ref, CEED_MEM_DEVICE, &d_x)); 2961dda9c1aSJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 2971dda9c1aSJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 298db2becc9SJeremy L Thompson if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 299db2becc9SJeremy L Thompson else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 3001dda9c1aSJeremy L Thompson 3011dda9c1aSJeremy L Thompson // Clear v for transpose operation 302db2becc9SJeremy L Thompson if (is_transpose && !apply_add) { 30319a04db8SJeremy L Thompson CeedInt num_comp, q_comp, num_nodes; 3041dda9c1aSJeremy L Thompson CeedSize length; 3051dda9c1aSJeremy L Thompson 30619a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 30719a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 30819a04db8SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes)); 30919a04db8SJeremy L Thompson length = 31019a04db8SJeremy L Thompson (CeedSize)num_elem * (CeedSize)num_comp * (t_mode == CEED_TRANSPOSE ? (CeedSize)num_nodes : ((CeedSize)max_num_points * (CeedSize)q_comp)); 3111dda9c1aSJeremy L Thompson CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar))); 3121dda9c1aSJeremy L Thompson } 3131dda9c1aSJeremy L Thompson 3141dda9c1aSJeremy L Thompson // Basis action 3151dda9c1aSJeremy L Thompson switch (eval_mode) { 3161dda9c1aSJeremy L Thompson case CEED_EVAL_INTERP: { 3179e1d4b82SJeremy L Thompson CeedInt P_1d, Q_1d; 3181dda9c1aSJeremy L Thompson 3199e1d4b82SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 3209e1d4b82SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 3219e1d4b82SJeremy L Thompson CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 3229e1d4b82SJeremy L Thompson 323*aa4002adSJeremy L Thompson void *interp_args[] = {(void *)&num_elem, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v}; 3249e1d4b82SJeremy L Thompson 3259e1d4b82SJeremy L Thompson if (dim == 1) { 326dac82721SJeremy L Thompson // avoid >512 total threads 327dac82721SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1)); 328a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3299e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 3309e1d4b82SJeremy L Thompson 3319e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d, 1, 3329e1d4b82SJeremy L Thompson elems_per_block, shared_mem, interp_args)); 3339e1d4b82SJeremy L Thompson } else if (dim == 2) { 3349e1d4b82SJeremy L Thompson const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8}; 3359e1d4b82SJeremy L Thompson // elems_per_block must be at least 1 3369e1d4b82SJeremy L Thompson CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1); 337a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3389e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 3399e1d4b82SJeremy L Thompson 3409e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d, 3419e1d4b82SJeremy L Thompson thread_1d, elems_per_block, shared_mem, interp_args)); 3429e1d4b82SJeremy L Thompson } else if (dim == 3) { 3439e1d4b82SJeremy L Thompson CeedInt elems_per_block = 1; 344a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3459e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 3469e1d4b82SJeremy L Thompson 3479e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d, 3489e1d4b82SJeremy L Thompson thread_1d, elems_per_block, shared_mem, interp_args)); 3499e1d4b82SJeremy L Thompson } 3501dda9c1aSJeremy L Thompson } break; 3511dda9c1aSJeremy L Thompson case CEED_EVAL_GRAD: { 3529e1d4b82SJeremy L Thompson CeedInt P_1d, Q_1d; 3531dda9c1aSJeremy L Thompson 3549e1d4b82SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 3559e1d4b82SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 3569e1d4b82SJeremy L Thompson CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 3579e1d4b82SJeremy L Thompson 3589e1d4b82SJeremy L Thompson void *grad_args[] = {(void *)&num_elem, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v}; 3599e1d4b82SJeremy L Thompson 3609e1d4b82SJeremy L Thompson if (dim == 1) { 361dac82721SJeremy L Thompson // avoid >512 total threads 362dac82721SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1)); 363a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3649e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 3659e1d4b82SJeremy L Thompson 3669e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, 1, 3679e1d4b82SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 3689e1d4b82SJeremy L Thompson } else if (dim == 2) { 3699e1d4b82SJeremy L Thompson const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8}; 3709e1d4b82SJeremy L Thompson // elems_per_block must be at least 1 3719e1d4b82SJeremy L Thompson CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1); 372a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3739e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 3749e1d4b82SJeremy L Thompson 3759e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, thread_1d, 3769e1d4b82SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 3779e1d4b82SJeremy L Thompson } else if (dim == 3) { 3789e1d4b82SJeremy L Thompson CeedInt elems_per_block = 1; 379a8d440fbSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 3809e1d4b82SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 3819e1d4b82SJeremy L Thompson 3829e1d4b82SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, thread_1d, 3839e1d4b82SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 3849e1d4b82SJeremy L Thompson } 3851dda9c1aSJeremy L Thompson } break; 3861dda9c1aSJeremy L Thompson case CEED_EVAL_WEIGHT: 3871dda9c1aSJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 3881dda9c1aSJeremy L Thompson break; 3891dda9c1aSJeremy L Thompson // LCOV_EXCL_START 3901dda9c1aSJeremy L Thompson case CEED_EVAL_DIV: 3911dda9c1aSJeremy L Thompson case CEED_EVAL_CURL: 3921dda9c1aSJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 3931dda9c1aSJeremy L Thompson // LCOV_EXCL_STOP 3941dda9c1aSJeremy L Thompson } 3951dda9c1aSJeremy L Thompson 3961dda9c1aSJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 3971dda9c1aSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArrayRead(x_ref, &d_x)); 3981dda9c1aSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 3991dda9c1aSJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 4001dda9c1aSJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 4019bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 4021dda9c1aSJeremy L Thompson return CEED_ERROR_SUCCESS; 4031dda9c1aSJeremy L Thompson } 4041dda9c1aSJeremy L Thompson 405db2becc9SJeremy L Thompson static int CeedBasisApplyAtPoints_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, false, 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 411db2becc9SJeremy L Thompson static int CeedBasisApplyAddAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, 412db2becc9SJeremy L Thompson CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) { 413db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda_shared(basis, true, num_elem, num_points, t_mode, eval_mode, x_ref, u, v)); 414db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 415db2becc9SJeremy L Thompson } 416db2becc9SJeremy L Thompson 4171dda9c1aSJeremy L Thompson //------------------------------------------------------------------------------ 4189ff05d55SJeremy L Thompson // Apply non-tensor basis 4199ff05d55SJeremy L Thompson //------------------------------------------------------------------------------ 4209ff05d55SJeremy L Thompson static int CeedBasisApplyNonTensorCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode, 4219ff05d55SJeremy L Thompson CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 4229ff05d55SJeremy L Thompson Ceed ceed; 4239ff05d55SJeremy L Thompson Ceed_Cuda *ceed_Cuda; 4249ff05d55SJeremy L Thompson CeedInt dim; 4259ff05d55SJeremy L Thompson const CeedScalar *d_u; 4269ff05d55SJeremy L Thompson CeedScalar *d_v; 4279ff05d55SJeremy L Thompson CeedBasis_Cuda_shared *data; 4289ff05d55SJeremy L Thompson 4299ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 4309ff05d55SJeremy L Thompson CeedCallBackend(CeedGetData(ceed, &ceed_Cuda)); 4319ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 4329ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 4339ff05d55SJeremy L Thompson 4349ff05d55SJeremy L Thompson // Get read/write access to u, v 4359ff05d55SJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 4369ff05d55SJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 4379ff05d55SJeremy L Thompson if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 4389ff05d55SJeremy L Thompson else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 4399ff05d55SJeremy L Thompson 4409ff05d55SJeremy L Thompson // Apply basis operation 4419ff05d55SJeremy L Thompson switch (eval_mode) { 4429ff05d55SJeremy L Thompson case CEED_EVAL_INTERP: { 4439ff05d55SJeremy L Thompson CeedInt P, Q; 4449ff05d55SJeremy L Thompson 4459ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes(basis, &P)); 4469ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q)); 4479ff05d55SJeremy L Thompson CeedInt thread = CeedIntMax(Q, P); 4489ff05d55SJeremy L Thompson 449*aa4002adSJeremy L Thompson void *interp_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_u, &d_v}; 4509ff05d55SJeremy L Thompson 4519ff05d55SJeremy L Thompson { 4529ff05d55SJeremy L Thompson // avoid >512 total threads 4539ff05d55SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1)); 4549ff05d55SJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 4559ff05d55SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread * sizeof(CeedScalar); 4569ff05d55SJeremy L Thompson 4579ff05d55SJeremy L Thompson if (t_mode == CEED_TRANSPOSE) { 4589ff05d55SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread, 1, 4599ff05d55SJeremy L Thompson elems_per_block, shared_mem, interp_args)); 4609ff05d55SJeremy L Thompson } else { 4619ff05d55SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread, 1, elems_per_block, shared_mem, interp_args)); 4629ff05d55SJeremy L Thompson } 4639ff05d55SJeremy L Thompson } 4649ff05d55SJeremy L Thompson } break; 4659ff05d55SJeremy L Thompson case CEED_EVAL_GRAD: { 4669ff05d55SJeremy L Thompson CeedInt P, Q; 4679ff05d55SJeremy L Thompson 4689ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes(basis, &P)); 4699ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q)); 4709ff05d55SJeremy L Thompson CeedInt thread = CeedIntMax(Q, P); 4719ff05d55SJeremy L Thompson 472*aa4002adSJeremy L Thompson void *grad_args[] = {(void *)&num_elem, &data->d_grad_1d, &d_u, &d_v}; 4739ff05d55SJeremy L Thompson 4749ff05d55SJeremy L Thompson { 4759ff05d55SJeremy L Thompson // avoid >512 total threads 4769ff05d55SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1)); 4779ff05d55SJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 4789ff05d55SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread * sizeof(CeedScalar); 4799ff05d55SJeremy L Thompson 4809ff05d55SJeremy L Thompson if (t_mode == CEED_TRANSPOSE) { 4819ff05d55SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread, 1, 4829ff05d55SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 4839ff05d55SJeremy L Thompson } else { 4849ff05d55SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread, 1, elems_per_block, shared_mem, grad_args)); 4859ff05d55SJeremy L Thompson } 4869ff05d55SJeremy L Thompson } 4879ff05d55SJeremy L Thompson } break; 4889ff05d55SJeremy L Thompson case CEED_EVAL_WEIGHT: { 4899ff05d55SJeremy L Thompson CeedInt Q; 4909ff05d55SJeremy L Thompson 4919ff05d55SJeremy L Thompson CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights_1d not set", CeedEvalModes[eval_mode]); 4929ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q)); 4939ff05d55SJeremy L Thompson void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v}; 4949ff05d55SJeremy L Thompson 4959ff05d55SJeremy L Thompson { 4969ff05d55SJeremy L Thompson // avoid >512 total threads 4979ff05d55SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / Q, 1)); 4989ff05d55SJeremy L Thompson CeedInt grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 4999ff05d55SJeremy L Thompson 5009ff05d55SJeremy L Thompson CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, Q, elems_per_block, 1, weight_args)); 5019ff05d55SJeremy L Thompson } 5029ff05d55SJeremy L Thompson } break; 5039ff05d55SJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 5049ff05d55SJeremy L Thompson break; 5059ff05d55SJeremy L Thompson // LCOV_EXCL_START 5069ff05d55SJeremy L Thompson case CEED_EVAL_DIV: 5079ff05d55SJeremy L Thompson case CEED_EVAL_CURL: 5089ff05d55SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 5099ff05d55SJeremy L Thompson // LCOV_EXCL_STOP 5109ff05d55SJeremy L Thompson } 5119ff05d55SJeremy L Thompson 5129ff05d55SJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 5139ff05d55SJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 5149ff05d55SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 5159ff05d55SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 5169ff05d55SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 5179ff05d55SJeremy L Thompson return CEED_ERROR_SUCCESS; 5189ff05d55SJeremy L Thompson } 5199ff05d55SJeremy L Thompson 5209ff05d55SJeremy L Thompson static int CeedBasisApplyNonTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, 5219ff05d55SJeremy L Thompson CeedVector u, CeedVector v) { 5229ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda_shared(basis, false, num_elem, t_mode, eval_mode, u, v)); 5239ff05d55SJeremy L Thompson return CEED_ERROR_SUCCESS; 5249ff05d55SJeremy L Thompson } 5259ff05d55SJeremy L Thompson 5269ff05d55SJeremy L Thompson static int CeedBasisApplyAddNonTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, 5279ff05d55SJeremy L Thompson CeedVector u, CeedVector v) { 5289ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda_shared(basis, true, num_elem, t_mode, eval_mode, u, v)); 5299ff05d55SJeremy L Thompson return CEED_ERROR_SUCCESS; 5309ff05d55SJeremy L Thompson } 5319ff05d55SJeremy L Thompson 5329ff05d55SJeremy L Thompson //------------------------------------------------------------------------------ 533ab213215SJeremy L Thompson // Destroy basis 534ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 535c532df63SYohann static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) { 536c532df63SYohann Ceed ceed; 537c532df63SYohann CeedBasis_Cuda_shared *data; 538ca735530SJeremy L Thompson 539ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 5402b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 5412b730f8bSJeremy L Thompson CeedCallCuda(ceed, cuModuleUnload(data->module)); 5421dda9c1aSJeremy L Thompson if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints)); 543097cc795SJames Wright if (data->d_q_weight_1d) CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d)); 544111870feSJeremy L Thompson CeedCallBackend(CeedFree(&data->h_points_per_elem)); 545111870feSJeremy L Thompson if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem)); 5462b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_interp_1d)); 5472b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_grad_1d)); 5482b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_collo_grad_1d)); 5491dda9c1aSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_chebyshev_interp_1d)); 5502b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&data)); 5519bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 552e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 553c532df63SYohann } 554c532df63SYohann 555ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 556ab213215SJeremy L Thompson // Create tensor basis 557ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 5582b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d, 5592b730f8bSJeremy L Thompson const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) { 560c532df63SYohann Ceed ceed; 561ca735530SJeremy L Thompson CeedInt num_comp; 562ca735530SJeremy L Thompson const CeedInt q_bytes = Q_1d * sizeof(CeedScalar); 563ca735530SJeremy L Thompson const CeedInt interp_bytes = q_bytes * P_1d; 564c532df63SYohann CeedBasis_Cuda_shared *data; 565ca735530SJeremy L Thompson 566ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 5672b730f8bSJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 568c532df63SYohann 569ab213215SJeremy L Thompson // Copy basis data to GPU 570097cc795SJames Wright if (q_weight_1d) { 5712b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes)); 5722b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice)); 573097cc795SJames Wright } 5742b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes)); 5752b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 5762b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes)); 5772b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice)); 578c532df63SYohann 579ab213215SJeremy L Thompson // Compute collocated gradient and copy to GPU 580437930d1SJeremy L Thompson data->d_collo_grad_1d = NULL; 5819e201c85SYohann bool has_collocated_grad = dim == 3 && Q_1d >= P_1d; 582ca735530SJeremy L Thompson 5839e201c85SYohann if (has_collocated_grad) { 584437930d1SJeremy L Thompson CeedScalar *collo_grad_1d; 585ca735530SJeremy L Thompson 5862b730f8bSJeremy L Thompson CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d)); 5872b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d)); 5882b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d)); 5892b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, cudaMemcpyHostToDevice)); 5902b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&collo_grad_1d)); 591ac421f39SYohann } 592ac421f39SYohann 593ab213215SJeremy L Thompson // Compile basis kernels 5949c25dd66SJeremy L Thompson const char basis_kernel_source[] = "// Tensor basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-tensor.h>\n"; 5959c25dd66SJeremy L Thompson 5962b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 597eb7e6cafSJeremy L Thompson CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D", 598eb7e6cafSJeremy L Thompson CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), 599eb7e6cafSJeremy L Thompson "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_HAS_COLLOCATED_GRAD", has_collocated_grad)); 600eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 601eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 602db2becc9SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd)); 603eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad)); 604eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose)); 605db2becc9SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd)); 606eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 607c532df63SYohann 6082b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 609ab213215SJeremy L Thompson 610ab213215SJeremy L Thompson // Register backend functions 6112b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared)); 612db2becc9SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddTensor_Cuda_shared)); 6131dda9c1aSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAtPoints", CeedBasisApplyAtPoints_Cuda_shared)); 614d474dafeSZach Atkins CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAddAtPoints", CeedBasisApplyAddAtPoints_Cuda_shared)); 6152b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared)); 6169bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 617e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 618c532df63SYohann } 6192a86cc9dSSebastian Grimberg 620ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 6219ff05d55SJeremy L Thompson // Create non-tensor basis 6229ff05d55SJeremy L Thompson //------------------------------------------------------------------------------ 6239ff05d55SJeremy L Thompson int CeedBasisCreateH1_Cuda_shared(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 6249ff05d55SJeremy L Thompson const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) { 6259ff05d55SJeremy L Thompson Ceed ceed; 6269ff05d55SJeremy L Thompson CeedInt num_comp, q_comp_interp, q_comp_grad; 6279ff05d55SJeremy L Thompson const CeedInt q_bytes = num_qpts * sizeof(CeedScalar); 6289ff05d55SJeremy L Thompson CeedBasis_Cuda_shared *data; 6299ff05d55SJeremy L Thompson 6309ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 6319ff05d55SJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 6329ff05d55SJeremy L Thompson 633*aa4002adSJeremy L Thompson // Check max sizes 634*aa4002adSJeremy L Thompson CeedCheck(dim <= 3, ceed, CEED_ERROR_BACKEND, "Backend does not implement nontensor bases with dim > 3"); 635*aa4002adSJeremy L Thompson CeedCheck(num_nodes * num_qpts * dim < 52 * 52 * 3, ceed, CEED_ERROR_BACKEND, "Backend does not implement nontensor bases with P * Q this large"); 636*aa4002adSJeremy L Thompson 6379ff05d55SJeremy L Thompson // Copy basis data to GPU 6389ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 6399ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad)); 6409ff05d55SJeremy L Thompson if (q_weight) { 6419ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes)); 6429ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight, q_bytes, cudaMemcpyHostToDevice)); 6439ff05d55SJeremy L Thompson } 6449ff05d55SJeremy L Thompson if (interp) { 6459ff05d55SJeremy L Thompson const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp; 6469ff05d55SJeremy L Thompson 6479ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes)); 6489ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp, interp_bytes, cudaMemcpyHostToDevice)); 6499ff05d55SJeremy L Thompson } 6509ff05d55SJeremy L Thompson if (grad) { 6519ff05d55SJeremy L Thompson const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad; 6529ff05d55SJeremy L Thompson 6539ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, grad_bytes)); 6549ff05d55SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad, grad_bytes, cudaMemcpyHostToDevice)); 6559ff05d55SJeremy L Thompson } 6569ff05d55SJeremy L Thompson 6579ff05d55SJeremy L Thompson // Compile basis kernels 6589ff05d55SJeremy L Thompson const char basis_kernel_source[] = "// Non-tensor basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-nontensor.h>\n"; 6599ff05d55SJeremy L Thompson 6609ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 6619ff05d55SJeremy L Thompson CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "T_1D", 6629ff05d55SJeremy L Thompson CeedIntMax(num_qpts, num_nodes), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp)); 6639ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 6649ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 6659ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd)); 6669ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad)); 6679ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose)); 6689ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd)); 6699ff05d55SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 6709ff05d55SJeremy L Thompson 6719ff05d55SJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 6729ff05d55SJeremy L Thompson 6739ff05d55SJeremy L Thompson // Register backend functions 6749ff05d55SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda_shared)); 6759ff05d55SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Cuda_shared)); 6769ff05d55SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared)); 6779ff05d55SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 6789ff05d55SJeremy L Thompson return CEED_ERROR_SUCCESS; 6799ff05d55SJeremy L Thompson } 6809ff05d55SJeremy L Thompson 6819ff05d55SJeremy L Thompson //------------------------------------------------------------------------------ 682