xref: /libCEED/backends/hip-shared/ceed-hip-shared-basis.c (revision c9d5affad74485f8d1e55e6be07e3d9f76bd4cae)
1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 #include <ceed.h>
9 #include <ceed/backend.h>
10 #include <ceed/jit-tools.h>
11 #include <stdbool.h>
12 #include <stddef.h>
13 #include <string.h>
14 #include <hip/hip_runtime.h>
15 
16 #include "../hip/ceed-hip-common.h"
17 #include "../hip/ceed-hip-compile.h"
18 #include "ceed-hip-shared.h"
19 
20 //------------------------------------------------------------------------------
21 // Compute a block size based on required minimum threads
22 //------------------------------------------------------------------------------
23 static CeedInt ComputeBlockSizeFromRequirement(const CeedInt required) {
24   CeedInt maxSize     = 1024;  // Max total threads per block
25   CeedInt currentSize = 64;    // Start with one group
26 
27   while (currentSize < maxSize) {
28     if (currentSize > required) break;
29     else currentSize = currentSize * 2;
30   }
31   return currentSize;
32 }
33 
34 //------------------------------------------------------------------------------
35 // Compute required thread block sizes for basis kernels given P, Q, dim, and
36 // num_comp (num_comp not currently used, but may be again in other basis
37 // parallelization options)
38 //------------------------------------------------------------------------------
39 static int ComputeBasisThreadBlockSizes(const CeedInt dim, const CeedInt P_1d, const CeedInt Q_1d, const CeedInt num_comp, CeedInt *block_sizes) {
40   // Note that this will use the same block sizes for all dimensions when compiling,
41   // but as each basis object is defined for a particular dimension, we will never
42   // call any kernels except the ones for the dimension for which we have computed the
43   // block sizes.
44   const CeedInt thread_1d = CeedIntMax(P_1d, Q_1d);
45 
46   switch (dim) {
47     case 1: {
48       // Interp kernels:
49       block_sizes[0] = 256;
50 
51       // Grad kernels:
52       block_sizes[1] = 256;
53 
54       // Weight kernels:
55       block_sizes[2] = 256;
56     } break;
57     case 2: {
58       // Interp kernels:
59       CeedInt required = thread_1d * thread_1d;
60 
61       block_sizes[0] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
62 
63       // Grad kernels: currently use same required minimum threads
64       block_sizes[1] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
65 
66       // Weight kernels:
67       required       = CeedIntMax(64, Q_1d * Q_1d);
68       block_sizes[2] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
69 
70     } break;
71     case 3: {
72       // Interp kernels:
73       CeedInt required = thread_1d * thread_1d;
74 
75       block_sizes[0] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
76 
77       // Grad kernels: currently use same required minimum threads
78       block_sizes[1] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
79 
80       // Weight kernels:
81       required       = Q_1d * Q_1d * Q_1d;
82       block_sizes[2] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
83     }
84   }
85   return CEED_ERROR_SUCCESS;
86 }
87 
88 //------------------------------------------------------------------------------
89 // Apply basis
90 //------------------------------------------------------------------------------
91 static int CeedBasisApplyTensorCore_Hip_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode,
92                                                CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
93   Ceed                  ceed;
94   Ceed_Hip             *ceed_Hip;
95   CeedInt               dim, num_comp;
96   const CeedScalar     *d_u;
97   CeedScalar           *d_v;
98   CeedBasis_Hip_shared *data;
99 
100   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
101   CeedCallBackend(CeedGetData(ceed, &ceed_Hip));
102   CeedCallBackend(CeedBasisGetData(basis, &data));
103   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
104   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
105 
106   // Get read/write access to u, v
107   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
108   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
109   if (apply_add) {
110     CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
111   } else {
112     CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
113   }
114 
115   // Apply basis operation
116   switch (eval_mode) {
117     case CEED_EVAL_INTERP: {
118       CeedInt P_1d, Q_1d;
119       CeedInt block_size = data->block_sizes[0];
120 
121       CeedCheck(data->d_interp_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; interp_1d not set", CeedEvalModes[eval_mode]);
122       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
123       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
124       CeedInt thread_1d     = CeedIntMax(Q_1d, P_1d);
125       void   *interp_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_u, &d_v};
126 
127       if (dim == 1) {
128         CeedInt elems_per_block = 64 * thread_1d > 256 ? 256 / thread_1d : 64;
129         elems_per_block         = elems_per_block > 0 ? elems_per_block : 1;
130         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
131         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
132 
133         if (t_mode == CEED_TRANSPOSE) {
134           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, 1,
135                                                      elems_per_block, shared_mem, interp_args));
136         } else {
137           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Interp, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
138         }
139       } else if (dim == 2) {
140         // Check if required threads is small enough to do multiple elems
141         const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1);
142         CeedInt       grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
143         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
144 
145         if (t_mode == CEED_TRANSPOSE) {
146           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, thread_1d,
147                                                      elems_per_block, shared_mem, interp_args));
148         } else {
149           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
150         }
151       } else if (dim == 3) {
152         const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1);
153         CeedInt       grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
154         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
155 
156         if (t_mode == CEED_TRANSPOSE) {
157           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread_1d, thread_1d,
158                                                      elems_per_block, shared_mem, interp_args));
159         } else {
160           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
161         }
162       }
163     } break;
164     case CEED_EVAL_GRAD: {
165       CeedInt P_1d, Q_1d;
166       CeedInt block_size = data->block_sizes[1];
167 
168       CeedCheck(data->d_grad_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; grad_1d not set", CeedEvalModes[eval_mode]);
169       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
170       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
171       CeedInt     thread_1d = CeedIntMax(Q_1d, P_1d);
172       CeedScalar *d_grad_1d = data->d_grad_1d;
173 
174       if (data->d_collo_grad_1d) {
175         d_grad_1d = data->d_collo_grad_1d;
176       }
177       void *grad_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_grad_1d, &d_u, &d_v};
178 
179       if (dim == 1) {
180         CeedInt elems_per_block = 64 * thread_1d > 256 ? 256 / thread_1d : 64;
181         elems_per_block         = elems_per_block > 0 ? elems_per_block : 1;
182         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
183         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
184 
185         if (t_mode == CEED_TRANSPOSE) {
186           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, 1,
187                                                      elems_per_block, shared_mem, grad_args));
188         } else {
189           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
190         }
191       } else if (dim == 2) {
192         // Check if required threads is small enough to do multiple elems
193         const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1);
194         CeedInt       grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
195         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
196 
197         if (t_mode == CEED_TRANSPOSE) {
198           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, thread_1d,
199                                                      elems_per_block, shared_mem, grad_args));
200         } else {
201           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
202         }
203       } else if (dim == 3) {
204         const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1);
205         CeedInt       grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
206         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
207 
208         if (t_mode == CEED_TRANSPOSE) {
209           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread_1d, thread_1d,
210                                                      elems_per_block, shared_mem, grad_args));
211         } else {
212           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
213         }
214       }
215     } break;
216     case CEED_EVAL_WEIGHT: {
217       CeedInt Q_1d;
218       CeedInt block_size = data->block_sizes[2];
219 
220       CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights_1d not set", CeedEvalModes[eval_mode]);
221       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
222       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
223 
224       if (dim == 1) {
225         const CeedInt opt_elems       = block_size / Q_1d;
226         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
227         const CeedInt grid_size       = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
228 
229         CeedCallBackend(CeedRunKernelDim_Hip(ceed, data->Weight, grid_size, Q_1d, elems_per_block, 1, weight_args));
230       } else if (dim == 2) {
231         const CeedInt opt_elems       = block_size / (Q_1d * Q_1d);
232         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
233         const CeedInt grid_size       = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
234 
235         CeedCallBackend(CeedRunKernelDim_Hip(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
236       } else if (dim == 3) {
237         const CeedInt opt_elems       = block_size / (Q_1d * Q_1d);
238         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
239         const CeedInt grid_size       = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
240 
241         CeedCallBackend(CeedRunKernelDim_Hip(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
242       }
243     } break;
244     case CEED_EVAL_NONE: /* handled separately below */
245       break;
246     // LCOV_EXCL_START
247     case CEED_EVAL_DIV:
248     case CEED_EVAL_CURL:
249       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
250       // LCOV_EXCL_STOP
251   }
252 
253   // Restore vectors, cover CEED_EVAL_NONE
254   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
255   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
256   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
257   CeedCallBackend(CeedDestroy(&ceed));
258   return CEED_ERROR_SUCCESS;
259 }
260 
261 int CeedBasisApplyTensor_Hip_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
262                                     CeedVector v) {
263   CeedCallBackend(CeedBasisApplyTensorCore_Hip_shared(basis, false, num_elem, t_mode, eval_mode, u, v));
264   return CEED_ERROR_SUCCESS;
265 }
266 
267 int CeedBasisApplyAddTensor_Hip_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
268                                        CeedVector v) {
269   CeedCallBackend(CeedBasisApplyTensorCore_Hip_shared(basis, true, num_elem, t_mode, eval_mode, u, v));
270   return CEED_ERROR_SUCCESS;
271 }
272 
273 //------------------------------------------------------------------------------
274 // Basis apply - tensor AtPoints
275 //------------------------------------------------------------------------------
276 static int CeedBasisApplyAtPointsCore_Hip_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, const CeedInt *num_points,
277                                                  CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
278   Ceed                  ceed;
279   CeedInt               Q_1d, dim, max_num_points = num_points[0];
280   const CeedInt         is_transpose = t_mode == CEED_TRANSPOSE;
281   const CeedScalar     *d_x, *d_u;
282   CeedScalar           *d_v;
283   CeedBasis_Hip_shared *data;
284 
285   CeedCallBackend(CeedBasisGetData(basis, &data));
286   CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
287   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
288 
289   // Weight handled separately
290   if (eval_mode == CEED_EVAL_WEIGHT) {
291     CeedCallBackend(CeedVectorSetValue(v, 1.0));
292     return CEED_ERROR_SUCCESS;
293   }
294 
295   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
296 
297   // Check padded to uniform number of points per elem
298   for (CeedInt i = 1; i < num_elem; i++) max_num_points = CeedIntMax(max_num_points, num_points[i]);
299   {
300     CeedInt  num_comp, q_comp;
301     CeedSize len, len_required;
302 
303     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
304     CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
305     CeedCallBackend(CeedVectorGetLength(is_transpose ? u : v, &len));
306     len_required = (CeedSize)num_comp * (CeedSize)q_comp * (CeedSize)num_elem * (CeedSize)max_num_points;
307     CeedCheck(len >= len_required, ceed, CEED_ERROR_BACKEND,
308               "Vector at points must be padded to the same number of points in each element for BasisApplyAtPoints on GPU backends."
309               " Found %" CeedSize_FMT ", Required %" CeedSize_FMT,
310               len, len_required);
311   }
312 
313   // Move num_points array to device
314   if (is_transpose) {
315     const CeedInt num_bytes = num_elem * sizeof(CeedInt);
316 
317     if (num_elem != data->num_elem_at_points) {
318       data->num_elem_at_points = num_elem;
319 
320       if (data->d_points_per_elem) CeedCallHip(ceed, hipFree(data->d_points_per_elem));
321       CeedCallHip(ceed, hipMalloc((void **)&data->d_points_per_elem, num_bytes));
322       CeedCallBackend(CeedFree(&data->h_points_per_elem));
323       CeedCallBackend(CeedCalloc(num_elem, &data->h_points_per_elem));
324     }
325     if (memcmp(data->h_points_per_elem, num_points, num_bytes)) {
326       memcpy(data->h_points_per_elem, num_points, num_bytes);
327       CeedCallHip(ceed, hipMemcpy(data->d_points_per_elem, num_points, num_bytes, hipMemcpyHostToDevice));
328     }
329   }
330 
331   // Build kernels if needed
332   if (data->num_points != max_num_points) {
333     CeedInt P_1d;
334 
335     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
336     data->num_points = max_num_points;
337 
338     // -- Create interp matrix to Chebyshev coefficients
339     if (!data->d_chebyshev_interp_1d) {
340       CeedSize    interp_bytes;
341       CeedScalar *chebyshev_interp_1d;
342 
343       interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
344       CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
345       CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
346       CeedCallHip(ceed, hipMalloc((void **)&data->d_chebyshev_interp_1d, interp_bytes));
347       CeedCallHip(ceed, hipMemcpy(data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, hipMemcpyHostToDevice));
348       CeedCallBackend(CeedFree(&chebyshev_interp_1d));
349     }
350 
351     // -- Compile kernels
352     const char basis_kernel_source[] = "// AtPoints basis source\n#include <ceed/jit-source/hip/hip-shared-basis-tensor-at-points.h>\n";
353     CeedInt    num_comp;
354 
355     if (data->moduleAtPoints) CeedCallHip(ceed, hipModuleUnload(data->moduleAtPoints));
356     CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
357     CeedCallBackend(CeedCompile_Hip(ceed, basis_kernel_source, &data->moduleAtPoints, 9, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D",
358                                     CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim),
359                                     "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_NUM_PTS", max_num_points, "BASIS_INTERP_BLOCK_SIZE",
360                                     data->block_sizes[0]));
361     CeedCallBackend(CeedGetKernel_Hip(ceed, data->moduleAtPoints, "InterpAtPoints", &data->InterpAtPoints));
362     CeedCallBackend(CeedGetKernel_Hip(ceed, data->moduleAtPoints, "InterpTransposeAtPoints", &data->InterpTransposeAtPoints));
363     CeedCallBackend(CeedGetKernel_Hip(ceed, data->moduleAtPoints, "GradAtPoints", &data->GradAtPoints));
364     CeedCallBackend(CeedGetKernel_Hip(ceed, data->moduleAtPoints, "GradTransposeAtPoints", &data->GradTransposeAtPoints));
365   }
366 
367   // Get read/write access to u, v
368   CeedCallBackend(CeedVectorGetArrayRead(x_ref, CEED_MEM_DEVICE, &d_x));
369   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
370   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
371   if (apply_add) {
372     CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
373   } else {
374     // Clear v for transpose operation
375     if (is_transpose) CeedCallBackend(CeedVectorSetValue(v, 0.0));
376     CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
377   }
378 
379   // Basis action
380   switch (eval_mode) {
381     case CEED_EVAL_INTERP: {
382       CeedInt P_1d, Q_1d;
383       CeedInt block_size = data->block_sizes[0];
384 
385       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
386       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
387       CeedInt thread_1d     = CeedIntMax(Q_1d, P_1d);
388       void   *interp_args[] = {(void *)&num_elem, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v};
389 
390       if (dim == 1) {
391         CeedInt elems_per_block = 64 * thread_1d > 256 ? 256 / thread_1d : 64;
392         elems_per_block         = elems_per_block > 0 ? elems_per_block : 1;
393         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
394         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
395 
396         CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d, 1,
397                                                    elems_per_block, shared_mem, interp_args));
398       } else if (dim == 2) {
399         // Check if required threads is small enough to do multiple elems
400         const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1);
401         CeedInt       grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
402         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
403 
404         CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d,
405                                                    thread_1d, elems_per_block, shared_mem, interp_args));
406       } else if (dim == 3) {
407         const CeedInt elems_per_block = 1;
408         CeedInt       grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
409         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
410 
411         CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, is_transpose ? data->InterpTransposeAtPoints : data->InterpAtPoints, grid, thread_1d,
412                                                    thread_1d, elems_per_block, shared_mem, interp_args));
413       }
414     } break;
415     case CEED_EVAL_GRAD: {
416       CeedInt P_1d, Q_1d;
417       CeedInt block_size = data->block_sizes[0];
418 
419       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
420       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
421       CeedInt thread_1d   = CeedIntMax(Q_1d, P_1d);
422       void   *grad_args[] = {(void *)&num_elem, &data->d_chebyshev_interp_1d, &data->d_points_per_elem, &d_x, &d_u, &d_v};
423 
424       if (dim == 1) {
425         CeedInt elems_per_block = 64 * thread_1d > 256 ? 256 / thread_1d : 64;
426         elems_per_block         = elems_per_block > 0 ? elems_per_block : 1;
427         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
428         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
429 
430         CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, 1,
431                                                    elems_per_block, shared_mem, grad_args));
432       } else if (dim == 2) {
433         // Check if required threads is small enough to do multiple elems
434         const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1);
435         CeedInt       grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
436         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
437 
438         CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, thread_1d,
439                                                    elems_per_block, shared_mem, grad_args));
440       } else if (dim == 3) {
441         const CeedInt elems_per_block = 1;
442         CeedInt       grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
443         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
444 
445         CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, is_transpose ? data->GradTransposeAtPoints : data->GradAtPoints, grid, thread_1d, thread_1d,
446                                                    elems_per_block, shared_mem, grad_args));
447       }
448     } break;
449     case CEED_EVAL_WEIGHT:
450     case CEED_EVAL_NONE: /* handled separately below */
451       break;
452     // LCOV_EXCL_START
453     case CEED_EVAL_DIV:
454     case CEED_EVAL_CURL:
455       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
456       // LCOV_EXCL_STOP
457   }
458 
459   // Restore vectors, cover CEED_EVAL_NONE
460   CeedCallBackend(CeedVectorRestoreArrayRead(x_ref, &d_x));
461   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
462   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
463   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
464   return CEED_ERROR_SUCCESS;
465 }
466 
467 static int CeedBasisApplyAtPoints_Hip_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
468                                              CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
469   CeedCallBackend(CeedBasisApplyAtPointsCore_Hip_shared(basis, false, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
470   return CEED_ERROR_SUCCESS;
471 }
472 
473 static int CeedBasisApplyAddAtPoints_Hip_shared(CeedBasis basis, const CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
474                                                 CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
475   CeedCallBackend(CeedBasisApplyAtPointsCore_Hip_shared(basis, true, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
476   return CEED_ERROR_SUCCESS;
477 }
478 
479 //------------------------------------------------------------------------------
480 // Apply basis
481 //------------------------------------------------------------------------------
482 static int CeedBasisApplyNonTensorCore_Hip_shared(CeedBasis basis, bool apply_add, const CeedInt num_elem, CeedTransposeMode t_mode,
483                                                   CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
484   Ceed                  ceed;
485   Ceed_Hip             *ceed_Hip;
486   CeedInt               dim, num_comp;
487   const CeedScalar     *d_u;
488   CeedScalar           *d_v;
489   CeedBasis_Hip_shared *data;
490 
491   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
492   CeedCallBackend(CeedGetData(ceed, &ceed_Hip));
493   CeedCallBackend(CeedBasisGetData(basis, &data));
494   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
495   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
496 
497   // Get read/write access to u, v
498   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
499   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
500   if (apply_add) {
501     CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
502   } else {
503     CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
504   }
505 
506   // Apply basis operation
507   switch (eval_mode) {
508     case CEED_EVAL_INTERP: {
509       CeedInt P, Q;
510 
511       CeedCheck(data->d_interp_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; interp not set", CeedEvalModes[eval_mode]);
512       CeedCallBackend(CeedBasisGetNumNodes(basis, &P));
513       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q));
514       CeedInt thread        = CeedIntMax(Q, P);
515       void   *interp_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_u, &d_v};
516 
517       {
518         CeedInt elems_per_block = 64 * thread > 256 ? 256 / thread : 64;
519         elems_per_block         = elems_per_block > 0 ? elems_per_block : 1;
520         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
521         CeedInt shared_mem      = elems_per_block * thread * sizeof(CeedScalar);
522 
523         if (t_mode == CEED_TRANSPOSE) {
524           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, apply_add ? data->InterpTransposeAdd : data->InterpTranspose, grid, thread, 1,
525                                                      elems_per_block, shared_mem, interp_args));
526         } else {
527           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Interp, grid, thread, 1, elems_per_block, shared_mem, interp_args));
528         }
529       }
530     } break;
531     case CEED_EVAL_GRAD: {
532       CeedInt P, Q;
533 
534       CeedCheck(data->d_grad_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; grad not set", CeedEvalModes[eval_mode]);
535       CeedCallBackend(CeedBasisGetNumNodes(basis, &P));
536       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q));
537       CeedInt thread      = CeedIntMax(Q, P);
538       void   *grad_args[] = {(void *)&num_elem, &data->d_grad_1d, &d_u, &d_v};
539 
540       {
541         CeedInt elems_per_block = 64 * thread > 256 ? 256 / thread : 64;
542         elems_per_block         = elems_per_block > 0 ? elems_per_block : 1;
543         CeedInt grid            = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
544         CeedInt shared_mem      = elems_per_block * thread * sizeof(CeedScalar);
545 
546         if (t_mode == CEED_TRANSPOSE) {
547           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, apply_add ? data->GradTransposeAdd : data->GradTranspose, grid, thread, 1, elems_per_block,
548                                                      shared_mem, grad_args));
549         } else {
550           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Grad, grid, thread, 1, elems_per_block, shared_mem, grad_args));
551         }
552       }
553     } break;
554     case CEED_EVAL_WEIGHT: {
555       CeedInt P, Q;
556 
557       CeedCheck(data->d_q_weight_1d, ceed, CEED_ERROR_BACKEND, "%s not supported; q_weights not set", CeedEvalModes[eval_mode]);
558       CeedCallBackend(CeedBasisGetNumNodes(basis, &P));
559       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &Q));
560       CeedInt thread        = CeedIntMax(Q, P);
561       void   *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
562 
563       {
564         CeedInt elems_per_block = 64 * thread > 256 ? 256 / thread : 64;
565         elems_per_block         = elems_per_block > 0 ? elems_per_block : 1;
566         const CeedInt grid_size = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
567 
568         CeedCallBackend(CeedRunKernelDim_Hip(ceed, data->Weight, grid_size, thread, elems_per_block, 1, weight_args));
569       }
570     } break;
571     case CEED_EVAL_NONE: /* handled separately below */
572       break;
573     // LCOV_EXCL_START
574     case CEED_EVAL_DIV:
575     case CEED_EVAL_CURL:
576       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
577       // LCOV_EXCL_STOP
578   }
579 
580   // Restore vectors, cover CEED_EVAL_NONE
581   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
582   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
583   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
584   CeedCallBackend(CeedDestroy(&ceed));
585   return CEED_ERROR_SUCCESS;
586 }
587 
588 int CeedBasisApplyNonTensor_Hip_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
589                                        CeedVector v) {
590   CeedCallBackend(CeedBasisApplyNonTensorCore_Hip_shared(basis, false, num_elem, t_mode, eval_mode, u, v));
591   return CEED_ERROR_SUCCESS;
592 }
593 
594 int CeedBasisApplyAddNonTensor_Hip_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
595                                           CeedVector v) {
596   CeedCallBackend(CeedBasisApplyNonTensorCore_Hip_shared(basis, true, num_elem, t_mode, eval_mode, u, v));
597   return CEED_ERROR_SUCCESS;
598 }
599 
600 //------------------------------------------------------------------------------
601 // Destroy basis
602 //------------------------------------------------------------------------------
603 static int CeedBasisDestroy_Hip_shared(CeedBasis basis) {
604   Ceed                  ceed;
605   CeedBasis_Hip_shared *data;
606 
607   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
608   CeedCallBackend(CeedBasisGetData(basis, &data));
609   CeedCallHip(ceed, hipModuleUnload(data->module));
610   if (data->moduleAtPoints) CeedCallHip(ceed, hipModuleUnload(data->moduleAtPoints));
611   if (data->d_q_weight_1d) CeedCallHip(ceed, hipFree(data->d_q_weight_1d));
612   CeedCallBackend(CeedFree(&data->h_points_per_elem));
613   if (data->d_points_per_elem) CeedCallHip(ceed, hipFree(data->d_points_per_elem));
614   CeedCallHip(ceed, hipFree(data->d_interp_1d));
615   CeedCallHip(ceed, hipFree(data->d_grad_1d));
616   CeedCallHip(ceed, hipFree(data->d_collo_grad_1d));
617   CeedCallHip(ceed, hipFree(data->d_chebyshev_interp_1d));
618   CeedCallBackend(CeedFree(&data));
619   return CEED_ERROR_SUCCESS;
620 }
621 
622 //------------------------------------------------------------------------------
623 // Create tensor basis
624 //------------------------------------------------------------------------------
625 int CeedBasisCreateTensorH1_Hip_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
626                                        const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
627   Ceed                  ceed;
628   CeedInt               num_comp;
629   const CeedInt         q_bytes      = Q_1d * sizeof(CeedScalar);
630   const CeedInt         interp_bytes = q_bytes * P_1d;
631   CeedBasis_Hip_shared *data;
632 
633   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
634   CeedCallBackend(CeedCalloc(1, &data));
635 
636   // Copy basis data to GPU
637   if (q_weight_1d) {
638     CeedCallHip(ceed, hipMalloc((void **)&data->d_q_weight_1d, q_bytes));
639     CeedCallHip(ceed, hipMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, hipMemcpyHostToDevice));
640   }
641   CeedCallHip(ceed, hipMalloc((void **)&data->d_interp_1d, interp_bytes));
642   CeedCallHip(ceed, hipMemcpy(data->d_interp_1d, interp_1d, interp_bytes, hipMemcpyHostToDevice));
643   CeedCallHip(ceed, hipMalloc((void **)&data->d_grad_1d, interp_bytes));
644   CeedCallHip(ceed, hipMemcpy(data->d_grad_1d, grad_1d, interp_bytes, hipMemcpyHostToDevice));
645 
646   // Compute collocated gradient and copy to GPU
647   data->d_collo_grad_1d    = NULL;
648   bool has_collocated_grad = dim == 3 && Q_1d >= P_1d;
649 
650   if (has_collocated_grad) {
651     CeedScalar *collo_grad_1d;
652 
653     CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d));
654     CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d));
655     CeedCallHip(ceed, hipMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d));
656     CeedCallHip(ceed, hipMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, hipMemcpyHostToDevice));
657     CeedCallBackend(CeedFree(&collo_grad_1d));
658   }
659 
660   // Set number of threads per block for basis kernels
661   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
662   CeedCallBackend(ComputeBasisThreadBlockSizes(dim, P_1d, Q_1d, num_comp, data->block_sizes));
663 
664   // Compile basis kernels
665   const char basis_kernel_source[] = "// Tensor basis source\n#include <ceed/jit-source/hip/hip-shared-basis-tensor.h>\n";
666 
667   CeedCallBackend(CeedCompile_Hip(ceed, basis_kernel_source, &data->module, 11, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D",
668                                   CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim),
669                                   "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_INTERP_BLOCK_SIZE", data->block_sizes[0], "BASIS_GRAD_BLOCK_SIZE",
670                                   data->block_sizes[1], "BASIS_WEIGHT_BLOCK_SIZE", data->block_sizes[2], "BASIS_HAS_COLLOCATED_GRAD",
671                                   has_collocated_grad));
672   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "Interp", &data->Interp));
673   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
674   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd));
675   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "Grad", &data->Grad));
676   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "GradTranspose", &data->GradTranspose));
677   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd));
678   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "Weight", &data->Weight));
679 
680   CeedCallBackend(CeedBasisSetData(basis, data));
681 
682   // Register backend functions
683   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Hip_shared));
684   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddTensor_Hip_shared));
685   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAtPoints", CeedBasisApplyAtPoints_Hip_shared));
686   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAddAtPoints", CeedBasisApplyAddAtPoints_Hip_shared));
687   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Hip_shared));
688   CeedCallBackend(CeedDestroy(&ceed));
689   return CEED_ERROR_SUCCESS;
690 }
691 
692 //------------------------------------------------------------------------------
693 // Create non-tensor basis
694 //------------------------------------------------------------------------------
695 int CeedBasisCreateH1_Hip_shared(CeedElemTopology topo, CeedInt dim, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
696                                  const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
697   Ceed                  ceed;
698   CeedInt               num_comp, q_comp_interp, q_comp_grad;
699   const CeedInt         q_bytes = num_qpts * sizeof(CeedScalar);
700   CeedBasis_Hip_shared *data;
701 
702   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
703 
704   // Check shared memory size
705   {
706     Ceed_Hip *hip_data;
707 
708     CeedCallBackend(CeedGetData(ceed, &hip_data));
709     if (((size_t)num_nodes * (size_t)num_qpts * (size_t)dim + (size_t)CeedIntMax(num_nodes, num_qpts)) * sizeof(CeedScalar) >
710         hip_data->device_prop.sharedMemPerBlock) {
711       CeedCallBackend(CeedBasisCreateH1Fallback(ceed, topo, dim, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis));
712       return CEED_ERROR_SUCCESS;
713     }
714   }
715 
716   CeedCallBackend(CeedCalloc(1, &data));
717 
718   // Copy basis data to GPU
719   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
720   CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad));
721   if (q_weight) {
722     CeedCallHip(ceed, hipMalloc((void **)&data->d_q_weight_1d, q_bytes));
723     CeedCallHip(ceed, hipMemcpy(data->d_q_weight_1d, q_weight, q_bytes, hipMemcpyHostToDevice));
724   }
725   if (interp) {
726     const CeedInt interp_bytes = q_bytes * num_nodes * q_comp_interp;
727 
728     CeedCallHip(ceed, hipMalloc((void **)&data->d_interp_1d, interp_bytes));
729     CeedCallHip(ceed, hipMemcpy(data->d_interp_1d, interp, interp_bytes, hipMemcpyHostToDevice));
730   }
731   if (grad) {
732     const CeedInt grad_bytes = q_bytes * num_nodes * q_comp_grad;
733 
734     CeedCallHip(ceed, hipMalloc((void **)&data->d_grad_1d, grad_bytes));
735     CeedCallHip(ceed, hipMemcpy(data->d_grad_1d, grad, grad_bytes, hipMemcpyHostToDevice));
736   }
737 
738   // Compile basis kernels
739   const char basis_kernel_source[] = "// Non-tensor basis source\n#include <ceed/jit-source/hip/hip-shared-basis-nontensor.h>\n";
740 
741   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
742   CeedCallBackend(ComputeBasisThreadBlockSizes(dim, num_nodes, num_qpts, num_comp, data->block_sizes));
743   CeedCallBackend(CeedCompile_Hip(ceed, basis_kernel_source, &data->module, 6, "BASIS_Q", num_qpts, "BASIS_P", num_nodes, "T_1D",
744                                   CeedIntMax(num_qpts, num_nodes), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_INTERP_BLOCK_SIZE",
745                                   data->block_sizes[0]));
746   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "Interp", &data->Interp));
747   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
748   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "InterpTransposeAdd", &data->InterpTransposeAdd));
749   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "Grad", &data->Grad));
750   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "GradTranspose", &data->GradTranspose));
751   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "GradTransposeAdd", &data->GradTransposeAdd));
752   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "Weight", &data->Weight));
753 
754   CeedCallBackend(CeedBasisSetData(basis, data));
755 
756   // Register backend functions
757   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyNonTensor_Hip_shared));
758   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "ApplyAdd", CeedBasisApplyAddNonTensor_Hip_shared));
759   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Hip_shared));
760   CeedCallBackend(CeedDestroy(&ceed));
761   return CEED_ERROR_SUCCESS;
762 }
763 
764 //------------------------------------------------------------------------------
765