xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-ref/ceed-cuda-ref-basis.c (revision 5aed82e4fa97acf4ba24a7f10a35f5303a6798e0)
1*5aed82e4SJeremy 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: {
63437930d1SJeremy L Thompson       void     *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
64b2165e7aSSebastian Grimberg       const int block_size_x  = Q_1d;
65b2165e7aSSebastian Grimberg       const int block_size_y  = dim >= 2 ? Q_1d : 1;
66b2165e7aSSebastian Grimberg 
67b2165e7aSSebastian Grimberg       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, num_elem, block_size_x, block_size_y, 1, weight_args));
680d0321e0SJeremy L Thompson     } break;
699ea2cfd9SJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
709ea2cfd9SJeremy L Thompson       break;
710d0321e0SJeremy L Thompson     // LCOV_EXCL_START
720d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
730d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
74bcbe1c99SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
750d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
760d0321e0SJeremy L Thompson   }
770d0321e0SJeremy L Thompson 
789ea2cfd9SJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
792b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
809ea2cfd9SJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
819ea2cfd9SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
820d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
830d0321e0SJeremy L Thompson }
840d0321e0SJeremy L Thompson 
850d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
860d0321e0SJeremy L Thompson // Basis apply - non-tensor
870d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
882b730f8bSJeremy L Thompson int CeedBasisApplyNonTensor_Cuda(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
892b730f8bSJeremy L Thompson                                  CeedVector v) {
900d0321e0SJeremy L Thompson   Ceed                     ceed;
91437930d1SJeremy L Thompson   CeedInt                  num_nodes, num_qpts;
927bbbfca3SJeremy L Thompson   const CeedInt            is_transpose    = t_mode == CEED_TRANSPOSE;
93d075f50bSSebastian Grimberg   const int                elems_per_block = 1;
94d075f50bSSebastian Grimberg   const int                grid            = CeedDivUpInt(num_elem, elems_per_block);
950d0321e0SJeremy L Thompson   const CeedScalar        *d_u;
960d0321e0SJeremy L Thompson   CeedScalar              *d_v;
97ca735530SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
98ca735530SJeremy L Thompson 
99ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
100ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
101ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
102ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumNodes(basis, &num_nodes));
103ca735530SJeremy L Thompson 
1049ea2cfd9SJeremy L Thompson   // Get read/write access to u, v
1059ea2cfd9SJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
1069ea2cfd9SJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
1072b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
1080d0321e0SJeremy L Thompson 
1090d0321e0SJeremy L Thompson   // Clear v for transpose operation
1107bbbfca3SJeremy L Thompson   if (is_transpose) {
1111f9221feSJeremy L Thompson     CeedSize length;
112ca735530SJeremy L Thompson 
1132b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(v, &length));
1142b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
1150d0321e0SJeremy L Thompson   }
1160d0321e0SJeremy L Thompson 
1170d0321e0SJeremy L Thompson   // Apply basis operation
118437930d1SJeremy L Thompson   switch (eval_mode) {
1190d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP: {
120d075f50bSSebastian Grimberg       void     *interp_args[] = {(void *)&num_elem, &data->d_interp, &d_u, &d_v};
1217bbbfca3SJeremy L Thompson       const int block_size_x  = is_transpose ? num_nodes : num_qpts;
122b2165e7aSSebastian Grimberg 
1237bbbfca3SJeremy L Thompson       if (is_transpose) {
124d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->InterpTranspose, grid, block_size_x, 1, elems_per_block, interp_args));
125d075f50bSSebastian Grimberg       } else {
126b2165e7aSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Interp, grid, block_size_x, 1, elems_per_block, interp_args));
127d075f50bSSebastian Grimberg       }
1280d0321e0SJeremy L Thompson     } break;
1290d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD: {
130d075f50bSSebastian Grimberg       void     *grad_args[]  = {(void *)&num_elem, &data->d_grad, &d_u, &d_v};
1317bbbfca3SJeremy L Thompson       const int block_size_x = is_transpose ? num_nodes : num_qpts;
132b2165e7aSSebastian Grimberg 
1337bbbfca3SJeremy L Thompson       if (is_transpose) {
134d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, grad_args));
135d075f50bSSebastian Grimberg       } else {
136d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, grad_args));
137d075f50bSSebastian Grimberg       }
138d075f50bSSebastian Grimberg     } break;
139d075f50bSSebastian Grimberg     case CEED_EVAL_DIV: {
140d075f50bSSebastian Grimberg       void     *div_args[]   = {(void *)&num_elem, &data->d_div, &d_u, &d_v};
1417bbbfca3SJeremy L Thompson       const int block_size_x = is_transpose ? num_nodes : num_qpts;
142d075f50bSSebastian Grimberg 
1437bbbfca3SJeremy L Thompson       if (is_transpose) {
144d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, div_args));
145d075f50bSSebastian Grimberg       } else {
146d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, div_args));
147d075f50bSSebastian Grimberg       }
148d075f50bSSebastian Grimberg     } break;
149d075f50bSSebastian Grimberg     case CEED_EVAL_CURL: {
150d075f50bSSebastian Grimberg       void     *curl_args[]  = {(void *)&num_elem, &data->d_curl, &d_u, &d_v};
1517bbbfca3SJeremy L Thompson       const int block_size_x = is_transpose ? num_nodes : num_qpts;
152d075f50bSSebastian Grimberg 
1537bbbfca3SJeremy L Thompson       if (is_transpose) {
154d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->DerivTranspose, grid, block_size_x, 1, elems_per_block, curl_args));
155d075f50bSSebastian Grimberg       } else {
156d075f50bSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Deriv, grid, block_size_x, 1, elems_per_block, curl_args));
157d075f50bSSebastian Grimberg       }
1580d0321e0SJeremy L Thompson     } break;
1590d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
160437930d1SJeremy L Thompson       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight, &d_v};
161b2165e7aSSebastian Grimberg 
162eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid, num_qpts, 1, elems_per_block, weight_args));
1630d0321e0SJeremy L Thompson     } break;
1649ea2cfd9SJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
1659ea2cfd9SJeremy L Thompson       break;
1660d0321e0SJeremy L Thompson   }
1670d0321e0SJeremy L Thompson 
1689ea2cfd9SJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
1692b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
1709ea2cfd9SJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
1719ea2cfd9SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
1720d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1730d0321e0SJeremy L Thompson }
1740d0321e0SJeremy L Thompson 
1750d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1760d0321e0SJeremy L Thompson // Destroy tensor basis
1770d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1780d0321e0SJeremy L Thompson static int CeedBasisDestroy_Cuda(CeedBasis basis) {
1790d0321e0SJeremy L Thompson   Ceed            ceed;
1800d0321e0SJeremy L Thompson   CeedBasis_Cuda *data;
181ca735530SJeremy L Thompson 
182ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
1832b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
1842b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
1852b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d));
1862b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp_1d));
1872b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad_1d));
1882b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
1890d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1900d0321e0SJeremy L Thompson }
1910d0321e0SJeremy L Thompson 
1920d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1930d0321e0SJeremy L Thompson // Destroy non-tensor basis
1940d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1950d0321e0SJeremy L Thompson static int CeedBasisDestroyNonTensor_Cuda(CeedBasis basis) {
1960d0321e0SJeremy L Thompson   Ceed                     ceed;
1970d0321e0SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
198ca735530SJeremy L Thompson 
199ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
2002b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
2012b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
2022b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_q_weight));
2032b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp));
2042b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad));
205d075f50bSSebastian Grimberg   CeedCallCuda(ceed, cudaFree(data->d_div));
206d075f50bSSebastian Grimberg   CeedCallCuda(ceed, cudaFree(data->d_curl));
2072b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
2080d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2090d0321e0SJeremy L Thompson }
2100d0321e0SJeremy L Thompson 
2110d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2120d0321e0SJeremy L Thompson // Create tensor
2130d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2142b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
2152b730f8bSJeremy L Thompson                                  const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
2160d0321e0SJeremy L Thompson   Ceed            ceed;
21722070f95SJeremy L Thompson   char           *basis_kernel_source;
21822070f95SJeremy L Thompson   const char     *basis_kernel_path;
219ca735530SJeremy L Thompson   CeedInt         num_comp;
220ca735530SJeremy L Thompson   const CeedInt   q_bytes      = Q_1d * sizeof(CeedScalar);
221ca735530SJeremy L Thompson   const CeedInt   interp_bytes = q_bytes * P_1d;
2220d0321e0SJeremy L Thompson   CeedBasis_Cuda *data;
223ca735530SJeremy L Thompson 
224ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
2252b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
2260d0321e0SJeremy L Thompson 
2270d0321e0SJeremy L Thompson   // Copy data to GPU
2282b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
2292b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice));
2302b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
2312b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice));
2322b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes));
2332b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice));
2340d0321e0SJeremy L Thompson 
235ecc88aebSJeremy L Thompson   // Compile basis kernels
2362b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
2372b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-tensor.h", &basis_kernel_path));
23823d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
2392b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
24023d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
241eb7e6cafSJeremy 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",
2422b730f8bSJeremy L Thompson                                    num_comp * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp,
2432b730f8bSJeremy L Thompson                                    "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim)));
244eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
245eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
246eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
2472b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_path));
2482b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_source));
249437930d1SJeremy L Thompson 
2502b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
2510d0321e0SJeremy L Thompson 
252d075f50bSSebastian Grimberg   // Register backend functions
2532b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApply_Cuda));
2542b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda));
2550d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2560d0321e0SJeremy L Thompson }
2570d0321e0SJeremy L Thompson 
2580d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
259d075f50bSSebastian Grimberg // Create non-tensor H^1
2600d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2612b730f8bSJeremy L Thompson int CeedBasisCreateH1_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *grad,
26251475c7cSJeremy L Thompson                            const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
2630d0321e0SJeremy L Thompson   Ceed                     ceed;
26422070f95SJeremy L Thompson   char                    *basis_kernel_source;
26522070f95SJeremy L Thompson   const char              *basis_kernel_path;
266d075f50bSSebastian Grimberg   CeedInt                  num_comp, q_comp_interp, q_comp_grad;
267ca735530SJeremy L Thompson   const CeedInt            q_bytes = num_qpts * sizeof(CeedScalar);
2680d0321e0SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
269ca735530SJeremy L Thompson 
270ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
2712b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
2720d0321e0SJeremy L Thompson 
2730d0321e0SJeremy L Thompson   // Copy basis data to GPU
274d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
275d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad));
2762b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes));
2772b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice));
278d075f50bSSebastian Grimberg   if (interp) {
279d075f50bSSebastian Grimberg     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
280d075f50bSSebastian Grimberg 
2812b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes));
2822b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice));
283d075f50bSSebastian Grimberg   }
284d075f50bSSebastian Grimberg   if (grad) {
285d075f50bSSebastian Grimberg     const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad;
286d075f50bSSebastian Grimberg 
2872b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad, grad_bytes));
2882b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_grad, grad, grad_bytes, cudaMemcpyHostToDevice));
289d075f50bSSebastian Grimberg   }
2900d0321e0SJeremy L Thompson 
2910d0321e0SJeremy L Thompson   // Compile basis kernels
2922b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
2932b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path));
29423d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
2952b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
29623d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
297d075f50bSSebastian Grimberg   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP",
298d075f50bSSebastian Grimberg                                    q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_grad, "BASIS_NUM_COMP", num_comp));
299d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
300d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
301d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv));
302d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose));
303d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
304d075f50bSSebastian Grimberg   CeedCallBackend(CeedFree(&basis_kernel_path));
305d075f50bSSebastian Grimberg   CeedCallBackend(CeedFree(&basis_kernel_source));
306d075f50bSSebastian Grimberg 
307d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisSetData(basis, data));
308d075f50bSSebastian Grimberg 
309d075f50bSSebastian Grimberg   // Register backend functions
310d075f50bSSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda));
311d075f50bSSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda));
312d075f50bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
313d075f50bSSebastian Grimberg }
314d075f50bSSebastian Grimberg 
315d075f50bSSebastian Grimberg //------------------------------------------------------------------------------
316d075f50bSSebastian Grimberg // Create non-tensor H(div)
317d075f50bSSebastian Grimberg //------------------------------------------------------------------------------
318d075f50bSSebastian Grimberg int CeedBasisCreateHdiv_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *div,
319d075f50bSSebastian Grimberg                              const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
320d075f50bSSebastian Grimberg   Ceed                     ceed;
32122070f95SJeremy L Thompson   char                    *basis_kernel_source;
32222070f95SJeremy L Thompson   const char              *basis_kernel_path;
323d075f50bSSebastian Grimberg   CeedInt                  num_comp, q_comp_interp, q_comp_div;
324d075f50bSSebastian Grimberg   const CeedInt            q_bytes = num_qpts * sizeof(CeedScalar);
325d075f50bSSebastian Grimberg   CeedBasisNonTensor_Cuda *data;
326d075f50bSSebastian Grimberg 
327d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
328d075f50bSSebastian Grimberg   CeedCallBackend(CeedCalloc(1, &data));
329d075f50bSSebastian Grimberg 
330d075f50bSSebastian Grimberg   // Copy basis data to GPU
331d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
332d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div));
333d075f50bSSebastian Grimberg   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes));
334d075f50bSSebastian Grimberg   CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice));
335d075f50bSSebastian Grimberg   if (interp) {
336d075f50bSSebastian Grimberg     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
337d075f50bSSebastian Grimberg 
338d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes));
339d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice));
340d075f50bSSebastian Grimberg   }
341d075f50bSSebastian Grimberg   if (div) {
342d075f50bSSebastian Grimberg     const CeedInt div_bytes = q_bytes * num_nodes * q_comp_div;
343d075f50bSSebastian Grimberg 
344d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_div, div_bytes));
345d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_div, div, div_bytes, cudaMemcpyHostToDevice));
346d075f50bSSebastian Grimberg   }
347d075f50bSSebastian Grimberg 
348d075f50bSSebastian Grimberg   // Compile basis kernels
349d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
350d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path));
351d075f50bSSebastian Grimberg   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
352d075f50bSSebastian Grimberg   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
353d075f50bSSebastian Grimberg   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
354d075f50bSSebastian Grimberg   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP",
355d075f50bSSebastian Grimberg                                    q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_div, "BASIS_NUM_COMP", num_comp));
356d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
357d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
358d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv));
359d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose));
360d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
361d075f50bSSebastian Grimberg   CeedCallBackend(CeedFree(&basis_kernel_path));
362d075f50bSSebastian Grimberg   CeedCallBackend(CeedFree(&basis_kernel_source));
363d075f50bSSebastian Grimberg 
364d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisSetData(basis, data));
365d075f50bSSebastian Grimberg 
366d075f50bSSebastian Grimberg   // Register backend functions
367d075f50bSSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda));
368d075f50bSSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda));
369d075f50bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
370d075f50bSSebastian Grimberg }
371d075f50bSSebastian Grimberg 
372d075f50bSSebastian Grimberg //------------------------------------------------------------------------------
373d075f50bSSebastian Grimberg // Create non-tensor H(curl)
374d075f50bSSebastian Grimberg //------------------------------------------------------------------------------
375d075f50bSSebastian Grimberg int CeedBasisCreateHcurl_Cuda(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
376d075f50bSSebastian Grimberg                               const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
377d075f50bSSebastian Grimberg   Ceed                     ceed;
37822070f95SJeremy L Thompson   char                    *basis_kernel_source;
37922070f95SJeremy L Thompson   const char              *basis_kernel_path;
380d075f50bSSebastian Grimberg   CeedInt                  num_comp, q_comp_interp, q_comp_curl;
381d075f50bSSebastian Grimberg   const CeedInt            q_bytes = num_qpts * sizeof(CeedScalar);
382d075f50bSSebastian Grimberg   CeedBasisNonTensor_Cuda *data;
383d075f50bSSebastian Grimberg 
384d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
385d075f50bSSebastian Grimberg   CeedCallBackend(CeedCalloc(1, &data));
386d075f50bSSebastian Grimberg 
387d075f50bSSebastian Grimberg   // Copy basis data to GPU
388d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
389d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl));
390d075f50bSSebastian Grimberg   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight, q_bytes));
391d075f50bSSebastian Grimberg   CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight, q_weight, q_bytes, cudaMemcpyHostToDevice));
392d075f50bSSebastian Grimberg   if (interp) {
393d075f50bSSebastian Grimberg     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
394d075f50bSSebastian Grimberg 
395d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp, interp_bytes));
396d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_interp, interp, interp_bytes, cudaMemcpyHostToDevice));
397d075f50bSSebastian Grimberg   }
398d075f50bSSebastian Grimberg   if (curl) {
399d075f50bSSebastian Grimberg     const CeedInt curl_bytes = q_bytes * num_nodes * q_comp_curl;
400d075f50bSSebastian Grimberg 
401d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_curl, curl_bytes));
402d075f50bSSebastian Grimberg     CeedCallCuda(ceed, cudaMemcpy(data->d_curl, curl, curl_bytes, cudaMemcpyHostToDevice));
403d075f50bSSebastian Grimberg   }
404d075f50bSSebastian Grimberg 
405d075f50bSSebastian Grimberg   // Compile basis kernels
406d075f50bSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
407d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-nontensor.h", &basis_kernel_path));
408d075f50bSSebastian Grimberg   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
409d075f50bSSebastian Grimberg   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
410d075f50bSSebastian Grimberg   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
411d075f50bSSebastian Grimberg   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 5, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "BASIS_Q_COMP_INTERP",
412d075f50bSSebastian Grimberg                                    q_comp_interp, "BASIS_Q_COMP_DERIV", q_comp_curl, "BASIS_NUM_COMP", num_comp));
413d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
414d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
415d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Deriv", &data->Deriv));
416d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "DerivTranspose", &data->DerivTranspose));
417d075f50bSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
4182b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_path));
4192b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_source));
4200d0321e0SJeremy L Thompson 
4212b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
4220d0321e0SJeremy L Thompson 
4230d0321e0SJeremy L Thompson   // Register backend functions
4242b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Cuda));
4252b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroyNonTensor_Cuda));
4260d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4270d0321e0SJeremy L Thompson }
4282a86cc9dSSebastian Grimberg 
4290d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
430