xref: /libCEED/backends/cuda-gen/ceed-cuda-gen-operator-build.cpp (revision 412e56839ae1393ae29619b6fd39ea10cdd89adf)
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 
2545a787f7SJeremy L Thompson struct FieldReuse_Cuda {
2645a787f7SJeremy L Thompson   CeedInt      index;
2745a787f7SJeremy L Thompson   bool         is_input;
2845a787f7SJeremy L Thompson   CeedEvalMode eval_mode;
2945a787f7SJeremy L Thompson };
3045a787f7SJeremy L Thompson 
31ab213215SJeremy L Thompson //------------------------------------------------------------------------------
324b3e95d5SJeremy L Thompson // Determine type of operator
334b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
344b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelData_Cuda_gen(Ceed ceed, CeedInt num_input_fields, CeedOperatorField *op_input_fields,
354b3e95d5SJeremy L Thompson                                                 CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, CeedOperatorField *op_output_fields,
36*412e5683SJeremy L Thompson                                                 CeedQFunctionField *qf_output_fields, CeedInt *max_P, CeedInt *max_P_1d, CeedInt *Q, CeedInt *Q_1d,
37*412e5683SJeremy L Thompson                                                 CeedInt *max_dim, bool *is_all_tensor, bool *use_3d_slices) {
38*412e5683SJeremy L Thompson   // Check if all are tensor
39*412e5683SJeremy L Thompson   *is_all_tensor = true;
404b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
414b3e95d5SJeremy L Thompson     CeedBasis basis;
424b3e95d5SJeremy L Thompson 
434b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
444b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
454b3e95d5SJeremy L Thompson       bool is_field_tensor;
464b3e95d5SJeremy L Thompson 
474b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
48*412e5683SJeremy L Thompson       *is_all_tensor = *is_all_tensor && is_field_tensor;
494b3e95d5SJeremy L Thompson     }
50681d0ea7SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis));
514b3e95d5SJeremy L Thompson   }
524b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
534b3e95d5SJeremy L Thompson     CeedBasis basis;
544b3e95d5SJeremy L Thompson 
554b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
564b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
574b3e95d5SJeremy L Thompson       bool is_field_tensor;
584b3e95d5SJeremy L Thompson 
594b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
60*412e5683SJeremy L Thompson       *is_all_tensor = *is_all_tensor && is_field_tensor;
61*412e5683SJeremy L Thompson     }
62*412e5683SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis));
63*412e5683SJeremy L Thompson   }
64*412e5683SJeremy L Thompson 
65*412e5683SJeremy L Thompson   // Find max_P, max_P_1d, Q, and Q_1d
66*412e5683SJeremy L Thompson   bool is_all_3d = true;
67*412e5683SJeremy L Thompson 
68*412e5683SJeremy L Thompson   *max_P    = 0;
69*412e5683SJeremy L Thompson   *max_P_1d = 0;
70*412e5683SJeremy L Thompson   *Q        = 0;
71*412e5683SJeremy L Thompson   *Q_1d     = 0;
72*412e5683SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
73*412e5683SJeremy L Thompson     CeedBasis basis;
74*412e5683SJeremy L Thompson 
75*412e5683SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
76*412e5683SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
77*412e5683SJeremy L Thompson       bool    is_field_tensor;
78*412e5683SJeremy L Thompson       CeedInt field_dim = 0, field_P = 0, field_P_1d = 0, field_Q = 0, field_Q_1d = 0;
79*412e5683SJeremy L Thompson 
80*412e5683SJeremy L Thompson       // Check if 3D
814b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
82*412e5683SJeremy L Thompson       is_all_3d = is_all_3d && (field_dim == 3);
83*412e5683SJeremy L Thompson       *max_dim  = CeedIntMax(*max_dim, field_dim);
84*412e5683SJeremy L Thompson 
85*412e5683SJeremy L Thompson       // Collect P, P_1d, Q, and Q_1d
86*412e5683SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P));
87*412e5683SJeremy L Thompson       *max_P = CeedIntMax(*max_P, field_P);
88*412e5683SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
89*412e5683SJeremy L Thompson       if (is_field_tensor) {
90*412e5683SJeremy L Thompson         CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
91*412e5683SJeremy L Thompson         *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
92*412e5683SJeremy L Thompson       }
93*412e5683SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q));
94*412e5683SJeremy L Thompson       CeedCheck(*Q == 0 || field_Q == *Q, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
95*412e5683SJeremy L Thompson       *Q = field_Q;
96*412e5683SJeremy L Thompson       if (is_field_tensor) {
97*412e5683SJeremy L Thompson         CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
984b3e95d5SJeremy L Thompson         CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
994b3e95d5SJeremy L Thompson         *Q_1d = field_Q_1d;
1004b3e95d5SJeremy L Thompson       }
101*412e5683SJeremy L Thompson     }
102*412e5683SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis));
103*412e5683SJeremy L Thompson   }
104*412e5683SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
105*412e5683SJeremy L Thompson     CeedBasis basis;
106*412e5683SJeremy L Thompson 
107*412e5683SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
108*412e5683SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
109*412e5683SJeremy L Thompson       bool    is_field_tensor;
110*412e5683SJeremy L Thompson       CeedInt field_dim = 0, field_P = 0, field_P_1d = 0, field_Q = 0, field_Q_1d = 0;
111*412e5683SJeremy L Thompson 
112*412e5683SJeremy L Thompson       // Check if 3D
113*412e5683SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
114*412e5683SJeremy L Thompson       is_all_3d = is_all_3d && (field_dim == 3);
115*412e5683SJeremy L Thompson       *max_dim  = CeedIntMax(*max_dim, field_dim);
116*412e5683SJeremy L Thompson 
117*412e5683SJeremy L Thompson       // Collect P, P_1d, Q, and Q_1d
118*412e5683SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P));
119*412e5683SJeremy L Thompson       *max_P = CeedIntMax(*max_P, field_P);
120*412e5683SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
121*412e5683SJeremy L Thompson       if (is_field_tensor) {
122*412e5683SJeremy L Thompson         CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
123*412e5683SJeremy L Thompson         *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
124*412e5683SJeremy L Thompson       }
125*412e5683SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q));
126*412e5683SJeremy L Thompson       CeedCheck(*Q == 0 || field_Q == *Q, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
127*412e5683SJeremy L Thompson       *Q = field_Q;
128*412e5683SJeremy L Thompson       if (is_field_tensor) {
129*412e5683SJeremy L Thompson         CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
130*412e5683SJeremy L Thompson         CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
131*412e5683SJeremy L Thompson         *Q_1d = field_Q_1d;
132*412e5683SJeremy L Thompson       }
133*412e5683SJeremy L Thompson     }
134681d0ea7SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis));
1354b3e95d5SJeremy L Thompson   }
1364b3e95d5SJeremy L Thompson 
1374b3e95d5SJeremy L Thompson   // Only use 3D collocated gradient parallelization strategy when gradient is computed
1384b3e95d5SJeremy L Thompson   *use_3d_slices = false;
139*412e5683SJeremy L Thompson   if (is_all_3d && *is_all_tensor) {
1404b3e95d5SJeremy L Thompson     bool was_grad_found = false;
1414b3e95d5SJeremy L Thompson 
1424b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
1434b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
1444b3e95d5SJeremy L Thompson 
1454b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1464b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
1474b3e95d5SJeremy L Thompson         CeedBasis_Cuda_shared *basis_data;
1484b3e95d5SJeremy L Thompson         CeedBasis              basis;
1494b3e95d5SJeremy L Thompson 
1504b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1514b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1524b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
1534b3e95d5SJeremy L Thompson         was_grad_found = true;
154681d0ea7SJeremy L Thompson         CeedCallBackend(CeedBasisDestroy(&basis));
1554b3e95d5SJeremy L Thompson       }
1564b3e95d5SJeremy L Thompson     }
1574b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
1584b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
1594b3e95d5SJeremy L Thompson 
1604b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1614b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
1624b3e95d5SJeremy L Thompson         CeedBasis_Cuda_shared *basis_data;
1634b3e95d5SJeremy L Thompson         CeedBasis              basis;
1644b3e95d5SJeremy L Thompson 
1654b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1664b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1674b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
1684b3e95d5SJeremy L Thompson         was_grad_found = true;
169681d0ea7SJeremy L Thompson         CeedCallBackend(CeedBasisDestroy(&basis));
1704b3e95d5SJeremy L Thompson       }
1714b3e95d5SJeremy L Thompson     }
1724b3e95d5SJeremy L Thompson   }
1734b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1744b3e95d5SJeremy L Thompson }
1754b3e95d5SJeremy L Thompson 
1764b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1774b3e95d5SJeremy L Thompson // Setup fields
1784b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1794b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelFieldData_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedOperatorField op_field,
18045a787f7SJeremy L Thompson                                                      CeedQFunctionField qf_field, FieldReuse_Cuda field_reuse, CeedInt Q_1d, bool is_input,
181*412e5683SJeremy L Thompson                                                      bool is_all_tensor, bool is_at_points, bool use_3d_slices) {
182*412e5683SJeremy L Thompson   bool      is_tensor = true;
183*412e5683SJeremy L Thompson   CeedBasis basis;
184*412e5683SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
185*412e5683SJeremy L Thompson   if (basis != CEED_BASIS_NONE) CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
186*412e5683SJeremy L Thompson 
18759fa3f92SJeremy L Thompson   const char            *field_name;
1884b3e95d5SJeremy L Thompson   std::string            var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
189dc007f05SJeremy L Thompson   std::string            P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q";
1904b3e95d5SJeremy L Thompson   std::string            option_name = (is_input ? "inputs" : "outputs");
1914b3e95d5SJeremy L Thompson   CeedEvalMode           eval_mode   = CEED_EVAL_NONE;
192*412e5683SJeremy L Thompson   CeedInt                elem_size = 0, num_comp = 0, dim = 1, P_1d = 0;
1934b3e95d5SJeremy L Thompson   CeedElemRestriction    elem_rstr;
1944b3e95d5SJeremy L Thompson   CeedBasis_Cuda_shared *basis_data;
1954b3e95d5SJeremy L Thompson 
1960a2a6492SJeremy L Thompson   // Field reuse info
19745a787f7SJeremy L Thompson   bool use_previous_field = field_reuse.index != -1;
1980a2a6492SJeremy L Thompson 
19959fa3f92SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetName(op_field, &field_name));
20059fa3f92SJeremy L Thompson   code << "  // -- " << (is_input ? "Input" : "Output") << " field " << i << ": " << field_name << "\n";
2014b3e95d5SJeremy L Thompson 
2024b3e95d5SJeremy L Thompson   // Get field data
2034b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
2044b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
2054b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
2064b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
2074b3e95d5SJeremy L Thompson   }
208681d0ea7SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
2094b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
2104b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetData(basis, &basis_data));
211*412e5683SJeremy L Thompson     CeedCallBackend(CeedBasisGetDimension(basis, &dim));
212dc007f05SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
213dc007f05SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
2144b3e95d5SJeremy L Thompson   }
2154b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
2164b3e95d5SJeremy L Thompson 
2174b3e95d5SJeremy L Thompson   // Set field constants
218*412e5683SJeremy L Thompson   code << "  const CeedInt dim" << var_suffix << " = " << dim << ";\n";
219*412e5683SJeremy L Thompson   if (is_tensor && !is_all_tensor) {
220*412e5683SJeremy L Thompson     CeedInt P = 0;
221*412e5683SJeremy L Thompson 
222*412e5683SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes(basis, &P));
223*412e5683SJeremy L Thompson     code << "  const CeedInt P" << var_suffix << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P) << ";\n";
224*412e5683SJeremy L Thompson   }
2254b3e95d5SJeremy L Thompson   code << "  const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n";
226343e3094SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
2274b3e95d5SJeremy L Thompson     code << "  const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n";
2284b3e95d5SJeremy L Thompson   }
2294b3e95d5SJeremy L Thompson 
2304b3e95d5SJeremy L Thompson   // Load basis data
2314b3e95d5SJeremy L Thompson   code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
2324b3e95d5SJeremy L Thompson   switch (eval_mode) {
2334b3e95d5SJeremy L Thompson     case CEED_EVAL_NONE:
2344b3e95d5SJeremy L Thompson       break;
2354b3e95d5SJeremy L Thompson     case CEED_EVAL_INTERP:
2368b97b69aSJeremy L Thompson       if (is_at_points) {
2378b97b69aSJeremy L Thompson         // AtPoints
2388b97b69aSJeremy L Thompson         if (!basis_data->d_chebyshev_interp_1d) {
2398b97b69aSJeremy L Thompson           CeedSize    interp_bytes;
2408b97b69aSJeremy L Thompson           CeedScalar *chebyshev_interp_1d;
2418b97b69aSJeremy L Thompson 
2428b97b69aSJeremy L Thompson           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
2438b97b69aSJeremy L Thompson           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
2448b97b69aSJeremy L Thompson           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
2458b97b69aSJeremy L Thompson           CeedCallCuda(CeedBasisReturnCeed(basis), cudaMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
2468b97b69aSJeremy L Thompson           CeedCallCuda(CeedBasisReturnCeed(basis),
2478b97b69aSJeremy L Thompson                        cudaMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
2488b97b69aSJeremy L Thompson           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
2498b97b69aSJeremy L Thompson         }
2508b97b69aSJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
2518b97b69aSJeremy L Thompson         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
2528b97b69aSJeremy L Thompson       } else {
2538b97b69aSJeremy L Thompson         // Standard quadrature
2544b3e95d5SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
2554b3e95d5SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_interp_1d;
2568b97b69aSJeremy L Thompson       }
2570a2a6492SJeremy L Thompson       if (use_previous_field) {
25845a787f7SJeremy L Thompson         std::string reuse_var = "s_B" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
2590a2a6492SJeremy L Thompson 
2600a2a6492SJeremy L Thompson         code << "  CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n";
2610a2a6492SJeremy L Thompson       } else {
262dc007f05SJeremy L Thompson         code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
263f815fac9SJeremy L Thompson         code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
2640a2a6492SJeremy L Thompson       }
2654b3e95d5SJeremy L Thompson       break;
2664b3e95d5SJeremy L Thompson     case CEED_EVAL_GRAD:
2678b97b69aSJeremy L Thompson       if (is_at_points) {
2688b97b69aSJeremy L Thompson         // AtPoints
2698b97b69aSJeremy L Thompson         if (!basis_data->d_chebyshev_interp_1d) {
2708b97b69aSJeremy L Thompson           CeedSize    interp_bytes;
2718b97b69aSJeremy L Thompson           CeedScalar *chebyshev_interp_1d;
2728b97b69aSJeremy L Thompson 
2738b97b69aSJeremy L Thompson           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
2748b97b69aSJeremy L Thompson           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
2758b97b69aSJeremy L Thompson           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
2768b97b69aSJeremy L Thompson           CeedCallCuda(CeedBasisReturnCeed(basis), cudaMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
2778b97b69aSJeremy L Thompson           CeedCallCuda(CeedBasisReturnCeed(basis),
2788b97b69aSJeremy L Thompson                        cudaMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
2798b97b69aSJeremy L Thompson           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
2808b97b69aSJeremy L Thompson         }
2818b97b69aSJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
2828b97b69aSJeremy L Thompson         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
2838b97b69aSJeremy L Thompson       } else {
2848b97b69aSJeremy L Thompson         // Standard quadrature
2854b3e95d5SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
2864b3e95d5SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_interp_1d;
2878b97b69aSJeremy L Thompson       }
288dc007f05SJeremy L Thompson       if (is_tensor) {
2890a2a6492SJeremy L Thompson         if (use_previous_field) {
29045a787f7SJeremy L Thompson           std::string reuse_var = "s_B" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
2910a2a6492SJeremy L Thompson 
2920a2a6492SJeremy L Thompson           code << "  CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n";
2930a2a6492SJeremy L Thompson         } else {
294dc007f05SJeremy L Thompson           code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
295f815fac9SJeremy L Thompson           code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
296dc007f05SJeremy L Thompson         }
2970a2a6492SJeremy L Thompson       }
2988b97b69aSJeremy L Thompson       if (is_at_points) break;  // No G mat for AtPoints
2994b3e95d5SJeremy L Thompson       if (use_3d_slices) {
3004b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d;
3014b3e95d5SJeremy L Thompson         else data->G.outputs[i] = basis_data->d_collo_grad_1d;
30245a787f7SJeremy L Thompson         if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD) {
30345a787f7SJeremy L Thompson           std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
3040a2a6492SJeremy L Thompson 
3050a2a6492SJeremy L Thompson           code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
3060a2a6492SJeremy L Thompson         } else {
307dc007f05SJeremy L Thompson           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
308f815fac9SJeremy L Thompson           code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
3090a2a6492SJeremy L Thompson         }
3104b3e95d5SJeremy L Thompson       } else {
3114b3e95d5SJeremy L Thompson         bool has_collo_grad = basis_data->d_collo_grad_1d;
3124b3e95d5SJeremy L Thompson 
3134b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
3144b3e95d5SJeremy L Thompson         else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
3154b3e95d5SJeremy L Thompson         if (has_collo_grad) {
31645a787f7SJeremy L Thompson           if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD) {
31745a787f7SJeremy L Thompson             std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
3180a2a6492SJeremy L Thompson 
3190a2a6492SJeremy L Thompson             code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
3200a2a6492SJeremy L Thompson           } else {
321dc007f05SJeremy L Thompson             code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
322f815fac9SJeremy L Thompson             code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
3230a2a6492SJeremy L Thompson           }
3240a2a6492SJeremy L Thompson         } else {
32545a787f7SJeremy L Thompson           if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD) {
32645a787f7SJeremy L Thompson             std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
3270a2a6492SJeremy L Thompson 
3280a2a6492SJeremy L Thompson             code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
3294b3e95d5SJeremy L Thompson           } else {
330*412e5683SJeremy L Thompson             code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << P_name << "*" << Q_name << (is_tensor ? "" : "*dim")
331*412e5683SJeremy L Thompson                  << (is_tensor ? "" : var_suffix) << "];\n";
332*412e5683SJeremy L Thompson             code << "  LoadMatrix<" << P_name << ", " << Q_name << (is_tensor ? "" : "*dim") << (is_tensor ? "" : var_suffix) << ">(data, G."
333*412e5683SJeremy L Thompson                  << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
3344b3e95d5SJeremy L Thompson           }
3354b3e95d5SJeremy L Thompson         }
3360a2a6492SJeremy L Thompson       }
3374b3e95d5SJeremy L Thompson       break;
3384b3e95d5SJeremy L Thompson     case CEED_EVAL_WEIGHT:
3394b3e95d5SJeremy L Thompson       break;  // No action
3404b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
3414b3e95d5SJeremy L Thompson     case CEED_EVAL_DIV:
3424b3e95d5SJeremy L Thompson     case CEED_EVAL_CURL:
3434b3e95d5SJeremy L Thompson       break;  // TODO: Not implemented
3444b3e95d5SJeremy L Thompson               // LCOV_EXCL_STOP
3454b3e95d5SJeremy L Thompson   }
3468b97b69aSJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
3474b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3484b3e95d5SJeremy L Thompson }
3494b3e95d5SJeremy L Thompson 
3504b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
3514b3e95d5SJeremy L Thompson // Restriction
3524b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
353*412e5683SJeremy L Thompson static int CeedOperatorBuildKernelRestriction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedInt max_dim,
354e93651e5SJeremy L Thompson                                                        CeedInt field_input_buffer[], CeedOperatorField op_field, CeedQFunctionField qf_field,
355*412e5683SJeremy L Thompson                                                        CeedInt Q_1d, bool is_input, bool is_all_tensor, bool is_at_points, bool use_3d_slices) {
3564b3e95d5SJeremy L Thompson   std::string               var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
357*412e5683SJeremy L Thompson   std::string               P_name     = (is_all_tensor ? "P_1d" : "P") + var_suffix;
3584b3e95d5SJeremy L Thompson   CeedEvalMode              eval_mode  = CEED_EVAL_NONE;
359*412e5683SJeremy L Thompson   CeedInt                   elem_size = 0, num_comp = 0;
3604b3e95d5SJeremy L Thompson   CeedSize                  l_size;
361f815fac9SJeremy L Thompson   CeedRestrictionType       rstr_type = CEED_RESTRICTION_STANDARD;
3624b3e95d5SJeremy L Thompson   CeedElemRestriction_Cuda *rstr_data;
3634b3e95d5SJeremy L Thompson   CeedElemRestriction       elem_rstr;
3644b3e95d5SJeremy L Thompson 
3654b3e95d5SJeremy L Thompson   // Get field data
3664b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
3674b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
368f815fac9SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
3694b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
3704b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
3714b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
3724b3e95d5SJeremy L Thompson   }
3734b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
3744b3e95d5SJeremy L Thompson 
3754b3e95d5SJeremy L Thompson   // Restriction
3764b3e95d5SJeremy L Thompson   if (is_input) {
3774b3e95d5SJeremy L Thompson     // Input
378e93651e5SJeremy L Thompson     if (field_input_buffer[i] != i) {
379e93651e5SJeremy L Thompson       std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]);
380e93651e5SJeremy L Thompson 
381e93651e5SJeremy L Thompson       // Restriction was already done for previous input
382e93651e5SJeremy L Thompson       code << "    CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n";
3838b97b69aSJeremy L Thompson     } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices && is_at_points)) {
3848b97b69aSJeremy L Thompson       if (eval_mode == CEED_EVAL_NONE && rstr_type != CEED_RESTRICTION_POINTS) {
385e93651e5SJeremy L Thompson         // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated
3864b3e95d5SJeremy L Thompson         code << "    CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n";
3878b97b69aSJeremy L Thompson       } else if (rstr_type != CEED_RESTRICTION_POINTS) {
388e93651e5SJeremy L Thompson         // Otherwise we're using the scratch space
389e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
390e93651e5SJeremy L Thompson       }
391f815fac9SJeremy L Thompson       switch (rstr_type) {
392f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STANDARD: {
3934b3e95d5SJeremy L Thompson           CeedInt comp_stride;
3944b3e95d5SJeremy L Thompson 
3954b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
3964b3e95d5SJeremy L Thompson           code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
3974b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3984b3e95d5SJeremy L Thompson           code << "    // CompStride: " << comp_stride << "\n";
3994b3e95d5SJeremy L Thompson           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
400*412e5683SJeremy L Thompson           code << "    ReadLVecStandard" << (is_all_tensor ? max_dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name
401dc007f05SJeremy L Thompson                << ">(data, l_size" << var_suffix << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n";
402f815fac9SJeremy L Thompson           break;
403f815fac9SJeremy L Thompson         }
404f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STRIDED: {
4054b3e95d5SJeremy L Thompson           bool    has_backend_strides;
4064b3e95d5SJeremy L Thompson           CeedInt num_elem;
4074b3e95d5SJeremy L Thompson 
4084b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
4094b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
4104b3e95d5SJeremy L Thompson           CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
4114b3e95d5SJeremy L Thompson 
4124b3e95d5SJeremy L Thompson           if (!has_backend_strides) {
4134b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
4144b3e95d5SJeremy L Thompson           }
4154b3e95d5SJeremy L Thompson           code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
416*412e5683SJeremy L Thompson           code << "    ReadLVecStrided" << (is_all_tensor ? max_dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", "
417dc007f05SJeremy L Thompson                << strides[1] << ", " << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n";
418f815fac9SJeremy L Thompson           break;
419f815fac9SJeremy L Thompson         }
4208b97b69aSJeremy L Thompson         case CEED_RESTRICTION_POINTS: {
4218b97b69aSJeremy L Thompson           CeedInt comp_stride;
4228b97b69aSJeremy L Thompson 
4238b97b69aSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
4248b97b69aSJeremy L Thompson           code << "    const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
4258b97b69aSJeremy L Thompson           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
4268b97b69aSJeremy L Thompson           break;
4278b97b69aSJeremy L Thompson         }
428f815fac9SJeremy L Thompson         // LCOV_EXCL_START
429f815fac9SJeremy L Thompson         case CEED_RESTRICTION_ORIENTED:
430f815fac9SJeremy L Thompson         case CEED_RESTRICTION_CURL_ORIENTED:
431f815fac9SJeremy L Thompson           break;  // TODO: Not implemented
432f815fac9SJeremy L Thompson                   // LCOV_EXCL_STOP
4334b3e95d5SJeremy L Thompson       }
4344b3e95d5SJeremy L Thompson     }
4354b3e95d5SJeremy L Thompson   } else {
4364b3e95d5SJeremy L Thompson     // Output
437f815fac9SJeremy L Thompson     switch (rstr_type) {
438f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STANDARD: {
4394b3e95d5SJeremy L Thompson         CeedInt comp_stride;
4404b3e95d5SJeremy L Thompson 
4414b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
4424b3e95d5SJeremy L Thompson         code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
4434b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
4444b3e95d5SJeremy L Thompson         code << "    // CompStride: " << comp_stride << "\n";
4454b3e95d5SJeremy L Thompson         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
446*412e5683SJeremy L Thompson         code << "    WriteLVecStandard" << (is_all_tensor ? max_dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name
447dc007f05SJeremy L Thompson              << ">(data, l_size" << var_suffix << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n";
448f815fac9SJeremy L Thompson         break;
449f815fac9SJeremy L Thompson       }
450f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STRIDED: {
4514b3e95d5SJeremy L Thompson         bool    has_backend_strides;
4524b3e95d5SJeremy L Thompson         CeedInt num_elem;
4534b3e95d5SJeremy L Thompson 
4544b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
4554b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
4564b3e95d5SJeremy L Thompson         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
4574b3e95d5SJeremy L Thompson 
4584b3e95d5SJeremy L Thompson         if (!has_backend_strides) {
4594b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
4604b3e95d5SJeremy L Thompson         }
4614b3e95d5SJeremy L Thompson         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
462*412e5683SJeremy L Thompson         code << "    WriteLVecStrided" << (is_all_tensor ? max_dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", "
463dc007f05SJeremy L Thompson              << strides[1] << ", " << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n";
464f815fac9SJeremy L Thompson         break;
465f815fac9SJeremy L Thompson       }
4668b97b69aSJeremy L Thompson       case CEED_RESTRICTION_POINTS:
4678b97b69aSJeremy L Thompson         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
4688b97b69aSJeremy L Thompson         break;
469f815fac9SJeremy L Thompson       // LCOV_EXCL_START
470f815fac9SJeremy L Thompson       case CEED_RESTRICTION_ORIENTED:
471f815fac9SJeremy L Thompson       case CEED_RESTRICTION_CURL_ORIENTED:
472f815fac9SJeremy L Thompson         break;  // TODO: Not implemented
473f815fac9SJeremy L Thompson                 // LCOV_EXCL_STOP
4744b3e95d5SJeremy L Thompson     }
4754b3e95d5SJeremy L Thompson   }
476681d0ea7SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
4774b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4784b3e95d5SJeremy L Thompson }
4794b3e95d5SJeremy L Thompson 
4804b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4814b3e95d5SJeremy L Thompson // Basis
4824b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
483*412e5683SJeremy L Thompson static int CeedOperatorBuildKernelBasis_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedInt max_dim,
484*412e5683SJeremy L Thompson                                                  CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input,
485*412e5683SJeremy L Thompson                                                  bool is_all_tensor, bool is_at_points, bool use_3d_slices) {
486*412e5683SJeremy L Thompson   bool      is_tensor = true;
487*412e5683SJeremy L Thompson   CeedBasis basis;
488*412e5683SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
489*412e5683SJeremy L Thompson   CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
490*412e5683SJeremy L Thompson 
4914b3e95d5SJeremy L Thompson   std::string         var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
492dc007f05SJeremy L Thompson   std::string         P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q";
4934b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
494*412e5683SJeremy L Thompson   CeedInt             dim = max_dim, elem_size = 0, num_comp = 0, P_1d = 0;
4954b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
4964b3e95d5SJeremy L Thompson 
4974b3e95d5SJeremy L Thompson   // Get field data
4984b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
4994b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
5004b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
5014b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
5024b3e95d5SJeremy L Thompson   }
503681d0ea7SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
5044b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
505*412e5683SJeremy L Thompson     CeedCallBackend(CeedBasisGetDimension(basis, &dim));
506dc007f05SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
507dc007f05SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
5084b3e95d5SJeremy L Thompson   }
5094b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
5104b3e95d5SJeremy L Thompson 
5114b3e95d5SJeremy L Thompson   // Basis
5124b3e95d5SJeremy L Thompson   code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
5134b3e95d5SJeremy L Thompson   if (is_input) {
5144b3e95d5SJeremy L Thompson     switch (eval_mode) {
5154b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
5168b97b69aSJeremy L Thompson         if (!use_3d_slices && !is_at_points) {
5174b3e95d5SJeremy L Thompson           code << "    CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n";
5184b3e95d5SJeremy L Thompson         }
5194b3e95d5SJeremy L Thompson         break;
5204b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
5218b97b69aSJeremy L Thompson         if (is_at_points) {
522dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
523dc007f05SJeremy L Thompson 
524dc007f05SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
52599421279SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_e" << var_suffix
52699421279SJeremy L Thompson                << ", s_B" << var_suffix << ", r_c" << var_suffix << ");\n";
5278b97b69aSJeremy L Thompson         } else {
528*412e5683SJeremy L Thompson           std::string function_name = is_tensor
529*412e5683SJeremy L Thompson                                           ? ((dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d" + (is_all_tensor ? "" : "Flattened"))
530*412e5683SJeremy L Thompson                                           : "InterpNonTensor";
531dc007f05SJeremy L Thompson 
532dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
533*412e5683SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", " << (P_1d > Q_1d ? P_name : Q_name)
534*412e5683SJeremy L Thompson                << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n";
5358b97b69aSJeremy L Thompson         }
5364b3e95d5SJeremy L Thompson         break;
5374b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
5388b97b69aSJeremy L Thompson         if (is_at_points) {
539dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
540dc007f05SJeremy L Thompson 
541dc007f05SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
54299421279SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_e" << var_suffix
54399421279SJeremy L Thompson                << ", s_B" << var_suffix << ", r_c" << var_suffix << ");\n";
5448b97b69aSJeremy L Thompson         } else if (use_3d_slices) {
545dc007f05SJeremy L Thompson           std::string function_name = (dim > 1 ? "InterpTensor" : "Interp") + std::to_string(dim) + "d";
546dc007f05SJeremy L Thompson 
5474b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
54899421279SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_e" << var_suffix
54999421279SJeremy L Thompson                << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n";
550dc007f05SJeremy L Thompson         } else if (is_tensor) {
551dc007f05SJeremy L Thompson           bool        is_collocated = dim == 3 && Q_1d >= P_1d;
552*412e5683SJeremy L Thompson           std::string function_name = (dim == 1 ? "Grad" : (is_collocated ? "GradTensorCollocated" : "GradTensor")) + std::to_string(dim) + "d" +
553*412e5683SJeremy L Thompson                                       (is_all_tensor ? "" : "Flattened");
554dc007f05SJeremy L Thompson 
555*412e5683SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "*" << (dim >= 3 ? Q_name : "1")
556*412e5683SJeremy L Thompson                << "];\n";
557*412e5683SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", " << (P_1d > Q_1d ? P_name : Q_name)
558*412e5683SJeremy L Thompson                << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
5594b3e95d5SJeremy L Thompson         } else {
560dc007f05SJeremy L Thompson           std::string function_name = "GradNonTensor";
561dc007f05SJeremy L Thompson 
562*412e5683SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "];\n";
563*412e5683SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", dim" << var_suffix << ", " << P_name << ", " << Q_name << ", "
564*412e5683SJeremy L Thompson                << (P_1d > Q_1d ? P_name : Q_name) << ">(data, r_e" << var_suffix << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
5654b3e95d5SJeremy L Thompson         }
5664b3e95d5SJeremy L Thompson         break;
5674b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
5688b97b69aSJeremy L Thompson         if (is_at_points) {
5698b97b69aSJeremy L Thompson           code << "    // Nothing to do AtPoints\n";
5708b97b69aSJeremy L Thompson         } else {
5714b3e95d5SJeremy L Thompson           CeedBasis_Cuda_shared *basis_data;
572*412e5683SJeremy L Thompson           std::string            function_name = is_tensor
573*412e5683SJeremy L Thompson                                                      ? ((dim == 1 ? "Weight" : "WeightTensor") + std::to_string(dim) + "d" + (is_all_tensor ? "" : "Flattened"))
574*412e5683SJeremy L Thompson                                                      : "WeightNonTensor";
5754b3e95d5SJeremy L Thompson 
576dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
5774b3e95d5SJeremy L Thompson           CeedCallBackend(CeedBasisGetData(basis, &basis_data));
5784b3e95d5SJeremy L Thompson           data->W = basis_data->d_q_weight_1d;
579343e3094SJeremy L Thompson           code << "    " << function_name << "<" << P_name << ", " << Q_name << ">(data, W, r_q" << var_suffix << ");\n";
5808b97b69aSJeremy L Thompson         }
5814b3e95d5SJeremy L Thompson         break;
5824b3e95d5SJeremy L Thompson       }
5834b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
5844b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
5854b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
5864b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
5874b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
5884b3e95d5SJeremy L Thompson     }
5894b3e95d5SJeremy L Thompson   } else {
5904b3e95d5SJeremy L Thompson     switch (eval_mode) {
5914b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
5924b3e95d5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n";
5934b3e95d5SJeremy L Thompson         break;  // No action
5944b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
595e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
5968b97b69aSJeremy L Thompson         if (is_at_points) {
597dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
598dc007f05SJeremy L Thompson 
59999421279SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_c" << var_suffix
60099421279SJeremy L Thompson                << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
6018b97b69aSJeremy L Thompson         } else {
602dc007f05SJeremy L Thompson           std::string function_name =
603*412e5683SJeremy L Thompson               is_tensor ? ((dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d" + (is_all_tensor ? "" : "Flattened"))
604*412e5683SJeremy L Thompson                         : "InterpTransposeNonTensor";
605dc007f05SJeremy L Thompson 
606*412e5683SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", " << (P_1d > Q_1d ? P_name : Q_name)
607*412e5683SJeremy L Thompson                << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
6088b97b69aSJeremy L Thompson         }
6094b3e95d5SJeremy L Thompson         break;
6104b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
611e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
6128b97b69aSJeremy L Thompson         if (is_at_points) {
613dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
614dc007f05SJeremy L Thompson 
61599421279SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_c" << var_suffix
61699421279SJeremy L Thompson                << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
6178b97b69aSJeremy L Thompson         } else if (use_3d_slices) {
618dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
619dc007f05SJeremy L Thompson 
62099421279SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_q" << var_suffix
62199421279SJeremy L Thompson                << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
622dc007f05SJeremy L Thompson         } else if (is_tensor) {
623dc007f05SJeremy L Thompson           bool        is_collocated = dim == 3 && Q_1d >= P_1d;
624*412e5683SJeremy L Thompson           std::string function_name = (dim == 1 ? "GradTranspose" : (is_collocated ? "GradTransposeTensorCollocated" : "GradTransposeTensor")) +
625*412e5683SJeremy L Thompson                                       std::to_string(dim) + "d" + (is_all_tensor ? "" : "Flattened");
626dc007f05SJeremy L Thompson 
627*412e5683SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", " << (P_1d > Q_1d ? P_name : Q_name)
628*412e5683SJeremy L Thompson                << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
6294b3e95d5SJeremy L Thompson         } else {
630dc007f05SJeremy L Thompson           std::string function_name = "GradTransposeNonTensor";
631dc007f05SJeremy L Thompson 
632*412e5683SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", dim" << var_suffix << ", " << P_name << ", " << Q_name << ", "
633*412e5683SJeremy L Thompson                << (P_1d > Q_1d ? P_name : Q_name) << ">(data, r_q" << var_suffix << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
6344b3e95d5SJeremy L Thompson         }
6354b3e95d5SJeremy L Thompson         break;
6364b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
6374b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT:
6384b3e95d5SJeremy L Thompson         break;  // Should not occur
6394b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
6404b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
6414b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
6424b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
6434b3e95d5SJeremy L Thompson     }
6444b3e95d5SJeremy L Thompson   }
645681d0ea7SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
6464b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6474b3e95d5SJeremy L Thompson }
6484b3e95d5SJeremy L Thompson 
6494b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
6504b3e95d5SJeremy L Thompson // QFunction
6514b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
652*412e5683SJeremy L Thompson static int CeedOperatorBuildKernelQFunction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt max_dim, CeedInt max_num_points,
6538b97b69aSJeremy L Thompson                                                      CeedInt num_input_fields, CeedOperatorField *op_input_fields,
6548b97b69aSJeremy L Thompson                                                      CeedQFunctionField *qf_input_fields, CeedInt num_output_fields,
6558b97b69aSJeremy L Thompson                                                      CeedOperatorField *op_output_fields, CeedQFunctionField *qf_output_fields,
656*412e5683SJeremy L Thompson                                                      std::string qfunction_name, CeedInt Q_1d, bool is_all_tensor, bool is_at_points,
657dc007f05SJeremy L Thompson                                                      bool use_3d_slices) {
658*412e5683SJeremy L Thompson   std::string         Q_name    = is_all_tensor ? "Q_1d" : "Q";
6594b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
6604b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
6614b3e95d5SJeremy L Thompson 
6628b97b69aSJeremy L Thompson   // Setup output arrays
6634b3e95d5SJeremy L Thompson   code << "\n    // -- Output field setup\n";
6644b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
66559fa3f92SJeremy L Thompson     const char *field_name;
6664b3e95d5SJeremy L Thompson     std::string var_suffix = "_out_" + std::to_string(i);
6674b3e95d5SJeremy L Thompson 
66859fa3f92SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
66959fa3f92SJeremy L Thompson     code << "    // ---- Output field " << i << ": " << field_name << "\n";
6704b3e95d5SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
6718b97b69aSJeremy L Thompson     switch (eval_mode) {
6728b97b69aSJeremy L Thompson       case CEED_EVAL_NONE:
6738b97b69aSJeremy L Thompson         if (is_at_points) {
6748b97b69aSJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "];\n";
6758b97b69aSJeremy L Thompson         } else {
676*412e5683SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_all_tensor && (max_dim >= 3) ? Q_name : "1")
677*412e5683SJeremy L Thompson                << "];\n";
6784b3e95d5SJeremy L Thompson         }
6798b97b69aSJeremy L Thompson         break;
6808b97b69aSJeremy L Thompson       case CEED_EVAL_INTERP:
6818b97b69aSJeremy L Thompson         if (is_at_points) {
6828b97b69aSJeremy L Thompson           // Accumulator for point data
683*412e5683SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (max_dim >= 3 ? Q_name : "1") << "];\n";
684*412e5683SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (max_dim >= 3 ? Q_name : "1") << "; i++) {\n";
6858b97b69aSJeremy L Thompson           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
6868b97b69aSJeremy L Thompson           code << "    }\n";
6878b97b69aSJeremy L Thompson         } else {
688*412e5683SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_all_tensor && (max_dim >= 3) ? Q_name : "1")
689*412e5683SJeremy L Thompson                << "];\n";
6908b97b69aSJeremy L Thompson         }
6918b97b69aSJeremy L Thompson         break;
6928b97b69aSJeremy L Thompson       case CEED_EVAL_GRAD:
6938b97b69aSJeremy L Thompson         if (is_at_points) {
6948b97b69aSJeremy L Thompson           // Accumulator for point data
695*412e5683SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (max_dim >= 3 ? Q_name : "1") << "*dim" << var_suffix
696*412e5683SJeremy L Thompson                << "];\n";
697*412e5683SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (max_dim >= 3 ? Q_name : "1") << "; i++) {\n";
6988b97b69aSJeremy L Thompson           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
6998b97b69aSJeremy L Thompson           code << "    }\n";
7008b97b69aSJeremy L Thompson         } else if (use_3d_slices) {
7014b3e95d5SJeremy L Thompson           // Accumulator for gradient slices
7024b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
7034b3e95d5SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n";
7044b3e95d5SJeremy L Thompson           code << "      r_q" << var_suffix << "[i] = 0.0;\n";
7054b3e95d5SJeremy L Thompson           code << "    }\n";
7064b3e95d5SJeremy L Thompson         } else {
707*412e5683SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "*"
708*412e5683SJeremy L Thompson                << (is_all_tensor && (max_dim >= 3) ? Q_name : "1") << "];\n";
7094b3e95d5SJeremy L Thompson         }
7108b97b69aSJeremy L Thompson         break;
7118b97b69aSJeremy L Thompson       case CEED_EVAL_WEIGHT:
7128b97b69aSJeremy L Thompson         break;
7138b97b69aSJeremy L Thompson         // LCOV_EXCL_START
7148b97b69aSJeremy L Thompson       case CEED_EVAL_DIV:
7158b97b69aSJeremy L Thompson       case CEED_EVAL_CURL:
7168b97b69aSJeremy L Thompson         break;  // TODO: Not implemented
7178b97b69aSJeremy L Thompson                 // LCOV_EXCL_STOP
7184b3e95d5SJeremy L Thompson     }
7194b3e95d5SJeremy L Thompson   }
7204b3e95d5SJeremy L Thompson 
7218b97b69aSJeremy L Thompson   if (is_at_points) {
7228b97b69aSJeremy L Thompson     // We need to handle batches of points
7238b97b69aSJeremy L Thompson     code << "\n    // Note: Using batches of points\n";
7248b97b69aSJeremy L Thompson     code << "    const CeedInt point_loop_bound = (blockDim.x * blockDim.y) * ceil(1.0 * max_num_points / (blockDim.x * blockDim.y));\n\n";
7258b97b69aSJeremy L Thompson     code << "    #pragma unroll\n";
7268b97b69aSJeremy L Thompson     code << "    for (CeedInt i = threadIdx.x + threadIdx.y * blockDim.x; i < point_loop_bound; i += blockDim.x * blockDim.y) {\n";
7278b97b69aSJeremy L Thompson     code << "      const CeedInt p = i % max_num_points;\n\n";
7288b97b69aSJeremy L Thompson 
7298b97b69aSJeremy L Thompson     code << "      // -- Coordinates\n";
730*412e5683SJeremy L Thompson     code << "      CeedScalar r_x[max_dim];\n";
731*412e5683SJeremy L Thompson     code << "      ReadPoint<max_dim, coords_comp_stride, max_num_points>(data, elem, p, max_num_points, points.indices, points.coords, r_x);\n\n";
7328b97b69aSJeremy L Thompson 
7338b97b69aSJeremy L Thompson     code << "      // -- Input fields\n";
7348b97b69aSJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
73559fa3f92SJeremy L Thompson       const char *field_name;
7368b97b69aSJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
737f725b54bSJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
7388b97b69aSJeremy L Thompson 
73959fa3f92SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name));
74059fa3f92SJeremy L Thompson       code << "      // ---- Input field " << i << ": " << field_name << "\n";
7418b97b69aSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
7428b97b69aSJeremy L Thompson       // Basis action
7438b97b69aSJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
7448b97b69aSJeremy L Thompson       switch (eval_mode) {
7458b97b69aSJeremy L Thompson         case CEED_EVAL_NONE:
7468b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7478b97b69aSJeremy L Thompson           code << "      ReadPoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
7488b97b69aSJeremy L Thompson                << ", max_num_points>(data, elem, p, max_num_points, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
7498b97b69aSJeremy L Thompson           break;
7508b97b69aSJeremy L Thompson         case CEED_EVAL_INTERP:
7518b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
752*412e5683SJeremy L Thompson           code << "      InterpAtPoints" << max_dim << "d<num_comp" << var_suffix << ", max_num_points, " << P_name << ", " << Q_name
753*412e5683SJeremy L Thompson                << ">(data, i, r_c" << var_suffix << ", r_x, r_s" << var_suffix << ");\n";
7548b97b69aSJeremy L Thompson           break;
7558b97b69aSJeremy L Thompson         case CEED_EVAL_GRAD:
756*412e5683SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "];\n";
757*412e5683SJeremy L Thompson           code << "      GradAtPoints" << max_dim << "d<num_comp" << var_suffix << ", max_num_points, " << P_name << ", " << Q_name
758*412e5683SJeremy L Thompson                << ">(data, i, r_c" << var_suffix << ", r_x, r_s" << var_suffix << ");\n";
7598b97b69aSJeremy L Thompson           break;
7608b97b69aSJeremy L Thompson         case CEED_EVAL_WEIGHT:
7618b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
7628b97b69aSJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = 1.0;\n";
7638b97b69aSJeremy L Thompson           break;
7648b97b69aSJeremy L Thompson           // LCOV_EXCL_START
7658b97b69aSJeremy L Thompson         case CEED_EVAL_DIV:
7668b97b69aSJeremy L Thompson         case CEED_EVAL_CURL:
7678b97b69aSJeremy L Thompson           break;  // TODO: Not implemented
7688b97b69aSJeremy L Thompson                   // LCOV_EXCL_STOP
7698b97b69aSJeremy L Thompson       }
7708b97b69aSJeremy L Thompson     }
7718b97b69aSJeremy L Thompson     code << "\n      // -- Output fields\n";
7728b97b69aSJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
77359fa3f92SJeremy L Thompson       const char *field_name;
7748b97b69aSJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
7758b97b69aSJeremy L Thompson 
77659fa3f92SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
77759fa3f92SJeremy L Thompson       code << "      // ---- Output field " << i << ": " << field_name << "\n";
7788b97b69aSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
7798b97b69aSJeremy L Thompson       // Basis action
7808b97b69aSJeremy L Thompson       switch (eval_mode) {
7818b97b69aSJeremy L Thompson         case CEED_EVAL_NONE:
7828b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7838b97b69aSJeremy L Thompson           break;
7848b97b69aSJeremy L Thompson         case CEED_EVAL_INTERP:
7858b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7868b97b69aSJeremy L Thompson           break;
7878b97b69aSJeremy L Thompson         case CEED_EVAL_GRAD:
788*412e5683SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "];\n";
7898b97b69aSJeremy L Thompson           break;
7908b97b69aSJeremy L Thompson           // LCOV_EXCL_START
7918b97b69aSJeremy L Thompson         case CEED_EVAL_WEIGHT:
7928b97b69aSJeremy L Thompson           break;  // Should not occur
7938b97b69aSJeremy L Thompson         case CEED_EVAL_DIV:
7948b97b69aSJeremy L Thompson         case CEED_EVAL_CURL:
7958b97b69aSJeremy L Thompson           break;  // TODO: Not implemented
7968b97b69aSJeremy L Thompson                   // LCOV_EXCL_STOP
7978b97b69aSJeremy L Thompson       }
7988b97b69aSJeremy L Thompson     }
7998b97b69aSJeremy L Thompson 
8008b97b69aSJeremy L Thompson   } else if (use_3d_slices) {
8014b3e95d5SJeremy L Thompson     // We treat quadrature points per slice in 3d to save registers
8024b3e95d5SJeremy L Thompson     code << "\n    // Note: Using planes of 3D elements\n";
8034b3e95d5SJeremy L Thompson     code << "    #pragma unroll\n";
8044b3e95d5SJeremy L Thompson     code << "    for (CeedInt q = 0; q < " << Q_name << "; q++) {\n";
8054b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
8064b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
80759fa3f92SJeremy L Thompson       const char *field_name;
8084b3e95d5SJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
8094b3e95d5SJeremy L Thompson 
81059fa3f92SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name));
81159fa3f92SJeremy L Thompson       code << "      // ---- Input field " << i << ": " << field_name << "\n";
8124b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
8134b3e95d5SJeremy L Thompson       // Basis action
8144b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
8154b3e95d5SJeremy L Thompson       switch (eval_mode) {
8164b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
8174b3e95d5SJeremy L Thompson           bool is_strided;
8184b3e95d5SJeremy L Thompson 
8194b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
8204b3e95d5SJeremy L Thompson 
8214b3e95d5SJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
8224b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
8234b3e95d5SJeremy L Thompson           if (is_strided) {
8244b3e95d5SJeremy L Thompson             bool    has_backend_strides;
8254b3e95d5SJeremy L Thompson             CeedInt num_elem, elem_size;
8264b3e95d5SJeremy L Thompson 
8274b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
8284b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
8294b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
8304b3e95d5SJeremy L Thompson             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
8314b3e95d5SJeremy L Thompson 
8324b3e95d5SJeremy L Thompson             if (!has_backend_strides) {
8334b3e95d5SJeremy L Thompson               CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
8344b3e95d5SJeremy L Thompson             }
8354b3e95d5SJeremy L Thompson             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
836f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStrided3d<num_comp" << var_suffix << ", " << Q_name << ", " << strides[0] << ", " << strides[1] << ", "
8374b3e95d5SJeremy L Thompson                  << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n";
8384b3e95d5SJeremy L Thompson           } else {
8394b3e95d5SJeremy L Thompson             CeedSize                  l_size = 0;
8404b3e95d5SJeremy L Thompson             CeedInt                   comp_stride;
8414b3e95d5SJeremy L Thompson             CeedElemRestriction_Cuda *rstr_data;
8424b3e95d5SJeremy L Thompson 
8434b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
8444b3e95d5SJeremy L Thompson             code << "      const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
8454b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
8464b3e95d5SJeremy L Thompson             code << "      // CompStride: " << comp_stride << "\n";
8474b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
8484b3e95d5SJeremy L Thompson             data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
849f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStandard3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix
8504b3e95d5SJeremy L Thompson                  << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
8514b3e95d5SJeremy L Thompson           }
852681d0ea7SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
8534b3e95d5SJeremy L Thompson           break;
8544b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
8554b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
8564b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n";
8574b3e95d5SJeremy L Thompson           code << "        r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n";
8584b3e95d5SJeremy L Thompson           code << "      }\n";
8594b3e95d5SJeremy L Thompson           break;
8604b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
861*412e5683SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "];\n";
86299421279SJeremy L Thompson           code << "      GradColloSlice3d<num_comp" << var_suffix << ", " << Q_name << ", OP_T_1D>(data, q, r_q" << var_suffix << ", s_G"
86399421279SJeremy L Thompson                << var_suffix << ", r_s" << var_suffix << ");\n";
8644b3e95d5SJeremy L Thompson           break;
8654b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
8664b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
8674b3e95d5SJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n";
8688b97b69aSJeremy L Thompson           break;
8694b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
8704b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
8714b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
8724b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
8734b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
8744b3e95d5SJeremy L Thompson       }
8754b3e95d5SJeremy L Thompson     }
8764b3e95d5SJeremy L Thompson     code << "\n      // -- Output fields\n";
8774b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
87859fa3f92SJeremy L Thompson       const char *field_name;
8794b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
8804b3e95d5SJeremy L Thompson 
88159fa3f92SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
88259fa3f92SJeremy L Thompson       code << "      // ---- Output field " << i << ": " << field_name << "\n";
8834b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
8844b3e95d5SJeremy L Thompson       // Basis action
8854b3e95d5SJeremy L Thompson       switch (eval_mode) {
8864b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
8874b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
8888b97b69aSJeremy L Thompson           break;
8894b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
8904b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
8914b3e95d5SJeremy L Thompson           break;
8924b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
893*412e5683SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "];\n";
8944b3e95d5SJeremy L Thompson           break;
8954b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
8964b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
8974b3e95d5SJeremy L Thompson           break;  // Should not occur
8984b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
8994b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
9004b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
9014b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
9024b3e95d5SJeremy L Thompson       }
9034b3e95d5SJeremy L Thompson     }
9044b3e95d5SJeremy L Thompson   } else {
9054b3e95d5SJeremy L Thompson     code << "\n    // Note: Using full elements\n";
9064b3e95d5SJeremy L Thompson     code << "    {\n";
9074b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
9084b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
90959fa3f92SJeremy L Thompson       const char *field_name;
91059fa3f92SJeremy L Thompson 
91159fa3f92SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name));
91259fa3f92SJeremy L Thompson       code << "      // ---- Input field " << i << ": " << field_name << "\n";
9134b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n";
9144b3e95d5SJeremy L Thompson     }
9154b3e95d5SJeremy L Thompson     code << "      // -- Output fields\n";
9164b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
91759fa3f92SJeremy L Thompson       const char *field_name;
91859fa3f92SJeremy L Thompson 
91959fa3f92SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
92059fa3f92SJeremy L Thompson       code << "      // ---- Output field " << i << ": " << field_name << "\n";
9214b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n";
9224b3e95d5SJeremy L Thompson     }
9234b3e95d5SJeremy L Thompson   }
9244b3e95d5SJeremy L Thompson 
9254b3e95d5SJeremy L Thompson   // Input and output buffers
9264b3e95d5SJeremy L Thompson   code << "\n      // -- QFunction inputs and outputs\n";
9274b3e95d5SJeremy L Thompson   code << "      // ---- Inputs\n";
9284b3e95d5SJeremy L Thompson   code << "      CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n";
9294b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
93059fa3f92SJeremy L Thompson     const char *field_name;
93159fa3f92SJeremy L Thompson 
93259fa3f92SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name));
93359fa3f92SJeremy L Thompson     code << "      // ------ Input field " << i << ": " << field_name << "\n";
9344b3e95d5SJeremy L Thompson     code << "      inputs[" << i << "] = r_s_in_" << i << ";\n";
9354b3e95d5SJeremy L Thompson   }
9364b3e95d5SJeremy L Thompson   code << "      // ---- Outputs\n";
9374b3e95d5SJeremy L Thompson   code << "      CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n";
9384b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
93959fa3f92SJeremy L Thompson     const char *field_name;
94059fa3f92SJeremy L Thompson 
94159fa3f92SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
94259fa3f92SJeremy L Thompson     code << "      // ------ Output field " << i << ": " << field_name << "\n";
9434b3e95d5SJeremy L Thompson     code << "      outputs[" << i << "] = r_s_out_" << i << ";\n";
9444b3e95d5SJeremy L Thompson   }
9454b3e95d5SJeremy L Thompson 
9464b3e95d5SJeremy L Thompson   // Apply QFunction
9474b3e95d5SJeremy L Thompson   code << "\n      // -- Apply QFunction\n";
9484b3e95d5SJeremy L Thompson   code << "      " << qfunction_name << "(ctx, ";
949*412e5683SJeremy L Thompson   if (max_dim != 3 || is_at_points || use_3d_slices || !is_all_tensor) {
9504b3e95d5SJeremy L Thompson     code << "1";
9514b3e95d5SJeremy L Thompson   } else {
952dc007f05SJeremy L Thompson     code << Q_name;
9534b3e95d5SJeremy L Thompson   }
9544b3e95d5SJeremy L Thompson   code << ", inputs, outputs);\n";
9554b3e95d5SJeremy L Thompson 
9568b97b69aSJeremy L Thompson   if (is_at_points) {
9578b97b69aSJeremy L Thompson     // Map back to coefficients
9588b97b69aSJeremy L Thompson     code << "\n      // -- Output fields\n";
9598b97b69aSJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
96059fa3f92SJeremy L Thompson       const char *field_name;
9618b97b69aSJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
9628b97b69aSJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
9638b97b69aSJeremy L Thompson 
96459fa3f92SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
96559fa3f92SJeremy L Thompson       code << "      // ---- Output field " << i << ": " << field_name << "\n";
9668b97b69aSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
9678b97b69aSJeremy L Thompson       // Basis action
9688b97b69aSJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
9698b97b69aSJeremy L Thompson       switch (eval_mode) {
9708b97b69aSJeremy L Thompson         case CEED_EVAL_NONE: {
9718b97b69aSJeremy L Thompson           CeedInt             comp_stride;
9728b97b69aSJeremy L Thompson           CeedElemRestriction elem_rstr;
9738b97b69aSJeremy L Thompson 
9748b97b69aSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
9758b97b69aSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
9768b97b69aSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
9778b97b69aSJeremy L Thompson           code << "      const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
9788b97b69aSJeremy L Thompson           code << "      WritePoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
9798b97b69aSJeremy L Thompson                << ", max_num_points>(data, elem, i, points.num_per_elem[elem], indices.outputs[" << i << "]"
9808b97b69aSJeremy L Thompson                << ", r_s" << var_suffix << ", d" << var_suffix << ");\n";
9818b97b69aSJeremy L Thompson           break;
9828b97b69aSJeremy L Thompson         }
9838b97b69aSJeremy L Thompson         case CEED_EVAL_INTERP:
9848b97b69aSJeremy L Thompson           code << "      if (i >= points.num_per_elem[elem]) {\n";
9858b97b69aSJeremy L Thompson           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) r_s" << var_suffix << "[j] = 0.0;\n";
9868b97b69aSJeremy L Thompson           code << "      }\n";
987*412e5683SJeremy L Thompson           code << "      InterpTransposeAtPoints" << max_dim << "d<num_comp" << var_suffix << ", max_num_points, " << P_name << ", " << Q_name
988f725b54bSJeremy L Thompson                << ">(data, i, r_s" << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
9898b97b69aSJeremy L Thompson           break;
9908b97b69aSJeremy L Thompson         case CEED_EVAL_GRAD:
9918b97b69aSJeremy L Thompson           code << "      if (i >= points.num_per_elem[elem]) {\n";
992*412e5683SJeremy L Thompson           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "*dim" << var_suffix << "; j++) r_s" << var_suffix << "[j] = 0.0;\n";
9938b97b69aSJeremy L Thompson           code << "      }\n";
994*412e5683SJeremy L Thompson           code << "      GradTransposeAtPoints" << max_dim << "d<num_comp" << var_suffix << ", max_num_points, " << P_name << ", " << Q_name
995f725b54bSJeremy L Thompson                << ">(data, i, r_s" << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
9968b97b69aSJeremy L Thompson           break;
9978b97b69aSJeremy L Thompson           // LCOV_EXCL_START
9988b97b69aSJeremy L Thompson         case CEED_EVAL_WEIGHT:
9998b97b69aSJeremy L Thompson           break;  // Should not occur
10008b97b69aSJeremy L Thompson         case CEED_EVAL_DIV:
10018b97b69aSJeremy L Thompson         case CEED_EVAL_CURL:
10028b97b69aSJeremy L Thompson           break;  // TODO: Not implemented
10038b97b69aSJeremy L Thompson                   // LCOV_EXCL_STOP
10048b97b69aSJeremy L Thompson       }
10058b97b69aSJeremy L Thompson     }
10068b97b69aSJeremy L Thompson   } else if (use_3d_slices) {
10074b3e95d5SJeremy L Thompson     // Copy or apply transpose grad, if needed
10088b97b69aSJeremy L Thompson     code << "\n      // -- Output fields\n";
10094b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
101059fa3f92SJeremy L Thompson       const char *field_name;
10114b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
10124b3e95d5SJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
10134b3e95d5SJeremy L Thompson 
101459fa3f92SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
101559fa3f92SJeremy L Thompson       code << "      // ---- Output field " << i << ": " << field_name << "\n";
10164b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
10174b3e95d5SJeremy L Thompson       // Basis action
10184b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
10194b3e95d5SJeremy L Thompson       switch (eval_mode) {
10204b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
10214b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
10224b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
10234b3e95d5SJeremy L Thompson           code << "      }\n";
10248b97b69aSJeremy L Thompson           break;
10254b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
10264b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
10274b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
10284b3e95d5SJeremy L Thompson           code << "      }\n";
10294b3e95d5SJeremy L Thompson           break;
10304b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
103199421279SJeremy L Thompson           code << "      GradColloSliceTranspose3d<num_comp" << var_suffix << ", " << Q_name << ", OP_T_1D>(data, q, r_s" << var_suffix << ", s_G"
1032f815fac9SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
10334b3e95d5SJeremy L Thompson           break;
10344b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
10354b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
10364b3e95d5SJeremy L Thompson           break;  // Should not occur
10374b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
10384b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
10394b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
10404b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
10414b3e95d5SJeremy L Thompson       }
10424b3e95d5SJeremy L Thompson     }
10434b3e95d5SJeremy L Thompson   }
10444b3e95d5SJeremy L Thompson   code << "    }\n";
10454b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10464b3e95d5SJeremy L Thompson }
10474b3e95d5SJeremy L Thompson 
10484b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1049b2165e7aSSebastian Grimberg // Build single operator kernel
1050ab213215SJeremy L Thompson //------------------------------------------------------------------------------
1051ddae5012SJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Cuda_gen(CeedOperator op, bool *is_good_build) {
1052*412e5683SJeremy L Thompson   bool                    is_all_tensor = true, is_all_nontensor = true, is_at_points = false, use_3d_slices = false;
1053241a4b83SYohann   Ceed                    ceed;
1054*412e5683SJeremy L Thompson   CeedInt                 Q, Q_1d, num_input_fields, num_output_fields, max_dim = 1, max_num_points = 0, coords_comp_stride = 0;
1055ca735530SJeremy L Thompson   CeedQFunctionField     *qf_input_fields, *qf_output_fields;
1056ca735530SJeremy L Thompson   CeedQFunction_Cuda_gen *qf_data;
1057ca735530SJeremy L Thompson   CeedQFunction           qf;
1058ca735530SJeremy L Thompson   CeedOperatorField      *op_input_fields, *op_output_fields;
1059ca735530SJeremy L Thompson   CeedOperator_Cuda_gen  *data;
10604b3e95d5SJeremy L Thompson   std::ostringstream      code;
10614b3e95d5SJeremy L Thompson 
1062ddae5012SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
10634b3e95d5SJeremy L Thompson   {
10644b3e95d5SJeremy L Thompson     bool is_setup_done;
1065ca735530SJeremy L Thompson 
1066ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
1067ddae5012SJeremy L Thompson     if (is_setup_done) {
1068ddae5012SJeremy L Thompson       *is_good_build = !data->use_fallback;
1069ddae5012SJeremy L Thompson       return CEED_ERROR_SUCCESS;
1070ddae5012SJeremy L Thompson     }
1071ddae5012SJeremy L Thompson   }
1072ddae5012SJeremy L Thompson 
1073ddae5012SJeremy L Thompson   // Check field compatibility
10748d12f40eSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
1075ddae5012SJeremy L Thompson   {
1076*412e5683SJeremy L Thompson     bool has_shared_bases = true;
1077ddae5012SJeremy L Thompson 
1078ddae5012SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
1079ddae5012SJeremy L Thompson       CeedBasis basis;
1080ddae5012SJeremy L Thompson 
1081ddae5012SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1082ddae5012SJeremy L Thompson       if (basis != CEED_BASIS_NONE) {
1083ddae5012SJeremy L Thompson         bool        is_tensor = true;
1084ddae5012SJeremy L Thompson         const char *resource;
1085ddae5012SJeremy L Thompson         char       *resource_root;
1086ddae5012SJeremy L Thompson         Ceed        basis_ceed;
1087ddae5012SJeremy L Thompson 
1088ddae5012SJeremy L Thompson         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
1089c9192acaSJeremy L Thompson         is_all_tensor    = is_all_tensor && is_tensor;
109045a787f7SJeremy L Thompson         is_all_nontensor = is_all_nontensor && !is_tensor;
1091ddae5012SJeremy L Thompson         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
1092ddae5012SJeremy L Thompson         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
1093ddae5012SJeremy L Thompson         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
1094c9192acaSJeremy L Thompson         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/cuda/shared");
1095ddae5012SJeremy L Thompson         CeedCallBackend(CeedFree(&resource_root));
1096ddae5012SJeremy L Thompson         CeedCallBackend(CeedDestroy(&basis_ceed));
1097ddae5012SJeremy L Thompson       }
1098ddae5012SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis));
10994b3e95d5SJeremy L Thompson     }
1100ca735530SJeremy L Thompson 
1101ddae5012SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
1102ddae5012SJeremy L Thompson       CeedBasis basis;
1103ddae5012SJeremy L Thompson 
1104ddae5012SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1105ddae5012SJeremy L Thompson       if (basis != CEED_BASIS_NONE) {
1106ddae5012SJeremy L Thompson         bool        is_tensor = true;
1107ddae5012SJeremy L Thompson         const char *resource;
1108ddae5012SJeremy L Thompson         char       *resource_root;
1109ddae5012SJeremy L Thompson         Ceed        basis_ceed;
1110ddae5012SJeremy L Thompson 
1111ddae5012SJeremy L Thompson         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
1112c9192acaSJeremy L Thompson         is_all_tensor    = is_all_tensor && is_tensor;
1113c9192acaSJeremy L Thompson         is_all_nontensor = is_all_nontensor && !is_tensor;
1114ddae5012SJeremy L Thompson 
1115ddae5012SJeremy L Thompson         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
1116ddae5012SJeremy L Thompson         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
1117ddae5012SJeremy L Thompson         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
1118c9192acaSJeremy L Thompson         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/cuda/shared");
1119ddae5012SJeremy L Thompson         CeedCallBackend(CeedFree(&resource_root));
1120ddae5012SJeremy L Thompson         CeedCallBackend(CeedDestroy(&basis_ceed));
1121ddae5012SJeremy L Thompson       }
1122ddae5012SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis));
1123ddae5012SJeremy L Thompson     }
1124ddae5012SJeremy L Thompson     // -- Fallback to ref if not all bases are shared
1125*412e5683SJeremy L Thompson     if (!has_shared_bases) {
1126ddae5012SJeremy L Thompson       *is_good_build = false;
1127ddae5012SJeremy L Thompson       return CEED_ERROR_SUCCESS;
1128ddae5012SJeremy L Thompson     }
1129ddae5012SJeremy L Thompson   }
1130ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1131ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1132ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
1133ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
1134241a4b83SYohann 
11354b3e95d5SJeremy L Thompson   // Get operator data
11368b97b69aSJeremy L Thompson   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
1137*412e5683SJeremy L Thompson   {
1138*412e5683SJeremy L Thompson     CeedInt max_P, max_P_1d;
1139*412e5683SJeremy L Thompson 
1140*412e5683SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelData_Cuda_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields,
1141*412e5683SJeremy L Thompson                                                          op_output_fields, qf_output_fields, &max_P, &max_P_1d, &Q, &Q_1d, &max_dim, &is_all_tensor,
1142*412e5683SJeremy L Thompson                                                          &use_3d_slices));
1143*412e5683SJeremy L Thompson     data->max_P_1d = is_all_tensor ? max_P_1d : max_P;
1144*412e5683SJeremy L Thompson   }
1145*412e5683SJeremy L Thompson   if (max_dim == 0) max_dim = 1;
1146*412e5683SJeremy L Thompson   data->dim = max_dim;
11478b97b69aSJeremy L Thompson   if (is_at_points) {
11488b97b69aSJeremy L Thompson     CeedElemRestriction_Cuda *rstr_data;
11498b97b69aSJeremy L Thompson     CeedElemRestriction       rstr_points = NULL;
11504b3e95d5SJeremy L Thompson 
11518b97b69aSJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
11528b97b69aSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points));
11538b97b69aSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetCompStride(rstr_points, &coords_comp_stride));
11548b97b69aSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(rstr_points, &rstr_data));
11558b97b69aSJeremy L Thompson     data->points.indices = (CeedInt *)rstr_data->d_offsets;
11568b97b69aSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
11578b97b69aSJeremy L Thompson   }
11588b97b69aSJeremy L Thompson   if (is_at_points) use_3d_slices = false;
11598b97b69aSJeremy L Thompson   if (Q_1d == 0) {
11608b97b69aSJeremy L Thompson     if (is_at_points) Q_1d = max_num_points;
11618b97b69aSJeremy L Thompson     else CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q_1d));
11624b3e95d5SJeremy L Thompson   }
11634b3e95d5SJeremy L Thompson   data->Q_1d = Q_1d;
11644b3e95d5SJeremy L Thompson 
11650b454692Sjeremylt   // Check for restriction only identity operator
11664b3e95d5SJeremy L Thompson   {
11674b3e95d5SJeremy L Thompson     bool is_identity_qf;
11684b3e95d5SJeremy L Thompson 
11692b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
11700b454692Sjeremylt     if (is_identity_qf) {
11719e201c85SYohann       CeedEvalMode eval_mode_in, eval_mode_out;
1172ca735530SJeremy L Thompson 
11732b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
11742b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
11756574a04fSJeremy L Thompson       CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
11766574a04fSJeremy L Thompson                 "Backend does not implement restriction only identity operators");
11770b454692Sjeremylt     }
11784b3e95d5SJeremy L Thompson   }
11790b454692Sjeremylt 
1180f1a13f77SYohann Dudouit   // Add atomicAdd function for old NVidia architectures
11814b3e95d5SJeremy L Thompson   {
11824b3e95d5SJeremy L Thompson     Ceed_Cuda            *ceed_data;
11834b3e95d5SJeremy L Thompson     struct cudaDeviceProp prop;
11844b3e95d5SJeremy L Thompson 
11852b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetData(ceed, &ceed_data));
11862b730f8bSJeremy L Thompson     CeedCallBackend(cudaGetDeviceProperties(&prop, ceed_data->device_id));
118780a9ef05SNatalie Beams     if ((prop.major < 6) && (CEED_SCALAR_TYPE != CEED_SCALAR_FP32)) {
11889c25dd66SJeremy L Thompson       code << "// AtomicAdd fallback source\n";
11899c25dd66SJeremy L Thompson       code << "#include <ceed/jit-source/cuda/cuda-atomic-add-fallback.h>\n\n";
1190f1a13f77SYohann Dudouit     }
11914b3e95d5SJeremy L Thompson   }
1192f1a13f77SYohann Dudouit 
11939e201c85SYohann   // Load basis source files
1194*412e5683SJeremy L Thompson   if (!is_all_nontensor) {
11959c25dd66SJeremy L Thompson     code << "// Tensor basis source\n";
11969c25dd66SJeremy L Thompson     code << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-templates.h>\n\n";
1197*412e5683SJeremy L Thompson   }
1198*412e5683SJeremy L Thompson   if (!is_all_tensor) {
1199dc007f05SJeremy L Thompson     code << "// Non-tensor basis source\n";
1200dc007f05SJeremy L Thompson     code << "#include <ceed/jit-source/cuda/cuda-shared-basis-nontensor-templates.h>\n\n";
1201dc007f05SJeremy L Thompson   }
1202dc007f05SJeremy L Thompson   if (is_at_points) {
12038b97b69aSJeremy L Thompson     code << "// AtPoints basis source\n";
12048b97b69aSJeremy L Thompson     code << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-at-points-templates.h>\n\n";
1205dc007f05SJeremy L Thompson   }
12069c25dd66SJeremy L Thompson   code << "// CodeGen operator source\n";
12079c25dd66SJeremy L Thompson   code << "#include <ceed/jit-source/cuda/cuda-gen-templates.h>\n\n";
1208241a4b83SYohann 
12094b3e95d5SJeremy L Thompson   // Get QFunction name
12104b3e95d5SJeremy L Thompson   std::string qfunction_name(qf_data->qfunction_name);
12114b3e95d5SJeremy L Thompson   std::string operator_name;
12124b3e95d5SJeremy L Thompson 
121309095acaSJeremy L Thompson   operator_name = "CeedKernelCudaGenOperator_" + qfunction_name;
12144d537eeaSYohann 
12159e201c85SYohann   // Define CEED_Q_VLA
12169e201c85SYohann   code << "\n#undef CEED_Q_VLA\n";
1217*412e5683SJeremy L Thompson   if (max_dim != 3 || is_at_points || use_3d_slices || !is_all_tensor) {
12189e201c85SYohann     code << "#define CEED_Q_VLA 1\n\n";
12199e201c85SYohann   } else {
12209e201c85SYohann     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
12219e201c85SYohann   }
12229e201c85SYohann 
12234b3e95d5SJeremy L Thompson   // Add user QFunction source
12244b3e95d5SJeremy L Thompson   {
12259c25dd66SJeremy L Thompson     const char *source_path;
12264b3e95d5SJeremy L Thompson 
12279c25dd66SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path));
12289c25dd66SJeremy L Thompson     CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/cuda/gen backend requires QFunction source code file");
12299c25dd66SJeremy L Thompson 
12309c25dd66SJeremy L Thompson     code << "// User QFunction source\n";
12319c25dd66SJeremy L Thompson     code << "#include \"" << source_path << "\"\n\n";
12324b3e95d5SJeremy L Thompson   }
1233241a4b83SYohann 
1234241a4b83SYohann   // Setup
1235818e0025SJeremy L Thompson   code << "\n// -----------------------------------------------------------------------------\n";
12364b3e95d5SJeremy L Thompson   code << "// Operator Kernel\n";
12374b3e95d5SJeremy L Thompson   code << "// \n";
12384b3e95d5SJeremy L Thompson   code << "// d_[in,out]_i:   CeedVector device array\n";
12394b3e95d5SJeremy L Thompson   code << "// r_[in,out]_e_i: Element vector register\n";
12404b3e95d5SJeremy L Thompson   code << "// r_[in,out]_q_i: Quadrature space vector register\n";
12419123fb08SJeremy L Thompson   code << "// r_[in,out]_c_i: AtPoints Chebyshev coefficients register\n";
12424b3e95d5SJeremy L Thompson   code << "// r_[in,out]_s_i: Quadrature space slice vector register\n";
12434b3e95d5SJeremy L Thompson   code << "// \n";
12444b3e95d5SJeremy L Thompson   code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n";
12454b3e95d5SJeremy L Thompson   code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n";
12464b3e95d5SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n";
12474b3e95d5SJeremy L Thompson   code << "extern \"C\" __global__ void " << operator_name
12488b97b69aSJeremy L Thompson        << "(CeedInt num_elem, void* ctx, FieldsInt_Cuda indices, Fields_Cuda fields, Fields_Cuda B, Fields_Cuda G, CeedScalar *W, Points_Cuda "
12498b97b69aSJeremy L Thompson           "points) {\n";
12504b3e95d5SJeremy L Thompson 
12514b3e95d5SJeremy L Thompson   // Scratch buffers
12529e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
12534b3e95d5SJeremy L Thompson     CeedEvalMode eval_mode;
12544b3e95d5SJeremy L Thompson 
12552b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
12569e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
1257826538b3SJeremy L Thompson       code << "  const CeedScalar *__restrict__ d_in_" << i << " = fields.inputs[" << i << "];\n";
1258241a4b83SYohann     }
1259241a4b83SYohann   }
12609e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1261826538b3SJeremy L Thompson     code << "  CeedScalar *__restrict__ d_out_" << i << " = fields.outputs[" << i << "];\n";
1262241a4b83SYohann   }
12631da99368SJeremy L Thompson 
1264*412e5683SJeremy L Thompson   code << "  const CeedInt max_dim = " << max_dim << ";\n";
1265*412e5683SJeremy L Thompson   if (!is_all_tensor) {
1266*412e5683SJeremy L Thompson     code << "  const CeedInt Q = " << Q << ";\n";
1267*412e5683SJeremy L Thompson   }
1268*412e5683SJeremy L Thompson   if (!is_all_nontensor) {
1269*412e5683SJeremy L Thompson     code << "  const CeedInt Q_1d = " << Q_1d << ";\n";
1270*412e5683SJeremy L Thompson   }
12718b97b69aSJeremy L Thompson   if (is_at_points) {
12728b97b69aSJeremy L Thompson     code << "  const CeedInt max_num_points = " << max_num_points << ";\n";
12738b97b69aSJeremy L Thompson     code << "  const CeedInt coords_comp_stride = " << coords_comp_stride << ";\n";
12748b97b69aSJeremy L Thompson   }
12751da99368SJeremy L Thompson 
12764b3e95d5SJeremy L Thompson   // Shared data
1277241a4b83SYohann   code << "  extern __shared__ CeedScalar slice[];\n";
12789e201c85SYohann   code << "  SharedData_Cuda data;\n";
12799e201c85SYohann   code << "  data.t_id_x = threadIdx.x;\n";
12809e201c85SYohann   code << "  data.t_id_y = threadIdx.y;\n";
12819e201c85SYohann   code << "  data.t_id_z = threadIdx.z;\n";
12829e201c85SYohann   code << "  data.t_id  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
1283*412e5683SJeremy L Thompson   code << "  data.slice = slice + data.t_id_z*OP_T_1D" << ((!is_all_tensor || max_dim == 1) ? "" : "*OP_T_1D") << ";\n";
1284920dcdc4Sjeremylt 
12850a2a6492SJeremy L Thompson   // -- Determine input mat reuse
128645a787f7SJeremy L Thompson   FieldReuse_Cuda input_matrix_reuse[CEED_FIELD_MAX];
12870a2a6492SJeremy L Thompson 
12880a2a6492SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
128945a787f7SJeremy L Thompson     input_matrix_reuse[i].index = -1;
12900a2a6492SJeremy L Thompson   }
12910a2a6492SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1292*412e5683SJeremy L Thompson     bool         is_tensor = true;
12930a2a6492SJeremy L Thompson     CeedEvalMode eval_mode_i;
12940a2a6492SJeremy L Thompson     CeedBasis    basis_i;
12950a2a6492SJeremy L Thompson 
12960a2a6492SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode_i));
12970a2a6492SJeremy L Thompson     if (eval_mode_i == CEED_EVAL_WEIGHT) continue;
12980a2a6492SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis_i));
1299*412e5683SJeremy L Thompson     CeedCallBackend(CeedBasisIsTensor(basis_i, &is_tensor));
130045a787f7SJeremy L Thompson     for (CeedInt j = 0; (input_matrix_reuse[i].index == -1) && (j < i); j++) {
13010a2a6492SJeremy L Thompson       CeedEvalMode eval_mode_j;
13020a2a6492SJeremy L Thompson       CeedBasis    basis_j;
13030a2a6492SJeremy L Thompson 
13040a2a6492SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
13050a2a6492SJeremy L Thompson       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
13060a2a6492SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
13070a2a6492SJeremy L Thompson       if (basis_i == basis_j) {
13080a2a6492SJeremy L Thompson         if (is_tensor) {
130945a787f7SJeremy L Thompson           input_matrix_reuse[i].index     = j;
131045a787f7SJeremy L Thompson           input_matrix_reuse[i].is_input  = true;
131145a787f7SJeremy L Thompson           input_matrix_reuse[i].eval_mode = eval_mode_j;
13120a2a6492SJeremy L Thompson         } else {
13130a2a6492SJeremy L Thompson           // For non-tensor can only re-use with the same eval mode
13140a2a6492SJeremy L Thompson           if (eval_mode_i == eval_mode_j) {
131545a787f7SJeremy L Thompson             input_matrix_reuse[i].index     = j;
131645a787f7SJeremy L Thompson             input_matrix_reuse[i].is_input  = true;
131745a787f7SJeremy L Thompson             input_matrix_reuse[i].eval_mode = eval_mode_j;
13180a2a6492SJeremy L Thompson           }
13190a2a6492SJeremy L Thompson         }
13200a2a6492SJeremy L Thompson       }
13210a2a6492SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis_j));
13220a2a6492SJeremy L Thompson     }
13230a2a6492SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis_i));
13240a2a6492SJeremy L Thompson   }
13250a2a6492SJeremy L Thompson 
13260a2a6492SJeremy L Thompson   // -- Determine output mat reuse
132745a787f7SJeremy L Thompson   FieldReuse_Cuda output_matrix_reuse[CEED_FIELD_MAX];
13280a2a6492SJeremy L Thompson 
13290a2a6492SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
133045a787f7SJeremy L Thompson     output_matrix_reuse[i].index = -1;
13310a2a6492SJeremy L Thompson   }
13320a2a6492SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1333*412e5683SJeremy L Thompson     bool         is_tensor = true;
13340a2a6492SJeremy L Thompson     CeedEvalMode eval_mode_i;
13350a2a6492SJeremy L Thompson     CeedBasis    basis_i;
13360a2a6492SJeremy L Thompson 
13370a2a6492SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode_i));
13380a2a6492SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis_i));
1339*412e5683SJeremy L Thompson     CeedCallBackend(CeedBasisIsTensor(basis_i, &is_tensor));
134045a787f7SJeremy L Thompson     for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < num_input_fields); j++) {
13410a2a6492SJeremy L Thompson       CeedEvalMode eval_mode_j;
13420a2a6492SJeremy L Thompson       CeedBasis    basis_j;
13430a2a6492SJeremy L Thompson 
13440a2a6492SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
13450a2a6492SJeremy L Thompson       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
13460a2a6492SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
13470a2a6492SJeremy L Thompson       if (basis_i == basis_j) {
13480a2a6492SJeremy L Thompson         if (is_tensor) {
134945a787f7SJeremy L Thompson           output_matrix_reuse[i].index     = j;
135045a787f7SJeremy L Thompson           output_matrix_reuse[i].is_input  = true;
135145a787f7SJeremy L Thompson           output_matrix_reuse[i].eval_mode = eval_mode_j;
13520a2a6492SJeremy L Thompson         } else {
13530a2a6492SJeremy L Thompson           // For non-tensor can only re-use with the same eval mode
13540a2a6492SJeremy L Thompson           if (eval_mode_i == eval_mode_j) {
135545a787f7SJeremy L Thompson             output_matrix_reuse[i].index     = j;
135645a787f7SJeremy L Thompson             output_matrix_reuse[i].is_input  = true;
135745a787f7SJeremy L Thompson             output_matrix_reuse[i].eval_mode = eval_mode_j;
13580a2a6492SJeremy L Thompson           }
13590a2a6492SJeremy L Thompson         }
13600a2a6492SJeremy L Thompson       }
13610a2a6492SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis_j));
13620a2a6492SJeremy L Thompson     }
136345a787f7SJeremy L Thompson     for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < i); j++) {
13640a2a6492SJeremy L Thompson       CeedEvalMode eval_mode_j;
13650a2a6492SJeremy L Thompson       CeedBasis    basis_j;
13660a2a6492SJeremy L Thompson 
13670a2a6492SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[j], &eval_mode_j));
13680a2a6492SJeremy L Thompson       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
13690a2a6492SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[j], &basis_j));
13700a2a6492SJeremy L Thompson       if (basis_i == basis_j) {
13710a2a6492SJeremy L Thompson         if (is_tensor) {
137245a787f7SJeremy L Thompson           output_matrix_reuse[i].index     = j;
137345a787f7SJeremy L Thompson           output_matrix_reuse[i].is_input  = false;
137445a787f7SJeremy L Thompson           output_matrix_reuse[i].eval_mode = eval_mode_j;
13750a2a6492SJeremy L Thompson         } else {
13760a2a6492SJeremy L Thompson           // For non-tensor can only re-use with the same eval mode
13770a2a6492SJeremy L Thompson           if (eval_mode_i == eval_mode_j) {
137845a787f7SJeremy L Thompson             output_matrix_reuse[i].index     = j;
137945a787f7SJeremy L Thompson             output_matrix_reuse[i].is_input  = false;
138045a787f7SJeremy L Thompson             output_matrix_reuse[i].eval_mode = eval_mode_j;
13810a2a6492SJeremy L Thompson           }
13820a2a6492SJeremy L Thompson         }
13830a2a6492SJeremy L Thompson       }
13840a2a6492SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis_j));
13850a2a6492SJeremy L Thompson     }
13860a2a6492SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis_i));
13870a2a6492SJeremy L Thompson   }
13880a2a6492SJeremy L Thompson 
1389ac421f39SYohann   // Initialize constants, and matrices B and G
13904b3e95d5SJeremy L Thompson   code << "\n  // Input field constants and basis data\n";
13919e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
13920a2a6492SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, i, op_input_fields[i], qf_input_fields[i], input_matrix_reuse[i], Q_1d,
1393*412e5683SJeremy L Thompson                                                               true, is_all_tensor, is_at_points, use_3d_slices));
1394920dcdc4Sjeremylt   }
13954b3e95d5SJeremy L Thompson   code << "\n  // Output field constants and basis data\n";
13969e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
13970a2a6492SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, i, op_output_fields[i], qf_output_fields[i], output_matrix_reuse[i], Q_1d,
1398*412e5683SJeremy L Thompson                                                               false, is_all_tensor, is_at_points, use_3d_slices));
13994b3e95d5SJeremy L Thompson   }
1400920dcdc4Sjeremylt 
14014b3e95d5SJeremy L Thompson   // Loop over all elements
14024b3e95d5SJeremy L Thompson   code << "\n  // Element loop\n";
1403ac421f39SYohann   code << "  __syncthreads();\n";
14049e201c85SYohann   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
14054b3e95d5SJeremy L Thompson 
1406e93651e5SJeremy L Thompson   // -- Compute minimum buffer space needed
14078b97b69aSJeremy L Thompson   CeedInt max_rstr_buffer_size = 1;
1408e93651e5SJeremy L Thompson 
14099e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1410e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
1411e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
1412e93651e5SJeremy L Thompson 
1413e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
1414e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1415e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1416*412e5683SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_all_tensor && (max_dim >= 3) ? elem_size : 1));
1417681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1418e93651e5SJeremy L Thompson   }
1419e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1420e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
1421e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
1422e93651e5SJeremy L Thompson 
1423e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
1424e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1425e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1426*412e5683SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_all_tensor && (max_dim >= 3) ? elem_size : 1));
1427681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1428e93651e5SJeremy L Thompson   }
1429e93651e5SJeremy L Thompson   code << "    // Scratch restriction buffer space\n";
1430e93651e5SJeremy L Thompson   code << "    CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n";
1431e93651e5SJeremy L Thompson 
1432e93651e5SJeremy L Thompson   // -- Determine best input field processing order
1433e93651e5SJeremy L Thompson   CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX];
1434e93651e5SJeremy L Thompson 
1435e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1436e93651e5SJeremy L Thompson     field_rstr_in_buffer[i] = -1;
1437e93651e5SJeremy L Thompson     input_field_order[i]    = -1;
1438e93651e5SJeremy L Thompson   }
1439e93651e5SJeremy L Thompson   {
1440e93651e5SJeremy L Thompson     bool    is_ordered[CEED_FIELD_MAX];
1441e93651e5SJeremy L Thompson     CeedInt curr_index = 0;
1442e93651e5SJeremy L Thompson 
1443e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
1444e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
1445e93651e5SJeremy L Thompson       CeedVector          vec_i;
1446e93651e5SJeremy L Thompson       CeedElemRestriction rstr_i;
1447e93651e5SJeremy L Thompson 
1448e93651e5SJeremy L Thompson       if (is_ordered[i]) continue;
1449e93651e5SJeremy L Thompson       field_rstr_in_buffer[i]       = i;
1450e93651e5SJeremy L Thompson       is_ordered[i]                 = true;
1451e93651e5SJeremy L Thompson       input_field_order[curr_index] = i;
1452e93651e5SJeremy L Thompson       curr_index++;
1453034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
1454e93651e5SJeremy L Thompson       if (vec_i == CEED_VECTOR_NONE) continue;  // CEED_EVAL_WEIGHT
1455e93651e5SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
1456e93651e5SJeremy L Thompson       for (CeedInt j = i + 1; j < num_input_fields; j++) {
1457e93651e5SJeremy L Thompson         CeedVector          vec_j;
1458e93651e5SJeremy L Thompson         CeedElemRestriction rstr_j;
1459e93651e5SJeremy L Thompson 
1460e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
1461e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
1462e93651e5SJeremy L Thompson         if (rstr_i == rstr_j && vec_i == vec_j) {
1463e93651e5SJeremy L Thompson           field_rstr_in_buffer[j]       = i;
1464e93651e5SJeremy L Thompson           is_ordered[j]                 = true;
1465e93651e5SJeremy L Thompson           input_field_order[curr_index] = j;
1466e93651e5SJeremy L Thompson           curr_index++;
1467e93651e5SJeremy L Thompson         }
1468681d0ea7SJeremy L Thompson         CeedCallBackend(CeedVectorDestroy(&vec_j));
1469681d0ea7SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
1470e93651e5SJeremy L Thompson       }
1471681d0ea7SJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec_i));
1472681d0ea7SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
1473e93651e5SJeremy L Thompson     }
1474e93651e5SJeremy L Thompson   }
1475e93651e5SJeremy L Thompson 
1476e93651e5SJeremy L Thompson   // -- Input restriction and basis
1477e93651e5SJeremy L Thompson   code << "\n    // -- Input field restrictions and basis actions\n";
1478e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
147959fa3f92SJeremy L Thompson     const char   *field_name;
148059fa3f92SJeremy L Thompson     const CeedInt f = input_field_order[i];
1481e93651e5SJeremy L Thompson 
148259fa3f92SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[f], &field_name));
148359fa3f92SJeremy L Thompson     code << "    // ---- Input field " << f << ": " << field_name << "\n";
1484920dcdc4Sjeremylt 
14854b3e95d5SJeremy L Thompson     // ---- Restriction
1486*412e5683SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, f, max_dim, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f],
1487*412e5683SJeremy L Thompson                                                                 Q_1d, true, is_all_tensor, is_at_points, use_3d_slices));
14889a2291e3SJeremy L Thompson 
14894b3e95d5SJeremy L Thompson     // ---- Basis action
1490*412e5683SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, f, max_dim, op_input_fields[f], qf_input_fields[f], Q_1d, true, is_all_tensor,
1491dc007f05SJeremy L Thompson                                                           is_at_points, use_3d_slices));
1492920dcdc4Sjeremylt   }
1493920dcdc4Sjeremylt 
14944b3e95d5SJeremy L Thompson   // -- Q function
1495*412e5683SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelQFunction_Cuda_gen(code, data, max_dim, max_num_points, num_input_fields, op_input_fields, qf_input_fields,
1496*412e5683SJeremy L Thompson                                                             num_output_fields, op_output_fields, qf_output_fields, qfunction_name, Q_1d,
1497*412e5683SJeremy L Thompson                                                             is_all_tensor, is_at_points, use_3d_slices));
1498ca735530SJeremy L Thompson 
14994b3e95d5SJeremy L Thompson   // -- Output basis and restriction
15004b3e95d5SJeremy L Thompson   code << "\n    // -- Output field basis action and restrictions\n";
15019e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
150259fa3f92SJeremy L Thompson     const char *field_name;
150359fa3f92SJeremy L Thompson 
150459fa3f92SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
150559fa3f92SJeremy L Thompson     code << "    // ---- Output field " << i << ": " << field_name << "\n";
1506ca735530SJeremy L Thompson 
15074b3e95d5SJeremy L Thompson     // ---- Basis action
1508*412e5683SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, i, max_dim, op_output_fields[i], qf_output_fields[i], Q_1d, false,
1509*412e5683SJeremy L Thompson                                                           is_all_tensor, is_at_points, use_3d_slices));
15109a2291e3SJeremy L Thompson 
15114b3e95d5SJeremy L Thompson     // ---- Restriction
1512*412e5683SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, i, max_dim, NULL, op_output_fields[i], qf_output_fields[i], Q_1d, false,
1513*412e5683SJeremy L Thompson                                                                 is_all_tensor, is_at_points, use_3d_slices));
1514ac421f39SYohann   }
1515241a4b83SYohann 
15164b3e95d5SJeremy L Thompson   // Close loop and function
1517241a4b83SYohann   code << "  }\n";
1518818e0025SJeremy L Thompson   code << "}\n";
1519818e0025SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n\n";
1520241a4b83SYohann 
15213a2968d6SJeremy L Thompson   // Compile
1522ddae5012SJeremy L Thompson   {
1523ddae5012SJeremy L Thompson     bool          is_compile_good = false;
1524*412e5683SJeremy L Thompson     const CeedInt T_1d            = CeedIntMax(is_all_tensor ? Q_1d : Q, data->max_P_1d);
1525ddae5012SJeremy L Thompson 
1526*412e5683SJeremy L Thompson     CeedCallBackend(CeedTryCompile_Cuda(ceed, code.str().c_str(), &is_compile_good, &data->module, 1, "OP_T_1D", T_1d));
1527ddae5012SJeremy L Thompson     if (is_compile_good) {
1528ddae5012SJeremy L Thompson       *is_good_build = true;
1529eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, operator_name.c_str(), &data->op));
1530ddae5012SJeremy L Thompson     } else {
1531ddae5012SJeremy L Thompson       *is_good_build     = false;
1532ddae5012SJeremy L Thompson       data->use_fallback = true;
1533ddae5012SJeremy L Thompson     }
1534ddae5012SJeremy L Thompson   }
15352b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
15369bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
1537c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionDestroy(&qf));
1538e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1539241a4b83SYohann }
15402a86cc9dSSebastian Grimberg 
1541ab213215SJeremy L Thompson //------------------------------------------------------------------------------
1542