xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-ref/ceed-cuda-ref-basis.c (revision f7c9815f6bdd2aac16e9c7ed574a2a61404669c3)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
30d0321e0SJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
50d0321e0SJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
70d0321e0SJeremy L Thompson 
849aac155SJeremy L Thompson #include <ceed.h>
90d0321e0SJeremy L Thompson #include <ceed/backend.h>
10437930d1SJeremy L Thompson #include <ceed/jit-tools.h>
110d0321e0SJeremy L Thompson #include <cuda.h>
120d0321e0SJeremy L Thompson #include <cuda_runtime.h>
132b730f8bSJeremy L Thompson 
1449aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
150d0321e0SJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
162b730f8bSJeremy L Thompson #include "ceed-cuda-ref.h"
170d0321e0SJeremy L Thompson 
180d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
190d0321e0SJeremy L Thompson // Basis apply - tensor
200d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
212b730f8bSJeremy L Thompson int CeedBasisApply_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
220d0321e0SJeremy L Thompson   Ceed              ceed;
23ca735530SJeremy L Thompson   CeedInt           Q_1d, dim;
247bbbfca3SJeremy L Thompson   const CeedInt     is_transpose   = t_mode == CEED_TRANSPOSE;
25437930d1SJeremy L Thompson   const int         max_block_size = 32;
260d0321e0SJeremy L Thompson   const CeedScalar *d_u;
270d0321e0SJeremy L Thompson   CeedScalar       *d_v;
28ca735530SJeremy L Thompson   CeedBasis_Cuda   *data;
29ca735530SJeremy L Thompson 
30ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
31ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
32ca735530SJeremy L Thompson 
339ea2cfd9SJeremy L Thompson   // Get read/write access to u, v
346574a04fSJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
356574a04fSJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
362b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
370d0321e0SJeremy L Thompson 
380d0321e0SJeremy L Thompson   // Clear v for transpose operation
397bbbfca3SJeremy L Thompson   if (is_transpose) {
401f9221feSJeremy L Thompson     CeedSize length;
41ca735530SJeremy L Thompson 
422b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(v, &length));
432b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
440d0321e0SJeremy L Thompson   }
452b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
462b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
470d0321e0SJeremy L Thompson 
480d0321e0SJeremy L Thompson   // Basis action
49437930d1SJeremy L Thompson   switch (eval_mode) {
500d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP: {
517bbbfca3SJeremy L Thompson       void         *interp_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_interp_1d, &d_u, &d_v};
52b2165e7aSSebastian Grimberg       const CeedInt block_size    = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
530d0321e0SJeremy L Thompson 
54eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->Interp, num_elem, block_size, interp_args));
550d0321e0SJeremy L Thompson     } break;
560d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD: {
577bbbfca3SJeremy L Thompson       void         *grad_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_interp_1d, &data->d_grad_1d, &d_u, &d_v};
58b2165e7aSSebastian Grimberg       const CeedInt block_size  = max_block_size;
590d0321e0SJeremy L Thompson 
60eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->Grad, num_elem, block_size, grad_args));
610d0321e0SJeremy L Thompson     } break;
620d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
63097cc795SJames Wright       CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights_1d not set", CeedEvalModes[eval_mode]);
64437930d1SJeremy L Thompson       void     *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
65b2165e7aSSebastian Grimberg       const int block_size_x  = Q_1d;
66b2165e7aSSebastian Grimberg       const int block_size_y  = dim >= 2 ? Q_1d : 1;
67b2165e7aSSebastian Grimberg 
68b2165e7aSSebastian Grimberg       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, num_elem, block_size_x, block_size_y, 1, weight_args));
690d0321e0SJeremy L Thompson     } break;
709ea2cfd9SJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
719ea2cfd9SJeremy L Thompson       break;
720d0321e0SJeremy L Thompson     // LCOV_EXCL_START
730d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
740d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
75bcbe1c99SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
760d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
770d0321e0SJeremy L Thompson   }
780d0321e0SJeremy L Thompson 
799ea2cfd9SJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
802b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
819ea2cfd9SJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
829ea2cfd9SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
830d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
840d0321e0SJeremy L Thompson }
850d0321e0SJeremy L Thompson 
860d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8734d14614SJeremy L Thompson // Basis apply - tensor AtPoints
8834d14614SJeremy L Thompson //------------------------------------------------------------------------------
8934d14614SJeremy L Thompson int CeedBasisApplyAtPoints_Cuda(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
9034d14614SJeremy L Thompson                                 CeedVector x_ref, CeedVector u, CeedVector v) {
9134d14614SJeremy L Thompson   Ceed              ceed;
9234d14614SJeremy L Thompson   CeedInt           Q_1d, dim, max_num_points = num_points[0];
9334d14614SJeremy L Thompson   const CeedInt     is_transpose   = t_mode == CEED_TRANSPOSE;
9434d14614SJeremy L Thompson   const int         max_block_size = 32;
9534d14614SJeremy L Thompson   const CeedScalar *d_x, *d_u;
9634d14614SJeremy L Thompson   CeedScalar       *d_v;
9734d14614SJeremy L Thompson   CeedBasis_Cuda   *data;
9834d14614SJeremy L Thompson 
9934d14614SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
10034d14614SJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
10134d14614SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
10234d14614SJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
10334d14614SJeremy L Thompson 
10434d14614SJeremy L Thompson   // Check uniform number of points per elem
10534d14614SJeremy L Thompson   for (CeedInt i = 1; i < num_elem; i++) {
10634d14614SJeremy L Thompson     CeedCheck(max_num_points == num_points[i], ceed, CEED_ERROR_BACKEND,
10734d14614SJeremy L Thompson               "BasisApplyAtPoints only supported for the same number of points in each element");
10834d14614SJeremy L Thompson   }
10934d14614SJeremy L Thompson 
11034d14614SJeremy L Thompson   // Weight handled separately
11134d14614SJeremy L Thompson   if (eval_mode == CEED_EVAL_WEIGHT) {
11234d14614SJeremy L Thompson     CeedCall(CeedVectorSetValue(v, 1.0));
11334d14614SJeremy L Thompson     return CEED_ERROR_SUCCESS;
11434d14614SJeremy L Thompson   }
11534d14614SJeremy L Thompson 
11634d14614SJeremy L Thompson   // Build kernels if needed
11734d14614SJeremy L Thompson   if (data->num_points != max_num_points) {
11834d14614SJeremy L Thompson     CeedInt P_1d;
11934d14614SJeremy L Thompson 
12034d14614SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
12134d14614SJeremy L Thompson     data->num_points = max_num_points;
12234d14614SJeremy L Thompson 
12334d14614SJeremy L Thompson     // -- Create interp matrix to Chebyshev coefficients
12434d14614SJeremy L Thompson     if (!data->d_chebyshev_interp_1d) {
12534d14614SJeremy L Thompson       CeedSize    interp_bytes;
12634d14614SJeremy L Thompson       CeedScalar *chebyshev_interp_1d;
12734d14614SJeremy L Thompson 
12834d14614SJeremy L Thompson       interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
12934d14614SJeremy L Thompson       CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
13034d14614SJeremy L Thompson       CeedCall(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
13134d14614SJeremy L Thompson       CeedCallCuda(ceed, cudaMalloc((void **)&data->d_chebyshev_interp_1d, interp_bytes));
13234d14614SJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy(data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
13334d14614SJeremy L Thompson       CeedCallBackend(CeedFree(&chebyshev_interp_1d));
13434d14614SJeremy L Thompson     }
13534d14614SJeremy L Thompson 
13634d14614SJeremy L Thompson     // -- Compile kernels
13734d14614SJeremy L Thompson     char       *basis_kernel_source;
13834d14614SJeremy L Thompson     const char *basis_kernel_path;
13934d14614SJeremy L Thompson     CeedInt     num_comp;
14034d14614SJeremy L Thompson 
14134d14614SJeremy L Thompson     if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
14234d14614SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
14334d14614SJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-tensor-at-points.h", &basis_kernel_path));
14434d14614SJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
14534d14614SJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
14634d14614SJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
14734d14614SJeremy L Thompson     CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->moduleAtPoints, 9, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "BASIS_BUF_LEN",
148*f7c9815fSJeremy L Thompson                                      Q_1d * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim - 1), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp,
14934d14614SJeremy L Thompson                                      "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_NUM_PTS",
150*f7c9815fSJeremy L Thompson                                      max_num_points, "POINTS_BUFF_LEN", CeedIntPow(Q_1d, dim - 1)));
15134d14614SJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpAtPoints", &data->InterpAtPoints));
15234d14614SJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradAtPoints", &data->GradAtPoints));
15334d14614SJeremy L Thompson     CeedCallBackend(CeedFree(&basis_kernel_path));
15434d14614SJeremy L Thompson     CeedCallBackend(CeedFree(&basis_kernel_source));
15534d14614SJeremy L Thompson   }
15634d14614SJeremy L Thompson 
15734d14614SJeremy L Thompson   // Get read/write access to u, v
15834d14614SJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayRead(x_ref, CEED_MEM_DEVICE, &d_x));
15934d14614SJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
16034d14614SJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
16134d14614SJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
16234d14614SJeremy L Thompson 
16334d14614SJeremy L Thompson   // Clear v for transpose operation
16434d14614SJeremy L Thompson   if (is_transpose) {
16534d14614SJeremy L Thompson     CeedSize length;
16634d14614SJeremy L Thompson 
16734d14614SJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(v, &length));
16834d14614SJeremy L Thompson     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
16934d14614SJeremy L Thompson   }
17034d14614SJeremy L Thompson 
17134d14614SJeremy L Thompson   // Basis action
17234d14614SJeremy L Thompson   switch (eval_mode) {
17334d14614SJeremy L Thompson     case CEED_EVAL_INTERP: {
17434d14614SJeremy L Thompson       void         *interp_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_chebyshev_interp_1d, &d_x, &d_u, &d_v};
17534d14614SJeremy L Thompson       const CeedInt block_size    = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
17634d14614SJeremy L Thompson 
17734d14614SJeremy L Thompson       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->InterpAtPoints, num_elem, block_size, interp_args));
17834d14614SJeremy L Thompson     } break;
17934d14614SJeremy L Thompson     case CEED_EVAL_GRAD: {
18034d14614SJeremy L Thompson       void         *grad_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_chebyshev_interp_1d, &d_x, &d_u, &d_v};
1812d10e82cSJeremy L Thompson       const CeedInt block_size  = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
18234d14614SJeremy L Thompson 
18334d14614SJeremy L Thompson       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->GradAtPoints, num_elem, block_size, grad_args));
18434d14614SJeremy L Thompson     } break;
18534d14614SJeremy L Thompson     case CEED_EVAL_WEIGHT:
18634d14614SJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
18734d14614SJeremy L Thompson       break;
18834d14614SJeremy L Thompson     // LCOV_EXCL_START
18934d14614SJeremy L Thompson     case CEED_EVAL_DIV:
19034d14614SJeremy L Thompson     case CEED_EVAL_CURL:
19134d14614SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
19234d14614SJeremy L Thompson       // LCOV_EXCL_STOP
19334d14614SJeremy L Thompson   }
19434d14614SJeremy L Thompson 
19534d14614SJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
19634d14614SJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArrayRead(x_ref, &d_x));
19734d14614SJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
19834d14614SJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
19934d14614SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
20034d14614SJeremy L Thompson   return CEED_ERROR_SUCCESS;
20134d14614SJeremy L Thompson }
20234d14614SJeremy L Thompson 
20334d14614SJeremy L Thompson //------------------------------------------------------------------------------
2040d0321e0SJeremy L Thompson // Basis apply - non-tensor
2050d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2062b730f8bSJeremy L Thompson int CeedBasisApplyNonTensor_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
2072b730f8bSJeremy L Thompson                                  CeedVector v) {
2080d0321e0SJeremy L Thompson   Ceed                     ceed;
209437930d1SJeremy L Thompson   CeedInt                  num_nodes, num_qpts;
2107bbbfca3SJeremy L Thompson   const CeedInt            is_transpose    = t_mode == CEED_TRANSPOSE;
211d075f50bSSebastian Grimberg   const int                elems_per_block = 1;
212d075f50bSSebastian Grimberg   const int                grid            = CeedDivUpInt(num_elem, elems_per_block);
2130d0321e0SJeremy L Thompson   const CeedScalar        *d_u;
2140d0321e0SJeremy L Thompson   CeedScalar              *d_v;
215ca735530SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
216ca735530SJeremy L Thompson 
217ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
218ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
219ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
220ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes));
221ca735530SJeremy L Thompson 
2229ea2cfd9SJeremy L Thompson   // Get read/write access to u, v
2239ea2cfd9SJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
2249ea2cfd9SJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
2252b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
2260d0321e0SJeremy L Thompson 
2270d0321e0SJeremy L Thompson   // Clear v for transpose operation
2287bbbfca3SJeremy L Thompson   if (is_transpose) {
2291f9221feSJeremy L Thompson     CeedSize length;
230ca735530SJeremy L Thompson 
2312b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(v, &length));
2322b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
2330d0321e0SJeremy L Thompson   }
2340d0321e0SJeremy L Thompson 
2350d0321e0SJeremy L Thompson   // Apply basis operation
236437930d1SJeremy L Thompson   switch (eval_mode) {
2370d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP: {
238d075f50bSSebastian Grimberg       void     *interp_args[] = {(void *)&num_elem, &data->d_interp, &d_u, &d_v};
2397bbbfca3SJeremy L Thompson       const int block_size_x  = is_transpose ? num_nodes : num_qpts;
240b2165e7aSSebastian Grimberg 
2417bbbfca3SJeremy L Thompson       if (is_transpose) {
242d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->InterpTranspose, grid, block_size_x, 1, elems_per_block, interp_args));
243d075f50bSSebastian Grimberg       } else {
244b2165e7aSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Interp, grid, block_size_x, 1, elems_per_block, interp_args));
245d075f50bSSebastian Grimberg       }
2460d0321e0SJeremy L Thompson     } break;
2470d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD: {
248d075f50bSSebastian Grimberg       void     *grad_args[]  = {(void *)&num_elem, &data->d_grad, &d_u, &d_v};
2497bbbfca3SJeremy L Thompson       const int block_size_x = is_transpose ? num_nodes : num_qpts;
250b2165e7aSSebastian Grimberg 
2517bbbfca3SJeremy L Thompson       if (is_transpose) {
252d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, grad_args));
253d075f50bSSebastian Grimberg       } else {
254d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, grad_args));
255d075f50bSSebastian Grimberg       }
256d075f50bSSebastian Grimberg     } break;
257d075f50bSSebastian Grimberg     case CEED_EVAL_DIV: {
258d075f50bSSebastian Grimberg       void     *div_args[]   = {(void *)&num_elem, &data->d_div, &d_u, &d_v};
2597bbbfca3SJeremy L Thompson       const int block_size_x = is_transpose ? num_nodes : num_qpts;
260d075f50bSSebastian Grimberg 
2617bbbfca3SJeremy L Thompson       if (is_transpose) {
262d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, div_args));
263d075f50bSSebastian Grimberg       } else {
264d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, div_args));
265d075f50bSSebastian Grimberg       }
266d075f50bSSebastian Grimberg     } break;
267d075f50bSSebastian Grimberg     case CEED_EVAL_CURL: {
268d075f50bSSebastian Grimberg       void     *curl_args[]  = {(void *)&num_elem, &data->d_curl, &d_u, &d_v};
2697bbbfca3SJeremy L Thompson       const int block_size_x = is_transpose ? num_nodes : num_qpts;
270d075f50bSSebastian Grimberg 
2717bbbfca3SJeremy L Thompson       if (is_transpose) {
272d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, curl_args));
273d075f50bSSebastian Grimberg       } else {
274d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, curl_args));
275d075f50bSSebastian Grimberg       }
2760d0321e0SJeremy L Thompson     } break;
2770d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
278097cc795SJames Wright       CeedCheck(data->d_q_weight, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights not set", CeedEvalModes[eval_mode]);
279437930d1SJeremy L Thompson       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight, &d_v};
280b2165e7aSSebastian Grimberg 
281eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, num_qpts, 1, elems_per_block, weight_args));
2820d0321e0SJeremy L Thompson     } break;
2839ea2cfd9SJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
2849ea2cfd9SJeremy L Thompson       break;
2850d0321e0SJeremy L Thompson   }
2860d0321e0SJeremy L Thompson 
2879ea2cfd9SJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
2882b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
2899ea2cfd9SJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
2909ea2cfd9SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
2910d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2920d0321e0SJeremy L Thompson }
2930d0321e0SJeremy L Thompson 
2940d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2950d0321e0SJeremy L Thompson // Destroy tensor basis
2960d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2970d0321e0SJeremy L Thompson static int CeedBasisDestroy_Cuda(CeedBasis basis) {
2980d0321e0SJeremy L Thompson   Ceed            ceed;
2990d0321e0SJeremy L Thompson   CeedBasis_Cuda *data;
300ca735530SJeremy L Thompson 
301ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
3022b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
3032b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
30434d14614SJeremy L Thompson   if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
305097cc795SJames Wright   if (data->d_q_weight_1d) CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d));
3062b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp_1d));
3072b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad_1d));
30834d14614SJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_chebyshev_interp_1d));
3092b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
3100d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3110d0321e0SJeremy L Thompson }
3120d0321e0SJeremy L Thompson 
3130d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3140d0321e0SJeremy L Thompson // Destroy non-tensor basis
3150d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3160d0321e0SJeremy L Thompson static int CeedBasisDestroyNonTensor_Cuda(CeedBasis basis) {
3170d0321e0SJeremy L Thompson   Ceed                     ceed;
3180d0321e0SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
319ca735530SJeremy L Thompson 
320ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
3212b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
3222b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
323097cc795SJames Wright   if (data->d_q_weight) CeedCallCuda(ceed, cudaFree(data->d_q_weight));
3242b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp));
3252b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad));
326d075f50bSSebastian Grimberg   CeedCallCuda(ceed, cudaFree(data->d_div));
327d075f50bSSebastian Grimberg   CeedCallCuda(ceed, cudaFree(data->d_curl));
3282b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
3290d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3300d0321e0SJeremy L Thompson }
3310d0321e0SJeremy L Thompson 
3320d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3330d0321e0SJeremy L Thompson // Create tensor
3340d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3352b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
3362b730f8bSJeremy L Thompson                                  const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
3370d0321e0SJeremy L Thompson   Ceed            ceed;
33822070f95SJeremy L Thompson   char           *basis_kernel_source;
33922070f95SJeremy L Thompson   const char     *basis_kernel_path;
340ca735530SJeremy L Thompson   CeedInt         num_comp;
341ca735530SJeremy L Thompson   const CeedInt   q_bytes      = Q_1d * sizeof(CeedScalar);
342ca735530SJeremy L Thompson   const CeedInt   interp_bytes = q_bytes * P_1d;
3430d0321e0SJeremy L Thompson   CeedBasis_Cuda *data;
344ca735530SJeremy L Thompson 
345ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
3462b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
3470d0321e0SJeremy L Thompson 
3480d0321e0SJeremy L Thompson   // Copy data to GPU
349097cc795SJames Wright   if (q_weight_1d) {
3502b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
3512b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice));
352097cc795SJames Wright   }
3532b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
3542b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice));
3552b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes));
3562b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice));
3570d0321e0SJeremy L Thompson 
358ecc88aebSJeremy L Thompson   // Compile basis kernels
3592b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
3602b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-tensor.h", &basis_kernel_path));
36123d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
3622b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
36323d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
364eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 7, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "BASIS_BUF_LEN",
365*f7c9815fSJeremy L Thompson                                    Q_1d * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim - 1), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp,
3662b730f8bSJeremy L Thompson                                    "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim)));
367eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
368eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
369eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
3702b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_path));
3712b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_source));
372437930d1SJeremy L Thompson 
3732b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
3740d0321e0SJeremy L Thompson 
375d075f50bSSebastian Grimberg   // Register backend functions
3762b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApply_Cuda));
37734d14614SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAtPoints", CeedBasisApplyAtPoints_Cuda));
3782b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda));
3790d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3800d0321e0SJeremy L Thompson }
3810d0321e0SJeremy L Thompson 
3820d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
383d075f50bSSebastian Grimberg // Create non-tensor H^1
3840d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3852b730f8bSJeremy L Thompson int CeedBasisCreateH1_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *grad,
38651475c7cSJeremy L Thompson                            const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
3870d0321e0SJeremy L Thompson   Ceed                     ceed;
38822070f95SJeremy L Thompson   char                    *basis_kernel_source;
38922070f95SJeremy L Thompson   const char              *basis_kernel_path;
390d075f50bSSebastian Grimberg   CeedInt                  num_comp, q_comp_interp, q_comp_grad;
391ca735530SJeremy L Thompson   const CeedInt            q_bytes = num_qpts * sizeof(CeedScalar);
3920d0321e0SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
393ca735530SJeremy L Thompson 
394ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
3952b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
3960d0321e0SJeremy L Thompson 
3970d0321e0SJeremy L Thompson   // Copy basis data to GPU
398d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
399d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad));
400097cc795SJames Wright   if (q_weight) {
4012b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes));
4022b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice));
403097cc795SJames Wright   }
404d075f50bSSebastian Grimberg   if (interp) {
405d075f50bSSebastian Grimberg     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
406d075f50bSSebastian Grimberg 
4072b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes));
4082b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice));
409d075f50bSSebastian Grimberg   }
410d075f50bSSebastian Grimberg   if (grad) {
411d075f50bSSebastian Grimberg     const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad;
412d075f50bSSebastian Grimberg 
4132b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad, grad_bytes));
4142b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_grad, grad, grad_bytes, cudaMemcpyHostToDevice));
415d075f50bSSebastian Grimberg   }
4160d0321e0SJeremy L Thompson 
4170d0321e0SJeremy L Thompson   // Compile basis kernels
4182b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
4192b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path));
42023d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
4212b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
42223d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
423d075f50bSSebastian Grimberg   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP",
424d075f50bSSebastian Grimberg                                    q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_grad, "BASIS_NUM_COMP", num_comp));
425d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
426d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
427d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv));
428d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose));
429d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
430d075f50bSSebastian Grimberg   CeedCallBackend(CeedFree(&basis_kernel_path));
431d075f50bSSebastian Grimberg   CeedCallBackend(CeedFree(&basis_kernel_source));
432d075f50bSSebastian Grimberg 
433d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisSetData(basis, data));
434d075f50bSSebastian Grimberg 
435d075f50bSSebastian Grimberg   // Register backend functions
436d075f50bSSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda));
437d075f50bSSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda));
438d075f50bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
439d075f50bSSebastian Grimberg }
440d075f50bSSebastian Grimberg 
441d075f50bSSebastian Grimberg //------------------------------------------------------------------------------
442d075f50bSSebastian Grimberg // Create non-tensor H(div)
443d075f50bSSebastian Grimberg //------------------------------------------------------------------------------
444d075f50bSSebastian Grimberg int CeedBasisCreateHdiv_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *div,
445d075f50bSSebastian Grimberg                              const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
446d075f50bSSebastian Grimberg   Ceed                     ceed;
44722070f95SJeremy L Thompson   char                    *basis_kernel_source;
44822070f95SJeremy L Thompson   const char              *basis_kernel_path;
449d075f50bSSebastian Grimberg   CeedInt                  num_comp, q_comp_interp, q_comp_div;
450d075f50bSSebastian Grimberg   const CeedInt            q_bytes = num_qpts * sizeof(CeedScalar);
451d075f50bSSebastian Grimberg   CeedBasisNonTensor_Cuda *data;
452d075f50bSSebastian Grimberg 
453d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
454d075f50bSSebastian Grimberg   CeedCallBackend(CeedCalloc(1, &data));
455d075f50bSSebastian Grimberg 
456d075f50bSSebastian Grimberg   // Copy basis data to GPU
457d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
458d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div));
459097cc795SJames Wright   if (q_weight) {
460d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes));
461d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice));
462097cc795SJames Wright   }
463d075f50bSSebastian Grimberg   if (interp) {
464d075f50bSSebastian Grimberg     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
465d075f50bSSebastian Grimberg 
466d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes));
467d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice));
468d075f50bSSebastian Grimberg   }
469d075f50bSSebastian Grimberg   if (div) {
470d075f50bSSebastian Grimberg     const CeedInt div_bytes = q_bytes * num_nodes * q_comp_div;
471d075f50bSSebastian Grimberg 
472d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_div, div_bytes));
473d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_div, div, div_bytes, cudaMemcpyHostToDevice));
474d075f50bSSebastian Grimberg   }
475d075f50bSSebastian Grimberg 
476d075f50bSSebastian Grimberg   // Compile basis kernels
477d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
478d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path));
479d075f50bSSebastian Grimberg   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
480d075f50bSSebastian Grimberg   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
481d075f50bSSebastian Grimberg   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
482d075f50bSSebastian Grimberg   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP",
483d075f50bSSebastian Grimberg                                    q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_div, "BASIS_NUM_COMP", num_comp));
484d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
485d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
486d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv));
487d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose));
488d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
489d075f50bSSebastian Grimberg   CeedCallBackend(CeedFree(&basis_kernel_path));
490d075f50bSSebastian Grimberg   CeedCallBackend(CeedFree(&basis_kernel_source));
491d075f50bSSebastian Grimberg 
492d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisSetData(basis, data));
493d075f50bSSebastian Grimberg 
494d075f50bSSebastian Grimberg   // Register backend functions
495d075f50bSSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda));
496d075f50bSSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda));
497d075f50bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
498d075f50bSSebastian Grimberg }
499d075f50bSSebastian Grimberg 
500d075f50bSSebastian Grimberg //------------------------------------------------------------------------------
501d075f50bSSebastian Grimberg // Create non-tensor H(curl)
502d075f50bSSebastian Grimberg //------------------------------------------------------------------------------
503d075f50bSSebastian Grimberg int CeedBasisCreateHcurl_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
504d075f50bSSebastian Grimberg                               const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
505d075f50bSSebastian Grimberg   Ceed                     ceed;
50622070f95SJeremy L Thompson   char                    *basis_kernel_source;
50722070f95SJeremy L Thompson   const char              *basis_kernel_path;
508d075f50bSSebastian Grimberg   CeedInt                  num_comp, q_comp_interp, q_comp_curl;
509d075f50bSSebastian Grimberg   const CeedInt            q_bytes = num_qpts * sizeof(CeedScalar);
510d075f50bSSebastian Grimberg   CeedBasisNonTensor_Cuda *data;
511d075f50bSSebastian Grimberg 
512d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
513d075f50bSSebastian Grimberg   CeedCallBackend(CeedCalloc(1, &data));
514d075f50bSSebastian Grimberg 
515d075f50bSSebastian Grimberg   // Copy basis data to GPU
516d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
517d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl));
518097cc795SJames Wright   if (q_weight) {
519d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes));
520d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice));
521097cc795SJames Wright   }
522d075f50bSSebastian Grimberg   if (interp) {
523d075f50bSSebastian Grimberg     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
524d075f50bSSebastian Grimberg 
525d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes));
526d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice));
527d075f50bSSebastian Grimberg   }
528d075f50bSSebastian Grimberg   if (curl) {
529d075f50bSSebastian Grimberg     const CeedInt curl_bytes = q_bytes * num_nodes * q_comp_curl;
530d075f50bSSebastian Grimberg 
531d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_curl, curl_bytes));
532d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_curl, curl, curl_bytes, cudaMemcpyHostToDevice));
533d075f50bSSebastian Grimberg   }
534d075f50bSSebastian Grimberg 
535d075f50bSSebastian Grimberg   // Compile basis kernels
536d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
537d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path));
538d075f50bSSebastian Grimberg   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
539d075f50bSSebastian Grimberg   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
540d075f50bSSebastian Grimberg   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
541d075f50bSSebastian Grimberg   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP",
542d075f50bSSebastian Grimberg                                    q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_curl, "BASIS_NUM_COMP", num_comp));
543d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
544d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
545d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv));
546d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose));
547d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
5482b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_path));
5492b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_source));
5500d0321e0SJeremy L Thompson 
5512b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
5520d0321e0SJeremy L Thompson 
5530d0321e0SJeremy L Thompson   // Register backend functions
5542b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda));
5552b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda));
5560d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5570d0321e0SJeremy L Thompson }
5582a86cc9dSSebastian Grimberg 
5590d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
560