xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-gen/ceed-cuda-gen-operator.c (revision 8b97b69a34af3dd9e1bc9b8fd8651212448dc9ec)
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>
11*8b97b69aSJeremy L Thompson #include <cuda.h>
12*8b97b69aSJeremy L Thompson #include <cuda_runtime.h>
133d576824SJeremy L Thompson #include <stddef.h>
142b730f8bSJeremy L Thompson 
1549aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
166d69246aSJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
172b730f8bSJeremy L Thompson #include "ceed-cuda-gen-operator-build.h"
182b730f8bSJeremy L Thompson #include "ceed-cuda-gen.h"
19241a4b83SYohann 
20ab213215SJeremy L Thompson //------------------------------------------------------------------------------
21ab213215SJeremy L Thompson // Destroy operator
22ab213215SJeremy L Thompson //------------------------------------------------------------------------------
23241a4b83SYohann static int CeedOperatorDestroy_Cuda_gen(CeedOperator op) {
24*8b97b69aSJeremy L Thompson   Ceed                   ceed;
25241a4b83SYohann   CeedOperator_Cuda_gen *impl;
26ca735530SJeremy L Thompson 
27*8b97b69aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
282b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
29*8b97b69aSJeremy L Thompson   if (impl->points.num_per_elem) CeedCallCuda(ceed, cudaFree((void **)impl->points.num_per_elem));
302b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl));
31*8b97b69aSJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
32e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
33241a4b83SYohann }
34241a4b83SYohann 
352b730f8bSJeremy L Thompson static int Waste(int threads_per_sm, int warp_size, int threads_per_elem, int elems_per_block) {
3639532cebSJed Brown   int useful_threads_per_block = threads_per_elem * elems_per_block;
3739532cebSJed Brown   // round up to nearest multiple of warp_size
38b2165e7aSSebastian Grimberg   int block_size    = CeedDivUpInt(useful_threads_per_block, warp_size) * warp_size;
3939532cebSJed Brown   int blocks_per_sm = threads_per_sm / block_size;
4039532cebSJed Brown   return threads_per_sm - useful_threads_per_block * blocks_per_sm;
4139532cebSJed Brown }
4239532cebSJed Brown 
43ea61e9acSJeremy L Thompson // Choose the least wasteful block size constrained by blocks_per_sm of max_threads_per_block.
4439532cebSJed Brown //
45ea61e9acSJeremy L Thompson // The x and y part of block[] contains per-element sizes (specified on input) while the z part is number of elements.
4639532cebSJed Brown //
47ea61e9acSJeremy L Thompson // Problem setting: we'd like to make occupancy high with relatively few inactive threads. CUDA (cuOccupancyMaxPotentialBlockSize) can tell us how
4839532cebSJed Brown // many threads can run.
4939532cebSJed Brown //
50ea61e9acSJeremy L Thompson // Note that full occupancy sometimes can't be achieved by one thread block.
51ea61e9acSJeremy L Thompson // For example, an SM might support 1536 threads in total, but only 1024 within a single thread block.
52ea61e9acSJeremy 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
53ea61e9acSJeremy 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
54ea61e9acSJeremy L Thompson // relatively low occupancy and smaller thread blocks, but we solve a reasonably general problem here. Empirically, we find that blocks bigger than
55ea61e9acSJeremy L Thompson // about 256 have higher latency and worse load balancing when the number of elements is modest.
5639532cebSJed Brown //
57ea61e9acSJeremy 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).
58ea61e9acSJeremy L Thompson // They also have a lot of __syncthreads(), which is another point against excessively large thread blocks.
59ea61e9acSJeremy L Thompson // Suppose I have elements with 7x7x7 quadrature points.
60ea61e9acSJeremy L Thompson // This will loop over the last dimension, so we have 7*7=49 threads per element.
61ea61e9acSJeremy L Thompson // Suppose we have two elements = 2*49=98 useful threads.
62ea61e9acSJeremy L Thompson // CUDA schedules in units of full warps (32 threads), so 128 CUDA hardware threads are effectively committed to that block.
63ea61e9acSJeremy L Thompson // Now suppose cuOccupancyMaxPotentialBlockSize returned 352.
64ea61e9acSJeremy 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
65ea61e9acSJeremy L Thompson // hardware threads).
6639532cebSJed Brown //
67ea61e9acSJeremy 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.
68ea61e9acSJeremy L Thompson // Alternatively, we could pack a single block of 7 elements (2*49=343 useful threads) into the 354 slots.
69ea61e9acSJeremy L Thompson // The latter has the least "waste", but __syncthreads() over-synchronizes and it might not pay off relative to smaller blocks.
702b730f8bSJeremy 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],
712b730f8bSJeremy L Thompson                               int *grid) {
7239532cebSJed Brown   const int threads_per_sm   = blocks_per_sm * max_threads_per_block;
7339532cebSJed Brown   const int threads_per_elem = block[0] * block[1];
7439532cebSJed Brown   int       elems_per_block  = 1;
7539532cebSJed Brown   int       waste            = Waste(threads_per_sm, warp_size, threads_per_elem, 1);
76ca735530SJeremy L Thompson 
772b730f8bSJeremy L Thompson   for (int i = 2; i <= CeedIntMin(max_threads_per_block / threads_per_elem, num_elem); i++) {
7839532cebSJed Brown     int i_waste = Waste(threads_per_sm, warp_size, threads_per_elem, i);
79ca735530SJeremy L Thompson 
80ea61e9acSJeremy 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
8139532cebSJed Brown     // waste as a smaller one, go ahead and prefer the smaller block.
8239532cebSJed Brown     if (i_waste < waste || (i_waste == waste && threads_per_elem * i <= 128)) {
8339532cebSJed Brown       elems_per_block = i;
8439532cebSJed Brown       waste           = i_waste;
8539532cebSJed Brown     }
8639532cebSJed Brown   }
87ea61e9acSJeremy 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
88ea61e9acSJeremy L Thompson   // check before setting the z-dimension size of the block.
8913516544Snbeams   block[2] = CeedIntMin(elems_per_block, max_threads_z);
90b2165e7aSSebastian Grimberg   *grid    = CeedDivUpInt(num_elem, elems_per_block);
9139532cebSJed Brown   return CEED_ERROR_SUCCESS;
9239532cebSJed Brown }
9339532cebSJed Brown 
94ea61e9acSJeremy L Thompson // callback for cuOccupancyMaxPotentialBlockSize, providing the amount of dynamic shared memory required for a thread block of size threads.
9539532cebSJed Brown static size_t dynamicSMemSize(int threads) { return threads * sizeof(CeedScalar); }
9639532cebSJed Brown 
97ab213215SJeremy L Thompson //------------------------------------------------------------------------------
98ab213215SJeremy L Thompson // Apply and add to output
99ab213215SJeremy L Thompson //------------------------------------------------------------------------------
1002b730f8bSJeremy L Thompson static int CeedOperatorApplyAdd_Cuda_gen(CeedOperator op, CeedVector input_vec, CeedVector output_vec, CeedRequest *request) {
101*8b97b69aSJeremy L Thompson   bool                    is_at_points;
102241a4b83SYohann   Ceed                    ceed;
10339532cebSJed Brown   Ceed_Cuda              *cuda_data;
104ca735530SJeremy L Thompson   CeedInt                 num_elem, num_input_fields, num_output_fields;
105ca735530SJeremy L Thompson   CeedEvalMode            eval_mode;
106ca735530SJeremy L Thompson   CeedVector              output_vecs[CEED_FIELD_MAX] = {NULL};
107ca735530SJeremy L Thompson   CeedQFunctionField     *qf_input_fields, *qf_output_fields;
108241a4b83SYohann   CeedQFunction_Cuda_gen *qf_data;
109ca735530SJeremy L Thompson   CeedQFunction           qf;
110ca735530SJeremy L Thompson   CeedOperatorField      *op_input_fields, *op_output_fields;
111ca735530SJeremy L Thompson   CeedOperator_Cuda_gen  *data;
112ca735530SJeremy L Thompson 
113f6eafd79SJeremy L Thompson   // Check for tensor-product bases
114f6eafd79SJeremy L Thompson   {
115f6eafd79SJeremy L Thompson     bool has_tensor_bases;
116f6eafd79SJeremy L Thompson 
117f6eafd79SJeremy L Thompson     CeedCallBackend(CeedOperatorHasTensorBases(op, &has_tensor_bases));
1184535e697SJeremy L Thompson     // -- Fallback to ref if not all bases are tensor-product
119f6eafd79SJeremy L Thompson     if (!has_tensor_bases) {
120f6eafd79SJeremy L Thompson       CeedOperator op_fallback;
121f6eafd79SJeremy L Thompson 
122c11e12f4SJeremy L Thompson       CeedDebug256(CeedOperatorReturnCeed(op), CEED_DEBUG_COLOR_SUCCESS, "Falling back to /gpu/cuda/ref CeedOperator due to non-tensor bases");
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));
135c11e12f4SJeremy 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   // Creation of the operator
139eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernel_Cuda_gen(op));
140241a4b83SYohann 
141241a4b83SYohann   // Input vectors
1429e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1432b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1449e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
1459e201c85SYohann       data->fields.inputs[i] = NULL;
146241a4b83SYohann     } else {
147681d0ea7SJeremy L Thompson       bool       is_active;
148681d0ea7SJeremy L Thompson       CeedVector vec;
149681d0ea7SJeremy L Thompson 
150241a4b83SYohann       // Get input vector
1512b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
152681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
153681d0ea7SJeremy L Thompson       if (is_active) vec = input_vec;
1542b730f8bSJeremy L Thompson       CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->fields.inputs[i]));
155681d0ea7SJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec));
156241a4b83SYohann     }
157241a4b83SYohann   }
158241a4b83SYohann 
159241a4b83SYohann   // Output vectors
1609e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1612b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1629e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
1639e201c85SYohann       data->fields.outputs[i] = NULL;
164241a4b83SYohann     } else {
165681d0ea7SJeremy L Thompson       bool       is_active;
166681d0ea7SJeremy L Thompson       CeedVector vec;
167681d0ea7SJeremy L Thompson 
168241a4b83SYohann       // Get output vector
1692b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
170681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
171681d0ea7SJeremy L Thompson       if (is_active) vec = output_vec;
1729e201c85SYohann       output_vecs[i] = vec;
1733b2939feSjeremylt       // Check for multiple output modes
1743b2939feSjeremylt       CeedInt index = -1;
175ca735530SJeremy L Thompson 
1763b2939feSjeremylt       for (CeedInt j = 0; j < i; j++) {
1779e201c85SYohann         if (vec == output_vecs[j]) {
1783b2939feSjeremylt           index = j;
1793b2939feSjeremylt           break;
1803b2939feSjeremylt         }
1813b2939feSjeremylt       }
1823b2939feSjeremylt       if (index == -1) {
1832b730f8bSJeremy L Thompson         CeedCallBackend(CeedVectorGetArray(vec, CEED_MEM_DEVICE, &data->fields.outputs[i]));
1843b2939feSjeremylt       } else {
1859e201c85SYohann         data->fields.outputs[i] = data->fields.outputs[index];
1863b2939feSjeremylt       }
187681d0ea7SJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec));
188241a4b83SYohann     }
189241a4b83SYohann   }
190241a4b83SYohann 
191*8b97b69aSJeremy L Thompson   // Point coordinates, if needed
192*8b97b69aSJeremy L Thompson   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
193*8b97b69aSJeremy L Thompson   if (is_at_points) {
194*8b97b69aSJeremy L Thompson     // Coords
195*8b97b69aSJeremy L Thompson     CeedVector vec;
196*8b97b69aSJeremy L Thompson 
197*8b97b69aSJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec));
198*8b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->points.coords));
199*8b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&vec));
200*8b97b69aSJeremy L Thompson 
201*8b97b69aSJeremy L Thompson     // Points per elem
202*8b97b69aSJeremy L Thompson     if (num_elem != data->points.num_elem) {
203*8b97b69aSJeremy L Thompson       CeedInt            *points_per_elem;
204*8b97b69aSJeremy L Thompson       const CeedInt       num_bytes   = num_elem * sizeof(CeedInt);
205*8b97b69aSJeremy L Thompson       CeedElemRestriction rstr_points = NULL;
206*8b97b69aSJeremy L Thompson 
207*8b97b69aSJeremy L Thompson       data->points.num_elem = num_elem;
208*8b97b69aSJeremy L Thompson       CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
209*8b97b69aSJeremy L Thompson       CeedCallBackend(CeedCalloc(num_elem, &points_per_elem));
210*8b97b69aSJeremy L Thompson       for (CeedInt e = 0; e < num_elem; e++) {
211*8b97b69aSJeremy L Thompson         CeedInt num_points_elem;
212*8b97b69aSJeremy L Thompson 
213*8b97b69aSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumPointsInElement(rstr_points, e, &num_points_elem));
214*8b97b69aSJeremy L Thompson         points_per_elem[e] = num_points_elem;
215*8b97b69aSJeremy L Thompson       }
216*8b97b69aSJeremy L Thompson       if (data->points.num_per_elem) CeedCallCuda(ceed, cudaFree((void **)data->points.num_per_elem));
217*8b97b69aSJeremy L Thompson       CeedCallCuda(ceed, cudaMalloc((void **)&data->points.num_per_elem, num_bytes));
218*8b97b69aSJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy((void *)data->points.num_per_elem, points_per_elem, num_bytes, cudaMemcpyHostToDevice));
219*8b97b69aSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
220*8b97b69aSJeremy L Thompson       CeedCallBackend(CeedFree(&points_per_elem));
221*8b97b69aSJeremy L Thompson     }
222*8b97b69aSJeremy L Thompson   }
223*8b97b69aSJeremy L Thompson 
224777ff853SJeremy L Thompson   // Get context data
2252b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetInnerContextData(qf, CEED_MEM_DEVICE, &qf_data->d_c));
226241a4b83SYohann 
227241a4b83SYohann   // Apply operator
228*8b97b69aSJeremy L Thompson   void         *opargs[]  = {(void *)&num_elem, &qf_data->d_c, &data->indices, &data->fields, &data->B, &data->G, &data->W, &data->points};
229241a4b83SYohann   const CeedInt dim       = data->dim;
2309e201c85SYohann   const CeedInt Q_1d      = data->Q_1d;
2319e201c85SYohann   const CeedInt P_1d      = data->max_P_1d;
2329e201c85SYohann   const CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
233*8b97b69aSJeremy L Thompson   int           max_threads_per_block, min_grid_size, grid;
234ca735530SJeremy L Thompson 
2352b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuOccupancyMaxPotentialBlockSize(&min_grid_size, &max_threads_per_block, data->op, dynamicSMemSize, 0, 0x10000));
236*8b97b69aSJeremy L Thompson   int block[3] = {thread_1d, dim < 2 ? 1 : thread_1d, -1};
237ca735530SJeremy L Thompson 
23893f4dbf1SSebastian Grimberg   CeedCallBackend(BlockGridCalculate(num_elem, min_grid_size / cuda_data->device_prop.multiProcessorCount, max_threads_per_block,
2392b730f8bSJeremy L Thompson                                      cuda_data->device_prop.maxThreadsDim[2], cuda_data->device_prop.warpSize, block, &grid));
24039532cebSJed Brown   CeedInt shared_mem = block[0] * block[1] * block[2] * sizeof(CeedScalar);
241ca735530SJeremy L Thompson 
242eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->op, grid, block[0], block[1], block[2], shared_mem, opargs));
243241a4b83SYohann 
244241a4b83SYohann   // Restore input arrays
2459e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
2462b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
2479e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
248241a4b83SYohann     } else {
249681d0ea7SJeremy L Thompson       bool       is_active;
250681d0ea7SJeremy L Thompson       CeedVector vec;
251681d0ea7SJeremy L Thompson 
2522b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
253681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
254681d0ea7SJeremy L Thompson       if (is_active) vec = input_vec;
2552b730f8bSJeremy L Thompson       CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->fields.inputs[i]));
256681d0ea7SJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec));
257241a4b83SYohann     }
258241a4b83SYohann   }
259241a4b83SYohann 
260241a4b83SYohann   // Restore output arrays
2619e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
2622b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
2639e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
264241a4b83SYohann     } else {
265681d0ea7SJeremy L Thompson       bool       is_active;
266681d0ea7SJeremy L Thompson       CeedVector vec;
267681d0ea7SJeremy L Thompson 
2682b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
269681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
270681d0ea7SJeremy L Thompson       if (is_active) vec = output_vec;
2713b2939feSjeremylt       // Check for multiple output modes
2723b2939feSjeremylt       CeedInt index = -1;
273*8b97b69aSJeremy L Thompson 
2743b2939feSjeremylt       for (CeedInt j = 0; j < i; j++) {
2759e201c85SYohann         if (vec == output_vecs[j]) {
2763b2939feSjeremylt           index = j;
2773b2939feSjeremylt           break;
2783b2939feSjeremylt         }
2793b2939feSjeremylt       }
2803b2939feSjeremylt       if (index == -1) {
2812b730f8bSJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArray(vec, &data->fields.outputs[i]));
282241a4b83SYohann       }
283681d0ea7SJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec));
284241a4b83SYohann     }
2853b2939feSjeremylt   }
286777ff853SJeremy L Thompson 
287*8b97b69aSJeremy L Thompson   // Restore point coordinates, if needed
288*8b97b69aSJeremy L Thompson   if (is_at_points) {
289*8b97b69aSJeremy L Thompson     CeedVector vec;
290*8b97b69aSJeremy L Thompson 
291*8b97b69aSJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec));
292*8b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->points.coords));
293*8b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&vec));
294*8b97b69aSJeremy L Thompson   }
295*8b97b69aSJeremy L Thompson 
296777ff853SJeremy L Thompson   // Restore context data
2972b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionRestoreInnerContextData(qf, &qf_data->d_c));
2989bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
299c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionDestroy(&qf));
300e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
301241a4b83SYohann }
302241a4b83SYohann 
303ab213215SJeremy L Thompson //------------------------------------------------------------------------------
304ab213215SJeremy L Thompson // Create operator
305ab213215SJeremy L Thompson //------------------------------------------------------------------------------
306241a4b83SYohann int CeedOperatorCreate_Cuda_gen(CeedOperator op) {
307241a4b83SYohann   Ceed                   ceed;
308241a4b83SYohann   CeedOperator_Cuda_gen *impl;
309241a4b83SYohann 
310ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
3112b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl));
3122b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetData(op, impl));
3132b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda_gen));
3142b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda_gen));
3159bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
316e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
317241a4b83SYohann }
3186aa95790SJeremy L Thompson 
319ab213215SJeremy L Thompson //------------------------------------------------------------------------------
320