xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-gen/ceed-cuda-gen-operator-build.cpp (revision 6e536b992ff6bc401c55631f1bc4464446496b52)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3241a4b83SYohann //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5241a4b83SYohann //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
73d576824SJeremy L Thompson 
8b8e71988SJeremy L Thompson #define CEED_DEBUG_COLOR 12
97df94212SJeremy L Thompson 
1049aac155SJeremy L Thompson #include <ceed.h>
11ec3da8bcSJed Brown #include <ceed/backend.h>
129e201c85SYohann #include <ceed/jit-tools.h>
133d576824SJeremy L Thompson #include <cuda_runtime.h>
142b730f8bSJeremy L Thompson 
15241a4b83SYohann #include <iostream>
16241a4b83SYohann #include <sstream>
17b2165e7aSSebastian Grimberg #include <string>
182b730f8bSJeremy L Thompson 
190d0321e0SJeremy L Thompson #include "../cuda-ref/ceed-cuda-ref.h"
20241a4b83SYohann #include "../cuda-shared/ceed-cuda-shared.h"
2149aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
222b730f8bSJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
232b730f8bSJeremy L Thompson #include "ceed-cuda-gen.h"
24241a4b83SYohann 
25ab213215SJeremy L Thompson //------------------------------------------------------------------------------
26b2165e7aSSebastian Grimberg // Build single operator kernel
27ab213215SJeremy L Thompson //------------------------------------------------------------------------------
28eb7e6cafSJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Cuda_gen(CeedOperator op) {
29241a4b83SYohann   using std::ostringstream;
30241a4b83SYohann   using std::string;
31ca735530SJeremy L Thompson 
32ca735530SJeremy L Thompson   bool                      is_setup_done, is_identity_qf;
33ca735530SJeremy L Thompson   struct cudaDeviceProp     prop;
34241a4b83SYohann   Ceed                      ceed;
35ca735530SJeremy L Thompson   Ceed_Cuda                *ceed_data;
36ca735530SJeremy L Thompson   CeedSize                  l_size;
372b730f8bSJeremy L Thompson   CeedInt                   Q, P_1d = 0, Q_1d = 0, elem_size, num_input_fields, num_output_fields, num_comp, dim = 1;
389e201c85SYohann   CeedEvalMode              eval_mode;
39edb2538eSJeremy L Thompson   CeedElemRestriction       elem_rstr;
40edb2538eSJeremy L Thompson   CeedElemRestriction_Cuda *rstr_data;
41241a4b83SYohann   CeedBasis                 basis;
42241a4b83SYohann   CeedBasis_Cuda_shared    *basis_data;
43ca735530SJeremy L Thompson   CeedQFunctionField       *qf_input_fields, *qf_output_fields;
44ca735530SJeremy L Thompson   CeedQFunction_Cuda_gen   *qf_data;
45ca735530SJeremy L Thompson   CeedQFunction             qf;
46ca735530SJeremy L Thompson   CeedOperatorField        *op_input_fields, *op_output_fields;
47ca735530SJeremy L Thompson   CeedOperator_Cuda_gen    *data;
48ca735530SJeremy L Thompson 
49ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
50ca735530SJeremy L Thompson   if (is_setup_done) return CEED_ERROR_SUCCESS;
51ca735530SJeremy L Thompson 
52ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
53ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
54ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
55ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
56ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
57ca735530SJeremy L Thompson   Q_1d = Q;
58ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
59ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
60241a4b83SYohann 
619e201c85SYohann   // TODO: put in a function?
620b454692Sjeremylt   // Check for restriction only identity operator
632b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
640b454692Sjeremylt   if (is_identity_qf) {
659e201c85SYohann     CeedEvalMode eval_mode_in, eval_mode_out;
66ca735530SJeremy L Thompson 
672b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
682b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
696574a04fSJeremy L Thompson     CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
706574a04fSJeremy L Thompson               "Backend does not implement restriction only identity operators");
710b454692Sjeremylt   }
720b454692Sjeremylt 
73241a4b83SYohann   ostringstream code;
74241a4b83SYohann 
759e201c85SYohann   // TODO: put in a function?
76f1a13f77SYohann Dudouit   // Add atomicAdd function for old NVidia architectures
772b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_data));
782b730f8bSJeremy L Thompson   CeedCallBackend(cudaGetDeviceProperties(&prop, ceed_data->device_id));
7980a9ef05SNatalie Beams   if ((prop.major < 6) && (CEED_SCALAR_TYPE != CEED_SCALAR_FP32)) {
8022070f95SJeremy L Thompson     char       *atomic_add_source;
8122070f95SJeremy L Thompson     const char *atomic_add_path;
82ca735530SJeremy L Thompson 
832b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-atomic-add-fallback.h", &atomic_add_path));
8423d4529eSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Atomic Add Source -----\n");
852b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, atomic_add_path, &atomic_add_source));
869e201c85SYohann     code << atomic_add_source;
872b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&atomic_add_path));
882b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&atomic_add_source));
89f1a13f77SYohann Dudouit   }
90f1a13f77SYohann Dudouit 
919e201c85SYohann   // Load basis source files
929e201c85SYohann   // TODO: generalize to accept different device functions?
939e201c85SYohann   {
9422070f95SJeremy L Thompson     char       *tensor_basis_kernel_source;
9522070f95SJeremy L Thompson     const char *tensor_basis_kernel_path;
96ca735530SJeremy L Thompson 
972b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-shared-basis-tensor-templates.h", &tensor_basis_kernel_path));
9823d4529eSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Tensor Basis Kernel Source -----\n");
992b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, tensor_basis_kernel_path, &tensor_basis_kernel_source));
1009e201c85SYohann     code << tensor_basis_kernel_source;
1012b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&tensor_basis_kernel_path));
1022b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&tensor_basis_kernel_source));
1039e201c85SYohann   }
1049e201c85SYohann   {
10522070f95SJeremy L Thompson     char       *cuda_gen_template_source;
10622070f95SJeremy L Thompson     const char *cuda_gen_template_path;
107ca735530SJeremy L Thompson 
1082b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-gen-templates.h", &cuda_gen_template_path));
10923d4529eSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Cuda-Gen Template Source -----\n");
1102b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, cuda_gen_template_path, &cuda_gen_template_source));
1119e201c85SYohann     code << cuda_gen_template_source;
1122b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&cuda_gen_template_path));
1132b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&cuda_gen_template_source));
1149e201c85SYohann   }
115241a4b83SYohann 
1169e201c85SYohann   // Get QFunction source and name
11709095acaSJeremy L Thompson   string qfunction_source(qf_data->qfunction_source);
11809095acaSJeremy L Thompson   string qfunction_name(qf_data->qfunction_name);
1199e201c85SYohann   string operator_name;
12009095acaSJeremy L Thompson   operator_name = "CeedKernelCudaGenOperator_" + qfunction_name;
1214d537eeaSYohann 
1229e201c85SYohann   // Find dim, P_1d, Q_1d
1239e201c85SYohann   data->max_P_1d = 0;
1249e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1252b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
126356036faSJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
127ca735530SJeremy L Thompson       bool is_tensor;
128ca735530SJeremy L Thompson 
1292b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1302b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1311da99368SJeremy L Thompson 
1329e201c85SYohann       // Collect dim, P_1d, and Q_1d
1332b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &dim));
134ca735530SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
135ca735530SJeremy L Thompson       CeedCheck(is_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
1362b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
1372b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
138ca735530SJeremy L Thompson       data->max_P_1d = CeedIntMax(data->max_P_1d, P_1d);
1391da99368SJeremy L Thompson     }
1401da99368SJeremy L Thompson   }
1419e201c85SYohann   // Check output bases for Q_1d, dim as well
142356036faSJeremy L Thompson   //   The only input basis might be CEED_BASIS_NONE
1439e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1442b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
145356036faSJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
146ca735530SJeremy L Thompson       bool is_tensor;
147ca735530SJeremy L Thompson 
1482b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1492b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1500f54b25eSjeremylt 
1519e201c85SYohann       // Collect Q_1d
1522b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &dim));
153ca735530SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
154ca735530SJeremy L Thompson       CeedCheck(is_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
1552b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
1561da99368SJeremy L Thompson     }
1571da99368SJeremy L Thompson   }
1581da99368SJeremy L Thompson   data->dim  = dim;
1599e201c85SYohann   data->Q_1d = Q_1d;
1601da99368SJeremy L Thompson 
1619e201c85SYohann   // Only use 3D collocated gradient parallelization strategy when gradient is computed
1629e201c85SYohann   // TODO: put in a function?
1639e201c85SYohann   bool use_collograd_parallelization = false;
164ca735530SJeremy L Thompson 
1659e201c85SYohann   if (dim == 3) {
1669e201c85SYohann     bool was_grad_found = false;
167ca735530SJeremy L Thompson 
1689e201c85SYohann     for (CeedInt i = 0; i < num_input_fields; i++) {
1692b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1709e201c85SYohann       if (eval_mode == CEED_EVAL_GRAD) {
1712b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1722b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1731c66c397SJeremy L Thompson         use_collograd_parallelization = basis_data->d_collo_grad_1d && (was_grad_found ? use_collograd_parallelization : true);
1749e201c85SYohann         was_grad_found                = true;
1759e201c85SYohann       }
1769e201c85SYohann     }
1779e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
1782b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1799e201c85SYohann       if (eval_mode == CEED_EVAL_GRAD) {
1802b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1812b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1821c66c397SJeremy L Thompson         use_collograd_parallelization = basis_data->d_collo_grad_1d && (was_grad_found ? use_collograd_parallelization : true);
1839e201c85SYohann         was_grad_found                = true;
1849e201c85SYohann       }
1859e201c85SYohann     }
1861da99368SJeremy L Thompson   }
1871da99368SJeremy L Thompson 
1889e201c85SYohann   // Define CEED_Q_VLA
1899e201c85SYohann   code << "\n#undef CEED_Q_VLA\n";
1909e201c85SYohann   if (dim != 3 || use_collograd_parallelization) {
1919e201c85SYohann     code << "#define CEED_Q_VLA 1\n\n";
1929e201c85SYohann   } else {
1939e201c85SYohann     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
1949e201c85SYohann   }
1959e201c85SYohann 
19609095acaSJeremy L Thompson   code << qfunction_source;
197241a4b83SYohann 
198241a4b83SYohann   // Setup
199818e0025SJeremy L Thompson   code << "\n// -----------------------------------------------------------------------------\n";
2002b730f8bSJeremy L Thompson   code << "\nextern \"C\" __global__ void " << operator_name
2012b730f8bSJeremy L Thompson        << "(CeedInt num_elem, void* ctx, FieldsInt_Cuda indices, Fields_Cuda fields, Fields_Cuda B, Fields_Cuda G, CeedScalar* W) {\n";
2029e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
2032b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
2049e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
2059e201c85SYohann       code << "  const CeedScalar* d_u_" << i << " = fields.inputs[" << i << "];\n";
206241a4b83SYohann     }
207241a4b83SYohann   }
208241a4b83SYohann 
2099e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
2109e201c85SYohann     code << "  CeedScalar* d_v_" << i << " = fields.outputs[" << i << "];\n";
211241a4b83SYohann   }
2121da99368SJeremy L Thompson 
2139e201c85SYohann   code << "  const CeedInt dim = " << dim << ";\n";
2149e201c85SYohann   code << "  const CeedInt Q_1d = " << Q_1d << ";\n";
2151da99368SJeremy L Thompson 
216241a4b83SYohann   code << "  extern __shared__ CeedScalar slice[];\n";
2179e201c85SYohann   // TODO put in a function? InitSharedData_Cuda?
2189e201c85SYohann   code << "  SharedData_Cuda data;\n";
2199e201c85SYohann   code << "  data.t_id_x = threadIdx.x;\n";
2209e201c85SYohann   code << "  data.t_id_y = threadIdx.y;\n";
2219e201c85SYohann   code << "  data.t_id_z = threadIdx.z;\n";
2229e201c85SYohann   code << "  data.t_id  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
2239e201c85SYohann   code << "  data.slice = slice+data.t_id_z*T_1D" << (dim > 1 ? "*T_1D" : "") << ";\n";
224920dcdc4Sjeremylt 
225818e0025SJeremy L Thompson   code << "\n  // -- Input field constants and basis data --\n";
2269e201c85SYohann   // TODO: Put in a function?
227ac421f39SYohann   // Initialize constants, and matrices B and G
2289e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
229818e0025SJeremy L Thompson     code << "  // ---- Input field " << i << " ----\n";
2309e201c85SYohann     // Get elem_size, eval_mode, num_comp
231edb2538eSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
232edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
2332b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
234edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
235920dcdc4Sjeremylt 
236920dcdc4Sjeremylt     // Set field constants
2379e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {
2382b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
239356036faSJeremy L Thompson       if (basis != CEED_BASIS_NONE) {
2402b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
2419e201c85SYohann         code << "  const CeedInt P_in_" << i << " = " << P_1d << ";\n";
242920dcdc4Sjeremylt       } else {
2439e201c85SYohann         code << "  const CeedInt P_in_" << i << " = " << Q_1d << ";\n";
244920dcdc4Sjeremylt       }
2459e201c85SYohann       code << "  const CeedInt num_comp_in_" << i << " = " << num_comp << ";\n";
246920dcdc4Sjeremylt     }
247920dcdc4Sjeremylt 
248920dcdc4Sjeremylt     // Load basis data
2499e201c85SYohann     code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
2509e201c85SYohann     switch (eval_mode) {
251241a4b83SYohann       case CEED_EVAL_NONE:
252241a4b83SYohann         break;
253241a4b83SYohann       case CEED_EVAL_INTERP:
2542b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
2559e201c85SYohann         data->B.inputs[i] = basis_data->d_interp_1d;
2569e201c85SYohann         code << "  __shared__ CeedScalar s_B_in_" << i << "[" << P_1d * Q_1d << "];\n";
2579e201c85SYohann         code << "  loadMatrix<P_in_" << i << ",Q_1d>(data, B.inputs[" << i << "], s_B_in_" << i << ");\n";
258241a4b83SYohann         break;
259241a4b83SYohann       case CEED_EVAL_GRAD:
2602b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
2619e201c85SYohann         data->B.inputs[i] = basis_data->d_interp_1d;
2629e201c85SYohann         code << "  __shared__ CeedScalar s_B_in_" << i << "[" << P_1d * Q_1d << "];\n";
2639e201c85SYohann         code << "  loadMatrix<P_in_" << i << ",Q_1d>(data, B.inputs[" << i << "], s_B_in_" << i << ");\n";
2649e201c85SYohann         if (use_collograd_parallelization) {
2659e201c85SYohann           data->G.inputs[i] = basis_data->d_collo_grad_1d;
2669e201c85SYohann           code << "  __shared__ CeedScalar s_G_in_" << i << "[" << Q_1d * Q_1d << "];\n";
2679e201c85SYohann           code << "  loadMatrix<Q_1d,Q_1d>(data, G.inputs[" << i << "], s_G_in_" << i << ");\n";
268ac421f39SYohann         } else {
2691c66c397SJeremy L Thompson           bool has_collo_grad = basis_data->d_collo_grad_1d;
2709e201c85SYohann           data->G.inputs[i]   = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2719e201c85SYohann           code << "  __shared__ CeedScalar s_G_in_" << i << "[" << Q_1d * (has_collo_grad ? Q_1d : P_1d) << "];\n";
2722b730f8bSJeremy 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
2732b730f8bSJeremy L Thompson                << ");\n";
274ac421f39SYohann         }
275241a4b83SYohann         break;
276241a4b83SYohann       case CEED_EVAL_WEIGHT:
277241a4b83SYohann         break;  // No action
278241a4b83SYohann       case CEED_EVAL_DIV:
279241a4b83SYohann         break;  // TODO: Not implemented
280241a4b83SYohann       case CEED_EVAL_CURL:
281241a4b83SYohann         break;  // TODO: Not implemented
282241a4b83SYohann     }
283241a4b83SYohann   }
284920dcdc4Sjeremylt 
285818e0025SJeremy L Thompson   code << "\n  // -- Output field constants and basis data --\n";
2869e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
287818e0025SJeremy L Thompson     code << "  // ---- Output field " << i << " ----\n";
2889e201c85SYohann     // Get elem_size, eval_mode, num_comp
289edb2538eSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
290edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
2912b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
292edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
293920dcdc4Sjeremylt 
294920dcdc4Sjeremylt     // Set field constants
2952b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
296356036faSJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
2972b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
2989e201c85SYohann       code << "  const CeedInt P_out_" << i << " = " << P_1d << ";\n";
299920dcdc4Sjeremylt     } else {
3009e201c85SYohann       code << "  const CeedInt P_out_" << i << " = " << Q_1d << ";\n";
301920dcdc4Sjeremylt     }
3029e201c85SYohann     code << "  const CeedInt num_comp_out_" << i << " = " << num_comp << ";\n";
303920dcdc4Sjeremylt 
304920dcdc4Sjeremylt     // Load basis data
3059e201c85SYohann     code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
3069e201c85SYohann     switch (eval_mode) {
307920dcdc4Sjeremylt       case CEED_EVAL_NONE:
308920dcdc4Sjeremylt         break;  // No action
309920dcdc4Sjeremylt       case CEED_EVAL_INTERP:
3102b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
3119e201c85SYohann         data->B.outputs[i] = basis_data->d_interp_1d;
3129e201c85SYohann         code << "  __shared__ CeedScalar s_B_out_" << i << "[" << P_1d * Q_1d << "];\n";
3139e201c85SYohann         code << "  loadMatrix<P_out_" << i << ",Q_1d>(data, B.outputs[" << i << "], s_B_out_" << i << ");\n";
314241a4b83SYohann         break;
315241a4b83SYohann       case CEED_EVAL_GRAD:
3162b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
3179e201c85SYohann         data->B.outputs[i] = basis_data->d_interp_1d;
3189e201c85SYohann         code << "  __shared__ CeedScalar s_B_out_" << i << "[" << P_1d * Q_1d << "];\n";
3199e201c85SYohann         code << "  loadMatrix<P_out_" << i << ",Q_1d>(data, B.outputs[" << i << "], s_B_out_" << i << ");\n";
3209e201c85SYohann         if (use_collograd_parallelization) {
3219e201c85SYohann           data->G.outputs[i] = basis_data->d_collo_grad_1d;
3229e201c85SYohann           code << "  __shared__ CeedScalar s_G_out_" << i << "[" << Q_1d * Q_1d << "];\n";
3239e201c85SYohann           code << "  loadMatrix<Q_1d,Q_1d>(data, G.outputs[" << i << "], s_G_out_" << i << ");\n";
324ac421f39SYohann         } else {
3251c66c397SJeremy L Thompson           bool has_collo_grad = basis_data->d_collo_grad_1d;
3269e201c85SYohann           data->G.outputs[i]  = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
3279e201c85SYohann           code << "  __shared__ CeedScalar s_G_out_" << i << "[" << Q_1d * (has_collo_grad ? Q_1d : P_1d) << "];\n";
3282b730f8bSJeremy L Thompson           code << "  loadMatrix<" << (has_collo_grad ? "Q_1d" : ("P_out_" + std::to_string(i))) << ",Q_1d>(data, G.outputs[" << i << "], s_G_out_"
3292b730f8bSJeremy L Thompson                << i << ");\n";
330ac421f39SYohann         }
331241a4b83SYohann         break;
332e9f4dca0SJeremy L Thompson       // LCOV_EXCL_START
333241a4b83SYohann       case CEED_EVAL_WEIGHT: {
334*6e536b99SJeremy L Thompson         return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
335241a4b83SYohann         break;  // Should not occur
336241a4b83SYohann       }
337241a4b83SYohann       case CEED_EVAL_DIV:
338bcbe1c99SJeremy L Thompson       case CEED_EVAL_CURL: {
339*6e536b99SJeremy L Thompson         return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
340bcbe1c99SJeremy L Thompson         break;  // Should not occur
341bcbe1c99SJeremy L Thompson       }
342e9f4dca0SJeremy L Thompson         // LCOV_EXCL_STOP
343241a4b83SYohann     }
344241a4b83SYohann   }
345818e0025SJeremy L Thompson   code << "\n  // -- Element loop --\n";
346ac421f39SYohann   code << "  __syncthreads();\n";
3479e201c85SYohann   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
348241a4b83SYohann   // Input basis apply if needed
349ac421f39SYohann   // Generate the correct eval mode code for each input
350818e0025SJeremy L Thompson   code << "    // -- Input field restrictions and basis actions --\n";
3519e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
352818e0025SJeremy L Thompson     code << "    // ---- Input field " << i << " ----\n";
3539e201c85SYohann     // Get elem_size, eval_mode, num_comp
354edb2538eSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
355edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
3562b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
357edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
358920dcdc4Sjeremylt 
3599e201c85SYohann     // TODO: put in a function?
360920dcdc4Sjeremylt     // Restriction
3612b730f8bSJeremy L Thompson     if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_collograd_parallelization)) {
3629e201c85SYohann       code << "    CeedScalar r_u_" << i << "[num_comp_in_" << i << "*P_in_" << i << "];\n";
3639a2291e3SJeremy L Thompson 
3649e201c85SYohann       bool is_strided;
365ca735530SJeremy L Thompson 
366edb2538eSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
3679e201c85SYohann       if (!is_strided) {
3689e201c85SYohann         CeedInt comp_stride;
369ca735530SJeremy L Thompson 
370edb2538eSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
371ca735530SJeremy L Thompson         code << "    const CeedInt l_size_in_" << i << " = " << l_size << ";\n";
372edb2538eSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3739e201c85SYohann         code << "    // CompStride: " << comp_stride << "\n";
374edb2538eSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
375edb2538eSJeremy L Thompson         data->indices.inputs[i] = rstr_data->d_ind;
376ca735530SJeremy L Thompson         code << "    readDofsOffset" << dim << "d<num_comp_in_" << i << ", " << comp_stride << ", P_in_" << i << ">(data, l_size_in_" << i
3772b730f8bSJeremy L Thompson              << ", elem, indices.inputs[" << i << "], d_u_" << i << ", r_u_" << i << ");\n";
378ccedf6b0Sjeremylt       } else {
379b2165e7aSSebastian Grimberg         bool    has_backend_strides;
3809e201c85SYohann         CeedInt num_elem;
381ca735530SJeremy L Thompson 
382edb2538eSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
383edb2538eSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
3849e201c85SYohann         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
385ca735530SJeremy L Thompson 
386b2165e7aSSebastian Grimberg         if (!has_backend_strides) {
38756c48462SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
388ccedf6b0Sjeremylt         }
389920dcdc4Sjeremylt         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
3902b730f8bSJeremy L Thompson         code << "    readDofsStrided" << dim << "d<num_comp_in_" << i << ",P_in_" << i << "," << strides[0] << "," << strides[1] << "," << strides[2]
3912b730f8bSJeremy L Thompson              << ">(data, elem, d_u_" << i << ", r_u_" << i << ");\n";
392920dcdc4Sjeremylt       }
393920dcdc4Sjeremylt     }
394920dcdc4Sjeremylt 
3959e201c85SYohann     // TODO: put in a function?
396920dcdc4Sjeremylt     // Basis action
3979e201c85SYohann     code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
3989e201c85SYohann     switch (eval_mode) {
399920dcdc4Sjeremylt       case CEED_EVAL_NONE:
4009e201c85SYohann         if (!use_collograd_parallelization) {
4019e201c85SYohann           code << "    CeedScalar* r_t_" << i << " = r_u_" << i << ";\n";
402920dcdc4Sjeremylt         }
403920dcdc4Sjeremylt         break;
404920dcdc4Sjeremylt       case CEED_EVAL_INTERP:
4059e201c85SYohann         code << "    CeedScalar r_t_" << i << "[num_comp_in_" << i << "*Q_1d];\n";
4062b730f8bSJeremy 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_"
4072b730f8bSJeremy L Thompson              << i << ", r_t_" << i << ");\n";
408241a4b83SYohann         break;
409241a4b83SYohann       case CEED_EVAL_GRAD:
4109e201c85SYohann         if (use_collograd_parallelization) {
4119e201c85SYohann           code << "    CeedScalar r_t_" << i << "[num_comp_in_" << i << "*Q_1d];\n";
4122b730f8bSJeremy L Thompson           code << "    Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_in_" << i << ",P_in_" << i << ",Q_1d>(data, r_u_" << i
4132b730f8bSJeremy L Thompson                << ", s_B_in_" << i << ", r_t_" << i << ");\n";
414ac421f39SYohann         } else {
4159e201c85SYohann           CeedInt P_1d;
416ca735530SJeremy L Thompson 
4172b730f8bSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
4182b730f8bSJeremy L Thompson           CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
4199e201c85SYohann           code << "    CeedScalar r_t_" << i << "[num_comp_in_" << i << "*dim*Q_1d];\n";
4202b730f8bSJeremy L Thompson           code << "    Grad" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp_in_" << i
4212b730f8bSJeremy L Thompson                << ",P_in_" << i << ",Q_1d>(data, r_u_" << i << ", s_B_in_" << i << ", s_G_in_" << i << ", r_t_" << i << ");\n";
422ac421f39SYohann         }
423241a4b83SYohann         break;
424241a4b83SYohann       case CEED_EVAL_WEIGHT:
4259e201c85SYohann         code << "    CeedScalar r_t_" << i << "[Q_1d];\n";
4262b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
4272b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
428437930d1SJeremy L Thompson         data->W = basis_data->d_q_weight_1d;
4299e201c85SYohann         code << "    Weight" << (dim > 1 ? "Tensor" : "") << dim << "d<Q_1d>(data, W, r_t_" << i << ");\n";
430241a4b83SYohann         break;  // No action
431241a4b83SYohann       case CEED_EVAL_DIV:
432241a4b83SYohann         break;  // TODO: Not implemented
433241a4b83SYohann       case CEED_EVAL_CURL:
434241a4b83SYohann         break;  // TODO: Not implemented
435241a4b83SYohann     }
436241a4b83SYohann   }
437ac421f39SYohann 
438ea61e9acSJeremy L Thompson   // TODO: put in a function + separate collograd logic
439241a4b83SYohann   // Q function
440818e0025SJeremy L Thompson   code << "\n    // -- Output field setup --\n";
4419e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
442818e0025SJeremy L Thompson     code << "\n    // ---- Output field " << i << " ----\n";
4432b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
4442b730f8bSJeremy L Thompson     if (eval_mode == CEED_EVAL_GRAD) {
4459e201c85SYohann       if (use_collograd_parallelization) {
446ac421f39SYohann         // Accumulator for gradient slices
4479e201c85SYohann         code << "    CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*Q_1d];\n";
4489e201c85SYohann         code << "    for (CeedInt i = 0; i < num_comp_out_" << i << "; i++) {\n";
4499e201c85SYohann         code << "      for (CeedInt j = 0; j < Q_1d; ++j) {\n";
4509e201c85SYohann         code << "        r_tt_" << i << "[j + i*Q_1d] = 0.0;\n";
451ac421f39SYohann         code << "      }\n";
452ac421f39SYohann         code << "    }\n";
453ac421f39SYohann       } else {
4549e201c85SYohann         code << "    CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*dim*Q_1d];\n";
455241a4b83SYohann       }
456ac421f39SYohann     }
4572b730f8bSJeremy L Thompson     if (eval_mode == CEED_EVAL_NONE || eval_mode == CEED_EVAL_INTERP) {
4589e201c85SYohann       code << "    CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*Q_1d];\n";
459241a4b83SYohann     }
460241a4b83SYohann   }
461ac421f39SYohann   // We treat quadrature points per slice in 3d to save registers
4629e201c85SYohann   if (use_collograd_parallelization) {
4639e201c85SYohann     code << "\n    // Note: Using planes of 3D elements\n";
464ac421f39SYohann     code << "#pragma unroll\n";
4659e201c85SYohann     code << "    for (CeedInt q = 0; q < Q_1d; q++) {\n";
466818e0025SJeremy L Thompson     code << "      // -- Input fields --\n";
4679e201c85SYohann     for (CeedInt i = 0; i < num_input_fields; i++) {
468818e0025SJeremy L Thompson       code << "      // ---- Input field " << i << " ----\n";
4699e201c85SYohann       // Get elem_size, eval_mode, num_comp
4702b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
471ac421f39SYohann       // Basis action
4729e201c85SYohann       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
4739e201c85SYohann       switch (eval_mode) {
474ac421f39SYohann         case CEED_EVAL_NONE:
475ca735530SJeremy L Thompson           bool is_strided;
476ca735530SJeremy L Thompson 
4779e201c85SYohann           code << "      CeedScalar r_q_" << i << "[num_comp_in_" << i << "];\n";
4789a2291e3SJeremy L Thompson 
479edb2538eSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
480edb2538eSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
4819e201c85SYohann           if (!is_strided) {
4829e201c85SYohann             CeedInt comp_stride;
483ca735530SJeremy L Thompson 
484edb2538eSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
485ca735530SJeremy L Thompson             code << "      const CeedInt l_size_in_" << i << " = " << l_size << ";\n";
486edb2538eSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
4879e201c85SYohann             code << "      // CompStride: " << comp_stride << "\n";
488edb2538eSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
489edb2538eSJeremy L Thompson             data->indices.inputs[i] = rstr_data->d_ind;
4902b730f8bSJeremy L Thompson             code << "      readSliceQuadsOffset"
491ca735530SJeremy L Thompson                  << "3d<num_comp_in_" << i << ", " << comp_stride << ", Q_1d>(data, l_size_in_" << i << ", elem, q, indices.inputs[" << i << "], d_u_"
4922b730f8bSJeremy L Thompson                  << i << ", r_q_" << i << ");\n";
493920dcdc4Sjeremylt           } else {
494b2165e7aSSebastian Grimberg             bool    has_backend_strides;
4959e201c85SYohann             CeedInt num_elem;
496ca735530SJeremy L Thompson 
497edb2538eSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
498edb2538eSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
499edb2538eSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
5009e201c85SYohann             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
501ca735530SJeremy L Thompson 
502b2165e7aSSebastian Grimberg             if (!has_backend_strides) {
50356c48462SJeremy L Thompson               CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
504920dcdc4Sjeremylt             }
505920dcdc4Sjeremylt             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
5062b730f8bSJeremy L Thompson             code << "      readSliceQuadsStrided"
5072b730f8bSJeremy L Thompson                  << "3d<num_comp_in_" << i
5082b730f8bSJeremy L Thompson                  << ",Q_1d"
5092b730f8bSJeremy L Thompson                     ","
5102b730f8bSJeremy L Thompson                  << strides[0] << "," << strides[1] << "," << strides[2] << ">(data, elem, q, d_u_" << i << ", r_q_" << i << ");\n";
511920dcdc4Sjeremylt           }
512ac421f39SYohann           break;
513ac421f39SYohann         case CEED_EVAL_INTERP:
5149e201c85SYohann           code << "      CeedScalar r_q_" << i << "[num_comp_in_" << i << "];\n";
5159e201c85SYohann           code << "      for (CeedInt j = 0; j < num_comp_in_" << i << " ; ++j) {\n";
5169e201c85SYohann           code << "        r_q_" << i << "[j] = r_t_" << i << "[q + j*Q_1d];\n";
517ac421f39SYohann           code << "      }\n";
518ac421f39SYohann           break;
519ac421f39SYohann         case CEED_EVAL_GRAD:
5209e201c85SYohann           code << "      CeedScalar r_q_" << i << "[num_comp_in_" << i << "*dim];\n";
5219e201c85SYohann           code << "      gradCollo3d<num_comp_in_" << i << ",Q_1d>(data, q, r_t_" << i << ", s_G_in_" << i << ", r_q_" << i << ");\n";
522ac421f39SYohann           break;
523ac421f39SYohann         case CEED_EVAL_WEIGHT:
5249e201c85SYohann           code << "      CeedScalar r_q_" << i << "[1];\n";
5259e201c85SYohann           code << "      r_q_" << i << "[0] = r_t_" << i << "[q];\n";
526ac421f39SYohann           break;  // No action
527ac421f39SYohann         case CEED_EVAL_DIV:
528ac421f39SYohann           break;  // TODO: Not implemented
529ac421f39SYohann         case CEED_EVAL_CURL:
530ac421f39SYohann           break;  // TODO: Not implemented
531ac421f39SYohann       }
532ac421f39SYohann     }
533818e0025SJeremy L Thompson     code << "\n      // -- Output fields --\n";
5349e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
535818e0025SJeremy L Thompson       code << "      // ---- Output field " << i << " ----\n";
5362b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
537ac421f39SYohann       // Basis action
5389e201c85SYohann       switch (eval_mode) {
539ac421f39SYohann         case CEED_EVAL_NONE:
5409e201c85SYohann           code << "      CeedScalar r_qq_" << i << "[num_comp_out_" << i << "];\n";
541ac421f39SYohann           break;  // No action
542ac421f39SYohann         case CEED_EVAL_INTERP:
5439e201c85SYohann           code << "      CeedScalar r_qq_" << i << "[num_comp_out_" << i << "];\n";
544ac421f39SYohann           break;
545ac421f39SYohann         case CEED_EVAL_GRAD:
5469e201c85SYohann           code << "      CeedScalar r_qq_" << i << "[num_comp_out_" << i << "*dim];\n";
547ac421f39SYohann           break;
548ac421f39SYohann         case CEED_EVAL_WEIGHT:
549ac421f39SYohann           break;  // Should not occur
550ac421f39SYohann         case CEED_EVAL_DIV:
551ac421f39SYohann           break;  // TODO: Not implemented
552ac421f39SYohann         case CEED_EVAL_CURL:
553ac421f39SYohann           break;  // TODO: Not implemented
554ac421f39SYohann       }
555ac421f39SYohann     }
556ac421f39SYohann   } else {
5579e201c85SYohann     code << "\n      // Note: Using full elements\n";
558818e0025SJeremy L Thompson     code << "      // -- Input fields --\n";
5599e201c85SYohann     for (CeedInt i = 0; i < num_input_fields; i++) {
560818e0025SJeremy L Thompson       code << "      // ---- Input field " << i << " ----\n";
5619e201c85SYohann       code << "      CeedScalar* r_q_" << i << " = r_t_" << i << ";\n";
562ac421f39SYohann     }
563818e0025SJeremy L Thompson     code << "      // -- Output fields --\n";
5649e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
565818e0025SJeremy L Thompson       code << "      // ---- Output field " << i << " ----\n";
5669e201c85SYohann       code << "      CeedScalar* r_qq_" << i << " = r_tt_" << i << ";\n";
567ac421f39SYohann     }
568ac421f39SYohann   }
569818e0025SJeremy L Thompson   code << "\n      // -- QFunction Inputs and outputs --\n";
5709e201c85SYohann   code << "      CeedScalar* in[" << num_input_fields << "];\n";
5719e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
572818e0025SJeremy L Thompson     code << "      // ---- Input field " << i << " ----\n";
5739e201c85SYohann     code << "      in[" << i << "] = r_q_" << i << ";\n";
5744d537eeaSYohann   }
5759e201c85SYohann   code << "      CeedScalar* out[" << num_output_fields << "];\n";
5769e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
577818e0025SJeremy L Thompson     code << "      // ---- Output field " << i << " ----\n";
5789e201c85SYohann     code << "      out[" << i << "] = r_qq_" << i << ";\n";
5794d537eeaSYohann   }
580818e0025SJeremy L Thompson   code << "\n      // -- Apply QFunction --\n";
58109095acaSJeremy L Thompson   code << "      " << qfunction_name << "(ctx, ";
5829e201c85SYohann   if (dim != 3 || use_collograd_parallelization) {
583ac421f39SYohann     code << "1";
584ac421f39SYohann   } else {
5859e201c85SYohann     code << "Q_1d";
586ac421f39SYohann   }
587ac421f39SYohann   code << ", in, out);\n";
5889e201c85SYohann   if (use_collograd_parallelization) {
589818e0025SJeremy L Thompson     code << "      // -- Output fields --\n";
5909e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
591818e0025SJeremy L Thompson       code << "      // ---- Output field " << i << " ----\n";
5922b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
593ac421f39SYohann       // Basis action
5949e201c85SYohann       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
5959e201c85SYohann       switch (eval_mode) {
596ac421f39SYohann         case CEED_EVAL_NONE:
5979e201c85SYohann           code << "      for (CeedInt j = 0; j < num_comp_out_" << i << " ; ++j) {\n";
5989e201c85SYohann           code << "        r_tt_" << i << "[q + j*Q_1d] = r_qq_" << i << "[j];\n";
599ac421f39SYohann           code << "      }\n";
600ac421f39SYohann           break;  // No action
601ac421f39SYohann         case CEED_EVAL_INTERP:
6029e201c85SYohann           code << "      for (CeedInt j = 0; j < num_comp_out_" << i << " ; ++j) {\n";
6039e201c85SYohann           code << "        r_tt_" << i << "[q + j*Q_1d] = r_qq_" << i << "[j];\n";
604ac421f39SYohann           code << "      }\n";
605ac421f39SYohann           break;
606ac421f39SYohann         case CEED_EVAL_GRAD:
6079e201c85SYohann           code << "      gradColloTranspose3d<num_comp_out_" << i << ",Q_1d>(data, q, r_qq_" << i << ", s_G_out_" << i << ", r_tt_" << i << ");\n";
608ac421f39SYohann           break;
609ac421f39SYohann         case CEED_EVAL_WEIGHT:
610ac421f39SYohann           break;  // Should not occur
611ac421f39SYohann         case CEED_EVAL_DIV:
612ac421f39SYohann           break;  // TODO: Not implemented
613ac421f39SYohann         case CEED_EVAL_CURL:
614ac421f39SYohann           break;  // TODO: Not implemented
615ac421f39SYohann       }
616ac421f39SYohann     }
617ac421f39SYohann     code << "    }\n";
618ac421f39SYohann   }
619241a4b83SYohann 
620241a4b83SYohann   // Output basis apply if needed
621ac421f39SYohann   // Generate the correct eval mode code for each output
622818e0025SJeremy L Thompson   code << "\n    // -- Output field basis action and restrictions --\n";
6239e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
624818e0025SJeremy L Thompson     code << "    // ---- Output field " << i << " ----\n";
6259e201c85SYohann     // Get elem_size, eval_mode, num_comp
626edb2538eSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
627edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
6282b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
629edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
6309e201c85SYohann     // TODO put in a function
631241a4b83SYohann     // Basis action
6329e201c85SYohann     code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
6339e201c85SYohann     switch (eval_mode) {
634241a4b83SYohann       case CEED_EVAL_NONE:
6359e201c85SYohann         code << "    CeedScalar* r_v_" << i << " = r_tt_" << i << ";\n";
636241a4b83SYohann         break;  // No action
637241a4b83SYohann       case CEED_EVAL_INTERP:
6389e201c85SYohann         code << "    CeedScalar r_v_" << i << "[num_comp_out_" << i << "*P_out_" << i << "];\n";
6392b730f8bSJeremy L Thompson         code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_out_" << i << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i
6402b730f8bSJeremy L Thompson              << ", s_B_out_" << i << ", r_v_" << i << ");\n";
641241a4b83SYohann         break;
642241a4b83SYohann       case CEED_EVAL_GRAD:
6439e201c85SYohann         code << "    CeedScalar r_v_" << i << "[num_comp_out_" << i << "*P_out_" << i << "];\n";
6449e201c85SYohann         if (use_collograd_parallelization) {
6452b730f8bSJeremy L Thompson           code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_out_" << i << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i
6462b730f8bSJeremy L Thompson                << ", s_B_out_" << i << ", r_v_" << i << ");\n";
647ac421f39SYohann         } else {
6489e201c85SYohann           CeedInt P_1d;
6492b730f8bSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
6502b730f8bSJeremy L Thompson           CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
6512b730f8bSJeremy L Thompson           code << "    GradTranspose" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp_out_" << i
6522b730f8bSJeremy L Thompson                << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i << ", s_B_out_" << i << ", s_G_out_" << i << ", r_v_" << i << ");\n";
653ac421f39SYohann         }
654241a4b83SYohann         break;
655e9f4dca0SJeremy L Thompson       // LCOV_EXCL_START
656241a4b83SYohann       case CEED_EVAL_WEIGHT: {
657*6e536b99SJeremy L Thompson         return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
658241a4b83SYohann         break;  // Should not occur
659241a4b83SYohann       }
660241a4b83SYohann       case CEED_EVAL_DIV:
661bcbe1c99SJeremy L Thompson       case CEED_EVAL_CURL: {
662*6e536b99SJeremy L Thompson         return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
663bcbe1c99SJeremy L Thompson         break;  // Should not occur
664bcbe1c99SJeremy L Thompson       }
665e9f4dca0SJeremy L Thompson         // LCOV_EXCL_STOP
666241a4b83SYohann     }
6679e201c85SYohann     // TODO put in a function
668920dcdc4Sjeremylt     // Restriction
6699e201c85SYohann     bool is_strided;
670edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
6719e201c85SYohann     if (!is_strided) {
6729e201c85SYohann       CeedInt comp_stride;
673ca735530SJeremy L Thompson 
674edb2538eSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
675ca735530SJeremy L Thompson       code << "    const CeedInt l_size_out_" << i << " = " << l_size << ";\n";
676edb2538eSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
6779e201c85SYohann       code << "    // CompStride: " << comp_stride << "\n";
678edb2538eSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
679edb2538eSJeremy L Thompson       data->indices.outputs[i] = rstr_data->d_ind;
680ca735530SJeremy L Thompson       code << "    writeDofsOffset" << dim << "d<num_comp_out_" << i << ", " << comp_stride << ", P_out_" << i << ">(data, l_size_out_" << i
6812b730f8bSJeremy L Thompson            << ", elem, indices.outputs[" << i << "], r_v_" << i << ", d_v_" << i << ");\n";
682920dcdc4Sjeremylt     } else {
6839e201c85SYohann       bool    has_backend_strides;
6849e201c85SYohann       CeedInt num_elem;
685ca735530SJeremy L Thompson 
686edb2538eSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
687edb2538eSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
6889e201c85SYohann       CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
689ca735530SJeremy L Thompson 
6909e201c85SYohann       if (!has_backend_strides) {
69156c48462SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
692920dcdc4Sjeremylt       }
693920dcdc4Sjeremylt       code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
6942b730f8bSJeremy L Thompson       code << "    writeDofsStrided" << dim << "d<num_comp_out_" << i << ",P_out_" << i << "," << strides[0] << "," << strides[1] << "," << strides[2]
6952b730f8bSJeremy L Thompson            << ">(data, elem, r_v_" << i << ", d_v_" << i << ");\n";
696920dcdc4Sjeremylt     }
697241a4b83SYohann   }
698241a4b83SYohann 
699241a4b83SYohann   code << "  }\n";
700818e0025SJeremy L Thompson   code << "}\n";
701818e0025SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n\n";
702241a4b83SYohann 
703ab213215SJeremy L Thompson   // View kernel for debugging
70423d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "Generated Operator Kernels:\n");
7053f21f6b1SJeremy L Thompson   CeedDebug(ceed, code.str().c_str());
706241a4b83SYohann 
707eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Cuda(ceed, code.str().c_str(), &data->module, 1, "T_1D", CeedIntMax(Q_1d, data->max_P_1d)));
708eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, operator_name.c_str(), &data->op));
709241a4b83SYohann 
7102b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
711e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
712241a4b83SYohann }
7132a86cc9dSSebastian Grimberg 
714ab213215SJeremy L Thompson //------------------------------------------------------------------------------
715