xref: /libCEED/rust/libceed-sys/c-src/backends/hip-gen/ceed-hip-gen-operator-build.cpp (revision 45a787f75d07c3f171308ee0d0be6daefd4754b0)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
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 
24*45a787f7SJeremy L Thompson struct FieldReuse_Hip {
25*45a787f7SJeremy L Thompson   CeedInt      index;
26*45a787f7SJeremy L Thompson   bool         is_input;
27*45a787f7SJeremy L Thompson   CeedEvalMode eval_mode;
28*45a787f7SJeremy L Thompson };
29*45a787f7SJeremy L Thompson 
30b3e1519bSnbeams //------------------------------------------------------------------------------
31b3e1519bSnbeams // Calculate the block size used for launching the operator kernel
32b3e1519bSnbeams //------------------------------------------------------------------------------
332b730f8bSJeremy 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) {
343a2968d6SJeremy L Thompson   const CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
35b3e1519bSnbeams   if (dim == 1) {
363a2968d6SJeremy L Thompson     CeedInt elems_per_block = 64 * thread_1d > 256 ? 256 / thread_1d : 64;
37b7453713SJeremy L Thompson 
389e201c85SYohann     elems_per_block = elems_per_block > 0 ? elems_per_block : 1;
393a2968d6SJeremy L Thompson     block_sizes[0]  = thread_1d;
40b3e1519bSnbeams     block_sizes[1]  = 1;
419e201c85SYohann     block_sizes[2]  = elems_per_block;
42b3e1519bSnbeams   } else if (dim == 2) {
433a2968d6SJeremy L Thompson     const CeedInt elems_per_block = thread_1d < 4 ? 16 : 2;
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   } else if (dim == 3) {
493a2968d6SJeremy L Thompson     const CeedInt elems_per_block = thread_1d < 6 ? 4 : (thread_1d < 8 ? 2 : 1);
50b7453713SJeremy L Thompson 
513a2968d6SJeremy L Thompson     block_sizes[0] = thread_1d;
523a2968d6SJeremy L Thompson     block_sizes[1] = thread_1d;
539e201c85SYohann     block_sizes[2] = elems_per_block;
54b3e1519bSnbeams   }
55b3e1519bSnbeams   return CEED_ERROR_SUCCESS;
56b3e1519bSnbeams }
57b3e1519bSnbeams 
587d8d0e25Snbeams //------------------------------------------------------------------------------
594b3e95d5SJeremy L Thompson // Determine type of operator
604b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
614b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelData_Hip_gen(Ceed ceed, CeedInt num_input_fields, CeedOperatorField *op_input_fields,
624b3e95d5SJeremy L Thompson                                                CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, CeedOperatorField *op_output_fields,
634b3e95d5SJeremy L Thompson                                                CeedQFunctionField *qf_output_fields, CeedInt *max_P_1d, CeedInt *Q_1d, CeedInt *dim, bool *is_tensor,
644b3e95d5SJeremy L Thompson                                                bool *use_3d_slices) {
654b3e95d5SJeremy L Thompson   // Find dim, P_1d, Q_1d
664b3e95d5SJeremy L Thompson   *max_P_1d  = 0;
674b3e95d5SJeremy L Thompson   *Q_1d      = 0;
684b3e95d5SJeremy L Thompson   *dim       = 0;
694b3e95d5SJeremy L Thompson   *is_tensor = true;
709123fb08SJeremy L Thompson 
714b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
724b3e95d5SJeremy L Thompson     CeedBasis basis;
734b3e95d5SJeremy L Thompson 
744b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
754b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
764b3e95d5SJeremy L Thompson       bool    is_field_tensor;
774b3e95d5SJeremy L Thompson       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
784b3e95d5SJeremy L Thompson 
794b3e95d5SJeremy L Thompson       // Collect dim, P_1d, and Q_1d
804b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
814b3e95d5SJeremy L Thompson       *is_tensor = *is_tensor && is_field_tensor;
829123fb08SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
839123fb08SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P_1d));
844b3e95d5SJeremy L Thompson       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
854b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
864b3e95d5SJeremy L Thompson       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
874b3e95d5SJeremy L Thompson       *dim = field_dim;
889123fb08SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
899123fb08SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q_1d));
904b3e95d5SJeremy L Thompson       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
914b3e95d5SJeremy L Thompson       *Q_1d = field_Q_1d;
924b3e95d5SJeremy L Thompson     }
933a2968d6SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis));
944b3e95d5SJeremy L Thompson   }
954b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
964b3e95d5SJeremy L Thompson     CeedBasis basis;
974b3e95d5SJeremy L Thompson 
984b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
994b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
1004b3e95d5SJeremy L Thompson       bool    is_field_tensor;
1014b3e95d5SJeremy L Thompson       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
1024b3e95d5SJeremy L Thompson 
1034b3e95d5SJeremy L Thompson       // Collect dim, P_1d, and Q_1d
1044b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
1054b3e95d5SJeremy L Thompson       *is_tensor = *is_tensor && is_field_tensor;
1069123fb08SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
1079123fb08SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P_1d));
1084b3e95d5SJeremy L Thompson       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
1094b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
1104b3e95d5SJeremy L Thompson       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
1114b3e95d5SJeremy L Thompson       *dim = field_dim;
1129123fb08SJeremy L Thompson       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
1139123fb08SJeremy L Thompson       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q_1d));
1144b3e95d5SJeremy L Thompson       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
1154b3e95d5SJeremy L Thompson       *Q_1d = field_Q_1d;
1164b3e95d5SJeremy L Thompson     }
1173a2968d6SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis));
1184b3e95d5SJeremy L Thompson   }
1194b3e95d5SJeremy L Thompson 
1204b3e95d5SJeremy L Thompson   // Only use 3D collocated gradient parallelization strategy when gradient is computed
1214b3e95d5SJeremy L Thompson   *use_3d_slices = false;
1224b3e95d5SJeremy L Thompson   if (*dim == 3) {
1234b3e95d5SJeremy L Thompson     bool was_grad_found = false;
1244b3e95d5SJeremy L Thompson 
1254b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
1264b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
1274b3e95d5SJeremy L Thompson 
1284b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1294b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
1304b3e95d5SJeremy L Thompson         CeedBasis_Hip_shared *basis_data;
1314b3e95d5SJeremy L Thompson         CeedBasis             basis;
1324b3e95d5SJeremy L Thompson 
1334b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1344b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1354b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
1364b3e95d5SJeremy L Thompson         was_grad_found = true;
1373a2968d6SJeremy L Thompson         CeedCallBackend(CeedBasisDestroy(&basis));
1384b3e95d5SJeremy L Thompson       }
1394b3e95d5SJeremy L Thompson     }
1404b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
1414b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
1424b3e95d5SJeremy L Thompson 
1434b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1444b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
1454b3e95d5SJeremy L Thompson         CeedBasis_Hip_shared *basis_data;
1464b3e95d5SJeremy L Thompson         CeedBasis             basis;
1474b3e95d5SJeremy L Thompson 
1484b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1494b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1504b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
1514b3e95d5SJeremy L Thompson         was_grad_found = true;
1523a2968d6SJeremy L Thompson         CeedCallBackend(CeedBasisDestroy(&basis));
1534b3e95d5SJeremy L Thompson       }
1544b3e95d5SJeremy L Thompson     }
1554b3e95d5SJeremy L Thompson   }
1564b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1574b3e95d5SJeremy L Thompson }
1584b3e95d5SJeremy L Thompson 
1594b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1604b3e95d5SJeremy L Thompson // Setup fields
1614b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1624b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelFieldData_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedOperatorField op_field,
163*45a787f7SJeremy L Thompson                                                     CeedQFunctionField qf_field, FieldReuse_Hip field_reuse, CeedInt Q_1d, bool is_input,
164*45a787f7SJeremy L Thompson                                                     bool is_tensor, bool is_at_points, bool use_3d_slices) {
1654b3e95d5SJeremy L Thompson   std::string           var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
1669123fb08SJeremy L Thompson   std::string           P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q";
1674b3e95d5SJeremy L Thompson   std::string           option_name = (is_input ? "inputs" : "outputs");
1684b3e95d5SJeremy L Thompson   CeedEvalMode          eval_mode   = CEED_EVAL_NONE;
1694b3e95d5SJeremy L Thompson   CeedInt               elem_size = 0, num_comp = 0, P_1d = 0;
1704b3e95d5SJeremy L Thompson   CeedElemRestriction   elem_rstr;
1714b3e95d5SJeremy L Thompson   CeedBasis_Hip_shared *basis_data;
1724b3e95d5SJeremy L Thompson   CeedBasis             basis;
1734b3e95d5SJeremy L Thompson 
1749ee499e5SJeremy L Thompson   // Field reuse info
175*45a787f7SJeremy L Thompson   bool use_previous_field = field_reuse.index != -1;
1769ee499e5SJeremy L Thompson 
1774b3e95d5SJeremy L Thompson   code << "  // -- " << (is_input ? "Input" : "Output") << " field " << i << "\n";
1784b3e95d5SJeremy L Thompson 
1794b3e95d5SJeremy L Thompson   // Get field data
1804b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
1814b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
1824b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1834b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1844b3e95d5SJeremy L Thompson   }
1853a2968d6SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1864b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
1874b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
1884b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1899123fb08SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
1909123fb08SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
1914b3e95d5SJeremy L Thompson   }
1924b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
1934b3e95d5SJeremy L Thompson 
1944b3e95d5SJeremy L Thompson   // Set field constants
1954b3e95d5SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
1964b3e95d5SJeremy L Thompson     code << "  const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n";
1974b3e95d5SJeremy L Thompson     code << "  const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n";
1984b3e95d5SJeremy L Thompson   }
1994b3e95d5SJeremy L Thompson 
2004b3e95d5SJeremy L Thompson   // Load basis data
2014b3e95d5SJeremy L Thompson   code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
2024b3e95d5SJeremy L Thompson   switch (eval_mode) {
2034b3e95d5SJeremy L Thompson     case CEED_EVAL_NONE:
2044b3e95d5SJeremy L Thompson       break;
2054b3e95d5SJeremy L Thompson     case CEED_EVAL_INTERP:
2063a2968d6SJeremy L Thompson       if (is_at_points) {
2073a2968d6SJeremy L Thompson         // AtPoints
2083a2968d6SJeremy L Thompson         if (!basis_data->d_chebyshev_interp_1d) {
2093a2968d6SJeremy L Thompson           CeedSize    interp_bytes;
2103a2968d6SJeremy L Thompson           CeedScalar *chebyshev_interp_1d;
2113a2968d6SJeremy L Thompson 
2123a2968d6SJeremy L Thompson           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
2133a2968d6SJeremy L Thompson           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
2143a2968d6SJeremy L Thompson           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
2153a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis), hipMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
2163a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis),
2173a2968d6SJeremy L Thompson                       hipMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, hipMemcpyHostToDevice));
2183a2968d6SJeremy L Thompson           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
2193a2968d6SJeremy L Thompson         }
2203a2968d6SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
2213a2968d6SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
2223a2968d6SJeremy L Thompson       } else {
2233a2968d6SJeremy L Thompson         // Standard quadrature
2244b3e95d5SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
2254b3e95d5SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_interp_1d;
2263a2968d6SJeremy L Thompson       }
2279ee499e5SJeremy L Thompson       if (use_previous_field) {
228*45a787f7SJeremy L Thompson         std::string reuse_var = "s_B" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
2299ee499e5SJeremy L Thompson 
2309ee499e5SJeremy L Thompson         code << "  CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n";
2319ee499e5SJeremy L Thompson       } else {
2329123fb08SJeremy L Thompson         code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
233f815fac9SJeremy L Thompson         code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
2349ee499e5SJeremy L Thompson       }
2354b3e95d5SJeremy L Thompson       break;
2364b3e95d5SJeremy L Thompson     case CEED_EVAL_GRAD:
2373a2968d6SJeremy L Thompson       if (is_at_points) {
2383a2968d6SJeremy L Thompson         // AtPoints
2393a2968d6SJeremy L Thompson         if (!basis_data->d_chebyshev_interp_1d) {
2403a2968d6SJeremy L Thompson           CeedSize    interp_bytes;
2413a2968d6SJeremy L Thompson           CeedScalar *chebyshev_interp_1d;
2423a2968d6SJeremy L Thompson 
2433a2968d6SJeremy L Thompson           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
2443a2968d6SJeremy L Thompson           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
2453a2968d6SJeremy L Thompson           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
2463a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis), hipMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
2473a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis),
2483a2968d6SJeremy L Thompson                       hipMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, hipMemcpyHostToDevice));
2493a2968d6SJeremy L Thompson           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
2503a2968d6SJeremy L Thompson         }
2513a2968d6SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
2523a2968d6SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
2533a2968d6SJeremy L Thompson       } else {
2543a2968d6SJeremy L Thompson         // Standard quadrature
2554b3e95d5SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
2564b3e95d5SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_interp_1d;
2573a2968d6SJeremy L Thompson       }
2589123fb08SJeremy L Thompson       if (is_tensor) {
2599ee499e5SJeremy L Thompson         if (use_previous_field) {
260*45a787f7SJeremy L Thompson           std::string reuse_var = "s_B" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
2619ee499e5SJeremy L Thompson 
2629ee499e5SJeremy L Thompson           code << "  CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n";
2639ee499e5SJeremy L Thompson         } else {
2649123fb08SJeremy L Thompson           code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
265f815fac9SJeremy L Thompson           code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
2669123fb08SJeremy L Thompson         }
2679ee499e5SJeremy L Thompson       }
2683a2968d6SJeremy L Thompson       if (is_at_points) break;  // No G mat for AtPoints
2694b3e95d5SJeremy L Thompson       if (use_3d_slices) {
2704b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d;
2714b3e95d5SJeremy L Thompson         else data->G.outputs[i] = basis_data->d_collo_grad_1d;
272*45a787f7SJeremy L Thompson         if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD) {
273*45a787f7SJeremy L Thompson           std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
2749ee499e5SJeremy L Thompson 
2759ee499e5SJeremy L Thompson           code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
2769ee499e5SJeremy L Thompson         } else {
2779123fb08SJeremy L Thompson           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
278f815fac9SJeremy L Thompson           code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
2799ee499e5SJeremy L Thompson         }
2804b3e95d5SJeremy L Thompson       } else {
2814b3e95d5SJeremy L Thompson         bool has_collo_grad = basis_data->d_collo_grad_1d;
2824b3e95d5SJeremy L Thompson 
2834b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2844b3e95d5SJeremy L Thompson         else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2854b3e95d5SJeremy L Thompson         if (has_collo_grad) {
286*45a787f7SJeremy L Thompson           if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD) {
287*45a787f7SJeremy L Thompson             std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
2889ee499e5SJeremy L Thompson 
2899ee499e5SJeremy L Thompson             code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
2909ee499e5SJeremy L Thompson           } else {
2919123fb08SJeremy L Thompson             code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
292f815fac9SJeremy L Thompson             code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
2939ee499e5SJeremy L Thompson           }
2949ee499e5SJeremy L Thompson         } else {
295*45a787f7SJeremy L Thompson           if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD) {
296*45a787f7SJeremy L Thompson             std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
2979ee499e5SJeremy L Thompson 
2989ee499e5SJeremy L Thompson             code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
2994b3e95d5SJeremy L Thompson           } else {
3009123fb08SJeremy L Thompson             code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << P_name << "*" << Q_name << (is_tensor ? "" : "*dim") << "];\n";
3019123fb08SJeremy L Thompson             code << "  LoadMatrix<" << P_name << ", " << Q_name << (is_tensor ? "" : "*dim") << ">(data, G." << option_name << "[" << i << "], s_G"
3029123fb08SJeremy L Thompson                  << var_suffix << ");\n";
3034b3e95d5SJeremy L Thompson           }
3044b3e95d5SJeremy L Thompson         }
3059ee499e5SJeremy L Thompson       }
3064b3e95d5SJeremy L Thompson       break;
3074b3e95d5SJeremy L Thompson     case CEED_EVAL_WEIGHT:
3084b3e95d5SJeremy L Thompson       break;  // No action
3094b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
3104b3e95d5SJeremy L Thompson     case CEED_EVAL_DIV:
3114b3e95d5SJeremy L Thompson     case CEED_EVAL_CURL:
3124b3e95d5SJeremy L Thompson       break;  // TODO: Not implemented
3134b3e95d5SJeremy L Thompson               // LCOV_EXCL_STOP
3144b3e95d5SJeremy L Thompson   }
3153a2968d6SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
3164b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3174b3e95d5SJeremy L Thompson }
3184b3e95d5SJeremy L Thompson 
3194b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
3204b3e95d5SJeremy L Thompson // Restriction
3214b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
3224b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelRestriction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim,
323e93651e5SJeremy L Thompson                                                       CeedInt field_input_buffer[], CeedOperatorField op_field, CeedQFunctionField qf_field,
3249123fb08SJeremy L Thompson                                                       CeedInt Q_1d, bool is_input, bool is_tensor, bool is_at_points, bool use_3d_slices) {
3254b3e95d5SJeremy L Thompson   std::string              var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
3269123fb08SJeremy L Thompson   std::string              P_name     = (is_tensor ? "P_1d" : "P") + var_suffix;
3274b3e95d5SJeremy L Thompson   CeedEvalMode             eval_mode  = CEED_EVAL_NONE;
3284b3e95d5SJeremy L Thompson   CeedInt                  elem_size = 0, num_comp = 0, P_1d = 0;
3294b3e95d5SJeremy L Thompson   CeedSize                 l_size;
330f815fac9SJeremy L Thompson   CeedRestrictionType      rstr_type = CEED_RESTRICTION_STANDARD;
3314b3e95d5SJeremy L Thompson   CeedElemRestriction_Hip *rstr_data;
3324b3e95d5SJeremy L Thompson   CeedElemRestriction      elem_rstr;
3334b3e95d5SJeremy L Thompson   CeedBasis                basis;
3344b3e95d5SJeremy L Thompson 
3354b3e95d5SJeremy L Thompson   // Get field data
3364b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
3374b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
338f815fac9SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
3394b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
3404b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
3414b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
3424b3e95d5SJeremy L Thompson   }
3434b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
3444b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
3459123fb08SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
3469123fb08SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
3474b3e95d5SJeremy L Thompson   }
3483a2968d6SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
3494b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
3504b3e95d5SJeremy L Thompson 
3514b3e95d5SJeremy L Thompson   // Restriction
3524b3e95d5SJeremy L Thompson   if (is_input) {
3534b3e95d5SJeremy L Thompson     // Input
354e93651e5SJeremy L Thompson     if (field_input_buffer[i] != i) {
355e93651e5SJeremy L Thompson       std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]);
356e93651e5SJeremy L Thompson 
357e93651e5SJeremy L Thompson       // Restriction was already done for previous input
358e93651e5SJeremy L Thompson       code << "    CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n";
3593a2968d6SJeremy L Thompson     } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices && is_at_points)) {
3603a2968d6SJeremy L Thompson       if (eval_mode == CEED_EVAL_NONE && rstr_type != CEED_RESTRICTION_POINTS) {
361e93651e5SJeremy L Thompson         // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated
3624b3e95d5SJeremy L Thompson         code << "    CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n";
3633a2968d6SJeremy L Thompson       } else if (rstr_type != CEED_RESTRICTION_POINTS) {
364e93651e5SJeremy L Thompson         // Otherwise we're using the scratch space
365e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
366e93651e5SJeremy L Thompson       }
367f815fac9SJeremy L Thompson       switch (rstr_type) {
368f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STANDARD: {
3694b3e95d5SJeremy L Thompson           CeedInt comp_stride;
3704b3e95d5SJeremy L Thompson 
3714b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
3724b3e95d5SJeremy L Thompson           code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
3734b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3744b3e95d5SJeremy L Thompson           code << "    // CompStride: " << comp_stride << "\n";
3754b3e95d5SJeremy L Thompson           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
3769123fb08SJeremy L Thompson           code << "    ReadLVecStandard" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name
3779123fb08SJeremy L Thompson                << ">(data, l_size" << var_suffix << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n";
378f815fac9SJeremy L Thompson           break;
379f815fac9SJeremy L Thompson         }
380f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STRIDED: {
3814b3e95d5SJeremy L Thompson           bool    has_backend_strides;
3824b3e95d5SJeremy L Thompson           CeedInt num_elem;
3834b3e95d5SJeremy L Thompson 
3844b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
3854b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
3864b3e95d5SJeremy L Thompson           CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
3874b3e95d5SJeremy L Thompson 
3884b3e95d5SJeremy L Thompson           if (!has_backend_strides) {
3894b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
3904b3e95d5SJeremy L Thompson           }
3914b3e95d5SJeremy L Thompson           code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
3929123fb08SJeremy L Thompson           code << "    ReadLVecStrided" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", "
3939123fb08SJeremy L Thompson                << strides[1] << ", " << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n";
394f815fac9SJeremy L Thompson           break;
395f815fac9SJeremy L Thompson         }
3963a2968d6SJeremy L Thompson         case CEED_RESTRICTION_POINTS: {
3973a2968d6SJeremy L Thompson           CeedInt comp_stride;
3983a2968d6SJeremy L Thompson 
3993a2968d6SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
4003a2968d6SJeremy L Thompson           code << "    const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
4013a2968d6SJeremy L Thompson           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
4023a2968d6SJeremy L Thompson           break;
4033a2968d6SJeremy L Thompson         }
404f815fac9SJeremy L Thompson         // LCOV_EXCL_START
405f815fac9SJeremy L Thompson         case CEED_RESTRICTION_ORIENTED:
406f815fac9SJeremy L Thompson         case CEED_RESTRICTION_CURL_ORIENTED:
407f815fac9SJeremy L Thompson           break;  // TODO: Not implemented
408f815fac9SJeremy L Thompson                   // LCOV_EXCL_STOP
4094b3e95d5SJeremy L Thompson       }
4104b3e95d5SJeremy L Thompson     }
4114b3e95d5SJeremy L Thompson   } else {
4124b3e95d5SJeremy L Thompson     // Output
413f815fac9SJeremy L Thompson     switch (rstr_type) {
414f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STANDARD: {
4154b3e95d5SJeremy L Thompson         CeedInt comp_stride;
4164b3e95d5SJeremy L Thompson 
4174b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
4184b3e95d5SJeremy L Thompson         code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
4194b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
4204b3e95d5SJeremy L Thompson         code << "    // CompStride: " << comp_stride << "\n";
4214b3e95d5SJeremy L Thompson         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
4229123fb08SJeremy L Thompson         code << "    WriteLVecStandard" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name
4239123fb08SJeremy L Thompson              << ">(data, l_size" << var_suffix << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n";
424f815fac9SJeremy L Thompson         break;
425f815fac9SJeremy L Thompson       }
426f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STRIDED: {
4274b3e95d5SJeremy L Thompson         bool    has_backend_strides;
4284b3e95d5SJeremy L Thompson         CeedInt num_elem;
4294b3e95d5SJeremy L Thompson 
4304b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
4314b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
4324b3e95d5SJeremy L Thompson         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
4334b3e95d5SJeremy L Thompson 
4344b3e95d5SJeremy L Thompson         if (!has_backend_strides) {
4354b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
4364b3e95d5SJeremy L Thompson         }
4374b3e95d5SJeremy L Thompson         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
4389123fb08SJeremy L Thompson         code << "    WriteLVecStrided" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", "
4399123fb08SJeremy L Thompson              << strides[1] << ", " << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n";
440f815fac9SJeremy L Thompson         break;
441f815fac9SJeremy L Thompson       }
4423a2968d6SJeremy L Thompson       case CEED_RESTRICTION_POINTS:
4433a2968d6SJeremy L Thompson         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
4443a2968d6SJeremy L Thompson         break;
445f815fac9SJeremy L Thompson       // LCOV_EXCL_START
446f815fac9SJeremy L Thompson       case CEED_RESTRICTION_ORIENTED:
447f815fac9SJeremy L Thompson       case CEED_RESTRICTION_CURL_ORIENTED:
448f815fac9SJeremy L Thompson         break;  // TODO: Not implemented
449f815fac9SJeremy L Thompson                 // LCOV_EXCL_STOP
4504b3e95d5SJeremy L Thompson     }
4514b3e95d5SJeremy L Thompson   }
4523a2968d6SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
4534b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4544b3e95d5SJeremy L Thompson }
4554b3e95d5SJeremy L Thompson 
4564b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4574b3e95d5SJeremy L Thompson // Basis
4584b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4594b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelBasis_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim,
4609123fb08SJeremy L Thompson                                                 CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool is_tensor,
4613a2968d6SJeremy L Thompson                                                 bool is_at_points, bool use_3d_slices) {
4624b3e95d5SJeremy L Thompson   std::string         var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
4639123fb08SJeremy L Thompson   std::string         P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q";
4644b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
4654b3e95d5SJeremy L Thompson   CeedInt             elem_size = 0, num_comp = 0, P_1d = 0;
4664b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
4674b3e95d5SJeremy L Thompson   CeedBasis           basis;
4684b3e95d5SJeremy L Thompson 
4694b3e95d5SJeremy L Thompson   // Get field data
4704b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
4714b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
4724b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
4734b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
4744b3e95d5SJeremy L Thompson   }
4753a2968d6SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
4764b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
4774b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
4789123fb08SJeremy L Thompson     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
4799123fb08SJeremy L Thompson     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
4804b3e95d5SJeremy L Thompson   }
4814b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
4824b3e95d5SJeremy L Thompson 
4834b3e95d5SJeremy L Thompson   // Basis
4844b3e95d5SJeremy L Thompson   code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
4854b3e95d5SJeremy L Thompson   if (is_input) {
4864b3e95d5SJeremy L Thompson     switch (eval_mode) {
4874b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
4883a2968d6SJeremy L Thompson         if (!use_3d_slices && !is_at_points) {
4894b3e95d5SJeremy L Thompson           code << "    CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n";
4904b3e95d5SJeremy L Thompson         }
4914b3e95d5SJeremy L Thompson         break;
4924b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
4933a2968d6SJeremy L Thompson         if (is_at_points) {
4949123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
4959123fb08SJeremy L Thompson 
4969123fb08SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
4979123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
4989123fb08SJeremy L Thompson                << var_suffix << ", r_c" << var_suffix << ");\n";
4993a2968d6SJeremy L Thompson         } else {
5009123fb08SJeremy L Thompson           std::string function_name = is_tensor ? ((dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d") : "InterpNonTensor";
5019123fb08SJeremy L Thompson 
5029123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
5039123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
5049123fb08SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
5053a2968d6SJeremy L Thompson         }
5064b3e95d5SJeremy L Thompson         break;
5074b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
5083a2968d6SJeremy L Thompson         if (is_at_points) {
5099123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
5109123fb08SJeremy L Thompson 
5119123fb08SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
5129123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
5139123fb08SJeremy L Thompson                << var_suffix << ", r_c" << var_suffix << ");\n";
5143a2968d6SJeremy L Thompson         } else if (use_3d_slices) {
5159123fb08SJeremy L Thompson           std::string function_name = (dim > 1 ? "InterpTensor" : "Interp") + std::to_string(dim) + "d";
5169123fb08SJeremy L Thompson 
5174b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
5189123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
5199123fb08SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
5209123fb08SJeremy L Thompson         } else if (is_tensor) {
5219123fb08SJeremy L Thompson           bool        is_collocated = dim == 3 && Q_1d >= P_1d;
5229123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "Grad" : (is_collocated ? "GradTensorCollocated" : "GradTensor")) + std::to_string(dim) + "d";
5239123fb08SJeremy L Thompson 
5249123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << (dim >= 3 ? Q_name : "1") << "];\n";
5259123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B"
5269123fb08SJeremy L Thompson                << var_suffix << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
5274b3e95d5SJeremy L Thompson         } else {
5289123fb08SJeremy L Thompson           std::string function_name = "GradNonTensor";
5299123fb08SJeremy L Thompson 
5309123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
5319123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", dim, " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix
5329123fb08SJeremy L Thompson                << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
5334b3e95d5SJeremy L Thompson         }
5344b3e95d5SJeremy L Thompson         break;
5354b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
5363a2968d6SJeremy L Thompson         if (is_at_points) {
5373a2968d6SJeremy L Thompson           code << "    // Nothing to do AtPoints\n";
5383a2968d6SJeremy L Thompson         } else {
5394b3e95d5SJeremy L Thompson           CeedBasis_Hip_shared *basis_data;
5409123fb08SJeremy L Thompson           std::string           function_name = is_tensor ? ((dim == 1 ? "Weight" : "WeightTensor") + std::to_string(dim) + "d") : "WeightNonTensor";
5414b3e95d5SJeremy L Thompson 
5429123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
5434b3e95d5SJeremy L Thompson           CeedCallBackend(CeedBasisGetData(basis, &basis_data));
5444b3e95d5SJeremy L Thompson           data->W = basis_data->d_q_weight_1d;
5459123fb08SJeremy L Thompson           code << "    " << function_name << "<" << Q_name << ">(data, W, r_q" << var_suffix << ");\n";
5463a2968d6SJeremy L Thompson         }
5474b3e95d5SJeremy L Thompson         break;
5484b3e95d5SJeremy L Thompson       }
5494b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
5504b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
5514b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
5524b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
5534b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
5544b3e95d5SJeremy L Thompson     }
5554b3e95d5SJeremy L Thompson   } else {
5564b3e95d5SJeremy L Thompson     switch (eval_mode) {
5574b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
5584b3e95d5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n";
5594b3e95d5SJeremy L Thompson         break;  // No action
5604b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
561e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
5623a2968d6SJeremy L Thompson         if (is_at_points) {
5639123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
5649123fb08SJeremy L Thompson 
5659123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_c" << var_suffix << ", s_B"
5669123fb08SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5673a2968d6SJeremy L Thompson         } else {
5689123fb08SJeremy L Thompson           std::string function_name =
5699123fb08SJeremy L Thompson               is_tensor ? ((dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d") : "InterpTransposeNonTensor";
5709123fb08SJeremy L Thompson 
5719123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B"
5729123fb08SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5733a2968d6SJeremy L Thompson         }
5744b3e95d5SJeremy L Thompson         break;
5754b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
576e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
5773a2968d6SJeremy L Thompson         if (is_at_points) {
5789123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
5799123fb08SJeremy L Thompson 
5809123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_c" << var_suffix << ", s_B"
5819123fb08SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5823a2968d6SJeremy L Thompson         } else if (use_3d_slices) {
5839123fb08SJeremy L Thompson           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
5849123fb08SJeremy L Thompson 
5859123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B"
5869123fb08SJeremy L Thompson                << var_suffix << ", r_e" << var_suffix << ");\n";
5879123fb08SJeremy L Thompson         } else if (is_tensor) {
5889123fb08SJeremy L Thompson           bool        is_collocated = dim == 3 && Q_1d >= P_1d;
5899123fb08SJeremy L Thompson           std::string function_name =
5909123fb08SJeremy L Thompson               (dim == 1 ? "GradTranspose" : (is_collocated ? "GradTransposeTensorCollocated" : "GradTransposeTensor")) + std::to_string(dim) + "d";
5919123fb08SJeremy L Thompson 
5929123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B"
5939123fb08SJeremy L Thompson                << var_suffix << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
5944b3e95d5SJeremy L Thompson         } else {
5959123fb08SJeremy L Thompson           std::string function_name = "GradTransposeNonTensor";
5969123fb08SJeremy L Thompson 
5979123fb08SJeremy L Thompson           code << "    " << function_name << "<num_comp" << var_suffix << ", dim, " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix
5989123fb08SJeremy L Thompson                << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
5994b3e95d5SJeremy L Thompson         }
6004b3e95d5SJeremy L Thompson         break;
6014b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
6024b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT:
6034b3e95d5SJeremy L Thompson         break;  // Should not occur
6044b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
6054b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
6064b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
6074b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
6084b3e95d5SJeremy L Thompson     }
6094b3e95d5SJeremy L Thompson   }
6103a2968d6SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
6114b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6124b3e95d5SJeremy L Thompson }
6134b3e95d5SJeremy L Thompson 
6144b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
6154b3e95d5SJeremy L Thompson // QFunction
6164b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
6173a2968d6SJeremy L Thompson static int CeedOperatorBuildKernelQFunction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt dim, CeedInt max_num_points,
6183a2968d6SJeremy L Thompson                                                     CeedInt num_input_fields, CeedOperatorField *op_input_fields, CeedQFunctionField *qf_input_fields,
6194b3e95d5SJeremy L Thompson                                                     CeedInt num_output_fields, CeedOperatorField *op_output_fields,
6209123fb08SJeremy L Thompson                                                     CeedQFunctionField *qf_output_fields, std::string qfunction_name, CeedInt Q_1d, bool is_tensor,
6219123fb08SJeremy L Thompson                                                     bool is_at_points, bool use_3d_slices) {
6229123fb08SJeremy L Thompson   std::string         Q_name    = is_tensor ? "Q_1d" : "Q";
6234b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
6244b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
6254b3e95d5SJeremy L Thompson 
6268b97b69aSJeremy L Thompson   // Setup output arrays
6274b3e95d5SJeremy L Thompson   code << "\n    // -- Output field setup\n";
6284b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
6294b3e95d5SJeremy L Thompson     std::string var_suffix = "_out_" + std::to_string(i);
6304b3e95d5SJeremy L Thompson 
6314b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
6324b3e95d5SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
6333a2968d6SJeremy L Thompson     switch (eval_mode) {
6343a2968d6SJeremy L Thompson       case CEED_EVAL_NONE:
6353a2968d6SJeremy L Thompson         if (is_at_points) {
6363a2968d6SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "];\n";
6373a2968d6SJeremy L Thompson         } else {
6389123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
6394b3e95d5SJeremy L Thompson         }
6403a2968d6SJeremy L Thompson         break;
6413a2968d6SJeremy L Thompson       case CEED_EVAL_INTERP:
6423a2968d6SJeremy L Thompson         if (is_at_points) {
6433a2968d6SJeremy L Thompson           // Accumulator for point data
6449123fb08SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
6459123fb08SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "; i++) {\n";
6463a2968d6SJeremy L Thompson           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
6473a2968d6SJeremy L Thompson           code << "    }\n";
6483a2968d6SJeremy L Thompson         } else {
6499123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
6503a2968d6SJeremy L Thompson         }
6513a2968d6SJeremy L Thompson         break;
6523a2968d6SJeremy L Thompson       case CEED_EVAL_GRAD:
6533a2968d6SJeremy L Thompson         if (is_at_points) {
6543a2968d6SJeremy L Thompson           // Accumulator for point data
6559123fb08SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "*dim];\n";
6569123fb08SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "; i++) {\n";
6573a2968d6SJeremy L Thompson           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
6583a2968d6SJeremy L Thompson           code << "    }\n";
6593a2968d6SJeremy L Thompson         } else if (use_3d_slices) {
6604b3e95d5SJeremy L Thompson           // Accumulator for gradient slices
6614b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
6624b3e95d5SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n";
6634b3e95d5SJeremy L Thompson           code << "      r_q" << var_suffix << "[i] = 0.0;\n";
6644b3e95d5SJeremy L Thompson           code << "    }\n";
6654b3e95d5SJeremy L Thompson         } else {
6669123fb08SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
6674b3e95d5SJeremy L Thompson         }
6683a2968d6SJeremy L Thompson         break;
6693a2968d6SJeremy L Thompson       case CEED_EVAL_WEIGHT:
6703a2968d6SJeremy L Thompson         break;
6713a2968d6SJeremy L Thompson         // LCOV_EXCL_START
6723a2968d6SJeremy L Thompson       case CEED_EVAL_DIV:
6733a2968d6SJeremy L Thompson       case CEED_EVAL_CURL:
6743a2968d6SJeremy L Thompson         break;  // TODO: Not implemented
6753a2968d6SJeremy L Thompson                 // LCOV_EXCL_STOP
6764b3e95d5SJeremy L Thompson     }
6774b3e95d5SJeremy L Thompson   }
6784b3e95d5SJeremy L Thompson 
6793a2968d6SJeremy L Thompson   if (is_at_points) {
6803a2968d6SJeremy L Thompson     // We need to handle batches of points
6813a2968d6SJeremy L Thompson     code << "\n    // Note: Using batches of points\n";
6823a2968d6SJeremy L Thompson     code << "    const CeedInt point_loop_bound = (blockDim.x * blockDim.y) * ceil(1.0 * max_num_points / (blockDim.x * blockDim.y));\n\n";
6833a2968d6SJeremy L Thompson     code << "    #pragma unroll\n";
6843a2968d6SJeremy L Thompson     code << "    for (CeedInt i = threadIdx.x + threadIdx.y * blockDim.x; i < point_loop_bound; i += blockDim.x * blockDim.y) {\n";
6853a2968d6SJeremy L Thompson     code << "      const CeedInt p = i % max_num_points;\n\n";
6863a2968d6SJeremy L Thompson 
6873a2968d6SJeremy L Thompson     code << "      // -- Coordinates\n";
6883a2968d6SJeremy L Thompson     code << "      CeedScalar r_x[dim];\n";
6893a2968d6SJeremy 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";
6903a2968d6SJeremy L Thompson 
6913a2968d6SJeremy L Thompson     code << "      // -- Input fields\n";
6923a2968d6SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
6933a2968d6SJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
6943a2968d6SJeremy L Thompson 
6953a2968d6SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
6963a2968d6SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
6973a2968d6SJeremy L Thompson       // Basis action
6983a2968d6SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
6993a2968d6SJeremy L Thompson       switch (eval_mode) {
7003a2968d6SJeremy L Thompson         case CEED_EVAL_NONE:
7013a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7023a2968d6SJeremy L Thompson           code << "      ReadPoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
7033a2968d6SJeremy L Thompson                << ", max_num_points>(data, elem, p, max_num_points, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
7043a2968d6SJeremy L Thompson           break;
7053a2968d6SJeremy L Thompson         case CEED_EVAL_INTERP:
7063a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7073a2968d6SJeremy L Thompson           code << "      InterpAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix
7083a2968d6SJeremy L Thompson                << ", r_x, r_s" << var_suffix << ");\n";
7093a2968d6SJeremy L Thompson           break;
7103a2968d6SJeremy L Thompson         case CEED_EVAL_GRAD:
7113a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
7123a2968d6SJeremy L Thompson           code << "      GradAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix
7133a2968d6SJeremy L Thompson                << ", r_x, r_s" << var_suffix << ");\n";
7143a2968d6SJeremy L Thompson           break;
7153a2968d6SJeremy L Thompson         case CEED_EVAL_WEIGHT:
7163a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
7173a2968d6SJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = 1.0;\n";
7183a2968d6SJeremy L Thompson           break;
7193a2968d6SJeremy L Thompson           // LCOV_EXCL_START
7203a2968d6SJeremy L Thompson         case CEED_EVAL_DIV:
7213a2968d6SJeremy L Thompson         case CEED_EVAL_CURL:
7223a2968d6SJeremy L Thompson           break;  // TODO: Not implemented
7233a2968d6SJeremy L Thompson                   // LCOV_EXCL_STOP
7243a2968d6SJeremy L Thompson       }
7253a2968d6SJeremy L Thompson     }
7263a2968d6SJeremy L Thompson     code << "\n      // -- Output fields\n";
7273a2968d6SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
7283a2968d6SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
7293a2968d6SJeremy L Thompson 
7303a2968d6SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
7313a2968d6SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
7323a2968d6SJeremy L Thompson       // Basis action
7333a2968d6SJeremy L Thompson       switch (eval_mode) {
7343a2968d6SJeremy L Thompson         case CEED_EVAL_NONE:
7353a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7363a2968d6SJeremy L Thompson           break;
7373a2968d6SJeremy L Thompson         case CEED_EVAL_INTERP:
7383a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7393a2968d6SJeremy L Thompson           break;
7403a2968d6SJeremy L Thompson         case CEED_EVAL_GRAD:
7413a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
7423a2968d6SJeremy L Thompson           break;
7433a2968d6SJeremy L Thompson           // LCOV_EXCL_START
7443a2968d6SJeremy L Thompson         case CEED_EVAL_WEIGHT:
7453a2968d6SJeremy L Thompson           break;  // Should not occur
7463a2968d6SJeremy L Thompson         case CEED_EVAL_DIV:
7473a2968d6SJeremy L Thompson         case CEED_EVAL_CURL:
7483a2968d6SJeremy L Thompson           break;  // TODO: Not implemented
7493a2968d6SJeremy L Thompson                   // LCOV_EXCL_STOP
7503a2968d6SJeremy L Thompson       }
7513a2968d6SJeremy L Thompson     }
7523a2968d6SJeremy L Thompson 
7533a2968d6SJeremy L Thompson   } else if (use_3d_slices) {
7544b3e95d5SJeremy L Thompson     // We treat quadrature points per slice in 3d to save registers
7554b3e95d5SJeremy L Thompson     code << "\n    // Note: Using planes of 3D elements\n";
7564b3e95d5SJeremy L Thompson     code << "    #pragma unroll\n";
7574b3e95d5SJeremy L Thompson     code << "    for (CeedInt q = 0; q < " << Q_name << "; q++) {\n";
7584b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
7594b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
7604b3e95d5SJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
7614b3e95d5SJeremy L Thompson 
7624b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
7634b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
7644b3e95d5SJeremy L Thompson       // Basis action
7654b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
7664b3e95d5SJeremy L Thompson       switch (eval_mode) {
7674b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
7684b3e95d5SJeremy L Thompson           bool is_strided;
7694b3e95d5SJeremy L Thompson 
7704b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7714b3e95d5SJeremy L Thompson 
7724b3e95d5SJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
7734b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
7744b3e95d5SJeremy L Thompson           if (is_strided) {
7754b3e95d5SJeremy L Thompson             bool    has_backend_strides;
7764b3e95d5SJeremy L Thompson             CeedInt num_elem, elem_size;
7774b3e95d5SJeremy L Thompson 
7784b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
7794b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
7804b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
7814b3e95d5SJeremy L Thompson             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
7824b3e95d5SJeremy L Thompson 
7834b3e95d5SJeremy L Thompson             if (!has_backend_strides) {
7844b3e95d5SJeremy L Thompson               CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
7854b3e95d5SJeremy L Thompson             }
7864b3e95d5SJeremy L Thompson             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
787f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStrided3d<num_comp" << var_suffix << ", " << Q_name << ", " << strides[0] << ", " << strides[1] << ", "
7884b3e95d5SJeremy L Thompson                  << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n";
7894b3e95d5SJeremy L Thompson           } else {
7904b3e95d5SJeremy L Thompson             CeedSize                 l_size = 0;
7914b3e95d5SJeremy L Thompson             CeedInt                  comp_stride;
7924b3e95d5SJeremy L Thompson             CeedElemRestriction_Hip *rstr_data;
7934b3e95d5SJeremy L Thompson 
7944b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
7954b3e95d5SJeremy L Thompson             code << "      const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
7964b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
7974b3e95d5SJeremy L Thompson             code << "      // CompStride: " << comp_stride << "\n";
7984b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
7994b3e95d5SJeremy L Thompson             data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
800f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStandard3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix
8014b3e95d5SJeremy L Thompson                  << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
8024b3e95d5SJeremy L Thompson           }
8039123fb08SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
8044b3e95d5SJeremy L Thompson           break;
8054b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
8064b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
8074b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n";
8084b3e95d5SJeremy L Thompson           code << "        r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n";
8094b3e95d5SJeremy L Thompson           code << "      }\n";
8104b3e95d5SJeremy L Thompson           break;
8114b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
8124b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
813f815fac9SJeremy L Thompson           code << "      GradColloSlice3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_q" << var_suffix << ", s_G" << var_suffix
814f815fac9SJeremy L Thompson                << ", r_s" << var_suffix << ");\n";
8154b3e95d5SJeremy L Thompson           break;
8164b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
8174b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
8184b3e95d5SJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n";
8193a2968d6SJeremy L Thompson           break;
8204b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
8214b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
8224b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
8234b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
8244b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
8254b3e95d5SJeremy L Thompson       }
8264b3e95d5SJeremy L Thompson     }
8274b3e95d5SJeremy L Thompson     code << "\n      // -- Output fields\n";
8284b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
8294b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
8304b3e95d5SJeremy L Thompson 
8314b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
8324b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
8334b3e95d5SJeremy L Thompson       // Basis action
8344b3e95d5SJeremy L Thompson       switch (eval_mode) {
8354b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
8364b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
8373a2968d6SJeremy L Thompson           break;
8384b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
8394b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
8404b3e95d5SJeremy L Thompson           break;
8414b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
8424b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
8434b3e95d5SJeremy L Thompson           break;
8444b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
8454b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
8464b3e95d5SJeremy L Thompson           break;  // Should not occur
8474b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
8484b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
8494b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
8504b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
8514b3e95d5SJeremy L Thompson       }
8524b3e95d5SJeremy L Thompson     }
8534b3e95d5SJeremy L Thompson   } else {
8544b3e95d5SJeremy L Thompson     code << "\n    // Note: Using full elements\n";
8554b3e95d5SJeremy L Thompson     code << "    {\n";
8564b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
8574b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
8584b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
8594b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n";
8604b3e95d5SJeremy L Thompson     }
8614b3e95d5SJeremy L Thompson     code << "      // -- Output fields\n";
8624b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
8634b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
8644b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n";
8654b3e95d5SJeremy L Thompson     }
8664b3e95d5SJeremy L Thompson   }
8674b3e95d5SJeremy L Thompson 
8684b3e95d5SJeremy L Thompson   // Input and output buffers
8694b3e95d5SJeremy L Thompson   code << "\n      // -- QFunction inputs and outputs\n";
8704b3e95d5SJeremy L Thompson   code << "      // ---- Inputs\n";
8714b3e95d5SJeremy L Thompson   code << "      CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n";
8724b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
8734b3e95d5SJeremy L Thompson     code << "      // ------ Input field " << i << "\n";
8744b3e95d5SJeremy L Thompson     code << "      inputs[" << i << "] = r_s_in_" << i << ";\n";
8754b3e95d5SJeremy L Thompson   }
8764b3e95d5SJeremy L Thompson   code << "      // ---- Outputs\n";
8774b3e95d5SJeremy L Thompson   code << "      CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n";
8784b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
8794b3e95d5SJeremy L Thompson     code << "      // ------ Output field " << i << "\n";
8804b3e95d5SJeremy L Thompson     code << "      outputs[" << i << "] = r_s_out_" << i << ";\n";
8814b3e95d5SJeremy L Thompson   }
8824b3e95d5SJeremy L Thompson 
8834b3e95d5SJeremy L Thompson   // Apply QFunction
8844b3e95d5SJeremy L Thompson   code << "\n      // -- Apply QFunction\n";
8854b3e95d5SJeremy L Thompson   code << "      " << qfunction_name << "(ctx, ";
8869123fb08SJeremy L Thompson   if (dim != 3 || is_at_points || use_3d_slices || !is_tensor) {
8874b3e95d5SJeremy L Thompson     code << "1";
8884b3e95d5SJeremy L Thompson   } else {
8899123fb08SJeremy L Thompson     code << Q_name;
8904b3e95d5SJeremy L Thompson   }
8914b3e95d5SJeremy L Thompson   code << ", inputs, outputs);\n";
8924b3e95d5SJeremy L Thompson 
8933a2968d6SJeremy L Thompson   if (is_at_points) {
8943a2968d6SJeremy L Thompson     // Map back to coefficients
8953a2968d6SJeremy L Thompson     code << "\n      // -- Output fields\n";
8963a2968d6SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
8973a2968d6SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
8983a2968d6SJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
8993a2968d6SJeremy L Thompson 
9003a2968d6SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
9013a2968d6SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
9023a2968d6SJeremy L Thompson       // Basis action
9033a2968d6SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
9043a2968d6SJeremy L Thompson       switch (eval_mode) {
9053a2968d6SJeremy L Thompson         case CEED_EVAL_NONE: {
9063a2968d6SJeremy L Thompson           CeedInt             comp_stride;
9073a2968d6SJeremy L Thompson           CeedElemRestriction elem_rstr;
9083a2968d6SJeremy L Thompson 
9093a2968d6SJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
9103a2968d6SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
9113a2968d6SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
9123a2968d6SJeremy L Thompson           code << "      const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
9133a2968d6SJeremy L Thompson           code << "      WritePoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
9143a2968d6SJeremy L Thompson                << ", max_num_points>(data, elem, i, points.num_per_elem[elem], indices.outputs[" << i << "]"
9153a2968d6SJeremy L Thompson                << ", r_s" << var_suffix << ", d" << var_suffix << ");\n";
9163a2968d6SJeremy L Thompson           break;
9173a2968d6SJeremy L Thompson         }
9183a2968d6SJeremy L Thompson         case CEED_EVAL_INTERP:
9193a2968d6SJeremy L Thompson           code << "      if (i >= points.num_per_elem[elem]) {\n";
9203a2968d6SJeremy L Thompson           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) r_s" << var_suffix << "[j] = 0.0;\n";
9213a2968d6SJeremy L Thompson           code << "      }\n";
9223a2968d6SJeremy L Thompson           code << "      InterpTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s"
9233a2968d6SJeremy L Thompson                << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
9243a2968d6SJeremy L Thompson           break;
9253a2968d6SJeremy L Thompson         case CEED_EVAL_GRAD:
9263a2968d6SJeremy L Thompson           code << "      if (i >= points.num_per_elem[elem]) {\n";
9273a2968d6SJeremy L Thompson           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "*dim; j++) r_s" << var_suffix << "[j] = 0.0;\n";
9283a2968d6SJeremy L Thompson           code << "      }\n";
9293a2968d6SJeremy L Thompson           code << "      GradTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s"
9303a2968d6SJeremy L Thompson                << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
9313a2968d6SJeremy L Thompson           break;
9323a2968d6SJeremy L Thompson           // LCOV_EXCL_START
9333a2968d6SJeremy L Thompson         case CEED_EVAL_WEIGHT:
9343a2968d6SJeremy L Thompson           break;  // Should not occur
9353a2968d6SJeremy L Thompson         case CEED_EVAL_DIV:
9363a2968d6SJeremy L Thompson         case CEED_EVAL_CURL:
9373a2968d6SJeremy L Thompson           break;  // TODO: Not implemented
9383a2968d6SJeremy L Thompson                   // LCOV_EXCL_STOP
9393a2968d6SJeremy L Thompson       }
9403a2968d6SJeremy L Thompson     }
9413a2968d6SJeremy L Thompson   } else if (use_3d_slices) {
9424b3e95d5SJeremy L Thompson     // Copy or apply transpose grad, if needed
9433a2968d6SJeremy L Thompson     code << "\n      // -- Output fields\n";
9444b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
9454b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
9464b3e95d5SJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
9474b3e95d5SJeremy L Thompson 
9484b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
9494b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
9504b3e95d5SJeremy L Thompson       // Basis action
9514b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
9524b3e95d5SJeremy L Thompson       switch (eval_mode) {
9534b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
9544b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
9554b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
9564b3e95d5SJeremy L Thompson           code << "      }\n";
9573a2968d6SJeremy L Thompson           break;
9584b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
9594b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
9604b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
9614b3e95d5SJeremy L Thompson           code << "      }\n";
9624b3e95d5SJeremy L Thompson           break;
9634b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
964f815fac9SJeremy L Thompson           code << "      GradColloSliceTranspose3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_s" << var_suffix << ", s_G"
965f815fac9SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
9664b3e95d5SJeremy L Thompson           break;
9674b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
9684b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
9694b3e95d5SJeremy L Thompson           break;  // Should not occur
9704b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
9714b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
9724b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
9734b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
9744b3e95d5SJeremy L Thompson       }
9754b3e95d5SJeremy L Thompson     }
9764b3e95d5SJeremy L Thompson   }
9774b3e95d5SJeremy L Thompson   code << "    }\n";
9784b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9794b3e95d5SJeremy L Thompson }
9804b3e95d5SJeremy L Thompson 
9814b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
9829e201c85SYohann // Build single operator kernel
9837d8d0e25Snbeams //------------------------------------------------------------------------------
9848d12f40eSJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Hip_gen(CeedOperator op, bool *is_good_build) {
9853a2968d6SJeremy L Thompson   bool                   is_tensor = true, is_at_points = false, use_3d_slices = false;
9867d8d0e25Snbeams   Ceed                   ceed;
9873a2968d6SJeremy L Thompson   CeedInt                Q_1d, num_input_fields, num_output_fields, dim = 1, max_num_points = 0, coords_comp_stride = 0;
988b7453713SJeremy L Thompson   CeedQFunctionField    *qf_input_fields, *qf_output_fields;
989b7453713SJeremy L Thompson   CeedQFunction_Hip_gen *qf_data;
990b7453713SJeremy L Thompson   CeedQFunction          qf;
991b7453713SJeremy L Thompson   CeedOperatorField     *op_input_fields, *op_output_fields;
992b7453713SJeremy L Thompson   CeedOperator_Hip_gen  *data;
9934b3e95d5SJeremy L Thompson   std::ostringstream     code;
9944b3e95d5SJeremy L Thompson 
9958d12f40eSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
9964b3e95d5SJeremy L Thompson   {
9974b3e95d5SJeremy L Thompson     bool is_setup_done;
998b7453713SJeremy L Thompson 
999b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
10008d12f40eSJeremy L Thompson     if (is_setup_done) {
10018d12f40eSJeremy L Thompson       *is_good_build = !data->use_fallback;
10028d12f40eSJeremy L Thompson       return CEED_ERROR_SUCCESS;
10038d12f40eSJeremy L Thompson     }
10044b3e95d5SJeremy L Thompson   }
1005b7453713SJeremy L Thompson 
10068d12f40eSJeremy L Thompson   // Check field compatibility
10078d12f40eSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
10088d12f40eSJeremy L Thompson   {
10098d12f40eSJeremy L Thompson     bool has_shared_bases = true, is_all_tensor = true, is_all_nontensor = true;
10108d12f40eSJeremy L Thompson 
10118d12f40eSJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
10128d12f40eSJeremy L Thompson       CeedBasis basis;
10138d12f40eSJeremy L Thompson 
10148d12f40eSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
10158d12f40eSJeremy L Thompson       if (basis != CEED_BASIS_NONE) {
10168d12f40eSJeremy L Thompson         bool        is_tensor = true;
10178d12f40eSJeremy L Thompson         const char *resource;
10188d12f40eSJeremy L Thompson         char       *resource_root;
10198d12f40eSJeremy L Thompson         Ceed        basis_ceed;
10208d12f40eSJeremy L Thompson 
10218d12f40eSJeremy L Thompson         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
1022c9192acaSJeremy L Thompson         is_all_tensor    = is_all_tensor && is_tensor;
1023c9192acaSJeremy L Thompson         is_all_nontensor = is_all_nontensor && !is_tensor;
10248d12f40eSJeremy L Thompson         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
10258d12f40eSJeremy L Thompson         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
10268d12f40eSJeremy L Thompson         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
1027c9192acaSJeremy L Thompson         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/hip/shared");
10288d12f40eSJeremy L Thompson         CeedCallBackend(CeedFree(&resource_root));
10298d12f40eSJeremy L Thompson         CeedCallBackend(CeedDestroy(&basis_ceed));
10308d12f40eSJeremy L Thompson       }
10318d12f40eSJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis));
10328d12f40eSJeremy L Thompson     }
10338d12f40eSJeremy L Thompson 
10348d12f40eSJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
10358d12f40eSJeremy L Thompson       CeedBasis basis;
10368d12f40eSJeremy L Thompson 
10378d12f40eSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
10388d12f40eSJeremy L Thompson       if (basis != CEED_BASIS_NONE) {
10398d12f40eSJeremy L Thompson         bool        is_tensor = true;
10408d12f40eSJeremy L Thompson         const char *resource;
10418d12f40eSJeremy L Thompson         char       *resource_root;
10428d12f40eSJeremy L Thompson         Ceed        basis_ceed;
10438d12f40eSJeremy L Thompson 
10448d12f40eSJeremy L Thompson         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
1045c9192acaSJeremy L Thompson         is_all_tensor    = is_all_tensor && is_tensor;
1046c9192acaSJeremy L Thompson         is_all_nontensor = is_all_nontensor && !is_tensor;
10478d12f40eSJeremy L Thompson 
10488d12f40eSJeremy L Thompson         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
10498d12f40eSJeremy L Thompson         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
10508d12f40eSJeremy L Thompson         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
1051c9192acaSJeremy L Thompson         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/hip/shared");
10528d12f40eSJeremy L Thompson         CeedCallBackend(CeedFree(&resource_root));
10538d12f40eSJeremy L Thompson         CeedCallBackend(CeedDestroy(&basis_ceed));
10548d12f40eSJeremy L Thompson       }
10558d12f40eSJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis));
10568d12f40eSJeremy L Thompson     }
10578d12f40eSJeremy L Thompson     // -- Fallback to ref if not all bases are shared
10588d12f40eSJeremy L Thompson     if (!has_shared_bases || (!is_all_tensor && !is_all_nontensor)) {
10598d12f40eSJeremy L Thompson       *is_good_build = false;
10608d12f40eSJeremy L Thompson       return CEED_ERROR_SUCCESS;
10618d12f40eSJeremy L Thompson     }
10628d12f40eSJeremy L Thompson   }
1063b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1064b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1065b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
1066b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
10677d8d0e25Snbeams 
10684b3e95d5SJeremy L Thompson   // Get operator data
10693a2968d6SJeremy L Thompson   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
10704b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelData_Hip_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, op_output_fields,
10714b3e95d5SJeremy L Thompson                                                       qf_output_fields, &data->max_P_1d, &Q_1d, &dim, &is_tensor, &use_3d_slices));
10724b3e95d5SJeremy L Thompson   if (dim == 0) dim = 1;
10734b3e95d5SJeremy L Thompson   data->dim = dim;
10743a2968d6SJeremy L Thompson   if (is_at_points) {
10753a2968d6SJeremy L Thompson     CeedElemRestriction_Hip *rstr_data;
10763a2968d6SJeremy L Thompson     CeedElemRestriction      rstr_points = NULL;
10774b3e95d5SJeremy L Thompson 
10783a2968d6SJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
10793a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points));
10803a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetCompStride(rstr_points, &coords_comp_stride));
10813a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(rstr_points, &rstr_data));
10823a2968d6SJeremy L Thompson     data->points.indices = (CeedInt *)rstr_data->d_offsets;
10833a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
10843a2968d6SJeremy L Thompson   }
10853a2968d6SJeremy L Thompson   if (is_at_points) use_3d_slices = false;
10863a2968d6SJeremy L Thompson   if (Q_1d == 0) {
10873a2968d6SJeremy L Thompson     if (is_at_points) Q_1d = max_num_points;
10883a2968d6SJeremy L Thompson     else CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q_1d));
10894b3e95d5SJeremy L Thompson   }
10904b3e95d5SJeremy L Thompson   data->Q_1d = Q_1d;
10914b3e95d5SJeremy L Thompson 
10920b454692Sjeremylt   // Check for restriction only identity operator
10934b3e95d5SJeremy L Thompson   {
10944b3e95d5SJeremy L Thompson     bool is_identity_qf;
10954b3e95d5SJeremy L Thompson 
10962b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
10970b454692Sjeremylt     if (is_identity_qf) {
10989e201c85SYohann       CeedEvalMode eval_mode_in, eval_mode_out;
1099b7453713SJeremy L Thompson 
11002b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
11012b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
11026574a04fSJeremy L Thompson       CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
11036574a04fSJeremy L Thompson                 "Backend does not implement restriction only identity operators");
11040b454692Sjeremylt     }
11054b3e95d5SJeremy L Thompson   }
1106b2165e7aSSebastian Grimberg 
1107b2165e7aSSebastian Grimberg   // Load basis source files
11089123fb08SJeremy L Thompson   if (is_tensor) {
11099c25dd66SJeremy L Thompson     code << "// Tensor basis source\n";
11109c25dd66SJeremy L Thompson     code << "#include <ceed/jit-source/hip/hip-shared-basis-tensor-templates.h>\n\n";
11119123fb08SJeremy L Thompson   } else {
11129123fb08SJeremy L Thompson     code << "// Non-tensor basis source\n";
11139123fb08SJeremy L Thompson     code << "#include <ceed/jit-source/hip/hip-shared-basis-nontensor-templates.h>\n\n";
11149123fb08SJeremy L Thompson   }
11159123fb08SJeremy L Thompson   if (is_at_points) {
11163a2968d6SJeremy L Thompson     code << "// AtPoints basis source\n";
11173a2968d6SJeremy L Thompson     code << "#include <ceed/jit-source/hip/hip-shared-basis-tensor-at-points-templates.h>\n\n";
11189123fb08SJeremy L Thompson   }
11199c25dd66SJeremy L Thompson   code << "// CodeGen operator source\n";
11209c25dd66SJeremy L Thompson   code << "#include <ceed/jit-source/hip/hip-gen-templates.h>\n\n";
11217d8d0e25Snbeams 
11224b3e95d5SJeremy L Thompson   // Get QFunction name
11234b3e95d5SJeremy L Thompson   std::string qfunction_name(qf_data->qfunction_name);
11244b3e95d5SJeremy L Thompson   std::string operator_name;
11254b3e95d5SJeremy L Thompson 
112609095acaSJeremy L Thompson   operator_name = "CeedKernelHipGenOperator_" + qfunction_name;
11277d8d0e25Snbeams 
11289e201c85SYohann   // Define CEED_Q_VLA
11299e201c85SYohann   code << "\n#undef CEED_Q_VLA\n";
11309123fb08SJeremy L Thompson   if (dim != 3 || is_at_points || use_3d_slices || !is_tensor) {
11319e201c85SYohann     code << "#define CEED_Q_VLA 1\n\n";
11329e201c85SYohann   } else {
11339e201c85SYohann     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
11349e201c85SYohann   }
11359e201c85SYohann 
11364b3e95d5SJeremy L Thompson   // Add user QFunction source
11374b3e95d5SJeremy L Thompson   {
11389c25dd66SJeremy L Thompson     const char *source_path;
11394b3e95d5SJeremy L Thompson 
11409c25dd66SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path));
11419c25dd66SJeremy L Thompson     CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/hip/gen backend requires QFunction source code file");
11429c25dd66SJeremy L Thompson 
11439c25dd66SJeremy L Thompson     code << "// User QFunction source\n";
11449c25dd66SJeremy L Thompson     code << "#include \"" << source_path << "\"\n\n";
11454b3e95d5SJeremy L Thompson   }
11467d8d0e25Snbeams 
11477d8d0e25Snbeams   // Setup
11487d8d0e25Snbeams   code << "\n// -----------------------------------------------------------------------------\n";
11494b3e95d5SJeremy L Thompson   code << "// Operator Kernel\n";
11504b3e95d5SJeremy L Thompson   code << "// \n";
11514b3e95d5SJeremy L Thompson   code << "// d_[in,out]_i:   CeedVector device array\n";
11524b3e95d5SJeremy L Thompson   code << "// r_[in,out]_e_i: Element vector register\n";
11534b3e95d5SJeremy L Thompson   code << "// r_[in,out]_q_i: Quadrature space vector register\n";
11549123fb08SJeremy L Thompson   code << "// r_[in,out]_c_i: AtPoints Chebyshev coefficients register\n";
11554b3e95d5SJeremy L Thompson   code << "// r_[in,out]_s_i: Quadrature space slice vector register\n";
11564b3e95d5SJeremy L Thompson   code << "// \n";
11574b3e95d5SJeremy L Thompson   code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n";
11584b3e95d5SJeremy L Thompson   code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n";
11594b3e95d5SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n";
1160b3e1519bSnbeams   code << "\nextern \"C\" __launch_bounds__(BLOCK_SIZE)\n";
11612b730f8bSJeremy L Thompson   code << "__global__ void " << operator_name
11623a2968d6SJeremy 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";
11634b3e95d5SJeremy L Thompson 
11644b3e95d5SJeremy L Thompson   // Scratch buffers
11659e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
11664b3e95d5SJeremy L Thompson     CeedEvalMode eval_mode;
11674b3e95d5SJeremy L Thompson 
11682b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
11699e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
11704b3e95d5SJeremy L Thompson       code << "  const CeedScalar *d_in_" << i << " = fields.inputs[" << i << "];\n";
11717d8d0e25Snbeams     }
11727d8d0e25Snbeams   }
11739e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
11744b3e95d5SJeremy L Thompson     code << "  CeedScalar *d_out_" << i << " = fields.outputs[" << i << "];\n";
11757d8d0e25Snbeams   }
11767d8d0e25Snbeams 
11779e201c85SYohann   code << "  const CeedInt dim = " << dim << ";\n";
11789123fb08SJeremy L Thompson   code << "  const CeedInt " << (is_tensor ? "Q_1d" : "Q") << " = " << Q_1d << ";\n";
11793a2968d6SJeremy L Thompson   if (is_at_points) {
11803a2968d6SJeremy L Thompson     code << "  const CeedInt max_num_points = " << max_num_points << ";\n";
11813a2968d6SJeremy L Thompson     code << "  const CeedInt coords_comp_stride = " << coords_comp_stride << ";\n";
11823a2968d6SJeremy L Thompson   }
11837d8d0e25Snbeams 
11844b3e95d5SJeremy L Thompson   // Shared data
11854b3e95d5SJeremy L Thompson   code << "  extern __shared__ CeedScalar slice[];\n";
11869e201c85SYohann   code << "  SharedData_Hip data;\n";
11879e201c85SYohann   code << "  data.t_id_x = threadIdx.x;\n";
11889e201c85SYohann   code << "  data.t_id_y = threadIdx.y;\n";
11899e201c85SYohann   code << "  data.t_id_z = threadIdx.z;\n";
11909e201c85SYohann   code << "  data.t_id  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
11919123fb08SJeremy L Thompson   code << "  data.slice = slice + data.t_id_z*T_1D" << ((!is_tensor || dim == 1) ? "" : "*T_1D") << ";\n";
11927d8d0e25Snbeams 
11939ee499e5SJeremy L Thompson   // -- Determine input mat reuse
1194*45a787f7SJeremy L Thompson   FieldReuse_Hip input_matrix_reuse[CEED_FIELD_MAX];
11959ee499e5SJeremy L Thompson 
11969ee499e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1197*45a787f7SJeremy L Thompson     input_matrix_reuse[i].index = -1;
11989ee499e5SJeremy L Thompson   }
11999ee499e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
12009ee499e5SJeremy L Thompson     CeedEvalMode eval_mode_i;
12019ee499e5SJeremy L Thompson     CeedBasis    basis_i;
12029ee499e5SJeremy L Thompson 
12039ee499e5SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode_i));
12049ee499e5SJeremy L Thompson     if (eval_mode_i == CEED_EVAL_WEIGHT) continue;
12059ee499e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis_i));
1206*45a787f7SJeremy L Thompson     for (CeedInt j = 0; (input_matrix_reuse[i].index == -1) && (j < i); j++) {
12079ee499e5SJeremy L Thompson       CeedEvalMode eval_mode_j;
12089ee499e5SJeremy L Thompson       CeedBasis    basis_j;
12099ee499e5SJeremy L Thompson 
12109ee499e5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
12119ee499e5SJeremy L Thompson       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
12129ee499e5SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
12139ee499e5SJeremy L Thompson       if (basis_i == basis_j) {
12149ee499e5SJeremy L Thompson         if (is_tensor) {
1215*45a787f7SJeremy L Thompson           input_matrix_reuse[i].index     = j;
1216*45a787f7SJeremy L Thompson           input_matrix_reuse[i].is_input  = true;
1217*45a787f7SJeremy L Thompson           input_matrix_reuse[i].eval_mode = eval_mode_j;
12189ee499e5SJeremy L Thompson         } else {
12199ee499e5SJeremy L Thompson           // For non-tensor can only re-use with the same eval mode
12209ee499e5SJeremy L Thompson           if (eval_mode_i == eval_mode_j) {
1221*45a787f7SJeremy L Thompson             input_matrix_reuse[i].index     = j;
1222*45a787f7SJeremy L Thompson             input_matrix_reuse[i].is_input  = true;
1223*45a787f7SJeremy L Thompson             input_matrix_reuse[i].eval_mode = eval_mode_j;
12249ee499e5SJeremy L Thompson           }
12259ee499e5SJeremy L Thompson         }
12269ee499e5SJeremy L Thompson       }
12279ee499e5SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis_j));
12289ee499e5SJeremy L Thompson     }
12299ee499e5SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis_i));
12309ee499e5SJeremy L Thompson   }
12319ee499e5SJeremy L Thompson 
12329ee499e5SJeremy L Thompson   // -- Determine output mat reuse
1233*45a787f7SJeremy L Thompson   FieldReuse_Hip output_matrix_reuse[CEED_FIELD_MAX];
12349ee499e5SJeremy L Thompson 
12359ee499e5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1236*45a787f7SJeremy L Thompson     output_matrix_reuse[i].index = -1;
12379ee499e5SJeremy L Thompson   }
12389ee499e5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
12399ee499e5SJeremy L Thompson     CeedEvalMode eval_mode_i;
12409ee499e5SJeremy L Thompson     CeedBasis    basis_i;
12419ee499e5SJeremy L Thompson 
12429ee499e5SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode_i));
12439ee499e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis_i));
1244*45a787f7SJeremy L Thompson     for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < num_input_fields); j++) {
12459ee499e5SJeremy L Thompson       CeedEvalMode eval_mode_j;
12469ee499e5SJeremy L Thompson       CeedBasis    basis_j;
12479ee499e5SJeremy L Thompson 
12489ee499e5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
12499ee499e5SJeremy L Thompson       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
12509ee499e5SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
12519ee499e5SJeremy L Thompson       if (basis_i == basis_j) {
12529ee499e5SJeremy L Thompson         if (is_tensor) {
1253*45a787f7SJeremy L Thompson           output_matrix_reuse[i].index     = j;
1254*45a787f7SJeremy L Thompson           output_matrix_reuse[i].is_input  = true;
1255*45a787f7SJeremy L Thompson           output_matrix_reuse[i].eval_mode = eval_mode_j;
12569ee499e5SJeremy L Thompson         } else {
12579ee499e5SJeremy L Thompson           // For non-tensor can only re-use with the same eval mode
12589ee499e5SJeremy L Thompson           if (eval_mode_i == eval_mode_j) {
1259*45a787f7SJeremy L Thompson             output_matrix_reuse[i].index     = j;
1260*45a787f7SJeremy L Thompson             output_matrix_reuse[i].is_input  = true;
1261*45a787f7SJeremy L Thompson             output_matrix_reuse[i].eval_mode = eval_mode_j;
12629ee499e5SJeremy L Thompson           }
12639ee499e5SJeremy L Thompson         }
12649ee499e5SJeremy L Thompson       }
12659ee499e5SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis_j));
12669ee499e5SJeremy L Thompson     }
1267*45a787f7SJeremy L Thompson     for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < i); j++) {
12689ee499e5SJeremy L Thompson       CeedEvalMode eval_mode_j;
12699ee499e5SJeremy L Thompson       CeedBasis    basis_j;
12709ee499e5SJeremy L Thompson 
12719ee499e5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[j], &eval_mode_j));
12729ee499e5SJeremy L Thompson       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
12739ee499e5SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[j], &basis_j));
12749ee499e5SJeremy L Thompson       if (basis_i == basis_j) {
12759ee499e5SJeremy L Thompson         if (is_tensor) {
1276*45a787f7SJeremy L Thompson           output_matrix_reuse[i].index     = j;
1277*45a787f7SJeremy L Thompson           output_matrix_reuse[i].is_input  = false;
1278*45a787f7SJeremy L Thompson           output_matrix_reuse[i].eval_mode = eval_mode_j;
12799ee499e5SJeremy L Thompson         } else {
12809ee499e5SJeremy L Thompson           // For non-tensor can only re-use with the same eval mode
12819ee499e5SJeremy L Thompson           if (eval_mode_i == eval_mode_j) {
1282*45a787f7SJeremy L Thompson             output_matrix_reuse[i].index     = j;
1283*45a787f7SJeremy L Thompson             output_matrix_reuse[i].is_input  = false;
1284*45a787f7SJeremy L Thompson             output_matrix_reuse[i].eval_mode = eval_mode_j;
12859ee499e5SJeremy L Thompson           }
12869ee499e5SJeremy L Thompson         }
12879ee499e5SJeremy L Thompson       }
12889ee499e5SJeremy L Thompson       CeedCallBackend(CeedBasisDestroy(&basis_j));
12899ee499e5SJeremy L Thompson     }
12909ee499e5SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis_i));
12919ee499e5SJeremy L Thompson   }
12929ee499e5SJeremy L Thompson 
12937d8d0e25Snbeams   // Initialize constants, and matrices B and G
12944b3e95d5SJeremy L Thompson   code << "\n  // Input field constants and basis data\n";
12959e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
12969ee499e5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_input_fields[i], qf_input_fields[i], input_matrix_reuse[i], Q_1d, true,
12979ee499e5SJeremy L Thompson                                                              is_tensor, is_at_points, use_3d_slices));
12987d8d0e25Snbeams   }
12994b3e95d5SJeremy L Thompson   code << "\n  // Output field constants and basis data\n";
13009e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
13019ee499e5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_output_fields[i], qf_output_fields[i], output_matrix_reuse[i], Q_1d,
13029ee499e5SJeremy L Thompson                                                              false, is_tensor, is_at_points, use_3d_slices));
13034b3e95d5SJeremy L Thompson   }
13047d8d0e25Snbeams 
13054b3e95d5SJeremy L Thompson   // Loop over all elements
13064b3e95d5SJeremy L Thompson   code << "\n  // Element loop\n";
13077d8d0e25Snbeams   code << "  __syncthreads();\n";
13089e201c85SYohann   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
13094b3e95d5SJeremy L Thompson 
1310e93651e5SJeremy L Thompson   // -- Compute minimum buffer space needed
13113a2968d6SJeremy L Thompson   CeedInt max_rstr_buffer_size = 1;
1312e93651e5SJeremy L Thompson 
1313e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1314e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
1315e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
1316e93651e5SJeremy L Thompson 
1317e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
1318e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1319e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
13209123fb08SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_tensor && (dim >= 3) ? elem_size : 1));
1321681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1322e93651e5SJeremy L Thompson   }
1323e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1324e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
1325e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
1326e93651e5SJeremy L Thompson 
1327e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
1328e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1329e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
13309123fb08SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_tensor && (dim >= 3) ? elem_size : 1));
1331681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1332e93651e5SJeremy L Thompson   }
1333e93651e5SJeremy L Thompson   code << "    // Scratch restriction buffer space\n";
1334e93651e5SJeremy L Thompson   code << "    CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n";
1335e93651e5SJeremy L Thompson 
1336e93651e5SJeremy L Thompson   // -- Determine best input field processing order
1337e93651e5SJeremy L Thompson   CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX];
1338e93651e5SJeremy L Thompson 
1339e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1340e93651e5SJeremy L Thompson     field_rstr_in_buffer[i] = -1;
1341e93651e5SJeremy L Thompson     input_field_order[i]    = -1;
1342e93651e5SJeremy L Thompson   }
1343e93651e5SJeremy L Thompson   {
1344e93651e5SJeremy L Thompson     bool    is_ordered[CEED_FIELD_MAX];
1345e93651e5SJeremy L Thompson     CeedInt curr_index = 0;
1346e93651e5SJeremy L Thompson 
1347e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
1348e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
1349e93651e5SJeremy L Thompson       CeedVector          vec_i;
1350e93651e5SJeremy L Thompson       CeedElemRestriction rstr_i;
1351e93651e5SJeremy L Thompson 
1352e93651e5SJeremy L Thompson       if (is_ordered[i]) continue;
1353e93651e5SJeremy L Thompson       field_rstr_in_buffer[i]       = i;
1354e93651e5SJeremy L Thompson       is_ordered[i]                 = true;
1355e93651e5SJeremy L Thompson       input_field_order[curr_index] = i;
1356e93651e5SJeremy L Thompson       curr_index++;
1357034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
1358e93651e5SJeremy L Thompson       if (vec_i == CEED_VECTOR_NONE) continue;  // CEED_EVAL_WEIGHT
1359e93651e5SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
1360e93651e5SJeremy L Thompson       for (CeedInt j = i + 1; j < num_input_fields; j++) {
1361e93651e5SJeremy L Thompson         CeedVector          vec_j;
1362e93651e5SJeremy L Thompson         CeedElemRestriction rstr_j;
1363e93651e5SJeremy L Thompson 
1364e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
1365e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
1366e93651e5SJeremy L Thompson         if (rstr_i == rstr_j && vec_i == vec_j) {
1367e93651e5SJeremy L Thompson           field_rstr_in_buffer[j]       = i;
1368e93651e5SJeremy L Thompson           is_ordered[j]                 = true;
1369e93651e5SJeremy L Thompson           input_field_order[curr_index] = j;
1370e93651e5SJeremy L Thompson           curr_index++;
1371e93651e5SJeremy L Thompson         }
13723a2968d6SJeremy L Thompson         CeedCallBackend(CeedVectorDestroy(&vec_j));
13733a2968d6SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
1374e93651e5SJeremy L Thompson       }
13753a2968d6SJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec_i));
13763a2968d6SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
1377e93651e5SJeremy L Thompson     }
1378e93651e5SJeremy L Thompson   }
1379e93651e5SJeremy L Thompson 
13804b3e95d5SJeremy L Thompson   // -- Input restriction and basis
13813a2968d6SJeremy L Thompson   code << "\n    // -- Input field restrictions and basis actions\n";
13829e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1383e93651e5SJeremy L Thompson     CeedInt f = input_field_order[i];
1384e93651e5SJeremy L Thompson 
1385e93651e5SJeremy L Thompson     code << "    // ---- Input field " << f << "\n";
13867d8d0e25Snbeams 
13874b3e95d5SJeremy L Thompson     // ---- Restriction
1388e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Hip_gen(code, data, f, dim, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f], Q_1d,
13899123fb08SJeremy L Thompson                                                                true, is_tensor, is_at_points, use_3d_slices));
1390b7453713SJeremy L Thompson 
13914b3e95d5SJeremy L Thompson     // ---- Basis action
13929123fb08SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, f, dim, op_input_fields[f], qf_input_fields[f], Q_1d, true, is_tensor,
13939123fb08SJeremy L Thompson                                                          is_at_points, use_3d_slices));
13947d8d0e25Snbeams   }
13957d8d0e25Snbeams 
13964b3e95d5SJeremy L Thompson   // -- Q function
13973a2968d6SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelQFunction_Hip_gen(code, data, dim, max_num_points, num_input_fields, op_input_fields, qf_input_fields,
13989123fb08SJeremy L Thompson                                                            num_output_fields, op_output_fields, qf_output_fields, qfunction_name, Q_1d, is_tensor,
13999123fb08SJeremy L Thompson                                                            is_at_points, use_3d_slices));
14007d8d0e25Snbeams 
14014b3e95d5SJeremy L Thompson   // -- Output basis and restriction
14024b3e95d5SJeremy L Thompson   code << "\n    // -- Output field basis action and restrictions\n";
14039e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
14044b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
1405b7453713SJeremy L Thompson 
14064b3e95d5SJeremy L Thompson     // ---- Basis action
14079123fb08SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, i, dim, op_output_fields[i], qf_output_fields[i], Q_1d, false, is_tensor,
14089123fb08SJeremy L Thompson                                                          is_at_points, use_3d_slices));
14097d8d0e25Snbeams 
14104b3e95d5SJeremy L Thompson     // ---- Restriction
14113a2968d6SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Hip_gen(code, data, i, dim, NULL, op_output_fields[i], qf_output_fields[i], Q_1d, false,
14129123fb08SJeremy L Thompson                                                                is_tensor, is_at_points, use_3d_slices));
14137d8d0e25Snbeams   }
14147d8d0e25Snbeams 
14154b3e95d5SJeremy L Thompson   // Close loop and function
14167d8d0e25Snbeams   code << "  }\n";
14177d8d0e25Snbeams   code << "}\n";
14187d8d0e25Snbeams   code << "// -----------------------------------------------------------------------------\n\n";
14197d8d0e25Snbeams 
1420539ec17dSJeremy L Thompson   CeedInt block_sizes[3] = {0, 0, 0};
14219e201c85SYohann   CeedInt num_elem;
1422b7453713SJeremy L Thompson 
14233a2968d6SJeremy L Thompson   // Compile
14242b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
14259123fb08SJeremy L Thompson   CeedCallBackend(BlockGridCalculate_Hip_gen(is_tensor ? dim : 1, num_elem, data->max_P_1d, Q_1d, block_sizes));
14268d12f40eSJeremy L Thompson   {
14278d12f40eSJeremy L Thompson     bool is_compile_good = false;
14288d12f40eSJeremy L Thompson 
14298d12f40eSJeremy L Thompson     CeedCallBackend(CeedTryCompile_Hip(ceed, code.str().c_str(), &is_compile_good, &data->module, 2, "T_1D", block_sizes[0], "BLOCK_SIZE",
14302b730f8bSJeremy L Thompson                                        block_sizes[0] * block_sizes[1] * block_sizes[2]));
14318d12f40eSJeremy L Thompson     if (is_compile_good) {
14328d12f40eSJeremy L Thompson       *is_good_build = true;
1433eb7e6cafSJeremy L Thompson       CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, operator_name.c_str(), &data->op));
14348d12f40eSJeremy L Thompson     } else {
14358d12f40eSJeremy L Thompson       *is_good_build     = false;
14368d12f40eSJeremy L Thompson       data->use_fallback = true;
14378d12f40eSJeremy L Thompson     }
14388d12f40eSJeremy L Thompson   }
14392b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
14409bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
1441c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionDestroy(&qf));
1442e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14437d8d0e25Snbeams }
14442a86cc9dSSebastian Grimberg 
14457d8d0e25Snbeams //------------------------------------------------------------------------------
1446