xref: /libCEED/backends/cuda-gen/ceed-cuda-gen-operator-build.cpp (revision 034f99fd28bbd688e3b7ed077b6d504c2f765b9c)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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 //------------------------------------------------------------------------------
264b3e95d5SJeremy L Thompson // Determine type of operator
274b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
284b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelData_Cuda_gen(Ceed ceed, CeedInt num_input_fields, CeedOperatorField *op_input_fields,
294b3e95d5SJeremy L Thompson                                                 CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, CeedOperatorField *op_output_fields,
304b3e95d5SJeremy L Thompson                                                 CeedQFunctionField *qf_output_fields, CeedInt *max_P_1d, CeedInt *Q_1d, CeedInt *dim, bool *is_tensor,
314b3e95d5SJeremy L Thompson                                                 bool *use_3d_slices) {
324b3e95d5SJeremy L Thompson   // Find dim, P_1d, Q_1d
334b3e95d5SJeremy L Thompson   *max_P_1d  = 0;
344b3e95d5SJeremy L Thompson   *Q_1d      = 0;
354b3e95d5SJeremy L Thompson   *dim       = 0;
364b3e95d5SJeremy L Thompson   *is_tensor = true;
374b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
384b3e95d5SJeremy L Thompson     CeedBasis basis;
394b3e95d5SJeremy L Thompson 
404b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
414b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
424b3e95d5SJeremy L Thompson       bool    is_field_tensor;
434b3e95d5SJeremy L Thompson       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
444b3e95d5SJeremy L Thompson 
454b3e95d5SJeremy L Thompson       // Collect dim, P_1d, and Q_1d
464b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
474b3e95d5SJeremy L Thompson       CeedCheck(is_field_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
484b3e95d5SJeremy L Thompson       *is_tensor = *is_tensor && is_field_tensor;
494b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
504b3e95d5SJeremy L Thompson       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
514b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
524b3e95d5SJeremy L Thompson       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
534b3e95d5SJeremy L Thompson       *dim = field_dim;
544b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
554b3e95d5SJeremy L Thompson       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
564b3e95d5SJeremy L Thompson       *Q_1d = field_Q_1d;
574b3e95d5SJeremy L Thompson     }
584b3e95d5SJeremy L Thompson   }
594b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
604b3e95d5SJeremy L Thompson     CeedBasis basis;
614b3e95d5SJeremy L Thompson 
624b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
634b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
644b3e95d5SJeremy L Thompson       bool    is_field_tensor;
654b3e95d5SJeremy L Thompson       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
664b3e95d5SJeremy L Thompson 
674b3e95d5SJeremy L Thompson       // Collect dim, P_1d, and Q_1d
684b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
694b3e95d5SJeremy L Thompson       CeedCheck(is_field_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
704b3e95d5SJeremy L Thompson       *is_tensor = *is_tensor && is_field_tensor;
714b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
724b3e95d5SJeremy L Thompson       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
734b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
744b3e95d5SJeremy L Thompson       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
754b3e95d5SJeremy L Thompson       *dim = field_dim;
764b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
774b3e95d5SJeremy L Thompson       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
784b3e95d5SJeremy L Thompson       *Q_1d = field_Q_1d;
794b3e95d5SJeremy L Thompson     }
804b3e95d5SJeremy L Thompson   }
814b3e95d5SJeremy L Thompson 
824b3e95d5SJeremy L Thompson   // Only use 3D collocated gradient parallelization strategy when gradient is computed
834b3e95d5SJeremy L Thompson   *use_3d_slices = false;
844b3e95d5SJeremy L Thompson   if (*dim == 3) {
854b3e95d5SJeremy L Thompson     bool was_grad_found = false;
864b3e95d5SJeremy L Thompson 
874b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
884b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
894b3e95d5SJeremy L Thompson 
904b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
914b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
924b3e95d5SJeremy L Thompson         CeedBasis_Cuda_shared *basis_data;
934b3e95d5SJeremy L Thompson         CeedBasis              basis;
944b3e95d5SJeremy L Thompson 
954b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
964b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
974b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
984b3e95d5SJeremy L Thompson         was_grad_found = true;
994b3e95d5SJeremy L Thompson       }
1004b3e95d5SJeremy L Thompson     }
1014b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
1024b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
1034b3e95d5SJeremy L Thompson 
1044b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1054b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
1064b3e95d5SJeremy L Thompson         CeedBasis_Cuda_shared *basis_data;
1074b3e95d5SJeremy L Thompson         CeedBasis              basis;
1084b3e95d5SJeremy L Thompson 
1094b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1104b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1114b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
1124b3e95d5SJeremy L Thompson         was_grad_found = true;
1134b3e95d5SJeremy L Thompson       }
1144b3e95d5SJeremy L Thompson     }
1154b3e95d5SJeremy L Thompson   }
1164b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1174b3e95d5SJeremy L Thompson }
1184b3e95d5SJeremy L Thompson 
1194b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1204b3e95d5SJeremy L Thompson // Setup fields
1214b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1224b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelFieldData_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedOperatorField op_field,
1234b3e95d5SJeremy L Thompson                                                      CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool use_3d_slices) {
1244b3e95d5SJeremy L Thompson   std::string            var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
1254b3e95d5SJeremy L Thompson   std::string            P_name = "P_1d" + var_suffix, Q_name = "Q_1d";
1264b3e95d5SJeremy L Thompson   std::string            option_name = (is_input ? "inputs" : "outputs");
1274b3e95d5SJeremy L Thompson   CeedEvalMode           eval_mode   = CEED_EVAL_NONE;
1284b3e95d5SJeremy L Thompson   CeedInt                elem_size = 0, num_comp = 0, P_1d = 0;
1294b3e95d5SJeremy L Thompson   CeedElemRestriction    elem_rstr;
1304b3e95d5SJeremy L Thompson   CeedBasis_Cuda_shared *basis_data;
1314b3e95d5SJeremy L Thompson   CeedBasis              basis;
1324b3e95d5SJeremy L Thompson 
1334b3e95d5SJeremy L Thompson   code << "  // -- " << (is_input ? "Input" : "Output") << " field " << i << "\n";
1344b3e95d5SJeremy L Thompson 
1354b3e95d5SJeremy L Thompson   // Get field data
1364b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
1374b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
1384b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1394b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1404b3e95d5SJeremy L Thompson   }
1414b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
1424b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
1434b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1444b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
1454b3e95d5SJeremy L Thompson   }
1464b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
1474b3e95d5SJeremy L Thompson 
1484b3e95d5SJeremy L Thompson   // Set field constants
1494b3e95d5SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
1504b3e95d5SJeremy L Thompson     code << "  const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n";
1514b3e95d5SJeremy L Thompson     code << "  const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n";
1524b3e95d5SJeremy L Thompson   }
1534b3e95d5SJeremy L Thompson 
1544b3e95d5SJeremy L Thompson   // Load basis data
1554b3e95d5SJeremy L Thompson   code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
1564b3e95d5SJeremy L Thompson   switch (eval_mode) {
1574b3e95d5SJeremy L Thompson     case CEED_EVAL_NONE:
1584b3e95d5SJeremy L Thompson       break;
1594b3e95d5SJeremy L Thompson     case CEED_EVAL_INTERP:
1604b3e95d5SJeremy L Thompson       if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
1614b3e95d5SJeremy L Thompson       else data->B.outputs[i] = basis_data->d_interp_1d;
1624b3e95d5SJeremy L Thompson       code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_1d * Q_1d << "];\n";
1634b3e95d5SJeremy L Thompson       code << "  loadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
1644b3e95d5SJeremy L Thompson       break;
1654b3e95d5SJeremy L Thompson     case CEED_EVAL_GRAD:
1664b3e95d5SJeremy L Thompson       if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
1674b3e95d5SJeremy L Thompson       else data->B.outputs[i] = basis_data->d_interp_1d;
1684b3e95d5SJeremy L Thompson       code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_1d * Q_1d << "];\n";
1694b3e95d5SJeremy L Thompson       code << "  loadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
1704b3e95d5SJeremy L Thompson       if (use_3d_slices) {
1714b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d;
1724b3e95d5SJeremy L Thompson         else data->G.outputs[i] = basis_data->d_collo_grad_1d;
1734b3e95d5SJeremy L Thompson         code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * Q_1d << "];\n";
1744b3e95d5SJeremy L Thompson         code << "  loadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
1754b3e95d5SJeremy L Thompson       } else {
1764b3e95d5SJeremy L Thompson         bool has_collo_grad = basis_data->d_collo_grad_1d;
1774b3e95d5SJeremy L Thompson 
1784b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
1794b3e95d5SJeremy L Thompson         else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
1804b3e95d5SJeremy L Thompson         if (has_collo_grad) {
1814b3e95d5SJeremy L Thompson           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * Q_1d << "];\n";
1824b3e95d5SJeremy L Thompson           code << "  loadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
1834b3e95d5SJeremy L Thompson         } else {
1844b3e95d5SJeremy L Thompson           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * P_1d << "];\n";
1854b3e95d5SJeremy L Thompson           code << "  loadMatrix<" << P_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
1864b3e95d5SJeremy L Thompson         }
1874b3e95d5SJeremy L Thompson       }
1884b3e95d5SJeremy L Thompson       break;
1894b3e95d5SJeremy L Thompson     case CEED_EVAL_WEIGHT:
1904b3e95d5SJeremy L Thompson       break;  // No action
1914b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
1924b3e95d5SJeremy L Thompson     case CEED_EVAL_DIV:
1934b3e95d5SJeremy L Thompson       break;  // TODO: Not implemented
1944b3e95d5SJeremy L Thompson     case CEED_EVAL_CURL:
1954b3e95d5SJeremy L Thompson       break;  // TODO: Not implemented
1964b3e95d5SJeremy L Thompson               // LCOV_EXCL_STOP
1974b3e95d5SJeremy L Thompson   }
1984b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1994b3e95d5SJeremy L Thompson }
2004b3e95d5SJeremy L Thompson 
2014b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
2024b3e95d5SJeremy L Thompson // Restriction
2034b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
2044b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelRestriction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedInt dim,
205e93651e5SJeremy L Thompson                                                        CeedInt field_input_buffer[], CeedOperatorField op_field, CeedQFunctionField qf_field,
206e93651e5SJeremy L Thompson                                                        CeedInt Q_1d, bool is_input, bool use_3d_slices) {
2074b3e95d5SJeremy L Thompson   std::string               var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
2084b3e95d5SJeremy L Thompson   std::string               P_name     = "P_1d" + var_suffix;
2094b3e95d5SJeremy L Thompson   CeedEvalMode              eval_mode  = CEED_EVAL_NONE;
2104b3e95d5SJeremy L Thompson   CeedInt                   elem_size = 0, num_comp = 0, P_1d = 0;
2114b3e95d5SJeremy L Thompson   CeedSize                  l_size;
2124b3e95d5SJeremy L Thompson   CeedElemRestriction_Cuda *rstr_data;
2134b3e95d5SJeremy L Thompson   CeedElemRestriction       elem_rstr;
2144b3e95d5SJeremy L Thompson   CeedBasis                 basis;
2154b3e95d5SJeremy L Thompson 
2164b3e95d5SJeremy L Thompson   // Get field data
2174b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
2184b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
2194b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
2204b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
2214b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
2224b3e95d5SJeremy L Thompson   }
2234b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
2244b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
2254b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
2264b3e95d5SJeremy L Thompson   }
2274b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
2284b3e95d5SJeremy L Thompson 
2294b3e95d5SJeremy L Thompson   // Restriction
2304b3e95d5SJeremy L Thompson   if (is_input) {
2314b3e95d5SJeremy L Thompson     // Input
232e93651e5SJeremy L Thompson     if (field_input_buffer[i] != i) {
233e93651e5SJeremy L Thompson       std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]);
234e93651e5SJeremy L Thompson 
235e93651e5SJeremy L Thompson       // Restriction was already done for previous input
236e93651e5SJeremy L Thompson       code << "    CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n";
237e93651e5SJeremy L Thompson     } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices)) {
2384b3e95d5SJeremy L Thompson       bool is_strided;
2394b3e95d5SJeremy L Thompson 
240e93651e5SJeremy L Thompson       if (eval_mode == CEED_EVAL_NONE) {
241e93651e5SJeremy L Thompson         // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated
2424b3e95d5SJeremy L Thompson         code << "    CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n";
243e93651e5SJeremy L Thompson       } else {
244e93651e5SJeremy L Thompson         // Otherwise we're using the scratch space
245e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
246e93651e5SJeremy L Thompson       }
2474b3e95d5SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
2484b3e95d5SJeremy L Thompson       if (!is_strided) {
2494b3e95d5SJeremy L Thompson         CeedInt comp_stride;
2504b3e95d5SJeremy L Thompson 
2514b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
2524b3e95d5SJeremy L Thompson         code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
2534b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
2544b3e95d5SJeremy L Thompson         code << "    // CompStride: " << comp_stride << "\n";
2554b3e95d5SJeremy L Thompson         data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
2564b3e95d5SJeremy L Thompson         code << "    readDofsOffset" << dim << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name << ">(data, l_size" << var_suffix
2574b3e95d5SJeremy L Thompson              << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n";
2584b3e95d5SJeremy L Thompson       } else {
2594b3e95d5SJeremy L Thompson         bool    has_backend_strides;
2604b3e95d5SJeremy L Thompson         CeedInt num_elem;
2614b3e95d5SJeremy L Thompson 
2624b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
2634b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
2644b3e95d5SJeremy L Thompson         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
2654b3e95d5SJeremy L Thompson 
2664b3e95d5SJeremy L Thompson         if (!has_backend_strides) {
2674b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
2684b3e95d5SJeremy L Thompson         }
2694b3e95d5SJeremy L Thompson         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
2704b3e95d5SJeremy L Thompson         code << "    readDofsStrided" << dim << "d<num_comp" << var_suffix << ", " << P_name << "," << strides[0] << "," << strides[1] << ","
2714b3e95d5SJeremy L Thompson              << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n";
2724b3e95d5SJeremy L Thompson       }
2734b3e95d5SJeremy L Thompson     }
2744b3e95d5SJeremy L Thompson   } else {
2754b3e95d5SJeremy L Thompson     // Output
2764b3e95d5SJeremy L Thompson     bool is_strided;
2774b3e95d5SJeremy L Thompson 
2784b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
2794b3e95d5SJeremy L Thompson     if (!is_strided) {
2804b3e95d5SJeremy L Thompson       CeedInt comp_stride;
2814b3e95d5SJeremy L Thompson 
2824b3e95d5SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
2834b3e95d5SJeremy L Thompson       code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
2844b3e95d5SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
2854b3e95d5SJeremy L Thompson       code << "    // CompStride: " << comp_stride << "\n";
2864b3e95d5SJeremy L Thompson       data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
2874b3e95d5SJeremy L Thompson       code << "    writeDofsOffset" << dim << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name << ">(data, l_size" << var_suffix
2884b3e95d5SJeremy L Thompson            << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n";
2894b3e95d5SJeremy L Thompson     } else {
2904b3e95d5SJeremy L Thompson       bool    has_backend_strides;
2914b3e95d5SJeremy L Thompson       CeedInt num_elem;
2924b3e95d5SJeremy L Thompson 
2934b3e95d5SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
2944b3e95d5SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
2954b3e95d5SJeremy L Thompson       CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
2964b3e95d5SJeremy L Thompson 
2974b3e95d5SJeremy L Thompson       if (!has_backend_strides) {
2984b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
2994b3e95d5SJeremy L Thompson       }
3004b3e95d5SJeremy L Thompson       code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
3014b3e95d5SJeremy L Thompson       code << "    writeDofsStrided" << dim << "d<num_comp" << var_suffix << ", " << P_name << "," << strides[0] << "," << strides[1] << ","
3024b3e95d5SJeremy L Thompson            << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n";
3034b3e95d5SJeremy L Thompson     }
3044b3e95d5SJeremy L Thompson   }
3054b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3064b3e95d5SJeremy L Thompson }
3074b3e95d5SJeremy L Thompson 
3084b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
3094b3e95d5SJeremy L Thompson // Basis
3104b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
3114b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelBasis_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedInt dim,
3124b3e95d5SJeremy L Thompson                                                  CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input,
3134b3e95d5SJeremy L Thompson                                                  bool use_3d_slices) {
3144b3e95d5SJeremy L Thompson   std::string         var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
3154b3e95d5SJeremy L Thompson   std::string         P_name = "P_1d" + var_suffix, Q_name = "Q_1d";
3164b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
3174b3e95d5SJeremy L Thompson   CeedInt             elem_size = 0, num_comp = 0, P_1d = 0;
3184b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
3194b3e95d5SJeremy L Thompson   CeedBasis           basis;
3204b3e95d5SJeremy L Thompson 
3214b3e95d5SJeremy L Thompson   // Get field data
3224b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
3234b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
3244b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
3254b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
3264b3e95d5SJeremy L Thompson   }
3274b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
3284b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
3294b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
3304b3e95d5SJeremy L Thompson   }
3314b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
3324b3e95d5SJeremy L Thompson 
3334b3e95d5SJeremy L Thompson   // Basis
3344b3e95d5SJeremy L Thompson   code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
3354b3e95d5SJeremy L Thompson   if (is_input) {
3364b3e95d5SJeremy L Thompson     switch (eval_mode) {
3374b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
3384b3e95d5SJeremy L Thompson         if (!use_3d_slices) {
3394b3e95d5SJeremy L Thompson           code << "    CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n";
3404b3e95d5SJeremy L Thompson         }
3414b3e95d5SJeremy L Thompson         break;
3424b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
3434b3e95d5SJeremy L Thompson         code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
3444b3e95d5SJeremy L Thompson         code << "    Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name
3454b3e95d5SJeremy L Thompson              << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n";
3464b3e95d5SJeremy L Thompson         break;
3474b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
3484b3e95d5SJeremy L Thompson         if (use_3d_slices) {
3494b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
3504b3e95d5SJeremy L Thompson           code << "    Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name
3514b3e95d5SJeremy L Thompson                << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n";
3524b3e95d5SJeremy L Thompson         } else {
3534b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << Q_name << "];\n";
3544b3e95d5SJeremy L Thompson           code << "    Grad" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp" << var_suffix
3554b3e95d5SJeremy L Thompson                << ", P_1d" << var_suffix << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix << ", r_q"
3564b3e95d5SJeremy L Thompson                << var_suffix << ");\n";
3574b3e95d5SJeremy L Thompson         }
3584b3e95d5SJeremy L Thompson         break;
3594b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
3604b3e95d5SJeremy L Thompson         CeedBasis_Cuda_shared *basis_data;
3614b3e95d5SJeremy L Thompson 
3624b3e95d5SJeremy L Thompson         code << "    CeedScalar r_q" << var_suffix << "[" << Q_name << "];\n";
3634b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
3644b3e95d5SJeremy L Thompson         data->W = basis_data->d_q_weight_1d;
3654b3e95d5SJeremy L Thompson         code << "    Weight" << (dim > 1 ? "Tensor" : "") << dim << "d<" << Q_name << ">(data, W, r_q" << var_suffix << ");\n";
3664b3e95d5SJeremy L Thompson         break;
3674b3e95d5SJeremy L Thompson       }
3684b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
3694b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
3704b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
3714b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
3724b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
3734b3e95d5SJeremy L Thompson     }
3744b3e95d5SJeremy L Thompson   } else {
3754b3e95d5SJeremy L Thompson     switch (eval_mode) {
3764b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
3774b3e95d5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n";
3784b3e95d5SJeremy L Thompson         break;  // No action
3794b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
380e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
3814b3e95d5SJeremy L Thompson         code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name
3824b3e95d5SJeremy L Thompson              << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
3834b3e95d5SJeremy L Thompson         break;
3844b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
385e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
3864b3e95d5SJeremy L Thompson         if (use_3d_slices) {
3874b3e95d5SJeremy L Thompson           code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name
3884b3e95d5SJeremy L Thompson                << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
3894b3e95d5SJeremy L Thompson         } else {
3904b3e95d5SJeremy L Thompson           code << "    GradTranspose" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp"
3914b3e95d5SJeremy L Thompson                << var_suffix << ", " << P_name << "," << Q_name << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix
3924b3e95d5SJeremy L Thompson                << ", r_e" << var_suffix << ");\n";
3934b3e95d5SJeremy L Thompson         }
3944b3e95d5SJeremy L Thompson         break;
3954b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
3964b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT:
3974b3e95d5SJeremy L Thompson         break;  // Should not occur
3984b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
3994b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
4004b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
4014b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
4024b3e95d5SJeremy L Thompson     }
4034b3e95d5SJeremy L Thompson   }
4044b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4054b3e95d5SJeremy L Thompson }
4064b3e95d5SJeremy L Thompson 
4074b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4084b3e95d5SJeremy L Thompson // QFunction
4094b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4104b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelQFunction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt dim, CeedInt num_input_fields,
4114b3e95d5SJeremy L Thompson                                                      CeedOperatorField *op_input_fields, CeedQFunctionField *qf_input_fields,
4124b3e95d5SJeremy L Thompson                                                      CeedInt num_output_fields, CeedOperatorField *op_output_fields,
4134b3e95d5SJeremy L Thompson                                                      CeedQFunctionField *qf_output_fields, std::string qfunction_name, CeedInt Q_1d,
4144b3e95d5SJeremy L Thompson                                                      bool use_3d_slices) {
4154b3e95d5SJeremy L Thompson   std::string         Q_name    = "Q_1d";
4164b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
4174b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
4184b3e95d5SJeremy L Thompson 
4194b3e95d5SJeremy L Thompson   // Setup output arays
4204b3e95d5SJeremy L Thompson   code << "\n    // -- Output field setup\n";
4214b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
4224b3e95d5SJeremy L Thompson     std::string var_suffix = "_out_" + std::to_string(i);
4234b3e95d5SJeremy L Thompson 
4244b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
4254b3e95d5SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
4264b3e95d5SJeremy L Thompson     if (eval_mode == CEED_EVAL_NONE || eval_mode == CEED_EVAL_INTERP) {
4274b3e95d5SJeremy L Thompson       code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
4284b3e95d5SJeremy L Thompson     }
4294b3e95d5SJeremy L Thompson     if (eval_mode == CEED_EVAL_GRAD) {
4304b3e95d5SJeremy L Thompson       if (use_3d_slices) {
4314b3e95d5SJeremy L Thompson         // Accumulator for gradient slices
4324b3e95d5SJeremy L Thompson         code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
4334b3e95d5SJeremy L Thompson         code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n";
4344b3e95d5SJeremy L Thompson         code << "      r_q" << var_suffix << "[i] = 0.0;\n";
4354b3e95d5SJeremy L Thompson         code << "    }\n";
4364b3e95d5SJeremy L Thompson       } else {
4374b3e95d5SJeremy L Thompson         code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << Q_name << "];\n";
4384b3e95d5SJeremy L Thompson       }
4394b3e95d5SJeremy L Thompson     }
4404b3e95d5SJeremy L Thompson   }
4414b3e95d5SJeremy L Thompson 
4424b3e95d5SJeremy L Thompson   // We treat quadrature points per slice in 3d to save registers
4434b3e95d5SJeremy L Thompson   if (use_3d_slices) {
4444b3e95d5SJeremy L Thompson     code << "\n    // Note: Using planes of 3D elements\n";
4454b3e95d5SJeremy L Thompson     code << "    #pragma unroll\n";
4464b3e95d5SJeremy L Thompson     code << "    for (CeedInt q = 0; q < " << Q_name << "; q++) {\n";
4474b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
4484b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
4494b3e95d5SJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
4504b3e95d5SJeremy L Thompson 
4514b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
4524b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
4534b3e95d5SJeremy L Thompson       // Basis action
4544b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
4554b3e95d5SJeremy L Thompson       switch (eval_mode) {
4564b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
4574b3e95d5SJeremy L Thompson           bool is_strided;
4584b3e95d5SJeremy L Thompson 
4594b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
4604b3e95d5SJeremy L Thompson 
4614b3e95d5SJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
4624b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
4634b3e95d5SJeremy L Thompson           if (is_strided) {
4644b3e95d5SJeremy L Thompson             bool    has_backend_strides;
4654b3e95d5SJeremy L Thompson             CeedInt num_elem, elem_size;
4664b3e95d5SJeremy L Thompson 
4674b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
4684b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
4694b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
4704b3e95d5SJeremy L Thompson             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
4714b3e95d5SJeremy L Thompson 
4724b3e95d5SJeremy L Thompson             if (!has_backend_strides) {
4734b3e95d5SJeremy L Thompson               CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
4744b3e95d5SJeremy L Thompson             }
4754b3e95d5SJeremy L Thompson             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
4764b3e95d5SJeremy L Thompson             code << "      readSliceQuadsStrided3d<num_comp" << var_suffix << ", " << Q_name << "," << strides[0] << "," << strides[1] << ","
4774b3e95d5SJeremy L Thompson                  << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n";
4784b3e95d5SJeremy L Thompson           } else {
4794b3e95d5SJeremy L Thompson             CeedSize                  l_size = 0;
4804b3e95d5SJeremy L Thompson             CeedInt                   comp_stride;
4814b3e95d5SJeremy L Thompson             CeedElemRestriction_Cuda *rstr_data;
4824b3e95d5SJeremy L Thompson 
4834b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
4844b3e95d5SJeremy L Thompson             code << "      const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
4854b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
4864b3e95d5SJeremy L Thompson             code << "      // CompStride: " << comp_stride << "\n";
4874b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
4884b3e95d5SJeremy L Thompson             data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
4894b3e95d5SJeremy L Thompson             code << "      readSliceQuadsOffset3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix
4904b3e95d5SJeremy L Thompson                  << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
4914b3e95d5SJeremy L Thompson           }
4924b3e95d5SJeremy L Thompson           break;
4934b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
4944b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
4954b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n";
4964b3e95d5SJeremy L Thompson           code << "        r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n";
4974b3e95d5SJeremy L Thompson           code << "      }\n";
4984b3e95d5SJeremy L Thompson           break;
4994b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
5004b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
5014b3e95d5SJeremy L Thompson           code << "      gradCollo3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_q" << var_suffix << ", s_G" << var_suffix << ", r_s"
5024b3e95d5SJeremy L Thompson                << var_suffix << ");\n";
5034b3e95d5SJeremy L Thompson           break;
5044b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
5054b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
5064b3e95d5SJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n";
5074b3e95d5SJeremy L Thompson           break;  // No action
5084b3e95d5SJeremy L Thompson                   // LCOV_EXCL_START
5094b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
5104b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
5114b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
5124b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
5134b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
5144b3e95d5SJeremy L Thompson       }
5154b3e95d5SJeremy L Thompson     }
5164b3e95d5SJeremy L Thompson     code << "\n      // -- Output fields\n";
5174b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
5184b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
5194b3e95d5SJeremy L Thompson 
5204b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
5214b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
5224b3e95d5SJeremy L Thompson       // Basis action
5234b3e95d5SJeremy L Thompson       switch (eval_mode) {
5244b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
5254b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
5264b3e95d5SJeremy L Thompson           break;  // No action
5274b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
5284b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
5294b3e95d5SJeremy L Thompson           break;
5304b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
5314b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
5324b3e95d5SJeremy L Thompson           break;
5334b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
5344b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
5354b3e95d5SJeremy L Thompson           break;  // Should not occur
5364b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
5374b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
5384b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
5394b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
5404b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
5414b3e95d5SJeremy L Thompson       }
5424b3e95d5SJeremy L Thompson     }
5434b3e95d5SJeremy L Thompson   } else {
5444b3e95d5SJeremy L Thompson     code << "\n    // Note: Using full elements\n";
5454b3e95d5SJeremy L Thompson     code << "    {\n";
5464b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
5474b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
5484b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
5494b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n";
5504b3e95d5SJeremy L Thompson     }
5514b3e95d5SJeremy L Thompson     code << "      // -- Output fields\n";
5524b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
5534b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
5544b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n";
5554b3e95d5SJeremy L Thompson     }
5564b3e95d5SJeremy L Thompson   }
5574b3e95d5SJeremy L Thompson 
5584b3e95d5SJeremy L Thompson   // Input and output buffers
5594b3e95d5SJeremy L Thompson   code << "\n      // -- QFunction inputs and outputs\n";
5604b3e95d5SJeremy L Thompson   code << "      // ---- Inputs\n";
5614b3e95d5SJeremy L Thompson   code << "      CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n";
5624b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
5634b3e95d5SJeremy L Thompson     code << "      // ------ Input field " << i << "\n";
5644b3e95d5SJeremy L Thompson     code << "      inputs[" << i << "] = r_s_in_" << i << ";\n";
5654b3e95d5SJeremy L Thompson   }
5664b3e95d5SJeremy L Thompson   code << "      // ---- Outputs\n";
5674b3e95d5SJeremy L Thompson   code << "      CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n";
5684b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
5694b3e95d5SJeremy L Thompson     code << "      // ------ Output field " << i << "\n";
5704b3e95d5SJeremy L Thompson     code << "      outputs[" << i << "] = r_s_out_" << i << ";\n";
5714b3e95d5SJeremy L Thompson   }
5724b3e95d5SJeremy L Thompson 
5734b3e95d5SJeremy L Thompson   // Apply QFunction
5744b3e95d5SJeremy L Thompson   code << "\n      // -- Apply QFunction\n";
5754b3e95d5SJeremy L Thompson   code << "      " << qfunction_name << "(ctx, ";
5764b3e95d5SJeremy L Thompson   if (dim != 3 || use_3d_slices) {
5774b3e95d5SJeremy L Thompson     code << "1";
5784b3e95d5SJeremy L Thompson   } else {
5794b3e95d5SJeremy L Thompson     code << "Q_1d";
5804b3e95d5SJeremy L Thompson   }
5814b3e95d5SJeremy L Thompson   code << ", inputs, outputs);\n";
5824b3e95d5SJeremy L Thompson 
5834b3e95d5SJeremy L Thompson   // Copy or apply transpose grad, if needed
5844b3e95d5SJeremy L Thompson   if (use_3d_slices) {
5854b3e95d5SJeremy L Thompson     code << "      // -- Output fields\n";
5864b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
5874b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
5884b3e95d5SJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
5894b3e95d5SJeremy L Thompson 
5904b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
5914b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
5924b3e95d5SJeremy L Thompson       // Basis action
5934b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
5944b3e95d5SJeremy L Thompson       switch (eval_mode) {
5954b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
5964b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
5974b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
5984b3e95d5SJeremy L Thompson           code << "      }\n";
5994b3e95d5SJeremy L Thompson           break;  // No action
6004b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
6014b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
6024b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
6034b3e95d5SJeremy L Thompson           code << "      }\n";
6044b3e95d5SJeremy L Thompson           break;
6054b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
6064b3e95d5SJeremy L Thompson           code << "      gradColloTranspose3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_s" << var_suffix << ", s_G" << var_suffix
6074b3e95d5SJeremy L Thompson                << ", r_q" << var_suffix << ");\n";
6084b3e95d5SJeremy L Thompson           break;
6094b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
6104b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
6114b3e95d5SJeremy L Thompson           break;  // Should not occur
6124b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
6134b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
6144b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
6154b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
6164b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
6174b3e95d5SJeremy L Thompson       }
6184b3e95d5SJeremy L Thompson     }
6194b3e95d5SJeremy L Thompson   }
6204b3e95d5SJeremy L Thompson   code << "    }\n";
6214b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6224b3e95d5SJeremy L Thompson }
6234b3e95d5SJeremy L Thompson 
6244b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
625b2165e7aSSebastian Grimberg // Build single operator kernel
626ab213215SJeremy L Thompson //------------------------------------------------------------------------------
627eb7e6cafSJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Cuda_gen(CeedOperator op) {
6284b3e95d5SJeremy L Thompson   bool                    is_tensor = true, use_3d_slices = false;
629241a4b83SYohann   Ceed                    ceed;
6304b3e95d5SJeremy L Thompson   CeedInt                 Q_1d, num_input_fields, num_output_fields, dim = 1;
631ca735530SJeremy L Thompson   CeedQFunctionField     *qf_input_fields, *qf_output_fields;
632ca735530SJeremy L Thompson   CeedQFunction_Cuda_gen *qf_data;
633ca735530SJeremy L Thompson   CeedQFunction           qf;
634ca735530SJeremy L Thompson   CeedOperatorField      *op_input_fields, *op_output_fields;
635ca735530SJeremy L Thompson   CeedOperator_Cuda_gen  *data;
6364b3e95d5SJeremy L Thompson   std::ostringstream      code;
6374b3e95d5SJeremy L Thompson 
6384b3e95d5SJeremy L Thompson   {
6394b3e95d5SJeremy L Thompson     bool is_setup_done;
640ca735530SJeremy L Thompson 
641ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
642ca735530SJeremy L Thompson     if (is_setup_done) return CEED_ERROR_SUCCESS;
6434b3e95d5SJeremy L Thompson   }
644ca735530SJeremy L Thompson 
645ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
646ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
647ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
648ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
649ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
650ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
651241a4b83SYohann 
6524b3e95d5SJeremy L Thompson   // Get operator data
6534b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelData_Cuda_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, op_output_fields,
6544b3e95d5SJeremy L Thompson                                                        qf_output_fields, &data->max_P_1d, &Q_1d, &dim, &is_tensor, &use_3d_slices));
6554b3e95d5SJeremy L Thompson   if (dim == 0) dim = 1;
6564b3e95d5SJeremy L Thompson   data->dim = dim;
6574b3e95d5SJeremy L Thompson   if (Q_1d == 0) {
6584b3e95d5SJeremy L Thompson     CeedInt Q;
6594b3e95d5SJeremy L Thompson 
6604b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
6614b3e95d5SJeremy L Thompson     Q_1d = Q;
6624b3e95d5SJeremy L Thompson   }
6634b3e95d5SJeremy L Thompson   data->Q_1d = Q_1d;
6644b3e95d5SJeremy L Thompson 
6650b454692Sjeremylt   // Check for restriction only identity operator
6664b3e95d5SJeremy L Thompson   {
6674b3e95d5SJeremy L Thompson     bool is_identity_qf;
6684b3e95d5SJeremy L Thompson 
6692b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
6700b454692Sjeremylt     if (is_identity_qf) {
6719e201c85SYohann       CeedEvalMode eval_mode_in, eval_mode_out;
672ca735530SJeremy L Thompson 
6732b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
6742b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
6756574a04fSJeremy L Thompson       CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
6766574a04fSJeremy L Thompson                 "Backend does not implement restriction only identity operators");
6770b454692Sjeremylt     }
6784b3e95d5SJeremy L Thompson   }
6790b454692Sjeremylt 
680f1a13f77SYohann Dudouit   // Add atomicAdd function for old NVidia architectures
6814b3e95d5SJeremy L Thompson   {
6824b3e95d5SJeremy L Thompson     Ceed_Cuda            *ceed_data;
6834b3e95d5SJeremy L Thompson     struct cudaDeviceProp prop;
6844b3e95d5SJeremy L Thompson 
6852b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetData(ceed, &ceed_data));
6862b730f8bSJeremy L Thompson     CeedCallBackend(cudaGetDeviceProperties(&prop, ceed_data->device_id));
68780a9ef05SNatalie Beams     if ((prop.major < 6) && (CEED_SCALAR_TYPE != CEED_SCALAR_FP32)) {
68822070f95SJeremy L Thompson       char       *atomic_add_source;
68922070f95SJeremy L Thompson       const char *atomic_add_path;
690ca735530SJeremy L Thompson 
6912b730f8bSJeremy L Thompson       CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-atomic-add-fallback.h", &atomic_add_path));
69223d4529eSJeremy L Thompson       CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Atomic Add Source -----\n");
6932b730f8bSJeremy L Thompson       CeedCallBackend(CeedLoadSourceToBuffer(ceed, atomic_add_path, &atomic_add_source));
6949e201c85SYohann       code << atomic_add_source;
6952b730f8bSJeremy L Thompson       CeedCallBackend(CeedFree(&atomic_add_path));
6962b730f8bSJeremy L Thompson       CeedCallBackend(CeedFree(&atomic_add_source));
697f1a13f77SYohann Dudouit     }
6984b3e95d5SJeremy L Thompson   }
699f1a13f77SYohann Dudouit 
7009e201c85SYohann   // Load basis source files
7014b3e95d5SJeremy L Thompson   // TODO: Add non-tensor, AtPoints
7029e201c85SYohann   {
70322070f95SJeremy L Thompson     char       *tensor_basis_kernel_source;
70422070f95SJeremy L Thompson     const char *tensor_basis_kernel_path;
705ca735530SJeremy L Thompson 
7062b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-shared-basis-tensor-templates.h", &tensor_basis_kernel_path));
70723d4529eSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Tensor Basis Kernel Source -----\n");
7082b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, tensor_basis_kernel_path, &tensor_basis_kernel_source));
7099e201c85SYohann     code << tensor_basis_kernel_source;
7102b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&tensor_basis_kernel_path));
7112b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&tensor_basis_kernel_source));
7129e201c85SYohann   }
7139e201c85SYohann   {
71422070f95SJeremy L Thompson     char       *cuda_gen_template_source;
71522070f95SJeremy L Thompson     const char *cuda_gen_template_path;
716ca735530SJeremy L Thompson 
7172b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-gen-templates.h", &cuda_gen_template_path));
71823d4529eSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Cuda-Gen Template Source -----\n");
7192b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, cuda_gen_template_path, &cuda_gen_template_source));
7209e201c85SYohann     code << cuda_gen_template_source;
7212b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&cuda_gen_template_path));
7222b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&cuda_gen_template_source));
7239e201c85SYohann   }
724241a4b83SYohann 
7254b3e95d5SJeremy L Thompson   // Get QFunction name
7264b3e95d5SJeremy L Thompson   std::string qfunction_name(qf_data->qfunction_name);
7274b3e95d5SJeremy L Thompson   std::string operator_name;
7284b3e95d5SJeremy L Thompson 
72909095acaSJeremy L Thompson   operator_name = "CeedKernelCudaGenOperator_" + qfunction_name;
7304d537eeaSYohann 
7319e201c85SYohann   // Define CEED_Q_VLA
7329e201c85SYohann   code << "\n#undef CEED_Q_VLA\n";
7334b3e95d5SJeremy L Thompson   if (dim != 3 || use_3d_slices) {
7349e201c85SYohann     code << "#define CEED_Q_VLA 1\n\n";
7359e201c85SYohann   } else {
7369e201c85SYohann     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
7379e201c85SYohann   }
7389e201c85SYohann 
7394b3e95d5SJeremy L Thompson   // Add user QFunction source
7404b3e95d5SJeremy L Thompson   {
7414b3e95d5SJeremy L Thompson     std::string qfunction_source(qf_data->qfunction_source);
7424b3e95d5SJeremy L Thompson 
74309095acaSJeremy L Thompson     code << qfunction_source;
7444b3e95d5SJeremy L Thompson   }
745241a4b83SYohann 
746241a4b83SYohann   // Setup
747818e0025SJeremy L Thompson   code << "\n// -----------------------------------------------------------------------------\n";
7484b3e95d5SJeremy L Thompson   code << "// Operator Kernel\n";
7494b3e95d5SJeremy L Thompson   code << "// \n";
7504b3e95d5SJeremy L Thompson   code << "// d_[in,out]_i:   CeedVector device array\n";
7514b3e95d5SJeremy L Thompson   code << "// r_[in,out]_e_i: Element vector register\n";
7524b3e95d5SJeremy L Thompson   code << "// r_[in,out]_q_i: Quadrature space vector register\n";
7534b3e95d5SJeremy L Thompson   code << "// r_[in,out]_s_i: Quadrature space slice  vector register\n";
7544b3e95d5SJeremy L Thompson   code << "// \n";
7554b3e95d5SJeremy L Thompson   code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n";
7564b3e95d5SJeremy L Thompson   code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n";
7574b3e95d5SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n";
7584b3e95d5SJeremy L Thompson   code << "extern \"C\" __global__ void " << operator_name
7592b730f8bSJeremy L Thompson        << "(CeedInt num_elem, void* ctx, FieldsInt_Cuda indices, Fields_Cuda fields, Fields_Cuda B, Fields_Cuda G, CeedScalar *W) {\n";
7604b3e95d5SJeremy L Thompson 
7614b3e95d5SJeremy L Thompson   // Scratch buffers
7629e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
7634b3e95d5SJeremy L Thompson     CeedEvalMode eval_mode;
7644b3e95d5SJeremy L Thompson 
7652b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
7669e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
7674b3e95d5SJeremy L Thompson       code << "  const CeedScalar *d_in_" << i << " = fields.inputs[" << i << "];\n";
768241a4b83SYohann     }
769241a4b83SYohann   }
7709e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
7714b3e95d5SJeremy L Thompson     code << "  CeedScalar *d_out_" << i << " = fields.outputs[" << i << "];\n";
772241a4b83SYohann   }
7731da99368SJeremy L Thompson 
7749e201c85SYohann   code << "  const CeedInt dim = " << dim << ";\n";
7759e201c85SYohann   code << "  const CeedInt Q_1d = " << Q_1d << ";\n";
7761da99368SJeremy L Thompson 
7774b3e95d5SJeremy L Thompson   // Shared data
778241a4b83SYohann   code << "  extern __shared__ CeedScalar slice[];\n";
7799e201c85SYohann   code << "  SharedData_Cuda data;\n";
7809e201c85SYohann   code << "  data.t_id_x = threadIdx.x;\n";
7819e201c85SYohann   code << "  data.t_id_y = threadIdx.y;\n";
7829e201c85SYohann   code << "  data.t_id_z = threadIdx.z;\n";
7839e201c85SYohann   code << "  data.t_id  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
7849e201c85SYohann   code << "  data.slice = slice + data.t_id_z*T_1D" << (dim > 1 ? "*T_1D" : "") << ";\n";
785920dcdc4Sjeremylt 
786ac421f39SYohann   // Initialize constants, and matrices B and G
7874b3e95d5SJeremy L Thompson   code << "\n  // Input field constants and basis data\n";
7889e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
7895a5594ffSJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, i, op_input_fields[i], qf_input_fields[i], Q_1d, true, use_3d_slices));
790920dcdc4Sjeremylt   }
7914b3e95d5SJeremy L Thompson   code << "\n  // Output field constants and basis data\n";
7929e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
7935a5594ffSJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, i, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices));
7944b3e95d5SJeremy L Thompson   }
795920dcdc4Sjeremylt 
7964b3e95d5SJeremy L Thompson   // Loop over all elements
7974b3e95d5SJeremy L Thompson   code << "\n  // Element loop\n";
798ac421f39SYohann   code << "  __syncthreads();\n";
7999e201c85SYohann   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
8004b3e95d5SJeremy L Thompson 
801e93651e5SJeremy L Thompson   // -- Compute minimum buffer space needed
802e93651e5SJeremy L Thompson   CeedInt max_rstr_buffer_size = 0;
803e93651e5SJeremy L Thompson 
8049e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
805e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
806e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
807e93651e5SJeremy L Thompson 
808e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
809e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
810e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
811e93651e5SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * elem_size);
812e93651e5SJeremy L Thompson   }
813e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
814e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
815e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
816e93651e5SJeremy L Thompson 
817e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
818e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
819e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
820e93651e5SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * elem_size);
821e93651e5SJeremy L Thompson   }
822e93651e5SJeremy L Thompson   code << "    // Scratch restriction buffer space\n";
823e93651e5SJeremy L Thompson   code << "    CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n";
824e93651e5SJeremy L Thompson 
825e93651e5SJeremy L Thompson   // -- Determine best input field processing order
826e93651e5SJeremy L Thompson   CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX];
827e93651e5SJeremy L Thompson 
828e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
829e93651e5SJeremy L Thompson     field_rstr_in_buffer[i] = -1;
830e93651e5SJeremy L Thompson     input_field_order[i]    = -1;
831e93651e5SJeremy L Thompson   }
832e93651e5SJeremy L Thompson   {
833e93651e5SJeremy L Thompson     bool    is_ordered[CEED_FIELD_MAX];
834e93651e5SJeremy L Thompson     CeedInt curr_index = 0;
835e93651e5SJeremy L Thompson 
836e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
837e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
838e93651e5SJeremy L Thompson       CeedVector          vec_i;
839e93651e5SJeremy L Thompson       CeedElemRestriction rstr_i;
840e93651e5SJeremy L Thompson 
841e93651e5SJeremy L Thompson       if (is_ordered[i]) continue;
842e93651e5SJeremy L Thompson       field_rstr_in_buffer[i]       = i;
843e93651e5SJeremy L Thompson       is_ordered[i]                 = true;
844e93651e5SJeremy L Thompson       input_field_order[curr_index] = i;
845e93651e5SJeremy L Thompson       curr_index++;
846*034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
847e93651e5SJeremy L Thompson       if (vec_i == CEED_VECTOR_NONE) continue;  // CEED_EVAL_WEIGHT
848e93651e5SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
849e93651e5SJeremy L Thompson       for (CeedInt j = i + 1; j < num_input_fields; j++) {
850e93651e5SJeremy L Thompson         CeedVector          vec_j;
851e93651e5SJeremy L Thompson         CeedElemRestriction rstr_j;
852e93651e5SJeremy L Thompson 
853e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
854e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
855e93651e5SJeremy L Thompson         if (rstr_i == rstr_j && vec_i == vec_j) {
856e93651e5SJeremy L Thompson           field_rstr_in_buffer[j]       = i;
857e93651e5SJeremy L Thompson           is_ordered[j]                 = true;
858e93651e5SJeremy L Thompson           input_field_order[curr_index] = j;
859e93651e5SJeremy L Thompson           curr_index++;
860e93651e5SJeremy L Thompson         }
861e93651e5SJeremy L Thompson       }
862e93651e5SJeremy L Thompson     }
863e93651e5SJeremy L Thompson   }
864e93651e5SJeremy L Thompson 
865e93651e5SJeremy L Thompson   // -- Input restriction and basis
866e93651e5SJeremy L Thompson   code << "\n    // -- Input field restrictions and basis actions\n";
867e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
868e93651e5SJeremy L Thompson     CeedInt f = input_field_order[i];
869e93651e5SJeremy L Thompson 
870e93651e5SJeremy L Thompson     code << "    // ---- Input field " << f << "\n";
871920dcdc4Sjeremylt 
8724b3e95d5SJeremy L Thompson     // ---- Restriction
873e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, f, dim, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f],
874e93651e5SJeremy L Thompson                                                                 Q_1d, true, use_3d_slices));
8759a2291e3SJeremy L Thompson 
8764b3e95d5SJeremy L Thompson     // ---- Basis action
877e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, f, dim, op_input_fields[f], qf_input_fields[f], Q_1d, true, use_3d_slices));
878920dcdc4Sjeremylt   }
879920dcdc4Sjeremylt 
8804b3e95d5SJeremy L Thompson   // -- Q function
8814b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelQFunction_Cuda_gen(code, data, dim, num_input_fields, op_input_fields, qf_input_fields, num_output_fields,
8824b3e95d5SJeremy L Thompson                                                             op_output_fields, qf_output_fields, qfunction_name, Q_1d, use_3d_slices));
883ca735530SJeremy L Thompson 
8844b3e95d5SJeremy L Thompson   // -- Output basis and restriction
8854b3e95d5SJeremy L Thompson   code << "\n    // -- Output field basis action and restrictions\n";
8869e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
8874b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
888ca735530SJeremy L Thompson 
8894b3e95d5SJeremy L Thompson     // ---- Basis action
8904b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, i, dim, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices));
8919a2291e3SJeremy L Thompson 
8924b3e95d5SJeremy L Thompson     // ---- Restriction
8934b3e95d5SJeremy L Thompson     CeedCallBackend(
894e93651e5SJeremy L Thompson         CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, i, dim, NULL, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices));
895ac421f39SYohann   }
896241a4b83SYohann 
8974b3e95d5SJeremy L Thompson   // Close loop and function
898241a4b83SYohann   code << "  }\n";
899818e0025SJeremy L Thompson   code << "}\n";
900818e0025SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n\n";
901241a4b83SYohann 
902ab213215SJeremy L Thompson   // View kernel for debugging
90323d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "Generated Operator Kernels:\n");
9043f21f6b1SJeremy L Thompson   CeedDebug(ceed, code.str().c_str());
905241a4b83SYohann 
906eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Cuda(ceed, code.str().c_str(), &data->module, 1, "T_1D", CeedIntMax(Q_1d, data->max_P_1d)));
907eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, operator_name.c_str(), &data->op));
908241a4b83SYohann 
9092b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
910e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
911241a4b83SYohann }
9122a86cc9dSSebastian Grimberg 
913ab213215SJeremy L Thompson //------------------------------------------------------------------------------
914