xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-shared/ceed-cuda-shared-basis.c (revision 9ff05d55386b4e6413be60b7231511258906fd9f)
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 //------------------------------------------------------------------------------
22ab213215SJeremy L Thompson // Device initalization
23ab213215SJeremy L Thompson //------------------------------------------------------------------------------
24eb7e6cafSJeremy L Thompson int CeedInit_CudaInterp(CeedScalar *d_B, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B);
25eb7e6cafSJeremy 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);
26eb7e6cafSJeremy 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);
27c532df63SYohann 
28ab213215SJeremy L Thompson //------------------------------------------------------------------------------
29*9ff05d55SJeremy L Thompson // Apply tensor basis
30ab213215SJeremy L Thompson //------------------------------------------------------------------------------
31db2becc9SJeremy L Thompson static int CeedBasisApplyTensorCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode,
32db2becc9SJeremy L Thompson                                                 CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
33c532df63SYohann   Ceed                   ceed;
346dbfb411Snbeams   Ceed_Cuda             *ceed_Cuda;
35437930d1SJeremy L Thompson   CeedInt                dim, num_comp;
36ca735530SJeremy L Thompson   const CeedScalar      *d_u;
37ca735530SJeremy L Thompson   CeedScalar            *d_v;
38ca735530SJeremy L Thompson   CeedBasis_Cuda_shared *data;
39ca735530SJeremy L Thompson 
40ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
41ca735530SJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
42ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
432b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
442b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
45c532df63SYohann 
469ea2cfd9SJeremy L Thompson   // Get read/write access to u, v
476574a04fSJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
486574a04fSJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
49db2becc9SJeremy L Thompson   if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
50db2becc9SJeremy L Thompson   else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
51c532df63SYohann 
52ab213215SJeremy L Thompson   // Apply basis operation
53437930d1SJeremy L Thompson   switch (eval_mode) {
54ab213215SJeremy L Thompson     case CEED_EVAL_INTERP: {
55437930d1SJeremy L Thompson       CeedInt P_1d, Q_1d;
56ca735530SJeremy L Thompson 
572b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
582b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
59437930d1SJeremy L Thompson       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
60ca735530SJeremy L Thompson 
61eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedInit_CudaInterp(data->d_interp_1d, P_1d, Q_1d, &data->c_B));
622b730f8bSJeremy L Thompson       void *interp_args[] = {(void *)&num_elem, &data->c_B, &d_u, &d_v};
63ca735530SJeremy L Thompson 
644d537eeaSYohann       if (dim == 1) {
65dac82721SJeremy L Thompson         // avoid >512 total threads
66dac82721SJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1));
67a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
68437930d1SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
69b2165e7aSSebastian Grimberg 
709e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
71db2becc9SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, 1,
72db2becc9SJeremy L Thompson                                                       elems_per_block, shared_mem, interp_args));
739e201c85SYohann         } else {
74eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
759e201c85SYohann         }
76074be161SYohann Dudouit       } else if (dim == 2) {
77437930d1SJeremy L Thompson         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
78437930d1SJeremy L Thompson         // elems_per_block must be at least 1
792b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
80a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
812b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
82b2165e7aSSebastian Grimberg 
839e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
84db2becc9SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, thread_1d,
85db2becc9SJeremy L Thompson                                                       elems_per_block, shared_mem, interp_args));
869e201c85SYohann         } else {
87eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
889e201c85SYohann         }
89074be161SYohann Dudouit       } else if (dim == 3) {
90437930d1SJeremy L Thompson         CeedInt elems_per_block = 1;
91a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
922b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
93b2165e7aSSebastian Grimberg 
949e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
95db2becc9SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, thread_1d,
96db2becc9SJeremy L Thompson                                                       elems_per_block, shared_mem, interp_args));
979e201c85SYohann         } else {
98eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
99074be161SYohann Dudouit         }
1009e201c85SYohann       }
101ab213215SJeremy L Thompson     } break;
102ab213215SJeremy L Thompson     case CEED_EVAL_GRAD: {
103437930d1SJeremy L Thompson       CeedInt P_1d, Q_1d;
104ca735530SJeremy L Thompson 
1052b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
1062b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
107437930d1SJeremy L Thompson       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
108ca735530SJeremy L Thompson 
1099e201c85SYohann       if (data->d_collo_grad_1d) {
110eb7e6cafSJeremy L Thompson         CeedCallBackend(CeedInit_CudaCollocatedGrad(data->d_interp_1d, data->d_collo_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G));
1119e201c85SYohann       } else {
112eb7e6cafSJeremy L Thompson         CeedCallBackend(CeedInit_CudaGrad(data->d_interp_1d, data->d_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G));
1139e201c85SYohann       }
1142b730f8bSJeremy L Thompson       void *grad_args[] = {(void *)&num_elem, &data->c_B, &data->c_G, &d_u, &d_v};
1154d537eeaSYohann       if (dim == 1) {
116dac82721SJeremy L Thompson         // avoid >512 total threads
117dac82721SJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1));
118a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
119437930d1SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
120b2165e7aSSebastian Grimberg 
1219e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
122db2becc9SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, 1,
123db2becc9SJeremy L Thompson                                                       elems_per_block, shared_mem, grad_args));
1249e201c85SYohann         } else {
125eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
1269e201c85SYohann         }
127074be161SYohann Dudouit       } else if (dim == 2) {
128437930d1SJeremy L Thompson         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
129437930d1SJeremy L Thompson         // elems_per_block must be at least 1
1302b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
131a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
1322b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
133b2165e7aSSebastian Grimberg 
1349e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
135db2becc9SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, thread_1d,
136db2becc9SJeremy L Thompson                                                       elems_per_block, shared_mem, grad_args));
1379e201c85SYohann         } else {
138eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
1399e201c85SYohann         }
140074be161SYohann Dudouit       } else if (dim == 3) {
141437930d1SJeremy L Thompson         CeedInt elems_per_block = 1;
142a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
1432b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
144b2165e7aSSebastian Grimberg 
1459e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
146db2becc9SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, thread_1d,
147db2becc9SJeremy L Thompson                                                       elems_per_block, shared_mem, grad_args));
1489e201c85SYohann         } else {
149eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
1509e201c85SYohann         }
151074be161SYohann Dudouit       }
152ab213215SJeremy L Thompson     } break;
153ab213215SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
154437930d1SJeremy L Thompson       CeedInt Q_1d;
155397164e9SSebastian Grimberg       CeedInt block_size = 32;
156ca735530SJeremy L Thompson 
157097cc795SJames Wright       CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights_1d not set", CeedEvalModes[eval_mode]);
1582b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
159437930d1SJeremy L Thompson       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
160074be161SYohann Dudouit       if (dim == 1) {
161397164e9SSebastian Grimberg         const CeedInt elems_per_block = block_size / Q_1d;
162a8d440fbSJeremy L Thompson         const CeedInt grid_size       = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
163b2165e7aSSebastian Grimberg 
164b2165e7aSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, elems_per_block, 1, weight_args));
165074be161SYohann Dudouit       } else if (dim == 2) {
166397164e9SSebastian Grimberg         const CeedInt opt_elems       = block_size / (Q_1d * Q_1d);
167437930d1SJeremy L Thompson         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
168a8d440fbSJeremy L Thompson         const CeedInt grid_size       = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
169b2165e7aSSebastian Grimberg 
170b2165e7aSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
171074be161SYohann Dudouit       } else if (dim == 3) {
172397164e9SSebastian Grimberg         const CeedInt opt_elems       = block_size / (Q_1d * Q_1d);
1739e201c85SYohann         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
174a8d440fbSJeremy L Thompson         const CeedInt grid_size       = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
175b2165e7aSSebastian Grimberg 
176b2165e7aSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
177074be161SYohann Dudouit       }
178ab213215SJeremy L Thompson     } break;
1799ea2cfd9SJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
1809ea2cfd9SJeremy L Thompson       break;
181ab213215SJeremy L Thompson     // LCOV_EXCL_START
182ab213215SJeremy L Thompson     case CEED_EVAL_DIV:
183ab213215SJeremy L Thompson     case CEED_EVAL_CURL:
184bcbe1c99SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
185ab213215SJeremy L Thompson       // LCOV_EXCL_STOP
186c532df63SYohann   }
187c532df63SYohann 
1889ea2cfd9SJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
1892b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
1909ea2cfd9SJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
1919ea2cfd9SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
1929bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
193e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
194c532df63SYohann }
195c532df63SYohann 
196db2becc9SJeremy L Thompson static int CeedBasisApplyTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
197db2becc9SJeremy L Thompson                                             CeedVector v) {
198db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyTensorCore_Cuda_shared(basis, false, num_elem, t_mode, eval_mode, u, v));
199db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
200db2becc9SJeremy L Thompson }
201db2becc9SJeremy L Thompson 
202db2becc9SJeremy L Thompson static int CeedBasisApplyAddTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
203db2becc9SJeremy L Thompson                                                CeedVector u, CeedVector v) {
204db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyTensorCore_Cuda_shared(basis, true, num_elem, t_mode, eval_mode, u, v));
205db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
206db2becc9SJeremy L Thompson }
207db2becc9SJeremy L Thompson 
208ab213215SJeremy L Thompson //------------------------------------------------------------------------------
2091dda9c1aSJeremy L Thompson // Basis apply - tensor AtPoints
2101dda9c1aSJeremy L Thompson //------------------------------------------------------------------------------
211db2becc9SJeremy L Thompson static int CeedBasisApplyAtPointsCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, const CeedInt *num_points,
212db2becc9SJeremy L Thompson                                                   CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
2131dda9c1aSJeremy L Thompson   Ceed                   ceed;
2149e1d4b82SJeremy L Thompson   Ceed_Cuda             *ceed_Cuda;
2159e1d4b82SJeremy L Thompson   CeedInt                Q_1d, dim, num_comp, max_num_points = num_points[0];
2161dda9c1aSJeremy L Thompson   const CeedInt          is_transpose = t_mode == CEED_TRANSPOSE;
2171dda9c1aSJeremy L Thompson   const CeedScalar      *d_x, *d_u;
2181dda9c1aSJeremy L Thompson   CeedScalar            *d_v;
2191dda9c1aSJeremy L Thompson   CeedBasis_Cuda_shared *data;
2201dda9c1aSJeremy L Thompson 
2211dda9c1aSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
2221dda9c1aSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
2231dda9c1aSJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
2249e1d4b82SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
2251dda9c1aSJeremy L Thompson 
2261dda9c1aSJeremy L Thompson   // Weight handled separately
2271dda9c1aSJeremy L Thompson   if (eval_mode == CEED_EVAL_WEIGHT) {
2285a5594ffSJeremy L Thompson     CeedCallBackend(CeedVectorSetValue(v, 1.0));
2291dda9c1aSJeremy L Thompson     return CEED_ERROR_SUCCESS;
2301dda9c1aSJeremy L Thompson   }
2311dda9c1aSJeremy L Thompson 
2329bc66399SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
2339e1d4b82SJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
2349bc66399SJeremy L Thompson 
235111870feSJeremy L Thompson   // Check padded to uniform number of points per elem
236111870feSJeremy L Thompson   for (CeedInt i = 1; i < num_elem; i++) max_num_points = CeedIntMax(max_num_points, num_points[i]);
237111870feSJeremy L Thompson   {
2389e1d4b82SJeremy L Thompson     CeedInt  q_comp;
239111870feSJeremy L Thompson     CeedSize len, len_required;
240111870feSJeremy L Thompson     CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
241111870feSJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(is_transpose ? u : v, &len));
242111870feSJeremy L Thompson     len_required = (CeedSize)num_comp * (CeedSize)q_comp * (CeedSize)num_elem * (CeedSize)max_num_points;
243111870feSJeremy L Thompson     CeedCheck(len >= len_required, ceed, CEED_ERROR_BACKEND,
244111870feSJeremy L Thompson               "Vector at points must be padded to the same number of points in each element for BasisApplyAtPoints on GPU backends."
245111870feSJeremy L Thompson               " Found %" CeedSize_FMT ", Required %" CeedSize_FMT,
246111870feSJeremy L Thompson               len, len_required);
247111870feSJeremy L Thompson   }
248111870feSJeremy L Thompson 
249111870feSJeremy L Thompson   // Move num_points array to device
250111870feSJeremy L Thompson   if (is_transpose) {
251111870feSJeremy L Thompson     const CeedInt num_bytes = num_elem * sizeof(CeedInt);
252111870feSJeremy L Thompson 
253111870feSJeremy L Thompson     if (num_elem != data->num_elem_at_points) {
254111870feSJeremy L Thompson       data->num_elem_at_points = num_elem;
255111870feSJeremy L Thompson 
256111870feSJeremy L Thompson       if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem));
257111870feSJeremy L Thompson       CeedCallCuda(ceed, cudaMalloc((void **)&data->d_points_per_elem, num_bytes));
258111870feSJeremy L Thompson       CeedCallBackend(CeedFree(&data->h_points_per_elem));
259111870feSJeremy L Thompson       CeedCallBackend(CeedCalloc(num_elem, &data->h_points_per_elem));
260111870feSJeremy L Thompson     }
2619e511c80SJeremy L Thompson     if (memcmp(data->h_points_per_elem, num_points, num_bytes)) {
262111870feSJeremy L Thompson       memcpy(data->h_points_per_elem, num_points, num_bytes);
263111870feSJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy(data->d_points_per_elem, num_points, num_bytes, cudaMemcpyHostToDevice));
264111870feSJeremy L Thompson     }
265111870feSJeremy L Thompson   }
266111870feSJeremy L Thompson 
2671dda9c1aSJeremy L Thompson   // Build kernels if needed
2681dda9c1aSJeremy L Thompson   if (data->num_points != max_num_points) {
2691dda9c1aSJeremy L Thompson     CeedInt P_1d;
2701dda9c1aSJeremy L Thompson 
2711dda9c1aSJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
2721dda9c1aSJeremy L Thompson     data->num_points = max_num_points;
2731dda9c1aSJeremy L Thompson 
2741dda9c1aSJeremy L Thompson     // -- Create interp matrix to Chebyshev coefficients
2751dda9c1aSJeremy L Thompson     if (!data->d_chebyshev_interp_1d) {
2761dda9c1aSJeremy L Thompson       CeedSize    interp_bytes;
2771dda9c1aSJeremy L Thompson       CeedScalar *chebyshev_interp_1d;
2781dda9c1aSJeremy L Thompson 
2791dda9c1aSJeremy L Thompson       interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
2801dda9c1aSJeremy L Thompson       CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
2815a5594ffSJeremy L Thompson       CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
2821dda9c1aSJeremy L Thompson       CeedCallCuda(ceed, cudaMalloc((void **)&data->d_chebyshev_interp_1d, interp_bytes));
2831dda9c1aSJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy(data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
2841dda9c1aSJeremy L Thompson       CeedCallBackend(CeedFree(&chebyshev_interp_1d));
2851dda9c1aSJeremy L Thompson     }
2861dda9c1aSJeremy L Thompson 
2871dda9c1aSJeremy L Thompson     // -- Compile kernels
2889e1d4b82SJeremy L Thompson     const char basis_kernel_source[] = "// AtPoints basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-at-points.h>\n";
2891dda9c1aSJeremy L Thompson     CeedInt    num_comp;
2901dda9c1aSJeremy L Thompson 
2911dda9c1aSJeremy L Thompson     if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
2921dda9c1aSJeremy L Thompson     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
2939e1d4b82SJeremy L Thompson     CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->moduleAtPoints, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D",
2949e1d4b82SJeremy L Thompson                                      CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim),
2959e1d4b82SJeremy L Thompson                                      "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_NUM_PTS", max_num_points));
2961dda9c1aSJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpAtPoints", &data->InterpAtPoints));
29781ae6159SJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpTransposeAtPoints", &data->InterpTransposeAtPoints));
2981dda9c1aSJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradAtPoints", &data->GradAtPoints));
29981ae6159SJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradTransposeAtPoints", &data->GradTransposeAtPoints));
3001dda9c1aSJeremy L Thompson   }
3011dda9c1aSJeremy L Thompson 
3021dda9c1aSJeremy L Thompson   // Get read/write access to u, v
3031dda9c1aSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayRead(x_ref, CEED_MEM_DEVICE, &d_x));
3041dda9c1aSJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
3051dda9c1aSJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
306db2becc9SJeremy L Thompson   if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
307db2becc9SJeremy L Thompson   else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
3081dda9c1aSJeremy L Thompson 
3091dda9c1aSJeremy L Thompson   // Clear v for transpose operation
310db2becc9SJeremy L Thompson   if (is_transpose && !apply_add) {
31119a04db8SJeremy L Thompson     CeedInt  num_comp, q_comp, num_nodes;
3121dda9c1aSJeremy L Thompson     CeedSize length;
3131dda9c1aSJeremy L Thompson 
31419a04db8SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
31519a04db8SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
31619a04db8SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes));
31719a04db8SJeremy L Thompson     length =
31819a04db8SJeremy L Thompson         (CeedSize)num_elem * (CeedSize)num_comp * (t_mode == CEED_TRANSPOSE ? (CeedSize)num_nodes : ((CeedSize)max_num_points * (CeedSize)q_comp));
3191dda9c1aSJeremy L Thompson     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
3201dda9c1aSJeremy L Thompson   }
3211dda9c1aSJeremy L Thompson 
3221dda9c1aSJeremy L Thompson   // Basis action
3231dda9c1aSJeremy L Thompson   switch (eval_mode) {
3241dda9c1aSJeremy L Thompson     case CEED_EVAL_INTERP: {
3259e1d4b82SJeremy L Thompson       CeedInt P_1d, Q_1d;
3261dda9c1aSJeremy L Thompson 
3279e1d4b82SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
3289e1d4b82SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
3299e1d4b82SJeremy L Thompson       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
3309e1d4b82SJeremy L Thompson 
3319e1d4b82SJeremy L Thompson       CeedCallBackend(CeedInit_CudaInterp(data->d_chebyshev_interp_1d, P_1d, Q_1d, &data->c_B));
3329e1d4b82SJeremy L Thompson       void *interp_args[] = {(void *)&num_elem, &data->c_B, &data->d_points_per_elem, &d_x, &d_u, &d_v};
3339e1d4b82SJeremy L Thompson 
3349e1d4b82SJeremy L Thompson       if (dim == 1) {
335dac82721SJeremy L Thompson         // avoid >512 total threads
336dac82721SJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 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 * sizeof(CeedScalar);
3399e1d4b82SJeremy L Thompson 
3409e1d4b82SJeremy L Thompson         CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d, 1,
3419e1d4b82SJeremy L Thompson                                                     elems_per_block, shared_mem, interp_args));
3429e1d4b82SJeremy L Thompson       } else if (dim == 2) {
3439e1d4b82SJeremy L Thompson         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
3449e1d4b82SJeremy L Thompson         // elems_per_block must be at least 1
3459e1d4b82SJeremy L Thompson         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
346a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
3479e1d4b82SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
3489e1d4b82SJeremy L Thompson 
3499e1d4b82SJeremy L Thompson         CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d,
3509e1d4b82SJeremy L Thompson                                                     thread_1d, elems_per_block, shared_mem, interp_args));
3519e1d4b82SJeremy L Thompson       } else if (dim == 3) {
3529e1d4b82SJeremy L Thompson         CeedInt elems_per_block = 1;
353a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
3549e1d4b82SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
3559e1d4b82SJeremy L Thompson 
3569e1d4b82SJeremy L Thompson         CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d,
3579e1d4b82SJeremy L Thompson                                                     thread_1d, elems_per_block, shared_mem, interp_args));
3589e1d4b82SJeremy L Thompson       }
3591dda9c1aSJeremy L Thompson     } break;
3601dda9c1aSJeremy L Thompson     case CEED_EVAL_GRAD: {
3619e1d4b82SJeremy L Thompson       CeedInt P_1d, Q_1d;
3621dda9c1aSJeremy L Thompson 
3639e1d4b82SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
3649e1d4b82SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
3659e1d4b82SJeremy L Thompson       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
3669e1d4b82SJeremy L Thompson 
3679e1d4b82SJeremy L Thompson       CeedCallBackend(CeedInit_CudaInterp(data->d_chebyshev_interp_1d, P_1d, Q_1d, &data->c_B));
3689e1d4b82SJeremy L Thompson       void *grad_args[] = {(void *)&num_elem, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v};
3699e1d4b82SJeremy L Thompson 
3709e1d4b82SJeremy L Thompson       if (dim == 1) {
371dac82721SJeremy L Thompson         // avoid >512 total threads
372dac82721SJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1));
373a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
3749e1d4b82SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
3759e1d4b82SJeremy L Thompson 
3769e1d4b82SJeremy L Thompson         CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, 1,
3779e1d4b82SJeremy L Thompson                                                     elems_per_block, shared_mem, grad_args));
3789e1d4b82SJeremy L Thompson       } else if (dim == 2) {
3799e1d4b82SJeremy L Thompson         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
3809e1d4b82SJeremy L Thompson         // elems_per_block must be at least 1
3819e1d4b82SJeremy L Thompson         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
382a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
3839e1d4b82SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
3849e1d4b82SJeremy L Thompson 
3859e1d4b82SJeremy L Thompson         CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, thread_1d,
3869e1d4b82SJeremy L Thompson                                                     elems_per_block, shared_mem, grad_args));
3879e1d4b82SJeremy L Thompson       } else if (dim == 3) {
3889e1d4b82SJeremy L Thompson         CeedInt elems_per_block = 1;
389a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
3909e1d4b82SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
3919e1d4b82SJeremy L Thompson 
3929e1d4b82SJeremy L Thompson         CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, thread_1d,
3939e1d4b82SJeremy L Thompson                                                     elems_per_block, shared_mem, grad_args));
3949e1d4b82SJeremy L Thompson       }
3951dda9c1aSJeremy L Thompson     } break;
3961dda9c1aSJeremy L Thompson     case CEED_EVAL_WEIGHT:
3971dda9c1aSJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
3981dda9c1aSJeremy L Thompson       break;
3991dda9c1aSJeremy L Thompson     // LCOV_EXCL_START
4001dda9c1aSJeremy L Thompson     case CEED_EVAL_DIV:
4011dda9c1aSJeremy L Thompson     case CEED_EVAL_CURL:
4021dda9c1aSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
4031dda9c1aSJeremy L Thompson       // LCOV_EXCL_STOP
4041dda9c1aSJeremy L Thompson   }
4051dda9c1aSJeremy L Thompson 
4061dda9c1aSJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
4071dda9c1aSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArrayRead(x_ref, &d_x));
4081dda9c1aSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
4091dda9c1aSJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
4101dda9c1aSJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
4119bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
4121dda9c1aSJeremy L Thompson   return CEED_ERROR_SUCCESS;
4131dda9c1aSJeremy L Thompson }
4141dda9c1aSJeremy L Thompson 
415db2becc9SJeremy L Thompson static int CeedBasisApplyAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
416db2becc9SJeremy L Thompson                                               CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
417db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda_shared(basis, false, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
418db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
419db2becc9SJeremy L Thompson }
420db2becc9SJeremy L Thompson 
421db2becc9SJeremy L Thompson static int CeedBasisApplyAddAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
422db2becc9SJeremy L Thompson                                                  CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
423db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda_shared(basis, true, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
424db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
425db2becc9SJeremy L Thompson }
426db2becc9SJeremy L Thompson 
4271dda9c1aSJeremy L Thompson //------------------------------------------------------------------------------
428*9ff05d55SJeremy L Thompson // Apply non-tensor basis
429*9ff05d55SJeremy L Thompson //------------------------------------------------------------------------------
430*9ff05d55SJeremy L Thompson static int CeedBasisApplyNonTensorCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode,
431*9ff05d55SJeremy L Thompson                                                    CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
432*9ff05d55SJeremy L Thompson   Ceed                   ceed;
433*9ff05d55SJeremy L Thompson   Ceed_Cuda             *ceed_Cuda;
434*9ff05d55SJeremy L Thompson   CeedInt                dim;
435*9ff05d55SJeremy L Thompson   const CeedScalar      *d_u;
436*9ff05d55SJeremy L Thompson   CeedScalar            *d_v;
437*9ff05d55SJeremy L Thompson   CeedBasis_Cuda_shared *data;
438*9ff05d55SJeremy L Thompson 
439*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
440*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
441*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
442*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
443*9ff05d55SJeremy L Thompson 
444*9ff05d55SJeremy L Thompson   // Get read/write access to u, v
445*9ff05d55SJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
446*9ff05d55SJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
447*9ff05d55SJeremy L Thompson   if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
448*9ff05d55SJeremy L Thompson   else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
449*9ff05d55SJeremy L Thompson 
450*9ff05d55SJeremy L Thompson   // Apply basis operation
451*9ff05d55SJeremy L Thompson   switch (eval_mode) {
452*9ff05d55SJeremy L Thompson     case CEED_EVAL_INTERP: {
453*9ff05d55SJeremy L Thompson       CeedInt P, Q;
454*9ff05d55SJeremy L Thompson 
455*9ff05d55SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes(basis, &P));
456*9ff05d55SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q));
457*9ff05d55SJeremy L Thompson       CeedInt thread = CeedIntMax(Q, P);
458*9ff05d55SJeremy L Thompson 
459*9ff05d55SJeremy L Thompson       CeedCallBackend(CeedInit_CudaInterp(data->d_interp_1d, P, Q, &data->c_B));
460*9ff05d55SJeremy L Thompson       void *interp_args[] = {(void *)&num_elem, &data->c_B, &d_u, &d_v};
461*9ff05d55SJeremy L Thompson 
462*9ff05d55SJeremy L Thompson       {
463*9ff05d55SJeremy L Thompson         // avoid >512 total threads
464*9ff05d55SJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1));
465*9ff05d55SJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
466*9ff05d55SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread * sizeof(CeedScalar);
467*9ff05d55SJeremy L Thompson 
468*9ff05d55SJeremy L Thompson         if (t_mode == CEED_TRANSPOSE) {
469*9ff05d55SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread, 1,
470*9ff05d55SJeremy L Thompson                                                       elems_per_block, shared_mem, interp_args));
471*9ff05d55SJeremy L Thompson         } else {
472*9ff05d55SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread, 1, elems_per_block, shared_mem, interp_args));
473*9ff05d55SJeremy L Thompson         }
474*9ff05d55SJeremy L Thompson       }
475*9ff05d55SJeremy L Thompson     } break;
476*9ff05d55SJeremy L Thompson     case CEED_EVAL_GRAD: {
477*9ff05d55SJeremy L Thompson       CeedInt P, Q;
478*9ff05d55SJeremy L Thompson 
479*9ff05d55SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes(basis, &P));
480*9ff05d55SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q));
481*9ff05d55SJeremy L Thompson       CeedInt thread = CeedIntMax(Q, P);
482*9ff05d55SJeremy L Thompson 
483*9ff05d55SJeremy L Thompson       CeedCallBackend(CeedInit_CudaInterp(data->d_grad_1d, P, Q * dim, &data->c_G));
484*9ff05d55SJeremy L Thompson       void *grad_args[] = {(void *)&num_elem, &data->c_G, &d_u, &d_v};
485*9ff05d55SJeremy L Thompson 
486*9ff05d55SJeremy L Thompson       {
487*9ff05d55SJeremy L Thompson         // avoid >512 total threads
488*9ff05d55SJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1));
489*9ff05d55SJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
490*9ff05d55SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread * sizeof(CeedScalar);
491*9ff05d55SJeremy L Thompson 
492*9ff05d55SJeremy L Thompson         if (t_mode == CEED_TRANSPOSE) {
493*9ff05d55SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread, 1,
494*9ff05d55SJeremy L Thompson                                                       elems_per_block, shared_mem, grad_args));
495*9ff05d55SJeremy L Thompson         } else {
496*9ff05d55SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread, 1, elems_per_block, shared_mem, grad_args));
497*9ff05d55SJeremy L Thompson         }
498*9ff05d55SJeremy L Thompson       }
499*9ff05d55SJeremy L Thompson     } break;
500*9ff05d55SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
501*9ff05d55SJeremy L Thompson       CeedInt Q;
502*9ff05d55SJeremy L Thompson 
503*9ff05d55SJeremy L Thompson       CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights_1d not set", CeedEvalModes[eval_mode]);
504*9ff05d55SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q));
505*9ff05d55SJeremy L Thompson       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
506*9ff05d55SJeremy L Thompson 
507*9ff05d55SJeremy L Thompson       {
508*9ff05d55SJeremy L Thompson         // avoid >512 total threads
509*9ff05d55SJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / Q, 1));
510*9ff05d55SJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
511*9ff05d55SJeremy L Thompson 
512*9ff05d55SJeremy L Thompson         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, Q, elems_per_block, 1, weight_args));
513*9ff05d55SJeremy L Thompson       }
514*9ff05d55SJeremy L Thompson     } break;
515*9ff05d55SJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
516*9ff05d55SJeremy L Thompson       break;
517*9ff05d55SJeremy L Thompson     // LCOV_EXCL_START
518*9ff05d55SJeremy L Thompson     case CEED_EVAL_DIV:
519*9ff05d55SJeremy L Thompson     case CEED_EVAL_CURL:
520*9ff05d55SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
521*9ff05d55SJeremy L Thompson       // LCOV_EXCL_STOP
522*9ff05d55SJeremy L Thompson   }
523*9ff05d55SJeremy L Thompson 
524*9ff05d55SJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
525*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
526*9ff05d55SJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
527*9ff05d55SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
528*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
529*9ff05d55SJeremy L Thompson   return CEED_ERROR_SUCCESS;
530*9ff05d55SJeremy L Thompson }
531*9ff05d55SJeremy L Thompson 
532*9ff05d55SJeremy L Thompson static int CeedBasisApplyNonTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
533*9ff05d55SJeremy L Thompson                                                CeedVector u, CeedVector v) {
534*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda_shared(basis, false, num_elem, t_mode, eval_mode, u, v));
535*9ff05d55SJeremy L Thompson   return CEED_ERROR_SUCCESS;
536*9ff05d55SJeremy L Thompson }
537*9ff05d55SJeremy L Thompson 
538*9ff05d55SJeremy L Thompson static int CeedBasisApplyAddNonTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
539*9ff05d55SJeremy L Thompson                                                   CeedVector u, CeedVector v) {
540*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda_shared(basis, true, num_elem, t_mode, eval_mode, u, v));
541*9ff05d55SJeremy L Thompson   return CEED_ERROR_SUCCESS;
542*9ff05d55SJeremy L Thompson }
543*9ff05d55SJeremy L Thompson 
544*9ff05d55SJeremy L Thompson //------------------------------------------------------------------------------
545ab213215SJeremy L Thompson // Destroy basis
546ab213215SJeremy L Thompson //------------------------------------------------------------------------------
547c532df63SYohann static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) {
548c532df63SYohann   Ceed                   ceed;
549c532df63SYohann   CeedBasis_Cuda_shared *data;
550ca735530SJeremy L Thompson 
551ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
5522b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
5532b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
5541dda9c1aSJeremy L Thompson   if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
555097cc795SJames Wright   if (data->d_q_weight_1d) CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d));
556111870feSJeremy L Thompson   CeedCallBackend(CeedFree(&data->h_points_per_elem));
557111870feSJeremy L Thompson   if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem));
5582b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp_1d));
5592b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad_1d));
5602b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_collo_grad_1d));
5611dda9c1aSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_chebyshev_interp_1d));
5622b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
5639bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
564e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
565c532df63SYohann }
566c532df63SYohann 
567ab213215SJeremy L Thompson //------------------------------------------------------------------------------
568ab213215SJeremy L Thompson // Create tensor basis
569ab213215SJeremy L Thompson //------------------------------------------------------------------------------
5702b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
5712b730f8bSJeremy L Thompson                                         const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
572c532df63SYohann   Ceed                   ceed;
573ca735530SJeremy L Thompson   CeedInt                num_comp;
574ca735530SJeremy L Thompson   const CeedInt          q_bytes      = Q_1d * sizeof(CeedScalar);
575ca735530SJeremy L Thompson   const CeedInt          interp_bytes = q_bytes * P_1d;
576c532df63SYohann   CeedBasis_Cuda_shared *data;
577ca735530SJeremy L Thompson 
578ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
5792b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
580c532df63SYohann 
581ab213215SJeremy L Thompson   // Copy basis data to GPU
582097cc795SJames Wright   if (q_weight_1d) {
5832b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
5842b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice));
585097cc795SJames Wright   }
5862b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
5872b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice));
5882b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes));
5892b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice));
590c532df63SYohann 
591ab213215SJeremy L Thompson   // Compute collocated gradient and copy to GPU
592437930d1SJeremy L Thompson   data->d_collo_grad_1d    = NULL;
5939e201c85SYohann   bool has_collocated_grad = dim == 3 && Q_1d >= P_1d;
594ca735530SJeremy L Thompson 
5959e201c85SYohann   if (has_collocated_grad) {
596437930d1SJeremy L Thompson     CeedScalar *collo_grad_1d;
597ca735530SJeremy L Thompson 
5982b730f8bSJeremy L Thompson     CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d));
5992b730f8bSJeremy L Thompson     CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d));
6002b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d));
6012b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, cudaMemcpyHostToDevice));
6022b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&collo_grad_1d));
603ac421f39SYohann   }
604ac421f39SYohann 
605ab213215SJeremy L Thompson   // Compile basis kernels
6069c25dd66SJeremy L Thompson   const char basis_kernel_source[] = "// Tensor basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-tensor.h>\n";
6079c25dd66SJeremy L Thompson 
6082b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
609eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D",
610eb7e6cafSJeremy L Thompson                                    CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim),
611eb7e6cafSJeremy L Thompson                                    "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_HAS_COLLOCATED_GRAD", has_collocated_grad));
612eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
613eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
614db2becc9SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd));
615eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
616eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose));
617db2becc9SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd));
618eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
619c532df63SYohann 
6202b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
621ab213215SJeremy L Thompson 
622ab213215SJeremy L Thompson   // Register backend functions
6232b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared));
624db2becc9SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddTensor_Cuda_shared));
6251dda9c1aSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAtPoints", CeedBasisApplyAtPoints_Cuda_shared));
626d474dafeSZach Atkins   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAddAtPoints", CeedBasisApplyAddAtPoints_Cuda_shared));
6272b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared));
6289bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
629e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
630c532df63SYohann }
6312a86cc9dSSebastian Grimberg 
632ab213215SJeremy L Thompson //------------------------------------------------------------------------------
633*9ff05d55SJeremy L Thompson // Create non-tensor basis
634*9ff05d55SJeremy L Thompson //------------------------------------------------------------------------------
635*9ff05d55SJeremy L Thompson int CeedBasisCreateH1_Cuda_shared(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
636*9ff05d55SJeremy L Thompson                                   const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
637*9ff05d55SJeremy L Thompson   Ceed                   ceed;
638*9ff05d55SJeremy L Thompson   CeedInt                num_comp, q_comp_interp, q_comp_grad;
639*9ff05d55SJeremy L Thompson   const CeedInt          q_bytes = num_qpts * sizeof(CeedScalar);
640*9ff05d55SJeremy L Thompson   CeedBasis_Cuda_shared *data;
641*9ff05d55SJeremy L Thompson 
642*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
643*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
644*9ff05d55SJeremy L Thompson 
645*9ff05d55SJeremy L Thompson   // Copy basis data to GPU
646*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
647*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad));
648*9ff05d55SJeremy L Thompson   if (q_weight) {
649*9ff05d55SJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
650*9ff05d55SJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight, q_bytes, cudaMemcpyHostToDevice));
651*9ff05d55SJeremy L Thompson   }
652*9ff05d55SJeremy L Thompson   if (interp) {
653*9ff05d55SJeremy L Thompson     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
654*9ff05d55SJeremy L Thompson 
655*9ff05d55SJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
656*9ff05d55SJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp, interp_bytes, cudaMemcpyHostToDevice));
657*9ff05d55SJeremy L Thompson   }
658*9ff05d55SJeremy L Thompson   if (grad) {
659*9ff05d55SJeremy L Thompson     const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad;
660*9ff05d55SJeremy L Thompson 
661*9ff05d55SJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, grad_bytes));
662*9ff05d55SJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad, grad_bytes, cudaMemcpyHostToDevice));
663*9ff05d55SJeremy L Thompson   }
664*9ff05d55SJeremy L Thompson 
665*9ff05d55SJeremy L Thompson   // Compile basis kernels
666*9ff05d55SJeremy L Thompson   const char basis_kernel_source[] = "// Non-tensor basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-nontensor.h>\n";
667*9ff05d55SJeremy L Thompson 
668*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
669*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "T_1D",
670*9ff05d55SJeremy L Thompson                                    CeedIntMax(num_qpts, num_nodes), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp));
671*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
672*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
673*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd));
674*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
675*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose));
676*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd));
677*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
678*9ff05d55SJeremy L Thompson 
679*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
680*9ff05d55SJeremy L Thompson 
681*9ff05d55SJeremy L Thompson   // Register backend functions
682*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda_shared));
683*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Cuda_shared));
684*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared));
685*9ff05d55SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
686*9ff05d55SJeremy L Thompson   return CEED_ERROR_SUCCESS;
687*9ff05d55SJeremy L Thompson }
688*9ff05d55SJeremy L Thompson 
689*9ff05d55SJeremy L Thompson //------------------------------------------------------------------------------
690