xref: /libCEED/rust/libceed-sys/c-src/backends/hip-gen/ceed-hip-gen-operator-build.cpp (revision 9ee499e57bead864c0606ad2f4caeb6962c54f83)
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;
649123fb08SJeremy 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;
769123fb08SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
779123fb08SJeremy 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;
829123fb08SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
839123fb08SJeremy 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;
1009123fb08SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
1019123fb08SJeremy 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;
1069123fb08SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
1079123fb08SJeremy 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*9ee499e5SJeremy L Thompson                                                     CeedQFunctionField qf_field, CeedInt field_reuse[3], CeedInt Q_1d, bool is_input, bool is_tensor,
158*9ee499e5SJeremy L Thompson                                                     bool is_at_points, bool use_3d_slices) {
1594b3e95d5SJeremy L Thompson   std::string           var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
1609123fb08SJeremy 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 
168*9ee499e5SJeremy L Thompson   // Field reuse info
169*9ee499e5SJeremy L Thompson   bool         use_previous_field = field_reuse[0] != -1;
170*9ee499e5SJeremy L Thompson   bool         reuse_input        = field_reuse[1];
171*9ee499e5SJeremy L Thompson   CeedInt      reuse_field        = field_reuse[0];
172*9ee499e5SJeremy L Thompson   CeedEvalMode reuse_mode         = (CeedEvalMode)field_reuse[2];
173*9ee499e5SJeremy L Thompson 
1744b3e95d5SJeremy L Thompson   code << "  // -- " << (is_input ? "Input" : "Output") << " field " << i << "\n";
1754b3e95d5SJeremy L Thompson 
1764b3e95d5SJeremy L Thompson   // Get field data
1774b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
1784b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
1794b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1804b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1814b3e95d5SJeremy L Thompson   }
1823a2968d6SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1834b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
1844b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
1854b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1869123fb08SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
1879123fb08SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
1884b3e95d5SJeremy L Thompson   }
1894b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
1904b3e95d5SJeremy L Thompson 
1914b3e95d5SJeremy L Thompson   // Set field constants
1924b3e95d5SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
1934b3e95d5SJeremy L Thompson     code << "  const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n";
1944b3e95d5SJeremy L Thompson     code << "  const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n";
1954b3e95d5SJeremy L Thompson   }
1964b3e95d5SJeremy L Thompson 
1974b3e95d5SJeremy L Thompson   // Load basis data
1984b3e95d5SJeremy L Thompson   code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
1994b3e95d5SJeremy L Thompson   switch (eval_mode) {
2004b3e95d5SJeremy L Thompson     case CEED_EVAL_NONE:
2014b3e95d5SJeremy L Thompson       break;
2024b3e95d5SJeremy L Thompson     case CEED_EVAL_INTERP:
2033a2968d6SJeremy L Thompson       if (is_at_points) {
2043a2968d6SJeremy L Thompson         // AtPoints
2053a2968d6SJeremy L Thompson         if (!basis_data->d_chebyshev_interp_1d) {
2063a2968d6SJeremy L Thompson           CeedSize    interp_bytes;
2073a2968d6SJeremy L Thompson           CeedScalar *chebyshev_interp_1d;
2083a2968d6SJeremy L Thompson 
2093a2968d6SJeremy L Thompson           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
2103a2968d6SJeremy L Thompson           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
2113a2968d6SJeremy L Thompson           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
2123a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis), hipMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
2133a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis),
2143a2968d6SJeremy L Thompson                       hipMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, hipMemcpyHostToDevice));
2153a2968d6SJeremy L Thompson           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
2163a2968d6SJeremy L Thompson         }
2173a2968d6SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
2183a2968d6SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
2193a2968d6SJeremy L Thompson       } else {
2203a2968d6SJeremy L Thompson         // Standard quadrature
2214b3e95d5SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
2224b3e95d5SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_interp_1d;
2233a2968d6SJeremy L Thompson       }
224*9ee499e5SJeremy L Thompson       if (use_previous_field) {
225*9ee499e5SJeremy L Thompson         std::string reuse_var = "s_B" + ((reuse_input ? "_in_" : "_out_") + std::to_string(reuse_field));
226*9ee499e5SJeremy L Thompson 
227*9ee499e5SJeremy L Thompson         code << "  CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n";
228*9ee499e5SJeremy L Thompson       } else {
2299123fb08SJeremy L Thompson         code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
230f815fac9SJeremy L Thompson         code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
231*9ee499e5SJeremy L Thompson       }
2324b3e95d5SJeremy L Thompson       break;
2334b3e95d5SJeremy L Thompson     case CEED_EVAL_GRAD:
2343a2968d6SJeremy L Thompson       if (is_at_points) {
2353a2968d6SJeremy L Thompson         // AtPoints
2363a2968d6SJeremy L Thompson         if (!basis_data->d_chebyshev_interp_1d) {
2373a2968d6SJeremy L Thompson           CeedSize    interp_bytes;
2383a2968d6SJeremy L Thompson           CeedScalar *chebyshev_interp_1d;
2393a2968d6SJeremy L Thompson 
2403a2968d6SJeremy L Thompson           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
2413a2968d6SJeremy L Thompson           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
2423a2968d6SJeremy L Thompson           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
2433a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis), hipMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
2443a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis),
2453a2968d6SJeremy L Thompson                       hipMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, hipMemcpyHostToDevice));
2463a2968d6SJeremy L Thompson           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
2473a2968d6SJeremy L Thompson         }
2483a2968d6SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
2493a2968d6SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
2503a2968d6SJeremy L Thompson       } else {
2513a2968d6SJeremy L Thompson         // Standard quadrature
2524b3e95d5SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
2534b3e95d5SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_interp_1d;
2543a2968d6SJeremy L Thompson       }
2559123fb08SJeremy L Thompson       if (is_tensor) {
256*9ee499e5SJeremy L Thompson         if (use_previous_field) {
257*9ee499e5SJeremy L Thompson           std::string reuse_var = "s_B" + ((reuse_input ? "_in_" : "_out_") + std::to_string(reuse_field));
258*9ee499e5SJeremy L Thompson 
259*9ee499e5SJeremy L Thompson           code << "  CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n";
260*9ee499e5SJeremy L Thompson         } else {
2619123fb08SJeremy L Thompson           code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
262f815fac9SJeremy L Thompson           code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
2639123fb08SJeremy L Thompson         }
264*9ee499e5SJeremy L Thompson       }
2653a2968d6SJeremy L Thompson       if (is_at_points) break;  // No G mat for AtPoints
2664b3e95d5SJeremy L Thompson       if (use_3d_slices) {
2674b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d;
2684b3e95d5SJeremy L Thompson         else data->G.outputs[i] = basis_data->d_collo_grad_1d;
269*9ee499e5SJeremy L Thompson         if (use_previous_field && reuse_mode == CEED_EVAL_GRAD) {
270*9ee499e5SJeremy L Thompson           std::string reuse_var = "s_G" + ((reuse_input ? "_in_" : "_out_") + std::to_string(reuse_field));
271*9ee499e5SJeremy L Thompson 
272*9ee499e5SJeremy L Thompson           code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
273*9ee499e5SJeremy L Thompson         } else {
2749123fb08SJeremy L Thompson           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
275f815fac9SJeremy L Thompson           code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
276*9ee499e5SJeremy L Thompson         }
2774b3e95d5SJeremy L Thompson       } else {
2784b3e95d5SJeremy L Thompson         bool has_collo_grad = basis_data->d_collo_grad_1d;
2794b3e95d5SJeremy L Thompson 
2804b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2814b3e95d5SJeremy L Thompson         else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2824b3e95d5SJeremy L Thompson         if (has_collo_grad) {
283*9ee499e5SJeremy L Thompson           if (use_previous_field && reuse_mode == CEED_EVAL_GRAD) {
284*9ee499e5SJeremy L Thompson             std::string reuse_var = "s_G" + ((reuse_input ? "_in_" : "_out_") + std::to_string(reuse_field));
285*9ee499e5SJeremy L Thompson 
286*9ee499e5SJeremy L Thompson             code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
287*9ee499e5SJeremy L Thompson           } else {
2889123fb08SJeremy L Thompson             code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
289f815fac9SJeremy L Thompson             code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
290*9ee499e5SJeremy L Thompson           }
291*9ee499e5SJeremy L Thompson         } else {
292*9ee499e5SJeremy L Thompson           if (use_previous_field && reuse_mode == CEED_EVAL_GRAD) {
293*9ee499e5SJeremy L Thompson             std::string reuse_var = "s_G" + ((reuse_input ? "_in_" : "_out_") + std::to_string(reuse_field));
294*9ee499e5SJeremy L Thompson 
295*9ee499e5SJeremy L Thompson             code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
2964b3e95d5SJeremy L Thompson           } else {
2979123fb08SJeremy L Thompson             code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << P_name << "*" << Q_name << (is_tensor ? "" : "*dim") << "];\n";
2989123fb08SJeremy L Thompson             code << "  LoadMatrix<" << P_name << ", " << Q_name << (is_tensor ? "" : "*dim") << ">(data, G." << option_name << "[" << i << "], s_G"
2999123fb08SJeremy L Thompson                  << var_suffix << ");\n";
3004b3e95d5SJeremy L Thompson           }
3014b3e95d5SJeremy L Thompson         }
302*9ee499e5SJeremy L Thompson       }
3034b3e95d5SJeremy L Thompson       break;
3044b3e95d5SJeremy L Thompson     case CEED_EVAL_WEIGHT:
3054b3e95d5SJeremy L Thompson       break;  // No action
3064b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
3074b3e95d5SJeremy L Thompson     case CEED_EVAL_DIV:
3084b3e95d5SJeremy L Thompson     case CEED_EVAL_CURL:
3094b3e95d5SJeremy L Thompson       break;  // TODO: Not implemented
3104b3e95d5SJeremy L Thompson               // LCOV_EXCL_STOP
3114b3e95d5SJeremy L Thompson   }
3123a2968d6SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
3134b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3144b3e95d5SJeremy L Thompson }
3154b3e95d5SJeremy L Thompson 
3164b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
3174b3e95d5SJeremy L Thompson // Restriction
3184b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
3194b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelRestriction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim,
320e93651e5SJeremy L Thompson                                                       CeedInt field_input_buffer[], CeedOperatorField op_field, CeedQFunctionField qf_field,
3219123fb08SJeremy L Thompson                                                       CeedInt Q_1d, bool is_input, bool is_tensor, bool is_at_points, bool use_3d_slices) {
3224b3e95d5SJeremy L Thompson   std::string              var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
3239123fb08SJeremy L Thompson   std::string              P_name     = (is_tensor ? "P_1d" : "P") + var_suffix;
3244b3e95d5SJeremy L Thompson   CeedEvalMode             eval_mode  = CEED_EVAL_NONE;
3254b3e95d5SJeremy L Thompson   CeedInt                  elem_size = 0, num_comp = 0, P_1d = 0;
3264b3e95d5SJeremy L Thompson   CeedSize                 l_size;
327f815fac9SJeremy L Thompson   CeedRestrictionType      rstr_type = CEED_RESTRICTION_STANDARD;
3284b3e95d5SJeremy L Thompson   CeedElemRestriction_Hip *rstr_data;
3294b3e95d5SJeremy L Thompson   CeedElemRestriction      elem_rstr;
3304b3e95d5SJeremy L Thompson   CeedBasis                basis;
3314b3e95d5SJeremy L Thompson 
3324b3e95d5SJeremy L Thompson   // Get field data
3334b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
3344b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
335f815fac9SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
3364b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
3374b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
3384b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
3394b3e95d5SJeremy L Thompson   }
3404b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
3414b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
3429123fb08SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
3439123fb08SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
3444b3e95d5SJeremy L Thompson   }
3453a2968d6SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
3464b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
3474b3e95d5SJeremy L Thompson 
3484b3e95d5SJeremy L Thompson   // Restriction
3494b3e95d5SJeremy L Thompson   if (is_input) {
3504b3e95d5SJeremy L Thompson     // Input
351e93651e5SJeremy L Thompson     if (field_input_buffer[i] != i) {
352e93651e5SJeremy L Thompson       std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]);
353e93651e5SJeremy L Thompson 
354e93651e5SJeremy L Thompson       // Restriction was already done for previous input
355e93651e5SJeremy L Thompson       code << "    CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n";
3563a2968d6SJeremy L Thompson     } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices && is_at_points)) {
3573a2968d6SJeremy L Thompson       if (eval_mode == CEED_EVAL_NONE && rstr_type != CEED_RESTRICTION_POINTS) {
358e93651e5SJeremy L Thompson         // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated
3594b3e95d5SJeremy L Thompson         code << "    CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n";
3603a2968d6SJeremy L Thompson       } else if (rstr_type != CEED_RESTRICTION_POINTS) {
361e93651e5SJeremy L Thompson         // Otherwise we're using the scratch space
362e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
363e93651e5SJeremy L Thompson       }
364f815fac9SJeremy L Thompson       switch (rstr_type) {
365f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STANDARD: {
3664b3e95d5SJeremy L Thompson           CeedInt comp_stride;
3674b3e95d5SJeremy L Thompson 
3684b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
3694b3e95d5SJeremy L Thompson           code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
3704b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3714b3e95d5SJeremy L Thompson           code << "    // CompStride: " << comp_stride << "\n";
3724b3e95d5SJeremy L Thompson           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
3739123fb08SJeremy L Thompson           code << "    ReadLVecStandard" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name
3749123fb08SJeremy L Thompson                << ">(data, l_size" << var_suffix << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n";
375f815fac9SJeremy L Thompson           break;
376f815fac9SJeremy L Thompson         }
377f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STRIDED: {
3784b3e95d5SJeremy L Thompson           bool    has_backend_strides;
3794b3e95d5SJeremy L Thompson           CeedInt num_elem;
3804b3e95d5SJeremy L Thompson 
3814b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
3824b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
3834b3e95d5SJeremy L Thompson           CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
3844b3e95d5SJeremy L Thompson 
3854b3e95d5SJeremy L Thompson           if (!has_backend_strides) {
3864b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
3874b3e95d5SJeremy L Thompson           }
3884b3e95d5SJeremy L Thompson           code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
3899123fb08SJeremy L Thompson           code << "    ReadLVecStrided" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", "
3909123fb08SJeremy L Thompson                << strides[1] << ", " << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n";
391f815fac9SJeremy L Thompson           break;
392f815fac9SJeremy L Thompson         }
3933a2968d6SJeremy L Thompson         case CEED_RESTRICTION_POINTS: {
3943a2968d6SJeremy L Thompson           CeedInt comp_stride;
3953a2968d6SJeremy L Thompson 
3963a2968d6SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3973a2968d6SJeremy L Thompson           code << "    const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
3983a2968d6SJeremy L Thompson           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
3993a2968d6SJeremy L Thompson           break;
4003a2968d6SJeremy L Thompson         }
401f815fac9SJeremy L Thompson         // LCOV_EXCL_START
402f815fac9SJeremy L Thompson         case CEED_RESTRICTION_ORIENTED:
403f815fac9SJeremy L Thompson         case CEED_RESTRICTION_CURL_ORIENTED:
404f815fac9SJeremy L Thompson           break;  // TODO: Not implemented
405f815fac9SJeremy L Thompson                   // LCOV_EXCL_STOP
4064b3e95d5SJeremy L Thompson       }
4074b3e95d5SJeremy L Thompson     }
4084b3e95d5SJeremy L Thompson   } else {
4094b3e95d5SJeremy L Thompson     // Output
410f815fac9SJeremy L Thompson     switch (rstr_type) {
411f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STANDARD: {
4124b3e95d5SJeremy L Thompson         CeedInt comp_stride;
4134b3e95d5SJeremy L Thompson 
4144b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
4154b3e95d5SJeremy L Thompson         code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
4164b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
4174b3e95d5SJeremy L Thompson         code << "    // CompStride: " << comp_stride << "\n";
4184b3e95d5SJeremy L Thompson         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
4199123fb08SJeremy L Thompson         code << "    WriteLVecStandard" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name
4209123fb08SJeremy L Thompson              << ">(data, l_size" << var_suffix << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n";
421f815fac9SJeremy L Thompson         break;
422f815fac9SJeremy L Thompson       }
423f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STRIDED: {
4244b3e95d5SJeremy L Thompson         bool    has_backend_strides;
4254b3e95d5SJeremy L Thompson         CeedInt num_elem;
4264b3e95d5SJeremy L Thompson 
4274b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
4284b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
4294b3e95d5SJeremy L Thompson         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
4304b3e95d5SJeremy L Thompson 
4314b3e95d5SJeremy L Thompson         if (!has_backend_strides) {
4324b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
4334b3e95d5SJeremy L Thompson         }
4344b3e95d5SJeremy L Thompson         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
4359123fb08SJeremy L Thompson         code << "    WriteLVecStrided" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", "
4369123fb08SJeremy L Thompson              << strides[1] << ", " << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n";
437f815fac9SJeremy L Thompson         break;
438f815fac9SJeremy L Thompson       }
4393a2968d6SJeremy L Thompson       case CEED_RESTRICTION_POINTS:
4403a2968d6SJeremy L Thompson         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
4413a2968d6SJeremy L Thompson         break;
442f815fac9SJeremy L Thompson       // LCOV_EXCL_START
443f815fac9SJeremy L Thompson       case CEED_RESTRICTION_ORIENTED:
444f815fac9SJeremy L Thompson       case CEED_RESTRICTION_CURL_ORIENTED:
445f815fac9SJeremy L Thompson         break;  // TODO: Not implemented
446f815fac9SJeremy L Thompson                 // LCOV_EXCL_STOP
4474b3e95d5SJeremy L Thompson     }
4484b3e95d5SJeremy L Thompson   }
4493a2968d6SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
4504b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4514b3e95d5SJeremy L Thompson }
4524b3e95d5SJeremy L Thompson 
4534b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4544b3e95d5SJeremy L Thompson // Basis
4554b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4564b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelBasis_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim,
4579123fb08SJeremy L Thompson                                                 CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool is_tensor,
4583a2968d6SJeremy L Thompson                                                 bool is_at_points, bool use_3d_slices) {
4594b3e95d5SJeremy L Thompson   std::string         var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
4609123fb08SJeremy L Thompson   std::string         P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q";
4614b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
4624b3e95d5SJeremy L Thompson   CeedInt             elem_size = 0, num_comp = 0, P_1d = 0;
4634b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
4644b3e95d5SJeremy L Thompson   CeedBasis           basis;
4654b3e95d5SJeremy L Thompson 
4664b3e95d5SJeremy L Thompson   // Get field data
4674b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
4684b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
4694b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
4704b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
4714b3e95d5SJeremy L Thompson   }
4723a2968d6SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
4734b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
4744b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
4759123fb08SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
4769123fb08SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
4774b3e95d5SJeremy L Thompson   }
4784b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
4794b3e95d5SJeremy L Thompson 
4804b3e95d5SJeremy L Thompson   // Basis
4814b3e95d5SJeremy L Thompson   code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
4824b3e95d5SJeremy L Thompson   if (is_input) {
4834b3e95d5SJeremy L Thompson     switch (eval_mode) {
4844b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
4853a2968d6SJeremy L Thompson         if (!use_3d_slices && !is_at_points) {
4864b3e95d5SJeremy L Thompson           code << "    CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n";
4874b3e95d5SJeremy L Thompson         }
4884b3e95d5SJeremy L Thompson         break;
4894b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
4903a2968d6SJeremy L Thompson         if (is_at_points) {
4919123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
4929123fb08SJeremy L Thompson 
4939123fb08SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
4949123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
4959123fb08SJeremy L Thompson                << var_suffix << ", r_c" << var_suffix << ");\n";
4963a2968d6SJeremy L Thompson         } else {
4979123fb08SJeremy L Thompson           std::string function_name = is_tensor ? ((dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d") : "InterpNonTensor";
4989123fb08SJeremy L Thompson 
4999123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
5009123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
5019123fb08SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
5023a2968d6SJeremy L Thompson         }
5034b3e95d5SJeremy L Thompson         break;
5044b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
5053a2968d6SJeremy L Thompson         if (is_at_points) {
5069123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
5079123fb08SJeremy L Thompson 
5089123fb08SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
5099123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
5109123fb08SJeremy L Thompson                << var_suffix << ", r_c" << var_suffix << ");\n";
5113a2968d6SJeremy L Thompson         } else if (use_3d_slices) {
5129123fb08SJeremy L Thompson           std::string function_name = (dim > 1 ? "InterpTensor" : "Interp") + std::to_string(dim) + "d";
5139123fb08SJeremy L Thompson 
5144b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
5159123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
5169123fb08SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
5179123fb08SJeremy L Thompson         } else if (is_tensor) {
5189123fb08SJeremy L Thompson           bool        is_collocated = dim == 3 && Q_1d >= P_1d;
5199123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "Grad" : (is_collocated ? "GradTensorCollocated" : "GradTensor")) + std::to_string(dim) + "d";
5209123fb08SJeremy L Thompson 
5219123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << (dim >= 3 ? Q_name : "1") << "];\n";
5229123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
5239123fb08SJeremy L Thompson                << var_suffix << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
5244b3e95d5SJeremy L Thompson         } else {
5259123fb08SJeremy L Thompson           std::string function_name = "GradNonTensor";
5269123fb08SJeremy L Thompson 
5279123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
5289123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", dim, " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix
5299123fb08SJeremy L Thompson                << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
5304b3e95d5SJeremy L Thompson         }
5314b3e95d5SJeremy L Thompson         break;
5324b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
5333a2968d6SJeremy L Thompson         if (is_at_points) {
5343a2968d6SJeremy L Thompson           code << "    // Nothing to do AtPoints\n";
5353a2968d6SJeremy L Thompson         } else {
5364b3e95d5SJeremy L Thompson           CeedBasis_Hip_shared *basis_data;
5379123fb08SJeremy L Thompson           std::string           function_name = is_tensor ? ((dim == 1 ? "Weight" : "WeightTensor") + std::to_string(dim) + "d") : "WeightNonTensor";
5384b3e95d5SJeremy L Thompson 
5399123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
5404b3e95d5SJeremy L Thompson           CeedCallBackend(CeedBasisGetData(basis, &basis_data));
5414b3e95d5SJeremy L Thompson           data->W = basis_data->d_q_weight_1d;
5429123fb08SJeremy L Thompson           code << "    " << function_name << "<" << Q_name << ">(data, W, r_q" << var_suffix << ");\n";
5433a2968d6SJeremy L Thompson         }
5444b3e95d5SJeremy L Thompson         break;
5454b3e95d5SJeremy L Thompson       }
5464b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
5474b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
5484b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
5494b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
5504b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
5514b3e95d5SJeremy L Thompson     }
5524b3e95d5SJeremy L Thompson   } else {
5534b3e95d5SJeremy L Thompson     switch (eval_mode) {
5544b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
5554b3e95d5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n";
5564b3e95d5SJeremy L Thompson         break;  // No action
5574b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
558e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
5593a2968d6SJeremy L Thompson         if (is_at_points) {
5609123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
5619123fb08SJeremy L Thompson 
5629123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_c" << var_suffix << ", s_B"
5639123fb08SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5643a2968d6SJeremy L Thompson         } else {
5659123fb08SJeremy L Thompson           std::string function_name =
5669123fb08SJeremy L Thompson               is_tensor ? ((dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d") : "InterpTransposeNonTensor";
5679123fb08SJeremy L Thompson 
5689123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B"
5699123fb08SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5703a2968d6SJeremy L Thompson         }
5714b3e95d5SJeremy L Thompson         break;
5724b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
573e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
5743a2968d6SJeremy L Thompson         if (is_at_points) {
5759123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
5769123fb08SJeremy L Thompson 
5779123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_c" << var_suffix << ", s_B"
5789123fb08SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5793a2968d6SJeremy L Thompson         } else if (use_3d_slices) {
5809123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
5819123fb08SJeremy L Thompson 
5829123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B"
5839123fb08SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5849123fb08SJeremy L Thompson         } else if (is_tensor) {
5859123fb08SJeremy L Thompson           bool        is_collocated = dim == 3 && Q_1d >= P_1d;
5869123fb08SJeremy L Thompson           std::string function_name =
5879123fb08SJeremy L Thompson               (dim == 1 ? "GradTranspose" : (is_collocated ? "GradTransposeTensorCollocated" : "GradTransposeTensor")) + std::to_string(dim) + "d";
5889123fb08SJeremy L Thompson 
5899123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B"
5909123fb08SJeremy L Thompson                << var_suffix << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
5914b3e95d5SJeremy L Thompson         } else {
5929123fb08SJeremy L Thompson           std::string function_name = "GradTransposeNonTensor";
5939123fb08SJeremy L Thompson 
5949123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", dim, " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix
5959123fb08SJeremy L Thompson                << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
5964b3e95d5SJeremy L Thompson         }
5974b3e95d5SJeremy L Thompson         break;
5984b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
5994b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT:
6004b3e95d5SJeremy L Thompson         break;  // Should not occur
6014b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
6024b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
6034b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
6044b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
6054b3e95d5SJeremy L Thompson     }
6064b3e95d5SJeremy L Thompson   }
6073a2968d6SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
6084b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6094b3e95d5SJeremy L Thompson }
6104b3e95d5SJeremy L Thompson 
6114b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
6124b3e95d5SJeremy L Thompson // QFunction
6134b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
6143a2968d6SJeremy L Thompson static int CeedOperatorBuildKernelQFunction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt dim, CeedInt max_num_points,
6153a2968d6SJeremy L Thompson                                                     CeedInt num_input_fields, CeedOperatorField *op_input_fields, CeedQFunctionField *qf_input_fields,
6164b3e95d5SJeremy L Thompson                                                     CeedInt num_output_fields, CeedOperatorField *op_output_fields,
6179123fb08SJeremy L Thompson                                                     CeedQFunctionField *qf_output_fields, std::string qfunction_name, CeedInt Q_1d, bool is_tensor,
6189123fb08SJeremy L Thompson                                                     bool is_at_points, bool use_3d_slices) {
6199123fb08SJeremy L Thompson   std::string         Q_name    = is_tensor ? "Q_1d" : "Q";
6204b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
6214b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
6224b3e95d5SJeremy L Thompson 
6238b97b69aSJeremy L Thompson   // Setup output arrays
6244b3e95d5SJeremy L Thompson   code << "\n    // -- Output field setup\n";
6254b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
6264b3e95d5SJeremy L Thompson     std::string var_suffix = "_out_" + std::to_string(i);
6274b3e95d5SJeremy L Thompson 
6284b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
6294b3e95d5SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
6303a2968d6SJeremy L Thompson     switch (eval_mode) {
6313a2968d6SJeremy L Thompson       case CEED_EVAL_NONE:
6323a2968d6SJeremy L Thompson         if (is_at_points) {
6333a2968d6SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "];\n";
6343a2968d6SJeremy L Thompson         } else {
6359123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
6364b3e95d5SJeremy L Thompson         }
6373a2968d6SJeremy L Thompson         break;
6383a2968d6SJeremy L Thompson       case CEED_EVAL_INTERP:
6393a2968d6SJeremy L Thompson         if (is_at_points) {
6403a2968d6SJeremy L Thompson           // Accumulator for point data
6419123fb08SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
6429123fb08SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "; i++) {\n";
6433a2968d6SJeremy L Thompson           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
6443a2968d6SJeremy L Thompson           code << "    }\n";
6453a2968d6SJeremy L Thompson         } else {
6469123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
6473a2968d6SJeremy L Thompson         }
6483a2968d6SJeremy L Thompson         break;
6493a2968d6SJeremy L Thompson       case CEED_EVAL_GRAD:
6503a2968d6SJeremy L Thompson         if (is_at_points) {
6513a2968d6SJeremy L Thompson           // Accumulator for point data
6529123fb08SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "*dim];\n";
6539123fb08SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "; i++) {\n";
6543a2968d6SJeremy L Thompson           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
6553a2968d6SJeremy L Thompson           code << "    }\n";
6563a2968d6SJeremy L Thompson         } else if (use_3d_slices) {
6574b3e95d5SJeremy L Thompson           // Accumulator for gradient slices
6584b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
6594b3e95d5SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n";
6604b3e95d5SJeremy L Thompson           code << "      r_q" << var_suffix << "[i] = 0.0;\n";
6614b3e95d5SJeremy L Thompson           code << "    }\n";
6624b3e95d5SJeremy L Thompson         } else {
6639123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
6644b3e95d5SJeremy L Thompson         }
6653a2968d6SJeremy L Thompson         break;
6663a2968d6SJeremy L Thompson       case CEED_EVAL_WEIGHT:
6673a2968d6SJeremy L Thompson         break;
6683a2968d6SJeremy L Thompson         // LCOV_EXCL_START
6693a2968d6SJeremy L Thompson       case CEED_EVAL_DIV:
6703a2968d6SJeremy L Thompson       case CEED_EVAL_CURL:
6713a2968d6SJeremy L Thompson         break;  // TODO: Not implemented
6723a2968d6SJeremy L Thompson                 // LCOV_EXCL_STOP
6734b3e95d5SJeremy L Thompson     }
6744b3e95d5SJeremy L Thompson   }
6754b3e95d5SJeremy L Thompson 
6763a2968d6SJeremy L Thompson   if (is_at_points) {
6773a2968d6SJeremy L Thompson     // We need to handle batches of points
6783a2968d6SJeremy L Thompson     code << "\n    // Note: Using batches of points\n";
6793a2968d6SJeremy L Thompson     code << "    const CeedInt point_loop_bound = (blockDim.x * blockDim.y) * ceil(1.0 * max_num_points / (blockDim.x * blockDim.y));\n\n";
6803a2968d6SJeremy L Thompson     code << "    #pragma unroll\n";
6813a2968d6SJeremy L Thompson     code << "    for (CeedInt i = threadIdx.x + threadIdx.y * blockDim.x; i < point_loop_bound; i += blockDim.x * blockDim.y) {\n";
6823a2968d6SJeremy L Thompson     code << "      const CeedInt p = i % max_num_points;\n\n";
6833a2968d6SJeremy L Thompson 
6843a2968d6SJeremy L Thompson     code << "      // -- Coordinates\n";
6853a2968d6SJeremy L Thompson     code << "      CeedScalar r_x[dim];\n";
6863a2968d6SJeremy 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";
6873a2968d6SJeremy L Thompson 
6883a2968d6SJeremy L Thompson     code << "      // -- Input fields\n";
6893a2968d6SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
6903a2968d6SJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
6913a2968d6SJeremy L Thompson 
6923a2968d6SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
6933a2968d6SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
6943a2968d6SJeremy L Thompson       // Basis action
6953a2968d6SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
6963a2968d6SJeremy L Thompson       switch (eval_mode) {
6973a2968d6SJeremy L Thompson         case CEED_EVAL_NONE:
6983a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
6993a2968d6SJeremy L Thompson           code << "      ReadPoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
7003a2968d6SJeremy L Thompson                << ", max_num_points>(data, elem, p, max_num_points, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
7013a2968d6SJeremy L Thompson           break;
7023a2968d6SJeremy L Thompson         case CEED_EVAL_INTERP:
7033a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7043a2968d6SJeremy L Thompson           code << "      InterpAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix
7053a2968d6SJeremy L Thompson                << ", r_x, r_s" << var_suffix << ");\n";
7063a2968d6SJeremy L Thompson           break;
7073a2968d6SJeremy L Thompson         case CEED_EVAL_GRAD:
7083a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
7093a2968d6SJeremy L Thompson           code << "      GradAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix
7103a2968d6SJeremy L Thompson                << ", r_x, r_s" << var_suffix << ");\n";
7113a2968d6SJeremy L Thompson           break;
7123a2968d6SJeremy L Thompson         case CEED_EVAL_WEIGHT:
7133a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
7143a2968d6SJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = 1.0;\n";
7153a2968d6SJeremy L Thompson           break;
7163a2968d6SJeremy L Thompson           // LCOV_EXCL_START
7173a2968d6SJeremy L Thompson         case CEED_EVAL_DIV:
7183a2968d6SJeremy L Thompson         case CEED_EVAL_CURL:
7193a2968d6SJeremy L Thompson           break;  // TODO: Not implemented
7203a2968d6SJeremy L Thompson                   // LCOV_EXCL_STOP
7213a2968d6SJeremy L Thompson       }
7223a2968d6SJeremy L Thompson     }
7233a2968d6SJeremy L Thompson     code << "\n      // -- Output fields\n";
7243a2968d6SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
7253a2968d6SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
7263a2968d6SJeremy L Thompson 
7273a2968d6SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
7283a2968d6SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
7293a2968d6SJeremy L Thompson       // Basis action
7303a2968d6SJeremy L Thompson       switch (eval_mode) {
7313a2968d6SJeremy L Thompson         case CEED_EVAL_NONE:
7323a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7333a2968d6SJeremy L Thompson           break;
7343a2968d6SJeremy L Thompson         case CEED_EVAL_INTERP:
7353a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7363a2968d6SJeremy L Thompson           break;
7373a2968d6SJeremy L Thompson         case CEED_EVAL_GRAD:
7383a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
7393a2968d6SJeremy L Thompson           break;
7403a2968d6SJeremy L Thompson           // LCOV_EXCL_START
7413a2968d6SJeremy L Thompson         case CEED_EVAL_WEIGHT:
7423a2968d6SJeremy L Thompson           break;  // Should not occur
7433a2968d6SJeremy L Thompson         case CEED_EVAL_DIV:
7443a2968d6SJeremy L Thompson         case CEED_EVAL_CURL:
7453a2968d6SJeremy L Thompson           break;  // TODO: Not implemented
7463a2968d6SJeremy L Thompson                   // LCOV_EXCL_STOP
7473a2968d6SJeremy L Thompson       }
7483a2968d6SJeremy L Thompson     }
7493a2968d6SJeremy L Thompson 
7503a2968d6SJeremy L Thompson   } else if (use_3d_slices) {
7514b3e95d5SJeremy L Thompson     // We treat quadrature points per slice in 3d to save registers
7524b3e95d5SJeremy L Thompson     code << "\n    // Note: Using planes of 3D elements\n";
7534b3e95d5SJeremy L Thompson     code << "    #pragma unroll\n";
7544b3e95d5SJeremy L Thompson     code << "    for (CeedInt q = 0; q < " << Q_name << "; q++) {\n";
7554b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
7564b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
7574b3e95d5SJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
7584b3e95d5SJeremy L Thompson 
7594b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
7604b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
7614b3e95d5SJeremy L Thompson       // Basis action
7624b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
7634b3e95d5SJeremy L Thompson       switch (eval_mode) {
7644b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
7654b3e95d5SJeremy L Thompson           bool is_strided;
7664b3e95d5SJeremy L Thompson 
7674b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7684b3e95d5SJeremy L Thompson 
7694b3e95d5SJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
7704b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
7714b3e95d5SJeremy L Thompson           if (is_strided) {
7724b3e95d5SJeremy L Thompson             bool    has_backend_strides;
7734b3e95d5SJeremy L Thompson             CeedInt num_elem, elem_size;
7744b3e95d5SJeremy L Thompson 
7754b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
7764b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
7774b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
7784b3e95d5SJeremy L Thompson             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
7794b3e95d5SJeremy L Thompson 
7804b3e95d5SJeremy L Thompson             if (!has_backend_strides) {
7814b3e95d5SJeremy L Thompson               CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
7824b3e95d5SJeremy L Thompson             }
7834b3e95d5SJeremy L Thompson             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
784f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStrided3d<num_comp" << var_suffix << ", " << Q_name << ", " << strides[0] << ", " << strides[1] << ", "
7854b3e95d5SJeremy L Thompson                  << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n";
7864b3e95d5SJeremy L Thompson           } else {
7874b3e95d5SJeremy L Thompson             CeedSize                 l_size = 0;
7884b3e95d5SJeremy L Thompson             CeedInt                  comp_stride;
7894b3e95d5SJeremy L Thompson             CeedElemRestriction_Hip *rstr_data;
7904b3e95d5SJeremy L Thompson 
7914b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
7924b3e95d5SJeremy L Thompson             code << "      const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
7934b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
7944b3e95d5SJeremy L Thompson             code << "      // CompStride: " << comp_stride << "\n";
7954b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
7964b3e95d5SJeremy L Thompson             data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
797f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStandard3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix
7984b3e95d5SJeremy L Thompson                  << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
7994b3e95d5SJeremy L Thompson           }
8009123fb08SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
8014b3e95d5SJeremy L Thompson           break;
8024b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
8034b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
8044b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n";
8054b3e95d5SJeremy L Thompson           code << "        r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n";
8064b3e95d5SJeremy L Thompson           code << "      }\n";
8074b3e95d5SJeremy L Thompson           break;
8084b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
8094b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
810f815fac9SJeremy L Thompson           code << "      GradColloSlice3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_q" << var_suffix << ", s_G" << var_suffix
811f815fac9SJeremy L Thompson                << ", r_s" << var_suffix << ");\n";
8124b3e95d5SJeremy L Thompson           break;
8134b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
8144b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
8154b3e95d5SJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n";
8163a2968d6SJeremy L Thompson           break;
8174b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
8184b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
8194b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
8204b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
8214b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
8224b3e95d5SJeremy L Thompson       }
8234b3e95d5SJeremy L Thompson     }
8244b3e95d5SJeremy L Thompson     code << "\n      // -- Output fields\n";
8254b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
8264b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
8274b3e95d5SJeremy L Thompson 
8284b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
8294b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
8304b3e95d5SJeremy L Thompson       // Basis action
8314b3e95d5SJeremy L Thompson       switch (eval_mode) {
8324b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
8334b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
8343a2968d6SJeremy L Thompson           break;
8354b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
8364b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
8374b3e95d5SJeremy L Thompson           break;
8384b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
8394b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
8404b3e95d5SJeremy L Thompson           break;
8414b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
8424b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
8434b3e95d5SJeremy L Thompson           break;  // Should not occur
8444b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
8454b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
8464b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
8474b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
8484b3e95d5SJeremy L Thompson       }
8494b3e95d5SJeremy L Thompson     }
8504b3e95d5SJeremy L Thompson   } else {
8514b3e95d5SJeremy L Thompson     code << "\n    // Note: Using full elements\n";
8524b3e95d5SJeremy L Thompson     code << "    {\n";
8534b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
8544b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
8554b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
8564b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n";
8574b3e95d5SJeremy L Thompson     }
8584b3e95d5SJeremy L Thompson     code << "      // -- Output fields\n";
8594b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
8604b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
8614b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n";
8624b3e95d5SJeremy L Thompson     }
8634b3e95d5SJeremy L Thompson   }
8644b3e95d5SJeremy L Thompson 
8654b3e95d5SJeremy L Thompson   // Input and output buffers
8664b3e95d5SJeremy L Thompson   code << "\n      // -- QFunction inputs and outputs\n";
8674b3e95d5SJeremy L Thompson   code << "      // ---- Inputs\n";
8684b3e95d5SJeremy L Thompson   code << "      CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n";
8694b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
8704b3e95d5SJeremy L Thompson     code << "      // ------ Input field " << i << "\n";
8714b3e95d5SJeremy L Thompson     code << "      inputs[" << i << "] = r_s_in_" << i << ";\n";
8724b3e95d5SJeremy L Thompson   }
8734b3e95d5SJeremy L Thompson   code << "      // ---- Outputs\n";
8744b3e95d5SJeremy L Thompson   code << "      CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n";
8754b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
8764b3e95d5SJeremy L Thompson     code << "      // ------ Output field " << i << "\n";
8774b3e95d5SJeremy L Thompson     code << "      outputs[" << i << "] = r_s_out_" << i << ";\n";
8784b3e95d5SJeremy L Thompson   }
8794b3e95d5SJeremy L Thompson 
8804b3e95d5SJeremy L Thompson   // Apply QFunction
8814b3e95d5SJeremy L Thompson   code << "\n      // -- Apply QFunction\n";
8824b3e95d5SJeremy L Thompson   code << "      " << qfunction_name << "(ctx, ";
8839123fb08SJeremy L Thompson   if (dim != 3 || is_at_points || use_3d_slices || !is_tensor) {
8844b3e95d5SJeremy L Thompson     code << "1";
8854b3e95d5SJeremy L Thompson   } else {
8869123fb08SJeremy L Thompson     code << Q_name;
8874b3e95d5SJeremy L Thompson   }
8884b3e95d5SJeremy L Thompson   code << ", inputs, outputs);\n";
8894b3e95d5SJeremy L Thompson 
8903a2968d6SJeremy L Thompson   if (is_at_points) {
8913a2968d6SJeremy L Thompson     // Map back to coefficients
8923a2968d6SJeremy L Thompson     code << "\n      // -- Output fields\n";
8933a2968d6SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
8943a2968d6SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
8953a2968d6SJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
8963a2968d6SJeremy L Thompson 
8973a2968d6SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
8983a2968d6SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
8993a2968d6SJeremy L Thompson       // Basis action
9003a2968d6SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
9013a2968d6SJeremy L Thompson       switch (eval_mode) {
9023a2968d6SJeremy L Thompson         case CEED_EVAL_NONE: {
9033a2968d6SJeremy L Thompson           CeedInt             comp_stride;
9043a2968d6SJeremy L Thompson           CeedElemRestriction elem_rstr;
9053a2968d6SJeremy L Thompson 
9063a2968d6SJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
9073a2968d6SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
9083a2968d6SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
9093a2968d6SJeremy L Thompson           code << "      const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
9103a2968d6SJeremy L Thompson           code << "      WritePoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
9113a2968d6SJeremy L Thompson                << ", max_num_points>(data, elem, i, points.num_per_elem[elem], indices.outputs[" << i << "]"
9123a2968d6SJeremy L Thompson                << ", r_s" << var_suffix << ", d" << var_suffix << ");\n";
9133a2968d6SJeremy L Thompson           break;
9143a2968d6SJeremy L Thompson         }
9153a2968d6SJeremy L Thompson         case CEED_EVAL_INTERP:
9163a2968d6SJeremy L Thompson           code << "      if (i >= points.num_per_elem[elem]) {\n";
9173a2968d6SJeremy L Thompson           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) r_s" << var_suffix << "[j] = 0.0;\n";
9183a2968d6SJeremy L Thompson           code << "      }\n";
9193a2968d6SJeremy L Thompson           code << "      InterpTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s"
9203a2968d6SJeremy L Thompson                << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
9213a2968d6SJeremy L Thompson           break;
9223a2968d6SJeremy L Thompson         case CEED_EVAL_GRAD:
9233a2968d6SJeremy L Thompson           code << "      if (i >= points.num_per_elem[elem]) {\n";
9243a2968d6SJeremy L Thompson           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "*dim; j++) r_s" << var_suffix << "[j] = 0.0;\n";
9253a2968d6SJeremy L Thompson           code << "      }\n";
9263a2968d6SJeremy L Thompson           code << "      GradTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s"
9273a2968d6SJeremy L Thompson                << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
9283a2968d6SJeremy L Thompson           break;
9293a2968d6SJeremy L Thompson           // LCOV_EXCL_START
9303a2968d6SJeremy L Thompson         case CEED_EVAL_WEIGHT:
9313a2968d6SJeremy L Thompson           break;  // Should not occur
9323a2968d6SJeremy L Thompson         case CEED_EVAL_DIV:
9333a2968d6SJeremy L Thompson         case CEED_EVAL_CURL:
9343a2968d6SJeremy L Thompson           break;  // TODO: Not implemented
9353a2968d6SJeremy L Thompson                   // LCOV_EXCL_STOP
9363a2968d6SJeremy L Thompson       }
9373a2968d6SJeremy L Thompson     }
9383a2968d6SJeremy L Thompson   } else if (use_3d_slices) {
9394b3e95d5SJeremy L Thompson     // Copy or apply transpose grad, if needed
9403a2968d6SJeremy L Thompson     code << "\n      // -- Output fields\n";
9414b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
9424b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
9434b3e95d5SJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
9444b3e95d5SJeremy L Thompson 
9454b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
9464b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
9474b3e95d5SJeremy L Thompson       // Basis action
9484b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
9494b3e95d5SJeremy L Thompson       switch (eval_mode) {
9504b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
9514b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
9524b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
9534b3e95d5SJeremy L Thompson           code << "      }\n";
9543a2968d6SJeremy L Thompson           break;
9554b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
9564b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
9574b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
9584b3e95d5SJeremy L Thompson           code << "      }\n";
9594b3e95d5SJeremy L Thompson           break;
9604b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
961f815fac9SJeremy L Thompson           code << "      GradColloSliceTranspose3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_s" << var_suffix << ", s_G"
962f815fac9SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
9634b3e95d5SJeremy L Thompson           break;
9644b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
9654b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
9664b3e95d5SJeremy L Thompson           break;  // Should not occur
9674b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
9684b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
9694b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
9704b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
9714b3e95d5SJeremy L Thompson       }
9724b3e95d5SJeremy L Thompson     }
9734b3e95d5SJeremy L Thompson   }
9744b3e95d5SJeremy L Thompson   code << "    }\n";
9754b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9764b3e95d5SJeremy L Thompson }
9774b3e95d5SJeremy L Thompson 
9784b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
9799e201c85SYohann // Build single operator kernel
9807d8d0e25Snbeams //------------------------------------------------------------------------------
9818d12f40eSJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Hip_gen(CeedOperator op, bool *is_good_build) {
9823a2968d6SJeremy L Thompson   bool                   is_tensor = true, is_at_points = false, use_3d_slices = false;
9837d8d0e25Snbeams   Ceed                   ceed;
9843a2968d6SJeremy L Thompson   CeedInt                Q_1d, num_input_fields, num_output_fields, dim = 1, max_num_points = 0, coords_comp_stride = 0;
985b7453713SJeremy L Thompson   CeedQFunctionField    *qf_input_fields, *qf_output_fields;
986b7453713SJeremy L Thompson   CeedQFunction_Hip_gen *qf_data;
987b7453713SJeremy L Thompson   CeedQFunction          qf;
988b7453713SJeremy L Thompson   CeedOperatorField     *op_input_fields, *op_output_fields;
989b7453713SJeremy L Thompson   CeedOperator_Hip_gen  *data;
9904b3e95d5SJeremy L Thompson   std::ostringstream     code;
9914b3e95d5SJeremy L Thompson 
9928d12f40eSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
9934b3e95d5SJeremy L Thompson   {
9944b3e95d5SJeremy L Thompson     bool is_setup_done;
995b7453713SJeremy L Thompson 
996b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
9978d12f40eSJeremy L Thompson     if (is_setup_done) {
9988d12f40eSJeremy L Thompson       *is_good_build = !data->use_fallback;
9998d12f40eSJeremy L Thompson       return CEED_ERROR_SUCCESS;
10008d12f40eSJeremy L Thompson     }
10014b3e95d5SJeremy L Thompson   }
1002b7453713SJeremy L Thompson 
10038d12f40eSJeremy L Thompson   // Check field compatibility
10048d12f40eSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
10058d12f40eSJeremy L Thompson   {
10068d12f40eSJeremy L Thompson     bool has_shared_bases = true, is_all_tensor = true, is_all_nontensor = true;
10078d12f40eSJeremy L Thompson 
10088d12f40eSJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
10098d12f40eSJeremy L Thompson       CeedBasis basis;
10108d12f40eSJeremy L Thompson 
10118d12f40eSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
10128d12f40eSJeremy L Thompson       if (basis != CEED_BASIS_NONE) {
10138d12f40eSJeremy L Thompson         bool        is_tensor = true;
10148d12f40eSJeremy L Thompson         const char *resource;
10158d12f40eSJeremy L Thompson         char       *resource_root;
10168d12f40eSJeremy L Thompson         Ceed        basis_ceed;
10178d12f40eSJeremy L Thompson 
10188d12f40eSJeremy L Thompson         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
1019c9192acaSJeremy L Thompson         is_all_tensor    = is_all_tensor && is_tensor;
1020c9192acaSJeremy L Thompson         is_all_nontensor = is_all_nontensor && !is_tensor;
10218d12f40eSJeremy L Thompson         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
10228d12f40eSJeremy L Thompson         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
10238d12f40eSJeremy L Thompson         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
1024c9192acaSJeremy L Thompson         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/hip/shared");
10258d12f40eSJeremy L Thompson         CeedCallBackend(CeedFree(&resource_root));
10268d12f40eSJeremy L Thompson         CeedCallBackend(CeedDestroy(&basis_ceed));
10278d12f40eSJeremy L Thompson       }
10288d12f40eSJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis));
10298d12f40eSJeremy L Thompson     }
10308d12f40eSJeremy L Thompson 
10318d12f40eSJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
10328d12f40eSJeremy L Thompson       CeedBasis basis;
10338d12f40eSJeremy L Thompson 
10348d12f40eSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
10358d12f40eSJeremy L Thompson       if (basis != CEED_BASIS_NONE) {
10368d12f40eSJeremy L Thompson         bool        is_tensor = true;
10378d12f40eSJeremy L Thompson         const char *resource;
10388d12f40eSJeremy L Thompson         char       *resource_root;
10398d12f40eSJeremy L Thompson         Ceed        basis_ceed;
10408d12f40eSJeremy L Thompson 
10418d12f40eSJeremy L Thompson         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
1042c9192acaSJeremy L Thompson         is_all_tensor    = is_all_tensor && is_tensor;
1043c9192acaSJeremy L Thompson         is_all_nontensor = is_all_nontensor && !is_tensor;
10448d12f40eSJeremy L Thompson 
10458d12f40eSJeremy L Thompson         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
10468d12f40eSJeremy L Thompson         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
10478d12f40eSJeremy L Thompson         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
1048c9192acaSJeremy L Thompson         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/hip/shared");
10498d12f40eSJeremy L Thompson         CeedCallBackend(CeedFree(&resource_root));
10508d12f40eSJeremy L Thompson         CeedCallBackend(CeedDestroy(&basis_ceed));
10518d12f40eSJeremy L Thompson       }
10528d12f40eSJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis));
10538d12f40eSJeremy L Thompson     }
10548d12f40eSJeremy L Thompson     // -- Fallback to ref if not all bases are shared
10558d12f40eSJeremy L Thompson     if (!has_shared_bases || (!is_all_tensor && !is_all_nontensor)) {
10568d12f40eSJeremy L Thompson       *is_good_build = false;
10578d12f40eSJeremy L Thompson       return CEED_ERROR_SUCCESS;
10588d12f40eSJeremy L Thompson     }
10598d12f40eSJeremy L Thompson   }
1060b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1061b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1062b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
1063b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
10647d8d0e25Snbeams 
10654b3e95d5SJeremy L Thompson   // Get operator data
10663a2968d6SJeremy L Thompson   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
10674b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelData_Hip_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, op_output_fields,
10684b3e95d5SJeremy L Thompson                                                       qf_output_fields, &data->max_P_1d, &Q_1d, &dim, &is_tensor, &use_3d_slices));
10694b3e95d5SJeremy L Thompson   if (dim == 0) dim = 1;
10704b3e95d5SJeremy L Thompson   data->dim = dim;
10713a2968d6SJeremy L Thompson   if (is_at_points) {
10723a2968d6SJeremy L Thompson     CeedElemRestriction_Hip *rstr_data;
10733a2968d6SJeremy L Thompson     CeedElemRestriction      rstr_points = NULL;
10744b3e95d5SJeremy L Thompson 
10753a2968d6SJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
10763a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points));
10773a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetCompStride(rstr_points, &coords_comp_stride));
10783a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(rstr_points, &rstr_data));
10793a2968d6SJeremy L Thompson     data->points.indices = (CeedInt *)rstr_data->d_offsets;
10803a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
10813a2968d6SJeremy L Thompson   }
10823a2968d6SJeremy L Thompson   if (is_at_points) use_3d_slices = false;
10833a2968d6SJeremy L Thompson   if (Q_1d == 0) {
10843a2968d6SJeremy L Thompson     if (is_at_points) Q_1d = max_num_points;
10853a2968d6SJeremy L Thompson     else CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q_1d));
10864b3e95d5SJeremy L Thompson   }
10874b3e95d5SJeremy L Thompson   data->Q_1d = Q_1d;
10884b3e95d5SJeremy L Thompson 
10890b454692Sjeremylt   // Check for restriction only identity operator
10904b3e95d5SJeremy L Thompson   {
10914b3e95d5SJeremy L Thompson     bool is_identity_qf;
10924b3e95d5SJeremy L Thompson 
10932b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
10940b454692Sjeremylt     if (is_identity_qf) {
10959e201c85SYohann       CeedEvalMode eval_mode_in, eval_mode_out;
1096b7453713SJeremy L Thompson 
10972b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
10982b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
10996574a04fSJeremy L Thompson       CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
11006574a04fSJeremy L Thompson                 "Backend does not implement restriction only identity operators");
11010b454692Sjeremylt     }
11024b3e95d5SJeremy L Thompson   }
1103b2165e7aSSebastian Grimberg 
1104b2165e7aSSebastian Grimberg   // Load basis source files
11059123fb08SJeremy L Thompson   if (is_tensor) {
11069c25dd66SJeremy L Thompson     code << "// Tensor basis source\n";
11079c25dd66SJeremy L Thompson     code << "#include <ceed/jit-source/hip/hip-shared-basis-tensor-templates.h>\n\n";
11089123fb08SJeremy L Thompson   } else {
11099123fb08SJeremy L Thompson     code << "// Non-tensor basis source\n";
11109123fb08SJeremy L Thompson     code << "#include <ceed/jit-source/hip/hip-shared-basis-nontensor-templates.h>\n\n";
11119123fb08SJeremy L Thompson   }
11129123fb08SJeremy L Thompson   if (is_at_points) {
11133a2968d6SJeremy L Thompson     code << "// AtPoints basis source\n";
11143a2968d6SJeremy L Thompson     code << "#include <ceed/jit-source/hip/hip-shared-basis-tensor-at-points-templates.h>\n\n";
11159123fb08SJeremy L Thompson   }
11169c25dd66SJeremy L Thompson   code << "// CodeGen operator source\n";
11179c25dd66SJeremy L Thompson   code << "#include <ceed/jit-source/hip/hip-gen-templates.h>\n\n";
11187d8d0e25Snbeams 
11194b3e95d5SJeremy L Thompson   // Get QFunction name
11204b3e95d5SJeremy L Thompson   std::string qfunction_name(qf_data->qfunction_name);
11214b3e95d5SJeremy L Thompson   std::string operator_name;
11224b3e95d5SJeremy L Thompson 
112309095acaSJeremy L Thompson   operator_name = "CeedKernelHipGenOperator_" + qfunction_name;
11247d8d0e25Snbeams 
11259e201c85SYohann   // Define CEED_Q_VLA
11269e201c85SYohann   code << "\n#undef CEED_Q_VLA\n";
11279123fb08SJeremy L Thompson   if (dim != 3 || is_at_points || use_3d_slices || !is_tensor) {
11289e201c85SYohann     code << "#define CEED_Q_VLA 1\n\n";
11299e201c85SYohann   } else {
11309e201c85SYohann     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
11319e201c85SYohann   }
11329e201c85SYohann 
11334b3e95d5SJeremy L Thompson   // Add user QFunction source
11344b3e95d5SJeremy L Thompson   {
11359c25dd66SJeremy L Thompson     const char *source_path;
11364b3e95d5SJeremy L Thompson 
11379c25dd66SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path));
11389c25dd66SJeremy L Thompson     CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/hip/gen backend requires QFunction source code file");
11399c25dd66SJeremy L Thompson 
11409c25dd66SJeremy L Thompson     code << "// User QFunction source\n";
11419c25dd66SJeremy L Thompson     code << "#include \"" << source_path << "\"\n\n";
11424b3e95d5SJeremy L Thompson   }
11437d8d0e25Snbeams 
11447d8d0e25Snbeams   // Setup
11457d8d0e25Snbeams   code << "\n// -----------------------------------------------------------------------------\n";
11464b3e95d5SJeremy L Thompson   code << "// Operator Kernel\n";
11474b3e95d5SJeremy L Thompson   code << "// \n";
11484b3e95d5SJeremy L Thompson   code << "// d_[in,out]_i:   CeedVector device array\n";
11494b3e95d5SJeremy L Thompson   code << "// r_[in,out]_e_i: Element vector register\n";
11504b3e95d5SJeremy L Thompson   code << "// r_[in,out]_q_i: Quadrature space vector register\n";
11519123fb08SJeremy L Thompson   code << "// r_[in,out]_c_i: AtPoints Chebyshev coefficients register\n";
11524b3e95d5SJeremy L Thompson   code << "// r_[in,out]_s_i: Quadrature space slice vector register\n";
11534b3e95d5SJeremy L Thompson   code << "// \n";
11544b3e95d5SJeremy L Thompson   code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n";
11554b3e95d5SJeremy L Thompson   code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n";
11564b3e95d5SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n";
1157b3e1519bSnbeams   code << "\nextern \"C\" __launch_bounds__(BLOCK_SIZE)\n";
11582b730f8bSJeremy L Thompson   code << "__global__ void " << operator_name
11593a2968d6SJeremy 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";
11604b3e95d5SJeremy L Thompson 
11614b3e95d5SJeremy L Thompson   // Scratch buffers
11629e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
11634b3e95d5SJeremy L Thompson     CeedEvalMode eval_mode;
11644b3e95d5SJeremy L Thompson 
11652b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
11669e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
11674b3e95d5SJeremy L Thompson       code << "  const CeedScalar *d_in_" << i << " = fields.inputs[" << i << "];\n";
11687d8d0e25Snbeams     }
11697d8d0e25Snbeams   }
11709e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
11714b3e95d5SJeremy L Thompson     code << "  CeedScalar *d_out_" << i << " = fields.outputs[" << i << "];\n";
11727d8d0e25Snbeams   }
11737d8d0e25Snbeams 
11749e201c85SYohann   code << "  const CeedInt dim = " << dim << ";\n";
11759123fb08SJeremy L Thompson   code << "  const CeedInt " << (is_tensor ? "Q_1d" : "Q") << " = " << Q_1d << ";\n";
11763a2968d6SJeremy L Thompson   if (is_at_points) {
11773a2968d6SJeremy L Thompson     code << "  const CeedInt max_num_points = " << max_num_points << ";\n";
11783a2968d6SJeremy L Thompson     code << "  const CeedInt coords_comp_stride = " << coords_comp_stride << ";\n";
11793a2968d6SJeremy L Thompson   }
11807d8d0e25Snbeams 
11814b3e95d5SJeremy L Thompson   // Shared data
11824b3e95d5SJeremy L Thompson   code << "  extern __shared__ CeedScalar slice[];\n";
11839e201c85SYohann   code << "  SharedData_Hip data;\n";
11849e201c85SYohann   code << "  data.t_id_x = threadIdx.x;\n";
11859e201c85SYohann   code << "  data.t_id_y = threadIdx.y;\n";
11869e201c85SYohann   code << "  data.t_id_z = threadIdx.z;\n";
11879e201c85SYohann   code << "  data.t_id  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
11889123fb08SJeremy L Thompson   code << "  data.slice = slice + data.t_id_z*T_1D" << ((!is_tensor || dim == 1) ? "" : "*T_1D") << ";\n";
11897d8d0e25Snbeams 
1190*9ee499e5SJeremy L Thompson   // -- Determine input mat reuse
1191*9ee499e5SJeremy L Thompson   CeedInt input_matrix_reuse[CEED_FIELD_MAX][3];  // field, is_input, eval_mode
1192*9ee499e5SJeremy L Thompson 
1193*9ee499e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1194*9ee499e5SJeremy L Thompson     input_matrix_reuse[i][0] = -1;
1195*9ee499e5SJeremy L Thompson   }
1196*9ee499e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1197*9ee499e5SJeremy L Thompson     CeedEvalMode eval_mode_i;
1198*9ee499e5SJeremy L Thompson     CeedBasis    basis_i;
1199*9ee499e5SJeremy L Thompson 
1200*9ee499e5SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode_i));
1201*9ee499e5SJeremy L Thompson     if (eval_mode_i == CEED_EVAL_WEIGHT) continue;
1202*9ee499e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis_i));
1203*9ee499e5SJeremy L Thompson     for (CeedInt j = 0; (input_matrix_reuse[i][0] == -1) && (j < i); j++) {
1204*9ee499e5SJeremy L Thompson       CeedEvalMode eval_mode_j;
1205*9ee499e5SJeremy L Thompson       CeedBasis    basis_j;
1206*9ee499e5SJeremy L Thompson 
1207*9ee499e5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
1208*9ee499e5SJeremy L Thompson       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
1209*9ee499e5SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
1210*9ee499e5SJeremy L Thompson       if (basis_i == basis_j) {
1211*9ee499e5SJeremy L Thompson         if (is_tensor) {
1212*9ee499e5SJeremy L Thompson           input_matrix_reuse[i][0] = j;
1213*9ee499e5SJeremy L Thompson           input_matrix_reuse[i][1] = true;
1214*9ee499e5SJeremy L Thompson           input_matrix_reuse[i][2] = eval_mode_j;
1215*9ee499e5SJeremy L Thompson         } else {
1216*9ee499e5SJeremy L Thompson           // For non-tensor can only re-use with the same eval mode
1217*9ee499e5SJeremy L Thompson           if (eval_mode_i == eval_mode_j) {
1218*9ee499e5SJeremy L Thompson             input_matrix_reuse[i][0] = j;
1219*9ee499e5SJeremy L Thompson             input_matrix_reuse[i][1] = true;
1220*9ee499e5SJeremy L Thompson             input_matrix_reuse[i][2] = eval_mode_j;
1221*9ee499e5SJeremy L Thompson           }
1222*9ee499e5SJeremy L Thompson         }
1223*9ee499e5SJeremy L Thompson       }
1224*9ee499e5SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis_j));
1225*9ee499e5SJeremy L Thompson     }
1226*9ee499e5SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis_i));
1227*9ee499e5SJeremy L Thompson   }
1228*9ee499e5SJeremy L Thompson 
1229*9ee499e5SJeremy L Thompson   // -- Determine output mat reuse
1230*9ee499e5SJeremy L Thompson   CeedInt output_matrix_reuse[CEED_FIELD_MAX][3];  // field, is_input, eval_mode
1231*9ee499e5SJeremy L Thompson 
1232*9ee499e5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1233*9ee499e5SJeremy L Thompson     output_matrix_reuse[i][0] = -1;
1234*9ee499e5SJeremy L Thompson   }
1235*9ee499e5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1236*9ee499e5SJeremy L Thompson     CeedEvalMode eval_mode_i;
1237*9ee499e5SJeremy L Thompson     CeedBasis    basis_i;
1238*9ee499e5SJeremy L Thompson 
1239*9ee499e5SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode_i));
1240*9ee499e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis_i));
1241*9ee499e5SJeremy L Thompson     for (CeedInt j = 0; (output_matrix_reuse[i][0] == -1) && (j < num_input_fields); j++) {
1242*9ee499e5SJeremy L Thompson       CeedEvalMode eval_mode_j;
1243*9ee499e5SJeremy L Thompson       CeedBasis    basis_j;
1244*9ee499e5SJeremy L Thompson 
1245*9ee499e5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
1246*9ee499e5SJeremy L Thompson       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
1247*9ee499e5SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
1248*9ee499e5SJeremy L Thompson       if (basis_i == basis_j) {
1249*9ee499e5SJeremy L Thompson         if (is_tensor) {
1250*9ee499e5SJeremy L Thompson           output_matrix_reuse[i][0] = j;
1251*9ee499e5SJeremy L Thompson           output_matrix_reuse[i][1] = true;
1252*9ee499e5SJeremy L Thompson           output_matrix_reuse[i][2] = eval_mode_j;
1253*9ee499e5SJeremy L Thompson         } else {
1254*9ee499e5SJeremy L Thompson           // For non-tensor can only re-use with the same eval mode
1255*9ee499e5SJeremy L Thompson           if (eval_mode_i == eval_mode_j) {
1256*9ee499e5SJeremy L Thompson             output_matrix_reuse[i][0] = j;
1257*9ee499e5SJeremy L Thompson             output_matrix_reuse[i][1] = true;
1258*9ee499e5SJeremy L Thompson             output_matrix_reuse[i][2] = eval_mode_j;
1259*9ee499e5SJeremy L Thompson           }
1260*9ee499e5SJeremy L Thompson         }
1261*9ee499e5SJeremy L Thompson       }
1262*9ee499e5SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis_j));
1263*9ee499e5SJeremy L Thompson     }
1264*9ee499e5SJeremy L Thompson     for (CeedInt j = 0; (output_matrix_reuse[i][0] == -1) && (j < i); j++) {
1265*9ee499e5SJeremy L Thompson       CeedEvalMode eval_mode_j;
1266*9ee499e5SJeremy L Thompson       CeedBasis    basis_j;
1267*9ee499e5SJeremy L Thompson 
1268*9ee499e5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[j], &eval_mode_j));
1269*9ee499e5SJeremy L Thompson       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
1270*9ee499e5SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[j], &basis_j));
1271*9ee499e5SJeremy L Thompson       if (basis_i == basis_j) {
1272*9ee499e5SJeremy L Thompson         if (is_tensor) {
1273*9ee499e5SJeremy L Thompson           output_matrix_reuse[i][0] = j;
1274*9ee499e5SJeremy L Thompson           output_matrix_reuse[i][1] = false;
1275*9ee499e5SJeremy L Thompson           output_matrix_reuse[i][2] = eval_mode_j;
1276*9ee499e5SJeremy L Thompson         } else {
1277*9ee499e5SJeremy L Thompson           // For non-tensor can only re-use with the same eval mode
1278*9ee499e5SJeremy L Thompson           if (eval_mode_i == eval_mode_j) {
1279*9ee499e5SJeremy L Thompson             output_matrix_reuse[i][0] = j;
1280*9ee499e5SJeremy L Thompson             output_matrix_reuse[i][1] = false;
1281*9ee499e5SJeremy L Thompson             output_matrix_reuse[i][2] = eval_mode_j;
1282*9ee499e5SJeremy L Thompson           }
1283*9ee499e5SJeremy L Thompson         }
1284*9ee499e5SJeremy L Thompson       }
1285*9ee499e5SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis_j));
1286*9ee499e5SJeremy L Thompson     }
1287*9ee499e5SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis_i));
1288*9ee499e5SJeremy L Thompson   }
1289*9ee499e5SJeremy L Thompson 
12907d8d0e25Snbeams   // Initialize constants, and matrices B and G
12914b3e95d5SJeremy L Thompson   code << "\n  // Input field constants and basis data\n";
12929e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1293*9ee499e5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_input_fields[i], qf_input_fields[i], input_matrix_reuse[i], Q_1d, true,
1294*9ee499e5SJeremy L Thompson                                                              is_tensor, is_at_points, use_3d_slices));
12957d8d0e25Snbeams   }
12964b3e95d5SJeremy L Thompson   code << "\n  // Output field constants and basis data\n";
12979e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1298*9ee499e5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_output_fields[i], qf_output_fields[i], output_matrix_reuse[i], Q_1d,
1299*9ee499e5SJeremy L Thompson                                                              false, is_tensor, is_at_points, use_3d_slices));
13004b3e95d5SJeremy L Thompson   }
13017d8d0e25Snbeams 
13024b3e95d5SJeremy L Thompson   // Loop over all elements
13034b3e95d5SJeremy L Thompson   code << "\n  // Element loop\n";
13047d8d0e25Snbeams   code << "  __syncthreads();\n";
13059e201c85SYohann   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
13064b3e95d5SJeremy L Thompson 
1307e93651e5SJeremy L Thompson   // -- Compute minimum buffer space needed
13083a2968d6SJeremy L Thompson   CeedInt max_rstr_buffer_size = 1;
1309e93651e5SJeremy L Thompson 
1310e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1311e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
1312e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
1313e93651e5SJeremy L Thompson 
1314e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
1315e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1316e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
13179123fb08SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_tensor && (dim >= 3) ? elem_size : 1));
1318681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1319e93651e5SJeremy L Thompson   }
1320e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1321e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
1322e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
1323e93651e5SJeremy L Thompson 
1324e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
1325e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1326e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
13279123fb08SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_tensor && (dim >= 3) ? elem_size : 1));
1328681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1329e93651e5SJeremy L Thompson   }
1330e93651e5SJeremy L Thompson   code << "    // Scratch restriction buffer space\n";
1331e93651e5SJeremy L Thompson   code << "    CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n";
1332e93651e5SJeremy L Thompson 
1333e93651e5SJeremy L Thompson   // -- Determine best input field processing order
1334e93651e5SJeremy L Thompson   CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX];
1335e93651e5SJeremy L Thompson 
1336e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1337e93651e5SJeremy L Thompson     field_rstr_in_buffer[i] = -1;
1338e93651e5SJeremy L Thompson     input_field_order[i]    = -1;
1339e93651e5SJeremy L Thompson   }
1340e93651e5SJeremy L Thompson   {
1341e93651e5SJeremy L Thompson     bool    is_ordered[CEED_FIELD_MAX];
1342e93651e5SJeremy L Thompson     CeedInt curr_index = 0;
1343e93651e5SJeremy L Thompson 
1344e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
1345e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
1346e93651e5SJeremy L Thompson       CeedVector          vec_i;
1347e93651e5SJeremy L Thompson       CeedElemRestriction rstr_i;
1348e93651e5SJeremy L Thompson 
1349e93651e5SJeremy L Thompson       if (is_ordered[i]) continue;
1350e93651e5SJeremy L Thompson       field_rstr_in_buffer[i]       = i;
1351e93651e5SJeremy L Thompson       is_ordered[i]                 = true;
1352e93651e5SJeremy L Thompson       input_field_order[curr_index] = i;
1353e93651e5SJeremy L Thompson       curr_index++;
1354034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
1355e93651e5SJeremy L Thompson       if (vec_i == CEED_VECTOR_NONE) continue;  // CEED_EVAL_WEIGHT
1356e93651e5SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
1357e93651e5SJeremy L Thompson       for (CeedInt j = i + 1; j < num_input_fields; j++) {
1358e93651e5SJeremy L Thompson         CeedVector          vec_j;
1359e93651e5SJeremy L Thompson         CeedElemRestriction rstr_j;
1360e93651e5SJeremy L Thompson 
1361e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
1362e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
1363e93651e5SJeremy L Thompson         if (rstr_i == rstr_j && vec_i == vec_j) {
1364e93651e5SJeremy L Thompson           field_rstr_in_buffer[j]       = i;
1365e93651e5SJeremy L Thompson           is_ordered[j]                 = true;
1366e93651e5SJeremy L Thompson           input_field_order[curr_index] = j;
1367e93651e5SJeremy L Thompson           curr_index++;
1368e93651e5SJeremy L Thompson         }
13693a2968d6SJeremy L Thompson         CeedCallBackend(CeedVectorDestroy(&vec_j));
13703a2968d6SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
1371e93651e5SJeremy L Thompson       }
13723a2968d6SJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec_i));
13733a2968d6SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
1374e93651e5SJeremy L Thompson     }
1375e93651e5SJeremy L Thompson   }
1376e93651e5SJeremy L Thompson 
13774b3e95d5SJeremy L Thompson   // -- Input restriction and basis
13783a2968d6SJeremy L Thompson   code << "\n    // -- Input field restrictions and basis actions\n";
13799e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1380e93651e5SJeremy L Thompson     CeedInt f = input_field_order[i];
1381e93651e5SJeremy L Thompson 
1382e93651e5SJeremy L Thompson     code << "    // ---- Input field " << f << "\n";
13837d8d0e25Snbeams 
13844b3e95d5SJeremy L Thompson     // ---- Restriction
1385e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Hip_gen(code, data, f, dim, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f], Q_1d,
13869123fb08SJeremy L Thompson                                                                true, is_tensor, is_at_points, use_3d_slices));
1387b7453713SJeremy L Thompson 
13884b3e95d5SJeremy L Thompson     // ---- Basis action
13899123fb08SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, f, dim, op_input_fields[f], qf_input_fields[f], Q_1d, true, is_tensor,
13909123fb08SJeremy L Thompson                                                          is_at_points, use_3d_slices));
13917d8d0e25Snbeams   }
13927d8d0e25Snbeams 
13934b3e95d5SJeremy L Thompson   // -- Q function
13943a2968d6SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelQFunction_Hip_gen(code, data, dim, max_num_points, num_input_fields, op_input_fields, qf_input_fields,
13959123fb08SJeremy L Thompson                                                            num_output_fields, op_output_fields, qf_output_fields, qfunction_name, Q_1d, is_tensor,
13969123fb08SJeremy L Thompson                                                            is_at_points, use_3d_slices));
13977d8d0e25Snbeams 
13984b3e95d5SJeremy L Thompson   // -- Output basis and restriction
13994b3e95d5SJeremy L Thompson   code << "\n    // -- Output field basis action and restrictions\n";
14009e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
14014b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
1402b7453713SJeremy L Thompson 
14034b3e95d5SJeremy L Thompson     // ---- Basis action
14049123fb08SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, i, dim, op_output_fields[i], qf_output_fields[i], Q_1d, false, is_tensor,
14059123fb08SJeremy L Thompson                                                          is_at_points, use_3d_slices));
14067d8d0e25Snbeams 
14074b3e95d5SJeremy L Thompson     // ---- Restriction
14083a2968d6SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Hip_gen(code, data, i, dim, NULL, op_output_fields[i], qf_output_fields[i], Q_1d, false,
14099123fb08SJeremy L Thompson                                                                is_tensor, is_at_points, use_3d_slices));
14107d8d0e25Snbeams   }
14117d8d0e25Snbeams 
14124b3e95d5SJeremy L Thompson   // Close loop and function
14137d8d0e25Snbeams   code << "  }\n";
14147d8d0e25Snbeams   code << "}\n";
14157d8d0e25Snbeams   code << "// -----------------------------------------------------------------------------\n\n";
14167d8d0e25Snbeams 
1417539ec17dSJeremy L Thompson   CeedInt block_sizes[3] = {0, 0, 0};
14189e201c85SYohann   CeedInt num_elem;
1419b7453713SJeremy L Thompson 
14203a2968d6SJeremy L Thompson   // Compile
14212b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
14229123fb08SJeremy L Thompson   CeedCallBackend(BlockGridCalculate_Hip_gen(is_tensor ? dim : 1, num_elem, data->max_P_1d, Q_1d, block_sizes));
14238d12f40eSJeremy L Thompson   {
14248d12f40eSJeremy L Thompson     bool is_compile_good = false;
14258d12f40eSJeremy L Thompson 
14268d12f40eSJeremy L Thompson     CeedCallBackend(CeedTryCompile_Hip(ceed, code.str().c_str(), &is_compile_good, &data->module, 2, "T_1D", block_sizes[0], "BLOCK_SIZE",
14272b730f8bSJeremy L Thompson                                        block_sizes[0] * block_sizes[1] * block_sizes[2]));
14288d12f40eSJeremy L Thompson     if (is_compile_good) {
14298d12f40eSJeremy L Thompson       *is_good_build = true;
1430eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, operator_name.c_str(), &data->op));
14318d12f40eSJeremy L Thompson     } else {
14328d12f40eSJeremy L Thompson       *is_good_build     = false;
14338d12f40eSJeremy L Thompson       data->use_fallback = true;
14348d12f40eSJeremy L Thompson     }
14358d12f40eSJeremy L Thompson   }
14362b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
14379bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
1438c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionDestroy(&qf));
1439e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14407d8d0e25Snbeams }
14412a86cc9dSSebastian Grimberg 
14427d8d0e25Snbeams //------------------------------------------------------------------------------
1443