xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-gen/ceed-cuda-gen-operator.c (revision e9c76bddc0f2a44f522e0176ed6b7e0c0aa1df73)
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.
3241a4b83SYohann //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5241a4b83SYohann //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7241a4b83SYohann 
849aac155SJeremy L Thompson #include <ceed.h>
9ec3da8bcSJed Brown #include <ceed/backend.h>
1049aac155SJeremy L Thompson #include <ceed/jit-source/cuda/cuda-types.h>
118b97b69aSJeremy L Thompson #include <cuda.h>
128b97b69aSJeremy L Thompson #include <cuda_runtime.h>
133d576824SJeremy L Thompson #include <stddef.h>
14dc007f05SJeremy L Thompson #include <string.h>
152b730f8bSJeremy L Thompson 
1649aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
176d69246aSJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
182b730f8bSJeremy L Thompson #include "ceed-cuda-gen-operator-build.h"
192b730f8bSJeremy L Thompson #include "ceed-cuda-gen.h"
20241a4b83SYohann 
21ab213215SJeremy L Thompson //------------------------------------------------------------------------------
22ab213215SJeremy L Thompson // Destroy operator
23ab213215SJeremy L Thompson //------------------------------------------------------------------------------
24241a4b83SYohann static int CeedOperatorDestroy_Cuda_gen(CeedOperator op) {
258b97b69aSJeremy L Thompson   Ceed                   ceed;
26241a4b83SYohann   CeedOperator_Cuda_gen *impl;
27ca735530SJeremy L Thompson 
288b97b69aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
292b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
308b97b69aSJeremy L Thompson   if (impl->points.num_per_elem) CeedCallCuda(ceed, cudaFree((void **)impl->points.num_per_elem));
312b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl));
328b97b69aSJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
33e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
34241a4b83SYohann }
35241a4b83SYohann 
362b730f8bSJeremy L Thompson static int Waste(int threads_per_sm, int warp_size, int threads_per_elem, int elems_per_block) {
3739532cebSJed Brown   int useful_threads_per_block = threads_per_elem * elems_per_block;
3839532cebSJed Brown   // round up to nearest multiple of warp_size
39b2165e7aSSebastian Grimberg   int block_size    = CeedDivUpInt(useful_threads_per_block, warp_size) * warp_size;
4039532cebSJed Brown   int blocks_per_sm = threads_per_sm / block_size;
4139532cebSJed Brown   return threads_per_sm - useful_threads_per_block * blocks_per_sm;
4239532cebSJed Brown }
4339532cebSJed Brown 
44ea61e9acSJeremy L Thompson // Choose the least wasteful block size constrained by blocks_per_sm of max_threads_per_block.
4539532cebSJed Brown //
46ea61e9acSJeremy L Thompson // The x and y part of block[] contains per-element sizes (specified on input) while the z part is number of elements.
4739532cebSJed Brown //
48ea61e9acSJeremy L Thompson // Problem setting: we'd like to make occupancy high with relatively few inactive threads. CUDA (cuOccupancyMaxPotentialBlockSize) can tell us how
4939532cebSJed Brown // many threads can run.
5039532cebSJed Brown //
51ea61e9acSJeremy L Thompson // Note that full occupancy sometimes can't be achieved by one thread block.
52ea61e9acSJeremy L Thompson // For example, an SM might support 1536 threads in total, but only 1024 within a single thread block.
53ea61e9acSJeremy L Thompson // So cuOccupancyMaxPotentialBlockSize may suggest a block size of 768 so that two blocks can run, versus one block of 1024 will prevent a second
54ea61e9acSJeremy L Thompson // block from running. The cuda-gen kernels are pretty heavy with lots of instruction-level parallelism (ILP) so we'll generally be okay with
55ea61e9acSJeremy L Thompson // relatively low occupancy and smaller thread blocks, but we solve a reasonably general problem here. Empirically, we find that blocks bigger than
56ea61e9acSJeremy L Thompson // about 256 have higher latency and worse load balancing when the number of elements is modest.
5739532cebSJed Brown //
58ea61e9acSJeremy L Thompson // cuda-gen can't choose block sizes arbitrarily; they need to be a multiple of the number of quadrature points (or number of basis functions).
59ea61e9acSJeremy L Thompson // They also have a lot of __syncthreads(), which is another point against excessively large thread blocks.
60ea61e9acSJeremy L Thompson // Suppose I have elements with 7x7x7 quadrature points.
61ea61e9acSJeremy L Thompson // This will loop over the last dimension, so we have 7*7=49 threads per element.
62ea61e9acSJeremy L Thompson // Suppose we have two elements = 2*49=98 useful threads.
63ea61e9acSJeremy L Thompson // CUDA schedules in units of full warps (32 threads), so 128 CUDA hardware threads are effectively committed to that block.
64ea61e9acSJeremy L Thompson // Now suppose cuOccupancyMaxPotentialBlockSize returned 352.
65ea61e9acSJeremy L Thompson // We can schedule 2 blocks of size 98 (196 useful threads using 256 hardware threads), but not a third block (which would need a total of 384
66ea61e9acSJeremy L Thompson // hardware threads).
6739532cebSJed Brown //
68ea61e9acSJeremy L Thompson // If instead, we had packed 3 elements, we'd have 3*49=147 useful threads occupying 160 slots, and could schedule two blocks.
69ea61e9acSJeremy L Thompson // Alternatively, we could pack a single block of 7 elements (2*49=343 useful threads) into the 354 slots.
70ea61e9acSJeremy L Thompson // The latter has the least "waste", but __syncthreads() over-synchronizes and it might not pay off relative to smaller blocks.
712b730f8bSJeremy L Thompson static int BlockGridCalculate(CeedInt num_elem, int blocks_per_sm, int max_threads_per_block, int max_threads_z, int warp_size, int block[3],
722b730f8bSJeremy L Thompson                               int *grid) {
7339532cebSJed Brown   const int threads_per_sm   = blocks_per_sm * max_threads_per_block;
7439532cebSJed Brown   const int threads_per_elem = block[0] * block[1];
7539532cebSJed Brown   int       elems_per_block  = 1;
7639532cebSJed Brown   int       waste            = Waste(threads_per_sm, warp_size, threads_per_elem, 1);
77ca735530SJeremy L Thompson 
782b730f8bSJeremy L Thompson   for (int i = 2; i <= CeedIntMin(max_threads_per_block / threads_per_elem, num_elem); i++) {
7939532cebSJed Brown     int i_waste = Waste(threads_per_sm, warp_size, threads_per_elem, i);
80ca735530SJeremy L Thompson 
81ea61e9acSJeremy L Thompson     // We want to minimize waste, but smaller kernels have lower latency and less __syncthreads() overhead so when a larger block size has the same
8239532cebSJed Brown     // waste as a smaller one, go ahead and prefer the smaller block.
8339532cebSJed Brown     if (i_waste < waste || (i_waste == waste && threads_per_elem * i <= 128)) {
8439532cebSJed Brown       elems_per_block = i;
8539532cebSJed Brown       waste           = i_waste;
8639532cebSJed Brown     }
8739532cebSJed Brown   }
88ea61e9acSJeremy L Thompson   // In low-order elements, threads_per_elem may be sufficiently low to give an elems_per_block greater than allowable for the device, so we must
89ea61e9acSJeremy L Thompson   // check before setting the z-dimension size of the block.
9013516544Snbeams   block[2] = CeedIntMin(elems_per_block, max_threads_z);
91b2165e7aSSebastian Grimberg   *grid    = CeedDivUpInt(num_elem, elems_per_block);
9239532cebSJed Brown   return CEED_ERROR_SUCCESS;
9339532cebSJed Brown }
9439532cebSJed Brown 
95ea61e9acSJeremy L Thompson // callback for cuOccupancyMaxPotentialBlockSize, providing the amount of dynamic shared memory required for a thread block of size threads.
9639532cebSJed Brown static size_t dynamicSMemSize(int threads) { return threads * sizeof(CeedScalar); }
9739532cebSJed Brown 
98ab213215SJeremy L Thompson //------------------------------------------------------------------------------
99ab213215SJeremy L Thompson // Apply and add to output
100ab213215SJeremy L Thompson //------------------------------------------------------------------------------
101*e9c76bddSJeremy L Thompson static int CeedOperatorApplyAddCore_Cuda_gen(CeedOperator op, CUstream stream, const CeedScalar *input_arr, CeedScalar *output_arr, bool *is_run_good,
102ea04d07fSJeremy L Thompson                                              CeedRequest *request) {
103ea04d07fSJeremy L Thompson   bool                    is_at_points, is_tensor;
104241a4b83SYohann   Ceed                    ceed;
10539532cebSJed Brown   Ceed_Cuda              *cuda_data;
106ca735530SJeremy L Thompson   CeedInt                 num_elem, num_input_fields, num_output_fields;
107ca735530SJeremy L Thompson   CeedEvalMode            eval_mode;
108ca735530SJeremy L Thompson   CeedQFunctionField     *qf_input_fields, *qf_output_fields;
109241a4b83SYohann   CeedQFunction_Cuda_gen *qf_data;
110ca735530SJeremy L Thompson   CeedQFunction           qf;
111ca735530SJeremy L Thompson   CeedOperatorField      *op_input_fields, *op_output_fields;
112ca735530SJeremy L Thompson   CeedOperator_Cuda_gen  *data;
113ca735530SJeremy L Thompson 
114ea04d07fSJeremy L Thompson   // Build the operator kernel
115ea04d07fSJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernel_Cuda_gen(op, is_run_good));
116ea04d07fSJeremy L Thompson   if (!(*is_run_good)) return CEED_ERROR_SUCCESS;
117f6eafd79SJeremy L Thompson 
118c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
119c11e12f4SJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &cuda_data));
120c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
121c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
122c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
123c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
124ddae5012SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
125c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
126c11e12f4SJeremy L Thompson 
127241a4b83SYohann   // Input vectors
1289e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1292b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1309e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
1319e201c85SYohann       data->fields.inputs[i] = NULL;
132241a4b83SYohann     } else {
133681d0ea7SJeremy L Thompson       bool       is_active;
134681d0ea7SJeremy L Thompson       CeedVector vec;
135681d0ea7SJeremy L Thompson 
136241a4b83SYohann       // Get input vector
1372b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
138681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
139ea04d07fSJeremy L Thompson       if (is_active) data->fields.inputs[i] = input_arr;
140ea04d07fSJeremy L Thompson       else CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->fields.inputs[i]));
141ea04d07fSJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec));
142241a4b83SYohann     }
143241a4b83SYohann   }
144241a4b83SYohann 
145241a4b83SYohann   // Output vectors
1469e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1472b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1489e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
1499e201c85SYohann       data->fields.outputs[i] = NULL;
150241a4b83SYohann     } else {
151681d0ea7SJeremy L Thompson       bool       is_active;
152681d0ea7SJeremy L Thompson       CeedVector vec;
153681d0ea7SJeremy L Thompson 
154241a4b83SYohann       // Get output vector
1552b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
156681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
157ea04d07fSJeremy L Thompson       if (is_active) data->fields.outputs[i] = output_arr;
158ea04d07fSJeremy L Thompson       else CeedCallBackend(CeedVectorGetArrayWrite(vec, CEED_MEM_DEVICE, &data->fields.outputs[i]));
159ea04d07fSJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec));
160241a4b83SYohann     }
161241a4b83SYohann   }
162241a4b83SYohann 
1638b97b69aSJeremy L Thompson   // Point coordinates, if needed
1648b97b69aSJeremy L Thompson   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
1658b97b69aSJeremy L Thompson   if (is_at_points) {
1668b97b69aSJeremy L Thompson     // Coords
1678b97b69aSJeremy L Thompson     CeedVector vec;
1688b97b69aSJeremy L Thompson 
1698b97b69aSJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec));
1708b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->points.coords));
1718b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&vec));
1728b97b69aSJeremy L Thompson 
1738b97b69aSJeremy L Thompson     // Points per elem
1748b97b69aSJeremy L Thompson     if (num_elem != data->points.num_elem) {
1758b97b69aSJeremy L Thompson       CeedInt            *points_per_elem;
1768b97b69aSJeremy L Thompson       const CeedInt       num_bytes   = num_elem * sizeof(CeedInt);
1778b97b69aSJeremy L Thompson       CeedElemRestriction rstr_points = NULL;
1788b97b69aSJeremy L Thompson 
1798b97b69aSJeremy L Thompson       data->points.num_elem = num_elem;
1808b97b69aSJeremy L Thompson       CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
1818b97b69aSJeremy L Thompson       CeedCallBackend(CeedCalloc(num_elem, &points_per_elem));
1828b97b69aSJeremy L Thompson       for (CeedInt e = 0; e < num_elem; e++) {
1838b97b69aSJeremy L Thompson         CeedInt num_points_elem;
1848b97b69aSJeremy L Thompson 
1858b97b69aSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumPointsInElement(rstr_points, e, &num_points_elem));
1868b97b69aSJeremy L Thompson         points_per_elem[e] = num_points_elem;
1878b97b69aSJeremy L Thompson       }
1888b97b69aSJeremy L Thompson       if (data->points.num_per_elem) CeedCallCuda(ceed, cudaFree((void **)data->points.num_per_elem));
1898b97b69aSJeremy L Thompson       CeedCallCuda(ceed, cudaMalloc((void **)&data->points.num_per_elem, num_bytes));
1908b97b69aSJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy((void *)data->points.num_per_elem, points_per_elem, num_bytes, cudaMemcpyHostToDevice));
1918b97b69aSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
1928b97b69aSJeremy L Thompson       CeedCallBackend(CeedFree(&points_per_elem));
1938b97b69aSJeremy L Thompson     }
1948b97b69aSJeremy L Thompson   }
1958b97b69aSJeremy L Thompson 
196777ff853SJeremy L Thompson   // Get context data
1972b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetInnerContextData(qf, CEED_MEM_DEVICE, &qf_data->d_c));
198241a4b83SYohann 
199241a4b83SYohann   // Apply operator
2008b97b69aSJeremy L Thompson   void         *opargs[]  = {(void *)&num_elem, &qf_data->d_c, &data->indices, &data->fields, &data->B, &data->G, &data->W, &data->points};
201241a4b83SYohann   const CeedInt dim       = data->dim;
2029e201c85SYohann   const CeedInt Q_1d      = data->Q_1d;
2039e201c85SYohann   const CeedInt P_1d      = data->max_P_1d;
2049e201c85SYohann   const CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
2058b97b69aSJeremy L Thompson   int           max_threads_per_block, min_grid_size, grid;
206ca735530SJeremy L Thompson 
207dc007f05SJeremy L Thompson   CeedCallBackend(CeedOperatorHasTensorBases(op, &is_tensor));
2082b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuOccupancyMaxPotentialBlockSize(&min_grid_size, &max_threads_per_block, data->op, dynamicSMemSize, 0, 0x10000));
209dc007f05SJeremy L Thompson   int block[3] = {thread_1d, ((!is_tensor || dim == 1) ? 1 : thread_1d), -1};
210ca735530SJeremy L Thompson 
211f82027a4SJeremy L Thompson   if (is_tensor) {
21293f4dbf1SSebastian Grimberg     CeedCallBackend(BlockGridCalculate(num_elem, min_grid_size / cuda_data->device_prop.multiProcessorCount, max_threads_per_block,
2132b730f8bSJeremy L Thompson                                        cuda_data->device_prop.maxThreadsDim[2], cuda_data->device_prop.warpSize, block, &grid));
214f82027a4SJeremy L Thompson   } else {
215f82027a4SJeremy L Thompson     CeedInt elems_per_block = CeedIntMin(cuda_data->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1));
216f82027a4SJeremy L Thompson 
217f82027a4SJeremy L Thompson     grid     = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
218f82027a4SJeremy L Thompson     block[2] = elems_per_block;
219f82027a4SJeremy L Thompson   }
22039532cebSJed Brown   CeedInt shared_mem = block[0] * block[1] * block[2] * sizeof(CeedScalar);
221ca735530SJeremy L Thompson 
222*e9c76bddSJeremy L Thompson   CeedCallBackend(CeedTryRunKernelDimShared_Cuda(ceed, data->op, stream, grid, block[0], block[1], block[2], shared_mem, is_run_good, opargs));
223241a4b83SYohann 
224241a4b83SYohann   // Restore input arrays
2259e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
2262b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
2279e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
228241a4b83SYohann     } else {
229681d0ea7SJeremy L Thompson       bool       is_active;
230681d0ea7SJeremy L Thompson       CeedVector vec;
231681d0ea7SJeremy L Thompson 
2322b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
233681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
234ea04d07fSJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->fields.inputs[i]));
235ea04d07fSJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec));
236241a4b83SYohann     }
237241a4b83SYohann   }
238241a4b83SYohann 
239241a4b83SYohann   // Restore output arrays
2409e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
2412b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
2429e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
243241a4b83SYohann     } else {
244681d0ea7SJeremy L Thompson       bool       is_active;
245681d0ea7SJeremy L Thompson       CeedVector vec;
246681d0ea7SJeremy L Thompson 
2472b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
248681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
249ea04d07fSJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorRestoreArray(vec, &data->fields.outputs[i]));
250ea04d07fSJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec));
251241a4b83SYohann     }
2523b2939feSjeremylt   }
253777ff853SJeremy L Thompson 
2548b97b69aSJeremy L Thompson   // Restore point coordinates, if needed
2558b97b69aSJeremy L Thompson   if (is_at_points) {
2568b97b69aSJeremy L Thompson     CeedVector vec;
2578b97b69aSJeremy L Thompson 
2588b97b69aSJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec));
2598b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->points.coords));
2608b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&vec));
2618b97b69aSJeremy L Thompson   }
2628b97b69aSJeremy L Thompson 
263777ff853SJeremy L Thompson   // Restore context data
2642b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionRestoreInnerContextData(qf, &qf_data->d_c));
265ddae5012SJeremy L Thompson 
266ddae5012SJeremy L Thompson   // Cleanup
2679bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
268c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionDestroy(&qf));
269ea04d07fSJeremy L Thompson   if (!(*is_run_good)) data->use_fallback = true;
270ea04d07fSJeremy L Thompson   return CEED_ERROR_SUCCESS;
271ea04d07fSJeremy L Thompson }
272ddae5012SJeremy L Thompson 
273ea04d07fSJeremy L Thompson static int CeedOperatorApplyAdd_Cuda_gen(CeedOperator op, CeedVector input_vec, CeedVector output_vec, CeedRequest *request) {
274ea04d07fSJeremy L Thompson   bool              is_run_good = false;
275ea04d07fSJeremy L Thompson   const CeedScalar *input_arr   = NULL;
276ea04d07fSJeremy L Thompson   CeedScalar       *output_arr  = NULL;
277ea04d07fSJeremy L Thompson 
278ea04d07fSJeremy L Thompson   // Try to run kernel
279ea04d07fSJeremy L Thompson   if (input_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(input_vec, CEED_MEM_DEVICE, &input_arr));
280ea04d07fSJeremy L Thompson   if (output_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArray(output_vec, CEED_MEM_DEVICE, &output_arr));
281*e9c76bddSJeremy L Thompson   CeedCallBackend(CeedOperatorApplyAddCore_Cuda_gen(op, NULL, input_arr, output_arr, &is_run_good, request));
282ea04d07fSJeremy L Thompson   if (input_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorRestoreArrayRead(input_vec, &input_arr));
283ea04d07fSJeremy L Thompson   if (output_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorRestoreArray(output_vec, &output_arr));
284ea04d07fSJeremy L Thompson 
285ea04d07fSJeremy L Thompson   // Fallback on unsuccessful run
2868d12f40eSJeremy L Thompson   if (!is_run_good) {
287ddae5012SJeremy L Thompson     CeedOperator op_fallback;
288ddae5012SJeremy L Thompson 
289ea04d07fSJeremy L Thompson     CeedDebug256(CeedOperatorReturnCeed(op), CEED_DEBUG_COLOR_SUCCESS, "Falling back to /gpu/cuda/ref CeedOperator");
290ddae5012SJeremy L Thompson     CeedCallBackend(CeedOperatorGetFallback(op, &op_fallback));
291ddae5012SJeremy L Thompson     CeedCallBackend(CeedOperatorApplyAdd(op_fallback, input_vec, output_vec, request));
292ddae5012SJeremy L Thompson   }
293e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
294241a4b83SYohann }
295241a4b83SYohann 
296ab213215SJeremy L Thompson //------------------------------------------------------------------------------
297ab213215SJeremy L Thompson // Create operator
298ab213215SJeremy L Thompson //------------------------------------------------------------------------------
299241a4b83SYohann int CeedOperatorCreate_Cuda_gen(CeedOperator op) {
300241a4b83SYohann   Ceed                   ceed;
301241a4b83SYohann   CeedOperator_Cuda_gen *impl;
302241a4b83SYohann 
303ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
3042b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl));
3052b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetData(op, impl));
3062b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda_gen));
3072b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda_gen));
3089bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
309e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
310241a4b83SYohann }
3116aa95790SJeremy L Thompson 
312ab213215SJeremy L Thompson //------------------------------------------------------------------------------
313