xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-gen/ceed-cuda-gen-operator-build.cpp (revision edb2538e3dd6743c029967fc4e89c6fcafedb8c2)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, 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
73d576824SJeremy L Thompson 
8b8e71988SJeremy L Thompson #define CEED_DEBUG_COLOR 12
97df94212SJeremy L Thompson 
1049aac155SJeremy L Thompson #include <ceed.h>
11ec3da8bcSJed Brown #include <ceed/backend.h>
129e201c85SYohann #include <ceed/jit-tools.h>
133d576824SJeremy L Thompson #include <cuda_runtime.h>
142b730f8bSJeremy L Thompson 
15241a4b83SYohann #include <iostream>
16241a4b83SYohann #include <sstream>
17b2165e7aSSebastian Grimberg #include <string>
182b730f8bSJeremy L Thompson 
190d0321e0SJeremy L Thompson #include "../cuda-ref/ceed-cuda-ref.h"
20241a4b83SYohann #include "../cuda-shared/ceed-cuda-shared.h"
2149aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
222b730f8bSJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
232b730f8bSJeremy L Thompson #include "ceed-cuda-gen.h"
24241a4b83SYohann 
25ab213215SJeremy L Thompson //------------------------------------------------------------------------------
26b2165e7aSSebastian Grimberg // Build single operator kernel
27ab213215SJeremy L Thompson //------------------------------------------------------------------------------
28eb7e6cafSJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Cuda_gen(CeedOperator op) {
29241a4b83SYohann   using std::ostringstream;
30241a4b83SYohann   using std::string;
31ca735530SJeremy L Thompson 
32ca735530SJeremy L Thompson   bool                      is_setup_done, is_identity_qf;
33ca735530SJeremy L Thompson   struct cudaDeviceProp     prop;
34241a4b83SYohann   Ceed                      ceed;
35ca735530SJeremy L Thompson   Ceed_Cuda                *ceed_data;
36ca735530SJeremy L Thompson   CeedSize                  l_size;
372b730f8bSJeremy L Thompson   CeedInt                   Q, P_1d = 0, Q_1d = 0, elem_size, num_input_fields, num_output_fields, num_comp, dim = 1;
389e201c85SYohann   CeedEvalMode              eval_mode;
39*edb2538eSJeremy L Thompson   CeedElemRestriction       elem_rstr;
40*edb2538eSJeremy L Thompson   CeedElemRestriction_Cuda *rstr_data;
41241a4b83SYohann   CeedBasis                 basis;
42241a4b83SYohann   CeedBasis_Cuda_shared    *basis_data;
43ca735530SJeremy L Thompson   CeedQFunctionField       *qf_input_fields, *qf_output_fields;
44ca735530SJeremy L Thompson   CeedQFunction_Cuda_gen   *qf_data;
45ca735530SJeremy L Thompson   CeedQFunction             qf;
46ca735530SJeremy L Thompson   CeedOperatorField        *op_input_fields, *op_output_fields;
47ca735530SJeremy L Thompson   CeedOperator_Cuda_gen    *data;
48ca735530SJeremy L Thompson 
49ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
50ca735530SJeremy L Thompson   if (is_setup_done) return CEED_ERROR_SUCCESS;
51ca735530SJeremy L Thompson 
52ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
53ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
54ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
55ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
56ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
57ca735530SJeremy L Thompson   Q_1d = Q;
58ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
59ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
60241a4b83SYohann 
619e201c85SYohann   // TODO: put in a function?
620b454692Sjeremylt   // Check for restriction only identity operator
632b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
640b454692Sjeremylt   if (is_identity_qf) {
659e201c85SYohann     CeedEvalMode eval_mode_in, eval_mode_out;
66ca735530SJeremy L Thompson 
672b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
682b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
696574a04fSJeremy L Thompson     CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
706574a04fSJeremy L Thompson               "Backend does not implement restriction only identity operators");
710b454692Sjeremylt   }
720b454692Sjeremylt 
73241a4b83SYohann   ostringstream code;
74241a4b83SYohann 
759e201c85SYohann   // TODO: put in a function?
76f1a13f77SYohann Dudouit   // Add atomicAdd function for old NVidia architectures
772b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_data));
782b730f8bSJeremy L Thompson   CeedCallBackend(cudaGetDeviceProperties(&prop, ceed_data->device_id));
7980a9ef05SNatalie Beams   if ((prop.major < 6) && (CEED_SCALAR_TYPE != CEED_SCALAR_FP32)) {
809e201c85SYohann     char *atomic_add_path, *atomic_add_source;
81ca735530SJeremy L Thompson 
822b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-atomic-add-fallback.h", &atomic_add_path));
8323d4529eSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Atomic Add Source -----\n");
842b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, atomic_add_path, &atomic_add_source));
859e201c85SYohann     code << atomic_add_source;
862b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&atomic_add_path));
872b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&atomic_add_source));
88f1a13f77SYohann Dudouit   }
89f1a13f77SYohann Dudouit 
909e201c85SYohann   // Load basis source files
919e201c85SYohann   // TODO: generalize to accept different device functions?
929e201c85SYohann   {
939e201c85SYohann     char *tensor_basis_kernel_path, *tensor_basis_kernel_source;
94ca735530SJeremy L Thompson 
952b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-shared-basis-tensor-templates.h", &tensor_basis_kernel_path));
9623d4529eSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Tensor Basis Kernel Source -----\n");
972b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, tensor_basis_kernel_path, &tensor_basis_kernel_source));
989e201c85SYohann     code << tensor_basis_kernel_source;
992b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&tensor_basis_kernel_path));
1002b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&tensor_basis_kernel_source));
1019e201c85SYohann   }
1029e201c85SYohann   {
1039e201c85SYohann     char *cuda_gen_template_path, *cuda_gen_template_source;
104ca735530SJeremy L Thompson 
1052b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-gen-templates.h", &cuda_gen_template_path));
10623d4529eSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Cuda-Gen Template Source -----\n");
1072b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, cuda_gen_template_path, &cuda_gen_template_source));
1089e201c85SYohann     code << cuda_gen_template_source;
1092b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&cuda_gen_template_path));
1102b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&cuda_gen_template_source));
1119e201c85SYohann   }
112241a4b83SYohann 
1139e201c85SYohann   // Get QFunction source and name
1149e201c85SYohann   string q_function_source(qf_data->q_function_source);
1159e201c85SYohann   string q_function_name(qf_data->q_function_name);
1169e201c85SYohann   string operator_name;
117204bfdd7SJeremy L Thompson   operator_name = "CeedKernelCudaGenOperator_" + q_function_name;
1184d537eeaSYohann 
1199e201c85SYohann   // Find dim, P_1d, Q_1d
1209e201c85SYohann   data->max_P_1d = 0;
1219e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1222b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
123356036faSJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
124ca735530SJeremy L Thompson       bool is_tensor;
125ca735530SJeremy L Thompson 
1262b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1272b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1281da99368SJeremy L Thompson 
1299e201c85SYohann       // Collect dim, P_1d, and Q_1d
1302b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &dim));
131ca735530SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
132ca735530SJeremy L Thompson       CeedCheck(is_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
1332b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
1342b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
135ca735530SJeremy L Thompson       data->max_P_1d = CeedIntMax(data->max_P_1d, P_1d);
1361da99368SJeremy L Thompson     }
1371da99368SJeremy L Thompson   }
1389e201c85SYohann   // Check output bases for Q_1d, dim as well
139356036faSJeremy L Thompson   //   The only input basis might be CEED_BASIS_NONE
1409e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1412b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
142356036faSJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
143ca735530SJeremy L Thompson       bool is_tensor;
144ca735530SJeremy L Thompson 
1452b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1462b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1470f54b25eSjeremylt 
1489e201c85SYohann       // Collect Q_1d
1492b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &dim));
150ca735530SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
151ca735530SJeremy L Thompson       CeedCheck(is_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
1522b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
1531da99368SJeremy L Thompson     }
1541da99368SJeremy L Thompson   }
1551da99368SJeremy L Thompson   data->dim  = dim;
1569e201c85SYohann   data->Q_1d = Q_1d;
1571da99368SJeremy L Thompson 
1589e201c85SYohann   // Only use 3D collocated gradient parallelization strategy when gradient is computed
1599e201c85SYohann   // TODO: put in a function?
1609e201c85SYohann   bool use_collograd_parallelization = false;
161ca735530SJeremy L Thompson 
1629e201c85SYohann   if (dim == 3) {
1639e201c85SYohann     bool was_grad_found = false;
164ca735530SJeremy L Thompson 
1659e201c85SYohann     for (CeedInt i = 0; i < num_input_fields; i++) {
1662b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1679e201c85SYohann       if (eval_mode == CEED_EVAL_GRAD) {
1682b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1692b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1701c66c397SJeremy L Thompson         use_collograd_parallelization = basis_data->d_collo_grad_1d && (was_grad_found ? use_collograd_parallelization : true);
1719e201c85SYohann         was_grad_found                = true;
1729e201c85SYohann       }
1739e201c85SYohann     }
1749e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
1752b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1769e201c85SYohann       if (eval_mode == CEED_EVAL_GRAD) {
1772b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1782b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1791c66c397SJeremy L Thompson         use_collograd_parallelization = basis_data->d_collo_grad_1d && (was_grad_found ? use_collograd_parallelization : true);
1809e201c85SYohann         was_grad_found                = true;
1819e201c85SYohann       }
1829e201c85SYohann     }
1831da99368SJeremy L Thompson   }
1841da99368SJeremy L Thompson 
1859e201c85SYohann   // Define CEED_Q_VLA
1869e201c85SYohann   code << "\n#undef CEED_Q_VLA\n";
1879e201c85SYohann   if (dim != 3 || use_collograd_parallelization) {
1889e201c85SYohann     code << "#define CEED_Q_VLA 1\n\n";
1899e201c85SYohann   } else {
1909e201c85SYohann     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
1919e201c85SYohann   }
1929e201c85SYohann 
1939e201c85SYohann   code << q_function_source;
194241a4b83SYohann 
195241a4b83SYohann   // Setup
196818e0025SJeremy L Thompson   code << "\n// -----------------------------------------------------------------------------\n";
1972b730f8bSJeremy L Thompson   code << "\nextern \"C\" __global__ void " << operator_name
1982b730f8bSJeremy L Thompson        << "(CeedInt num_elem, void* ctx, FieldsInt_Cuda indices, Fields_Cuda fields, Fields_Cuda B, Fields_Cuda G, CeedScalar* W) {\n";
1999e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
2002b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
2019e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
2029e201c85SYohann       code << "  const CeedScalar* d_u_" << i << " = fields.inputs[" << i << "];\n";
203241a4b83SYohann     }
204241a4b83SYohann   }
205241a4b83SYohann 
2069e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
2079e201c85SYohann     code << "  CeedScalar* d_v_" << i << " = fields.outputs[" << i << "];\n";
208241a4b83SYohann   }
2091da99368SJeremy L Thompson 
2109e201c85SYohann   code << "  const CeedInt dim = " << dim << ";\n";
2119e201c85SYohann   code << "  const CeedInt Q_1d = " << Q_1d << ";\n";
2121da99368SJeremy L Thompson 
213241a4b83SYohann   code << "  extern __shared__ CeedScalar slice[];\n";
2149e201c85SYohann   // TODO put in a function? InitSharedData_Cuda?
2159e201c85SYohann   code << "  SharedData_Cuda data;\n";
2169e201c85SYohann   code << "  data.t_id_x = threadIdx.x;\n";
2179e201c85SYohann   code << "  data.t_id_y = threadIdx.y;\n";
2189e201c85SYohann   code << "  data.t_id_z = threadIdx.z;\n";
2199e201c85SYohann   code << "  data.t_id  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
2209e201c85SYohann   code << "  data.slice = slice+data.t_id_z*T_1D" << (dim > 1 ? "*T_1D" : "") << ";\n";
221920dcdc4Sjeremylt 
222818e0025SJeremy L Thompson   code << "\n  // -- Input field constants and basis data --\n";
2239e201c85SYohann   // TODO: Put in a function?
224ac421f39SYohann   // Initialize constants, and matrices B and G
2259e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
226818e0025SJeremy L Thompson     code << "  // ---- Input field " << i << " ----\n";
2279e201c85SYohann     // Get elem_size, eval_mode, num_comp
228*edb2538eSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
229*edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
2302b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
231*edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
232920dcdc4Sjeremylt 
233920dcdc4Sjeremylt     // Set field constants
2349e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {
2352b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
236356036faSJeremy L Thompson       if (basis != CEED_BASIS_NONE) {
2372b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
2389e201c85SYohann         code << "  const CeedInt P_in_" << i << " = " << P_1d << ";\n";
239920dcdc4Sjeremylt       } else {
2409e201c85SYohann         code << "  const CeedInt P_in_" << i << " = " << Q_1d << ";\n";
241920dcdc4Sjeremylt       }
2429e201c85SYohann       code << "  const CeedInt num_comp_in_" << i << " = " << num_comp << ";\n";
243920dcdc4Sjeremylt     }
244920dcdc4Sjeremylt 
245920dcdc4Sjeremylt     // Load basis data
2469e201c85SYohann     code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
2479e201c85SYohann     switch (eval_mode) {
248241a4b83SYohann       case CEED_EVAL_NONE:
249241a4b83SYohann         break;
250241a4b83SYohann       case CEED_EVAL_INTERP:
2512b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
2529e201c85SYohann         data->B.inputs[i] = basis_data->d_interp_1d;
2539e201c85SYohann         code << "  __shared__ CeedScalar s_B_in_" << i << "[" << P_1d * Q_1d << "];\n";
2549e201c85SYohann         code << "  loadMatrix<P_in_" << i << ",Q_1d>(data, B.inputs[" << i << "], s_B_in_" << i << ");\n";
255241a4b83SYohann         break;
256241a4b83SYohann       case CEED_EVAL_GRAD:
2572b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
2589e201c85SYohann         data->B.inputs[i] = basis_data->d_interp_1d;
2599e201c85SYohann         code << "  __shared__ CeedScalar s_B_in_" << i << "[" << P_1d * Q_1d << "];\n";
2609e201c85SYohann         code << "  loadMatrix<P_in_" << i << ",Q_1d>(data, B.inputs[" << i << "], s_B_in_" << i << ");\n";
2619e201c85SYohann         if (use_collograd_parallelization) {
2629e201c85SYohann           data->G.inputs[i] = basis_data->d_collo_grad_1d;
2639e201c85SYohann           code << "  __shared__ CeedScalar s_G_in_" << i << "[" << Q_1d * Q_1d << "];\n";
2649e201c85SYohann           code << "  loadMatrix<Q_1d,Q_1d>(data, G.inputs[" << i << "], s_G_in_" << i << ");\n";
265ac421f39SYohann         } else {
2661c66c397SJeremy L Thompson           bool has_collo_grad = basis_data->d_collo_grad_1d;
2679e201c85SYohann           data->G.inputs[i]   = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2689e201c85SYohann           code << "  __shared__ CeedScalar s_G_in_" << i << "[" << Q_1d * (has_collo_grad ? Q_1d : P_1d) << "];\n";
2692b730f8bSJeremy L Thompson           code << "  loadMatrix<" << (has_collo_grad ? "Q_1d" : ("P_in_" + std::to_string(i))) << ",Q_1d>(data, G.inputs[" << i << "], s_G_in_" << i
2702b730f8bSJeremy L Thompson                << ");\n";
271ac421f39SYohann         }
272241a4b83SYohann         break;
273241a4b83SYohann       case CEED_EVAL_WEIGHT:
274241a4b83SYohann         break;  // No action
275241a4b83SYohann       case CEED_EVAL_DIV:
276241a4b83SYohann         break;  // TODO: Not implemented
277241a4b83SYohann       case CEED_EVAL_CURL:
278241a4b83SYohann         break;  // TODO: Not implemented
279241a4b83SYohann     }
280241a4b83SYohann   }
281920dcdc4Sjeremylt 
282818e0025SJeremy L Thompson   code << "\n  // -- Output field constants and basis data --\n";
2839e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
284818e0025SJeremy L Thompson     code << "  // ---- Output field " << i << " ----\n";
2859e201c85SYohann     // Get elem_size, eval_mode, num_comp
286*edb2538eSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
287*edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
2882b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
289*edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
290920dcdc4Sjeremylt 
291920dcdc4Sjeremylt     // Set field constants
2922b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
293356036faSJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
2942b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
2959e201c85SYohann       code << "  const CeedInt P_out_" << i << " = " << P_1d << ";\n";
296920dcdc4Sjeremylt     } else {
2979e201c85SYohann       code << "  const CeedInt P_out_" << i << " = " << Q_1d << ";\n";
298920dcdc4Sjeremylt     }
2999e201c85SYohann     code << "  const CeedInt num_comp_out_" << i << " = " << num_comp << ";\n";
300920dcdc4Sjeremylt 
301920dcdc4Sjeremylt     // Load basis data
3029e201c85SYohann     code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
3039e201c85SYohann     switch (eval_mode) {
304920dcdc4Sjeremylt       case CEED_EVAL_NONE:
305920dcdc4Sjeremylt         break;  // No action
306920dcdc4Sjeremylt       case CEED_EVAL_INTERP:
3072b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
3089e201c85SYohann         data->B.outputs[i] = basis_data->d_interp_1d;
3099e201c85SYohann         code << "  __shared__ CeedScalar s_B_out_" << i << "[" << P_1d * Q_1d << "];\n";
3109e201c85SYohann         code << "  loadMatrix<P_out_" << i << ",Q_1d>(data, B.outputs[" << i << "], s_B_out_" << i << ");\n";
311241a4b83SYohann         break;
312241a4b83SYohann       case CEED_EVAL_GRAD:
3132b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
3149e201c85SYohann         data->B.outputs[i] = basis_data->d_interp_1d;
3159e201c85SYohann         code << "  __shared__ CeedScalar s_B_out_" << i << "[" << P_1d * Q_1d << "];\n";
3169e201c85SYohann         code << "  loadMatrix<P_out_" << i << ",Q_1d>(data, B.outputs[" << i << "], s_B_out_" << i << ");\n";
3179e201c85SYohann         if (use_collograd_parallelization) {
3189e201c85SYohann           data->G.outputs[i] = basis_data->d_collo_grad_1d;
3199e201c85SYohann           code << "  __shared__ CeedScalar s_G_out_" << i << "[" << Q_1d * Q_1d << "];\n";
3209e201c85SYohann           code << "  loadMatrix<Q_1d,Q_1d>(data, G.outputs[" << i << "], s_G_out_" << i << ");\n";
321ac421f39SYohann         } else {
3221c66c397SJeremy L Thompson           bool has_collo_grad = basis_data->d_collo_grad_1d;
3239e201c85SYohann           data->G.outputs[i]  = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
3249e201c85SYohann           code << "  __shared__ CeedScalar s_G_out_" << i << "[" << Q_1d * (has_collo_grad ? Q_1d : P_1d) << "];\n";
3252b730f8bSJeremy L Thompson           code << "  loadMatrix<" << (has_collo_grad ? "Q_1d" : ("P_out_" + std::to_string(i))) << ",Q_1d>(data, G.outputs[" << i << "], s_G_out_"
3262b730f8bSJeremy L Thompson                << i << ");\n";
327ac421f39SYohann         }
328241a4b83SYohann         break;
329e9f4dca0SJeremy L Thompson       // LCOV_EXCL_START
330241a4b83SYohann       case CEED_EVAL_WEIGHT: {
331241a4b83SYohann         Ceed ceed;
332ca735530SJeremy L Thompson 
3332b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
3342b730f8bSJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
335241a4b83SYohann         break;  // Should not occur
336241a4b83SYohann       }
337241a4b83SYohann       case CEED_EVAL_DIV:
338241a4b83SYohann         break;  // TODO: Not implemented
339241a4b83SYohann       case CEED_EVAL_CURL:
340241a4b83SYohann         break;  // TODO: Not implemented
341e9f4dca0SJeremy L Thompson                 // LCOV_EXCL_STOP
342241a4b83SYohann     }
343241a4b83SYohann   }
344818e0025SJeremy L Thompson   code << "\n  // -- Element loop --\n";
345ac421f39SYohann   code << "  __syncthreads();\n";
3469e201c85SYohann   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
347241a4b83SYohann   // Input basis apply if needed
348ac421f39SYohann   // Generate the correct eval mode code for each input
349818e0025SJeremy L Thompson   code << "    // -- Input field restrictions and basis actions --\n";
3509e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
351818e0025SJeremy L Thompson     code << "    // ---- Input field " << i << " ----\n";
3529e201c85SYohann     // Get elem_size, eval_mode, num_comp
353*edb2538eSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
354*edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
3552b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
356*edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
357920dcdc4Sjeremylt 
3589e201c85SYohann     // TODO: put in a function?
359920dcdc4Sjeremylt     // Restriction
3602b730f8bSJeremy L Thompson     if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_collograd_parallelization)) {
3619e201c85SYohann       code << "    CeedScalar r_u_" << i << "[num_comp_in_" << i << "*P_in_" << i << "];\n";
3629a2291e3SJeremy L Thompson 
3639e201c85SYohann       bool is_strided;
364ca735530SJeremy L Thompson 
365*edb2538eSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
3669e201c85SYohann       if (!is_strided) {
3679e201c85SYohann         CeedInt comp_stride;
368ca735530SJeremy L Thompson 
369*edb2538eSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
370ca735530SJeremy L Thompson         code << "    const CeedInt l_size_in_" << i << " = " << l_size << ";\n";
371*edb2538eSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3729e201c85SYohann         code << "    // CompStride: " << comp_stride << "\n";
373*edb2538eSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
374*edb2538eSJeremy L Thompson         data->indices.inputs[i] = rstr_data->d_ind;
375ca735530SJeremy L Thompson         code << "    readDofsOffset" << dim << "d<num_comp_in_" << i << ", " << comp_stride << ", P_in_" << i << ">(data, l_size_in_" << i
3762b730f8bSJeremy L Thompson              << ", elem, indices.inputs[" << i << "], d_u_" << i << ", r_u_" << i << ");\n";
377ccedf6b0Sjeremylt       } else {
378b2165e7aSSebastian Grimberg         bool    has_backend_strides;
3799e201c85SYohann         CeedInt num_elem;
380ca735530SJeremy L Thompson 
381*edb2538eSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
382*edb2538eSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
3839e201c85SYohann         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
384ca735530SJeremy L Thompson 
385b2165e7aSSebastian Grimberg         if (!has_backend_strides) {
386*edb2538eSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, &strides));
387ccedf6b0Sjeremylt         }
388920dcdc4Sjeremylt         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
3892b730f8bSJeremy L Thompson         code << "    readDofsStrided" << dim << "d<num_comp_in_" << i << ",P_in_" << i << "," << strides[0] << "," << strides[1] << "," << strides[2]
3902b730f8bSJeremy L Thompson              << ">(data, elem, d_u_" << i << ", r_u_" << i << ");\n";
391920dcdc4Sjeremylt       }
392920dcdc4Sjeremylt     }
393920dcdc4Sjeremylt 
3949e201c85SYohann     // TODO: put in a function?
395920dcdc4Sjeremylt     // Basis action
3969e201c85SYohann     code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
3979e201c85SYohann     switch (eval_mode) {
398920dcdc4Sjeremylt       case CEED_EVAL_NONE:
3999e201c85SYohann         if (!use_collograd_parallelization) {
4009e201c85SYohann           code << "    CeedScalar* r_t_" << i << " = r_u_" << i << ";\n";
401920dcdc4Sjeremylt         }
402920dcdc4Sjeremylt         break;
403920dcdc4Sjeremylt       case CEED_EVAL_INTERP:
4049e201c85SYohann         code << "    CeedScalar r_t_" << i << "[num_comp_in_" << i << "*Q_1d];\n";
4052b730f8bSJeremy L Thompson         code << "    Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_in_" << i << ",P_in_" << i << ",Q_1d>(data, r_u_" << i << ", s_B_in_"
4062b730f8bSJeremy L Thompson              << i << ", r_t_" << i << ");\n";
407241a4b83SYohann         break;
408241a4b83SYohann       case CEED_EVAL_GRAD:
4099e201c85SYohann         if (use_collograd_parallelization) {
4109e201c85SYohann           code << "    CeedScalar r_t_" << i << "[num_comp_in_" << i << "*Q_1d];\n";
4112b730f8bSJeremy L Thompson           code << "    Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_in_" << i << ",P_in_" << i << ",Q_1d>(data, r_u_" << i
4122b730f8bSJeremy L Thompson                << ", s_B_in_" << i << ", r_t_" << i << ");\n";
413ac421f39SYohann         } else {
4149e201c85SYohann           CeedInt P_1d;
415ca735530SJeremy L Thompson 
4162b730f8bSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
4172b730f8bSJeremy L Thompson           CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
4189e201c85SYohann           code << "    CeedScalar r_t_" << i << "[num_comp_in_" << i << "*dim*Q_1d];\n";
4192b730f8bSJeremy L Thompson           code << "    Grad" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp_in_" << i
4202b730f8bSJeremy L Thompson                << ",P_in_" << i << ",Q_1d>(data, r_u_" << i << ", s_B_in_" << i << ", s_G_in_" << i << ", r_t_" << i << ");\n";
421ac421f39SYohann         }
422241a4b83SYohann         break;
423241a4b83SYohann       case CEED_EVAL_WEIGHT:
4249e201c85SYohann         code << "    CeedScalar r_t_" << i << "[Q_1d];\n";
4252b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
4262b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
427437930d1SJeremy L Thompson         data->W = basis_data->d_q_weight_1d;
4289e201c85SYohann         code << "    Weight" << (dim > 1 ? "Tensor" : "") << dim << "d<Q_1d>(data, W, r_t_" << i << ");\n";
429241a4b83SYohann         break;  // No action
430241a4b83SYohann       case CEED_EVAL_DIV:
431241a4b83SYohann         break;  // TODO: Not implemented
432241a4b83SYohann       case CEED_EVAL_CURL:
433241a4b83SYohann         break;  // TODO: Not implemented
434241a4b83SYohann     }
435241a4b83SYohann   }
436ac421f39SYohann 
437ea61e9acSJeremy L Thompson   // TODO: put in a function + separate collograd logic
438241a4b83SYohann   // Q function
439818e0025SJeremy L Thompson   code << "\n    // -- Output field setup --\n";
4409e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
441818e0025SJeremy L Thompson     code << "\n    // ---- Output field " << i << " ----\n";
4422b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
4432b730f8bSJeremy L Thompson     if (eval_mode == CEED_EVAL_GRAD) {
4449e201c85SYohann       if (use_collograd_parallelization) {
445ac421f39SYohann         // Accumulator for gradient slices
4469e201c85SYohann         code << "    CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*Q_1d];\n";
4479e201c85SYohann         code << "    for (CeedInt i = 0; i < num_comp_out_" << i << "; i++) {\n";
4489e201c85SYohann         code << "      for (CeedInt j = 0; j < Q_1d; ++j) {\n";
4499e201c85SYohann         code << "        r_tt_" << i << "[j + i*Q_1d] = 0.0;\n";
450ac421f39SYohann         code << "      }\n";
451ac421f39SYohann         code << "    }\n";
452ac421f39SYohann       } else {
4539e201c85SYohann         code << "    CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*dim*Q_1d];\n";
454241a4b83SYohann       }
455ac421f39SYohann     }
4562b730f8bSJeremy L Thompson     if (eval_mode == CEED_EVAL_NONE || eval_mode == CEED_EVAL_INTERP) {
4579e201c85SYohann       code << "    CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*Q_1d];\n";
458241a4b83SYohann     }
459241a4b83SYohann   }
460ac421f39SYohann   // We treat quadrature points per slice in 3d to save registers
4619e201c85SYohann   if (use_collograd_parallelization) {
4629e201c85SYohann     code << "\n    // Note: Using planes of 3D elements\n";
463ac421f39SYohann     code << "#pragma unroll\n";
4649e201c85SYohann     code << "    for (CeedInt q = 0; q < Q_1d; q++) {\n";
465818e0025SJeremy L Thompson     code << "      // -- Input fields --\n";
4669e201c85SYohann     for (CeedInt i = 0; i < num_input_fields; i++) {
467818e0025SJeremy L Thompson       code << "      // ---- Input field " << i << " ----\n";
4689e201c85SYohann       // Get elem_size, eval_mode, num_comp
4692b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
470ac421f39SYohann       // Basis action
4719e201c85SYohann       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
4729e201c85SYohann       switch (eval_mode) {
473ac421f39SYohann         case CEED_EVAL_NONE:
474ca735530SJeremy L Thompson           bool is_strided;
475ca735530SJeremy L Thompson 
4769e201c85SYohann           code << "      CeedScalar r_q_" << i << "[num_comp_in_" << i << "];\n";
4779a2291e3SJeremy L Thompson 
478*edb2538eSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
479*edb2538eSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
4809e201c85SYohann           if (!is_strided) {
4819e201c85SYohann             CeedInt comp_stride;
482ca735530SJeremy L Thompson 
483*edb2538eSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
484ca735530SJeremy L Thompson             code << "      const CeedInt l_size_in_" << i << " = " << l_size << ";\n";
485*edb2538eSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
4869e201c85SYohann             code << "      // CompStride: " << comp_stride << "\n";
487*edb2538eSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
488*edb2538eSJeremy L Thompson             data->indices.inputs[i] = rstr_data->d_ind;
4892b730f8bSJeremy L Thompson             code << "      readSliceQuadsOffset"
490ca735530SJeremy L Thompson                  << "3d<num_comp_in_" << i << ", " << comp_stride << ", Q_1d>(data, l_size_in_" << i << ", elem, q, indices.inputs[" << i << "], d_u_"
4912b730f8bSJeremy L Thompson                  << i << ", r_q_" << i << ");\n";
492920dcdc4Sjeremylt           } else {
493b2165e7aSSebastian Grimberg             bool    has_backend_strides;
4949e201c85SYohann             CeedInt num_elem;
495ca735530SJeremy L Thompson 
496*edb2538eSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
497*edb2538eSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
498*edb2538eSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
4999e201c85SYohann             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
500ca735530SJeremy L Thompson 
501b2165e7aSSebastian Grimberg             if (!has_backend_strides) {
502*edb2538eSJeremy L Thompson               CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, &strides));
503920dcdc4Sjeremylt             }
504920dcdc4Sjeremylt             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
5052b730f8bSJeremy L Thompson             code << "      readSliceQuadsStrided"
5062b730f8bSJeremy L Thompson                  << "3d<num_comp_in_" << i
5072b730f8bSJeremy L Thompson                  << ",Q_1d"
5082b730f8bSJeremy L Thompson                     ","
5092b730f8bSJeremy L Thompson                  << strides[0] << "," << strides[1] << "," << strides[2] << ">(data, elem, q, d_u_" << i << ", r_q_" << i << ");\n";
510920dcdc4Sjeremylt           }
511ac421f39SYohann           break;
512ac421f39SYohann         case CEED_EVAL_INTERP:
5139e201c85SYohann           code << "      CeedScalar r_q_" << i << "[num_comp_in_" << i << "];\n";
5149e201c85SYohann           code << "      for (CeedInt j = 0; j < num_comp_in_" << i << " ; ++j) {\n";
5159e201c85SYohann           code << "        r_q_" << i << "[j] = r_t_" << i << "[q + j*Q_1d];\n";
516ac421f39SYohann           code << "      }\n";
517ac421f39SYohann           break;
518ac421f39SYohann         case CEED_EVAL_GRAD:
5199e201c85SYohann           code << "      CeedScalar r_q_" << i << "[num_comp_in_" << i << "*dim];\n";
5209e201c85SYohann           code << "      gradCollo3d<num_comp_in_" << i << ",Q_1d>(data, q, r_t_" << i << ", s_G_in_" << i << ", r_q_" << i << ");\n";
521ac421f39SYohann           break;
522ac421f39SYohann         case CEED_EVAL_WEIGHT:
5239e201c85SYohann           code << "      CeedScalar r_q_" << i << "[1];\n";
5249e201c85SYohann           code << "      r_q_" << i << "[0] = r_t_" << i << "[q];\n";
525ac421f39SYohann           break;  // No action
526ac421f39SYohann         case CEED_EVAL_DIV:
527ac421f39SYohann           break;  // TODO: Not implemented
528ac421f39SYohann         case CEED_EVAL_CURL:
529ac421f39SYohann           break;  // TODO: Not implemented
530ac421f39SYohann       }
531ac421f39SYohann     }
532818e0025SJeremy L Thompson     code << "\n      // -- Output fields --\n";
5339e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
534818e0025SJeremy L Thompson       code << "      // ---- Output field " << i << " ----\n";
5352b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
536ac421f39SYohann       // Basis action
5379e201c85SYohann       switch (eval_mode) {
538ac421f39SYohann         case CEED_EVAL_NONE:
5399e201c85SYohann           code << "      CeedScalar r_qq_" << i << "[num_comp_out_" << i << "];\n";
540ac421f39SYohann           break;  // No action
541ac421f39SYohann         case CEED_EVAL_INTERP:
5429e201c85SYohann           code << "      CeedScalar r_qq_" << i << "[num_comp_out_" << i << "];\n";
543ac421f39SYohann           break;
544ac421f39SYohann         case CEED_EVAL_GRAD:
5459e201c85SYohann           code << "      CeedScalar r_qq_" << i << "[num_comp_out_" << i << "*dim];\n";
546ac421f39SYohann           break;
547ac421f39SYohann         case CEED_EVAL_WEIGHT:
548ac421f39SYohann           break;  // Should not occur
549ac421f39SYohann         case CEED_EVAL_DIV:
550ac421f39SYohann           break;  // TODO: Not implemented
551ac421f39SYohann         case CEED_EVAL_CURL:
552ac421f39SYohann           break;  // TODO: Not implemented
553ac421f39SYohann       }
554ac421f39SYohann     }
555ac421f39SYohann   } else {
5569e201c85SYohann     code << "\n      // Note: Using full elements\n";
557818e0025SJeremy L Thompson     code << "      // -- Input fields --\n";
5589e201c85SYohann     for (CeedInt i = 0; i < num_input_fields; i++) {
559818e0025SJeremy L Thompson       code << "      // ---- Input field " << i << " ----\n";
5609e201c85SYohann       code << "      CeedScalar* r_q_" << i << " = r_t_" << i << ";\n";
561ac421f39SYohann     }
562818e0025SJeremy L Thompson     code << "      // -- Output fields --\n";
5639e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
564818e0025SJeremy L Thompson       code << "      // ---- Output field " << i << " ----\n";
5659e201c85SYohann       code << "      CeedScalar* r_qq_" << i << " = r_tt_" << i << ";\n";
566ac421f39SYohann     }
567ac421f39SYohann   }
568818e0025SJeremy L Thompson   code << "\n      // -- QFunction Inputs and outputs --\n";
5699e201c85SYohann   code << "      CeedScalar* in[" << num_input_fields << "];\n";
5709e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
571818e0025SJeremy L Thompson     code << "      // ---- Input field " << i << " ----\n";
5729e201c85SYohann     code << "      in[" << i << "] = r_q_" << i << ";\n";
5734d537eeaSYohann   }
5749e201c85SYohann   code << "      CeedScalar* out[" << num_output_fields << "];\n";
5759e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
576818e0025SJeremy L Thompson     code << "      // ---- Output field " << i << " ----\n";
5779e201c85SYohann     code << "      out[" << i << "] = r_qq_" << i << ";\n";
5784d537eeaSYohann   }
579818e0025SJeremy L Thompson   code << "\n      // -- Apply QFunction --\n";
5809e201c85SYohann   code << "      " << q_function_name << "(ctx, ";
5819e201c85SYohann   if (dim != 3 || use_collograd_parallelization) {
582ac421f39SYohann     code << "1";
583ac421f39SYohann   } else {
5849e201c85SYohann     code << "Q_1d";
585ac421f39SYohann   }
586ac421f39SYohann   code << ", in, out);\n";
5879e201c85SYohann   if (use_collograd_parallelization) {
588818e0025SJeremy L Thompson     code << "      // -- Output fields --\n";
5899e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
590818e0025SJeremy L Thompson       code << "      // ---- Output field " << i << " ----\n";
5912b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
592ac421f39SYohann       // Basis action
5939e201c85SYohann       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
5949e201c85SYohann       switch (eval_mode) {
595ac421f39SYohann         case CEED_EVAL_NONE:
5969e201c85SYohann           code << "      for (CeedInt j = 0; j < num_comp_out_" << i << " ; ++j) {\n";
5979e201c85SYohann           code << "        r_tt_" << i << "[q + j*Q_1d] = r_qq_" << i << "[j];\n";
598ac421f39SYohann           code << "      }\n";
599ac421f39SYohann           break;  // No action
600ac421f39SYohann         case CEED_EVAL_INTERP:
6019e201c85SYohann           code << "      for (CeedInt j = 0; j < num_comp_out_" << i << " ; ++j) {\n";
6029e201c85SYohann           code << "        r_tt_" << i << "[q + j*Q_1d] = r_qq_" << i << "[j];\n";
603ac421f39SYohann           code << "      }\n";
604ac421f39SYohann           break;
605ac421f39SYohann         case CEED_EVAL_GRAD:
6069e201c85SYohann           code << "      gradColloTranspose3d<num_comp_out_" << i << ",Q_1d>(data, q, r_qq_" << i << ", s_G_out_" << i << ", r_tt_" << i << ");\n";
607ac421f39SYohann           break;
608ac421f39SYohann         case CEED_EVAL_WEIGHT:
609ac421f39SYohann           break;  // Should not occur
610ac421f39SYohann         case CEED_EVAL_DIV:
611ac421f39SYohann           break;  // TODO: Not implemented
612ac421f39SYohann         case CEED_EVAL_CURL:
613ac421f39SYohann           break;  // TODO: Not implemented
614ac421f39SYohann       }
615ac421f39SYohann     }
616ac421f39SYohann     code << "    }\n";
617ac421f39SYohann   }
618241a4b83SYohann 
619241a4b83SYohann   // Output basis apply if needed
620ac421f39SYohann   // Generate the correct eval mode code for each output
621818e0025SJeremy L Thompson   code << "\n    // -- Output field basis action and restrictions --\n";
6229e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
623818e0025SJeremy L Thompson     code << "    // ---- Output field " << i << " ----\n";
6249e201c85SYohann     // Get elem_size, eval_mode, num_comp
625*edb2538eSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
626*edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
6272b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
628*edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
6299e201c85SYohann     // TODO put in a function
630241a4b83SYohann     // Basis action
6319e201c85SYohann     code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
6329e201c85SYohann     switch (eval_mode) {
633241a4b83SYohann       case CEED_EVAL_NONE:
6349e201c85SYohann         code << "    CeedScalar* r_v_" << i << " = r_tt_" << i << ";\n";
635241a4b83SYohann         break;  // No action
636241a4b83SYohann       case CEED_EVAL_INTERP:
6379e201c85SYohann         code << "    CeedScalar r_v_" << i << "[num_comp_out_" << i << "*P_out_" << i << "];\n";
6382b730f8bSJeremy L Thompson         code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_out_" << i << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i
6392b730f8bSJeremy L Thompson              << ", s_B_out_" << i << ", r_v_" << i << ");\n";
640241a4b83SYohann         break;
641241a4b83SYohann       case CEED_EVAL_GRAD:
6429e201c85SYohann         code << "    CeedScalar r_v_" << i << "[num_comp_out_" << i << "*P_out_" << i << "];\n";
6439e201c85SYohann         if (use_collograd_parallelization) {
6442b730f8bSJeremy L Thompson           code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_out_" << i << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i
6452b730f8bSJeremy L Thompson                << ", s_B_out_" << i << ", r_v_" << i << ");\n";
646ac421f39SYohann         } else {
6479e201c85SYohann           CeedInt P_1d;
6482b730f8bSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
6492b730f8bSJeremy L Thompson           CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
6502b730f8bSJeremy L Thompson           code << "    GradTranspose" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp_out_" << i
6512b730f8bSJeremy L Thompson                << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i << ", s_B_out_" << i << ", s_G_out_" << i << ", r_v_" << i << ");\n";
652ac421f39SYohann         }
653241a4b83SYohann         break;
654e9f4dca0SJeremy L Thompson       // LCOV_EXCL_START
655241a4b83SYohann       case CEED_EVAL_WEIGHT: {
656241a4b83SYohann         Ceed ceed;
6572b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
6582b730f8bSJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
659241a4b83SYohann         break;  // Should not occur
660241a4b83SYohann       }
661241a4b83SYohann       case CEED_EVAL_DIV:
662241a4b83SYohann         break;  // TODO: Not implemented
663241a4b83SYohann       case CEED_EVAL_CURL:
664241a4b83SYohann         break;  // TODO: Not implemented
665e9f4dca0SJeremy L Thompson                 // LCOV_EXCL_STOP
666241a4b83SYohann     }
6679e201c85SYohann     // TODO put in a function
668920dcdc4Sjeremylt     // Restriction
6699e201c85SYohann     bool is_strided;
670*edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
6719e201c85SYohann     if (!is_strided) {
6729e201c85SYohann       CeedInt comp_stride;
673ca735530SJeremy L Thompson 
674*edb2538eSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
675ca735530SJeremy L Thompson       code << "    const CeedInt l_size_out_" << i << " = " << l_size << ";\n";
676*edb2538eSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
6779e201c85SYohann       code << "    // CompStride: " << comp_stride << "\n";
678*edb2538eSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
679*edb2538eSJeremy L Thompson       data->indices.outputs[i] = rstr_data->d_ind;
680ca735530SJeremy L Thompson       code << "    writeDofsOffset" << dim << "d<num_comp_out_" << i << ", " << comp_stride << ", P_out_" << i << ">(data, l_size_out_" << i
6812b730f8bSJeremy L Thompson            << ", elem, indices.outputs[" << i << "], r_v_" << i << ", d_v_" << i << ");\n";
682920dcdc4Sjeremylt     } else {
6839e201c85SYohann       bool    has_backend_strides;
6849e201c85SYohann       CeedInt num_elem;
685ca735530SJeremy L Thompson 
686*edb2538eSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
687*edb2538eSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
6889e201c85SYohann       CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
689ca735530SJeremy L Thompson 
6909e201c85SYohann       if (!has_backend_strides) {
691*edb2538eSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, &strides));
692920dcdc4Sjeremylt       }
693920dcdc4Sjeremylt       code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
6942b730f8bSJeremy L Thompson       code << "    writeDofsStrided" << dim << "d<num_comp_out_" << i << ",P_out_" << i << "," << strides[0] << "," << strides[1] << "," << strides[2]
6952b730f8bSJeremy L Thompson            << ">(data, elem, r_v_" << i << ", d_v_" << i << ");\n";
696920dcdc4Sjeremylt     }
697241a4b83SYohann   }
698241a4b83SYohann 
699241a4b83SYohann   code << "  }\n";
700818e0025SJeremy L Thompson   code << "}\n";
701818e0025SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n\n";
702241a4b83SYohann 
703ab213215SJeremy L Thompson   // View kernel for debugging
70423d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "Generated Operator Kernels:\n");
7053f21f6b1SJeremy L Thompson   CeedDebug(ceed, code.str().c_str());
706241a4b83SYohann 
707eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Cuda(ceed, code.str().c_str(), &data->module, 1, "T_1D", CeedIntMax(Q_1d, data->max_P_1d)));
708eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, operator_name.c_str(), &data->op));
709241a4b83SYohann 
7102b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
711e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
712241a4b83SYohann }
7132a86cc9dSSebastian Grimberg 
714ab213215SJeremy L Thompson //------------------------------------------------------------------------------
715