xref: /libCEED/backends/cuda-gen/ceed-cuda-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.
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>
1249aac155SJeremy L Thompson #include <ceed/jit-source/cuda/cuda-types.h>
139e201c85SYohann #include <ceed/jit-tools.h>
143d576824SJeremy L Thompson #include <cuda_runtime.h>
152b730f8bSJeremy L Thompson 
16241a4b83SYohann #include <iostream>
17241a4b83SYohann #include <sstream>
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 //------------------------------------------------------------------------------
26ab213215SJeremy L Thompson // Build singe operator kernel
27ab213215SJeremy L Thompson //------------------------------------------------------------------------------
28*eb7e6cafSJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Cuda_gen(CeedOperator op) {
29241a4b83SYohann   using std::ostringstream;
30241a4b83SYohann   using std::string;
319e201c85SYohann   bool is_setup_done;
322b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
339e201c85SYohann   if (is_setup_done) return CEED_ERROR_SUCCESS;
34241a4b83SYohann   Ceed ceed;
352b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
36241a4b83SYohann   CeedOperator_Cuda_gen *data;
372b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
38241a4b83SYohann   CeedQFunction           qf;
39241a4b83SYohann   CeedQFunction_Cuda_gen *qf_data;
402b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
412b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
42e79b91d9SJeremy L Thompson   CeedSize lsize;
432b730f8bSJeremy L Thompson   CeedInt  Q, P_1d = 0, Q_1d = 0, elem_size, num_input_fields, num_output_fields, num_comp, dim = 1;
442b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
459e201c85SYohann   Q_1d = Q;
469e201c85SYohann   CeedOperatorField *op_input_fields, *op_output_fields;
472b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
489e201c85SYohann   CeedQFunctionField *qf_input_fields, *qf_output_fields;
492b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
509e201c85SYohann   CeedEvalMode              eval_mode;
51241a4b83SYohann   CeedBasis                 basis;
52241a4b83SYohann   CeedBasis_Cuda_shared    *basis_data;
53241a4b83SYohann   CeedElemRestriction       Erestrict;
54461525f5SNatalie Beams   CeedElemRestriction_Cuda *restr_data;
55241a4b83SYohann 
569e201c85SYohann   // TODO: put in a function?
570b454692Sjeremylt   // Check for restriction only identity operator
580b454692Sjeremylt   bool is_identity_qf;
592b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
600b454692Sjeremylt   if (is_identity_qf) {
619e201c85SYohann     CeedEvalMode eval_mode_in, eval_mode_out;
622b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
632b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
646574a04fSJeremy L Thompson     CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
656574a04fSJeremy L Thompson               "Backend does not implement restriction only identity operators");
660b454692Sjeremylt   }
670b454692Sjeremylt 
68241a4b83SYohann   ostringstream code;
69241a4b83SYohann 
709e201c85SYohann   // TODO: put in a function?
71f1a13f77SYohann Dudouit   // Add atomicAdd function for old NVidia architectures
72f1a13f77SYohann Dudouit   struct cudaDeviceProp prop;
73f1a13f77SYohann Dudouit   Ceed_Cuda            *ceed_data;
742b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_data));
752b730f8bSJeremy L Thompson   CeedCallBackend(cudaGetDeviceProperties(&prop, ceed_data->device_id));
7680a9ef05SNatalie Beams   if ((prop.major < 6) && (CEED_SCALAR_TYPE != CEED_SCALAR_FP32)) {
779e201c85SYohann     char *atomic_add_path, *atomic_add_source;
782b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-atomic-add-fallback.h", &atomic_add_path));
799e201c85SYohann     CeedDebug256(ceed, 2, "----- Loading Atomic Add Source -----\n");
802b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, atomic_add_path, &atomic_add_source));
819e201c85SYohann     code << atomic_add_source;
822b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&atomic_add_path));
832b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&atomic_add_source));
84f1a13f77SYohann Dudouit   }
85f1a13f77SYohann Dudouit 
869e201c85SYohann   // Load basis source files
879e201c85SYohann   // TODO: generalize to accept different device functions?
889e201c85SYohann   {
899e201c85SYohann     char *tensor_basis_kernel_path, *tensor_basis_kernel_source;
902b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-shared-basis-tensor-templates.h", &tensor_basis_kernel_path));
919e201c85SYohann     CeedDebug256(ceed, 2, "----- Loading Tensor Basis Kernel Source -----\n");
922b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, tensor_basis_kernel_path, &tensor_basis_kernel_source));
939e201c85SYohann     code << tensor_basis_kernel_source;
942b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&tensor_basis_kernel_path));
952b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&tensor_basis_kernel_source));
969e201c85SYohann   }
979e201c85SYohann   {
989e201c85SYohann     char *cuda_gen_template_path, *cuda_gen_template_source;
992b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-gen-templates.h", &cuda_gen_template_path));
1009e201c85SYohann     CeedDebug256(ceed, 2, "----- Loading Cuda-Gen Template Source -----\n");
1012b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, cuda_gen_template_path, &cuda_gen_template_source));
1029e201c85SYohann     code << cuda_gen_template_source;
1032b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&cuda_gen_template_path));
1042b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&cuda_gen_template_source));
1059e201c85SYohann   }
106241a4b83SYohann 
1079e201c85SYohann   // Get QFunction source and name
1089e201c85SYohann   string q_function_source(qf_data->q_function_source);
1099e201c85SYohann   string q_function_name(qf_data->q_function_name);
1109e201c85SYohann   string operator_name;
111204bfdd7SJeremy L Thompson   operator_name = "CeedKernelCudaGenOperator_" + q_function_name;
1124d537eeaSYohann 
1139e201c85SYohann   // Find dim, P_1d, Q_1d
1149e201c85SYohann   data->max_P_1d = 0;
1159e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1162b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1171da99368SJeremy L Thompson     if (basis != CEED_BASIS_COLLOCATED) {
1182b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1192b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1201da99368SJeremy L Thompson 
1219e201c85SYohann       // Collect dim, P_1d, and Q_1d
1222b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &dim));
1231da99368SJeremy L Thompson       bool isTensor;
1242b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &isTensor));
1256574a04fSJeremy L Thompson       CeedCheck(isTensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
1262b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
1272b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
1289e201c85SYohann       if (P_1d > data->max_P_1d) data->max_P_1d = P_1d;
1291da99368SJeremy L Thompson     }
1301da99368SJeremy L Thompson   }
1319e201c85SYohann   // Check output bases for Q_1d, dim as well
132dc64899eSYohann   //   The only input basis might be CEED_BASIS_COLLOCATED
1339e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1342b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1350f54b25eSjeremylt 
1361da99368SJeremy L Thompson     if (basis != CEED_BASIS_COLLOCATED) {
1372b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1382b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1390f54b25eSjeremylt 
1409e201c85SYohann       // Collect Q_1d
1412b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &dim));
1421da99368SJeremy L Thompson       bool isTensor;
1432b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &isTensor));
1446574a04fSJeremy L Thompson       CeedCheck(isTensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
1452b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
1461da99368SJeremy L Thompson     }
1471da99368SJeremy L Thompson   }
1481da99368SJeremy L Thompson   data->dim  = dim;
1499e201c85SYohann   data->Q_1d = Q_1d;
1501da99368SJeremy L Thompson 
1519e201c85SYohann   // Only use 3D collocated gradient parallelization strategy when gradient is computed
1529e201c85SYohann   // TODO: put in a function?
1539e201c85SYohann   bool use_collograd_parallelization = false;
1549e201c85SYohann   if (dim == 3) {
1559e201c85SYohann     bool was_grad_found = false;
1569e201c85SYohann     for (CeedInt i = 0; i < num_input_fields; i++) {
1572b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1589e201c85SYohann       if (eval_mode == CEED_EVAL_GRAD) {
1592b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1602b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1619e201c85SYohann         use_collograd_parallelization = !!basis_data->d_collo_grad_1d && (was_grad_found ? use_collograd_parallelization : true);
1629e201c85SYohann         was_grad_found                = true;
1639e201c85SYohann       }
1649e201c85SYohann     }
1659e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
1662b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1679e201c85SYohann       if (eval_mode == CEED_EVAL_GRAD) {
1682b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1692b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1709e201c85SYohann         use_collograd_parallelization = !!basis_data->d_collo_grad_1d && (was_grad_found ? use_collograd_parallelization : true);
1719e201c85SYohann         was_grad_found                = true;
1729e201c85SYohann       }
1739e201c85SYohann     }
1741da99368SJeremy L Thompson   }
1751da99368SJeremy L Thompson 
1769e201c85SYohann   // Define CEED_Q_VLA
1779e201c85SYohann   code << "\n#undef CEED_Q_VLA\n";
1789e201c85SYohann   if (dim != 3 || use_collograd_parallelization) {
1799e201c85SYohann     code << "#define CEED_Q_VLA 1\n\n";
1809e201c85SYohann   } else {
1819e201c85SYohann     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
1829e201c85SYohann   }
1839e201c85SYohann 
1849e201c85SYohann   code << q_function_source;
185241a4b83SYohann 
186241a4b83SYohann   // Setup
187818e0025SJeremy L Thompson   code << "\n// -----------------------------------------------------------------------------\n";
1882b730f8bSJeremy L Thompson   code << "\nextern \"C\" __global__ void " << operator_name
1892b730f8bSJeremy L Thompson        << "(CeedInt num_elem, void* ctx, FieldsInt_Cuda indices, Fields_Cuda fields, Fields_Cuda B, Fields_Cuda G, CeedScalar* W) {\n";
1909e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1912b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1929e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
1939e201c85SYohann       code << "  const CeedScalar* d_u_" << i << " = fields.inputs[" << i << "];\n";
194241a4b83SYohann     }
195241a4b83SYohann   }
196241a4b83SYohann 
1979e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1989e201c85SYohann     code << "  CeedScalar* d_v_" << i << " = fields.outputs[" << i << "];\n";
199241a4b83SYohann   }
2001da99368SJeremy L Thompson 
2019e201c85SYohann   code << "  const CeedInt dim = " << dim << ";\n";
2029e201c85SYohann   code << "  const CeedInt Q_1d = " << Q_1d << ";\n";
2031da99368SJeremy L Thompson 
204241a4b83SYohann   code << "  extern __shared__ CeedScalar slice[];\n";
2059e201c85SYohann   // TODO put in a function? InitSharedData_Cuda?
2069e201c85SYohann   code << "  SharedData_Cuda data;\n";
2079e201c85SYohann   code << "  data.t_id_x = threadIdx.x;\n";
2089e201c85SYohann   code << "  data.t_id_y = threadIdx.y;\n";
2099e201c85SYohann   code << "  data.t_id_z = threadIdx.z;\n";
2109e201c85SYohann   code << "  data.t_id  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
2119e201c85SYohann   code << "  data.slice = slice+data.t_id_z*T_1D" << (dim > 1 ? "*T_1D" : "") << ";\n";
212920dcdc4Sjeremylt 
213818e0025SJeremy L Thompson   code << "\n  // -- Input field constants and basis data --\n";
2149e201c85SYohann   // TODO: Put in a function?
215ac421f39SYohann   // Initialize constants, and matrices B and G
2169e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
217818e0025SJeremy L Thompson     code << "  // ---- Input field " << i << " ----\n";
2189e201c85SYohann     // Get elem_size, eval_mode, num_comp
2192b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &Erestrict));
2202b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size));
2212b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
2222b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp));
223920dcdc4Sjeremylt 
224920dcdc4Sjeremylt     // Set field constants
2259e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {
2262b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
227920dcdc4Sjeremylt       if (basis != CEED_BASIS_COLLOCATED) {
2282b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
2299e201c85SYohann         code << "  const CeedInt P_in_" << i << " = " << P_1d << ";\n";
230920dcdc4Sjeremylt       } else {
2319e201c85SYohann         code << "  const CeedInt P_in_" << i << " = " << Q_1d << ";\n";
232920dcdc4Sjeremylt       }
2339e201c85SYohann       code << "  const CeedInt num_comp_in_" << i << " = " << num_comp << ";\n";
234920dcdc4Sjeremylt     }
235920dcdc4Sjeremylt 
236920dcdc4Sjeremylt     // Load basis data
2379e201c85SYohann     code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
2389e201c85SYohann     switch (eval_mode) {
239241a4b83SYohann       case CEED_EVAL_NONE:
240241a4b83SYohann         break;
241241a4b83SYohann       case CEED_EVAL_INTERP:
2422b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
2439e201c85SYohann         data->B.inputs[i] = basis_data->d_interp_1d;
2449e201c85SYohann         code << "  __shared__ CeedScalar s_B_in_" << i << "[" << P_1d * Q_1d << "];\n";
2459e201c85SYohann         code << "  loadMatrix<P_in_" << i << ",Q_1d>(data, B.inputs[" << i << "], s_B_in_" << i << ");\n";
246241a4b83SYohann         break;
247241a4b83SYohann       case CEED_EVAL_GRAD:
2482b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
2499e201c85SYohann         data->B.inputs[i] = basis_data->d_interp_1d;
2509e201c85SYohann         code << "  __shared__ CeedScalar s_B_in_" << i << "[" << P_1d * Q_1d << "];\n";
2519e201c85SYohann         code << "  loadMatrix<P_in_" << i << ",Q_1d>(data, B.inputs[" << i << "], s_B_in_" << i << ");\n";
2529e201c85SYohann         if (use_collograd_parallelization) {
2539e201c85SYohann           data->G.inputs[i] = basis_data->d_collo_grad_1d;
2549e201c85SYohann           code << "  __shared__ CeedScalar s_G_in_" << i << "[" << Q_1d * Q_1d << "];\n";
2559e201c85SYohann           code << "  loadMatrix<Q_1d,Q_1d>(data, G.inputs[" << i << "], s_G_in_" << i << ");\n";
256ac421f39SYohann         } else {
2579e201c85SYohann           bool has_collo_grad = !!basis_data->d_collo_grad_1d;
2589e201c85SYohann           data->G.inputs[i]   = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2599e201c85SYohann           code << "  __shared__ CeedScalar s_G_in_" << i << "[" << Q_1d * (has_collo_grad ? Q_1d : P_1d) << "];\n";
2602b730f8bSJeremy 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
2612b730f8bSJeremy L Thompson                << ");\n";
262ac421f39SYohann         }
263241a4b83SYohann         break;
264241a4b83SYohann       case CEED_EVAL_WEIGHT:
265241a4b83SYohann         break;  // No action
266241a4b83SYohann       case CEED_EVAL_DIV:
267241a4b83SYohann         break;  // TODO: Not implemented
268241a4b83SYohann       case CEED_EVAL_CURL:
269241a4b83SYohann         break;  // TODO: Not implemented
270241a4b83SYohann     }
271241a4b83SYohann   }
272920dcdc4Sjeremylt 
273818e0025SJeremy L Thompson   code << "\n  // -- Output field constants and basis data --\n";
2749e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
275818e0025SJeremy L Thompson     code << "  // ---- Output field " << i << " ----\n";
2769e201c85SYohann     // Get elem_size, eval_mode, num_comp
2772b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &Erestrict));
2782b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size));
2792b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
2802b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp));
281920dcdc4Sjeremylt 
282920dcdc4Sjeremylt     // Set field constants
2832b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
284920dcdc4Sjeremylt     if (basis != CEED_BASIS_COLLOCATED) {
2852b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
2869e201c85SYohann       code << "  const CeedInt P_out_" << i << " = " << P_1d << ";\n";
287920dcdc4Sjeremylt     } else {
2889e201c85SYohann       code << "  const CeedInt P_out_" << i << " = " << Q_1d << ";\n";
289920dcdc4Sjeremylt     }
2909e201c85SYohann     code << "  const CeedInt num_comp_out_" << i << " = " << num_comp << ";\n";
291920dcdc4Sjeremylt 
292920dcdc4Sjeremylt     // Load basis data
2939e201c85SYohann     code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
2949e201c85SYohann     switch (eval_mode) {
295920dcdc4Sjeremylt       case CEED_EVAL_NONE:
296920dcdc4Sjeremylt         break;  // No action
297920dcdc4Sjeremylt       case CEED_EVAL_INTERP:
2982b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
2999e201c85SYohann         data->B.outputs[i] = basis_data->d_interp_1d;
3009e201c85SYohann         code << "  __shared__ CeedScalar s_B_out_" << i << "[" << P_1d * Q_1d << "];\n";
3019e201c85SYohann         code << "  loadMatrix<P_out_" << i << ",Q_1d>(data, B.outputs[" << i << "], s_B_out_" << i << ");\n";
302241a4b83SYohann         break;
303241a4b83SYohann       case CEED_EVAL_GRAD:
3042b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
3059e201c85SYohann         data->B.outputs[i] = basis_data->d_interp_1d;
3069e201c85SYohann         code << "  __shared__ CeedScalar s_B_out_" << i << "[" << P_1d * Q_1d << "];\n";
3079e201c85SYohann         code << "  loadMatrix<P_out_" << i << ",Q_1d>(data, B.outputs[" << i << "], s_B_out_" << i << ");\n";
3089e201c85SYohann         if (use_collograd_parallelization) {
3099e201c85SYohann           data->G.outputs[i] = basis_data->d_collo_grad_1d;
3109e201c85SYohann           code << "  __shared__ CeedScalar s_G_out_" << i << "[" << Q_1d * Q_1d << "];\n";
3119e201c85SYohann           code << "  loadMatrix<Q_1d,Q_1d>(data, G.outputs[" << i << "], s_G_out_" << i << ");\n";
312ac421f39SYohann         } else {
3139e201c85SYohann           bool has_collo_grad = !!basis_data->d_collo_grad_1d;
3149e201c85SYohann           data->G.outputs[i]  = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
3159e201c85SYohann           code << "  __shared__ CeedScalar s_G_out_" << i << "[" << Q_1d * (has_collo_grad ? Q_1d : P_1d) << "];\n";
3162b730f8bSJeremy L Thompson           code << "  loadMatrix<" << (has_collo_grad ? "Q_1d" : ("P_out_" + std::to_string(i))) << ",Q_1d>(data, G.outputs[" << i << "], s_G_out_"
3172b730f8bSJeremy L Thompson                << i << ");\n";
318ac421f39SYohann         }
319241a4b83SYohann         break;
320e9f4dca0SJeremy L Thompson       // LCOV_EXCL_START
321241a4b83SYohann       case CEED_EVAL_WEIGHT: {
322241a4b83SYohann         Ceed ceed;
3232b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
3242b730f8bSJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
325241a4b83SYohann         break;  // Should not occur
326241a4b83SYohann       }
327241a4b83SYohann       case CEED_EVAL_DIV:
328241a4b83SYohann         break;  // TODO: Not implemented
329241a4b83SYohann       case CEED_EVAL_CURL:
330241a4b83SYohann         break;  // TODO: Not implemented
331e9f4dca0SJeremy L Thompson                 // LCOV_EXCL_STOP
332241a4b83SYohann     }
333241a4b83SYohann   }
334818e0025SJeremy L Thompson   code << "\n  // -- Element loop --\n";
335ac421f39SYohann   code << "  __syncthreads();\n";
3369e201c85SYohann   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
337241a4b83SYohann   // Input basis apply if needed
338ac421f39SYohann   // Generate the correct eval mode code for each input
339818e0025SJeremy L Thompson   code << "    // -- Input field restrictions and basis actions --\n";
3409e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
341818e0025SJeremy L Thompson     code << "    // ---- Input field " << i << " ----\n";
3429e201c85SYohann     // Get elem_size, eval_mode, num_comp
3432b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &Erestrict));
3442b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size));
3452b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
3462b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp));
347920dcdc4Sjeremylt 
3489e201c85SYohann     // TODO: put in a function?
349920dcdc4Sjeremylt     // Restriction
3502b730f8bSJeremy L Thompson     if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_collograd_parallelization)) {
3519e201c85SYohann       code << "    CeedScalar r_u_" << i << "[num_comp_in_" << i << "*P_in_" << i << "];\n";
3529a2291e3SJeremy L Thompson 
3539e201c85SYohann       bool is_strided;
3542b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionIsStrided(Erestrict, &is_strided));
3559e201c85SYohann       if (!is_strided) {
3562b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetLVectorSize(Erestrict, &lsize));
3575c7b696cSjeremylt         code << "    const CeedInt lsize_in_" << i << " = " << lsize << ";\n";
3589e201c85SYohann         CeedInt comp_stride;
3592b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetCompStride(Erestrict, &comp_stride));
3609e201c85SYohann         code << "    // CompStride: " << comp_stride << "\n";
3612b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetData(Erestrict, &restr_data));
3629e201c85SYohann         data->indices.inputs[i] = restr_data->d_ind;
3632b730f8bSJeremy L Thompson         code << "    readDofsOffset" << dim << "d<num_comp_in_" << i << ", " << comp_stride << ", P_in_" << i << ">(data, lsize_in_" << i
3642b730f8bSJeremy L Thompson              << ", elem, indices.inputs[" << i << "], d_u_" << i << ", r_u_" << i << ");\n";
365ccedf6b0Sjeremylt       } else {
366920dcdc4Sjeremylt         bool backendstrides;
3672b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionHasBackendStrides(Erestrict, &backendstrides));
3689e201c85SYohann         CeedInt num_elem;
3692b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumElements(Erestrict, &num_elem));
3709e201c85SYohann         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
371920dcdc4Sjeremylt         if (!backendstrides) {
3722b730f8bSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetStrides(Erestrict, &strides));
373ccedf6b0Sjeremylt         }
374920dcdc4Sjeremylt         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
3752b730f8bSJeremy L Thompson         code << "    readDofsStrided" << dim << "d<num_comp_in_" << i << ",P_in_" << i << "," << strides[0] << "," << strides[1] << "," << strides[2]
3762b730f8bSJeremy L Thompson              << ">(data, elem, d_u_" << i << ", r_u_" << i << ");\n";
377920dcdc4Sjeremylt       }
378920dcdc4Sjeremylt     }
379920dcdc4Sjeremylt 
3809e201c85SYohann     // TODO: put in a function?
381920dcdc4Sjeremylt     // Basis action
3829e201c85SYohann     code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
3839e201c85SYohann     switch (eval_mode) {
384920dcdc4Sjeremylt       case CEED_EVAL_NONE:
3859e201c85SYohann         if (!use_collograd_parallelization) {
3869e201c85SYohann           code << "    CeedScalar* r_t_" << i << " = r_u_" << i << ";\n";
387920dcdc4Sjeremylt         }
388920dcdc4Sjeremylt         break;
389920dcdc4Sjeremylt       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";
393241a4b83SYohann         break;
394241a4b83SYohann       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";
399ac421f39SYohann         } 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";
406ac421f39SYohann         }
407241a4b83SYohann         break;
408241a4b83SYohann       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";
414241a4b83SYohann         break;  // No action
415241a4b83SYohann       case CEED_EVAL_DIV:
416241a4b83SYohann         break;  // TODO: Not implemented
417241a4b83SYohann       case CEED_EVAL_CURL:
418241a4b83SYohann         break;  // TODO: Not implemented
419241a4b83SYohann     }
420241a4b83SYohann   }
421ac421f39SYohann 
422ea61e9acSJeremy L Thompson   // TODO: put in a function + separate collograd logic
423241a4b83SYohann   // Q function
424818e0025SJeremy L Thompson   code << "\n    // -- Output field setup --\n";
4259e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
426818e0025SJeremy L Thompson     code << "\n    // ---- Output field " << i << " ----\n";
4272b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
4282b730f8bSJeremy L Thompson     if (eval_mode == CEED_EVAL_GRAD) {
4299e201c85SYohann       if (use_collograd_parallelization) {
430ac421f39SYohann         // Accumulator for gradient slices
4319e201c85SYohann         code << "    CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*Q_1d];\n";
4329e201c85SYohann         code << "    for (CeedInt i = 0; i < num_comp_out_" << i << "; i++) {\n";
4339e201c85SYohann         code << "      for (CeedInt j = 0; j < Q_1d; ++j) {\n";
4349e201c85SYohann         code << "        r_tt_" << i << "[j + i*Q_1d] = 0.0;\n";
435ac421f39SYohann         code << "      }\n";
436ac421f39SYohann         code << "    }\n";
437ac421f39SYohann       } else {
4389e201c85SYohann         code << "    CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*dim*Q_1d];\n";
439241a4b83SYohann       }
440ac421f39SYohann     }
4412b730f8bSJeremy L Thompson     if (eval_mode == CEED_EVAL_NONE || eval_mode == CEED_EVAL_INTERP) {
4429e201c85SYohann       code << "    CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*Q_1d];\n";
443241a4b83SYohann     }
444241a4b83SYohann   }
445ac421f39SYohann   // We treat quadrature points per slice in 3d to save registers
4469e201c85SYohann   if (use_collograd_parallelization) {
4479e201c85SYohann     code << "\n    // Note: Using planes of 3D elements\n";
448ac421f39SYohann     code << "#pragma unroll\n";
4499e201c85SYohann     code << "    for (CeedInt q = 0; q < Q_1d; q++) {\n";
450818e0025SJeremy L Thompson     code << "      // -- Input fields --\n";
4519e201c85SYohann     for (CeedInt i = 0; i < num_input_fields; i++) {
452818e0025SJeremy L Thompson       code << "      // ---- Input field " << i << " ----\n";
4539e201c85SYohann       // Get elem_size, eval_mode, num_comp
4542b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
455ac421f39SYohann       // Basis action
4569e201c85SYohann       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
4579e201c85SYohann       switch (eval_mode) {
458ac421f39SYohann         case CEED_EVAL_NONE:
4599e201c85SYohann           code << "      CeedScalar r_q_" << i << "[num_comp_in_" << i << "];\n";
4609a2291e3SJeremy L Thompson 
4619e201c85SYohann           bool is_strided;
4622b730f8bSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &Erestrict));
4632b730f8bSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionIsStrided(Erestrict, &is_strided));
4649e201c85SYohann           if (!is_strided) {
4652b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetLVectorSize(Erestrict, &lsize));
4665c7b696cSjeremylt             code << "      const CeedInt lsize_in_" << i << " = " << lsize << ";\n";
4679e201c85SYohann             CeedInt comp_stride;
4682b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetCompStride(Erestrict, &comp_stride));
4699e201c85SYohann             code << "      // CompStride: " << comp_stride << "\n";
4702b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetData(Erestrict, &restr_data));
4719e201c85SYohann             data->indices.inputs[i] = restr_data->d_ind;
4722b730f8bSJeremy L Thompson             code << "      readSliceQuadsOffset"
4732b730f8bSJeremy L Thompson                  << "3d<num_comp_in_" << i << ", " << comp_stride << ", Q_1d>(data, lsize_in_" << i << ", elem, q, indices.inputs[" << i << "], d_u_"
4742b730f8bSJeremy L Thompson                  << i << ", r_q_" << i << ");\n";
475920dcdc4Sjeremylt           } else {
4762b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size));
477920dcdc4Sjeremylt             bool backendstrides;
4782b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionHasBackendStrides(Erestrict, &backendstrides));
4799e201c85SYohann             CeedInt num_elem;
4802b730f8bSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetNumElements(Erestrict, &num_elem));
4819e201c85SYohann             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
482920dcdc4Sjeremylt             if (!backendstrides) {
4832b730f8bSJeremy L Thompson               CeedCallBackend(CeedElemRestrictionGetStrides(Erestrict, &strides));
484920dcdc4Sjeremylt             }
485920dcdc4Sjeremylt             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
4862b730f8bSJeremy L Thompson             code << "      readSliceQuadsStrided"
4872b730f8bSJeremy L Thompson                  << "3d<num_comp_in_" << i
4882b730f8bSJeremy L Thompson                  << ",Q_1d"
4892b730f8bSJeremy L Thompson                     ","
4902b730f8bSJeremy L Thompson                  << strides[0] << "," << strides[1] << "," << strides[2] << ">(data, elem, q, d_u_" << i << ", r_q_" << i << ");\n";
491920dcdc4Sjeremylt           }
492ac421f39SYohann           break;
493ac421f39SYohann         case CEED_EVAL_INTERP:
4949e201c85SYohann           code << "      CeedScalar r_q_" << i << "[num_comp_in_" << i << "];\n";
4959e201c85SYohann           code << "      for (CeedInt j = 0; j < num_comp_in_" << i << " ; ++j) {\n";
4969e201c85SYohann           code << "        r_q_" << i << "[j] = r_t_" << i << "[q + j*Q_1d];\n";
497ac421f39SYohann           code << "      }\n";
498ac421f39SYohann           break;
499ac421f39SYohann         case CEED_EVAL_GRAD:
5009e201c85SYohann           code << "      CeedScalar r_q_" << i << "[num_comp_in_" << i << "*dim];\n";
5019e201c85SYohann           code << "      gradCollo3d<num_comp_in_" << i << ",Q_1d>(data, q, r_t_" << i << ", s_G_in_" << i << ", r_q_" << i << ");\n";
502ac421f39SYohann           break;
503ac421f39SYohann         case CEED_EVAL_WEIGHT:
5049e201c85SYohann           code << "      CeedScalar r_q_" << i << "[1];\n";
5059e201c85SYohann           code << "      r_q_" << i << "[0] = r_t_" << i << "[q];\n";
506ac421f39SYohann           break;  // No action
507ac421f39SYohann         case CEED_EVAL_DIV:
508ac421f39SYohann           break;  // TODO: Not implemented
509ac421f39SYohann         case CEED_EVAL_CURL:
510ac421f39SYohann           break;  // TODO: Not implemented
511ac421f39SYohann       }
512ac421f39SYohann     }
513818e0025SJeremy L Thompson     code << "\n      // -- Output fields --\n";
5149e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
515818e0025SJeremy L Thompson       code << "      // ---- Output field " << i << " ----\n";
5162b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
517ac421f39SYohann       // Basis action
5189e201c85SYohann       switch (eval_mode) {
519ac421f39SYohann         case CEED_EVAL_NONE:
5209e201c85SYohann           code << "      CeedScalar r_qq_" << i << "[num_comp_out_" << i << "];\n";
521ac421f39SYohann           break;  // No action
522ac421f39SYohann         case CEED_EVAL_INTERP:
5239e201c85SYohann           code << "      CeedScalar r_qq_" << i << "[num_comp_out_" << i << "];\n";
524ac421f39SYohann           break;
525ac421f39SYohann         case CEED_EVAL_GRAD:
5269e201c85SYohann           code << "      CeedScalar r_qq_" << i << "[num_comp_out_" << i << "*dim];\n";
527ac421f39SYohann           break;
528ac421f39SYohann         case CEED_EVAL_WEIGHT:
529ac421f39SYohann           break;  // Should not occur
530ac421f39SYohann         case CEED_EVAL_DIV:
531ac421f39SYohann           break;  // TODO: Not implemented
532ac421f39SYohann         case CEED_EVAL_CURL:
533ac421f39SYohann           break;  // TODO: Not implemented
534ac421f39SYohann       }
535ac421f39SYohann     }
536ac421f39SYohann   } else {
5379e201c85SYohann     code << "\n      // Note: Using full elements\n";
538818e0025SJeremy L Thompson     code << "      // -- Input fields --\n";
5399e201c85SYohann     for (CeedInt i = 0; i < num_input_fields; i++) {
540818e0025SJeremy L Thompson       code << "      // ---- Input field " << i << " ----\n";
5419e201c85SYohann       code << "      CeedScalar* r_q_" << i << " = r_t_" << i << ";\n";
542ac421f39SYohann     }
543818e0025SJeremy L Thompson     code << "      // -- Output fields --\n";
5449e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
545818e0025SJeremy L Thompson       code << "      // ---- Output field " << i << " ----\n";
5469e201c85SYohann       code << "      CeedScalar* r_qq_" << i << " = r_tt_" << i << ";\n";
547ac421f39SYohann     }
548ac421f39SYohann   }
549818e0025SJeremy L Thompson   code << "\n      // -- QFunction Inputs and outputs --\n";
5509e201c85SYohann   code << "      CeedScalar* in[" << num_input_fields << "];\n";
5519e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
552818e0025SJeremy L Thompson     code << "      // ---- Input field " << i << " ----\n";
5539e201c85SYohann     code << "      in[" << i << "] = r_q_" << i << ";\n";
5544d537eeaSYohann   }
5559e201c85SYohann   code << "      CeedScalar* out[" << num_output_fields << "];\n";
5569e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
557818e0025SJeremy L Thompson     code << "      // ---- Output field " << i << " ----\n";
5589e201c85SYohann     code << "      out[" << i << "] = r_qq_" << i << ";\n";
5594d537eeaSYohann   }
560818e0025SJeremy L Thompson   code << "\n      // -- Apply QFunction --\n";
5619e201c85SYohann   code << "      " << q_function_name << "(ctx, ";
5629e201c85SYohann   if (dim != 3 || use_collograd_parallelization) {
563ac421f39SYohann     code << "1";
564ac421f39SYohann   } else {
5659e201c85SYohann     code << "Q_1d";
566ac421f39SYohann   }
567ac421f39SYohann   code << ", in, out);\n";
5689e201c85SYohann   if (use_collograd_parallelization) {
569818e0025SJeremy L Thompson     code << "      // -- Output fields --\n";
5709e201c85SYohann     for (CeedInt i = 0; i < num_output_fields; i++) {
571818e0025SJeremy L Thompson       code << "      // ---- Output field " << i << " ----\n";
5722b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
573ac421f39SYohann       // Basis action
5749e201c85SYohann       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
5759e201c85SYohann       switch (eval_mode) {
576ac421f39SYohann         case CEED_EVAL_NONE:
5779e201c85SYohann           code << "      for (CeedInt j = 0; j < num_comp_out_" << i << " ; ++j) {\n";
5789e201c85SYohann           code << "        r_tt_" << i << "[q + j*Q_1d] = r_qq_" << i << "[j];\n";
579ac421f39SYohann           code << "      }\n";
580ac421f39SYohann           break;  // No action
581ac421f39SYohann         case CEED_EVAL_INTERP:
5829e201c85SYohann           code << "      for (CeedInt j = 0; j < num_comp_out_" << i << " ; ++j) {\n";
5839e201c85SYohann           code << "        r_tt_" << i << "[q + j*Q_1d] = r_qq_" << i << "[j];\n";
584ac421f39SYohann           code << "      }\n";
585ac421f39SYohann           break;
586ac421f39SYohann         case CEED_EVAL_GRAD:
5879e201c85SYohann           code << "      gradColloTranspose3d<num_comp_out_" << i << ",Q_1d>(data, q, r_qq_" << i << ", s_G_out_" << i << ", r_tt_" << i << ");\n";
588ac421f39SYohann           break;
589ac421f39SYohann         case CEED_EVAL_WEIGHT:
590ac421f39SYohann           break;  // Should not occur
591ac421f39SYohann         case CEED_EVAL_DIV:
592ac421f39SYohann           break;  // TODO: Not implemented
593ac421f39SYohann         case CEED_EVAL_CURL:
594ac421f39SYohann           break;  // TODO: Not implemented
595ac421f39SYohann       }
596ac421f39SYohann     }
597ac421f39SYohann     code << "    }\n";
598ac421f39SYohann   }
599241a4b83SYohann 
600241a4b83SYohann   // Output basis apply if needed
601ac421f39SYohann   // Generate the correct eval mode code for each output
602818e0025SJeremy L Thompson   code << "\n    // -- Output field basis action and restrictions --\n";
6039e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
604818e0025SJeremy L Thompson     code << "    // ---- Output field " << i << " ----\n";
6059e201c85SYohann     // Get elem_size, eval_mode, num_comp
6062b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &Erestrict));
6072b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size));
6082b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
6092b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp));
6109e201c85SYohann     // TODO put in a function
611241a4b83SYohann     // Basis action
6129e201c85SYohann     code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
6139e201c85SYohann     switch (eval_mode) {
614241a4b83SYohann       case CEED_EVAL_NONE:
6159e201c85SYohann         code << "    CeedScalar* r_v_" << i << " = r_tt_" << i << ";\n";
616241a4b83SYohann         break;  // No action
617241a4b83SYohann       case CEED_EVAL_INTERP:
6189e201c85SYohann         code << "    CeedScalar r_v_" << i << "[num_comp_out_" << i << "*P_out_" << i << "];\n";
6192b730f8bSJeremy L Thompson         code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_out_" << i << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i
6202b730f8bSJeremy L Thompson              << ", s_B_out_" << i << ", r_v_" << i << ");\n";
621241a4b83SYohann         break;
622241a4b83SYohann       case CEED_EVAL_GRAD:
6239e201c85SYohann         code << "    CeedScalar r_v_" << i << "[num_comp_out_" << i << "*P_out_" << i << "];\n";
6249e201c85SYohann         if (use_collograd_parallelization) {
6252b730f8bSJeremy L Thompson           code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_out_" << i << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i
6262b730f8bSJeremy L Thompson                << ", s_B_out_" << i << ", r_v_" << i << ");\n";
627ac421f39SYohann         } else {
6289e201c85SYohann           CeedInt P_1d;
6292b730f8bSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
6302b730f8bSJeremy L Thompson           CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
6312b730f8bSJeremy L Thompson           code << "    GradTranspose" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp_out_" << i
6322b730f8bSJeremy L Thompson                << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i << ", s_B_out_" << i << ", s_G_out_" << i << ", r_v_" << i << ");\n";
633ac421f39SYohann         }
634241a4b83SYohann         break;
635e9f4dca0SJeremy L Thompson       // LCOV_EXCL_START
636241a4b83SYohann       case CEED_EVAL_WEIGHT: {
637241a4b83SYohann         Ceed ceed;
6382b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
6392b730f8bSJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
640241a4b83SYohann         break;  // Should not occur
641241a4b83SYohann       }
642241a4b83SYohann       case CEED_EVAL_DIV:
643241a4b83SYohann         break;  // TODO: Not implemented
644241a4b83SYohann       case CEED_EVAL_CURL:
645241a4b83SYohann         break;  // TODO: Not implemented
646e9f4dca0SJeremy L Thompson                 // LCOV_EXCL_STOP
647241a4b83SYohann     }
6489e201c85SYohann     // TODO put in a function
649920dcdc4Sjeremylt     // Restriction
6509e201c85SYohann     bool is_strided;
6512b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionIsStrided(Erestrict, &is_strided));
6529e201c85SYohann     if (!is_strided) {
6532b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetLVectorSize(Erestrict, &lsize));
6545c7b696cSjeremylt       code << "    const CeedInt lsize_out_" << i << " = " << lsize << ";\n";
6559e201c85SYohann       CeedInt comp_stride;
6562b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetCompStride(Erestrict, &comp_stride));
6579e201c85SYohann       code << "    // CompStride: " << comp_stride << "\n";
6582b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetData(Erestrict, &restr_data));
6599e201c85SYohann       data->indices.outputs[i] = restr_data->d_ind;
6602b730f8bSJeremy L Thompson       code << "    writeDofsOffset" << dim << "d<num_comp_out_" << i << ", " << comp_stride << ", P_out_" << i << ">(data, lsize_out_" << i
6612b730f8bSJeremy L Thompson            << ", elem, indices.outputs[" << i << "], r_v_" << i << ", d_v_" << i << ");\n";
662920dcdc4Sjeremylt     } else {
6639e201c85SYohann       bool has_backend_strides;
6642b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionHasBackendStrides(Erestrict, &has_backend_strides));
6659e201c85SYohann       CeedInt num_elem;
6662b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetNumElements(Erestrict, &num_elem));
6679e201c85SYohann       CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
6689e201c85SYohann       if (!has_backend_strides) {
6692b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetStrides(Erestrict, &strides));
670920dcdc4Sjeremylt       }
671920dcdc4Sjeremylt       code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
6722b730f8bSJeremy L Thompson       code << "    writeDofsStrided" << dim << "d<num_comp_out_" << i << ",P_out_" << i << "," << strides[0] << "," << strides[1] << "," << strides[2]
6732b730f8bSJeremy L Thompson            << ">(data, elem, r_v_" << i << ", d_v_" << i << ");\n";
674920dcdc4Sjeremylt     }
675241a4b83SYohann   }
676241a4b83SYohann 
677241a4b83SYohann   code << "  }\n";
678818e0025SJeremy L Thompson   code << "}\n";
679818e0025SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n\n";
680241a4b83SYohann 
681ab213215SJeremy L Thompson   // View kernel for debugging
68246dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "Generated Operator Kernels:\n");
6833f21f6b1SJeremy L Thompson   CeedDebug(ceed, code.str().c_str());
684241a4b83SYohann 
685*eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Cuda(ceed, code.str().c_str(), &data->module, 1, "T_1D", CeedIntMax(Q_1d, data->max_P_1d)));
686*eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, operator_name.c_str(), &data->op));
687241a4b83SYohann 
6882b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
689e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
690241a4b83SYohann }
6912a86cc9dSSebastian Grimberg 
692ab213215SJeremy L Thompson //------------------------------------------------------------------------------
693