xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-gen/ceed-cuda-gen-operator.c (revision 915834c9f1e582e3fdfc87db6b4fa4e010d293bb)
1d275d636SJeremy L Thompson // Copyright (c) 2017-2025, 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));
308b7d3340SJeremy L Thompson   if (impl->module) CeedCallCuda(ceed, cuModuleUnload(impl->module));
310183ed61SJeremy L Thompson   if (impl->module_assemble_full) CeedCallCuda(ceed, cuModuleUnload(impl->module_assemble_full));
320183ed61SJeremy L Thompson   if (impl->module_assemble_diagonal) CeedCallCuda(ceed, cuModuleUnload(impl->module_assemble_diagonal));
338b97b69aSJeremy L Thompson   if (impl->points.num_per_elem) CeedCallCuda(ceed, cudaFree((void **)impl->points.num_per_elem));
342b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl));
358b97b69aSJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
36e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
37241a4b83SYohann }
38241a4b83SYohann 
392b730f8bSJeremy L Thompson static int Waste(int threads_per_sm, int warp_size, int threads_per_elem, int elems_per_block) {
4039532cebSJed Brown   int useful_threads_per_block = threads_per_elem * elems_per_block;
4139532cebSJed Brown   // round up to nearest multiple of warp_size
42b2165e7aSSebastian Grimberg   int block_size    = CeedDivUpInt(useful_threads_per_block, warp_size) * warp_size;
4339532cebSJed Brown   int blocks_per_sm = threads_per_sm / block_size;
4439532cebSJed Brown   return threads_per_sm - useful_threads_per_block * blocks_per_sm;
4539532cebSJed Brown }
4639532cebSJed Brown 
47ea61e9acSJeremy L Thompson // Choose the least wasteful block size constrained by blocks_per_sm of max_threads_per_block.
4839532cebSJed Brown //
49ea61e9acSJeremy L Thompson // The x and y part of block[] contains per-element sizes (specified on input) while the z part is number of elements.
5039532cebSJed Brown //
51ea61e9acSJeremy L Thompson // Problem setting: we'd like to make occupancy high with relatively few inactive threads. CUDA (cuOccupancyMaxPotentialBlockSize) can tell us how
5239532cebSJed Brown // many threads can run.
5339532cebSJed Brown //
54ea61e9acSJeremy L Thompson // Note that full occupancy sometimes can't be achieved by one thread block.
55ea61e9acSJeremy L Thompson // For example, an SM might support 1536 threads in total, but only 1024 within a single thread block.
56ea61e9acSJeremy 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
57ea61e9acSJeremy 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
58ea61e9acSJeremy L Thompson // relatively low occupancy and smaller thread blocks, but we solve a reasonably general problem here. Empirically, we find that blocks bigger than
59ea61e9acSJeremy L Thompson // about 256 have higher latency and worse load balancing when the number of elements is modest.
6039532cebSJed Brown //
61ea61e9acSJeremy 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).
62ea61e9acSJeremy L Thompson // They also have a lot of __syncthreads(), which is another point against excessively large thread blocks.
63ea61e9acSJeremy L Thompson // Suppose I have elements with 7x7x7 quadrature points.
64ea61e9acSJeremy L Thompson // This will loop over the last dimension, so we have 7*7=49 threads per element.
65ea61e9acSJeremy L Thompson // Suppose we have two elements = 2*49=98 useful threads.
66ea61e9acSJeremy L Thompson // CUDA schedules in units of full warps (32 threads), so 128 CUDA hardware threads are effectively committed to that block.
67ea61e9acSJeremy L Thompson // Now suppose cuOccupancyMaxPotentialBlockSize returned 352.
68ea61e9acSJeremy 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
69ea61e9acSJeremy L Thompson // hardware threads).
7039532cebSJed Brown //
71ea61e9acSJeremy 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.
72ea61e9acSJeremy L Thompson // Alternatively, we could pack a single block of 7 elements (2*49=343 useful threads) into the 354 slots.
73ea61e9acSJeremy L Thompson // The latter has the least "waste", but __syncthreads() over-synchronizes and it might not pay off relative to smaller blocks.
742b730f8bSJeremy 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],
752b730f8bSJeremy L Thompson                               int *grid) {
7639532cebSJed Brown   const int threads_per_sm   = blocks_per_sm * max_threads_per_block;
7739532cebSJed Brown   const int threads_per_elem = block[0] * block[1];
7839532cebSJed Brown   int       elems_per_block  = 1;
7939532cebSJed Brown   int       waste            = Waste(threads_per_sm, warp_size, threads_per_elem, 1);
80ca735530SJeremy L Thompson 
812b730f8bSJeremy L Thompson   for (int i = 2; i <= CeedIntMin(max_threads_per_block / threads_per_elem, num_elem); i++) {
8239532cebSJed Brown     int i_waste = Waste(threads_per_sm, warp_size, threads_per_elem, i);
83ca735530SJeremy L Thompson 
84ea61e9acSJeremy 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
8539532cebSJed Brown     // waste as a smaller one, go ahead and prefer the smaller block.
8639532cebSJed Brown     if (i_waste < waste || (i_waste == waste && threads_per_elem * i <= 128)) {
8739532cebSJed Brown       elems_per_block = i;
8839532cebSJed Brown       waste           = i_waste;
8939532cebSJed Brown     }
9039532cebSJed Brown   }
91ea61e9acSJeremy 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
92ea61e9acSJeremy L Thompson   // check before setting the z-dimension size of the block.
9313516544Snbeams   block[2] = CeedIntMin(elems_per_block, max_threads_z);
94b2165e7aSSebastian Grimberg   *grid    = CeedDivUpInt(num_elem, elems_per_block);
9539532cebSJed Brown   return CEED_ERROR_SUCCESS;
9639532cebSJed Brown }
9739532cebSJed Brown 
98ea61e9acSJeremy L Thompson // callback for cuOccupancyMaxPotentialBlockSize, providing the amount of dynamic shared memory required for a thread block of size threads.
9939532cebSJed Brown static size_t dynamicSMemSize(int threads) { return threads * sizeof(CeedScalar); }
10039532cebSJed Brown 
101ab213215SJeremy L Thompson //------------------------------------------------------------------------------
102ab213215SJeremy L Thompson // Apply and add to output
103ab213215SJeremy L Thompson //------------------------------------------------------------------------------
104e9c76bddSJeremy L Thompson static int CeedOperatorApplyAddCore_Cuda_gen(CeedOperator op, CUstream stream, const CeedScalar *input_arr, CeedScalar *output_arr, bool *is_run_good,
105ea04d07fSJeremy L Thompson                                              CeedRequest *request) {
106ea04d07fSJeremy L Thompson   bool                    is_at_points, is_tensor;
107241a4b83SYohann   Ceed                    ceed;
10839532cebSJed Brown   Ceed_Cuda              *cuda_data;
109ca735530SJeremy L Thompson   CeedInt                 num_elem, num_input_fields, num_output_fields;
110ca735530SJeremy L Thompson   CeedEvalMode            eval_mode;
111ca735530SJeremy L Thompson   CeedQFunctionField     *qf_input_fields, *qf_output_fields;
112241a4b83SYohann   CeedQFunction_Cuda_gen *qf_data;
113ca735530SJeremy L Thompson   CeedQFunction           qf;
114ca735530SJeremy L Thompson   CeedOperatorField      *op_input_fields, *op_output_fields;
115ca735530SJeremy L Thompson   CeedOperator_Cuda_gen  *data;
116ca735530SJeremy L Thompson 
117ea04d07fSJeremy L Thompson   // Build the operator kernel
118ea04d07fSJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernel_Cuda_gen(op, is_run_good));
119ea04d07fSJeremy L Thompson   if (!(*is_run_good)) return CEED_ERROR_SUCCESS;
120f6eafd79SJeremy L Thompson 
121c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
122c11e12f4SJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &cuda_data));
123c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
124c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
125c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
126c11e12f4SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
127ddae5012SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
128c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
129c11e12f4SJeremy L Thompson 
130241a4b83SYohann   // Input vectors
1319e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1322b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1339e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
1349e201c85SYohann       data->fields.inputs[i] = NULL;
135241a4b83SYohann     } else {
136681d0ea7SJeremy L Thompson       bool       is_active;
137681d0ea7SJeremy L Thompson       CeedVector vec;
138681d0ea7SJeremy L Thompson 
139241a4b83SYohann       // Get input vector
1402b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
141681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
142ea04d07fSJeremy L Thompson       if (is_active) data->fields.inputs[i] = input_arr;
143ea04d07fSJeremy L Thompson       else CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->fields.inputs[i]));
144ea04d07fSJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec));
145241a4b83SYohann     }
146241a4b83SYohann   }
147241a4b83SYohann 
148241a4b83SYohann   // Output vectors
1499e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1502b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1519e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
1529e201c85SYohann       data->fields.outputs[i] = NULL;
153241a4b83SYohann     } else {
154681d0ea7SJeremy L Thompson       bool       is_active;
155681d0ea7SJeremy L Thompson       CeedVector vec;
156681d0ea7SJeremy L Thompson 
157241a4b83SYohann       // Get output vector
1582b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
159681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
160ea04d07fSJeremy L Thompson       if (is_active) data->fields.outputs[i] = output_arr;
1610c8fbeedSJeremy L Thompson       else CeedCallBackend(CeedVectorGetArray(vec, CEED_MEM_DEVICE, &data->fields.outputs[i]));
162ea04d07fSJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec));
163241a4b83SYohann     }
164241a4b83SYohann   }
165241a4b83SYohann 
1668b97b69aSJeremy L Thompson   // Point coordinates, if needed
1678b97b69aSJeremy L Thompson   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
1688b97b69aSJeremy L Thompson   if (is_at_points) {
1698b97b69aSJeremy L Thompson     // Coords
1708b97b69aSJeremy L Thompson     CeedVector vec;
1718b97b69aSJeremy L Thompson 
1728b97b69aSJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec));
1738b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->points.coords));
1748b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&vec));
1758b97b69aSJeremy L Thompson 
1768b97b69aSJeremy L Thompson     // Points per elem
1778b97b69aSJeremy L Thompson     if (num_elem != data->points.num_elem) {
1788b97b69aSJeremy L Thompson       CeedInt            *points_per_elem;
1798b97b69aSJeremy L Thompson       const CeedInt       num_bytes   = num_elem * sizeof(CeedInt);
1808b97b69aSJeremy L Thompson       CeedElemRestriction rstr_points = NULL;
1818b97b69aSJeremy L Thompson 
1828b97b69aSJeremy L Thompson       data->points.num_elem = num_elem;
1838b97b69aSJeremy L Thompson       CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
1848b97b69aSJeremy L Thompson       CeedCallBackend(CeedCalloc(num_elem, &points_per_elem));
1858b97b69aSJeremy L Thompson       for (CeedInt e = 0; e < num_elem; e++) {
1868b97b69aSJeremy L Thompson         CeedInt num_points_elem;
1878b97b69aSJeremy L Thompson 
1888b97b69aSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumPointsInElement(rstr_points, e, &num_points_elem));
1898b97b69aSJeremy L Thompson         points_per_elem[e] = num_points_elem;
1908b97b69aSJeremy L Thompson       }
1918b97b69aSJeremy L Thompson       if (data->points.num_per_elem) CeedCallCuda(ceed, cudaFree((void **)data->points.num_per_elem));
1928b97b69aSJeremy L Thompson       CeedCallCuda(ceed, cudaMalloc((void **)&data->points.num_per_elem, num_bytes));
1938b97b69aSJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy((void *)data->points.num_per_elem, points_per_elem, num_bytes, cudaMemcpyHostToDevice));
1948b97b69aSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
1958b97b69aSJeremy L Thompson       CeedCallBackend(CeedFree(&points_per_elem));
1968b97b69aSJeremy L Thompson     }
1978b97b69aSJeremy L Thompson   }
1988b97b69aSJeremy L Thompson 
199777ff853SJeremy L Thompson   // Get context data
2002b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetInnerContextData(qf, CEED_MEM_DEVICE, &qf_data->d_c));
201241a4b83SYohann 
202241a4b83SYohann   // Apply operator
2038b97b69aSJeremy L Thompson   void *opargs[] = {(void *)&num_elem, &qf_data->d_c, &data->indices, &data->fields, &data->B, &data->G, &data->W, &data->points};
2048b97b69aSJeremy L Thompson   int   max_threads_per_block, min_grid_size, grid;
205ca735530SJeremy L Thompson 
206dc007f05SJeremy L Thompson   CeedCallBackend(CeedOperatorHasTensorBases(op, &is_tensor));
2072b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuOccupancyMaxPotentialBlockSize(&min_grid_size, &max_threads_per_block, data->op, dynamicSMemSize, 0, 0x10000));
208a61b1c91SJeremy L Thompson   int block[3] = {data->thread_1d, ((!is_tensor || data->dim == 1) ? 1 : data->thread_1d), -1};
209ca735530SJeremy L Thompson 
210f82027a4SJeremy L Thompson   if (is_tensor) {
21190c30374SJeremy L Thompson     CeedCallBackend(BlockGridCalculate(num_elem, min_grid_size / cuda_data->device_prop.multiProcessorCount, is_at_points ? 1 : max_threads_per_block,
2122b730f8bSJeremy L Thompson                                        cuda_data->device_prop.maxThreadsDim[2], cuda_data->device_prop.warpSize, block, &grid));
213f82027a4SJeremy L Thompson   } else {
214a61b1c91SJeremy L Thompson     CeedInt elems_per_block = CeedIntMin(cuda_data->device_prop.maxThreadsDim[2], CeedIntMax(512 / data->thread_1d, 1));
215f82027a4SJeremy L Thompson 
216f82027a4SJeremy L Thompson     grid     = num_elem / elems_per_block + (num_elem % elems_per_block > 0);
217f82027a4SJeremy L Thompson     block[2] = elems_per_block;
218f82027a4SJeremy L Thompson   }
21939532cebSJed Brown   CeedInt shared_mem = block[0] * block[1] * block[2] * sizeof(CeedScalar);
220ca735530SJeremy L Thompson 
221e9c76bddSJeremy L Thompson   CeedCallBackend(CeedTryRunKernelDimShared_Cuda(ceed, data->op, stream, grid, block[0], block[1], block[2], shared_mem, is_run_good, opargs));
222241a4b83SYohann 
223241a4b83SYohann   // Restore input arrays
2249e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
2252b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
2269e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
227241a4b83SYohann     } else {
228681d0ea7SJeremy L Thompson       bool       is_active;
229681d0ea7SJeremy L Thompson       CeedVector vec;
230681d0ea7SJeremy L Thompson 
2312b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
232681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
233ea04d07fSJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->fields.inputs[i]));
234ea04d07fSJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec));
235241a4b83SYohann     }
236241a4b83SYohann   }
237241a4b83SYohann 
238241a4b83SYohann   // Restore output arrays
2399e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
2402b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
2419e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
242241a4b83SYohann     } else {
243681d0ea7SJeremy L Thompson       bool       is_active;
244681d0ea7SJeremy L Thompson       CeedVector vec;
245681d0ea7SJeremy L Thompson 
2462b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
247681d0ea7SJeremy L Thompson       is_active = vec == CEED_VECTOR_ACTIVE;
248ea04d07fSJeremy L Thompson       if (!is_active) CeedCallBackend(CeedVectorRestoreArray(vec, &data->fields.outputs[i]));
249ea04d07fSJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec));
250241a4b83SYohann     }
2513b2939feSjeremylt   }
252777ff853SJeremy L Thompson 
2538b97b69aSJeremy L Thompson   // Restore point coordinates, if needed
2548b97b69aSJeremy L Thompson   if (is_at_points) {
2558b97b69aSJeremy L Thompson     CeedVector vec;
2568b97b69aSJeremy L Thompson 
2578b97b69aSJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec));
2588b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->points.coords));
2598b97b69aSJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&vec));
2608b97b69aSJeremy L Thompson   }
2618b97b69aSJeremy L Thompson 
262777ff853SJeremy L Thompson   // Restore context data
2632b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionRestoreInnerContextData(qf, &qf_data->d_c));
264ddae5012SJeremy L Thompson 
265ddae5012SJeremy L Thompson   // Cleanup
2669bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
267c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionDestroy(&qf));
268ea04d07fSJeremy L Thompson   if (!(*is_run_good)) data->use_fallback = true;
269ea04d07fSJeremy L Thompson   return CEED_ERROR_SUCCESS;
270ea04d07fSJeremy L Thompson }
271ddae5012SJeremy L Thompson 
272ea04d07fSJeremy L Thompson static int CeedOperatorApplyAdd_Cuda_gen(CeedOperator op, CeedVector input_vec, CeedVector output_vec, CeedRequest *request) {
273ea04d07fSJeremy L Thompson   bool              is_run_good = false;
274ea04d07fSJeremy L Thompson   const CeedScalar *input_arr   = NULL;
275ea04d07fSJeremy L Thompson   CeedScalar       *output_arr  = NULL;
276ea04d07fSJeremy L Thompson 
277ea04d07fSJeremy L Thompson   // Try to run kernel
278ea04d07fSJeremy L Thompson   if (input_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(input_vec, CEED_MEM_DEVICE, &input_arr));
279ea04d07fSJeremy L Thompson   if (output_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArray(output_vec, CEED_MEM_DEVICE, &output_arr));
280e9c76bddSJeremy L Thompson   CeedCallBackend(CeedOperatorApplyAddCore_Cuda_gen(op, NULL, input_arr, output_arr, &is_run_good, request));
281ea04d07fSJeremy L Thompson   if (input_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorRestoreArrayRead(input_vec, &input_arr));
282ea04d07fSJeremy L Thompson   if (output_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorRestoreArray(output_vec, &output_arr));
283ea04d07fSJeremy L Thompson 
284ea04d07fSJeremy L Thompson   // Fallback on unsuccessful run
2858d12f40eSJeremy L Thompson   if (!is_run_good) {
286ddae5012SJeremy L Thompson     CeedOperator op_fallback;
287ddae5012SJeremy L Thompson 
288ea04d07fSJeremy L Thompson     CeedDebug256(CeedOperatorReturnCeed(op), CEED_DEBUG_COLOR_SUCCESS, "Falling back to /gpu/cuda/ref CeedOperator");
289ddae5012SJeremy L Thompson     CeedCallBackend(CeedOperatorGetFallback(op, &op_fallback));
290ddae5012SJeremy L Thompson     CeedCallBackend(CeedOperatorApplyAdd(op_fallback, input_vec, output_vec, request));
291ddae5012SJeremy L Thompson   }
292e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
293241a4b83SYohann }
294241a4b83SYohann 
295c99afcd8SJeremy L Thompson static int CeedOperatorApplyAddComposite_Cuda_gen(CeedOperator op, CeedVector input_vec, CeedVector output_vec, CeedRequest *request) {
296c99afcd8SJeremy L Thompson   bool              is_run_good[CEED_COMPOSITE_MAX] = {false};
297c99afcd8SJeremy L Thompson   CeedInt           num_suboperators;
298c99afcd8SJeremy L Thompson   const CeedScalar *input_arr  = NULL;
299c99afcd8SJeremy L Thompson   CeedScalar       *output_arr = NULL;
300087855afSJeremy L Thompson   Ceed              ceed;
301c99afcd8SJeremy L Thompson   CeedOperator     *sub_operators;
302c99afcd8SJeremy L Thompson 
303087855afSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
304c99afcd8SJeremy L Thompson   CeedCall(CeedCompositeOperatorGetNumSub(op, &num_suboperators));
305c99afcd8SJeremy L Thompson   CeedCall(CeedCompositeOperatorGetSubList(op, &sub_operators));
306c99afcd8SJeremy L Thompson   if (input_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(input_vec, CEED_MEM_DEVICE, &input_arr));
307c99afcd8SJeremy L Thompson   if (output_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArray(output_vec, CEED_MEM_DEVICE, &output_arr));
308c99afcd8SJeremy L Thompson   for (CeedInt i = 0; i < num_suboperators; i++) {
309c99afcd8SJeremy L Thompson     CeedInt num_elem = 0;
310c99afcd8SJeremy L Thompson 
311c99afcd8SJeremy L Thompson     CeedCall(CeedOperatorGetNumElements(sub_operators[i], &num_elem));
312c99afcd8SJeremy L Thompson     if (num_elem > 0) {
313087855afSJeremy L Thompson       cudaStream_t stream = NULL;
314087855afSJeremy L Thompson 
315087855afSJeremy L Thompson       CeedCallCuda(ceed, cudaStreamCreate(&stream));
316087855afSJeremy L Thompson       CeedCallBackend(CeedOperatorApplyAddCore_Cuda_gen(sub_operators[i], stream, input_arr, output_arr, &is_run_good[i], request));
317087855afSJeremy L Thompson       CeedCallCuda(ceed, cudaStreamDestroy(stream));
318c99afcd8SJeremy L Thompson     }
319c99afcd8SJeremy L Thompson   }
320c99afcd8SJeremy L Thompson   if (input_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorRestoreArrayRead(input_vec, &input_arr));
321c99afcd8SJeremy L Thompson   if (output_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorRestoreArray(output_vec, &output_arr));
322087855afSJeremy L Thompson   CeedCallCuda(ceed, cudaDeviceSynchronize());
323c99afcd8SJeremy L Thompson 
324c99afcd8SJeremy L Thompson   // Fallback on unsuccessful run
325c99afcd8SJeremy L Thompson   for (CeedInt i = 0; i < num_suboperators; i++) {
326c99afcd8SJeremy L Thompson     if (!is_run_good[i]) {
327c99afcd8SJeremy L Thompson       CeedOperator op_fallback;
328c99afcd8SJeremy L Thompson 
329087855afSJeremy L Thompson       CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "Falling back to /gpu/cuda/ref CeedOperator");
330c99afcd8SJeremy L Thompson       CeedCallBackend(CeedOperatorGetFallback(sub_operators[i], &op_fallback));
331c99afcd8SJeremy L Thompson       CeedCallBackend(CeedOperatorApplyAdd(op_fallback, input_vec, output_vec, request));
332c99afcd8SJeremy L Thompson     }
333c99afcd8SJeremy L Thompson   }
334087855afSJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
335c99afcd8SJeremy L Thompson   return CEED_ERROR_SUCCESS;
336c99afcd8SJeremy L Thompson }
337c99afcd8SJeremy L Thompson 
338ab213215SJeremy L Thompson //------------------------------------------------------------------------------
3390183ed61SJeremy L Thompson // AtPoints diagonal assembly
3400183ed61SJeremy L Thompson //------------------------------------------------------------------------------
3410183ed61SJeremy L Thompson static int CeedOperatorLinearAssembleAddDiagonalAtPoints_Cuda_gen(CeedOperator op, CeedVector assembled, CeedRequest *request) {
3420183ed61SJeremy L Thompson   Ceed                   ceed;
3430183ed61SJeremy L Thompson   CeedOperator_Cuda_gen *data;
3440183ed61SJeremy L Thompson 
3450183ed61SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
3460183ed61SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
3470183ed61SJeremy L Thompson 
3480183ed61SJeremy L Thompson   // Build the assembly kernel
3490183ed61SJeremy L Thompson   if (!data->assemble_diagonal && !data->use_assembly_fallback) {
3500183ed61SJeremy L Thompson     bool                     is_build_good = false;
3510183ed61SJeremy L Thompson     CeedInt                  num_active_bases_in, num_active_bases_out;
3520183ed61SJeremy L Thompson     CeedOperatorAssemblyData assembly_data;
3530183ed61SJeremy L Thompson 
3540183ed61SJeremy L Thompson     CeedCallBackend(CeedOperatorGetOperatorAssemblyData(op, &assembly_data));
3550183ed61SJeremy L Thompson     CeedCallBackend(
3560183ed61SJeremy L Thompson         CeedOperatorAssemblyDataGetEvalModes(assembly_data, &num_active_bases_in, NULL, NULL, NULL, &num_active_bases_out, NULL, NULL, NULL, NULL));
3570183ed61SJeremy L Thompson     if (num_active_bases_in == num_active_bases_out) {
3580183ed61SJeremy L Thompson       CeedCallBackend(CeedOperatorBuildKernel_Cuda_gen(op, &is_build_good));
3590183ed61SJeremy L Thompson       if (is_build_good) CeedCallBackend(CeedOperatorBuildKernelDiagonalAssemblyAtPoints_Cuda_gen(op, &is_build_good));
3600183ed61SJeremy L Thompson     }
3610183ed61SJeremy L Thompson     if (!is_build_good) data->use_assembly_fallback = true;
3620183ed61SJeremy L Thompson   }
3630183ed61SJeremy L Thompson 
3640183ed61SJeremy L Thompson   // Try assembly
3650183ed61SJeremy L Thompson   if (!data->use_assembly_fallback) {
3660183ed61SJeremy L Thompson     bool                    is_run_good = true;
3670183ed61SJeremy L Thompson     Ceed_Cuda              *cuda_data;
3680183ed61SJeremy L Thompson     CeedInt                 num_elem, num_input_fields, num_output_fields;
3690183ed61SJeremy L Thompson     CeedEvalMode            eval_mode;
3700183ed61SJeremy L Thompson     CeedScalar             *assembled_array;
3710183ed61SJeremy L Thompson     CeedQFunctionField     *qf_input_fields, *qf_output_fields;
3720183ed61SJeremy L Thompson     CeedQFunction_Cuda_gen *qf_data;
3730183ed61SJeremy L Thompson     CeedQFunction           qf;
3740183ed61SJeremy L Thompson     CeedOperatorField      *op_input_fields, *op_output_fields;
3750183ed61SJeremy L Thompson 
3760183ed61SJeremy L Thompson     CeedCallBackend(CeedGetData(ceed, &cuda_data));
3770183ed61SJeremy L Thompson     CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
3780183ed61SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
3790183ed61SJeremy L Thompson     CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
3800183ed61SJeremy L Thompson     CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
3810183ed61SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
3820183ed61SJeremy L Thompson 
3830183ed61SJeremy L Thompson     // Input vectors
3840183ed61SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
3850183ed61SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
3860183ed61SJeremy L Thompson       if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
3870183ed61SJeremy L Thompson         data->fields.inputs[i] = NULL;
3880183ed61SJeremy L Thompson       } else {
3890183ed61SJeremy L Thompson         bool       is_active;
3900183ed61SJeremy L Thompson         CeedVector vec;
3910183ed61SJeremy L Thompson 
3920183ed61SJeremy L Thompson         // Get input vector
3930183ed61SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
3940183ed61SJeremy L Thompson         is_active = vec == CEED_VECTOR_ACTIVE;
3950183ed61SJeremy L Thompson         if (is_active) data->fields.inputs[i] = NULL;
3960183ed61SJeremy L Thompson         else CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->fields.inputs[i]));
3970183ed61SJeremy L Thompson         CeedCallBackend(CeedVectorDestroy(&vec));
3980183ed61SJeremy L Thompson       }
3990183ed61SJeremy L Thompson     }
4000183ed61SJeremy L Thompson 
4010183ed61SJeremy L Thompson     // Point coordinates
4020183ed61SJeremy L Thompson     {
4030183ed61SJeremy L Thompson       CeedVector vec;
4040183ed61SJeremy L Thompson 
4050183ed61SJeremy L Thompson       CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec));
4060183ed61SJeremy L Thompson       CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->points.coords));
4070183ed61SJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec));
4080183ed61SJeremy L Thompson 
4090183ed61SJeremy L Thompson       // Points per elem
4100183ed61SJeremy L Thompson       if (num_elem != data->points.num_elem) {
4110183ed61SJeremy L Thompson         CeedInt            *points_per_elem;
4120183ed61SJeremy L Thompson         const CeedInt       num_bytes   = num_elem * sizeof(CeedInt);
4130183ed61SJeremy L Thompson         CeedElemRestriction rstr_points = NULL;
4140183ed61SJeremy L Thompson 
4150183ed61SJeremy L Thompson         data->points.num_elem = num_elem;
4160183ed61SJeremy L Thompson         CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
4170183ed61SJeremy L Thompson         CeedCallBackend(CeedCalloc(num_elem, &points_per_elem));
4180183ed61SJeremy L Thompson         for (CeedInt e = 0; e < num_elem; e++) {
4190183ed61SJeremy L Thompson           CeedInt num_points_elem;
4200183ed61SJeremy L Thompson 
4210183ed61SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetNumPointsInElement(rstr_points, e, &num_points_elem));
4220183ed61SJeremy L Thompson           points_per_elem[e] = num_points_elem;
4230183ed61SJeremy L Thompson         }
4240183ed61SJeremy L Thompson         if (data->points.num_per_elem) CeedCallCuda(ceed, cudaFree((void **)data->points.num_per_elem));
4250183ed61SJeremy L Thompson         CeedCallCuda(ceed, cudaMalloc((void **)&data->points.num_per_elem, num_bytes));
4260183ed61SJeremy L Thompson         CeedCallCuda(ceed, cudaMemcpy((void *)data->points.num_per_elem, points_per_elem, num_bytes, cudaMemcpyHostToDevice));
4270183ed61SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
4280183ed61SJeremy L Thompson         CeedCallBackend(CeedFree(&points_per_elem));
4290183ed61SJeremy L Thompson       }
4300183ed61SJeremy L Thompson     }
4310183ed61SJeremy L Thompson 
4320183ed61SJeremy L Thompson     // Get context data
4330183ed61SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetInnerContextData(qf, CEED_MEM_DEVICE, &qf_data->d_c));
4340183ed61SJeremy L Thompson 
4350183ed61SJeremy L Thompson     // Assembly array
4360183ed61SJeremy L Thompson     CeedCallBackend(CeedVectorGetArray(assembled, CEED_MEM_DEVICE, &assembled_array));
4370183ed61SJeremy L Thompson 
4380183ed61SJeremy L Thompson     // Assemble diagonal
4390183ed61SJeremy L Thompson     void *opargs[] = {(void *)&num_elem, &qf_data->d_c, &data->indices, &data->fields, &data->B, &data->G, &data->W, &data->points, &assembled_array};
4400183ed61SJeremy L Thompson     int   max_threads_per_block, min_grid_size, grid;
4410183ed61SJeremy L Thompson 
4420183ed61SJeremy L Thompson     CeedCallCuda(ceed, cuOccupancyMaxPotentialBlockSize(&min_grid_size, &max_threads_per_block, data->op, dynamicSMemSize, 0, 0x10000));
4430183ed61SJeremy L Thompson     int block[3] = {data->thread_1d, (data->dim == 1 ? 1 : data->thread_1d), -1};
4440183ed61SJeremy L Thompson 
4450183ed61SJeremy L Thompson     CeedCallBackend(BlockGridCalculate(num_elem, min_grid_size / cuda_data->device_prop.multiProcessorCount, 1,
4460183ed61SJeremy L Thompson                                        cuda_data->device_prop.maxThreadsDim[2], cuda_data->device_prop.warpSize, block, &grid));
4470183ed61SJeremy L Thompson     CeedInt shared_mem = block[0] * block[1] * block[2] * sizeof(CeedScalar);
4480183ed61SJeremy L Thompson 
4490183ed61SJeremy L Thompson     CeedCallBackend(
4500183ed61SJeremy L Thompson         CeedTryRunKernelDimShared_Cuda(ceed, data->assemble_diagonal, NULL, grid, block[0], block[1], block[2], shared_mem, &is_run_good, opargs));
451*915834c9SZach Atkins     CeedCallCuda(ceed, cudaDeviceSynchronize());
4520183ed61SJeremy L Thompson 
4530183ed61SJeremy L Thompson     // Restore input arrays
4540183ed61SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
4550183ed61SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
4560183ed61SJeremy L Thompson       if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
4570183ed61SJeremy L Thompson       } else {
4580183ed61SJeremy L Thompson         bool       is_active;
4590183ed61SJeremy L Thompson         CeedVector vec;
4600183ed61SJeremy L Thompson 
4610183ed61SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
4620183ed61SJeremy L Thompson         is_active = vec == CEED_VECTOR_ACTIVE;
4630183ed61SJeremy L Thompson         if (!is_active) CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->fields.inputs[i]));
4640183ed61SJeremy L Thompson         CeedCallBackend(CeedVectorDestroy(&vec));
4650183ed61SJeremy L Thompson       }
4660183ed61SJeremy L Thompson     }
4670183ed61SJeremy L Thompson 
4680183ed61SJeremy L Thompson     // Restore point coordinates
4690183ed61SJeremy L Thompson     {
4700183ed61SJeremy L Thompson       CeedVector vec;
4710183ed61SJeremy L Thompson 
4720183ed61SJeremy L Thompson       CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec));
4730183ed61SJeremy L Thompson       CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->points.coords));
4740183ed61SJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec));
4750183ed61SJeremy L Thompson     }
4760183ed61SJeremy L Thompson 
4770183ed61SJeremy L Thompson     // Restore context data
4780183ed61SJeremy L Thompson     CeedCallBackend(CeedQFunctionRestoreInnerContextData(qf, &qf_data->d_c));
4790183ed61SJeremy L Thompson 
4800183ed61SJeremy L Thompson     // Restore assembly array
4810183ed61SJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArray(assembled, &assembled_array));
4820183ed61SJeremy L Thompson 
4830183ed61SJeremy L Thompson     // Cleanup
4840183ed61SJeremy L Thompson     CeedCallBackend(CeedQFunctionDestroy(&qf));
4850183ed61SJeremy L Thompson     if (!is_run_good) data->use_assembly_fallback = true;
4860183ed61SJeremy L Thompson   }
4870183ed61SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
4880183ed61SJeremy L Thompson 
4890183ed61SJeremy L Thompson   // Fallback, if needed
4900183ed61SJeremy L Thompson   if (data->use_assembly_fallback) {
4910183ed61SJeremy L Thompson     CeedOperator op_fallback;
4920183ed61SJeremy L Thompson 
4930183ed61SJeremy L Thompson     CeedDebug256(CeedOperatorReturnCeed(op), CEED_DEBUG_COLOR_SUCCESS, "Falling back to /gpu/cuda/ref CeedOperator");
4940183ed61SJeremy L Thompson     CeedCallBackend(CeedOperatorGetFallback(op, &op_fallback));
4950183ed61SJeremy L Thompson     CeedCallBackend(CeedOperatorLinearAssembleAddDiagonal(op_fallback, assembled, request));
4960183ed61SJeremy L Thompson     return CEED_ERROR_SUCCESS;
4970183ed61SJeremy L Thompson   }
4980183ed61SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4990183ed61SJeremy L Thompson }
5000183ed61SJeremy L Thompson 
5010183ed61SJeremy L Thompson //------------------------------------------------------------------------------
502*915834c9SZach Atkins // AtPoints full assembly
503*915834c9SZach Atkins //------------------------------------------------------------------------------
504*915834c9SZach Atkins static int CeedSingleOperatorAssembleAtPoints_Cuda_gen(CeedOperator op, CeedInt offset, CeedVector assembled) {
505*915834c9SZach Atkins   Ceed                   ceed;
506*915834c9SZach Atkins   CeedOperator_Cuda_gen *data;
507*915834c9SZach Atkins 
508*915834c9SZach Atkins   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
509*915834c9SZach Atkins   CeedCallBackend(CeedOperatorGetData(op, &data));
510*915834c9SZach Atkins 
511*915834c9SZach Atkins   // Build the assembly kernel
512*915834c9SZach Atkins   if (!data->assemble_full && !data->use_assembly_fallback) {
513*915834c9SZach Atkins     bool                     is_build_good = false;
514*915834c9SZach Atkins     CeedInt                  num_active_bases_in, num_active_bases_out;
515*915834c9SZach Atkins     CeedOperatorAssemblyData assembly_data;
516*915834c9SZach Atkins 
517*915834c9SZach Atkins     CeedCallBackend(CeedOperatorGetOperatorAssemblyData(op, &assembly_data));
518*915834c9SZach Atkins     CeedCallBackend(
519*915834c9SZach Atkins         CeedOperatorAssemblyDataGetEvalModes(assembly_data, &num_active_bases_in, NULL, NULL, NULL, &num_active_bases_out, NULL, NULL, NULL, NULL));
520*915834c9SZach Atkins     if (num_active_bases_in == num_active_bases_out) {
521*915834c9SZach Atkins       CeedCallBackend(CeedOperatorBuildKernel_Cuda_gen(op, &is_build_good));
522*915834c9SZach Atkins       if (is_build_good) CeedCallBackend(CeedOperatorBuildKernelFullAssemblyAtPoints_Cuda_gen(op, &is_build_good));
523*915834c9SZach Atkins     }
524*915834c9SZach Atkins     if (!is_build_good) data->use_assembly_fallback = true;
525*915834c9SZach Atkins   }
526*915834c9SZach Atkins 
527*915834c9SZach Atkins   // Try assembly
528*915834c9SZach Atkins   if (!data->use_assembly_fallback) {
529*915834c9SZach Atkins     bool                    is_run_good = true;
530*915834c9SZach Atkins     Ceed_Cuda              *cuda_data;
531*915834c9SZach Atkins     CeedInt                 num_elem, num_input_fields, num_output_fields;
532*915834c9SZach Atkins     CeedEvalMode            eval_mode;
533*915834c9SZach Atkins     CeedScalar             *assembled_array;
534*915834c9SZach Atkins     CeedQFunctionField     *qf_input_fields, *qf_output_fields;
535*915834c9SZach Atkins     CeedQFunction_Cuda_gen *qf_data;
536*915834c9SZach Atkins     CeedQFunction           qf;
537*915834c9SZach Atkins     CeedOperatorField      *op_input_fields, *op_output_fields;
538*915834c9SZach Atkins 
539*915834c9SZach Atkins     CeedCallBackend(CeedGetData(ceed, &cuda_data));
540*915834c9SZach Atkins     CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
541*915834c9SZach Atkins     CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
542*915834c9SZach Atkins     CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
543*915834c9SZach Atkins     CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
544*915834c9SZach Atkins     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
545*915834c9SZach Atkins 
546*915834c9SZach Atkins     // Input vectors
547*915834c9SZach Atkins     for (CeedInt i = 0; i < num_input_fields; i++) {
548*915834c9SZach Atkins       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
549*915834c9SZach Atkins       if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
550*915834c9SZach Atkins         data->fields.inputs[i] = NULL;
551*915834c9SZach Atkins       } else {
552*915834c9SZach Atkins         bool       is_active;
553*915834c9SZach Atkins         CeedVector vec;
554*915834c9SZach Atkins 
555*915834c9SZach Atkins         // Get input vector
556*915834c9SZach Atkins         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
557*915834c9SZach Atkins         is_active = vec == CEED_VECTOR_ACTIVE;
558*915834c9SZach Atkins         if (is_active) data->fields.inputs[i] = NULL;
559*915834c9SZach Atkins         else CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->fields.inputs[i]));
560*915834c9SZach Atkins         CeedCallBackend(CeedVectorDestroy(&vec));
561*915834c9SZach Atkins       }
562*915834c9SZach Atkins     }
563*915834c9SZach Atkins 
564*915834c9SZach Atkins     // Point coordinates
565*915834c9SZach Atkins     {
566*915834c9SZach Atkins       CeedVector vec;
567*915834c9SZach Atkins 
568*915834c9SZach Atkins       CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec));
569*915834c9SZach Atkins       CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->points.coords));
570*915834c9SZach Atkins       CeedCallBackend(CeedVectorDestroy(&vec));
571*915834c9SZach Atkins 
572*915834c9SZach Atkins       // Points per elem
573*915834c9SZach Atkins       if (num_elem != data->points.num_elem) {
574*915834c9SZach Atkins         CeedInt            *points_per_elem;
575*915834c9SZach Atkins         const CeedInt       num_bytes   = num_elem * sizeof(CeedInt);
576*915834c9SZach Atkins         CeedElemRestriction rstr_points = NULL;
577*915834c9SZach Atkins 
578*915834c9SZach Atkins         data->points.num_elem = num_elem;
579*915834c9SZach Atkins         CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
580*915834c9SZach Atkins         CeedCallBackend(CeedCalloc(num_elem, &points_per_elem));
581*915834c9SZach Atkins         for (CeedInt e = 0; e < num_elem; e++) {
582*915834c9SZach Atkins           CeedInt num_points_elem;
583*915834c9SZach Atkins 
584*915834c9SZach Atkins           CeedCallBackend(CeedElemRestrictionGetNumPointsInElement(rstr_points, e, &num_points_elem));
585*915834c9SZach Atkins           points_per_elem[e] = num_points_elem;
586*915834c9SZach Atkins         }
587*915834c9SZach Atkins         if (data->points.num_per_elem) CeedCallCuda(ceed, cudaFree((void **)data->points.num_per_elem));
588*915834c9SZach Atkins         CeedCallCuda(ceed, cudaMalloc((void **)&data->points.num_per_elem, num_bytes));
589*915834c9SZach Atkins         CeedCallCuda(ceed, cudaMemcpy((void *)data->points.num_per_elem, points_per_elem, num_bytes, cudaMemcpyHostToDevice));
590*915834c9SZach Atkins         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
591*915834c9SZach Atkins         CeedCallBackend(CeedFree(&points_per_elem));
592*915834c9SZach Atkins       }
593*915834c9SZach Atkins     }
594*915834c9SZach Atkins 
595*915834c9SZach Atkins     // Get context data
596*915834c9SZach Atkins     CeedCallBackend(CeedQFunctionGetInnerContextData(qf, CEED_MEM_DEVICE, &qf_data->d_c));
597*915834c9SZach Atkins 
598*915834c9SZach Atkins     // Assembly array
599*915834c9SZach Atkins     CeedCallBackend(CeedVectorGetArray(assembled, CEED_MEM_DEVICE, &assembled_array));
600*915834c9SZach Atkins     CeedScalar *assembled_offset_array = &assembled_array[offset];
601*915834c9SZach Atkins 
602*915834c9SZach Atkins     // Assemble diagonal
603*915834c9SZach Atkins     void *opargs[] = {(void *)&num_elem, &qf_data->d_c, &data->indices, &data->fields,          &data->B,
604*915834c9SZach Atkins                       &data->G,          &data->W,      &data->points,  &assembled_offset_array};
605*915834c9SZach Atkins     int   max_threads_per_block, min_grid_size, grid;
606*915834c9SZach Atkins 
607*915834c9SZach Atkins     CeedCallCuda(ceed, cuOccupancyMaxPotentialBlockSize(&min_grid_size, &max_threads_per_block, data->op, dynamicSMemSize, 0, 0x10000));
608*915834c9SZach Atkins     int block[3] = {data->thread_1d, (data->dim == 1 ? 1 : data->thread_1d), -1};
609*915834c9SZach Atkins 
610*915834c9SZach Atkins     CeedCallBackend(BlockGridCalculate(num_elem, min_grid_size / cuda_data->device_prop.multiProcessorCount, 1,
611*915834c9SZach Atkins                                        cuda_data->device_prop.maxThreadsDim[2], cuda_data->device_prop.warpSize, block, &grid));
612*915834c9SZach Atkins     CeedInt shared_mem = block[0] * block[1] * block[2] * sizeof(CeedScalar);
613*915834c9SZach Atkins 
614*915834c9SZach Atkins     CeedCallBackend(
615*915834c9SZach Atkins         CeedTryRunKernelDimShared_Cuda(ceed, data->assemble_full, NULL, grid, block[0], block[1], block[2], shared_mem, &is_run_good, opargs));
616*915834c9SZach Atkins     CeedCallCuda(ceed, cudaDeviceSynchronize());
617*915834c9SZach Atkins 
618*915834c9SZach Atkins     // Restore input arrays
619*915834c9SZach Atkins     for (CeedInt i = 0; i < num_input_fields; i++) {
620*915834c9SZach Atkins       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
621*915834c9SZach Atkins       if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
622*915834c9SZach Atkins       } else {
623*915834c9SZach Atkins         bool       is_active;
624*915834c9SZach Atkins         CeedVector vec;
625*915834c9SZach Atkins 
626*915834c9SZach Atkins         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
627*915834c9SZach Atkins         is_active = vec == CEED_VECTOR_ACTIVE;
628*915834c9SZach Atkins         if (!is_active) CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->fields.inputs[i]));
629*915834c9SZach Atkins         CeedCallBackend(CeedVectorDestroy(&vec));
630*915834c9SZach Atkins       }
631*915834c9SZach Atkins     }
632*915834c9SZach Atkins 
633*915834c9SZach Atkins     // Restore point coordinates
634*915834c9SZach Atkins     {
635*915834c9SZach Atkins       CeedVector vec;
636*915834c9SZach Atkins 
637*915834c9SZach Atkins       CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec));
638*915834c9SZach Atkins       CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->points.coords));
639*915834c9SZach Atkins       CeedCallBackend(CeedVectorDestroy(&vec));
640*915834c9SZach Atkins     }
641*915834c9SZach Atkins 
642*915834c9SZach Atkins     // Restore context data
643*915834c9SZach Atkins     CeedCallBackend(CeedQFunctionRestoreInnerContextData(qf, &qf_data->d_c));
644*915834c9SZach Atkins 
645*915834c9SZach Atkins     // Restore assembly array
646*915834c9SZach Atkins     CeedCallBackend(CeedVectorRestoreArray(assembled, &assembled_array));
647*915834c9SZach Atkins 
648*915834c9SZach Atkins     // Cleanup
649*915834c9SZach Atkins     CeedCallBackend(CeedQFunctionDestroy(&qf));
650*915834c9SZach Atkins     if (!is_run_good) data->use_assembly_fallback = true;
651*915834c9SZach Atkins   }
652*915834c9SZach Atkins   CeedCallBackend(CeedDestroy(&ceed));
653*915834c9SZach Atkins 
654*915834c9SZach Atkins   // Fallback, if needed
655*915834c9SZach Atkins   if (data->use_assembly_fallback) {
656*915834c9SZach Atkins     CeedOperator op_fallback;
657*915834c9SZach Atkins 
658*915834c9SZach Atkins     CeedDebug256(CeedOperatorReturnCeed(op), CEED_DEBUG_COLOR_SUCCESS, "Falling back to /gpu/cuda/ref CeedOperator");
659*915834c9SZach Atkins     CeedCallBackend(CeedOperatorGetFallback(op, &op_fallback));
660*915834c9SZach Atkins     CeedCallBackend(CeedSingleOperatorAssemble(op_fallback, offset, assembled));
661*915834c9SZach Atkins     return CEED_ERROR_SUCCESS;
662*915834c9SZach Atkins   }
663*915834c9SZach Atkins   return CEED_ERROR_SUCCESS;
664*915834c9SZach Atkins }
665*915834c9SZach Atkins 
666*915834c9SZach Atkins //------------------------------------------------------------------------------
667ab213215SJeremy L Thompson // Create operator
668ab213215SJeremy L Thompson //------------------------------------------------------------------------------
669241a4b83SYohann int CeedOperatorCreate_Cuda_gen(CeedOperator op) {
6700183ed61SJeremy L Thompson   bool                   is_composite, is_at_points;
671241a4b83SYohann   Ceed                   ceed;
672241a4b83SYohann   CeedOperator_Cuda_gen *impl;
673241a4b83SYohann 
674ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
6752b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl));
6762b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetData(op, impl));
677c99afcd8SJeremy L Thompson   CeedCall(CeedOperatorIsComposite(op, &is_composite));
678c99afcd8SJeremy L Thompson   if (is_composite) {
679c99afcd8SJeremy L Thompson     CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAddComposite", CeedOperatorApplyAddComposite_Cuda_gen));
680c99afcd8SJeremy L Thompson   } else {
6812b730f8bSJeremy L Thompson     CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda_gen));
682c99afcd8SJeremy L Thompson   }
6830183ed61SJeremy L Thompson   CeedCall(CeedOperatorIsAtPoints(op, &is_at_points));
6840183ed61SJeremy L Thompson   if (is_at_points) {
6850183ed61SJeremy L Thompson     CeedCallBackend(
6860183ed61SJeremy L Thompson         CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonalAtPoints_Cuda_gen));
687*915834c9SZach Atkins     CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleSingle", CeedSingleOperatorAssembleAtPoints_Cuda_gen));
6880183ed61SJeremy L Thompson   }
6892b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda_gen));
6909bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
691e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
692241a4b83SYohann }
6936aa95790SJeremy L Thompson 
694ab213215SJeremy L Thompson //------------------------------------------------------------------------------
695