xref: /libCEED/backends/cuda-shared/ceed-cuda-shared-basis.c (revision e9c76bddc0f2a44f522e0176ed6b7e0c0aa1df73)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3c532df63SYohann //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5c532df63SYohann //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7c532df63SYohann 
849aac155SJeremy L Thompson #include <ceed.h>
9ec3da8bcSJed Brown #include <ceed/backend.h>
10437930d1SJeremy L Thompson #include <ceed/jit-tools.h>
113d576824SJeremy L Thompson #include <cuda.h>
123d576824SJeremy L Thompson #include <cuda_runtime.h>
1349aac155SJeremy L Thompson #include <stdbool.h>
143d576824SJeremy L Thompson #include <stddef.h>
15111870feSJeremy L Thompson #include <string.h>
16c532df63SYohann 
1749aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
182b730f8bSJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
192b730f8bSJeremy L Thompson #include "ceed-cuda-shared.h"
20c532df63SYohann 
21ab213215SJeremy L Thompson //------------------------------------------------------------------------------
229ff05d55SJeremy L Thompson // Apply tensor basis
23ab213215SJeremy L Thompson //------------------------------------------------------------------------------
24db2becc9SJeremy L Thompson static int CeedBasisApplyTensorCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode,
25db2becc9SJeremy L Thompson                                                 CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
26c532df63SYohann   Ceed                   ceed;
276dbfb411Snbeams   Ceed_Cuda             *ceed_Cuda;
28437930d1SJeremy L Thompson   CeedInt                dim, num_comp;
29ca735530SJeremy L Thompson   const CeedScalar      *d_u;
30ca735530SJeremy L Thompson   CeedScalar            *d_v;
31ca735530SJeremy L Thompson   CeedBasis_Cuda_shared *data;
32ca735530SJeremy L Thompson 
33ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
34ca735530SJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
35ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
362b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
372b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
38c532df63SYohann 
399ea2cfd9SJeremy L Thompson   // Get read/write access to u, v
406574a04fSJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
416574a04fSJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
42759e0bc3SJeremy L Thompson   if (apply_add) {
43759e0bc3SJeremy L Thompson     CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
44759e0bc3SJeremy L Thompson   } else {
45759e0bc3SJeremy L Thompson     CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
46759e0bc3SJeremy L Thompson   }
47c532df63SYohann 
48ab213215SJeremy L Thompson   // Apply basis operation
49437930d1SJeremy L Thompson   switch (eval_mode) {
50ab213215SJeremy L Thompson     case CEED_EVAL_INTERP: {
51437930d1SJeremy L Thompson       CeedInt P_1d, Q_1d;
52ca735530SJeremy L Thompson 
534cbc44e0SJeremy L Thompson       CeedCheck(data->d_interp_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; interp_1d not set", CeedEvalModes[eval_mode]);
542b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
552b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
56437930d1SJeremy L Thompson       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
57ca735530SJeremy L Thompson 
58aa4002adSJeremy L Thompson       void *interp_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_u, &d_v};
59ca735530SJeremy L Thompson 
604d537eeaSYohann       if (dim == 1) {
61dac82721SJeremy L Thompson         // avoid >512 total threads
62dac82721SJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1));
63a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
64437930d1SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
65b2165e7aSSebastian Grimberg 
669e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
67*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, NULL, grid, thread_1d, 1,
68db2becc9SJeremy L Thompson                                                       elems_per_block, shared_mem, interp_args));
699e201c85SYohann         } else {
70*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, NULL, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
719e201c85SYohann         }
72074be161SYohann Dudouit       } else if (dim == 2) {
73437930d1SJeremy L Thompson         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
74437930d1SJeremy L Thompson         // elems_per_block must be at least 1
752b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
76a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
772b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
78b2165e7aSSebastian Grimberg 
799e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
80*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, NULL, grid, thread_1d,
81*e9c76bddSJeremy L Thompson                                                       thread_1d, elems_per_block, shared_mem, interp_args));
829e201c85SYohann         } else {
83*e9c76bddSJeremy L Thompson           CeedCallBackend(
84*e9c76bddSJeremy L Thompson               CeedRunKernelDimShared_Cuda(ceed, data->Interp, NULL, 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;
88a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
892b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
90b2165e7aSSebastian Grimberg 
919e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
92*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, NULL, grid, thread_1d,
93*e9c76bddSJeremy L Thompson                                                       thread_1d, elems_per_block, shared_mem, interp_args));
949e201c85SYohann         } else {
95*e9c76bddSJeremy L Thompson           CeedCallBackend(
96*e9c76bddSJeremy L Thompson               CeedRunKernelDimShared_Cuda(ceed, data->Interp, NULL, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
97074be161SYohann Dudouit         }
989e201c85SYohann       }
99ab213215SJeremy L Thompson     } break;
100ab213215SJeremy L Thompson     case CEED_EVAL_GRAD: {
101437930d1SJeremy L Thompson       CeedInt P_1d, Q_1d;
102ca735530SJeremy L Thompson 
1034cbc44e0SJeremy L Thompson       CeedCheck(data->d_grad_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; grad_1d not set", CeedEvalModes[eval_mode]);
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);
107aa4002adSJeremy L Thompson       CeedScalar *d_grad_1d = data->d_grad_1d;
108ca735530SJeremy L Thompson 
1099e201c85SYohann       if (data->d_collo_grad_1d) {
110aa4002adSJeremy L Thompson         d_grad_1d = data->d_collo_grad_1d;
1119e201c85SYohann       }
112aa4002adSJeremy L Thompson       void *grad_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_grad_1d, &d_u, &d_v};
113aa4002adSJeremy L Thompson 
1144d537eeaSYohann       if (dim == 1) {
115dac82721SJeremy L Thompson         // avoid >512 total threads
116dac82721SJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1));
117a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
118437930d1SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
119b2165e7aSSebastian Grimberg 
1209e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
121*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, NULL, grid, thread_1d, 1,
122db2becc9SJeremy L Thompson                                                       elems_per_block, shared_mem, grad_args));
1239e201c85SYohann         } else {
124*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, NULL, 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);
130a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 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*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, NULL, grid, thread_1d,
135*e9c76bddSJeremy L Thompson                                                       thread_1d, elems_per_block, shared_mem, grad_args));
1369e201c85SYohann         } else {
137*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, NULL, 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;
141a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 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*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, NULL, grid, thread_1d,
146*e9c76bddSJeremy L Thompson                                                       thread_1d, elems_per_block, shared_mem, grad_args));
1479e201c85SYohann         } else {
148*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, NULL, 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;
161a8d440fbSJeremy L Thompson         const CeedInt grid_size       = num_elem / elems_per_block + (num_elem % elems_per_block > 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;
167a8d440fbSJeremy L Thompson         const CeedInt grid_size       = num_elem / elems_per_block + (num_elem % elems_per_block > 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;
173a8d440fbSJeremy L Thompson         const CeedInt grid_size       = num_elem / elems_per_block + (num_elem % elems_per_block > 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));
1919bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
192e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
193c532df63SYohann }
194c532df63SYohann 
195db2becc9SJeremy L Thompson static int CeedBasisApplyTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
196db2becc9SJeremy L Thompson                                             CeedVector v) {
197db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyTensorCore_Cuda_shared(basis, false, num_elem, t_mode, eval_mode, u, v));
198db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
199db2becc9SJeremy L Thompson }
200db2becc9SJeremy L Thompson 
201db2becc9SJeremy L Thompson static int CeedBasisApplyAddTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
202db2becc9SJeremy L Thompson                                                CeedVector u, CeedVector v) {
203db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyTensorCore_Cuda_shared(basis, true, num_elem, t_mode, eval_mode, u, v));
204db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
205db2becc9SJeremy L Thompson }
206db2becc9SJeremy L Thompson 
207ab213215SJeremy L Thompson //------------------------------------------------------------------------------
2081dda9c1aSJeremy L Thompson // Basis apply - tensor AtPoints
2091dda9c1aSJeremy L Thompson //------------------------------------------------------------------------------
210db2becc9SJeremy L Thompson static int CeedBasisApplyAtPointsCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, const CeedInt *num_points,
211db2becc9SJeremy L Thompson                                                   CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
2121dda9c1aSJeremy L Thompson   Ceed                   ceed;
2139e1d4b82SJeremy L Thompson   Ceed_Cuda             *ceed_Cuda;
2149e1d4b82SJeremy L Thompson   CeedInt                Q_1d, dim, num_comp, max_num_points = num_points[0];
2151dda9c1aSJeremy L Thompson   const CeedInt          is_transpose = t_mode == CEED_TRANSPOSE;
2161dda9c1aSJeremy L Thompson   const CeedScalar      *d_x, *d_u;
2171dda9c1aSJeremy L Thompson   CeedScalar            *d_v;
2181dda9c1aSJeremy L Thompson   CeedBasis_Cuda_shared *data;
2191dda9c1aSJeremy L Thompson 
2201dda9c1aSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
2211dda9c1aSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
2221dda9c1aSJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
2239e1d4b82SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
2241dda9c1aSJeremy L Thompson 
2251dda9c1aSJeremy L Thompson   // Weight handled separately
2261dda9c1aSJeremy L Thompson   if (eval_mode == CEED_EVAL_WEIGHT) {
2275a5594ffSJeremy L Thompson     CeedCallBackend(CeedVectorSetValue(v, 1.0));
2281dda9c1aSJeremy L Thompson     return CEED_ERROR_SUCCESS;
2291dda9c1aSJeremy L Thompson   }
2301dda9c1aSJeremy L Thompson 
2319bc66399SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
2329e1d4b82SJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
2339bc66399SJeremy L Thompson 
234111870feSJeremy L Thompson   // Check padded to uniform number of points per elem
235111870feSJeremy L Thompson   for (CeedInt i = 1; i < num_elem; i++) max_num_points = CeedIntMax(max_num_points, num_points[i]);
236111870feSJeremy L Thompson   {
2379e1d4b82SJeremy L Thompson     CeedInt  q_comp;
238111870feSJeremy L Thompson     CeedSize len, len_required;
239111870feSJeremy L Thompson     CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
240111870feSJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(is_transpose ? u : v, &len));
241111870feSJeremy L Thompson     len_required = (CeedSize)num_comp * (CeedSize)q_comp * (CeedSize)num_elem * (CeedSize)max_num_points;
242111870feSJeremy L Thompson     CeedCheck(len >= len_required, ceed, CEED_ERROR_BACKEND,
243111870feSJeremy L Thompson               "Vector at points must be padded to the same number of points in each element for BasisApplyAtPoints on GPU backends."
244111870feSJeremy L Thompson               " Found %" CeedSize_FMT ", Required %" CeedSize_FMT,
245111870feSJeremy L Thompson               len, len_required);
246111870feSJeremy L Thompson   }
247111870feSJeremy L Thompson 
248111870feSJeremy L Thompson   // Move num_points array to device
249111870feSJeremy L Thompson   if (is_transpose) {
250111870feSJeremy L Thompson     const CeedInt num_bytes = num_elem * sizeof(CeedInt);
251111870feSJeremy L Thompson 
252111870feSJeremy L Thompson     if (num_elem != data->num_elem_at_points) {
253111870feSJeremy L Thompson       data->num_elem_at_points = num_elem;
254111870feSJeremy L Thompson 
255111870feSJeremy L Thompson       if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem));
256111870feSJeremy L Thompson       CeedCallCuda(ceed, cudaMalloc((void **)&data->d_points_per_elem, num_bytes));
257111870feSJeremy L Thompson       CeedCallBackend(CeedFree(&data->h_points_per_elem));
258111870feSJeremy L Thompson       CeedCallBackend(CeedCalloc(num_elem, &data->h_points_per_elem));
259111870feSJeremy L Thompson     }
2609e511c80SJeremy L Thompson     if (memcmp(data->h_points_per_elem, num_points, num_bytes)) {
261111870feSJeremy L Thompson       memcpy(data->h_points_per_elem, num_points, num_bytes);
262111870feSJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy(data->d_points_per_elem, num_points, num_bytes, cudaMemcpyHostToDevice));
263111870feSJeremy L Thompson     }
264111870feSJeremy L Thompson   }
265111870feSJeremy L Thompson 
2661dda9c1aSJeremy L Thompson   // Build kernels if needed
2671dda9c1aSJeremy L Thompson   if (data->num_points != max_num_points) {
2681dda9c1aSJeremy L Thompson     CeedInt P_1d;
2691dda9c1aSJeremy L Thompson 
2701dda9c1aSJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
2711dda9c1aSJeremy L Thompson     data->num_points = max_num_points;
2721dda9c1aSJeremy L Thompson 
2731dda9c1aSJeremy L Thompson     // -- Create interp matrix to Chebyshev coefficients
2741dda9c1aSJeremy L Thompson     if (!data->d_chebyshev_interp_1d) {
2751dda9c1aSJeremy L Thompson       CeedSize    interp_bytes;
2761dda9c1aSJeremy L Thompson       CeedScalar *chebyshev_interp_1d;
2771dda9c1aSJeremy L Thompson 
2781dda9c1aSJeremy L Thompson       interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
2791dda9c1aSJeremy L Thompson       CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
2805a5594ffSJeremy L Thompson       CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
2811dda9c1aSJeremy L Thompson       CeedCallCuda(ceed, cudaMalloc((void **)&data->d_chebyshev_interp_1d, interp_bytes));
2821dda9c1aSJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy(data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
2831dda9c1aSJeremy L Thompson       CeedCallBackend(CeedFree(&chebyshev_interp_1d));
2841dda9c1aSJeremy L Thompson     }
2851dda9c1aSJeremy L Thompson 
2861dda9c1aSJeremy L Thompson     // -- Compile kernels
2879e1d4b82SJeremy L Thompson     const char basis_kernel_source[] = "// AtPoints basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-at-points.h>\n";
2881dda9c1aSJeremy L Thompson     CeedInt    num_comp;
2891dda9c1aSJeremy L Thompson 
2901dda9c1aSJeremy L Thompson     if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
2911dda9c1aSJeremy L Thompson     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
2929e1d4b82SJeremy L Thompson     CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->moduleAtPoints, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D",
2939e1d4b82SJeremy L Thompson                                      CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim),
2949e1d4b82SJeremy L Thompson                                      "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_NUM_PTS", max_num_points));
2951dda9c1aSJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpAtPoints", &data->InterpAtPoints));
29681ae6159SJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpTransposeAtPoints", &data->InterpTransposeAtPoints));
297af0e6e89SJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpTransposeAddAtPoints", &data->InterpTransposeAddAtPoints));
2981dda9c1aSJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradAtPoints", &data->GradAtPoints));
29981ae6159SJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradTransposeAtPoints", &data->GradTransposeAtPoints));
300af0e6e89SJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradTransposeAddAtPoints", &data->GradTransposeAddAtPoints));
3011dda9c1aSJeremy L Thompson   }
3021dda9c1aSJeremy L Thompson 
3031dda9c1aSJeremy L Thompson   // Get read/write access to u, v
3041dda9c1aSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayRead(x_ref, CEED_MEM_DEVICE, &d_x));
3051dda9c1aSJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
3061dda9c1aSJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
3079dafd6dfSZach Atkins   if (apply_add) {
3089dafd6dfSZach Atkins     CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
3099dafd6dfSZach Atkins   } else {
3109dafd6dfSZach Atkins     CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
3111dda9c1aSJeremy L Thompson   }
3121dda9c1aSJeremy L Thompson 
3131dda9c1aSJeremy L Thompson   // Basis action
3141dda9c1aSJeremy L Thompson   switch (eval_mode) {
3151dda9c1aSJeremy L Thompson     case CEED_EVAL_INTERP: {
3169e1d4b82SJeremy L Thompson       CeedInt P_1d, Q_1d;
3171dda9c1aSJeremy L Thompson 
3189e1d4b82SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
3199e1d4b82SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
3209e1d4b82SJeremy L Thompson       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
3219e1d4b82SJeremy L Thompson 
322aa4002adSJeremy L Thompson       void *interp_args[] = {(void *)&num_elem, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v};
3239e1d4b82SJeremy L Thompson 
3249e1d4b82SJeremy L Thompson       if (dim == 1) {
325dac82721SJeremy L Thompson         // avoid >512 total threads
326dac82721SJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1));
327a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
3289e1d4b82SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
3299e1d4b82SJeremy L Thompson 
330af0e6e89SJeremy L Thompson         if (is_transpose) {
331*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAddAtPoints : data->InterpTransposeAtPoints, NULL, grid,
332af0e6e89SJeremy L Thompson                                                       thread_1d, 1, elems_per_block, shared_mem, interp_args));
333af0e6e89SJeremy L Thompson         } else {
334*e9c76bddSJeremy L Thompson           CeedCallBackend(
335*e9c76bddSJeremy L Thompson               CeedRunKernelDimShared_Cuda(ceed, data->InterpAtPoints, NULL, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
336af0e6e89SJeremy L Thompson         }
3379e1d4b82SJeremy L Thompson       } else if (dim == 2) {
3389e1d4b82SJeremy L Thompson         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
3399e1d4b82SJeremy L Thompson         // elems_per_block must be at least 1
3409e1d4b82SJeremy L Thompson         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
341a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
3429e1d4b82SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
3439e1d4b82SJeremy L Thompson 
344af0e6e89SJeremy L Thompson         if (is_transpose) {
345*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAddAtPoints : data->InterpTransposeAtPoints, NULL, grid,
346af0e6e89SJeremy L Thompson                                                       thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
347af0e6e89SJeremy L Thompson         } else {
348af0e6e89SJeremy L Thompson           CeedCallBackend(
349*e9c76bddSJeremy L Thompson               CeedRunKernelDimShared_Cuda(ceed, data->InterpAtPoints, NULL, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
350af0e6e89SJeremy L Thompson         }
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 
356af0e6e89SJeremy L Thompson         if (is_transpose) {
357*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAddAtPoints : data->InterpTransposeAtPoints, NULL, grid,
358af0e6e89SJeremy L Thompson                                                       thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
359af0e6e89SJeremy L Thompson         } else {
360af0e6e89SJeremy L Thompson           CeedCallBackend(
361*e9c76bddSJeremy L Thompson               CeedRunKernelDimShared_Cuda(ceed, data->InterpAtPoints, NULL, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
362af0e6e89SJeremy L Thompson         }
3639e1d4b82SJeremy L Thompson       }
3641dda9c1aSJeremy L Thompson     } break;
3651dda9c1aSJeremy L Thompson     case CEED_EVAL_GRAD: {
3669e1d4b82SJeremy L Thompson       CeedInt P_1d, Q_1d;
3671dda9c1aSJeremy L Thompson 
3689e1d4b82SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
3699e1d4b82SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
3709e1d4b82SJeremy L Thompson       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
3719e1d4b82SJeremy L Thompson 
3729e1d4b82SJeremy L Thompson       void *grad_args[] = {(void *)&num_elem, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v};
3739e1d4b82SJeremy L Thompson 
3749e1d4b82SJeremy L Thompson       if (dim == 1) {
375dac82721SJeremy L Thompson         // avoid >512 total threads
376dac82721SJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1));
377a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
3789e1d4b82SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
3799e1d4b82SJeremy L Thompson 
380af0e6e89SJeremy L Thompson         if (is_transpose) {
381*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAddAtPoints : data->GradTransposeAtPoints, NULL, grid,
382*e9c76bddSJeremy L Thompson                                                       thread_1d, 1, elems_per_block, shared_mem, grad_args));
383af0e6e89SJeremy L Thompson         } else {
384*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->GradAtPoints, NULL, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
385af0e6e89SJeremy L Thompson         }
3869e1d4b82SJeremy L Thompson       } else if (dim == 2) {
3879e1d4b82SJeremy L Thompson         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
3889e1d4b82SJeremy L Thompson         // elems_per_block must be at least 1
3899e1d4b82SJeremy L Thompson         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
390a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
3919e1d4b82SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
3929e1d4b82SJeremy L Thompson 
393af0e6e89SJeremy L Thompson         if (is_transpose) {
394*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAddAtPoints : data->GradTransposeAtPoints, NULL, grid,
395*e9c76bddSJeremy L Thompson                                                       thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
396af0e6e89SJeremy L Thompson         } else {
397*e9c76bddSJeremy L Thompson           CeedCallBackend(
398*e9c76bddSJeremy L Thompson               CeedRunKernelDimShared_Cuda(ceed, data->GradAtPoints, NULL, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
399af0e6e89SJeremy L Thompson         }
4009e1d4b82SJeremy L Thompson       } else if (dim == 3) {
4019e1d4b82SJeremy L Thompson         CeedInt elems_per_block = 1;
402a8d440fbSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
4039e1d4b82SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
4049e1d4b82SJeremy L Thompson 
405af0e6e89SJeremy L Thompson         if (is_transpose) {
406*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAddAtPoints : data->GradTransposeAtPoints, NULL, grid,
407*e9c76bddSJeremy L Thompson                                                       thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
408af0e6e89SJeremy L Thompson         } else {
409*e9c76bddSJeremy L Thompson           CeedCallBackend(
410*e9c76bddSJeremy L Thompson               CeedRunKernelDimShared_Cuda(ceed, data->GradAtPoints, NULL, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
411af0e6e89SJeremy L Thompson         }
4129e1d4b82SJeremy L Thompson       }
4131dda9c1aSJeremy L Thompson     } break;
4141dda9c1aSJeremy L Thompson     case CEED_EVAL_WEIGHT:
4151dda9c1aSJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
4161dda9c1aSJeremy L Thompson       break;
4171dda9c1aSJeremy L Thompson     // LCOV_EXCL_START
4181dda9c1aSJeremy L Thompson     case CEED_EVAL_DIV:
4191dda9c1aSJeremy L Thompson     case CEED_EVAL_CURL:
4201dda9c1aSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
4211dda9c1aSJeremy L Thompson       // LCOV_EXCL_STOP
4221dda9c1aSJeremy L Thompson   }
4231dda9c1aSJeremy L Thompson 
4241dda9c1aSJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
4251dda9c1aSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArrayRead(x_ref, &d_x));
4261dda9c1aSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
4271dda9c1aSJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
4281dda9c1aSJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
4299bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
4301dda9c1aSJeremy L Thompson   return CEED_ERROR_SUCCESS;
4311dda9c1aSJeremy L Thompson }
4321dda9c1aSJeremy L Thompson 
433db2becc9SJeremy L Thompson static int CeedBasisApplyAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
434db2becc9SJeremy L Thompson                                               CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
435db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda_shared(basis, false, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
436db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
437db2becc9SJeremy L Thompson }
438db2becc9SJeremy L Thompson 
439db2becc9SJeremy L Thompson static int CeedBasisApplyAddAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
440db2becc9SJeremy L Thompson                                                  CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
441db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda_shared(basis, true, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
442db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
443db2becc9SJeremy L Thompson }
444db2becc9SJeremy L Thompson 
4451dda9c1aSJeremy L Thompson //------------------------------------------------------------------------------
4469ff05d55SJeremy L Thompson // Apply non-tensor basis
4479ff05d55SJeremy L Thompson //------------------------------------------------------------------------------
4489ff05d55SJeremy L Thompson static int CeedBasisApplyNonTensorCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode,
4499ff05d55SJeremy L Thompson                                                    CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
4509ff05d55SJeremy L Thompson   Ceed                   ceed;
4519ff05d55SJeremy L Thompson   Ceed_Cuda             *ceed_Cuda;
4529ff05d55SJeremy L Thompson   CeedInt                dim;
4539ff05d55SJeremy L Thompson   const CeedScalar      *d_u;
4549ff05d55SJeremy L Thompson   CeedScalar            *d_v;
4559ff05d55SJeremy L Thompson   CeedBasis_Cuda_shared *data;
4569ff05d55SJeremy L Thompson 
4579ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
4589ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
4599ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
4609ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
4619ff05d55SJeremy L Thompson 
4629ff05d55SJeremy L Thompson   // Get read/write access to u, v
4639ff05d55SJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
4649ff05d55SJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
465759e0bc3SJeremy L Thompson   if (apply_add) {
466759e0bc3SJeremy L Thompson     CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
467759e0bc3SJeremy L Thompson   } else {
468759e0bc3SJeremy L Thompson     CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
469759e0bc3SJeremy L Thompson   }
4709ff05d55SJeremy L Thompson 
4719ff05d55SJeremy L Thompson   // Apply basis operation
4729ff05d55SJeremy L Thompson   switch (eval_mode) {
4739ff05d55SJeremy L Thompson     case CEED_EVAL_INTERP: {
4749ff05d55SJeremy L Thompson       CeedInt P, Q;
4759ff05d55SJeremy L Thompson 
4764cbc44e0SJeremy L Thompson       CeedCheck(data->d_interp_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; interp not set", CeedEvalModes[eval_mode]);
4779ff05d55SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes(basis, &P));
4789ff05d55SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q));
4799ff05d55SJeremy L Thompson       CeedInt thread = CeedIntMax(Q, P);
4809ff05d55SJeremy L Thompson 
481aa4002adSJeremy L Thompson       void *interp_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_u, &d_v};
4829ff05d55SJeremy L Thompson 
4839ff05d55SJeremy L Thompson       {
4849ff05d55SJeremy L Thompson         // avoid >512 total threads
4859ff05d55SJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1));
4869ff05d55SJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
4879ff05d55SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread * sizeof(CeedScalar);
4889ff05d55SJeremy L Thompson 
4899ff05d55SJeremy L Thompson         if (t_mode == CEED_TRANSPOSE) {
490*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, NULL, grid, thread, 1,
4919ff05d55SJeremy L Thompson                                                       elems_per_block, shared_mem, interp_args));
4929ff05d55SJeremy L Thompson         } else {
493*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, NULL, grid, thread, 1, elems_per_block, shared_mem, interp_args));
4949ff05d55SJeremy L Thompson         }
4959ff05d55SJeremy L Thompson       }
4969ff05d55SJeremy L Thompson     } break;
4979ff05d55SJeremy L Thompson     case CEED_EVAL_GRAD: {
4989ff05d55SJeremy L Thompson       CeedInt P, Q;
4999ff05d55SJeremy L Thompson 
5004cbc44e0SJeremy L Thompson       CeedCheck(data->d_grad_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; grad not set", CeedEvalModes[eval_mode]);
5019ff05d55SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes(basis, &P));
5029ff05d55SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q));
5039ff05d55SJeremy L Thompson       CeedInt thread = CeedIntMax(Q, P);
5049ff05d55SJeremy L Thompson 
505aa4002adSJeremy L Thompson       void *grad_args[] = {(void *)&num_elem, &data->d_grad_1d, &d_u, &d_v};
5069ff05d55SJeremy L Thompson 
5079ff05d55SJeremy L Thompson       {
5089ff05d55SJeremy L Thompson         // avoid >512 total threads
5099ff05d55SJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1));
5109ff05d55SJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
5119ff05d55SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread * sizeof(CeedScalar);
5129ff05d55SJeremy L Thompson 
5139ff05d55SJeremy L Thompson         if (t_mode == CEED_TRANSPOSE) {
514*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, NULL, grid, thread, 1,
5159ff05d55SJeremy L Thompson                                                       elems_per_block, shared_mem, grad_args));
5169ff05d55SJeremy L Thompson         } else {
517*e9c76bddSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, NULL, grid, thread, 1, elems_per_block, shared_mem, grad_args));
5189ff05d55SJeremy L Thompson         }
5199ff05d55SJeremy L Thompson       }
5209ff05d55SJeremy L Thompson     } break;
5219ff05d55SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
52297011eabSZach Atkins       CeedInt P, Q;
5239ff05d55SJeremy L Thompson 
5244cbc44e0SJeremy L Thompson       CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights not set", CeedEvalModes[eval_mode]);
52597011eabSZach Atkins       CeedCallBackend(CeedBasisGetNumNodes(basis, &P));
5269ff05d55SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q));
52797011eabSZach Atkins       CeedInt thread = CeedIntMax(Q, P);
52897011eabSZach Atkins 
5299ff05d55SJeremy L Thompson       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
5309ff05d55SJeremy L Thompson 
5319ff05d55SJeremy L Thompson       {
5329ff05d55SJeremy L Thompson         // avoid >512 total threads
53397011eabSZach Atkins         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1));
5349ff05d55SJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
5359ff05d55SJeremy L Thompson 
53697011eabSZach Atkins         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, thread, elems_per_block, 1, weight_args));
5379ff05d55SJeremy L Thompson       }
5389ff05d55SJeremy L Thompson     } break;
5399ff05d55SJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
5409ff05d55SJeremy L Thompson       break;
5419ff05d55SJeremy L Thompson     // LCOV_EXCL_START
5429ff05d55SJeremy L Thompson     case CEED_EVAL_DIV:
5439ff05d55SJeremy L Thompson     case CEED_EVAL_CURL:
5449ff05d55SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
5459ff05d55SJeremy L Thompson       // LCOV_EXCL_STOP
5469ff05d55SJeremy L Thompson   }
5479ff05d55SJeremy L Thompson 
5489ff05d55SJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
5499ff05d55SJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
5509ff05d55SJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
5519ff05d55SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
5529ff05d55SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
5539ff05d55SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5549ff05d55SJeremy L Thompson }
5559ff05d55SJeremy L Thompson 
5569ff05d55SJeremy L Thompson static int CeedBasisApplyNonTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
5579ff05d55SJeremy L Thompson                                                CeedVector u, CeedVector v) {
5589ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda_shared(basis, false, num_elem, t_mode, eval_mode, u, v));
5599ff05d55SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5609ff05d55SJeremy L Thompson }
5619ff05d55SJeremy L Thompson 
5629ff05d55SJeremy L Thompson static int CeedBasisApplyAddNonTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
5639ff05d55SJeremy L Thompson                                                   CeedVector u, CeedVector v) {
5649ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda_shared(basis, true, num_elem, t_mode, eval_mode, u, v));
5659ff05d55SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5669ff05d55SJeremy L Thompson }
5679ff05d55SJeremy L Thompson 
5689ff05d55SJeremy L Thompson //------------------------------------------------------------------------------
569ab213215SJeremy L Thompson // Destroy basis
570ab213215SJeremy L Thompson //------------------------------------------------------------------------------
571c532df63SYohann static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) {
572c532df63SYohann   Ceed                   ceed;
573c532df63SYohann   CeedBasis_Cuda_shared *data;
574ca735530SJeremy L Thompson 
575ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
5762b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
5772b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
5781dda9c1aSJeremy L Thompson   if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
579097cc795SJames Wright   if (data->d_q_weight_1d) CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d));
580111870feSJeremy L Thompson   CeedCallBackend(CeedFree(&data->h_points_per_elem));
581111870feSJeremy L Thompson   if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem));
5822b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp_1d));
5832b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad_1d));
5842b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_collo_grad_1d));
5851dda9c1aSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_chebyshev_interp_1d));
5862b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
5879bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
588e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
589c532df63SYohann }
590c532df63SYohann 
591ab213215SJeremy L Thompson //------------------------------------------------------------------------------
592ab213215SJeremy L Thompson // Create tensor basis
593ab213215SJeremy L Thompson //------------------------------------------------------------------------------
5942b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
5952b730f8bSJeremy L Thompson                                         const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
596c532df63SYohann   Ceed                   ceed;
597ca735530SJeremy L Thompson   CeedInt                num_comp;
598ca735530SJeremy L Thompson   const CeedInt          q_bytes      = Q_1d * sizeof(CeedScalar);
599ca735530SJeremy L Thompson   const CeedInt          interp_bytes = q_bytes * P_1d;
600c532df63SYohann   CeedBasis_Cuda_shared *data;
601ca735530SJeremy L Thompson 
602ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
6032b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
604c532df63SYohann 
605ab213215SJeremy L Thompson   // Copy basis data to GPU
606097cc795SJames Wright   if (q_weight_1d) {
6072b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
6082b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice));
609097cc795SJames Wright   }
6102b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
6112b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice));
6122b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes));
6132b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice));
614c532df63SYohann 
615ab213215SJeremy L Thompson   // Compute collocated gradient and copy to GPU
616437930d1SJeremy L Thompson   data->d_collo_grad_1d    = NULL;
6179e201c85SYohann   bool has_collocated_grad = dim == 3 && Q_1d >= P_1d;
618ca735530SJeremy L Thompson 
6199e201c85SYohann   if (has_collocated_grad) {
620437930d1SJeremy L Thompson     CeedScalar *collo_grad_1d;
621ca735530SJeremy L Thompson 
6222b730f8bSJeremy L Thompson     CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d));
6232b730f8bSJeremy L Thompson     CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d));
6242b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d));
6252b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, cudaMemcpyHostToDevice));
6262b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&collo_grad_1d));
627ac421f39SYohann   }
628ac421f39SYohann 
629ab213215SJeremy L Thompson   // Compile basis kernels
6309c25dd66SJeremy L Thompson   const char basis_kernel_source[] = "// Tensor basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-tensor.h>\n";
6319c25dd66SJeremy L Thompson 
6322b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
633eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D",
634eb7e6cafSJeremy L Thompson                                    CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim),
635eb7e6cafSJeremy L Thompson                                    "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_HAS_COLLOCATED_GRAD", has_collocated_grad));
636eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
637eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
638db2becc9SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd));
639eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
640eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose));
641db2becc9SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd));
642eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
643c532df63SYohann 
6442b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
645ab213215SJeremy L Thompson 
646ab213215SJeremy L Thompson   // Register backend functions
6472b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared));
648db2becc9SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddTensor_Cuda_shared));
6491dda9c1aSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAtPoints", CeedBasisApplyAtPoints_Cuda_shared));
650d474dafeSZach Atkins   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAddAtPoints", CeedBasisApplyAddAtPoints_Cuda_shared));
6512b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared));
6529bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
653e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
654c532df63SYohann }
6552a86cc9dSSebastian Grimberg 
656ab213215SJeremy L Thompson //------------------------------------------------------------------------------
6579ff05d55SJeremy L Thompson // Create non-tensor basis
6589ff05d55SJeremy L Thompson //------------------------------------------------------------------------------
6599ff05d55SJeremy L Thompson int CeedBasisCreateH1_Cuda_shared(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
6609ff05d55SJeremy L Thompson                                   const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
6619ff05d55SJeremy L Thompson   Ceed                   ceed;
6629ff05d55SJeremy L Thompson   CeedInt                num_comp, q_comp_interp, q_comp_grad;
6639ff05d55SJeremy L Thompson   const CeedInt          q_bytes = num_qpts * sizeof(CeedScalar);
6649ff05d55SJeremy L Thompson   CeedBasis_Cuda_shared *data;
6659ff05d55SJeremy L Thompson 
6669ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
6679ff05d55SJeremy L Thompson 
668fda26546SJeremy L Thompson   // Check shared memory size
669fda26546SJeremy L Thompson   {
670fda26546SJeremy L Thompson     Ceed_Cuda *cuda_data;
671fda26546SJeremy L Thompson 
672fda26546SJeremy L Thompson     CeedCallBackend(CeedGetData(ceed, &cuda_data));
673fda26546SJeremy L Thompson     if (((size_t)num_nodes * (size_t)num_qpts * (size_t)dim + (size_t)CeedIntMax(num_nodes, num_qpts)) * sizeof(CeedScalar) >
674fda26546SJeremy L Thompson         cuda_data->device_prop.sharedMemPerBlock) {
675fda26546SJeremy L Thompson       CeedCallBackend(CeedBasisCreateH1Fallback(ceed, topo, dim, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis));
676fda26546SJeremy L Thompson       CeedCallBackend(CeedDestroy(&ceed));
677fda26546SJeremy L Thompson       return CEED_ERROR_SUCCESS;
678fda26546SJeremy L Thompson     }
679fda26546SJeremy L Thompson   }
680fda26546SJeremy L Thompson 
681fda26546SJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
682aa4002adSJeremy L Thompson 
6839ff05d55SJeremy L Thompson   // Copy basis data to GPU
6849ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
6859ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad));
6869ff05d55SJeremy L Thompson   if (q_weight) {
6879ff05d55SJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
6889ff05d55SJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight, q_bytes, cudaMemcpyHostToDevice));
6899ff05d55SJeremy L Thompson   }
6909ff05d55SJeremy L Thompson   if (interp) {
6919ff05d55SJeremy L Thompson     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
6929ff05d55SJeremy L Thompson 
6939ff05d55SJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
6949ff05d55SJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp, interp_bytes, cudaMemcpyHostToDevice));
6959ff05d55SJeremy L Thompson   }
6969ff05d55SJeremy L Thompson   if (grad) {
6979ff05d55SJeremy L Thompson     const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad;
6989ff05d55SJeremy L Thompson 
6999ff05d55SJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, grad_bytes));
7009ff05d55SJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad, grad_bytes, cudaMemcpyHostToDevice));
7019ff05d55SJeremy L Thompson   }
7029ff05d55SJeremy L Thompson 
7039ff05d55SJeremy L Thompson   // Compile basis kernels
7049ff05d55SJeremy L Thompson   const char basis_kernel_source[] = "// Non-tensor basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-nontensor.h>\n";
7059ff05d55SJeremy L Thompson 
7069ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
7079ff05d55SJeremy L Thompson   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "T_1D",
7089ff05d55SJeremy L Thompson                                    CeedIntMax(num_qpts, num_nodes), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp));
7099ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
7109ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
7119ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd));
7129ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
7139ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose));
7149ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd));
7159ff05d55SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
7169ff05d55SJeremy L Thompson 
7179ff05d55SJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
7189ff05d55SJeremy L Thompson 
7199ff05d55SJeremy L Thompson   // Register backend functions
7209ff05d55SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda_shared));
7219ff05d55SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Cuda_shared));
7229ff05d55SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared));
7239ff05d55SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
7249ff05d55SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7259ff05d55SJeremy L Thompson }
7269ff05d55SJeremy L Thompson 
7279ff05d55SJeremy L Thompson //------------------------------------------------------------------------------
728