xref: /libCEED/rust/libceed-sys/c-src/backends/hip-gen/ceed-hip-gen-operator-build.cpp (revision 356036fa84f714fa73ef64c9a80ce2028dde816f)
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.
37d8d0e25Snbeams //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
57d8d0e25Snbeams //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
73d576824SJeremy L Thompson 
87d8d0e25Snbeams #define CEED_DEBUG_COLOR 12
97d8d0e25Snbeams 
1049aac155SJeremy L Thompson #include <ceed.h>
11ec3da8bcSJed Brown #include <ceed/backend.h>
129e201c85SYohann #include <ceed/jit-tools.h>
132b730f8bSJeremy L Thompson 
147d8d0e25Snbeams #include <iostream>
157d8d0e25Snbeams #include <sstream>
162b730f8bSJeremy L Thompson #include <string>
172b730f8bSJeremy L Thompson 
180d0321e0SJeremy L Thompson #include "../hip-ref/ceed-hip-ref.h"
197d8d0e25Snbeams #include "../hip-shared/ceed-hip-shared.h"
20b2165e7aSSebastian Grimberg #include "../hip/ceed-hip-common.h"
217d8d0e25Snbeams #include "../hip/ceed-hip-compile.h"
222b730f8bSJeremy L Thompson #include "ceed-hip-gen.h"
23b3e1519bSnbeams 
24b3e1519bSnbeams //------------------------------------------------------------------------------
25b3e1519bSnbeams // Calculate the block size used for launching the operator kernel
26b3e1519bSnbeams //------------------------------------------------------------------------------
272b730f8bSJeremy L Thompson extern "C" int BlockGridCalculate_Hip_gen(const CeedInt dim, const CeedInt num_elem, const CeedInt P_1d, const CeedInt Q_1d, CeedInt *block_sizes) {
289e201c85SYohann   const CeedInt thread1d = CeedIntMax(Q_1d, P_1d);
29b3e1519bSnbeams   if (dim == 1) {
309e201c85SYohann     CeedInt elems_per_block = 64 * thread1d > 256 ? 256 / thread1d : 64;
319e201c85SYohann     elems_per_block         = elems_per_block > 0 ? elems_per_block : 1;
32b3e1519bSnbeams     block_sizes[0]          = thread1d;
33b3e1519bSnbeams     block_sizes[1]          = 1;
349e201c85SYohann     block_sizes[2]          = elems_per_block;
35b3e1519bSnbeams   } else if (dim == 2) {
369e201c85SYohann     const CeedInt elems_per_block = thread1d < 4 ? 16 : 2;
37b3e1519bSnbeams     block_sizes[0]                = thread1d;
38b3e1519bSnbeams     block_sizes[1]                = thread1d;
399e201c85SYohann     block_sizes[2]                = elems_per_block;
40b3e1519bSnbeams   } else if (dim == 3) {
419e201c85SYohann     const CeedInt elems_per_block = thread1d < 6 ? 4 : (thread1d < 8 ? 2 : 1);
42b3e1519bSnbeams     block_sizes[0]                = thread1d;
43b3e1519bSnbeams     block_sizes[1]                = thread1d;
449e201c85SYohann     block_sizes[2]                = elems_per_block;
45b3e1519bSnbeams   }
46b3e1519bSnbeams   return CEED_ERROR_SUCCESS;
47b3e1519bSnbeams }
48b3e1519bSnbeams 
497d8d0e25Snbeams //------------------------------------------------------------------------------
509e201c85SYohann // Build single operator kernel
517d8d0e25Snbeams //------------------------------------------------------------------------------
52eb7e6cafSJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Hip_gen(CeedOperator op) {
537d8d0e25Snbeams   using std::ostringstream;
547d8d0e25Snbeams   using std::string;
559e201c85SYohann   bool is_setup_done;
562b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
579e201c85SYohann   if (is_setup_done) return CEED_ERROR_SUCCESS;
587d8d0e25Snbeams   Ceed ceed;
592b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
607d8d0e25Snbeams   CeedOperator_Hip_gen *data;
612b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
627d8d0e25Snbeams   CeedQFunction          qf;
637d8d0e25Snbeams   CeedQFunction_Hip_gen *qf_data;
642b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
652b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
66e79b91d9SJeremy L Thompson   CeedSize lsize;
672b730f8bSJeremy L Thompson   CeedInt  Q, P_1d = 0, Q_1d = 0, elem_size, num_input_fields, num_output_fields, num_comp, dim = 1;
682b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
699e201c85SYohann   Q_1d = Q;
709e201c85SYohann   CeedOperatorField *op_input_fields, *op_output_fields;
712b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
729e201c85SYohann   CeedQFunctionField *qf_input_fields, *qf_output_fields;
732b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
749e201c85SYohann   CeedEvalMode             eval_mode;
757d8d0e25Snbeams   CeedBasis                basis;
767d8d0e25Snbeams   CeedBasis_Hip_shared    *basis_data;
777d8d0e25Snbeams   CeedElemRestriction      Erestrict;
787d8d0e25Snbeams   CeedElemRestriction_Hip *restr_data;
797d8d0e25Snbeams 
80b2165e7aSSebastian Grimberg   // TODO: put in a function?
810b454692Sjeremylt   // Check for restriction only identity operator
820b454692Sjeremylt   bool is_identity_qf;
832b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
840b454692Sjeremylt   if (is_identity_qf) {
859e201c85SYohann     CeedEvalMode eval_mode_in, eval_mode_out;
862b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
872b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
886574a04fSJeremy L Thompson     CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
896574a04fSJeremy L Thompson               "Backend does not implement restriction only identity operators");
900b454692Sjeremylt   }
910b454692Sjeremylt 
927d8d0e25Snbeams   ostringstream code;
93b2165e7aSSebastian Grimberg 
94b2165e7aSSebastian Grimberg   // Load basis source files
959e201c85SYohann   // TODO: generalize to accept different device functions?
969e201c85SYohann   {
979e201c85SYohann     char *tensor_basis_kernel_path, *tensor_basis_kernel_source;
982b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-shared-basis-tensor-templates.h", &tensor_basis_kernel_path));
9923d4529eSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Tensor Basis Kernel Source -----\n");
1002b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, tensor_basis_kernel_path, &tensor_basis_kernel_source));
1019e201c85SYohann     code << tensor_basis_kernel_source;
1022b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&tensor_basis_kernel_path));
1032b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&tensor_basis_kernel_source));
1049e201c85SYohann   }
1059e201c85SYohann   {
1069e201c85SYohann     char *hip_gen_template_path, *hip_gen_template_source;
1072b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-gen-templates.h", &hip_gen_template_path));
10823d4529eSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Hip-Gen Template Source -----\n");
1092b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, hip_gen_template_path, &hip_gen_template_source));
1109e201c85SYohann     code << hip_gen_template_source;
1112b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&hip_gen_template_path));
1122b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&hip_gen_template_source));
1139e201c85SYohann   }
1147d8d0e25Snbeams 
115b2165e7aSSebastian Grimberg   // Get QFunction source and name
1169e201c85SYohann   string q_function_source(qf_data->q_function_source);
1179e201c85SYohann   string q_function_name(qf_data->q_function_name);
1189e201c85SYohann   string operator_name;
119204bfdd7SJeremy L Thompson   operator_name = "CeedKernelHipGenOperator_" + q_function_name;
1207d8d0e25Snbeams 
1219e201c85SYohann   // Find dim, P_1d, Q_1d
1229e201c85SYohann   data->max_P_1d = 0;
1239e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1242b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
125*356036faSJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
1262b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1272b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1287d8d0e25Snbeams 
1299e201c85SYohann       // Collect dim, P_1d, and Q_1d
1302b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &dim));
1317d8d0e25Snbeams       bool isTensor;
1322b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &isTensor));
1336574a04fSJeremy L Thompson       CeedCheck(isTensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
1342b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
1352b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
1369e201c85SYohann       if (P_1d > data->max_P_1d) data->max_P_1d = P_1d;
1377d8d0e25Snbeams     }
1387d8d0e25Snbeams   }
1399e201c85SYohann   // Check output bases for Q_1d, dim as well
140*356036faSJeremy L Thompson   //   The only input basis might be CEED_BASIS_NONE
1419e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1422b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1437d8d0e25Snbeams 
144*356036faSJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
1452b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1462b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1477d8d0e25Snbeams 
1489e201c85SYohann       // Collect Q_1d
1492b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &dim));
1507d8d0e25Snbeams       bool isTensor;
1512b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &isTensor));
1526574a04fSJeremy L Thompson       CeedCheck(isTensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
1532b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
1547d8d0e25Snbeams     }
1557d8d0e25Snbeams   }
1567d8d0e25Snbeams   data->dim  = dim;
1579e201c85SYohann   data->Q_1d = Q_1d;
1587d8d0e25Snbeams 
1599e201c85SYohann   // Only use 3D collocated gradient parallelization strategy when gradient is computed
1609e201c85SYohann   // TODO: put in a function?
1619e201c85SYohann   bool use_collograd_parallelization = false;
1629e201c85SYohann   if (dim == 3) {
1639e201c85SYohann     bool was_grad_found = false;
1649e201c85SYohann     for (CeedInt i = 0; i < num_input_fields; i++) {
1652b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1669e201c85SYohann       if (eval_mode == CEED_EVAL_GRAD) {
1672b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1682b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1691c66c397SJeremy L Thompson         use_collograd_parallelization = basis_data->d_collo_grad_1d && (was_grad_found ? use_collograd_parallelization : true);
1709e201c85SYohann         was_grad_found                = true;
1719e201c85SYohann       }
1729e201c85SYohann     }
1739e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
1742b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1759e201c85SYohann       if (eval_mode == CEED_EVAL_GRAD) {
1762b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1772b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1781c66c397SJeremy L Thompson         use_collograd_parallelization = basis_data->d_collo_grad_1d && (was_grad_found ? use_collograd_parallelization : true);
1799e201c85SYohann         was_grad_found                = true;
1809e201c85SYohann       }
1819e201c85SYohann     }
1827d8d0e25Snbeams   }
1837d8d0e25Snbeams 
1849e201c85SYohann   // Define CEED_Q_VLA
1859e201c85SYohann   code << "\n#undef CEED_Q_VLA\n";
1869e201c85SYohann   if (dim != 3 || use_collograd_parallelization) {
1879e201c85SYohann     code << "#define CEED_Q_VLA 1\n\n";
1889e201c85SYohann   } else {
1899e201c85SYohann     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
1909e201c85SYohann   }
1919e201c85SYohann 
1929e201c85SYohann   code << q_function_source;
1937d8d0e25Snbeams 
1947d8d0e25Snbeams   // Setup
1957d8d0e25Snbeams   code << "\n// -----------------------------------------------------------------------------\n";
196b3e1519bSnbeams   code << "\nextern \"C\" __launch_bounds__(BLOCK_SIZE)\n";
1972b730f8bSJeremy L Thompson   code << "__global__ void " << operator_name
1982b730f8bSJeremy L Thompson        << "(CeedInt num_elem, void* ctx, FieldsInt_Hip indices, Fields_Hip fields, Fields_Hip B, Fields_Hip 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";
2037d8d0e25Snbeams     }
2047d8d0e25Snbeams   }
2057d8d0e25Snbeams 
2069e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
2079e201c85SYohann     code << "  CeedScalar* d_v_" << i << " = fields.outputs[" << i << "];\n";
2087d8d0e25Snbeams   }
2097d8d0e25Snbeams 
2109e201c85SYohann   code << "  const CeedInt dim = " << dim << ";\n";
2119e201c85SYohann   code << "  const CeedInt Q_1d = " << Q_1d << ";\n";
2127d8d0e25Snbeams 
2137d8d0e25Snbeams   code << "  HIP_DYNAMIC_SHARED( CeedScalar, slice)\n";
214b2165e7aSSebastian Grimberg   // TODO put in a function? InitSharedData_Hip?
2159e201c85SYohann   code << "  SharedData_Hip 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";
2217d8d0e25Snbeams 
2227d8d0e25Snbeams   code << "\n  // -- Input field constants and basis data --\n";
223b2165e7aSSebastian Grimberg   // TODO: Put in a function?
2247d8d0e25Snbeams   // Initialize constants, and matrices B and G
2259e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
2267d8d0e25Snbeams     code << "  // ---- Input field " << i << " ----\n";
2279e201c85SYohann     // Get elem_size, eval_mode, num_comp
2282b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &Erestrict));
2292b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size));
2302b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
2312b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp));
2327d8d0e25Snbeams 
2337d8d0e25Snbeams     // Set field constants
2349e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {
2352b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
236*356036faSJeremy 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";
2397d8d0e25Snbeams       } else {
2409e201c85SYohann         code << "  const CeedInt P_in_" << i << " = " << Q_1d << ";\n";
2417d8d0e25Snbeams       }
2429e201c85SYohann       code << "  const CeedInt num_comp_in_" << i << " = " << num_comp << ";\n";
2437d8d0e25Snbeams     }
2447d8d0e25Snbeams 
2457d8d0e25Snbeams     // Load basis data
2469e201c85SYohann     code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
2479e201c85SYohann     switch (eval_mode) {
2487d8d0e25Snbeams       case CEED_EVAL_NONE:
2497d8d0e25Snbeams         break;
2507d8d0e25Snbeams       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";
2557d8d0e25Snbeams         break;
2567d8d0e25Snbeams       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";
2657d8d0e25Snbeams         } 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";
2717d8d0e25Snbeams         }
2727d8d0e25Snbeams         break;
2737d8d0e25Snbeams       case CEED_EVAL_WEIGHT:
2747d8d0e25Snbeams         break;  // No action
2757d8d0e25Snbeams       case CEED_EVAL_DIV:
2767d8d0e25Snbeams         break;  // TODO: Not implemented
2777d8d0e25Snbeams       case CEED_EVAL_CURL:
2787d8d0e25Snbeams         break;  // TODO: Not implemented
2797d8d0e25Snbeams     }
2807d8d0e25Snbeams   }
2817d8d0e25Snbeams 
2827d8d0e25Snbeams   code << "\n  // -- Output field constants and basis data --\n";
2839e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
2847d8d0e25Snbeams     code << "  // ---- Output field " << i << " ----\n";
2859e201c85SYohann     // Get elem_size, eval_mode, num_comp
2862b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &Erestrict));
2872b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size));
2882b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
2892b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp));
2907d8d0e25Snbeams 
2917d8d0e25Snbeams     // Set field constants
2922b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
293*356036faSJeremy 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";
2967d8d0e25Snbeams     } else {
2979e201c85SYohann       code << "  const CeedInt P_out_" << i << " = " << Q_1d << ";\n";
2987d8d0e25Snbeams     }
2999e201c85SYohann     code << "  const CeedInt num_comp_out_" << i << " = " << num_comp << ";\n";
3007d8d0e25Snbeams 
3017d8d0e25Snbeams     // Load basis data
3029e201c85SYohann     code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
3039e201c85SYohann     switch (eval_mode) {
3047d8d0e25Snbeams       case CEED_EVAL_NONE:
3057d8d0e25Snbeams         break;  // No action
3067d8d0e25Snbeams       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";
3117d8d0e25Snbeams         break;
3127d8d0e25Snbeams       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";
3217d8d0e25Snbeams         } 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";
3277d8d0e25Snbeams         }
3287d8d0e25Snbeams         break;
3297d8d0e25Snbeams       // LCOV_EXCL_START
3307d8d0e25Snbeams       case CEED_EVAL_WEIGHT: {
3317d8d0e25Snbeams         Ceed ceed;
3322b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
3332b730f8bSJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
3347d8d0e25Snbeams         break;  // Should not occur
3357d8d0e25Snbeams       }
3367d8d0e25Snbeams       case CEED_EVAL_DIV:
3377d8d0e25Snbeams         break;  // TODO: Not implemented
3387d8d0e25Snbeams       case CEED_EVAL_CURL:
3397d8d0e25Snbeams         break;  // TODO: Not implemented
3407d8d0e25Snbeams                 // LCOV_EXCL_STOP
3417d8d0e25Snbeams     }
3427d8d0e25Snbeams   }
3437d8d0e25Snbeams   code << "\n  // -- Element loop --\n";
3447d8d0e25Snbeams   code << "  __syncthreads();\n";
3459e201c85SYohann   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
3467d8d0e25Snbeams   // Input basis apply if needed
3477d8d0e25Snbeams   // Generate the correct eval mode code for each input
3487d8d0e25Snbeams   code << "    // -- Input field restrictions and basis actions --\n";
3499e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
3507d8d0e25Snbeams     code << "    // ---- Input field " << i << " ----\n";
3519e201c85SYohann     // Get elem_size, eval_mode, num_comp
3522b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &Erestrict));
3532b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size));
3542b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
3552b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp));
3567d8d0e25Snbeams 
3577d8d0e25Snbeams     // Restriction
3582b730f8bSJeremy L Thompson     if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_collograd_parallelization)) {
3599e201c85SYohann       code << "    CeedScalar r_u_" << i << "[num_comp_in_" << i << "*P_in_" << i << "];\n";
3607d8d0e25Snbeams 
3619e201c85SYohann       bool is_strided;
3622b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionIsStrided(Erestrict, &is_strided));
3639e201c85SYohann       if (!is_strided) {
3642b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetLVectorSize(Erestrict, &lsize));
3657d8d0e25Snbeams         code << "    const CeedInt lsize_in_" << i << " = " << lsize << ";\n";
3669e201c85SYohann         CeedInt comp_stride;
3672b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetCompStride(Erestrict, &comp_stride));
3689e201c85SYohann         code << "    // CompStride: " << comp_stride << "\n";
3692b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetData(Erestrict, &restr_data));
3709e201c85SYohann         data->indices.inputs[i] = restr_data->d_ind;
3712b730f8bSJeremy L Thompson         code << "    readDofsOffset" << dim << "d<num_comp_in_" << i << ", " << comp_stride << ", P_in_" << i << ">(data, lsize_in_" << i
3722b730f8bSJeremy L Thompson              << ", elem, indices.inputs[" << i << "], d_u_" << i << ", r_u_" << i << ");\n";
3737d8d0e25Snbeams       } else {
3749e201c85SYohann         bool has_backend_strides;
3752b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionHasBackendStrides(Erestrict, &has_backend_strides));
3769e201c85SYohann         CeedInt num_elem;
3772b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumElements(Erestrict, &num_elem));
3789e201c85SYohann         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
3799e201c85SYohann         if (!has_backend_strides) {
3802b730f8bSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetStrides(Erestrict, &strides));
3817d8d0e25Snbeams         }
3827d8d0e25Snbeams         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
3832b730f8bSJeremy L Thompson         code << "    readDofsStrided" << dim << "d<num_comp_in_" << i << ",P_in_" << i << "," << strides[0] << "," << strides[1] << "," << strides[2]
3842b730f8bSJeremy L Thompson              << ">(data, elem, d_u_" << i << ", r_u_" << i << ");\n";
3857d8d0e25Snbeams       }
3867d8d0e25Snbeams     }
3877d8d0e25Snbeams 
388b2165e7aSSebastian Grimberg     // TODO: put in a function?
3897d8d0e25Snbeams     // Basis action
3909e201c85SYohann     code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
3919e201c85SYohann     switch (eval_mode) {
3927d8d0e25Snbeams       case CEED_EVAL_NONE:
3939e201c85SYohann         if (!use_collograd_parallelization) {
3949e201c85SYohann           code << "    CeedScalar* r_t_" << i << " = r_u_" << i << ";\n";
3957d8d0e25Snbeams         }
3967d8d0e25Snbeams         break;
3977d8d0e25Snbeams       case CEED_EVAL_INTERP:
3989e201c85SYohann         code << "    CeedScalar r_t_" << i << "[num_comp_in_" << i << "*Q_1d];\n";
3992b730f8bSJeremy 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_"
4002b730f8bSJeremy L Thompson              << i << ", r_t_" << i << ");\n";
4017d8d0e25Snbeams         break;
4027d8d0e25Snbeams       case CEED_EVAL_GRAD:
4039e201c85SYohann         if (use_collograd_parallelization) {
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
4062b730f8bSJeremy L Thompson                << ", s_B_in_" << i << ", r_t_" << i << ");\n";
4077d8d0e25Snbeams         } else {
4089e201c85SYohann           CeedInt P_1d;
4092b730f8bSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
4102b730f8bSJeremy L Thompson           CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
4119e201c85SYohann           code << "    CeedScalar r_t_" << i << "[num_comp_in_" << i << "*dim*Q_1d];\n";
4122b730f8bSJeremy L Thompson           code << "    Grad" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp_in_" << i
4132b730f8bSJeremy L Thompson                << ",P_in_" << i << ",Q_1d>(data, r_u_" << i << ", s_B_in_" << i << ", s_G_in_" << i << ", r_t_" << i << ");\n";
4147d8d0e25Snbeams         }
4157d8d0e25Snbeams         break;
4167d8d0e25Snbeams       case CEED_EVAL_WEIGHT:
4179e201c85SYohann         code << "    CeedScalar r_t_" << i << "[Q_1d];\n";
4182b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
4192b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
420437930d1SJeremy L Thompson         data->W = basis_data->d_q_weight_1d;
4219e201c85SYohann         code << "    Weight" << (dim > 1 ? "Tensor" : "") << dim << "d<Q_1d>(data, W, r_t_" << i << ");\n";
4227d8d0e25Snbeams         break;  // No action
4237d8d0e25Snbeams       case CEED_EVAL_DIV:
4247d8d0e25Snbeams         break;  // TODO: Not implemented
4257d8d0e25Snbeams       case CEED_EVAL_CURL:
4267d8d0e25Snbeams         break;  // TODO: Not implemented
4277d8d0e25Snbeams     }
4287d8d0e25Snbeams   }
4297d8d0e25Snbeams 
430b2165e7aSSebastian Grimberg   // TODO: put in a function + separate collograd logic
4317d8d0e25Snbeams   // Q function
4327d8d0e25Snbeams   code << "\n    // -- Output field setup --\n";
4339e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
4347d8d0e25Snbeams     code << "\n    // ---- Output field " << i << " ----\n";
4352b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
4362b730f8bSJeremy L Thompson     if (eval_mode == CEED_EVAL_GRAD) {
4379e201c85SYohann       if (use_collograd_parallelization) {
4387d8d0e25Snbeams         // Accumulator for gradient slices
4399e201c85SYohann         code << "    CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*Q_1d];\n";
4409e201c85SYohann         code << "    for (CeedInt i = 0; i < num_comp_out_" << i << "; i++) {\n";
4419e201c85SYohann         code << "      for (CeedInt j = 0; j < Q_1d; ++j) {\n";
4429e201c85SYohann         code << "        r_tt_" << i << "[j + i*Q_1d] = 0.0;\n";
4437d8d0e25Snbeams         code << "      }\n";
4447d8d0e25Snbeams         code << "    }\n";
4457d8d0e25Snbeams       } else {
4469e201c85SYohann         code << "    CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*dim*Q_1d];\n";
4477d8d0e25Snbeams       }
4487d8d0e25Snbeams     }
4492b730f8bSJeremy L Thompson     if (eval_mode == CEED_EVAL_NONE || eval_mode == CEED_EVAL_INTERP) {
4509e201c85SYohann       code << "    CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*Q_1d];\n";
4517d8d0e25Snbeams     }
4527d8d0e25Snbeams   }
4537d8d0e25Snbeams   // We treat quadrature points per slice in 3d to save registers
4549e201c85SYohann   if (use_collograd_parallelization) {
4559e201c85SYohann     code << "\n    // Note: Using planes of 3D elements\n";
4567d8d0e25Snbeams     code << "#pragma unroll\n";
4579e201c85SYohann     code << "    for (CeedInt q = 0; q < Q_1d; q++) {\n";
4587d8d0e25Snbeams     code << "      // -- Input fields --\n";
4599e201c85SYohann     for (CeedInt i = 0; i < num_input_fields; i++) {
4607d8d0e25Snbeams       code << "      // ---- Input field " << i << " ----\n";
4619e201c85SYohann       // Get elem_size, eval_mode, num_comp
4622b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
4637d8d0e25Snbeams       // Basis action
4649e201c85SYohann       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
4659e201c85SYohann       switch (eval_mode) {
4667d8d0e25Snbeams         case CEED_EVAL_NONE:
4679e201c85SYohann           code << "      CeedScalar r_q_" << i << "[num_comp_in_" << i << "];\n";
4687d8d0e25Snbeams 
4699e201c85SYohann           bool is_strided;
4702b730f8bSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &Erestrict));
4712b730f8bSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionIsStrided(Erestrict, &is_strided));
4729e201c85SYohann           if (!is_strided) {
4732b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetLVectorSize(Erestrict, &lsize));
4747d8d0e25Snbeams             code << "      const CeedInt lsize_in_" << i << " = " << lsize << ";\n";
4759e201c85SYohann             CeedInt comp_stride;
4762b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetCompStride(Erestrict, &comp_stride));
4779e201c85SYohann             code << "      // CompStride: " << comp_stride << "\n";
4782b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetData(Erestrict, &restr_data));
4799e201c85SYohann             data->indices.inputs[i] = restr_data->d_ind;
4802b730f8bSJeremy L Thompson             code << "      readSliceQuadsOffset"
4812b730f8bSJeremy L Thompson                  << "3d<num_comp_in_" << i << ", " << comp_stride << ", Q_1d>(data, lsize_in_" << i << ", elem, q, indices.inputs[" << i << "], d_u_"
4822b730f8bSJeremy L Thompson                  << i << ", r_q_" << i << ");\n";
4837d8d0e25Snbeams           } else {
4842b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size));
4859e201c85SYohann             bool has_backend_strides;
4862b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionHasBackendStrides(Erestrict, &has_backend_strides));
4879e201c85SYohann             CeedInt num_elem;
4882b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetNumElements(Erestrict, &num_elem));
4899e201c85SYohann             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
4909e201c85SYohann             if (!has_backend_strides) {
4912b730f8bSJeremy L Thompson               CeedCallBackend(CeedElemRestrictionGetStrides(Erestrict, &strides));
4927d8d0e25Snbeams             }
4937d8d0e25Snbeams             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
4942b730f8bSJeremy L Thompson             code << "      readSliceQuadsStrided"
4952b730f8bSJeremy L Thompson                  << "3d<num_comp_in_" << i
4962b730f8bSJeremy L Thompson                  << ",Q_1d"
4972b730f8bSJeremy L Thompson                     ","
4982b730f8bSJeremy L Thompson                  << strides[0] << "," << strides[1] << "," << strides[2] << ">(data, elem, q, d_u_" << i << ", r_q_" << i << ");\n";
4997d8d0e25Snbeams           }
5007d8d0e25Snbeams           break;
5017d8d0e25Snbeams         case CEED_EVAL_INTERP:
5029e201c85SYohann           code << "      CeedScalar r_q_" << i << "[num_comp_in_" << i << "];\n";
5039e201c85SYohann           code << "      for (CeedInt j = 0; j < num_comp_in_" << i << " ; ++j) {\n";
5049e201c85SYohann           code << "        r_q_" << i << "[j] = r_t_" << i << "[q + j*Q_1d];\n";
5057d8d0e25Snbeams           code << "      }\n";
5067d8d0e25Snbeams           break;
5077d8d0e25Snbeams         case CEED_EVAL_GRAD:
5089e201c85SYohann           code << "      CeedScalar r_q_" << i << "[num_comp_in_" << i << "*dim];\n";
5099e201c85SYohann           code << "      gradCollo3d<num_comp_in_" << i << ",Q_1d>(data, q, r_t_" << i << ", s_G_in_" << i << ", r_q_" << i << ");\n";
5107d8d0e25Snbeams           break;
5117d8d0e25Snbeams         case CEED_EVAL_WEIGHT:
5129e201c85SYohann           code << "      CeedScalar r_q_" << i << "[1];\n";
5139e201c85SYohann           code << "      r_q_" << i << "[0] = r_t_" << i << "[q];\n";
5147d8d0e25Snbeams           break;  // No action
5157d8d0e25Snbeams         case CEED_EVAL_DIV:
5167d8d0e25Snbeams           break;  // TODO: Not implemented
5177d8d0e25Snbeams         case CEED_EVAL_CURL:
5187d8d0e25Snbeams           break;  // TODO: Not implemented
5197d8d0e25Snbeams       }
5207d8d0e25Snbeams     }
5217d8d0e25Snbeams     code << "\n      // -- Output fields --\n";
5229e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
5237d8d0e25Snbeams       code << "      // ---- Output field " << i << " ----\n";
5242b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
5257d8d0e25Snbeams       // Basis action
5269e201c85SYohann       switch (eval_mode) {
5277d8d0e25Snbeams         case CEED_EVAL_NONE:
5289e201c85SYohann           code << "      CeedScalar r_qq_" << i << "[num_comp_out_" << i << "];\n";
5297d8d0e25Snbeams           break;  // No action
5307d8d0e25Snbeams         case CEED_EVAL_INTERP:
5319e201c85SYohann           code << "      CeedScalar r_qq_" << i << "[num_comp_out_" << i << "];\n";
5327d8d0e25Snbeams           break;
5337d8d0e25Snbeams         case CEED_EVAL_GRAD:
5349e201c85SYohann           code << "      CeedScalar r_qq_" << i << "[num_comp_out_" << i << "*dim];\n";
5357d8d0e25Snbeams           break;
5367d8d0e25Snbeams         case CEED_EVAL_WEIGHT:
5377d8d0e25Snbeams           break;  // Should not occur
5387d8d0e25Snbeams         case CEED_EVAL_DIV:
5397d8d0e25Snbeams           break;  // TODO: Not implemented
5407d8d0e25Snbeams         case CEED_EVAL_CURL:
5417d8d0e25Snbeams           break;  // TODO: Not implemented
5427d8d0e25Snbeams       }
5437d8d0e25Snbeams     }
5447d8d0e25Snbeams   } else {
5459e201c85SYohann     code << "\n      // Note: Using full elements\n";
5467d8d0e25Snbeams     code << "      // -- Input fields --\n";
5479e201c85SYohann     for (CeedInt i = 0; i < num_input_fields; i++) {
5487d8d0e25Snbeams       code << "      // ---- Input field " << i << " ----\n";
5499e201c85SYohann       code << "      CeedScalar* r_q_" << i << " = r_t_" << i << ";\n";
5507d8d0e25Snbeams     }
5517d8d0e25Snbeams     code << "      // -- Output fields --\n";
5529e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
5537d8d0e25Snbeams       code << "      // ---- Output field " << i << " ----\n";
5549e201c85SYohann       code << "      CeedScalar* r_qq_" << i << " = r_tt_" << i << ";\n";
5557d8d0e25Snbeams     }
5567d8d0e25Snbeams   }
5577d8d0e25Snbeams   code << "\n      // -- QFunction Inputs and outputs --\n";
5589e201c85SYohann   code << "      CeedScalar* in[" << num_input_fields << "];\n";
5599e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
5607d8d0e25Snbeams     code << "      // ---- Input field " << i << " ----\n";
5619e201c85SYohann     code << "      in[" << i << "] = r_q_" << i << ";\n";
5627d8d0e25Snbeams   }
5639e201c85SYohann   code << "      CeedScalar* out[" << num_output_fields << "];\n";
5649e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
5657d8d0e25Snbeams     code << "      // ---- Output field " << i << " ----\n";
5669e201c85SYohann     code << "      out[" << i << "] = r_qq_" << i << ";\n";
5677d8d0e25Snbeams   }
5687d8d0e25Snbeams   code << "\n      // -- Apply QFunction --\n";
5699e201c85SYohann   code << "      " << q_function_name << "(ctx, ";
5709e201c85SYohann   if (dim != 3 || use_collograd_parallelization) {
5717d8d0e25Snbeams     code << "1";
5727d8d0e25Snbeams   } else {
5739e201c85SYohann     code << "Q_1d";
5747d8d0e25Snbeams   }
5757d8d0e25Snbeams   code << ", in, out);\n";
5769e201c85SYohann   if (use_collograd_parallelization) {
5777d8d0e25Snbeams     code << "      // -- Output fields --\n";
5789e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
5797d8d0e25Snbeams       code << "      // ---- Output field " << i << " ----\n";
5802b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
5817d8d0e25Snbeams       // Basis action
5829e201c85SYohann       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
5839e201c85SYohann       switch (eval_mode) {
5847d8d0e25Snbeams         case CEED_EVAL_NONE:
5859e201c85SYohann           code << "      for (CeedInt j = 0; j < num_comp_out_" << i << " ; ++j) {\n";
5869e201c85SYohann           code << "        r_tt_" << i << "[q + j*Q_1d] = r_qq_" << i << "[j];\n";
5877d8d0e25Snbeams           code << "      }\n";
5887d8d0e25Snbeams           break;  // No action
5897d8d0e25Snbeams         case CEED_EVAL_INTERP:
5909e201c85SYohann           code << "      for (CeedInt j = 0; j < num_comp_out_" << i << " ; ++j) {\n";
5919e201c85SYohann           code << "        r_tt_" << i << "[q + j*Q_1d] = r_qq_" << i << "[j];\n";
5927d8d0e25Snbeams           code << "      }\n";
5937d8d0e25Snbeams           break;
5947d8d0e25Snbeams         case CEED_EVAL_GRAD:
5959e201c85SYohann           code << "      gradColloTranspose3d<num_comp_out_" << i << ",Q_1d>(data, q, r_qq_" << i << ", s_G_out_" << i << ", r_tt_" << i << ");\n";
5967d8d0e25Snbeams           break;
5977d8d0e25Snbeams         case CEED_EVAL_WEIGHT:
5987d8d0e25Snbeams           break;  // Should not occur
5997d8d0e25Snbeams         case CEED_EVAL_DIV:
6007d8d0e25Snbeams           break;  // TODO: Not implemented
6017d8d0e25Snbeams         case CEED_EVAL_CURL:
6027d8d0e25Snbeams           break;  // TODO: Not implemented
6037d8d0e25Snbeams       }
6047d8d0e25Snbeams     }
6057d8d0e25Snbeams     code << "    }\n";
6067d8d0e25Snbeams   }
6077d8d0e25Snbeams 
6087d8d0e25Snbeams   // Output basis apply if needed
6097d8d0e25Snbeams   // Generate the correct eval mode code for each output
6107d8d0e25Snbeams   code << "\n    // -- Output field basis action and restrictions --\n";
6119e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
6127d8d0e25Snbeams     code << "    // ---- Output field " << i << " ----\n";
6139e201c85SYohann     // Get elem_size, eval_mode, num_comp
6142b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &Erestrict));
6152b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size));
6162b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
6172b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp));
618b2165e7aSSebastian Grimberg     // TODO put in a function
6197d8d0e25Snbeams     // Basis action
6209e201c85SYohann     code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
6219e201c85SYohann     switch (eval_mode) {
6227d8d0e25Snbeams       case CEED_EVAL_NONE:
6239e201c85SYohann         code << "    CeedScalar* r_v_" << i << " = r_tt_" << i << ";\n";
6247d8d0e25Snbeams         break;  // No action
6257d8d0e25Snbeams       case CEED_EVAL_INTERP:
6269e201c85SYohann         code << "    CeedScalar r_v_" << i << "[num_comp_out_" << i << "*P_out_" << i << "];\n";
6272b730f8bSJeremy L Thompson         code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_out_" << i << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i
6282b730f8bSJeremy L Thompson              << ", s_B_out_" << i << ", r_v_" << i << ");\n";
6297d8d0e25Snbeams         break;
6307d8d0e25Snbeams       case CEED_EVAL_GRAD:
6319e201c85SYohann         code << "    CeedScalar r_v_" << i << "[num_comp_out_" << i << "*P_out_" << i << "];\n";
6329e201c85SYohann         if (use_collograd_parallelization) {
6332b730f8bSJeremy L Thompson           code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_out_" << i << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i
6342b730f8bSJeremy L Thompson                << ", s_B_out_" << i << ", r_v_" << i << ");\n";
6357d8d0e25Snbeams         } else {
6369e201c85SYohann           CeedInt P_1d;
6372b730f8bSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
6382b730f8bSJeremy L Thompson           CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
6392b730f8bSJeremy L Thompson           code << "    GradTranspose" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp_out_" << i
6402b730f8bSJeremy L Thompson                << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i << ", s_B_out_" << i << ", s_G_out_" << i << ", r_v_" << i << ");\n";
6417d8d0e25Snbeams         }
6427d8d0e25Snbeams         break;
6437d8d0e25Snbeams       // LCOV_EXCL_START
6447d8d0e25Snbeams       case CEED_EVAL_WEIGHT: {
6457d8d0e25Snbeams         Ceed ceed;
6462b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
6472b730f8bSJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
6487d8d0e25Snbeams         break;  // Should not occur
6497d8d0e25Snbeams       }
6507d8d0e25Snbeams       case CEED_EVAL_DIV:
6517d8d0e25Snbeams         break;  // TODO: Not implemented
6527d8d0e25Snbeams       case CEED_EVAL_CURL:
6537d8d0e25Snbeams         break;  // TODO: Not implemented
6547d8d0e25Snbeams                 // LCOV_EXCL_STOP
6557d8d0e25Snbeams     }
656b2165e7aSSebastian Grimberg     // TODO put in a function
6577d8d0e25Snbeams     // Restriction
6589e201c85SYohann     bool is_strided;
6592b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionIsStrided(Erestrict, &is_strided));
6609e201c85SYohann     if (!is_strided) {
6612b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetLVectorSize(Erestrict, &lsize));
6627d8d0e25Snbeams       code << "    const CeedInt lsize_out_" << i << " = " << lsize << ";\n";
6639e201c85SYohann       CeedInt comp_stride;
6642b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetCompStride(Erestrict, &comp_stride));
6659e201c85SYohann       code << "    // CompStride: " << comp_stride << "\n";
6662b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetData(Erestrict, &restr_data));
6679e201c85SYohann       data->indices.outputs[i] = restr_data->d_ind;
6682b730f8bSJeremy L Thompson       code << "    writeDofsOffset" << dim << "d<num_comp_out_" << i << ", " << comp_stride << ", P_out_" << i << ">(data, lsize_out_" << i
6692b730f8bSJeremy L Thompson            << ", elem, indices.outputs[" << i << "], r_v_" << i << ", d_v_" << i << ");\n";
6707d8d0e25Snbeams     } else {
6719e201c85SYohann       bool has_backend_strides;
6722b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionHasBackendStrides(Erestrict, &has_backend_strides));
6739e201c85SYohann       CeedInt num_elem;
6742b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetNumElements(Erestrict, &num_elem));
6759e201c85SYohann       CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
6769e201c85SYohann       if (!has_backend_strides) {
6772b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetStrides(Erestrict, &strides));
6787d8d0e25Snbeams       }
6797d8d0e25Snbeams       code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
6802b730f8bSJeremy L Thompson       code << "    writeDofsStrided" << dim << "d<num_comp_out_" << i << ",P_out_" << i << "," << strides[0] << "," << strides[1] << "," << strides[2]
6812b730f8bSJeremy L Thompson            << ">(data, elem, r_v_" << i << ", d_v_" << i << ");\n";
6827d8d0e25Snbeams     }
6837d8d0e25Snbeams   }
6847d8d0e25Snbeams 
6857d8d0e25Snbeams   code << "  }\n";
6867d8d0e25Snbeams   code << "}\n";
6877d8d0e25Snbeams   code << "// -----------------------------------------------------------------------------\n\n";
6887d8d0e25Snbeams 
6897d8d0e25Snbeams   // View kernel for debugging
69023d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "Generated Operator Kernels:\n");
6913f21f6b1SJeremy L Thompson   CeedDebug(ceed, code.str().c_str());
6927d8d0e25Snbeams 
693539ec17dSJeremy L Thompson   CeedInt block_sizes[3] = {0, 0, 0};
6949e201c85SYohann   CeedInt num_elem;
6952b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
6962b730f8bSJeremy L Thompson   CeedCallBackend(BlockGridCalculate_Hip_gen(dim, num_elem, data->max_P_1d, Q_1d, block_sizes));
697eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Hip(ceed, code.str().c_str(), &data->module, 2, "T_1D", block_sizes[0], "BLOCK_SIZE",
6982b730f8bSJeremy L Thompson                                   block_sizes[0] * block_sizes[1] * block_sizes[2]));
699eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, operator_name.c_str(), &data->op));
7007d8d0e25Snbeams 
7012b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
702e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7037d8d0e25Snbeams }
7042a86cc9dSSebastian Grimberg 
7057d8d0e25Snbeams //------------------------------------------------------------------------------
706