xref: /libCEED/backends/cuda-ref/ceed-cuda-ref-basis.c (revision 3f08121c3514eef2c33780965530b5ca922abfd3)
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 
14 #include "../cuda/ceed-cuda-common.h"
15 #include "../cuda/ceed-cuda-compile.h"
16 #include "ceed-cuda-ref.h"
17 
18 //------------------------------------------------------------------------------
19 // Basis apply - tensor
20 //------------------------------------------------------------------------------
21 int CeedBasisApply_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
22   Ceed              ceed;
23   CeedInt           Q_1d, dim;
24   const CeedInt     is_transpose   = t_mode == CEED_TRANSPOSE;
25   const int         max_block_size = 32;
26   const CeedScalar *d_u;
27   CeedScalar       *d_v;
28   CeedBasis_Cuda   *data;
29 
30   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
31   CeedCallBackend(CeedBasisGetData(basis, &data));
32 
33   // Get read/write access to u, v
34   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
35   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
36   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
37 
38   // Clear v for transpose operation
39   if (is_transpose) {
40     CeedSize length;
41 
42     CeedCallBackend(CeedVectorGetLength(v, &length));
43     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
44   }
45   CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
46   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
47 
48   // Basis action
49   switch (eval_mode) {
50     case CEED_EVAL_INTERP: {
51       void         *interp_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_interp_1d, &d_u, &d_v};
52       const CeedInt block_size    = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
53 
54       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->Interp, num_elem, block_size, interp_args));
55     } break;
56     case CEED_EVAL_GRAD: {
57       void         *grad_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_interp_1d, &data->d_grad_1d, &d_u, &d_v};
58       const CeedInt block_size  = max_block_size;
59 
60       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->Grad, num_elem, block_size, grad_args));
61     } break;
62     case CEED_EVAL_WEIGHT: {
63       CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights_1d not set", CeedEvalModes[eval_mode]);
64       void     *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
65       const int block_size_x  = Q_1d;
66       const int block_size_y  = dim >= 2 ? Q_1d : 1;
67 
68       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, num_elem, block_size_x, block_size_y, 1, weight_args));
69     } break;
70     case CEED_EVAL_NONE: /* handled separately below */
71       break;
72     // LCOV_EXCL_START
73     case CEED_EVAL_DIV:
74     case CEED_EVAL_CURL:
75       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
76       // LCOV_EXCL_STOP
77   }
78 
79   // Restore vectors, cover CEED_EVAL_NONE
80   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
81   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
82   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
83   return CEED_ERROR_SUCCESS;
84 }
85 
86 //------------------------------------------------------------------------------
87 // Basis apply - tensor AtPoints
88 //------------------------------------------------------------------------------
89 int CeedBasisApplyAtPoints_Cuda(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
90                                 CeedVector x_ref, CeedVector u, CeedVector v) {
91   Ceed              ceed;
92   CeedInt           Q_1d, dim, max_num_points = num_points[0];
93   const CeedInt     is_transpose   = t_mode == CEED_TRANSPOSE;
94   const int         max_block_size = 32;
95   const CeedScalar *d_x, *d_u;
96   CeedScalar       *d_v;
97   CeedBasis_Cuda   *data;
98 
99   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
100   CeedCallBackend(CeedBasisGetData(basis, &data));
101   CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
102   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
103 
104   // Check uniform number of points per elem
105   for (CeedInt i = 1; i < num_elem; i++) {
106     CeedCheck(max_num_points == num_points[i], ceed, CEED_ERROR_BACKEND,
107               "BasisApplyAtPoints only supported for the same number of points in each element");
108   }
109 
110   // Weight handled separately
111   if (eval_mode == CEED_EVAL_WEIGHT) {
112     CeedCall(CeedVectorSetValue(v, 1.0));
113     return CEED_ERROR_SUCCESS;
114   }
115 
116   // Build kernels if needed
117   if (data->num_points != max_num_points) {
118     CeedInt P_1d;
119 
120     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
121     data->num_points = max_num_points;
122 
123     // -- Create interp matrix to Chebyshev coefficients
124     if (!data->d_chebyshev_interp_1d) {
125       CeedSize    interp_bytes;
126       CeedScalar *chebyshev_interp_1d;
127 
128       interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
129       CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
130       CeedCall(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
131       CeedCallCuda(ceed, cudaMalloc((void **)&data->d_chebyshev_interp_1d, interp_bytes));
132       CeedCallCuda(ceed, cudaMemcpy(data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
133       CeedCallBackend(CeedFree(&chebyshev_interp_1d));
134     }
135 
136     // -- Compile kernels
137     char       *basis_kernel_source;
138     const char *basis_kernel_path;
139     CeedInt     num_comp;
140 
141     if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
142     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
143     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-tensor-at-points.h", &basis_kernel_path));
144     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
145     CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
146     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
147     CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->moduleAtPoints, 9, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "BASIS_BUF_LEN",
148                                      Q_1d * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim - 1), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp,
149                                      "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_NUM_PTS",
150                                      max_num_points, "POINTS_BUFF_LEN", CeedIntPow(Q_1d, dim - 1)));
151     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpAtPoints", &data->InterpAtPoints));
152     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradAtPoints", &data->GradAtPoints));
153     CeedCallBackend(CeedFree(&basis_kernel_path));
154     CeedCallBackend(CeedFree(&basis_kernel_source));
155   }
156 
157   // Get read/write access to u, v
158   CeedCallBackend(CeedVectorGetArrayRead(x_ref, CEED_MEM_DEVICE, &d_x));
159   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
160   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
161   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
162 
163   // Clear v for transpose operation
164   if (is_transpose) {
165     CeedSize length;
166 
167     CeedCallBackend(CeedVectorGetLength(v, &length));
168     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
169   }
170 
171   // Basis action
172   switch (eval_mode) {
173     case CEED_EVAL_INTERP: {
174       void         *interp_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_chebyshev_interp_1d, &d_x, &d_u, &d_v};
175       const CeedInt block_size    = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
176 
177       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->InterpAtPoints, num_elem, block_size, interp_args));
178     } break;
179     case CEED_EVAL_GRAD: {
180       void         *grad_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_chebyshev_interp_1d, &d_x, &d_u, &d_v};
181       const CeedInt block_size  = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
182 
183       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->GradAtPoints, num_elem, block_size, grad_args));
184     } break;
185     case CEED_EVAL_WEIGHT:
186     case CEED_EVAL_NONE: /* handled separately below */
187       break;
188     // LCOV_EXCL_START
189     case CEED_EVAL_DIV:
190     case CEED_EVAL_CURL:
191       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
192       // LCOV_EXCL_STOP
193   }
194 
195   // Restore vectors, cover CEED_EVAL_NONE
196   CeedCallBackend(CeedVectorRestoreArrayRead(x_ref, &d_x));
197   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
198   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
199   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
200   return CEED_ERROR_SUCCESS;
201 }
202 
203 //------------------------------------------------------------------------------
204 // Basis apply - non-tensor
205 //------------------------------------------------------------------------------
206 int CeedBasisApplyNonTensor_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
207                                  CeedVector v) {
208   Ceed                     ceed;
209   CeedInt                  num_nodes, num_qpts;
210   const CeedInt            is_transpose    = t_mode == CEED_TRANSPOSE;
211   const int                elems_per_block = 1;
212   const int                grid            = CeedDivUpInt(num_elem, elems_per_block);
213   const CeedScalar        *d_u;
214   CeedScalar              *d_v;
215   CeedBasisNonTensor_Cuda *data;
216 
217   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
218   CeedCallBackend(CeedBasisGetData(basis, &data));
219   CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
220   CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes));
221 
222   // Get read/write access to u, v
223   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
224   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
225   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
226 
227   // Clear v for transpose operation
228   if (is_transpose) {
229     CeedSize length;
230 
231     CeedCallBackend(CeedVectorGetLength(v, &length));
232     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
233   }
234 
235   // Apply basis operation
236   switch (eval_mode) {
237     case CEED_EVAL_INTERP: {
238       void     *interp_args[] = {(void *)&num_elem, &data->d_interp, &d_u, &d_v};
239       const int block_size_x  = is_transpose ? num_nodes : num_qpts;
240 
241       if (is_transpose) {
242         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->InterpTranspose, grid, block_size_x, 1, elems_per_block, interp_args));
243       } else {
244         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Interp, grid, block_size_x, 1, elems_per_block, interp_args));
245       }
246     } break;
247     case CEED_EVAL_GRAD: {
248       void     *grad_args[]  = {(void *)&num_elem, &data->d_grad, &d_u, &d_v};
249       const int block_size_x = is_transpose ? num_nodes : num_qpts;
250 
251       if (is_transpose) {
252         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, grad_args));
253       } else {
254         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, grad_args));
255       }
256     } break;
257     case CEED_EVAL_DIV: {
258       void     *div_args[]   = {(void *)&num_elem, &data->d_div, &d_u, &d_v};
259       const int block_size_x = is_transpose ? num_nodes : num_qpts;
260 
261       if (is_transpose) {
262         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, div_args));
263       } else {
264         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, div_args));
265       }
266     } break;
267     case CEED_EVAL_CURL: {
268       void     *curl_args[]  = {(void *)&num_elem, &data->d_curl, &d_u, &d_v};
269       const int block_size_x = is_transpose ? num_nodes : num_qpts;
270 
271       if (is_transpose) {
272         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, curl_args));
273       } else {
274         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, curl_args));
275       }
276     } break;
277     case CEED_EVAL_WEIGHT: {
278       CeedCheck(data->d_q_weight, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights not set", CeedEvalModes[eval_mode]);
279       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight, &d_v};
280 
281       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, num_qpts, 1, elems_per_block, weight_args));
282     } break;
283     case CEED_EVAL_NONE: /* handled separately below */
284       break;
285   }
286 
287   // Restore vectors, cover CEED_EVAL_NONE
288   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
289   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
290   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
291   return CEED_ERROR_SUCCESS;
292 }
293 
294 //------------------------------------------------------------------------------
295 // Destroy tensor basis
296 //------------------------------------------------------------------------------
297 static int CeedBasisDestroy_Cuda(CeedBasis basis) {
298   Ceed            ceed;
299   CeedBasis_Cuda *data;
300 
301   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
302   CeedCallBackend(CeedBasisGetData(basis, &data));
303   CeedCallCuda(ceed, cuModuleUnload(data->module));
304   if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
305   if (data->d_q_weight_1d) CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d));
306   CeedCallCuda(ceed, cudaFree(data->d_interp_1d));
307   CeedCallCuda(ceed, cudaFree(data->d_grad_1d));
308   CeedCallCuda(ceed, cudaFree(data->d_chebyshev_interp_1d));
309   CeedCallBackend(CeedFree(&data));
310   return CEED_ERROR_SUCCESS;
311 }
312 
313 //------------------------------------------------------------------------------
314 // Destroy non-tensor basis
315 //------------------------------------------------------------------------------
316 static int CeedBasisDestroyNonTensor_Cuda(CeedBasis basis) {
317   Ceed                     ceed;
318   CeedBasisNonTensor_Cuda *data;
319 
320   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
321   CeedCallBackend(CeedBasisGetData(basis, &data));
322   CeedCallCuda(ceed, cuModuleUnload(data->module));
323   if (data->d_q_weight) CeedCallCuda(ceed, cudaFree(data->d_q_weight));
324   CeedCallCuda(ceed, cudaFree(data->d_interp));
325   CeedCallCuda(ceed, cudaFree(data->d_grad));
326   CeedCallCuda(ceed, cudaFree(data->d_div));
327   CeedCallCuda(ceed, cudaFree(data->d_curl));
328   CeedCallBackend(CeedFree(&data));
329   return CEED_ERROR_SUCCESS;
330 }
331 
332 //------------------------------------------------------------------------------
333 // Create tensor
334 //------------------------------------------------------------------------------
335 int CeedBasisCreateTensorH1_Cuda(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
336                                  const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
337   Ceed            ceed;
338   char           *basis_kernel_source;
339   const char     *basis_kernel_path;
340   CeedInt         num_comp;
341   const CeedInt   q_bytes      = Q_1d * sizeof(CeedScalar);
342   const CeedInt   interp_bytes = q_bytes * P_1d;
343   CeedBasis_Cuda *data;
344 
345   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
346   CeedCallBackend(CeedCalloc(1, &data));
347 
348   // Copy data to GPU
349   if (q_weight_1d) {
350     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
351     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice));
352   }
353   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
354   CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice));
355   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes));
356   CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice));
357 
358   // Compile basis kernels
359   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
360   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-tensor.h", &basis_kernel_path));
361   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
362   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
363   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
364   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 7, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "BASIS_BUF_LEN",
365                                    Q_1d * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim - 1), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp,
366                                    "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim)));
367   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
368   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
369   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
370   CeedCallBackend(CeedFree(&basis_kernel_path));
371   CeedCallBackend(CeedFree(&basis_kernel_source));
372 
373   CeedCallBackend(CeedBasisSetData(basis, data));
374 
375   // Register backend functions
376   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApply_Cuda));
377   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAtPoints", CeedBasisApplyAtPoints_Cuda));
378   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda));
379   return CEED_ERROR_SUCCESS;
380 }
381 
382 //------------------------------------------------------------------------------
383 // Create non-tensor H^1
384 //------------------------------------------------------------------------------
385 int CeedBasisCreateH1_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *grad,
386                            const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
387   Ceed                     ceed;
388   char                    *basis_kernel_source;
389   const char              *basis_kernel_path;
390   CeedInt                  num_comp, q_comp_interp, q_comp_grad;
391   const CeedInt            q_bytes = num_qpts * sizeof(CeedScalar);
392   CeedBasisNonTensor_Cuda *data;
393 
394   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
395   CeedCallBackend(CeedCalloc(1, &data));
396 
397   // Copy basis data to GPU
398   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
399   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad));
400   if (q_weight) {
401     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes));
402     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice));
403   }
404   if (interp) {
405     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
406 
407     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes));
408     CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice));
409   }
410   if (grad) {
411     const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad;
412 
413     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad, grad_bytes));
414     CeedCallCuda(ceed, cudaMemcpy(data->d_grad, grad, grad_bytes, cudaMemcpyHostToDevice));
415   }
416 
417   // Compile basis kernels
418   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
419   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path));
420   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
421   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
422   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
423   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP",
424                                    q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_grad, "BASIS_NUM_COMP", num_comp));
425   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
426   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
427   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv));
428   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose));
429   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
430   CeedCallBackend(CeedFree(&basis_kernel_path));
431   CeedCallBackend(CeedFree(&basis_kernel_source));
432 
433   CeedCallBackend(CeedBasisSetData(basis, data));
434 
435   // Register backend functions
436   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda));
437   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda));
438   return CEED_ERROR_SUCCESS;
439 }
440 
441 //------------------------------------------------------------------------------
442 // Create non-tensor H(div)
443 //------------------------------------------------------------------------------
444 int CeedBasisCreateHdiv_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *div,
445                              const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
446   Ceed                     ceed;
447   char                    *basis_kernel_source;
448   const char              *basis_kernel_path;
449   CeedInt                  num_comp, q_comp_interp, q_comp_div;
450   const CeedInt            q_bytes = num_qpts * sizeof(CeedScalar);
451   CeedBasisNonTensor_Cuda *data;
452 
453   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
454   CeedCallBackend(CeedCalloc(1, &data));
455 
456   // Copy basis data to GPU
457   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
458   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div));
459   if (q_weight) {
460     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes));
461     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice));
462   }
463   if (interp) {
464     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
465 
466     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes));
467     CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice));
468   }
469   if (div) {
470     const CeedInt div_bytes = q_bytes * num_nodes * q_comp_div;
471 
472     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_div, div_bytes));
473     CeedCallCuda(ceed, cudaMemcpy(data->d_div, div, div_bytes, cudaMemcpyHostToDevice));
474   }
475 
476   // Compile basis kernels
477   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
478   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path));
479   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
480   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
481   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
482   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP",
483                                    q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_div, "BASIS_NUM_COMP", num_comp));
484   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
485   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
486   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv));
487   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose));
488   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
489   CeedCallBackend(CeedFree(&basis_kernel_path));
490   CeedCallBackend(CeedFree(&basis_kernel_source));
491 
492   CeedCallBackend(CeedBasisSetData(basis, data));
493 
494   // Register backend functions
495   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda));
496   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda));
497   return CEED_ERROR_SUCCESS;
498 }
499 
500 //------------------------------------------------------------------------------
501 // Create non-tensor H(curl)
502 //------------------------------------------------------------------------------
503 int CeedBasisCreateHcurl_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
504                               const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
505   Ceed                     ceed;
506   char                    *basis_kernel_source;
507   const char              *basis_kernel_path;
508   CeedInt                  num_comp, q_comp_interp, q_comp_curl;
509   const CeedInt            q_bytes = num_qpts * sizeof(CeedScalar);
510   CeedBasisNonTensor_Cuda *data;
511 
512   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
513   CeedCallBackend(CeedCalloc(1, &data));
514 
515   // Copy basis data to GPU
516   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
517   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl));
518   if (q_weight) {
519     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes));
520     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice));
521   }
522   if (interp) {
523     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
524 
525     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes));
526     CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice));
527   }
528   if (curl) {
529     const CeedInt curl_bytes = q_bytes * num_nodes * q_comp_curl;
530 
531     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_curl, curl_bytes));
532     CeedCallCuda(ceed, cudaMemcpy(data->d_curl, curl, curl_bytes, cudaMemcpyHostToDevice));
533   }
534 
535   // Compile basis kernels
536   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
537   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path));
538   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
539   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
540   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
541   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP",
542                                    q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_curl, "BASIS_NUM_COMP", num_comp));
543   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
544   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
545   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv));
546   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose));
547   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
548   CeedCallBackend(CeedFree(&basis_kernel_path));
549   CeedCallBackend(CeedFree(&basis_kernel_source));
550 
551   CeedCallBackend(CeedBasisSetData(basis, data));
552 
553   // Register backend functions
554   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda));
555   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda));
556   return CEED_ERROR_SUCCESS;
557 }
558 
559 //------------------------------------------------------------------------------
560