xref: /libCEED/backends/cuda-shared/ceed-cuda-shared-basis.c (revision f7c9815f6bdd2aac16e9c7ed574a2a61404669c3)
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 //------------------------------------------------------------------------------
302b730f8bSJeremy L Thompson int CeedBasisApplyTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
317f823360Sjeremylt                                      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");
482b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
49c532df63SYohann 
50ab213215SJeremy L Thompson   // Apply basis operation
51437930d1SJeremy L Thompson   switch (eval_mode) {
52ab213215SJeremy L Thompson     case CEED_EVAL_INTERP: {
53437930d1SJeremy L Thompson       CeedInt P_1d, Q_1d;
54ca735530SJeremy L Thompson 
552b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
562b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
57437930d1SJeremy L Thompson       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
58ca735530SJeremy L Thompson 
59eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedInit_CudaInterp(data->d_interp_1d, P_1d, Q_1d, &data->c_B));
602b730f8bSJeremy L Thompson       void *interp_args[] = {(void *)&num_elem, &data->c_B, &d_u, &d_v};
61ca735530SJeremy L Thompson 
624d537eeaSYohann       if (dim == 1) {
632b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d,
64e6f67ff7SJed Brown                                                                                                  1));  // avoid >512 total threads
652b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
66437930d1SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
67b2165e7aSSebastian Grimberg 
689e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
69eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->InterpTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
709e201c85SYohann         } else {
71eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
729e201c85SYohann         }
73074be161SYohann Dudouit       } else if (dim == 2) {
74437930d1SJeremy L Thompson         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
75437930d1SJeremy L Thompson         // elems_per_block must be at least 1
762b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
772b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
782b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
79b2165e7aSSebastian Grimberg 
809e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
812b730f8bSJeremy L Thompson           CeedCallBackend(
82eb7e6cafSJeremy L Thompson               CeedRunKernelDimShared_Cuda(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
839e201c85SYohann         } else {
84eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
859e201c85SYohann         }
86074be161SYohann Dudouit       } else if (dim == 3) {
87437930d1SJeremy L Thompson         CeedInt elems_per_block = 1;
882b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
892b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
90b2165e7aSSebastian Grimberg 
919e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
922b730f8bSJeremy L Thompson           CeedCallBackend(
93eb7e6cafSJeremy L Thompson               CeedRunKernelDimShared_Cuda(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
949e201c85SYohann         } else {
95eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
96074be161SYohann Dudouit         }
979e201c85SYohann       }
98ab213215SJeremy L Thompson     } break;
99ab213215SJeremy L Thompson     case CEED_EVAL_GRAD: {
100437930d1SJeremy L Thompson       CeedInt P_1d, Q_1d;
101ca735530SJeremy L Thompson 
1022b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
1032b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
104437930d1SJeremy L Thompson       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
105ca735530SJeremy L Thompson 
1069e201c85SYohann       if (data->d_collo_grad_1d) {
107eb7e6cafSJeremy L Thompson         CeedCallBackend(CeedInit_CudaCollocatedGrad(data->d_interp_1d, data->d_collo_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G));
1089e201c85SYohann       } else {
109eb7e6cafSJeremy L Thompson         CeedCallBackend(CeedInit_CudaGrad(data->d_interp_1d, data->d_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G));
1109e201c85SYohann       }
1112b730f8bSJeremy L Thompson       void *grad_args[] = {(void *)&num_elem, &data->c_B, &data->c_G, &d_u, &d_v};
1124d537eeaSYohann       if (dim == 1) {
1132b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d,
114e6f67ff7SJed Brown                                                                                                  1));  // avoid >512 total threads
1152b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
116437930d1SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
117b2165e7aSSebastian Grimberg 
1189e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
119eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->GradTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
1209e201c85SYohann         } else {
121eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
1229e201c85SYohann         }
123074be161SYohann Dudouit       } else if (dim == 2) {
124437930d1SJeremy L Thompson         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
125437930d1SJeremy L Thompson         // elems_per_block must be at least 1
1262b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
1272b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
1282b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
129b2165e7aSSebastian Grimberg 
1309e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
131eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
1329e201c85SYohann         } else {
133eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
1349e201c85SYohann         }
135074be161SYohann Dudouit       } else if (dim == 3) {
136437930d1SJeremy L Thompson         CeedInt elems_per_block = 1;
1372b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
1382b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
139b2165e7aSSebastian Grimberg 
1409e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
141eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->GradTranspose, grid, thread_1d, thread_1d, 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;
156b2165e7aSSebastian Grimberg         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 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;
162b2165e7aSSebastian Grimberg         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 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;
168b2165e7aSSebastian Grimberg         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 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));
186e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
187c532df63SYohann }
188c532df63SYohann 
189ab213215SJeremy L Thompson //------------------------------------------------------------------------------
1901dda9c1aSJeremy L Thompson // Basis apply - tensor AtPoints
1911dda9c1aSJeremy L Thompson //------------------------------------------------------------------------------
1921dda9c1aSJeremy L Thompson int CeedBasisApplyAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
1931dda9c1aSJeremy L Thompson                                        CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
1941dda9c1aSJeremy L Thompson   Ceed                   ceed;
1951dda9c1aSJeremy L Thompson   CeedInt                Q_1d, dim, max_num_points = num_points[0];
1961dda9c1aSJeremy L Thompson   const CeedInt          is_transpose   = t_mode == CEED_TRANSPOSE;
1971dda9c1aSJeremy L Thompson   const int              max_block_size = 32;
1981dda9c1aSJeremy L Thompson   const CeedScalar      *d_x, *d_u;
1991dda9c1aSJeremy L Thompson   CeedScalar            *d_v;
2001dda9c1aSJeremy L Thompson   CeedBasis_Cuda_shared *data;
2011dda9c1aSJeremy L Thompson 
2021dda9c1aSJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
2031dda9c1aSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
2041dda9c1aSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
2051dda9c1aSJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
2061dda9c1aSJeremy L Thompson 
2071dda9c1aSJeremy L Thompson   // Check uniform number of points per elem
2081dda9c1aSJeremy L Thompson   for (CeedInt i = 1; i < num_elem; i++) {
2091dda9c1aSJeremy L Thompson     CeedCheck(max_num_points == num_points[i], ceed, CEED_ERROR_BACKEND,
2101dda9c1aSJeremy L Thompson               "BasisApplyAtPoints only supported for the same number of points in each element");
2111dda9c1aSJeremy L Thompson   }
2121dda9c1aSJeremy L Thompson 
2131dda9c1aSJeremy L Thompson   // Weight handled separately
2141dda9c1aSJeremy L Thompson   if (eval_mode == CEED_EVAL_WEIGHT) {
2151dda9c1aSJeremy L Thompson     CeedCall(CeedVectorSetValue(v, 1.0));
2161dda9c1aSJeremy L Thompson     return CEED_ERROR_SUCCESS;
2171dda9c1aSJeremy L Thompson   }
2181dda9c1aSJeremy L Thompson 
2191dda9c1aSJeremy L Thompson   // Build kernels if needed
2201dda9c1aSJeremy L Thompson   if (data->num_points != max_num_points) {
2211dda9c1aSJeremy L Thompson     CeedInt P_1d;
2221dda9c1aSJeremy L Thompson 
2231dda9c1aSJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
2241dda9c1aSJeremy L Thompson     data->num_points = max_num_points;
2251dda9c1aSJeremy L Thompson 
2261dda9c1aSJeremy L Thompson     // -- Create interp matrix to Chebyshev coefficients
2271dda9c1aSJeremy L Thompson     if (!data->d_chebyshev_interp_1d) {
2281dda9c1aSJeremy L Thompson       CeedSize    interp_bytes;
2291dda9c1aSJeremy L Thompson       CeedScalar *chebyshev_interp_1d;
2301dda9c1aSJeremy L Thompson 
2311dda9c1aSJeremy L Thompson       interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
2321dda9c1aSJeremy L Thompson       CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
2331dda9c1aSJeremy L Thompson       CeedCall(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
2341dda9c1aSJeremy L Thompson       CeedCallCuda(ceed, cudaMalloc((void **)&data->d_chebyshev_interp_1d, interp_bytes));
2351dda9c1aSJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy(data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
2361dda9c1aSJeremy L Thompson       CeedCallBackend(CeedFree(&chebyshev_interp_1d));
2371dda9c1aSJeremy L Thompson     }
2381dda9c1aSJeremy L Thompson 
2391dda9c1aSJeremy L Thompson     // -- Compile kernels
2401dda9c1aSJeremy L Thompson     char       *basis_kernel_source;
2411dda9c1aSJeremy L Thompson     const char *basis_kernel_path;
2421dda9c1aSJeremy L Thompson     CeedInt     num_comp;
2431dda9c1aSJeremy L Thompson 
2441dda9c1aSJeremy L Thompson     if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
2451dda9c1aSJeremy L Thompson     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
2461dda9c1aSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-tensor-at-points.h", &basis_kernel_path));
2471dda9c1aSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
2481dda9c1aSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
2491dda9c1aSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
2501dda9c1aSJeremy 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",
251*f7c9815fSJeremy L Thompson                                      Q_1d * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim - 1), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp,
2521dda9c1aSJeremy L Thompson                                      "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_NUM_PTS",
253*f7c9815fSJeremy L Thompson                                      max_num_points, "POINTS_BUFF_LEN", CeedIntPow(Q_1d, dim - 1)));
2541dda9c1aSJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpAtPoints", &data->InterpAtPoints));
2551dda9c1aSJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradAtPoints", &data->GradAtPoints));
2561dda9c1aSJeremy L Thompson     CeedCallBackend(CeedFree(&basis_kernel_path));
2571dda9c1aSJeremy L Thompson     CeedCallBackend(CeedFree(&basis_kernel_source));
2581dda9c1aSJeremy L Thompson   }
2591dda9c1aSJeremy L Thompson 
2601dda9c1aSJeremy L Thompson   // Get read/write access to u, v
2611dda9c1aSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayRead(x_ref, CEED_MEM_DEVICE, &d_x));
2621dda9c1aSJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
2631dda9c1aSJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
2641dda9c1aSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
2651dda9c1aSJeremy L Thompson 
2661dda9c1aSJeremy L Thompson   // Clear v for transpose operation
2671dda9c1aSJeremy L Thompson   if (is_transpose) {
2681dda9c1aSJeremy L Thompson     CeedSize length;
2691dda9c1aSJeremy L Thompson 
2701dda9c1aSJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(v, &length));
2711dda9c1aSJeremy L Thompson     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
2721dda9c1aSJeremy L Thompson   }
2731dda9c1aSJeremy L Thompson 
2741dda9c1aSJeremy L Thompson   // Basis action
2751dda9c1aSJeremy L Thompson   switch (eval_mode) {
2761dda9c1aSJeremy L Thompson     case CEED_EVAL_INTERP: {
2771dda9c1aSJeremy L Thompson       void         *interp_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_chebyshev_interp_1d, &d_x, &d_u, &d_v};
2781dda9c1aSJeremy L Thompson       const CeedInt block_size    = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
2791dda9c1aSJeremy L Thompson 
2801dda9c1aSJeremy L Thompson       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->InterpAtPoints, num_elem, block_size, interp_args));
2811dda9c1aSJeremy L Thompson     } break;
2821dda9c1aSJeremy L Thompson     case CEED_EVAL_GRAD: {
2831dda9c1aSJeremy L Thompson       void         *grad_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_chebyshev_interp_1d, &d_x, &d_u, &d_v};
2842d10e82cSJeremy L Thompson       const CeedInt block_size  = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
2851dda9c1aSJeremy L Thompson 
2861dda9c1aSJeremy L Thompson       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->GradAtPoints, num_elem, block_size, grad_args));
2871dda9c1aSJeremy L Thompson     } break;
2881dda9c1aSJeremy L Thompson     case CEED_EVAL_WEIGHT:
2891dda9c1aSJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
2901dda9c1aSJeremy L Thompson       break;
2911dda9c1aSJeremy L Thompson     // LCOV_EXCL_START
2921dda9c1aSJeremy L Thompson     case CEED_EVAL_DIV:
2931dda9c1aSJeremy L Thompson     case CEED_EVAL_CURL:
2941dda9c1aSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
2951dda9c1aSJeremy L Thompson       // LCOV_EXCL_STOP
2961dda9c1aSJeremy L Thompson   }
2971dda9c1aSJeremy L Thompson 
2981dda9c1aSJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
2991dda9c1aSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArrayRead(x_ref, &d_x));
3001dda9c1aSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
3011dda9c1aSJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
3021dda9c1aSJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
3031dda9c1aSJeremy L Thompson   return CEED_ERROR_SUCCESS;
3041dda9c1aSJeremy L Thompson }
3051dda9c1aSJeremy L Thompson 
3061dda9c1aSJeremy L Thompson //------------------------------------------------------------------------------
307ab213215SJeremy L Thompson // Destroy basis
308ab213215SJeremy L Thompson //------------------------------------------------------------------------------
309c532df63SYohann static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) {
310c532df63SYohann   Ceed                   ceed;
311c532df63SYohann   CeedBasis_Cuda_shared *data;
312ca735530SJeremy L Thompson 
313ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
3142b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
3152b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
3161dda9c1aSJeremy L Thompson   if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
317097cc795SJames Wright   if (data->d_q_weight_1d) CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d));
3182b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp_1d));
3192b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad_1d));
3202b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_collo_grad_1d));
3211dda9c1aSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_chebyshev_interp_1d));
3222b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
323e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
324c532df63SYohann }
325c532df63SYohann 
326ab213215SJeremy L Thompson //------------------------------------------------------------------------------
327ab213215SJeremy L Thompson // Create tensor basis
328ab213215SJeremy L Thompson //------------------------------------------------------------------------------
3292b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
3302b730f8bSJeremy L Thompson                                         const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
331c532df63SYohann   Ceed                   ceed;
33222070f95SJeremy L Thompson   char                  *basis_kernel_source;
33322070f95SJeremy L Thompson   const char            *basis_kernel_path;
334ca735530SJeremy L Thompson   CeedInt                num_comp;
335ca735530SJeremy L Thompson   const CeedInt          q_bytes      = Q_1d * sizeof(CeedScalar);
336ca735530SJeremy L Thompson   const CeedInt          interp_bytes = q_bytes * P_1d;
337c532df63SYohann   CeedBasis_Cuda_shared *data;
338ca735530SJeremy L Thompson 
339ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
3402b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
341c532df63SYohann 
342ab213215SJeremy L Thompson   // Copy basis data to GPU
343097cc795SJames Wright   if (q_weight_1d) {
3442b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
3452b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice));
346097cc795SJames Wright   }
3472b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
3482b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice));
3492b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes));
3502b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice));
351c532df63SYohann 
352ab213215SJeremy L Thompson   // Compute collocated gradient and copy to GPU
353437930d1SJeremy L Thompson   data->d_collo_grad_1d    = NULL;
3549e201c85SYohann   bool has_collocated_grad = dim == 3 && Q_1d >= P_1d;
355ca735530SJeremy L Thompson 
3569e201c85SYohann   if (has_collocated_grad) {
357437930d1SJeremy L Thompson     CeedScalar *collo_grad_1d;
358ca735530SJeremy L Thompson 
3592b730f8bSJeremy L Thompson     CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d));
3602b730f8bSJeremy L Thompson     CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d));
3612b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d));
3622b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, cudaMemcpyHostToDevice));
3632b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&collo_grad_1d));
364ac421f39SYohann   }
365ac421f39SYohann 
366ab213215SJeremy L Thompson   // Compile basis kernels
3672b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
3682b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-shared-basis-tensor.h", &basis_kernel_path));
36923d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
3702b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
37123d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete -----\n");
372eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D",
373eb7e6cafSJeremy L Thompson                                    CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim),
374eb7e6cafSJeremy L Thompson                                    "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_HAS_COLLOCATED_GRAD", has_collocated_grad));
375eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
376eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
377eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
378eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose));
379eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
3802b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_path));
3812b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_source));
382c532df63SYohann 
3832b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
384ab213215SJeremy L Thompson 
385ab213215SJeremy L Thompson   // Register backend functions
3862b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared));
3871dda9c1aSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAtPoints", CeedBasisApplyAtPoints_Cuda_shared));
3882b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared));
389e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
390c532df63SYohann }
3912a86cc9dSSebastian Grimberg 
392ab213215SJeremy L Thompson //------------------------------------------------------------------------------
393