xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-gen/ceed-cuda-gen-operator-build.cpp (revision 0a2a64927a7a47782a31ad89eee1373d25546e0c)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3241a4b83SYohann //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5241a4b83SYohann //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
73d576824SJeremy L Thompson 
8b8e71988SJeremy L Thompson #define CEED_DEBUG_COLOR 12
97df94212SJeremy L Thompson 
1049aac155SJeremy L Thompson #include <ceed.h>
11ec3da8bcSJed Brown #include <ceed/backend.h>
129e201c85SYohann #include <ceed/jit-tools.h>
133d576824SJeremy L Thompson #include <cuda_runtime.h>
142b730f8bSJeremy L Thompson 
15241a4b83SYohann #include <iostream>
16241a4b83SYohann #include <sstream>
17b2165e7aSSebastian Grimberg #include <string>
182b730f8bSJeremy L Thompson 
190d0321e0SJeremy L Thompson #include "../cuda-ref/ceed-cuda-ref.h"
20241a4b83SYohann #include "../cuda-shared/ceed-cuda-shared.h"
2149aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
222b730f8bSJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
232b730f8bSJeremy L Thompson #include "ceed-cuda-gen.h"
24241a4b83SYohann 
25ab213215SJeremy L Thompson //------------------------------------------------------------------------------
264b3e95d5SJeremy L Thompson // Determine type of operator
274b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
284b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelData_Cuda_gen(Ceed ceed, CeedInt num_input_fields, CeedOperatorField *op_input_fields,
294b3e95d5SJeremy L Thompson                                                 CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, CeedOperatorField *op_output_fields,
304b3e95d5SJeremy L Thompson                                                 CeedQFunctionField *qf_output_fields, CeedInt *max_P_1d, CeedInt *Q_1d, CeedInt *dim, bool *is_tensor,
314b3e95d5SJeremy L Thompson                                                 bool *use_3d_slices) {
324b3e95d5SJeremy L Thompson   // Find dim, P_1d, Q_1d
334b3e95d5SJeremy L Thompson   *max_P_1d  = 0;
344b3e95d5SJeremy L Thompson   *Q_1d      = 0;
354b3e95d5SJeremy L Thompson   *dim       = 0;
364b3e95d5SJeremy L Thompson   *is_tensor = true;
37dc007f05SJeremy L Thompson 
384b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
394b3e95d5SJeremy L Thompson     CeedBasis basis;
404b3e95d5SJeremy L Thompson 
414b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
424b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
434b3e95d5SJeremy L Thompson       bool    is_field_tensor;
444b3e95d5SJeremy L Thompson       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
454b3e95d5SJeremy L Thompson 
464b3e95d5SJeremy L Thompson       // Collect dim, P_1d, and Q_1d
474b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
484b3e95d5SJeremy L Thompson       *is_tensor = *is_tensor && is_field_tensor;
49dc007f05SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
50dc007f05SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P_1d));
514b3e95d5SJeremy L Thompson       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
524b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
534b3e95d5SJeremy L Thompson       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
544b3e95d5SJeremy L Thompson       *dim = field_dim;
55dc007f05SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
56dc007f05SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q_1d));
574b3e95d5SJeremy L Thompson       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
584b3e95d5SJeremy L Thompson       *Q_1d = field_Q_1d;
594b3e95d5SJeremy L Thompson     }
60681d0ea7SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis));
614b3e95d5SJeremy L Thompson   }
624b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
634b3e95d5SJeremy L Thompson     CeedBasis basis;
644b3e95d5SJeremy L Thompson 
654b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
664b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
674b3e95d5SJeremy L Thompson       bool    is_field_tensor;
684b3e95d5SJeremy L Thompson       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
694b3e95d5SJeremy L Thompson 
704b3e95d5SJeremy L Thompson       // Collect dim, P_1d, and Q_1d
714b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
724b3e95d5SJeremy L Thompson       *is_tensor = *is_tensor && is_field_tensor;
73dc007f05SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
74dc007f05SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P_1d));
754b3e95d5SJeremy L Thompson       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
764b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
774b3e95d5SJeremy L Thompson       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
784b3e95d5SJeremy L Thompson       *dim = field_dim;
79dc007f05SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
80dc007f05SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q_1d));
814b3e95d5SJeremy L Thompson       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
824b3e95d5SJeremy L Thompson       *Q_1d = field_Q_1d;
834b3e95d5SJeremy L Thompson     }
84681d0ea7SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis));
854b3e95d5SJeremy L Thompson   }
864b3e95d5SJeremy L Thompson 
874b3e95d5SJeremy L Thompson   // Only use 3D collocated gradient parallelization strategy when gradient is computed
884b3e95d5SJeremy L Thompson   *use_3d_slices = false;
894b3e95d5SJeremy L Thompson   if (*dim == 3) {
904b3e95d5SJeremy L Thompson     bool was_grad_found = false;
914b3e95d5SJeremy L Thompson 
924b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
934b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
944b3e95d5SJeremy L Thompson 
954b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
964b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
974b3e95d5SJeremy L Thompson         CeedBasis_Cuda_shared *basis_data;
984b3e95d5SJeremy L Thompson         CeedBasis              basis;
994b3e95d5SJeremy L Thompson 
1004b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1014b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1024b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
1034b3e95d5SJeremy L Thompson         was_grad_found = true;
104681d0ea7SJeremy L Thompson         CeedCallBackend(CeedBasisDestroy(&basis));
1054b3e95d5SJeremy L Thompson       }
1064b3e95d5SJeremy L Thompson     }
1074b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
1084b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
1094b3e95d5SJeremy L Thompson 
1104b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1114b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
1124b3e95d5SJeremy L Thompson         CeedBasis_Cuda_shared *basis_data;
1134b3e95d5SJeremy L Thompson         CeedBasis              basis;
1144b3e95d5SJeremy L Thompson 
1154b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1164b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1174b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
1184b3e95d5SJeremy L Thompson         was_grad_found = true;
119681d0ea7SJeremy L Thompson         CeedCallBackend(CeedBasisDestroy(&basis));
1204b3e95d5SJeremy L Thompson       }
1214b3e95d5SJeremy L Thompson     }
1224b3e95d5SJeremy L Thompson   }
1234b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1244b3e95d5SJeremy L Thompson }
1254b3e95d5SJeremy L Thompson 
1264b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1274b3e95d5SJeremy L Thompson // Setup fields
1284b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1294b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelFieldData_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedOperatorField op_field,
130*0a2a6492SJeremy L Thompson                                                      CeedQFunctionField qf_field, CeedInt field_reuse[3], CeedInt Q_1d, bool is_input, bool is_tensor,
131*0a2a6492SJeremy L Thompson                                                      bool is_at_points, bool use_3d_slices) {
1324b3e95d5SJeremy L Thompson   std::string            var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
133dc007f05SJeremy L Thompson   std::string            P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q";
1344b3e95d5SJeremy L Thompson   std::string            option_name = (is_input ? "inputs" : "outputs");
1354b3e95d5SJeremy L Thompson   CeedEvalMode           eval_mode   = CEED_EVAL_NONE;
1364b3e95d5SJeremy L Thompson   CeedInt                elem_size = 0, num_comp = 0, P_1d = 0;
1374b3e95d5SJeremy L Thompson   CeedElemRestriction    elem_rstr;
1384b3e95d5SJeremy L Thompson   CeedBasis_Cuda_shared *basis_data;
1394b3e95d5SJeremy L Thompson   CeedBasis              basis;
1404b3e95d5SJeremy L Thompson 
141*0a2a6492SJeremy L Thompson   // Field reuse info
142*0a2a6492SJeremy L Thompson   bool         use_previous_field = field_reuse[0] != -1;
143*0a2a6492SJeremy L Thompson   bool         reuse_input        = field_reuse[1];
144*0a2a6492SJeremy L Thompson   CeedInt      reuse_field        = field_reuse[0];
145*0a2a6492SJeremy L Thompson   CeedEvalMode reuse_mode         = (CeedEvalMode)field_reuse[2];
146*0a2a6492SJeremy L Thompson 
1474b3e95d5SJeremy L Thompson   code << "  // -- " << (is_input ? "Input" : "Output") << " field " << i << "\n";
1484b3e95d5SJeremy L Thompson 
1494b3e95d5SJeremy L Thompson   // Get field data
1504b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
1514b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
1524b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1534b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1544b3e95d5SJeremy L Thompson   }
155681d0ea7SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1564b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
1574b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
1584b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetData(basis, &basis_data));
159dc007f05SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
160dc007f05SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
1614b3e95d5SJeremy L Thompson   }
1624b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
1634b3e95d5SJeremy L Thompson 
1644b3e95d5SJeremy L Thompson   // Set field constants
1654b3e95d5SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
1664b3e95d5SJeremy L Thompson     code << "  const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n";
1674b3e95d5SJeremy L Thompson     code << "  const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n";
1684b3e95d5SJeremy L Thompson   }
1694b3e95d5SJeremy L Thompson 
1704b3e95d5SJeremy L Thompson   // Load basis data
1714b3e95d5SJeremy L Thompson   code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
1724b3e95d5SJeremy L Thompson   switch (eval_mode) {
1734b3e95d5SJeremy L Thompson     case CEED_EVAL_NONE:
1744b3e95d5SJeremy L Thompson       break;
1754b3e95d5SJeremy L Thompson     case CEED_EVAL_INTERP:
1768b97b69aSJeremy L Thompson       if (is_at_points) {
1778b97b69aSJeremy L Thompson         // AtPoints
1788b97b69aSJeremy L Thompson         if (!basis_data->d_chebyshev_interp_1d) {
1798b97b69aSJeremy L Thompson           CeedSize    interp_bytes;
1808b97b69aSJeremy L Thompson           CeedScalar *chebyshev_interp_1d;
1818b97b69aSJeremy L Thompson 
1828b97b69aSJeremy L Thompson           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
1838b97b69aSJeremy L Thompson           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
1848b97b69aSJeremy L Thompson           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
1858b97b69aSJeremy L Thompson           CeedCallCuda(CeedBasisReturnCeed(basis), cudaMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
1868b97b69aSJeremy L Thompson           CeedCallCuda(CeedBasisReturnCeed(basis),
1878b97b69aSJeremy L Thompson                        cudaMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
1888b97b69aSJeremy L Thompson           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
1898b97b69aSJeremy L Thompson         }
1908b97b69aSJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
1918b97b69aSJeremy L Thompson         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
1928b97b69aSJeremy L Thompson       } else {
1938b97b69aSJeremy L Thompson         // Standard quadrature
1944b3e95d5SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
1954b3e95d5SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_interp_1d;
1968b97b69aSJeremy L Thompson       }
197*0a2a6492SJeremy L Thompson       if (use_previous_field) {
198*0a2a6492SJeremy L Thompson         std::string reuse_var = "s_B" + ((reuse_input ? "_in_" : "_out_") + std::to_string(reuse_field));
199*0a2a6492SJeremy L Thompson 
200*0a2a6492SJeremy L Thompson         code << "  CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n";
201*0a2a6492SJeremy L Thompson       } else {
202dc007f05SJeremy L Thompson         code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
203f815fac9SJeremy L Thompson         code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
204*0a2a6492SJeremy L Thompson       }
2054b3e95d5SJeremy L Thompson       break;
2064b3e95d5SJeremy L Thompson     case CEED_EVAL_GRAD:
2078b97b69aSJeremy L Thompson       if (is_at_points) {
2088b97b69aSJeremy L Thompson         // AtPoints
2098b97b69aSJeremy L Thompson         if (!basis_data->d_chebyshev_interp_1d) {
2108b97b69aSJeremy L Thompson           CeedSize    interp_bytes;
2118b97b69aSJeremy L Thompson           CeedScalar *chebyshev_interp_1d;
2128b97b69aSJeremy L Thompson 
2138b97b69aSJeremy L Thompson           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
2148b97b69aSJeremy L Thompson           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
2158b97b69aSJeremy L Thompson           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
2168b97b69aSJeremy L Thompson           CeedCallCuda(CeedBasisReturnCeed(basis), cudaMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
2178b97b69aSJeremy L Thompson           CeedCallCuda(CeedBasisReturnCeed(basis),
2188b97b69aSJeremy L Thompson                        cudaMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
2198b97b69aSJeremy L Thompson           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
2208b97b69aSJeremy L Thompson         }
2218b97b69aSJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
2228b97b69aSJeremy L Thompson         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
2238b97b69aSJeremy L Thompson       } else {
2248b97b69aSJeremy L Thompson         // Standard quadrature
2254b3e95d5SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
2264b3e95d5SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_interp_1d;
2278b97b69aSJeremy L Thompson       }
228dc007f05SJeremy L Thompson       if (is_tensor) {
229*0a2a6492SJeremy L Thompson         if (use_previous_field) {
230*0a2a6492SJeremy L Thompson           std::string reuse_var = "s_B" + ((reuse_input ? "_in_" : "_out_") + std::to_string(reuse_field));
231*0a2a6492SJeremy L Thompson 
232*0a2a6492SJeremy L Thompson           code << "  CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n";
233*0a2a6492SJeremy L Thompson         } else {
234dc007f05SJeremy L Thompson           code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
235f815fac9SJeremy L Thompson           code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
236dc007f05SJeremy L Thompson         }
237*0a2a6492SJeremy L Thompson       }
2388b97b69aSJeremy L Thompson       if (is_at_points) break;  // No G mat for AtPoints
2394b3e95d5SJeremy L Thompson       if (use_3d_slices) {
2404b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d;
2414b3e95d5SJeremy L Thompson         else data->G.outputs[i] = basis_data->d_collo_grad_1d;
242*0a2a6492SJeremy L Thompson         if (use_previous_field && reuse_mode == CEED_EVAL_GRAD) {
243*0a2a6492SJeremy L Thompson           std::string reuse_var = "s_G" + ((reuse_input ? "_in_" : "_out_") + std::to_string(reuse_field));
244*0a2a6492SJeremy L Thompson 
245*0a2a6492SJeremy L Thompson           code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
246*0a2a6492SJeremy L Thompson         } else {
247dc007f05SJeremy L Thompson           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
248f815fac9SJeremy L Thompson           code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
249*0a2a6492SJeremy L Thompson         }
2504b3e95d5SJeremy L Thompson       } else {
2514b3e95d5SJeremy L Thompson         bool has_collo_grad = basis_data->d_collo_grad_1d;
2524b3e95d5SJeremy L Thompson 
2534b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2544b3e95d5SJeremy L Thompson         else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2554b3e95d5SJeremy L Thompson         if (has_collo_grad) {
256*0a2a6492SJeremy L Thompson           if (use_previous_field && reuse_mode == CEED_EVAL_GRAD) {
257*0a2a6492SJeremy L Thompson             std::string reuse_var = "s_G" + ((reuse_input ? "_in_" : "_out_") + std::to_string(reuse_field));
258*0a2a6492SJeremy L Thompson 
259*0a2a6492SJeremy L Thompson             code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
260*0a2a6492SJeremy L Thompson           } else {
261dc007f05SJeremy L Thompson             code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
262f815fac9SJeremy L Thompson             code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
263*0a2a6492SJeremy L Thompson           }
264*0a2a6492SJeremy L Thompson         } else {
265*0a2a6492SJeremy L Thompson           if (use_previous_field && reuse_mode == CEED_EVAL_GRAD) {
266*0a2a6492SJeremy L Thompson             std::string reuse_var = "s_G" + ((reuse_input ? "_in_" : "_out_") + std::to_string(reuse_field));
267*0a2a6492SJeremy L Thompson 
268*0a2a6492SJeremy L Thompson             code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
2694b3e95d5SJeremy L Thompson           } else {
270dc007f05SJeremy L Thompson             code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << P_name << "*" << Q_name << (is_tensor ? "" : "*dim") << "];\n";
271dc007f05SJeremy L Thompson             code << "  LoadMatrix<" << P_name << ", " << Q_name << (is_tensor ? "" : "*dim") << ">(data, G." << option_name << "[" << i << "], s_G"
272dc007f05SJeremy L Thompson                  << var_suffix << ");\n";
2734b3e95d5SJeremy L Thompson           }
2744b3e95d5SJeremy L Thompson         }
275*0a2a6492SJeremy L Thompson       }
2764b3e95d5SJeremy L Thompson       break;
2774b3e95d5SJeremy L Thompson     case CEED_EVAL_WEIGHT:
2784b3e95d5SJeremy L Thompson       break;  // No action
2794b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
2804b3e95d5SJeremy L Thompson     case CEED_EVAL_DIV:
2814b3e95d5SJeremy L Thompson     case CEED_EVAL_CURL:
2824b3e95d5SJeremy L Thompson       break;  // TODO: Not implemented
2834b3e95d5SJeremy L Thompson               // LCOV_EXCL_STOP
2844b3e95d5SJeremy L Thompson   }
2858b97b69aSJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
2864b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2874b3e95d5SJeremy L Thompson }
2884b3e95d5SJeremy L Thompson 
2894b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
2904b3e95d5SJeremy L Thompson // Restriction
2914b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
2924b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelRestriction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedInt dim,
293e93651e5SJeremy L Thompson                                                        CeedInt field_input_buffer[], CeedOperatorField op_field, CeedQFunctionField qf_field,
294dc007f05SJeremy L Thompson                                                        CeedInt Q_1d, bool is_input, bool is_tensor, bool is_at_points, bool use_3d_slices) {
2954b3e95d5SJeremy L Thompson   std::string               var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
296dc007f05SJeremy L Thompson   std::string               P_name     = (is_tensor ? "P_1d" : "P") + var_suffix;
2974b3e95d5SJeremy L Thompson   CeedEvalMode              eval_mode  = CEED_EVAL_NONE;
2984b3e95d5SJeremy L Thompson   CeedInt                   elem_size = 0, num_comp = 0, P_1d = 0;
2994b3e95d5SJeremy L Thompson   CeedSize                  l_size;
300f815fac9SJeremy L Thompson   CeedRestrictionType       rstr_type = CEED_RESTRICTION_STANDARD;
3014b3e95d5SJeremy L Thompson   CeedElemRestriction_Cuda *rstr_data;
3024b3e95d5SJeremy L Thompson   CeedElemRestriction       elem_rstr;
3034b3e95d5SJeremy L Thompson   CeedBasis                 basis;
3044b3e95d5SJeremy L Thompson 
3054b3e95d5SJeremy L Thompson   // Get field data
3064b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
3074b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
308f815fac9SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
3094b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
3104b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
3114b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
3124b3e95d5SJeremy L Thompson   }
3134b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
3144b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
315dc007f05SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
316dc007f05SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
3174b3e95d5SJeremy L Thompson   }
318681d0ea7SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
3194b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
3204b3e95d5SJeremy L Thompson 
3214b3e95d5SJeremy L Thompson   // Restriction
3224b3e95d5SJeremy L Thompson   if (is_input) {
3234b3e95d5SJeremy L Thompson     // Input
324e93651e5SJeremy L Thompson     if (field_input_buffer[i] != i) {
325e93651e5SJeremy L Thompson       std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]);
326e93651e5SJeremy L Thompson 
327e93651e5SJeremy L Thompson       // Restriction was already done for previous input
328e93651e5SJeremy L Thompson       code << "    CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n";
3298b97b69aSJeremy L Thompson     } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices && is_at_points)) {
3308b97b69aSJeremy L Thompson       if (eval_mode == CEED_EVAL_NONE && rstr_type != CEED_RESTRICTION_POINTS) {
331e93651e5SJeremy L Thompson         // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated
3324b3e95d5SJeremy L Thompson         code << "    CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n";
3338b97b69aSJeremy L Thompson       } else if (rstr_type != CEED_RESTRICTION_POINTS) {
334e93651e5SJeremy L Thompson         // Otherwise we're using the scratch space
335e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
336e93651e5SJeremy L Thompson       }
337f815fac9SJeremy L Thompson       switch (rstr_type) {
338f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STANDARD: {
3394b3e95d5SJeremy L Thompson           CeedInt comp_stride;
3404b3e95d5SJeremy L Thompson 
3414b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
3424b3e95d5SJeremy L Thompson           code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
3434b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3444b3e95d5SJeremy L Thompson           code << "    // CompStride: " << comp_stride << "\n";
3454b3e95d5SJeremy L Thompson           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
346dc007f05SJeremy L Thompson           code << "    ReadLVecStandard" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name
347dc007f05SJeremy L Thompson                << ">(data, l_size" << var_suffix << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n";
348f815fac9SJeremy L Thompson           break;
349f815fac9SJeremy L Thompson         }
350f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STRIDED: {
3514b3e95d5SJeremy L Thompson           bool    has_backend_strides;
3524b3e95d5SJeremy L Thompson           CeedInt num_elem;
3534b3e95d5SJeremy L Thompson 
3544b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
3554b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
3564b3e95d5SJeremy L Thompson           CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
3574b3e95d5SJeremy L Thompson 
3584b3e95d5SJeremy L Thompson           if (!has_backend_strides) {
3594b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
3604b3e95d5SJeremy L Thompson           }
3614b3e95d5SJeremy L Thompson           code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
362dc007f05SJeremy L Thompson           code << "    ReadLVecStrided" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", "
363dc007f05SJeremy L Thompson                << strides[1] << ", " << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n";
364f815fac9SJeremy L Thompson           break;
365f815fac9SJeremy L Thompson         }
3668b97b69aSJeremy L Thompson         case CEED_RESTRICTION_POINTS: {
3678b97b69aSJeremy L Thompson           CeedInt comp_stride;
3688b97b69aSJeremy L Thompson 
3698b97b69aSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3708b97b69aSJeremy L Thompson           code << "    const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
3718b97b69aSJeremy L Thompson           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
3728b97b69aSJeremy L Thompson           break;
3738b97b69aSJeremy L Thompson         }
374f815fac9SJeremy L Thompson         // LCOV_EXCL_START
375f815fac9SJeremy L Thompson         case CEED_RESTRICTION_ORIENTED:
376f815fac9SJeremy L Thompson         case CEED_RESTRICTION_CURL_ORIENTED:
377f815fac9SJeremy L Thompson           break;  // TODO: Not implemented
378f815fac9SJeremy L Thompson                   // LCOV_EXCL_STOP
3794b3e95d5SJeremy L Thompson       }
3804b3e95d5SJeremy L Thompson     }
3814b3e95d5SJeremy L Thompson   } else {
3824b3e95d5SJeremy L Thompson     // Output
383f815fac9SJeremy L Thompson     switch (rstr_type) {
384f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STANDARD: {
3854b3e95d5SJeremy L Thompson         CeedInt comp_stride;
3864b3e95d5SJeremy L Thompson 
3874b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
3884b3e95d5SJeremy L Thompson         code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
3894b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3904b3e95d5SJeremy L Thompson         code << "    // CompStride: " << comp_stride << "\n";
3914b3e95d5SJeremy L Thompson         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
392dc007f05SJeremy L Thompson         code << "    WriteLVecStandard" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name
393dc007f05SJeremy L Thompson              << ">(data, l_size" << var_suffix << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n";
394f815fac9SJeremy L Thompson         break;
395f815fac9SJeremy L Thompson       }
396f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STRIDED: {
3974b3e95d5SJeremy L Thompson         bool    has_backend_strides;
3984b3e95d5SJeremy L Thompson         CeedInt num_elem;
3994b3e95d5SJeremy L Thompson 
4004b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
4014b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
4024b3e95d5SJeremy L Thompson         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
4034b3e95d5SJeremy L Thompson 
4044b3e95d5SJeremy L Thompson         if (!has_backend_strides) {
4054b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
4064b3e95d5SJeremy L Thompson         }
4074b3e95d5SJeremy L Thompson         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
408dc007f05SJeremy L Thompson         code << "    WriteLVecStrided" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", "
409dc007f05SJeremy L Thompson              << strides[1] << ", " << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n";
410f815fac9SJeremy L Thompson         break;
411f815fac9SJeremy L Thompson       }
4128b97b69aSJeremy L Thompson       case CEED_RESTRICTION_POINTS:
4138b97b69aSJeremy L Thompson         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
4148b97b69aSJeremy L Thompson         break;
415f815fac9SJeremy L Thompson       // LCOV_EXCL_START
416f815fac9SJeremy L Thompson       case CEED_RESTRICTION_ORIENTED:
417f815fac9SJeremy L Thompson       case CEED_RESTRICTION_CURL_ORIENTED:
418f815fac9SJeremy L Thompson         break;  // TODO: Not implemented
419f815fac9SJeremy L Thompson                 // LCOV_EXCL_STOP
4204b3e95d5SJeremy L Thompson     }
4214b3e95d5SJeremy L Thompson   }
422681d0ea7SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
4234b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4244b3e95d5SJeremy L Thompson }
4254b3e95d5SJeremy L Thompson 
4264b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4274b3e95d5SJeremy L Thompson // Basis
4284b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4294b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelBasis_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedInt dim,
430dc007f05SJeremy L Thompson                                                  CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool is_tensor,
4318b97b69aSJeremy L Thompson                                                  bool is_at_points, bool use_3d_slices) {
4324b3e95d5SJeremy L Thompson   std::string         var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
433dc007f05SJeremy L Thompson   std::string         P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q";
4344b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
4354b3e95d5SJeremy L Thompson   CeedInt             elem_size = 0, num_comp = 0, P_1d = 0;
4364b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
4374b3e95d5SJeremy L Thompson   CeedBasis           basis;
4384b3e95d5SJeremy L Thompson 
4394b3e95d5SJeremy L Thompson   // Get field data
4404b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
4414b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
4424b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
4434b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
4444b3e95d5SJeremy L Thompson   }
445681d0ea7SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
4464b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
4474b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
448dc007f05SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
449dc007f05SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
4504b3e95d5SJeremy L Thompson   }
4514b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
4524b3e95d5SJeremy L Thompson 
4534b3e95d5SJeremy L Thompson   // Basis
4544b3e95d5SJeremy L Thompson   code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
4554b3e95d5SJeremy L Thompson   if (is_input) {
4564b3e95d5SJeremy L Thompson     switch (eval_mode) {
4574b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
4588b97b69aSJeremy L Thompson         if (!use_3d_slices && !is_at_points) {
4594b3e95d5SJeremy L Thompson           code << "    CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n";
4604b3e95d5SJeremy L Thompson         }
4614b3e95d5SJeremy L Thompson         break;
4624b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
4638b97b69aSJeremy L Thompson         if (is_at_points) {
464dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
465dc007f05SJeremy L Thompson 
466dc007f05SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
467dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
468dc007f05SJeremy L Thompson                << var_suffix << ", r_c" << var_suffix << ");\n";
4698b97b69aSJeremy L Thompson         } else {
470dc007f05SJeremy L Thompson           std::string function_name = is_tensor ? ((dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d") : "InterpNonTensor";
471dc007f05SJeremy L Thompson 
472dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
473dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
474dc007f05SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
4758b97b69aSJeremy L Thompson         }
4764b3e95d5SJeremy L Thompson         break;
4774b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
4788b97b69aSJeremy L Thompson         if (is_at_points) {
479dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
480dc007f05SJeremy L Thompson 
481dc007f05SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
482dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
483dc007f05SJeremy L Thompson                << var_suffix << ", r_c" << var_suffix << ");\n";
4848b97b69aSJeremy L Thompson         } else if (use_3d_slices) {
485dc007f05SJeremy L Thompson           std::string function_name = (dim > 1 ? "InterpTensor" : "Interp") + std::to_string(dim) + "d";
486dc007f05SJeremy L Thompson 
4874b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
488dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
489dc007f05SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
490dc007f05SJeremy L Thompson         } else if (is_tensor) {
491dc007f05SJeremy L Thompson           bool        is_collocated = dim == 3 && Q_1d >= P_1d;
492dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "Grad" : (is_collocated ? "GradTensorCollocated" : "GradTensor")) + std::to_string(dim) + "d";
493dc007f05SJeremy L Thompson 
494dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << (dim >= 3 ? Q_name : "1") << "];\n";
495dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
496dc007f05SJeremy L Thompson                << var_suffix << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
4974b3e95d5SJeremy L Thompson         } else {
498dc007f05SJeremy L Thompson           std::string function_name = "GradNonTensor";
499dc007f05SJeremy L Thompson 
500dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
501dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", dim, " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix
502dc007f05SJeremy L Thompson                << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
5034b3e95d5SJeremy L Thompson         }
5044b3e95d5SJeremy L Thompson         break;
5054b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
5068b97b69aSJeremy L Thompson         if (is_at_points) {
5078b97b69aSJeremy L Thompson           code << "    // Nothing to do AtPoints\n";
5088b97b69aSJeremy L Thompson         } else {
5094b3e95d5SJeremy L Thompson           CeedBasis_Cuda_shared *basis_data;
510dc007f05SJeremy L Thompson           std::string            function_name = is_tensor ? ((dim == 1 ? "Weight" : "WeightTensor") + std::to_string(dim) + "d") : "WeightNonTensor";
5114b3e95d5SJeremy L Thompson 
512dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
5134b3e95d5SJeremy L Thompson           CeedCallBackend(CeedBasisGetData(basis, &basis_data));
5144b3e95d5SJeremy L Thompson           data->W = basis_data->d_q_weight_1d;
515dc007f05SJeremy L Thompson           code << "    " << function_name << "<" << Q_name << ">(data, W, r_q" << var_suffix << ");\n";
5168b97b69aSJeremy L Thompson         }
5174b3e95d5SJeremy L Thompson         break;
5184b3e95d5SJeremy L Thompson       }
5194b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
5204b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
5214b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
5224b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
5234b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
5244b3e95d5SJeremy L Thompson     }
5254b3e95d5SJeremy L Thompson   } else {
5264b3e95d5SJeremy L Thompson     switch (eval_mode) {
5274b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
5284b3e95d5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n";
5294b3e95d5SJeremy L Thompson         break;  // No action
5304b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
531e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
5328b97b69aSJeremy L Thompson         if (is_at_points) {
533dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
534dc007f05SJeremy L Thompson 
535dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_c" << var_suffix << ", s_B"
536dc007f05SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5378b97b69aSJeremy L Thompson         } else {
538dc007f05SJeremy L Thompson           std::string function_name =
539dc007f05SJeremy L Thompson               is_tensor ? ((dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d") : "InterpTransposeNonTensor";
540dc007f05SJeremy L Thompson 
541dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B"
542dc007f05SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5438b97b69aSJeremy L Thompson         }
5444b3e95d5SJeremy L Thompson         break;
5454b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
546e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
5478b97b69aSJeremy L Thompson         if (is_at_points) {
548dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
549dc007f05SJeremy L Thompson 
550dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_c" << var_suffix << ", s_B"
551dc007f05SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5528b97b69aSJeremy L Thompson         } else if (use_3d_slices) {
553dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
554dc007f05SJeremy L Thompson 
555dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B"
556dc007f05SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
557dc007f05SJeremy L Thompson         } else if (is_tensor) {
558dc007f05SJeremy L Thompson           bool        is_collocated = dim == 3 && Q_1d >= P_1d;
559dc007f05SJeremy L Thompson           std::string function_name =
560dc007f05SJeremy L Thompson               (dim == 1 ? "GradTranspose" : (is_collocated ? "GradTransposeTensorCollocated" : "GradTransposeTensor")) + std::to_string(dim) + "d";
561dc007f05SJeremy L Thompson 
562dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B"
563dc007f05SJeremy L Thompson                << var_suffix << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
5644b3e95d5SJeremy L Thompson         } else {
565dc007f05SJeremy L Thompson           std::string function_name = "GradTransposeNonTensor";
566dc007f05SJeremy L Thompson 
567dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", dim, " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix
568dc007f05SJeremy L Thompson                << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
5694b3e95d5SJeremy L Thompson         }
5704b3e95d5SJeremy L Thompson         break;
5714b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
5724b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT:
5734b3e95d5SJeremy L Thompson         break;  // Should not occur
5744b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
5754b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
5764b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
5774b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
5784b3e95d5SJeremy L Thompson     }
5794b3e95d5SJeremy L Thompson   }
580681d0ea7SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
5814b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5824b3e95d5SJeremy L Thompson }
5834b3e95d5SJeremy L Thompson 
5844b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
5854b3e95d5SJeremy L Thompson // QFunction
5864b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
5878b97b69aSJeremy L Thompson static int CeedOperatorBuildKernelQFunction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt dim, CeedInt max_num_points,
5888b97b69aSJeremy L Thompson                                                      CeedInt num_input_fields, CeedOperatorField *op_input_fields,
5898b97b69aSJeremy L Thompson                                                      CeedQFunctionField *qf_input_fields, CeedInt num_output_fields,
5908b97b69aSJeremy L Thompson                                                      CeedOperatorField *op_output_fields, CeedQFunctionField *qf_output_fields,
591dc007f05SJeremy L Thompson                                                      std::string qfunction_name, CeedInt Q_1d, bool is_tensor, bool is_at_points,
592dc007f05SJeremy L Thompson                                                      bool use_3d_slices) {
593dc007f05SJeremy L Thompson   std::string         Q_name    = is_tensor ? "Q_1d" : "Q";
5944b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
5954b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
5964b3e95d5SJeremy L Thompson 
5978b97b69aSJeremy L Thompson   // Setup output arrays
5984b3e95d5SJeremy L Thompson   code << "\n    // -- Output field setup\n";
5994b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
6004b3e95d5SJeremy L Thompson     std::string var_suffix = "_out_" + std::to_string(i);
6014b3e95d5SJeremy L Thompson 
6024b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
6034b3e95d5SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
6048b97b69aSJeremy L Thompson     switch (eval_mode) {
6058b97b69aSJeremy L Thompson       case CEED_EVAL_NONE:
6068b97b69aSJeremy L Thompson         if (is_at_points) {
6078b97b69aSJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "];\n";
6088b97b69aSJeremy L Thompson         } else {
609dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
6104b3e95d5SJeremy L Thompson         }
6118b97b69aSJeremy L Thompson         break;
6128b97b69aSJeremy L Thompson       case CEED_EVAL_INTERP:
6138b97b69aSJeremy L Thompson         if (is_at_points) {
6148b97b69aSJeremy L Thompson           // Accumulator for point data
615dc007f05SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
616dc007f05SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "; i++) {\n";
6178b97b69aSJeremy L Thompson           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
6188b97b69aSJeremy L Thompson           code << "    }\n";
6198b97b69aSJeremy L Thompson         } else {
620dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
6218b97b69aSJeremy L Thompson         }
6228b97b69aSJeremy L Thompson         break;
6238b97b69aSJeremy L Thompson       case CEED_EVAL_GRAD:
6248b97b69aSJeremy L Thompson         if (is_at_points) {
6258b97b69aSJeremy L Thompson           // Accumulator for point data
626dc007f05SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "*dim];\n";
627dc007f05SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "; i++) {\n";
6288b97b69aSJeremy L Thompson           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
6298b97b69aSJeremy L Thompson           code << "    }\n";
6308b97b69aSJeremy L Thompson         } else if (use_3d_slices) {
6314b3e95d5SJeremy L Thompson           // Accumulator for gradient slices
6324b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
6334b3e95d5SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n";
6344b3e95d5SJeremy L Thompson           code << "      r_q" << var_suffix << "[i] = 0.0;\n";
6354b3e95d5SJeremy L Thompson           code << "    }\n";
6364b3e95d5SJeremy L Thompson         } else {
637dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
6384b3e95d5SJeremy L Thompson         }
6398b97b69aSJeremy L Thompson         break;
6408b97b69aSJeremy L Thompson       case CEED_EVAL_WEIGHT:
6418b97b69aSJeremy L Thompson         break;
6428b97b69aSJeremy L Thompson         // LCOV_EXCL_START
6438b97b69aSJeremy L Thompson       case CEED_EVAL_DIV:
6448b97b69aSJeremy L Thompson       case CEED_EVAL_CURL:
6458b97b69aSJeremy L Thompson         break;  // TODO: Not implemented
6468b97b69aSJeremy L Thompson                 // LCOV_EXCL_STOP
6474b3e95d5SJeremy L Thompson     }
6484b3e95d5SJeremy L Thompson   }
6494b3e95d5SJeremy L Thompson 
6508b97b69aSJeremy L Thompson   if (is_at_points) {
6518b97b69aSJeremy L Thompson     // We need to handle batches of points
6528b97b69aSJeremy L Thompson     code << "\n    // Note: Using batches of points\n";
6538b97b69aSJeremy L Thompson     code << "    const CeedInt point_loop_bound = (blockDim.x * blockDim.y) * ceil(1.0 * max_num_points / (blockDim.x * blockDim.y));\n\n";
6548b97b69aSJeremy L Thompson     code << "    #pragma unroll\n";
6558b97b69aSJeremy L Thompson     code << "    for (CeedInt i = threadIdx.x + threadIdx.y * blockDim.x; i < point_loop_bound; i += blockDim.x * blockDim.y) {\n";
6568b97b69aSJeremy L Thompson     code << "      const CeedInt p = i % max_num_points;\n\n";
6578b97b69aSJeremy L Thompson 
6588b97b69aSJeremy L Thompson     code << "      // -- Coordinates\n";
6598b97b69aSJeremy L Thompson     code << "      CeedScalar r_x[dim];\n";
6608b97b69aSJeremy L Thompson     code << "      ReadPoint<dim, coords_comp_stride, max_num_points>(data, elem, p, max_num_points, points.indices, points.coords, r_x);\n\n";
6618b97b69aSJeremy L Thompson 
6628b97b69aSJeremy L Thompson     code << "      // -- Input fields\n";
6638b97b69aSJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
6648b97b69aSJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
6658b97b69aSJeremy L Thompson 
6668b97b69aSJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
6678b97b69aSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
6688b97b69aSJeremy L Thompson       // Basis action
6698b97b69aSJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
6708b97b69aSJeremy L Thompson       switch (eval_mode) {
6718b97b69aSJeremy L Thompson         case CEED_EVAL_NONE:
6728b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
6738b97b69aSJeremy L Thompson           code << "      ReadPoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
6748b97b69aSJeremy L Thompson                << ", max_num_points>(data, elem, p, max_num_points, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
6758b97b69aSJeremy L Thompson           break;
6768b97b69aSJeremy L Thompson         case CEED_EVAL_INTERP:
6778b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
6788b97b69aSJeremy L Thompson           code << "      InterpAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix
6798b97b69aSJeremy L Thompson                << ", r_x, r_s" << var_suffix << ");\n";
6808b97b69aSJeremy L Thompson           break;
6818b97b69aSJeremy L Thompson         case CEED_EVAL_GRAD:
6828b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
6838b97b69aSJeremy L Thompson           code << "      GradAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix
6848b97b69aSJeremy L Thompson                << ", r_x, r_s" << var_suffix << ");\n";
6858b97b69aSJeremy L Thompson           break;
6868b97b69aSJeremy L Thompson         case CEED_EVAL_WEIGHT:
6878b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
6888b97b69aSJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = 1.0;\n";
6898b97b69aSJeremy L Thompson           break;
6908b97b69aSJeremy L Thompson           // LCOV_EXCL_START
6918b97b69aSJeremy L Thompson         case CEED_EVAL_DIV:
6928b97b69aSJeremy L Thompson         case CEED_EVAL_CURL:
6938b97b69aSJeremy L Thompson           break;  // TODO: Not implemented
6948b97b69aSJeremy L Thompson                   // LCOV_EXCL_STOP
6958b97b69aSJeremy L Thompson       }
6968b97b69aSJeremy L Thompson     }
6978b97b69aSJeremy L Thompson     code << "\n      // -- Output fields\n";
6988b97b69aSJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
6998b97b69aSJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
7008b97b69aSJeremy L Thompson 
7018b97b69aSJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
7028b97b69aSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
7038b97b69aSJeremy L Thompson       // Basis action
7048b97b69aSJeremy L Thompson       switch (eval_mode) {
7058b97b69aSJeremy L Thompson         case CEED_EVAL_NONE:
7068b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7078b97b69aSJeremy L Thompson           break;
7088b97b69aSJeremy L Thompson         case CEED_EVAL_INTERP:
7098b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7108b97b69aSJeremy L Thompson           break;
7118b97b69aSJeremy L Thompson         case CEED_EVAL_GRAD:
7128b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
7138b97b69aSJeremy L Thompson           break;
7148b97b69aSJeremy L Thompson           // LCOV_EXCL_START
7158b97b69aSJeremy L Thompson         case CEED_EVAL_WEIGHT:
7168b97b69aSJeremy L Thompson           break;  // Should not occur
7178b97b69aSJeremy L Thompson         case CEED_EVAL_DIV:
7188b97b69aSJeremy L Thompson         case CEED_EVAL_CURL:
7198b97b69aSJeremy L Thompson           break;  // TODO: Not implemented
7208b97b69aSJeremy L Thompson                   // LCOV_EXCL_STOP
7218b97b69aSJeremy L Thompson       }
7228b97b69aSJeremy L Thompson     }
7238b97b69aSJeremy L Thompson 
7248b97b69aSJeremy L Thompson   } else if (use_3d_slices) {
7254b3e95d5SJeremy L Thompson     // We treat quadrature points per slice in 3d to save registers
7264b3e95d5SJeremy L Thompson     code << "\n    // Note: Using planes of 3D elements\n";
7274b3e95d5SJeremy L Thompson     code << "    #pragma unroll\n";
7284b3e95d5SJeremy L Thompson     code << "    for (CeedInt q = 0; q < " << Q_name << "; q++) {\n";
7294b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
7304b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
7314b3e95d5SJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
7324b3e95d5SJeremy L Thompson 
7334b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
7344b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
7354b3e95d5SJeremy L Thompson       // Basis action
7364b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
7374b3e95d5SJeremy L Thompson       switch (eval_mode) {
7384b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
7394b3e95d5SJeremy L Thompson           bool is_strided;
7404b3e95d5SJeremy L Thompson 
7414b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7424b3e95d5SJeremy L Thompson 
7434b3e95d5SJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
7444b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
7454b3e95d5SJeremy L Thompson           if (is_strided) {
7464b3e95d5SJeremy L Thompson             bool    has_backend_strides;
7474b3e95d5SJeremy L Thompson             CeedInt num_elem, elem_size;
7484b3e95d5SJeremy L Thompson 
7494b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
7504b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
7514b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
7524b3e95d5SJeremy L Thompson             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
7534b3e95d5SJeremy L Thompson 
7544b3e95d5SJeremy L Thompson             if (!has_backend_strides) {
7554b3e95d5SJeremy L Thompson               CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
7564b3e95d5SJeremy L Thompson             }
7574b3e95d5SJeremy L Thompson             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
758f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStrided3d<num_comp" << var_suffix << ", " << Q_name << ", " << strides[0] << ", " << strides[1] << ", "
7594b3e95d5SJeremy L Thompson                  << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n";
7604b3e95d5SJeremy L Thompson           } else {
7614b3e95d5SJeremy L Thompson             CeedSize                  l_size = 0;
7624b3e95d5SJeremy L Thompson             CeedInt                   comp_stride;
7634b3e95d5SJeremy L Thompson             CeedElemRestriction_Cuda *rstr_data;
7644b3e95d5SJeremy L Thompson 
7654b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
7664b3e95d5SJeremy L Thompson             code << "      const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
7674b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
7684b3e95d5SJeremy L Thompson             code << "      // CompStride: " << comp_stride << "\n";
7694b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
7704b3e95d5SJeremy L Thompson             data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
771f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStandard3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix
7724b3e95d5SJeremy L Thompson                  << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
7734b3e95d5SJeremy L Thompson           }
774681d0ea7SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
7754b3e95d5SJeremy L Thompson           break;
7764b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
7774b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7784b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n";
7794b3e95d5SJeremy L Thompson           code << "        r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n";
7804b3e95d5SJeremy L Thompson           code << "      }\n";
7814b3e95d5SJeremy L Thompson           break;
7824b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
7834b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
784f815fac9SJeremy L Thompson           code << "      GradColloSlice3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_q" << var_suffix << ", s_G" << var_suffix
785f815fac9SJeremy L Thompson                << ", r_s" << var_suffix << ");\n";
7864b3e95d5SJeremy L Thompson           break;
7874b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
7884b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
7894b3e95d5SJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n";
7908b97b69aSJeremy L Thompson           break;
7914b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
7924b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
7934b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
7944b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
7954b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
7964b3e95d5SJeremy L Thompson       }
7974b3e95d5SJeremy L Thompson     }
7984b3e95d5SJeremy L Thompson     code << "\n      // -- Output fields\n";
7994b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
8004b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
8014b3e95d5SJeremy L Thompson 
8024b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
8034b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
8044b3e95d5SJeremy L Thompson       // Basis action
8054b3e95d5SJeremy L Thompson       switch (eval_mode) {
8064b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
8074b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
8088b97b69aSJeremy L Thompson           break;
8094b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
8104b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
8114b3e95d5SJeremy L Thompson           break;
8124b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
8134b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
8144b3e95d5SJeremy L Thompson           break;
8154b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
8164b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
8174b3e95d5SJeremy L Thompson           break;  // Should not occur
8184b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
8194b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
8204b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
8214b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
8224b3e95d5SJeremy L Thompson       }
8234b3e95d5SJeremy L Thompson     }
8244b3e95d5SJeremy L Thompson   } else {
8254b3e95d5SJeremy L Thompson     code << "\n    // Note: Using full elements\n";
8264b3e95d5SJeremy L Thompson     code << "    {\n";
8274b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
8284b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
8294b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
8304b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n";
8314b3e95d5SJeremy L Thompson     }
8324b3e95d5SJeremy L Thompson     code << "      // -- Output fields\n";
8334b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
8344b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
8354b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n";
8364b3e95d5SJeremy L Thompson     }
8374b3e95d5SJeremy L Thompson   }
8384b3e95d5SJeremy L Thompson 
8394b3e95d5SJeremy L Thompson   // Input and output buffers
8404b3e95d5SJeremy L Thompson   code << "\n      // -- QFunction inputs and outputs\n";
8414b3e95d5SJeremy L Thompson   code << "      // ---- Inputs\n";
8424b3e95d5SJeremy L Thompson   code << "      CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n";
8434b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
8444b3e95d5SJeremy L Thompson     code << "      // ------ Input field " << i << "\n";
8454b3e95d5SJeremy L Thompson     code << "      inputs[" << i << "] = r_s_in_" << i << ";\n";
8464b3e95d5SJeremy L Thompson   }
8474b3e95d5SJeremy L Thompson   code << "      // ---- Outputs\n";
8484b3e95d5SJeremy L Thompson   code << "      CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n";
8494b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
8504b3e95d5SJeremy L Thompson     code << "      // ------ Output field " << i << "\n";
8514b3e95d5SJeremy L Thompson     code << "      outputs[" << i << "] = r_s_out_" << i << ";\n";
8524b3e95d5SJeremy L Thompson   }
8534b3e95d5SJeremy L Thompson 
8544b3e95d5SJeremy L Thompson   // Apply QFunction
8554b3e95d5SJeremy L Thompson   code << "\n      // -- Apply QFunction\n";
8564b3e95d5SJeremy L Thompson   code << "      " << qfunction_name << "(ctx, ";
857dc007f05SJeremy L Thompson   if (dim != 3 || is_at_points || use_3d_slices || !is_tensor) {
8584b3e95d5SJeremy L Thompson     code << "1";
8594b3e95d5SJeremy L Thompson   } else {
860dc007f05SJeremy L Thompson     code << Q_name;
8614b3e95d5SJeremy L Thompson   }
8624b3e95d5SJeremy L Thompson   code << ", inputs, outputs);\n";
8634b3e95d5SJeremy L Thompson 
8648b97b69aSJeremy L Thompson   if (is_at_points) {
8658b97b69aSJeremy L Thompson     // Map back to coefficients
8668b97b69aSJeremy L Thompson     code << "\n      // -- Output fields\n";
8678b97b69aSJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
8688b97b69aSJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
8698b97b69aSJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
8708b97b69aSJeremy L Thompson 
8718b97b69aSJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
8728b97b69aSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
8738b97b69aSJeremy L Thompson       // Basis action
8748b97b69aSJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
8758b97b69aSJeremy L Thompson       switch (eval_mode) {
8768b97b69aSJeremy L Thompson         case CEED_EVAL_NONE: {
8778b97b69aSJeremy L Thompson           CeedInt             comp_stride;
8788b97b69aSJeremy L Thompson           CeedElemRestriction elem_rstr;
8798b97b69aSJeremy L Thompson 
8808b97b69aSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
8818b97b69aSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
8828b97b69aSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
8838b97b69aSJeremy L Thompson           code << "      const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
8848b97b69aSJeremy L Thompson           code << "      WritePoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
8858b97b69aSJeremy L Thompson                << ", max_num_points>(data, elem, i, points.num_per_elem[elem], indices.outputs[" << i << "]"
8868b97b69aSJeremy L Thompson                << ", r_s" << var_suffix << ", d" << var_suffix << ");\n";
8878b97b69aSJeremy L Thompson           break;
8888b97b69aSJeremy L Thompson         }
8898b97b69aSJeremy L Thompson         case CEED_EVAL_INTERP:
8908b97b69aSJeremy L Thompson           code << "      if (i >= points.num_per_elem[elem]) {\n";
8918b97b69aSJeremy L Thompson           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) r_s" << var_suffix << "[j] = 0.0;\n";
8928b97b69aSJeremy L Thompson           code << "      }\n";
8938b97b69aSJeremy L Thompson           code << "      InterpTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s"
8948b97b69aSJeremy L Thompson                << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
8958b97b69aSJeremy L Thompson           break;
8968b97b69aSJeremy L Thompson         case CEED_EVAL_GRAD:
8978b97b69aSJeremy L Thompson           code << "      if (i >= points.num_per_elem[elem]) {\n";
8988b97b69aSJeremy L Thompson           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "*dim; j++) r_s" << var_suffix << "[j] = 0.0;\n";
8998b97b69aSJeremy L Thompson           code << "      }\n";
9008b97b69aSJeremy L Thompson           code << "      GradTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s"
9018b97b69aSJeremy L Thompson                << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
9028b97b69aSJeremy L Thompson           break;
9038b97b69aSJeremy L Thompson           // LCOV_EXCL_START
9048b97b69aSJeremy L Thompson         case CEED_EVAL_WEIGHT:
9058b97b69aSJeremy L Thompson           break;  // Should not occur
9068b97b69aSJeremy L Thompson         case CEED_EVAL_DIV:
9078b97b69aSJeremy L Thompson         case CEED_EVAL_CURL:
9088b97b69aSJeremy L Thompson           break;  // TODO: Not implemented
9098b97b69aSJeremy L Thompson                   // LCOV_EXCL_STOP
9108b97b69aSJeremy L Thompson       }
9118b97b69aSJeremy L Thompson     }
9128b97b69aSJeremy L Thompson   } else if (use_3d_slices) {
9134b3e95d5SJeremy L Thompson     // Copy or apply transpose grad, if needed
9148b97b69aSJeremy L Thompson     code << "\n      // -- Output fields\n";
9154b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
9164b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
9174b3e95d5SJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
9184b3e95d5SJeremy L Thompson 
9194b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
9204b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
9214b3e95d5SJeremy L Thompson       // Basis action
9224b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
9234b3e95d5SJeremy L Thompson       switch (eval_mode) {
9244b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
9254b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
9264b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
9274b3e95d5SJeremy L Thompson           code << "      }\n";
9288b97b69aSJeremy L Thompson           break;
9294b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
9304b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
9314b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
9324b3e95d5SJeremy L Thompson           code << "      }\n";
9334b3e95d5SJeremy L Thompson           break;
9344b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
935f815fac9SJeremy L Thompson           code << "      GradColloSliceTranspose3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_s" << var_suffix << ", s_G"
936f815fac9SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
9374b3e95d5SJeremy L Thompson           break;
9384b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
9394b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
9404b3e95d5SJeremy L Thompson           break;  // Should not occur
9414b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
9424b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
9434b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
9444b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
9454b3e95d5SJeremy L Thompson       }
9464b3e95d5SJeremy L Thompson     }
9474b3e95d5SJeremy L Thompson   }
9484b3e95d5SJeremy L Thompson   code << "    }\n";
9494b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9504b3e95d5SJeremy L Thompson }
9514b3e95d5SJeremy L Thompson 
9524b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
953b2165e7aSSebastian Grimberg // Build single operator kernel
954ab213215SJeremy L Thompson //------------------------------------------------------------------------------
955ddae5012SJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Cuda_gen(CeedOperator op, bool *is_good_build) {
9568b97b69aSJeremy L Thompson   bool                    is_tensor = true, is_at_points = false, use_3d_slices = false;
957241a4b83SYohann   Ceed                    ceed;
9588b97b69aSJeremy L Thompson   CeedInt                 Q_1d, num_input_fields, num_output_fields, dim = 1, max_num_points = 0, coords_comp_stride = 0;
959ca735530SJeremy L Thompson   CeedQFunctionField     *qf_input_fields, *qf_output_fields;
960ca735530SJeremy L Thompson   CeedQFunction_Cuda_gen *qf_data;
961ca735530SJeremy L Thompson   CeedQFunction           qf;
962ca735530SJeremy L Thompson   CeedOperatorField      *op_input_fields, *op_output_fields;
963ca735530SJeremy L Thompson   CeedOperator_Cuda_gen  *data;
9644b3e95d5SJeremy L Thompson   std::ostringstream      code;
9654b3e95d5SJeremy L Thompson 
966ddae5012SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
9674b3e95d5SJeremy L Thompson   {
9684b3e95d5SJeremy L Thompson     bool is_setup_done;
969ca735530SJeremy L Thompson 
970ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
971ddae5012SJeremy L Thompson     if (is_setup_done) {
972ddae5012SJeremy L Thompson       *is_good_build = !data->use_fallback;
973ddae5012SJeremy L Thompson       return CEED_ERROR_SUCCESS;
974ddae5012SJeremy L Thompson     }
975ddae5012SJeremy L Thompson   }
976ddae5012SJeremy L Thompson 
977ddae5012SJeremy L Thompson   // Check field compatibility
9788d12f40eSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
979ddae5012SJeremy L Thompson   {
980ddae5012SJeremy L Thompson     bool has_shared_bases = true, is_all_tensor = true, is_all_nontensor = true;
981ddae5012SJeremy L Thompson 
982ddae5012SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
983ddae5012SJeremy L Thompson       CeedBasis basis;
984ddae5012SJeremy L Thompson 
985ddae5012SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
986ddae5012SJeremy L Thompson       if (basis != CEED_BASIS_NONE) {
987ddae5012SJeremy L Thompson         bool        is_tensor = true;
988ddae5012SJeremy L Thompson         const char *resource;
989ddae5012SJeremy L Thompson         char       *resource_root;
990ddae5012SJeremy L Thompson         Ceed        basis_ceed;
991ddae5012SJeremy L Thompson 
992ddae5012SJeremy L Thompson         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
993c9192acaSJeremy L Thompson         is_all_tensor    = is_all_tensor && is_tensor;
994c9192acaSJeremy L Thompson         is_all_nontensor = is_all_not_tensor && !is_tensor;
995ddae5012SJeremy L Thompson         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
996ddae5012SJeremy L Thompson         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
997ddae5012SJeremy L Thompson         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
998c9192acaSJeremy L Thompson         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/cuda/shared");
999ddae5012SJeremy L Thompson         CeedCallBackend(CeedFree(&resource_root));
1000ddae5012SJeremy L Thompson         CeedCallBackend(CeedDestroy(&basis_ceed));
1001ddae5012SJeremy L Thompson       }
1002ddae5012SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis));
10034b3e95d5SJeremy L Thompson     }
1004ca735530SJeremy L Thompson 
1005ddae5012SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
1006ddae5012SJeremy L Thompson       CeedBasis basis;
1007ddae5012SJeremy L Thompson 
1008ddae5012SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1009ddae5012SJeremy L Thompson       if (basis != CEED_BASIS_NONE) {
1010ddae5012SJeremy L Thompson         bool        is_tensor = true;
1011ddae5012SJeremy L Thompson         const char *resource;
1012ddae5012SJeremy L Thompson         char       *resource_root;
1013ddae5012SJeremy L Thompson         Ceed        basis_ceed;
1014ddae5012SJeremy L Thompson 
1015ddae5012SJeremy L Thompson         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
1016c9192acaSJeremy L Thompson         is_all_tensor    = is_all_tensor && is_tensor;
1017c9192acaSJeremy L Thompson         is_all_nontensor = is_all_nontensor && !is_tensor;
1018ddae5012SJeremy L Thompson 
1019ddae5012SJeremy L Thompson         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
1020ddae5012SJeremy L Thompson         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
1021ddae5012SJeremy L Thompson         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
1022c9192acaSJeremy L Thompson         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/cuda/shared");
1023ddae5012SJeremy L Thompson         CeedCallBackend(CeedFree(&resource_root));
1024ddae5012SJeremy L Thompson         CeedCallBackend(CeedDestroy(&basis_ceed));
1025ddae5012SJeremy L Thompson       }
1026ddae5012SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis));
1027ddae5012SJeremy L Thompson     }
1028ddae5012SJeremy L Thompson     // -- Fallback to ref if not all bases are shared
1029ddae5012SJeremy L Thompson     if (!has_shared_bases || (!is_all_tensor && !is_all_nontensor)) {
1030ddae5012SJeremy L Thompson       *is_good_build = false;
1031ddae5012SJeremy L Thompson       return CEED_ERROR_SUCCESS;
1032ddae5012SJeremy L Thompson     }
1033ddae5012SJeremy L Thompson   }
1034ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1035ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1036ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
1037ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
1038241a4b83SYohann 
10394b3e95d5SJeremy L Thompson   // Get operator data
10408b97b69aSJeremy L Thompson   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
10414b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelData_Cuda_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, op_output_fields,
10424b3e95d5SJeremy L Thompson                                                        qf_output_fields, &data->max_P_1d, &Q_1d, &dim, &is_tensor, &use_3d_slices));
10434b3e95d5SJeremy L Thompson   if (dim == 0) dim = 1;
10444b3e95d5SJeremy L Thompson   data->dim = dim;
10458b97b69aSJeremy L Thompson   if (is_at_points) {
10468b97b69aSJeremy L Thompson     CeedElemRestriction_Cuda *rstr_data;
10478b97b69aSJeremy L Thompson     CeedElemRestriction       rstr_points = NULL;
10484b3e95d5SJeremy L Thompson 
10498b97b69aSJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
10508b97b69aSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points));
10518b97b69aSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetCompStride(rstr_points, &coords_comp_stride));
10528b97b69aSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(rstr_points, &rstr_data));
10538b97b69aSJeremy L Thompson     data->points.indices = (CeedInt *)rstr_data->d_offsets;
10548b97b69aSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
10558b97b69aSJeremy L Thompson   }
10568b97b69aSJeremy L Thompson   if (is_at_points) use_3d_slices = false;
10578b97b69aSJeremy L Thompson   if (Q_1d == 0) {
10588b97b69aSJeremy L Thompson     if (is_at_points) Q_1d = max_num_points;
10598b97b69aSJeremy L Thompson     else CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q_1d));
10604b3e95d5SJeremy L Thompson   }
10614b3e95d5SJeremy L Thompson   data->Q_1d = Q_1d;
10624b3e95d5SJeremy L Thompson 
10630b454692Sjeremylt   // Check for restriction only identity operator
10644b3e95d5SJeremy L Thompson   {
10654b3e95d5SJeremy L Thompson     bool is_identity_qf;
10664b3e95d5SJeremy L Thompson 
10672b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
10680b454692Sjeremylt     if (is_identity_qf) {
10699e201c85SYohann       CeedEvalMode eval_mode_in, eval_mode_out;
1070ca735530SJeremy L Thompson 
10712b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
10722b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
10736574a04fSJeremy L Thompson       CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
10746574a04fSJeremy L Thompson                 "Backend does not implement restriction only identity operators");
10750b454692Sjeremylt     }
10764b3e95d5SJeremy L Thompson   }
10770b454692Sjeremylt 
1078f1a13f77SYohann Dudouit   // Add atomicAdd function for old NVidia architectures
10794b3e95d5SJeremy L Thompson   {
10804b3e95d5SJeremy L Thompson     Ceed_Cuda            *ceed_data;
10814b3e95d5SJeremy L Thompson     struct cudaDeviceProp prop;
10824b3e95d5SJeremy L Thompson 
10832b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetData(ceed, &ceed_data));
10842b730f8bSJeremy L Thompson     CeedCallBackend(cudaGetDeviceProperties(&prop, ceed_data->device_id));
108580a9ef05SNatalie Beams     if ((prop.major < 6) && (CEED_SCALAR_TYPE != CEED_SCALAR_FP32)) {
10869c25dd66SJeremy L Thompson       code << "// AtomicAdd fallback source\n";
10879c25dd66SJeremy L Thompson       code << "#include <ceed/jit-source/cuda/cuda-atomic-add-fallback.h>\n\n";
1088f1a13f77SYohann Dudouit     }
10894b3e95d5SJeremy L Thompson   }
1090f1a13f77SYohann Dudouit 
10919e201c85SYohann   // Load basis source files
1092dc007f05SJeremy L Thompson   if (is_tensor) {
10939c25dd66SJeremy L Thompson     code << "// Tensor basis source\n";
10949c25dd66SJeremy L Thompson     code << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-templates.h>\n\n";
1095dc007f05SJeremy L Thompson   } else {
1096dc007f05SJeremy L Thompson     code << "// Non-tensor basis source\n";
1097dc007f05SJeremy L Thompson     code << "#include <ceed/jit-source/cuda/cuda-shared-basis-nontensor-templates.h>\n\n";
1098dc007f05SJeremy L Thompson   }
1099dc007f05SJeremy L Thompson   if (is_at_points) {
11008b97b69aSJeremy L Thompson     code << "// AtPoints basis source\n";
11018b97b69aSJeremy L Thompson     code << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-at-points-templates.h>\n\n";
1102dc007f05SJeremy L Thompson   }
11039c25dd66SJeremy L Thompson   code << "// CodeGen operator source\n";
11049c25dd66SJeremy L Thompson   code << "#include <ceed/jit-source/cuda/cuda-gen-templates.h>\n\n";
1105241a4b83SYohann 
11064b3e95d5SJeremy L Thompson   // Get QFunction name
11074b3e95d5SJeremy L Thompson   std::string qfunction_name(qf_data->qfunction_name);
11084b3e95d5SJeremy L Thompson   std::string operator_name;
11094b3e95d5SJeremy L Thompson 
111009095acaSJeremy L Thompson   operator_name = "CeedKernelCudaGenOperator_" + qfunction_name;
11114d537eeaSYohann 
11129e201c85SYohann   // Define CEED_Q_VLA
11139e201c85SYohann   code << "\n#undef CEED_Q_VLA\n";
1114dc007f05SJeremy L Thompson   if (dim != 3 || is_at_points || use_3d_slices || !is_tensor) {
11159e201c85SYohann     code << "#define CEED_Q_VLA 1\n\n";
11169e201c85SYohann   } else {
11179e201c85SYohann     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
11189e201c85SYohann   }
11199e201c85SYohann 
11204b3e95d5SJeremy L Thompson   // Add user QFunction source
11214b3e95d5SJeremy L Thompson   {
11229c25dd66SJeremy L Thompson     const char *source_path;
11234b3e95d5SJeremy L Thompson 
11249c25dd66SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path));
11259c25dd66SJeremy L Thompson     CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/cuda/gen backend requires QFunction source code file");
11269c25dd66SJeremy L Thompson 
11279c25dd66SJeremy L Thompson     code << "// User QFunction source\n";
11289c25dd66SJeremy L Thompson     code << "#include \"" << source_path << "\"\n\n";
11294b3e95d5SJeremy L Thompson   }
1130241a4b83SYohann 
1131241a4b83SYohann   // Setup
1132818e0025SJeremy L Thompson   code << "\n// -----------------------------------------------------------------------------\n";
11334b3e95d5SJeremy L Thompson   code << "// Operator Kernel\n";
11344b3e95d5SJeremy L Thompson   code << "// \n";
11354b3e95d5SJeremy L Thompson   code << "// d_[in,out]_i:   CeedVector device array\n";
11364b3e95d5SJeremy L Thompson   code << "// r_[in,out]_e_i: Element vector register\n";
11374b3e95d5SJeremy L Thompson   code << "// r_[in,out]_q_i: Quadrature space vector register\n";
11389123fb08SJeremy L Thompson   code << "// r_[in,out]_c_i: AtPoints Chebyshev coefficients register\n";
11394b3e95d5SJeremy L Thompson   code << "// r_[in,out]_s_i: Quadrature space slice vector register\n";
11404b3e95d5SJeremy L Thompson   code << "// \n";
11414b3e95d5SJeremy L Thompson   code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n";
11424b3e95d5SJeremy L Thompson   code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n";
11434b3e95d5SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n";
11444b3e95d5SJeremy L Thompson   code << "extern \"C\" __global__ void " << operator_name
11458b97b69aSJeremy L Thompson        << "(CeedInt num_elem, void* ctx, FieldsInt_Cuda indices, Fields_Cuda fields, Fields_Cuda B, Fields_Cuda G, CeedScalar *W, Points_Cuda "
11468b97b69aSJeremy L Thompson           "points) {\n";
11474b3e95d5SJeremy L Thompson 
11484b3e95d5SJeremy L Thompson   // Scratch buffers
11499e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
11504b3e95d5SJeremy L Thompson     CeedEvalMode eval_mode;
11514b3e95d5SJeremy L Thompson 
11522b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
11539e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
11544b3e95d5SJeremy L Thompson       code << "  const CeedScalar *d_in_" << i << " = fields.inputs[" << i << "];\n";
1155241a4b83SYohann     }
1156241a4b83SYohann   }
11579e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
11584b3e95d5SJeremy L Thompson     code << "  CeedScalar *d_out_" << i << " = fields.outputs[" << i << "];\n";
1159241a4b83SYohann   }
11601da99368SJeremy L Thompson 
11619e201c85SYohann   code << "  const CeedInt dim = " << dim << ";\n";
1162dc007f05SJeremy L Thompson   code << "  const CeedInt " << (is_tensor ? "Q_1d" : "Q") << " = " << Q_1d << ";\n";
11638b97b69aSJeremy L Thompson   if (is_at_points) {
11648b97b69aSJeremy L Thompson     code << "  const CeedInt max_num_points = " << max_num_points << ";\n";
11658b97b69aSJeremy L Thompson     code << "  const CeedInt coords_comp_stride = " << coords_comp_stride << ";\n";
11668b97b69aSJeremy L Thompson   }
11671da99368SJeremy L Thompson 
11684b3e95d5SJeremy L Thompson   // Shared data
1169241a4b83SYohann   code << "  extern __shared__ CeedScalar slice[];\n";
11709e201c85SYohann   code << "  SharedData_Cuda data;\n";
11719e201c85SYohann   code << "  data.t_id_x = threadIdx.x;\n";
11729e201c85SYohann   code << "  data.t_id_y = threadIdx.y;\n";
11739e201c85SYohann   code << "  data.t_id_z = threadIdx.z;\n";
11749e201c85SYohann   code << "  data.t_id  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
1175dc007f05SJeremy L Thompson   code << "  data.slice = slice + data.t_id_z*T_1D" << ((!is_tensor || dim == 1) ? "" : "*T_1D") << ";\n";
1176920dcdc4Sjeremylt 
1177*0a2a6492SJeremy L Thompson   // -- Determine input mat reuse
1178*0a2a6492SJeremy L Thompson   CeedInt input_matrix_reuse[CEED_FIELD_MAX][3];  // field, is_input, eval_mode
1179*0a2a6492SJeremy L Thompson 
1180*0a2a6492SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1181*0a2a6492SJeremy L Thompson     input_matrix_reuse[i][0] = -1;
1182*0a2a6492SJeremy L Thompson   }
1183*0a2a6492SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1184*0a2a6492SJeremy L Thompson     CeedEvalMode eval_mode_i;
1185*0a2a6492SJeremy L Thompson     CeedBasis    basis_i;
1186*0a2a6492SJeremy L Thompson 
1187*0a2a6492SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode_i));
1188*0a2a6492SJeremy L Thompson     if (eval_mode_i == CEED_EVAL_WEIGHT) continue;
1189*0a2a6492SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis_i));
1190*0a2a6492SJeremy L Thompson     for (CeedInt j = 0; (input_matrix_reuse[i][0] == -1) && (j < i); j++) {
1191*0a2a6492SJeremy L Thompson       CeedEvalMode eval_mode_j;
1192*0a2a6492SJeremy L Thompson       CeedBasis    basis_j;
1193*0a2a6492SJeremy L Thompson 
1194*0a2a6492SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
1195*0a2a6492SJeremy L Thompson       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
1196*0a2a6492SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
1197*0a2a6492SJeremy L Thompson       if (basis_i == basis_j) {
1198*0a2a6492SJeremy L Thompson         if (is_tensor) {
1199*0a2a6492SJeremy L Thompson           input_matrix_reuse[i][0] = j;
1200*0a2a6492SJeremy L Thompson           input_matrix_reuse[i][1] = true;
1201*0a2a6492SJeremy L Thompson           input_matrix_reuse[i][2] = eval_mode_j;
1202*0a2a6492SJeremy L Thompson         } else {
1203*0a2a6492SJeremy L Thompson           // For non-tensor can only re-use with the same eval mode
1204*0a2a6492SJeremy L Thompson           if (eval_mode_i == eval_mode_j) {
1205*0a2a6492SJeremy L Thompson             input_matrix_reuse[i][0] = j;
1206*0a2a6492SJeremy L Thompson             input_matrix_reuse[i][1] = true;
1207*0a2a6492SJeremy L Thompson             input_matrix_reuse[i][2] = eval_mode_j;
1208*0a2a6492SJeremy L Thompson           }
1209*0a2a6492SJeremy L Thompson         }
1210*0a2a6492SJeremy L Thompson       }
1211*0a2a6492SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis_j));
1212*0a2a6492SJeremy L Thompson     }
1213*0a2a6492SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis_i));
1214*0a2a6492SJeremy L Thompson   }
1215*0a2a6492SJeremy L Thompson 
1216*0a2a6492SJeremy L Thompson   // -- Determine output mat reuse
1217*0a2a6492SJeremy L Thompson   CeedInt output_matrix_reuse[CEED_FIELD_MAX][3];  // field, is_input, eval_mode
1218*0a2a6492SJeremy L Thompson 
1219*0a2a6492SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1220*0a2a6492SJeremy L Thompson     output_matrix_reuse[i][0] = -1;
1221*0a2a6492SJeremy L Thompson   }
1222*0a2a6492SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1223*0a2a6492SJeremy L Thompson     CeedEvalMode eval_mode_i;
1224*0a2a6492SJeremy L Thompson     CeedBasis    basis_i;
1225*0a2a6492SJeremy L Thompson 
1226*0a2a6492SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode_i));
1227*0a2a6492SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis_i));
1228*0a2a6492SJeremy L Thompson     for (CeedInt j = 0; (output_matrix_reuse[i][0] == -1) && (j < num_input_fields); j++) {
1229*0a2a6492SJeremy L Thompson       CeedEvalMode eval_mode_j;
1230*0a2a6492SJeremy L Thompson       CeedBasis    basis_j;
1231*0a2a6492SJeremy L Thompson 
1232*0a2a6492SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
1233*0a2a6492SJeremy L Thompson       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
1234*0a2a6492SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
1235*0a2a6492SJeremy L Thompson       if (basis_i == basis_j) {
1236*0a2a6492SJeremy L Thompson         if (is_tensor) {
1237*0a2a6492SJeremy L Thompson           output_matrix_reuse[i][0] = j;
1238*0a2a6492SJeremy L Thompson           output_matrix_reuse[i][1] = true;
1239*0a2a6492SJeremy L Thompson           output_matrix_reuse[i][2] = eval_mode_j;
1240*0a2a6492SJeremy L Thompson         } else {
1241*0a2a6492SJeremy L Thompson           // For non-tensor can only re-use with the same eval mode
1242*0a2a6492SJeremy L Thompson           if (eval_mode_i == eval_mode_j) {
1243*0a2a6492SJeremy L Thompson             output_matrix_reuse[i][0] = j;
1244*0a2a6492SJeremy L Thompson             output_matrix_reuse[i][1] = true;
1245*0a2a6492SJeremy L Thompson             output_matrix_reuse[i][2] = eval_mode_j;
1246*0a2a6492SJeremy L Thompson           }
1247*0a2a6492SJeremy L Thompson         }
1248*0a2a6492SJeremy L Thompson       }
1249*0a2a6492SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis_j));
1250*0a2a6492SJeremy L Thompson     }
1251*0a2a6492SJeremy L Thompson     for (CeedInt j = 0; (output_matrix_reuse[i][0] == -1) && (j < i); j++) {
1252*0a2a6492SJeremy L Thompson       CeedEvalMode eval_mode_j;
1253*0a2a6492SJeremy L Thompson       CeedBasis    basis_j;
1254*0a2a6492SJeremy L Thompson 
1255*0a2a6492SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[j], &eval_mode_j));
1256*0a2a6492SJeremy L Thompson       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
1257*0a2a6492SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[j], &basis_j));
1258*0a2a6492SJeremy L Thompson       if (basis_i == basis_j) {
1259*0a2a6492SJeremy L Thompson         if (is_tensor) {
1260*0a2a6492SJeremy L Thompson           output_matrix_reuse[i][0] = j;
1261*0a2a6492SJeremy L Thompson           output_matrix_reuse[i][1] = false;
1262*0a2a6492SJeremy L Thompson           output_matrix_reuse[i][2] = eval_mode_j;
1263*0a2a6492SJeremy L Thompson         } else {
1264*0a2a6492SJeremy L Thompson           // For non-tensor can only re-use with the same eval mode
1265*0a2a6492SJeremy L Thompson           if (eval_mode_i == eval_mode_j) {
1266*0a2a6492SJeremy L Thompson             output_matrix_reuse[i][0] = j;
1267*0a2a6492SJeremy L Thompson             output_matrix_reuse[i][1] = false;
1268*0a2a6492SJeremy L Thompson             output_matrix_reuse[i][2] = eval_mode_j;
1269*0a2a6492SJeremy L Thompson           }
1270*0a2a6492SJeremy L Thompson         }
1271*0a2a6492SJeremy L Thompson       }
1272*0a2a6492SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis_j));
1273*0a2a6492SJeremy L Thompson     }
1274*0a2a6492SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis_i));
1275*0a2a6492SJeremy L Thompson   }
1276*0a2a6492SJeremy L Thompson 
1277ac421f39SYohann   // Initialize constants, and matrices B and G
12784b3e95d5SJeremy L Thompson   code << "\n  // Input field constants and basis data\n";
12799e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1280*0a2a6492SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, i, op_input_fields[i], qf_input_fields[i], input_matrix_reuse[i], Q_1d,
1281*0a2a6492SJeremy L Thompson                                                               true, is_tensor, is_at_points, use_3d_slices));
1282920dcdc4Sjeremylt   }
12834b3e95d5SJeremy L Thompson   code << "\n  // Output field constants and basis data\n";
12849e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1285*0a2a6492SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, i, op_output_fields[i], qf_output_fields[i], output_matrix_reuse[i], Q_1d,
1286*0a2a6492SJeremy L Thompson                                                               false, is_tensor, is_at_points, use_3d_slices));
12874b3e95d5SJeremy L Thompson   }
1288920dcdc4Sjeremylt 
12894b3e95d5SJeremy L Thompson   // Loop over all elements
12904b3e95d5SJeremy L Thompson   code << "\n  // Element loop\n";
1291ac421f39SYohann   code << "  __syncthreads();\n";
12929e201c85SYohann   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
12934b3e95d5SJeremy L Thompson 
1294e93651e5SJeremy L Thompson   // -- Compute minimum buffer space needed
12958b97b69aSJeremy L Thompson   CeedInt max_rstr_buffer_size = 1;
1296e93651e5SJeremy L Thompson 
12979e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1298e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
1299e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
1300e93651e5SJeremy L Thompson 
1301e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
1302e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1303e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1304dc007f05SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_tensor && (dim >= 3) ? elem_size : 1));
1305681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1306e93651e5SJeremy L Thompson   }
1307e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1308e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
1309e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
1310e93651e5SJeremy L Thompson 
1311e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
1312e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1313e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1314dc007f05SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_tensor && (dim >= 3) ? elem_size : 1));
1315681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1316e93651e5SJeremy L Thompson   }
1317e93651e5SJeremy L Thompson   code << "    // Scratch restriction buffer space\n";
1318e93651e5SJeremy L Thompson   code << "    CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n";
1319e93651e5SJeremy L Thompson 
1320e93651e5SJeremy L Thompson   // -- Determine best input field processing order
1321e93651e5SJeremy L Thompson   CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX];
1322e93651e5SJeremy L Thompson 
1323e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1324e93651e5SJeremy L Thompson     field_rstr_in_buffer[i] = -1;
1325e93651e5SJeremy L Thompson     input_field_order[i]    = -1;
1326e93651e5SJeremy L Thompson   }
1327e93651e5SJeremy L Thompson   {
1328e93651e5SJeremy L Thompson     bool    is_ordered[CEED_FIELD_MAX];
1329e93651e5SJeremy L Thompson     CeedInt curr_index = 0;
1330e93651e5SJeremy L Thompson 
1331e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
1332e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
1333e93651e5SJeremy L Thompson       CeedVector          vec_i;
1334e93651e5SJeremy L Thompson       CeedElemRestriction rstr_i;
1335e93651e5SJeremy L Thompson 
1336e93651e5SJeremy L Thompson       if (is_ordered[i]) continue;
1337e93651e5SJeremy L Thompson       field_rstr_in_buffer[i]       = i;
1338e93651e5SJeremy L Thompson       is_ordered[i]                 = true;
1339e93651e5SJeremy L Thompson       input_field_order[curr_index] = i;
1340e93651e5SJeremy L Thompson       curr_index++;
1341034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
1342e93651e5SJeremy L Thompson       if (vec_i == CEED_VECTOR_NONE) continue;  // CEED_EVAL_WEIGHT
1343e93651e5SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
1344e93651e5SJeremy L Thompson       for (CeedInt j = i + 1; j < num_input_fields; j++) {
1345e93651e5SJeremy L Thompson         CeedVector          vec_j;
1346e93651e5SJeremy L Thompson         CeedElemRestriction rstr_j;
1347e93651e5SJeremy L Thompson 
1348e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
1349e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
1350e93651e5SJeremy L Thompson         if (rstr_i == rstr_j && vec_i == vec_j) {
1351e93651e5SJeremy L Thompson           field_rstr_in_buffer[j]       = i;
1352e93651e5SJeremy L Thompson           is_ordered[j]                 = true;
1353e93651e5SJeremy L Thompson           input_field_order[curr_index] = j;
1354e93651e5SJeremy L Thompson           curr_index++;
1355e93651e5SJeremy L Thompson         }
1356681d0ea7SJeremy L Thompson         CeedCallBackend(CeedVectorDestroy(&vec_j));
1357681d0ea7SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
1358e93651e5SJeremy L Thompson       }
1359681d0ea7SJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec_i));
1360681d0ea7SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
1361e93651e5SJeremy L Thompson     }
1362e93651e5SJeremy L Thompson   }
1363e93651e5SJeremy L Thompson 
1364e93651e5SJeremy L Thompson   // -- Input restriction and basis
1365e93651e5SJeremy L Thompson   code << "\n    // -- Input field restrictions and basis actions\n";
1366e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1367e93651e5SJeremy L Thompson     CeedInt f = input_field_order[i];
1368e93651e5SJeremy L Thompson 
1369e93651e5SJeremy L Thompson     code << "    // ---- Input field " << f << "\n";
1370920dcdc4Sjeremylt 
13714b3e95d5SJeremy L Thompson     // ---- Restriction
1372e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, f, dim, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f],
1373dc007f05SJeremy L Thompson                                                                 Q_1d, true, is_tensor, is_at_points, use_3d_slices));
13749a2291e3SJeremy L Thompson 
13754b3e95d5SJeremy L Thompson     // ---- Basis action
1376dc007f05SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, f, dim, op_input_fields[f], qf_input_fields[f], Q_1d, true, is_tensor,
1377dc007f05SJeremy L Thompson                                                           is_at_points, use_3d_slices));
1378920dcdc4Sjeremylt   }
1379920dcdc4Sjeremylt 
13804b3e95d5SJeremy L Thompson   // -- Q function
13818b97b69aSJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelQFunction_Cuda_gen(code, data, dim, max_num_points, num_input_fields, op_input_fields, qf_input_fields,
1382dc007f05SJeremy L Thompson                                                             num_output_fields, op_output_fields, qf_output_fields, qfunction_name, Q_1d, is_tensor,
1383dc007f05SJeremy L Thompson                                                             is_at_points, use_3d_slices));
1384ca735530SJeremy L Thompson 
13854b3e95d5SJeremy L Thompson   // -- Output basis and restriction
13864b3e95d5SJeremy L Thompson   code << "\n    // -- Output field basis action and restrictions\n";
13879e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
13884b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
1389ca735530SJeremy L Thompson 
13904b3e95d5SJeremy L Thompson     // ---- Basis action
1391dc007f05SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, i, dim, op_output_fields[i], qf_output_fields[i], Q_1d, false, is_tensor,
1392dc007f05SJeremy L Thompson                                                           is_at_points, use_3d_slices));
13939a2291e3SJeremy L Thompson 
13944b3e95d5SJeremy L Thompson     // ---- Restriction
13958b97b69aSJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, i, dim, NULL, op_output_fields[i], qf_output_fields[i], Q_1d, false,
1396dc007f05SJeremy L Thompson                                                                 is_tensor, is_at_points, use_3d_slices));
1397ac421f39SYohann   }
1398241a4b83SYohann 
13994b3e95d5SJeremy L Thompson   // Close loop and function
1400241a4b83SYohann   code << "  }\n";
1401818e0025SJeremy L Thompson   code << "}\n";
1402818e0025SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n\n";
1403241a4b83SYohann 
14043a2968d6SJeremy L Thompson   // Compile
1405ddae5012SJeremy L Thompson   {
1406ddae5012SJeremy L Thompson     bool is_compile_good = false;
1407ddae5012SJeremy L Thompson 
1408ddae5012SJeremy L Thompson     CeedCallBackend(CeedTryCompile_Cuda(ceed, code.str().c_str(), &is_compile_good, &data->module, 1, "T_1D", CeedIntMax(Q_1d, data->max_P_1d)));
1409ddae5012SJeremy L Thompson     if (is_compile_good) {
1410ddae5012SJeremy L Thompson       *is_good_build = true;
1411eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, operator_name.c_str(), &data->op));
1412ddae5012SJeremy L Thompson     } else {
1413ddae5012SJeremy L Thompson       *is_good_build     = false;
1414ddae5012SJeremy L Thompson       data->use_fallback = true;
1415ddae5012SJeremy L Thompson     }
1416ddae5012SJeremy L Thompson   }
14172b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
14189bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
1419c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionDestroy(&qf));
1420e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1421241a4b83SYohann }
14222a86cc9dSSebastian Grimberg 
1423ab213215SJeremy L Thompson //------------------------------------------------------------------------------
1424