xref: /libCEED/backends/cuda-shared/ceed-cuda-shared-basis.c (revision 9e511c80210b4447ceabc14bc1a44d55934e0fbf)
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.
3c532df63SYohann //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5c532df63SYohann //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7c532df63SYohann 
849aac155SJeremy L Thompson #include <ceed.h>
9ec3da8bcSJed Brown #include <ceed/backend.h>
10437930d1SJeremy L Thompson #include <ceed/jit-tools.h>
113d576824SJeremy L Thompson #include <cuda.h>
123d576824SJeremy L Thompson #include <cuda_runtime.h>
1349aac155SJeremy L Thompson #include <stdbool.h>
143d576824SJeremy L Thompson #include <stddef.h>
15111870feSJeremy L Thompson #include <string.h>
16c532df63SYohann 
1749aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
182b730f8bSJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
192b730f8bSJeremy L Thompson #include "ceed-cuda-shared.h"
20c532df63SYohann 
21ab213215SJeremy L Thompson //------------------------------------------------------------------------------
22ab213215SJeremy L Thompson // Device initalization
23ab213215SJeremy L Thompson //------------------------------------------------------------------------------
24eb7e6cafSJeremy L Thompson int CeedInit_CudaInterp(CeedScalar *d_B, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B);
25eb7e6cafSJeremy L Thompson int CeedInit_CudaGrad(CeedScalar *d_B, CeedScalar *d_G, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B_ptr, CeedScalar **c_G_ptr);
26eb7e6cafSJeremy L Thompson int CeedInit_CudaCollocatedGrad(CeedScalar *d_B, CeedScalar *d_G, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B_ptr, CeedScalar **c_G_ptr);
27c532df63SYohann 
28ab213215SJeremy L Thompson //------------------------------------------------------------------------------
29ab213215SJeremy L Thompson // Apply basis
30ab213215SJeremy L Thompson //------------------------------------------------------------------------------
31db2becc9SJeremy L Thompson static int CeedBasisApplyTensorCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode,
32db2becc9SJeremy L Thompson                                                 CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
33c532df63SYohann   Ceed                   ceed;
346dbfb411Snbeams   Ceed_Cuda             *ceed_Cuda;
35437930d1SJeremy L Thompson   CeedInt                dim, num_comp;
36ca735530SJeremy L Thompson   const CeedScalar      *d_u;
37ca735530SJeremy L Thompson   CeedScalar            *d_v;
38ca735530SJeremy L Thompson   CeedBasis_Cuda_shared *data;
39ca735530SJeremy L Thompson 
40ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
41ca735530SJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
42ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
432b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
442b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
45c532df63SYohann 
469ea2cfd9SJeremy L Thompson   // Get read/write access to u, v
476574a04fSJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
486574a04fSJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
49db2becc9SJeremy L Thompson   if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
50db2becc9SJeremy L Thompson   else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
51c532df63SYohann 
52ab213215SJeremy L Thompson   // Apply basis operation
53437930d1SJeremy L Thompson   switch (eval_mode) {
54ab213215SJeremy L Thompson     case CEED_EVAL_INTERP: {
55437930d1SJeremy L Thompson       CeedInt P_1d, Q_1d;
56ca735530SJeremy L Thompson 
572b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
582b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
59437930d1SJeremy L Thompson       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
60ca735530SJeremy L Thompson 
61eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedInit_CudaInterp(data->d_interp_1d, P_1d, Q_1d, &data->c_B));
622b730f8bSJeremy L Thompson       void *interp_args[] = {(void *)&num_elem, &data->c_B, &d_u, &d_v};
63ca735530SJeremy L Thompson 
644d537eeaSYohann       if (dim == 1) {
652b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d,
66e6f67ff7SJed Brown                                                                                                  1));  // avoid >512 total threads
672b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
68437930d1SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
69b2165e7aSSebastian Grimberg 
709e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
71db2becc9SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, 1,
72db2becc9SJeremy L Thompson                                                       elems_per_block, shared_mem, interp_args));
739e201c85SYohann         } else {
74eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
759e201c85SYohann         }
76074be161SYohann Dudouit       } else if (dim == 2) {
77437930d1SJeremy L Thompson         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
78437930d1SJeremy L Thompson         // elems_per_block must be at least 1
792b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
802b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
812b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
82b2165e7aSSebastian Grimberg 
839e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
84db2becc9SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, thread_1d,
85db2becc9SJeremy L Thompson                                                       elems_per_block, shared_mem, interp_args));
869e201c85SYohann         } else {
87eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
889e201c85SYohann         }
89074be161SYohann Dudouit       } else if (dim == 3) {
90437930d1SJeremy L Thompson         CeedInt elems_per_block = 1;
912b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
922b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
93b2165e7aSSebastian Grimberg 
949e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
95db2becc9SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, thread_1d,
96db2becc9SJeremy L Thompson                                                       elems_per_block, shared_mem, interp_args));
979e201c85SYohann         } else {
98eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
99074be161SYohann Dudouit         }
1009e201c85SYohann       }
101ab213215SJeremy L Thompson     } break;
102ab213215SJeremy L Thompson     case CEED_EVAL_GRAD: {
103437930d1SJeremy L Thompson       CeedInt P_1d, Q_1d;
104ca735530SJeremy L Thompson 
1052b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
1062b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
107437930d1SJeremy L Thompson       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
108ca735530SJeremy L Thompson 
1099e201c85SYohann       if (data->d_collo_grad_1d) {
110eb7e6cafSJeremy L Thompson         CeedCallBackend(CeedInit_CudaCollocatedGrad(data->d_interp_1d, data->d_collo_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G));
1119e201c85SYohann       } else {
112eb7e6cafSJeremy L Thompson         CeedCallBackend(CeedInit_CudaGrad(data->d_interp_1d, data->d_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G));
1139e201c85SYohann       }
1142b730f8bSJeremy L Thompson       void *grad_args[] = {(void *)&num_elem, &data->c_B, &data->c_G, &d_u, &d_v};
1154d537eeaSYohann       if (dim == 1) {
1162b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d,
117e6f67ff7SJed Brown                                                                                                  1));  // avoid >512 total threads
1182b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
119437930d1SJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
120b2165e7aSSebastian Grimberg 
1219e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
122db2becc9SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, 1,
123db2becc9SJeremy L Thompson                                                       elems_per_block, shared_mem, grad_args));
1249e201c85SYohann         } else {
125eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
1269e201c85SYohann         }
127074be161SYohann Dudouit       } else if (dim == 2) {
128437930d1SJeremy L Thompson         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
129437930d1SJeremy L Thompson         // elems_per_block must be at least 1
1302b730f8bSJeremy L Thompson         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
1312b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
1322b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
133b2165e7aSSebastian Grimberg 
1349e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
135db2becc9SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, thread_1d,
136db2becc9SJeremy L Thompson                                                       elems_per_block, shared_mem, grad_args));
1379e201c85SYohann         } else {
138eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
1399e201c85SYohann         }
140074be161SYohann Dudouit       } else if (dim == 3) {
141437930d1SJeremy L Thompson         CeedInt elems_per_block = 1;
1422b730f8bSJeremy L Thompson         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
1432b730f8bSJeremy L Thompson         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
144b2165e7aSSebastian Grimberg 
1459e201c85SYohann         if (t_mode == CEED_TRANSPOSE) {
146db2becc9SJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, thread_1d,
147db2becc9SJeremy L Thompson                                                       elems_per_block, shared_mem, grad_args));
1489e201c85SYohann         } else {
149eb7e6cafSJeremy L Thompson           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
1509e201c85SYohann         }
151074be161SYohann Dudouit       }
152ab213215SJeremy L Thompson     } break;
153ab213215SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
154437930d1SJeremy L Thompson       CeedInt Q_1d;
155397164e9SSebastian Grimberg       CeedInt block_size = 32;
156ca735530SJeremy L Thompson 
157097cc795SJames Wright       CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights_1d not set", CeedEvalModes[eval_mode]);
1582b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
159437930d1SJeremy L Thompson       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
160074be161SYohann Dudouit       if (dim == 1) {
161397164e9SSebastian Grimberg         const CeedInt elems_per_block = block_size / Q_1d;
162b2165e7aSSebastian Grimberg         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
163b2165e7aSSebastian Grimberg 
164b2165e7aSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, elems_per_block, 1, weight_args));
165074be161SYohann Dudouit       } else if (dim == 2) {
166397164e9SSebastian Grimberg         const CeedInt opt_elems       = block_size / (Q_1d * Q_1d);
167437930d1SJeremy L Thompson         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
168b2165e7aSSebastian Grimberg         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
169b2165e7aSSebastian Grimberg 
170b2165e7aSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
171074be161SYohann Dudouit       } else if (dim == 3) {
172397164e9SSebastian Grimberg         const CeedInt opt_elems       = block_size / (Q_1d * Q_1d);
1739e201c85SYohann         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
174b2165e7aSSebastian Grimberg         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
175b2165e7aSSebastian Grimberg 
176b2165e7aSSebastian Grimberg         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
177074be161SYohann Dudouit       }
178ab213215SJeremy L Thompson     } break;
1799ea2cfd9SJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
1809ea2cfd9SJeremy L Thompson       break;
181ab213215SJeremy L Thompson     // LCOV_EXCL_START
182ab213215SJeremy L Thompson     case CEED_EVAL_DIV:
183ab213215SJeremy L Thompson     case CEED_EVAL_CURL:
184bcbe1c99SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
185ab213215SJeremy L Thompson       // LCOV_EXCL_STOP
186c532df63SYohann   }
187c532df63SYohann 
1889ea2cfd9SJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
1892b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
1909ea2cfd9SJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
1919ea2cfd9SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
192e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
193c532df63SYohann }
194c532df63SYohann 
195db2becc9SJeremy L Thompson static int CeedBasisApplyTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
196db2becc9SJeremy L Thompson                                             CeedVector v) {
197db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyTensorCore_Cuda_shared(basis, false, num_elem, t_mode, eval_mode, u, v));
198db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
199db2becc9SJeremy L Thompson }
200db2becc9SJeremy L Thompson 
201db2becc9SJeremy L Thompson static int CeedBasisApplyAddTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
202db2becc9SJeremy L Thompson                                                CeedVector u, CeedVector v) {
203db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyTensorCore_Cuda_shared(basis, true, num_elem, t_mode, eval_mode, u, v));
204db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
205db2becc9SJeremy L Thompson }
206db2becc9SJeremy L Thompson 
207ab213215SJeremy L Thompson //------------------------------------------------------------------------------
2081dda9c1aSJeremy L Thompson // Basis apply - tensor AtPoints
2091dda9c1aSJeremy L Thompson //------------------------------------------------------------------------------
210db2becc9SJeremy L Thompson static int CeedBasisApplyAtPointsCore_Cuda_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, const CeedInt *num_points,
211db2becc9SJeremy L Thompson                                                   CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
2121dda9c1aSJeremy L Thompson   Ceed                   ceed;
2131dda9c1aSJeremy L Thompson   CeedInt                Q_1d, dim, max_num_points = num_points[0];
2141dda9c1aSJeremy L Thompson   const CeedInt          is_transpose   = t_mode == CEED_TRANSPOSE;
2151dda9c1aSJeremy L Thompson   const int              max_block_size = 32;
2161dda9c1aSJeremy L Thompson   const CeedScalar      *d_x, *d_u;
2171dda9c1aSJeremy L Thompson   CeedScalar            *d_v;
2181dda9c1aSJeremy L Thompson   CeedBasis_Cuda_shared *data;
2191dda9c1aSJeremy L Thompson 
2201dda9c1aSJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
2211dda9c1aSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
2221dda9c1aSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
2231dda9c1aSJeremy L Thompson   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
2241dda9c1aSJeremy L Thompson 
2251dda9c1aSJeremy L Thompson   // Weight handled separately
2261dda9c1aSJeremy L Thompson   if (eval_mode == CEED_EVAL_WEIGHT) {
2275a5594ffSJeremy L Thompson     CeedCallBackend(CeedVectorSetValue(v, 1.0));
2281dda9c1aSJeremy L Thompson     return CEED_ERROR_SUCCESS;
2291dda9c1aSJeremy L Thompson   }
2301dda9c1aSJeremy L Thompson 
231111870feSJeremy L Thompson   // Check padded to uniform number of points per elem
232111870feSJeremy L Thompson   for (CeedInt i = 1; i < num_elem; i++) max_num_points = CeedIntMax(max_num_points, num_points[i]);
233111870feSJeremy L Thompson   {
234111870feSJeremy L Thompson     CeedInt  num_comp, q_comp;
235111870feSJeremy L Thompson     CeedSize len, len_required;
236111870feSJeremy L Thompson 
237111870feSJeremy L Thompson     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
238111870feSJeremy L Thompson     CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
239111870feSJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(is_transpose ? u : v, &len));
240111870feSJeremy L Thompson     len_required = (CeedSize)num_comp * (CeedSize)q_comp * (CeedSize)num_elem * (CeedSize)max_num_points;
241111870feSJeremy L Thompson     CeedCheck(len >= len_required, ceed, CEED_ERROR_BACKEND,
242111870feSJeremy L Thompson               "Vector at points must be padded to the same number of points in each element for BasisApplyAtPoints on GPU backends."
243111870feSJeremy L Thompson               " Found %" CeedSize_FMT ", Required %" CeedSize_FMT,
244111870feSJeremy L Thompson               len, len_required);
245111870feSJeremy L Thompson   }
246111870feSJeremy L Thompson 
247111870feSJeremy L Thompson   // Move num_points array to device
248111870feSJeremy L Thompson   if (is_transpose) {
249111870feSJeremy L Thompson     const CeedInt num_bytes = num_elem * sizeof(CeedInt);
250111870feSJeremy L Thompson 
251111870feSJeremy L Thompson     if (num_elem != data->num_elem_at_points) {
252111870feSJeremy L Thompson       data->num_elem_at_points = num_elem;
253111870feSJeremy L Thompson 
254111870feSJeremy L Thompson       if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem));
255111870feSJeremy L Thompson       CeedCallCuda(ceed, cudaMalloc((void **)&data->d_points_per_elem, num_bytes));
256111870feSJeremy L Thompson       CeedCallBackend(CeedFree(&data->h_points_per_elem));
257111870feSJeremy L Thompson       CeedCallBackend(CeedCalloc(num_elem, &data->h_points_per_elem));
258111870feSJeremy L Thompson     }
259*9e511c80SJeremy L Thompson     if (memcmp(data->h_points_per_elem, num_points, num_bytes)) {
260111870feSJeremy L Thompson       memcpy(data->h_points_per_elem, num_points, num_bytes);
261111870feSJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy(data->d_points_per_elem, num_points, num_bytes, cudaMemcpyHostToDevice));
262111870feSJeremy L Thompson     }
263111870feSJeremy L Thompson   }
264111870feSJeremy L Thompson 
2651dda9c1aSJeremy L Thompson   // Build kernels if needed
2661dda9c1aSJeremy L Thompson   if (data->num_points != max_num_points) {
2671dda9c1aSJeremy L Thompson     CeedInt P_1d;
2681dda9c1aSJeremy L Thompson 
2691dda9c1aSJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
2701dda9c1aSJeremy L Thompson     data->num_points = max_num_points;
2711dda9c1aSJeremy L Thompson 
2721dda9c1aSJeremy L Thompson     // -- Create interp matrix to Chebyshev coefficients
2731dda9c1aSJeremy L Thompson     if (!data->d_chebyshev_interp_1d) {
2741dda9c1aSJeremy L Thompson       CeedSize    interp_bytes;
2751dda9c1aSJeremy L Thompson       CeedScalar *chebyshev_interp_1d;
2761dda9c1aSJeremy L Thompson 
2771dda9c1aSJeremy L Thompson       interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
2781dda9c1aSJeremy L Thompson       CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
2795a5594ffSJeremy L Thompson       CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
2801dda9c1aSJeremy L Thompson       CeedCallCuda(ceed, cudaMalloc((void **)&data->d_chebyshev_interp_1d, interp_bytes));
2811dda9c1aSJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy(data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
2821dda9c1aSJeremy L Thompson       CeedCallBackend(CeedFree(&chebyshev_interp_1d));
2831dda9c1aSJeremy L Thompson     }
2841dda9c1aSJeremy L Thompson 
2851dda9c1aSJeremy L Thompson     // -- Compile kernels
2861dda9c1aSJeremy L Thompson     char       *basis_kernel_source;
2871dda9c1aSJeremy L Thompson     const char *basis_kernel_path;
2881dda9c1aSJeremy L Thompson     CeedInt     num_comp;
2891dda9c1aSJeremy L Thompson 
2901dda9c1aSJeremy L Thompson     if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
2911dda9c1aSJeremy L Thompson     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
2921dda9c1aSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-basis-tensor-at-points.h", &basis_kernel_path));
2931dda9c1aSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
2941dda9c1aSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
2951dda9c1aSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
2961dda9c1aSJeremy 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",
297f7c9815fSJeremy L Thompson                                      Q_1d * CeedIntPow(Q_1d > P_1d ? Q_1d : P_1d, dim - 1), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp,
2981dda9c1aSJeremy L Thompson                                      "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_NUM_PTS",
299f7c9815fSJeremy L Thompson                                      max_num_points, "POINTS_BUFF_LEN", CeedIntPow(Q_1d, dim - 1)));
3001dda9c1aSJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "InterpAtPoints", &data->InterpAtPoints));
3011dda9c1aSJeremy L Thompson     CeedCallBackend(CeedGetKernel_Cuda(ceed, data->moduleAtPoints, "GradAtPoints", &data->GradAtPoints));
3021dda9c1aSJeremy L Thompson     CeedCallBackend(CeedFree(&basis_kernel_path));
3031dda9c1aSJeremy L Thompson     CeedCallBackend(CeedFree(&basis_kernel_source));
3041dda9c1aSJeremy L Thompson   }
3051dda9c1aSJeremy L Thompson 
3061dda9c1aSJeremy L Thompson   // Get read/write access to u, v
3071dda9c1aSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayRead(x_ref, CEED_MEM_DEVICE, &d_x));
3081dda9c1aSJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
3091dda9c1aSJeremy L Thompson   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
310db2becc9SJeremy L Thompson   if (apply_add) CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
311db2becc9SJeremy L Thompson   else CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
3121dda9c1aSJeremy L Thompson 
3131dda9c1aSJeremy L Thompson   // Clear v for transpose operation
314db2becc9SJeremy L Thompson   if (is_transpose && !apply_add) {
3151dda9c1aSJeremy L Thompson     CeedSize length;
3161dda9c1aSJeremy L Thompson 
3171dda9c1aSJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(v, &length));
3181dda9c1aSJeremy L Thompson     CeedCallCuda(ceed, cudaMemset(d_v, 0, length * sizeof(CeedScalar)));
3191dda9c1aSJeremy L Thompson   }
3201dda9c1aSJeremy L Thompson 
3211dda9c1aSJeremy L Thompson   // Basis action
3221dda9c1aSJeremy L Thompson   switch (eval_mode) {
3231dda9c1aSJeremy L Thompson     case CEED_EVAL_INTERP: {
324111870feSJeremy 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};
3251dda9c1aSJeremy L Thompson       const CeedInt block_size = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
3261dda9c1aSJeremy L Thompson 
3271dda9c1aSJeremy L Thompson       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->InterpAtPoints, num_elem, block_size, interp_args));
3281dda9c1aSJeremy L Thompson     } break;
3291dda9c1aSJeremy L Thompson     case CEED_EVAL_GRAD: {
330111870feSJeremy 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};
3312d10e82cSJeremy L Thompson       const CeedInt block_size = CeedIntMin(CeedIntPow(Q_1d, dim), max_block_size);
3321dda9c1aSJeremy L Thompson 
3331dda9c1aSJeremy L Thompson       CeedCallBackend(CeedRunKernel_Cuda(ceed, data->GradAtPoints, num_elem, block_size, grad_args));
3341dda9c1aSJeremy L Thompson     } break;
3351dda9c1aSJeremy L Thompson     case CEED_EVAL_WEIGHT:
3361dda9c1aSJeremy L Thompson     case CEED_EVAL_NONE: /* handled separately below */
3371dda9c1aSJeremy L Thompson       break;
3381dda9c1aSJeremy L Thompson     // LCOV_EXCL_START
3391dda9c1aSJeremy L Thompson     case CEED_EVAL_DIV:
3401dda9c1aSJeremy L Thompson     case CEED_EVAL_CURL:
3411dda9c1aSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
3421dda9c1aSJeremy L Thompson       // LCOV_EXCL_STOP
3431dda9c1aSJeremy L Thompson   }
3441dda9c1aSJeremy L Thompson 
3451dda9c1aSJeremy L Thompson   // Restore vectors, cover CEED_EVAL_NONE
3461dda9c1aSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArrayRead(x_ref, &d_x));
3471dda9c1aSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
3481dda9c1aSJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
3491dda9c1aSJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
3501dda9c1aSJeremy L Thompson   return CEED_ERROR_SUCCESS;
3511dda9c1aSJeremy L Thompson }
3521dda9c1aSJeremy L Thompson 
353db2becc9SJeremy L Thompson static int CeedBasisApplyAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
354db2becc9SJeremy L Thompson                                               CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
355db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda_shared(basis, false, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
356db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
357db2becc9SJeremy L Thompson }
358db2becc9SJeremy L Thompson 
359db2becc9SJeremy L Thompson static int CeedBasisApplyAddAtPoints_Cuda_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
360db2becc9SJeremy L Thompson                                                  CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
361db2becc9SJeremy L Thompson   CeedCallBackend(CeedBasisApplyAtPointsCore_Cuda_shared(basis, true, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
362db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
363db2becc9SJeremy L Thompson }
364db2becc9SJeremy L Thompson 
3651dda9c1aSJeremy L Thompson //------------------------------------------------------------------------------
366ab213215SJeremy L Thompson // Destroy basis
367ab213215SJeremy L Thompson //------------------------------------------------------------------------------
368c532df63SYohann static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) {
369c532df63SYohann   Ceed                   ceed;
370c532df63SYohann   CeedBasis_Cuda_shared *data;
371ca735530SJeremy L Thompson 
372ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
3732b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetData(basis, &data));
3742b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuModuleUnload(data->module));
3751dda9c1aSJeremy L Thompson   if (data->moduleAtPoints) CeedCallCuda(ceed, cuModuleUnload(data->moduleAtPoints));
376097cc795SJames Wright   if (data->d_q_weight_1d) CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d));
377111870feSJeremy L Thompson   CeedCallBackend(CeedFree(&data->h_points_per_elem));
378111870feSJeremy L Thompson   if (data->d_points_per_elem) CeedCallCuda(ceed, cudaFree(data->d_points_per_elem));
3792b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_interp_1d));
3802b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_grad_1d));
3812b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_collo_grad_1d));
3821dda9c1aSJeremy L Thompson   CeedCallCuda(ceed, cudaFree(data->d_chebyshev_interp_1d));
3832b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data));
384e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
385c532df63SYohann }
386c532df63SYohann 
387ab213215SJeremy L Thompson //------------------------------------------------------------------------------
388ab213215SJeremy L Thompson // Create tensor basis
389ab213215SJeremy L Thompson //------------------------------------------------------------------------------
3902b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
3912b730f8bSJeremy L Thompson                                         const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
392c532df63SYohann   Ceed                   ceed;
39322070f95SJeremy L Thompson   char                  *basis_kernel_source;
39422070f95SJeremy L Thompson   const char            *basis_kernel_path;
395ca735530SJeremy L Thompson   CeedInt                num_comp;
396ca735530SJeremy L Thompson   const CeedInt          q_bytes      = Q_1d * sizeof(CeedScalar);
397ca735530SJeremy L Thompson   const CeedInt          interp_bytes = q_bytes * P_1d;
398c532df63SYohann   CeedBasis_Cuda_shared *data;
399ca735530SJeremy L Thompson 
400ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
4012b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &data));
402c532df63SYohann 
403ab213215SJeremy L Thompson   // Copy basis data to GPU
404097cc795SJames Wright   if (q_weight_1d) {
4052b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
4062b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice));
407097cc795SJames Wright   }
4082b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
4092b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice));
4102b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes));
4112b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice));
412c532df63SYohann 
413ab213215SJeremy L Thompson   // Compute collocated gradient and copy to GPU
414437930d1SJeremy L Thompson   data->d_collo_grad_1d    = NULL;
4159e201c85SYohann   bool has_collocated_grad = dim == 3 && Q_1d >= P_1d;
416ca735530SJeremy L Thompson 
4179e201c85SYohann   if (has_collocated_grad) {
418437930d1SJeremy L Thompson     CeedScalar *collo_grad_1d;
419ca735530SJeremy L Thompson 
4202b730f8bSJeremy L Thompson     CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d));
4212b730f8bSJeremy L Thompson     CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d));
4222b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d));
4232b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, cudaMemcpyHostToDevice));
4242b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&collo_grad_1d));
425ac421f39SYohann   }
426ac421f39SYohann 
427ab213215SJeremy L Thompson   // Compile basis kernels
4282b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
4292b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-shared-basis-tensor.h", &basis_kernel_path));
43023d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
4312b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
43223d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete -----\n");
433eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D",
434eb7e6cafSJeremy L Thompson                                    CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim),
435eb7e6cafSJeremy L Thompson                                    "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_HAS_COLLOCATED_GRAD", has_collocated_grad));
436eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
437eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
438db2becc9SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd));
439eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
440eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose));
441db2becc9SJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd));
442eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
4432b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_path));
4442b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&basis_kernel_source));
445c532df63SYohann 
4462b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisSetData(basis, data));
447ab213215SJeremy L Thompson 
448ab213215SJeremy L Thompson   // Register backend functions
4492b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared));
450db2becc9SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddTensor_Cuda_shared));
4511dda9c1aSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAtPoints", CeedBasisApplyAtPoints_Cuda_shared));
4522b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared));
453e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
454c532df63SYohann }
4552a86cc9dSSebastian Grimberg 
456ab213215SJeremy L Thompson //------------------------------------------------------------------------------
457