xref: /libCEED/backends/cuda-ref/ceed-cuda-ref-basis.c (revision d7d111ec11fb843d4be607afc8fa89b13b3a5edf)
10d0321e0SJeremy L Thompson // Copyright (c) 2017-2018, Lawrence Livermore National Security, LLC.
20d0321e0SJeremy L Thompson // Produced at the Lawrence Livermore National Laboratory. LLNL-CODE-734707.
30d0321e0SJeremy L Thompson // All Rights reserved. See files LICENSE and NOTICE for details.
40d0321e0SJeremy L Thompson //
50d0321e0SJeremy L Thompson // This file is part of CEED, a collection of benchmarks, miniapps, software
60d0321e0SJeremy L Thompson // libraries and APIs for efficient high-order finite element and spectral
70d0321e0SJeremy L Thompson // element discretizations for exascale applications. For more information and
80d0321e0SJeremy L Thompson // source code availability see http://github.com/ceed.
90d0321e0SJeremy L Thompson //
100d0321e0SJeremy L Thompson // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
110d0321e0SJeremy L Thompson // a collaborative effort of two U.S. Department of Energy organizations (Office
120d0321e0SJeremy L Thompson // of Science and the National Nuclear Security Administration) responsible for
130d0321e0SJeremy L Thompson // the planning and preparation of a capable exascale ecosystem, including
140d0321e0SJeremy L Thompson // software, applications, hardware, advanced system engineering and early
150d0321e0SJeremy L Thompson // testbed platforms, in support of the nation's exascale computing imperative.
160d0321e0SJeremy L Thompson 
170d0321e0SJeremy L Thompson #include <ceed/ceed.h>
180d0321e0SJeremy L Thompson #include <ceed/backend.h>
19437930d1SJeremy L Thompson #include <ceed/jit-tools.h>
200d0321e0SJeremy L Thompson #include <cuda.h>
210d0321e0SJeremy L Thompson #include <cuda_runtime.h>
220d0321e0SJeremy L Thompson #include "ceed-cuda-ref.h"
230d0321e0SJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
240d0321e0SJeremy L Thompson 
250d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
260d0321e0SJeremy L Thompson // Basis apply - tensor
270d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
28437930d1SJeremy L Thompson int CeedBasisApply_Cuda(CeedBasis basis, const CeedInt num_elem,
29437930d1SJeremy L Thompson                         CeedTransposeMode t_mode, CeedEvalMode eval_mode,
30437930d1SJeremy L Thompson                         CeedVector u, CeedVector v) {
310d0321e0SJeremy L Thompson   int ierr;
320d0321e0SJeremy L Thompson   Ceed ceed;
330d0321e0SJeremy L Thompson   ierr = CeedBasisGetCeed(basis, &ceed); CeedChkBackend(ierr);
340d0321e0SJeremy L Thompson   Ceed_Cuda *ceed_Cuda;
350d0321e0SJeremy L Thompson   ierr = CeedGetData(ceed, &ceed_Cuda); CeedChkBackend(ierr);
360d0321e0SJeremy L Thompson   CeedBasis_Cuda *data;
370d0321e0SJeremy L Thompson   ierr = CeedBasisGetData(basis, &data); CeedChkBackend(ierr);
38437930d1SJeremy L Thompson   const CeedInt transpose = t_mode == CEED_TRANSPOSE;
39437930d1SJeremy L Thompson   const int max_block_size = 32;
400d0321e0SJeremy L Thompson 
410d0321e0SJeremy L Thompson   // Read vectors
420d0321e0SJeremy L Thompson   const CeedScalar *d_u;
430d0321e0SJeremy L Thompson   CeedScalar *d_v;
44437930d1SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
450d0321e0SJeremy L Thompson     ierr = CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u); CeedChkBackend(ierr);
460d0321e0SJeremy L Thompson   }
470d0321e0SJeremy L Thompson   ierr = CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v); CeedChkBackend(ierr);
480d0321e0SJeremy L Thompson 
490d0321e0SJeremy L Thompson   // Clear v for transpose operation
50437930d1SJeremy L Thompson   if (t_mode == CEED_TRANSPOSE) {
510d0321e0SJeremy L Thompson     CeedInt length;
520d0321e0SJeremy L Thompson     ierr = CeedVectorGetLength(v, &length); CeedChkBackend(ierr);
530d0321e0SJeremy L Thompson     ierr = cudaMemset(d_v, 0, length * sizeof(CeedScalar));
540d0321e0SJeremy L Thompson     CeedChk_Cu(ceed, ierr);
550d0321e0SJeremy L Thompson   }
56437930d1SJeremy L Thompson   CeedInt Q_1d, dim;
57437930d1SJeremy L Thompson   ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d); CeedChkBackend(ierr);
580d0321e0SJeremy L Thompson   ierr = CeedBasisGetDimension(basis, &dim); CeedChkBackend(ierr);
590d0321e0SJeremy L Thompson 
600d0321e0SJeremy L Thompson   // Basis action
61437930d1SJeremy L Thompson   switch (eval_mode) {
620d0321e0SJeremy L Thompson   case CEED_EVAL_INTERP: {
63437930d1SJeremy L Thompson     void *interp_args[] = {(void *) &num_elem, (void *) &transpose,
64437930d1SJeremy L Thompson                            &data->d_interp_1d, &d_u, &d_v
650d0321e0SJeremy L Thompson                           };
66437930d1SJeremy L Thompson     CeedInt block_size = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
670d0321e0SJeremy L Thompson 
68437930d1SJeremy L Thompson     ierr = CeedRunKernelCuda(ceed, data->Interp, num_elem, block_size, interp_args);
690d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
700d0321e0SJeremy L Thompson   } break;
710d0321e0SJeremy L Thompson   case CEED_EVAL_GRAD: {
72437930d1SJeremy L Thompson     void *grad_args[] = {(void *) &num_elem, (void *) &transpose, &data->d_interp_1d,
73437930d1SJeremy L Thompson                          &data->d_grad_1d, &d_u, &d_v
740d0321e0SJeremy L Thompson                         };
75437930d1SJeremy L Thompson     CeedInt block_size = max_block_size;
760d0321e0SJeremy L Thompson 
77437930d1SJeremy L Thompson     ierr = CeedRunKernelCuda(ceed, data->Grad, num_elem, block_size, grad_args);
780d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
790d0321e0SJeremy L Thompson   } break;
800d0321e0SJeremy L Thompson   case CEED_EVAL_WEIGHT: {
81437930d1SJeremy L Thompson     void *weight_args[] = {(void *) &num_elem, (void *) &data->d_q_weight_1d, &d_v};
82437930d1SJeremy L Thompson     const int grid_size = num_elem;
83437930d1SJeremy L Thompson     ierr = CeedRunKernelDimCuda(ceed, data->Weight, grid_size,
84437930d1SJeremy L Thompson                                 Q_1d, dim >= 2 ? Q_1d : 1, 1,
85437930d1SJeremy L Thompson                                 weight_args); CeedChkBackend(ierr);
860d0321e0SJeremy L Thompson   } break;
870d0321e0SJeremy L Thompson   // LCOV_EXCL_START
880d0321e0SJeremy L Thompson   // Evaluate the divergence to/from the quadrature points
890d0321e0SJeremy L Thompson   case CEED_EVAL_DIV:
900d0321e0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_DIV not supported");
910d0321e0SJeremy L Thompson   // Evaluate the curl to/from the quadrature points
920d0321e0SJeremy L Thompson   case CEED_EVAL_CURL:
930d0321e0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_CURL not supported");
940d0321e0SJeremy L Thompson   // Take no action, BasisApply should not have been called
950d0321e0SJeremy L Thompson   case CEED_EVAL_NONE:
960d0321e0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_BACKEND,
970d0321e0SJeremy L Thompson                      "CEED_EVAL_NONE does not make sense in this context");
980d0321e0SJeremy L Thompson     // LCOV_EXCL_STOP
990d0321e0SJeremy L Thompson   }
1000d0321e0SJeremy L Thompson 
1010d0321e0SJeremy L Thompson   // Restore vectors
102437930d1SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
1030d0321e0SJeremy L Thompson     ierr = CeedVectorRestoreArrayRead(u, &d_u); CeedChkBackend(ierr);
1040d0321e0SJeremy L Thompson   }
1050d0321e0SJeremy L Thompson   ierr = CeedVectorRestoreArray(v, &d_v); CeedChkBackend(ierr);
1060d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1070d0321e0SJeremy L Thompson }
1080d0321e0SJeremy L Thompson 
1090d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1100d0321e0SJeremy L Thompson // Basis apply - non-tensor
1110d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
112437930d1SJeremy L Thompson int CeedBasisApplyNonTensor_Cuda(CeedBasis basis, const CeedInt num_elem,
113437930d1SJeremy L Thompson                                  CeedTransposeMode t_mode, CeedEvalMode eval_mode,
1140d0321e0SJeremy L Thompson                                  CeedVector u, CeedVector v) {
1150d0321e0SJeremy L Thompson   int ierr;
1160d0321e0SJeremy L Thompson   Ceed ceed;
1170d0321e0SJeremy L Thompson   ierr = CeedBasisGetCeed(basis, &ceed); CeedChkBackend(ierr);
1180d0321e0SJeremy L Thompson   Ceed_Cuda *ceed_Cuda;
1190d0321e0SJeremy L Thompson   ierr = CeedGetData(ceed, &ceed_Cuda); CeedChkBackend(ierr);
1200d0321e0SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
1210d0321e0SJeremy L Thompson   ierr = CeedBasisGetData(basis, &data); CeedChkBackend(ierr);
122437930d1SJeremy L Thompson   CeedInt num_nodes, num_qpts;
123437930d1SJeremy L Thompson   ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChkBackend(ierr);
124437930d1SJeremy L Thompson   ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChkBackend(ierr);
125437930d1SJeremy L Thompson   const CeedInt transpose = t_mode == CEED_TRANSPOSE;
126437930d1SJeremy L Thompson   int elems_per_block = 1;
127437930d1SJeremy L Thompson   int grid = num_elem / elems_per_block +
128437930d1SJeremy L Thompson              ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
1290d0321e0SJeremy L Thompson 
1300d0321e0SJeremy L Thompson   // Read vectors
1310d0321e0SJeremy L Thompson   const CeedScalar *d_u;
1320d0321e0SJeremy L Thompson   CeedScalar *d_v;
133437930d1SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
1340d0321e0SJeremy L Thompson     ierr = CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u); CeedChkBackend(ierr);
1350d0321e0SJeremy L Thompson   }
1360d0321e0SJeremy L Thompson   ierr = CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v); CeedChkBackend(ierr);
1370d0321e0SJeremy L Thompson 
1380d0321e0SJeremy L Thompson   // Clear v for transpose operation
139437930d1SJeremy L Thompson   if (t_mode == CEED_TRANSPOSE) {
1400d0321e0SJeremy L Thompson     CeedInt length;
1410d0321e0SJeremy L Thompson     ierr = CeedVectorGetLength(v, &length); CeedChkBackend(ierr);
1420d0321e0SJeremy L Thompson     ierr = cudaMemset(d_v, 0, length * sizeof(CeedScalar));
1430d0321e0SJeremy L Thompson     CeedChk_Cu(ceed, ierr);
1440d0321e0SJeremy L Thompson   }
1450d0321e0SJeremy L Thompson 
1460d0321e0SJeremy L Thompson   // Apply basis operation
147437930d1SJeremy L Thompson   switch (eval_mode) {
1480d0321e0SJeremy L Thompson   case CEED_EVAL_INTERP: {
149437930d1SJeremy L Thompson     void *interp_args[] = {(void *) &num_elem, (void *) &transpose,
1500d0321e0SJeremy L Thompson                            &data->d_interp, &d_u, &d_v
1510d0321e0SJeremy L Thompson                           };
152437930d1SJeremy L Thompson     if (transpose) {
153437930d1SJeremy L Thompson       ierr = CeedRunKernelDimCuda(ceed, data->Interp, grid, num_nodes, 1,
154437930d1SJeremy L Thompson                                   elems_per_block, interp_args); CeedChkBackend(ierr);
1550d0321e0SJeremy L Thompson     } else {
156437930d1SJeremy L Thompson       ierr = CeedRunKernelDimCuda(ceed, data->Interp, grid, num_qpts, 1,
157437930d1SJeremy L Thompson                                   elems_per_block, interp_args); CeedChkBackend(ierr);
1580d0321e0SJeremy L Thompson     }
1590d0321e0SJeremy L Thompson   } break;
1600d0321e0SJeremy L Thompson   case CEED_EVAL_GRAD: {
161437930d1SJeremy L Thompson     void *grad_args[] = {(void *) &num_elem, (void *) &transpose, &data->d_grad,
1620d0321e0SJeremy L Thompson                          &d_u, &d_v
1630d0321e0SJeremy L Thompson                         };
164437930d1SJeremy L Thompson     if (transpose) {
165437930d1SJeremy L Thompson       ierr = CeedRunKernelDimCuda(ceed, data->Grad, grid, num_nodes, 1,
166437930d1SJeremy L Thompson                                   elems_per_block, grad_args); CeedChkBackend(ierr);
1670d0321e0SJeremy L Thompson     } else {
168437930d1SJeremy L Thompson       ierr = CeedRunKernelDimCuda(ceed, data->Grad, grid, num_qpts, 1,
169437930d1SJeremy L Thompson                                   elems_per_block, grad_args); CeedChkBackend(ierr);
1700d0321e0SJeremy L Thompson     }
1710d0321e0SJeremy L Thompson   } break;
1720d0321e0SJeremy L Thompson   case CEED_EVAL_WEIGHT: {
173437930d1SJeremy L Thompson     void *weight_args[] = {(void *) &num_elem, (void *) &data->d_q_weight, &d_v};
174437930d1SJeremy L Thompson     ierr = CeedRunKernelDimCuda(ceed, data->Weight, grid, num_qpts, 1,
175437930d1SJeremy L Thompson                                 elems_per_block, weight_args); CeedChkBackend(ierr);
1760d0321e0SJeremy L Thompson   } break;
1770d0321e0SJeremy L Thompson   // LCOV_EXCL_START
1780d0321e0SJeremy L Thompson   // Evaluate the divergence to/from the quadrature points
1790d0321e0SJeremy L Thompson   case CEED_EVAL_DIV:
1800d0321e0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_DIV not supported");
1810d0321e0SJeremy L Thompson   // Evaluate the curl to/from the quadrature points
1820d0321e0SJeremy L Thompson   case CEED_EVAL_CURL:
1830d0321e0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_CURL not supported");
1840d0321e0SJeremy L Thompson   // Take no action, BasisApply should not have been called
1850d0321e0SJeremy L Thompson   case CEED_EVAL_NONE:
1860d0321e0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_BACKEND,
1870d0321e0SJeremy L Thompson                      "CEED_EVAL_NONE does not make sense in this context");
1880d0321e0SJeremy L Thompson     // LCOV_EXCL_STOP
1890d0321e0SJeremy L Thompson   }
1900d0321e0SJeremy L Thompson 
1910d0321e0SJeremy L Thompson   // Restore vectors
192437930d1SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
1930d0321e0SJeremy L Thompson     ierr = CeedVectorRestoreArrayRead(u, &d_u); CeedChkBackend(ierr);
1940d0321e0SJeremy L Thompson   }
1950d0321e0SJeremy L Thompson   ierr = CeedVectorRestoreArray(v, &d_v); CeedChkBackend(ierr);
1960d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1970d0321e0SJeremy L Thompson }
1980d0321e0SJeremy L Thompson 
1990d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2000d0321e0SJeremy L Thompson // Destroy tensor basis
2010d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2020d0321e0SJeremy L Thompson static int CeedBasisDestroy_Cuda(CeedBasis basis) {
2030d0321e0SJeremy L Thompson   int ierr;
2040d0321e0SJeremy L Thompson   Ceed ceed;
2050d0321e0SJeremy L Thompson   ierr = CeedBasisGetCeed(basis, &ceed); CeedChkBackend(ierr);
2060d0321e0SJeremy L Thompson 
2070d0321e0SJeremy L Thompson   CeedBasis_Cuda *data;
2080d0321e0SJeremy L Thompson   ierr = CeedBasisGetData(basis, &data); CeedChkBackend(ierr);
2090d0321e0SJeremy L Thompson 
2100d0321e0SJeremy L Thompson   CeedChk_Cu(ceed, cuModuleUnload(data->module));
2110d0321e0SJeremy L Thompson 
212437930d1SJeremy L Thompson   ierr = cudaFree(data->d_q_weight_1d); CeedChk_Cu(ceed, ierr);
213437930d1SJeremy L Thompson   ierr = cudaFree(data->d_interp_1d); CeedChk_Cu(ceed, ierr);
214437930d1SJeremy L Thompson   ierr = cudaFree(data->d_grad_1d); CeedChk_Cu(ceed, ierr);
2150d0321e0SJeremy L Thompson   ierr = CeedFree(&data); CeedChkBackend(ierr);
216437930d1SJeremy L Thompson 
2170d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2180d0321e0SJeremy L Thompson }
2190d0321e0SJeremy L Thompson 
2200d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2210d0321e0SJeremy L Thompson // Destroy non-tensor basis
2220d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2230d0321e0SJeremy L Thompson static int CeedBasisDestroyNonTensor_Cuda(CeedBasis basis) {
2240d0321e0SJeremy L Thompson   int ierr;
2250d0321e0SJeremy L Thompson   Ceed ceed;
2260d0321e0SJeremy L Thompson   ierr = CeedBasisGetCeed(basis, &ceed); CeedChkBackend(ierr);
2270d0321e0SJeremy L Thompson 
2280d0321e0SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
2290d0321e0SJeremy L Thompson   ierr = CeedBasisGetData(basis, &data); CeedChkBackend(ierr);
2300d0321e0SJeremy L Thompson 
2310d0321e0SJeremy L Thompson   CeedChk_Cu(ceed, cuModuleUnload(data->module));
2320d0321e0SJeremy L Thompson 
233437930d1SJeremy L Thompson   ierr = cudaFree(data->d_q_weight); CeedChk_Cu(ceed, ierr);
2340d0321e0SJeremy L Thompson   ierr = cudaFree(data->d_interp); CeedChk_Cu(ceed, ierr);
2350d0321e0SJeremy L Thompson   ierr = cudaFree(data->d_grad); CeedChk_Cu(ceed, ierr);
2360d0321e0SJeremy L Thompson   ierr = CeedFree(&data); CeedChkBackend(ierr);
237437930d1SJeremy L Thompson 
2380d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2390d0321e0SJeremy L Thompson }
2400d0321e0SJeremy L Thompson 
2410d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2420d0321e0SJeremy L Thompson // Create tensor
2430d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
244437930d1SJeremy L Thompson int CeedBasisCreateTensorH1_Cuda(CeedInt dim, CeedInt P_1d, CeedInt Q_1d,
245437930d1SJeremy L Thompson                                  const CeedScalar *interp_1d,
246437930d1SJeremy L Thompson                                  const CeedScalar *grad_1d,
247437930d1SJeremy L Thompson                                  const CeedScalar *q_ref_1d,
248437930d1SJeremy L Thompson                                  const CeedScalar *q_weight_1d,
2490d0321e0SJeremy L Thompson                                  CeedBasis basis) {
2500d0321e0SJeremy L Thompson   int ierr;
2510d0321e0SJeremy L Thompson   Ceed ceed;
2520d0321e0SJeremy L Thompson   ierr = CeedBasisGetCeed(basis, &ceed); CeedChkBackend(ierr);
2530d0321e0SJeremy L Thompson   CeedBasis_Cuda *data;
2540d0321e0SJeremy L Thompson   ierr = CeedCalloc(1, &data); CeedChkBackend(ierr);
2550d0321e0SJeremy L Thompson 
2560d0321e0SJeremy L Thompson   // Copy data to GPU
257437930d1SJeremy L Thompson   const CeedInt q_bytes = Q_1d * sizeof(CeedScalar);
258437930d1SJeremy L Thompson   ierr = cudaMalloc((void **)&data->d_q_weight_1d, q_bytes);
259437930d1SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
260437930d1SJeremy L Thompson   ierr = cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes,
2610d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
2620d0321e0SJeremy L Thompson 
263437930d1SJeremy L Thompson   const CeedInt interp_bytes = q_bytes * P_1d;
264437930d1SJeremy L Thompson   ierr = cudaMalloc((void **)&data->d_interp_1d, interp_bytes);
265437930d1SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
266437930d1SJeremy L Thompson   ierr = cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes,
2670d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
2680d0321e0SJeremy L Thompson 
269437930d1SJeremy L Thompson   ierr = cudaMalloc((void **)&data->d_grad_1d, interp_bytes);
270437930d1SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
271437930d1SJeremy L Thompson   ierr = cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes,
2720d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
2730d0321e0SJeremy L Thompson 
2740d0321e0SJeremy L Thompson   // Complie basis kernels
275437930d1SJeremy L Thompson   CeedInt num_comp;
276437930d1SJeremy L Thompson   ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChkBackend(ierr);
277437930d1SJeremy L Thompson   char *basis_kernel_path, *basis_kernel_source;
278437930d1SJeremy L Thompson   ierr = CeedPathConcatenate(ceed, __FILE__, "kernels/cuda-ref-basis-tensor.h",
279437930d1SJeremy L Thompson                              &basis_kernel_path); CeedChkBackend(ierr);
28046dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading Basis Kernel Source -----\n");
281437930d1SJeremy L Thompson   ierr = CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source);
282437930d1SJeremy L Thompson   CeedChkBackend(ierr);
28346dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading Basis Kernel Source Complete! -----\n");
284437930d1SJeremy L Thompson   ierr = CeedCompileCuda(ceed, basis_kernel_source, &data->module, 7,
285*d7d111ecSJeremy L Thompson                          "BASIS_Q_1D", Q_1d,
286*d7d111ecSJeremy L Thompson                          "BASIS_P_1D", P_1d,
287437930d1SJeremy L Thompson                          "BASIS_BUF_LEN", num_comp * CeedIntPow(Q_1d > P_1d ?
288437930d1SJeremy L Thompson                              Q_1d : P_1d, dim),
2890d0321e0SJeremy L Thompson                          "BASIS_DIM", dim,
290*d7d111ecSJeremy L Thompson                          "BASIS_NUM_COMP", num_comp,
291*d7d111ecSJeremy L Thompson                          "BASIS_NUM_NODES", CeedIntPow(P_1d, dim),
292*d7d111ecSJeremy L Thompson                          "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim)
2930d0321e0SJeremy L Thompson                         ); CeedChkBackend(ierr);
294437930d1SJeremy L Thompson   ierr = CeedGetKernelCuda(ceed, data->module, "Interp", &data->Interp);
2950d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
296437930d1SJeremy L Thompson   ierr = CeedGetKernelCuda(ceed, data->module, "Grad", &data->Grad);
2970d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
298437930d1SJeremy L Thompson   ierr = CeedGetKernelCuda(ceed, data->module, "Weight", &data->Weight);
2990d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
300437930d1SJeremy L Thompson   ierr = CeedFree(&basis_kernel_path); CeedChkBackend(ierr);
301437930d1SJeremy L Thompson   ierr = CeedFree(&basis_kernel_source); CeedChkBackend(ierr);
302437930d1SJeremy L Thompson 
3030d0321e0SJeremy L Thompson   ierr = CeedBasisSetData(basis, data); CeedChkBackend(ierr);
3040d0321e0SJeremy L Thompson 
3050d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Basis", basis, "Apply",
3060d0321e0SJeremy L Thompson                                 CeedBasisApply_Cuda); CeedChkBackend(ierr);
3070d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Basis", basis, "Destroy",
3080d0321e0SJeremy L Thompson                                 CeedBasisDestroy_Cuda); CeedChkBackend(ierr);
3090d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3100d0321e0SJeremy L Thompson }
3110d0321e0SJeremy L Thompson 
3120d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3130d0321e0SJeremy L Thompson // Create non-tensor
3140d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
315437930d1SJeremy L Thompson int CeedBasisCreateH1_Cuda(CeedElemTopology topo, CeedInt dim,
316437930d1SJeremy L Thompson                            CeedInt num_nodes,
317437930d1SJeremy L Thompson                            CeedInt num_qpts, const CeedScalar *interp,
3180d0321e0SJeremy L Thompson                            const CeedScalar *grad, const CeedScalar *qref,
319437930d1SJeremy L Thompson                            const CeedScalar *q_weight, CeedBasis basis) {
3200d0321e0SJeremy L Thompson   int ierr;
3210d0321e0SJeremy L Thompson   Ceed ceed;
3220d0321e0SJeremy L Thompson   ierr = CeedBasisGetCeed(basis, &ceed); CeedChkBackend(ierr);
3230d0321e0SJeremy L Thompson   CeedBasisNonTensor_Cuda *data;
3240d0321e0SJeremy L Thompson   ierr = CeedCalloc(1, &data); CeedChkBackend(ierr);
3250d0321e0SJeremy L Thompson 
3260d0321e0SJeremy L Thompson   // Copy basis data to GPU
327437930d1SJeremy L Thompson   const CeedInt q_bytes = num_qpts * sizeof(CeedScalar);
328437930d1SJeremy L Thompson   ierr = cudaMalloc((void **)&data->d_q_weight, q_bytes); CeedChk_Cu(ceed, ierr);
329437930d1SJeremy L Thompson   ierr = cudaMemcpy(data->d_q_weight, q_weight, q_bytes,
3300d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
3310d0321e0SJeremy L Thompson 
332437930d1SJeremy L Thompson   const CeedInt interp_bytes = q_bytes * num_nodes;
333437930d1SJeremy L Thompson   ierr = cudaMalloc((void **)&data->d_interp, interp_bytes);
334437930d1SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
335437930d1SJeremy L Thompson   ierr = cudaMemcpy(data->d_interp, interp, interp_bytes,
3360d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
3370d0321e0SJeremy L Thompson 
338437930d1SJeremy L Thompson   const CeedInt grad_bytes = q_bytes * num_nodes * dim;
339437930d1SJeremy L Thompson   ierr = cudaMalloc((void **)&data->d_grad, grad_bytes); CeedChk_Cu(ceed, ierr);
340437930d1SJeremy L Thompson   ierr = cudaMemcpy(data->d_grad, grad, grad_bytes,
3410d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
3420d0321e0SJeremy L Thompson 
3430d0321e0SJeremy L Thompson   // Compile basis kernels
344437930d1SJeremy L Thompson   CeedInt num_comp;
345437930d1SJeremy L Thompson   ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChkBackend(ierr);
346437930d1SJeremy L Thompson   char *basis_kernel_path, *basis_kernel_source;
347437930d1SJeremy L Thompson   ierr = CeedPathConcatenate(ceed, __FILE__, "kernels/cuda-ref-basis-nontensor.h",
348437930d1SJeremy L Thompson                              &basis_kernel_path); CeedChkBackend(ierr);
34946dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading Basis Kernel Source -----\n");
350437930d1SJeremy L Thompson   ierr = CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source);
351437930d1SJeremy L Thompson   CeedChkBackend(ierr);
35246dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading Basis Kernel Source Complete! -----\n");
353437930d1SJeremy L Thompson   ierr = CeedCompileCuda(ceed, basis_kernel_source, &data->module, 4,
354*d7d111ecSJeremy L Thompson                          "BASIS_Q", num_qpts,
355*d7d111ecSJeremy L Thompson                          "BASIS_P", num_nodes,
3560d0321e0SJeremy L Thompson                          "BASIS_DIM", dim,
357*d7d111ecSJeremy L Thompson                          "BASIS_NUM_COMP", num_comp
3580d0321e0SJeremy L Thompson                         ); CeedChk_Cu(ceed, ierr);
359437930d1SJeremy L Thompson   ierr = CeedGetKernelCuda(ceed, data->module, "Interp", &data->Interp);
3600d0321e0SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
361437930d1SJeremy L Thompson   ierr = CeedGetKernelCuda(ceed, data->module, "Grad", &data->Grad);
3620d0321e0SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
363437930d1SJeremy L Thompson   ierr = CeedGetKernelCuda(ceed, data->module, "Weight", &data->Weight);
3640d0321e0SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
365437930d1SJeremy L Thompson   ierr = CeedFree(&basis_kernel_path); CeedChkBackend(ierr);
366437930d1SJeremy L Thompson   ierr = CeedFree(&basis_kernel_source); CeedChkBackend(ierr);
3670d0321e0SJeremy L Thompson 
3680d0321e0SJeremy L Thompson   ierr = CeedBasisSetData(basis, data); CeedChkBackend(ierr);
3690d0321e0SJeremy L Thompson 
3700d0321e0SJeremy L Thompson   // Register backend functions
3710d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Basis", basis, "Apply",
3720d0321e0SJeremy L Thompson                                 CeedBasisApplyNonTensor_Cuda); CeedChkBackend(ierr);
3730d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Basis", basis, "Destroy",
3740d0321e0SJeremy L Thompson                                 CeedBasisDestroyNonTensor_Cuda); CeedChkBackend(ierr);
3750d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3760d0321e0SJeremy L Thompson }
3770d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
378