xref: /libCEED/rust/libceed-sys/c-src/backends/hip-gen/ceed-hip-gen-operator-build.cpp (revision 9123fb08d52f01bdd0d1f3a790ba84e4ab900e9f)
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.
37d8d0e25Snbeams //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
57d8d0e25Snbeams //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
73d576824SJeremy L Thompson 
87d8d0e25Snbeams #define CEED_DEBUG_COLOR 12
97d8d0e25Snbeams 
1049aac155SJeremy L Thompson #include <ceed.h>
11ec3da8bcSJed Brown #include <ceed/backend.h>
129e201c85SYohann #include <ceed/jit-tools.h>
132b730f8bSJeremy L Thompson 
147d8d0e25Snbeams #include <iostream>
157d8d0e25Snbeams #include <sstream>
162b730f8bSJeremy L Thompson #include <string>
172b730f8bSJeremy L Thompson 
180d0321e0SJeremy L Thompson #include "../hip-ref/ceed-hip-ref.h"
197d8d0e25Snbeams #include "../hip-shared/ceed-hip-shared.h"
20b2165e7aSSebastian Grimberg #include "../hip/ceed-hip-common.h"
217d8d0e25Snbeams #include "../hip/ceed-hip-compile.h"
222b730f8bSJeremy L Thompson #include "ceed-hip-gen.h"
23b3e1519bSnbeams 
24b3e1519bSnbeams //------------------------------------------------------------------------------
25b3e1519bSnbeams // Calculate the block size used for launching the operator kernel
26b3e1519bSnbeams //------------------------------------------------------------------------------
272b730f8bSJeremy L Thompson extern "C" int BlockGridCalculate_Hip_gen(const CeedInt dim, const CeedInt num_elem, const CeedInt P_1d, const CeedInt Q_1d, CeedInt *block_sizes) {
283a2968d6SJeremy L Thompson   const CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
29b3e1519bSnbeams   if (dim == 1) {
303a2968d6SJeremy L Thompson     CeedInt elems_per_block = 64 * thread_1d > 256 ? 256 / thread_1d : 64;
31b7453713SJeremy L Thompson 
329e201c85SYohann     elems_per_block = elems_per_block > 0 ? elems_per_block : 1;
333a2968d6SJeremy L Thompson     block_sizes[0]  = thread_1d;
34b3e1519bSnbeams     block_sizes[1]  = 1;
359e201c85SYohann     block_sizes[2]  = elems_per_block;
36b3e1519bSnbeams   } else if (dim == 2) {
373a2968d6SJeremy L Thompson     const CeedInt elems_per_block = thread_1d < 4 ? 16 : 2;
38b7453713SJeremy L Thompson 
393a2968d6SJeremy L Thompson     block_sizes[0] = thread_1d;
403a2968d6SJeremy L Thompson     block_sizes[1] = thread_1d;
419e201c85SYohann     block_sizes[2] = elems_per_block;
42b3e1519bSnbeams   } else if (dim == 3) {
433a2968d6SJeremy L Thompson     const CeedInt elems_per_block = thread_1d < 6 ? 4 : (thread_1d < 8 ? 2 : 1);
44b7453713SJeremy L Thompson 
453a2968d6SJeremy L Thompson     block_sizes[0] = thread_1d;
463a2968d6SJeremy L Thompson     block_sizes[1] = thread_1d;
479e201c85SYohann     block_sizes[2] = elems_per_block;
48b3e1519bSnbeams   }
49b3e1519bSnbeams   return CEED_ERROR_SUCCESS;
50b3e1519bSnbeams }
51b3e1519bSnbeams 
527d8d0e25Snbeams //------------------------------------------------------------------------------
534b3e95d5SJeremy L Thompson // Determine type of operator
544b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
554b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelData_Hip_gen(Ceed ceed, CeedInt num_input_fields, CeedOperatorField *op_input_fields,
564b3e95d5SJeremy L Thompson                                                CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, CeedOperatorField *op_output_fields,
574b3e95d5SJeremy L Thompson                                                CeedQFunctionField *qf_output_fields, CeedInt *max_P_1d, CeedInt *Q_1d, CeedInt *dim, bool *is_tensor,
584b3e95d5SJeremy L Thompson                                                bool *use_3d_slices) {
594b3e95d5SJeremy L Thompson   // Find dim, P_1d, Q_1d
604b3e95d5SJeremy L Thompson   *max_P_1d  = 0;
614b3e95d5SJeremy L Thompson   *Q_1d      = 0;
624b3e95d5SJeremy L Thompson   *dim       = 0;
634b3e95d5SJeremy L Thompson   *is_tensor = true;
64*9123fb08SJeremy L Thompson 
654b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
664b3e95d5SJeremy L Thompson     CeedBasis basis;
674b3e95d5SJeremy L Thompson 
684b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
694b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
704b3e95d5SJeremy L Thompson       bool    is_field_tensor;
714b3e95d5SJeremy L Thompson       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
724b3e95d5SJeremy L Thompson 
734b3e95d5SJeremy L Thompson       // Collect dim, P_1d, and Q_1d
744b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
754b3e95d5SJeremy L Thompson       *is_tensor = *is_tensor && is_field_tensor;
76*9123fb08SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
77*9123fb08SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P_1d));
784b3e95d5SJeremy L Thompson       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
794b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
804b3e95d5SJeremy L Thompson       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
814b3e95d5SJeremy L Thompson       *dim = field_dim;
82*9123fb08SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
83*9123fb08SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q_1d));
844b3e95d5SJeremy L Thompson       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
854b3e95d5SJeremy L Thompson       *Q_1d = field_Q_1d;
864b3e95d5SJeremy L Thompson     }
873a2968d6SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis));
884b3e95d5SJeremy L Thompson   }
894b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
904b3e95d5SJeremy L Thompson     CeedBasis basis;
914b3e95d5SJeremy L Thompson 
924b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
934b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
944b3e95d5SJeremy L Thompson       bool    is_field_tensor;
954b3e95d5SJeremy L Thompson       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
964b3e95d5SJeremy L Thompson 
974b3e95d5SJeremy L Thompson       // Collect dim, P_1d, and Q_1d
984b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
994b3e95d5SJeremy L Thompson       *is_tensor = *is_tensor && is_field_tensor;
100*9123fb08SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
101*9123fb08SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P_1d));
1024b3e95d5SJeremy L Thompson       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
1034b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
1044b3e95d5SJeremy L Thompson       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
1054b3e95d5SJeremy L Thompson       *dim = field_dim;
106*9123fb08SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
107*9123fb08SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q_1d));
1084b3e95d5SJeremy L Thompson       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
1094b3e95d5SJeremy L Thompson       *Q_1d = field_Q_1d;
1104b3e95d5SJeremy L Thompson     }
1113a2968d6SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis));
1124b3e95d5SJeremy L Thompson   }
1134b3e95d5SJeremy L Thompson 
1144b3e95d5SJeremy L Thompson   // Only use 3D collocated gradient parallelization strategy when gradient is computed
1154b3e95d5SJeremy L Thompson   *use_3d_slices = false;
1164b3e95d5SJeremy L Thompson   if (*dim == 3) {
1174b3e95d5SJeremy L Thompson     bool was_grad_found = false;
1184b3e95d5SJeremy L Thompson 
1194b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
1204b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
1214b3e95d5SJeremy L Thompson 
1224b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1234b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
1244b3e95d5SJeremy L Thompson         CeedBasis_Hip_shared *basis_data;
1254b3e95d5SJeremy L Thompson         CeedBasis             basis;
1264b3e95d5SJeremy L Thompson 
1274b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1284b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1294b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
1304b3e95d5SJeremy L Thompson         was_grad_found = true;
1313a2968d6SJeremy L Thompson         CeedCallBackend(CeedBasisDestroy(&basis));
1324b3e95d5SJeremy L Thompson       }
1334b3e95d5SJeremy L Thompson     }
1344b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
1354b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
1364b3e95d5SJeremy L Thompson 
1374b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1384b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
1394b3e95d5SJeremy L Thompson         CeedBasis_Hip_shared *basis_data;
1404b3e95d5SJeremy L Thompson         CeedBasis             basis;
1414b3e95d5SJeremy L Thompson 
1424b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1434b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1444b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
1454b3e95d5SJeremy L Thompson         was_grad_found = true;
1463a2968d6SJeremy L Thompson         CeedCallBackend(CeedBasisDestroy(&basis));
1474b3e95d5SJeremy L Thompson       }
1484b3e95d5SJeremy L Thompson     }
1494b3e95d5SJeremy L Thompson   }
1504b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1514b3e95d5SJeremy L Thompson }
1524b3e95d5SJeremy L Thompson 
1534b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1544b3e95d5SJeremy L Thompson // Setup fields
1554b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1564b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelFieldData_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedOperatorField op_field,
157*9123fb08SJeremy L Thompson                                                     CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool is_tensor, bool is_at_points,
158*9123fb08SJeremy L Thompson                                                     bool use_3d_slices) {
1594b3e95d5SJeremy L Thompson   std::string           var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
160*9123fb08SJeremy L Thompson   std::string           P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q";
1614b3e95d5SJeremy L Thompson   std::string           option_name = (is_input ? "inputs" : "outputs");
1624b3e95d5SJeremy L Thompson   CeedEvalMode          eval_mode   = CEED_EVAL_NONE;
1634b3e95d5SJeremy L Thompson   CeedInt               elem_size = 0, num_comp = 0, P_1d = 0;
1644b3e95d5SJeremy L Thompson   CeedElemRestriction   elem_rstr;
1654b3e95d5SJeremy L Thompson   CeedBasis_Hip_shared *basis_data;
1664b3e95d5SJeremy L Thompson   CeedBasis             basis;
1674b3e95d5SJeremy L Thompson 
1684b3e95d5SJeremy L Thompson   code << "  // -- " << (is_input ? "Input" : "Output") << " field " << i << "\n";
1694b3e95d5SJeremy L Thompson 
1704b3e95d5SJeremy L Thompson   // Get field data
1714b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
1724b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
1734b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1744b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1754b3e95d5SJeremy L Thompson   }
1763a2968d6SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1774b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
1784b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
1794b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetData(basis, &basis_data));
180*9123fb08SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
181*9123fb08SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
1824b3e95d5SJeremy L Thompson   }
1834b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
1844b3e95d5SJeremy L Thompson 
1854b3e95d5SJeremy L Thompson   // Set field constants
1864b3e95d5SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
1874b3e95d5SJeremy L Thompson     code << "  const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n";
1884b3e95d5SJeremy L Thompson     code << "  const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n";
1894b3e95d5SJeremy L Thompson   }
1904b3e95d5SJeremy L Thompson 
1914b3e95d5SJeremy L Thompson   // Load basis data
1924b3e95d5SJeremy L Thompson   code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
1934b3e95d5SJeremy L Thompson   switch (eval_mode) {
1944b3e95d5SJeremy L Thompson     case CEED_EVAL_NONE:
1954b3e95d5SJeremy L Thompson       break;
1964b3e95d5SJeremy L Thompson     case CEED_EVAL_INTERP:
1973a2968d6SJeremy L Thompson       if (is_at_points) {
1983a2968d6SJeremy L Thompson         // AtPoints
1993a2968d6SJeremy L Thompson         if (!basis_data->d_chebyshev_interp_1d) {
2003a2968d6SJeremy L Thompson           CeedSize    interp_bytes;
2013a2968d6SJeremy L Thompson           CeedScalar *chebyshev_interp_1d;
2023a2968d6SJeremy L Thompson 
2033a2968d6SJeremy L Thompson           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
2043a2968d6SJeremy L Thompson           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
2053a2968d6SJeremy L Thompson           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
2063a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis), hipMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
2073a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis),
2083a2968d6SJeremy L Thompson                       hipMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, hipMemcpyHostToDevice));
2093a2968d6SJeremy L Thompson           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
2103a2968d6SJeremy L Thompson         }
2113a2968d6SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
2123a2968d6SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
2133a2968d6SJeremy L Thompson       } else {
2143a2968d6SJeremy L Thompson         // Standard quadrature
2154b3e95d5SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
2164b3e95d5SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_interp_1d;
2173a2968d6SJeremy L Thompson       }
218*9123fb08SJeremy L Thompson       code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
219f815fac9SJeremy L Thompson       code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
2204b3e95d5SJeremy L Thompson       break;
2214b3e95d5SJeremy L Thompson     case CEED_EVAL_GRAD:
2223a2968d6SJeremy L Thompson       if (is_at_points) {
2233a2968d6SJeremy L Thompson         // AtPoints
2243a2968d6SJeremy L Thompson         if (!basis_data->d_chebyshev_interp_1d) {
2253a2968d6SJeremy L Thompson           CeedSize    interp_bytes;
2263a2968d6SJeremy L Thompson           CeedScalar *chebyshev_interp_1d;
2273a2968d6SJeremy L Thompson 
2283a2968d6SJeremy L Thompson           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
2293a2968d6SJeremy L Thompson           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
2303a2968d6SJeremy L Thompson           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
2313a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis), hipMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
2323a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis),
2333a2968d6SJeremy L Thompson                       hipMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, hipMemcpyHostToDevice));
2343a2968d6SJeremy L Thompson           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
2353a2968d6SJeremy L Thompson         }
2363a2968d6SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
2373a2968d6SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
2383a2968d6SJeremy L Thompson       } else {
2393a2968d6SJeremy L Thompson         // Standard quadrature
2404b3e95d5SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
2414b3e95d5SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_interp_1d;
2423a2968d6SJeremy L Thompson       }
243*9123fb08SJeremy L Thompson       if (is_tensor) {
244*9123fb08SJeremy L Thompson         code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
245f815fac9SJeremy L Thompson         code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
246*9123fb08SJeremy L Thompson       }
2473a2968d6SJeremy L Thompson       if (is_at_points) break;  // No G mat for AtPoints
2484b3e95d5SJeremy L Thompson       if (use_3d_slices) {
2494b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d;
2504b3e95d5SJeremy L Thompson         else data->G.outputs[i] = basis_data->d_collo_grad_1d;
251*9123fb08SJeremy L Thompson         code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
252f815fac9SJeremy L Thompson         code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
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*9123fb08SJeremy L Thompson           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
260f815fac9SJeremy L Thompson           code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
2614b3e95d5SJeremy L Thompson         } else {
262*9123fb08SJeremy L Thompson           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << P_name << "*" << Q_name << (is_tensor ? "" : "*dim") << "];\n";
263*9123fb08SJeremy L Thompson           code << "  LoadMatrix<" << P_name << ", " << Q_name << (is_tensor ? "" : "*dim") << ">(data, G." << option_name << "[" << i << "], s_G"
264*9123fb08SJeremy L Thompson                << var_suffix << ");\n";
2654b3e95d5SJeremy L Thompson         }
2664b3e95d5SJeremy L Thompson       }
2674b3e95d5SJeremy L Thompson       break;
2684b3e95d5SJeremy L Thompson     case CEED_EVAL_WEIGHT:
2694b3e95d5SJeremy L Thompson       break;  // No action
2704b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
2714b3e95d5SJeremy L Thompson     case CEED_EVAL_DIV:
2724b3e95d5SJeremy L Thompson     case CEED_EVAL_CURL:
2734b3e95d5SJeremy L Thompson       break;  // TODO: Not implemented
2744b3e95d5SJeremy L Thompson               // LCOV_EXCL_STOP
2754b3e95d5SJeremy L Thompson   }
2763a2968d6SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
2774b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2784b3e95d5SJeremy L Thompson }
2794b3e95d5SJeremy L Thompson 
2804b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
2814b3e95d5SJeremy L Thompson // Restriction
2824b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
2834b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelRestriction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim,
284e93651e5SJeremy L Thompson                                                       CeedInt field_input_buffer[], CeedOperatorField op_field, CeedQFunctionField qf_field,
285*9123fb08SJeremy L Thompson                                                       CeedInt Q_1d, bool is_input, bool is_tensor, bool is_at_points, bool use_3d_slices) {
2864b3e95d5SJeremy L Thompson   std::string              var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
287*9123fb08SJeremy L Thompson   std::string              P_name     = (is_tensor ? "P_1d" : "P") + var_suffix;
2884b3e95d5SJeremy L Thompson   CeedEvalMode             eval_mode  = CEED_EVAL_NONE;
2894b3e95d5SJeremy L Thompson   CeedInt                  elem_size = 0, num_comp = 0, P_1d = 0;
2904b3e95d5SJeremy L Thompson   CeedSize                 l_size;
291f815fac9SJeremy L Thompson   CeedRestrictionType      rstr_type = CEED_RESTRICTION_STANDARD;
2924b3e95d5SJeremy L Thompson   CeedElemRestriction_Hip *rstr_data;
2934b3e95d5SJeremy L Thompson   CeedElemRestriction      elem_rstr;
2944b3e95d5SJeremy L Thompson   CeedBasis                basis;
2954b3e95d5SJeremy L Thompson 
2964b3e95d5SJeremy L Thompson   // Get field data
2974b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
2984b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
299f815fac9SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
3004b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
3014b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
3024b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
3034b3e95d5SJeremy L Thompson   }
3044b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
3054b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
306*9123fb08SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
307*9123fb08SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
3084b3e95d5SJeremy L Thompson   }
3093a2968d6SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
3104b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
3114b3e95d5SJeremy L Thompson 
3124b3e95d5SJeremy L Thompson   // Restriction
3134b3e95d5SJeremy L Thompson   if (is_input) {
3144b3e95d5SJeremy L Thompson     // Input
315e93651e5SJeremy L Thompson     if (field_input_buffer[i] != i) {
316e93651e5SJeremy L Thompson       std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]);
317e93651e5SJeremy L Thompson 
318e93651e5SJeremy L Thompson       // Restriction was already done for previous input
319e93651e5SJeremy L Thompson       code << "    CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n";
3203a2968d6SJeremy L Thompson     } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices && is_at_points)) {
3213a2968d6SJeremy L Thompson       if (eval_mode == CEED_EVAL_NONE && rstr_type != CEED_RESTRICTION_POINTS) {
322e93651e5SJeremy L Thompson         // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated
3234b3e95d5SJeremy L Thompson         code << "    CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n";
3243a2968d6SJeremy L Thompson       } else if (rstr_type != CEED_RESTRICTION_POINTS) {
325e93651e5SJeremy L Thompson         // Otherwise we're using the scratch space
326e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
327e93651e5SJeremy L Thompson       }
328f815fac9SJeremy L Thompson       switch (rstr_type) {
329f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STANDARD: {
3304b3e95d5SJeremy L Thompson           CeedInt comp_stride;
3314b3e95d5SJeremy L Thompson 
3324b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
3334b3e95d5SJeremy L Thompson           code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
3344b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3354b3e95d5SJeremy L Thompson           code << "    // CompStride: " << comp_stride << "\n";
3364b3e95d5SJeremy L Thompson           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
337*9123fb08SJeremy L Thompson           code << "    ReadLVecStandard" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name
338*9123fb08SJeremy L Thompson                << ">(data, l_size" << var_suffix << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n";
339f815fac9SJeremy L Thompson           break;
340f815fac9SJeremy L Thompson         }
341f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STRIDED: {
3424b3e95d5SJeremy L Thompson           bool    has_backend_strides;
3434b3e95d5SJeremy L Thompson           CeedInt num_elem;
3444b3e95d5SJeremy L Thompson 
3454b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
3464b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
3474b3e95d5SJeremy L Thompson           CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
3484b3e95d5SJeremy L Thompson 
3494b3e95d5SJeremy L Thompson           if (!has_backend_strides) {
3504b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
3514b3e95d5SJeremy L Thompson           }
3524b3e95d5SJeremy L Thompson           code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
353*9123fb08SJeremy L Thompson           code << "    ReadLVecStrided" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", "
354*9123fb08SJeremy L Thompson                << strides[1] << ", " << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n";
355f815fac9SJeremy L Thompson           break;
356f815fac9SJeremy L Thompson         }
3573a2968d6SJeremy L Thompson         case CEED_RESTRICTION_POINTS: {
3583a2968d6SJeremy L Thompson           CeedInt comp_stride;
3593a2968d6SJeremy L Thompson 
3603a2968d6SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3613a2968d6SJeremy L Thompson           code << "    const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
3623a2968d6SJeremy L Thompson           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
3633a2968d6SJeremy L Thompson           break;
3643a2968d6SJeremy L Thompson         }
365f815fac9SJeremy L Thompson         // LCOV_EXCL_START
366f815fac9SJeremy L Thompson         case CEED_RESTRICTION_ORIENTED:
367f815fac9SJeremy L Thompson         case CEED_RESTRICTION_CURL_ORIENTED:
368f815fac9SJeremy L Thompson           break;  // TODO: Not implemented
369f815fac9SJeremy L Thompson                   // LCOV_EXCL_STOP
3704b3e95d5SJeremy L Thompson       }
3714b3e95d5SJeremy L Thompson     }
3724b3e95d5SJeremy L Thompson   } else {
3734b3e95d5SJeremy L Thompson     // Output
374f815fac9SJeremy L Thompson     switch (rstr_type) {
375f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STANDARD: {
3764b3e95d5SJeremy L Thompson         CeedInt comp_stride;
3774b3e95d5SJeremy L Thompson 
3784b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
3794b3e95d5SJeremy L Thompson         code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
3804b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3814b3e95d5SJeremy L Thompson         code << "    // CompStride: " << comp_stride << "\n";
3824b3e95d5SJeremy L Thompson         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
383*9123fb08SJeremy L Thompson         code << "    WriteLVecStandard" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name
384*9123fb08SJeremy L Thompson              << ">(data, l_size" << var_suffix << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n";
385f815fac9SJeremy L Thompson         break;
386f815fac9SJeremy L Thompson       }
387f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STRIDED: {
3884b3e95d5SJeremy L Thompson         bool    has_backend_strides;
3894b3e95d5SJeremy L Thompson         CeedInt num_elem;
3904b3e95d5SJeremy L Thompson 
3914b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
3924b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
3934b3e95d5SJeremy L Thompson         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
3944b3e95d5SJeremy L Thompson 
3954b3e95d5SJeremy L Thompson         if (!has_backend_strides) {
3964b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
3974b3e95d5SJeremy L Thompson         }
3984b3e95d5SJeremy L Thompson         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
399*9123fb08SJeremy L Thompson         code << "    WriteLVecStrided" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", "
400*9123fb08SJeremy L Thompson              << strides[1] << ", " << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n";
401f815fac9SJeremy L Thompson         break;
402f815fac9SJeremy L Thompson       }
4033a2968d6SJeremy L Thompson       case CEED_RESTRICTION_POINTS:
4043a2968d6SJeremy L Thompson         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
4053a2968d6SJeremy L Thompson         break;
406f815fac9SJeremy L Thompson       // LCOV_EXCL_START
407f815fac9SJeremy L Thompson       case CEED_RESTRICTION_ORIENTED:
408f815fac9SJeremy L Thompson       case CEED_RESTRICTION_CURL_ORIENTED:
409f815fac9SJeremy L Thompson         break;  // TODO: Not implemented
410f815fac9SJeremy L Thompson                 // LCOV_EXCL_STOP
4114b3e95d5SJeremy L Thompson     }
4124b3e95d5SJeremy L Thompson   }
4133a2968d6SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
4144b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4154b3e95d5SJeremy L Thompson }
4164b3e95d5SJeremy L Thompson 
4174b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4184b3e95d5SJeremy L Thompson // Basis
4194b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4204b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelBasis_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim,
421*9123fb08SJeremy L Thompson                                                 CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool is_tensor,
4223a2968d6SJeremy L Thompson                                                 bool is_at_points, bool use_3d_slices) {
4234b3e95d5SJeremy L Thompson   std::string         var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
424*9123fb08SJeremy L Thompson   std::string         P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q";
4254b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
4264b3e95d5SJeremy L Thompson   CeedInt             elem_size = 0, num_comp = 0, P_1d = 0;
4274b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
4284b3e95d5SJeremy L Thompson   CeedBasis           basis;
4294b3e95d5SJeremy L Thompson 
4304b3e95d5SJeremy L Thompson   // Get field data
4314b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
4324b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
4334b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
4344b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
4354b3e95d5SJeremy L Thompson   }
4363a2968d6SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
4374b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
4384b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
439*9123fb08SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
440*9123fb08SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
4414b3e95d5SJeremy L Thompson   }
4424b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
4434b3e95d5SJeremy L Thompson 
4444b3e95d5SJeremy L Thompson   // Basis
4454b3e95d5SJeremy L Thompson   code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
4464b3e95d5SJeremy L Thompson   if (is_input) {
4474b3e95d5SJeremy L Thompson     switch (eval_mode) {
4484b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
4493a2968d6SJeremy L Thompson         if (!use_3d_slices && !is_at_points) {
4504b3e95d5SJeremy L Thompson           code << "    CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n";
4514b3e95d5SJeremy L Thompson         }
4524b3e95d5SJeremy L Thompson         break;
4534b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
4543a2968d6SJeremy L Thompson         if (is_at_points) {
455*9123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
456*9123fb08SJeremy L Thompson 
457*9123fb08SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
458*9123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
459*9123fb08SJeremy L Thompson                << var_suffix << ", r_c" << var_suffix << ");\n";
4603a2968d6SJeremy L Thompson         } else {
461*9123fb08SJeremy L Thompson           std::string function_name = is_tensor ? ((dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d") : "InterpNonTensor";
462*9123fb08SJeremy L Thompson 
463*9123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
464*9123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
465*9123fb08SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
4663a2968d6SJeremy L Thompson         }
4674b3e95d5SJeremy L Thompson         break;
4684b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
4693a2968d6SJeremy L Thompson         if (is_at_points) {
470*9123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
471*9123fb08SJeremy L Thompson 
472*9123fb08SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
473*9123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
474*9123fb08SJeremy L Thompson                << var_suffix << ", r_c" << var_suffix << ");\n";
4753a2968d6SJeremy L Thompson         } else if (use_3d_slices) {
476*9123fb08SJeremy L Thompson           std::string function_name = (dim > 1 ? "InterpTensor" : "Interp") + std::to_string(dim) + "d";
477*9123fb08SJeremy L Thompson 
4784b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
479*9123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
480*9123fb08SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
481*9123fb08SJeremy L Thompson         } else if (is_tensor) {
482*9123fb08SJeremy L Thompson           bool        is_collocated = dim == 3 && Q_1d >= P_1d;
483*9123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "Grad" : (is_collocated ? "GradTensorCollocated" : "GradTensor")) + std::to_string(dim) + "d";
484*9123fb08SJeremy L Thompson 
485*9123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << (dim >= 3 ? Q_name : "1") << "];\n";
486*9123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
487*9123fb08SJeremy L Thompson                << var_suffix << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
4884b3e95d5SJeremy L Thompson         } else {
489*9123fb08SJeremy L Thompson           std::string function_name = "GradNonTensor";
490*9123fb08SJeremy L Thompson 
491*9123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
492*9123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", dim, " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix
493*9123fb08SJeremy L Thompson                << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
4944b3e95d5SJeremy L Thompson         }
4954b3e95d5SJeremy L Thompson         break;
4964b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
4973a2968d6SJeremy L Thompson         if (is_at_points) {
4983a2968d6SJeremy L Thompson           code << "    // Nothing to do AtPoints\n";
4993a2968d6SJeremy L Thompson         } else {
5004b3e95d5SJeremy L Thompson           CeedBasis_Hip_shared *basis_data;
501*9123fb08SJeremy L Thompson           std::string           function_name = is_tensor ? ((dim == 1 ? "Weight" : "WeightTensor") + std::to_string(dim) + "d") : "WeightNonTensor";
5024b3e95d5SJeremy L Thompson 
503*9123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
5044b3e95d5SJeremy L Thompson           CeedCallBackend(CeedBasisGetData(basis, &basis_data));
5054b3e95d5SJeremy L Thompson           data->W = basis_data->d_q_weight_1d;
506*9123fb08SJeremy L Thompson           code << "    " << function_name << "<" << Q_name << ">(data, W, r_q" << var_suffix << ");\n";
5073a2968d6SJeremy L Thompson         }
5084b3e95d5SJeremy L Thompson         break;
5094b3e95d5SJeremy L Thompson       }
5104b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
5114b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
5124b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
5134b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
5144b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
5154b3e95d5SJeremy L Thompson     }
5164b3e95d5SJeremy L Thompson   } else {
5174b3e95d5SJeremy L Thompson     switch (eval_mode) {
5184b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
5194b3e95d5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n";
5204b3e95d5SJeremy L Thompson         break;  // No action
5214b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
522e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
5233a2968d6SJeremy L Thompson         if (is_at_points) {
524*9123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
525*9123fb08SJeremy L Thompson 
526*9123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_c" << var_suffix << ", s_B"
527*9123fb08SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5283a2968d6SJeremy L Thompson         } else {
529*9123fb08SJeremy L Thompson           std::string function_name =
530*9123fb08SJeremy L Thompson               is_tensor ? ((dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d") : "InterpTransposeNonTensor";
531*9123fb08SJeremy L Thompson 
532*9123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B"
533*9123fb08SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5343a2968d6SJeremy L Thompson         }
5354b3e95d5SJeremy L Thompson         break;
5364b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
537e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
5383a2968d6SJeremy L Thompson         if (is_at_points) {
539*9123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
540*9123fb08SJeremy L Thompson 
541*9123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_c" << var_suffix << ", s_B"
542*9123fb08SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5433a2968d6SJeremy L Thompson         } else if (use_3d_slices) {
544*9123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
545*9123fb08SJeremy L Thompson 
546*9123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B"
547*9123fb08SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
548*9123fb08SJeremy L Thompson         } else if (is_tensor) {
549*9123fb08SJeremy L Thompson           bool        is_collocated = dim == 3 && Q_1d >= P_1d;
550*9123fb08SJeremy L Thompson           std::string function_name =
551*9123fb08SJeremy L Thompson               (dim == 1 ? "GradTranspose" : (is_collocated ? "GradTransposeTensorCollocated" : "GradTransposeTensor")) + std::to_string(dim) + "d";
552*9123fb08SJeremy L Thompson 
553*9123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B"
554*9123fb08SJeremy L Thompson                << var_suffix << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
5554b3e95d5SJeremy L Thompson         } else {
556*9123fb08SJeremy L Thompson           std::string function_name = "GradTransposeNonTensor";
557*9123fb08SJeremy L Thompson 
558*9123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", dim, " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix
559*9123fb08SJeremy L Thompson                << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
5604b3e95d5SJeremy L Thompson         }
5614b3e95d5SJeremy L Thompson         break;
5624b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
5634b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT:
5644b3e95d5SJeremy L Thompson         break;  // Should not occur
5654b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
5664b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
5674b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
5684b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
5694b3e95d5SJeremy L Thompson     }
5704b3e95d5SJeremy L Thompson   }
5713a2968d6SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
5724b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5734b3e95d5SJeremy L Thompson }
5744b3e95d5SJeremy L Thompson 
5754b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
5764b3e95d5SJeremy L Thompson // QFunction
5774b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
5783a2968d6SJeremy L Thompson static int CeedOperatorBuildKernelQFunction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt dim, CeedInt max_num_points,
5793a2968d6SJeremy L Thompson                                                     CeedInt num_input_fields, CeedOperatorField *op_input_fields, CeedQFunctionField *qf_input_fields,
5804b3e95d5SJeremy L Thompson                                                     CeedInt num_output_fields, CeedOperatorField *op_output_fields,
581*9123fb08SJeremy L Thompson                                                     CeedQFunctionField *qf_output_fields, std::string qfunction_name, CeedInt Q_1d, bool is_tensor,
582*9123fb08SJeremy L Thompson                                                     bool is_at_points, bool use_3d_slices) {
583*9123fb08SJeremy L Thompson   std::string         Q_name    = is_tensor ? "Q_1d" : "Q";
5844b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
5854b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
5864b3e95d5SJeremy L Thompson 
5878b97b69aSJeremy L Thompson   // Setup output arrays
5884b3e95d5SJeremy L Thompson   code << "\n    // -- Output field setup\n";
5894b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
5904b3e95d5SJeremy L Thompson     std::string var_suffix = "_out_" + std::to_string(i);
5914b3e95d5SJeremy L Thompson 
5924b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
5934b3e95d5SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
5943a2968d6SJeremy L Thompson     switch (eval_mode) {
5953a2968d6SJeremy L Thompson       case CEED_EVAL_NONE:
5963a2968d6SJeremy L Thompson         if (is_at_points) {
5973a2968d6SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "];\n";
5983a2968d6SJeremy L Thompson         } else {
599*9123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
6004b3e95d5SJeremy L Thompson         }
6013a2968d6SJeremy L Thompson         break;
6023a2968d6SJeremy L Thompson       case CEED_EVAL_INTERP:
6033a2968d6SJeremy L Thompson         if (is_at_points) {
6043a2968d6SJeremy L Thompson           // Accumulator for point data
605*9123fb08SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
606*9123fb08SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "; i++) {\n";
6073a2968d6SJeremy L Thompson           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
6083a2968d6SJeremy L Thompson           code << "    }\n";
6093a2968d6SJeremy L Thompson         } else {
610*9123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
6113a2968d6SJeremy L Thompson         }
6123a2968d6SJeremy L Thompson         break;
6133a2968d6SJeremy L Thompson       case CEED_EVAL_GRAD:
6143a2968d6SJeremy L Thompson         if (is_at_points) {
6153a2968d6SJeremy L Thompson           // Accumulator for point data
616*9123fb08SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "*dim];\n";
617*9123fb08SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "; i++) {\n";
6183a2968d6SJeremy L Thompson           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
6193a2968d6SJeremy L Thompson           code << "    }\n";
6203a2968d6SJeremy L Thompson         } else if (use_3d_slices) {
6214b3e95d5SJeremy L Thompson           // Accumulator for gradient slices
6224b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
6234b3e95d5SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n";
6244b3e95d5SJeremy L Thompson           code << "      r_q" << var_suffix << "[i] = 0.0;\n";
6254b3e95d5SJeremy L Thompson           code << "    }\n";
6264b3e95d5SJeremy L Thompson         } else {
627*9123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
6284b3e95d5SJeremy L Thompson         }
6293a2968d6SJeremy L Thompson         break;
6303a2968d6SJeremy L Thompson       case CEED_EVAL_WEIGHT:
6313a2968d6SJeremy L Thompson         break;
6323a2968d6SJeremy L Thompson         // LCOV_EXCL_START
6333a2968d6SJeremy L Thompson       case CEED_EVAL_DIV:
6343a2968d6SJeremy L Thompson       case CEED_EVAL_CURL:
6353a2968d6SJeremy L Thompson         break;  // TODO: Not implemented
6363a2968d6SJeremy L Thompson                 // LCOV_EXCL_STOP
6374b3e95d5SJeremy L Thompson     }
6384b3e95d5SJeremy L Thompson   }
6394b3e95d5SJeremy L Thompson 
6403a2968d6SJeremy L Thompson   if (is_at_points) {
6413a2968d6SJeremy L Thompson     // We need to handle batches of points
6423a2968d6SJeremy L Thompson     code << "\n    // Note: Using batches of points\n";
6433a2968d6SJeremy L Thompson     code << "    const CeedInt point_loop_bound = (blockDim.x * blockDim.y) * ceil(1.0 * max_num_points / (blockDim.x * blockDim.y));\n\n";
6443a2968d6SJeremy L Thompson     code << "    #pragma unroll\n";
6453a2968d6SJeremy L Thompson     code << "    for (CeedInt i = threadIdx.x + threadIdx.y * blockDim.x; i < point_loop_bound; i += blockDim.x * blockDim.y) {\n";
6463a2968d6SJeremy L Thompson     code << "      const CeedInt p = i % max_num_points;\n\n";
6473a2968d6SJeremy L Thompson 
6483a2968d6SJeremy L Thompson     code << "      // -- Coordinates\n";
6493a2968d6SJeremy L Thompson     code << "      CeedScalar r_x[dim];\n";
6503a2968d6SJeremy 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";
6513a2968d6SJeremy L Thompson 
6523a2968d6SJeremy L Thompson     code << "      // -- Input fields\n";
6533a2968d6SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
6543a2968d6SJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
6553a2968d6SJeremy L Thompson 
6563a2968d6SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
6573a2968d6SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
6583a2968d6SJeremy L Thompson       // Basis action
6593a2968d6SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
6603a2968d6SJeremy L Thompson       switch (eval_mode) {
6613a2968d6SJeremy L Thompson         case CEED_EVAL_NONE:
6623a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
6633a2968d6SJeremy L Thompson           code << "      ReadPoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
6643a2968d6SJeremy L Thompson                << ", max_num_points>(data, elem, p, max_num_points, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
6653a2968d6SJeremy L Thompson           break;
6663a2968d6SJeremy L Thompson         case CEED_EVAL_INTERP:
6673a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
6683a2968d6SJeremy L Thompson           code << "      InterpAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix
6693a2968d6SJeremy L Thompson                << ", r_x, r_s" << var_suffix << ");\n";
6703a2968d6SJeremy L Thompson           break;
6713a2968d6SJeremy L Thompson         case CEED_EVAL_GRAD:
6723a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
6733a2968d6SJeremy L Thompson           code << "      GradAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix
6743a2968d6SJeremy L Thompson                << ", r_x, r_s" << var_suffix << ");\n";
6753a2968d6SJeremy L Thompson           break;
6763a2968d6SJeremy L Thompson         case CEED_EVAL_WEIGHT:
6773a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
6783a2968d6SJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = 1.0;\n";
6793a2968d6SJeremy L Thompson           break;
6803a2968d6SJeremy L Thompson           // LCOV_EXCL_START
6813a2968d6SJeremy L Thompson         case CEED_EVAL_DIV:
6823a2968d6SJeremy L Thompson         case CEED_EVAL_CURL:
6833a2968d6SJeremy L Thompson           break;  // TODO: Not implemented
6843a2968d6SJeremy L Thompson                   // LCOV_EXCL_STOP
6853a2968d6SJeremy L Thompson       }
6863a2968d6SJeremy L Thompson     }
6873a2968d6SJeremy L Thompson     code << "\n      // -- Output fields\n";
6883a2968d6SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
6893a2968d6SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
6903a2968d6SJeremy L Thompson 
6913a2968d6SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
6923a2968d6SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
6933a2968d6SJeremy L Thompson       // Basis action
6943a2968d6SJeremy L Thompson       switch (eval_mode) {
6953a2968d6SJeremy L Thompson         case CEED_EVAL_NONE:
6963a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
6973a2968d6SJeremy L Thompson           break;
6983a2968d6SJeremy L Thompson         case CEED_EVAL_INTERP:
6993a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7003a2968d6SJeremy L Thompson           break;
7013a2968d6SJeremy L Thompson         case CEED_EVAL_GRAD:
7023a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
7033a2968d6SJeremy L Thompson           break;
7043a2968d6SJeremy L Thompson           // LCOV_EXCL_START
7053a2968d6SJeremy L Thompson         case CEED_EVAL_WEIGHT:
7063a2968d6SJeremy L Thompson           break;  // Should not occur
7073a2968d6SJeremy L Thompson         case CEED_EVAL_DIV:
7083a2968d6SJeremy L Thompson         case CEED_EVAL_CURL:
7093a2968d6SJeremy L Thompson           break;  // TODO: Not implemented
7103a2968d6SJeremy L Thompson                   // LCOV_EXCL_STOP
7113a2968d6SJeremy L Thompson       }
7123a2968d6SJeremy L Thompson     }
7133a2968d6SJeremy L Thompson 
7143a2968d6SJeremy L Thompson   } else if (use_3d_slices) {
7154b3e95d5SJeremy L Thompson     // We treat quadrature points per slice in 3d to save registers
7164b3e95d5SJeremy L Thompson     code << "\n    // Note: Using planes of 3D elements\n";
7174b3e95d5SJeremy L Thompson     code << "    #pragma unroll\n";
7184b3e95d5SJeremy L Thompson     code << "    for (CeedInt q = 0; q < " << Q_name << "; q++) {\n";
7194b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
7204b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
7214b3e95d5SJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
7224b3e95d5SJeremy L Thompson 
7234b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
7244b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
7254b3e95d5SJeremy L Thompson       // Basis action
7264b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
7274b3e95d5SJeremy L Thompson       switch (eval_mode) {
7284b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
7294b3e95d5SJeremy L Thompson           bool is_strided;
7304b3e95d5SJeremy L Thompson 
7314b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7324b3e95d5SJeremy L Thompson 
7334b3e95d5SJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
7344b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
7354b3e95d5SJeremy L Thompson           if (is_strided) {
7364b3e95d5SJeremy L Thompson             bool    has_backend_strides;
7374b3e95d5SJeremy L Thompson             CeedInt num_elem, elem_size;
7384b3e95d5SJeremy L Thompson 
7394b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
7404b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
7414b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
7424b3e95d5SJeremy L Thompson             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
7434b3e95d5SJeremy L Thompson 
7444b3e95d5SJeremy L Thompson             if (!has_backend_strides) {
7454b3e95d5SJeremy L Thompson               CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
7464b3e95d5SJeremy L Thompson             }
7474b3e95d5SJeremy L Thompson             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
748f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStrided3d<num_comp" << var_suffix << ", " << Q_name << ", " << strides[0] << ", " << strides[1] << ", "
7494b3e95d5SJeremy L Thompson                  << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n";
7504b3e95d5SJeremy L Thompson           } else {
7514b3e95d5SJeremy L Thompson             CeedSize                 l_size = 0;
7524b3e95d5SJeremy L Thompson             CeedInt                  comp_stride;
7534b3e95d5SJeremy L Thompson             CeedElemRestriction_Hip *rstr_data;
7544b3e95d5SJeremy L Thompson 
7554b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
7564b3e95d5SJeremy L Thompson             code << "      const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
7574b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
7584b3e95d5SJeremy L Thompson             code << "      // CompStride: " << comp_stride << "\n";
7594b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
7604b3e95d5SJeremy L Thompson             data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
761f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStandard3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix
7624b3e95d5SJeremy L Thompson                  << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
7634b3e95d5SJeremy L Thompson           }
764*9123fb08SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
7654b3e95d5SJeremy L Thompson           break;
7664b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
7674b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7684b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n";
7694b3e95d5SJeremy L Thompson           code << "        r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n";
7704b3e95d5SJeremy L Thompson           code << "      }\n";
7714b3e95d5SJeremy L Thompson           break;
7724b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
7734b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
774f815fac9SJeremy L Thompson           code << "      GradColloSlice3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_q" << var_suffix << ", s_G" << var_suffix
775f815fac9SJeremy L Thompson                << ", r_s" << var_suffix << ");\n";
7764b3e95d5SJeremy L Thompson           break;
7774b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
7784b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
7794b3e95d5SJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n";
7803a2968d6SJeremy L Thompson           break;
7814b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
7824b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
7834b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
7844b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
7854b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
7864b3e95d5SJeremy L Thompson       }
7874b3e95d5SJeremy L Thompson     }
7884b3e95d5SJeremy L Thompson     code << "\n      // -- Output fields\n";
7894b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
7904b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
7914b3e95d5SJeremy L Thompson 
7924b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
7934b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
7944b3e95d5SJeremy L Thompson       // Basis action
7954b3e95d5SJeremy L Thompson       switch (eval_mode) {
7964b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
7974b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7983a2968d6SJeremy L Thompson           break;
7994b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
8004b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
8014b3e95d5SJeremy L Thompson           break;
8024b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
8034b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
8044b3e95d5SJeremy L Thompson           break;
8054b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
8064b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
8074b3e95d5SJeremy L Thompson           break;  // Should not occur
8084b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
8094b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
8104b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
8114b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
8124b3e95d5SJeremy L Thompson       }
8134b3e95d5SJeremy L Thompson     }
8144b3e95d5SJeremy L Thompson   } else {
8154b3e95d5SJeremy L Thompson     code << "\n    // Note: Using full elements\n";
8164b3e95d5SJeremy L Thompson     code << "    {\n";
8174b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
8184b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
8194b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
8204b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n";
8214b3e95d5SJeremy L Thompson     }
8224b3e95d5SJeremy L Thompson     code << "      // -- Output fields\n";
8234b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
8244b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
8254b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n";
8264b3e95d5SJeremy L Thompson     }
8274b3e95d5SJeremy L Thompson   }
8284b3e95d5SJeremy L Thompson 
8294b3e95d5SJeremy L Thompson   // Input and output buffers
8304b3e95d5SJeremy L Thompson   code << "\n      // -- QFunction inputs and outputs\n";
8314b3e95d5SJeremy L Thompson   code << "      // ---- Inputs\n";
8324b3e95d5SJeremy L Thompson   code << "      CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n";
8334b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
8344b3e95d5SJeremy L Thompson     code << "      // ------ Input field " << i << "\n";
8354b3e95d5SJeremy L Thompson     code << "      inputs[" << i << "] = r_s_in_" << i << ";\n";
8364b3e95d5SJeremy L Thompson   }
8374b3e95d5SJeremy L Thompson   code << "      // ---- Outputs\n";
8384b3e95d5SJeremy L Thompson   code << "      CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n";
8394b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
8404b3e95d5SJeremy L Thompson     code << "      // ------ Output field " << i << "\n";
8414b3e95d5SJeremy L Thompson     code << "      outputs[" << i << "] = r_s_out_" << i << ";\n";
8424b3e95d5SJeremy L Thompson   }
8434b3e95d5SJeremy L Thompson 
8444b3e95d5SJeremy L Thompson   // Apply QFunction
8454b3e95d5SJeremy L Thompson   code << "\n      // -- Apply QFunction\n";
8464b3e95d5SJeremy L Thompson   code << "      " << qfunction_name << "(ctx, ";
847*9123fb08SJeremy L Thompson   if (dim != 3 || is_at_points || use_3d_slices || !is_tensor) {
8484b3e95d5SJeremy L Thompson     code << "1";
8494b3e95d5SJeremy L Thompson   } else {
850*9123fb08SJeremy L Thompson     code << Q_name;
8514b3e95d5SJeremy L Thompson   }
8524b3e95d5SJeremy L Thompson   code << ", inputs, outputs);\n";
8534b3e95d5SJeremy L Thompson 
8543a2968d6SJeremy L Thompson   if (is_at_points) {
8553a2968d6SJeremy L Thompson     // Map back to coefficients
8563a2968d6SJeremy L Thompson     code << "\n      // -- Output fields\n";
8573a2968d6SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
8583a2968d6SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
8593a2968d6SJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
8603a2968d6SJeremy L Thompson 
8613a2968d6SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
8623a2968d6SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
8633a2968d6SJeremy L Thompson       // Basis action
8643a2968d6SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
8653a2968d6SJeremy L Thompson       switch (eval_mode) {
8663a2968d6SJeremy L Thompson         case CEED_EVAL_NONE: {
8673a2968d6SJeremy L Thompson           CeedInt             comp_stride;
8683a2968d6SJeremy L Thompson           CeedElemRestriction elem_rstr;
8693a2968d6SJeremy L Thompson 
8703a2968d6SJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
8713a2968d6SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
8723a2968d6SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
8733a2968d6SJeremy L Thompson           code << "      const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
8743a2968d6SJeremy L Thompson           code << "      WritePoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
8753a2968d6SJeremy L Thompson                << ", max_num_points>(data, elem, i, points.num_per_elem[elem], indices.outputs[" << i << "]"
8763a2968d6SJeremy L Thompson                << ", r_s" << var_suffix << ", d" << var_suffix << ");\n";
8773a2968d6SJeremy L Thompson           break;
8783a2968d6SJeremy L Thompson         }
8793a2968d6SJeremy L Thompson         case CEED_EVAL_INTERP:
8803a2968d6SJeremy L Thompson           code << "      if (i >= points.num_per_elem[elem]) {\n";
8813a2968d6SJeremy L Thompson           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) r_s" << var_suffix << "[j] = 0.0;\n";
8823a2968d6SJeremy L Thompson           code << "      }\n";
8833a2968d6SJeremy L Thompson           code << "      InterpTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s"
8843a2968d6SJeremy L Thompson                << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
8853a2968d6SJeremy L Thompson           break;
8863a2968d6SJeremy L Thompson         case CEED_EVAL_GRAD:
8873a2968d6SJeremy L Thompson           code << "      if (i >= points.num_per_elem[elem]) {\n";
8883a2968d6SJeremy L Thompson           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "*dim; j++) r_s" << var_suffix << "[j] = 0.0;\n";
8893a2968d6SJeremy L Thompson           code << "      }\n";
8903a2968d6SJeremy L Thompson           code << "      GradTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s"
8913a2968d6SJeremy L Thompson                << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
8923a2968d6SJeremy L Thompson           break;
8933a2968d6SJeremy L Thompson           // LCOV_EXCL_START
8943a2968d6SJeremy L Thompson         case CEED_EVAL_WEIGHT:
8953a2968d6SJeremy L Thompson           break;  // Should not occur
8963a2968d6SJeremy L Thompson         case CEED_EVAL_DIV:
8973a2968d6SJeremy L Thompson         case CEED_EVAL_CURL:
8983a2968d6SJeremy L Thompson           break;  // TODO: Not implemented
8993a2968d6SJeremy L Thompson                   // LCOV_EXCL_STOP
9003a2968d6SJeremy L Thompson       }
9013a2968d6SJeremy L Thompson     }
9023a2968d6SJeremy L Thompson   } else if (use_3d_slices) {
9034b3e95d5SJeremy L Thompson     // Copy or apply transpose grad, if needed
9043a2968d6SJeremy L Thompson     code << "\n      // -- Output fields\n";
9054b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
9064b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
9074b3e95d5SJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
9084b3e95d5SJeremy L Thompson 
9094b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
9104b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
9114b3e95d5SJeremy L Thompson       // Basis action
9124b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
9134b3e95d5SJeremy L Thompson       switch (eval_mode) {
9144b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
9154b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
9164b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
9174b3e95d5SJeremy L Thompson           code << "      }\n";
9183a2968d6SJeremy L Thompson           break;
9194b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
9204b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
9214b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
9224b3e95d5SJeremy L Thompson           code << "      }\n";
9234b3e95d5SJeremy L Thompson           break;
9244b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
925f815fac9SJeremy L Thompson           code << "      GradColloSliceTranspose3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_s" << var_suffix << ", s_G"
926f815fac9SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
9274b3e95d5SJeremy L Thompson           break;
9284b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
9294b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
9304b3e95d5SJeremy L Thompson           break;  // Should not occur
9314b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
9324b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
9334b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
9344b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
9354b3e95d5SJeremy L Thompson       }
9364b3e95d5SJeremy L Thompson     }
9374b3e95d5SJeremy L Thompson   }
9384b3e95d5SJeremy L Thompson   code << "    }\n";
9394b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9404b3e95d5SJeremy L Thompson }
9414b3e95d5SJeremy L Thompson 
9424b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
9439e201c85SYohann // Build single operator kernel
9447d8d0e25Snbeams //------------------------------------------------------------------------------
945eb7e6cafSJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Hip_gen(CeedOperator op) {
9463a2968d6SJeremy L Thompson   bool                   is_tensor = true, is_at_points = false, use_3d_slices = false;
9477d8d0e25Snbeams   Ceed                   ceed;
9483a2968d6SJeremy L Thompson   CeedInt                Q_1d, num_input_fields, num_output_fields, dim = 1, max_num_points = 0, coords_comp_stride = 0;
949b7453713SJeremy L Thompson   CeedQFunctionField    *qf_input_fields, *qf_output_fields;
950b7453713SJeremy L Thompson   CeedQFunction_Hip_gen *qf_data;
951b7453713SJeremy L Thompson   CeedQFunction          qf;
952b7453713SJeremy L Thompson   CeedOperatorField     *op_input_fields, *op_output_fields;
953b7453713SJeremy L Thompson   CeedOperator_Hip_gen  *data;
9544b3e95d5SJeremy L Thompson   std::ostringstream     code;
9554b3e95d5SJeremy L Thompson 
9564b3e95d5SJeremy L Thompson   {
9574b3e95d5SJeremy L Thompson     bool is_setup_done;
958b7453713SJeremy L Thompson 
959b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
960b7453713SJeremy L Thompson     if (is_setup_done) return CEED_ERROR_SUCCESS;
9614b3e95d5SJeremy L Thompson   }
962b7453713SJeremy L Thompson 
963b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
964b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
965b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
966b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
967b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
968b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
9697d8d0e25Snbeams 
9704b3e95d5SJeremy L Thompson   // Get operator data
9713a2968d6SJeremy L Thompson   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
9724b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelData_Hip_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, op_output_fields,
9734b3e95d5SJeremy L Thompson                                                       qf_output_fields, &data->max_P_1d, &Q_1d, &dim, &is_tensor, &use_3d_slices));
9744b3e95d5SJeremy L Thompson   if (dim == 0) dim = 1;
9754b3e95d5SJeremy L Thompson   data->dim = dim;
9763a2968d6SJeremy L Thompson   if (is_at_points) {
9773a2968d6SJeremy L Thompson     CeedElemRestriction_Hip *rstr_data;
9783a2968d6SJeremy L Thompson     CeedElemRestriction      rstr_points = NULL;
9794b3e95d5SJeremy L Thompson 
9803a2968d6SJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
9813a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points));
9823a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetCompStride(rstr_points, &coords_comp_stride));
9833a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(rstr_points, &rstr_data));
9843a2968d6SJeremy L Thompson     data->points.indices = (CeedInt *)rstr_data->d_offsets;
9853a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
9863a2968d6SJeremy L Thompson   }
9873a2968d6SJeremy L Thompson   if (is_at_points) use_3d_slices = false;
9883a2968d6SJeremy L Thompson   if (Q_1d == 0) {
9893a2968d6SJeremy L Thompson     if (is_at_points) Q_1d = max_num_points;
9903a2968d6SJeremy L Thompson     else CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q_1d));
9914b3e95d5SJeremy L Thompson   }
9924b3e95d5SJeremy L Thompson   data->Q_1d = Q_1d;
9934b3e95d5SJeremy L Thompson 
9940b454692Sjeremylt   // Check for restriction only identity operator
9954b3e95d5SJeremy L Thompson   {
9964b3e95d5SJeremy L Thompson     bool is_identity_qf;
9974b3e95d5SJeremy L Thompson 
9982b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
9990b454692Sjeremylt     if (is_identity_qf) {
10009e201c85SYohann       CeedEvalMode eval_mode_in, eval_mode_out;
1001b7453713SJeremy L Thompson 
10022b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
10032b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
10046574a04fSJeremy L Thompson       CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
10056574a04fSJeremy L Thompson                 "Backend does not implement restriction only identity operators");
10060b454692Sjeremylt     }
10074b3e95d5SJeremy L Thompson   }
1008b2165e7aSSebastian Grimberg 
1009b2165e7aSSebastian Grimberg   // Load basis source files
1010*9123fb08SJeremy L Thompson   if (is_tensor) {
10119c25dd66SJeremy L Thompson     code << "// Tensor basis source\n";
10129c25dd66SJeremy L Thompson     code << "#include <ceed/jit-source/hip/hip-shared-basis-tensor-templates.h>\n\n";
1013*9123fb08SJeremy L Thompson   } else {
1014*9123fb08SJeremy L Thompson     code << "// Non-tensor basis source\n";
1015*9123fb08SJeremy L Thompson     code << "#include <ceed/jit-source/hip/hip-shared-basis-nontensor-templates.h>\n\n";
1016*9123fb08SJeremy L Thompson   }
1017*9123fb08SJeremy L Thompson   if (is_at_points) {
10183a2968d6SJeremy L Thompson     code << "// AtPoints basis source\n";
10193a2968d6SJeremy L Thompson     code << "#include <ceed/jit-source/hip/hip-shared-basis-tensor-at-points-templates.h>\n\n";
1020*9123fb08SJeremy L Thompson   }
10219c25dd66SJeremy L Thompson   code << "// CodeGen operator source\n";
10229c25dd66SJeremy L Thompson   code << "#include <ceed/jit-source/hip/hip-gen-templates.h>\n\n";
10237d8d0e25Snbeams 
10244b3e95d5SJeremy L Thompson   // Get QFunction name
10254b3e95d5SJeremy L Thompson   std::string qfunction_name(qf_data->qfunction_name);
10264b3e95d5SJeremy L Thompson   std::string operator_name;
10274b3e95d5SJeremy L Thompson 
102809095acaSJeremy L Thompson   operator_name = "CeedKernelHipGenOperator_" + qfunction_name;
10297d8d0e25Snbeams 
10309e201c85SYohann   // Define CEED_Q_VLA
10319e201c85SYohann   code << "\n#undef CEED_Q_VLA\n";
1032*9123fb08SJeremy L Thompson   if (dim != 3 || is_at_points || use_3d_slices || !is_tensor) {
10339e201c85SYohann     code << "#define CEED_Q_VLA 1\n\n";
10349e201c85SYohann   } else {
10359e201c85SYohann     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
10369e201c85SYohann   }
10379e201c85SYohann 
10384b3e95d5SJeremy L Thompson   // Add user QFunction source
10394b3e95d5SJeremy L Thompson   {
10409c25dd66SJeremy L Thompson     const char *source_path;
10414b3e95d5SJeremy L Thompson 
10429c25dd66SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path));
10439c25dd66SJeremy L Thompson     CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/hip/gen backend requires QFunction source code file");
10449c25dd66SJeremy L Thompson 
10459c25dd66SJeremy L Thompson     code << "// User QFunction source\n";
10469c25dd66SJeremy L Thompson     code << "#include \"" << source_path << "\"\n\n";
10474b3e95d5SJeremy L Thompson   }
10487d8d0e25Snbeams 
10497d8d0e25Snbeams   // Setup
10507d8d0e25Snbeams   code << "\n// -----------------------------------------------------------------------------\n";
10514b3e95d5SJeremy L Thompson   code << "// Operator Kernel\n";
10524b3e95d5SJeremy L Thompson   code << "// \n";
10534b3e95d5SJeremy L Thompson   code << "// d_[in,out]_i:   CeedVector device array\n";
10544b3e95d5SJeremy L Thompson   code << "// r_[in,out]_e_i: Element vector register\n";
10554b3e95d5SJeremy L Thompson   code << "// r_[in,out]_q_i: Quadrature space vector register\n";
1056*9123fb08SJeremy L Thompson   code << "// r_[in,out]_c_i: AtPoints Chebyshev coefficients register\n";
10574b3e95d5SJeremy L Thompson   code << "// r_[in,out]_s_i: Quadrature space slice vector register\n";
10584b3e95d5SJeremy L Thompson   code << "// \n";
10594b3e95d5SJeremy L Thompson   code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n";
10604b3e95d5SJeremy L Thompson   code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n";
10614b3e95d5SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n";
1062b3e1519bSnbeams   code << "\nextern \"C\" __launch_bounds__(BLOCK_SIZE)\n";
10632b730f8bSJeremy L Thompson   code << "__global__ void " << operator_name
10643a2968d6SJeremy L Thompson        << "(CeedInt num_elem, void* ctx, FieldsInt_Hip indices, Fields_Hip fields, Fields_Hip B, Fields_Hip G, CeedScalar* W, Points_Hip points) {\n";
10654b3e95d5SJeremy L Thompson 
10664b3e95d5SJeremy L Thompson   // Scratch buffers
10679e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
10684b3e95d5SJeremy L Thompson     CeedEvalMode eval_mode;
10694b3e95d5SJeremy L Thompson 
10702b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
10719e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
10724b3e95d5SJeremy L Thompson       code << "  const CeedScalar *d_in_" << i << " = fields.inputs[" << i << "];\n";
10737d8d0e25Snbeams     }
10747d8d0e25Snbeams   }
10759e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
10764b3e95d5SJeremy L Thompson     code << "  CeedScalar *d_out_" << i << " = fields.outputs[" << i << "];\n";
10777d8d0e25Snbeams   }
10787d8d0e25Snbeams 
10799e201c85SYohann   code << "  const CeedInt dim = " << dim << ";\n";
1080*9123fb08SJeremy L Thompson   code << "  const CeedInt " << (is_tensor ? "Q_1d" : "Q") << " = " << Q_1d << ";\n";
10813a2968d6SJeremy L Thompson   if (is_at_points) {
10823a2968d6SJeremy L Thompson     code << "  const CeedInt max_num_points = " << max_num_points << ";\n";
10833a2968d6SJeremy L Thompson     code << "  const CeedInt coords_comp_stride = " << coords_comp_stride << ";\n";
10843a2968d6SJeremy L Thompson   }
10857d8d0e25Snbeams 
10864b3e95d5SJeremy L Thompson   // Shared data
10874b3e95d5SJeremy L Thompson   code << "  extern __shared__ CeedScalar slice[];\n";
10889e201c85SYohann   code << "  SharedData_Hip data;\n";
10899e201c85SYohann   code << "  data.t_id_x = threadIdx.x;\n";
10909e201c85SYohann   code << "  data.t_id_y = threadIdx.y;\n";
10919e201c85SYohann   code << "  data.t_id_z = threadIdx.z;\n";
10929e201c85SYohann   code << "  data.t_id  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
1093*9123fb08SJeremy L Thompson   code << "  data.slice = slice + data.t_id_z*T_1D" << ((!is_tensor || dim == 1) ? "" : "*T_1D") << ";\n";
10947d8d0e25Snbeams 
10957d8d0e25Snbeams   // Initialize constants, and matrices B and G
10964b3e95d5SJeremy L Thompson   code << "\n  // Input field constants and basis data\n";
10979e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1098*9123fb08SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_input_fields[i], qf_input_fields[i], Q_1d, true, is_tensor,
1099*9123fb08SJeremy L Thompson                                                              is_at_points, use_3d_slices));
11007d8d0e25Snbeams   }
11014b3e95d5SJeremy L Thompson   code << "\n  // Output field constants and basis data\n";
11029e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1103*9123fb08SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_output_fields[i], qf_output_fields[i], Q_1d, false, is_tensor,
1104*9123fb08SJeremy L Thompson                                                              is_at_points, use_3d_slices));
11054b3e95d5SJeremy L Thompson   }
11067d8d0e25Snbeams 
11074b3e95d5SJeremy L Thompson   // Loop over all elements
11084b3e95d5SJeremy L Thompson   code << "\n  // Element loop\n";
11097d8d0e25Snbeams   code << "  __syncthreads();\n";
11109e201c85SYohann   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
11114b3e95d5SJeremy L Thompson 
1112e93651e5SJeremy L Thompson   // -- Compute minimum buffer space needed
11133a2968d6SJeremy L Thompson   CeedInt max_rstr_buffer_size = 1;
1114e93651e5SJeremy L Thompson 
1115e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1116e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
1117e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
1118e93651e5SJeremy L Thompson 
1119e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
1120e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1121e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1122*9123fb08SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_tensor && (dim >= 3) ? elem_size : 1));
1123681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1124e93651e5SJeremy L Thompson   }
1125e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1126e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
1127e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
1128e93651e5SJeremy L Thompson 
1129e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
1130e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1131e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1132*9123fb08SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_tensor && (dim >= 3) ? elem_size : 1));
1133681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1134e93651e5SJeremy L Thompson   }
1135e93651e5SJeremy L Thompson   code << "    // Scratch restriction buffer space\n";
1136e93651e5SJeremy L Thompson   code << "    CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n";
1137e93651e5SJeremy L Thompson 
1138e93651e5SJeremy L Thompson   // -- Determine best input field processing order
1139e93651e5SJeremy L Thompson   CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX];
1140e93651e5SJeremy L Thompson 
1141e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1142e93651e5SJeremy L Thompson     field_rstr_in_buffer[i] = -1;
1143e93651e5SJeremy L Thompson     input_field_order[i]    = -1;
1144e93651e5SJeremy L Thompson   }
1145e93651e5SJeremy L Thompson   {
1146e93651e5SJeremy L Thompson     bool    is_ordered[CEED_FIELD_MAX];
1147e93651e5SJeremy L Thompson     CeedInt curr_index = 0;
1148e93651e5SJeremy L Thompson 
1149e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
1150e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
1151e93651e5SJeremy L Thompson       CeedVector          vec_i;
1152e93651e5SJeremy L Thompson       CeedElemRestriction rstr_i;
1153e93651e5SJeremy L Thompson 
1154e93651e5SJeremy L Thompson       if (is_ordered[i]) continue;
1155e93651e5SJeremy L Thompson       field_rstr_in_buffer[i]       = i;
1156e93651e5SJeremy L Thompson       is_ordered[i]                 = true;
1157e93651e5SJeremy L Thompson       input_field_order[curr_index] = i;
1158e93651e5SJeremy L Thompson       curr_index++;
1159034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
1160e93651e5SJeremy L Thompson       if (vec_i == CEED_VECTOR_NONE) continue;  // CEED_EVAL_WEIGHT
1161e93651e5SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
1162e93651e5SJeremy L Thompson       for (CeedInt j = i + 1; j < num_input_fields; j++) {
1163e93651e5SJeremy L Thompson         CeedVector          vec_j;
1164e93651e5SJeremy L Thompson         CeedElemRestriction rstr_j;
1165e93651e5SJeremy L Thompson 
1166e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
1167e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
1168e93651e5SJeremy L Thompson         if (rstr_i == rstr_j && vec_i == vec_j) {
1169e93651e5SJeremy L Thompson           field_rstr_in_buffer[j]       = i;
1170e93651e5SJeremy L Thompson           is_ordered[j]                 = true;
1171e93651e5SJeremy L Thompson           input_field_order[curr_index] = j;
1172e93651e5SJeremy L Thompson           curr_index++;
1173e93651e5SJeremy L Thompson         }
11743a2968d6SJeremy L Thompson         CeedCallBackend(CeedVectorDestroy(&vec_j));
11753a2968d6SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
1176e93651e5SJeremy L Thompson       }
11773a2968d6SJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec_i));
11783a2968d6SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
1179e93651e5SJeremy L Thompson     }
1180e93651e5SJeremy L Thompson   }
1181e93651e5SJeremy L Thompson 
11824b3e95d5SJeremy L Thompson   // -- Input restriction and basis
11833a2968d6SJeremy L Thompson   code << "\n    // -- Input field restrictions and basis actions\n";
11849e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1185e93651e5SJeremy L Thompson     CeedInt f = input_field_order[i];
1186e93651e5SJeremy L Thompson 
1187e93651e5SJeremy L Thompson     code << "    // ---- Input field " << f << "\n";
11887d8d0e25Snbeams 
11894b3e95d5SJeremy L Thompson     // ---- Restriction
1190e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Hip_gen(code, data, f, dim, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f], Q_1d,
1191*9123fb08SJeremy L Thompson                                                                true, is_tensor, is_at_points, use_3d_slices));
1192b7453713SJeremy L Thompson 
11934b3e95d5SJeremy L Thompson     // ---- Basis action
1194*9123fb08SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, f, dim, op_input_fields[f], qf_input_fields[f], Q_1d, true, is_tensor,
1195*9123fb08SJeremy L Thompson                                                          is_at_points, use_3d_slices));
11967d8d0e25Snbeams   }
11977d8d0e25Snbeams 
11984b3e95d5SJeremy L Thompson   // -- Q function
11993a2968d6SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelQFunction_Hip_gen(code, data, dim, max_num_points, num_input_fields, op_input_fields, qf_input_fields,
1200*9123fb08SJeremy L Thompson                                                            num_output_fields, op_output_fields, qf_output_fields, qfunction_name, Q_1d, is_tensor,
1201*9123fb08SJeremy L Thompson                                                            is_at_points, use_3d_slices));
12027d8d0e25Snbeams 
12034b3e95d5SJeremy L Thompson   // -- Output basis and restriction
12044b3e95d5SJeremy L Thompson   code << "\n    // -- Output field basis action and restrictions\n";
12059e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
12064b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
1207b7453713SJeremy L Thompson 
12084b3e95d5SJeremy L Thompson     // ---- Basis action
1209*9123fb08SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, i, dim, op_output_fields[i], qf_output_fields[i], Q_1d, false, is_tensor,
1210*9123fb08SJeremy L Thompson                                                          is_at_points, use_3d_slices));
12117d8d0e25Snbeams 
12124b3e95d5SJeremy L Thompson     // ---- Restriction
12133a2968d6SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Hip_gen(code, data, i, dim, NULL, op_output_fields[i], qf_output_fields[i], Q_1d, false,
1214*9123fb08SJeremy L Thompson                                                                is_tensor, is_at_points, use_3d_slices));
12157d8d0e25Snbeams   }
12167d8d0e25Snbeams 
12174b3e95d5SJeremy L Thompson   // Close loop and function
12187d8d0e25Snbeams   code << "  }\n";
12197d8d0e25Snbeams   code << "}\n";
12207d8d0e25Snbeams   code << "// -----------------------------------------------------------------------------\n\n";
12217d8d0e25Snbeams 
1222539ec17dSJeremy L Thompson   CeedInt block_sizes[3] = {0, 0, 0};
12239e201c85SYohann   CeedInt num_elem;
1224b7453713SJeremy L Thompson 
12253a2968d6SJeremy L Thompson   // Compile
12262b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
1227*9123fb08SJeremy L Thompson   CeedCallBackend(BlockGridCalculate_Hip_gen(is_tensor ? dim : 1, num_elem, data->max_P_1d, Q_1d, block_sizes));
1228eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Hip(ceed, code.str().c_str(), &data->module, 2, "T_1D", block_sizes[0], "BLOCK_SIZE",
12292b730f8bSJeremy L Thompson                                   block_sizes[0] * block_sizes[1] * block_sizes[2]));
1230eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, operator_name.c_str(), &data->op));
12312b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
12329bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
1233c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionDestroy(&qf));
1234e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12357d8d0e25Snbeams }
12362a86cc9dSSebastian Grimberg 
12377d8d0e25Snbeams //------------------------------------------------------------------------------
1238