xref: /libCEED/rust/libceed-sys/c-src/backends/hip-gen/ceed-hip-gen-operator-build.cpp (revision eb7e6cafeb5d6cc94d59355f95e7bc9ae3fc1c25)
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"
207d8d0e25Snbeams #include "../hip/ceed-hip-compile.h"
212b730f8bSJeremy L Thompson #include "ceed-hip-gen.h"
22b3e1519bSnbeams 
23b3e1519bSnbeams //------------------------------------------------------------------------------
24b3e1519bSnbeams // Calculate the block size used for launching the operator kernel
25b3e1519bSnbeams //------------------------------------------------------------------------------
262b730f8bSJeremy 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) {
279e201c85SYohann   const CeedInt thread1d = CeedIntMax(Q_1d, P_1d);
28b3e1519bSnbeams   if (dim == 1) {
299e201c85SYohann     CeedInt elems_per_block = 64 * thread1d > 256 ? 256 / thread1d : 64;
309e201c85SYohann     elems_per_block         = elems_per_block > 0 ? elems_per_block : 1;
31b3e1519bSnbeams     block_sizes[0]          = thread1d;
32b3e1519bSnbeams     block_sizes[1]          = 1;
339e201c85SYohann     block_sizes[2]          = elems_per_block;
34b3e1519bSnbeams   } else if (dim == 2) {
359e201c85SYohann     const CeedInt elems_per_block = thread1d < 4 ? 16 : 2;
36b3e1519bSnbeams     block_sizes[0]                = thread1d;
37b3e1519bSnbeams     block_sizes[1]                = thread1d;
389e201c85SYohann     block_sizes[2]                = elems_per_block;
39b3e1519bSnbeams   } else if (dim == 3) {
409e201c85SYohann     const CeedInt elems_per_block = thread1d < 6 ? 4 : (thread1d < 8 ? 2 : 1);
41b3e1519bSnbeams     block_sizes[0]                = thread1d;
42b3e1519bSnbeams     block_sizes[1]                = thread1d;
439e201c85SYohann     block_sizes[2]                = elems_per_block;
44b3e1519bSnbeams   }
45b3e1519bSnbeams   return CEED_ERROR_SUCCESS;
46b3e1519bSnbeams }
47b3e1519bSnbeams 
487d8d0e25Snbeams //------------------------------------------------------------------------------
499e201c85SYohann // Build single operator kernel
507d8d0e25Snbeams //------------------------------------------------------------------------------
51*eb7e6cafSJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Hip_gen(CeedOperator op) {
527d8d0e25Snbeams   using std::ostringstream;
537d8d0e25Snbeams   using std::string;
549e201c85SYohann   bool is_setup_done;
552b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
569e201c85SYohann   if (is_setup_done) return CEED_ERROR_SUCCESS;
577d8d0e25Snbeams   Ceed ceed;
582b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
597d8d0e25Snbeams   CeedOperator_Hip_gen *data;
602b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
617d8d0e25Snbeams   CeedQFunction          qf;
627d8d0e25Snbeams   CeedQFunction_Hip_gen *qf_data;
632b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
642b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
65e79b91d9SJeremy L Thompson   CeedSize lsize;
662b730f8bSJeremy L Thompson   CeedInt  Q, P_1d = 0, Q_1d = 0, elem_size, num_input_fields, num_output_fields, num_comp, dim = 1;
672b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
689e201c85SYohann   Q_1d = Q;
699e201c85SYohann   CeedOperatorField *op_input_fields, *op_output_fields;
702b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
719e201c85SYohann   CeedQFunctionField *qf_input_fields, *qf_output_fields;
722b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
739e201c85SYohann   CeedEvalMode             eval_mode;
747d8d0e25Snbeams   CeedBasis                basis;
757d8d0e25Snbeams   CeedBasis_Hip_shared    *basis_data;
767d8d0e25Snbeams   CeedElemRestriction      Erestrict;
777d8d0e25Snbeams   CeedElemRestriction_Hip *restr_data;
787d8d0e25Snbeams 
790b454692Sjeremylt   // Check for restriction only identity operator
800b454692Sjeremylt   bool is_identity_qf;
812b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
820b454692Sjeremylt   if (is_identity_qf) {
839e201c85SYohann     CeedEvalMode eval_mode_in, eval_mode_out;
842b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
852b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
866574a04fSJeremy L Thompson     CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
876574a04fSJeremy L Thompson               "Backend does not implement restriction only identity operators");
880b454692Sjeremylt   }
890b454692Sjeremylt 
907d8d0e25Snbeams   ostringstream code;
919e201c85SYohann   // TODO: generalize to accept different device functions?
929e201c85SYohann   {
939e201c85SYohann     char *tensor_basis_kernel_path, *tensor_basis_kernel_source;
942b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-shared-basis-tensor-templates.h", &tensor_basis_kernel_path));
959e201c85SYohann     CeedDebug256(ceed, 2, "----- Loading Tensor Basis Kernel Source -----\n");
962b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, tensor_basis_kernel_path, &tensor_basis_kernel_source));
979e201c85SYohann     code << tensor_basis_kernel_source;
982b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&tensor_basis_kernel_path));
992b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&tensor_basis_kernel_source));
1009e201c85SYohann   }
1019e201c85SYohann   {
1029e201c85SYohann     char *hip_gen_template_path, *hip_gen_template_source;
1032b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-gen-templates.h", &hip_gen_template_path));
1049e201c85SYohann     CeedDebug256(ceed, 2, "----- Loading Hip-Gen Template Source -----\n");
1052b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, hip_gen_template_path, &hip_gen_template_source));
1069e201c85SYohann     code << hip_gen_template_source;
1072b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&hip_gen_template_path));
1082b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&hip_gen_template_source));
1099e201c85SYohann   }
1107d8d0e25Snbeams 
1119e201c85SYohann   string q_function_source(qf_data->q_function_source);
1129e201c85SYohann   string q_function_name(qf_data->q_function_name);
1139e201c85SYohann   string operator_name;
114204bfdd7SJeremy L Thompson   operator_name = "CeedKernelHipGenOperator_" + q_function_name;
1157d8d0e25Snbeams 
1169e201c85SYohann   // Find dim, P_1d, Q_1d
1179e201c85SYohann   data->max_P_1d = 0;
1189e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1192b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1207d8d0e25Snbeams     if (basis != CEED_BASIS_COLLOCATED) {
1212b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1222b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1237d8d0e25Snbeams 
1249e201c85SYohann       // Collect dim, P_1d, and Q_1d
1252b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &dim));
1267d8d0e25Snbeams       bool isTensor;
1272b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &isTensor));
1286574a04fSJeremy L Thompson       CeedCheck(isTensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
1292b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
1302b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
1319e201c85SYohann       if (P_1d > data->max_P_1d) data->max_P_1d = P_1d;
1327d8d0e25Snbeams     }
1337d8d0e25Snbeams   }
1349e201c85SYohann   // Check output bases for Q_1d, dim as well
135dc64899eSYohann   //   The only input basis might be CEED_BASIS_COLLOCATED
1369e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1372b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1387d8d0e25Snbeams 
1397d8d0e25Snbeams     if (basis != CEED_BASIS_COLLOCATED) {
1402b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1412b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1427d8d0e25Snbeams 
1439e201c85SYohann       // Collect Q_1d
1442b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &dim));
1457d8d0e25Snbeams       bool isTensor;
1462b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &isTensor));
1476574a04fSJeremy L Thompson       CeedCheck(isTensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
1482b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
1497d8d0e25Snbeams     }
1507d8d0e25Snbeams   }
1517d8d0e25Snbeams   data->dim  = dim;
1529e201c85SYohann   data->Q_1d = Q_1d;
1537d8d0e25Snbeams 
1549e201c85SYohann   // Only use 3D collocated gradient parallelization strategy when gradient is computed
1559e201c85SYohann   // TODO: put in a function?
1569e201c85SYohann   bool use_collograd_parallelization = false;
1579e201c85SYohann   if (dim == 3) {
1589e201c85SYohann     bool was_grad_found = false;
1599e201c85SYohann     for (CeedInt i = 0; i < num_input_fields; i++) {
1602b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1619e201c85SYohann       if (eval_mode == CEED_EVAL_GRAD) {
1622b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1632b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1649e201c85SYohann         use_collograd_parallelization = !!basis_data->d_collo_grad_1d && (was_grad_found ? use_collograd_parallelization : true);
1659e201c85SYohann         was_grad_found                = true;
1669e201c85SYohann       }
1679e201c85SYohann     }
1689e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
1692b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1709e201c85SYohann       if (eval_mode == CEED_EVAL_GRAD) {
1712b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1722b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1739e201c85SYohann         use_collograd_parallelization = !!basis_data->d_collo_grad_1d && (was_grad_found ? use_collograd_parallelization : true);
1749e201c85SYohann         was_grad_found                = true;
1759e201c85SYohann       }
1769e201c85SYohann     }
1777d8d0e25Snbeams   }
1787d8d0e25Snbeams 
1799e201c85SYohann   // Define CEED_Q_VLA
1809e201c85SYohann   code << "\n#undef CEED_Q_VLA\n";
1819e201c85SYohann   if (dim != 3 || use_collograd_parallelization) {
1829e201c85SYohann     code << "#define CEED_Q_VLA 1\n\n";
1839e201c85SYohann   } else {
1849e201c85SYohann     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
1859e201c85SYohann   }
1869e201c85SYohann 
1879e201c85SYohann   code << q_function_source;
1887d8d0e25Snbeams 
1897d8d0e25Snbeams   // Setup
1907d8d0e25Snbeams   code << "\n// -----------------------------------------------------------------------------\n";
191b3e1519bSnbeams   code << "\nextern \"C\" __launch_bounds__(BLOCK_SIZE)\n";
1922b730f8bSJeremy L Thompson   code << "__global__ void " << operator_name
1932b730f8bSJeremy L Thompson        << "(CeedInt num_elem, void* ctx, FieldsInt_Hip indices, Fields_Hip fields, Fields_Hip B, Fields_Hip G, CeedScalar* W) {\n";
1949e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1952b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1969e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
1979e201c85SYohann       code << "  const CeedScalar* d_u_" << i << " = fields.inputs[" << i << "];\n";
1987d8d0e25Snbeams     }
1997d8d0e25Snbeams   }
2007d8d0e25Snbeams 
2019e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
2029e201c85SYohann     code << "  CeedScalar* d_v_" << i << " = fields.outputs[" << i << "];\n";
2037d8d0e25Snbeams   }
2047d8d0e25Snbeams 
2059e201c85SYohann   code << "  const CeedInt dim = " << dim << ";\n";
2069e201c85SYohann   code << "  const CeedInt Q_1d = " << Q_1d << ";\n";
2077d8d0e25Snbeams 
2087d8d0e25Snbeams   code << "  HIP_DYNAMIC_SHARED( CeedScalar, slice)\n";
2099e201c85SYohann   code << "  SharedData_Hip data;\n";
2109e201c85SYohann   code << "  data.t_id_x = threadIdx.x;\n";
2119e201c85SYohann   code << "  data.t_id_y = threadIdx.y;\n";
2129e201c85SYohann   code << "  data.t_id_z = threadIdx.z;\n";
2139e201c85SYohann   code << "  data.t_id  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
2149e201c85SYohann   code << "  data.slice = slice+data.t_id_z*T_1D" << (dim > 1 ? "*T_1D" : "") << ";\n";
2157d8d0e25Snbeams 
2167d8d0e25Snbeams   code << "\n  // -- Input field constants and basis data --\n";
2177d8d0e25Snbeams   // Initialize constants, and matrices B and G
2189e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
2197d8d0e25Snbeams     code << "  // ---- Input field " << i << " ----\n";
2209e201c85SYohann     // Get elem_size, eval_mode, num_comp
2212b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &Erestrict));
2222b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size));
2232b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
2242b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp));
2257d8d0e25Snbeams 
2267d8d0e25Snbeams     // Set field constants
2279e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {
2282b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
2297d8d0e25Snbeams       if (basis != CEED_BASIS_COLLOCATED) {
2302b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
2319e201c85SYohann         code << "  const CeedInt P_in_" << i << " = " << P_1d << ";\n";
2327d8d0e25Snbeams       } else {
2339e201c85SYohann         code << "  const CeedInt P_in_" << i << " = " << Q_1d << ";\n";
2347d8d0e25Snbeams       }
2359e201c85SYohann       code << "  const CeedInt num_comp_in_" << i << " = " << num_comp << ";\n";
2367d8d0e25Snbeams     }
2377d8d0e25Snbeams 
2387d8d0e25Snbeams     // Load basis data
2399e201c85SYohann     code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
2409e201c85SYohann     switch (eval_mode) {
2417d8d0e25Snbeams       case CEED_EVAL_NONE:
2427d8d0e25Snbeams         break;
2437d8d0e25Snbeams       case CEED_EVAL_INTERP:
2442b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
2459e201c85SYohann         data->B.inputs[i] = basis_data->d_interp_1d;
2469e201c85SYohann         code << "  __shared__ CeedScalar s_B_in_" << i << "[" << P_1d * Q_1d << "];\n";
2479e201c85SYohann         code << "  loadMatrix<P_in_" << i << ",Q_1d>(data, B.inputs[" << i << "], s_B_in_" << i << ");\n";
2487d8d0e25Snbeams         break;
2497d8d0e25Snbeams       case CEED_EVAL_GRAD:
2502b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
2519e201c85SYohann         data->B.inputs[i] = basis_data->d_interp_1d;
2529e201c85SYohann         code << "  __shared__ CeedScalar s_B_in_" << i << "[" << P_1d * Q_1d << "];\n";
2539e201c85SYohann         code << "  loadMatrix<P_in_" << i << ",Q_1d>(data, B.inputs[" << i << "], s_B_in_" << i << ");\n";
2549e201c85SYohann         if (use_collograd_parallelization) {
2559e201c85SYohann           data->G.inputs[i] = basis_data->d_collo_grad_1d;
2569e201c85SYohann           code << "  __shared__ CeedScalar s_G_in_" << i << "[" << Q_1d * Q_1d << "];\n";
2579e201c85SYohann           code << "  loadMatrix<Q_1d,Q_1d>(data, G.inputs[" << i << "], s_G_in_" << i << ");\n";
2587d8d0e25Snbeams         } else {
2599e201c85SYohann           bool has_collo_grad = !!basis_data->d_collo_grad_1d;
2609e201c85SYohann           data->G.inputs[i]   = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2619e201c85SYohann           code << "  __shared__ CeedScalar s_G_in_" << i << "[" << Q_1d * (has_collo_grad ? Q_1d : P_1d) << "];\n";
2622b730f8bSJeremy 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
2632b730f8bSJeremy L Thompson                << ");\n";
2647d8d0e25Snbeams         }
2657d8d0e25Snbeams         break;
2667d8d0e25Snbeams       case CEED_EVAL_WEIGHT:
2677d8d0e25Snbeams         break;  // No action
2687d8d0e25Snbeams       case CEED_EVAL_DIV:
2697d8d0e25Snbeams         break;  // TODO: Not implemented
2707d8d0e25Snbeams       case CEED_EVAL_CURL:
2717d8d0e25Snbeams         break;  // TODO: Not implemented
2727d8d0e25Snbeams     }
2737d8d0e25Snbeams   }
2747d8d0e25Snbeams 
2757d8d0e25Snbeams   code << "\n  // -- Output field constants and basis data --\n";
2769e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
2777d8d0e25Snbeams     code << "  // ---- Output field " << i << " ----\n";
2789e201c85SYohann     // Get elem_size, eval_mode, num_comp
2792b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &Erestrict));
2802b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size));
2812b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
2822b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp));
2837d8d0e25Snbeams 
2847d8d0e25Snbeams     // Set field constants
2852b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
2867d8d0e25Snbeams     if (basis != CEED_BASIS_COLLOCATED) {
2872b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
2889e201c85SYohann       code << "  const CeedInt P_out_" << i << " = " << P_1d << ";\n";
2897d8d0e25Snbeams     } else {
2909e201c85SYohann       code << "  const CeedInt P_out_" << i << " = " << Q_1d << ";\n";
2917d8d0e25Snbeams     }
2929e201c85SYohann     code << "  const CeedInt num_comp_out_" << i << " = " << num_comp << ";\n";
2937d8d0e25Snbeams 
2947d8d0e25Snbeams     // Load basis data
2959e201c85SYohann     code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
2969e201c85SYohann     switch (eval_mode) {
2977d8d0e25Snbeams       case CEED_EVAL_NONE:
2987d8d0e25Snbeams         break;  // No action
2997d8d0e25Snbeams       case CEED_EVAL_INTERP:
3002b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
3019e201c85SYohann         data->B.outputs[i] = basis_data->d_interp_1d;
3029e201c85SYohann         code << "  __shared__ CeedScalar s_B_out_" << i << "[" << P_1d * Q_1d << "];\n";
3039e201c85SYohann         code << "  loadMatrix<P_out_" << i << ",Q_1d>(data, B.outputs[" << i << "], s_B_out_" << i << ");\n";
3047d8d0e25Snbeams         break;
3057d8d0e25Snbeams       case CEED_EVAL_GRAD:
3062b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
3079e201c85SYohann         data->B.outputs[i] = basis_data->d_interp_1d;
3089e201c85SYohann         code << "  __shared__ CeedScalar s_B_out_" << i << "[" << P_1d * Q_1d << "];\n";
3099e201c85SYohann         code << "  loadMatrix<P_out_" << i << ",Q_1d>(data, B.outputs[" << i << "], s_B_out_" << i << ");\n";
3109e201c85SYohann         if (use_collograd_parallelization) {
3119e201c85SYohann           data->G.outputs[i] = basis_data->d_collo_grad_1d;
3129e201c85SYohann           code << "  __shared__ CeedScalar s_G_out_" << i << "[" << Q_1d * Q_1d << "];\n";
3139e201c85SYohann           code << "  loadMatrix<Q_1d,Q_1d>(data, G.outputs[" << i << "], s_G_out_" << i << ");\n";
3147d8d0e25Snbeams         } else {
3159e201c85SYohann           bool has_collo_grad = !!basis_data->d_collo_grad_1d;
3169e201c85SYohann           data->G.outputs[i]  = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
3179e201c85SYohann           code << "  __shared__ CeedScalar s_G_out_" << i << "[" << Q_1d * (has_collo_grad ? Q_1d : P_1d) << "];\n";
3182b730f8bSJeremy L Thompson           code << "  loadMatrix<" << (has_collo_grad ? "Q_1d" : ("P_out_" + std::to_string(i))) << ",Q_1d>(data, G.outputs[" << i << "], s_G_out_"
3192b730f8bSJeremy L Thompson                << i << ");\n";
3207d8d0e25Snbeams         }
3217d8d0e25Snbeams         break;
3227d8d0e25Snbeams       // LCOV_EXCL_START
3237d8d0e25Snbeams       case CEED_EVAL_WEIGHT: {
3247d8d0e25Snbeams         Ceed ceed;
3252b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
3262b730f8bSJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
3277d8d0e25Snbeams         break;  // Should not occur
3287d8d0e25Snbeams       }
3297d8d0e25Snbeams       case CEED_EVAL_DIV:
3307d8d0e25Snbeams         break;  // TODO: Not implemented
3317d8d0e25Snbeams       case CEED_EVAL_CURL:
3327d8d0e25Snbeams         break;  // TODO: Not implemented
3337d8d0e25Snbeams                 // LCOV_EXCL_STOP
3347d8d0e25Snbeams     }
3357d8d0e25Snbeams   }
3367d8d0e25Snbeams   code << "\n  // -- Element loop --\n";
3377d8d0e25Snbeams   code << "  __syncthreads();\n";
3389e201c85SYohann   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
3397d8d0e25Snbeams   // Input basis apply if needed
3407d8d0e25Snbeams   // Generate the correct eval mode code for each input
3417d8d0e25Snbeams   code << "    // -- Input field restrictions and basis actions --\n";
3429e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
3437d8d0e25Snbeams     code << "    // ---- Input field " << i << " ----\n";
3449e201c85SYohann     // Get elem_size, eval_mode, num_comp
3452b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &Erestrict));
3462b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size));
3472b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
3482b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp));
3497d8d0e25Snbeams 
3507d8d0e25Snbeams     // Restriction
3512b730f8bSJeremy L Thompson     if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_collograd_parallelization)) {
3529e201c85SYohann       code << "    CeedScalar r_u_" << i << "[num_comp_in_" << i << "*P_in_" << i << "];\n";
3537d8d0e25Snbeams 
3549e201c85SYohann       bool is_strided;
3552b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionIsStrided(Erestrict, &is_strided));
3569e201c85SYohann       if (!is_strided) {
3572b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetLVectorSize(Erestrict, &lsize));
3587d8d0e25Snbeams         code << "    const CeedInt lsize_in_" << i << " = " << lsize << ";\n";
3599e201c85SYohann         CeedInt comp_stride;
3602b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetCompStride(Erestrict, &comp_stride));
3619e201c85SYohann         code << "    // CompStride: " << comp_stride << "\n";
3622b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetData(Erestrict, &restr_data));
3639e201c85SYohann         data->indices.inputs[i] = restr_data->d_ind;
3642b730f8bSJeremy L Thompson         code << "    readDofsOffset" << dim << "d<num_comp_in_" << i << ", " << comp_stride << ", P_in_" << i << ">(data, lsize_in_" << i
3652b730f8bSJeremy L Thompson              << ", elem, indices.inputs[" << i << "], d_u_" << i << ", r_u_" << i << ");\n";
3667d8d0e25Snbeams       } else {
3679e201c85SYohann         bool has_backend_strides;
3682b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionHasBackendStrides(Erestrict, &has_backend_strides));
3699e201c85SYohann         CeedInt num_elem;
3702b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumElements(Erestrict, &num_elem));
3719e201c85SYohann         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
3729e201c85SYohann         if (!has_backend_strides) {
3732b730f8bSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetStrides(Erestrict, &strides));
3747d8d0e25Snbeams         }
3757d8d0e25Snbeams         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
3762b730f8bSJeremy L Thompson         code << "    readDofsStrided" << dim << "d<num_comp_in_" << i << ",P_in_" << i << "," << strides[0] << "," << strides[1] << "," << strides[2]
3772b730f8bSJeremy L Thompson              << ">(data, elem, d_u_" << i << ", r_u_" << i << ");\n";
3787d8d0e25Snbeams       }
3797d8d0e25Snbeams     }
3807d8d0e25Snbeams 
3817d8d0e25Snbeams     // Basis action
3829e201c85SYohann     code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
3839e201c85SYohann     switch (eval_mode) {
3847d8d0e25Snbeams       case CEED_EVAL_NONE:
3859e201c85SYohann         if (!use_collograd_parallelization) {
3869e201c85SYohann           code << "    CeedScalar* r_t_" << i << " = r_u_" << i << ";\n";
3877d8d0e25Snbeams         }
3887d8d0e25Snbeams         break;
3897d8d0e25Snbeams       case CEED_EVAL_INTERP:
3909e201c85SYohann         code << "    CeedScalar r_t_" << i << "[num_comp_in_" << i << "*Q_1d];\n";
3912b730f8bSJeremy 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_"
3922b730f8bSJeremy L Thompson              << i << ", r_t_" << i << ");\n";
3937d8d0e25Snbeams         break;
3947d8d0e25Snbeams       case CEED_EVAL_GRAD:
3959e201c85SYohann         if (use_collograd_parallelization) {
3969e201c85SYohann           code << "    CeedScalar r_t_" << i << "[num_comp_in_" << i << "*Q_1d];\n";
3972b730f8bSJeremy L Thompson           code << "    Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_in_" << i << ",P_in_" << i << ",Q_1d>(data, r_u_" << i
3982b730f8bSJeremy L Thompson                << ", s_B_in_" << i << ", r_t_" << i << ");\n";
3997d8d0e25Snbeams         } else {
4009e201c85SYohann           CeedInt P_1d;
4012b730f8bSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
4022b730f8bSJeremy L Thompson           CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
4039e201c85SYohann           code << "    CeedScalar r_t_" << i << "[num_comp_in_" << i << "*dim*Q_1d];\n";
4042b730f8bSJeremy L Thompson           code << "    Grad" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp_in_" << i
4052b730f8bSJeremy L Thompson                << ",P_in_" << i << ",Q_1d>(data, r_u_" << i << ", s_B_in_" << i << ", s_G_in_" << i << ", r_t_" << i << ");\n";
4067d8d0e25Snbeams         }
4077d8d0e25Snbeams         break;
4087d8d0e25Snbeams       case CEED_EVAL_WEIGHT:
4099e201c85SYohann         code << "    CeedScalar r_t_" << i << "[Q_1d];\n";
4102b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
4112b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
412437930d1SJeremy L Thompson         data->W = basis_data->d_q_weight_1d;
4139e201c85SYohann         code << "    Weight" << (dim > 1 ? "Tensor" : "") << dim << "d<Q_1d>(data, W, r_t_" << i << ");\n";
4147d8d0e25Snbeams         break;  // No action
4157d8d0e25Snbeams       case CEED_EVAL_DIV:
4167d8d0e25Snbeams         break;  // TODO: Not implemented
4177d8d0e25Snbeams       case CEED_EVAL_CURL:
4187d8d0e25Snbeams         break;  // TODO: Not implemented
4197d8d0e25Snbeams     }
4207d8d0e25Snbeams   }
4217d8d0e25Snbeams 
4227d8d0e25Snbeams   // Q function
4237d8d0e25Snbeams   code << "\n    // -- Output field setup --\n";
4249e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
4257d8d0e25Snbeams     code << "\n    // ---- Output field " << i << " ----\n";
4262b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
4272b730f8bSJeremy L Thompson     if (eval_mode == CEED_EVAL_GRAD) {
4289e201c85SYohann       if (use_collograd_parallelization) {
4297d8d0e25Snbeams         // Accumulator for gradient slices
4309e201c85SYohann         code << "    CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*Q_1d];\n";
4319e201c85SYohann         code << "    for (CeedInt i = 0; i < num_comp_out_" << i << "; i++) {\n";
4329e201c85SYohann         code << "      for (CeedInt j = 0; j < Q_1d; ++j) {\n";
4339e201c85SYohann         code << "        r_tt_" << i << "[j + i*Q_1d] = 0.0;\n";
4347d8d0e25Snbeams         code << "      }\n";
4357d8d0e25Snbeams         code << "    }\n";
4367d8d0e25Snbeams       } else {
4379e201c85SYohann         code << "    CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*dim*Q_1d];\n";
4387d8d0e25Snbeams       }
4397d8d0e25Snbeams     }
4402b730f8bSJeremy L Thompson     if (eval_mode == CEED_EVAL_NONE || eval_mode == CEED_EVAL_INTERP) {
4419e201c85SYohann       code << "    CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*Q_1d];\n";
4427d8d0e25Snbeams     }
4437d8d0e25Snbeams   }
4447d8d0e25Snbeams   // We treat quadrature points per slice in 3d to save registers
4459e201c85SYohann   if (use_collograd_parallelization) {
4469e201c85SYohann     code << "\n    // Note: Using planes of 3D elements\n";
4477d8d0e25Snbeams     code << "#pragma unroll\n";
4489e201c85SYohann     code << "    for (CeedInt q = 0; q < Q_1d; q++) {\n";
4497d8d0e25Snbeams     code << "      // -- Input fields --\n";
4509e201c85SYohann     for (CeedInt i = 0; i < num_input_fields; i++) {
4517d8d0e25Snbeams       code << "      // ---- Input field " << i << " ----\n";
4529e201c85SYohann       // Get elem_size, eval_mode, num_comp
4532b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
4547d8d0e25Snbeams       // Basis action
4559e201c85SYohann       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
4569e201c85SYohann       switch (eval_mode) {
4577d8d0e25Snbeams         case CEED_EVAL_NONE:
4589e201c85SYohann           code << "      CeedScalar r_q_" << i << "[num_comp_in_" << i << "];\n";
4597d8d0e25Snbeams 
4609e201c85SYohann           bool is_strided;
4612b730f8bSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &Erestrict));
4622b730f8bSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionIsStrided(Erestrict, &is_strided));
4639e201c85SYohann           if (!is_strided) {
4642b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetLVectorSize(Erestrict, &lsize));
4657d8d0e25Snbeams             code << "      const CeedInt lsize_in_" << i << " = " << lsize << ";\n";
4669e201c85SYohann             CeedInt comp_stride;
4672b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetCompStride(Erestrict, &comp_stride));
4689e201c85SYohann             code << "      // CompStride: " << comp_stride << "\n";
4692b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetData(Erestrict, &restr_data));
4709e201c85SYohann             data->indices.inputs[i] = restr_data->d_ind;
4712b730f8bSJeremy L Thompson             code << "      readSliceQuadsOffset"
4722b730f8bSJeremy L Thompson                  << "3d<num_comp_in_" << i << ", " << comp_stride << ", Q_1d>(data, lsize_in_" << i << ", elem, q, indices.inputs[" << i << "], d_u_"
4732b730f8bSJeremy L Thompson                  << i << ", r_q_" << i << ");\n";
4747d8d0e25Snbeams           } else {
4752b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size));
4769e201c85SYohann             bool has_backend_strides;
4772b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionHasBackendStrides(Erestrict, &has_backend_strides));
4789e201c85SYohann             CeedInt num_elem;
4792b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetNumElements(Erestrict, &num_elem));
4809e201c85SYohann             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
4819e201c85SYohann             if (!has_backend_strides) {
4822b730f8bSJeremy L Thompson               CeedCallBackend(CeedElemRestrictionGetStrides(Erestrict, &strides));
4837d8d0e25Snbeams             }
4847d8d0e25Snbeams             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
4852b730f8bSJeremy L Thompson             code << "      readSliceQuadsStrided"
4862b730f8bSJeremy L Thompson                  << "3d<num_comp_in_" << i
4872b730f8bSJeremy L Thompson                  << ",Q_1d"
4882b730f8bSJeremy L Thompson                     ","
4892b730f8bSJeremy L Thompson                  << strides[0] << "," << strides[1] << "," << strides[2] << ">(data, elem, q, d_u_" << i << ", r_q_" << i << ");\n";
4907d8d0e25Snbeams           }
4917d8d0e25Snbeams           break;
4927d8d0e25Snbeams         case CEED_EVAL_INTERP:
4939e201c85SYohann           code << "      CeedScalar r_q_" << i << "[num_comp_in_" << i << "];\n";
4949e201c85SYohann           code << "      for (CeedInt j = 0; j < num_comp_in_" << i << " ; ++j) {\n";
4959e201c85SYohann           code << "        r_q_" << i << "[j] = r_t_" << i << "[q + j*Q_1d];\n";
4967d8d0e25Snbeams           code << "      }\n";
4977d8d0e25Snbeams           break;
4987d8d0e25Snbeams         case CEED_EVAL_GRAD:
4999e201c85SYohann           code << "      CeedScalar r_q_" << i << "[num_comp_in_" << i << "*dim];\n";
5009e201c85SYohann           code << "      gradCollo3d<num_comp_in_" << i << ",Q_1d>(data, q, r_t_" << i << ", s_G_in_" << i << ", r_q_" << i << ");\n";
5017d8d0e25Snbeams           break;
5027d8d0e25Snbeams         case CEED_EVAL_WEIGHT:
5039e201c85SYohann           code << "      CeedScalar r_q_" << i << "[1];\n";
5049e201c85SYohann           code << "      r_q_" << i << "[0] = r_t_" << i << "[q];\n";
5057d8d0e25Snbeams           break;  // No action
5067d8d0e25Snbeams         case CEED_EVAL_DIV:
5077d8d0e25Snbeams           break;  // TODO: Not implemented
5087d8d0e25Snbeams         case CEED_EVAL_CURL:
5097d8d0e25Snbeams           break;  // TODO: Not implemented
5107d8d0e25Snbeams       }
5117d8d0e25Snbeams     }
5127d8d0e25Snbeams     code << "\n      // -- Output fields --\n";
5139e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
5147d8d0e25Snbeams       code << "      // ---- Output field " << i << " ----\n";
5152b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
5167d8d0e25Snbeams       // Basis action
5179e201c85SYohann       switch (eval_mode) {
5187d8d0e25Snbeams         case CEED_EVAL_NONE:
5199e201c85SYohann           code << "      CeedScalar r_qq_" << i << "[num_comp_out_" << i << "];\n";
5207d8d0e25Snbeams           break;  // No action
5217d8d0e25Snbeams         case CEED_EVAL_INTERP:
5229e201c85SYohann           code << "      CeedScalar r_qq_" << i << "[num_comp_out_" << i << "];\n";
5237d8d0e25Snbeams           break;
5247d8d0e25Snbeams         case CEED_EVAL_GRAD:
5259e201c85SYohann           code << "      CeedScalar r_qq_" << i << "[num_comp_out_" << i << "*dim];\n";
5267d8d0e25Snbeams           break;
5277d8d0e25Snbeams         case CEED_EVAL_WEIGHT:
5287d8d0e25Snbeams           break;  // Should not occur
5297d8d0e25Snbeams         case CEED_EVAL_DIV:
5307d8d0e25Snbeams           break;  // TODO: Not implemented
5317d8d0e25Snbeams         case CEED_EVAL_CURL:
5327d8d0e25Snbeams           break;  // TODO: Not implemented
5337d8d0e25Snbeams       }
5347d8d0e25Snbeams     }
5357d8d0e25Snbeams   } else {
5369e201c85SYohann     code << "\n      // Note: Using full elements\n";
5377d8d0e25Snbeams     code << "      // -- Input fields --\n";
5389e201c85SYohann     for (CeedInt i = 0; i < num_input_fields; i++) {
5397d8d0e25Snbeams       code << "      // ---- Input field " << i << " ----\n";
5409e201c85SYohann       code << "      CeedScalar* r_q_" << i << " = r_t_" << i << ";\n";
5417d8d0e25Snbeams     }
5427d8d0e25Snbeams     code << "      // -- Output fields --\n";
5439e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
5447d8d0e25Snbeams       code << "      // ---- Output field " << i << " ----\n";
5459e201c85SYohann       code << "      CeedScalar* r_qq_" << i << " = r_tt_" << i << ";\n";
5467d8d0e25Snbeams     }
5477d8d0e25Snbeams   }
5487d8d0e25Snbeams   code << "\n      // -- QFunction Inputs and outputs --\n";
5499e201c85SYohann   code << "      CeedScalar* in[" << num_input_fields << "];\n";
5509e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
5517d8d0e25Snbeams     code << "      // ---- Input field " << i << " ----\n";
5529e201c85SYohann     code << "      in[" << i << "] = r_q_" << i << ";\n";
5537d8d0e25Snbeams   }
5549e201c85SYohann   code << "      CeedScalar* out[" << num_output_fields << "];\n";
5559e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
5567d8d0e25Snbeams     code << "      // ---- Output field " << i << " ----\n";
5579e201c85SYohann     code << "      out[" << i << "] = r_qq_" << i << ";\n";
5587d8d0e25Snbeams   }
5597d8d0e25Snbeams   code << "\n      // -- Apply QFunction --\n";
5609e201c85SYohann   code << "      " << q_function_name << "(ctx, ";
5619e201c85SYohann   if (dim != 3 || use_collograd_parallelization) {
5627d8d0e25Snbeams     code << "1";
5637d8d0e25Snbeams   } else {
5649e201c85SYohann     code << "Q_1d";
5657d8d0e25Snbeams   }
5667d8d0e25Snbeams   code << ", in, out);\n";
5679e201c85SYohann   if (use_collograd_parallelization) {
5687d8d0e25Snbeams     code << "      // -- Output fields --\n";
5699e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
5707d8d0e25Snbeams       code << "      // ---- Output field " << i << " ----\n";
5712b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
5727d8d0e25Snbeams       // Basis action
5739e201c85SYohann       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
5749e201c85SYohann       switch (eval_mode) {
5757d8d0e25Snbeams         case CEED_EVAL_NONE:
5769e201c85SYohann           code << "      for (CeedInt j = 0; j < num_comp_out_" << i << " ; ++j) {\n";
5779e201c85SYohann           code << "        r_tt_" << i << "[q + j*Q_1d] = r_qq_" << i << "[j];\n";
5787d8d0e25Snbeams           code << "      }\n";
5797d8d0e25Snbeams           break;  // No action
5807d8d0e25Snbeams         case CEED_EVAL_INTERP:
5819e201c85SYohann           code << "      for (CeedInt j = 0; j < num_comp_out_" << i << " ; ++j) {\n";
5829e201c85SYohann           code << "        r_tt_" << i << "[q + j*Q_1d] = r_qq_" << i << "[j];\n";
5837d8d0e25Snbeams           code << "      }\n";
5847d8d0e25Snbeams           break;
5857d8d0e25Snbeams         case CEED_EVAL_GRAD:
5869e201c85SYohann           code << "      gradColloTranspose3d<num_comp_out_" << i << ",Q_1d>(data, q, r_qq_" << i << ", s_G_out_" << i << ", r_tt_" << i << ");\n";
5877d8d0e25Snbeams           break;
5887d8d0e25Snbeams         case CEED_EVAL_WEIGHT:
5897d8d0e25Snbeams           break;  // Should not occur
5907d8d0e25Snbeams         case CEED_EVAL_DIV:
5917d8d0e25Snbeams           break;  // TODO: Not implemented
5927d8d0e25Snbeams         case CEED_EVAL_CURL:
5937d8d0e25Snbeams           break;  // TODO: Not implemented
5947d8d0e25Snbeams       }
5957d8d0e25Snbeams     }
5967d8d0e25Snbeams     code << "    }\n";
5977d8d0e25Snbeams   }
5987d8d0e25Snbeams 
5997d8d0e25Snbeams   // Output basis apply if needed
6007d8d0e25Snbeams   // Generate the correct eval mode code for each output
6017d8d0e25Snbeams   code << "\n    // -- Output field basis action and restrictions --\n";
6029e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
6037d8d0e25Snbeams     code << "    // ---- Output field " << i << " ----\n";
6049e201c85SYohann     // Get elem_size, eval_mode, num_comp
6052b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &Erestrict));
6062b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size));
6072b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
6082b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp));
6097d8d0e25Snbeams     // Basis action
6109e201c85SYohann     code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
6119e201c85SYohann     switch (eval_mode) {
6127d8d0e25Snbeams       case CEED_EVAL_NONE:
6139e201c85SYohann         code << "    CeedScalar* r_v_" << i << " = r_tt_" << i << ";\n";
6147d8d0e25Snbeams         break;  // No action
6157d8d0e25Snbeams       case CEED_EVAL_INTERP:
6169e201c85SYohann         code << "    CeedScalar r_v_" << i << "[num_comp_out_" << i << "*P_out_" << i << "];\n";
6172b730f8bSJeremy L Thompson         code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_out_" << i << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i
6182b730f8bSJeremy L Thompson              << ", s_B_out_" << i << ", r_v_" << i << ");\n";
6197d8d0e25Snbeams         break;
6207d8d0e25Snbeams       case CEED_EVAL_GRAD:
6219e201c85SYohann         code << "    CeedScalar r_v_" << i << "[num_comp_out_" << i << "*P_out_" << i << "];\n";
6229e201c85SYohann         if (use_collograd_parallelization) {
6232b730f8bSJeremy L Thompson           code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_out_" << i << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i
6242b730f8bSJeremy L Thompson                << ", s_B_out_" << i << ", r_v_" << i << ");\n";
6257d8d0e25Snbeams         } else {
6269e201c85SYohann           CeedInt P_1d;
6272b730f8bSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
6282b730f8bSJeremy L Thompson           CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
6292b730f8bSJeremy L Thompson           code << "    GradTranspose" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp_out_" << i
6302b730f8bSJeremy L Thompson                << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i << ", s_B_out_" << i << ", s_G_out_" << i << ", r_v_" << i << ");\n";
6317d8d0e25Snbeams         }
6327d8d0e25Snbeams         break;
6337d8d0e25Snbeams       // LCOV_EXCL_START
6347d8d0e25Snbeams       case CEED_EVAL_WEIGHT: {
6357d8d0e25Snbeams         Ceed ceed;
6362b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
6372b730f8bSJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
6387d8d0e25Snbeams         break;  // Should not occur
6397d8d0e25Snbeams       }
6407d8d0e25Snbeams       case CEED_EVAL_DIV:
6417d8d0e25Snbeams         break;  // TODO: Not implemented
6427d8d0e25Snbeams       case CEED_EVAL_CURL:
6437d8d0e25Snbeams         break;  // TODO: Not implemented
6447d8d0e25Snbeams                 // LCOV_EXCL_STOP
6457d8d0e25Snbeams     }
6467d8d0e25Snbeams     // Restriction
6479e201c85SYohann     bool is_strided;
6482b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionIsStrided(Erestrict, &is_strided));
6499e201c85SYohann     if (!is_strided) {
6502b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetLVectorSize(Erestrict, &lsize));
6517d8d0e25Snbeams       code << "    const CeedInt lsize_out_" << i << " = " << lsize << ";\n";
6529e201c85SYohann       CeedInt comp_stride;
6532b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetCompStride(Erestrict, &comp_stride));
6549e201c85SYohann       code << "    // CompStride: " << comp_stride << "\n";
6552b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetData(Erestrict, &restr_data));
6569e201c85SYohann       data->indices.outputs[i] = restr_data->d_ind;
6572b730f8bSJeremy L Thompson       code << "    writeDofsOffset" << dim << "d<num_comp_out_" << i << ", " << comp_stride << ", P_out_" << i << ">(data, lsize_out_" << i
6582b730f8bSJeremy L Thompson            << ", elem, indices.outputs[" << i << "], r_v_" << i << ", d_v_" << i << ");\n";
6597d8d0e25Snbeams     } else {
6609e201c85SYohann       bool has_backend_strides;
6612b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionHasBackendStrides(Erestrict, &has_backend_strides));
6629e201c85SYohann       CeedInt num_elem;
6632b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetNumElements(Erestrict, &num_elem));
6649e201c85SYohann       CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
6659e201c85SYohann       if (!has_backend_strides) {
6662b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetStrides(Erestrict, &strides));
6677d8d0e25Snbeams       }
6687d8d0e25Snbeams       code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
6692b730f8bSJeremy L Thompson       code << "    writeDofsStrided" << dim << "d<num_comp_out_" << i << ",P_out_" << i << "," << strides[0] << "," << strides[1] << "," << strides[2]
6702b730f8bSJeremy L Thompson            << ">(data, elem, r_v_" << i << ", d_v_" << i << ");\n";
6717d8d0e25Snbeams     }
6727d8d0e25Snbeams   }
6737d8d0e25Snbeams 
6747d8d0e25Snbeams   code << "  }\n";
6757d8d0e25Snbeams   code << "}\n";
6767d8d0e25Snbeams   code << "// -----------------------------------------------------------------------------\n\n";
6777d8d0e25Snbeams 
6787d8d0e25Snbeams   // View kernel for debugging
67946dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "Generated Operator Kernels:\n");
6803f21f6b1SJeremy L Thompson   CeedDebug(ceed, code.str().c_str());
6817d8d0e25Snbeams 
682539ec17dSJeremy L Thompson   CeedInt block_sizes[3] = {0, 0, 0};
6839e201c85SYohann   CeedInt num_elem;
6842b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
6852b730f8bSJeremy L Thompson   CeedCallBackend(BlockGridCalculate_Hip_gen(dim, num_elem, data->max_P_1d, Q_1d, block_sizes));
686*eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Hip(ceed, code.str().c_str(), &data->module, 2, "T_1D", block_sizes[0], "BLOCK_SIZE",
6872b730f8bSJeremy L Thompson                                   block_sizes[0] * block_sizes[1] * block_sizes[2]));
688*eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, operator_name.c_str(), &data->op));
6897d8d0e25Snbeams 
6902b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
691e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6927d8d0e25Snbeams }
6932a86cc9dSSebastian Grimberg 
6947d8d0e25Snbeams //------------------------------------------------------------------------------
695