xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-gen/ceed-cuda-gen-operator.c (revision 8d12f40e0e187f71c4a1a78742076f931e72da09)
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 //------------------------------------------------------------------------------
1012b730f8bSJeremy L Thompson static int CeedOperatorApplyAdd_Cuda_gen(CeedOperator op, CeedVector input_vec, CeedVector output_vec, CeedRequest *request) {
102*8d12f40eSJeremy L Thompson   bool                    is_at_points, is_tensor, is_run_good = true;
103241a4b83SYohann   Ceed                    ceed;
10439532cebSJed Brown   Ceed_Cuda              *cuda_data;
105ca735530SJeremy L Thompson   CeedInt                 num_elem, num_input_fields, num_output_fields;
106ca735530SJeremy L Thompson   CeedEvalMode            eval_mode;
107ca735530SJeremy L Thompson   CeedVector              output_vecs[CEED_FIELD_MAX] = {NULL};
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 
114ddae5012SJeremy L Thompson   // Creation of the operator
115f6eafd79SJeremy L Thompson   {
116ddae5012SJeremy L Thompson     bool is_good_build = false;
117f6eafd79SJeremy L Thompson 
118ddae5012SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernel_Cuda_gen(op, &is_good_build));
119ddae5012SJeremy L Thompson     if (!is_good_build) {
120f6eafd79SJeremy L Thompson       CeedOperator op_fallback;
121f6eafd79SJeremy L Thompson 
122ddae5012SJeremy L Thompson       CeedDebug256(CeedOperatorReturnCeed(op), CEED_DEBUG_COLOR_SUCCESS, "Falling back to /gpu/cuda/ref CeedOperator due to code generation issue");
123f6eafd79SJeremy L Thompson       CeedCallBackend(CeedOperatorGetFallback(op, &op_fallback));
124f6eafd79SJeremy L Thompson       CeedCallBackend(CeedOperatorApplyAdd(op_fallback, input_vec, output_vec, request));
125f6eafd79SJeremy L Thompson       return CEED_ERROR_SUCCESS;
126f6eafd79SJeremy L Thompson     }
127f6eafd79SJeremy L Thompson   }
128f6eafd79SJeremy L Thompson 
129c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
130c11e12f4SJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &cuda_data));
131c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
132c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
133c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
134c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
135ddae5012SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
136c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
137c11e12f4SJeremy L Thompson 
138241a4b83SYohann   // Input vectors
1399e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1402b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1419e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
1429e201c85SYohann       data->fields.inputs[i] = NULL;
143241a4b83SYohann     } else {
144681d0ea7SJeremy L Thompson       bool       is_active;
145681d0ea7SJeremy L Thompson       CeedVector vec;
146681d0ea7SJeremy L Thompson 
147241a4b83SYohann       // Get input vector
1482b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
149681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
150681d0ea7SJeremy L Thompson       if (is_active) vec = input_vec;
1512b730f8bSJeremy L Thompson       CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->fields.inputs[i]));
152681d0ea7SJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec));
153241a4b83SYohann     }
154241a4b83SYohann   }
155241a4b83SYohann 
156241a4b83SYohann   // Output vectors
1579e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1582b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1599e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
1609e201c85SYohann       data->fields.outputs[i] = NULL;
161241a4b83SYohann     } else {
162681d0ea7SJeremy L Thompson       bool       is_active;
163681d0ea7SJeremy L Thompson       CeedVector vec;
164681d0ea7SJeremy L Thompson 
165241a4b83SYohann       // Get output vector
1662b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
167681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
168681d0ea7SJeremy L Thompson       if (is_active) vec = output_vec;
1699e201c85SYohann       output_vecs[i] = vec;
1703b2939feSjeremylt       // Check for multiple output modes
1713b2939feSjeremylt       CeedInt index = -1;
172ca735530SJeremy L Thompson 
1733b2939feSjeremylt       for (CeedInt j = 0; j < i; j++) {
1749e201c85SYohann         if (vec == output_vecs[j]) {
1753b2939feSjeremylt           index = j;
1763b2939feSjeremylt           break;
1773b2939feSjeremylt         }
1783b2939feSjeremylt       }
1793b2939feSjeremylt       if (index == -1) {
1802b730f8bSJeremy L Thompson         CeedCallBackend(CeedVectorGetArray(vec, CEED_MEM_DEVICE, &data->fields.outputs[i]));
1813b2939feSjeremylt       } else {
1829e201c85SYohann         data->fields.outputs[i] = data->fields.outputs[index];
1833b2939feSjeremylt       }
184681d0ea7SJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec));
185241a4b83SYohann     }
186241a4b83SYohann   }
187241a4b83SYohann 
1888b97b69aSJeremy L Thompson   // Point coordinates, if needed
1898b97b69aSJeremy L Thompson   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
1908b97b69aSJeremy L Thompson   if (is_at_points) {
1918b97b69aSJeremy L Thompson     // Coords
1928b97b69aSJeremy L Thompson     CeedVector vec;
1938b97b69aSJeremy L Thompson 
1948b97b69aSJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec));
1958b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->points.coords));
1968b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&vec));
1978b97b69aSJeremy L Thompson 
1988b97b69aSJeremy L Thompson     // Points per elem
1998b97b69aSJeremy L Thompson     if (num_elem != data->points.num_elem) {
2008b97b69aSJeremy L Thompson       CeedInt            *points_per_elem;
2018b97b69aSJeremy L Thompson       const CeedInt       num_bytes   = num_elem * sizeof(CeedInt);
2028b97b69aSJeremy L Thompson       CeedElemRestriction rstr_points = NULL;
2038b97b69aSJeremy L Thompson 
2048b97b69aSJeremy L Thompson       data->points.num_elem = num_elem;
2058b97b69aSJeremy L Thompson       CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
2068b97b69aSJeremy L Thompson       CeedCallBackend(CeedCalloc(num_elem, &points_per_elem));
2078b97b69aSJeremy L Thompson       for (CeedInt e = 0; e < num_elem; e++) {
2088b97b69aSJeremy L Thompson         CeedInt num_points_elem;
2098b97b69aSJeremy L Thompson 
2108b97b69aSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumPointsInElement(rstr_points, e, &num_points_elem));
2118b97b69aSJeremy L Thompson         points_per_elem[e] = num_points_elem;
2128b97b69aSJeremy L Thompson       }
2138b97b69aSJeremy L Thompson       if (data->points.num_per_elem) CeedCallCuda(ceed, cudaFree((void **)data->points.num_per_elem));
2148b97b69aSJeremy L Thompson       CeedCallCuda(ceed, cudaMalloc((void **)&data->points.num_per_elem, num_bytes));
2158b97b69aSJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy((void *)data->points.num_per_elem, points_per_elem, num_bytes, cudaMemcpyHostToDevice));
2168b97b69aSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
2178b97b69aSJeremy L Thompson       CeedCallBackend(CeedFree(&points_per_elem));
2188b97b69aSJeremy L Thompson     }
2198b97b69aSJeremy L Thompson   }
2208b97b69aSJeremy L Thompson 
221777ff853SJeremy L Thompson   // Get context data
2222b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetInnerContextData(qf, CEED_MEM_DEVICE, &qf_data->d_c));
223241a4b83SYohann 
224241a4b83SYohann   // Apply operator
2258b97b69aSJeremy L Thompson   void         *opargs[]  = {(void *)&num_elem, &qf_data->d_c, &data->indices, &data->fields, &data->B, &data->G, &data->W, &data->points};
226241a4b83SYohann   const CeedInt dim       = data->dim;
2279e201c85SYohann   const CeedInt Q_1d      = data->Q_1d;
2289e201c85SYohann   const CeedInt P_1d      = data->max_P_1d;
2299e201c85SYohann   const CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
2308b97b69aSJeremy L Thompson   int           max_threads_per_block, min_grid_size, grid;
231ca735530SJeremy L Thompson 
232dc007f05SJeremy L Thompson   CeedCallBackend(CeedOperatorHasTensorBases(op, &is_tensor));
2332b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuOccupancyMaxPotentialBlockSize(&min_grid_size, &max_threads_per_block, data->op, dynamicSMemSize, 0, 0x10000));
234dc007f05SJeremy L Thompson   int block[3] = {thread_1d, ((!is_tensor || dim == 1) ? 1 : thread_1d), -1};
235ca735530SJeremy L Thompson 
236f82027a4SJeremy L Thompson   if (is_tensor) {
23793f4dbf1SSebastian Grimberg     CeedCallBackend(BlockGridCalculate(num_elem, min_grid_size / cuda_data->device_prop.multiProcessorCount, max_threads_per_block,
2382b730f8bSJeremy L Thompson                                        cuda_data->device_prop.maxThreadsDim[2], cuda_data->device_prop.warpSize, block, &grid));
239f82027a4SJeremy L Thompson   } else {
240f82027a4SJeremy L Thompson     CeedInt elems_per_block = CeedIntMin(cuda_data->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1));
241f82027a4SJeremy L Thompson 
242f82027a4SJeremy L Thompson     grid     = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
243f82027a4SJeremy L Thompson     block[2] = elems_per_block;
244f82027a4SJeremy L Thompson   }
24539532cebSJed Brown   CeedInt shared_mem = block[0] * block[1] * block[2] * sizeof(CeedScalar);
246ca735530SJeremy L Thompson 
247*8d12f40eSJeremy L Thompson   CeedCallBackend(CeedTryRunKernelDimShared_Cuda(ceed, data->op, grid, block[0], block[1], block[2], shared_mem, &is_run_good, opargs));
248241a4b83SYohann 
249241a4b83SYohann   // Restore input arrays
2509e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
2512b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
2529e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
253241a4b83SYohann     } else {
254681d0ea7SJeremy L Thompson       bool       is_active;
255681d0ea7SJeremy L Thompson       CeedVector vec;
256681d0ea7SJeremy L Thompson 
2572b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
258681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
259681d0ea7SJeremy L Thompson       if (is_active) vec = input_vec;
2602b730f8bSJeremy L Thompson       CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->fields.inputs[i]));
261681d0ea7SJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec));
262241a4b83SYohann     }
263241a4b83SYohann   }
264241a4b83SYohann 
265241a4b83SYohann   // Restore output arrays
2669e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
2672b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
2689e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
269241a4b83SYohann     } else {
270681d0ea7SJeremy L Thompson       bool       is_active;
271681d0ea7SJeremy L Thompson       CeedVector vec;
272681d0ea7SJeremy L Thompson 
2732b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
274681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
275681d0ea7SJeremy L Thompson       if (is_active) vec = output_vec;
2763b2939feSjeremylt       // Check for multiple output modes
2773b2939feSjeremylt       CeedInt index = -1;
2788b97b69aSJeremy L Thompson 
2793b2939feSjeremylt       for (CeedInt j = 0; j < i; j++) {
2809e201c85SYohann         if (vec == output_vecs[j]) {
2813b2939feSjeremylt           index = j;
2823b2939feSjeremylt           break;
2833b2939feSjeremylt         }
2843b2939feSjeremylt       }
2853b2939feSjeremylt       if (index == -1) {
2862b730f8bSJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArray(vec, &data->fields.outputs[i]));
287241a4b83SYohann       }
288681d0ea7SJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec));
289241a4b83SYohann     }
2903b2939feSjeremylt   }
291777ff853SJeremy L Thompson 
2928b97b69aSJeremy L Thompson   // Restore point coordinates, if needed
2938b97b69aSJeremy L Thompson   if (is_at_points) {
2948b97b69aSJeremy L Thompson     CeedVector vec;
2958b97b69aSJeremy L Thompson 
2968b97b69aSJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec));
2978b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->points.coords));
2988b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&vec));
2998b97b69aSJeremy L Thompson   }
3008b97b69aSJeremy L Thompson 
301777ff853SJeremy L Thompson   // Restore context data
3022b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionRestoreInnerContextData(qf, &qf_data->d_c));
303ddae5012SJeremy L Thompson 
304ddae5012SJeremy L Thompson   // Cleanup
3059bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
306c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionDestroy(&qf));
307ddae5012SJeremy L Thompson 
308ddae5012SJeremy L Thompson   // Fallback if run was bad (out of resources)
309*8d12f40eSJeremy L Thompson   if (!is_run_good) {
310ddae5012SJeremy L Thompson     CeedOperator op_fallback;
311ddae5012SJeremy L Thompson 
312ddae5012SJeremy L Thompson     data->use_fallback = true;
313ddae5012SJeremy L Thompson     CeedDebug256(CeedOperatorReturnCeed(op), CEED_DEBUG_COLOR_SUCCESS, "Falling back to /gpu/cuda/ref CeedOperator due to kernel execution issue");
314ddae5012SJeremy L Thompson     CeedCallBackend(CeedOperatorGetFallback(op, &op_fallback));
315ddae5012SJeremy L Thompson     CeedCallBackend(CeedOperatorApplyAdd(op_fallback, input_vec, output_vec, request));
316ddae5012SJeremy L Thompson     return CEED_ERROR_SUCCESS;
317ddae5012SJeremy L Thompson   }
318e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
319241a4b83SYohann }
320241a4b83SYohann 
321ab213215SJeremy L Thompson //------------------------------------------------------------------------------
322ab213215SJeremy L Thompson // Create operator
323ab213215SJeremy L Thompson //------------------------------------------------------------------------------
324241a4b83SYohann int CeedOperatorCreate_Cuda_gen(CeedOperator op) {
325241a4b83SYohann   Ceed                   ceed;
326241a4b83SYohann   CeedOperator_Cuda_gen *impl;
327241a4b83SYohann 
328ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
3292b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl));
3302b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetData(op, impl));
3312b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda_gen));
3322b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda_gen));
3339bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
334e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
335241a4b83SYohann }
3366aa95790SJeremy L Thompson 
337ab213215SJeremy L Thompson //------------------------------------------------------------------------------
338