xref: /libCEED/backends/cuda-gen/ceed-cuda-gen-operator-build.cpp (revision 45a787f75d07c3f171308ee0d0be6daefd4754b0)
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 
25*45a787f7SJeremy L Thompson struct FieldReuse_Cuda {
26*45a787f7SJeremy L Thompson   CeedInt      index;
27*45a787f7SJeremy L Thompson   bool         is_input;
28*45a787f7SJeremy L Thompson   CeedEvalMode eval_mode;
29*45a787f7SJeremy L Thompson };
30*45a787f7SJeremy L Thompson 
31ab213215SJeremy L Thompson //------------------------------------------------------------------------------
324b3e95d5SJeremy L Thompson // Determine type of operator
334b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
344b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelData_Cuda_gen(Ceed ceed, CeedInt num_input_fields, CeedOperatorField *op_input_fields,
354b3e95d5SJeremy L Thompson                                                 CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, CeedOperatorField *op_output_fields,
364b3e95d5SJeremy L Thompson                                                 CeedQFunctionField *qf_output_fields, CeedInt *max_P_1d, CeedInt *Q_1d, CeedInt *dim, bool *is_tensor,
374b3e95d5SJeremy L Thompson                                                 bool *use_3d_slices) {
384b3e95d5SJeremy L Thompson   // Find dim, P_1d, Q_1d
394b3e95d5SJeremy L Thompson   *max_P_1d  = 0;
404b3e95d5SJeremy L Thompson   *Q_1d      = 0;
414b3e95d5SJeremy L Thompson   *dim       = 0;
424b3e95d5SJeremy L Thompson   *is_tensor = true;
43dc007f05SJeremy L Thompson 
444b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
454b3e95d5SJeremy L Thompson     CeedBasis basis;
464b3e95d5SJeremy L Thompson 
474b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
484b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
494b3e95d5SJeremy L Thompson       bool    is_field_tensor;
504b3e95d5SJeremy L Thompson       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
514b3e95d5SJeremy L Thompson 
524b3e95d5SJeremy L Thompson       // Collect dim, P_1d, and Q_1d
534b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
544b3e95d5SJeremy L Thompson       *is_tensor = *is_tensor && is_field_tensor;
55dc007f05SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
56dc007f05SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P_1d));
574b3e95d5SJeremy L Thompson       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
584b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
594b3e95d5SJeremy L Thompson       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
604b3e95d5SJeremy L Thompson       *dim = field_dim;
61dc007f05SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
62dc007f05SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q_1d));
634b3e95d5SJeremy L Thompson       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
644b3e95d5SJeremy L Thompson       *Q_1d = field_Q_1d;
654b3e95d5SJeremy L Thompson     }
66681d0ea7SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis));
674b3e95d5SJeremy L Thompson   }
684b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
694b3e95d5SJeremy L Thompson     CeedBasis basis;
704b3e95d5SJeremy L Thompson 
714b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
724b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
734b3e95d5SJeremy L Thompson       bool    is_field_tensor;
744b3e95d5SJeremy L Thompson       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
754b3e95d5SJeremy L Thompson 
764b3e95d5SJeremy L Thompson       // Collect dim, P_1d, and Q_1d
774b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
784b3e95d5SJeremy L Thompson       *is_tensor = *is_tensor && is_field_tensor;
79dc007f05SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
80dc007f05SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P_1d));
814b3e95d5SJeremy L Thompson       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
824b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
834b3e95d5SJeremy L Thompson       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
844b3e95d5SJeremy L Thompson       *dim = field_dim;
85dc007f05SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
86dc007f05SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q_1d));
874b3e95d5SJeremy L Thompson       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
884b3e95d5SJeremy L Thompson       *Q_1d = field_Q_1d;
894b3e95d5SJeremy L Thompson     }
90681d0ea7SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis));
914b3e95d5SJeremy L Thompson   }
924b3e95d5SJeremy L Thompson 
934b3e95d5SJeremy L Thompson   // Only use 3D collocated gradient parallelization strategy when gradient is computed
944b3e95d5SJeremy L Thompson   *use_3d_slices = false;
954b3e95d5SJeremy L Thompson   if (*dim == 3) {
964b3e95d5SJeremy L Thompson     bool was_grad_found = false;
974b3e95d5SJeremy L Thompson 
984b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
994b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
1004b3e95d5SJeremy L Thompson 
1014b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1024b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
1034b3e95d5SJeremy L Thompson         CeedBasis_Cuda_shared *basis_data;
1044b3e95d5SJeremy L Thompson         CeedBasis              basis;
1054b3e95d5SJeremy L Thompson 
1064b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1074b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1084b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
1094b3e95d5SJeremy L Thompson         was_grad_found = true;
110681d0ea7SJeremy L Thompson         CeedCallBackend(CeedBasisDestroy(&basis));
1114b3e95d5SJeremy L Thompson       }
1124b3e95d5SJeremy L Thompson     }
1134b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
1144b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
1154b3e95d5SJeremy L Thompson 
1164b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1174b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
1184b3e95d5SJeremy L Thompson         CeedBasis_Cuda_shared *basis_data;
1194b3e95d5SJeremy L Thompson         CeedBasis              basis;
1204b3e95d5SJeremy L Thompson 
1214b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1224b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1234b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
1244b3e95d5SJeremy L Thompson         was_grad_found = true;
125681d0ea7SJeremy L Thompson         CeedCallBackend(CeedBasisDestroy(&basis));
1264b3e95d5SJeremy L Thompson       }
1274b3e95d5SJeremy L Thompson     }
1284b3e95d5SJeremy L Thompson   }
1294b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1304b3e95d5SJeremy L Thompson }
1314b3e95d5SJeremy L Thompson 
1324b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1334b3e95d5SJeremy L Thompson // Setup fields
1344b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1354b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelFieldData_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedOperatorField op_field,
136*45a787f7SJeremy L Thompson                                                      CeedQFunctionField qf_field, FieldReuse_Cuda field_reuse, CeedInt Q_1d, bool is_input,
137*45a787f7SJeremy L Thompson                                                      bool is_tensor, bool is_at_points, bool use_3d_slices) {
1384b3e95d5SJeremy L Thompson   std::string            var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
139dc007f05SJeremy L Thompson   std::string            P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q";
1404b3e95d5SJeremy L Thompson   std::string            option_name = (is_input ? "inputs" : "outputs");
1414b3e95d5SJeremy L Thompson   CeedEvalMode           eval_mode   = CEED_EVAL_NONE;
1424b3e95d5SJeremy L Thompson   CeedInt                elem_size = 0, num_comp = 0, P_1d = 0;
1434b3e95d5SJeremy L Thompson   CeedElemRestriction    elem_rstr;
1444b3e95d5SJeremy L Thompson   CeedBasis_Cuda_shared *basis_data;
1454b3e95d5SJeremy L Thompson   CeedBasis              basis;
1464b3e95d5SJeremy L Thompson 
1470a2a6492SJeremy L Thompson   // Field reuse info
148*45a787f7SJeremy L Thompson   bool use_previous_field = field_reuse.index != -1;
1490a2a6492SJeremy L Thompson 
1504b3e95d5SJeremy L Thompson   code << "  // -- " << (is_input ? "Input" : "Output") << " field " << i << "\n";
1514b3e95d5SJeremy L Thompson 
1524b3e95d5SJeremy L Thompson   // Get field data
1534b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
1544b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
1554b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1564b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1574b3e95d5SJeremy L Thompson   }
158681d0ea7SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1594b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
1604b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
1614b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetData(basis, &basis_data));
162dc007f05SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
163dc007f05SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
1644b3e95d5SJeremy L Thompson   }
1654b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
1664b3e95d5SJeremy L Thompson 
1674b3e95d5SJeremy L Thompson   // Set field constants
1684b3e95d5SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
1694b3e95d5SJeremy L Thompson     code << "  const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n";
1704b3e95d5SJeremy L Thompson     code << "  const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n";
1714b3e95d5SJeremy L Thompson   }
1724b3e95d5SJeremy L Thompson 
1734b3e95d5SJeremy L Thompson   // Load basis data
1744b3e95d5SJeremy L Thompson   code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
1754b3e95d5SJeremy L Thompson   switch (eval_mode) {
1764b3e95d5SJeremy L Thompson     case CEED_EVAL_NONE:
1774b3e95d5SJeremy L Thompson       break;
1784b3e95d5SJeremy L Thompson     case CEED_EVAL_INTERP:
1798b97b69aSJeremy L Thompson       if (is_at_points) {
1808b97b69aSJeremy L Thompson         // AtPoints
1818b97b69aSJeremy L Thompson         if (!basis_data->d_chebyshev_interp_1d) {
1828b97b69aSJeremy L Thompson           CeedSize    interp_bytes;
1838b97b69aSJeremy L Thompson           CeedScalar *chebyshev_interp_1d;
1848b97b69aSJeremy L Thompson 
1858b97b69aSJeremy L Thompson           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
1868b97b69aSJeremy L Thompson           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
1878b97b69aSJeremy L Thompson           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
1888b97b69aSJeremy L Thompson           CeedCallCuda(CeedBasisReturnCeed(basis), cudaMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
1898b97b69aSJeremy L Thompson           CeedCallCuda(CeedBasisReturnCeed(basis),
1908b97b69aSJeremy L Thompson                        cudaMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
1918b97b69aSJeremy L Thompson           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
1928b97b69aSJeremy L Thompson         }
1938b97b69aSJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
1948b97b69aSJeremy L Thompson         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
1958b97b69aSJeremy L Thompson       } else {
1968b97b69aSJeremy L Thompson         // Standard quadrature
1974b3e95d5SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
1984b3e95d5SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_interp_1d;
1998b97b69aSJeremy L Thompson       }
2000a2a6492SJeremy L Thompson       if (use_previous_field) {
201*45a787f7SJeremy L Thompson         std::string reuse_var = "s_B" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
2020a2a6492SJeremy L Thompson 
2030a2a6492SJeremy L Thompson         code << "  CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n";
2040a2a6492SJeremy L Thompson       } else {
205dc007f05SJeremy L Thompson         code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
206f815fac9SJeremy L Thompson         code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
2070a2a6492SJeremy L Thompson       }
2084b3e95d5SJeremy L Thompson       break;
2094b3e95d5SJeremy L Thompson     case CEED_EVAL_GRAD:
2108b97b69aSJeremy L Thompson       if (is_at_points) {
2118b97b69aSJeremy L Thompson         // AtPoints
2128b97b69aSJeremy L Thompson         if (!basis_data->d_chebyshev_interp_1d) {
2138b97b69aSJeremy L Thompson           CeedSize    interp_bytes;
2148b97b69aSJeremy L Thompson           CeedScalar *chebyshev_interp_1d;
2158b97b69aSJeremy L Thompson 
2168b97b69aSJeremy L Thompson           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
2178b97b69aSJeremy L Thompson           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
2188b97b69aSJeremy L Thompson           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
2198b97b69aSJeremy L Thompson           CeedCallCuda(CeedBasisReturnCeed(basis), cudaMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
2208b97b69aSJeremy L Thompson           CeedCallCuda(CeedBasisReturnCeed(basis),
2218b97b69aSJeremy L Thompson                        cudaMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
2228b97b69aSJeremy L Thompson           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
2238b97b69aSJeremy L Thompson         }
2248b97b69aSJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
2258b97b69aSJeremy L Thompson         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
2268b97b69aSJeremy L Thompson       } else {
2278b97b69aSJeremy L Thompson         // Standard quadrature
2284b3e95d5SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
2294b3e95d5SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_interp_1d;
2308b97b69aSJeremy L Thompson       }
231dc007f05SJeremy L Thompson       if (is_tensor) {
2320a2a6492SJeremy L Thompson         if (use_previous_field) {
233*45a787f7SJeremy L Thompson           std::string reuse_var = "s_B" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
2340a2a6492SJeremy L Thompson 
2350a2a6492SJeremy L Thompson           code << "  CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n";
2360a2a6492SJeremy L Thompson         } else {
237dc007f05SJeremy L Thompson           code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
238f815fac9SJeremy L Thompson           code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
239dc007f05SJeremy L Thompson         }
2400a2a6492SJeremy L Thompson       }
2418b97b69aSJeremy L Thompson       if (is_at_points) break;  // No G mat for AtPoints
2424b3e95d5SJeremy L Thompson       if (use_3d_slices) {
2434b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d;
2444b3e95d5SJeremy L Thompson         else data->G.outputs[i] = basis_data->d_collo_grad_1d;
245*45a787f7SJeremy L Thompson         if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD) {
246*45a787f7SJeremy L Thompson           std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
2470a2a6492SJeremy L Thompson 
2480a2a6492SJeremy L Thompson           code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
2490a2a6492SJeremy L Thompson         } else {
250dc007f05SJeremy L Thompson           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
251f815fac9SJeremy L Thompson           code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
2520a2a6492SJeremy L Thompson         }
2534b3e95d5SJeremy L Thompson       } else {
2544b3e95d5SJeremy L Thompson         bool has_collo_grad = basis_data->d_collo_grad_1d;
2554b3e95d5SJeremy L Thompson 
2564b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2574b3e95d5SJeremy L Thompson         else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2584b3e95d5SJeremy L Thompson         if (has_collo_grad) {
259*45a787f7SJeremy L Thompson           if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD) {
260*45a787f7SJeremy L Thompson             std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
2610a2a6492SJeremy L Thompson 
2620a2a6492SJeremy L Thompson             code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
2630a2a6492SJeremy L Thompson           } else {
264dc007f05SJeremy L Thompson             code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
265f815fac9SJeremy L Thompson             code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
2660a2a6492SJeremy L Thompson           }
2670a2a6492SJeremy L Thompson         } else {
268*45a787f7SJeremy L Thompson           if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD) {
269*45a787f7SJeremy L Thompson             std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
2700a2a6492SJeremy L Thompson 
2710a2a6492SJeremy L Thompson             code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
2724b3e95d5SJeremy L Thompson           } else {
273dc007f05SJeremy L Thompson             code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << P_name << "*" << Q_name << (is_tensor ? "" : "*dim") << "];\n";
274dc007f05SJeremy L Thompson             code << "  LoadMatrix<" << P_name << ", " << Q_name << (is_tensor ? "" : "*dim") << ">(data, G." << option_name << "[" << i << "], s_G"
275dc007f05SJeremy L Thompson                  << var_suffix << ");\n";
2764b3e95d5SJeremy L Thompson           }
2774b3e95d5SJeremy L Thompson         }
2780a2a6492SJeremy L Thompson       }
2794b3e95d5SJeremy L Thompson       break;
2804b3e95d5SJeremy L Thompson     case CEED_EVAL_WEIGHT:
2814b3e95d5SJeremy L Thompson       break;  // No action
2824b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
2834b3e95d5SJeremy L Thompson     case CEED_EVAL_DIV:
2844b3e95d5SJeremy L Thompson     case CEED_EVAL_CURL:
2854b3e95d5SJeremy L Thompson       break;  // TODO: Not implemented
2864b3e95d5SJeremy L Thompson               // LCOV_EXCL_STOP
2874b3e95d5SJeremy L Thompson   }
2888b97b69aSJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
2894b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2904b3e95d5SJeremy L Thompson }
2914b3e95d5SJeremy L Thompson 
2924b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
2934b3e95d5SJeremy L Thompson // Restriction
2944b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
2954b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelRestriction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedInt dim,
296e93651e5SJeremy L Thompson                                                        CeedInt field_input_buffer[], CeedOperatorField op_field, CeedQFunctionField qf_field,
297dc007f05SJeremy L Thompson                                                        CeedInt Q_1d, bool is_input, bool is_tensor, bool is_at_points, bool use_3d_slices) {
2984b3e95d5SJeremy L Thompson   std::string               var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
299dc007f05SJeremy L Thompson   std::string               P_name     = (is_tensor ? "P_1d" : "P") + var_suffix;
3004b3e95d5SJeremy L Thompson   CeedEvalMode              eval_mode  = CEED_EVAL_NONE;
3014b3e95d5SJeremy L Thompson   CeedInt                   elem_size = 0, num_comp = 0, P_1d = 0;
3024b3e95d5SJeremy L Thompson   CeedSize                  l_size;
303f815fac9SJeremy L Thompson   CeedRestrictionType       rstr_type = CEED_RESTRICTION_STANDARD;
3044b3e95d5SJeremy L Thompson   CeedElemRestriction_Cuda *rstr_data;
3054b3e95d5SJeremy L Thompson   CeedElemRestriction       elem_rstr;
3064b3e95d5SJeremy L Thompson   CeedBasis                 basis;
3074b3e95d5SJeremy L Thompson 
3084b3e95d5SJeremy L Thompson   // Get field data
3094b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
3104b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
311f815fac9SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
3124b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
3134b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
3144b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
3154b3e95d5SJeremy L Thompson   }
3164b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
3174b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
318dc007f05SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
319dc007f05SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
3204b3e95d5SJeremy L Thompson   }
321681d0ea7SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
3224b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
3234b3e95d5SJeremy L Thompson 
3244b3e95d5SJeremy L Thompson   // Restriction
3254b3e95d5SJeremy L Thompson   if (is_input) {
3264b3e95d5SJeremy L Thompson     // Input
327e93651e5SJeremy L Thompson     if (field_input_buffer[i] != i) {
328e93651e5SJeremy L Thompson       std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]);
329e93651e5SJeremy L Thompson 
330e93651e5SJeremy L Thompson       // Restriction was already done for previous input
331e93651e5SJeremy L Thompson       code << "    CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n";
3328b97b69aSJeremy L Thompson     } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices && is_at_points)) {
3338b97b69aSJeremy L Thompson       if (eval_mode == CEED_EVAL_NONE && rstr_type != CEED_RESTRICTION_POINTS) {
334e93651e5SJeremy L Thompson         // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated
3354b3e95d5SJeremy L Thompson         code << "    CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n";
3368b97b69aSJeremy L Thompson       } else if (rstr_type != CEED_RESTRICTION_POINTS) {
337e93651e5SJeremy L Thompson         // Otherwise we're using the scratch space
338e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
339e93651e5SJeremy L Thompson       }
340f815fac9SJeremy L Thompson       switch (rstr_type) {
341f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STANDARD: {
3424b3e95d5SJeremy L Thompson           CeedInt comp_stride;
3434b3e95d5SJeremy L Thompson 
3444b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
3454b3e95d5SJeremy L Thompson           code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
3464b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3474b3e95d5SJeremy L Thompson           code << "    // CompStride: " << comp_stride << "\n";
3484b3e95d5SJeremy L Thompson           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
349dc007f05SJeremy L Thompson           code << "    ReadLVecStandard" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name
350dc007f05SJeremy L Thompson                << ">(data, l_size" << var_suffix << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n";
351f815fac9SJeremy L Thompson           break;
352f815fac9SJeremy L Thompson         }
353f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STRIDED: {
3544b3e95d5SJeremy L Thompson           bool    has_backend_strides;
3554b3e95d5SJeremy L Thompson           CeedInt num_elem;
3564b3e95d5SJeremy L Thompson 
3574b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
3584b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
3594b3e95d5SJeremy L Thompson           CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
3604b3e95d5SJeremy L Thompson 
3614b3e95d5SJeremy L Thompson           if (!has_backend_strides) {
3624b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
3634b3e95d5SJeremy L Thompson           }
3644b3e95d5SJeremy L Thompson           code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
365dc007f05SJeremy L Thompson           code << "    ReadLVecStrided" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", "
366dc007f05SJeremy L Thompson                << strides[1] << ", " << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n";
367f815fac9SJeremy L Thompson           break;
368f815fac9SJeremy L Thompson         }
3698b97b69aSJeremy L Thompson         case CEED_RESTRICTION_POINTS: {
3708b97b69aSJeremy L Thompson           CeedInt comp_stride;
3718b97b69aSJeremy L Thompson 
3728b97b69aSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3738b97b69aSJeremy L Thompson           code << "    const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
3748b97b69aSJeremy L Thompson           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
3758b97b69aSJeremy L Thompson           break;
3768b97b69aSJeremy L Thompson         }
377f815fac9SJeremy L Thompson         // LCOV_EXCL_START
378f815fac9SJeremy L Thompson         case CEED_RESTRICTION_ORIENTED:
379f815fac9SJeremy L Thompson         case CEED_RESTRICTION_CURL_ORIENTED:
380f815fac9SJeremy L Thompson           break;  // TODO: Not implemented
381f815fac9SJeremy L Thompson                   // LCOV_EXCL_STOP
3824b3e95d5SJeremy L Thompson       }
3834b3e95d5SJeremy L Thompson     }
3844b3e95d5SJeremy L Thompson   } else {
3854b3e95d5SJeremy L Thompson     // Output
386f815fac9SJeremy L Thompson     switch (rstr_type) {
387f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STANDARD: {
3884b3e95d5SJeremy L Thompson         CeedInt comp_stride;
3894b3e95d5SJeremy L Thompson 
3904b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
3914b3e95d5SJeremy L Thompson         code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
3924b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3934b3e95d5SJeremy L Thompson         code << "    // CompStride: " << comp_stride << "\n";
3944b3e95d5SJeremy L Thompson         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
395dc007f05SJeremy L Thompson         code << "    WriteLVecStandard" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name
396dc007f05SJeremy L Thompson              << ">(data, l_size" << var_suffix << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n";
397f815fac9SJeremy L Thompson         break;
398f815fac9SJeremy L Thompson       }
399f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STRIDED: {
4004b3e95d5SJeremy L Thompson         bool    has_backend_strides;
4014b3e95d5SJeremy L Thompson         CeedInt num_elem;
4024b3e95d5SJeremy L Thompson 
4034b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
4044b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
4054b3e95d5SJeremy L Thompson         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
4064b3e95d5SJeremy L Thompson 
4074b3e95d5SJeremy L Thompson         if (!has_backend_strides) {
4084b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
4094b3e95d5SJeremy L Thompson         }
4104b3e95d5SJeremy L Thompson         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
411dc007f05SJeremy L Thompson         code << "    WriteLVecStrided" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", "
412dc007f05SJeremy L Thompson              << strides[1] << ", " << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n";
413f815fac9SJeremy L Thompson         break;
414f815fac9SJeremy L Thompson       }
4158b97b69aSJeremy L Thompson       case CEED_RESTRICTION_POINTS:
4168b97b69aSJeremy L Thompson         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
4178b97b69aSJeremy L Thompson         break;
418f815fac9SJeremy L Thompson       // LCOV_EXCL_START
419f815fac9SJeremy L Thompson       case CEED_RESTRICTION_ORIENTED:
420f815fac9SJeremy L Thompson       case CEED_RESTRICTION_CURL_ORIENTED:
421f815fac9SJeremy L Thompson         break;  // TODO: Not implemented
422f815fac9SJeremy L Thompson                 // LCOV_EXCL_STOP
4234b3e95d5SJeremy L Thompson     }
4244b3e95d5SJeremy L Thompson   }
425681d0ea7SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
4264b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4274b3e95d5SJeremy L Thompson }
4284b3e95d5SJeremy L Thompson 
4294b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4304b3e95d5SJeremy L Thompson // Basis
4314b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4324b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelBasis_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedInt dim,
433dc007f05SJeremy L Thompson                                                  CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool is_tensor,
4348b97b69aSJeremy L Thompson                                                  bool is_at_points, bool use_3d_slices) {
4354b3e95d5SJeremy L Thompson   std::string         var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
436dc007f05SJeremy L Thompson   std::string         P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q";
4374b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
4384b3e95d5SJeremy L Thompson   CeedInt             elem_size = 0, num_comp = 0, P_1d = 0;
4394b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
4404b3e95d5SJeremy L Thompson   CeedBasis           basis;
4414b3e95d5SJeremy L Thompson 
4424b3e95d5SJeremy L Thompson   // Get field data
4434b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
4444b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
4454b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
4464b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
4474b3e95d5SJeremy L Thompson   }
448681d0ea7SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
4494b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
4504b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
451dc007f05SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
452dc007f05SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
4534b3e95d5SJeremy L Thompson   }
4544b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
4554b3e95d5SJeremy L Thompson 
4564b3e95d5SJeremy L Thompson   // Basis
4574b3e95d5SJeremy L Thompson   code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
4584b3e95d5SJeremy L Thompson   if (is_input) {
4594b3e95d5SJeremy L Thompson     switch (eval_mode) {
4604b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
4618b97b69aSJeremy L Thompson         if (!use_3d_slices && !is_at_points) {
4624b3e95d5SJeremy L Thompson           code << "    CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n";
4634b3e95d5SJeremy L Thompson         }
4644b3e95d5SJeremy L Thompson         break;
4654b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
4668b97b69aSJeremy L Thompson         if (is_at_points) {
467dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
468dc007f05SJeremy L Thompson 
469dc007f05SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
470dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
471dc007f05SJeremy L Thompson                << var_suffix << ", r_c" << var_suffix << ");\n";
4728b97b69aSJeremy L Thompson         } else {
473dc007f05SJeremy L Thompson           std::string function_name = is_tensor ? ((dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d") : "InterpNonTensor";
474dc007f05SJeremy L Thompson 
475dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
476dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
477dc007f05SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
4788b97b69aSJeremy L Thompson         }
4794b3e95d5SJeremy L Thompson         break;
4804b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
4818b97b69aSJeremy L Thompson         if (is_at_points) {
482dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
483dc007f05SJeremy L Thompson 
484dc007f05SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
485dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
486dc007f05SJeremy L Thompson                << var_suffix << ", r_c" << var_suffix << ");\n";
4878b97b69aSJeremy L Thompson         } else if (use_3d_slices) {
488dc007f05SJeremy L Thompson           std::string function_name = (dim > 1 ? "InterpTensor" : "Interp") + std::to_string(dim) + "d";
489dc007f05SJeremy L Thompson 
4904b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
491dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
492dc007f05SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
493dc007f05SJeremy L Thompson         } else if (is_tensor) {
494dc007f05SJeremy L Thompson           bool        is_collocated = dim == 3 && Q_1d >= P_1d;
495dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "Grad" : (is_collocated ? "GradTensorCollocated" : "GradTensor")) + std::to_string(dim) + "d";
496dc007f05SJeremy L Thompson 
497dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << (dim >= 3 ? Q_name : "1") << "];\n";
498dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
499dc007f05SJeremy L Thompson                << var_suffix << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
5004b3e95d5SJeremy L Thompson         } else {
501dc007f05SJeremy L Thompson           std::string function_name = "GradNonTensor";
502dc007f05SJeremy L Thompson 
503dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
504dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", dim, " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix
505dc007f05SJeremy L Thompson                << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
5064b3e95d5SJeremy L Thompson         }
5074b3e95d5SJeremy L Thompson         break;
5084b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
5098b97b69aSJeremy L Thompson         if (is_at_points) {
5108b97b69aSJeremy L Thompson           code << "    // Nothing to do AtPoints\n";
5118b97b69aSJeremy L Thompson         } else {
5124b3e95d5SJeremy L Thompson           CeedBasis_Cuda_shared *basis_data;
513dc007f05SJeremy L Thompson           std::string            function_name = is_tensor ? ((dim == 1 ? "Weight" : "WeightTensor") + std::to_string(dim) + "d") : "WeightNonTensor";
5144b3e95d5SJeremy L Thompson 
515dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
5164b3e95d5SJeremy L Thompson           CeedCallBackend(CeedBasisGetData(basis, &basis_data));
5174b3e95d5SJeremy L Thompson           data->W = basis_data->d_q_weight_1d;
518dc007f05SJeremy L Thompson           code << "    " << function_name << "<" << Q_name << ">(data, W, r_q" << var_suffix << ");\n";
5198b97b69aSJeremy L Thompson         }
5204b3e95d5SJeremy L Thompson         break;
5214b3e95d5SJeremy L Thompson       }
5224b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
5234b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
5244b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
5254b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
5264b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
5274b3e95d5SJeremy L Thompson     }
5284b3e95d5SJeremy L Thompson   } else {
5294b3e95d5SJeremy L Thompson     switch (eval_mode) {
5304b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
5314b3e95d5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n";
5324b3e95d5SJeremy L Thompson         break;  // No action
5334b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
534e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
5358b97b69aSJeremy L Thompson         if (is_at_points) {
536dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
537dc007f05SJeremy L Thompson 
538dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_c" << var_suffix << ", s_B"
539dc007f05SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5408b97b69aSJeremy L Thompson         } else {
541dc007f05SJeremy L Thompson           std::string function_name =
542dc007f05SJeremy L Thompson               is_tensor ? ((dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d") : "InterpTransposeNonTensor";
543dc007f05SJeremy L Thompson 
544dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B"
545dc007f05SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5468b97b69aSJeremy L Thompson         }
5474b3e95d5SJeremy L Thompson         break;
5484b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
549e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
5508b97b69aSJeremy L Thompson         if (is_at_points) {
551dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
552dc007f05SJeremy L Thompson 
553dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_c" << var_suffix << ", s_B"
554dc007f05SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5558b97b69aSJeremy L Thompson         } else if (use_3d_slices) {
556dc007f05SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
557dc007f05SJeremy L Thompson 
558dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B"
559dc007f05SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
560dc007f05SJeremy L Thompson         } else if (is_tensor) {
561dc007f05SJeremy L Thompson           bool        is_collocated = dim == 3 && Q_1d >= P_1d;
562dc007f05SJeremy L Thompson           std::string function_name =
563dc007f05SJeremy L Thompson               (dim == 1 ? "GradTranspose" : (is_collocated ? "GradTransposeTensorCollocated" : "GradTransposeTensor")) + std::to_string(dim) + "d";
564dc007f05SJeremy L Thompson 
565dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B"
566dc007f05SJeremy L Thompson                << var_suffix << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
5674b3e95d5SJeremy L Thompson         } else {
568dc007f05SJeremy L Thompson           std::string function_name = "GradTransposeNonTensor";
569dc007f05SJeremy L Thompson 
570dc007f05SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", dim, " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix
571dc007f05SJeremy L Thompson                << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
5724b3e95d5SJeremy L Thompson         }
5734b3e95d5SJeremy L Thompson         break;
5744b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
5754b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT:
5764b3e95d5SJeremy L Thompson         break;  // Should not occur
5774b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
5784b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
5794b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
5804b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
5814b3e95d5SJeremy L Thompson     }
5824b3e95d5SJeremy L Thompson   }
583681d0ea7SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
5844b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5854b3e95d5SJeremy L Thompson }
5864b3e95d5SJeremy L Thompson 
5874b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
5884b3e95d5SJeremy L Thompson // QFunction
5894b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
5908b97b69aSJeremy L Thompson static int CeedOperatorBuildKernelQFunction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt dim, CeedInt max_num_points,
5918b97b69aSJeremy L Thompson                                                      CeedInt num_input_fields, CeedOperatorField *op_input_fields,
5928b97b69aSJeremy L Thompson                                                      CeedQFunctionField *qf_input_fields, CeedInt num_output_fields,
5938b97b69aSJeremy L Thompson                                                      CeedOperatorField *op_output_fields, CeedQFunctionField *qf_output_fields,
594dc007f05SJeremy L Thompson                                                      std::string qfunction_name, CeedInt Q_1d, bool is_tensor, bool is_at_points,
595dc007f05SJeremy L Thompson                                                      bool use_3d_slices) {
596dc007f05SJeremy L Thompson   std::string         Q_name    = is_tensor ? "Q_1d" : "Q";
5974b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
5984b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
5994b3e95d5SJeremy L Thompson 
6008b97b69aSJeremy L Thompson   // Setup output arrays
6014b3e95d5SJeremy L Thompson   code << "\n    // -- Output field setup\n";
6024b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
6034b3e95d5SJeremy L Thompson     std::string var_suffix = "_out_" + std::to_string(i);
6044b3e95d5SJeremy L Thompson 
6054b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
6064b3e95d5SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
6078b97b69aSJeremy L Thompson     switch (eval_mode) {
6088b97b69aSJeremy L Thompson       case CEED_EVAL_NONE:
6098b97b69aSJeremy L Thompson         if (is_at_points) {
6108b97b69aSJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "];\n";
6118b97b69aSJeremy L Thompson         } else {
612dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
6134b3e95d5SJeremy L Thompson         }
6148b97b69aSJeremy L Thompson         break;
6158b97b69aSJeremy L Thompson       case CEED_EVAL_INTERP:
6168b97b69aSJeremy L Thompson         if (is_at_points) {
6178b97b69aSJeremy L Thompson           // Accumulator for point data
618dc007f05SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
619dc007f05SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "; i++) {\n";
6208b97b69aSJeremy L Thompson           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
6218b97b69aSJeremy L Thompson           code << "    }\n";
6228b97b69aSJeremy L Thompson         } else {
623dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
6248b97b69aSJeremy L Thompson         }
6258b97b69aSJeremy L Thompson         break;
6268b97b69aSJeremy L Thompson       case CEED_EVAL_GRAD:
6278b97b69aSJeremy L Thompson         if (is_at_points) {
6288b97b69aSJeremy L Thompson           // Accumulator for point data
629dc007f05SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "*dim];\n";
630dc007f05SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "; i++) {\n";
6318b97b69aSJeremy L Thompson           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
6328b97b69aSJeremy L Thompson           code << "    }\n";
6338b97b69aSJeremy L Thompson         } else if (use_3d_slices) {
6344b3e95d5SJeremy L Thompson           // Accumulator for gradient slices
6354b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
6364b3e95d5SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n";
6374b3e95d5SJeremy L Thompson           code << "      r_q" << var_suffix << "[i] = 0.0;\n";
6384b3e95d5SJeremy L Thompson           code << "    }\n";
6394b3e95d5SJeremy L Thompson         } else {
640dc007f05SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
6414b3e95d5SJeremy L Thompson         }
6428b97b69aSJeremy L Thompson         break;
6438b97b69aSJeremy L Thompson       case CEED_EVAL_WEIGHT:
6448b97b69aSJeremy L Thompson         break;
6458b97b69aSJeremy L Thompson         // LCOV_EXCL_START
6468b97b69aSJeremy L Thompson       case CEED_EVAL_DIV:
6478b97b69aSJeremy L Thompson       case CEED_EVAL_CURL:
6488b97b69aSJeremy L Thompson         break;  // TODO: Not implemented
6498b97b69aSJeremy L Thompson                 // LCOV_EXCL_STOP
6504b3e95d5SJeremy L Thompson     }
6514b3e95d5SJeremy L Thompson   }
6524b3e95d5SJeremy L Thompson 
6538b97b69aSJeremy L Thompson   if (is_at_points) {
6548b97b69aSJeremy L Thompson     // We need to handle batches of points
6558b97b69aSJeremy L Thompson     code << "\n    // Note: Using batches of points\n";
6568b97b69aSJeremy L Thompson     code << "    const CeedInt point_loop_bound = (blockDim.x * blockDim.y) * ceil(1.0 * max_num_points / (blockDim.x * blockDim.y));\n\n";
6578b97b69aSJeremy L Thompson     code << "    #pragma unroll\n";
6588b97b69aSJeremy L Thompson     code << "    for (CeedInt i = threadIdx.x + threadIdx.y * blockDim.x; i < point_loop_bound; i += blockDim.x * blockDim.y) {\n";
6598b97b69aSJeremy L Thompson     code << "      const CeedInt p = i % max_num_points;\n\n";
6608b97b69aSJeremy L Thompson 
6618b97b69aSJeremy L Thompson     code << "      // -- Coordinates\n";
6628b97b69aSJeremy L Thompson     code << "      CeedScalar r_x[dim];\n";
6638b97b69aSJeremy 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";
6648b97b69aSJeremy L Thompson 
6658b97b69aSJeremy L Thompson     code << "      // -- Input fields\n";
6668b97b69aSJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
6678b97b69aSJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
6688b97b69aSJeremy L Thompson 
6698b97b69aSJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
6708b97b69aSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
6718b97b69aSJeremy L Thompson       // Basis action
6728b97b69aSJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
6738b97b69aSJeremy L Thompson       switch (eval_mode) {
6748b97b69aSJeremy L Thompson         case CEED_EVAL_NONE:
6758b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
6768b97b69aSJeremy L Thompson           code << "      ReadPoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
6778b97b69aSJeremy L Thompson                << ", max_num_points>(data, elem, p, max_num_points, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
6788b97b69aSJeremy L Thompson           break;
6798b97b69aSJeremy L Thompson         case CEED_EVAL_INTERP:
6808b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
6818b97b69aSJeremy L Thompson           code << "      InterpAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix
6828b97b69aSJeremy L Thompson                << ", r_x, r_s" << var_suffix << ");\n";
6838b97b69aSJeremy L Thompson           break;
6848b97b69aSJeremy L Thompson         case CEED_EVAL_GRAD:
6858b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
6868b97b69aSJeremy L Thompson           code << "      GradAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix
6878b97b69aSJeremy L Thompson                << ", r_x, r_s" << var_suffix << ");\n";
6888b97b69aSJeremy L Thompson           break;
6898b97b69aSJeremy L Thompson         case CEED_EVAL_WEIGHT:
6908b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
6918b97b69aSJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = 1.0;\n";
6928b97b69aSJeremy L Thompson           break;
6938b97b69aSJeremy L Thompson           // LCOV_EXCL_START
6948b97b69aSJeremy L Thompson         case CEED_EVAL_DIV:
6958b97b69aSJeremy L Thompson         case CEED_EVAL_CURL:
6968b97b69aSJeremy L Thompson           break;  // TODO: Not implemented
6978b97b69aSJeremy L Thompson                   // LCOV_EXCL_STOP
6988b97b69aSJeremy L Thompson       }
6998b97b69aSJeremy L Thompson     }
7008b97b69aSJeremy L Thompson     code << "\n      // -- Output fields\n";
7018b97b69aSJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
7028b97b69aSJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
7038b97b69aSJeremy L Thompson 
7048b97b69aSJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
7058b97b69aSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
7068b97b69aSJeremy L Thompson       // Basis action
7078b97b69aSJeremy L Thompson       switch (eval_mode) {
7088b97b69aSJeremy L Thompson         case CEED_EVAL_NONE:
7098b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7108b97b69aSJeremy L Thompson           break;
7118b97b69aSJeremy L Thompson         case CEED_EVAL_INTERP:
7128b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7138b97b69aSJeremy L Thompson           break;
7148b97b69aSJeremy L Thompson         case CEED_EVAL_GRAD:
7158b97b69aSJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
7168b97b69aSJeremy L Thompson           break;
7178b97b69aSJeremy L Thompson           // LCOV_EXCL_START
7188b97b69aSJeremy L Thompson         case CEED_EVAL_WEIGHT:
7198b97b69aSJeremy L Thompson           break;  // Should not occur
7208b97b69aSJeremy L Thompson         case CEED_EVAL_DIV:
7218b97b69aSJeremy L Thompson         case CEED_EVAL_CURL:
7228b97b69aSJeremy L Thompson           break;  // TODO: Not implemented
7238b97b69aSJeremy L Thompson                   // LCOV_EXCL_STOP
7248b97b69aSJeremy L Thompson       }
7258b97b69aSJeremy L Thompson     }
7268b97b69aSJeremy L Thompson 
7278b97b69aSJeremy L Thompson   } else if (use_3d_slices) {
7284b3e95d5SJeremy L Thompson     // We treat quadrature points per slice in 3d to save registers
7294b3e95d5SJeremy L Thompson     code << "\n    // Note: Using planes of 3D elements\n";
7304b3e95d5SJeremy L Thompson     code << "    #pragma unroll\n";
7314b3e95d5SJeremy L Thompson     code << "    for (CeedInt q = 0; q < " << Q_name << "; q++) {\n";
7324b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
7334b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
7344b3e95d5SJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
7354b3e95d5SJeremy L Thompson 
7364b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
7374b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
7384b3e95d5SJeremy L Thompson       // Basis action
7394b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
7404b3e95d5SJeremy L Thompson       switch (eval_mode) {
7414b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
7424b3e95d5SJeremy L Thompson           bool is_strided;
7434b3e95d5SJeremy L Thompson 
7444b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7454b3e95d5SJeremy L Thompson 
7464b3e95d5SJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
7474b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
7484b3e95d5SJeremy L Thompson           if (is_strided) {
7494b3e95d5SJeremy L Thompson             bool    has_backend_strides;
7504b3e95d5SJeremy L Thompson             CeedInt num_elem, elem_size;
7514b3e95d5SJeremy L Thompson 
7524b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
7534b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
7544b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
7554b3e95d5SJeremy L Thompson             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
7564b3e95d5SJeremy L Thompson 
7574b3e95d5SJeremy L Thompson             if (!has_backend_strides) {
7584b3e95d5SJeremy L Thompson               CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
7594b3e95d5SJeremy L Thompson             }
7604b3e95d5SJeremy L Thompson             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
761f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStrided3d<num_comp" << var_suffix << ", " << Q_name << ", " << strides[0] << ", " << strides[1] << ", "
7624b3e95d5SJeremy L Thompson                  << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n";
7634b3e95d5SJeremy L Thompson           } else {
7644b3e95d5SJeremy L Thompson             CeedSize                  l_size = 0;
7654b3e95d5SJeremy L Thompson             CeedInt                   comp_stride;
7664b3e95d5SJeremy L Thompson             CeedElemRestriction_Cuda *rstr_data;
7674b3e95d5SJeremy L Thompson 
7684b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
7694b3e95d5SJeremy L Thompson             code << "      const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
7704b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
7714b3e95d5SJeremy L Thompson             code << "      // CompStride: " << comp_stride << "\n";
7724b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
7734b3e95d5SJeremy L Thompson             data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
774f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStandard3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix
7754b3e95d5SJeremy L Thompson                  << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
7764b3e95d5SJeremy L Thompson           }
777681d0ea7SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
7784b3e95d5SJeremy L Thompson           break;
7794b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
7804b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7814b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n";
7824b3e95d5SJeremy L Thompson           code << "        r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n";
7834b3e95d5SJeremy L Thompson           code << "      }\n";
7844b3e95d5SJeremy L Thompson           break;
7854b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
7864b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
787f815fac9SJeremy L Thompson           code << "      GradColloSlice3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_q" << var_suffix << ", s_G" << var_suffix
788f815fac9SJeremy L Thompson                << ", r_s" << var_suffix << ");\n";
7894b3e95d5SJeremy L Thompson           break;
7904b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
7914b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
7924b3e95d5SJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n";
7938b97b69aSJeremy L Thompson           break;
7944b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
7954b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
7964b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
7974b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
7984b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
7994b3e95d5SJeremy L Thompson       }
8004b3e95d5SJeremy L Thompson     }
8014b3e95d5SJeremy L Thompson     code << "\n      // -- Output fields\n";
8024b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
8034b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
8044b3e95d5SJeremy L Thompson 
8054b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
8064b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
8074b3e95d5SJeremy L Thompson       // Basis action
8084b3e95d5SJeremy L Thompson       switch (eval_mode) {
8094b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
8104b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
8118b97b69aSJeremy L Thompson           break;
8124b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
8134b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
8144b3e95d5SJeremy L Thompson           break;
8154b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
8164b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
8174b3e95d5SJeremy L Thompson           break;
8184b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
8194b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
8204b3e95d5SJeremy L Thompson           break;  // Should not occur
8214b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
8224b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
8234b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
8244b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
8254b3e95d5SJeremy L Thompson       }
8264b3e95d5SJeremy L Thompson     }
8274b3e95d5SJeremy L Thompson   } else {
8284b3e95d5SJeremy L Thompson     code << "\n    // Note: Using full elements\n";
8294b3e95d5SJeremy L Thompson     code << "    {\n";
8304b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
8314b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
8324b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
8334b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n";
8344b3e95d5SJeremy L Thompson     }
8354b3e95d5SJeremy L Thompson     code << "      // -- Output fields\n";
8364b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
8374b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
8384b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n";
8394b3e95d5SJeremy L Thompson     }
8404b3e95d5SJeremy L Thompson   }
8414b3e95d5SJeremy L Thompson 
8424b3e95d5SJeremy L Thompson   // Input and output buffers
8434b3e95d5SJeremy L Thompson   code << "\n      // -- QFunction inputs and outputs\n";
8444b3e95d5SJeremy L Thompson   code << "      // ---- Inputs\n";
8454b3e95d5SJeremy L Thompson   code << "      CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n";
8464b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
8474b3e95d5SJeremy L Thompson     code << "      // ------ Input field " << i << "\n";
8484b3e95d5SJeremy L Thompson     code << "      inputs[" << i << "] = r_s_in_" << i << ";\n";
8494b3e95d5SJeremy L Thompson   }
8504b3e95d5SJeremy L Thompson   code << "      // ---- Outputs\n";
8514b3e95d5SJeremy L Thompson   code << "      CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n";
8524b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
8534b3e95d5SJeremy L Thompson     code << "      // ------ Output field " << i << "\n";
8544b3e95d5SJeremy L Thompson     code << "      outputs[" << i << "] = r_s_out_" << i << ";\n";
8554b3e95d5SJeremy L Thompson   }
8564b3e95d5SJeremy L Thompson 
8574b3e95d5SJeremy L Thompson   // Apply QFunction
8584b3e95d5SJeremy L Thompson   code << "\n      // -- Apply QFunction\n";
8594b3e95d5SJeremy L Thompson   code << "      " << qfunction_name << "(ctx, ";
860dc007f05SJeremy L Thompson   if (dim != 3 || is_at_points || use_3d_slices || !is_tensor) {
8614b3e95d5SJeremy L Thompson     code << "1";
8624b3e95d5SJeremy L Thompson   } else {
863dc007f05SJeremy L Thompson     code << Q_name;
8644b3e95d5SJeremy L Thompson   }
8654b3e95d5SJeremy L Thompson   code << ", inputs, outputs);\n";
8664b3e95d5SJeremy L Thompson 
8678b97b69aSJeremy L Thompson   if (is_at_points) {
8688b97b69aSJeremy L Thompson     // Map back to coefficients
8698b97b69aSJeremy L Thompson     code << "\n      // -- Output fields\n";
8708b97b69aSJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
8718b97b69aSJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
8728b97b69aSJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
8738b97b69aSJeremy L Thompson 
8748b97b69aSJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
8758b97b69aSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
8768b97b69aSJeremy L Thompson       // Basis action
8778b97b69aSJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
8788b97b69aSJeremy L Thompson       switch (eval_mode) {
8798b97b69aSJeremy L Thompson         case CEED_EVAL_NONE: {
8808b97b69aSJeremy L Thompson           CeedInt             comp_stride;
8818b97b69aSJeremy L Thompson           CeedElemRestriction elem_rstr;
8828b97b69aSJeremy L Thompson 
8838b97b69aSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
8848b97b69aSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
8858b97b69aSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
8868b97b69aSJeremy L Thompson           code << "      const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
8878b97b69aSJeremy L Thompson           code << "      WritePoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
8888b97b69aSJeremy L Thompson                << ", max_num_points>(data, elem, i, points.num_per_elem[elem], indices.outputs[" << i << "]"
8898b97b69aSJeremy L Thompson                << ", r_s" << var_suffix << ", d" << var_suffix << ");\n";
8908b97b69aSJeremy L Thompson           break;
8918b97b69aSJeremy L Thompson         }
8928b97b69aSJeremy L Thompson         case CEED_EVAL_INTERP:
8938b97b69aSJeremy L Thompson           code << "      if (i >= points.num_per_elem[elem]) {\n";
8948b97b69aSJeremy L Thompson           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) r_s" << var_suffix << "[j] = 0.0;\n";
8958b97b69aSJeremy L Thompson           code << "      }\n";
8968b97b69aSJeremy L Thompson           code << "      InterpTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s"
8978b97b69aSJeremy L Thompson                << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
8988b97b69aSJeremy L Thompson           break;
8998b97b69aSJeremy L Thompson         case CEED_EVAL_GRAD:
9008b97b69aSJeremy L Thompson           code << "      if (i >= points.num_per_elem[elem]) {\n";
9018b97b69aSJeremy L Thompson           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "*dim; j++) r_s" << var_suffix << "[j] = 0.0;\n";
9028b97b69aSJeremy L Thompson           code << "      }\n";
9038b97b69aSJeremy L Thompson           code << "      GradTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s"
9048b97b69aSJeremy L Thompson                << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
9058b97b69aSJeremy L Thompson           break;
9068b97b69aSJeremy L Thompson           // LCOV_EXCL_START
9078b97b69aSJeremy L Thompson         case CEED_EVAL_WEIGHT:
9088b97b69aSJeremy L Thompson           break;  // Should not occur
9098b97b69aSJeremy L Thompson         case CEED_EVAL_DIV:
9108b97b69aSJeremy L Thompson         case CEED_EVAL_CURL:
9118b97b69aSJeremy L Thompson           break;  // TODO: Not implemented
9128b97b69aSJeremy L Thompson                   // LCOV_EXCL_STOP
9138b97b69aSJeremy L Thompson       }
9148b97b69aSJeremy L Thompson     }
9158b97b69aSJeremy L Thompson   } else if (use_3d_slices) {
9164b3e95d5SJeremy L Thompson     // Copy or apply transpose grad, if needed
9178b97b69aSJeremy L Thompson     code << "\n      // -- Output fields\n";
9184b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
9194b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
9204b3e95d5SJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
9214b3e95d5SJeremy L Thompson 
9224b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
9234b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
9244b3e95d5SJeremy L Thompson       // Basis action
9254b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
9264b3e95d5SJeremy L Thompson       switch (eval_mode) {
9274b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
9284b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
9294b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
9304b3e95d5SJeremy L Thompson           code << "      }\n";
9318b97b69aSJeremy L Thompson           break;
9324b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
9334b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
9344b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
9354b3e95d5SJeremy L Thompson           code << "      }\n";
9364b3e95d5SJeremy L Thompson           break;
9374b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
938f815fac9SJeremy L Thompson           code << "      GradColloSliceTranspose3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_s" << var_suffix << ", s_G"
939f815fac9SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
9404b3e95d5SJeremy L Thompson           break;
9414b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
9424b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
9434b3e95d5SJeremy L Thompson           break;  // Should not occur
9444b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
9454b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
9464b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
9474b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
9484b3e95d5SJeremy L Thompson       }
9494b3e95d5SJeremy L Thompson     }
9504b3e95d5SJeremy L Thompson   }
9514b3e95d5SJeremy L Thompson   code << "    }\n";
9524b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9534b3e95d5SJeremy L Thompson }
9544b3e95d5SJeremy L Thompson 
9554b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
956b2165e7aSSebastian Grimberg // Build single operator kernel
957ab213215SJeremy L Thompson //------------------------------------------------------------------------------
958ddae5012SJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Cuda_gen(CeedOperator op, bool *is_good_build) {
9598b97b69aSJeremy L Thompson   bool                    is_tensor = true, is_at_points = false, use_3d_slices = false;
960241a4b83SYohann   Ceed                    ceed;
9618b97b69aSJeremy L Thompson   CeedInt                 Q_1d, num_input_fields, num_output_fields, dim = 1, max_num_points = 0, coords_comp_stride = 0;
962ca735530SJeremy L Thompson   CeedQFunctionField     *qf_input_fields, *qf_output_fields;
963ca735530SJeremy L Thompson   CeedQFunction_Cuda_gen *qf_data;
964ca735530SJeremy L Thompson   CeedQFunction           qf;
965ca735530SJeremy L Thompson   CeedOperatorField      *op_input_fields, *op_output_fields;
966ca735530SJeremy L Thompson   CeedOperator_Cuda_gen  *data;
9674b3e95d5SJeremy L Thompson   std::ostringstream      code;
9684b3e95d5SJeremy L Thompson 
969ddae5012SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
9704b3e95d5SJeremy L Thompson   {
9714b3e95d5SJeremy L Thompson     bool is_setup_done;
972ca735530SJeremy L Thompson 
973ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
974ddae5012SJeremy L Thompson     if (is_setup_done) {
975ddae5012SJeremy L Thompson       *is_good_build = !data->use_fallback;
976ddae5012SJeremy L Thompson       return CEED_ERROR_SUCCESS;
977ddae5012SJeremy L Thompson     }
978ddae5012SJeremy L Thompson   }
979ddae5012SJeremy L Thompson 
980ddae5012SJeremy L Thompson   // Check field compatibility
9818d12f40eSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
982ddae5012SJeremy L Thompson   {
983ddae5012SJeremy L Thompson     bool has_shared_bases = true, is_all_tensor = true, is_all_nontensor = true;
984ddae5012SJeremy L Thompson 
985ddae5012SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
986ddae5012SJeremy L Thompson       CeedBasis basis;
987ddae5012SJeremy L Thompson 
988ddae5012SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
989ddae5012SJeremy L Thompson       if (basis != CEED_BASIS_NONE) {
990ddae5012SJeremy L Thompson         bool        is_tensor = true;
991ddae5012SJeremy L Thompson         const char *resource;
992ddae5012SJeremy L Thompson         char       *resource_root;
993ddae5012SJeremy L Thompson         Ceed        basis_ceed;
994ddae5012SJeremy L Thompson 
995ddae5012SJeremy L Thompson         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
996c9192acaSJeremy L Thompson         is_all_tensor    = is_all_tensor && is_tensor;
997*45a787f7SJeremy L Thompson         is_all_nontensor = is_all_nontensor && !is_tensor;
998ddae5012SJeremy L Thompson         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
999ddae5012SJeremy L Thompson         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
1000ddae5012SJeremy L Thompson         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
1001c9192acaSJeremy L Thompson         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/cuda/shared");
1002ddae5012SJeremy L Thompson         CeedCallBackend(CeedFree(&resource_root));
1003ddae5012SJeremy L Thompson         CeedCallBackend(CeedDestroy(&basis_ceed));
1004ddae5012SJeremy L Thompson       }
1005ddae5012SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis));
10064b3e95d5SJeremy L Thompson     }
1007ca735530SJeremy L Thompson 
1008ddae5012SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
1009ddae5012SJeremy L Thompson       CeedBasis basis;
1010ddae5012SJeremy L Thompson 
1011ddae5012SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1012ddae5012SJeremy L Thompson       if (basis != CEED_BASIS_NONE) {
1013ddae5012SJeremy L Thompson         bool        is_tensor = true;
1014ddae5012SJeremy L Thompson         const char *resource;
1015ddae5012SJeremy L Thompson         char       *resource_root;
1016ddae5012SJeremy L Thompson         Ceed        basis_ceed;
1017ddae5012SJeremy L Thompson 
1018ddae5012SJeremy L Thompson         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
1019c9192acaSJeremy L Thompson         is_all_tensor    = is_all_tensor && is_tensor;
1020c9192acaSJeremy L Thompson         is_all_nontensor = is_all_nontensor && !is_tensor;
1021ddae5012SJeremy L Thompson 
1022ddae5012SJeremy L Thompson         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
1023ddae5012SJeremy L Thompson         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
1024ddae5012SJeremy L Thompson         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
1025c9192acaSJeremy L Thompson         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/cuda/shared");
1026ddae5012SJeremy L Thompson         CeedCallBackend(CeedFree(&resource_root));
1027ddae5012SJeremy L Thompson         CeedCallBackend(CeedDestroy(&basis_ceed));
1028ddae5012SJeremy L Thompson       }
1029ddae5012SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis));
1030ddae5012SJeremy L Thompson     }
1031ddae5012SJeremy L Thompson     // -- Fallback to ref if not all bases are shared
1032ddae5012SJeremy L Thompson     if (!has_shared_bases || (!is_all_tensor && !is_all_nontensor)) {
1033ddae5012SJeremy L Thompson       *is_good_build = false;
1034ddae5012SJeremy L Thompson       return CEED_ERROR_SUCCESS;
1035ddae5012SJeremy L Thompson     }
1036ddae5012SJeremy L Thompson   }
1037ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1038ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1039ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
1040ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
1041241a4b83SYohann 
10424b3e95d5SJeremy L Thompson   // Get operator data
10438b97b69aSJeremy L Thompson   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
10444b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelData_Cuda_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, op_output_fields,
10454b3e95d5SJeremy L Thompson                                                        qf_output_fields, &data->max_P_1d, &Q_1d, &dim, &is_tensor, &use_3d_slices));
10464b3e95d5SJeremy L Thompson   if (dim == 0) dim = 1;
10474b3e95d5SJeremy L Thompson   data->dim = dim;
10488b97b69aSJeremy L Thompson   if (is_at_points) {
10498b97b69aSJeremy L Thompson     CeedElemRestriction_Cuda *rstr_data;
10508b97b69aSJeremy L Thompson     CeedElemRestriction       rstr_points = NULL;
10514b3e95d5SJeremy L Thompson 
10528b97b69aSJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
10538b97b69aSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points));
10548b97b69aSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetCompStride(rstr_points, &coords_comp_stride));
10558b97b69aSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(rstr_points, &rstr_data));
10568b97b69aSJeremy L Thompson     data->points.indices = (CeedInt *)rstr_data->d_offsets;
10578b97b69aSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
10588b97b69aSJeremy L Thompson   }
10598b97b69aSJeremy L Thompson   if (is_at_points) use_3d_slices = false;
10608b97b69aSJeremy L Thompson   if (Q_1d == 0) {
10618b97b69aSJeremy L Thompson     if (is_at_points) Q_1d = max_num_points;
10628b97b69aSJeremy L Thompson     else CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q_1d));
10634b3e95d5SJeremy L Thompson   }
10644b3e95d5SJeremy L Thompson   data->Q_1d = Q_1d;
10654b3e95d5SJeremy L Thompson 
10660b454692Sjeremylt   // Check for restriction only identity operator
10674b3e95d5SJeremy L Thompson   {
10684b3e95d5SJeremy L Thompson     bool is_identity_qf;
10694b3e95d5SJeremy L Thompson 
10702b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
10710b454692Sjeremylt     if (is_identity_qf) {
10729e201c85SYohann       CeedEvalMode eval_mode_in, eval_mode_out;
1073ca735530SJeremy L Thompson 
10742b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
10752b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
10766574a04fSJeremy L Thompson       CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
10776574a04fSJeremy L Thompson                 "Backend does not implement restriction only identity operators");
10780b454692Sjeremylt     }
10794b3e95d5SJeremy L Thompson   }
10800b454692Sjeremylt 
1081f1a13f77SYohann Dudouit   // Add atomicAdd function for old NVidia architectures
10824b3e95d5SJeremy L Thompson   {
10834b3e95d5SJeremy L Thompson     Ceed_Cuda            *ceed_data;
10844b3e95d5SJeremy L Thompson     struct cudaDeviceProp prop;
10854b3e95d5SJeremy L Thompson 
10862b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetData(ceed, &ceed_data));
10872b730f8bSJeremy L Thompson     CeedCallBackend(cudaGetDeviceProperties(&prop, ceed_data->device_id));
108880a9ef05SNatalie Beams     if ((prop.major < 6) && (CEED_SCALAR_TYPE != CEED_SCALAR_FP32)) {
10899c25dd66SJeremy L Thompson       code << "// AtomicAdd fallback source\n";
10909c25dd66SJeremy L Thompson       code << "#include <ceed/jit-source/cuda/cuda-atomic-add-fallback.h>\n\n";
1091f1a13f77SYohann Dudouit     }
10924b3e95d5SJeremy L Thompson   }
1093f1a13f77SYohann Dudouit 
10949e201c85SYohann   // Load basis source files
1095dc007f05SJeremy L Thompson   if (is_tensor) {
10969c25dd66SJeremy L Thompson     code << "// Tensor basis source\n";
10979c25dd66SJeremy L Thompson     code << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-templates.h>\n\n";
1098dc007f05SJeremy L Thompson   } else {
1099dc007f05SJeremy L Thompson     code << "// Non-tensor basis source\n";
1100dc007f05SJeremy L Thompson     code << "#include <ceed/jit-source/cuda/cuda-shared-basis-nontensor-templates.h>\n\n";
1101dc007f05SJeremy L Thompson   }
1102dc007f05SJeremy L Thompson   if (is_at_points) {
11038b97b69aSJeremy L Thompson     code << "// AtPoints basis source\n";
11048b97b69aSJeremy L Thompson     code << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-at-points-templates.h>\n\n";
1105dc007f05SJeremy L Thompson   }
11069c25dd66SJeremy L Thompson   code << "// CodeGen operator source\n";
11079c25dd66SJeremy L Thompson   code << "#include <ceed/jit-source/cuda/cuda-gen-templates.h>\n\n";
1108241a4b83SYohann 
11094b3e95d5SJeremy L Thompson   // Get QFunction name
11104b3e95d5SJeremy L Thompson   std::string qfunction_name(qf_data->qfunction_name);
11114b3e95d5SJeremy L Thompson   std::string operator_name;
11124b3e95d5SJeremy L Thompson 
111309095acaSJeremy L Thompson   operator_name = "CeedKernelCudaGenOperator_" + qfunction_name;
11144d537eeaSYohann 
11159e201c85SYohann   // Define CEED_Q_VLA
11169e201c85SYohann   code << "\n#undef CEED_Q_VLA\n";
1117dc007f05SJeremy L Thompson   if (dim != 3 || is_at_points || use_3d_slices || !is_tensor) {
11189e201c85SYohann     code << "#define CEED_Q_VLA 1\n\n";
11199e201c85SYohann   } else {
11209e201c85SYohann     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
11219e201c85SYohann   }
11229e201c85SYohann 
11234b3e95d5SJeremy L Thompson   // Add user QFunction source
11244b3e95d5SJeremy L Thompson   {
11259c25dd66SJeremy L Thompson     const char *source_path;
11264b3e95d5SJeremy L Thompson 
11279c25dd66SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path));
11289c25dd66SJeremy L Thompson     CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/cuda/gen backend requires QFunction source code file");
11299c25dd66SJeremy L Thompson 
11309c25dd66SJeremy L Thompson     code << "// User QFunction source\n";
11319c25dd66SJeremy L Thompson     code << "#include \"" << source_path << "\"\n\n";
11324b3e95d5SJeremy L Thompson   }
1133241a4b83SYohann 
1134241a4b83SYohann   // Setup
1135818e0025SJeremy L Thompson   code << "\n// -----------------------------------------------------------------------------\n";
11364b3e95d5SJeremy L Thompson   code << "// Operator Kernel\n";
11374b3e95d5SJeremy L Thompson   code << "// \n";
11384b3e95d5SJeremy L Thompson   code << "// d_[in,out]_i:   CeedVector device array\n";
11394b3e95d5SJeremy L Thompson   code << "// r_[in,out]_e_i: Element vector register\n";
11404b3e95d5SJeremy L Thompson   code << "// r_[in,out]_q_i: Quadrature space vector register\n";
11419123fb08SJeremy L Thompson   code << "// r_[in,out]_c_i: AtPoints Chebyshev coefficients register\n";
11424b3e95d5SJeremy L Thompson   code << "// r_[in,out]_s_i: Quadrature space slice vector register\n";
11434b3e95d5SJeremy L Thompson   code << "// \n";
11444b3e95d5SJeremy L Thompson   code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n";
11454b3e95d5SJeremy L Thompson   code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n";
11464b3e95d5SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n";
11474b3e95d5SJeremy L Thompson   code << "extern \"C\" __global__ void " << operator_name
11488b97b69aSJeremy L Thompson        << "(CeedInt num_elem, void* ctx, FieldsInt_Cuda indices, Fields_Cuda fields, Fields_Cuda B, Fields_Cuda G, CeedScalar *W, Points_Cuda "
11498b97b69aSJeremy L Thompson           "points) {\n";
11504b3e95d5SJeremy L Thompson 
11514b3e95d5SJeremy L Thompson   // Scratch buffers
11529e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
11534b3e95d5SJeremy L Thompson     CeedEvalMode eval_mode;
11544b3e95d5SJeremy L Thompson 
11552b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
11569e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
11574b3e95d5SJeremy L Thompson       code << "  const CeedScalar *d_in_" << i << " = fields.inputs[" << i << "];\n";
1158241a4b83SYohann     }
1159241a4b83SYohann   }
11609e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
11614b3e95d5SJeremy L Thompson     code << "  CeedScalar *d_out_" << i << " = fields.outputs[" << i << "];\n";
1162241a4b83SYohann   }
11631da99368SJeremy L Thompson 
11649e201c85SYohann   code << "  const CeedInt dim = " << dim << ";\n";
1165dc007f05SJeremy L Thompson   code << "  const CeedInt " << (is_tensor ? "Q_1d" : "Q") << " = " << Q_1d << ";\n";
11668b97b69aSJeremy L Thompson   if (is_at_points) {
11678b97b69aSJeremy L Thompson     code << "  const CeedInt max_num_points = " << max_num_points << ";\n";
11688b97b69aSJeremy L Thompson     code << "  const CeedInt coords_comp_stride = " << coords_comp_stride << ";\n";
11698b97b69aSJeremy L Thompson   }
11701da99368SJeremy L Thompson 
11714b3e95d5SJeremy L Thompson   // Shared data
1172241a4b83SYohann   code << "  extern __shared__ CeedScalar slice[];\n";
11739e201c85SYohann   code << "  SharedData_Cuda data;\n";
11749e201c85SYohann   code << "  data.t_id_x = threadIdx.x;\n";
11759e201c85SYohann   code << "  data.t_id_y = threadIdx.y;\n";
11769e201c85SYohann   code << "  data.t_id_z = threadIdx.z;\n";
11779e201c85SYohann   code << "  data.t_id  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
1178dc007f05SJeremy L Thompson   code << "  data.slice = slice + data.t_id_z*T_1D" << ((!is_tensor || dim == 1) ? "" : "*T_1D") << ";\n";
1179920dcdc4Sjeremylt 
11800a2a6492SJeremy L Thompson   // -- Determine input mat reuse
1181*45a787f7SJeremy L Thompson   FieldReuse_Cuda input_matrix_reuse[CEED_FIELD_MAX];
11820a2a6492SJeremy L Thompson 
11830a2a6492SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1184*45a787f7SJeremy L Thompson     input_matrix_reuse[i].index = -1;
11850a2a6492SJeremy L Thompson   }
11860a2a6492SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
11870a2a6492SJeremy L Thompson     CeedEvalMode eval_mode_i;
11880a2a6492SJeremy L Thompson     CeedBasis    basis_i;
11890a2a6492SJeremy L Thompson 
11900a2a6492SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode_i));
11910a2a6492SJeremy L Thompson     if (eval_mode_i == CEED_EVAL_WEIGHT) continue;
11920a2a6492SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis_i));
1193*45a787f7SJeremy L Thompson     for (CeedInt j = 0; (input_matrix_reuse[i].index == -1) && (j < i); j++) {
11940a2a6492SJeremy L Thompson       CeedEvalMode eval_mode_j;
11950a2a6492SJeremy L Thompson       CeedBasis    basis_j;
11960a2a6492SJeremy L Thompson 
11970a2a6492SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
11980a2a6492SJeremy L Thompson       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
11990a2a6492SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
12000a2a6492SJeremy L Thompson       if (basis_i == basis_j) {
12010a2a6492SJeremy L Thompson         if (is_tensor) {
1202*45a787f7SJeremy L Thompson           input_matrix_reuse[i].index     = j;
1203*45a787f7SJeremy L Thompson           input_matrix_reuse[i].is_input  = true;
1204*45a787f7SJeremy L Thompson           input_matrix_reuse[i].eval_mode = eval_mode_j;
12050a2a6492SJeremy L Thompson         } else {
12060a2a6492SJeremy L Thompson           // For non-tensor can only re-use with the same eval mode
12070a2a6492SJeremy L Thompson           if (eval_mode_i == eval_mode_j) {
1208*45a787f7SJeremy L Thompson             input_matrix_reuse[i].index     = j;
1209*45a787f7SJeremy L Thompson             input_matrix_reuse[i].is_input  = true;
1210*45a787f7SJeremy L Thompson             input_matrix_reuse[i].eval_mode = eval_mode_j;
12110a2a6492SJeremy L Thompson           }
12120a2a6492SJeremy L Thompson         }
12130a2a6492SJeremy L Thompson       }
12140a2a6492SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis_j));
12150a2a6492SJeremy L Thompson     }
12160a2a6492SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis_i));
12170a2a6492SJeremy L Thompson   }
12180a2a6492SJeremy L Thompson 
12190a2a6492SJeremy L Thompson   // -- Determine output mat reuse
1220*45a787f7SJeremy L Thompson   FieldReuse_Cuda output_matrix_reuse[CEED_FIELD_MAX];
12210a2a6492SJeremy L Thompson 
12220a2a6492SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1223*45a787f7SJeremy L Thompson     output_matrix_reuse[i].index = -1;
12240a2a6492SJeremy L Thompson   }
12250a2a6492SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
12260a2a6492SJeremy L Thompson     CeedEvalMode eval_mode_i;
12270a2a6492SJeremy L Thompson     CeedBasis    basis_i;
12280a2a6492SJeremy L Thompson 
12290a2a6492SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode_i));
12300a2a6492SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis_i));
1231*45a787f7SJeremy L Thompson     for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < num_input_fields); j++) {
12320a2a6492SJeremy L Thompson       CeedEvalMode eval_mode_j;
12330a2a6492SJeremy L Thompson       CeedBasis    basis_j;
12340a2a6492SJeremy L Thompson 
12350a2a6492SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
12360a2a6492SJeremy L Thompson       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
12370a2a6492SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
12380a2a6492SJeremy L Thompson       if (basis_i == basis_j) {
12390a2a6492SJeremy L Thompson         if (is_tensor) {
1240*45a787f7SJeremy L Thompson           output_matrix_reuse[i].index     = j;
1241*45a787f7SJeremy L Thompson           output_matrix_reuse[i].is_input  = true;
1242*45a787f7SJeremy L Thompson           output_matrix_reuse[i].eval_mode = eval_mode_j;
12430a2a6492SJeremy L Thompson         } else {
12440a2a6492SJeremy L Thompson           // For non-tensor can only re-use with the same eval mode
12450a2a6492SJeremy L Thompson           if (eval_mode_i == eval_mode_j) {
1246*45a787f7SJeremy L Thompson             output_matrix_reuse[i].index     = j;
1247*45a787f7SJeremy L Thompson             output_matrix_reuse[i].is_input  = true;
1248*45a787f7SJeremy L Thompson             output_matrix_reuse[i].eval_mode = eval_mode_j;
12490a2a6492SJeremy L Thompson           }
12500a2a6492SJeremy L Thompson         }
12510a2a6492SJeremy L Thompson       }
12520a2a6492SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis_j));
12530a2a6492SJeremy L Thompson     }
1254*45a787f7SJeremy L Thompson     for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < i); j++) {
12550a2a6492SJeremy L Thompson       CeedEvalMode eval_mode_j;
12560a2a6492SJeremy L Thompson       CeedBasis    basis_j;
12570a2a6492SJeremy L Thompson 
12580a2a6492SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[j], &eval_mode_j));
12590a2a6492SJeremy L Thompson       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
12600a2a6492SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[j], &basis_j));
12610a2a6492SJeremy L Thompson       if (basis_i == basis_j) {
12620a2a6492SJeremy L Thompson         if (is_tensor) {
1263*45a787f7SJeremy L Thompson           output_matrix_reuse[i].index     = j;
1264*45a787f7SJeremy L Thompson           output_matrix_reuse[i].is_input  = false;
1265*45a787f7SJeremy L Thompson           output_matrix_reuse[i].eval_mode = eval_mode_j;
12660a2a6492SJeremy L Thompson         } else {
12670a2a6492SJeremy L Thompson           // For non-tensor can only re-use with the same eval mode
12680a2a6492SJeremy L Thompson           if (eval_mode_i == eval_mode_j) {
1269*45a787f7SJeremy L Thompson             output_matrix_reuse[i].index     = j;
1270*45a787f7SJeremy L Thompson             output_matrix_reuse[i].is_input  = false;
1271*45a787f7SJeremy L Thompson             output_matrix_reuse[i].eval_mode = eval_mode_j;
12720a2a6492SJeremy L Thompson           }
12730a2a6492SJeremy L Thompson         }
12740a2a6492SJeremy L Thompson       }
12750a2a6492SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis_j));
12760a2a6492SJeremy L Thompson     }
12770a2a6492SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis_i));
12780a2a6492SJeremy L Thompson   }
12790a2a6492SJeremy L Thompson 
1280ac421f39SYohann   // Initialize constants, and matrices B and G
12814b3e95d5SJeremy L Thompson   code << "\n  // Input field constants and basis data\n";
12829e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
12830a2a6492SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, i, op_input_fields[i], qf_input_fields[i], input_matrix_reuse[i], Q_1d,
12840a2a6492SJeremy L Thompson                                                               true, is_tensor, is_at_points, use_3d_slices));
1285920dcdc4Sjeremylt   }
12864b3e95d5SJeremy L Thompson   code << "\n  // Output field constants and basis data\n";
12879e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
12880a2a6492SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, i, op_output_fields[i], qf_output_fields[i], output_matrix_reuse[i], Q_1d,
12890a2a6492SJeremy L Thompson                                                               false, is_tensor, is_at_points, use_3d_slices));
12904b3e95d5SJeremy L Thompson   }
1291920dcdc4Sjeremylt 
12924b3e95d5SJeremy L Thompson   // Loop over all elements
12934b3e95d5SJeremy L Thompson   code << "\n  // Element loop\n";
1294ac421f39SYohann   code << "  __syncthreads();\n";
12959e201c85SYohann   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
12964b3e95d5SJeremy L Thompson 
1297e93651e5SJeremy L Thompson   // -- Compute minimum buffer space needed
12988b97b69aSJeremy L Thompson   CeedInt max_rstr_buffer_size = 1;
1299e93651e5SJeremy L Thompson 
13009e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1301e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
1302e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
1303e93651e5SJeremy L Thompson 
1304e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
1305e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1306e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1307dc007f05SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_tensor && (dim >= 3) ? elem_size : 1));
1308681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1309e93651e5SJeremy L Thompson   }
1310e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1311e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
1312e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
1313e93651e5SJeremy L Thompson 
1314e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
1315e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1316e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1317dc007f05SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_tensor && (dim >= 3) ? elem_size : 1));
1318681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1319e93651e5SJeremy L Thompson   }
1320e93651e5SJeremy L Thompson   code << "    // Scratch restriction buffer space\n";
1321e93651e5SJeremy L Thompson   code << "    CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n";
1322e93651e5SJeremy L Thompson 
1323e93651e5SJeremy L Thompson   // -- Determine best input field processing order
1324e93651e5SJeremy L Thompson   CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX];
1325e93651e5SJeremy L Thompson 
1326e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1327e93651e5SJeremy L Thompson     field_rstr_in_buffer[i] = -1;
1328e93651e5SJeremy L Thompson     input_field_order[i]    = -1;
1329e93651e5SJeremy L Thompson   }
1330e93651e5SJeremy L Thompson   {
1331e93651e5SJeremy L Thompson     bool    is_ordered[CEED_FIELD_MAX];
1332e93651e5SJeremy L Thompson     CeedInt curr_index = 0;
1333e93651e5SJeremy L Thompson 
1334e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
1335e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
1336e93651e5SJeremy L Thompson       CeedVector          vec_i;
1337e93651e5SJeremy L Thompson       CeedElemRestriction rstr_i;
1338e93651e5SJeremy L Thompson 
1339e93651e5SJeremy L Thompson       if (is_ordered[i]) continue;
1340e93651e5SJeremy L Thompson       field_rstr_in_buffer[i]       = i;
1341e93651e5SJeremy L Thompson       is_ordered[i]                 = true;
1342e93651e5SJeremy L Thompson       input_field_order[curr_index] = i;
1343e93651e5SJeremy L Thompson       curr_index++;
1344034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
1345e93651e5SJeremy L Thompson       if (vec_i == CEED_VECTOR_NONE) continue;  // CEED_EVAL_WEIGHT
1346e93651e5SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
1347e93651e5SJeremy L Thompson       for (CeedInt j = i + 1; j < num_input_fields; j++) {
1348e93651e5SJeremy L Thompson         CeedVector          vec_j;
1349e93651e5SJeremy L Thompson         CeedElemRestriction rstr_j;
1350e93651e5SJeremy L Thompson 
1351e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
1352e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
1353e93651e5SJeremy L Thompson         if (rstr_i == rstr_j && vec_i == vec_j) {
1354e93651e5SJeremy L Thompson           field_rstr_in_buffer[j]       = i;
1355e93651e5SJeremy L Thompson           is_ordered[j]                 = true;
1356e93651e5SJeremy L Thompson           input_field_order[curr_index] = j;
1357e93651e5SJeremy L Thompson           curr_index++;
1358e93651e5SJeremy L Thompson         }
1359681d0ea7SJeremy L Thompson         CeedCallBackend(CeedVectorDestroy(&vec_j));
1360681d0ea7SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
1361e93651e5SJeremy L Thompson       }
1362681d0ea7SJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec_i));
1363681d0ea7SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
1364e93651e5SJeremy L Thompson     }
1365e93651e5SJeremy L Thompson   }
1366e93651e5SJeremy L Thompson 
1367e93651e5SJeremy L Thompson   // -- Input restriction and basis
1368e93651e5SJeremy L Thompson   code << "\n    // -- Input field restrictions and basis actions\n";
1369e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1370e93651e5SJeremy L Thompson     CeedInt f = input_field_order[i];
1371e93651e5SJeremy L Thompson 
1372e93651e5SJeremy L Thompson     code << "    // ---- Input field " << f << "\n";
1373920dcdc4Sjeremylt 
13744b3e95d5SJeremy L Thompson     // ---- Restriction
1375e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, f, dim, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f],
1376dc007f05SJeremy L Thompson                                                                 Q_1d, true, is_tensor, is_at_points, use_3d_slices));
13779a2291e3SJeremy L Thompson 
13784b3e95d5SJeremy L Thompson     // ---- Basis action
1379dc007f05SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, f, dim, op_input_fields[f], qf_input_fields[f], Q_1d, true, is_tensor,
1380dc007f05SJeremy L Thompson                                                           is_at_points, use_3d_slices));
1381920dcdc4Sjeremylt   }
1382920dcdc4Sjeremylt 
13834b3e95d5SJeremy L Thompson   // -- Q function
13848b97b69aSJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelQFunction_Cuda_gen(code, data, dim, max_num_points, num_input_fields, op_input_fields, qf_input_fields,
1385dc007f05SJeremy L Thompson                                                             num_output_fields, op_output_fields, qf_output_fields, qfunction_name, Q_1d, is_tensor,
1386dc007f05SJeremy L Thompson                                                             is_at_points, use_3d_slices));
1387ca735530SJeremy L Thompson 
13884b3e95d5SJeremy L Thompson   // -- Output basis and restriction
13894b3e95d5SJeremy L Thompson   code << "\n    // -- Output field basis action and restrictions\n";
13909e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
13914b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
1392ca735530SJeremy L Thompson 
13934b3e95d5SJeremy L Thompson     // ---- Basis action
1394dc007f05SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, i, dim, op_output_fields[i], qf_output_fields[i], Q_1d, false, is_tensor,
1395dc007f05SJeremy L Thompson                                                           is_at_points, use_3d_slices));
13969a2291e3SJeremy L Thompson 
13974b3e95d5SJeremy L Thompson     // ---- Restriction
13988b97b69aSJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, i, dim, NULL, op_output_fields[i], qf_output_fields[i], Q_1d, false,
1399dc007f05SJeremy L Thompson                                                                 is_tensor, is_at_points, use_3d_slices));
1400ac421f39SYohann   }
1401241a4b83SYohann 
14024b3e95d5SJeremy L Thompson   // Close loop and function
1403241a4b83SYohann   code << "  }\n";
1404818e0025SJeremy L Thompson   code << "}\n";
1405818e0025SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n\n";
1406241a4b83SYohann 
14073a2968d6SJeremy L Thompson   // Compile
1408ddae5012SJeremy L Thompson   {
1409ddae5012SJeremy L Thompson     bool is_compile_good = false;
1410ddae5012SJeremy L Thompson 
1411ddae5012SJeremy 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)));
1412ddae5012SJeremy L Thompson     if (is_compile_good) {
1413ddae5012SJeremy L Thompson       *is_good_build = true;
1414eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, operator_name.c_str(), &data->op));
1415ddae5012SJeremy L Thompson     } else {
1416ddae5012SJeremy L Thompson       *is_good_build     = false;
1417ddae5012SJeremy L Thompson       data->use_fallback = true;
1418ddae5012SJeremy L Thompson     }
1419ddae5012SJeremy L Thompson   }
14202b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
14219bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
1422c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionDestroy(&qf));
1423e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1424241a4b83SYohann }
14252a86cc9dSSebastian Grimberg 
1426ab213215SJeremy L Thompson //------------------------------------------------------------------------------
1427