xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-gen/ceed-cuda-gen-operator.c (revision dc007f05648c670dfdc3e42fab8d6c1219c0afbb)
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>
14*dc007f05SJeremy 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*dc007f05SJeremy L Thompson   bool                    is_at_points, is_tensor;
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 
114*dc007f05SJeremy L Thompson   // Check for shared bases
115*dc007f05SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
116f6eafd79SJeremy L Thompson   {
117*dc007f05SJeremy L Thompson     bool has_shared_bases = true, is_all_tensor = true, is_all_nontensor = true;
118f6eafd79SJeremy L Thompson 
119*dc007f05SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
120*dc007f05SJeremy L Thompson       CeedBasis basis;
121*dc007f05SJeremy L Thompson 
122*dc007f05SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
123*dc007f05SJeremy L Thompson       if (basis != CEED_BASIS_NONE) {
124*dc007f05SJeremy L Thompson         bool        is_tensor = true;
125*dc007f05SJeremy L Thompson         const char *resource;
126*dc007f05SJeremy L Thompson         char       *resource_root;
127*dc007f05SJeremy L Thompson         Ceed        basis_ceed;
128*dc007f05SJeremy L Thompson 
129*dc007f05SJeremy L Thompson         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
130*dc007f05SJeremy L Thompson         is_all_tensor    &= is_tensor;
131*dc007f05SJeremy L Thompson         is_all_nontensor &= !is_tensor;
132*dc007f05SJeremy L Thompson         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
133*dc007f05SJeremy L Thompson         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
134*dc007f05SJeremy L Thompson         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
135*dc007f05SJeremy L Thompson         has_shared_bases &= !strcmp(resource_root, "/gpu/cuda/shared");
136*dc007f05SJeremy L Thompson         CeedCallBackend(CeedFree(&resource_root));
137*dc007f05SJeremy L Thompson         CeedCallBackend(CeedDestroy(&basis_ceed));
138*dc007f05SJeremy L Thompson       }
139*dc007f05SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis));
140*dc007f05SJeremy L Thompson     }
141*dc007f05SJeremy L Thompson 
142*dc007f05SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
143*dc007f05SJeremy L Thompson       CeedBasis basis;
144*dc007f05SJeremy L Thompson 
145*dc007f05SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
146*dc007f05SJeremy L Thompson       if (basis != CEED_BASIS_NONE) {
147*dc007f05SJeremy L Thompson         bool        is_tensor = true;
148*dc007f05SJeremy L Thompson         const char *resource;
149*dc007f05SJeremy L Thompson         char       *resource_root;
150*dc007f05SJeremy L Thompson         Ceed        basis_ceed;
151*dc007f05SJeremy L Thompson 
152*dc007f05SJeremy L Thompson         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
153*dc007f05SJeremy L Thompson         is_all_tensor    &= is_tensor;
154*dc007f05SJeremy L Thompson         is_all_nontensor &= !is_tensor;
155*dc007f05SJeremy L Thompson 
156*dc007f05SJeremy L Thompson         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
157*dc007f05SJeremy L Thompson         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
158*dc007f05SJeremy L Thompson         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
159*dc007f05SJeremy L Thompson         has_shared_bases &= !strcmp(resource_root, "/gpu/cuda/shared");
160*dc007f05SJeremy L Thompson         CeedCallBackend(CeedFree(&resource_root));
161*dc007f05SJeremy L Thompson         CeedCallBackend(CeedDestroy(&basis_ceed));
162*dc007f05SJeremy L Thompson       }
163*dc007f05SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis));
164*dc007f05SJeremy L Thompson     }
165*dc007f05SJeremy L Thompson     // -- Fallback to ref if not all bases are shared
166*dc007f05SJeremy L Thompson     if (!has_shared_bases || (!is_all_tensor && !is_all_nontensor)) {
167f6eafd79SJeremy L Thompson       CeedOperator op_fallback;
168f6eafd79SJeremy L Thompson 
169*dc007f05SJeremy L Thompson       CeedDebug256(CeedOperatorReturnCeed(op), CEED_DEBUG_COLOR_SUCCESS, "Falling back to /gpu/cuda/ref CeedOperator due to large non-tensor bases");
170f6eafd79SJeremy L Thompson       CeedCallBackend(CeedOperatorGetFallback(op, &op_fallback));
171f6eafd79SJeremy L Thompson       CeedCallBackend(CeedOperatorApplyAdd(op_fallback, input_vec, output_vec, request));
172f6eafd79SJeremy L Thompson       return CEED_ERROR_SUCCESS;
173f6eafd79SJeremy L Thompson     }
174f6eafd79SJeremy L Thompson   }
175f6eafd79SJeremy L Thompson 
176c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
177c11e12f4SJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &cuda_data));
178c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
179c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
180c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
181c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
182c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
183c11e12f4SJeremy L Thompson 
184241a4b83SYohann   // Creation of the operator
185eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernel_Cuda_gen(op));
186241a4b83SYohann 
187241a4b83SYohann   // Input vectors
1889e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1892b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1909e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
1919e201c85SYohann       data->fields.inputs[i] = NULL;
192241a4b83SYohann     } else {
193681d0ea7SJeremy L Thompson       bool       is_active;
194681d0ea7SJeremy L Thompson       CeedVector vec;
195681d0ea7SJeremy L Thompson 
196241a4b83SYohann       // Get input vector
1972b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
198681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
199681d0ea7SJeremy L Thompson       if (is_active) vec = input_vec;
2002b730f8bSJeremy L Thompson       CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->fields.inputs[i]));
201681d0ea7SJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec));
202241a4b83SYohann     }
203241a4b83SYohann   }
204241a4b83SYohann 
205241a4b83SYohann   // Output vectors
2069e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
2072b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
2089e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
2099e201c85SYohann       data->fields.outputs[i] = NULL;
210241a4b83SYohann     } else {
211681d0ea7SJeremy L Thompson       bool       is_active;
212681d0ea7SJeremy L Thompson       CeedVector vec;
213681d0ea7SJeremy L Thompson 
214241a4b83SYohann       // Get output vector
2152b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
216681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
217681d0ea7SJeremy L Thompson       if (is_active) vec = output_vec;
2189e201c85SYohann       output_vecs[i] = vec;
2193b2939feSjeremylt       // Check for multiple output modes
2203b2939feSjeremylt       CeedInt index = -1;
221ca735530SJeremy L Thompson 
2223b2939feSjeremylt       for (CeedInt j = 0; j < i; j++) {
2239e201c85SYohann         if (vec == output_vecs[j]) {
2243b2939feSjeremylt           index = j;
2253b2939feSjeremylt           break;
2263b2939feSjeremylt         }
2273b2939feSjeremylt       }
2283b2939feSjeremylt       if (index == -1) {
2292b730f8bSJeremy L Thompson         CeedCallBackend(CeedVectorGetArray(vec, CEED_MEM_DEVICE, &data->fields.outputs[i]));
2303b2939feSjeremylt       } else {
2319e201c85SYohann         data->fields.outputs[i] = data->fields.outputs[index];
2323b2939feSjeremylt       }
233681d0ea7SJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec));
234241a4b83SYohann     }
235241a4b83SYohann   }
236241a4b83SYohann 
2378b97b69aSJeremy L Thompson   // Point coordinates, if needed
2388b97b69aSJeremy L Thompson   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
2398b97b69aSJeremy L Thompson   if (is_at_points) {
2408b97b69aSJeremy L Thompson     // Coords
2418b97b69aSJeremy L Thompson     CeedVector vec;
2428b97b69aSJeremy L Thompson 
2438b97b69aSJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec));
2448b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->points.coords));
2458b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&vec));
2468b97b69aSJeremy L Thompson 
2478b97b69aSJeremy L Thompson     // Points per elem
2488b97b69aSJeremy L Thompson     if (num_elem != data->points.num_elem) {
2498b97b69aSJeremy L Thompson       CeedInt            *points_per_elem;
2508b97b69aSJeremy L Thompson       const CeedInt       num_bytes   = num_elem * sizeof(CeedInt);
2518b97b69aSJeremy L Thompson       CeedElemRestriction rstr_points = NULL;
2528b97b69aSJeremy L Thompson 
2538b97b69aSJeremy L Thompson       data->points.num_elem = num_elem;
2548b97b69aSJeremy L Thompson       CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
2558b97b69aSJeremy L Thompson       CeedCallBackend(CeedCalloc(num_elem, &points_per_elem));
2568b97b69aSJeremy L Thompson       for (CeedInt e = 0; e < num_elem; e++) {
2578b97b69aSJeremy L Thompson         CeedInt num_points_elem;
2588b97b69aSJeremy L Thompson 
2598b97b69aSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumPointsInElement(rstr_points, e, &num_points_elem));
2608b97b69aSJeremy L Thompson         points_per_elem[e] = num_points_elem;
2618b97b69aSJeremy L Thompson       }
2628b97b69aSJeremy L Thompson       if (data->points.num_per_elem) CeedCallCuda(ceed, cudaFree((void **)data->points.num_per_elem));
2638b97b69aSJeremy L Thompson       CeedCallCuda(ceed, cudaMalloc((void **)&data->points.num_per_elem, num_bytes));
2648b97b69aSJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy((void *)data->points.num_per_elem, points_per_elem, num_bytes, cudaMemcpyHostToDevice));
2658b97b69aSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
2668b97b69aSJeremy L Thompson       CeedCallBackend(CeedFree(&points_per_elem));
2678b97b69aSJeremy L Thompson     }
2688b97b69aSJeremy L Thompson   }
2698b97b69aSJeremy L Thompson 
270777ff853SJeremy L Thompson   // Get context data
2712b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetInnerContextData(qf, CEED_MEM_DEVICE, &qf_data->d_c));
272241a4b83SYohann 
273241a4b83SYohann   // Apply operator
2748b97b69aSJeremy L Thompson   void         *opargs[]  = {(void *)&num_elem, &qf_data->d_c, &data->indices, &data->fields, &data->B, &data->G, &data->W, &data->points};
275241a4b83SYohann   const CeedInt dim       = data->dim;
2769e201c85SYohann   const CeedInt Q_1d      = data->Q_1d;
2779e201c85SYohann   const CeedInt P_1d      = data->max_P_1d;
2789e201c85SYohann   const CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
2798b97b69aSJeremy L Thompson   int           max_threads_per_block, min_grid_size, grid;
280ca735530SJeremy L Thompson 
281*dc007f05SJeremy L Thompson   CeedCallBackend(CeedOperatorHasTensorBases(op, &is_tensor));
2822b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuOccupancyMaxPotentialBlockSize(&min_grid_size, &max_threads_per_block, data->op, dynamicSMemSize, 0, 0x10000));
283*dc007f05SJeremy L Thompson   int block[3] = {thread_1d, ((!is_tensor || dim == 1) ? 1 : thread_1d), -1};
284ca735530SJeremy L Thompson 
28593f4dbf1SSebastian Grimberg   CeedCallBackend(BlockGridCalculate(num_elem, min_grid_size / cuda_data->device_prop.multiProcessorCount, max_threads_per_block,
2862b730f8bSJeremy L Thompson                                      cuda_data->device_prop.maxThreadsDim[2], cuda_data->device_prop.warpSize, block, &grid));
28739532cebSJed Brown   CeedInt shared_mem = block[0] * block[1] * block[2] * sizeof(CeedScalar);
288ca735530SJeremy L Thompson 
289eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->op, grid, block[0], block[1], block[2], shared_mem, opargs));
290241a4b83SYohann 
291241a4b83SYohann   // Restore input arrays
2929e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
2932b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
2949e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
295241a4b83SYohann     } else {
296681d0ea7SJeremy L Thompson       bool       is_active;
297681d0ea7SJeremy L Thompson       CeedVector vec;
298681d0ea7SJeremy L Thompson 
2992b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
300681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
301681d0ea7SJeremy L Thompson       if (is_active) vec = input_vec;
3022b730f8bSJeremy L Thompson       CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->fields.inputs[i]));
303681d0ea7SJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec));
304241a4b83SYohann     }
305241a4b83SYohann   }
306241a4b83SYohann 
307241a4b83SYohann   // Restore output arrays
3089e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
3092b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
3109e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
311241a4b83SYohann     } else {
312681d0ea7SJeremy L Thompson       bool       is_active;
313681d0ea7SJeremy L Thompson       CeedVector vec;
314681d0ea7SJeremy L Thompson 
3152b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
316681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
317681d0ea7SJeremy L Thompson       if (is_active) vec = output_vec;
3183b2939feSjeremylt       // Check for multiple output modes
3193b2939feSjeremylt       CeedInt index = -1;
3208b97b69aSJeremy L Thompson 
3213b2939feSjeremylt       for (CeedInt j = 0; j < i; j++) {
3229e201c85SYohann         if (vec == output_vecs[j]) {
3233b2939feSjeremylt           index = j;
3243b2939feSjeremylt           break;
3253b2939feSjeremylt         }
3263b2939feSjeremylt       }
3273b2939feSjeremylt       if (index == -1) {
3282b730f8bSJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArray(vec, &data->fields.outputs[i]));
329241a4b83SYohann       }
330681d0ea7SJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec));
331241a4b83SYohann     }
3323b2939feSjeremylt   }
333777ff853SJeremy L Thompson 
3348b97b69aSJeremy L Thompson   // Restore point coordinates, if needed
3358b97b69aSJeremy L Thompson   if (is_at_points) {
3368b97b69aSJeremy L Thompson     CeedVector vec;
3378b97b69aSJeremy L Thompson 
3388b97b69aSJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec));
3398b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->points.coords));
3408b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&vec));
3418b97b69aSJeremy L Thompson   }
3428b97b69aSJeremy L Thompson 
343777ff853SJeremy L Thompson   // Restore context data
3442b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionRestoreInnerContextData(qf, &qf_data->d_c));
3459bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
346c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionDestroy(&qf));
347e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
348241a4b83SYohann }
349241a4b83SYohann 
350ab213215SJeremy L Thompson //------------------------------------------------------------------------------
351ab213215SJeremy L Thompson // Create operator
352ab213215SJeremy L Thompson //------------------------------------------------------------------------------
353241a4b83SYohann int CeedOperatorCreate_Cuda_gen(CeedOperator op) {
354241a4b83SYohann   Ceed                   ceed;
355241a4b83SYohann   CeedOperator_Cuda_gen *impl;
356241a4b83SYohann 
357ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
3582b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl));
3592b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetData(op, impl));
3602b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda_gen));
3612b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda_gen));
3629bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
363e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
364241a4b83SYohann }
3656aa95790SJeremy L Thompson 
366ab213215SJeremy L Thompson //------------------------------------------------------------------------------
367