xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-shared/ceed-cuda-shared-basis.c (revision 39577a105a4a62dcc151f47b13c7c1be81f50807)
1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 #include <ceed.h>
9 #include <ceed/backend.h>
10 #include <ceed/jit-tools.h>
11 #include <cuda.h>
12 #include <cuda_runtime.h>
13 #include <stdbool.h>
14 #include <stddef.h>
15 #include <string.h>
16 
17 #include "../cuda/ceed-cuda-common.h"
18 #include "../cuda/ceed-cuda-compile.h"
19 #include "ceed-cuda-shared.h"
20 
21 //------------------------------------------------------------------------------
22 // Apply tensor basis
23 //------------------------------------------------------------------------------
24 static int CeedBasisApplyTensorCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode,
25                                                 CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
26   Ceed                   ceed;
27   Ceed_Cuda             *ceed_Cuda;
28   CeedInt                dim, num_comp;
29   const CeedScalar      *d_u;
30   CeedScalar            *d_v;
31   CeedBasis_Cuda_shared *data;
32 
33   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
34   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
35   CeedCallBackend(CeedBasisGetData(basis, &data));
36   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
37   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
38 
39   // Get read/write access to u, v
40   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
41   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
42   if (apply_add) {
43     CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
44   } else {
45     CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
46   }
47 
48   // Apply basis operation
49   switch (eval_mode) {
50     case CEED_EVAL_INTERP: {
51       CeedInt P_1d, Q_1d;
52 
53       CeedCheck(data->d_interp_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; interp_1d not set", CeedEvalModes[eval_mode]);
54       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
55       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
56       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
57 
58       void *interp_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_u, &d_v};
59 
60       if (dim == 1) {
61         // avoid >512 total threads
62         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1));
63         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
64         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
65 
66         if (t_mode == CEED_TRANSPOSE) {
67           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, 1,
68                                                       elems_per_block, shared_mem, interp_args));
69         } else {
70           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
71         }
72       } else if (dim == 2) {
73         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
74         // elems_per_block must be at least 1
75         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
76         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
77         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
78 
79         if (t_mode == CEED_TRANSPOSE) {
80           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, thread_1d,
81                                                       elems_per_block, shared_mem, interp_args));
82         } else {
83           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
84         }
85       } else if (dim == 3) {
86         CeedInt elems_per_block = 1;
87         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
88         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
89 
90         if (t_mode == CEED_TRANSPOSE) {
91           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, thread_1d,
92                                                       elems_per_block, shared_mem, interp_args));
93         } else {
94           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
95         }
96       }
97     } break;
98     case CEED_EVAL_GRAD: {
99       CeedInt P_1d, Q_1d;
100 
101       CeedCheck(data->d_grad_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; grad_1d not set", CeedEvalModes[eval_mode]);
102       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
103       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
104       CeedInt     thread_1d = CeedIntMax(Q_1d, P_1d);
105       CeedScalar *d_grad_1d = data->d_grad_1d;
106 
107       if (data->d_collo_grad_1d) {
108         d_grad_1d = data->d_collo_grad_1d;
109       }
110       void *grad_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_grad_1d, &d_u, &d_v};
111 
112       if (dim == 1) {
113         // avoid >512 total threads
114         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1));
115         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
116         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
117 
118         if (t_mode == CEED_TRANSPOSE) {
119           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, 1,
120                                                       elems_per_block, shared_mem, grad_args));
121         } else {
122           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
123         }
124       } else if (dim == 2) {
125         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
126         // elems_per_block must be at least 1
127         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
128         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
129         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
130 
131         if (t_mode == CEED_TRANSPOSE) {
132           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, thread_1d,
133                                                       elems_per_block, shared_mem, grad_args));
134         } else {
135           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
136         }
137       } else if (dim == 3) {
138         CeedInt elems_per_block = 1;
139         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
140         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
141 
142         if (t_mode == CEED_TRANSPOSE) {
143           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, thread_1d,
144                                                       elems_per_block, shared_mem, grad_args));
145         } else {
146           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
147         }
148       }
149     } break;
150     case CEED_EVAL_WEIGHT: {
151       CeedInt Q_1d;
152       CeedInt block_size = 32;
153 
154       CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights_1d not set", CeedEvalModes[eval_mode]);
155       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
156       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
157       if (dim == 1) {
158         const CeedInt elems_per_block = block_size / Q_1d;
159         const CeedInt grid_size       = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
160 
161         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, elems_per_block, 1, weight_args));
162       } else if (dim == 2) {
163         const CeedInt opt_elems       = block_size / (Q_1d * Q_1d);
164         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
165         const CeedInt grid_size       = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
166 
167         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
168       } else if (dim == 3) {
169         const CeedInt opt_elems       = block_size / (Q_1d * Q_1d);
170         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
171         const CeedInt grid_size       = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
172 
173         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
174       }
175     } break;
176     case CEED_EVAL_NONE: /* handled separately below */
177       break;
178     // LCOV_EXCL_START
179     case CEED_EVAL_DIV:
180     case CEED_EVAL_CURL:
181       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
182       // LCOV_EXCL_STOP
183   }
184 
185   // Restore vectors, cover CEED_EVAL_NONE
186   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
187   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
188   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
189   CeedCallBackend(CeedDestroy(&ceed));
190   return CEED_ERROR_SUCCESS;
191 }
192 
193 static int CeedBasisApplyTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
194                                             CeedVector v) {
195   CeedCallBackend(CeedBasisApplyTensorCore_Cuda_shared(basis, false, num_elem, t_mode, eval_mode, u, v));
196   return CEED_ERROR_SUCCESS;
197 }
198 
199 static int CeedBasisApplyAddTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
200                                                CeedVector u, CeedVector v) {
201   CeedCallBackend(CeedBasisApplyTensorCore_Cuda_shared(basis, true, num_elem, t_mode, eval_mode, u, v));
202   return CEED_ERROR_SUCCESS;
203 }
204 
205 //------------------------------------------------------------------------------
206 // Basis apply - tensor AtPoints
207 //------------------------------------------------------------------------------
208 static int CeedBasisApplyAtPointsCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, const CeedInt *num_points,
209                                                   CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
210   Ceed                   ceed;
211   Ceed_Cuda             *ceed_Cuda;
212   CeedInt                Q_1d, dim, num_comp, max_num_points = num_points[0];
213   const CeedInt          is_transpose = t_mode == CEED_TRANSPOSE;
214   const CeedScalar      *d_x, *d_u;
215   CeedScalar            *d_v;
216   CeedBasis_Cuda_shared *data;
217 
218   CeedCallBackend(CeedBasisGetData(basis, &data));
219   CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
220   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
221   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
222 
223   // Weight handled separately
224   if (eval_mode == CEED_EVAL_WEIGHT) {
225     CeedCallBackend(CeedVectorSetValue(v, 1.0));
226     return CEED_ERROR_SUCCESS;
227   }
228 
229   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
230   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
231 
232   // Check padded to uniform number of points per elem
233   for (CeedInt i = 1; i < num_elem; i++) max_num_points = CeedIntMax(max_num_points, num_points[i]);
234   {
235     CeedInt  q_comp;
236     CeedSize len, len_required;
237     CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
238     CeedCallBackend(CeedVectorGetLength(is_transpose ? u : v, &len));
239     len_required = (CeedSize)num_comp * (CeedSize)q_comp * (CeedSize)num_elem * (CeedSize)max_num_points;
240     CeedCheck(len >= len_required, ceed, CEED_ERROR_BACKEND,
241               "Vector at points must be padded to the same number of points in each element for BasisApplyAtPoints on GPU backends."
242               " Found %" CeedSize_FMT ", Required %" CeedSize_FMT,
243               len, len_required);
244   }
245 
246   // Move num_points array to device
247   if (is_transpose) {
248     const CeedInt num_bytes = num_elem * sizeof(CeedInt);
249 
250     if (num_elem != data->num_elem_at_points) {
251       data->num_elem_at_points = num_elem;
252 
253       if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem));
254       CeedCallCuda(ceed, cudaMalloc((void **)&data->d_points_per_elem, num_bytes));
255       CeedCallBackend(CeedFree(&data->h_points_per_elem));
256       CeedCallBackend(CeedCalloc(num_elem, &data->h_points_per_elem));
257     }
258     if (memcmp(data->h_points_per_elem, num_points, num_bytes)) {
259       memcpy(data->h_points_per_elem, num_points, num_bytes);
260       CeedCallCuda(ceed, cudaMemcpy(data->d_points_per_elem, num_points, num_bytes, cudaMemcpyHostToDevice));
261     }
262   }
263 
264   // Build kernels if needed
265   if (data->num_points != max_num_points) {
266     CeedInt P_1d;
267 
268     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
269     data->num_points = max_num_points;
270 
271     // -- Create interp matrix to Chebyshev coefficients
272     if (!data->d_chebyshev_interp_1d) {
273       CeedSize    interp_bytes;
274       CeedScalar *chebyshev_interp_1d;
275 
276       interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
277       CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
278       CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
279       CeedCallCuda(ceed, cudaMalloc((void **)&data->d_chebyshev_interp_1d, interp_bytes));
280       CeedCallCuda(ceed, cudaMemcpy(data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
281       CeedCallBackend(CeedFree(&chebyshev_interp_1d));
282     }
283 
284     // -- Compile kernels
285     const char basis_kernel_source[] = "// AtPoints basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-at-points.h>\n";
286     CeedInt    num_comp;
287 
288     if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
289     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
290     CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->moduleAtPoints, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D",
291                                      CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim),
292                                      "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_NUM_PTS", max_num_points));
293     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpAtPoints", &data->InterpAtPoints));
294     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpTransposeAtPoints", &data->InterpTransposeAtPoints));
295     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradAtPoints", &data->GradAtPoints));
296     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradTransposeAtPoints", &data->GradTransposeAtPoints));
297   }
298 
299   // Get read/write access to u, v
300   CeedCallBackend(CeedVectorGetArrayRead(x_ref, CEED_MEM_DEVICE, &d_x));
301   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
302   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
303   if (apply_add) {
304     CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
305   } else {
306     // Clear v for transpose operation
307     if (is_transpose) CeedCallBackend(CeedVectorSetValue(v, 0.0));
308     CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
309   }
310 
311   // Basis action
312   switch (eval_mode) {
313     case CEED_EVAL_INTERP: {
314       CeedInt P_1d, Q_1d;
315 
316       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
317       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
318       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
319 
320       void *interp_args[] = {(void *)&num_elem, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v};
321 
322       if (dim == 1) {
323         // avoid >512 total threads
324         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1));
325         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
326         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
327 
328         CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d, 1,
329                                                     elems_per_block, shared_mem, interp_args));
330       } else if (dim == 2) {
331         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
332         // elems_per_block must be at least 1
333         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
334         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
335         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
336 
337         CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d,
338                                                     thread_1d, elems_per_block, shared_mem, interp_args));
339       } else if (dim == 3) {
340         CeedInt elems_per_block = 1;
341         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
342         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
343 
344         CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d,
345                                                     thread_1d, elems_per_block, shared_mem, interp_args));
346       }
347     } break;
348     case CEED_EVAL_GRAD: {
349       CeedInt P_1d, Q_1d;
350 
351       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
352       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
353       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
354 
355       void *grad_args[] = {(void *)&num_elem, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v};
356 
357       if (dim == 1) {
358         // avoid >512 total threads
359         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1));
360         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
361         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
362 
363         CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, 1,
364                                                     elems_per_block, shared_mem, grad_args));
365       } else if (dim == 2) {
366         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
367         // elems_per_block must be at least 1
368         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
369         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
370         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
371 
372         CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, thread_1d,
373                                                     elems_per_block, shared_mem, grad_args));
374       } else if (dim == 3) {
375         CeedInt elems_per_block = 1;
376         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
377         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
378 
379         CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, thread_1d,
380                                                     elems_per_block, shared_mem, grad_args));
381       }
382     } break;
383     case CEED_EVAL_WEIGHT:
384     case CEED_EVAL_NONE: /* handled separately below */
385       break;
386     // LCOV_EXCL_START
387     case CEED_EVAL_DIV:
388     case CEED_EVAL_CURL:
389       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
390       // LCOV_EXCL_STOP
391   }
392 
393   // Restore vectors, cover CEED_EVAL_NONE
394   CeedCallBackend(CeedVectorRestoreArrayRead(x_ref, &d_x));
395   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
396   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
397   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
398   CeedCallBackend(CeedDestroy(&ceed));
399   return CEED_ERROR_SUCCESS;
400 }
401 
402 static int CeedBasisApplyAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
403                                               CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
404   CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda_shared(basis, false, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
405   return CEED_ERROR_SUCCESS;
406 }
407 
408 static int CeedBasisApplyAddAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
409                                                  CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
410   CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda_shared(basis, true, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
411   return CEED_ERROR_SUCCESS;
412 }
413 
414 //------------------------------------------------------------------------------
415 // Apply non-tensor basis
416 //------------------------------------------------------------------------------
417 static int CeedBasisApplyNonTensorCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode,
418                                                    CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
419   Ceed                   ceed;
420   Ceed_Cuda             *ceed_Cuda;
421   CeedInt                dim;
422   const CeedScalar      *d_u;
423   CeedScalar            *d_v;
424   CeedBasis_Cuda_shared *data;
425 
426   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
427   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
428   CeedCallBackend(CeedBasisGetData(basis, &data));
429   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
430 
431   // Get read/write access to u, v
432   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
433   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
434   if (apply_add) {
435     CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
436   } else {
437     CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
438   }
439 
440   // Apply basis operation
441   switch (eval_mode) {
442     case CEED_EVAL_INTERP: {
443       CeedInt P, Q;
444 
445       CeedCheck(data->d_interp_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; interp not set", CeedEvalModes[eval_mode]);
446       CeedCallBackend(CeedBasisGetNumNodes(basis, &P));
447       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q));
448       CeedInt thread = CeedIntMax(Q, P);
449 
450       void *interp_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_u, &d_v};
451 
452       {
453         // avoid >512 total threads
454         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1));
455         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
456         CeedInt shared_mem      = elems_per_block * thread * sizeof(CeedScalar);
457 
458         if (t_mode == CEED_TRANSPOSE) {
459           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread, 1,
460                                                       elems_per_block, shared_mem, interp_args));
461         } else {
462           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread, 1, elems_per_block, shared_mem, interp_args));
463         }
464       }
465     } break;
466     case CEED_EVAL_GRAD: {
467       CeedInt P, Q;
468 
469       CeedCheck(data->d_grad_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; grad not set", CeedEvalModes[eval_mode]);
470       CeedCallBackend(CeedBasisGetNumNodes(basis, &P));
471       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q));
472       CeedInt thread = CeedIntMax(Q, P);
473 
474       void *grad_args[] = {(void *)&num_elem, &data->d_grad_1d, &d_u, &d_v};
475 
476       {
477         // avoid >512 total threads
478         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1));
479         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
480         CeedInt shared_mem      = elems_per_block * thread * sizeof(CeedScalar);
481 
482         if (t_mode == CEED_TRANSPOSE) {
483           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread, 1,
484                                                       elems_per_block, shared_mem, grad_args));
485         } else {
486           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread, 1, elems_per_block, shared_mem, grad_args));
487         }
488       }
489     } break;
490     case CEED_EVAL_WEIGHT: {
491       CeedInt P, Q;
492 
493       CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights not set", CeedEvalModes[eval_mode]);
494       CeedCallBackend(CeedBasisGetNumNodes(basis, &P));
495       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q));
496       CeedInt thread = CeedIntMax(Q, P);
497 
498       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
499 
500       {
501         // avoid >512 total threads
502         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread, 1));
503         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
504 
505         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, thread, elems_per_block, 1, weight_args));
506       }
507     } break;
508     case CEED_EVAL_NONE: /* handled separately below */
509       break;
510     // LCOV_EXCL_START
511     case CEED_EVAL_DIV:
512     case CEED_EVAL_CURL:
513       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
514       // LCOV_EXCL_STOP
515   }
516 
517   // Restore vectors, cover CEED_EVAL_NONE
518   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
519   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
520   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
521   CeedCallBackend(CeedDestroy(&ceed));
522   return CEED_ERROR_SUCCESS;
523 }
524 
525 static int CeedBasisApplyNonTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
526                                                CeedVector u, CeedVector v) {
527   CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda_shared(basis, false, num_elem, t_mode, eval_mode, u, v));
528   return CEED_ERROR_SUCCESS;
529 }
530 
531 static int CeedBasisApplyAddNonTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
532                                                   CeedVector u, CeedVector v) {
533   CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda_shared(basis, true, num_elem, t_mode, eval_mode, u, v));
534   return CEED_ERROR_SUCCESS;
535 }
536 
537 //------------------------------------------------------------------------------
538 // Destroy basis
539 //------------------------------------------------------------------------------
540 static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) {
541   Ceed                   ceed;
542   CeedBasis_Cuda_shared *data;
543 
544   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
545   CeedCallBackend(CeedBasisGetData(basis, &data));
546   CeedCallCuda(ceed, cuModuleUnload(data->module));
547   if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
548   if (data->d_q_weight_1d) CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d));
549   CeedCallBackend(CeedFree(&data->h_points_per_elem));
550   if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem));
551   CeedCallCuda(ceed, cudaFree(data->d_interp_1d));
552   CeedCallCuda(ceed, cudaFree(data->d_grad_1d));
553   CeedCallCuda(ceed, cudaFree(data->d_collo_grad_1d));
554   CeedCallCuda(ceed, cudaFree(data->d_chebyshev_interp_1d));
555   CeedCallBackend(CeedFree(&data));
556   CeedCallBackend(CeedDestroy(&ceed));
557   return CEED_ERROR_SUCCESS;
558 }
559 
560 //------------------------------------------------------------------------------
561 // Create tensor basis
562 //------------------------------------------------------------------------------
563 int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
564                                         const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
565   Ceed                   ceed;
566   CeedInt                num_comp;
567   const CeedInt          q_bytes      = Q_1d * sizeof(CeedScalar);
568   const CeedInt          interp_bytes = q_bytes * P_1d;
569   CeedBasis_Cuda_shared *data;
570 
571   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
572   CeedCallBackend(CeedCalloc(1, &data));
573 
574   // Copy basis data to GPU
575   if (q_weight_1d) {
576     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
577     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice));
578   }
579   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
580   CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice));
581   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes));
582   CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice));
583 
584   // Compute collocated gradient and copy to GPU
585   data->d_collo_grad_1d    = NULL;
586   bool has_collocated_grad = dim == 3 && Q_1d >= P_1d;
587 
588   if (has_collocated_grad) {
589     CeedScalar *collo_grad_1d;
590 
591     CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d));
592     CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d));
593     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d));
594     CeedCallCuda(ceed, cudaMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, cudaMemcpyHostToDevice));
595     CeedCallBackend(CeedFree(&collo_grad_1d));
596   }
597 
598   // Compile basis kernels
599   const char basis_kernel_source[] = "// Tensor basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-tensor.h>\n";
600 
601   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
602   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D",
603                                    CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim),
604                                    "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_HAS_COLLOCATED_GRAD", has_collocated_grad));
605   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
606   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
607   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd));
608   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
609   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose));
610   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd));
611   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
612 
613   CeedCallBackend(CeedBasisSetData(basis, data));
614 
615   // Register backend functions
616   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared));
617   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddTensor_Cuda_shared));
618   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAtPoints", CeedBasisApplyAtPoints_Cuda_shared));
619   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAddAtPoints", CeedBasisApplyAddAtPoints_Cuda_shared));
620   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared));
621   CeedCallBackend(CeedDestroy(&ceed));
622   return CEED_ERROR_SUCCESS;
623 }
624 
625 //------------------------------------------------------------------------------
626 // Create non-tensor basis
627 //------------------------------------------------------------------------------
628 int CeedBasisCreateH1_Cuda_shared(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
629                                   const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
630   Ceed                   ceed;
631   CeedInt                num_comp, q_comp_interp, q_comp_grad;
632   const CeedInt          q_bytes = num_qpts * sizeof(CeedScalar);
633   CeedBasis_Cuda_shared *data;
634 
635   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
636 
637   // Check shared memory size
638   {
639     Ceed_Cuda *cuda_data;
640 
641     CeedCallBackend(CeedGetData(ceed, &cuda_data));
642     if (((size_t)num_nodes * (size_t)num_qpts * (size_t)dim + (size_t)CeedIntMax(num_nodes, num_qpts)) * sizeof(CeedScalar) >
643         cuda_data->device_prop.sharedMemPerBlock) {
644       CeedCallBackend(CeedBasisCreateH1Fallback(ceed, topo, dim, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis));
645       CeedCallBackend(CeedDestroy(&ceed));
646       return CEED_ERROR_SUCCESS;
647     }
648   }
649 
650   CeedCallBackend(CeedCalloc(1, &data));
651 
652   // Copy basis data to GPU
653   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
654   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad));
655   if (q_weight) {
656     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
657     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight, q_bytes, cudaMemcpyHostToDevice));
658   }
659   if (interp) {
660     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
661 
662     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
663     CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp, interp_bytes, cudaMemcpyHostToDevice));
664   }
665   if (grad) {
666     const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad;
667 
668     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, grad_bytes));
669     CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad, grad_bytes, cudaMemcpyHostToDevice));
670   }
671 
672   // Compile basis kernels
673   const char basis_kernel_source[] = "// Non-tensor basis source\n#include <ceed/jit-source/cuda/cuda-shared-basis-nontensor.h>\n";
674 
675   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
676   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "T_1D",
677                                    CeedIntMax(num_qpts, num_nodes), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp));
678   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
679   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
680   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd));
681   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
682   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose));
683   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd));
684   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
685 
686   CeedCallBackend(CeedBasisSetData(basis, data));
687 
688   // Register backend functions
689   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda_shared));
690   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Cuda_shared));
691   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared));
692   CeedCallBackend(CeedDestroy(&ceed));
693   return CEED_ERROR_SUCCESS;
694 }
695 
696 //------------------------------------------------------------------------------
697