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> 15c532df63SYohann 1649aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h" 172b730f8bSJeremy L Thompson #include "../cuda/ceed-cuda-compile.h" 182b730f8bSJeremy L Thompson #include "ceed-cuda-shared.h" 19c532df63SYohann 20ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 21ab213215SJeremy L Thompson // Device initalization 22ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 23eb7e6cafSJeremy L Thompson int CeedInit_CudaInterp(CeedScalar *d_B, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B); 24eb7e6cafSJeremy L Thompson int CeedInit_CudaGrad(CeedScalar *d_B, CeedScalar *d_G, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B_ptr, CeedScalar **c_G_ptr); 25eb7e6cafSJeremy L Thompson int CeedInit_CudaCollocatedGrad(CeedScalar *d_B, CeedScalar *d_G, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B_ptr, CeedScalar **c_G_ptr); 26c532df63SYohann 27ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 28ab213215SJeremy L Thompson // Apply basis 29ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 30*db2becc9SJeremy L Thompson static int CeedBasisApplyTensorCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode, 31*db2becc9SJeremy L Thompson CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 32c532df63SYohann Ceed ceed; 336dbfb411Snbeams Ceed_Cuda *ceed_Cuda; 34437930d1SJeremy L Thompson CeedInt dim, num_comp; 35ca735530SJeremy L Thompson const CeedScalar *d_u; 36ca735530SJeremy L Thompson CeedScalar *d_v; 37ca735530SJeremy L Thompson CeedBasis_Cuda_shared *data; 38ca735530SJeremy L Thompson 39ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 40ca735530SJeremy L Thompson CeedCallBackend(CeedGetData(ceed, &ceed_Cuda)); 41ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 422b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 432b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 44c532df63SYohann 459ea2cfd9SJeremy L Thompson // Get read/write access to u, v 466574a04fSJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 476574a04fSJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 48*db2becc9SJeremy L Thompson if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 49*db2becc9SJeremy L Thompson else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 50c532df63SYohann 51ab213215SJeremy L Thompson // Apply basis operation 52437930d1SJeremy L Thompson switch (eval_mode) { 53ab213215SJeremy L Thompson case CEED_EVAL_INTERP: { 54437930d1SJeremy L Thompson CeedInt P_1d, Q_1d; 55ca735530SJeremy L Thompson 562b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 572b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 58437930d1SJeremy L Thompson CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 59ca735530SJeremy L Thompson 60eb7e6cafSJeremy L Thompson CeedCallBackend(CeedInit_CudaInterp(data->d_interp_1d, P_1d, Q_1d, &data->c_B)); 612b730f8bSJeremy L Thompson void *interp_args[] = {(void *)&num_elem, &data->c_B, &d_u, &d_v}; 62ca735530SJeremy L Thompson 634d537eeaSYohann if (dim == 1) { 642b730f8bSJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 65e6f67ff7SJed Brown 1)); // avoid >512 total threads 662b730f8bSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 67437930d1SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 68b2165e7aSSebastian Grimberg 699e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 70*db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, 1, 71*db2becc9SJeremy L Thompson elems_per_block, shared_mem, interp_args)); 729e201c85SYohann } else { 73eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args)); 749e201c85SYohann } 75074be161SYohann Dudouit } else if (dim == 2) { 76437930d1SJeremy L Thompson const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8}; 77437930d1SJeremy L Thompson // elems_per_block must be at least 1 782b730f8bSJeremy L Thompson CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1); 792b730f8bSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 802b730f8bSJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 81b2165e7aSSebastian Grimberg 829e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 83*db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, thread_1d, 84*db2becc9SJeremy L Thompson elems_per_block, shared_mem, interp_args)); 859e201c85SYohann } else { 86eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 879e201c85SYohann } 88074be161SYohann Dudouit } else if (dim == 3) { 89437930d1SJeremy L Thompson CeedInt elems_per_block = 1; 902b730f8bSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 912b730f8bSJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 92b2165e7aSSebastian Grimberg 939e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 94*db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, thread_1d, 95*db2becc9SJeremy L Thompson elems_per_block, shared_mem, interp_args)); 969e201c85SYohann } else { 97eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 98074be161SYohann Dudouit } 999e201c85SYohann } 100ab213215SJeremy L Thompson } break; 101ab213215SJeremy L Thompson case CEED_EVAL_GRAD: { 102437930d1SJeremy L Thompson CeedInt P_1d, Q_1d; 103ca735530SJeremy L Thompson 1042b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 1052b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 106437930d1SJeremy L Thompson CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 107ca735530SJeremy L Thompson 1089e201c85SYohann if (data->d_collo_grad_1d) { 109eb7e6cafSJeremy L Thompson CeedCallBackend(CeedInit_CudaCollocatedGrad(data->d_interp_1d, data->d_collo_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G)); 1109e201c85SYohann } else { 111eb7e6cafSJeremy L Thompson CeedCallBackend(CeedInit_CudaGrad(data->d_interp_1d, data->d_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G)); 1129e201c85SYohann } 1132b730f8bSJeremy L Thompson void *grad_args[] = {(void *)&num_elem, &data->c_B, &data->c_G, &d_u, &d_v}; 1144d537eeaSYohann if (dim == 1) { 1152b730f8bSJeremy L Thompson CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 116e6f67ff7SJed Brown 1)); // avoid >512 total threads 1172b730f8bSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 118437930d1SJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 119b2165e7aSSebastian Grimberg 1209e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 121*db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, 1, 122*db2becc9SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 1239e201c85SYohann } else { 124eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args)); 1259e201c85SYohann } 126074be161SYohann Dudouit } else if (dim == 2) { 127437930d1SJeremy L Thompson const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8}; 128437930d1SJeremy L Thompson // elems_per_block must be at least 1 1292b730f8bSJeremy L Thompson CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1); 1302b730f8bSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 1312b730f8bSJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 132b2165e7aSSebastian Grimberg 1339e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 134*db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, thread_1d, 135*db2becc9SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 1369e201c85SYohann } else { 137eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 1389e201c85SYohann } 139074be161SYohann Dudouit } else if (dim == 3) { 140437930d1SJeremy L Thompson CeedInt elems_per_block = 1; 1412b730f8bSJeremy L Thompson CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 1422b730f8bSJeremy L Thompson CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 143b2165e7aSSebastian Grimberg 1449e201c85SYohann if (t_mode == CEED_TRANSPOSE) { 145*db2becc9SJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, thread_1d, 146*db2becc9SJeremy L Thompson elems_per_block, shared_mem, grad_args)); 1479e201c85SYohann } else { 148eb7e6cafSJeremy L Thompson CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 1499e201c85SYohann } 150074be161SYohann Dudouit } 151ab213215SJeremy L Thompson } break; 152ab213215SJeremy L Thompson case CEED_EVAL_WEIGHT: { 153437930d1SJeremy L Thompson CeedInt Q_1d; 154397164e9SSebastian Grimberg CeedInt block_size = 32; 155ca735530SJeremy L Thompson 156097cc795SJames Wright CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights_1d not set", CeedEvalModes[eval_mode]); 1572b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 158437930d1SJeremy L Thompson void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v}; 159074be161SYohann Dudouit if (dim == 1) { 160397164e9SSebastian Grimberg const CeedInt elems_per_block = block_size / Q_1d; 161b2165e7aSSebastian Grimberg const CeedInt grid_size = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 162b2165e7aSSebastian Grimberg 163b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, elems_per_block, 1, weight_args)); 164074be161SYohann Dudouit } else if (dim == 2) { 165397164e9SSebastian Grimberg const CeedInt opt_elems = block_size / (Q_1d * Q_1d); 166437930d1SJeremy L Thompson const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1; 167b2165e7aSSebastian Grimberg const CeedInt grid_size = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 168b2165e7aSSebastian Grimberg 169b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args)); 170074be161SYohann Dudouit } else if (dim == 3) { 171397164e9SSebastian Grimberg const CeedInt opt_elems = block_size / (Q_1d * Q_1d); 1729e201c85SYohann const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1; 173b2165e7aSSebastian Grimberg const CeedInt grid_size = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 174b2165e7aSSebastian Grimberg 175b2165e7aSSebastian Grimberg CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args)); 176074be161SYohann Dudouit } 177ab213215SJeremy L Thompson } break; 1789ea2cfd9SJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 1799ea2cfd9SJeremy L Thompson break; 180ab213215SJeremy L Thompson // LCOV_EXCL_START 181ab213215SJeremy L Thompson case CEED_EVAL_DIV: 182ab213215SJeremy L Thompson case CEED_EVAL_CURL: 183bcbe1c99SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 184ab213215SJeremy L Thompson // LCOV_EXCL_STOP 185c532df63SYohann } 186c532df63SYohann 1879ea2cfd9SJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 1882b730f8bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 1899ea2cfd9SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 1909ea2cfd9SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 191e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 192c532df63SYohann } 193c532df63SYohann 194*db2becc9SJeremy L Thompson static int CeedBasisApplyTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, 195*db2becc9SJeremy L Thompson CeedVector v) { 196*db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyTensorCore_Cuda_shared(basis, false, num_elem, t_mode, eval_mode, u, v)); 197*db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 198*db2becc9SJeremy L Thompson } 199*db2becc9SJeremy L Thompson 200*db2becc9SJeremy L Thompson static int CeedBasisApplyAddTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, 201*db2becc9SJeremy L Thompson CeedVector u, CeedVector v) { 202*db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyTensorCore_Cuda_shared(basis, true, num_elem, t_mode, eval_mode, u, v)); 203*db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 204*db2becc9SJeremy L Thompson } 205*db2becc9SJeremy L Thompson 206ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 2071dda9c1aSJeremy L Thompson // Basis apply - tensor AtPoints 2081dda9c1aSJeremy L Thompson //------------------------------------------------------------------------------ 209*db2becc9SJeremy L Thompson static int CeedBasisApplyAtPointsCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, const CeedInt *num_points, 210*db2becc9SJeremy L Thompson CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) { 2111dda9c1aSJeremy L Thompson Ceed ceed; 2121dda9c1aSJeremy L Thompson CeedInt Q_1d, dim, max_num_points = num_points[0]; 2131dda9c1aSJeremy L Thompson const CeedInt is_transpose = t_mode == CEED_TRANSPOSE; 2141dda9c1aSJeremy L Thompson const int max_block_size = 32; 2151dda9c1aSJeremy L Thompson const CeedScalar *d_x, *d_u; 2161dda9c1aSJeremy L Thompson CeedScalar *d_v; 2171dda9c1aSJeremy L Thompson CeedBasis_Cuda_shared *data; 2181dda9c1aSJeremy L Thompson 2191dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 2201dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 2211dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 2221dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 2231dda9c1aSJeremy L Thompson 2241dda9c1aSJeremy L Thompson // Check uniform number of points per elem 2251dda9c1aSJeremy L Thompson for (CeedInt i = 1; i < num_elem; i++) { 2261dda9c1aSJeremy L Thompson CeedCheck(max_num_points == num_points[i], ceed, CEED_ERROR_BACKEND, 2271dda9c1aSJeremy L Thompson "BasisApplyAtPoints only supported for the same number of points in each element"); 2281dda9c1aSJeremy L Thompson } 2291dda9c1aSJeremy L Thompson 2301dda9c1aSJeremy L Thompson // Weight handled separately 2311dda9c1aSJeremy L Thompson if (eval_mode == CEED_EVAL_WEIGHT) { 2321dda9c1aSJeremy L Thompson CeedCall(CeedVectorSetValue(v, 1.0)); 2331dda9c1aSJeremy L Thompson return CEED_ERROR_SUCCESS; 2341dda9c1aSJeremy L Thompson } 2351dda9c1aSJeremy L Thompson 2361dda9c1aSJeremy L Thompson // Build kernels if needed 2371dda9c1aSJeremy L Thompson if (data->num_points != max_num_points) { 2381dda9c1aSJeremy L Thompson CeedInt P_1d; 2391dda9c1aSJeremy L Thompson 2401dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 2411dda9c1aSJeremy L Thompson data->num_points = max_num_points; 2421dda9c1aSJeremy L Thompson 2431dda9c1aSJeremy L Thompson // -- Create interp matrix to Chebyshev coefficients 2441dda9c1aSJeremy L Thompson if (!data->d_chebyshev_interp_1d) { 2451dda9c1aSJeremy L Thompson CeedSize interp_bytes; 2461dda9c1aSJeremy L Thompson CeedScalar *chebyshev_interp_1d; 2471dda9c1aSJeremy L Thompson 2481dda9c1aSJeremy L Thompson interp_bytes = P_1d * Q_1d * sizeof(CeedScalar); 2491dda9c1aSJeremy L Thompson CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d)); 2501dda9c1aSJeremy L Thompson CeedCall(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d)); 2511dda9c1aSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_chebyshev_interp_1d, interp_bytes)); 2521dda9c1aSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 2531dda9c1aSJeremy L Thompson CeedCallBackend(CeedFree(&chebyshev_interp_1d)); 2541dda9c1aSJeremy L Thompson } 2551dda9c1aSJeremy L Thompson 2561dda9c1aSJeremy L Thompson // -- Compile kernels 2571dda9c1aSJeremy L Thompson char *basis_kernel_source; 2581dda9c1aSJeremy L Thompson const char *basis_kernel_path; 2591dda9c1aSJeremy L Thompson CeedInt num_comp; 2601dda9c1aSJeremy L Thompson 2611dda9c1aSJeremy L Thompson if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints)); 2621dda9c1aSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 2631dda9c1aSJeremy L Thompson CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-tensor-at-points.h", &basis_kernel_path)); 2641dda9c1aSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 2651dda9c1aSJeremy L Thompson CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 2661dda9c1aSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n"); 2671dda9c1aSJeremy L Thompson CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->moduleAtPoints, 9, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "BASIS_BUF_LEN", 268f7c9815fSJeremy L Thompson Q_1d * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim - 1), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, 2691dda9c1aSJeremy L Thompson "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_NUM_PTS", 270f7c9815fSJeremy L Thompson max_num_points, "POINTS_BUFF_LEN", CeedIntPow(Q_1d, dim - 1))); 2711dda9c1aSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpAtPoints", &data->InterpAtPoints)); 2721dda9c1aSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradAtPoints", &data->GradAtPoints)); 2731dda9c1aSJeremy L Thompson CeedCallBackend(CeedFree(&basis_kernel_path)); 2741dda9c1aSJeremy L Thompson CeedCallBackend(CeedFree(&basis_kernel_source)); 2751dda9c1aSJeremy L Thompson } 2761dda9c1aSJeremy L Thompson 2771dda9c1aSJeremy L Thompson // Get read/write access to u, v 2781dda9c1aSJeremy L Thompson CeedCallBackend(CeedVectorGetArrayRead(x_ref, CEED_MEM_DEVICE, &d_x)); 2791dda9c1aSJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 2801dda9c1aSJeremy L Thompson else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 281*db2becc9SJeremy L Thompson if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 282*db2becc9SJeremy L Thompson else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 2831dda9c1aSJeremy L Thompson 2841dda9c1aSJeremy L Thompson // Clear v for transpose operation 285*db2becc9SJeremy L Thompson if (is_transpose && !apply_add) { 2861dda9c1aSJeremy L Thompson CeedSize length; 2871dda9c1aSJeremy L Thompson 2881dda9c1aSJeremy L Thompson CeedCallBackend(CeedVectorGetLength(v, &length)); 2891dda9c1aSJeremy L Thompson CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar))); 2901dda9c1aSJeremy L Thompson } 2911dda9c1aSJeremy L Thompson 2921dda9c1aSJeremy L Thompson // Basis action 2931dda9c1aSJeremy L Thompson switch (eval_mode) { 2941dda9c1aSJeremy L Thompson case CEED_EVAL_INTERP: { 2951dda9c1aSJeremy L Thompson void *interp_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_chebyshev_interp_1d, &d_x, &d_u, &d_v}; 2961dda9c1aSJeremy L Thompson const CeedInt block_size = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size); 2971dda9c1aSJeremy L Thompson 2981dda9c1aSJeremy L Thompson CeedCallBackend(CeedRunKernel_Cuda(ceed, data->InterpAtPoints, num_elem, block_size, interp_args)); 2991dda9c1aSJeremy L Thompson } break; 3001dda9c1aSJeremy L Thompson case CEED_EVAL_GRAD: { 3011dda9c1aSJeremy L Thompson void *grad_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_chebyshev_interp_1d, &d_x, &d_u, &d_v}; 3022d10e82cSJeremy L Thompson const CeedInt block_size = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size); 3031dda9c1aSJeremy L Thompson 3041dda9c1aSJeremy L Thompson CeedCallBackend(CeedRunKernel_Cuda(ceed, data->GradAtPoints, num_elem, block_size, grad_args)); 3051dda9c1aSJeremy L Thompson } break; 3061dda9c1aSJeremy L Thompson case CEED_EVAL_WEIGHT: 3071dda9c1aSJeremy L Thompson case CEED_EVAL_NONE: /* handled separately below */ 3081dda9c1aSJeremy L Thompson break; 3091dda9c1aSJeremy L Thompson // LCOV_EXCL_START 3101dda9c1aSJeremy L Thompson case CEED_EVAL_DIV: 3111dda9c1aSJeremy L Thompson case CEED_EVAL_CURL: 3121dda9c1aSJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 3131dda9c1aSJeremy L Thompson // LCOV_EXCL_STOP 3141dda9c1aSJeremy L Thompson } 3151dda9c1aSJeremy L Thompson 3161dda9c1aSJeremy L Thompson // Restore vectors, cover CEED_EVAL_NONE 3171dda9c1aSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArrayRead(x_ref, &d_x)); 3181dda9c1aSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 3191dda9c1aSJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u)); 3201dda9c1aSJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 3211dda9c1aSJeremy L Thompson return CEED_ERROR_SUCCESS; 3221dda9c1aSJeremy L Thompson } 3231dda9c1aSJeremy L Thompson 324*db2becc9SJeremy L Thompson static int CeedBasisApplyAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, 325*db2becc9SJeremy L Thompson CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) { 326*db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda_shared(basis, false, num_elem, num_points, t_mode, eval_mode, x_ref, u, v)); 327*db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 328*db2becc9SJeremy L Thompson } 329*db2becc9SJeremy L Thompson 330*db2becc9SJeremy L Thompson static int CeedBasisApplyAddAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, 331*db2becc9SJeremy L Thompson CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) { 332*db2becc9SJeremy L Thompson CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda_shared(basis, true, num_elem, num_points, t_mode, eval_mode, x_ref, u, v)); 333*db2becc9SJeremy L Thompson return CEED_ERROR_SUCCESS; 334*db2becc9SJeremy L Thompson } 335*db2becc9SJeremy L Thompson 3361dda9c1aSJeremy L Thompson //------------------------------------------------------------------------------ 337ab213215SJeremy L Thompson // Destroy basis 338ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 339c532df63SYohann static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) { 340c532df63SYohann Ceed ceed; 341c532df63SYohann CeedBasis_Cuda_shared *data; 342ca735530SJeremy L Thompson 343ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 3442b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &data)); 3452b730f8bSJeremy L Thompson CeedCallCuda(ceed, cuModuleUnload(data->module)); 3461dda9c1aSJeremy L Thompson if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints)); 347097cc795SJames Wright if (data->d_q_weight_1d) CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d)); 3482b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_interp_1d)); 3492b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_grad_1d)); 3502b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_collo_grad_1d)); 3511dda9c1aSJeremy L Thompson CeedCallCuda(ceed, cudaFree(data->d_chebyshev_interp_1d)); 3522b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&data)); 353e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 354c532df63SYohann } 355c532df63SYohann 356ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 357ab213215SJeremy L Thompson // Create tensor basis 358ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 3592b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d, 3602b730f8bSJeremy L Thompson const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) { 361c532df63SYohann Ceed ceed; 36222070f95SJeremy L Thompson char *basis_kernel_source; 36322070f95SJeremy L Thompson const char *basis_kernel_path; 364ca735530SJeremy L Thompson CeedInt num_comp; 365ca735530SJeremy L Thompson const CeedInt q_bytes = Q_1d * sizeof(CeedScalar); 366ca735530SJeremy L Thompson const CeedInt interp_bytes = q_bytes * P_1d; 367c532df63SYohann CeedBasis_Cuda_shared *data; 368ca735530SJeremy L Thompson 369ca735530SJeremy L Thompson CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 3702b730f8bSJeremy L Thompson CeedCallBackend(CeedCalloc(1, &data)); 371c532df63SYohann 372ab213215SJeremy L Thompson // Copy basis data to GPU 373097cc795SJames Wright if (q_weight_1d) { 3742b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes)); 3752b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice)); 376097cc795SJames Wright } 3772b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes)); 3782b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 3792b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes)); 3802b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice)); 381c532df63SYohann 382ab213215SJeremy L Thompson // Compute collocated gradient and copy to GPU 383437930d1SJeremy L Thompson data->d_collo_grad_1d = NULL; 3849e201c85SYohann bool has_collocated_grad = dim == 3 && Q_1d >= P_1d; 385ca735530SJeremy L Thompson 3869e201c85SYohann if (has_collocated_grad) { 387437930d1SJeremy L Thompson CeedScalar *collo_grad_1d; 388ca735530SJeremy L Thompson 3892b730f8bSJeremy L Thompson CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d)); 3902b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d)); 3912b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d)); 3922b730f8bSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, cudaMemcpyHostToDevice)); 3932b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&collo_grad_1d)); 394ac421f39SYohann } 395ac421f39SYohann 396ab213215SJeremy L Thompson // Compile basis kernels 3972b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 3982b730f8bSJeremy L Thompson CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-shared-basis-tensor.h", &basis_kernel_path)); 39923d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 4002b730f8bSJeremy L Thompson CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 40123d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete -----\n"); 402eb7e6cafSJeremy L Thompson CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D", 403eb7e6cafSJeremy L Thompson CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), 404eb7e6cafSJeremy L Thompson "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_HAS_COLLOCATED_GRAD", has_collocated_grad)); 405eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 406eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 407*db2becc9SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd)); 408eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad)); 409eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose)); 410*db2becc9SJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd)); 411eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 4122b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&basis_kernel_path)); 4132b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&basis_kernel_source)); 414c532df63SYohann 4152b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisSetData(basis, data)); 416ab213215SJeremy L Thompson 417ab213215SJeremy L Thompson // Register backend functions 4182b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared)); 419*db2becc9SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddTensor_Cuda_shared)); 4201dda9c1aSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAtPoints", CeedBasisApplyAtPoints_Cuda_shared)); 4212b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared)); 422e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 423c532df63SYohann } 4242a86cc9dSSebastian Grimberg 425ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 426