xref: /libCEED/backends/cuda-ref/ceed-cuda-ref-basis.c (revision 9bc663991d6482bcb1d60b1f116148f11db83fa1)
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>
13111870feSJeremy L Thompson #include <string.h>
142b730f8bSJeremy L Thompson 
1549aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
160d0321e0SJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
172b730f8bSJeremy L Thompson #include "ceed-cuda-ref.h"
180d0321e0SJeremy L Thompson 
190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
200d0321e0SJeremy L Thompson // Basis apply - tensor
210d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
22db2becc9SJeremy L Thompson static int CeedBasisApplyCore_Cuda(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
23db2becc9SJeremy L Thompson                                    CeedVector u, CeedVector v) {
240d0321e0SJeremy L Thompson   Ceed              ceed;
25ca735530SJeremy L Thompson   CeedInt           Q_1d, dim;
267bbbfca3SJeremy L Thompson   const CeedInt     is_transpose   = t_mode == CEED_TRANSPOSE;
27437930d1SJeremy L Thompson   const int         max_block_size = 32;
280d0321e0SJeremy L Thompson   const CeedScalar *d_u;
290d0321e0SJeremy L Thompson   CeedScalar       *d_v;
30ca735530SJeremy L Thompson   CeedBasis_Cuda   *data;
31ca735530SJeremy L Thompson 
32ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
33ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
34ca735530SJeremy L Thompson 
359ea2cfd9SJeremy L Thompson   // Get read/write access to u, v
366574a04fSJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
376574a04fSJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
38db2becc9SJeremy L Thompson   if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
39db2becc9SJeremy L Thompson   else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
400d0321e0SJeremy L Thompson 
410d0321e0SJeremy L Thompson   // Clear v for transpose operation
42db2becc9SJeremy L Thompson   if (is_transpose && !apply_add) {
4319a04db8SJeremy L Thompson     CeedInt  num_comp, q_comp, num_nodes, num_qpts;
441f9221feSJeremy L Thompson     CeedSize length;
45ca735530SJeremy L Thompson 
4619a04db8SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
4719a04db8SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
4819a04db8SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes));
4919a04db8SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
5019a04db8SJeremy L Thompson     length = (CeedSize)num_elem * (CeedSize)num_comp * (t_mode == CEED_TRANSPOSE ? (CeedSize)num_nodes : ((CeedSize)num_qpts * (CeedSize)q_comp));
512b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
520d0321e0SJeremy L Thompson   }
532b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
542b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
550d0321e0SJeremy L Thompson 
560d0321e0SJeremy L Thompson   // Basis action
57437930d1SJeremy L Thompson   switch (eval_mode) {
580d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP: {
597bbbfca3SJeremy L Thompson       void         *interp_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_interp_1d, &d_u, &d_v};
60b2165e7aSSebastian Grimberg       const CeedInt block_size    = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
610d0321e0SJeremy L Thompson 
62eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->Interp, num_elem, block_size, interp_args));
630d0321e0SJeremy L Thompson     } break;
640d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD: {
657bbbfca3SJeremy L Thompson       void         *grad_args[] = {(void *)&num_elem, (void *)&is_transpose, &data->d_interp_1d, &data->d_grad_1d, &d_u, &d_v};
66b2165e7aSSebastian Grimberg       const CeedInt block_size  = max_block_size;
670d0321e0SJeremy L Thompson 
68eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->Grad, num_elem, block_size, grad_args));
690d0321e0SJeremy L Thompson     } break;
700d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
71097cc795SJames Wright       CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights_1d not set", CeedEvalModes[eval_mode]);
72437930d1SJeremy L Thompson       void     *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
73b2165e7aSSebastian Grimberg       const int block_size_x  = Q_1d;
74b2165e7aSSebastian Grimberg       const int block_size_y  = dim >= 2 ? Q_1d : 1;
75b2165e7aSSebastian Grimberg 
76b2165e7aSSebastian Grimberg       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, num_elem, block_size_x, block_size_y, 1, weight_args));
770d0321e0SJeremy L Thompson     } break;
789ea2cfd9SJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
799ea2cfd9SJeremy L Thompson       break;
800d0321e0SJeremy L Thompson     // LCOV_EXCL_START
810d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
820d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
83bcbe1c99SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
840d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
850d0321e0SJeremy L Thompson   }
860d0321e0SJeremy L Thompson 
879ea2cfd9SJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
882b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
899ea2cfd9SJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
909ea2cfd9SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
91*9bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
920d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
930d0321e0SJeremy L Thompson }
940d0321e0SJeremy L Thompson 
95db2becc9SJeremy L Thompson static int CeedBasisApply_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
96db2becc9SJeremy L Thompson                                CeedVector v) {
97db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyCore_Cuda(basis, false, num_elem, t_mode, eval_mode, u, v));
98db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
99db2becc9SJeremy L Thompson }
100db2becc9SJeremy L Thompson 
101db2becc9SJeremy L Thompson static int CeedBasisApplyAdd_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
102db2becc9SJeremy L Thompson                                   CeedVector v) {
103db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyCore_Cuda(basis, true, num_elem, t_mode, eval_mode, u, v));
104db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
105db2becc9SJeremy L Thompson }
106db2becc9SJeremy L Thompson 
1070d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
10834d14614SJeremy L Thompson // Basis apply - tensor AtPoints
10934d14614SJeremy L Thompson //------------------------------------------------------------------------------
110db2becc9SJeremy L Thompson static int CeedBasisApplyAtPointsCore_Cuda(CeedBasis basis, bool apply_add, const CeedInt num_elem, const CeedInt *num_points,
111db2becc9SJeremy L Thompson                                            CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
11234d14614SJeremy L Thompson   Ceed              ceed;
11334d14614SJeremy L Thompson   CeedInt           Q_1d, dim, max_num_points = num_points[0];
11434d14614SJeremy L Thompson   const CeedInt     is_transpose   = t_mode == CEED_TRANSPOSE;
11534d14614SJeremy L Thompson   const int         max_block_size = 32;
11634d14614SJeremy L Thompson   const CeedScalar *d_x, *d_u;
11734d14614SJeremy L Thompson   CeedScalar       *d_v;
11834d14614SJeremy L Thompson   CeedBasis_Cuda   *data;
11934d14614SJeremy L Thompson 
12034d14614SJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
12134d14614SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
12234d14614SJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
12334d14614SJeremy L Thompson 
12434d14614SJeremy L Thompson   // Weight handled separately
12534d14614SJeremy L Thompson   if (eval_mode == CEED_EVAL_WEIGHT) {
1265a5594ffSJeremy L Thompson     CeedCallBackend(CeedVectorSetValue(v, 1.0));
12734d14614SJeremy L Thompson     return CEED_ERROR_SUCCESS;
12834d14614SJeremy L Thompson   }
12934d14614SJeremy L Thompson 
130*9bc66399SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
131*9bc66399SJeremy L Thompson 
132111870feSJeremy L Thompson   // Check padded to uniform number of points per elem
133111870feSJeremy L Thompson   for (CeedInt i = 1; i < num_elem; i++) max_num_points = CeedIntMax(max_num_points, num_points[i]);
134111870feSJeremy L Thompson   {
135111870feSJeremy L Thompson     CeedInt  num_comp, q_comp;
136111870feSJeremy L Thompson     CeedSize len, len_required;
137111870feSJeremy L Thompson 
138111870feSJeremy L Thompson     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
139111870feSJeremy L Thompson     CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
140111870feSJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(is_transpose ? u : v, &len));
141111870feSJeremy L Thompson     len_required = (CeedSize)num_comp * (CeedSize)q_comp * (CeedSize)num_elem * (CeedSize)max_num_points;
142111870feSJeremy L Thompson     CeedCheck(len >= len_required, ceed, CEED_ERROR_BACKEND,
143111870feSJeremy L Thompson               "Vector at points must be padded to the same number of points in each element for BasisApplyAtPoints on GPU backends."
144111870feSJeremy L Thompson               " Found %" CeedSize_FMT ", Required %" CeedSize_FMT,
145111870feSJeremy L Thompson               len, len_required);
146111870feSJeremy L Thompson   }
147111870feSJeremy L Thompson 
148111870feSJeremy L Thompson   // Move num_points array to device
149111870feSJeremy L Thompson   if (is_transpose) {
150111870feSJeremy L Thompson     const CeedInt num_bytes = num_elem * sizeof(CeedInt);
151111870feSJeremy L Thompson 
152111870feSJeremy L Thompson     if (num_elem != data->num_elem_at_points) {
153111870feSJeremy L Thompson       data->num_elem_at_points = num_elem;
154111870feSJeremy L Thompson 
155111870feSJeremy L Thompson       if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem));
156111870feSJeremy L Thompson       CeedCallCuda(ceed, cudaMalloc((void **)&data->d_points_per_elem, num_bytes));
157111870feSJeremy L Thompson       CeedCallBackend(CeedFree(&data->h_points_per_elem));
158111870feSJeremy L Thompson       CeedCallBackend(CeedCalloc(num_elem, &data->h_points_per_elem));
159111870feSJeremy L Thompson     }
1609e511c80SJeremy L Thompson     if (memcmp(data->h_points_per_elem, num_points, num_bytes)) {
161111870feSJeremy L Thompson       memcpy(data->h_points_per_elem, num_points, num_bytes);
162111870feSJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy(data->d_points_per_elem, num_points, num_bytes, cudaMemcpyHostToDevice));
163111870feSJeremy L Thompson     }
164111870feSJeremy L Thompson   }
165111870feSJeremy L Thompson 
16634d14614SJeremy L Thompson   // Build kernels if needed
16734d14614SJeremy L Thompson   if (data->num_points != max_num_points) {
16834d14614SJeremy L Thompson     CeedInt P_1d;
16934d14614SJeremy L Thompson 
17034d14614SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
17134d14614SJeremy L Thompson     data->num_points = max_num_points;
17234d14614SJeremy L Thompson 
17334d14614SJeremy L Thompson     // -- Create interp matrix to Chebyshev coefficients
17434d14614SJeremy L Thompson     if (!data->d_chebyshev_interp_1d) {
17534d14614SJeremy L Thompson       CeedSize    interp_bytes;
17634d14614SJeremy L Thompson       CeedScalar *chebyshev_interp_1d;
17734d14614SJeremy L Thompson 
17834d14614SJeremy L Thompson       interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
17934d14614SJeremy L Thompson       CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
1805a5594ffSJeremy L Thompson       CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
18134d14614SJeremy L Thompson       CeedCallCuda(ceed, cudaMalloc((void **)&data->d_chebyshev_interp_1d, interp_bytes));
18234d14614SJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy(data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
18334d14614SJeremy L Thompson       CeedCallBackend(CeedFree(&chebyshev_interp_1d));
18434d14614SJeremy L Thompson     }
18534d14614SJeremy L Thompson 
18634d14614SJeremy L Thompson     // -- Compile kernels
1879c25dd66SJeremy L Thompson     const char basis_kernel_source[] = "// AtPoints basis source\n#include <ceed/jit-source/cuda/cuda-ref-basis-tensor-at-points.h>\n";
18834d14614SJeremy L Thompson     CeedInt    num_comp;
18934d14614SJeremy L Thompson 
19034d14614SJeremy L Thompson     if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
19134d14614SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
19234d14614SJeremy 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",
193f7c9815fSJeremy L Thompson                                      Q_1d * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim - 1), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp,
19434d14614SJeremy L Thompson                                      "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_NUM_PTS",
195f7c9815fSJeremy L Thompson                                      max_num_points, "POINTS_BUFF_LEN", CeedIntPow(Q_1d, dim - 1)));
19634d14614SJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpAtPoints", &data->InterpAtPoints));
19734d14614SJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradAtPoints", &data->GradAtPoints));
19834d14614SJeremy L Thompson   }
19934d14614SJeremy L Thompson 
20034d14614SJeremy L Thompson   // Get read/write access to u, v
20134d14614SJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayRead(x_ref, CEED_MEM_DEVICE, &d_x));
20234d14614SJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
20334d14614SJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
204db2becc9SJeremy L Thompson   if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
205db2becc9SJeremy L Thompson   else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
20634d14614SJeremy L Thompson 
20734d14614SJeremy L Thompson   // Clear v for transpose operation
208db2becc9SJeremy L Thompson   if (is_transpose && !apply_add) {
20919a04db8SJeremy L Thompson     CeedInt  num_comp, q_comp, num_nodes;
21034d14614SJeremy L Thompson     CeedSize length;
21134d14614SJeremy L Thompson 
21219a04db8SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
21319a04db8SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
21419a04db8SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes));
21519a04db8SJeremy L Thompson     length =
21619a04db8SJeremy L Thompson         (CeedSize)num_elem * (CeedSize)num_comp * (t_mode == CEED_TRANSPOSE ? (CeedSize)num_nodes : ((CeedSize)max_num_points * (CeedSize)q_comp));
21734d14614SJeremy L Thompson     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
21834d14614SJeremy L Thompson   }
21934d14614SJeremy L Thompson 
22034d14614SJeremy L Thompson   // Basis action
22134d14614SJeremy L Thompson   switch (eval_mode) {
22234d14614SJeremy L Thompson     case CEED_EVAL_INTERP: {
223111870feSJeremy L Thompson       void *interp_args[]      = {(void *)&num_elem, (void *)&is_transpose, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v};
22434d14614SJeremy L Thompson       const CeedInt block_size = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
22534d14614SJeremy L Thompson 
22634d14614SJeremy L Thompson       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->InterpAtPoints, num_elem, block_size, interp_args));
22734d14614SJeremy L Thompson     } break;
22834d14614SJeremy L Thompson     case CEED_EVAL_GRAD: {
229111870feSJeremy L Thompson       void *grad_args[]        = {(void *)&num_elem, (void *)&is_transpose, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v};
2302d10e82cSJeremy L Thompson       const CeedInt block_size = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
23134d14614SJeremy L Thompson 
23234d14614SJeremy L Thompson       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->GradAtPoints, num_elem, block_size, grad_args));
23334d14614SJeremy L Thompson     } break;
23434d14614SJeremy L Thompson     case CEED_EVAL_WEIGHT:
23534d14614SJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
23634d14614SJeremy L Thompson       break;
23734d14614SJeremy L Thompson     // LCOV_EXCL_START
23834d14614SJeremy L Thompson     case CEED_EVAL_DIV:
23934d14614SJeremy L Thompson     case CEED_EVAL_CURL:
24034d14614SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
24134d14614SJeremy L Thompson       // LCOV_EXCL_STOP
24234d14614SJeremy L Thompson   }
24334d14614SJeremy L Thompson 
24434d14614SJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
24534d14614SJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArrayRead(x_ref, &d_x));
24634d14614SJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
24734d14614SJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
24834d14614SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
249*9bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
25034d14614SJeremy L Thompson   return CEED_ERROR_SUCCESS;
25134d14614SJeremy L Thompson }
25234d14614SJeremy L Thompson 
253db2becc9SJeremy L Thompson static int CeedBasisApplyAtPoints_Cuda(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
254db2becc9SJeremy L Thompson                                        CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
255db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda(basis, false, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
256db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
257db2becc9SJeremy L Thompson }
258db2becc9SJeremy L Thompson 
259db2becc9SJeremy L Thompson static int CeedBasisApplyAddAtPoints_Cuda(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
260db2becc9SJeremy L Thompson                                           CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
261db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda(basis, true, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
262db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
263db2becc9SJeremy L Thompson }
264db2becc9SJeremy L Thompson 
26534d14614SJeremy L Thompson //------------------------------------------------------------------------------
2660d0321e0SJeremy L Thompson // Basis apply - non-tensor
2670d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
268db2becc9SJeremy L Thompson static int CeedBasisApplyNonTensorCore_Cuda(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
269db2becc9SJeremy L Thompson                                             CeedVector u, CeedVector v) {
2700d0321e0SJeremy L Thompson   Ceed                     ceed;
271437930d1SJeremy L Thompson   CeedInt                  num_nodes, num_qpts;
2727bbbfca3SJeremy L Thompson   const CeedInt            is_transpose    = t_mode == CEED_TRANSPOSE;
273d075f50bSSebastian Grimberg   const int                elems_per_block = 1;
274d075f50bSSebastian Grimberg   const int                grid            = CeedDivUpInt(num_elem, elems_per_block);
2750d0321e0SJeremy L Thompson   const CeedScalar        *d_u;
2760d0321e0SJeremy L Thompson   CeedScalar              *d_v;
277ca735530SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
278ca735530SJeremy L Thompson 
279ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
280ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
281ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
282ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes));
283ca735530SJeremy L Thompson 
2849ea2cfd9SJeremy L Thompson   // Get read/write access to u, v
2859ea2cfd9SJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
2869ea2cfd9SJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
287db2becc9SJeremy L Thompson   if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
288db2becc9SJeremy L Thompson   else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
2890d0321e0SJeremy L Thompson 
2900d0321e0SJeremy L Thompson   // Clear v for transpose operation
291db2becc9SJeremy L Thompson   if (is_transpose && !apply_add) {
29219a04db8SJeremy L Thompson     CeedInt  num_comp, q_comp;
2931f9221feSJeremy L Thompson     CeedSize length;
294ca735530SJeremy L Thompson 
29519a04db8SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
29619a04db8SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
29719a04db8SJeremy L Thompson     length = (CeedSize)num_elem * (CeedSize)num_comp * (t_mode == CEED_TRANSPOSE ? (CeedSize)num_nodes : ((CeedSize)num_qpts * (CeedSize)q_comp));
2982b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
2990d0321e0SJeremy L Thompson   }
3000d0321e0SJeremy L Thompson 
3010d0321e0SJeremy L Thompson   // Apply basis operation
302437930d1SJeremy L Thompson   switch (eval_mode) {
3030d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP: {
304d075f50bSSebastian Grimberg       void     *interp_args[] = {(void *)&num_elem, &data->d_interp, &d_u, &d_v};
3057bbbfca3SJeremy L Thompson       const int block_size_x  = is_transpose ? num_nodes : num_qpts;
306b2165e7aSSebastian Grimberg 
3077bbbfca3SJeremy L Thompson       if (is_transpose) {
308d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->InterpTranspose, grid, block_size_x, 1, elems_per_block, interp_args));
309d075f50bSSebastian Grimberg       } else {
310b2165e7aSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Interp, grid, block_size_x, 1, elems_per_block, interp_args));
311d075f50bSSebastian Grimberg       }
3120d0321e0SJeremy L Thompson     } break;
3130d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD: {
314d075f50bSSebastian Grimberg       void     *grad_args[]  = {(void *)&num_elem, &data->d_grad, &d_u, &d_v};
3157bbbfca3SJeremy L Thompson       const int block_size_x = is_transpose ? num_nodes : num_qpts;
316b2165e7aSSebastian Grimberg 
3177bbbfca3SJeremy L Thompson       if (is_transpose) {
318d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, grad_args));
319d075f50bSSebastian Grimberg       } else {
320d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, grad_args));
321d075f50bSSebastian Grimberg       }
322d075f50bSSebastian Grimberg     } break;
323d075f50bSSebastian Grimberg     case CEED_EVAL_DIV: {
324d075f50bSSebastian Grimberg       void     *div_args[]   = {(void *)&num_elem, &data->d_div, &d_u, &d_v};
3257bbbfca3SJeremy L Thompson       const int block_size_x = is_transpose ? num_nodes : num_qpts;
326d075f50bSSebastian Grimberg 
3277bbbfca3SJeremy L Thompson       if (is_transpose) {
328d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, div_args));
329d075f50bSSebastian Grimberg       } else {
330d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, div_args));
331d075f50bSSebastian Grimberg       }
332d075f50bSSebastian Grimberg     } break;
333d075f50bSSebastian Grimberg     case CEED_EVAL_CURL: {
334d075f50bSSebastian Grimberg       void     *curl_args[]  = {(void *)&num_elem, &data->d_curl, &d_u, &d_v};
3357bbbfca3SJeremy L Thompson       const int block_size_x = is_transpose ? num_nodes : num_qpts;
336d075f50bSSebastian Grimberg 
3377bbbfca3SJeremy L Thompson       if (is_transpose) {
338d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, curl_args));
339d075f50bSSebastian Grimberg       } else {
340d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, curl_args));
341d075f50bSSebastian Grimberg       }
3420d0321e0SJeremy L Thompson     } break;
3430d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
344097cc795SJames Wright       CeedCheck(data->d_q_weight, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights not set", CeedEvalModes[eval_mode]);
345437930d1SJeremy L Thompson       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight, &d_v};
346b2165e7aSSebastian Grimberg 
347eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, num_qpts, 1, elems_per_block, weight_args));
3480d0321e0SJeremy L Thompson     } break;
3499ea2cfd9SJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
3509ea2cfd9SJeremy L Thompson       break;
3510d0321e0SJeremy L Thompson   }
3520d0321e0SJeremy L Thompson 
3539ea2cfd9SJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
3542b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
3559ea2cfd9SJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
3569ea2cfd9SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
357*9bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
3580d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3590d0321e0SJeremy L Thompson }
3600d0321e0SJeremy L Thompson 
361db2becc9SJeremy L Thompson static int CeedBasisApplyNonTensor_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
362db2becc9SJeremy L Thompson                                         CeedVector v) {
363db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda(basis, false, num_elem, t_mode, eval_mode, u, v));
364db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
365db2becc9SJeremy L Thompson }
366db2becc9SJeremy L Thompson 
367db2becc9SJeremy L Thompson static int CeedBasisApplyAddNonTensor_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
368db2becc9SJeremy L Thompson                                            CeedVector v) {
369db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyNonTensorCore_Cuda(basis, true, num_elem, t_mode, eval_mode, u, v));
370db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
371db2becc9SJeremy L Thompson }
372db2becc9SJeremy L Thompson 
3730d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3740d0321e0SJeremy L Thompson // Destroy tensor basis
3750d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3760d0321e0SJeremy L Thompson static int CeedBasisDestroy_Cuda(CeedBasis basis) {
3770d0321e0SJeremy L Thompson   Ceed            ceed;
3780d0321e0SJeremy L Thompson   CeedBasis_Cuda *data;
379ca735530SJeremy L Thompson 
380ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
3812b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
3822b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
38334d14614SJeremy L Thompson   if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
384097cc795SJames Wright   if (data->d_q_weight_1d) CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d));
385111870feSJeremy L Thompson   CeedCallBackend(CeedFree(&data->h_points_per_elem));
386111870feSJeremy L Thompson   if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem));
3872b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp_1d));
3882b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad_1d));
38934d14614SJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_chebyshev_interp_1d));
3902b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
391*9bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
3920d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3930d0321e0SJeremy L Thompson }
3940d0321e0SJeremy L Thompson 
3950d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3960d0321e0SJeremy L Thompson // Destroy non-tensor basis
3970d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3980d0321e0SJeremy L Thompson static int CeedBasisDestroyNonTensor_Cuda(CeedBasis basis) {
3990d0321e0SJeremy L Thompson   Ceed                     ceed;
4000d0321e0SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
401ca735530SJeremy L Thompson 
402ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
4032b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
4042b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
405097cc795SJames Wright   if (data->d_q_weight) CeedCallCuda(ceed, cudaFree(data->d_q_weight));
4062b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp));
4072b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad));
408d075f50bSSebastian Grimberg   CeedCallCuda(ceed, cudaFree(data->d_div));
409d075f50bSSebastian Grimberg   CeedCallCuda(ceed, cudaFree(data->d_curl));
4102b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
411*9bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
4120d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4130d0321e0SJeremy L Thompson }
4140d0321e0SJeremy L Thompson 
4150d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
4160d0321e0SJeremy L Thompson // Create tensor
4170d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
4182b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
4192b730f8bSJeremy L Thompson                                  const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
4200d0321e0SJeremy L Thompson   Ceed            ceed;
421ca735530SJeremy L Thompson   CeedInt         num_comp;
422ca735530SJeremy L Thompson   const CeedInt   q_bytes      = Q_1d * sizeof(CeedScalar);
423ca735530SJeremy L Thompson   const CeedInt   interp_bytes = q_bytes * P_1d;
4240d0321e0SJeremy L Thompson   CeedBasis_Cuda *data;
425ca735530SJeremy L Thompson 
426ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
4272b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
4280d0321e0SJeremy L Thompson 
4290d0321e0SJeremy L Thompson   // Copy data to GPU
430097cc795SJames Wright   if (q_weight_1d) {
4312b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
4322b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice));
433097cc795SJames Wright   }
4342b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
4352b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice));
4362b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes));
4372b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice));
4380d0321e0SJeremy L Thompson 
439ecc88aebSJeremy L Thompson   // Compile basis kernels
4409c25dd66SJeremy L Thompson   const char basis_kernel_source[] = "// Tensor basis source\n#include <ceed/jit-source/cuda/cuda-ref-basis-tensor.h>\n";
4419c25dd66SJeremy L Thompson 
4422b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
443eb7e6cafSJeremy 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",
444f7c9815fSJeremy L Thompson                                    Q_1d * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim - 1), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp,
4452b730f8bSJeremy L Thompson                                    "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim)));
446eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
447eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
448eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
449437930d1SJeremy L Thompson 
4502b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
4510d0321e0SJeremy L Thompson 
452d075f50bSSebastian Grimberg   // Register backend functions
4532b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApply_Cuda));
454db2becc9SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAdd_Cuda));
45534d14614SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAtPoints", CeedBasisApplyAtPoints_Cuda));
456db2becc9SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAddAtPoints", CeedBasisApplyAddAtPoints_Cuda));
4572b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda));
458*9bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
4590d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4600d0321e0SJeremy L Thompson }
4610d0321e0SJeremy L Thompson 
4620d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
463d075f50bSSebastian Grimberg // Create non-tensor H^1
4640d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
4652b730f8bSJeremy L Thompson int CeedBasisCreateH1_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *grad,
46651475c7cSJeremy L Thompson                            const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
4670d0321e0SJeremy L Thompson   Ceed                     ceed;
468d075f50bSSebastian Grimberg   CeedInt                  num_comp, q_comp_interp, q_comp_grad;
469ca735530SJeremy L Thompson   const CeedInt            q_bytes = num_qpts * sizeof(CeedScalar);
4700d0321e0SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
471ca735530SJeremy L Thompson 
472ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
4732b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
4740d0321e0SJeremy L Thompson 
4750d0321e0SJeremy L Thompson   // Copy basis data to GPU
476d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
477d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad));
478097cc795SJames Wright   if (q_weight) {
4792b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes));
4802b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice));
481097cc795SJames Wright   }
482d075f50bSSebastian Grimberg   if (interp) {
483d075f50bSSebastian Grimberg     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
484d075f50bSSebastian Grimberg 
4852b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes));
4862b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice));
487d075f50bSSebastian Grimberg   }
488d075f50bSSebastian Grimberg   if (grad) {
489d075f50bSSebastian Grimberg     const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad;
490d075f50bSSebastian Grimberg 
4912b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad, grad_bytes));
4922b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_grad, grad, grad_bytes, cudaMemcpyHostToDevice));
493d075f50bSSebastian Grimberg   }
4940d0321e0SJeremy L Thompson 
4950d0321e0SJeremy L Thompson   // Compile basis kernels
4969c25dd66SJeremy L Thompson   const char basis_kernel_source[] = "// Nontensor basis source\n#include <ceed/jit-source/cuda/cuda-ref-basis-nontensor.h>\n";
4979c25dd66SJeremy L Thompson 
4982b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
499d075f50bSSebastian Grimberg   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP",
500d075f50bSSebastian Grimberg                                    q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_grad, "BASIS_NUM_COMP", num_comp));
501d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
502d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
503d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv));
504d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose));
505d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
506d075f50bSSebastian Grimberg 
507d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisSetData(basis, data));
508d075f50bSSebastian Grimberg 
509d075f50bSSebastian Grimberg   // Register backend functions
510d075f50bSSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda));
511db2becc9SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Cuda));
512d075f50bSSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda));
513*9bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
514d075f50bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
515d075f50bSSebastian Grimberg }
516d075f50bSSebastian Grimberg 
517d075f50bSSebastian Grimberg //------------------------------------------------------------------------------
518d075f50bSSebastian Grimberg // Create non-tensor H(div)
519d075f50bSSebastian Grimberg //------------------------------------------------------------------------------
520d075f50bSSebastian Grimberg int CeedBasisCreateHdiv_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *div,
521d075f50bSSebastian Grimberg                              const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
522d075f50bSSebastian Grimberg   Ceed                     ceed;
523d075f50bSSebastian Grimberg   CeedInt                  num_comp, q_comp_interp, q_comp_div;
524d075f50bSSebastian Grimberg   const CeedInt            q_bytes = num_qpts * sizeof(CeedScalar);
525d075f50bSSebastian Grimberg   CeedBasisNonTensor_Cuda *data;
526d075f50bSSebastian Grimberg 
527d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
528d075f50bSSebastian Grimberg   CeedCallBackend(CeedCalloc(1, &data));
529d075f50bSSebastian Grimberg 
530d075f50bSSebastian Grimberg   // Copy basis data to GPU
531d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
532d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div));
533097cc795SJames Wright   if (q_weight) {
534d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes));
535d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice));
536097cc795SJames Wright   }
537d075f50bSSebastian Grimberg   if (interp) {
538d075f50bSSebastian Grimberg     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
539d075f50bSSebastian Grimberg 
540d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes));
541d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice));
542d075f50bSSebastian Grimberg   }
543d075f50bSSebastian Grimberg   if (div) {
544d075f50bSSebastian Grimberg     const CeedInt div_bytes = q_bytes * num_nodes * q_comp_div;
545d075f50bSSebastian Grimberg 
546d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_div, div_bytes));
547d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_div, div, div_bytes, cudaMemcpyHostToDevice));
548d075f50bSSebastian Grimberg   }
549d075f50bSSebastian Grimberg 
550d075f50bSSebastian Grimberg   // Compile basis kernels
5519c25dd66SJeremy L Thompson   const char basis_kernel_source[] = "// Nontensor basis source\n#include <ceed/jit-source/cuda/cuda-ref-basis-nontensor.h>\n";
5529c25dd66SJeremy L Thompson 
553d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
554d075f50bSSebastian Grimberg   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP",
555d075f50bSSebastian Grimberg                                    q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_div, "BASIS_NUM_COMP", num_comp));
556d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
557d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
558d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv));
559d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose));
560d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
561d075f50bSSebastian Grimberg 
562d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisSetData(basis, data));
563d075f50bSSebastian Grimberg 
564d075f50bSSebastian Grimberg   // Register backend functions
565d075f50bSSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda));
566db2becc9SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Cuda));
567d075f50bSSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda));
568*9bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
569d075f50bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
570d075f50bSSebastian Grimberg }
571d075f50bSSebastian Grimberg 
572d075f50bSSebastian Grimberg //------------------------------------------------------------------------------
573d075f50bSSebastian Grimberg // Create non-tensor H(curl)
574d075f50bSSebastian Grimberg //------------------------------------------------------------------------------
575d075f50bSSebastian Grimberg int CeedBasisCreateHcurl_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
576d075f50bSSebastian Grimberg                               const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
577d075f50bSSebastian Grimberg   Ceed                     ceed;
578d075f50bSSebastian Grimberg   CeedInt                  num_comp, q_comp_interp, q_comp_curl;
579d075f50bSSebastian Grimberg   const CeedInt            q_bytes = num_qpts * sizeof(CeedScalar);
580d075f50bSSebastian Grimberg   CeedBasisNonTensor_Cuda *data;
581d075f50bSSebastian Grimberg 
582d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
583d075f50bSSebastian Grimberg   CeedCallBackend(CeedCalloc(1, &data));
584d075f50bSSebastian Grimberg 
585d075f50bSSebastian Grimberg   // Copy basis data to GPU
586d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
587d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl));
588097cc795SJames Wright   if (q_weight) {
589d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes));
590d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice));
591097cc795SJames Wright   }
592d075f50bSSebastian Grimberg   if (interp) {
593d075f50bSSebastian Grimberg     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
594d075f50bSSebastian Grimberg 
595d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes));
596d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice));
597d075f50bSSebastian Grimberg   }
598d075f50bSSebastian Grimberg   if (curl) {
599d075f50bSSebastian Grimberg     const CeedInt curl_bytes = q_bytes * num_nodes * q_comp_curl;
600d075f50bSSebastian Grimberg 
601d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_curl, curl_bytes));
602d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_curl, curl, curl_bytes, cudaMemcpyHostToDevice));
603d075f50bSSebastian Grimberg   }
604d075f50bSSebastian Grimberg 
605d075f50bSSebastian Grimberg   // Compile basis kernels
6069c25dd66SJeremy L Thompson   const char basis_kernel_source[] = "// Nontensor basis source\n#include <ceed/jit-source/cuda/cuda-ref-basis-nontensor.h>\n";
6079c25dd66SJeremy L Thompson 
608d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
609d075f50bSSebastian Grimberg   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP",
610d075f50bSSebastian Grimberg                                    q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_curl, "BASIS_NUM_COMP", num_comp));
611d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
612d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
613d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv));
614d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose));
615d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
6160d0321e0SJeremy L Thompson 
6172b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
6180d0321e0SJeremy L Thompson 
6190d0321e0SJeremy L Thompson   // Register backend functions
6202b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda));
621db2becc9SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Cuda));
6222b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda));
623*9bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
6240d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6250d0321e0SJeremy L Thompson }
6262a86cc9dSSebastian Grimberg 
6270d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
628